StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StTofMaker.cxx
1 /***************************************************************************
2  *
3  * $Id: StTofMaker.cxx,v 1.20 2018/02/26 23:26:50 smirnovd Exp $
4  *
5  * Author: W.J. Llope / Wei-Ming Zhang / Frank Geurts
6  *
7  ***************************************************************************
8  *
9  * Description: TOF offline software: maker to insert TOFp raw data onto DST
10  *
11  **************************************************************************/
13 
36 #include "StTofMaker.h"
37 #include <stdlib.h>
38 #include "StMessMgr.h"
39 #include "StEventTypes.h"
40 #include "StEvent/StTofData.h"
41 #include "StEvent/StTofRawData.h"
42 #include "StTofUtil/StTofGeometry.h"
43 #include "StTofUtil/StTofDataCollection.h"
44 #include "StTofUtil/StTofRawDataCollection.h" //RunV
45 #include "StDaqLib/GENERIC/EventReader.hh"
46 #include "StDaqLib/TOF/TOF_Reader.hh"
47 #include "StDAQMaker/StDAQReader.h"
48 #include "tables/St_TofTag_Table.h"
49 #include "TFile.h"
50 #include "TH1.h"
51 
52 
53 
55 StTofMaker::StTofMaker(const char *name):StMaker(name) {
56  drawinit=kFALSE;
57  mTofGeom = 0;
58  }
59 
62 
63 
64 
67  // create histograms
68  nAdcHitHisto = new TH1S("tof_nadchit","Crude No. ADC Hits/Event",181,-0.5,180.5);
69  nTdcHitHisto = new TH1S("tof_ntdchit","Crude No. TDC Hits/Event",185,-0.5,184.5);
70  return StMaker::Init();
71 }
72 
73 
75 Int_t StTofMaker::InitRun(int runnumber){
76  LOG_INFO << "StTofMaker::InitRun -- initializing TofGeometry --" << endm;
77  mTofGeom = new StTofGeometry();
78  mTofGeom->init(this);
79  return kStOK;
80 }
81 
82 
83 
85 Int_t StTofMaker::FinishRun(int runnumber){
86  LOG_INFO << "StTofMaker::FinishRun -- cleaning up TofGeometry --" << endm;
87  if (mTofGeom) delete mTofGeom;
88  mTofGeom=0;
89  return 0;
90 }
91 
92 
93 
94 //_________________________________________________________________________
97  //
98  LOG_INFO << "StTofMaker::Make() -- All Your Base Are Belong To Us --" << endm;
99  StTofDataCollection myDataCollection;
100  mDataCollection = &myDataCollection;
101  mTofCollectionPresent = 0;
102  mDataCollectionPresent = 0;
103 
104  // RunV
105  StTofRawDataCollection myRawDataCollection;
106  mRawDataCollection = &myRawDataCollection;
107  mRawDataCollectionPresent = 0;
108 
109  //--- check for existence of collections.................
110  //
111  tofTag = -99;
112  mEvent = (StEvent *) GetInputDS("StEvent");
113  if (mEvent) {
114  LOG_INFO << "StTofMaker Checking for Tof and TofData collections ..." << endm;
115  mTheTofCollection = mEvent->tofCollection();
116  if(mTheTofCollection) {
117  mTofCollectionPresent = 1;
118  LOG_INFO << " + tofCollection Exists" << endm;
119 //--
120  if(mTheTofCollection->dataPresent()) {
121  LOG_INFO << " + tofDataCollection Exists" << endm;
122  mDataCollectionPresent = 1;
123  } else if(mTheTofCollection->rawdataPresent()) {
124  mRawDataCollectionPresent = 1; //RunV
125  }
126  else LOG_INFO << " - tofDataCollection & tofRawDataCollection does not exist" << endm;
127  //--
128  } else {
129  LOG_INFO << " - tofCollection does not exist" << endm;
130  LOG_INFO << " - tofDataCollection & tofRawDataCollection does not exist" << endm; //RunV
131  }
132  //-- create tofcollection...
133  if (!mTheTofCollection){
134  LOG_INFO << "StTofMaker Creating tofCollection in StEvent ... ";
135  mTheTofCollection = new StTofCollection();
136  mEvent->setTofCollection(mTheTofCollection);
137  mTheTofCollection = mEvent->tofCollection();
138  if (mTheTofCollection) LOG_INFO << "OK"<< endm;
139  else LOG_INFO << "FAILED" << endm;
140  } // end !mTheTofCollection check
141  } // end mEvent check
142  else LOG_INFO << "StTofMaker No StEvent structures detected ... not good!" << endm;
143  //
144  //--- end check for collections....
145 
146  // check if slatCollection not already available, get data if not
147  if(!mDataCollectionPresent) {
148  if (m_Mode ==0) { // default DAQ-reader setting
149  LOG_INFO << "StTofMaker Will get real data from TOFpDAQ Reader ... " << endm;
150  mTheTofData = GetDataSet("StDAQReader");
151  if (!mTheTofData) {
152  LOG_INFO << "StTofMaker No StDAQReader dataset. Event skipped" << endm;
153  return kStWarn;
154  }
155  mTheDataReader = (StDAQReader*)(mTheTofData->GetObject());
156  if (!mTheDataReader) {
157  LOG_INFO << "StTofMaker No StDAQReader object. Event skipped" << endm;
158  return kStWarn;
159  }
160  if (!(mTheDataReader->TOFPresent())) {
161  LOG_INFO << "StTofMaker TOF is not in datastream. Event skipped" << endm;
162  return kStWarn;
163  }
164  //
165  mTheTofReader = mTheDataReader->getTOFReader();
166  if (!mTheTofReader) {
167  LOG_INFO << "StTofMaker Failed to getTofReader()...." << endm;
168  return kStWarn;
169  }
170 
171 
172 
173 #ifdef TOFP_DEBUG
174  TOF_Reader* MyTest = dynamic_cast<TOF_Reader*>(mTheDataReader->getTOFReader());
175  if (MyTest) MyTest->printRawData();
176  if (MyTest->year2Data()) LOG_INFO << "StTofMaker: year2 data - TOFp+pVPD" << endm;
177  if (MyTest->year3Data()) LOG_INFO << "StTofMaker: year3 data - TOFp+pVPD+TOFr" << endm;
178  if (MyTest->year4Data()) LOG_INFO << "StTofMaker: year4 data - TOFp+pVPD+TOFrprime" << endm;
179  if (MyTest->year5Data()) LOG_INFO << "StTofMaker: year5 data - TOFr5+pVPD " << endm;
180 #endif
181 
182  //
183  //--- copy daq data to collection. (TOFp and pVPD)
184  //
185  if (mTheTofReader->year2Data()||mTheTofReader->year3Data()||mTheTofReader->year4Data()){
186  int nadchit=0;
187  int ntdchit=0;
188  int iStrobe = 0;
189 
190  // set limits on basic event statistics
191  int tdcHitMax=1600;
192  int pvpdStrobeMin=1500;
193  if (mTheTofReader->year3Data()||mTheTofReader->year4Data()){
194  tdcHitMax=880;
195  pvpdStrobeMin=950;
196  }
197 
198  for (int i=0;i<48;i++){
199  unsigned short slatid = mTofGeom->daqToSlatId(i);
200  unsigned short rawAdc = mTheTofReader->GetAdcFromSlat(i);
201  unsigned short rawTdc = mTheTofReader->GetTdcFromSlat(i);
202  short rawTc = 0;
203  if (i<32) rawTc = mTheTofReader->GetTc(i);
204  unsigned short rawSc = 0;
205  if (i<12) rawSc = mTheTofReader->GetSc(i);
206  // new StTofData
207  StTofData *rawTofData = new StTofData(slatid,rawAdc,rawTdc,rawTc,rawSc,0,0);
208  mDataCollection->push_back(rawTofData); // local collection
209  if(i<41) {
210  if (rawAdc> 50) nadchit++;
211  if (rawTdc<tdcHitMax) ntdchit++;
212  }
213  if(i>41) {
214  if (rawTdc>pvpdStrobeMin) iStrobe++;
215  }
216  }
217  if (mTheTofReader->year3Data()){
218  for (int i=48;i<132;i++){
219  //arbitrary id, start from at 100.
220  unsigned short id = (100-48)+i;
221  unsigned short rawAdc = mTheTofReader->GetAdc(i);
222  unsigned short rawTdc = 0;
223  if (i<120) rawTdc = mTheTofReader->GetTdc(i);
224  // new StTofData
225  StTofData *rawTofData = new StTofData(id,rawAdc,rawTdc,0,0,0,0);
226  mDataCollection->push_back(rawTofData);
227  if ((i>59) && (rawAdc>50)) nadchit++;
228  if ((i<120) && (rawTdc<tdcHitMax)) ntdchit++;
229  }
230  }
231 
232  if (mTheTofReader->year4Data()){
233  for (int i=48;i<184;i++){
234  //arbitrary id, start from at 100.
235  unsigned short id = (100-48)+i;
236  unsigned short rawAdc = 0;
237  if(i<180) {
238  rawAdc = mTheTofReader->GetAdc(i);
239  }
240 
241  unsigned short rawTdc = 0;
242  rawTdc = mTheTofReader->GetTdc(i);
243  // new StTofData
244  StTofData *rawTofData = new StTofData(id,rawAdc,rawTdc,0,0,0,0);
245  mDataCollection->push_back(rawTofData);
246  if ((i>59) && (rawAdc>50)) nadchit++;
247  if ((rawTdc<tdcHitMax)) ntdchit++;
248  }
249  }
250 
251  //--- end loop over data words...
252  //
253  //---
254  if (iStrobe>4) {
255  LOG_INFO << "StTofMaker ... Strobe Event:" << iStrobe << endm;
256  tofTag = -1; // set tag to show this was a strobe event
257  nAdcHitHisto->Fill(-1.);
258  nTdcHitHisto->Fill(-1.);
259  } else {
260  LOG_INFO << "StTofMaker ... Beam Event: ~" << nadchit << " ADC and ~"
261  << ntdchit << " TDC Raw Hits in TRAY" << endm;
262  tofTag = nadchit; // set tag to show this was physics event (tofTag>=0)
263  nAdcHitHisto->Fill(nadchit);
264  nTdcHitHisto->Fill(ntdchit);
265  }
266  //----
267  //
268  } // end year2-year4
269 
270  //------ RunV Raw Data Haidong Liu
271  if (mTheTofReader->year5Data()){
272  for (int i=0;i<(int)mTheTofReader->GetNLeadingHits();i++){ // Leading Edge
273  unsigned short flag; // 1 - leading; 2 - trailing
274  unsigned short chn=mTheTofReader->GetLeadingGlobalTdcChan(i);
275  unsigned int tdc=mTheTofReader->GetLeadingTdc(i);
276  flag=1;
277  unsigned short quality = 1;
278  StTofRawData *rawTofData1 = new StTofRawData(flag,chn,tdc,quality);
279  mRawDataCollection->push_back(rawTofData1);
280  }
281  for (int i=0;i<(int)mTheTofReader->GetNTrailingHits();i++){ // Trailing Edge
282  unsigned short flag; // 1 - leading; 2 - trailing
283  unsigned short chn=mTheTofReader->GetTrailingGlobalTdcChan(i);
284  unsigned int tdc=mTheTofReader->GetTrailingTdc(i);
285  flag=2;
286  unsigned short quality = 1;
287  StTofRawData *rawTofData1 = new StTofRawData(flag,chn,tdc,quality);
288  mRawDataCollection->push_back(rawTofData1);
289  }
290  for (int i=0;i<TOF_MAX_TDC_CHANNELS;i++){ //TOF_MAX_TDC_CHANNELS=198;
291  unsigned int rawLdTdc=mTheTofReader->GetLdTdc(i);
292  unsigned int rawTrTdc=mTheTofReader->GetTrTdc(i);
293  StTofData *rawTofData = new StTofData(i,0,0,0,0,rawLdTdc,rawTrTdc);
294  mDataCollection->push_back(rawTofData);
295  }
296 
297  } // end year5
298 
299  } else if(m_Mode == 1) // no action
300  LOG_INFO << "StTofMaker No special action required for DataCollection"<<endm;
301  else {
302  LOG_INFO << "WRONG SetMode (" << m_Mode << ")! "
303  << "... should be either 0 (DAQ reader) or 1(no action, SIM)" << endm;
304  return kStOK;
305  }
306  //--- end m_mode check....
307 
308  } // end dataPresent
309 
310  //--- end section if !mDataCollectionPresent
311 
312  LOG_INFO << "StTofMaker tofTag = " << tofTag << endm;
313  if (mEvent) storeTag();
314 
315  //--- fill StEvent, clean-up, and return...
316  //
317  if (mEvent) this->fillStEvent();
318  mDataCollection=0;
319  LOG_INFO << "StTofMaker::Make() -- see you next event --" << endm;
320 
321  return kStOK;
322 }
323 
324 
325 //_________________________________________________________________________
327 void StTofMaker::fillStEvent() {
328 
329  LOG_INFO << "StTofMaker::fillStEvent() Starting..." << endm;
330 
331 //--- make sure we have a tofcollection
332  if(!mTheTofCollection){
333  mTheTofCollection = new StTofCollection();
334  mEvent->setTofCollection(mTheTofCollection);
335  }
336 
337  //--- fill mTheTofCollection with ALL RAW DATA
338  if(mDataCollectionPresent != 1) {
339  LOG_INFO << "StTofMaker::fillStEvent() Size of mDataCollection = " << mDataCollection->size() << endm;
340  for(size_t jj = 0; jj < mDataCollection->size(); jj++)
341  mTheTofCollection->addData(mDataCollection->getData(jj));
342  }
343  if(mRawDataCollectionPresent != 1) {
344  LOG_INFO << "StTofMaker::fillStEvent() Size of mRawDataCollection = " << mRawDataCollection->size() << endm;
345  for(size_t jj = 0; jj < mRawDataCollection->size(); jj++)
346  mTheTofCollection->addRawData(mRawDataCollection->getRawData(jj));
347  }
348 
349  //--- verify existence of tofCollection in StEvent (mEvent)
350  //
351  // also check contents of StTofRaws filled in TofCollection. This kind of check
352  // has never been done for StTofMCSlat too! (We only verified in local
353  // collections) Here, we would just make sure that raw data is really in
354  // TofCollection. WMZ
355  //
356  LOG_INFO << "StTofMaker::fillStEvent() Verifying TOF StEvent data ..." << endm;
357  StTofCollection* mmTheTofCollection = mEvent->tofCollection();
358  if(mmTheTofCollection){
359  LOG_INFO << " + StEvent tofCollection Exists" << endm;
360  if(mmTheTofCollection->dataPresent()) {
361  LOG_INFO << " + StEvent TofDataCollection Exists" << endm;
362  StSPtrVecTofData& rawTofVec = mmTheTofCollection->tofData();
363 #ifdef TOFP_DEBUG
364  for (size_t i = 0; i < rawTofVec.size(); i++) {
365  StTofData* p = rawTofVec[i];
366  LOG_INFO << *p;
367  }
368 #endif
369  LOG_INFO << " StEvent TofDataCollection has " << rawTofVec.size() << " entries..." << endm;
370  }
371  else LOG_INFO << " - StEvent TofDataCollection does not Exist" << endm;
372  }
373  else {
374  LOG_INFO << " - StEvent tofCollection does not Exist" << endm;
375  LOG_INFO << " - StEvent tofDataCollection does not Exist" << endm;
376  }
377 }
378 
379 //_____________________________________________________________________________
381 void StTofMaker::storeTag(){
382  // instantiate new TofTag class
383  St_TofTag *tagTable = new St_TofTag("TofTag",1);
384  if (!tagTable) return;
385  tagTable->SetNRows(1);
386 
387  // add table to root .data directory
388  AddData(tagTable,".data");
389  TofTag_st *pTofTag = tagTable->GetTable();
390 
391  if (pTofTag) {
392  pTofTag->tofEventType = tofTag;
393  LOG_INFO << "StTofMaker::storeTag() tofTag stored" << endm;
394  }
395  else
396  LOG_INFO << "StTofMaker::storeTag() unable to store tofTag" << endm;
397 
398 }
399 
400 //_____________________________________________________________________________
403  // nope
404  return kStOK;
405 }
406 
407 /***************************************************************************
408  *
409  * $Log: StTofMaker.cxx,v $
410  * Revision 1.20 2018/02/26 23:26:50 smirnovd
411  * StTof: Remove outdated ClassImp macro
412  *
413  * Revision 1.19 2018/02/26 23:13:20 smirnovd
414  * Move embedded CVS log messages to the end of file
415  *
416  * Revision 1.18 2007/04/17 23:00:41 dongx
417  * replaced with standard STAR Loggers
418  *
419  * Revision 1.17 2005/04/12 17:33:18 dongx
420  * update for year 5 new data format. Store into TofRawData from now on.
421  *
422  * Revision 1.16 2004/09/19 00:09:02 perev
423  * Small Walgrind leak fixed
424  *
425  * Revision 1.15 2004/09/10 22:09:21 perev
426  * more defence agains corrupted DAQ data
427  *
428  * Revision 1.14 2004/01/27 23:17:01 dongx
429  * change for year4 run (pVPD+TOFp+TOFr')
430  * - Additional TOFr' ADC and TDC channels put in
431  * - Add TOTs of TOFr' in
432  *
433  *
434  * Revision 1.12 2003/09/17 19:54:28 geurts
435  * zeroed geometry pointer
436  *
437  * Revision 1.11 2003/08/08 00:20:41 geurts
438  * moved local collection code to StTofUtil, changed StTofMaker accordingly
439  *
440  * Revision 1.10 2003/07/18 18:31:49 perev
441  * test for nonexistance of XXXReader added
442  *
443  * Revision 1.9 2003/02/06 05:02:05 geurts
444  * Added TOFr and extra pVPD-ADC channels to the datastream:
445  * StTofMaker is now aware of year2 (TOFp+pVPD) and year3 (TOFp+pVPD+TOFr) raw data.
446  *
447  * Revision 1.8 2002/01/22 06:50:34 geurts
448  * modifications for STAR dBase access and doxygenized
449  *
450  * Revision 1.7 2001/10/09 03:06:38 geurts
451  * TofTag introduced
452  *
453  * Revision 1.6 2001/10/07 19:01:46 geurts
454  * change default operation mode to DAQ-reader
455  *
456  * Revision 1.5 2001/10/05 21:08:40 geurts
457  * clean-up histograms and private root file
458  *
459  * Revision 1.4 2001/09/28 18:40:03 llope
460  * first release
461  *
462  */
Int_t m_Mode
counters
Definition: StMaker.h:81
virtual Int_t Make()
Make method, check for collections; create and fill them according to m_Mode.
Definition: StTofMaker.cxx:96
TH1S * nAdcHitHisto
method for testing classes of StTofPidMaker
Definition: StTofMaker.h:98
virtual void AddData(TDataSet *data, const char *dir=".data")
User methods.
Definition: StMaker.cxx:332
virtual Int_t Finish()
default Finish method (empty)
Definition: StTofMaker.cxx:402
virtual ~StTofMaker()
default empty destructor
Definition: StTofMaker.cxx:61
Int_t InitRun(int)
InitRun method, (re)initialize TOFp data from STAR dBase.
Definition: StTofMaker.cxx:75
Int_t FinishRun(int)
FinishRun method, clean up TOFp dBase entries.
Definition: StTofMaker.cxx:85
Time-of-Flight Geometry Utilities.
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
Definition: Stypes.h:42
Definition: Stypes.h:40
StTofMaker(const char *name="tof")
default constructor
Definition: StTofMaker.cxx:55
virtual Int_t Init()
Init method, book histograms.
Definition: StTofMaker.cxx:66
void init()
initialize geometry class from XDF file and set-up DAQ/Slat mappings