StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFlowEvent.h
1 //
3 // $Id: StFlowEvent.h,v 1.56 2009/11/24 19:23:03 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 with flow functions
12 //
14 
15 #ifndef StFlowEvent_h
16 #define StFlowEvent_h
17 #include "StObject.h"
18 #include "StFlowTrackCollection.h"
19 #include "StTrackTopologyMap.h"
20 #include "StThreeVectorF.hh"
21 #include "StEnumerations.h"
22 #include "Rtypes.h"
23 #include "TVector2.h"
24 #include "TComplex.h"
25 class StFlowSelection;
26 
27 class StFlowEvent : public StObject {
28 
29 public:
30 
31  StFlowEvent();
32  virtual ~StFlowEvent();
33 
34  Double_t PhiWeight(Int_t selN, Int_t harN, StFlowTrack* pFlowTrack) const;
35  Double_t PhiWeightRaw(Int_t selN, Int_t harN, StFlowTrack* pFlowTrack) const;
36  Double_t Weight(Int_t selN, Int_t harN, StFlowTrack* pFlowTrack) const;
37  Double_t ZDCSMD_PsiWgtEast();
38  Double_t ZDCSMD_PsiWgtWest();
39  Double_t ZDCSMD_PsiWgtFull();
40  Int_t EventID() const;
41  Int_t RunID() const;
42  Double_t CenterOfMassEnergy() const;
43  Double_t MagneticField() const;
44  Short_t BeamMassNumberEast() const;
45  Short_t BeamMassNumberWest() const;
46  UInt_t OrigMult() const;
47  UInt_t L0TriggerWord() const;
48  UInt_t UncorrNegMult() const;
49  UInt_t UncorrPosMult() const;
50  UInt_t MultEta() const;
51  UInt_t FlowEventMult() const;
52  UInt_t Centrality() const;
53  StThreeVectorF VertexPos() const;
54  UInt_t Mult(StFlowSelection*);
55  UInt_t MultPart(StFlowSelection*);
56  TVector2 Q(StFlowSelection*);
57  TVector2 NormQ(StFlowSelection* pFlowSelect);
58  TVector2 QPart(StFlowSelection*);
59  TVector2 ReCentEPPar(StFlowSelection*, const char*); // for ana
60  TVector2 ReCentPar(StFlowSelection*, const char*); // for LYZ
61  TVector2 ReCent(Int_t selN, Int_t harN, StFlowTrack* pFlowTrack) const;
62  TVector2 ReCentEP(Int_t selN, Int_t harN, StFlowTrack* pFlowTrack) const;
63  Float_t q(StFlowSelection*);
64  Float_t MeanPt(StFlowSelection*);
65  Float_t Qtheta(StFlowSelection*, Float_t theta);
66  TComplex Grtheta(StFlowSelection*, Float_t r, Float_t theta);
67  TComplex GV1r0theta(StFlowSelection*, Float_t r, Float_t theta1, Float_t theta);
68  TComplex Gder_r0theta(StFlowSelection*, Float_t r, Float_t theta);
69  Float_t Psi(StFlowSelection*);
70  Float_t ZDCSMD_PsiCorr();
71  Float_t ZDCSMD_PsiEst();
72  Float_t ZDCSMD_PsiWst();
73  Float_t ZDCSMD_GetPosition(int eastwest,int verthori,int strip);
74  Double_t G_New(StFlowSelection* pFlowSelect, Double_t Zx, Double_t Zy);
75  Double_t G_Mix(StFlowSelection* pFlowSelect, Double_t Z1x, Double_t Z1y, Double_t Z2x, Double_t Z2y);
76  Double_t SumWeightSquare(StFlowSelection* pFlowSelect);
77  Double_t PtAbsWgtValue(Double_t pt) const;
78  Double_t EtaAbsWgtValue(Double_t eta) const;
79  Float_t CTB() const;
80  Float_t ZDCe() const;
81  Float_t ZDCw() const;
82  Float_t ZDCSMD(int eastwest,int verthori,int strip) const;
83  Float_t PtWgtSaturation() const;
84  Bool_t PtWgt() const;
85  Bool_t EtaWgt() const;
86  Bool_t FirstLastPhiWgt() const;
87  Bool_t FirstLastPoints() const;
88  Bool_t UseZDCSMD() const;
89  Char_t* Pid();
90  Bool_t EtaSubs() const;
91  Bool_t RanSubs() const;
92 
93  Float_t V1TPCDetctWgtG_Mix(Int_t selN) const;
94  Float_t V1FtpcEastDetctWgtG_Mix(Int_t selN) const;
95  Float_t V1FtpcWestDetctWgtG_Mix(Int_t selN) const;
96 
97  Float_t V2TPCDetctWgtG_Mix(Int_t selN) const;
98  Float_t V2FtpcEastDetctWgtG_Mix(Int_t selN) const;
99  Float_t V2FtpcWestDetctWgtG_Mix(Int_t selN) const;
100 
101 
102  StFlowTrackCollection* TrackCollection() const;
103 
104  void SetSelections();
105  void SetPid(const Char_t*);
106  void SetPidsDeviant();
107  void SetPidsProb();
108  void PrintSelectionList();
109  void MakeSubEvents();
110  void MakeEtaSubEvents();
111  void SetEventID(const Int_t&);
112  void SetRunID(const Int_t&);
113  void SetCenterOfMassEnergy(const Double_t&);
114  void SetMagneticField(const Double_t&);
115  void SetBeamMassNumberEast(const Short_t&);
116  void SetBeamMassNumberWest(const Short_t&);
117  void SetOrigMult(const UInt_t&);
118  void SetL0TriggerWord(const UInt_t&);
119  void SetUncorrPosMult(const UInt_t&);
120  void SetUncorrNegMult(const UInt_t&);
121  void SetMultEta(const UInt_t&);
122  void SetCentrality();
123  void SetVertexPos(const StThreeVectorF&);
124  void SetCTB(const Float_t ctb);
125  void SetZDCe(const Float_t zdce);
126  void SetZDCw(const Float_t zdcw);
127  void SetZDCSMD(int eastwest,int verthori,int strip,const Float_t zdcsmd);
128 #ifndef __CINT__
129  void SetPhiWeightFarEast(const Flow::PhiWgt_t &pPhiWgt);
130  void SetPhiWeightEast(const Flow::PhiWgt_t &pPhiWgt);
131  void SetPhiWeightWest(const Flow::PhiWgt_t &pPhiWgt);
132  void SetPhiWeightFarWest(const Flow::PhiWgt_t &pPhiWgt);
133  void SetPhiWeightFtpcFarEast(const Flow::PhiWgtFtpc_t &pPhiWgt);
134  void SetPhiWeightFtpcEast(const Flow::PhiWgtFtpc_t &pPhiWgt);
135  void SetPhiWeightFtpcWest(const Flow::PhiWgtFtpc_t &pPhiWgt);
136  void SetPhiWeightFtpcFarWest(const Flow::PhiWgtFtpc_t &pPhiWgt);
137  void SetZDCSMD_PsiWeightWest(const Flow::ZDCSMD_PsiWgt_t& ZDCSMD_PsiWgtWest);
138  void SetZDCSMD_PsiWeightEast(const Flow::ZDCSMD_PsiWgt_t& ZDCSMD_PsiWgtEast);
139  void SetZDCSMD_PsiWeightFull(const Flow::ZDCSMD_PsiWgt_t& ZDCSMD_PsiWgtFull);
140  void SetZDCSMD_BeamCenter(Double_t ex,Double_t ey,Double_t wx,Double_t wy);
141  void SetReCentX(const Flow::ReCent_t &pReCentX);
142  void SetReCentY(const Flow::ReCent_t &pReCentY);
143 #endif
144  static void SetEtaTpcCut(Float_t lo, Float_t hi, Int_t harN, Int_t selN);
145  static void SetPtTpcCut(Float_t lo, Float_t hi, Int_t harN, Int_t selN);
146  static void SetEtaFtpcCut(Float_t lo_neg, Float_t hi_neg,
147  Float_t lo_pos, Float_t hi_pos,
148  Int_t harN, Int_t selN);
149  static void SetPtFtpcCut(Float_t lo, Float_t hi, Int_t harN, Int_t selN);
150  static void SetPiPlusCut(Float_t lo, Float_t hi);
151  static void SetPiMinusCut(Float_t lo, Float_t hi);
152  static void SetProtonCut(Float_t lo, Float_t hi);
153  static void SetAntiProtonCut(Float_t lo, Float_t hi);
154  static void SetKPlusCut(Float_t lo, Float_t hi);
155  static void SetKMinusCut(Float_t lo, Float_t hi);
156  static void SetElectronCut(Float_t lo, Float_t hi);
157  static void SetPositronCut(Float_t lo, Float_t hi);
158  static void SetDeuteronCut(Float_t lo, Float_t hi);
159  static void SetAntiDeuteronCut(Float_t lo, Float_t hi);
160  static void SetDcaGlobalTpcCut(Float_t lo, Float_t hi);
161  static void SetDcaGlobalFtpcCut(Float_t lo, Float_t hi);
162  static void SetProbPid();
163  static void SetEtaSubs();
164  static void SetRanSubs();
165  static void SetPtWgtSaturation(Float_t);
166  static void SetPtWgt(Bool_t);
167  static void SetEtaWgt(Bool_t);
168  static void SetFirstLastPhiWgt();
169  static void SetFirstLastPoints();
170  static void SetUseZDCSMD(Bool_t);
171 
172  static void SetV1TPCDetctWgtG_Mix(Float_t val, Int_t selN);
173  static void SetV1FtpcEastDetctWgtG_Mix(Float_t val, Int_t selN);
174  static void SetV1FtpcWestDetctWgtG_Mix(Float_t val, Int_t selN);
175 
176 
177  static void SetV2TPCDetctWgtG_Mix(Float_t val, Int_t selN);
178  static void SetV2FtpcEastDetctWgtG_Mix(Float_t val, Int_t selN);
179  static void SetV2FtpcWestDetctWgtG_Mix(Float_t val, Int_t selN);
180  static Bool_t ProbPid();
181 
182 
183 private:
184 
185  Int_t mEventID; // ID of the event
186  Int_t mRunID; // ID of the run
187  Double_t mMagneticField; //
188  Double_t mCenterOfMassEnergy; //
189  Short_t mBeamMassNumberEast; //
190  Short_t mBeamMassNumberWest; //
191  UInt_t mOrigMult; // number of tracks
192  UInt_t mL0TriggerWord; // L0 trigger word
193  UInt_t mUncorrNegMult; // number of h-
194  UInt_t mUncorrPosMult; // number of h+
195  UInt_t mMultEta; // number of tracks
196  // with pos. flag in 1.5 unit of eta
197  UInt_t mCentrality; // centrality bin
198  StThreeVectorF mVertexPos; // primary vertex position
199  Float_t mCTB; // CTB value sum
200  Float_t mZDCe; // ZDC east
201  Float_t mZDCw; // ZDC west
202  Float_t mZDCSMD[2][2][8]; // ZDCSMD
203  static Double_t mZDCSMDCenterex,mZDCSMDCenterey; // ZDCSMD Beam Center
204  static Double_t mZDCSMDCenterwx,mZDCSMDCenterwy; // ZDCSMD Beam Center
205  static Float_t mEtaTpcCuts[2][2][Flow::nSels]; // range absolute values
206  static Float_t mEtaFtpcCuts[4][2][Flow::nSels]; // range values
207  static Float_t mPtTpcCuts[2][2][Flow::nSels]; // range
208  static Float_t mPtFtpcCuts[2][2][Flow::nSels]; // range
209 
210  static Float_t mV1TPCDetctWgtG_Mix[Flow::nSels]; // detector wgt for G_Mix in v1{3} calc.
211  static Float_t mV1FtpcEastDetctWgtG_Mix[Flow::nSels];
212  static Float_t mV1FtpcWestDetctWgtG_Mix[Flow::nSels];
213  static Float_t mV2TPCDetctWgtG_Mix[Flow::nSels];
214  static Float_t mV2FtpcEastDetctWgtG_Mix[Flow::nSels];
215  static Float_t mV2FtpcWestDetctWgtG_Mix[Flow::nSels];
216 
217  Flow::PhiWgt_t mPhiWgtFarEast;
218  Flow::PhiWgt_t mPhiWgtEast;
219  Flow::PhiWgt_t mPhiWgtWest;
220  Flow::PhiWgt_t mPhiWgtFarWest;
221  Flow::PhiWgtFtpc_t mPhiWgtFtpcFarEast;
222  Flow::PhiWgtFtpc_t mPhiWgtFtpcEast;
223  Flow::PhiWgtFtpc_t mPhiWgtFtpcWest;
224  Flow::PhiWgtFtpc_t mPhiWgtFtpcFarWest;
225  Flow::ZDCSMD_PsiWgt_t mZDCSMD_PsiWgtWest;
226  Flow::ZDCSMD_PsiWgt_t mZDCSMD_PsiWgtEast;
227  Flow::ZDCSMD_PsiWgt_t mZDCSMD_PsiWgtFull;
228  Flow::ReCent_t mReCentX;
229  Flow::ReCent_t mReCentY;
230 
231  static Float_t mPiPlusCuts[2]; // PID cuts
232  static Float_t mPtWgtSaturation; // saturation value for pt weighting
233  static Bool_t mPtWgt; // flag for pt weighting
234  static Bool_t mEtaWgt; // flag for y weighting for odd harmonics
235  static Char_t mPid[10]; // h+, h-, pi-, pi+, pi, k+, k-, k, pr+, pr-, pr, e+, e-, e
236  static Bool_t mProbPid; // flag for probability pid
237  static Bool_t mEtaSubs; // flag for eta subevents
238  static Bool_t mRanSubs; // flag for random subevents
239  static Bool_t mFirstLastPhiWgt; // flag for using z of first and last points for reading phi weights
240  static Bool_t mFirstLastPoints; // flag for using z of first and last points for generating phi weights
241  static Bool_t mUseZDCSMD; // flag for using ZDC SMD for RP
242  static Float_t mPiMinusCuts[2];
243  static Float_t mProtonCuts[2];
244  static Float_t mKMinusCuts[2];
245  static Float_t mKPlusCuts[2];
246  static Float_t mAntiProtonCuts[2];
247  static Float_t mDeuteronCuts[2];
248  static Float_t mAntiDeuteronCuts[2];
249  static Float_t mElectronCuts[2];
250  static Float_t mPositronCuts[2];
251  static Float_t mDcaGlobalTpcCuts[2];
252  static Float_t mDcaGlobalFtpcCuts[2];
253 
254  StFlowTrackCollection* pTrackCollection;
255 
256  ClassDef(StFlowEvent,1) // macro for rootcint
257 
258 };
259 
260 inline StFlowTrackCollection* StFlowEvent::TrackCollection() const {
261  return pTrackCollection; }
262 
263 inline Int_t StFlowEvent::EventID() const { return mEventID; }
264 
265 inline Int_t StFlowEvent::RunID() const { return mRunID; }
266 
267 inline Double_t StFlowEvent::CenterOfMassEnergy() const { return mCenterOfMassEnergy; }
268 
269 inline Double_t StFlowEvent::MagneticField() const { return mMagneticField; }
270 
271 inline Short_t StFlowEvent::BeamMassNumberEast() const { return mBeamMassNumberEast; }
272 
273 inline Short_t StFlowEvent::BeamMassNumberWest() const { return mBeamMassNumberWest; }
274 
275 inline UInt_t StFlowEvent::OrigMult() const { return mOrigMult; }
276 
277 inline UInt_t StFlowEvent::L0TriggerWord() const { return mL0TriggerWord; }
278 
279 inline UInt_t StFlowEvent::UncorrNegMult() const { return mUncorrNegMult; }
280 
281 inline UInt_t StFlowEvent::UncorrPosMult() const { return mUncorrPosMult; }
282 
283 inline UInt_t StFlowEvent::MultEta() const { return mMultEta; }
284 
285 inline UInt_t StFlowEvent::FlowEventMult() const { return pTrackCollection->size(); }
286 
287 inline UInt_t StFlowEvent::Centrality() const { return mCentrality; }
288 
289 inline StThreeVectorF StFlowEvent::VertexPos() const { return mVertexPos; }
290 
291 inline Float_t StFlowEvent::CTB() const { return mCTB; }
292 
293 inline Float_t StFlowEvent::ZDCe() const { return mZDCe; }
294 
295 inline Float_t StFlowEvent::ZDCw() const { return mZDCw; }
296 
297 inline Float_t StFlowEvent::ZDCSMD(int eastwest,int verthori,int strip) const {return mZDCSMD[eastwest][verthori][strip-1];}
298 
299 inline Float_t StFlowEvent::PtWgtSaturation() const { return mPtWgtSaturation; }
300 
301 inline Bool_t StFlowEvent::PtWgt() const { return mPtWgt; }
302 
303 inline Bool_t StFlowEvent::EtaWgt() const { return mEtaWgt; }
304 
305 inline Bool_t StFlowEvent::FirstLastPhiWgt() const { return mFirstLastPhiWgt; }
306 
307 inline Bool_t StFlowEvent::FirstLastPoints() const { return mFirstLastPoints; }
308 
309 inline Char_t* StFlowEvent::Pid() { return mPid; }
310 
311 inline Bool_t StFlowEvent::ProbPid() { return mProbPid; }
312 
313 inline Bool_t StFlowEvent::EtaSubs() const { return mEtaSubs; }
314 
315 inline Bool_t StFlowEvent::RanSubs() const { return mRanSubs; }
316 
317 inline Bool_t StFlowEvent::UseZDCSMD() const { return mUseZDCSMD;}
318 
319 inline Float_t StFlowEvent::V1TPCDetctWgtG_Mix(Int_t selN) const {
320  return mV1TPCDetctWgtG_Mix[selN]; }
321 inline Float_t StFlowEvent::V1FtpcEastDetctWgtG_Mix(Int_t selN) const {
322  return mV1FtpcEastDetctWgtG_Mix[selN]; }
323 inline Float_t StFlowEvent::V1FtpcWestDetctWgtG_Mix(Int_t selN) const {
324  return mV1FtpcWestDetctWgtG_Mix[selN]; }
325 
326 inline Float_t StFlowEvent::V2TPCDetctWgtG_Mix(Int_t selN) const {
327  return mV2TPCDetctWgtG_Mix[selN]; }
328 inline Float_t StFlowEvent::V2FtpcEastDetctWgtG_Mix(Int_t selN) const {
329  return mV2FtpcEastDetctWgtG_Mix[selN]; }
330 inline Float_t StFlowEvent::V2FtpcWestDetctWgtG_Mix(Int_t selN) const {
331  return mV2FtpcWestDetctWgtG_Mix[selN]; }
332 
333 #ifndef __CINT__
334 inline void StFlowEvent::SetPhiWeightFarEast(const Flow::PhiWgt_t& pPhiWgtFarEast) {
335  memcpy (mPhiWgtFarEast, pPhiWgtFarEast, sizeof(Flow::PhiWgt_t)); }
336 
337 inline void StFlowEvent::SetPhiWeightEast(const Flow::PhiWgt_t& pPhiWgtEast) {
338  memcpy (mPhiWgtEast, pPhiWgtEast, sizeof(Flow::PhiWgt_t)); }
339 
340 inline void StFlowEvent::SetPhiWeightWest(const Flow::PhiWgt_t& pPhiWgtWest) {
341  memcpy (mPhiWgtWest, pPhiWgtWest, sizeof(Flow::PhiWgt_t)); }
342 
343 inline void StFlowEvent::SetPhiWeightFarWest(const Flow::PhiWgt_t& pPhiWgtFarWest) {
344  memcpy (mPhiWgtFarWest, pPhiWgtFarWest, sizeof(Flow::PhiWgt_t)); }
345 
346 inline void StFlowEvent::SetPhiWeightFtpcFarEast(const Flow::PhiWgtFtpc_t& pPhiWgtFtpcFarEast) {
347  memcpy (mPhiWgtFtpcFarEast, pPhiWgtFtpcFarEast, sizeof(Flow::PhiWgtFtpc_t)); }
348 
349 inline void StFlowEvent::SetPhiWeightFtpcEast(const Flow::PhiWgtFtpc_t& pPhiWgtFtpcEast) {
350  memcpy (mPhiWgtFtpcEast, pPhiWgtFtpcEast, sizeof(Flow::PhiWgtFtpc_t)); }
351 
352 inline void StFlowEvent::SetPhiWeightFtpcWest(const Flow::PhiWgtFtpc_t& pPhiWgtFtpcWest) {
353  memcpy (mPhiWgtFtpcWest, pPhiWgtFtpcWest, sizeof(Flow::PhiWgtFtpc_t)); }
354 
355 inline void StFlowEvent::SetPhiWeightFtpcFarWest(const Flow::PhiWgtFtpc_t& pPhiWgtFtpcFarWest) {
356  memcpy (mPhiWgtFtpcFarWest, pPhiWgtFtpcFarWest, sizeof(Flow::PhiWgtFtpc_t)); }
357 
358 inline void StFlowEvent::SetZDCSMD_PsiWeightEast(const Flow::ZDCSMD_PsiWgt_t& ZDCSMD_PsiWgtEast) {
359  memcpy (mZDCSMD_PsiWgtEast, ZDCSMD_PsiWgtEast, sizeof(Flow::ZDCSMD_PsiWgt_t)); }
360 
361 inline void StFlowEvent::SetZDCSMD_PsiWeightWest(const Flow::ZDCSMD_PsiWgt_t& ZDCSMD_PsiWgtWest) {
362  memcpy (mZDCSMD_PsiWgtWest, ZDCSMD_PsiWgtWest, sizeof(Flow::ZDCSMD_PsiWgt_t)); }
363 
364 inline void StFlowEvent::SetZDCSMD_PsiWeightFull(const Flow::ZDCSMD_PsiWgt_t& ZDCSMD_PsiWgtFull) {
365  memcpy (mZDCSMD_PsiWgtFull, ZDCSMD_PsiWgtFull, sizeof(Flow::ZDCSMD_PsiWgt_t)); }
366 
367 inline void StFlowEvent::SetZDCSMD_BeamCenter(Double_t ex,Double_t ey,Double_t wx,Double_t wy) {
368  mZDCSMDCenterex = ex; mZDCSMDCenterey = ey; mZDCSMDCenterwx = wx; mZDCSMDCenterwy = wy;}
369 
370 inline void StFlowEvent::SetReCentX(const Flow::ReCent_t& pReCentX) {
371  memcpy (mReCentX, pReCentX, sizeof(Flow::ReCent_t)); }
372 
373 inline void StFlowEvent::SetReCentY(const Flow::ReCent_t& pReCentY) {
374  memcpy (mReCentY, pReCentY, sizeof(Flow::ReCent_t)); }
375 #endif
376 
377 inline void StFlowEvent::SetEtaTpcCut(Float_t lo, Float_t hi, Int_t harN, Int_t selN)
378 { mEtaTpcCuts[0][harN][selN] = lo; mEtaTpcCuts[1][harN][selN] = hi; }
379 
380 inline void StFlowEvent::SetPtTpcCut(Float_t lo, Float_t hi, Int_t harN, Int_t selN)
381 { mPtTpcCuts[0][harN][selN] = lo; mPtTpcCuts[1][harN][selN] = hi; }
382 
383 inline void StFlowEvent::SetEtaFtpcCut(Float_t lo_neg, Float_t hi_neg,
384  Float_t lo_pos, Float_t hi_pos,
385  Int_t harN, Int_t selN)
386 { mEtaFtpcCuts[0][harN][selN] = lo_neg; mEtaFtpcCuts[1][harN][selN] = hi_neg;
387  mEtaFtpcCuts[2][harN][selN] = lo_pos; mEtaFtpcCuts[3][harN][selN] = hi_pos; }
388 
389 inline void StFlowEvent::SetPtFtpcCut(Float_t lo, Float_t hi, Int_t harN, Int_t selN)
390 { mPtFtpcCuts[0][harN][selN] = lo; mPtFtpcCuts[1][harN][selN] = hi; }
391 
392 inline void StFlowEvent::SetEventID(const Int_t& id) { mEventID = id; }
393 
394 inline void StFlowEvent::SetRunID(const Int_t& id) { mRunID = id; }
395 
396 inline void StFlowEvent::SetMagneticField(const Double_t& mf) { mMagneticField = mf; }
397 
398 inline void StFlowEvent::SetCenterOfMassEnergy(const Double_t& cms) { mCenterOfMassEnergy = cms; }
399 
400 inline void StFlowEvent::SetBeamMassNumberEast(const Short_t& bme) { mBeamMassNumberEast = bme; }
401 
402 inline void StFlowEvent::SetBeamMassNumberWest(const Short_t& bmw) { mBeamMassNumberWest = bmw; }
403 
404 inline void StFlowEvent::SetOrigMult(const UInt_t& tracks) {
405  mOrigMult = tracks; }
406 
407 inline void StFlowEvent::SetL0TriggerWord(const UInt_t& trigger) {
408  mL0TriggerWord = trigger; }
409 
410 inline void StFlowEvent::SetUncorrNegMult(const UInt_t& negtracks) {
411  mUncorrNegMult = negtracks; }
412 
413 inline void StFlowEvent::SetUncorrPosMult(const UInt_t& postracks) {
414  mUncorrPosMult = postracks; }
415 
416 inline void StFlowEvent::SetMultEta(const UInt_t& goodtracks) {
417  mMultEta = goodtracks; }
418 
419 inline void StFlowEvent::SetVertexPos(const StThreeVectorF& vertexPos) {
420  mVertexPos = vertexPos; }
421 
422 inline void StFlowEvent::SetPiPlusCut(Float_t lo, Float_t hi) {
423  mPiPlusCuts[0] = lo; mPiPlusCuts[1] = hi; }
424 
425 inline void StFlowEvent::SetPiMinusCut(Float_t lo, Float_t hi) {
426  mPiMinusCuts[0] = lo; mPiMinusCuts[1] = hi; }
427 
428 inline void StFlowEvent::SetProtonCut(Float_t lo, Float_t hi) {
429  mProtonCuts[0] = lo; mProtonCuts[1] = hi; }
430 
431 inline void StFlowEvent::SetKMinusCut(Float_t lo, Float_t hi) {
432  mKMinusCuts[0] = lo; mKMinusCuts[1] = hi; }
433 
434 inline void StFlowEvent::SetKPlusCut(Float_t lo, Float_t hi) {
435  mKPlusCuts[0] = lo; mKPlusCuts[1] = hi; }
436 
437 inline void StFlowEvent::SetAntiProtonCut(Float_t lo, Float_t hi) {
438  mAntiProtonCuts[0] = lo; mAntiProtonCuts[1] = hi; }
439 
440 inline void StFlowEvent::SetDeuteronCut(Float_t lo, Float_t hi) {
441  mDeuteronCuts[0] = lo; mDeuteronCuts[1] = hi; }
442 
443 inline void StFlowEvent::SetAntiDeuteronCut(Float_t lo, Float_t hi) {
444  mAntiDeuteronCuts[0] = lo; mAntiDeuteronCuts[1] = hi; }
445 
446 inline void StFlowEvent::SetElectronCut(Float_t lo, Float_t hi) {
447  mElectronCuts[0] = lo; mElectronCuts[1] = hi; }
448 
449 inline void StFlowEvent::SetPositronCut(Float_t lo, Float_t hi) {
450  mPositronCuts[0] = lo; mPositronCuts[1] = hi; }
451 
452 inline void StFlowEvent::SetDcaGlobalTpcCut(Float_t lo, Float_t hi) {
453  mDcaGlobalTpcCuts[0] = lo; mDcaGlobalTpcCuts[1] = hi; }
454 
455 inline void StFlowEvent::SetDcaGlobalFtpcCut(Float_t lo, Float_t hi) {
456  mDcaGlobalFtpcCuts[0] = lo; mDcaGlobalFtpcCuts[1] = hi; }
457 
458 inline void StFlowEvent::SetCTB(const Float_t ctb) { mCTB = ctb; }
459 
460 inline void StFlowEvent::SetZDCe(const Float_t zdce) { mZDCe = zdce; }
461 
462 inline void StFlowEvent::SetZDCw(const Float_t zdcw) { mZDCw = zdcw; }
463 
464 inline void StFlowEvent::SetZDCSMD(int eastwest,int verthori,int strip,const Float_t zdcsmd) {
465  mZDCSMD[eastwest][verthori][strip-1] = (zdcsmd >0.)? zdcsmd:0.;}
466 
467 inline void StFlowEvent::SetPid(const Char_t* pid) {
468  strncpy(mPid, pid, 9); mPid[9] = '\0'; }
469 
470 inline void StFlowEvent::SetProbPid() { mProbPid = kTRUE; }
471 
472 inline void StFlowEvent::SetEtaSubs() { mEtaSubs = kTRUE; }
473 
474 inline void StFlowEvent::SetRanSubs() { mRanSubs = kTRUE; }
475 
476 inline void StFlowEvent::SetFirstLastPhiWgt() { mFirstLastPhiWgt = kTRUE; }
477 
478 inline void StFlowEvent::SetFirstLastPoints() { mFirstLastPoints = kTRUE; }
479 
480 inline void StFlowEvent::SetUseZDCSMD(Bool_t UseZDCSMD) { mUseZDCSMD = UseZDCSMD; }
481 
482 inline void StFlowEvent::SetPtWgtSaturation(Float_t PtWgtSaturation) { mPtWgtSaturation = PtWgtSaturation; }
483 
484 inline void StFlowEvent::SetPtWgt(Bool_t PtWgt) { mPtWgt = PtWgt; }
485 
486 inline void StFlowEvent::SetEtaWgt(Bool_t EtaWgt) { mEtaWgt = EtaWgt; }
487 
488 inline void StFlowEvent::SetV1TPCDetctWgtG_Mix(Float_t val, Int_t selN){
489  mV1TPCDetctWgtG_Mix[selN]=val;}
490 inline void StFlowEvent::SetV1FtpcEastDetctWgtG_Mix(Float_t val, Int_t selN){
491  mV1FtpcEastDetctWgtG_Mix[selN]=val;}
492 inline void StFlowEvent::SetV1FtpcWestDetctWgtG_Mix(Float_t val, Int_t selN){
493  mV1FtpcWestDetctWgtG_Mix[selN]=val;}
494 
495 inline void StFlowEvent::SetV2TPCDetctWgtG_Mix(Float_t val, Int_t selN){
496  mV2TPCDetctWgtG_Mix[selN]=val;}
497 inline void StFlowEvent::SetV2FtpcEastDetctWgtG_Mix(Float_t val, Int_t selN){
498  mV2FtpcEastDetctWgtG_Mix[selN]=val;}
499 inline void StFlowEvent::SetV2FtpcWestDetctWgtG_Mix(Float_t val, Int_t selN){
500  mV2FtpcWestDetctWgtG_Mix[selN]=val;}
501 
502 #endif
503 
505 //
506 // $Log: StFlowEvent.h,v $
507 // Revision 1.56 2009/11/24 19:23:03 posk
508 // Added reCenter option to remove acceptance correlations instead of phiWgt.
509 //
510 // Revision 1.55 2007/02/06 18:57:54 posk
511 // In Lee Yang Zeros method, introduced recentering of Q vector.
512 // Reactivated eta symmetry cut.
513 //
514 // Revision 1.54 2006/07/06 16:56:01 posk
515 // Calculation of v1 for selection=2 is done with mixed harmonics.
516 //
517 // Revision 1.53 2006/02/22 19:29:16 posk
518 // Additions needed for the StFlowLeeYangZerosMaker
519 //
520 // Revision 1.52 2005/02/10 21:04:57 aihong
521 // test mProbPid of StFlowEvent before launch calculation pid on fly
522 //
523 // Revision 1.51 2004/12/17 22:33:08 aihong
524 // add in full Psi weight for ZDC SMD and fix a few bugs, done by Gang
525 //
526 // Revision 1.50 2004/12/17 15:50:08 aihong
527 // check in v1{3} code
528 //
529 // Revision 1.49 2004/12/07 17:04:46 posk
530 // Eliminated the very old mOnePhiWgt, which used one phiWgt histogram for flttening
531 // instead of four.
532 //
533 // Revision 1.48 2004/11/16 21:22:22 aihong
534 // removed old cumulant method
535 //
536 // Revision 1.47 2004/05/05 21:13:45 aihong
537 // Gang's code for ZDC-SMD added
538 //
539 // Revision 1.46 2004/03/11 17:58:42 posk
540 // Added Random Subs analysis method.
541 //
542 // Revision 1.45 2003/07/30 22:00:40 oldi
543 // Eta cuts for event plane selection separated for FTPC east and west.
544 // PtWgtSaturation parameter introduced (default set to 2. -> no change of default behavior).
545 //
546 // Revision 1.44 2003/06/18 17:00:59 posk
547 // Event plane cuts now only odd and even, instead of different for each harmonic.
548 //
549 // Revision 1.43 2003/04/01 00:27:07 posk
550 // Little q is now unweighted by pt or eta. Big Q is unaffected.
551 //
552 // Revision 1.42 2003/01/10 16:42:15 oldi
553 // Several changes to comply with FTPC tracks:
554 // - Switch to include/exclude FTPC tracks introduced.
555 // The same switch changes the range of the eta histograms.
556 // - Eta symmetry plots for FTPC tracks added and separated from TPC plots.
557 // - PhiWgts and related histograms for FTPC tracks split in FarEast, East,
558 // West, FarWest (depending on vertex.z()).
559 // - Psi_Diff plots for 2 different selections and the first 2 harmonics added.
560 // - Cut to exclude mu-events with no primary vertex introduced.
561 // (This is possible for UPC events and FTPC tracks.)
562 // - Global DCA cut for FTPC tracks added.
563 // - Global DCA cuts for event plane selection separated for TPC and FTPC tracks.
564 // - Charge cut for FTPC tracks added.
565 //
566 // Revision 1.41 2003/01/08 19:26:48 posk
567 // PhiWgt hists sorted on sign of z of first and last points.
568 // Version 6 of pico file.
569 //
570 // Revision 1.40 2002/05/23 18:54:11 posk
571 // Moved centrality cuts into StFlowConstants
572 //
573 // Revision 1.39 2002/03/15 16:43:22 snelling
574 // Added a method to recalculate the centrality in StFlowPicoEvent
575 //
576 // Revision 1.38 2002/03/14 18:51:50 snelling
577 // Added new centralities
578 //
579 // Revision 1.37 2002/02/13 22:29:28 posk
580 // Pt Weight now also weights Phi Weights. Added Eta Weight, default=FALSE.
581 //
582 // Revision 1.36 2002/02/01 23:06:26 snelling
583 // Added entries for header information in flowPico (not everthing is available yet)
584 //
585 // Revision 1.35 2002/01/31 01:04:47 posk
586 // *** empty log message ***
587 //
588 // Revision 1.34 2001/12/18 19:22:09 posk
589 // "proton" and "antiproton" changed to "pr+" and "pr-".
590 // Compiles on Solaris.
591 //
592 // Revision 1.33 2001/12/11 21:33:50 posk
593 // Went from one to four sets of histograms for making the event plane isotropic.
594 // StFlowEvent::PhiWeight() has changed arguments and return value.
595 // The ptWgt saturates above 2 GeV/c.
596 //
597 // Revision 1.32 2001/11/09 21:10:42 posk
598 // Switched from CERNLIB to TMath. Little q is now normalized.
599 //
600 // Revision 1.31 2001/11/02 04:49:52 aihong
601 // add func. for cumulant maker
602 //
603 // Revision 1.30 2001/08/01 19:39:35 snelling
604 // Added the trigger word
605 //
606 // Revision 1.29 2001/07/27 01:26:14 snelling
607 // Added and changed variables for picoEvent. Changed trackCut class to StTrack
608 //
609 // Revision 1.28 2001/07/24 22:29:12 snelling
610 // First attempt to get a standard root pico file again, added variables
611 //
612 // Revision 1.27 2001/06/07 20:06:20 posk
613 // Global Dca cut for event plane particles.
614 // Removed SetPtWgt().
615 //
616 // Revision 1.26 2001/06/06 13:02:58 rcwells
617 // Added SetPtWgt(Bool_t) function to StFlowEvent
618 //
619 // Revision 1.25 2001/06/04 18:57:05 rcwells
620 // Adding filling from HbtEvents
621 //
622 // Revision 1.24 2001/05/23 18:11:09 posk
623 // Removed SetPids().
624 //
625 // Revision 1.23 2001/05/22 20:17:30 posk
626 // Now can do pseudorapidity subevents.
627 //
628 // Revision 1.22 2001/04/03 17:47:23 oldi
629 // Bug fix that excluded FTPC tracks from the determination of the reaction plane.
630 //
631 // Revision 1.21 2000/12/12 20:22:05 posk
632 // Put log comments at end of files.
633 // Deleted persistent StFlowEvent (old micro DST).
634 //
635 // Revision 1.20 2000/12/10 02:01:13 oldi
636 // A new member (StTrackTopologyMap mTopology) was added to StFlowPicoTrack.
637 // The evaluation of either a track originates from the FTPC or not is
638 // unambiguous now. The evaluation itself is easily extendible for other
639 // detectors (e.g. SVT+TPC). Old flowpicoevent.root files are treated as if
640 // they contain TPC tracks only (backward compatibility).
641 //
642 // Revision 1.19 2000/12/08 17:03:38 oldi
643 // Phi weights for both FTPCs included.
644 //
645 // Revision 1.18 2000/10/12 22:46:36 snelling
646 // Added support for the new pDST's and the probability pid method
647 //
648 // Revision 1.16 2000/09/15 22:51:30 posk
649 // Added pt weighting for event plane calcualtion.
650 //
651 // Revision 1.15 2000/09/05 16:11:32 snelling
652 // Added global DCA, electron and positron
653 //
654 // Revision 1.14 2000/08/31 18:58:22 posk
655 // For picoDST, added version number, runID, and multEta for centrality.
656 // Added centrality cut when reading picoDST.
657 // Added pt and eta selections for particles corr. wrt event plane.
658 //
659 // Revision 1.13 2000/08/09 21:38:23 snelling
660 // PID added
661 //
662 // Revision 1.12 2000/08/05 22:07:18 fisyak
663 // less restrictive selection for ROOTCINT
664 //
665 // Revision 1.11 2000/08/05 21:21:33 fisyak
666 // hide from CINT inline functions
667 //
668 // Revision 1.10 2000/08/04 21:03:45 perev
669 // Leaks + Clear() cleanup
670 //
671 // Revision 1.9 2000/06/30 14:48:32 posk
672 // Using MessageMgr, changed Eta Symmetry cut.
673 //
674 // Revision 1.8 2000/06/20 16:34:25 snelling
675 // fixed cout/streamer problem for mPhiWgt under Solaris
676 //
677 // Revision 1.7 2000/05/26 21:29:27 posk
678 // Protected Track data members from overflow.
679 //
680 // Revision 1.5 2000/05/16 20:59:30 posk
681 // Voloshin's flownanoevent.root added.
682 //
683 // Revision 1.4 2000/05/12 22:42:04 snelling
684 // Additions for persistency and minor fix
685 //
686 // Revision 1.2 2000/03/15 23:28:51 posk
687 // Added StFlowSelection.
688 //
689 // Revision 1.1 2000/03/02 23:02:50 posk
690 // Changed extensions from .hh and .cc to .h and .cxx .
691 //
692 // Revision 1.10 2000/02/29 22:00:54 posk
693 // Made SetPhiWeight inline, changed ImpactPar to Dca, etc.
694 //
695 // Revision 1.9 2000/02/18 22:49:55 posk
696 // Added PID and centrality.
697 //
698 // Revision 1.5 1999/12/15 22:01:26 posk
699 // Added StFlowConstants.hh
700 //
701 // Revision 1.4 1999/12/04 00:10:33 posk
702 // Works with the new StEvent
703 //
704 // Revision 1.3 1999/11/30 18:52:52 snelling
705 // First modification for the new StEvent
706 //
707 // Revision 1.2 1999/11/24 18:17:14 posk
708 // Put the methods which act on the data in with the data in StFlowEvent.
709 //
710 // Revision 1.1 1999/11/04 19:02:06 snelling
711 // First check in of StFlowMaker. It contains the common code from
712 // StFlowTagMaker and StFlowAnalysisMaker.
713 //