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