/[PAMELA software]/DarthVader/ToFLevel2/inc/ToFLevel2.h
ViewVC logotype

Contents of /DarthVader/ToFLevel2/inc/ToFLevel2.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.24 - (show annotations) (download)
Mon Nov 23 09:50:48 2009 UTC (15 years ago) by mocchiut
Branch: MAIN
Changes since 1.23: +7 -42 lines
File MIME type: text/plain
Cleanup of both C++ and F77 routines

1 /**
2 * \file ToFLevel2.h
3 * \author Gianfranca DeRosa / Wolfgang Menn / Rita Carbone with E. M. supervision
4 */
5
6 #ifndef ToFLevel2_h
7 #define ToFLevel2_h
8 //
9 #include <TObject.h>
10 #include <TArrayI.h>
11 #include <TArrayF.h>
12 #include <TClonesArray.h>
13
14 #include <math.h> // EMILIANO
15 #include <iostream> // from ToFLevel2.cpp
16 #include <fstream> // Emiliano
17 #include <sstream> // Emiliano
18 #include <string> // Emiliano
19
20
21 #include <ToFStruct.h>
22
23 #include <TrkLevel2.h> // Emiliano
24 #include <TrigLevel2.h> // Emiliano
25 #include <GLTables.h> // Emiliano
26 #include <OrbitalInfo.h> // Emiliano
27 #include <ToFCore.h> // Emiliano
28 #include <physics/tof/TofEvent.h>
29
30 //
31 // Declaration of the core fortran routines
32 //
33 #define tofl2com tofl2com_
34 extern "C" int tofl2com();
35 #define toftrk toftrk_
36 extern "C" int toftrk();
37 #define rdtofcal rdtofcal_
38 extern "C" int rdtofcal(char [], int *);
39
40 //
41 // class which contains track related variables
42 //
43 #define ZTOF11 53.74
44 #define ZTOF12 53.04
45 #define ZTOF21 23.94
46 #define ZTOF22 23.44
47 #define ZTOF31 -23.49
48 #define ZTOF32 -24.34
49
50
51 class ToFGeom : public TObject {
52
53 private:
54 TArrayI ePlane, eXY;
55
56 public:
57 ToFGeom() {
58 int plane[24] = {
59 0, 0, 0, 0, 0, 0, 0, 0,
60 1, 1, 1, 1, 1, 1,
61 2, 2,
62 3, 3,
63 4, 4, 4,
64 5, 5, 5
65 };
66 int plXY[6]= { 2, 1, 1, 2, 2, 1 }; // X==1, Y==2 */
67 ePlane.Set(24,plane);
68 eXY.Set(6,plXY);
69 }
70
71 int GetPad( int idpmt) { return (int)((idpmt+0.5)/2.); }
72 int GetPlane( int idpmt) { return ePlane[ GetPad(idpmt) ]; }
73 int GetXY( int idpmt) { return eXY[ GetPlane(idpmt) ]; }
74
75 ClassDef(ToFGeom,1);
76
77 };
78
79
80 /**
81 * \brief Class which contains the PMT data
82 *
83 * If there is a valid ADC or a TDC value (value<4095) for a PMT, both ADC and TDC data
84 * are stored in the PMT class.
85 * Look in the ToFLevel2Ex.cxx example in the repository how to read the PMT class.
86 */
87 class ToFPMT : public TObject {
88
89 private:
90
91 public:
92 Int_t pmt_id; ///<the identification number of the PMT from 0 to 47
93 Float_t adc; ///<raw ADC value for this PMT
94 Float_t tdc; ///<raw TDC value for this PMT
95 Float_t tdc_tw; ///<time-walk corrected TDC value for this PMT
96 //
97 ToFPMT();
98 ToFPMT(const ToFPMT&);
99 //
100 ToFPMT* GetToFPMT(){return this;};
101 void Clear(Option_t *t="");
102
103
104
105 ClassDef(ToFPMT,2);
106 };
107
108 /**
109 * \brief Class used to calibrate adc to dEdx for each PMT
110 *
111 * Class used to calibrate adc to dEdx for each PMT
112 */
113 class ToFdEdx : public TObject {
114
115 private:
116 //
117 ToFGeom eGeom; // ToF geometry
118 //
119 Float_t adc_he;
120 TArrayF eDEDXpmt; // 0-47 pmt dEdx
121 // parameters:
122 TArrayF PMTsat; // 0-47 saturation parameters
123 Float_t adc[48];
124 //
125
126 TArrayF parAtt[48]; // 48 x 6
127 TArrayF parPos[48]; // 48 x 4
128 TArrayF parDesatBB[48]; // 48 x 3
129 TArrayF parBBneg[48]; // 48 x 3
130 TArrayF parBBpos; // 48 x 1
131
132 double f_adcPC( float x );
133 double f_BB( TArrayF &p, float x );
134 double f_BB5B( float x );
135 double f_att( TArrayF &p, float x ) ;
136 double f_att5B( float x );
137 double f_desatBB( TArrayF &p, float x );
138 double f_desatBB5B( float x );
139 double f_pos( TArrayF &p, float x );
140 double f_pos5B( float x );
141 float Get_adc_he( int id, float pl_x[6], float pl_y[6]);
142
143 Bool_t conn[12];
144
145 UInt_t ts[12];
146 UInt_t te[12];
147
148
149 public:
150 ToFdEdx(); // class constructor
151 ~ToFdEdx(){ Delete(); }; // class distructor
152 //
153 void Clear(Option_t *option="");
154 void Delete(Option_t *option="") { Clear(); }
155
156 void Init(pamela::tof::TofEvent *tofl0 ); // init parameters
157 void Init(Int_t i, Int_t j, Float_t adce);
158 void Define_PMTsat();
159
160 void ReadParAtt( const char *fname );
161 void ReadParPos( const char *fname );
162 void ReadParBBneg( const char *fname );
163 void ReadParBBpos( const char *fname );
164 void ReadParDesatBB( const char *fname );
165
166 void CheckConnectors(UInt_t atime, GL_PARAM *glparam, TSQLServer *dbc);
167
168 void Process( UInt_t atime, Float_t betamean, Float_t *xtr_tof, Float_t *ytr_tof); //
169 void Print(Option_t *option="");
170
171 Float_t GetdEdx_pmt(Int_t ipmt) { return eDEDXpmt[ipmt]; } // 0-47 dEdx for each PMT for tracked events
172 //
173 ToFdEdx* GetToFdEdx(){return this;};
174 ClassDef(ToFdEdx,2);
175 };
176
177
178 /**
179 * \brief Class which contains the tracker related variables
180 *
181 * We can use the ToF standalone to find hitted paddles, calculate beta, etc..
182 * These results are then stored with the "trkseqno" = -1.
183 * If we use the track from the tracker, then the penetration points in the
184 * scintillators are calculated, which defines the hitted paddles. For these paddles
185 * we calculate then all the output.
186 * Note: The artificial ADC values are stored as dEdx in the output, the dEdx will be
187 * by definition = 1. However, the artificial TDC values are just used internally
188 * and not stored in the output. But one can see in both cases which PMT has artificial
189 * values using "adcflag" and "tdcflag".
190 * Look in the ToFLevel2Ex.cxx example in the repository how to read the tracker related
191 * variables.
192 */
193 class ToFTrkVar : public TObject {
194
195 private:
196
197 public:
198 //
199 Int_t trkseqno; ///< tracker sequ. number: -1=ToF standalone, 0=first Tracker track, ...
200 //
201 Int_t npmttdc; ///<number of the TDC measurements used to evaluate beta
202 TArrayI pmttdc; ///<contains the ID (0..47) for the PMT used to evaluate beta
203 TArrayI tdcflag; ///<flag for artificial TDC, "0" if normal TDC value
204
205 /**
206 * \brief beta, 12 measurements for the 12 combinations, beta[13] is modified weighted mean
207 *
208 * The 12 measurements are S11-S31, S11-S32, S12-S31, S12-S32, and then analogue for
209 * S2-S3 and S1-S2.
210 * The calculation of beta[13] is now modified:
211 * We check the individual weights for artificial TDC values, then calculate
212 * am mean beta for the first time. In a second step we loop again through
213 * the single measurements, checking for the residual from the mean
214 * The cut on the residual reject measurements > "x"-sigma. A chi2 value is
215 * calculated, furthermore a "quality" value by adding the weights which
216 * are finally used. If all measurements are taken, "quality" will be = 505.
217 * A chi2 cut around 3-4 and a quality-cut > 400 is needed for clean beta
218 * The Level2 beta[12] which is derived in the fortran routines uses: 10.,200.,20.
219 * This is not a very high quality measurement. One can re-calculate a new beta[13]
220 * using the L2-method "CalcBeta"
221 */
222 Float_t beta[13];
223 //
224 Int_t npmtadc; ///<number of the ADC measurements used for dEdx evaluation
225 TArrayI pmtadc; ///<contains the ID (0..47) for the PMT used to evaluate dEdx
226 TArrayI adcflag; ///<flag for artificial ADCs, "0" if normal ADC value
227 TArrayF dedx; ///<energy loss for this PMT in mip
228 //
229 Float_t xtofpos[3]; ///<x-measurement using the TDC values and the calibration from S12, S21, S32
230 Float_t ytofpos[3]; ///<x-measurement using the TDC values and the calibration from S11, S22, S31
231 //
232 Float_t xtr_tof[6]; ///<x-measurement in the ToF layers from tracker
233 Float_t ytr_tof[6]; ///<x-measurement in the ToF layers from tracker
234 //
235 ToFTrkVar();
236 ToFTrkVar(const ToFTrkVar&);
237
238 ToFTrkVar* GetToFTrkVar(){return this;};
239 void Clear(Option_t *t="");
240
241 ClassDef(ToFTrkVar,1);
242 //
243 };
244
245 /**
246 * \brief Class to describe ToF LEVEL2 data
247 *
248 */
249
250 class ToFLevel2 : public TObject {
251 private:
252
253 public:
254 //
255 TClonesArray *PMT; ///<class needed to store PMT hit informations
256 TClonesArray *ToFTrk; ///<track related variable class
257 Int_t tof_j_flag[6]; ///<number of hitted paddle(s) for each ToF layer: flag = flag + 2**(paddlenumber-1)
258 //
259 Int_t unpackError;///< zero if no error presente
260 Int_t default_calib; ///< one if the default calibration has been used to process the data, zero otherwise
261 //
262 Float_t GetdEdx(Int_t notrack, Int_t plane, Int_t adcfl); // gf Apr 07
263
264 Float_t CalcBeta(Int_t notrack, Float_t resmax, Float_t qualitycut, Float_t chi2cut); // wm feb 08
265
266 //
267 // Float_t CalcBeta(Int_t notrack, Float_t resmax, Float_t chi2cut, Float_t qualitycut); // wm feb 08
268 //
269 // methods to make life simplier during the analysis, returns a pointer to the ToFTrkVar class containing track related variables
270 //
271 Int_t ntrk(){return ToFTrk->GetEntries();};
272 Int_t npmt(){return PMT->GetEntries();};
273
274 //
275 void GetLevel2Struct(cToFLevel2 *) const;
276 //
277 ToFTrkVar *GetToFTrkVar(Int_t notrack);
278 ToFPMT *GetToFPMT(Int_t nohit);
279 Int_t GetPMTid(Int_t gg, Int_t hh);
280 TString GetPMTName(Int_t ind);
281
282 Int_t GetPlaneIndex(Int_t pmt_id);
283 void GetMatrix(Int_t notrack, Float_t adc[4][12], Float_t tdc[4][12]);
284 void GetPMTIndex(Int_t pmt_id, Int_t &gg, Int_t &hh);
285
286 // gf Apr 07
287 void GetdEdxPaddle(Int_t notrack, Int_t paddleid, Int_t adcfl, Float_t &PadEdx, Int_t &SatWarning); // gf Apr 07
288 TString GetPMTName(Int_t ind, Int_t &iplane, Int_t &ipaddle,Int_t &ipmt);
289 Int_t GetPaddleIdOfTrack(Float_t xtr, Float_t ytr, Int_t plane); // gf Apr 07
290 Int_t GetPaddleIdOfTrack(Float_t xtr, Float_t ytr, Int_t plane, Float_t margin); // wm jun 2008
291 void GetPMTPaddle(Int_t pmt_id, Int_t &plane, Int_t &paddle); // gf Apr 07
292 void GetPaddlePMT(Int_t paddle, Int_t &pmtleft, Int_t &pmtright); // gf Apr 07
293 void GetPaddleGeometry(Int_t plane, Int_t paddle, Float_t &xleft, Float_t &xright, Float_t &yleft, Float_t &yright); // gf Apr 07
294 Int_t GetPaddleid(Int_t plane, Int_t paddle);
295 void GetPaddlePlane(Int_t padid, Int_t &plane, Int_t &paddle);
296 Int_t GetNPaddle(Int_t plane);
297 //
298 //
299 //
300 Int_t Process(TrkLevel2 *trk, TrigLevel2 *trg, GL_RUN *run, OrbitalInfo *orb, Bool_t force); // Emiliano
301
302 //
303 // constructor
304 //
305 ToFLevel2();
306 ~ToFLevel2(){Delete();}; //ELENA
307 void Delete(Option_t *t=""); //ELENA
308 void Set();//ELENA
309 //
310 //
311 ToFLevel2* GetToFLevel2(){return this;};
312
313 /**
314 * Method to get the z-position of the 6 TOF layers from the plane ID
315 * @param plane_id Plane ID (11 12 21 22 31 32)
316 */
317 Float_t GetZTOF(Int_t plane_id){
318 switch(plane_id){
319 case 11: return ZTOF11;
320 case 12: return ZTOF12;
321 case 21: return ZTOF21;
322 case 22: return ZTOF22;
323 case 31: return ZTOF31;
324 case 32: return ZTOF32;
325 default: return 0.;
326 };
327 };
328
329 //
330 // Paddles position
331 //
332 /*
333 S11 8 paddles 33.0 x 5.1 cm
334 S12 6 paddles 40.8 x 5.5 cm
335 S21 2 paddles 18.0 x 7.5 cm
336 S22 2 paddles 15.0 x 9.0 cm
337 S31 3 paddles 15.0 x 6.0 cm
338 S32 3 paddles 18.0 x 5.0 cm
339 */
340
341 Int_t GetToFPlaneID(Int_t ip);
342 Int_t GetToFPlaneIndex(Int_t plane_id);
343 Bool_t HitPaddle(Int_t ,Int_t);
344 Int_t GetNHitPaddles(Int_t plane);
345 void Clear(Option_t *t="");
346 //
347 ClassDef(ToFLevel2,4);
348 };
349
350 #endif
351

  ViewVC Help
Powered by ViewVC 1.1.23