StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StHbtExample_McEvent.C
1 
2 class StChain;
3 StChain *chain=0;
4 
5 // keep pointers to Correlation Functions global, so you can have access to them...
6 class QinvCorrFctn;
7 QinvCorrFctn* QinvCF;
8 class QvecCorrFctn;
9 QvecCorrFctn* QvecCF;
10 class MinvCorrFctn;
11 MinvCorrFctn* MinvCF;
12 
13 void StHbtExample_McEvent(Int_t nevents=1,
14  const char *MainFile="/disk00001/star/auau200/venus412/default/b0_3/year_1b/hadronic_on/tfs_4/set0364_01_35evts.geant.root")
15 {
16 
17  // Dynamically link needed shared libs
18  gSystem->Load("St_base");
19  gSystem->Load("StChain");
20  gSystem->Load("St_Tables");
21  gSystem->Load("StUtilities"); // new addition 22jul99
22  gSystem->Load("StAnalysisUtilities"); // needed by V0dstMaker
23  gSystem->Load("StMagF");
24  gSystem->Load("StIOMaker");
25  gSystem->Load("StarClassLibrary");
26  gSystem->Load("StEvent");
27  gSystem->Load("StEventMaker");
28  gSystem->Load("StMcEvent");
29  gSystem->Load("StMcEventMaker");
30  gSystem->Load("StHbtMaker");
31  gSystem->Load("StV0MiniDstMaker");
32 
33  cout << "Dynamic loading done" << endl;
34 
35  chain = new StChain("StChain");
36  chain->SetDebug();
37 
38 
39  // Now we add Makers to the chain...
40 
41  StIOMaker* ioMaker = new StIOMaker("IO","r",MainFile,"bfcTree");
42  ioMaker->SetDebug();
43 
44  ioMaker->SetIOMode("r");
45  ioMaker->SetDebug();
46  ioMaker->SetBranch("*",0,"0"); //deactivate all branches
47  // ioMaker->SetBranch("dstBranch",0,"r"); //activate EventBranch
48  ioMaker->SetBranch("geantBranch",0,"r"); //activate EventBranch
49 
50 
51  // StEventMaker* eventMaker = new StEventMaker("events","title");
52  // cout << "Just instantiated StEventMaker... lets go StHbtMaker!" << endl;
53  StMcEventMaker* mcEventMaker = new StMcEventMaker("mcEvents","title");
54  cout << "Just instantiated StMcEventMaker... lets go StHbtMaker!" << endl;
55 
56  // UNCOMMENT THIS NEXT PART OUT IF YOU WANT V0's
57  //StV0MiniDstMaker* v0dst = new StV0MiniDstMaker("v0dst");
58  //cout << "Just instantiated StV0MiniDstMaker... lets go StHbt!" << endl;
59  //v0dst.SetV0VertexType(); //Set v0MiniDstMaker to find v0s not Xis
60  //v0dst.SetOutputFile("muv0dst.root"); // Set V0MiniDStMaker output file
61 
62 
63 
64  StHbtMaker* hbtMaker = new StHbtMaker("HBT","title");
65  cout << "StHbtMaker instantiated"<<endl;
66 
67 
68 
69  /* -------------- set up of hbt stuff ----- */
70  cout << "StHbtMaker::Init - setting up Reader and Analyses..." << endl;
71 
72  StHbtManager* TheManager = hbtMaker->HbtManager();
73 
74  // here, we instantiate the appropriate StHbtEventReader
75  // for STAR analyses in root4star, we instantiate StStandardHbtEventReader
76  //StStandardHbtEventReader* Reader = new StStandardHbtEventReader;
77  //Reader->SetTheEventMaker(eventMaker); // gotta tell the reader where it should read from
79  Reader->SetTheMcEventMaker(mcEventMaker); // gotta tell the reader where it should read from
80 
81  // UNCOMMENT THIS NEXT LINE OUT IF YOU WANT V0's
82  // Reader->SetTheV0Maker(v0dst); //Gotta tell the reader where to read the v0 stuff from
83 
84  // here would be the palce to plug in any "front-loaded" Event or Particle Cuts...
85  TheManager->SetEventReader(Reader);
86 
87  cout << "READER SET UP.... " << endl;
88 
89  // Hey kids! Let's make a microDST!
90  // in StHbt we do this by instantiating and plugging in a StHbtEventReader as a writer!
91  // the particular StHbtEventReader that we will use will write (and read) ASCII files
92  //
93  // StHbtAsciiReader* Writer = new StHbtAsciiReader;
94  // Writer->SetFileName("FirstMicroDst.asc");
95  // TheManager->SetEventWriter(Writer);
96  // cout << "WRITER SET UP.... " << endl;
97 
98  // 0) now define an analysis...
99  StHbtAnalysis* anal = new StHbtAnalysis;
100  // 1) set the Event cuts for the analysis
101  mikesEventCut* evcut = new mikesEventCut; // use "mike's" event cut object
102  evcut->SetEventMult(0,100000); // selected multiplicity range
103  evcut->SetVertZPos(-35.0,35.0); // selected range of vertex z-position
104  anal->SetEventCut(evcut); // this is the event cut object for this analsys
105  // 2) set the Track (particle) cuts for the analysis
106  mikesTrackCut* trkcut = new mikesTrackCut; // use "mike's" particle cut object
107  trkcut->SetNSigmaPion(-0.5,0.5); // number of Sigma in TPC dEdx away from nominal pion dEdx
108  trkcut->SetNSigmaKaon(-1000.0,1000.0); // number of Sigma in TPC dEdx away from nominal kaon dEdx
109  trkcut->SetNSigmaProton(-1000.0,1000.0); // number of Sigma in TPC dEdx away from nominal proton dEdx
110  trkcut->SetNHits(5,50); // range on number of TPC hits on the track
111  trkcut->SetPt(0.1,1.0); // range in Pt
112  trkcut->SetRapidity(-1.0,1.0); // range in rapidity
113  trkcut->SetDCA(0.0,0.5); // range in Distance of Closest Approach to primary vertex
114  trkcut->SetCharge(-1); // want negative pions
115  trkcut->SetMass(0.139); // pion mass
116  anal->SetFirstParticleCut(trkcut); // this is the track cut for the "first" particle
117  anal->SetSecondParticleCut(trkcut); // NOTE - it is also for the "second" particle -- i.e. identical particle HBT
118  // 3) set the Pair cuts for the analysis
119  mikesPairCut* paircut = new mikesPairCut; // use "mike's" pair cut object
120  anal->SetPairCut(paircut); // this is the pair cut for this analysis
121  // 4) set the number of events to mix (per event)
122  anal->SetNumEventsToMix(5);
123  // 5) now set up the correlation functions that this analysis will make
124  // this particular analysis will have two: the first is a Q-invariant correlation function
125  QinvCF = new QinvCorrFctn("mikesQinvCF",50,0.0,0.2); // defines a Qinv correlation function
126  anal->AddCorrFctn(QinvCF); // adds the just-defined correlation function to the analysis
127  // for this analysis, we will also (simultaneously) build a Q-vector correlation function
128  QvecCF = new QvecCorrFctn("randysQvecCF",50,0.0,0.2);
129  anal->AddCorrFctn(QvecCF); // adds the just-defined correlation function to the analysis
130 
131  // now add as many more correlation functions to the Analysis as you like..
132 
133  // 6) add the Analysis to the AnalysisCollection
134  TheManager->AddAnalysis(anal);
135 
136  // now, we define another analysis that runs simultaneously with the previous one.
137  // this one looks at K+K- correlations (so NONidentical particles) in invariant mass
138 
139  /* ****************************************** */
140  /* * franks phi analysis - by Frank Laue, OSU */
141  /* ****************************************** */
142  // 0) now define an analysis...
143  StHbtAnalysis* phiAnal = new StHbtAnalysis;
144  // 1) set the Event cuts for the analysis
145  mikesEventCut* phiEvcut = new mikesEventCut; // use "mike's" event cut object
146  phiEvcut->SetEventMult(0,100000); // selected multiplicity range
147  phiEvcut->SetVertZPos(-35.0,35.0); // selected range of vertex z-position
148  phiAnal->SetEventCut(phiEvcut); // this is the event cut object for this analsys
149  // 2) set the Track (particle) cuts for the analysis
150  mikesTrackCut* kaonTrkcut = new mikesTrackCut; // use "mike's" particle cut object
151  kaonTrkcut->SetNSigmaPion(-1000,1000.0); // number of Sigma in TPC dEdx away from nominal pion dEdx
152  kaonTrkcut->SetNSigmaKaon(-3.0,3.0); // number of Sigma in TPC dEdx away from nominal kaon dEdx
153  kaonTrkcut->SetNSigmaProton(-1000.,1000.0); // number of Sigma in TPC dEdx away from nominal proton dEdx
154  kaonTrkcut->SetNHits(0,50); // range on number of TPC hits on the track
155  kaonTrkcut->SetPt(0.1,2.0); // range in Pt
156  kaonTrkcut->SetRapidity(-2.0,2.0); // range in rapidity
157  kaonTrkcut->SetDCA(0.0,0.5); // range in Distance of Closest Approach to primary vertex
158  kaonTrkcut->SetCharge(+1); // want positive kaons
159  kaonTrkcut->SetMass(0.494); // kaon mass
160  phiAnal->SetFirstParticleCut(kaonTrkcut); // this is the track cut for the "first" particle
161  mikesTrackCut* antikaonTrkcut = new mikesTrackCut; // use "mike's" particle cut object
162  antikaonTrkcut->SetNSigmaPion(-1000.0,1000.0); // number of Sigma in TPC dEdx away from nominal pion dEdx
163  antikaonTrkcut->SetNSigmaKaon(-3.0,3.0); // number of Sigma in TPC dEdx away from nominal kaon dEdx
164  antikaonTrkcut->SetNSigmaProton(-1000.0,1000.0); // number of Sigma in TPC dEdx away from nominal proton dEdx
165  antikaonTrkcut->SetNHits(0,50); // range on number of TPC hits on the track
166  antikaonTrkcut->SetPt(0.1,2.0); // range in Pt
167  antikaonTrkcut->SetRapidity(-2.0,2.0); // range in rapidity
168  antikaonTrkcut->SetDCA(0.0,0.5); // range in Distance of Closest Approach to primary vertex
169  antikaonTrkcut->SetCharge(-1); // want negative kaons
170  antikaonTrkcut->SetMass(0.494); // kaon mass
171  phiAnal->SetSecondParticleCut(antikaonTrkcut); // this is the track cut for the "second" particle
172  // 3) set the Pair cuts for the analysis
173  mikesPairCut* phiPaircut = new mikesPairCut; // use "mike's" pair cut object
174  phiAnal->SetPairCut(phiPaircut); // this is the pair cut for this analysis
175  // 4) set the number of events to mix (per event)
176  phiAnal->SetNumEventsToMix(5);
177  // 5) now set up the correlation functions that this analysis will make
178  MinvCF = new MinvCorrFctn("franksMinvCF",100,0.98,1.18); // defines a Minv correlation function
179  phiAnal->AddCorrFctn(MinvCF); // adds the just-defined correlation function to the analysis
180  // now add as many more correlation functions to the Analysis as you like..
181  // 6) add the Analysis to the AnalysisCollection
182  TheManager->AddAnalysis(phiAnal);
183 
184 
185 
186  /* ------------------ end of setting up hbt stuff ------------------ */
187 
188 
189  // now execute the chain member functions
190 
191  if (chain->Init()){ // This should call the Init() method in ALL makers
192  cout << "Initialization failed \n";
193  goto TheEnd;
194  }
195  chain->PrintInfo();
196 
197 
198  // Event loop
199  int istat=0,iev=1;
200  EventLoop: if (iev <= nevents && !istat) {
201  cout << "StHbtExample -- Working on eventNumber " << iev << " of " << nevents << endl;
202  chain->Clear();
203  istat = chain->Make(iev);
204  if (istat) {cout << "Last event processed. Status = " << istat << endl;}
205  iev++; goto EventLoop;
206  }
207 
208 // good old Cint can't even handle a for-loop
209 // for (Int_t iev=1;iev<=nevents; iev++) {
210 // chain->Clear();
211 // int iret = chain->Make(iev); // This should call the Make() method in ALL makers
212 // if (iret) {
213 // cout << "StHbtExample.C -- chain returned nonzero value " << iret
214 // << " on event " << iev << endl;
215 // break;
216 // }
217 // } // Event Loop
218 
219 
220 
221  cout << "StHbtExample -- Done with event loop" << endl;
222 
223  chain->Finish(); // This should call the Finish() method in ALL makers
224  TheEnd:
225 }
virtual void SetIOMode(Option_t *iomode="w")
number of transactions
Definition: StIOInterFace.h:35
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StChain.cxx:77
virtual Int_t Finish()
Definition: StChain.cxx:85
Filling of all StMcEvent classes from g2t tables Transform all the data in the g2t tables into the co...
virtual Int_t Make()
Definition: StChain.cxx:110