StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StEmcMicroDstMaker.cxx
1 #include "StEmcMicroDstMaker.h"
2 #include "StEvent.h"
3 #include "StEventTypes.h"
4 #include "StEmcMicroEvent.h"
5 #include "StIOMaker/StIOMaker.h"
6 #include "StEmcUtil/filters/StEmcFilter.h"
7 #include "StEmcMicroUtil.h"
8 #include "StChain/StEvtHddr.h"
9 #include "StTreeMaker/StTreeMaker.h"
10 #include "StTree.h"
11 
12 ClassImp(StEmcMicroDstMaker);
16 {
17 
18  mGFilter = new StEmcFilter();
19  mPFilter = new StEmcFilter();
20 
21  mPFilter->setEmcPresent(kTRUE);
22  mPFilter->setPrintLog(kTRUE);
23  mPFilter->setHaveVertex(kTRUE);
24  mPFilter->setHavePrimaries(kTRUE);
25  mPFilter->setZVertexCut(100);
26  mPFilter->setMinMult(0);
27  mPFilter->setMaxMult(500000);
28 
29  mGFilter->setPtCut(0.0);
30  mGFilter->setFitPointsCut(0);
31  mGFilter->setMustProjectEmc(kFALSE);
32 
33  mPFilter->setPtCut(0.0);
34  mPFilter->setFitPointsCut(0);
35  mPFilter->setDCACut(3);
36  mPFilter->setEtaMin(-4.0);
37  mPFilter->setEtaMax( 4.0);
38  mPFilter->setMustProjectEmc(kFALSE);
39  // end cuts
40 
41  mEventDir="./";
42  mEventFileOld="";
43  mEventFile="";
44  mMicroEventChain=NULL;
45  mDoRead=kFALSE;
46  mDoCreateStEvent=kTRUE;
47 
48  mDoSavePrimaries = kTRUE;
49  mDoSaveGlobals = kTRUE;
50  mDoSaveEmc = kTRUE;
51  mDoSaveFpd = kTRUE;
52  mDoSaveV0 = kTRUE;
53 
54  mOldMaker = NULL;
55 
56  mStart = 0;
57 
58 }
59 StEmcMicroDstMaker::~StEmcMicroDstMaker()
60 {
61  if(mMicroUtil) delete mMicroUtil;
62  if(mGFilter) delete mGFilter;
63  if(mPFilter) delete mPFilter;
64 }
65 //------------------------------------------------------------------------
69 {
70  mAcc=new TH1F("macc","Events accepted and rejected",2,-0.5,1.5);
71  mMicroUtil = new StEmcMicroUtil();
72  mMicroUtil->setPrimaryFilter(mPFilter);
73  mMicroUtil->setGlobalFilter(mGFilter);
74 
75  mMicroUtil->setSavePrimaries(mDoSavePrimaries);
76  mMicroUtil->setSaveGlobals(mDoSaveGlobals);
77  mMicroUtil->setSaveEmc(mDoSaveEmc);
78  mMicroUtil->setSaveFpd(mDoSaveFpd);
79  mMicroUtil->setSaveV0(mDoSaveV0);
80 
81  mCurMicroEvent=0;
82  mFileCounter = 0;
83  mAccEv = 0;
84  return StMaker::Init();
85 }
86 //------------------------------------------------------------------------
90 {
91  cout <<"***** StEmcMicroDstMaker::Make() *******************************\n";
92  if(mMicroEvent) delete mMicroEvent;
93  mMicroEvent = 0;
94  mStEvent = 0;
95 
96  if(!mDoRead)
97  {
98  mStEvent=NULL;
99  mStEvent = (StEvent*) GetInputDS("StEvent");
100  if(!mStEvent) return kStWarn;
101 
102  Float_t BF=0.5;
103  StEventSummary* summary = mStEvent->summary();
104  if(summary) BF = summary->magneticField()/10.;
105 
106  mPFilter->setBField(BF);
107  mGFilter->setBField(BF);
108 
109  //event filter
110  if (!mPFilter->accept(mStEvent))
111  {
112  cout << "StEmcMicroDstMaker::Make(): event not accepted." << endl;
113  mAcc->Fill(0);
114  return kStOK;
115  }
116  mAcc->Fill(1);
117  mAccEv++;
118 
119  // at this point, the event has been accepted.
120  StIOMaker* IO=(StIOMaker*)GetMaker("IOMaker");
121  if(IO)
122  {
123  const char* chp = strrchr(IO->GetFile(),'/');
124  if(chp) mEventFile = strrchr(IO->GetFile(),'/')+1;
125  else mEventFile = IO->GetFile();
126  }
127  // char *strrchr(const char *s, int c);
128  // DESCRIPTION
129  // The strrchr() function returns a pointer to the last
130  // occurrence of the character c in the string s.
131  else
132  {
133  //try to get from StTreeMaker
134  StTreeMaker *tree=(StTreeMaker*)GetMaker("outputStream");
135  if(tree)
136  {
137  StTree *t = tree->GetTree();;
138  mEventFile=t->GetBaseName();
139  }
140  else
141  {
142  if(mAccEv%100==0) mFileCounter++;
143  char f[100];
144  sprintf(f,"EmcEvent.file%04d",mFileCounter);
145  mEventFile = f;
146  }
147  }
148 
149  if(mOldMaker) mEventFile = mOldMaker->getCurrentFile();
150  // if the input event file has changed,
151  // close the Micro dst file and creates another one
152  if(mEventFile!=mEventFileOld)
153  {
154  if(mMicroDstFile) if(mMicroDstFile->IsOpen())
155  {
156  mMicroDstFile->Write(0,TObject::kOverwrite);
157  mMicroDstFile->Close();
158  }
159  if (mMicroEvent) delete mMicroEvent;
160  if (mMicroDstFile) delete mMicroDstFile;
161  mMicroEvent=NULL;
162  mMicroDstFile=NULL;
164  mEventFileOld=mEventFile;
165  }
166  else cout <<"Filename being written: "<<mEventFile.Data()<<endl;
167  mMicroEvent= mMicroUtil->getMicroEvent(mStEvent);
168  if(mOldMaker)
169  {
170  // stuff that are not converted back in StEvent...
171  StEmcMicroEvent* old=mOldMaker->getMicroEvent();
172  mMicroEvent->setCTB(old->getCTB());
173  mMicroEvent->setZDCe(old->getZDCe());
174  mMicroEvent->setZDCw(old->getZDCw());
175  mMicroEvent->setZVertexZDC(old->getZVertexZDC());
176  mMicroEvent->setBBCe(old->getBBCe());
177  mMicroEvent->setBBCw(old->getBBCw());
178  mMicroEvent->setBBCNHits(old->getBBCNHits());
179  mMicroEvent->setZVertexBBC(old->getZVertexBBC());
180  mMicroEvent->setBunchCrossing7bit(old->getBunchCrossing7bit());
181  mMicroEvent->setBunchCrossing(old->getBunchCrossing());
182  mMicroEvent->setSpinBits(old->getSpinBits());
184  StFpdMicroCollection *oldfpd = old->getFpd();
185  if(oldfpd)
186  {
187  fpd->setToken(oldfpd->getToken());
188  fpd->setSumAdcNorth(oldfpd->getSumAdcNorth());
189  fpd->setSumAdcSouth(oldfpd->getSumAdcSouth());
190  fpd->setSumAdcTop(oldfpd->getSumAdcTop());
191  fpd->setSumAdcBottom(oldfpd->getSumAdcBottom());
192  fpd->setSumAdcPreShower1(oldfpd->getSumAdcPreShower1());
193  fpd->setSumAdcPreShower2(oldfpd->getSumAdcPreShower2());
194  fpd->setSumAdcSmdX(oldfpd->getSumAdcSmdX());
195  fpd->setSumAdcSmdY( oldfpd->getSumAdcSmdY());
196  mMicroEvent->setFpd(fpd);
197  }
198  }
199  cout << "before Fill()" << endl;
200  mEmcTree->Fill();
201  cout << "after Fill()" << endl;
202  }
203  else //read Mode
204  {
205  if(mMicroEventChain)
206  {
207  if(mNMicroEvents>0)
208  {
209  if(mCurMicroEvent<mNMicroEvents)
210  {
211  mMicroEvent = new StEmcMicroEvent();
212  mStEvent = 0;
213  mMicroEventChain->SetBranchAddress("MicroEvent",&mMicroEvent);
214  mMicroEventChain->GetEntry(mStart+mCurMicroEvent++);
215 
216  UInt_t GMTTime = (UInt_t)mMicroEvent->getEventTime();
217  StEvtHddr *hd = (StEvtHddr*)GetDataSet("EvtHddr");
218  if(!hd) { hd = new StEvtHddr(); AddData(hd); }
219  hd->SetGMTime(GMTTime);
220  mEventFile = strrchr(mMicroEventChain->GetFile()->GetName(),'/')+1;
221  if(mEventFile.EndsWith(".emcEvent.root"))
222  {
223  Int_t size=mEventFile.Sizeof();
224  mEventFile.Remove(size-15,14);
225  }
226  if(mDoCreateStEvent)
227  {
228  mStEvent = mMicroUtil->getStEvent(mMicroEvent);
229  AddData(mStEvent);
230  }
231  }
232  else return kStWarn;
233  }
234  else return kStWarn;
235  }
236  else return kStWarn;
237  }
238 
239  return kStOK;
240 }
241 //-----------------------------------------------------------------------
245 {
246  if(!mDoRead)
247  {
248  if (mMicroDstFile)
249  if(mMicroDstFile->IsOpen())
250  {
251  mMicroDstFile->Write(0,TObject::kOverwrite);
252  mMicroDstFile->Close();
253  }
254  }
255  return StMaker::Finish();
256 }
257 //----------------------------------------------------------------------------
262 {
263  Int_t split = 99; // by default split Event into sub branches
264  Int_t comp = 1; // by default file is compressed
265  Int_t bufsize = 256000;
266  if (split) bufsize /= 16;
267  // creat a Microevent and an output file
268  mMicroEvent = new StEmcMicroEvent();
269 
270  cout <<"Input file = "<<mEventFile.Data()<<endl;
271  TString* filestring = new TString(mEventDir);
272  filestring->Append(mEventFile);
273 
274  if(filestring->EndsWith(".event.root"))
275  {
276  Int_t size=filestring->Sizeof();
277  filestring->Remove(size-12,11);
278  }
279 
280  filestring->Append(".emcEvent.root");
281  cout <<"Output file = "<<filestring->Data()<<endl;
282 
283  mMicroDstFile = new TFile(filestring->Data(),"RECREATE","EMC Micro DST file");
284  if (!mMicroDstFile)
285  {
286  cout << "##### EmcMicroEventMaker: Warning: no MicroEvents file = " << filestring->Data() << endl;
287  return kStFatal;
288  }
289  mMicroDstFile->SetCompressionLevel(comp);
290  mEmcTree = new TTree("EmcTree", "EMC Micro Tree");
291  if (!mEmcTree)
292  {
293  cout << "##### EmcMicroEventMaker: Warning: No EmcMicroTree" << endl;
294  return kStFatal;
295  }
296  mEmcTree->SetAutoSave(1000000); // autosave when 1 Mbyte written
297  mEmcTree->Branch("MicroEvent", "StEmcMicroEvent", &mMicroEvent,bufsize,split);
298 
299  delete filestring;
300  return kStOK;
301 }
302 //----------------------------------------------------------------------------
306 {
307  if(!mMicroEventChain)
308  {
309  mMicroEventChain = new TChain("EmcTree","EMC Micro chain");
310  //mMicroEvent = new StEmcMicroEvent();
311  //mMicroEventChain->SetBranchAddress("MicroEvent",&mMicroEvent);
312  //mMicroEventChain->GetBranch("MicroEvent")->SetAutoDelete(kTRUE);
313  }
314  mMicroEventChain->Add(file);
315  mNMicroEvents=(Int_t)mMicroEventChain->GetEntries();
316  cout <<"Total number of events in the chain = "<<mNMicroEvents<<endl;
317  mDoRead=kTRUE;
318 }
319 
320 
321 
322 
323 
324 
325 
326 
void setSavePrimaries(Bool_t a)
Save or don&#39;t primary tracks. Default is kTRUE. Track selection is done by PrimaryFilter.
Int_t getZDCe() const
Return ZDC east.
void setMinMult(Float_t a)
minimum multiplicity. Default is 0.
Definition: StEmcFilter.h:237
Definition: IO.h:20
Int_t getBBCe() const
Return BBC east.
Int_t getSumAdcPreShower1() const
Return FPD Sum ADC Pre Shower 1.
void setSaveEmc(Bool_t a)
Save or don&#39;t EMC data. Default is kTRUE.
virtual void AddData(TDataSet *data, const char *dir=".data")
User methods.
Definition: StMaker.cxx:332
void setPrintLog(Bool_t a)
Print log if event is rejected.
Definition: StEmcFilter.h:230
Int_t getSumAdcSouth() const
Return FPD Sum ADC South.
StEmcMicroDstMaker(const Char_t *name="EmcMicroDst")
Int_t getEventTime() const
Return event time (unix time)
Int_t getSumAdcSmdX() const
Return FPD Sum ADC Shower Max X.
StEmcMicroEvent * getMicroEvent()
Return current StEmcMicroEvent.
void setEmcPresent(Bool_t a)
requires EMC data on event (kTRUE, kFALSE). Default value is kTRUE.
Definition: StEmcFilter.h:231
StFpdMicroCollection * getFpd() const
Return FPD Information.
void setHavePrimaries(Bool_t a)
requires at least one primary track (kTRUE, kFALSE). Default is kTRUE.
Definition: StEmcFilter.h:233
const char * getCurrentFile()
Return currect .event.root file.
void setFitPointsCut(Int_t a)
number of fit points in the track. Default is 0.
Definition: StEmcFilter.h:249
Float_t getZVertexZDC() const
Return Z vertex from ZDC.
Int_t getZDCw() const
Return ZDC west.
void setSaveV0(Bool_t a)
Save or don&#39;t V0 data. Default is kFALSE.
void setHaveVertex(Bool_t a)
requires a valid vertex (kTRUE, kFALSE). Default is kTRUE.
Definition: StEmcFilter.h:232
void setMustProjectEmc(Bool_t a)
requires that the track is projected on a valid EMC tower (kTRUE, kFALSE). Default is kTRUE...
Definition: StEmcFilter.h:250
StEvent * getStEvent(StEmcMicroEvent *)
Returns StEvent.
Int_t getSumAdcNorth() const
Return FPD Sum ADC North.
void setSaveGlobals(Bool_t a)
Save or don&#39;t global tracks. Default is kTRUE. Track selection is done by GlobalFilter.
Int_t getSpinBits() const
Return Spin bits.
Int_t getBBCNHits() const
Return BBC number of hits.
void setBField(Float_t a)
just set the BField value for track projection in tesla. Default is 0.5 tesla.
Definition: StEmcFilter.h:239
void setPtCut(Float_t a)
transverse momentm cut (GeV/c). Default is 0.
Definition: StEmcFilter.h:245
Int_t getBBCw() const
Return BBC west.
Int_t getSumAdcPreShower2() const
Return FPD Sum ADC Pre Shower 2.
Bool_t accept(StEvent *)
Returns kTRUE/kFALSE if StEvent is accepted.
Int_t getSumAdcBottom() const
Return FPD Sum ADC Bottom.
Int_t getSumAdcTop() const
Return FPD Sum ADC Top.
Definition: Stypes.h:42
Definition: Stypes.h:40
void setDCACut(Float_t a)
DCA cut on track in cm. Default is 300000.
Definition: StEmcFilter.h:244
void setZVertexCut(Float_t a)
Z Vertex threshold in cm. Default is 20 cm.
Definition: StEmcFilter.h:234
void setSaveFpd(Bool_t a)
Save or don&#39;t FPD data. Default is kFALSE.
void setEtaMin(Float_t a)
minimum eta of the track (momentum). Default is -10000.
Definition: StEmcFilter.h:247
UInt_t getBunchCrossing7bit() const
Return Bunch Crossing Id 7 bits.
void setEtaMax(Float_t a)
maximum eta of the track (momentum). Default is +10000.
Definition: StEmcFilter.h:248
StEmcMicroEvent * getMicroEvent(StEvent *)
Returns StEmcMicroEvent.
virtual Int_t Finish()
Definition: StMaker.cxx:776
void setMaxMult(Float_t a)
maximum multiplicity. Default is 100000.
Definition: StEmcFilter.h:238
Int_t getSumAdcSmdY() const
Return FPD Sum ADC Shower Max Y.
Float_t getZVertexBBC() const
Return BBC Z Vertex.
Definition: StTree.h:84
void addMicroEventFile(char *)
Add EMC micro DST file to read.
UInt_t getBunchCrossing() const
Return Bunch Crossing.
Int_t getCTB() const
Return CTB sum.
Int_t getToken() const
Return FPD token.