StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StMuDstMaker.cxx
1 /***************************************************************************
2  *
3  * $Id: StMuDstMaker.cxx,v 1.129 2019/02/21 13:32:54 jdb Exp $
4  * Author: Frank Laue, BNL, laue@bnl.gov
5  *
6  **************************************************************************/
7 #include "TRegexp.h"
8 #include "Stiostream.h"
9 #include "Stsstream.h"
10 #include "StChain.h"
11 #include "THack.h"
12 #include "StEvent/StEvent.h"
13 #include "StEvent/StTrack.h"
14 #include "StEvent/StTrackNode.h"
15 #include "StEvent/StRichSpectra.h"
16 #include "StEvent/StDetectorState.h"
17 #include "StEvent/StEventTypes.h"
18 #include "StEvent/StRunInfo.h"
19 #include "StEvent/StEventInfo.h"
20 #include "StEvent/StDcaGeometry.h"
21 #include "StEvent/StFgtCollection.h"
22 #include "StEvent/StFgtStrip.h"
23 #include "StEvent/StFgtHit.h"
24 #include "StEvent/StEnumerations.h"
25 #include "StFgtUtil/StFgtConsts.h"
26 #include "StEventUtilities/StuRefMult.hh"
27 #include "StEventUtilities/StuProbabilityPidAlgorithm.h"
28 
29 #include "StarClassLibrary/StPhysicalHelixD.hh"
30 #include "StarClassLibrary/StTimer.hh"
31 #include "StarClassLibrary/StMatrixF.hh"
32 
33 #include "StIOMaker/StIOMaker.h"
34 #include "StTreeMaker/StTreeMaker.h"
35 #ifndef __NO_STRANGE_MUDST__
37 #include "StStrangeMuDstMaker/StStrangeEvMuDst.hh"
38 #include "StStrangeMuDstMaker/StV0MuDst.hh"
39 #include "StStrangeMuDstMaker/StV0Mc.hh"
40 #include "StStrangeMuDstMaker/StXiMuDst.hh"
41 #include "StStrangeMuDstMaker/StXiMc.hh"
42 #include "StStrangeMuDstMaker/StKinkMuDst.hh"
43 #include "StStrangeMuDstMaker/StKinkMc.hh"
44 #include "StStrangeMuDstMaker/StStrangeCuts.hh"
45 #endif
46 
47 #include "StMuException.hh"
48 #include "StMuEvent.h"
49 #include "StMuPrimaryVertex.h"
50 #include "StMuRpsCollection.h"
51 #include "StMuMtdCollection.h"
52 #include "StMuMtdRawHit.h"
53 #include "StMuMtdHit.h"
54 #include "StMuMtdHeader.h"
55 #include "StMuTrack.h"
56 #include "StMuDebug.h"
57 #include "StMuCut.h"
58 #include "StMuFilter.h"
59 #include "StMuL3Filter.h"
60 #include "StMuChainMaker.h"
61 #include "StMuEmcCollection.h"
62 #include "StMuEmcUtil.h"
63 #include "StMuFmsCollection.h"
64 #include "StMuFmsUtil.h"
65 #include "StMuFmsHit.h"
66 #include "StMuRHICfCollection.h"
67 #include "StMuRHICfUtil.h"
68 #include "StMuRHICfHit.h"
69 #include "StMuRHICfPoint.h"
70 #include "StMuFcsCollection.h"
71 #include "StMuFcsUtil.h"
72 #include "StMuFcsHit.h"
73 #include "StMuFttCollection.h"
74 #include "StMuFttUtil.h"
75 #include "StMuFttRawHit.h"
76 #include "StMuFstCollection.h"
77 #include "StMuFwdTrack.h"
78 #include "StMuFwdTrackUtil.h"
79 #include "StMuFwdTrackCollection.h"
80 #include "StMuFstUtil.h"
81 #include "StMuFstRawHit.h"
82 #include "StMuFstHit.h"
83 #include "StMuEpdHit.h" // MALisa
84 #include "StMuEpdHitCollection.h" // MALisa
85 #include "StEvent/StEpdCollection.h" // MALisa
86 #include "StMuEpdUtil.h" // MALisa
87 #include "StMuETofCollection.h" // fseck
88 #include "StMuETofHeader.h" // fseck
89 #include "StMuETofDigi.h" // fseck
90 #include "StMuETofHit.h" // fseck
91 #include "StMuFgtStrip.h"
92 #include "StMuFgtCluster.h"
93 #include "StMuFgtStripAssociation.h"
94 #include "StMuFgtAdc.h"
95 #include "StMuPmdCollection.h"
96 #include "StMuPmdUtil.h"
97 #include "StMuPmdHit.h"
98 #include "StMuPmdCluster.h"
99 #include "StMuTofHit.h"
100 #include "StMuTofHitCollection.h"
101 #include "StMuTofUtil.h"
103 #include "StEvent/StBTofCollection.h"
104 #include "StEvent/StBTofRawHit.h"
105 #include "StEvent/StBTofHeader.h"
106 #include "StEvent/StETofCollection.h" //fseck
107 
108 #include "StMuBTofHit.h"
109 #include "StMuBTofHitCollection.h"
110 #include "StMuBTofUtil.h"
111 #include "StMuPrimaryTrackCovariance.h"
112 #include "StMuEzTree.h"
113 #include "EztEventHeader.h"
114 #include "EztEmcRawData.h"
115 #include "EztTrigBlob.h"
116 #include "EztFpdBlob.h"
117 
118 #include "StMuDstMaker.h"
119 #include "StMuDst.h"
120 
121 #include "TFile.h"
122 #include "TTree.h"
123 #include "TClass.h"
124 #include "TChain.h"
125 #include "TStreamerInfo.h"
126 #include "TClonesArray.h"
127 #include "TEventList.h"
128 
129 #include "THack.h"
130 #include "StMuMcVertex.h"
131 #include "StMuMcTrack.h"
132 #include "StG2TrackVertexMap.h"
133 
134 class StEpdHit; // MALisa
135 
136 ClassImp(StMuDstMaker)
137 
138 #if !(ST_NO_NAMESPACES)
139  using namespace units;
140 #endif
141 
142 //-----------------------------------------------------------------------
143 //-----------------------------------------------------------------------
144 //-----------------------------------------------------------------------
154 StMuDstMaker::StMuDstMaker(const char* name) : StIOInterFace(name),
155  mStEvent(0), mStMuDst(0),
156 #ifndef __NO_STRANGE_MUDST__
157  mStStrangeMuDstMaker(0),
158 #endif
159  mIOMaker(0), mTreeMaker(0),
160  mIoMode(1), mIoNameMode((int)ioTreeMaker), mEventList(0),
161  mTrackType(256), mReadTracks(1),
162  mReadV0s(1), mReadXis(1), mReadKinks(1), mFinish(0),
163  mTrackFilter(0), mL3TrackFilter(0),
164  mCurrentFile(0),
165  mChain (0), mTTree(0),
166  mSplit(99), mCompression(9), mBufferSize(65536*4), mVtxList(100),
167  mProbabilityPidAlgorithm(0), mEmcCollectionArray(0), mEmcCollection(0),
168  mFmsCollection(0),mRHICfCollection(0),mFcsCollection(0),mFttCollection(0),mFstCollection(0),mFwdTrackCollection(0), mPmdCollectionArray(0), mPmdCollection(0)
169 
170 {
171  assignArrays();
172 
173  mDirName="./";
174  mFileName="";
175  streamerOff();
176  zeroArrays();
177  if (mIoMode==ioRead) openRead();
178  if (mIoMode==ioWrite) mProbabilityPidAlgorithm = new StuProbabilityPidAlgorithm();
179 
180  mEventCounter=0;
181  mStMuDst = new StMuDst();
182  mEmcUtil = new StMuEmcUtil();
183  mFmsUtil = new StMuFmsUtil();
184  mRHICfUtil = new StMuRHICfUtil();
185  mFcsUtil = new StMuFcsUtil();
186  mFttUtil = new StMuFttUtil();
187  mFstUtil = new StMuFstUtil();
188  mPmdUtil = new StMuPmdUtil();
189  mTofUtil = new StMuTofUtil();
190  mBTofUtil = new StMuBTofUtil();
191  mEpdUtil = new StMuEpdUtil();
192  mEzTree = new StMuEzTree();
193  if ( ! mStMuDst || ! mEmcUtil || ! mFmsUtil || ! mRHICfUtil || !mFcsUtil || ! mPmdUtil || ! mTofUtil || ! mBTofUtil || ! mEpdUtil || ! mEzTree )
194  throw StMuExceptionNullPointer("StMuDstMaker:: constructor. Something went horribly wrong, cannot allocate pointers",__PRETTYF__);
195 
196 
197  createArrays();
198 
199 
201  StMuL3Filter* l3Filter = new StMuL3Filter(); setL3TrackFilter(l3Filter);
202  StMuFilter* filter = new StMuFilter(); setTrackFilter(filter);
203  FORCEDDEBUGMESSAGE("ATTENTION: use standard MuFilter");
204  FORCEDDEBUGMESSAGE("ATTENTION: use standard l3 MuFilter");
205 
206 
207 }
208 
217 {
218  mArrays = mAArrays + 0;
219 #ifndef __NO_STRANGE_MUDST__
220  mStrangeArrays = mArrays + __NARRAYS__;
221  mMCArrays = mStrangeArrays + __NSTRANGEARRAYS__;
222 #else
223  mMCArrays = mArrays + __NARRAYS__;
224 #endif
225  mEmcArrays = mMCArrays + __NMCARRAYS__;
226  mPmdArrays = mEmcArrays + __NEMCARRAYS__;
227  mFmsArrays = mPmdArrays + __NPMDARRAYS__;
228  mRHICfArrays = mFmsArrays + __NFMSARRAYS__;
229  mFcsArrays = mRHICfArrays + __NRHICFARRAYS__;
230  mFttArrays = mFcsArrays + __NFCSARRAYS__;
231  mFstArrays = mFttArrays + __NFTTARRAYS__;
232  mFwdTrackArrays = mFstArrays + __NFSTARRAYS__;
233  mTofArrays = mFwdTrackArrays + __NFWDTRACKARRAYS__;
234  mBTofArrays = mTofArrays + __NTOFARRAYS__;
235  mETofArrays = mBTofArrays + __NBTOFARRAYS__;
236  mEpdArrays = mETofArrays + __NETOFARRAYS__;
237  mMtdArrays = mEpdArrays + __NEPDARRAYS__;
238  mFgtArrays = mMtdArrays + __NMTDARRAYS__;
239  mEztArrays = mFgtArrays + __NFGTARRAYS__;
240 }
241 
243 {
244  const int ezIndex=__NARRAYS__+
245 #ifndef __NO_STRANGE_MUDST__
246  __NSTRANGEARRAYS__+
247 #endif
248  __NMCARRAYS__+
249  __NEMCARRAYS__+
250  __NPMDARRAYS__+
251  __NFMSARRAYS__+
252  __NRHICFARRAYS__+
253  __NFCSARRAYS__+
254  __NFTTARRAYS__+
255  __NFSTARRAYS__+
256  __NFWDTRACKARRAYS__+
257  __NTOFARRAYS__+
258  __NBTOFARRAYS__+
259  __NETOFARRAYS__+ //jdb
260  __NEPDARRAYS__+ // MALisa
261  __NMTDARRAYS__+__NFGTARRAYS__;
262 
263  for ( int i=0; i<ezIndex; i++) {
264  mAArrays[i]->Clear("C");
266  }
267  // ezTree classes need Delete, because of TArrayS
268  for ( int i=ezIndex; i<__NALLARRAYS__; i++) {
269  mAArrays[i]->Delete();
271  }
272 }
273 
275 {
276  memset(mAArrays,0,sizeof(void*)*__NALLARRAYS__);
277  memset(mStatusArrays,(char)1,sizeof(mStatusArrays) ); //default all ON
278  // ezt arrays switched off
279  memset(&mStatusArrays[__NARRAYS__+
280 #ifndef __NO_STRANGE_MUDST__
281  __NSTRANGEARRAYS__+
282 #endif
283  __NMCARRAYS__+
284  __NEMCARRAYS__+
285  __NPMDARRAYS__+
286  __NFMSARRAYS__+
287  __NRHICFARRAYS__+
288  __NFCSARRAYS__+
289  __NFTTARRAYS__+
290  __NFSTARRAYS__+
291  __NFWDTRACKARRAYS__+
292  __NTOFARRAYS__+
293  __NBTOFARRAYS__+
294  __NETOFARRAYS__+ //jdb
295  __NEPDARRAYS__+ // MALisa
296  __NMTDARRAYS__+
297  __NFGTARRAYS__],(char)0,__NEZTARRAYS__);
298 
299 }
300 //-----------------------------------------------------------------------
301 //-----------------------------------------------------------------------
302 //-----------------------------------------------------------------------
334 void StMuDstMaker::SetStatus(const char *arrType,int status)
335 {
336 #ifndef __NO_STRANGE_MUDST__
337  static const char *specNames[]={"MuEventAll","StrangeAll","MCAll","EmcAll","PmdAll","FMSAll","RHICfAll","FcsAll","FttAll","FstAll","TofAll","BTofAll","ETofAll","EpdAll","MTDAll","FgtAll","EztAll",0};
338 #else
339  static const char *specNames[]={"MuEventAll", "MCAll","EmcAll","PmdAll","FMSAll","RHICfAll","FcsAll","FttAll","FstAll","TofAll","BTofAll","ETofAll","EpdAll","MTDAll","FgtAll","EztAll",0};
340 #endif
341  static const int specIndex[]={
342  0, __NARRAYS__,
343  #ifndef __NO_STRANGE_MUDST__
344  __NSTRANGEARRAYS__,
345  #endif
346  __NMCARRAYS__,__NEMCARRAYS__,__NPMDARRAYS__,__NFMSARRAYS__,__NRHICFARRAYS__,__NFCSARRAYS__,__NFTTARRAYS__,__NFSTARRAYS__,__NFWDTRACKARRAYS__,__NTOFARRAYS__,__NBTOFARRAYS__,__NETOFARRAYS__,__NEPDARRAYS__,__NMTDARRAYS__,__NFGTARRAYS__,__NEZTARRAYS__,-1};
347 
348  // jdb fixed with new implementation,
349  // this method was broken for several years
350  bool found = false;
351  int species_start = 0;
352  int num = 0;
353  for (int i=0;specNames[i];i++) {
354  species_start += specIndex[i];
355  num = specIndex[i+1];
356  int cm = strcmp(arrType,specNames[i]);
357  if ( 0 == cm ) {
358  found = true;
359  break; // matches
360  }
361  }
362 
363  // found a species name
364  if ( found ){
365  char *sta=mStatusArrays+species_start;
366  memset(sta,status,num);
367  LOG_INFO << "StMuDstMaker::SetStatus \"" << arrType << "\" to status=" << status << " array indices (" << species_start << " -> " << species_start + num << ")" << endm;
368  if (mIoMode==ioRead)
369  setBranchAddresses(mChain);
370  return;
371  }
372 
373  if ( strncmp(arrType,"Strange",7) != 0 && strncmp(arrType,"St",2)==0) arrType+=2; //Ignore first "St"
374 
375  TRegexp re(arrType,1);
376  for (int i=0;i<__NALLARRAYS__;i++) {
377  Ssiz_t len;
378  if (re.Index(StMuArrays::arrayNames[i],&len) < 0) continue;
379  LOG_INFO << "StMuDstMaker::SetStatus " << status << " to " << StMuArrays::arrayNames[i] << endm;
380  mStatusArrays[i]=status;
381  }
382  if (mIoMode==ioRead)
383  setBranchAddresses(mChain);
384 }
385 //-----------------------------------------------------------------------
386 //-----------------------------------------------------------------------
387 //-----------------------------------------------------------------------
388 StMuDstMaker::StMuDstMaker(int mode, int nameMode, const char* dirName, const char* fileName, const char* filter, int maxFiles, const char* name) :
389  StIOInterFace(name),
390  mStEvent(0), mStMuDst(0),
391 #ifndef __NO_STRANGE_MUDST__
392  mStStrangeMuDstMaker(0),
393 #endif
394  mIOMaker(0), mTreeMaker(0),
395  mIoMode(mode), mIoNameMode(nameMode),
396  mDirName(dirName), mFileName(fileName), mFilter(filter),
397  mMaxFiles(maxFiles), mEventList(0),
398  mTrackType(256), mReadTracks(1),
399  mReadV0s(1), mReadXis(1), mReadKinks(1), mFinish(0),
400  mTrackFilter(0), mL3TrackFilter(0), mCurrentFile(0),
401  mSplit(99), mCompression(9), mBufferSize(65536*4),
402  mProbabilityPidAlgorithm(0), mEmcCollectionArray(0), mEmcCollection(0),
403  mFmsCollection(0), mRHICfCollection(0), mFcsCollection(0), mFttCollection(0), mFstCollection(0), mFwdTrackCollection(0), mPmdCollectionArray(0), mPmdCollection(0)
404 {
405  assignArrays();
406  streamerOff();
407  zeroArrays();
408  createArrays();
409  if (mIoMode==ioRead) openRead();
410  if (mIoMode==ioWrite) mProbabilityPidAlgorithm = new StuProbabilityPidAlgorithm();
411 
413 
414  mEventCounter=0;
415  mStMuDst = new StMuDst();
416  mEmcUtil = new StMuEmcUtil();
417  mFmsUtil = new StMuFmsUtil();
418  mRHICfUtil = new StMuRHICfUtil();
419  mFcsUtil = new StMuFcsUtil();
420  mFttUtil = new StMuFttUtil();
421  mFstUtil = new StMuFstUtil();
422  mPmdUtil = new StMuPmdUtil();
423  mTofUtil = new StMuTofUtil();
424  mBTofUtil= new StMuBTofUtil();
425  mEpdUtil = new StMuEpdUtil(); // jdb, fix rt 3338
426  mEzTree = new StMuEzTree();
427 }
428 //-----------------------------------------------------------------------
429 //-----------------------------------------------------------------------
430 //-----------------------------------------------------------------------
431 /*StMuDstMaker::StMuDstMaker(int mode, int nameMode, const char* dirName, int year) : StIOInterFace("MuDst"),
432  mStEvent(0), mStMuDst(0), mStStrangeMuDstMaker(0),
433  mIOMaker(0), mTreeMaker(0),
434  mIoMode(mode), mIoNameMode(nameMode),
435  mDirName(dirName), mFileName(""), mFilter("."),
436  mMaxFiles(10), mEventList(0),
437  mTrackType(256), mReadTracks(1),
438  mReadV0s(1), mReadXis(1), mReadKinks(1), mFinish(0),
439  mTrackFilter(0), mL3TrackFilter(0), mCurrentFile(0),
440  mSplit(99), mCompression(9), mBufferSize(65536*4),mStTriggerYear(year),
441  mProbabilityPidAlgorithm(0), mEmcCollectionArray(0), mEmcCollection(0),
442 mFmsCollection(0), mRHICfCollection(0), mPmdCollectionArray(0), mPmdCollection(0)
443 {
444  assignArrays();
445  streamerOff();
446  zeroArrays();
447  createArrays();
448  if (mIoMode==ioRead) openRead();
449  if (mIoMode==ioWrite) mProbabilityPidAlgorithm = new StuProbabilityPidAlgorithm();
450 
451  setProbabilityPidFile();
452 
453  mEventCounter=0;
454  mStMuDst = new StMuDst();
455  mEmcUtil = new StMuEmcUtil();
456  mFmsUtil = new StMuFmsUtil();
457  mRHICfUtil = new StMuRHICfUtil();
458  mPmdUtil = new StMuPmdUtil();
459  mTofUtil = new StMuTofUtil();
460  mBTofUtil= new StMuBTofUtil(); /// dongx
461  mEzTree = new StMuEzTree();
462 }*/
463 //-----------------------------------------------------------------------
464 //-----------------------------------------------------------------------
465 //-----------------------------------------------------------------------
467  DEBUGMESSAGE1("");
468  //clear(999);
469  if (mStMuDst && mStMuDst->TestBit(kCanDelete)) SafeDelete(mStMuDst);
470  SafeDelete(mTofUtil);
471  SafeDelete(mBTofUtil);
472  SafeDelete(mEpdUtil);
473  DEBUGMESSAGE3("after arrays");
474  saveDelete(mProbabilityPidAlgorithm);
475  saveDelete(mTrackFilter);
476  saveDelete(mL3TrackFilter);
477  DEBUGMESSAGE3("after filter");
478  if (mIoMode== ioWrite ) closeWrite();
479  if (mIoMode== ioRead ) closeRead();
480  DEBUGMESSAGE3("after close");
481  saveDelete(mChain);
482 //VP saveDelete(mTTree);
483  SafeDelete(mEmcCollectionArray);
484  SafeDelete(mPmdCollectionArray);
485  DEBUGMESSAGE3("out");
486 }
487 //-----------------------------------------------------------------------
488 //-----------------------------------------------------------------------
489 //-----------------------------------------------------------------------
495  StMuEvent::Class()->IgnoreTObjectStreamer();
496  StMuL3EventSummary::Class()->IgnoreTObjectStreamer();
497 #ifndef __NO_STRANGE_MUDST__
498  StStrangeMuDst::Class()->IgnoreTObjectStreamer();
499  StStrangeAssoc::Class()->IgnoreTObjectStreamer();
500  StV0MuDst::Class()->IgnoreTObjectStreamer();
501  StXiMuDst::Class()->IgnoreTObjectStreamer();
502  StKinkMuDst::Class()->IgnoreTObjectStreamer();
503  StV0Mc::Class()->IgnoreTObjectStreamer();
504  StXiMc::Class()->IgnoreTObjectStreamer();
505  StKinkMc::Class()->IgnoreTObjectStreamer();
506 #endif
507  StMuMcVertex::Class()->IgnoreTObjectStreamer();
508  StMuMcTrack::Class()->IgnoreTObjectStreamer();
509  StMuTrack::Class()->IgnoreTObjectStreamer();
510  StMuPrimaryVertex::Class()->IgnoreTObjectStreamer();
511  // StDcaGeometry::Class()->IgnoreTObjectStreamer();
512  StMuPrimaryTrackCovariance::Class()->IgnoreTObjectStreamer();
513  StMuHelix::Class()->IgnoreTObjectStreamer();
514  StMuEmcHit::Class()->IgnoreTObjectStreamer();
515  StMuEmcTowerData::Class()->IgnoreTObjectStreamer();
516  StMuPmdHit::Class()->IgnoreTObjectStreamer();
517  StMuPmdCluster::Class()->IgnoreTObjectStreamer();
518  EztEventHeader::Class()->IgnoreTObjectStreamer();
519  EztTrigBlob::Class()->IgnoreTObjectStreamer();
520  EztFpdBlob::Class()->IgnoreTObjectStreamer();
521  EztEmcRawData::Class()->IgnoreTObjectStreamer();
522  StMuFgtStrip::Class()->IgnoreTObjectStreamer();
523  StMuFgtCluster::Class()->IgnoreTObjectStreamer();
524  StMuFgtStripAssociation::Class()->IgnoreTObjectStreamer();
525  StMuFgtAdc::Class()->IgnoreTObjectStreamer();
526 }
527 //-----------------------------------------------------------------------
528 //-----------------------------------------------------------------------
529 //-----------------------------------------------------------------------
532  for ( int i=0; i<__NALLARRAYS__; i++) {
533  DEBUGVALUE2(mAArrays[i]);
535  DEBUGVALUE2(mAArrays[i]);
536  }
537  mStMuDst->set(this);
538  // commented to include tof again (subhasis)
539  // mStMuDst->set(mArrays,mStrangeArrays,mEmcArrays,mPmdArrays);
540 }
541 //-----------------------------------------------------------------------
542 //-----------------------------------------------------------------------
543 //-----------------------------------------------------------------------
544 TClonesArray* StMuDstMaker::clonesArray(TClonesArray*& p, const char* type, int size, int& counter) {
545  DEBUGMESSAGE2("");
546  if (p) return p;
547  DEBUGVALUE2(type);
548  p = new TClonesArray(type, size);
549  counter=0;
550  return p;
551 }
552 //-----------------------------------------------------------------------
553 //-----------------------------------------------------------------------
554 //-----------------------------------------------------------------------
564  DEBUGMESSAGE2("");
565  mIOMaker = (StIOMaker*)GetMaker("IOMaker");
566  mTreeMaker = (StTreeMaker*)GetMaker("outputStream");
567 #ifndef __NO_STRANGE_MUDST__
568  mStStrangeMuDstMaker = (StStrangeMuDstMaker*)GetMaker("strangeMuDst");
569 #endif
570  TDataSet *muDstSet = AddObj(mStMuDst,".const");
571  if (muDstSet ) muDstSet ->SetName("muDst");
572 
573  return 0;
574 }
575 //-----------------------------------------------------------------------
576 //-----------------------------------------------------------------------
577 //-----------------------------------------------------------------------
578 void StMuDstMaker::Clear(const char *){
579  DEBUGMESSAGE2("");
580  // We need to clear the StMuFmsCluster array even when reading events, in
581  // order to clear the TRefArrays of hits and photons it holds. Otherwise
582  // the references from the previous event may be re-used.
583  const int fmsClusterIndex = __NARRAYS__ +
584 #ifndef __NO_STRANGE_MUDST__
585  __NSTRANGEARRAYS__ +
586 #endif
587  // FMS arrays follow PMD arrays. Hits = 0th, clusters = 1st, photons = 2nd
588  __NMCARRAYS__ + __NEMCARRAYS__ + __NPMDARRAYS__ + 1;
589  mAArrays[fmsClusterIndex]->Clear("C"); // "C" calls StMuFmsCluster::Clear
590  if (mIoMode==ioRead)
591  return;
592  clearArrays();
593  if(mStMuDst->event()) mStMuDst->event()->fmsTriggerDetector().clearFlag(); //Clear flag for every event
594 
595  DEBUGMESSAGE3("out");
596 }
597 //-----------------------------------------------------------------------
598 //-----------------------------------------------------------------------
599 //-----------------------------------------------------------------------
607 
608  DEBUGMESSAGE2("");
609  int returnStarCode = kStOK;
610  StTimer timer;
611  timer.start();
612  if (mIoMode == ioWrite) returnStarCode = MakeWrite();
613  else if (mIoMode == ioRead) returnStarCode = MakeRead();
614  DEBUGVALUE2(timer.elapsedTime());
615  return returnStarCode;
616 
617 
618 }
619 //-----------------------------------------------------------------------
620 //-----------------------------------------------------------------------
621 //-----------------------------------------------------------------------
622 Int_t StMuDstMaker::MakeRead(const StUKey &RunEvent)
623 {
624  return MakeRead();
625 }
626 //-----------------------------------------------------------------------
627 //-----------------------------------------------------------------------
628 //-----------------------------------------------------------------------
629 Int_t StMuDstMaker::MakeRead()
630 {
631  int returnStarCode = kStOK;
632  if (mIoMode == ioRead) {
633  try {
634  read();
635  }
636  catch(StMuExceptionEOF e) {
637  e.print();
638  returnStarCode = kStEOF;
639  }
640  catch(StMuException e) {
641  e.print();
642  returnStarCode = kStERR;
643  }
644  }
645  return returnStarCode;
646 }
647 
648 //-----------------------------------------------------------------------
649 //-----------------------------------------------------------------------
650 //-----------------------------------------------------------------------
651 Int_t StMuDstMaker::MakeWrite(){
652  int returnStarCode = kStOK;
653  if (mIoMode == ioWrite) {
654  try {
655  write();
656  }
657  catch(StMuExceptionEOF e) {
658  e.print();
659  returnStarCode = kStEOF;
660  }
661  catch(StMuException e) {
662  e.print();
663  returnStarCode = kStERR;
664  }
665  }
666  return returnStarCode;
667 }
668 //-----------------------------------------------------------------------
669 //-----------------------------------------------------------------------
670 //-----------------------------------------------------------------------
672  DEBUGMESSAGE2("");
673  mStEvent = (StEvent*) GetInputDS("StEvent");
674  if (!mStEvent) {
675  DEBUGMESSAGE2("no StEvent");
676  throw StMuExceptionNullPointer("no StEvent",__PRETTYF__);
677  }
680  if (mProbabilityPidAlgorithm) SafeDelete(mProbabilityPidAlgorithm);
681  mProbabilityPidAlgorithm = new StuProbabilityPidAlgorithm(*mStEvent);
682  StMuTrack::setProbabilityPidAlgorithm(mProbabilityPidAlgorithm);
683  StMuTrack::setProbabilityPidCentrality(uncorrectedNumberOfNegativePrimaries(*mStEvent));
684  try {
685  fillTrees(mStEvent);
686  }
687  catch(StMuException e) {
688  e.print();
689  throw e;
690  }
691 }
692 //-----------------------------------------------------------------------
693 //-----------------------------------------------------------------------
694 //-----------------------------------------------------------------------
695 void StMuDstMaker::write(){
696  DEBUGMESSAGE2("");
697  try {
698  fill();
699  }
700  catch (StMuException e) {
701  e.print();
702  return;
703  }
704 
705  string ioMakerFileName;
706  string theFileName("/dev/null");
707  DEBUGVALUE2(mIoNameMode);
708  switch (mIoNameMode) {
709  case ioFix:
710  DEBUGMESSAGE2("===> ioFix\n");
711  theFileName = buildFileName( mDirName+"/", basename(mFileName),".MuDst.root");
712  break;
713  case ioIOMaker:
714  DEBUGMESSAGE2("===> ioIOMaker\n");
715  ioMakerFileName = string(mIOMaker->GetFile());
716  DEBUGVALUE2(ioMakerFileName);
717  theFileName = buildFileName( mDirName+"/", basename(ioMakerFileName),".MuDst.root");
718  break;
719  case ioTreeMaker:
720  // ioMakerFileName = mTreeMaker->GetTree()->GetBaseName();
721  ioMakerFileName = mTreeMaker->GetTree()->GetBaseName();
722  theFileName = buildFileName(dirname(ioMakerFileName),basename(ioMakerFileName),".MuDst.root");
723  break;
724  default:
725  DEBUGMESSAGE("do not know where to get the filename from");
726  }
727 
728  DEBUGVALUE2(theFileName.c_str());
729 
730  if (theFileName != mCurrentFileName) {
731  closeWrite();
732  openWrite(theFileName);
733  mCurrentFileName = theFileName;
734  }
735 
736  DEBUGMESSAGE2("now fill tree");
737  mTTree->Fill(); THack::IsTreeWritable(mTTree);
738  DEBUGMESSAGE2("tree filled");
739 
740  return;
741 }
742 //-----------------------------------------------------------------------
743 //-----------------------------------------------------------------------
744 //-----------------------------------------------------------------------
746  DEBUGMESSAGE2("");
747  if (mFinish) {
748  for ( int i=0; i<10; i++) {
749  LOG_INFO << "why are you calling the Finish() again ???????" << endl;
750  LOG_INFO << "are you the stupid chain destructor ???????????" << endl;
751  }
752  }
753  else {
754  if (mIoMode== ioWrite ) closeWrite();
755  if (mIoMode== ioRead ) closeRead();
756  mFinish = true;
757  }
758  DEBUGMESSAGE3("out");
759  return 0;
760 }
761 //-----------------------------------------------------------------------
762 const char* StMuDstMaker::GetFile() const {
763  if (mIoMode== ioWrite ) return mCurrentFileName.c_str();
764  if (mIoMode== ioRead && mChain && mChain->GetFile())
765  return mChain->GetFile()->GetName();
766  return 0;
767 }
768 //-----------------------------------------------------------------------
769 //-----------------------------------------------------------------------
770 void StMuDstMaker::setBranchAddresses() {
771  setBranchAddresses(mChain);
772 }
773 
774 void StMuDstMaker::setBranchAddresses(TChain* chain) {
775  // all stuff
776  if (!chain) return;
777  chain->SetBranchStatus("*",0);
778  TString ts;
779  Int_t emc_oldformat=0;
780  Int_t pmd_oldformat=0;
781  chain->BranchRef(); // Activate autoloading of TRef-referenced objects
782  for ( int i=0; i<__NALLARRAYS__; i++) {
783  if (mStatusArrays[i]==0) continue;
784  const char *bname=StMuArrays::arrayNames[i];
785  TBranch *tb = chain->GetBranch(bname);
786  if(!tb) {
787 #ifndef __NO_STRANGE_MUDST__
788  if (i >= __NARRAYS__+__NSTRANGEARRAYS__+__NMCARRAYS__ &&
789  i < __NARRAYS__+__NSTRANGEARRAYS__+__NMCARRAYS__+__NEMCARRAYS__) {
790  emc_oldformat=1;
791  continue;
792  }
793 
794  if (i >= __NARRAYS__+__NSTRANGEARRAYS__+__NMCARRAYS__+__NEMCARRAYS__ &&
795  i < __NARRAYS__+__NSTRANGEARRAYS__+__NMCARRAYS__+__NEMCARRAYS__+__NPMDARRAYS__) {
796  pmd_oldformat=1;
797  continue;
798  }
799 #else
800  if (i >= __NARRAYS__+__NMCARRAYS__ &&
801  i < __NARRAYS__+__NMCARRAYS__+__NEMCARRAYS__) {
802  emc_oldformat=1;
803  continue;
804  }
805 
806  if (i >= __NARRAYS__+__NMCARRAYS__+__NEMCARRAYS__ &&
807  i < __NARRAYS__+__NMCARRAYS__+__NEMCARRAYS__+__NPMDARRAYS__) {
808  pmd_oldformat=1;
809  continue;
810  }
811 #endif
812  Warning("setBranchAddresses","Branch name %s does not exist",bname);
813  continue;
814  }
815  ts = bname; ts +="*";
816  chain->SetBranchStatus (ts,1);
817  if (strstr("MuEvent",bname) && mChain->GetBranch("MuEvent.mQA.fX")) {
818  // Need to manually switch off Q-vector branches to avoid root warnings
819  // Note: the Q-vectors are only present in SL07b
820  mChain->SetBranchStatus("MuEvent.mQA*",0);
821  mChain->SetBranchStatus("MuEvent.mQB*",0);
822  mChain->SetBranchStatus("MuEvent.mQNegEastA*",0);
823  mChain->SetBranchStatus("MuEvent.mQNegEastB*",0);
824  mChain->SetBranchStatus("MuEvent.mQPosEastA*",0);
825  mChain->SetBranchStatus("MuEvent.mQPosEastB*",0);
826  mChain->SetBranchStatus("MuEvent.mQNegWestA*",0);
827  mChain->SetBranchStatus("MuEvent.mQNegWestB*",0);
828  mChain->SetBranchStatus("MuEvent.mQPosWestA*",0);
829  mChain->SetBranchStatus("MuEvent.mQPosWestB*",0);
830  }
831  chain->SetBranchAddress(bname,mAArrays+i);
832  assert(tb->GetAddress() == (char*)(mAArrays+i));
833  }
834  if (emc_oldformat) {
835  TBranch *branch=chain->GetBranch("EmcCollection");
836  if (branch) {
837  Warning("setBranchAddresses","Using backward compatibility mode for EMC");
838  if (!mEmcCollectionArray) {
839  mEmcCollectionArray=new TClonesArray("StMuEmcCollection",1);
840  }
841  chain->SetBranchStatus("EmcCollection*",1);
842  chain->SetBranchAddress("EmcCollection",&mEmcCollectionArray);
843  StMuEmcHit::Class()->IgnoreTObjectStreamer(0);
844  mStMuDst->set(this);
845  }
846  }
847  else if (!mEmcCollection) {
848  mEmcCollection=new StMuEmcCollection();
850  mStMuDst->set(this);
851  }
852 
853  if (!mFmsCollection) {
854  mFmsCollection=new StMuFmsCollection();
855  connectFmsCollection();
856  mStMuDst->set(this);
857  }
858 
859  if (!mRHICfCollection) {
860  mRHICfCollection=new StMuRHICfCollection();
861  connectRHICfCollection();
862  mStMuDst->set(this);
863  }
864 
865  if (!mFcsCollection) {
866  mFcsCollection=new StMuFcsCollection();
867  connectFcsCollection();
868  mStMuDst->set(this);
869  }
870 
871  if (!mFttCollection) {
872  mFttCollection=new StMuFttCollection();
873  connectFttCollection();
874  mStMuDst->set(this);
875  }
876 
877  if (!mFstCollection) {
878  mFstCollection=new StMuFstCollection();
879  connectFstCollection();
880  mStMuDst->set(this);
881  }
882 
883  if (!mFwdTrackCollection) {
884  mFwdTrackCollection=new StMuFwdTrackCollection();
885  connectFwdTrackCollection();
886  mStMuDst->set(this);
887  }
888 
889 
890  if (pmd_oldformat) {
891  TBranch *branch=chain->GetBranch("PmdCollection");
892  if (branch) {
893  Warning("setBranchAddresses","Using backward compatibility mode for PMD");
894  if (!mPmdCollectionArray) {
895  mPmdCollectionArray=new TClonesArray("StMuPmdCollection",1);
896  }
897  chain->SetBranchStatus("PmdCollection*",1);
898  chain->SetBranchAddress("PmdCollection",&mPmdCollectionArray);
899  StMuPmdCluster::Class()->IgnoreTObjectStreamer(0);
900  mStMuDst->set(this);
901  }
902  }
903  else if (!mPmdCollection) {
904  mPmdCollection=new StMuPmdCollection();
905  connectPmdCollection();
906  mStMuDst->set(this);
907  }
908  mTTree = mChain->GetTree();
909 }
910 //-----------------------------------------------------------------------
911 //-----------------------------------------------------------------------
912 //-----------------------------------------------------------------------
913 int StMuDstMaker::openRead() {
914  DEBUGVALUE2(mDirName.c_str());
915  DEBUGVALUE2(mFileName.c_str());
916  DEBUGVALUE2(mFilter.c_str());
917 
918  StMuChainMaker chainMaker("MuDst");
919  mChain = chainMaker.make(mDirName, mFileName, mFilter, mMaxFiles);
920 
921  if (mChain){
922  DEBUGVALUE3(mChain);
923  setBranchAddresses(mChain);
924 
925  mStMuDst->set(this);
926  }
927 
928  return 0;
929 }
930 //-----------------------------------------------------------------------
931 //-----------------------------------------------------------------------
932 //-----------------------------------------------------------------------
933 void StMuDstMaker::read(){
934  if (!mChain){
935  DEBUGMESSAGE2("ATTENTION: No StMuChain ... results won't be exciting (nothing to do)");
936  throw StMuExceptionNullPointer("No input files",__PRETTYF__);
937  return;
938  }
939 
940  DEBUGMESSAGE2("");
941  if (mChain->GetCurrentFile()) {
942  DEBUGVALUE2(mChain->GetCurrentFile()->GetName());
943  }
944 
945  if ( !mEventList ) {
946  int bytes = mChain->GetEntry(mEventCounter++);
947  while (bytes<=0 ) {
948  DEBUGVALUE3(mEventCounter);
949  if ( mEventCounter >= mChain->GetEntriesFast() ) throw StMuExceptionEOF("end of input",__PRETTYF__);
950  bytes = mChain->GetEntry(mEventCounter++);
951  DEBUGVALUE3(bytes);
952  }
953  }
954  else {
955  int bytes = mChain->GetEntry( mEventList->GetEntry( mEventCounter++ ) );
956  while ( bytes<=0 ) {
957  DEBUGVALUE3(mEventCounter);
958  if ( mEventCounter >= mEventList->GetN() ) throw StMuExceptionEOF("end of event list",__PRETTYF__);
959  bytes = mChain->GetEntry( mEventList->GetEntry( mEventCounter++ ) );
960  DEBUGVALUE3(bytes);
961  }
962  }
963  if (GetDebug()>1) printArrays();
964  mStMuDst->set(this);
965  fillHddr();
966  mStMuDst->setVertexIndex(0);
967  mStMuDst->collectVertexTracks(); // Make temp list of tracks for current prim vtx
968  return;
969 }
970 //-----------------------------------------------------------------------
971 //-----------------------------------------------------------------------
972 //-----------------------------------------------------------------------
973 void StMuDstMaker::closeRead(){
974  DEBUGMESSAGE2("");
975  if (mChain) mChain->Delete();
976  mChain = 0;
977  }
978 //-----------------------------------------------------------------------
979 //-----------------------------------------------------------------------
980 //-----------------------------------------------------------------------
981 void StMuDstMaker::openWrite(string fileName) {
982  DEBUGVALUE2(fileName.c_str());
983  // creat a Picoevent and and output file
984  DEBUGMESSAGE2("now create file");
985  mCurrentFile = new TFile(fileName.c_str(),"RECREATE","StMuDst");
986 
987  if (mCurrentFile->IsZombie()) throw StMuExceptionNullPointer("no file openend",__PRETTYF__);
988 
989  mCurrentFile->SetCompressionLevel(mCompression);
990 
991  // Create a ROOT Tree and one superbranch
992  DEBUGMESSAGE2("now create trees and branches");
993 
994  int bufsize = mBufferSize;
995  if (mSplit) bufsize /= 4;
996  // all stuff
997  mTTree = new TTree("MuDst", "StMuDst",mSplit);
998  mTTree->BranchRef(); // Activate autoloading of TRef-referenced objects
999 #if ROOT_VERSION_CODE < ROOT_VERSION(5,26,0)
1000  Long64_t MAXLONG=100000000000LL; // 100 GB
1001  LOG_INFO << "Tree size MAX will be " << (float) MAXLONG/1000/1000/1000 << " GB " << endm;
1002  mTTree->SetMaxTreeSize(MAXLONG); // limited to 1.9 GB - set to maximum
1003 #endif
1004  // mTTree->SetAutoSave(1000000); // autosave when 1 Mbyte written
1005  DEBUGMESSAGE2("all arrays");
1006  for ( int i=0; i<__NALLARRAYS__; i++) {
1007  if (mStatusArrays[i]==0) continue;
1008  mTTree->Branch(StMuArrays::arrayNames[i],&mAArrays[i], bufsize, mSplit);
1009  }
1010  mCurrentFileName = fileName;
1011 }
1012 //-----------------------------------------------------------------------
1013 //-----------------------------------------------------------------------
1014 //-----------------------------------------------------------------------
1015 void StMuDstMaker::closeWrite(){
1016  DEBUGMESSAGE(__PRETTYF__);
1017  if (mTTree && mCurrentFile) {
1018  LOG_INFO << " ##### " << __PRETTYF__ << endm;
1019  LOG_INFO << " ##### File=" << mCurrentFile->GetName() << " ";
1020  LOG_INFO << " NumberOfEvents= " << mTTree->GetEntries() << " ";
1021  LOG_INFO << " ##### " << endm;
1022  }
1023  //if (mTTree) mTTree->AutoSave();
1024 
1025  if (mCurrentFile) {
1026  mCurrentFile->Write();
1027  mCurrentFile->Close();
1028  }
1029  mTTree = 0;
1030  mCurrentFile = 0;
1031 }
1032 //-----------------------------------------------------------------------
1033 //-----------------------------------------------------------------------
1034 //-----------------------------------------------------------------------
1035 void StMuDstMaker::fillTrees(StEvent* ev, StMuCut* cut){
1036 
1037  DEBUGMESSAGE2("");
1038  try {
1039  fillMC();
1040  }
1041  catch(StMuException e) {
1042  e.print();
1043  throw e;
1044  }
1045 
1046  try {
1047  fillEvent(ev);
1048  fillL3AlgorithmInfo(ev);
1049  fillDetectorStates(ev);
1050  fillEmc(ev);
1051  fillPmd(ev);
1052  fillFms(ev);
1053  fillRHICf(ev);
1054  fillFcs(ev);
1055  fillFtt(ev);
1056  fillFst(ev);
1057  fillFwdTrack(ev);
1058  fillTof(ev);
1059  fillBTof(ev);
1060  fillETof(ev);
1061  fillEpd(ev); // MALisa
1062  fillMtd(ev);
1063  fillFgt(ev);
1064  fillEzt(ev);
1065  }
1066  catch(StMuException e) {
1067  e.print();
1068  throw e;
1069  }
1070 
1071  try {
1072  fillVertices(ev);
1073  }
1074  catch(StMuException e) {
1075  e.print();
1076  throw e;
1077  }
1078 
1079  try {
1080  fillpp2pp(ev);
1081  }
1082  catch(StMuException e) {
1083  e.print();
1084  throw e;
1085  }
1086 
1087  try {
1088  fillTracks(ev,mTrackFilter);
1089  }
1090  catch(StMuException e) {
1091  e.print();
1092  throw e;
1093  }
1094 
1095  try {
1096  fillL3Tracks(ev, mL3TrackFilter);
1097  }
1098  catch(StMuException e) {
1099  e.print();
1100  throw e;
1101  }
1102 
1103 #ifndef __NO_STRANGE_MUDST__
1104  if (mStStrangeMuDstMaker) {
1105  try {
1106  fillStrange(mStStrangeMuDstMaker);
1107  }
1108  catch(StMuException e) {
1109  e.print();
1110  throw e;
1111  }
1112  }
1113 #endif
1114 
1115  //catch(StMuException e) {
1116  // e.print();
1117  // throw e;
1118  //}
1119 
1120  mStMuDst->set(this);
1121  mStMuDst->fixTofTrackIndices();
1122  mStMuDst->fixETofTrackIndices();
1123  mStMuDst->fixMtdTrackIndices();
1124  mStMuDst->fixTrackIndicesG(mStMuDst->numberOfPrimaryVertices());
1125 }
1126 
1127 
1128 
1129 //-----------------------------------------------------------------------
1130 //-----------------------------------------------------------------------
1131 //-----------------------------------------------------------------------
1132 void StMuDstMaker::fillEvent(StEvent* ev, StMuCut* cut) {
1133  DEBUGMESSAGE2("");
1134  StMuEvent *typeOfEvent=0;
1135  if (!ev) throw StMuExceptionNullPointer("no StEvent",__PRETTYF__);
1136  StTimer timer;
1137  timer.start();
1138  if (!cut || cut->pass(ev)) {
1139  DEBUGMESSAGE3("");
1140  addType(mArrays[muEvent],ev,typeOfEvent);
1141  }
1142  timer.stop();
1143  DEBUGVALUE2(timer.elapsedTime());
1144 }
1145 //-----------------------------------------------------------------------
1146 //-----------------------------------------------------------------------
1147 //-----------------------------------------------------------------------
1148 void StMuDstMaker::fillEmc(StEvent* ev) {
1149  DEBUGMESSAGE2("");
1150  StEmcCollection* emccol=(StEmcCollection*)ev->emcCollection();
1151  if (!emccol) return; //throw StMuExceptionNullPointer("no StEmcCollection",__PRETTYF__);
1152  StTimer timer;
1153  timer.start();
1154 
1155  TClonesArray *tca = mEmcArrays[muEmcTow];
1156  new((*tca)[0]) StMuEmcTowerData();
1157  if (!mEmcCollection) {
1158  mEmcCollection=new StMuEmcCollection();
1160  mStMuDst->set(this);
1161  }
1162  mEmcUtil->fillMuEmc(mEmcCollection,emccol);
1163 
1164  timer.stop();
1165  DEBUGVALUE2(timer.elapsedTime());
1166 }
1167 //-----------------------------------------------------------------------
1168 //-----------------------------------------------------------------------
1169 void StMuDstMaker::fillFms(StEvent* ev) {
1170  DEBUGMESSAGE2("");
1171  StFmsCollection* fmscol=(StFmsCollection*)ev->fmsCollection();
1172  if (!fmscol) return; //throw StMuExceptionNullPointer("no StFmsCollection",__PRETTYF__);
1173  StTimer timer;
1174  timer.start();
1175 
1176  if (!mFmsCollection) {
1177  mFmsCollection=new StMuFmsCollection();
1178  connectFmsCollection();
1179  mStMuDst->set(this);
1180  }
1181  LOG_DEBUG << "StMuDSTMaker filling StMuFmsCollection from StEvent" << endm;
1182  mFmsUtil->fillMuFms(mFmsCollection,fmscol);
1183 
1184  timer.stop();
1185  DEBUGVALUE2(timer.elapsedTime());
1186 }
1187 //-----------------------------------------------------------------------
1188 //-----------------------------------------------------------------------
1189 void StMuDstMaker::fillRHICf(StEvent* ev) {
1190  DEBUGMESSAGE2("");
1191  StRHICfCollection* rhicfcol=(StRHICfCollection*)ev->rhicfCollection();
1192  if (!rhicfcol) return; //throw StMuExceptionNullPointer("no StRHICfCollection",__PRETTYF__);
1193  StTimer timer;
1194  timer.start();
1195 
1196  if (!mRHICfCollection) {
1197  LOG_INFO << "no find a RHICfCollection !!! " << endm;
1198  mRHICfCollection=new StMuRHICfCollection();
1199  connectRHICfCollection();
1200  mStMuDst->set(this);
1201  }
1202  LOG_DEBUG << "StMuDSTMaker filling StMuRHICfCollection from StEvent" << endm;
1203  mRHICfUtil->fillMuRHICf(mRHICfCollection,rhicfcol);
1204 
1205  timer.stop();
1206  DEBUGVALUE2(timer.elapsedTime());
1207 }
1208 //-----------------------------------------------------------------------
1209 //-----------------------------------------------------------------------
1210 void StMuDstMaker::fillFcs(StEvent* ev) {
1211  DEBUGMESSAGE2("");
1212  StFcsCollection* fcscol=(StFcsCollection*)ev->fcsCollection();
1213  if (!fcscol) return; //throw StMuExceptionNullPointer("no StFcsCollection",__PRETTYF__);
1214  StTimer timer;
1215  timer.start();
1216 
1217 
1218  if (!mFcsCollection) {
1219  mFcsCollection=new StMuFcsCollection();
1220  connectFcsCollection();
1221  mStMuDst->set(this);
1222  }
1223 
1224  mFcsUtil->fillMuFcs(mFcsCollection,fcscol);
1225 
1226  timer.stop();
1227  DEBUGVALUE2(timer.elapsedTime());
1228 }
1229 //-----------------------------------------------------------------------
1230 //-----------------------------------------------------------------------
1231 void StMuDstMaker::fillFtt(StEvent* ev) {
1232  DEBUGMESSAGE2("");
1233  StFttCollection* fttcol=(StFttCollection*)ev->fttCollection();
1234  if (!fttcol) return; //throw StMuExceptionNullPointer("no StFttCollection",__PRETTYF__);
1235  StTimer timer;
1236  timer.start();
1237 
1238  if (!mFttCollection) {
1239  mFttCollection=new StMuFttCollection();
1240  connectFttCollection();
1241  mStMuDst->set(this);
1242  }
1243  mFttUtil->fillMuFtt(mFttCollection,fttcol);
1244 
1245  timer.stop();
1246  DEBUGVALUE2(timer.elapsedTime());
1247 }
1248 //-----------------------------------------------------------------------
1249 //-----------------------------------------------------------------------
1250 void StMuDstMaker::fillFst(StEvent* ev) {
1251  DEBUGMESSAGE2("");
1252  StFstHitCollection* fstcol=(StFstHitCollection*)ev->fstHitCollection();
1253  if (!fstcol) return; //throw StMuExceptionNullPointer("no StFstHitCollection",__PRETTYF__);
1254  StTimer timer;
1255  timer.start();
1256 
1257  if (!mFstCollection) {
1258  mFstCollection=new StMuFstCollection();
1259  connectFstCollection();
1260  mStMuDst->set(this);
1261  }
1262 
1263  //raw hit data input
1264  StFstEvtCollection *fstevtcol = 0;
1265  if (IAttr("fstMuRawHit")){//True for storing FST Raw hits
1266  fstevtcol=(StFstEvtCollection*)ev->fstEvtCollection();
1267  }
1268 
1269  mFstUtil->fillMuFst(mFstCollection,fstcol,fstevtcol);
1270 
1271  timer.stop();
1272  DEBUGVALUE2(timer.elapsedTime());
1273 }
1274 //-----------------------------------------------------------------------
1275 //-----------------------------------------------------------------------
1276 void StMuDstMaker::fillFwdTrack(StEvent* ev) {
1277  DEBUGMESSAGE2("");
1278  LOG_INFO << "StMuDstMaker::fillFwdTrack(StEvent* ev)" << endm;
1279  StFwdTrackCollection* fwdcol=(StFwdTrackCollection*)ev->fwdTrackCollection();
1280  if (!fwdcol) return; //throw StMuExceptionNullPointer("no StFstHitCollection",__PRETTYF__);
1281  StTimer timer;
1282  timer.start();
1283 
1284  LOG_INFO << "StMuDSTMaker filling StMuFwdTrackCollection from StEvent" << endm;
1285 
1286  if (!mFwdTrackCollection) {
1287  mFwdTrackCollection=new StMuFwdTrackCollection();
1288  LOG_INFO << "Connecting StMuFwdTrackCollection" << endm;
1289  connectFwdTrackCollection();
1290  mStMuDst->set(this);
1291  }
1292 
1293  mFwdTrackUtil->fillMuFwdTrack(mFwdTrackCollection,fwdcol,mFcsUtil);
1294 
1295  timer.stop();
1296  DEBUGVALUE2(timer.elapsedTime());
1297 }
1298 //-----------------------------------------------------------------------
1299 //-----------------------------------------------------------------------
1300 void StMuDstMaker::fillPmd(StEvent* ev) {
1301  DEBUGMESSAGE2("");
1302  StPhmdCollection* phmdColl=(StPhmdCollection*)ev->phmdCollection();
1303  if (!phmdColl) return; //throw StMuExceptionNullPointer("no StPhmdCollection",__PRETTYF__);
1304  StTimer timer;
1305  timer.start();
1306 
1307  if (!mPmdCollection) {
1308  mPmdCollection=new StMuPmdCollection();
1309  connectPmdCollection();
1310  mStMuDst->set(this);
1311  }
1312  mPmdUtil->fillMuPmd(phmdColl,mPmdCollection);
1313 
1314  timer.stop();
1315  DEBUGVALUE2(timer.elapsedTime());
1316 }
1317 
1318 //-----------------------------------------------------------------------
1319 //-----------------------------------------------------------------------
1320 
1321 void StMuDstMaker::fillTof(StEvent* ev) {
1322  DEBUGMESSAGE2("");
1323  StTofCollection *tofcol = ev->tofCollection();
1324  // run 5 - dongx
1325  if( !ev || !tofcol || (!tofcol->dataPresent()&&!tofcol->rawdataPresent()) )
1326  return; //throw StMuExceptionNullPointer("no StTofDataCollection",__PRETTYF__);
1327  StTimer timer;
1328  timer.start();
1329 
1330  // fill tofHit
1331  StMuTofHitCollection muTofHitColl;
1332  mTofUtil->fillMuTofHit(&muTofHitColl, tofcol);
1333  for(size_t i=0; i < muTofHitColl.size(); i++) {
1334  StMuTofHit* tofMuHit = (StMuTofHit *)muTofHitColl.getHit(i);
1335  addType( mTofArrays[muTofHit], *tofMuHit );
1336  }
1337 
1338  // fill tofData
1339  StSPtrVecTofData &tofData = tofcol->tofData();
1340  for(size_t i=0; i < tofData.size(); i++) {
1341  addType( mTofArrays[muTofData], *tofData[i] );
1342  }
1343 
1344  // run 5 - dongx
1345  // fill tofRawData
1346  StSPtrVecTofRawData &tofRawData = tofcol->tofRawData();
1347  for(size_t i=0; i < tofRawData.size(); i++) {
1348  addType( mTofArrays[muTofRawData], *tofRawData[i] );
1349  }
1350 
1351  timer.stop();
1352  DEBUGVALUE2(timer.elapsedTime());
1353 }
1354 //-----------------------------------------------------------------------
1355 //-----------------------------------------------------------------------
1356 //-----------------------------------------------------------------------
1359  DEBUGMESSAGE2("");
1360  StBTofCollection *btofcol = ev->btofCollection();
1361  if( !ev || !btofcol || !btofcol->rawHitsPresent() )
1362  return; //throw StMuExceptionNullPointer("no StBTofRawHitCollection",__PRETTYF__);
1363  StTimer timer;
1364  timer.start();
1365 
1366  // fill btofHit
1367  StMuBTofHitCollection muBTofHitColl;
1368  mBTofUtil->fillMuBTofHit(&muBTofHitColl, btofcol);
1369  for(size_t i=0; i < muBTofHitColl.size(); i++) {
1370  StMuBTofHit* btofMuHit = (StMuBTofHit *)muBTofHitColl.getHit(i);
1371  addType( mBTofArrays[muBTofHit], *btofMuHit );
1372  }
1373 
1374  // fill btofRawHit
1375  StSPtrVecBTofRawHit &btofRawHits = btofcol->tofRawHits();
1376  for(size_t i=0; i < btofRawHits.size(); i++) {
1377  addType( mBTofArrays[muBTofRawHit], *btofRawHits[i] );
1378  }
1379 
1380  // fill btofHeader
1381  StBTofHeader *btofHeader = btofcol->tofHeader();
1382  addType( mBTofArrays[muBTofHeader], *btofHeader);
1383 
1384  timer.stop();
1385  DEBUGVALUE2(timer.elapsedTime());
1386 }
1387 //-----------------------------------------------------------------------
1388 //-----------------------------------------------------------------------
1389 //-----------------------------------------------------------------------
1392  DEBUGMESSAGE2( "" );
1393 
1394  const StETofCollection* eTofColl = ev->etofCollection();
1395  if( !eTofColl )
1396  return;
1397 
1398  StTimer timer;
1399  timer.start();
1400 
1401  StMuETofCollection muETofColl( eTofColl );
1402 
1403  // fill etof digis from StMuETofCollection
1404  for(size_t i=0; i < (size_t) muETofColl.digisPresent(); i++ ) {
1405  StMuETofDigi* pDigi = ( StMuETofDigi* ) muETofColl.etofDigi( i );
1406  if( !pDigi ) continue;
1407  addType( mETofArrays[ muETofDigi ], *pDigi );
1408  }
1409 
1410  // fill etof hits from StMuETofCollection
1411  for(size_t i=0; i < (size_t) muETofColl.hitsPresent(); i++ ) {
1412  StMuETofHit* pHit = ( StMuETofHit* ) muETofColl.etofHit( i );
1413  if( !pHit ) continue;
1414  addType( mETofArrays[ muETofHit ], *pHit );
1415  }
1416 
1417  //fill etof header from StMuETofCollection
1418  StMuETofHeader* pHeader = muETofColl.etofHeader();
1419  if( pHeader ) {
1420  addType( mETofArrays[ muETofHeader ], *pHeader );
1421  }
1422 
1423  timer.stop();
1424  DEBUGVALUE2( timer.elapsedTime() );
1425 }
1426 // -------------------------------------------------------------------------
1427 // -------------------------------------------------------------------------
1428 void StMuDstMaker::fillEpd(StEvent* ev){ // MALisa
1429  DEBUGMESSAGE2("");
1430  StEpdCollection* epdcol = (StEpdCollection*)ev->epdCollection();
1431  if (!epdcol) return; //throw StMuExceptionNullPointer("no StEpdCollection"__PRETTYF__);
1432  StTimer timer;
1433  timer.start();
1434 
1435  // fill epdHit
1436  StMuEpdHitCollection muEpdHitColl;
1437  mEpdUtil->fillMuEpdHit(&muEpdHitColl, epdcol);
1438  for(size_t i=0; i < muEpdHitColl.size(); i++){
1439  StMuEpdHit* epdMuHit = (StMuEpdHit*)muEpdHitColl.getHit(i);
1440  addType( mEpdArrays[muEpdHit], *epdMuHit);
1441  }
1442  timer.stop();
1443  DEBUGVALUE2(timer.elapsedTime());
1444  /*
1445  if (!mMuEpdHitCollection){
1446  mMuEpdHitCollection = new StMuEpdHitCollection();
1447  mStMuDst->set(this);
1448  }
1449  LOG_DEBUG << "StMuDSTMaker filling StMuEpdHitCollection from StEvent" << endm;
1450  mMuEpdHitCollection->Fill(epdHitCollection);
1451 
1452  cout << "\n\n\n\n I just filled the muEpdHitCollection and it has " << mMuEpdHitCollection->nHits() << " hits\n\n\n\n";
1453  */
1454 }
1455 //-----------------------------------------------------------------------
1456 //-----------------------------------------------------------------------
1457 //-----------------------------------------------------------------------
1458 void StMuDstMaker::fillFgt(StEvent* ev) {
1459 
1460  const StFgtCollection* fgtCollPtr = 0;
1461  if( ev )
1462  fgtCollPtr = ev->fgtCollection();
1463 
1464  if( fgtCollPtr ){
1465  // assert existance of the arrays
1466  assert( mFgtArrays );
1467  assert( mFgtArrays[muFgtStrips] );
1468  assert( mFgtArrays[muFgtClusters] );
1469  assert( mFgtArrays[muFgtStripAssociations] );
1470 
1471  // need pointer types to enforce conversions in StMuDstMaker::addType
1472  StMuFgtCluster* clusterClassType = 0;
1473 
1474  // need a map to keep track of the index of each strip.
1475  // Key = geoId, value = idx in muFgtStrips TClonesArray
1476  std::map< Int_t, Int_t > stripGeoIdIdxMap;
1477  std::map< Int_t, Int_t >::iterator stripGeoIdIdxIter;
1478 
1479  // iterate over discs
1480  for( UShort_t discIdx = 0; discIdx < kFgtNumDiscs; ++discIdx ){
1481 
1482  // copy strips
1483  const StSPtrVecFgtStrip& stripVec = fgtCollPtr->mStripCollection[ discIdx ].getStripVec();
1484  for( const_StFgtStripIterator stripIter = stripVec.begin(); stripIter != stripVec.end(); ++stripIter ){
1485  Short_t stripType = (*stripIter)->getClusterSeedType();
1486 
1487  // save only the strips with a certain type
1488  if( stripType == kFgtSeedType1 ||
1489  stripType == kFgtSeedType2 ||
1490  stripType == kFgtSeedType3 ||
1491  stripType == kFgtSeedType4 ||
1492  stripType == kFgtSeedType5 ||
1493  stripType == kFgtSeedTypeMax ||
1494  stripType == kFgtClusterPart ||
1495  stripType == kFgtNextToDeadGuy ||
1496  stripType == kFgtClusterEndUp ||
1497  stripType == kFgtClusterEndDown ||
1498  stripType == kFgtStripShared ||
1499  stripType == kFgtClusterTooBig ||
1500  stripType == kFgtClusterSeedInSeaOfNoise ||
1501  stripType == kFgtNextToCluster ||
1502  stripType == kFgtKeepStrip){
1503 
1504  // make the strip using the "casting" constructor
1505  StMuFgtStrip fgtStrip( *(*stripIter) );
1506 
1507  // Set the range for the time bins to save.
1508  // Note: for now, use the obvious. Update the
1509  // following two lines later.`1
1510  Int_t nTBstart = 0;
1511  Int_t nTBend = nTBstart + fgtCollPtr->getNumTimeBins();
1512 
1513  // to do: comment about ensuring tb in same range as StFgtStrip class in StEvent
1514 
1515  for( Int_t tb = nTBstart; tb < nTBend; ++tb ){
1516  // make an StMuFgtAdc object
1517  StMuFgtAdc fgtAdc( (*stripIter)->getAdc( tb ), tb );
1518 
1519  Int_t adcIdx = addType( mFgtArrays[muFgtAdcs], fgtAdc );
1520 
1521  // just the first time
1522  if( tb == nTBstart ){
1523 
1524  // set the strip to know where the first ADC is and how many there are
1525  fgtStrip.setAdcInfo( adcIdx, nTBend-nTBstart );
1526  };
1527  };
1528 
1529  // add strip to the TClonesArray
1530  Int_t idx = addType( mFgtArrays[muFgtStrips], fgtStrip );
1531 
1532  // add to the map
1533  stripGeoIdIdxMap[ (*stripIter)->getGeoId() ] = idx;
1534 
1535  };
1536  };
1537 
1538  // faster to declare iterator outside of the loop
1539  stripWeightMap_t::const_iterator wIter;
1540 
1541  // copy clusters
1542  const StSPtrVecFgtHit& hitVec = fgtCollPtr->mHitCollection[ discIdx ].getHitVec();
1543  for( const_StFgtHitIterator hitIter = hitVec.begin(); hitIter != hitVec.end(); ++hitIter ){
1544  if( (*hitIter)->charge() > 0 ){
1545  // add the cluster to the array
1546  Int_t clusIdx = addType( mFgtArrays[muFgtClusters], *(*hitIter), clusterClassType );
1547 
1548  // get the map of associations
1549  const stripWeightMap_t& wMap = (*hitIter)->getStripWeightMap();
1550 
1551  // iterate over associated strips and add the associations
1552  Bool_t isFirst = 1;
1553  for( wIter = wMap.begin(); wIter != wMap.end(); ++wIter ){
1554 
1555  Int_t geoId = wIter->first->getGeoId(); // first is a StFgtStrip, which knows its own geoId
1556  Float_t weight = wIter->second; // second is the weight
1557 
1558  // determine the stripIdx
1559  stripGeoIdIdxIter = stripGeoIdIdxMap.find( geoId );
1560 
1561  // make sure the strip was saved
1562  assert( stripGeoIdIdxIter != stripGeoIdIdxMap.end() );
1563 
1564  Int_t stripIdx = stripGeoIdIdxIter->second;
1565 
1566  // make the association
1567  StMuFgtStripAssociation association( clusIdx, stripIdx, weight );
1568 
1569  // add it to the array
1570  Int_t associationIdx = addType( mFgtArrays[muFgtStripAssociations], association );
1571 
1572  if( isFirst ){
1573  // change the flag state
1574  isFirst = 0;
1575 
1576  // set the index in the cluster
1577  StMuFgtCluster* clusPtr = (StMuFgtCluster*)mFgtArrays[muFgtClusters]->UncheckedAt( clusIdx );
1578  clusPtr->setFirstStripAssociationIndex( associationIdx );
1579  };
1580  };
1581  };
1582  };
1583  };
1584  };
1585 };
1586 
1587 void StMuDstMaker::fillMtd(StEvent* ev) {
1588  DEBUGMESSAGE2("");
1589  StTimer timer;
1590  timer.start();
1591 
1592 // StMuMtdCollection *typeOfMtd=0;
1593 
1594  const StMtdCollection *mtd=ev->mtdCollection();
1595  if(!mtd) return;
1596  StMuMtdCollection mMTD(*mtd);
1597 
1598  for(size_t i=0; i < (size_t)mMTD.hitsPresent(); i++) {
1599  StMuMtdHit* mtdHit = (StMuMtdHit*)mMTD.MtdHit(i);
1600  addType( mMtdArrays[muMTDHit], *mtdHit );
1601  }
1602 
1603  for(size_t i=0; i < (size_t)mMTD.rawHitsPresent(); i++) {
1604  StMtdRawHit* mtdHit = (StMtdRawHit*)mMTD.RawMtdHit(i);
1605  addType( mMtdArrays[muMTDRawHit], *mtdHit );
1606  }
1607  StMuMtdHeader *mtdHead = mMTD.mtdHeader();
1608  if(mtdHead) addType(mMtdArrays[muMTDHeader],*mtdHead);
1609 
1610  //if (mtd) addType( mArrays[muMtd], *mtd, typeOfMtd );
1611  timer.stop();
1612  DEBUGVALUE2(timer.elapsedTime());
1613 }
1614 
1615 //-----------------------------------------------------------------------
1616 //-----------------------------------------------------------------------
1617 //-----------------------------------------------------------------------
1619  if (ev==0) return;
1620  char *eztArrayStatus=&mStatusArrays[__NARRAYS__+
1621 #ifndef __NO_STRANGE_MUDST__
1622  __NSTRANGEARRAYS__+
1623 #endif
1624  __NMCARRAYS__+__NEMCARRAYS__+__NPMDARRAYS__+__NFMSARRAYS__+__NRHICFARRAYS__+__NFCSARRAYS__+__NFTTARRAYS__+__NFSTARRAYS__+__NFWDTRACKARRAYS__+
1625  __NTOFARRAYS__+__NBTOFARRAYS__+__NETOFARRAYS__+__NEPDARRAYS__+__NMTDARRAYS__+__NFGTARRAYS__];
1626  if(eztArrayStatus[muEztHead]){
1627  EztEventHeader* header = mEzTree->copyHeader(ev);
1628  addType(mEztArrays[muEztHead], *header);
1629  }
1630 
1631  if(eztArrayStatus[muEztTrig]) {
1632  EztTrigBlob* trig = mEzTree->copyTrig(ev);
1633  if (trig)
1634  addType(mEztArrays[muEztTrig], *trig);
1635  }
1636 
1637  if(eztArrayStatus[muEztFpd]) {
1638  EztFpdBlob* fpd = mEzTree->copyFpd(ev);
1639  addType(mEztArrays[muEztFpd], *fpd);
1640  }
1641 
1642  if(eztArrayStatus[muEztETow] || eztArrayStatus[muEztESmd]) {
1643  StEmcCollection* emcCol=(StEmcCollection*)ev->emcCollection();
1644  if(emcCol==0){
1645  gMessMgr->Message("","W") << GetName()<<"::fillEzt(), missing StEmcCollection, EEMC raw data NOT saved in muDst" <<endm;
1646  } else { //........... EMC-Collection in StEvent exist
1647  StEmcRawData *eeRaw=emcCol->eemcRawData();
1648  if (eeRaw) {
1649  if(eztArrayStatus[muEztETow]) {
1650  EztEmcRawData* ETow = mEzTree->copyETow(eeRaw);
1651  addType(mEztArrays[muEztETow], *ETow);
1652  }
1653 
1654  if(eztArrayStatus[muEztESmd]) {
1655  EztEmcRawData* ESmd = mEzTree->copyESmd(eeRaw);
1656  addType(mEztArrays[muEztESmd], *ESmd);
1657  }
1658  }
1659  }
1660  }
1661 }
1662 //-----------------------------------------------------------------------
1663 //-----------------------------------------------------------------------
1664 //-----------------------------------------------------------------------
1665 void StMuDstMaker::fillL3AlgorithmInfo(StEvent* ev) {
1666  DEBUGMESSAGE2("");
1667  if ( !ev->l3Trigger() ) return;
1668  if ( !ev->l3Trigger()->l3EventSummary()) return;
1669 
1670  StTimer timer;
1671  timer.start();
1672  StL3EventSummary* l3 = ev->l3Trigger()->l3EventSummary();
1673  int n = l3->numberOfAlgorithms();
1674  for (int i=0; i<n; i++) {
1675  if (l3->algorithms()[i]->accept())
1676  addType( mArrays[muAccept], *l3->algorithms()[i] );
1677  else
1678  addType( mArrays[muReject], *l3->algorithms()[i] );
1679  }
1680  timer.stop();
1681  DEBUGVALUE2(timer.elapsedTime());
1682 }
1683 //-----------------------------------------------------------------------
1684 //-----------------------------------------------------------------------
1685 //-----------------------------------------------------------------------
1686 void StMuDstMaker::fillVertices(StEvent* ev) {
1687  DEBUGMESSAGE2("");
1688  StTimer timer;
1689  timer.start();
1690 
1691  StMuPrimaryVertex *typeOfVertex=0;
1692  Int_t n_vtx = ev->numberOfPrimaryVertices();
1693  DEBUGVALUE2(n_vtx);
1694  mVtxList.Clear();
1695  for (Int_t i_vtx=0; i_vtx < n_vtx; i_vtx++) {
1696  const StPrimaryVertex *vtx=ev->primaryVertex(i_vtx);
1697  addType( mArrays[muPrimaryVertex], vtx, typeOfVertex );
1698  mVtxList.AddAtAndExpand(ev->primaryVertex(i_vtx),i_vtx);
1699  }
1700  timer.stop();
1701  DEBUGVALUE2(timer.elapsedTime());
1702 }
1703 
1704 
1705 void StMuDstMaker::fillpp2pp(StEvent* ev) {
1706  DEBUGMESSAGE2("");
1707  StTimer timer;
1708  timer.start();
1709 
1710  StMuRpsCollection *typeOfRps=0;
1711 
1712  const StRpsCollection *rps=ev->rpsCollection();
1713  if (rps) addType( mArrays[mupp2pp], *rps, typeOfRps );
1714  timer.stop();
1715  DEBUGVALUE2(timer.elapsedTime());
1716 }
1717 
1718 
1719 //-----------------------------------------------------------------------
1720 //-----------------------------------------------------------------------
1721 //-----------------------------------------------------------------------
1722 void StMuDstMaker::fillTracks(StEvent* ev, StMuCut* cut) {
1723  DEBUGMESSAGE2("");
1724  StTimer timer;
1725  timer.start();
1726 
1727  StSPtrVecTrackNode& nodes= ev->trackNodes();
1728  DEBUGVALUE2(nodes.size());
1729  for (StSPtrVecTrackNodeConstIterator iter=nodes.begin(); iter!=nodes.end(); iter++) {
1730  addTrackNode(ev, *iter, cut, mArrays[muGlobal], mArrays[muPrimary], mArrays[muOther], mArrays[muCovGlobTrack], mArrays[muCovPrimTrack], false);
1731  }
1732  timer.stop();
1733  DEBUGVALUE2(timer.elapsedTime());
1734 }
1735 //-----------------------------------------------------------------------
1736 //-----------------------------------------------------------------------
1737 //-----------------------------------------------------------------------
1738 void StMuDstMaker::fillL3Tracks(StEvent* ev, StMuCut* cut) {
1739  DEBUGMESSAGE2("");
1740  if (!ev->l3Trigger()) return;
1741 
1742  StTimer timer;
1743  timer.start();
1744  StSPtrVecTrackNode& nodes= ev->l3Trigger()->trackNodes();
1745  DEBUGVALUE2(nodes.size());
1746  for (StSPtrVecTrackNodeConstIterator iter=nodes.begin(); iter!=nodes.end(); iter++) {
1747  addTrackNode(ev, *iter, cut, mArrays[muL3], 0, 0, 0, 0, true );
1748  }
1749  timer.stop();
1750  DEBUGVALUE2(timer.elapsedTime());
1751 }
1752 //-----------------------------------------------------------------------
1753 //-----------------------------------------------------------------------
1754 //-----------------------------------------------------------------------
1755 void StMuDstMaker::fillDetectorStates(StEvent* ev) {
1756  DEBUGMESSAGE2("");
1757  StTimer timer;
1758  timer.start();
1759  for (int i=0; i<StMuArrays::arraySizes[muState]; i++) {
1760  StDetectorState* state = ev->detectorState((StDetectorId) i);
1761  if (state) {
1762  addType( mArrays[muState], *state );
1763  }
1764  }
1765  timer.stop();
1766  DEBUGVALUE2(timer.elapsedTime());
1767 }
1768 //-----------------------------------------------------------------------
1769 //-----------------------------------------------------------------------
1770 //-----------------------------------------------------------------------
1771 void StMuDstMaker::addTrackNode(const StEvent* ev, const StTrackNode* node, StMuCut* cut,
1772  TClonesArray* gTCA, TClonesArray* pTCA, TClonesArray* oTCA, TClonesArray* covgTCA, TClonesArray* covpTCA, bool l3) {
1773  DEBUGMESSAGE3("");
1774  const StTrack* tr=0;
1775 
1777  int index2Global =-1;
1778  if (gTCA) {
1779  const StTrack *pr_tr = node->track(primary);
1780  const StVertex *vtx = 0;
1781  if (pr_tr)
1782  vtx = pr_tr->vertex();
1783  if (vtx==0)
1784  vtx = ev->primaryVertex();
1785 
1786  tr= node->track(global);
1787  if (tr && !tr->bad()) index2Global = addTrack(gTCA, ev, tr, vtx, cut, -1, l3, covgTCA, covpTCA);
1788  }
1790  if (pTCA) {
1791  tr = node->track(primary);
1792  if (tr && !tr->bad()) addTrack(pTCA, ev, tr, tr->vertex(), cut, index2Global, l3, covgTCA, covpTCA);
1793  }
1795  if (oTCA) {
1796  size_t nEntries = node->entries();
1797  for (size_t j=0; j<nEntries; j++) {
1798  tr = node->track(j);
1799  if (tr && !tr->bad() && (tr->type()!=global) && (tr->type()!=primary) ) {
1800  addTrack(oTCA, ev, tr, tr->vertex(), cut, index2Global, l3);
1801  }
1802  }
1803  }
1804 }
1805 //---------------------------------------------------------------------
1806 //---------------------------------------------------------------------
1807 //---------------------------------------------------------------------
1808 int StMuDstMaker::addTrack(TClonesArray* tca, const StEvent*event, const StTrack* track, const StVertex *vtx, StMuCut* cut, int index2Global, bool l3,
1809  TClonesArray* covgTCA, TClonesArray* covpTCA) {
1810  DEBUGMESSAGE3("");
1811  StRichSpectra typeOfStRichSpectra;
1812  int index = -1;
1813  int index2RichSpectra=-1;
1815  int counter = tca->GetEntries();
1816  try{
1817  if (cut && !cut->pass(track)) throw StMuExceptionBadValue("failed track cut",__PRETTYF__);
1818  // add StRichSpectra if StRichPidTraits are found
1819  // we have to do this more elegant
1820  StRichSpectra* rich = richSpectra(track);
1821  if (rich) {
1822  index2RichSpectra = addType( mArrays[muRich], *rich );
1823  }
1824  StMuTrack *muTrack = new((*tca)[counter]) StMuTrack(event, track, vtx, index2Global, index2RichSpectra, l3, &mVtxList);
1825  if (track->type() == primary) {
1826  if (covpTCA) {
1827  Int_t countCOVPTCA = covpTCA->GetEntries();
1828 #if 0
1829  const StMatrixF covMatrix = track->fitTraits().covariantMatrix();
1830  new((*covpTCA)[countCOVPTCA]) StMuPrimaryTrackCovariance(covMatrix);
1831 #else
1832  // cout << track->fitTraits().covariantMatrix() << endl;
1833  const Float_t* cov = track->fitTraits().covariance();
1834  new((*covpTCA)[countCOVPTCA]) StMuPrimaryTrackCovariance(cov);
1835 #endif
1836  muTrack->setIndex2Cov(countCOVPTCA);
1837  }
1838  }
1839  else {
1840  if (track->type() == global) {
1841  if (covgTCA) {
1842  Int_t countCOVGTCA = covgTCA->GetEntries();
1843  const StDcaGeometry *dcaGeometry = ((StGlobalTrack *)track)->dcaGeometry();
1844  if (dcaGeometry) {
1845  new((*covgTCA)[countCOVGTCA]) StDcaGeometry(*dcaGeometry);
1846  muTrack->setIndex2Cov(countCOVGTCA);
1847  }
1848  }
1849  }
1850  }
1851  index = counter;
1852  }
1853  catch (StMuException e) {
1854  IFDEBUG3(e.print());
1855  }
1856  return index;
1857 }
1858 //-----------------------------------------------------------------------
1859 //-----------------------------------------------------------------------
1860 //-----------------------------------------------------------------------
1861 StRichSpectra* StMuDstMaker::richSpectra(const StTrack* track) {
1862  DEBUGMESSAGE3("");
1863  const StPtrVecTrackPidTraits& traits = track->pidTraits(kRichId);
1864  for (StPtrVecTrackPidTraitsConstIterator traitIter=traits.begin();traitIter!=traits.end();++traitIter) {
1865  StRichPidTraits* pid = dynamic_cast<StRichPidTraits*>(*traitIter);
1866  if (pid) return pid->getRichSpectra();
1867  }
1868  return 0;
1869 }
1870 #ifndef __NO_STRANGE_MUDST__
1872  DEBUGMESSAGE2("");
1874  if (!maker) throw StMuExceptionNullPointer("no StrangeMuDstMaker",__PRETTYF__);
1875 
1876  StStrangeEvMuDst *ev=0;
1877  StV0MuDst *v0=0;
1878  StStrangeAssoc *assoc=0;
1879  StXiMuDst *xi=0;
1880  StKinkMuDst *kink=0;
1881  StV0Mc *v0Mc=0;
1882  StXiMc *xiMc=0;
1883  StKinkMc *kinkMc=0;
1884  TCut *strangeCut=0;
1885 
1886  addType(maker->GetEvClonesArray(), mStrangeArrays[0],ev);
1887  addType(maker->GetEvMcArray(), mStrangeArrays[1],ev);
1888 
1889  addType(maker->GetV0ClonesArray(), mStrangeArrays[2],v0);
1890  addType(maker->GetV0McArray(), mStrangeArrays[3],v0Mc);
1891  addType(maker->GetV0AssocArray(), mStrangeArrays[4],assoc);
1892 
1893  addType(maker->GetXiClonesArray(), mStrangeArrays[5],xi);
1894  addType(maker->GetXiMcArray(), mStrangeArrays[6],xiMc);
1895  addType(maker->GetXiAssocArray(), mStrangeArrays[7],assoc);
1896 
1897  addType(maker->GetKinkClonesArray(),mStrangeArrays[8],kink);
1898  addType(maker->GetKinkMcArray(), mStrangeArrays[9],kinkMc);
1899  addType(maker->GetKinkAssocArray(), mStrangeArrays[10],assoc);
1900 
1901  addType(maker->GetCutsArray(), mStrangeArrays[11],strangeCut);
1902 
1903 }
1904 #endif
1905 //-----------------------------------------------------------------------
1906 void StMuDstMaker::fillMC() {
1907  St_g2t_track *g2t_track = (St_g2t_track *) GetDataSet("geant/g2t_track"); if (!g2t_track) return;
1908  St_g2t_vertex *g2t_vertex = (St_g2t_vertex *) GetDataSet("geant/g2t_vertex"); if (!g2t_vertex) return;
1909  StG2TrackVertexMap::instance(g2t_track,g2t_vertex);
1910  mStEvent->setIdTruth();
1911  StMuMcVertex *mcvx = 0;
1912  StMuMcTrack *mctr = 0;
1913  g2t_vertex_st *vertex = g2t_vertex->GetTable();
1914  UInt_t NV = g2t_vertex->GetNRows();
1915  for (UInt_t i = 0; i < NV; i++) addType(mMCArrays[MCVertex], vertex[i], mcvx);
1916  g2t_track_st *track = g2t_track->GetTable();
1917  UInt_t NT = g2t_track->GetNRows();
1918  for (UInt_t i = 0; i < NT; i++) {
1919  if (track[i].pt<=1e-3) track[i].pt = -999;
1920  addType(mMCArrays[MCTrack], track[i], mctr);
1921  }
1922 
1923 }
1924 //-----------------------------------------------------------------------
1925 //-----------------------------------------------------------------------
1926 //-----------------------------------------------------------------------
1927 template <class T>
1928 int StMuDstMaker::addType(TClonesArray* tcaFrom, TClonesArray* &tcaTo ,T *t) {
1929  if (tcaFrom && tcaTo) {
1930  int n = tcaFrom->GetEntries();
1931  int counter = tcaTo->GetEntries();
1932  for (int i=0; i<n;i++) {
1933  // old new((*tcaTo)[counter++]) T( (T&)*tcaFrom->UncheckedAt(i) );
1934  new((*tcaTo)[counter++]) T( *(T*)(void*)tcaFrom->UncheckedAt(i) );
1935  }
1936  }
1937  return 0;
1938 }
1939 //-----------------------------------------------------------------------
1940 //-----------------------------------------------------------------------
1941 //-----------------------------------------------------------------------
1942 // int StMuDstMaker::addType(TClonesArray* tcaTo , StMuEmcCollection *t) {
1943 // int counter =-1;
1944 // if (tcaTo) {
1945 // counter = tcaTo->GetEntries();
1946 
1947 // new((*tcaTo)[counter]) StMuEmcCollection();
1948 // }
1949 // return counter;
1950 // }
1951 //-----------------------------------------------------------------------
1952 //-----------------------------------------------------------------------
1953 //-----------------------------------------------------------------------
1954 template <class T>
1955 int StMuDstMaker::addType(TClonesArray* tcaTo ,T &t) {
1956  int counter =-1;
1957  if (tcaTo) {
1958  counter = tcaTo->GetEntries();
1959  new((*tcaTo)[counter]) T( t );
1960  }
1961  return counter;
1962 }
1963 //-----------------------------------------------------------------------
1964 //-----------------------------------------------------------------------
1965 //-----------------------------------------------------------------------
1966 template <class T, class U>
1967 int StMuDstMaker::addType(TClonesArray* tcaTo ,U &u,T *t) {
1968  int counter =-1;
1969  if (tcaTo) {
1970  counter = tcaTo->GetEntries();
1971  DEBUGMESSAGE2("");
1972  new((*tcaTo)[counter]) T(u);
1973  }
1974  return counter;
1975 }
1976 //-----------------------------------------------------------------------
1977 //-----------------------------------------------------------------------
1978 //-----------------------------------------------------------------------
1979 string StMuDstMaker::buildFileName(string dir, string fileName, string extention){
1980  DEBUGMESSAGE3(dir.c_str());
1981  DEBUGMESSAGE3(fileName.c_str());
1982  DEBUGMESSAGE3(extention.c_str());
1983  fileName = dir + fileName + extention;
1984  while (fileName.find("//")!=string::npos) {
1985  int pos = fileName.find("//");
1986  fileName.erase(pos,1);
1987  }
1988  return fileName;
1989 }
1990 //-----------------------------------------------------------------------
1991 //-----------------------------------------------------------------------
1992 //-----------------------------------------------------------------------
1993 string StMuDstMaker::basename(string s){
1994  DEBUGVALUE3(s.c_str());
1995  string name(s);
1996  DEBUGVALUE3(name.c_str());
1997  size_t pos;
1998  pos = name.find_last_of("/");
1999  if (pos!=string::npos ) name.erase(0, pos );
2000  DEBUGVALUE3(name.c_str());
2001  pos = name.find_last_of(".");
2002  if (pos!=string::npos ) name.erase(pos,name.length()-pos );
2003  DEBUGVALUE3(name.c_str());
2004  pos = name.find_last_of(".");
2005  if (pos!=string::npos ) name.erase(pos,name.length()-pos );
2006  DEBUGVALUE3(name.c_str());
2007  return name;
2008 }
2009 //-----------------------------------------------------------------------
2010 //-----------------------------------------------------------------------
2011 //-----------------------------------------------------------------------
2012 string StMuDstMaker::dirname(string s){
2013  string name(s);
2014  DEBUGVALUE3(name.c_str());
2015  size_t pos;
2016  pos = name.find_last_of("/");
2017 
2018  if (pos != string::npos ) name.erase(pos, name.length());
2019  if (name == s) name =".";
2020 
2021  name=name+"/";
2022  DEBUGVALUE3(name);
2023  return name;
2024 }
2025 
2026 void StMuDstMaker::setProbabilityPidFile(const char* file) {
2027  std::ostringstream flnm;
2028 
2029  if ( ! file ){
2030  const char *PIDtable="PIDTableP01gl.root";
2031 
2032  flnm << getenv("STAR") << "/StarDb/dEdxModel/" << PIDtable << ends;
2033  FORCEDDEBUGMESSAGE("ATTENTION: pid table hardwired to " << flnm.str() );
2034 
2035  } else {
2036  flnm << file << ends;
2037  FORCEDDEBUGMESSAGE("Using pid table to user value " << flnm.str() );
2038  }
2039 
2040 
2041  if (mProbabilityPidAlgorithm)
2042  mProbabilityPidAlgorithm->readParametersFromFile(flnm.str());
2043 }
2044 //-----------------------------------------------------------------------
2045 //-----------------------------------------------------------------------
2046 //-----------------------------------------------------------------------
2047 void StMuDstMaker::printArrays()
2048 {
2049 // all stuff
2050  TClonesArray *tcl;
2051  for ( int i=0; i<__NALLARRAYS__; i++) {
2052  if (mStatusArrays[i]==0) continue;
2053  tcl = mAArrays[i];
2054  printf(" Array %s\t = %s::%s(%d)\n",
2055  StMuArrays::arrayNames[i],
2056  tcl->ClassName(),tcl->GetName(),tcl->GetEntriesFast());
2057  }
2058 
2059 }
2060 //-----------------------------------------------------------------------
2061 void StMuDstMaker::fillHddr()
2062 {
2063  StMuEvent *me = mStMuDst->event();
2064  if (me==0)
2065  return;
2066  StEventInfo &ei = me->eventInfo();
2067  StRunInfo &ri = me->runInfo();
2068  StEvtHddr *hd = GetEvtHddr();
2069 
2070  hd->SetRunNumber(ei.runId()) ;
2071  hd->SetEventType(ei.type().Data());
2072  hd->SetTriggerMask(ei.triggerMask()) ;
2073 //hd->SetInputTriggerMask(???);
2074  hd->SetBunchCrossingNumber(ei.bunchCrossingNumber(0),ei.bunchCrossingNumber(1));
2075  hd->SetCenterOfMassEnergy(ri.centerOfMassEnergy());
2076 //hd->SetBImpact (float b) ;
2077 //hd->SetPhiImpact(float p) ;
2078 //hd->SetPhImpact (float p) ;
2079 //hd->SetAEast(int a) ;
2080 //hd->SetZEast(int z) ;
2081 //hd->SetAWest(int a) ;
2082 //hd->SetZWest(int z) ;
2083 //hd->SetLuminosity(float lu) ;
2084  hd->SetGMTime((unsigned int)ei.time());
2085  hd->SetProdDateTime(ri.productionTime());
2086 //hd->SetIventNumber(int iv) ;
2087  hd->SetEventSize(ei.eventSize());
2088  hd->SetEventNumber(ei.id()) ;
2089 //hd->SetGenerType(int g);
2090 }
2091 
2092 //-----------------------------------------------------------------------
2094  mEmcCollection->setTowerData((StMuEmcTowerData*)(*(mEmcArrays[muEmcTow]))[0]);
2095  mEmcCollection->setPrsArray(2,mEmcArrays[muEmcPrs]);
2096  mEmcCollection->setSmdArray(3,mEmcArrays[muEmcSmde]);
2097  mEmcCollection->setSmdArray(4,mEmcArrays[muEmcSmdp]);
2098  mEmcCollection->setPrsArray(6,mEmcArrays[muEEmcPrs]);
2099  mEmcCollection->setSmdArray(7,mEmcArrays[muEEmcSmdu]);
2100  mEmcCollection->setSmdArray(8,mEmcArrays[muEEmcSmdv]);
2101 }
2102 //-----------------------------------------------------------------------
2103 void StMuDstMaker::connectFmsCollection() {
2104  mFmsCollection->setFmsHitArray(mFmsArrays[muFmsHit]);
2105  mFmsCollection->setFmsClusterArray(mFmsArrays[muFmsCluster]);
2106  mFmsCollection->setFmsPointArray(mFmsArrays[muFmsPoint]);
2107  mFmsCollection->setFmsInfoArray(mFmsArrays[muFmsInfo]);
2108 }
2109 //-----------------------------------------------------------------------
2110 void StMuDstMaker::connectRHICfCollection() {
2111  mRHICfCollection->setRHICfRawHitArray(mRHICfArrays[muRHICfRawHit]);
2112  mRHICfCollection->setRHICfHitArray(mRHICfArrays[muRHICfHit]);
2113  mRHICfCollection->setRHICfPointArray(mRHICfArrays[muRHICfPoint]);
2114 }
2115 //-----------------------------------------------------------------------
2116 void StMuDstMaker::connectFcsCollection() {
2117  LOG_INFO << "Setting Fcs arrays" << endm;
2118  mFcsCollection->setFcsHitArray(mFcsArrays[muFcsHit]);
2119  mFcsCollection->setFcsClusterArray(mFcsArrays[muFcsCluster]);
2120  mFcsCollection->setFcsPointArray(mFcsArrays[muFcsPoint]);
2121  mFcsCollection->setFcsInfoArray(mFcsArrays[muFcsInfo]);
2122 }
2123 //-----------------------------------------------------------------------
2124 void StMuDstMaker::connectFttCollection() {
2125  LOG_INFO << "Setting Ftt arrays" << endm;
2126  mFttCollection->setFttHitArray(mFttArrays[muFttRawHit]);
2127  mFttCollection->setFttClusterArray(mFttArrays[muFttCluster]);
2128  mFttCollection->setFttPointArray(mFttArrays[muFttPoint]);
2129 }
2130 //-----------------------------------------------------------------------
2131 void StMuDstMaker::connectFstCollection() {
2132  LOG_INFO << "Setting Fst arrays" << endm;
2133  mFstCollection->setFstRawHitArray(mFstArrays[muFstRawHit]);
2134  mFstCollection->setFstHitArray(mFstArrays[muFstHit]);
2135 }
2136 //-----------------------------------------------------------------------
2137 void StMuDstMaker::connectFwdTrackCollection() {
2138  LOG_INFO << "Setting FwdTrack arrays" << endm;
2139  mFwdTrackCollection->setFwdTrackArray(mFwdTrackArrays[muFwdTrack]);
2140 }
2141 //-----------------------------------------------------------------------
2142 void StMuDstMaker::connectPmdCollection() {
2143  mPmdCollection->setPmdHitArray(mPmdArrays[muPmdHit]);
2144  mPmdCollection->setCpvHitArray(mPmdArrays[muCpvHit]);
2145  mPmdCollection->setPmdClusterArray(mPmdArrays[muPmdCluster]);
2146  mPmdCollection->setCpvClusterArray(mPmdArrays[muCpvCluster]);
2147 }
2148 
2149 /***************************************************************************
2150  *
2151  * $Log: StMuDstMaker.cxx,v $
2152  * Revision 1.129 2019/02/21 13:32:54 jdb
2153  * Inclusion of ETOF MuDst code. This code adds support for the full set of ETOF data which includes EtofDigi, EtofHit, EtofHeader. The code essentially copies similar structures from StEvent and additionally rebuilds the maps between Digis and Hits. Accessor methods are added based on the pattern from BTOF to provide access to data at various levels. The code for accessing the PID traits provided by ETOF is also provided
2154  *
2155  * Revision 1.128 2018/03/01 20:42:38 jdb
2156  * mEpdUtil was not initialized in the second StMuDstMaker ctor, causing seg fault when deleted - added initialization to second ctor solves issue
2157  *
2158  * Revision 1.127 2018/02/27 04:15:02 jdb
2159  * added EPD support and fixed long standing bug in SetStatus, cleaned up
2160  *
2161  * Revision 1.126 2016/05/04 19:24:07 smirnovd
2162  * Addressed compiler warning by removing set but never used variables
2163  *
2164  * Revision 1.125 2015/11/06 17:47:16 jdb
2165  * Added StMuFmsInfo.{h,cxx} as a new branch for storing event-by-event FMS paramters
2166  *
2167  * Revision 1.124 2015/08/31 20:01:10 jdb
2168  * Added pT check back to StMuDstMaker::fillMC() - its removal was causing the test failure reported in [www.star.bnl.gov #3136] Assertion failed in StMuMcTrack.cxx for year 2009 MC test in DEV
2169  *
2170  * Revision 1.123 2015/08/28 18:36:03 jdb
2171  * Added Akios FMS codes
2172  *
2173  * Revision 1.118 2013/12/04 19:56:32 jdb
2174  * Added StMuMtdPidTraits.{cxx, h} added Mtd items to StMuMtdHit.h, StMuDst.{cxx,h}, StMuDstMaker.cxx, StMuTrack.{cxx,h}
2175  *
2176  * Revision 1.117 2013/07/23 11:02:59 jeromel
2177  * Undo changes (KF and other)
2178  *
2179  * Revision 1.115 2013/04/10 19:28:35 jeromel
2180  * Step back to 04/04 version (van aware) - previous changes may be recoverred
2181  *
2182  * Revision 1.113 2013/01/08 22:57:33 sangalin
2183  * Merged in FGT changes allowing for a variable number of timebins to be read out for each strip.
2184  *
2185  * Revision 1.112 2012/12/12 00:36:03 sangalin
2186  * Merged in updated StMuFgtCluster class format from Anselm Vossen.
2187  *
2188  * Revision 1.111 2012/11/30 20:29:41 sangalin
2189  * Removed cout debug statements in fillMtd().
2190  *
2191  * Revision 1.110 2012/11/16 12:31:37 jeromel
2192  * Fix catch without try
2193  *
2194  * Revision 1.109 2012/11/15 22:26:13 sangalin
2195  * Added the FGT. Fixed bugs in array offsets for the MTD.
2196  *
2197  * Revision 1.108 2012/10/04 18:57:59 fisyak
2198  * Add protection for empty emc raw data
2199  *
2200  * Revision 1.107 2012/09/28 22:38:05 tone421
2201  * Changed array stucture of MTD upon request of the TOF group. MTD arrays now on top level, rather than within __NARRAYS__
2202  *
2203  * Revision 1.106 2012/05/07 14:47:06 fisyak
2204  * Add handles for track to fast detector matching
2205  *
2206  * Revision 1.105 2011/10/17 00:19:13 fisyak
2207  * Active handing of IdTruth
2208  *
2209  * Revision 1.104 2011/08/18 18:41:36 fisyak
2210  * set max. tree size = 100 GB
2211  *
2212  * Revision 1.103 2011/05/04 19:51:32 tone421
2213  * Added MTD infomation
2214  *
2215  * Revision 1.102 2011/04/19 22:50:08 fisyak
2216  * Use default size of TTree (100 GB) for ROOT >= 5.26.0
2217  *
2218  * Revision 1.101 2011/04/08 01:25:50 fisyak
2219  * Add branches for MC track and vertex information, add IdTruth to tracks and vertices, reserve a possiblity to remove Strange MuDst
2220  *
2221  * Revision 1.100 2010/05/28 19:47:51 tone421
2222  * Removed a cout needed for test purposes in StMuDstMaker. Made sure StTriggerData objects copied into the MuDst have a debug value of 0..
2223  *
2224  * Revision 1.99 2010/05/26 04:25:50 tone421
2225  * Added StTriggerData arrays in muevent and fixed an issue with PMD arrays being read....
2226  *
2227  * Revision 1.98 2010/03/08 19:06:51 tone421
2228  * Two things. Global tracks how are filled with an index to primary at birth. Added StMuDst::fixTrackIndicesG(), which is used for matching the primary track indices to global tracks. Previously, this was quite slow - see this post:
2229  *
2230  * http://www.star.bnl.gov/HyperNews-star/protected/get/starsoft/8092/1/1/1.html
2231  *
2232  * for more details.
2233  *
2234  * Revision 1.97 2010/01/25 18:46:16 fine
2235  * RT #1826. Add protection against of the zero rps pointer and add the safe version of the StMuRpsCollection(const StRpsCollection & rps) ctor
2236  *
2237  * Revision 1.96 2010/01/25 03:57:39 tone421
2238  * Added FMS and Roman pot arrays
2239  *
2240  * Revision 1.95 2010/01/21 02:08:17 fine
2241  * RT #1803: Restore the broken MakeRead/MakeWrite interface to fix Skip event method
2242  *
2243  * Revision 1.94 2009/12/01 21:56:35 tone421
2244  * Implemented changes as per http://www.star.bnl.gov/rt2/Ticket/Display.html?id=1734
2245  *
2246  * Revision 1.93 2009/05/22 23:48:18 fine
2247  * Test I/O errors after filling the TTree
2248  *
2249  * Revision 1.92 2009/05/22 22:25:31 fine
2250  * Add the Zombue test for TFile ctors
2251  *
2252  * Revision 1.91 2009/04/27 20:50:43 perev
2253  * Change type of DataSet
2254  *
2255  * Revision 1.90 2009/03/10 23:43:53 jeromel
2256  * Set tree size to max size
2257  *
2258  * Revision 1.89 2009/02/20 16:37:44 tone421
2259  * *** empty log message ***
2260  *
2261  * Revision 1.87 2008/10/03 17:50:42 tone421
2262  * Added mVtxList(100); see http://www.star.bnl.gov/HyperNews-star/protected/get/starsoft/7529.html
2263  *
2264  * Revision 1.86 2008/04/14 21:32:12 fisyak
2265  * Remove stripping TObject from StDcaGeometry, because StDcaGeometry is inherit from StObject and this stripping brakes schema evolution
2266  *
2267  * Revision 1.85 2008/03/19 14:51:03 fisyak
2268  * Add two clone arrays for global and primary track covariance matrices, remove mSigmaDcaD and mSigmaDcaZ
2269  *
2270  * Revision 1.84 2007/08/31 01:55:05 mvl
2271  * Added protection against corrupted files by checking for return code -1 from TTree:GetEntry(). StMuDstMaker will silently skip these events; StMuIOMaker returns kStWarn.
2272  *
2273  * Revision 1.83 2007/08/02 20:46:47 mvl
2274  * Switch off Q-vector branhces in StMuDstMaker and increase version number in StMuEvent.
2275  * This is to avoid wranings when reading P07ib data which has Q-vector information stored with more recent libraries.
2276  *
2277  * Revision 1.82 2007/05/16 18:50:48 mvl
2278  * Cleanup of output. Replaced cout with LOG_INFO etc.
2279  *
2280  * Revision 1.81 2007/04/27 17:07:01 mvl
2281  * Added protection against StEvent::triggerData() == 0 in EZTREE.
2282  *
2283  * Revision 1.80 2007/04/20 06:26:00 mvl
2284  * Removed Q-vector calculation. Will implement utility class instead.
2285  *
2286  * Revision 1.78 2007/02/07 07:53:09 mvl
2287  * Added SetEventList function to read only pre-selected events (by J. Webb)
2288  *
2289  * Revision 1.77 2006/12/20 21:53:15 mvl
2290  * Added warning when file list not found (read mode)
2291  *
2292  * Revision 1.76 2006/09/10 00:58:43 mvl
2293  * Roll-back of previous changes for ROOT 5.12. Problem has been resolved inside ROOT.
2294  *
2295  * Revision 1.75 2006/08/18 20:09:51 fine
2296  * ROOT 5.12 bug workaround See: STAR Bug 741
2297  *
2298  * Revision 1.74 2006/07/28 18:25:11 mvl
2299  * Added call to StMuDst::setVertexIndex(0) to StMuDstMaker::read() to reset the current vertex index to 0 for every event
2300  *
2301  * Revision 1.73 2006/02/08 23:35:36 mvl
2302  * Added overloaded version for StIOInterface::GetFile() to return name
2303  * of current input or output file (depending on read or write mode)
2304  * StIOInterface::GetFileName() is an alias for this function.
2305  *
2306  * Revision 1.72 2005/10/18 17:55:43 mvl
2307  * Fixed initialisation problem of mCurrentFile, leading to potential segvio when creating MuDst
2308  *
2309  * Revision 1.71 2005/10/06 01:30:30 mvl
2310  * Changed some of the logic in StMuChainMaker: Now files are no longer opened
2311  * and checked at the start of the job, but simply added to the TChain. TChain
2312  * automatically skips corrupted files (this is a new feature).
2313  *
2314  * Revision 1.70 2005/08/19 19:46:05 mvl
2315  * Further updates for multiple vertices. The main changes are:
2316  * 1) StMudst::primaryTracks() now returns a list (TObjArray*) of tracks
2317  * belonging to the 'current' primary vertex. The index number of the
2318  * 'current' vertex can be set using StMuDst::setCurrentVertex().
2319  * This also affects StMuDst::primaryTracks(int i) and
2320  * StMuDst::numberOfprimaryTracks().
2321  * 2) refMult is now stored for all vertices, in StMuPrimaryVertex. The
2322  * obvious way to access these numbers is from the StMuprimaryVertex structures,
2323  * but for ebakcward compatibility a function is provided in StMuEvent as well
2324  * (this is the only function taht works for existing MuDst)
2325  *
2326  * As an aside, I've also changes the internals of StMuDst::createStEvent and
2327  * StMuDst::fixTrackIndices() to be able to deal with a larger range of index numbers for tracks as generated by Ittf.
2328  *
2329  * BIG FAT WARNING: StMudst2StEventMaker and StMuDstFilterMaker
2330  * do not fully support the multiple vertex functionality yet.
2331  *
2332  * Revision 1.69 2005/07/15 21:45:08 mvl
2333  * Added support for multiple primary vertices (StMuPrimaryVertex). Track Dcas are now calculated with repect to the first vertex in the list (highest rank), but another vertex number can be specified. Tarcks also store the index of the vertex they belong to (StMuTrack::vertexIndex())
2334  *
2335  * Revision 1.68 2005/04/12 21:56:29 mvl
2336  * Changes by Xin Dong for year-5 TOF data format: extra TClonesArray and routines to fill it from StEvent (StTofRawData).
2337  *
2338  * Revision 1.67 2004/11/29 15:53:22 mvl
2339  * Additions by Jan for Fpd ezTree
2340  *
2341  * Revision 1.66 2004/11/15 18:20:25 mvl
2342  * Added call to StMuDst::set() for V0-event-pointers in read()
2343  *
2344  * Revision 1.65 2004/10/31 23:43:21 mvl
2345  * Removed some warnings for files without EMC, PMD info.
2346  * Prevent filling of empty event when no stevent pointer.
2347  *
2348  * Revision 1.64 2004/10/28 00:11:33 mvl
2349  * Added stuff to support ezTree mode of MuDstMaker.
2350  * This is a special mode for fast-online processing of fast-detector data.
2351  *
2352  * Revision 1.63 2004/10/21 02:58:17 mvl
2353  * Removed some code from Make() (backward compatible EMc mode), to fix StMuIOMaker
2354  *
2355  * Revision 1.62 2004/10/19 01:42:29 mvl
2356  * Changes for splitting Emc and Pmd collections. Emc clusters and points dropped
2357  *
2358  * Revision 1.61 2004/09/18 01:28:11 jeromel
2359  * *** empty log message ***
2360  *
2361  * Revision 1.60 2004/05/04 13:26:23 jeromel
2362  * Oops .. Conflict resolution fixed.
2363  *
2364  * Revision 1.59 2004/05/04 13:17:11 jeromel
2365  * Changed to the documentation in doxygen format
2366  *
2367  * Revision 1.58 2004/05/04 00:09:23 perev
2368  *
2369  * // Selecting SetBranchStatus for particular MuDst branches
2370  * // Special names:
2371  * // MuEventAll - all branches related to StMuEvent
2372  * // StrangeAll - all branches related to StrangeMuDst
2373  * // EmcAll - all branches related to Emc
2374  * // PmdAll - all branches related to Pmd
2375  * // TofAll - all branches related to Tof
2376  * // By default all branches of MuDst are read. If user wants to read only some of
2377  * // them, then:
2378  * // SetStatus("*",0) // all branches off
2379  * // SetStatus("MuEventAll",1) // all standard MuEvent branches ON
2380  * // SetStatus("StrangeAll",1) // all standard Strange branches ON
2381  * // SetStatus("EmcAll" ,1) // all standard Emc branches ON
2382  * // SetStatus("PmdAll" ,1) // all standard Pmd branches ON
2383  * // SetStatus("TofAll" ,1) // all standard Tof branches ON
2384  * //
2385  * // SetStatus("XiAssoc" ,1) // Strange branch "XiAssoc" is ON
2386  * // Names of branches look StMuArrays::arrayTypes[]
2387  * // It allows to speed up reading MuDst significantly
2388  *
2389  * Revision 1.57 2004/04/26 00:13:28 perev
2390  * Cleanup+simplification
2391  *
2392  * Revision 1.56 2004/04/20 18:42:47 perev
2393  * remove redundant arrays
2394  *
2395  * Revision 1.55 2004/04/15 00:25:49 perev
2396  * fillHddr() added to fill time stamp ...
2397  *
2398  * Revision 1.54 2004/04/14 17:15:56 subhasis
2399  * Xin's TOF reinclusion
2400  *
2401  * Revision 1.53 2004/04/09 22:03:50 subhasis
2402  * after tof createevent fix by Xin
2403  *
2404  * Revision 1.52 2004/04/09 03:36:14 jeromel
2405  * Removed TOF support entirely for now as we need a working version ... Will
2406  * revisit later.
2407  *
2408  * Revision 1.51 2004/04/06 01:48:09 perev
2409  * Small leak + incorrect filing StMuTofHitCollection
2410  *
2411  * Revision 1.50 2004/04/02 03:24:54 jeromel
2412  * Changes implements PMD and TOF. TOF is clearly incomplete.
2413  *
2414  * Revision 1.49 2004/02/17 05:05:35 jeromel
2415  * One more hidden one found post-commit
2416  *
2417  * Revision 1.48 2004/02/17 04:56:36 jeromel
2418  * Extended help, added crs support, restored __GNUC__ for PRETTY_FUNCTION(checked once
2419  * more and yes, it is ONLY defined in GCC and so is __FUCTION__), use of a consistent
2420  * internal __PRETTYF__, return NULL if no case selected (+message) and protected against
2421  * NULL mChain.
2422  *
2423  * Revision 1.47 2003/11/17 22:16:55 perev
2424  * THack::DeleteClonesArray used for deleting, to avoid ROOT bad features
2425  *
2426  * Revision 1.46 2003/11/10 04:07:47 perev
2427  * again clear improved to avoid leaks
2428  *
2429  * Revision 1.45 2003/11/09 01:02:59 perev
2430  * more sofisticated clear() to fix leaks
2431  *
2432  * Revision 1.44 2003/11/03 22:24:45 perev
2433  * TClones::Clear added into StMuDstMaker::clear to avoid empty ebjects writing
2434  *
2435  * Revision 1.43 2003/10/30 20:08:13 perev
2436  * Check of quality added
2437  *
2438  * Revision 1.42 2003/10/28 00:03:46 perev
2439  * remove some debug lines
2440  *
2441  * Revision 1.41 2003/10/27 23:54:33 perev
2442  * weird template bug fized and templates simplified
2443  *
2444  * Revision 1.40 2003/10/23 04:08:29 perev
2445  * use SetBranchStatus fixed
2446  *
2447  * Revision 1.39 2003/10/20 19:50:13 perev
2448  * workaround added for TClonesArray::Delete + some cleanup of MuEmc
2449  *
2450  * Revision 1.38 2003/10/15 17:34:16 laue
2451  * StMuDstMaker: Reading fixed. Delete() changed back to Clear()
2452  * StMuEmcCollection: Re-implemented the DeleteThis() function,
2453  * This hoopefully fixed the memory leak when
2454  * writing MuDst again.
2455  * StMuTimer: ClassDef/ClassImp
2456  *
2457  * Revision 1.37 2003/10/12 03:43:56 perev
2458  * LeakOff TClonesArray::Clear replaced to Delete
2459  *
2460  * Revision 1.36 2003/10/08 21:17:15 laue
2461  * StMuEmcUtil updates from Alex Suaide
2462  * StMuDst and StMuDstMaker fixes to take the double inheritance of the
2463  * StKinkMuDsts into account. A void* had to be introduced when casting
2464  * TObject* to StKinkMuDst*.
2465  *
2466  * Revision 1.35 2003/10/03 15:26:07 laue
2467  * some moe arrays initialized
2468  *
2469  * Revision 1.34 2003/09/28 21:10:59 jeromel
2470  * More data members zeroed (would cause a crash on exit)
2471  *
2472  * Revision 1.33 2003/09/19 01:45:18 jeromel
2473  * A few problems hopefully fixed i.e. one constructor lacked zeroing
2474  * emcArrays were not zeroed, mStMuDst not zeroed.
2475  * For maintainability zeroArrays() added.
2476  *
2477  * Revision 1.32 2003/09/07 03:49:03 perev
2478  * gcc 3.2 + WarnOff
2479  *
2480  * Revision 1.31 2003/09/02 17:58:44 perev
2481  * gcc 3.2 updates + WarnOff
2482  *
2483  * Revision 1.30 2003/08/04 14:38:10 laue
2484  * Alex Suaide's updated for the EMC. Now EEMC is included.
2485  *
2486  * Revision 1.29 2003/03/25 19:08:06 laue
2487  * added StMuDst into TDataSet so that Valeri can pick it up for his
2488  * StEventDisplayMaker
2489  *
2490  * Revision 1.28 2003/03/06 01:34:18 laue
2491  * StAddRunInfoMaker is a make helper maker to add the StRunInfo for the
2492  * only year1 Au+Au 130GeV data
2493  *
2494  * Revision 1.27 2003/02/20 15:29:42 laue
2495  * StMuTriggerIdCollection added
2496  *
2497  * Revision 1.26 2003/02/19 15:38:10 jeromel
2498  * Modifications made to account for the new location of the PIDTable file.
2499  * The setProbabilityPidFile() method has been modified to take care of a default
2500  * file loading if unspecified. Messages will be displayed appropriatly.
2501  * Macros mdoofied to not call the method (leave it handled through the default
2502  * file).
2503  *
2504  * Revision 1.25 2003/02/07 23:47:53 laue
2505  * New EMC code. TObject arrays replaced by TClonesArrays (thanks to Alex)
2506  *
2507  * Revision 1.23 2003/02/05 22:10:00 laue
2508  * delete emc collection after being copied (when creating mudst)
2509  *
2510  * Revision 1.21 2003/01/29 03:04:57 laue
2511  * !!DIRTY FIX FOR StMuEmcCollection
2512  * !! Was memor leaking. Leak fixed, but slow and dirty.
2513  * !! Propose to change the structure as soon as possible.
2514  *
2515  * Revision 1.20 2003/01/09 18:59:45 laue
2516  * initial check in of new EMC classes and the changes required
2517  *
2518  * Revision 1.19 2002/11/27 15:07:31 laue
2519  * fix to run with standard root
2520  *
2521  * Revision 1.18 2002/11/07 17:12:22 laue
2522  * Comment changed.
2523  *
2524  * Revision 1.17 2002/08/27 19:05:56 laue
2525  * Minor updates to make the muDst from simulation work
2526  *
2527  * Revision 1.16 2002/08/20 19:55:49 laue
2528  * Doxygen comments added
2529  *
2530  * Revision 1.15 2002/05/20 18:57:18 laue
2531  * update for Christof
2532  *
2533  * Revision 1.14 2002/05/20 17:23:31 laue
2534  * StStrangeCuts added
2535  *
2536  * Revision 1.13 2002/05/04 23:56:30 laue
2537  * some documentation added
2538  *
2539  * Revision 1.12 2002/04/26 21:02:56 jeromel
2540  * Bug fix in dirname(). Still cannot get the arg3 bla/test.root mechanism to work
2541  * (but it does neither for everything else). Will come back to it.
2542  *
2543  * Revision 1.11 2002/04/23 21:35:32 laue
2544  * Changed name of StStraMuDstMaker to 'strangeMuDst' so that it can get picked
2545  * from the bfc.
2546  *
2547  * Revision 1.9 2002/04/11 14:19:30 laue
2548  * - update for RH 7.2
2549  * - decrease default arrays sizes
2550  * - add data base readerfor number of events in a file
2551  *
2552  * Revision 1.8 2002/04/01 22:42:30 laue
2553  * improved chain filter options
2554  *
2555  * Revision 1.7 2002/03/28 05:10:34 laue
2556  * update for running in the production
2557  *
2558  * Revision 1.6 2002/03/27 03:47:27 laue
2559  * better filter options
2560  *
2561  * Revision 1.5 2002/03/27 00:50:11 laue
2562  * bux fix from earlier check in
2563  *
2564  * Revision 1.3 2002/03/20 16:04:11 laue
2565  * minor changes, mostly added access functions
2566  *
2567  * Revision 1.2 2002/03/08 20:04:31 laue
2568  * change from two trees to 1 tree per file
2569  *
2570  * Revision 1.1 2002/03/08 17:04:17 laue
2571  * initial revision
2572  *
2573  *
2574  **************************************************************************/
2575 
2576 
2577 
2578 
2579 
2580 
2581 
2582 
2583 
2584 
2585 
2586 
2587 
2588 
2589 
2590 
2591 
2592 
void fillEzt(StEvent *ev)
static void fixMtdTrackIndices(TClonesArray *mtdHit, TClonesArray *primary, TClonesArray *global)
Definition: StMuDst.cxx:565
static void set(StMuDstMaker *maker)
set the pointers to the TClonesArrays
Definition: StMuDst.cxx:146
void connectEmcCollection()
routine to set up connection between mEmcCollection and Emc arrays
virtual int Make()
void fillStrange(StStrangeMuDstMaker *)
static const char * arrayTypes[__NALLARRAYS__]
&lt; names of the classes, the TClonesArrays are arrays of this type
Definition: StMuArrays.h:120
StMuDstMaker(const char *name="MuDst")
Default constructor.
void streamerOff()
void setProbabilityPidFile(const char *file=NULL)
Set the file from where the PID probability tables should be read.
virtual int Init()
virtual void Clear(Option_t *option="")
User defined functions.
void clearArrays()
virtual const char * GetFile() const
Returns name of current input or output file, depending on mode (GetFileName does the same...
static void collectVertexTracks()
Helper function to collect tracks for the current prim vertex.
Definition: StMuDst.cxx:260
static void fixTrackIndicesG(int mult=1)
Definition: StMuDst.cxx:338
void setL3TrackFilter(StMuCut *c)
Set the track filter used for L3 tracks when creating muDsts from StEvent and writing to disk...
Definition: StMuDstMaker.h:429
static void fixETofTrackIndices(TClonesArray *btofHit, TClonesArray *primary, TClonesArray *global)
Definition: StMuDst.cxx:485
void assignArrays()
void createArrays()
bool pass(const StEvent *)
called by user code, returns true if argument passes cuts, else false
Definition: StMuCut.h:68
void setTrackFilter(StMuCut *c)
Set the track filter used for all tracks (except the L3 tracks) when creating muDsts from StEvent and...
Definition: StMuDstMaker.h:428
Definition: Stypes.h:43
StMuBTofUtil * mBTofUtil
dongx
Definition: StMuDstMaker.h:248
static void fixTofTrackIndices(TClonesArray *btofHit, TClonesArray *primary, TClonesArray *global)
dongx
Definition: StMuDst.cxx:409
Definition: StV0Mc.hh:17
TClonesArray ** mBTofArrays
dongx
Definition: StMuDstMaker.h:401
static StMuEvent * event()
returns pointer to current StMuEvent (class holding the event wise information, e.g. event number, run number)
Definition: StMuDst.h:320
Stores information for tiles in STAR Event Plane Detector.
Definition: StEpdHit.h:43
static int arrayCounters[__NALLARRAYS__]
&lt; number of entries in current event, currently not used
Definition: StMuArrays.h:164
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:40
static void setProbabilityPidCentrality(Double_t cent)
Sets the centrality for calculating Aihong&#39;s pid.
Definition: StMuTrack.h:268
virtual void openWrite(string fileName)
protected:
Definition: StXiMc.hh:17
void SetStatus(const char *arrType, int status)
static void setProbabilityPidAlgorithm(StuProbabilityPidAlgorithm *)
Sets the StuProbabilityPidAlgorithm. Important in order to calculate Aihong&#39;s pids.
Definition: StMuTrack.h:267
virtual int Finish()
void fillBTof(StEvent *ev)
dongx
Definition: StFileI.h:13
static int arraySizes[__NALLARRAYS__]
&lt; maximum sizes of the TClonesArrays
Definition: StMuArrays.h:142
void addTrackNode(const StEvent *ev, const StTrackNode *node, StMuCut *cut, TClonesArray *gTCA=0, TClonesArray *pTCA=0, TClonesArray *oTCA=0, TClonesArray *covgTCA=0, TClonesArray *covpTCA=0, bool l3=false)
void setIndex2Cov(Int_t i)
Bingchu.
Definition: StMuTrack.h:148
void fillETof(StEvent *ev)
jdb &amp; fseck
static void setVertexIndex(Int_t vtx_id)
Set the index number of the current primary vertex (used by both primaryTracks() functions and for St...
Definition: StMuDst.cxx:273
int addTrack(TClonesArray *tca, const StEvent *event, const StTrack *track, const StVertex *vtx, StMuCut *cut, int index2Global, bool l3=false, TClonesArray *covgTCA=0, TClonesArray *covpTCA=0)
Definition: Stypes.h:45