StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StTofrMatchMaker.h
1 /*******************************************************************
2  *
3  * $Id: StTofrMatchMaker.h,v 1.17 2014/08/06 11:43:49 jeromel Exp $
4  *
5  * Author: Xin Dong
6  *****************************************************************
7  *
8  * Description: Tofr Match Maker to do the matching between the
9  * fired celles and TPC tracks ( similar to Frank's
10  * TofpMatchMaker )
11  *
12  *****************************************************************
13  *
14  * $Log: StTofrMatchMaker.h,v $
15  * Revision 1.17 2014/08/06 11:43:49 jeromel
16  * Suffix on literals need to be space (later gcc compiler makes it an error) - first wave of fixes
17  *
18  * Revision 1.16 2009/07/24 22:33:33 fine
19  * Make the code C++ compliant
20  *
21  * Revision 1.15 2009/06/09 19:45:35 jeromel
22  * Changes for BT#1428
23  *
24  * Revision 1.14 2008/07/23 19:22:03 dongx
25  * New track quality cuts for Run8
26  *
27  * Revision 1.13 2008/05/06 18:41:40 dongx
28  * - Fixed bug in ouput histogram filename switch
29  * - Added switch for tpc track tree output
30  *
31  * Revision 1.12 2008/03/27 18:12:12 dongx
32  * Update the HPTDC bin width to full precision
33  *
34  * Revision 1.11 2008/03/27 00:16:03 dongx
35  * update for Run8 finished.
36  *
37  * Revision 1.10 2007/11/29 22:43:12 dongx
38  * changed vpd trayId definition to 121 (East) and 122 (West)
39  *
40  * Revision 1.9 2007/11/22 00:22:37 dongx
41  * update for run8 - first version
42  *
43  * Revision 1.8 2007/02/28 23:32:00 dongx
44  * completion for Run V matching
45  * - trailing tdc and leading tdc stored as adc and tdc in StTofCell
46  * - multi-hit association cases use hit position priority
47  *
48  * Revision 1.7 2005/07/06 23:10:24 fisyak
49  * Use template StThreeVectorD
50  *
51  * Revision 1.6 2005/04/12 17:31:56 dongx
52  * update for year 5 data - not completed, leave as empty at present
53  *
54  * Revision 1.5 2004/05/03 23:08:50 dongx
55  * change according to the update of StTofrGeometry, save CPU time by 100 times
56  *
57  * Revision 1.4 2004/03/16 22:30:51 dongx
58  * fix the warning message when compiling
59  *
60  * Revision 1.3 2004/03/11 22:30:34 dongx
61  * -move m_Mode control to Init()
62  * -clear up
63  *
64  * Revision 1.2 2004/03/09 17:44:56 dongx
65  * first release
66  *
67  *
68  *******************************************************************/
69 #ifndef STTOFRMATCHMAKER_HH
70 #define STTOFRMATCHMAKER_HH
71 #include "StMaker.h"
72 #include "StThreeVectorD.hh"
73 
74 #define VHRBIN2PS 24.4140625 // Very High resolution mode, pico-second per bin
75  // 1000*25/1024 (ps/chn)
76 #define HRBIN2PS 97.65625 // High resolution mode, pico-second per bin
77  // 97.65625= 1000*100/1024 (ps/chn)
78 
79 #include <string>
80 #include <vector>
81 #ifndef ST_NO_NAMESPACES
82 using std::string;
83 using std::vector;
84 #endif
85 
86 class StEvent;
87 class StTrack;
88 class StGlobalTrack;
89 class StHelix;
90 #include "StThreeVectorD.hh"
91 class StTrackGeometry;
92 class StDcaGeometry;
93 #include "StTofUtil/StSortTofRawData.h"
94 class StTofINLCorr;
95 class StTofrGeometry;
96 class StTofCollection;
100 class StSPtrVecTofData;
101 class StTofrDaqMap;
102 class TH1D;
103 class TH2D;
104 class TTree;
105 
106 #if !defined(ST_NO_TEMPLATE_DEF_ARGS) || defined(__CINT__)
107 typedef vector<Int_t> IntVec;
108 typedef vector<Double_t> DoubleVec;
109 #else
110 typedef vector<Int_t, allocator<Int_t>> IntVec;
111 typedef vector<Double_t, allocator<Double_t>> DoubleVec;
112 #endif
113 
114 class StTofrMatchMaker : public StMaker {
115 public:
116  StTofrMatchMaker(const Char_t *name="tofrMatch");
117  ~StTofrMatchMaker();
118 
119  // void Clear(Option_t *option="");
120  Int_t Init();
121  Int_t InitRun(Int_t);
122  Int_t FinishRun(Int_t);
123  Int_t Make();
124  Int_t Finish();
125 
126  void setCreateHistoFlag(Bool_t histos=kTRUE);
127  void setCreateTreeFlag(Bool_t tree=kTRUE);
128  void setOuterTrackGeometry();
129  void setStandardTrackGeometry();
130  void setValidAdcRange(Int_t, Int_t);
131  void setValidTdcRange(Int_t, Int_t);
132  void setMinHitsPerTrack(Int_t);
133  void setMinFitPointsPerTrack(Int_t);
134  void setMinFitPointsOverMax(Float_t);
135  void setMaxDCA(Float_t);
136  void setHistoFileName(const Char_t*);
137  void setNtupleFileName(Char_t*);
138  void setSaveGeometry(Bool_t geomSave=kFALSE);
139 
140  Int_t processEventYear2to4();
141  Int_t processEventYear5();
142  Int_t processEventYear8();
143 
144 private:
145  StTrackGeometry* trackGeometry(StTrack*);
146  Int_t getTofData(StTofCollection*); // check, remap and fill local arrays with tof and pvpd data
147  Int_t storeMatchData(StTofCellCollection*, StTofCollection*);
148 
149  Bool_t strobeEvent(StSPtrVecTofData&);// check pVPD data for strobe event
150  void bookHistograms();
151  void writeHistogramsToFile();
152 
153  Bool_t validAdc(Float_t);
154  Bool_t validTdc(Float_t);
155  Bool_t validEvent(StEvent *);
156  Bool_t validTrack(StTrack*);
157  Bool_t validTrackRun8(StGlobalTrack*);
158  Bool_t validTofTrack(StTrack*);
159 
160  // year5 moved to calibration maker
161  // float GetINLcorr(int edgeid,int tempchan,int bin);
162  //
163 
164  //y8++ INL Corr moved here
165  float getINLCorr(int boardId, int tdcChanId);
166 
167 public:
168  Bool_t doPrintMemoryInfo;
169  Bool_t doPrintCpuInfo;
170 
171 private:
172  static const Int_t mDAQOVERFLOW;
173  static const Int_t mNTOFP;
174  static const Int_t mNPVPD;
175  static const Int_t mNTOFR;
176  static const Int_t mNTOFR5;
177  // static const Float_t mWidthPad = 3.45;
178 
179  static const Int_t mNTOF; // 192 for tof in Run 8++
180  static const Int_t mNModule; // 32 for tofr5++
181  static const Int_t mNCell;
182  static const Int_t mNVPD; //
183 
184  static const Int_t mEastVpdTrayId;
185  static const Int_t mWestVpdTrayId;
186 
187  static const Int_t mNValidTrays_Run3;
188  static const Int_t mNValidTrays_Run4;
189  static const Int_t mNValidTrays_Run5;
190  static const Int_t mNValidTrays_Run6;
191  static const Int_t mNValidTrays_Run7;
192  static const Int_t mNValidTrays_Run8;
193 
194  //year 5
195  static const Int_t mTdigBoard;
196  static const Int_t mTdcOnBoard;
197  static const Int_t mTdcChannel;
198 
199  Float_t mWidthPad;
200  vector<Float_t> mTofrAdc;
201  vector<Float_t> mTofrTdc;
202  vector<Float_t> mPvpdAdc;
203  vector<Float_t> mPvpdAdcLoRes;
204  vector<Float_t> mPvpdTdc;
205  //year 5
206  vector<Float_t> mPvpdToT;
207  IntVec mTofr5LdChan;
208  IntVec mTofr5LdTdc;
209  IntVec mTofr5TrChan;
210  IntVec mTofr5TrTdc;
211  vector<Float_t> mTofr5Tdc;
212  vector<Float_t> mTofr5ToT; // ToT as adc
213 
214  // moved to calibration maker
215  // Float_t mINLtable[mTdigBoard][mTdcOnBoard][mTdcChannel];
216  //
217 
218  IntVec mStrobeTdcMin;
219  IntVec mStrobeTdcMax;
220  DoubleVec mPedTOFr;
221 
222  StEvent *mEvent;
223  StTofrGeometry *mTofrGeom;
224  StTofrDaqMap *mDaqMap;
225  StSortTofRawData *mSortTofRawData; // sorted TOFr5 raw data
226  // tofr8++
227  StTofINLCorr *mTofINLCorr; // INL Correction
228 
229  Bool_t mHisto;
230  Bool_t mSaveTree;
231 
232  Bool_t mYear2;
233  Bool_t mYear3;
234  Bool_t mYear4;
235  Bool_t mYear5;
236  Bool_t mYear8;
237  Bool_t mYearX;
238 
239  Bool_t mOuterTrackGeometry;
240  Bool_t mGeometrySave;
241 
242  string mHistoFileName;
243 
244  // event counters
245  Int_t mEventCounter;
246  Int_t mAcceptedEventCounter;
247  Int_t mTofEventCounter;
248  Int_t mTofStrobeEventCounter;
249  Int_t mAcceptAndStrobe;
250  Int_t mAcceptAndBeam;
251 
252  // various cut-offs and ranges
253  Float_t mMinValidTdc;
254  Float_t mMaxValidTdc;
255  Float_t mMinValidAdc;
256  Float_t mMaxValidAdc;
257  unsigned int mMinHitsPerTrack;
258  unsigned int mMinFitPointsPerTrack;
259  Float_t mMinFitPointsOverMax;
260  Float_t mMaxDCA;
261 
262  //
263 
264  // TOFr histograms
265  TH2D* mADCTDCCorelation;
266 
267  TH1D* mEventCounterHisto;
268  TH1D* mCellsMultInEvent;
269  TH1D* mHitsMultInEvent;
270  TH1D* mHitsMultPerTrack;
271  TH1D* mDaqOccupancy;
272  TH1D* mDaqOccupancyValid;
273  TH1D* mDaqOccupancyProj;
274  TH2D* mHitsPosition;
275 
276  TH1D* mDaqOccupancyValidAll;
277  TH1D* mDaqOccupancyProjAll;
278  TH1D* mDaqOccupancyVpd;
279  TH1D* mDaqOccupancyValidVpd;
280 
281  vector<TH2D*> mHitCorr;
282  vector<TH2D*> mHitCorrModule;
283  TH2D* mHitCorrAll;
284 
285 
286  TH2D* mTrackPtEta;
287  TH2D* mTrackPtPhi;
288  TH1D* mTrackNFitPts;
289  TH2D* mTrackdEdxvsp;
290  TH2D* mNSigmaPivsPt;
291 
292  TTree* mTrackTree;
293 
294  TH1D* mCellsPerEventMatch1;
295  TH1D* mHitsPerEventMatch1;
296  TH1D* mCellsPerTrackMatch1;
297  TH1D* mTracksPerCellMatch1;
298  TH1D* mDaqOccupancyMatch1;
299  TH2D* mDeltaHitMatch1;
300 
301  TH1D* mCellsPerEventMatch2;
302  TH1D* mHitsPerEventMatch2;
303  TH1D* mCellsPerTrackMatch2;
304  TH1D* mTracksPerCellMatch2;
305  TH1D* mDaqOccupancyMatch2;
306  TH2D* mDeltaHitMatch2;
307 
308  TH1D* mCellsPerEventMatch3;
309  TH1D* mHitsPerEventMatch3;
310  TH1D* mCellsPerTrackMatch3;
311  TH1D* mTracksPerCellMatch3;
312  TH1D* mDaqOccupancyMatch3;
313  TH2D* mDeltaHitMatch3;
314 
315 #ifndef ST_NO_TEMPLATE_DEF_ARGS
316  typedef vector<Int_t> idVector;
317 #else
318  typedef vector<Int_t,allocator<Int_t>> idVector;
319 #endif
320  typedef idVector::iterator idVectorIter;
321 
322  struct StructCellHit{
323  Int_t channel;
324  Int_t tray;
325  Int_t module;
326  Int_t cell;
327  StThreeVectorD hitPosition;
328  idVector trackIdVec;
329  Int_t matchFlag;
330  Float_t zhit;
331  Float_t yhit;
332  };
333 
334  struct TRACKTREE{
335  Float_t pt;
336  Float_t eta;
337  Float_t phi;
338  Int_t nfitpts;
339  Float_t dEdx;
340  Int_t ndEdxpts;
341  Int_t charge;
342  Int_t projTrayId;
343  Int_t projCellChan;
344  Float_t projY;
345  Float_t projZ;
346  };
347  TRACKTREE trackTree;
348 
349 #ifndef ST_NO_TEMPLATE_DEF_ARGS
350  typedef vector<StructCellHit> tofCellHitVector;
351 #else
352  typedef vector<StructCellHit,allocator<StructCellHit>> tofCellHitVector;
353 #endif
354  typedef vector<StructCellHit>::iterator tofCellHitVectorIter;
355 
356 
357  virtual const char *GetCVS() const
358  {static const char cvs[]="Tag $Name: $ $Id: StTofrMatchMaker.h,v 1.17 2014/08/06 11:43:49 jeromel Exp $ built " __DATE__ " " __TIME__ ; return cvs;}
359 
360  ClassDef(StTofrMatchMaker,1)
361 };
362 
363 inline void StTofrMatchMaker::setValidAdcRange(Int_t min, Int_t max){
364  mMinValidAdc=min;
365  mMaxValidAdc=max;
366 }
367 
368 inline void StTofrMatchMaker::setValidTdcRange(Int_t min, Int_t max){
369  mMinValidTdc=min;
370  mMaxValidTdc=max;
371 }
372 
373 inline void StTofrMatchMaker::setOuterTrackGeometry(){mOuterTrackGeometry=true;}
374 
375 inline void StTofrMatchMaker::setStandardTrackGeometry(){mOuterTrackGeometry=false;}
376 
377 inline void StTofrMatchMaker::setMinHitsPerTrack(Int_t nhits){mMinHitsPerTrack=nhits;}
378 
379 inline void StTofrMatchMaker::setMinFitPointsPerTrack(Int_t nfitpnts){mMinFitPointsPerTrack=nfitpnts;}
380 
381 inline void StTofrMatchMaker::setMinFitPointsOverMax(Float_t ratio) {mMinFitPointsOverMax=ratio;}
382 
383 inline void StTofrMatchMaker::setMaxDCA(Float_t maxdca){mMaxDCA=maxdca;}
384 
385 inline void StTofrMatchMaker::setHistoFileName(const Char_t* filename){mHistoFileName=filename;}
386 
387 inline void StTofrMatchMaker::setCreateHistoFlag(Bool_t histos){mHisto = histos;}
388 
389 inline void StTofrMatchMaker::setCreateTreeFlag(Bool_t tree){mSaveTree = tree;}
390 
391 inline void StTofrMatchMaker::setSaveGeometry(Bool_t geomSave){mGeometrySave = geomSave; }
392 
393 inline Bool_t StTofrMatchMaker::validAdc(Float_t adc){return((adc>=mMinValidAdc) && (adc<=mMaxValidAdc));}
394 
395 inline Bool_t StTofrMatchMaker::validTdc(Float_t tdc){return((tdc>=mMinValidTdc) && (tdc<=mMaxValidTdc));}
396 
397 #endif