StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFlowEvent.cxx
1 //
3 // $Id: StFlowEvent.cxx,v 1.63 2010/09/30 19:30:26 posk Exp $
4 //
5 // Author: Raimond Snellings and Art Poskanzer
6 // FTPC added by Markus Oldenburg, MPI, Dec 2000
7 // Cumulants added by Aihong Tang, KSU, Nov 2001
8 //
10 //
11 // Description: A subset of StEvent
12 //
14 
15 #include <Stiostream.h>
16 #include <stdlib.h>
17 #include <math.h>
18 #include "StFlowEvent.h"
19 #include "StFlowTrackCollection.h"
20 #include "StFlowSelection.h"
21 #include "StFlowConstants.h"
22 #include "PhysicalConstants.h"
23 #include "SystemOfUnits.h"
24 #include "StEnumerations.h"
25 #include "TVector2.h"
26 #include "TComplex.h"
27 #include "TH1.h"
28 #define PR(x) cout << "##### FlowEvent: " << (#x) << " = " << (x) << endl;
29 
30 ClassImp(StFlowEvent)
31 Double_t StFlowEvent::mZDCSMDCenterex = Flow::zdcsmd_ex0;
32 Double_t StFlowEvent::mZDCSMDCenterey = Flow::zdcsmd_ey0;
33 Double_t StFlowEvent::mZDCSMDCenterwx = Flow::zdcsmd_wx0;
34 Double_t StFlowEvent::mZDCSMDCenterwy = Flow::zdcsmd_wy0;
35 
36 // CumalentMaker:
37 // Main TPC particles do not participate in G Mix calc. for the 1st Har. sel2
38 // for sel 1, C{3} is ( har1(all) + har1(all) + har2(all) )
39 // differential v1{3} is ( har1(all) + har1(all) + har2(all) )
40 //
41 // for sel 2, C{3} is ( har1(FTPC)+ har1(FTPC)+ har2(TPC) )
42 // differential v1{3} is ( har1(FTPC)+ har1(FTPC)+ har2(TPC) )
43 //
44 // Please note that for the mixed har analysis, only the 1st har. make sense anyway.
45 
46 Float_t StFlowEvent::mV1TPCDetctWgtG_Mix[Flow::nSels] ={1., 0.};
47 Float_t StFlowEvent::mV1FtpcWestDetctWgtG_Mix[Flow::nSels] ={1., 1.};
48 Float_t StFlowEvent::mV1FtpcEastDetctWgtG_Mix[Flow::nSels] ={1., 1.};
49 
50 Float_t StFlowEvent::mV2TPCDetctWgtG_Mix[Flow::nSels] ={1., 1.};
51 Float_t StFlowEvent::mV2FtpcWestDetctWgtG_Mix[Flow::nSels] ={1., 0.};
52 Float_t StFlowEvent::mV2FtpcEastDetctWgtG_Mix[Flow::nSels] ={1., 0.};
53 
54 Float_t StFlowEvent::mEtaTpcCuts[2][2][Flow::nSels] = {{{0.,0.5},
55  {0.,0.} },
56  {{1.0,2.},
57  {1.0,1.} }};
58 // Float_t StFlowEvent::mEtaTpcCuts[2][2][Flow::nSels] = {{{0.,0.},
59 // {0.,0.} },
60 // {{1.0,1.},
61 // {1.0,1.} }};
62 Float_t StFlowEvent::mEtaFtpcCuts[4][2][Flow::nSels] = {{{-4.0,-4.0},
63  {-4.0,-4.0}},
64  {{-2.7,-2.7},
65  {-2.7,-2.7} },
66  {{2.7,2.7},
67  {2.7,2.7} },
68  {{4.0,4.0},
69  {4.0,4.0} }};
70 
71 Float_t StFlowEvent::mPtTpcCuts[2][2][Flow::nSels] = {{{0.15,0.15},
72  {0.15,0.15} },
73  {{2.,2.},
74  {2.,2.} }};
75 Float_t StFlowEvent::mPtFtpcCuts[2][2][Flow::nSels] = {{{0.15,0.15},
76  {0.15,0.15} },
77  {{2.,2.},
78  {2.,2.} }};
79 
80 Float_t StFlowEvent::mPiPlusCuts[2] = {-3., 3.};
81 Float_t StFlowEvent::mPiMinusCuts[2] = {-3., 3.};
82 Float_t StFlowEvent::mProtonCuts[2] = {-3., 3.};
83 Float_t StFlowEvent::mAntiProtonCuts[2] = {-3., 3.};
84 Float_t StFlowEvent::mKMinusCuts[2] = {-3., 3.};
85 Float_t StFlowEvent::mKPlusCuts[2] = {-3., 3.};
86 Float_t StFlowEvent::mDeuteronCuts[2] = {-3., 3.};
87 Float_t StFlowEvent::mAntiDeuteronCuts[2] = {-3., 3.};
88 Float_t StFlowEvent::mElectronCuts[2] = {-3., 3.};
89 Float_t StFlowEvent::mPositronCuts[2] = {-3., 3.};
90 Float_t StFlowEvent::mDcaGlobalTpcCuts[2] = { 0., 0.};
91 Float_t StFlowEvent::mDcaGlobalFtpcCuts[2] = { 0., 0.};
92 Float_t StFlowEvent::mPtWgtSaturation = 2.;
93 Bool_t StFlowEvent::mPtWgt = kTRUE;
94 Bool_t StFlowEvent::mEtaWgt = kTRUE;
95 Bool_t StFlowEvent::mProbPid = kFALSE;
96 Bool_t StFlowEvent::mEtaSubs = kFALSE;
97 Bool_t StFlowEvent::mRanSubs = kFALSE;
98 Bool_t StFlowEvent::mFirstLastPhiWgt = kFALSE;
99 Bool_t StFlowEvent::mFirstLastPoints = kFALSE;
100 Bool_t StFlowEvent::mUseZDCSMD = kFALSE;
101 Char_t StFlowEvent::mPid[10] = {'\0'};
102 
103 //-----------------------------------------------------------
104 
105 StFlowEvent::StFlowEvent() {
106  // Make a new track collection
107 
108  pTrackCollection = new StFlowTrackCollection;
109 
110 }
111 
112 //-----------------------------------------------------------
113 
114 StFlowEvent::~StFlowEvent() {
115 
116  delete pTrackCollection;
117 
118 }
119 
120 //-------------------------------------------------------------
121 
122 Double_t StFlowEvent::PhiWeightRaw(Int_t selN, Int_t harN, StFlowTrack*
123  pFlowTrack) const {
124  // Weight for making the event plane isotropic in the lab.
125 
126  bool oddHar = (harN+1) % 2; // Choose odd or even
127  harN = oddHar ? 0 : 1;
128 
129  StTrackTopologyMap map = pFlowTrack->TopologyMap();
130  float phi = pFlowTrack->Phi();
131  if (phi < 0.) phi += twopi;
132  float eta = pFlowTrack->Eta();
133 
134  Double_t phiWgt = 0.;
135  int n = 0;
136 
137  float vertexZ = mVertexPos.z();
138  if (map.hasHitInDetector(kTpcId) || (map.data(0) == 0 && map.data(1) == 0)) {
139  // Tpc track, or no topologyMap
140  n = (int)((phi/twopi)*Flow::nPhiBins);
141  if (mFirstLastPhiWgt) {
142  float zFirstPoint = pFlowTrack->ZFirstPoint();
143  float zLastPoint = pFlowTrack->ZLastPoint();
144  if (zFirstPoint > 0. && zLastPoint > 0.) {
145  phiWgt = mPhiWgtFarWest[selN][harN][n];
146  } else if (zFirstPoint > 0. && zLastPoint < 0.) {
147  phiWgt = mPhiWgtWest[selN][harN][n];
148  } else if (zFirstPoint < 0. && zLastPoint > 0.) {
149  phiWgt = mPhiWgtEast[selN][harN][n];
150  } else {
151  phiWgt = mPhiWgtFarEast[selN][harN][n];
152  }
153  } else {
154  if (eta > 0. && vertexZ > 0.) {
155  phiWgt = mPhiWgtFarWest[selN][harN][n];
156  } else if (eta > 0. && vertexZ < 0.) {
157  phiWgt = mPhiWgtWest[selN][harN][n];
158  } else if (eta < 0. && vertexZ > 0.) {
159  phiWgt = mPhiWgtEast[selN][harN][n];
160  } else {
161  phiWgt = mPhiWgtFarEast[selN][harN][n];
162  }
163  }
164  }
165 
166  else if (map.trackFtpcEast()) { // Ftpc east track == eta < 0.
167  n = (int)((phi/twopi)*Flow::nPhiBinsFtpc);
168  if (vertexZ < 0.) {
169  phiWgt = mPhiWgtFtpcFarEast[selN][harN][n];
170  }
171  else { // vertexZ > 0.
172  phiWgt = mPhiWgtFtpcEast[selN][harN][n];
173  }
174  }
175 
176  else if (map.trackFtpcWest()) { // Ftpc west track == eta > 0.
177  n = (int)((phi/twopi)*Flow::nPhiBinsFtpc);
178  if (vertexZ > 0.) {
179  phiWgt = mPhiWgtFtpcFarWest[selN][harN][n];
180  }
181  else { // vertexZ < 0.
182  phiWgt = mPhiWgtFtpcWest[selN][harN][n];
183  }
184  }
185 
186  return phiWgt;
187 }
188 
189 //-------------------------------------------------------------
190 
191 Double_t StFlowEvent::Weight(Int_t selN, Int_t harN, StFlowTrack*
192  pFlowTrack) const {
193  // Weight for enhancing the resolution.
194 
195  //bool oddHar = (harN+1) % 2;
196  Double_t wgt = 1.;
197  wgt *= PtAbsWgtValue(pFlowTrack->Pt());
198  float eta = pFlowTrack->Eta();
199  wgt *= EtaAbsWgtValue(eta);
200  //if (oddHar && eta < 0.) { wgt *= -1. ; }
201  if (harN == 0 && eta < 0.) { wgt *= -1. ; } // reverse only for 1st har
202 
203  return wgt;
204 }
205 
206 //-------------------------------------------------------------
207 
208 Double_t StFlowEvent::PtAbsWgtValue(Double_t pt) const {
209 
210  return ((mPtWgt) ? ((pt < mPtWgtSaturation) ? pt : mPtWgtSaturation) : 1.);
211 }
212 
213 //-------------------------------------------------------------
214 
215 Double_t StFlowEvent::EtaAbsWgtValue(Double_t eta) const {
216 
217  return ((mEtaWgt) ? ((fabs(eta)<0.005) ? 0.005 : fabs(eta)) : 1.);
218 }
219 
220 //-----------------------------------------------------------------------
221 
222 Double_t StFlowEvent::PhiWeight(Int_t selN, Int_t harN, StFlowTrack*
223  pFlowTrack) const {
224  // Weight for making the event plane isotropic in the lab and
225  // enhancing the resolution.
226 
227  Double_t phiWgtRaw = PhiWeightRaw(selN, harN, pFlowTrack);
228  Double_t weight = Weight(selN, harN, pFlowTrack);
229 
230  return phiWgtRaw * weight;
231 }
232 
233 //-------------------------------------------------------------
234 
235 Double_t StFlowEvent::ZDCSMD_PsiWgtEast() {
236  // Get psi weight for east ZDCSMD
237 
238  TH1F *mZDCSMD_PsiWgt = new TH1F("ZDCSMD_PsiWgt","ZDCSMD_PsiWgt",
239  Flow::zdcsmd_nPsiBins,-twopi/2.,twopi/2.);
240  Int_t n = mZDCSMD_PsiWgt->FindBin(ZDCSMD_PsiEst());
241  mZDCSMD_PsiWgt->Delete();
242 
243  return mZDCSMD_PsiWgtEast[n-1];
244 }
245 
246 //-------------------------------------------------------------
247 
248 Double_t StFlowEvent::ZDCSMD_PsiWgtWest() {
249  // Get psi weight for west ZDCSMD
250 
251  TH1F *mZDCSMD_PsiWgt = new TH1F("ZDCSMD_PsiWgt","ZDCSMD_PsiWgt",
252  Flow::zdcsmd_nPsiBins,-twopi/2.,twopi/2.);
253  Int_t n = mZDCSMD_PsiWgt->FindBin(ZDCSMD_PsiWst());
254  mZDCSMD_PsiWgt->Delete();
255 
256  return mZDCSMD_PsiWgtWest[n-1];
257 }
258 
259 //-------------------------------------------------------------
260 Double_t StFlowEvent::ZDCSMD_PsiWgtFull() {
261 
262  TH1F *mZDCSMD_PsiWgt=new TH1F("ZDCSMD_PsiWgt","ZDCSMD_PsiWgt",Flow::zdcsmd_nPsiBins,0.,twopi);
263  StFlowSelection* mFlowSelect;
264  Int_t n =mZDCSMD_PsiWgt->FindBin(Q(mFlowSelect).Phi());
265  mZDCSMD_PsiWgt->Delete();
266 
267  return mZDCSMD_PsiWgtFull[n-1];
268 }
269 
270 //-------------------------------------------------------------
271 
272 UInt_t StFlowEvent::Mult(StFlowSelection* pFlowSelect) {
273  // Multiplicity of tracks selected for the event plane
274 
275  UInt_t mult = 0;
276 
277  StFlowTrackIterator itr;
278  for (itr = TrackCollection()->begin();
279  itr != TrackCollection()->end(); itr++) {
280  StFlowTrack* pFlowTrack = *itr;
281  if (pFlowSelect->Select(pFlowTrack)) mult++;
282  }
283 
284  return mult;
285 }
286 
287 //-------------------------------------------------------------
288 
289 UInt_t StFlowEvent::MultPart(StFlowSelection* pFlowSelect) {
290  // Multiplicity of tracks selected for correlation with the event plane
291 
292  UInt_t mult = 0;
293 
294  StFlowTrackIterator itr;
295  for (itr = TrackCollection()->begin();
296  itr != TrackCollection()->end(); itr++) {
297  StFlowTrack* pFlowTrack = *itr;
298  if (pFlowSelect->SelectPart(pFlowTrack)) mult++;
299  }
300 
301  return mult;
302 }
303 
304 //-------------------------------------------------------------
305 
306 Float_t StFlowEvent::MeanPt(StFlowSelection* pFlowSelect) {
307  // Mean pt of tracks selected for the event plane
308 
309  float sumPt = 0.;
310  UInt_t mult = 0;
311 
312  StFlowTrackIterator itr;
313  for (itr = TrackCollection()->begin();
314  itr != TrackCollection()->end(); itr++) {
315  StFlowTrack* pFlowTrack = *itr;
316  if (pFlowSelect->Select(pFlowTrack)) {
317  sumPt += pFlowTrack->Pt();
318  mult++;
319  }
320  }
321 
322  return (mult) ? sumPt/(float)mult : 0.;
323 }
324 
325 //-------------------------------------------------------------
326 
327 TVector2 StFlowEvent::Q(StFlowSelection* pFlowSelect) {
328  // Event plane vector
329 
330  TVector2 mQ, reCent;
331  Double_t mQx=0., mQy=0.;
332 
333  if (mUseZDCSMD) { // pFlowSelect is disabled; only 1st order Q generated
334  Float_t eXsum=0., eYsum=0., wXsum=0., wYsum=0.;
335  Float_t eXWgt=0., eYWgt=0., wXWgt=0., wYWgt=0.;
336  for (int v=1; v<8; v++) {
337  eXsum += ZDCSMD_GetPosition(0,0,v)*ZDCSMD(0,0,v);
338  wXsum += ZDCSMD_GetPosition(1,0,v)*ZDCSMD(1,0,v);
339  eXWgt += ZDCSMD(0,0,v);
340  wXWgt += ZDCSMD(1,0,v);
341  } //v
342  for (int h=1;h<9;h++) {
343  eYsum += ZDCSMD_GetPosition(0,1,h)*ZDCSMD(0,1,h);
344  wYsum += ZDCSMD_GetPosition(1,1,h)*ZDCSMD(1,1,h);
345  eYWgt += ZDCSMD(0,1,h);
346  wYWgt += ZDCSMD(1,1,h);
347  }
348  mQx= (eXWgt>0. && wXWgt>0. && eYWgt>0. && wYWgt>0.) ?
349  eXsum/eXWgt - wXsum/wXWgt:0.;
350  mQy= (eXWgt>0. && wXWgt>0. && eYWgt>0. && wYWgt>0.) ?
351  eYsum/eYWgt - wYsum/wYWgt:0.;
352  }//ZDCSMD
353  else { // ana
354  int selN = pFlowSelect->Sel();
355  int harN = pFlowSelect->Har();
356  double order = (double)(harN + 1);
357 
358  StFlowTrackIterator itr;
359  for (itr = TrackCollection()->begin();
360  itr != TrackCollection()->end(); itr++) {
361  StFlowTrack* pFlowTrack = *itr;
362  if (pFlowSelect->Select(pFlowTrack)) {
363  double phiWgt = PhiWeight(selN, harN, pFlowTrack);
364  float phi = pFlowTrack->Phi();
365  reCent = ReCentEP(selN, harN, pFlowTrack);
366  mQx += (phiWgt * cos(order * phi) - reCent.X());
367  mQy += (phiWgt * sin(order * phi) - reCent.Y());
368  }
369  }
370  }//ana
371 
372  mQ.Set(mQx, mQy);
373 
374  return mQ;
375 }
376 
377 //-------------------------------------------------------------
378 
379 TVector2 StFlowEvent::ReCentPar(StFlowSelection* pFlowSelect, const char* TPC) {
380  // For LYZ
381  // Calculate weighted recentering vector per particle for each TPC
382  // TPC can be "TPC", "TPCE", or "TPCW"
383  // for all particles that could be correlated with the event plane
384 
385  TVector2 mQ;
386  Double_t mQx=0., mQy=0., SumOfWeight=0.;
387 
388  int selN = pFlowSelect->Sel();
389  int harN = pFlowSelect->Har();
390  double order = (double)(harN + 1);
391  StTrackTopologyMap map;
392 
393  StFlowTrackIterator itr;
394  for (itr = TrackCollection()->begin();
395  itr != TrackCollection()->end(); itr++) {
396  StFlowTrack* pFlowTrack = *itr;
397  if (pFlowSelect->SelectPart(pFlowTrack)) {
398  map = pFlowTrack->TopologyMap();
399  if ((!strcmp(TPC,"TPCE") && map.trackFtpcEast()) ||
400  (!strcmp(TPC,"TPCW") && map.trackFtpcWest()) ||
401  (!strcmp(TPC,"TPC") && map.hasHitInDetector(kTpcId))) {
402  float phi = pFlowTrack->Phi();
403  double wgt = fabs(Weight(selN, harN, pFlowTrack));
404  mQx += wgt * cos(order * phi);
405  mQy += wgt * sin(order * phi);
406  SumOfWeight += wgt;
407  }
408  }
409  }
410 
411  if (SumOfWeight)
412  mQ.Set(mQx/SumOfWeight, mQy/SumOfWeight);
413  else mQ.Set(0.,0.);
414 
415  return mQ;
416 }
417 
418 //-------------------------------------------------------------
419 
420 TVector2 StFlowEvent::ReCentEPPar(StFlowSelection* pFlowSelect, const char* TPC) {
421  // For ana
422  // Calculate the recentering vector per particle for each TPC
423  // TPC can be "TPCE", "TPCW", "FTPCE", or "FTPCW"
424  // for particles that are used for the event plane
425 
426  TVector2 mQ;
427  Double_t mult=0., mQx=0., mQy=0.;
428 
429  int selN = pFlowSelect->Sel();
430  int harN = pFlowSelect->Har();
431  double order = (double)(harN + 1);
432  StTrackTopologyMap map;
433 
434  StFlowTrackIterator itr;
435  for (itr = TrackCollection()->begin();
436  itr != TrackCollection()->end(); itr++) {
437  StFlowTrack* pFlowTrack = *itr;
438  if (pFlowSelect->Select(pFlowTrack)) {
439  map = pFlowTrack->TopologyMap();
440  float eta = pFlowTrack->Eta();
441  if ((!strcmp(TPC,"FTPCE") && map.trackFtpcEast()) ||
442  (!strcmp(TPC,"FTPCW") && map.trackFtpcWest()) ||
443  (!strcmp(TPC,"TPCE") && eta < 0. && map.hasHitInDetector(kTpcId)) ||
444  (!strcmp(TPC,"TPCW") && eta > 0. && map.hasHitInDetector(kTpcId)) ) {
445  float phi = pFlowTrack->Phi();
446  double phiWgt = PhiWeight(selN, harN, pFlowTrack);
447  mQx += phiWgt * cos(order * phi);
448  mQy += phiWgt * sin(order * phi);
449  mult++;
450  }
451  }
452  }
453 
454  if (mult) { mQ.Set(mQx/mult, mQy/mult); }
455  else { mQ.Set(0.,0.); }
456 
457  return mQ;
458 }
459 
460 //-------------------------------------------------------------
461 
462 TVector2 StFlowEvent::ReCent(Int_t selN, Int_t harN, StFlowTrack* pFlowTrack) const {
463  // Get TVector2 for recentering in LYZ makers.
464 
465  TVector2 reCent;
466  Double_t reCentX, reCentY;
467  StTrackTopologyMap map = pFlowTrack->TopologyMap();
468 
469  if (map.hasHitInDetector(kTpcId)) {
470  reCentX = mReCentX[selN][harN][2];
471  reCentY = mReCentY[selN][harN][2];
472  } else if (map.trackFtpcEast()) {
473  reCentX = mReCentX[selN][harN][0];
474  reCentY = mReCentY[selN][harN][0];
475  } else if (map.trackFtpcWest()) {
476  reCentX = mReCentX[selN][harN][1];
477  reCentY = mReCentY[selN][harN][1];
478  } else {
479  reCentX = 0.;
480  reCentY = 0.;
481  }
482 
483  reCent.Set(reCentX, reCentY);
484  return reCent;
485 }
486 
487 //-------------------------------------------------------------
488 
489 TVector2 StFlowEvent::ReCentEP(Int_t selN, Int_t harN, StFlowTrack* pFlowTrack) const {
490  // Get TVector2 for recentering in ana makers.
491 
492  TVector2 reCent;
493  Double_t reCentX, reCentY;
494  StTrackTopologyMap map = pFlowTrack->TopologyMap();
495 
496  if (map.hasHitInDetector(kTpcId)) {
497  float eta = pFlowTrack->Eta();
498  if (eta > 0.) { // TPCW
499  reCentX = mReCentX[selN][harN][3];
500  reCentY = mReCentY[selN][harN][3];
501  } else { // TPCE
502  reCentX = mReCentX[selN][harN][2];
503  reCentY = mReCentY[selN][harN][2];
504  }
505  } else if (map.trackFtpcEast()) { // FTPCE
506  reCentX = mReCentX[selN][harN][0];
507  reCentY = mReCentY[selN][harN][0];
508  } else if (map.trackFtpcWest()) { // FTPCW
509  reCentX = mReCentX[selN][harN][1];
510  reCentY = mReCentY[selN][harN][1];
511  } else {
512  reCentX = 0.;
513  reCentY = 0.;
514  }
515 
516  reCent.Set(reCentX, reCentY);
517  return reCent;
518 }
519 
520 //-------------------------------------------------------------
521 
522 TVector2 StFlowEvent::QPart(StFlowSelection* pFlowSelect) {
523  // Event plane vector for LYZ method for all particles that could be correlated with the event plane
524 
525  TVector2 reCent, mQ;
526  Double_t mQx=0., mQy=0.;
527 
528  int selN = pFlowSelect->Sel();
529  int harN = pFlowSelect->Har();
530  double order = (double)(harN + 1);
531 
532  StFlowTrackIterator itr;
533  for (itr = TrackCollection()->begin();
534  itr != TrackCollection()->end(); itr++) {
535  StFlowTrack* pFlowTrack = *itr;
536  if (pFlowSelect->SelectPart(pFlowTrack)) {
537  double phiWgt = PhiWeight(selN, harN, pFlowTrack);
538  float phi = pFlowTrack->Phi();
539  reCent = ReCent(selN, harN, pFlowTrack);
540  mQx += phiWgt * (cos(order * phi) - reCent.X());
541  mQy += phiWgt * (sin(order * phi) - reCent.Y());
542  }
543  }
544 
545  mQ.Set(mQx, mQy);
546  return mQ;
547 }
548 
549 //-------------------------------------------------------------
550 
551 TVector2 StFlowEvent::NormQ(StFlowSelection* pFlowSelect) {
552  // Return normalized Q = Q / ::sqrt(sum of weights**2)
553  // Where is this used?
554 
555  TVector2 mQ, reCent;
556  Double_t mQx=0., mQy=0.;
557  int selN = pFlowSelect->Sel();
558  int harN = pFlowSelect->Har();
559  double order = (double)(harN + 1);
560  double SumOfWeightSqr = 0;
561 
562 
563  StFlowTrackIterator itr;
564  for (itr = TrackCollection()->begin();
565  itr != TrackCollection()->end(); itr++) {
566  StFlowTrack* pFlowTrack = *itr;
567  if (pFlowSelect->Select(pFlowTrack)) {
568 
569  double phiWgt = PhiWeight(selN, harN, pFlowTrack);
570  SumOfWeightSqr += phiWgt*phiWgt;
571 
572  float phi = pFlowTrack->Phi();
573  reCent = ReCent(selN, harN, pFlowTrack);
574  mQx += phiWgt * (cos(order * phi) - reCent.X());
575  mQy += phiWgt * (sin(order * phi) - reCent.Y());
576  }
577  }
578 
579  if (SumOfWeightSqr)
580  mQ.Set(mQx/::sqrt(SumOfWeightSqr), mQy/::sqrt(SumOfWeightSqr));
581  else mQ.Set(0.,0.);
582 
583  return mQ;
584 }
585 
586 //-------------------------------------------------------------
587 
588 Float_t StFlowEvent::q(StFlowSelection* pFlowSelect) {
589  // Magnitude of normalized Q vector without pt or eta weighting
590 
591  TVector2 mQ, reCent;
592  Double_t mQx=0., mQy=0.;
593  int selN = pFlowSelect->Sel();
594  int harN = pFlowSelect->Har();
595  double order = (double)(harN + 1);
596  double SumOfWeightSqr = 0;
597 
598 
599  StFlowTrackIterator itr;
600  for (itr = TrackCollection()->begin();
601  itr != TrackCollection()->end(); itr++) {
602  StFlowTrack* pFlowTrack = *itr;
603  if (pFlowSelect->Select(pFlowTrack)) { // use event plane selections
604 
605  double phiWgt = PhiWeightRaw(selN, harN, pFlowTrack); // Raw
606  SumOfWeightSqr += phiWgt*phiWgt;
607 
608  float phi = pFlowTrack->Phi();
609  reCent = ReCent(selN, harN, pFlowTrack);
610  mQx += phiWgt * (cos(order * phi) - reCent.X());
611  mQy += phiWgt * (sin(order * phi) - reCent.Y());
612  }
613  }
614 
615  if (SumOfWeightSqr)
616  mQ.Set(mQx/::sqrt(SumOfWeightSqr), mQy/::sqrt(SumOfWeightSqr));
617  else mQ.Set(0.,0.);
618 
619  return mQ.Mod();
620 }
621 
622 //-------------------------------------------------------------
623 
624 Double_t StFlowEvent::SumWeightSquare(StFlowSelection* pFlowSelect) {
625 
626  // Return sum of weights**2
627 
628  // This method is called only by cumulant method. It won't return the true
629  // SumWeightSquare if called in standard method, in which the raw phi weight
630  // is applied. -TAH
631 
632  int selN = pFlowSelect->Sel();
633  int harN = pFlowSelect->Har();
634  Double_t SumOfWeightSqr = 0;
635 
636  StFlowTrackIterator itr;
637  for (itr = TrackCollection()->begin();
638  itr != TrackCollection()->end(); itr++) {
639  StFlowTrack* pFlowTrack = *itr;
640  if (pFlowSelect->Select(pFlowTrack)) {
641 
642  double phiWgt = Weight(selN, harN, pFlowTrack);
643  SumOfWeightSqr += phiWgt*phiWgt;
644  }
645  }
646 
647  if (SumOfWeightSqr < 0.) return Mult(pFlowSelect);
648 
649  return SumOfWeightSqr;
650 }
651 
652 //-------------------------------------------------------------
653 
654 Float_t StFlowEvent::Qtheta(StFlowSelection* pFlowSelect, Float_t theta) {
655  // Q^{\theta} for LeeYangZeros method for all particles that could be correlated with the event plane
656  // BP Eq. 3 (Nucl. Phys. A 727, 373 (2003))
657 
658  Float_t Qtheta = 0.;
659 
660  int selN = pFlowSelect->Sel();
661  int harN = pFlowSelect->Har();
662  double order = (double)(harN + 1);
663  double reCentTheta;
664  TVector2 reCent;
665 
666  StFlowTrackIterator itr;
667  for (itr = TrackCollection()->begin();
668  itr != TrackCollection()->end(); itr++) {
669  StFlowTrack* pFlowTrack = *itr;
670  if (pFlowSelect->SelectPart(pFlowTrack)) {
671  double wgt = Weight(selN, harN, pFlowTrack);
672  float phi = pFlowTrack->Phi();
673  reCent = ReCent(selN, harN, pFlowTrack);
674  reCentTheta = reCent.X() * cos(order*theta) + reCent.Y() * sin(order*theta);
675  Qtheta += wgt * (cos(order * (phi - theta)) - reCentTheta);
676  }
677  }
678 
679  return Qtheta;
680 }
681 
682 //-------------------------------------------------------------
683 
684 TComplex StFlowEvent::Grtheta(StFlowSelection* pFlowSelect, Float_t r, Float_t theta) {
685  // Product Generating Function for LeeYangZeros method
686  // PG Eq. 3 (J. Phys. G Nucl. Part. Phys 30 S1213 (2004))
687 
688  TComplex G = TComplex::One();
689  int selN = pFlowSelect->Sel();
690  int harN = pFlowSelect->Har();
691  double order = (double)(harN + 1);
692  double reCentTheta;
693  TVector2 reCent;
694 
695  StFlowTrackIterator itr;
696  for (itr = TrackCollection()->begin();
697  itr != TrackCollection()->end(); itr++) {
698  StFlowTrack* pFlowTrack = *itr;
699  if (pFlowSelect->SelectPart(pFlowTrack)) {
700  double wgt = Weight(selN, harN, pFlowTrack);
701  float phi = pFlowTrack->Phi();
702  reCent = ReCent(selN, harN, pFlowTrack);
703  reCentTheta = reCent.X() * cos(order*theta) + reCent.Y() * sin(order*theta);
704  double Gim = r * wgt * (cos(order * (phi - theta)) - reCentTheta);
705  TComplex G_i(1., Gim);
706  G *= G_i;
707  }
708  }
709 
710  return G;
711 }
712 
713 //-------------------------------------------------------------
714 
715 
716 TComplex StFlowEvent::GV1r0theta(StFlowSelection* pFlowSelect, Float_t r0, Float_t theta1, Float_t theta) {
717  // Product Generating Function for LeeYangZeros method for v1 mixed harmonics
718  // DF Eq. 1
719 
720  TComplex G = TComplex::One();
721  double reCentTheta;
722  TVector2 reCent;
723 
724  StFlowTrackIterator itr;
725  for (itr = TrackCollection()->begin();
726  itr != TrackCollection()->end(); itr++) {
727  StFlowTrack* pFlowTrack = *itr;
728  if (pFlowSelect->SelectPart(pFlowTrack)) {
729  double wgt1 = Weight(1, 0, pFlowTrack); // selection 2, harmonic 1
730  double wgt2 = Weight(1, 1, pFlowTrack); // selection 2, harmonic 2
731  float phi = pFlowTrack->Phi();
732  //reCent = ReCent(1, 0, pFlowTrack); // selection 2, harmonic 1
733  reCent.Set(0.,0.);
734  reCentTheta = reCent.X() * cos(1.*theta1) + reCent.Y() * sin(1.*theta1);
735  double Gim1 = r0 * Flow::epsV1 * wgt1 * (cos(phi - theta1) - reCentTheta);
736  //reCent = ReCent(1, 1, pFlowTrack); // selection 2, harmonic 2
737  reCent.Set(0.,0.); // NOT recentered
738  reCentTheta = reCent.X() * cos(2.*theta) + reCent.Y() * sin(2.*theta);
739  double Gim2 = r0 * wgt2 * (cos(2*(phi - theta)) - reCentTheta);
740  TComplex G_i(1., Gim1+Gim2);
741  G *= G_i;
742  }
743  }
744 
745  return G;
746 }
747 
748 //-------------------------------------------------------------
749 TComplex StFlowEvent::Gder_r0theta(StFlowSelection* pFlowSelect, Float_t r0, Float_t theta) {
750  // Sum for the denominator for diff. flow for the Product Generating Function for LeeYangZeros method
751  // PG Eq. 9 (J. Phys. G Nucl. Part. Phys 30 S1213 (2004))
752  // Also for v1 mixed harmonics: DF Eq. 5
753  // It is the deriverative of Grtheta at r0 divided by Grtheta at r0
754 
755  TComplex Gder(0.,0.);
756  int selN = pFlowSelect->Sel();
757  int harN = pFlowSelect->Har();
758  double order = (double)(harN + 1);
759  double reCentTheta;
760  TVector2 reCent;
761 
762  StFlowTrackIterator itr;
763  for (itr = TrackCollection()->begin();
764  itr != TrackCollection()->end(); itr++) {
765  StFlowTrack* pFlowTrack = *itr;
766  if (pFlowSelect->SelectPart(pFlowTrack)) {
767  double wgt = Weight(selN, harN, pFlowTrack);
768  float phi = pFlowTrack->Phi();
769  reCent = ReCent(selN, harN, pFlowTrack);
770  reCentTheta = reCent.X() * cos(order*theta) + reCent.Y() * sin(order*theta);
771  double cosTerm = wgt * (cos(order * (phi - theta)) - reCentTheta);
772  TComplex denom(1., r0*cosTerm);
773  Gder += (cosTerm / denom);
774  }
775  }
776 
777  return Gder;
778 }
779 
780 //-------------------------------------------------------------
781 
782 Float_t StFlowEvent::Psi(StFlowSelection* pFlowSelect) {
783  // Event plane angle
784 
785  int harN = pFlowSelect->Har();
786  float order = (float)(harN + 1);
787  Float_t psi = 0.;
788 
789  TVector2 mQ = Q(pFlowSelect);
790 
791  if (mQ.Mod()) {
792  psi= mQ.Phi() / order;
793  if (psi < 0.) { psi += twopi / order; }
794  }
795 
796  return psi;
797 }
798 
799 //-------------------------------------------------------------
800 
801 Float_t StFlowEvent::ZDCSMD_PsiCorr() {
802  //difference between psi's from east and west ZDCSMD
803 
804  Float_t psi_e = ZDCSMD_PsiEst();
805  Float_t psi_w = ZDCSMD_PsiWst();
806  Float_t psi_w_shift = psi_w + twopi/2.;
807 
808  if(psi_w_shift > twopi/2.) psi_w_shift -= twopi;
809  if(psi_w_shift < -twopi/2.) psi_w_shift += twopi;
810 
811  Float_t psi_delta = psi_e - psi_w_shift;
812 
813  if(psi_delta > twopi/2.) psi_delta -= twopi;
814  if(psi_delta < -twopi/2.) psi_delta += twopi;
815 
816  return psi_delta;
817 }
818 
819 //----------------------------------------------------------------------
820 
821 Float_t StFlowEvent::ZDCSMD_PsiEst() {
822  //psi angle from east ZDCSMD
823 
824  Float_t eXsum=0., eYsum=0., eXWgt=0., eYWgt=0.;
825 
826  for(int v=1; v<8; v++) {
827  eXsum += ZDCSMD_GetPosition(0,0,v)*ZDCSMD(0,0,v);
828  eXWgt += ZDCSMD(0,0,v);
829  }//v
830  for(int h=1; h<9; h++) {
831  eYsum += ZDCSMD_GetPosition(0,1,h)*ZDCSMD(0,1,h);
832  eYWgt += ZDCSMD(0,1,h);
833  }//h
834 
835  Float_t psi_e = atan2((eYWgt>0.) ? eYsum/eYWgt:0.,(eXWgt>0.) ? eXsum/eXWgt:0.);
836 
837  return psi_e;
838 }
839 
840 //----------------------------------------------------------------------
841 
842 Float_t StFlowEvent::ZDCSMD_PsiWst() {
843  //psi angle from west ZDCSMD
844 
845  Float_t wXsum=0.,wYsum=0.,wXWgt=0.,wYWgt=0.;
846 
847  for(int v=1;v<8;v++) {
848  wXsum += ZDCSMD_GetPosition(1,0,v)*ZDCSMD(1,0,v);
849  wXWgt += ZDCSMD(1,0,v);
850  }//v
851  for(int h=1;h<9;h++) {
852  wYsum += ZDCSMD_GetPosition(1,1,h)*ZDCSMD(1,1,h);
853  wYWgt += ZDCSMD(1,1,h);
854  }//h
855 
856  Float_t psi_w = atan2((wYWgt>0.) ? wYsum/wYWgt:0.,(wXWgt>0.) ? wXsum/wXWgt:0.);
857 
858  return psi_w;
859 }
860 
861 //-----------------------------------------------------------------------
862 
863 Float_t StFlowEvent::ZDCSMD_GetPosition(int eastwest,int verthori,int strip) {
864  // Get position of each slat;strip starts from 1
865 
866  Float_t zdcsmd_x[7] = {0.5,2,3.5,5,6.5,8,9.5};
867  Float_t zdcsmd_y[8] = {1.25,3.25,5.25,7.25,9.25,11.25,13.25,15.25};
868 
869  if(eastwest==0 && verthori==0) return zdcsmd_x[strip-1]-mZDCSMDCenterex;
870  if(eastwest==1 && verthori==0) return mZDCSMDCenterwx-zdcsmd_x[strip-1];
871  if(eastwest==0 && verthori==1) return zdcsmd_y[strip-1]/sqrt(2.)-mZDCSMDCenterey;
872  if(eastwest==1 && verthori==1) return zdcsmd_y[strip-1]/sqrt(2.)-mZDCSMDCenterwy;
873 
874  return 0;
875  }
876 
877 //-----------------------------------------------------------------------
878 
879 Double_t StFlowEvent::G_New(StFlowSelection* pFlowSelect, Double_t Zx, Double_t Zy) {
880  // Generating function for the new cumulant method. Eq. 3 in the Practical Guide.
881 
882  int selN = pFlowSelect->Sel();
883  int harN = pFlowSelect->Har();
884  double order = (double)(harN + 1);
885 
886  double mult = (double)Mult(pFlowSelect);
887  Double_t theG = 1.;
888 
889  StFlowTrackIterator itr;
890  for (itr = TrackCollection()->begin();
891  itr != TrackCollection()->end(); itr++) {
892  StFlowTrack* pFlowTrack = *itr;
893  if (pFlowSelect->Select(pFlowTrack)) {
894 
895  double phiWgt = Weight(selN, harN, pFlowTrack);
896  float phi = pFlowTrack->Phi();
897  theG *= (1. + (phiWgt/mult) * (2.* Zx * cos(order * phi) +
898  2.* Zy * sin(order * phi) ) );
899  }
900  }
901 
902  return theG;
903 }
904 
905 //-----------------------------------------------------------------------
906 
907 Double_t StFlowEvent::G_Mix(StFlowSelection* pFlowSelect, Double_t Z1x, Double_t Z1y, Double_t Z2x, Double_t Z2y) {
908 
909  //only make sense in v1{3} calculation.
910  //which means, one should call this function with StFlowSelection with harN==0 only !
911 
912  int selN = pFlowSelect->Sel();
913  int harN = pFlowSelect->Har();
914  double order = (double)(harN + 1); //should be 1 always in this func.
915 
916  bool oddHar = (harN+1) % 2;
917  double etaWgt=1.;
918  double ptWgt=1.;
919 
920  double mult = (double)Mult(pFlowSelect);
921  Double_t theG = 1.;
922 
923  StFlowTrackIterator itr;
924  for (itr = TrackCollection()->begin();
925  itr != TrackCollection()->end(); itr++) {
926  StFlowTrack* pFlowTrack = *itr;
927  if (pFlowSelect->Select(pFlowTrack)) {
928 
929  double detectorV1Wgt = 1.;
930 
931  if (pFlowTrack->TopologyMap().hasHitInDetector(kTpcId) ||
932  (pFlowTrack->TopologyMap().data(0) == 0 &&
933  pFlowTrack->TopologyMap().data(1) == 0)) {
934  // hits in Tpc or TopologyMap not available
935  detectorV1Wgt =mV1TPCDetctWgtG_Mix[selN];
936  } else if (pFlowTrack->TopologyMap().trackFtpcEast() ) {
937  detectorV1Wgt =mV1FtpcEastDetctWgtG_Mix[selN];
938  } else if (pFlowTrack->TopologyMap().trackFtpcWest() ) {
939  detectorV1Wgt =mV1FtpcWestDetctWgtG_Mix[selN];
940  }
941 
942  double detectorV2Wgt = 1.;
943 
944  if (pFlowTrack->TopologyMap().hasHitInDetector(kTpcId) ||
945  (pFlowTrack->TopologyMap().data(0) == 0 &&
946  pFlowTrack->TopologyMap().data(1) == 0)) {
947  // hits in Tpc or TopologyMap not available
948  detectorV2Wgt =mV2TPCDetctWgtG_Mix[selN];
949  } else if (pFlowTrack->TopologyMap().trackFtpcEast() ) {
950  detectorV2Wgt =mV2FtpcEastDetctWgtG_Mix[selN];
951  } else if (pFlowTrack->TopologyMap().trackFtpcWest() ) {
952  detectorV2Wgt =mV2FtpcWestDetctWgtG_Mix[selN];
953  }
954 
955  double phiWgt = 1.;//no need raw phi weight for cumulant method.
956  float phi = pFlowTrack->Phi();
957  double pt = pFlowTrack->Pt();
958  double eta = pFlowTrack->Eta();
959 
960  if (oddHar) etaWgt =( (eta>0) ? (EtaAbsWgtValue(eta)) : (-1.*EtaAbsWgtValue(eta)) );
961 
962  ptWgt = PtAbsWgtValue(pt);
963 
964  theG *=
965  (1. + (phiWgt*etaWgt*detectorV1Wgt/mult) * (2.* Z1x * cos(order * phi) +
966  2.* Z1y * sin(order * phi) )
967  + (phiWgt*ptWgt*detectorV2Wgt/mult) * (2.* Z2x * cos(phi * order*2.) +
968  2.* Z2y * sin(phi * order*2.) ) );
969 
970  }
971  }
972 
973  return theG;
974 }
975 
976 //-------------------------------------------------------------
977 void StFlowEvent::SetSelections() {
978  // for particles used for the event plane
979 
980  StFlowTrackIterator itr;
981  for (itr = TrackCollection()->begin();
982  itr != TrackCollection()->end(); itr++) {
983  StFlowTrack* pFlowTrack = *itr;
984  double eta = (double)(pFlowTrack->Eta());
985  float Pt = pFlowTrack->Pt();
986  StTrackTopologyMap map = pFlowTrack->TopologyMap();
987 
988  // PID
989  if (mPid[0] != '\0') {
990  if (strstr(mPid, "h")!=0) {
991  int charge = pFlowTrack->Charge();
992  if (strcmp("h+", mPid)==0 && charge != 1) continue;
993  if (strcmp("h-", mPid)==0 && charge != -1) continue;
994  } else {
995  Char_t pid[10];
996  strcpy(pid, pFlowTrack->Pid());
997  if (strstr(pid, mPid)==0) continue;
998  }
999  }
1000 
1001  // Global DCA
1002  float gDca = pFlowTrack->DcaGlobal();
1003 
1004  if (map.hasHitInDetector(kTpcId) || (map.data(0) == 0 && map.data(1) == 0)) {
1005  // Tpc track or TopologyMap not available)
1006 
1007  if (mDcaGlobalTpcCuts[1] > mDcaGlobalTpcCuts[0] &&
1008  (gDca < mDcaGlobalTpcCuts[0] || gDca >= mDcaGlobalTpcCuts[1])) continue;
1009  }
1010  else if (map.trackFtpcEast() || map.trackFtpcWest()) {
1011  // Ftpc track
1012 
1013  if (mDcaGlobalFtpcCuts[1] > mDcaGlobalFtpcCuts[0] &&
1014  (gDca < mDcaGlobalFtpcCuts[0] || gDca >= mDcaGlobalFtpcCuts[1])) continue;
1015  }
1016 
1017  for (int selN = 0; selN < Flow::nSels; selN++) {
1018  for (int harN = 0; harN < Flow::nHars; harN++) {
1019 
1020  if (map.hasHitInDetector(kTpcId) || (map.data(0) == 0 && map.data(1) == 0)) {
1021  // Tpc track or TopologyMap not available
1022 
1023  // Eta
1024 // if (mEtaTpcCuts[1][harN%2][selN] > mEtaTpcCuts[0][harN%2][selN] &&
1025 // (fabs(eta) < mEtaTpcCuts[0][harN%2][selN] ||
1026 // fabs(eta) >= mEtaTpcCuts[1][harN%2][selN])) continue;
1027  int etaCut = harN ? 1 : 0;
1028  if (mEtaTpcCuts[1][etaCut][selN] > mEtaTpcCuts[0][etaCut][selN] &&
1029  (fabs(eta) < mEtaTpcCuts[0][etaCut][selN] ||
1030  fabs(eta) >= mEtaTpcCuts[1][etaCut][selN])) continue;
1031 
1032  // Pt
1033  if (mPtTpcCuts[1][harN%2][selN] > mPtTpcCuts[0][harN%2][selN] &&
1034  (Pt < mPtTpcCuts[0][harN%2][selN] ||
1035  Pt >= mPtTpcCuts[1][harN%2][selN])) continue;
1036  }
1037  else if (map.trackFtpcEast() || map.trackFtpcWest()) {
1038  // Ftpc track
1039 
1040  // Eta
1041  if (eta < 0.) {
1042  if (mEtaFtpcCuts[1][harN%2][selN] > mEtaFtpcCuts[0][harN%2][selN] &&
1043  (eta < mEtaFtpcCuts[0][harN%2][selN] ||
1044  eta >= mEtaFtpcCuts[1][harN%2][selN])) continue;
1045  }
1046  else { // eta > 0.
1047  if (mEtaFtpcCuts[3][harN%2][selN] > mEtaFtpcCuts[2][harN%2][selN] &&
1048  (eta < mEtaFtpcCuts[2][harN%2][selN] ||
1049  eta >= mEtaFtpcCuts[3][harN%2][selN])) continue;
1050  }
1051 
1052  // Pt
1053  if (mPtFtpcCuts[1][harN%2][selN] > mPtFtpcCuts[0][harN%2][selN] &&
1054  (Pt < mPtFtpcCuts[0][harN%2][selN] ||
1055  Pt >= mPtFtpcCuts[1][harN%2][selN])) continue;
1056  }
1057 
1058  pFlowTrack->SetSelect(harN, selN);
1059  }
1060  }
1061  }
1062 }
1063 
1064 //-------------------------------------------------------------
1065 
1066 void StFlowEvent::MakeSubEvents() {
1067  // Make random subevents
1068 
1069  StFlowTrackIterator itr;
1070  int eventMult[Flow::nHars][Flow::nSels] = {{0}};
1071  int harN, selN, subN = 0;
1072 
1073  // loop to count the total number of tracks for each selection
1074  for (itr = TrackCollection()->begin();
1075  itr != TrackCollection()->end(); itr++) {
1076  StFlowTrack* pFlowTrack = *itr;
1077  for (selN = 0; selN < Flow::nSels; selN++) {
1078  for (harN = 0; harN < Flow::nHars; harN++) {
1079  if (pFlowTrack->Select(harN, selN)) {
1080  eventMult[harN][selN]++;
1081  }
1082  }
1083  }
1084  }
1085 
1086  // loop to set the SubEvent member variable
1087  for (selN = 0; selN < Flow::nSels; selN++) {
1088  for (harN = 0; harN < Flow::nHars; harN++) {
1089  int subEventMult = eventMult[harN][selN] / Flow::nSubs;
1090  if (subEventMult) {
1091  subN = 0;
1092  int countN = 0;
1093  for (itr = TrackCollection()->begin();
1094  itr != TrackCollection()->end(); itr++) {
1095  StFlowTrack* pFlowTrack = *itr;
1096  if (pFlowTrack->Select(harN, selN)) {
1097  pFlowTrack->SetSubevent(harN, selN, subN);
1098  countN++;
1099  if (countN % subEventMult == 0.) subN++;
1100  }
1101  }
1102  }
1103  }
1104  }
1105 
1106 }
1107 
1108 //-------------------------------------------------------------
1109 
1110 void StFlowEvent::MakeEtaSubEvents() {
1111  // Make subevents with positive and negative eta
1112 
1113  StFlowTrackIterator itr;
1114  int harN, selN = 0;
1115 
1116  // loop to set the SubEvent member variable
1117  for (selN = 0; selN < Flow::nSels; selN++) {
1118  for (harN = 0; harN < Flow::nHars; harN++) {
1119  for (itr = TrackCollection()->begin();
1120  itr != TrackCollection()->end(); itr++) {
1121  StFlowTrack* pFlowTrack = *itr;
1122  if (pFlowTrack->Select(harN, selN)) {
1123  float eta = pFlowTrack->Eta();
1124  if (eta > 0.) {
1125  pFlowTrack->SetSubevent(harN, selN, 0);
1126  } else {
1127  pFlowTrack->SetSubevent(harN, selN, 1);
1128  }
1129  }
1130  }
1131  }
1132  }
1133 
1134 }
1135 
1136 //-----------------------------------------------------------------------
1137 
1138 void StFlowEvent::SetPidsDeviant() {
1139  // Set PID with PID Deviant method
1140 
1141  StFlowTrackIterator itr;
1142 
1143  for (itr = TrackCollection()->begin();
1144  itr != TrackCollection()->end(); itr++) {
1145 
1146  StFlowTrack* pFlowTrack = *itr;
1147  Char_t pid[10] = "NA";
1148  Short_t charge = pFlowTrack->Charge();
1149 
1150  bool bPiPlus = kFALSE;
1151  bool bPiMinus = kFALSE;
1152  bool bProton = kFALSE;
1153  bool bAntiProton = kFALSE;
1154  bool bKplus = kFALSE;
1155  bool bKminus = kFALSE;
1156  bool bDeuteron = kFALSE;
1157  bool bAntiDeuteron = kFALSE;
1158  bool bElectron = kFALSE;
1159  bool bPositron = kFALSE;
1160 
1161  if (charge == 1) {
1162  float piPlus = pFlowTrack->PidPiPlus();
1163  float proton = pFlowTrack->PidProton();
1164  float kPlus = pFlowTrack->PidKaonPlus();
1165  float deuteron = pFlowTrack->PidDeuteron();
1166  float positron = pFlowTrack->PidPositron();
1167  if (piPlus > mPiPlusCuts[0] &&
1168  piPlus < mPiPlusCuts[1]) {
1169  bPiPlus = kTRUE;
1170  }
1171  if ( proton > mProtonCuts[0] &&
1172  proton < mProtonCuts[1]) {
1173  bProton = kTRUE;
1174  }
1175  if ( kPlus > mKPlusCuts[0] &&
1176  kPlus < mKPlusCuts[1]) {
1177  bKplus = kTRUE;
1178  }
1179  if ( deuteron > mDeuteronCuts[0] &&
1180  deuteron < mDeuteronCuts[1]) {
1181  bDeuteron = kTRUE;
1182  }
1183  if ( positron > mPositronCuts[0] &&
1184  positron < mPositronCuts[1]) {
1185  bPositron = kTRUE;
1186  }
1187  } else if (charge == -1) {
1188  float piMinus = pFlowTrack->PidPiMinus();
1189  float antiProton = pFlowTrack->PidAntiProton();
1190  float kMinus = pFlowTrack->PidKaonMinus();
1191  float antiDeuteron = pFlowTrack->PidAntiDeuteron();
1192  float electron = pFlowTrack->PidElectron();
1193  if (piMinus > mPiMinusCuts[0] &&
1194  piMinus < mPiMinusCuts[1]) {
1195  bPiMinus = kTRUE;
1196  }
1197  if ( antiProton > mAntiProtonCuts[0] &&
1198  antiProton < mAntiProtonCuts[1]) {
1199  bAntiProton = kTRUE;
1200  }
1201  if ( kMinus > mKMinusCuts[0] &&
1202  kMinus < mKMinusCuts[1]) {
1203  bKminus = kTRUE;
1204  }
1205  if ( antiDeuteron > mAntiDeuteronCuts[0] &&
1206  antiDeuteron < mAntiDeuteronCuts[1]) {
1207  bAntiDeuteron = kTRUE;
1208  }
1209  if ( electron > mElectronCuts[0] &&
1210  electron < mElectronCuts[1]) {
1211  bElectron = kTRUE;
1212  }
1213  }
1214 
1215  if (bPiPlus && !bPiMinus && !bProton && !bAntiProton &&
1216  !bKplus && !bKminus && !bDeuteron && !bAntiDeuteron &&
1217  !bElectron && !bPositron) { strcpy(pid, "pi+"); }
1218  if (!bPiPlus && bPiMinus && !bProton && !bAntiProton &&
1219  !bKplus && !bKminus && !bDeuteron && !bAntiDeuteron &&
1220  !bElectron && !bPositron) { strcpy(pid, "pi-"); }
1221  if (!bPiPlus && !bPiMinus && bProton && !bAntiProton &&
1222  !bKplus && !bKminus && !bDeuteron && !bAntiDeuteron &&
1223  !bElectron && !bPositron) { strcpy(pid, "pr+"); }
1224  if (!bPiPlus && !bPiMinus && !bProton && bAntiProton &&
1225  !bKplus && !bKminus && !bDeuteron && !bAntiDeuteron &&
1226  !bElectron && !bPositron) { strcpy(pid, "pr-"); }
1227  if (!bPiPlus && !bPiMinus && !bProton && !bAntiProton &&
1228  bKplus && !bKminus && !bDeuteron && !bAntiDeuteron &&
1229  !bElectron && !bPositron) { strcpy(pid, "k+"); }
1230  if (!bPiPlus && !bPiMinus && !bProton && !bAntiProton &&
1231  !bKplus && bKminus && !bDeuteron && !bAntiDeuteron &&
1232  !bElectron && !bPositron) { strcpy(pid, "k-"); }
1233  if (!bPiPlus && !bPiMinus && !bProton && !bAntiProton &&
1234  !bKplus && !bKminus && bDeuteron && !bAntiDeuteron &&
1235  !bElectron && !bPositron) { strcpy(pid, "d+"); }
1236  if (!bPiPlus && !bPiMinus && !bProton && !bAntiProton &&
1237  !bKplus && !bKminus && !bDeuteron && bAntiDeuteron &&
1238  !bElectron && !bPositron) { strcpy(pid, "d-"); }
1239  if (!bPiPlus && !bPiMinus && !bProton && !bAntiProton &&
1240  !bKplus && !bKminus && !bDeuteron && !bAntiDeuteron &&
1241  bElectron && !bPositron) { strcpy(pid, "e-"); }
1242  if (!bPiPlus && !bPiMinus && !bProton && !bAntiProton &&
1243  !bKplus && !bKminus && !bDeuteron && !bAntiDeuteron &&
1244  !bElectron && bPositron) { strcpy(pid, "e+"); }
1245 
1246  pFlowTrack->SetPid(pid);
1247 
1248  }
1249 }
1250 
1251 //-----------------------------------------------------------------------
1252 
1253 void StFlowEvent::SetPidsProb() {
1254  // Set PID with PID Probability method
1255 
1256  StFlowTrackIterator itr;
1257 
1258  for (itr = TrackCollection()->begin();
1259  itr != TrackCollection()->end(); itr++) {
1260 
1261  StFlowTrack* pFlowTrack = *itr;
1262  Char_t pid[10] = "NA";
1263 
1264  if (pFlowTrack->MostLikelihoodPID() == 8 &&
1265  pFlowTrack->MostLikelihoodProb() > 0.9)
1266  { strcpy(pid, "pi+"); }
1267  if (pFlowTrack->MostLikelihoodPID() == 9 &&
1268  pFlowTrack->MostLikelihoodProb() > 0.9)
1269  { strcpy(pid, "pi-"); }
1270  if (pFlowTrack->MostLikelihoodPID() == 14 &&
1271  pFlowTrack->MostLikelihoodProb() > 0.9)
1272  { strcpy(pid, "pr+"); }
1273  if (pFlowTrack->MostLikelihoodPID() == 15 &&
1274  pFlowTrack->MostLikelihoodProb() > 0.9)
1275  { strcpy(pid, "pr-"); }
1276  if (pFlowTrack->MostLikelihoodPID() == 11 &&
1277  pFlowTrack->MostLikelihoodProb() > 0.9)
1278  { strcpy(pid, "k+"); }
1279  if (pFlowTrack->MostLikelihoodPID() == 12 &&
1280  pFlowTrack->MostLikelihoodProb() > 0.9)
1281  { strcpy(pid, "k-"); }
1282  if (pFlowTrack->MostLikelihoodPID() == 45 &&
1283  pFlowTrack->MostLikelihoodProb() > 0.9)
1284  { strcpy(pid, "d+"); }
1285 // if (pFlowTrack->MostLikelihoodPID() == &&
1286 // pFlowTrack->MostLikelihoodProb() > 0.9)
1287 // { strcpy(pid, "d-"); }
1288  if (pFlowTrack->MostLikelihoodPID() == 3 &&
1289  pFlowTrack->MostLikelihoodProb() > 0.9)
1290  { strcpy(pid, "e-"); }
1291  if (pFlowTrack->MostLikelihoodPID() == 2 &&
1292  pFlowTrack->MostLikelihoodProb() > 0.9)
1293  { strcpy(pid, "e+"); }
1294 
1295  pFlowTrack->SetPid(pid);
1296 
1297  }
1298 }
1299 
1300 //-----------------------------------------------------------------------
1301 
1302 void StFlowEvent::SetCentrality() {
1303  // Centrality=0 is not retrieveable
1304 
1305  Int_t* cent = 0;
1306  Int_t tracks = mMultEta; // converts UInt_t to Int_t
1307 
1308  if (mRunID > 8000000) { // year 7
1309  cent = Flow::cent200Year7;
1310  } else if (mRunID > 4000000) { // year 4
1311  if (mCenterOfMassEnergy >= 199.) {
1312  if (fabs(mMagneticField) >= 4.) { // year=4, Au+Au, Full Field
1313  cent = Flow::cent200Year4Full;
1314  } else { // year=4, Au+Au, Half Field
1315  cent = Flow::cent200Year4Half;
1316  }
1317  } else if (mCenterOfMassEnergy >60. && mCenterOfMassEnergy < 65.) { // 62 GeV
1318  cent = Flow::cent62;
1319  }
1320  } else if (mCenterOfMassEnergy >= 199.) {
1321  if (fabs(mMagneticField) >= 4.) { // year=2, Au+Au, Full Field
1322  cent = Flow::cent200Full;
1323  } else { // year=2, Au+Au, Half Field
1324  cent = Flow::cent200Half;
1325  }
1326  } else if (mCenterOfMassEnergy == 0.) { // year=1
1327  cent = Flow::cent130;
1328  } else if (mCenterOfMassEnergy <= 25.){ // year=2, 22 GeV
1329  cent = Flow::cent22;
1330  }
1331 
1332  if (tracks < cent[0]) { mCentrality = 0; }
1333  else if (tracks < cent[1]) { mCentrality = 1; }
1334  else if (tracks < cent[2]) { mCentrality = 2; }
1335  else if (tracks < cent[3]) { mCentrality = 3; }
1336  else if (tracks < cent[4]) { mCentrality = 4; }
1337  else if (tracks < cent[5]) { mCentrality = 5; }
1338  else if (tracks < cent[6]) { mCentrality = 6; }
1339  else if (tracks < cent[7]) { mCentrality = 7; }
1340  else if (tracks < cent[8]) { mCentrality = 8; }
1341  else { mCentrality = 9; }
1342 
1343 }
1344 
1345 //-----------------------------------------------------------------------
1346 
1347 void StFlowEvent::PrintSelectionList() {
1348  // Prints the list of selection cuts
1349  // Call in Finish
1350 
1351  cout << "#######################################################" << endl;
1352  cout << "# Weighting and Striping:" << endl;
1353  if (mPtWgt) {
1354  cout << "# PtWgt= TRUE, also for output of PhiWgt file" << endl;
1355  cout << "# PtWgt Saturation= " << mPtWgtSaturation << endl;
1356  } else {
1357  cout << "# PtWgt= FALSE" << endl;
1358  }
1359  if (mEtaWgt) {
1360  cout << "# EtaWgt= TRUE, also for output of PhiWgt file for 1st harmonic" << endl;
1361  } else {
1362  cout << "# EtaWgt= FALSE" << endl;
1363  }
1364  if (mEtaSubs) {
1365  cout << "# EtaSubs= TRUE" << endl;
1366  } else {
1367  cout << "# EtaSubs= FALSE" << endl;
1368  }
1369  if (mRanSubs) {
1370  cout << "# RanSubs= TRUE" << endl;
1371  } else {
1372  cout << "# RanSubs= FALSE" << endl;
1373  }
1374  cout << "#######################################################" << endl;
1375  cout << "# Pid Deviant Cuts:" << endl;
1376  cout << "# PiPlus cuts= " << mPiPlusCuts[0] << ", "
1377  << mPiPlusCuts[1] << endl;
1378  cout << "# PiMinus cuts= " << mPiMinusCuts[0] << ", "
1379  << mPiMinusCuts[1] << endl;
1380  cout << "# Proton cuts= " << mProtonCuts[0] << ", "
1381  << mProtonCuts[1] << endl;
1382  cout << "# Anti Proton cuts= " << mAntiProtonCuts[0] << ", "
1383  << mAntiProtonCuts[1] << endl;
1384  cout << "# Deuteron cuts= " << mDeuteronCuts[0] << ", "
1385  << mDeuteronCuts[1] << endl;
1386  cout << "# Anti Deuteron cuts= " << mAntiDeuteronCuts[0] << ", "
1387  << mAntiDeuteronCuts[1] << endl;
1388  cout << "# K- cuts= " << mKMinusCuts[0] << ", "
1389  << mKMinusCuts[1] << endl;
1390  cout << "# K+ cuts= " << mKPlusCuts[0] << ", "
1391  << mKPlusCuts[1] << endl;
1392  cout << "# Electron cuts= " << mElectronCuts[0] << ", "
1393  << mElectronCuts[1] << endl;
1394  cout << "# Positron cuts= " << mPositronCuts[0] << ", "
1395  << mPositronCuts[1] << endl;
1396  cout << "#######################################################" << endl;
1397  cout << "# Tracks used for the event plane:" << endl;
1398  cout << "# Particle ID= " << mPid << endl;
1399  cout << "# Global Dca Tpc cuts= " << mDcaGlobalTpcCuts[0] << ", "
1400  << mDcaGlobalTpcCuts[1] << endl;
1401  cout << "# Global Dca Ftpc cuts= " << mDcaGlobalFtpcCuts[0] << ", "
1402  << mDcaGlobalFtpcCuts[1] << endl;
1403  for (int k = 0; k < Flow::nSels; k++) {
1404  for (int j = 0; j < 2; j++) {
1405  cout << "# selection= " << k+1 << " harmonic= "
1406  << j+1 << endl;
1407  cout << "# abs(Eta) Tpc cuts= " << mEtaTpcCuts[0][j][k] << ", "
1408  << mEtaTpcCuts[1][j][k] << endl;
1409  cout << "# Eta Ftpc cuts= " << mEtaFtpcCuts[0][j][k] << ", "
1410  << mEtaFtpcCuts[1][j][k] << ", " << mEtaFtpcCuts[2][j][k] << ", "
1411  << mEtaFtpcCuts[3][j][k] << endl;
1412  cout << "# Pt Tpc cuts= " << mPtTpcCuts[0][j][k] << ", "
1413  << mPtTpcCuts[1][j][k] << endl;
1414  cout << "# Pt Ftpc cuts= " << mPtFtpcCuts[0][j][k] << ", "
1415  << mPtFtpcCuts[1][j][k] << endl;
1416  }
1417  }
1418  cout << "#######################################################" << endl;
1419 
1420 }
1421 
1423 //
1424 // $Log: StFlowEvent.cxx,v $
1425 // Revision 1.63 2010/09/30 19:30:26 posk
1426 // Instead of reversing the weight for negative pseudrapidity for odd harmonics,
1427 // it is now done only for the first harmonic.
1428 // Recentering is now done for all harmonics.
1429 //
1430 // Revision 1.62 2009/11/24 19:23:01 posk
1431 // Added reCenter option to remove acceptance correlations instead of phiWgt.
1432 //
1433 // Revision 1.61 2009/08/04 23:00:28 posk
1434 // Reads year 7 MuDsts.
1435 //
1436 // Revision 1.60 2007/02/06 18:57:52 posk
1437 // In Lee Yang Zeros method, introduced recentering of Q vector.
1438 // Reactivated eta symmetry cut.
1439 //
1440 // Revision 1.59 2006/07/06 16:56:00 posk
1441 // Calculation of v1 for selection=2 is done with mixed harmonics.
1442 //
1443 // Revision 1.58 2006/02/22 19:29:14 posk
1444 // Additions needed for the StFlowLeeYangZerosMaker
1445 //
1446 // Revision 1.57 2005/10/27 17:53:19 aihong
1447 // changes made so that the cumulant method won't pick up the raw phi weight
1448 //
1449 // Revision 1.56 2005/02/11 23:24:29 posk
1450 // SetCentrality works for year4.
1451 //
1452 // Revision 1.55 2005/02/08 20:57:36 psoren
1453 // trigger and centrality selections were updated for all runs after run 4 to be compatible with trigger collections. Added TriggersFound() and GetFlowTriggerBitMap() functions.
1454 //
1455 // Revision 1.54 2004/12/17 22:33:10 aihong
1456 // add in full Psi weight for ZDC SMD and fix a few bugs, done by Gang
1457 //
1458 // Revision 1.53 2004/12/17 15:50:00 aihong
1459 // check in v1{3} code
1460 //
1461 // Revision 1.52 2004/12/09 23:43:35 posk
1462 // Minor changes in code formatting.
1463 //
1464 // Revision 1.51 2004/12/07 23:08:11 posk
1465 // Only odd and even phiWgt hists. If the old phiWgt file contains more than
1466 // two harmonics, only the first two are read. Now writes only the first two.
1467 //
1468 // Revision 1.50 2004/12/07 17:04:44 posk
1469 // Eliminated the very old mOnePhiWgt, which used one phiWgt histogram for flttening
1470 // instead of four.
1471 //
1472 // Revision 1.49 2004/11/16 21:22:21 aihong
1473 // removed old cumulant method
1474 //
1475 // Revision 1.48 2004/08/24 20:24:33 oldi
1476 // Minor modifications to avoid compiler warnings.
1477 // Small bug fix (didn't affect anyone yet).
1478 //
1479 // Revision 1.47 2004/05/31 20:09:35 oldi
1480 // PicoDst format changed (Version 7) to hold ZDC SMD information.
1481 // Trigger cut modified to comply with TriggerCollections.
1482 // Centrality definition for 62 GeV data introduced.
1483 // Minor bug fixes.
1484 //
1485 // Revision 1.46 2004/05/05 21:13:45 aihong
1486 // Gang's code for ZDC-SMD added
1487 //
1488 // Revision 1.45 2004/03/11 17:58:40 posk
1489 // Added Random Subs analysis method.
1490 //
1491 // Revision 1.44 2003/09/02 17:58:11 perev
1492 // gcc 3.2 updates + WarnOff
1493 //
1494 // Revision 1.43 2003/07/30 22:00:39 oldi
1495 // Eta cuts for event plane selection separated for FTPC east and west.
1496 // PtWgtSaturation parameter introduced (default set to 2. -> no change of default behavior).
1497 //
1498 // Revision 1.42 2003/07/15 18:34:26 oldi
1499 // Printout for upper FTPC pt cut for the event plane determination was showing
1500 // upper TPC pt cut always. Changed.
1501 //
1502 // Revision 1.41 2003/06/18 17:00:49 posk
1503 // Event plane cuts now only odd and even, instead of different for each harmonic.
1504 //
1505 // Revision 1.40 2003/05/15 06:08:41 aihong
1506 // default PID is changed from none to NA, SetDedxPtsPart() added
1507 //
1508 // Revision 1.39 2003/05/02 21:09:41 posk
1509 // Reduced the number of harmonics from 3 to 2.
1510 //
1511 // Revision 1.38 2003/04/01 00:27:05 posk
1512 // Little q is now unweighted by pt or eta. Big Q is unaffected.
1513 //
1514 // Revision 1.37 2003/02/25 19:28:40 posk
1515 // Changed a few unimportant default cuts.
1516 //
1517 // Revision 1.36 2003/01/10 16:42:09 oldi
1518 // Several changes to comply with FTPC tracks:
1519 // - Switch to include/exclude FTPC tracks introduced.
1520 // The same switch changes the range of the eta histograms.
1521 // - Eta symmetry plots for FTPC tracks added and separated from TPC plots.
1522 // - PhiWgts and related histograms for FTPC tracks split in FarEast, East,
1523 // West, FarWest (depending on vertex.z()).
1524 // - Psi_Diff plots for 2 different selections and the first 2 harmonics added.
1525 // - Cut to exclude mu-events with no primary vertex introduced.
1526 // (This is possible for UPC events and FTPC tracks.)
1527 // - Global DCA cut for FTPC tracks added.
1528 // - Global DCA cuts for event plane selection separated for TPC and FTPC tracks.
1529 // - Charge cut for FTPC tracks added.
1530 //
1531 // Revision 1.35 2003/01/08 19:26:46 posk
1532 // PhiWgt hists sorted on sign of z of first and last points.
1533 // Version 6 of pico file.
1534 //
1535 // Revision 1.34 2002/06/10 22:50:59 posk
1536 // pt and eta weighting now default.
1537 // DcaGlobalPart default now 0 to 1 cm.
1538 // Event cut order changed.
1539 //
1540 // Revision 1.33 2002/05/23 18:54:10 posk
1541 // Moved centrality cuts into StFlowConstants
1542 //
1543 // Revision 1.32 2002/03/15 16:43:21 snelling
1544 // Added a method to recalculate the centrality in StFlowPicoEvent
1545 //
1546 // Revision 1.31 2002/03/14 18:51:49 snelling
1547 // Added new centralities
1548 //
1549 // Revision 1.30 2002/02/13 22:29:21 posk
1550 // Pt Weight now also weights Phi Weights. Added Eta Weight, default=FALSE.
1551 //
1552 // Revision 1.29 2002/01/31 01:04:43 posk
1553 // *** empty log message ***
1554 //
1555 // Revision 1.28 2001/12/18 19:22:02 posk
1556 // "proton" and "antiproton" changed to "pr+" and "pr-".
1557 // Compiles on Solaris.
1558 //
1559 // Revision 1.27 2001/12/11 21:33:43 posk
1560 // Went from one to four sets of histograms for making the event plane isotropic.
1561 // StFlowEvent::PhiWeight() has changed arguments and return value.
1562 // The ptWgt saturates above 2 GeV/c.
1563 //
1564 // Revision 1.26 2001/11/09 21:10:37 posk
1565 // Switched from CERNLIB to TMath. Little q is now normalized.
1566 //
1567 // Revision 1.25 2001/11/02 04:49:52 aihong
1568 // add func. for cumulant maker
1569 //
1570 // Revision 1.24 2001/08/08 10:35:07 oldi
1571 // Typo in output statement of cut lists fixed (mEtaTpcCuts[1][j][k] -> mEtaFtpcCuts[1][j][k]).
1572 //
1573 // Revision 1.23 2001/06/07 20:06:16 posk
1574 // Global Dca cut for event plane particles.
1575 // Removed SetPtWgt().
1576 //
1577 // Revision 1.22 2001/05/22 20:17:26 posk
1578 // Now can do pseudorapidity subevents.
1579 //
1580 // Revision 1.21 2001/04/03 17:47:17 oldi
1581 // Bug fix that excluded FTPC tracks from the determination of the reaction plane.
1582 //
1583 // Revision 1.20 2000/12/12 20:22:05 posk
1584 // Put log comments at end of files.
1585 // Deleted persistent StFlowEvent (old micro DST).
1586 //
1587 // Revision 1.19 2000/12/10 02:01:13 oldi
1588 // A new member (StTrackTopologyMap mTopology) was added to StFlowPicoTrack.
1589 // The evaluation of either a track originates from the FTPC or not is
1590 // unambiguous now. The evaluation itself is easily extendible for other
1591 // detectors (e.g. SVT+TPC). Old flowpicoevent.root files are treated as if
1592 // they contain TPC tracks only (backward compatibility).
1593 //
1594 // Revision 1.18 2000/12/08 17:03:38 oldi
1595 // Phi weights for both FTPCs included.
1596 //
1597 // Revision 1.17 2000/10/12 22:46:35 snelling
1598 // Added support for the new pDST's and the probability pid method
1599 //
1600 // Revision 1.16 2000/09/26 20:51:37 posk
1601 // Updated documentation.
1602 //
1603 // Revision 1.14 2000/09/15 22:51:28 posk
1604 // Added pt weighting for event plane calcualtion.
1605 //
1606 // Revision 1.13 2000/09/12 01:30:23 snelling
1607 // Changed PID selection
1608 //
1609 // Revision 1.12 2000/09/05 16:11:31 snelling
1610 // Added global DCA, electron and positron
1611 //
1612 // Revision 1.11 2000/08/31 18:58:21 posk
1613 // For picoDST, added version number, runID, and multEta for centrality.
1614 // Added centrality cut when reading picoDST.
1615 // Added pt and eta selections for particles corr. wrt event plane.
1616 //
1617 // Revision 1.10 2000/08/12 20:22:19 posk
1618 // Recalculate centrality in read from pico.
1619 //
1620 // Revision 1.9 2000/08/10 23:00:21 posk
1621 // New centralities. pt and eta cuts.
1622 //
1623 // Revision 1.8 2000/08/09 21:38:23 snelling
1624 // PID added
1625 //
1626 // Revision 1.7 2000/06/01 18:26:35 posk
1627 // Increased precision of Track integer data members.
1628 //
1629 // Revision 1.5 2000/05/16 20:59:29 posk
1630 // Voloshin's flownanoevent.root added.
1631 //
1632 // Revision 1.4 2000/05/12 22:42:04 snelling
1633 // Additions for persistency and minor fix
1634 //
1635 // Revision 1.2 2000/03/15 23:28:50 posk
1636 // Added StFlowSelection.
1637 //
1638 // Revision 1.1 2000/03/02 23:02:48 posk
1639 // Changed extensions from .hh and .cc to .h and .cxx .
1640 //
1641 // Revision 1.16 2000/02/29 22:00:53 posk
1642 // Made SetPhiWeight inline, changed ImpactPar to Dca, etc.
1643 //
1644 // Revision 1.15 2000/02/29 01:26:11 snelling
1645 // removed static const int& nxxx = Flow::nxxx;
1646 //
1647 // Revision 1.14 2000/02/18 22:49:54 posk
1648 // Added PID and centrality.
1649 //
1650 // Revision 1.9 1999/12/21 17:31:50 posk
1651 // Fixed random_shuffle in making the sub events.
1652 //
1653 // Revision 1.6 1999/12/15 22:01:25 posk
1654 // Added StFlowConstants.hh
1655 //
1656 // Revision 1.4 1999/12/04 00:10:32 posk
1657 // Works with the new StEvent
1658 //
1659 // Revision 1.3 1999/11/30 18:52:51 snelling
1660 // First modification for the new StEvent
1661 //
1662 // Revision 1.2 1999/11/24 18:17:13 posk
1663 // Put the methods which act on the data in with the data in StFlowEvent.
1664 //
1665 // Revision 1.1 1999/11/04 19:02:04 snelling
1666 // First check in of StFlowMaker. It contains the common code from
1667 // StFlowTagMaker and StFlowAnalysisMaker.
1668 //