/[PAMELA software]/calo/flight/CaloProfile/inc/CaloProfile.h
ViewVC logotype

Contents of /calo/flight/CaloProfile/inc/CaloProfile.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.12 - (show annotations) (download)
Thu Sep 10 12:53:56 2009 UTC (15 years, 4 months ago) by mocchiut
Branch: MAIN
Changes since 1.11: +12 -1 lines
File MIME type: text/plain
Bugs fixed and new methods implemented

1 /**
2 * \file CaloProfile.h
3 * \author Emiliano Mocchiutti
4 */
5 #ifndef caloprofile_h
6 #define caloprofile_h
7
8 #define EMPTY -1
9 #define MAX 50
10
11 #include <iostream>
12 #include <stdio.h>
13 #include <string.h>
14 #include <ctype.h>
15 #include <TTree.h>
16 #include <TFriendElement.h>
17 #include <TChain.h>
18 #include <TCanvas.h>
19 #include <TGraph.h>
20 #include <TH1F.h>
21 #include <TH2F.h>
22 #include <TFile.h>
23 #include <TPolyLine.h>
24 #include <TMath.h>
25 #include <TStyle.h>
26 #include <TList.h>
27 #include <TKey.h>
28 #include <TSystemFile.h>
29 #include <TSystemDirectory.h>
30 #include <TSQLServer.h>
31 #include <TF1.h>
32 #include <TGraphErrors.h>
33
34 #include <PamLevel2.h>
35
36 using namespace std;
37
38 struct stack{
39 char data[MAX];
40 int top;
41 };
42
43 /**
44 *
45 */
46 class CaloLat : public TObject {
47
48 private:
49 //
50 PamLevel2 *L2;
51 Bool_t debug;
52 //
53 // needed to avoid reprocessing the same event over and over to obtain the variables
54 //
55 UInt_t OBT;
56 UInt_t PKT;
57 UInt_t atime;
58 //
59 Float_t estrip[2][22][96];
60 TString suf;
61 //
62
63 public:
64 //
65 //
66 void Draw();
67 void Draw(Int_t,Int_t);
68 //
69 CaloLat();
70 CaloLat(PamLevel2 *L2);
71 ~CaloLat(){ Delete(); };
72 //
73 void SetDebug(Bool_t d){ debug=d; };
74 //
75 void Clear();
76 void Clear(Option_t *option){Clear();};
77 void Delete();
78 void Delete(Option_t *option){Delete();};
79 //
80 void Process(); ///< Process data
81 void Print();
82 void Print(Option_t *option){Print();};
83 //
84 void SetSuffix(TString suffix){ suf = suffix;};
85 //
86 ClassDef(CaloLat,2);
87 };
88
89 /**
90 *
91 */
92 class CaloLong : public TObject {
93
94 private:
95 //
96 PamLevel2 *L2;
97 Bool_t debug;
98 //
99 // needed to avoid reprocessing the same event over and over to obtain the variables
100 //
101 UInt_t OBT;
102 UInt_t PKT;
103 UInt_t atime;
104 UInt_t fOBT;
105 UInt_t fPKT;
106 UInt_t fatime;
107 //
108 Int_t N;
109 Int_t NC;
110 Bool_t sel;
111 Bool_t cont;
112 Int_t mask18b;
113 //
114 Float_t chi2;
115 Float_t ndf;
116 Float_t E0;
117 Float_t a;
118 Float_t b;
119 Float_t errE0;
120 Float_t erra;
121 Float_t errb;
122 Float_t etmax;
123 Float_t asymm;
124 Float_t X0pl;
125 Float_t defE0;
126 Float_t umax;
127 Float_t lmax;
128 TString sumax;
129 TString slmax;
130 Int_t fitresult;
131 //
132 Bool_t no18x;
133 Bool_t maskXE;
134 Bool_t maskYE;
135 Bool_t maskXO;
136 Bool_t maskYO;
137 //
138 Bool_t xyaverage;
139 //
140 Bool_t heavytail;
141 Float_t letmax;
142 Float_t lmipth;
143 //
144 Float_t eplane[2][22];
145 //
146 CaloLevel2 *clp;
147 //
148 Float_t Evaluate(TString s, Float_t tmax, Float_t X0pl ); // expression must be of the form "tmax+2.*X0pl", "5*tmax"."tmax+10","tmax-(4*tmax)+3.*X0pl"...
149 //
150 TString suf;
151
152 public:
153 //
154 //
155 void Fit();
156 void Fit(Bool_t draw);
157 //
158 // Double_t ccurve(Double_t *t, Double_t *par);
159 //
160 void SetCaloLevel2Pointer(CaloLevel2 *cp){ clp = cp;};
161 //
162 Float_t Get_E0(){this->Fit(); return E0;};
163 Float_t Get_defE0(){this->Fit(); return defE0;};
164 Float_t Get_a(){this->Fit(); return a;};
165 Float_t Get_b(){this->Fit(); return b;};
166 Float_t Get_errE0(){this->Fit(); return errE0;};
167 Float_t Get_erra(){this->Fit(); return erra;};
168 Float_t Get_errb(){this->Fit(); return errb;};
169 Float_t Get_chi2(){this->Fit(); return chi2;};
170 Float_t Get_ndf(){this->Fit(); return ndf;};
171 Float_t Get_nchi2(){this->Fit(); if ( ndf > 0 ) return (chi2/ndf); return 0;};
172 Float_t Get_tmax(){this->Fit(); if ( b != 0 ) return ((a-1.)/b); return 0;};
173 Float_t Get_asymm(){this->Fit(); return asymm;};
174 Float_t Get_exptmax(){this->Process(); return etmax;};
175 Float_t Get_X0pl(){this->Process(); return X0pl;};
176 //
177 Float_t Get_letmax(){ return letmax;};
178 Float_t Get_lmipth(){ return lmipth;};
179 Int_t Get_fitresult(){this->Fit(); return fitresult;};
180 //
181 void ForceNextFit(){atime=0;fatime=0;};
182 void Draw();
183 void Draw(Int_t);
184 //
185 CaloLong();
186 CaloLong(PamLevel2 *L2);
187 ~CaloLong(){ Delete(); };
188 //
189 void SetDebug(Bool_t d){ debug=d; };
190 void UsePlane18X(){ no18x=false; };
191 //
192 void UseAverage(){ xyaverage = true;};
193 void UseAllMeas(){ xyaverage = false;};
194 //
195 void MaskSection(TString);
196 void UnMaskSection(TString);
197 void UnMaskSections();
198 void Selection(){sel = true; cont = false;}; ///< Set selection mode: planes from 1 to 22-N are used, plane 18 - N is masked if "emulate18" variable is true (DEFAULT);
199 void Contamination(){sel = false; cont = true;}; ///< Set contamination mode: planes from N to 22 are used.
200 void SetNoWpreSampler(Int_t n);
201 void SetNoWcalo(Int_t n);
202 void SplitInto(Int_t NoWpreSampler, Int_t NoWcalo);
203 Int_t GetNoWpreSampler(){return N;}; ///< Get the number of W planes used as presampler.
204 Int_t GetNoWcalo(){return NC;}; ///< Get the number of W planes used as calorimeter.
205 void SetEnergies(Float_t myene[][22]);
206 //
207 void SetLowerLimit(Float_t l){ lmax = l; };
208 void SetUpperLimit(Float_t u){ umax = u; };
209 void SetLowerLimit(TString sl){ slmax = sl; };// expression must be of the form "5*t"."t+10","t-(4*t)"... where t will be replaced by the fitted maximum (X0)
210 void SetUpperLimit(TString su){ sumax = su; };// expression must be of the form "5*t"."t+10","t-(4*t)"... where t will be replaced by the fitted maximum (X0)
211 //
212 void Setletmax(Float_t l){ letmax = l;};
213 void Setlmipth(Float_t l){ lmipth = l;};
214 void HeavyTail(Bool_t b){ heavytail=b;};
215 //
216 Float_t GetLowerLimit(){ return lmax;};
217 Float_t GetUpperLimit(){ return umax;};
218 //
219 void SetSuffix(TString suffix){ suf = suffix;};
220 //
221 void Clear();
222 void Clear(Option_t *option){Clear();};
223 void Delete();
224 void Delete(Option_t *option){Delete();};
225 //
226 void Process(); ///< Process data
227 void Print();
228 void Print(Option_t *option){Print();};
229 //
230 ClassDef(CaloLong,3);
231 };
232
233 /**
234 *
235 */
236 class Calo2D : public TObject {
237
238 private:
239 //
240 PamLevel2 *L2;
241 Bool_t debug;
242 //
243 // needed to avoid reprocessing the same event over and over to obtain the variables
244 //
245 UInt_t OBT;
246 UInt_t PKT;
247 UInt_t atime;
248 //
249 Float_t estrip[23][96][96];
250 Int_t smax[23];
251 Int_t smay[23];
252 //
253 TString suf;
254
255 public:
256 //
257 //
258 void Draw();
259 void Draw(Int_t);
260 //
261 Calo2D();
262 Calo2D(PamLevel2 *L2);
263 ~Calo2D(){ Delete(); };
264 //
265 void SetDebug(Bool_t d){ debug=d; };
266 //
267 void Clear();
268 void Clear(Option_t *option){Clear();};
269 void Delete();
270 void Delete(Option_t *option){Delete();};
271 //
272 void SetSuffix(TString suffix){ suf = suffix;};
273 //
274 void Process(); ///< Process data
275 void Print();
276 void Print(Option_t *option){Print();};
277 //
278 ClassDef(Calo2D,2);
279 };
280
281 #endif
282

  ViewVC Help
Powered by ViewVC 1.1.23