ROOT logo
/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
* See cxx source for full Copyright notice */
/* $Id$ */

/*****************************************************************
  AliFlowEventSimple: A simple event 
  for flow analysis                  
                                     
  origin: Naomi van der Kolk (kolk@nikhef.nl)           
          Ante Bilandzic     (anteb@nikhef.nl)         
          Raimond Snellings  (Raimond.Snellings@nikhef.nl)    
  mods:   Mikolaj Krzewicki  (mikolaj.krzewicki@cern.ch)
*****************************************************************/

#ifndef ALIFLOWEVENTSIMPLE_H
#define ALIFLOWEVENTSIMPLE_H

#include "TObject.h"
#include "TParameter.h"
#include "TMath.h"
class TTree;
class TF1;
class AliFlowVector;
class AliFlowTrackSimple;
class AliFlowTrackSimpleCuts;

class AliFlowEventSimple: public TObject {

 public:

  enum ConstructionMethod {kEmpty,kGenerate};

  AliFlowEventSimple();
  AliFlowEventSimple( Int_t nParticles,
                      ConstructionMethod m=kEmpty,
                      TF1* ptDist=NULL,
                      Double_t phiMin=0.0,
                      Double_t phiMax=TMath::TwoPi(),
                      Double_t etaMin=-1.0,
                      Double_t etaMax= 1.0 );
  AliFlowEventSimple(TTree* anInput, const AliFlowTrackSimpleCuts* rpCuts, const AliFlowTrackSimpleCuts* poiCuts);
  AliFlowEventSimple(const AliFlowEventSimple& anEvent);
  AliFlowEventSimple& operator=(const AliFlowEventSimple& anEvent);
  virtual  ~AliFlowEventSimple();

  Bool_t  IsFolder() const {return kTRUE;};
  void    Browse(TBrowser *b); 
  void    Print(Option_t* option = "") const;      //method to print stats
  
  Int_t    NumberOfTracks() const                   { return fNumberOfTracks; }
  Int_t    GetReferenceMultiplicity() const         { return fReferenceMultiplicity; }
  void     SetReferenceMultiplicity( Int_t m )      { fReferenceMultiplicity = m; }
  Int_t    GetEventNSelTracksRP() const             { return GetNumberOfPOIs(0); } 
  void     SetEventNSelTracksRP(Int_t nr)           { SetNumberOfPOIs(nr,0); }  
  Int_t    GetEventNSelTracksPOI() const            { return GetNumberOfPOIs(1); } 
  void     SetEventNSelTracksPOI(Int_t np)          { SetNumberOfPOIs(np,1); }  
  Int_t    GetNumberOfRPs() const                   { return GetNumberOfPOIs(0); }
  void     SetNumberOfRPs( Int_t nr )               { SetNumberOfPOIs(nr,0); }
  Int_t    GetNumberOfPOIs(Int_t i=1) const         { return (i<fNumberOfPOItypes)?fNumberOfPOIs[i]:0; }
  void     SetNumberOfPOIs( Int_t nubmerOfPOIs, Int_t poiType=1 );
  void     IncrementNumberOfPOIs(Int_t poiType=1);

  void     SetUseGlauberMCSymmetryPlanes()          { fUseGlauberMCSymmetryPlanes = kTRUE; }
  void     SetUseExternalSymmetryPlanes(TF1 *gPsi1Psi3 = 0x0,
					TF1 *gPsi2Psi4 = 0x0,
					TF1 *gPsi3Psi5 = 0x0);
  void     SetPsi1(Double_t gPsi1)                  { fPsi1 = gPsi1; }
  void     SetPsi2(Double_t gPsi2)                  { fPsi2 = gPsi2; }
  void     SetPsi3(Double_t gPsi3)                  { fPsi3 = gPsi3; }
  void     SetPsi4(Double_t gPsi4)                  { fPsi4 = gPsi4; }
  void     SetPsi5(Double_t gPsi5)                  { fPsi5 = gPsi5; }
  Double_t GetPsi1() const                          { return fPsi1; }
  Double_t GetPsi2() const                          { return fPsi2; }
  Double_t GetPsi3() const                          { return fPsi3; }
  Double_t GetPsi4() const                          { return fPsi4; }
  Double_t GetPsi5() const                          { return fPsi5; }

  Double_t GetMCReactionPlaneAngle() const          { return fMCReactionPlaneAngle; }
  void     SetMCReactionPlaneAngle(Double_t fPhiRP) { fMCReactionPlaneAngle=fPhiRP; fMCReactionPlaneAngleIsSet=kTRUE; }
  Bool_t   IsSetMCReactionPlaneAngle() const        { return fMCReactionPlaneAngleIsSet; }
  void     SetAfterBurnerPrecision(Double_t p)      { fAfterBurnerPrecision=p; }
  Double_t GetAfterBurnerPrecision() const          { return fAfterBurnerPrecision; }
  void     SetUserModified(Bool_t s=kTRUE)          { fUserModified=s; }
  Bool_t   IsUserModified() const                   { return fUserModified; }
  void     SetShuffleTracks(Bool_t b)               {fShuffleTracks=b;}
  void     ShuffleTracks();

  void ResolutionPt(Double_t res);
  void TagSubeventsInEta(Double_t etaMinA, Double_t etaMaxA, Double_t etaMinB, Double_t etaMaxB );
  void TagSubeventsByCharge();
  void TagRP(const AliFlowTrackSimpleCuts* cuts );
  void TagPOI(const AliFlowTrackSimpleCuts* cuts, Int_t poiType=1);
  void TagTracks(const AliFlowTrackSimpleCuts* cutsRP, const AliFlowTrackSimpleCuts* cutsPOI);
  void CloneTracks(Int_t n);
  void AddV1( Double_t v1 );
  void AddV2( Double_t v2 );
  void AddV3( Double_t v3 );
  void AddV4( Double_t v4 );
  void AddV5( Double_t v5 );
  void AddFlow( Double_t v1, Double_t v2, Double_t v3, Double_t v4, Double_t v5 );
  void AddFlow(Double_t v1, Double_t v2, Double_t v3, Double_t v4, Double_t v5,
               Double_t rp1, Double_t rp2, Double_t rp3, Double_t rp4, Double_t rp5 );
  void AddV2( TF1* ptDepV2 );
  void DefineDeadZone( Double_t etaMin, Double_t etaMax, Double_t phiMin, Double_t phiMax );
  Int_t CleanUpDeadTracks();
  virtual void ClearFast();
 
  static TF1* SimplePtSpectrum();
  static TF1* SimplePtDepV2();

  AliFlowTrackSimple* GetTrack(Int_t i);
  void AddTrack( AliFlowTrackSimple* track ); 
  void TrackAdded();
  AliFlowTrackSimple* MakeNewTrack();
 
  virtual AliFlowVector GetQ(Int_t n=2, TList *weightsList=NULL, Bool_t usePhiWeights=kFALSE, Bool_t usePtWeights=kFALSE, Bool_t useEtaWeights=kFALSE);
  virtual void Get2Qsub(AliFlowVector* Qarray, Int_t n=2, TList *weightsList=NULL, Bool_t usePhiWeights=kFALSE, Bool_t usePtWeights=kFALSE, Bool_t useEtaWeights=kFALSE);  

  void SetCentrality(Double_t c) {fCentrality=c;}
  Double_t GetCentrality() const {return fCentrality;}

 protected:
  virtual void Generate( Int_t nParticles,
                         TF1* ptDist=NULL,
                         Double_t phiMin=0.0,
                         Double_t phiMax=TMath::TwoPi(),
                         Double_t etaMin=-1.0,
                         Double_t etaMax= 1.0 );

  //data members
  TObjArray*              fTrackCollection;           //-> collection of tracks
  Int_t                   fReferenceMultiplicity;     // reference multiplicity
  Int_t                   fNumberOfTracks;            // number of tracks
  Bool_t                  fUseGlauberMCSymmetryPlanes;// Use symmetry planes (Glauber MC)
  Bool_t                  fUseExternalSymmetryPlanes; // Use symmetry planes (external)
  Double_t                fPsi1;                      // Psi_1
  Double_t                fPsi2;                      // Psi_2
  Double_t                fPsi3;                      // Psi_3
  Double_t                fPsi4;                      // Psi_4
  Double_t                fPsi5;                      // Psi_5
  TF1*                    fPsi1Psi3;                  // Correlation between Psi_1 and Psi_3
  TF1*                    fPsi2Psi4;                  // Correlation between Psi_2 and Psi_4
  TF1*                    fPsi3Psi5;                  // Correlation between Psi_3 and Psi_5
  Double_t                fMCReactionPlaneAngle;      // the angle of the reaction plane from the MC truth
  Bool_t                  fMCReactionPlaneAngleIsSet; // did we set it from MC?
  Double_t                fAfterBurnerPrecision;      // iteration precision in afterburner
  Bool_t                  fUserModified;              // did we modify the event in any way (afterburner etc) ?
  TParameter<Int_t>*      fNumberOfTracksWrap;        //! number of tracks in TBrowser
  TParameter<Int_t>*      fNumberOfRPsWrap;           //! number of tracks that have passed the RP selection in TBrowser
  TParameter<Int_t>*      fNumberOfPOIsWrap;          //! number of tracks that have passed the POI selection in TBrowser
  TParameter<Double_t>*   fMCReactionPlaneAngleWrap;  //! the angle of the reaction plane from the MC truth in TBrowser
  Int_t*                  fShuffledIndexes;           //! placeholder for randomized indexes
  Bool_t                  fShuffleTracks;             // do we shuffle tracks on get?
  TObjArray*              fMothersCollection;         //!cache the particles with daughters
  Double_t                fCentrality;                //centrality
 
 private:
  Int_t                   fNumberOfPOItypes;    // how many different flow particle types do we have? (RP,POI,POI_2,...)
  Int_t*                  fNumberOfPOIs;          //[fNumberOfPOItypes] number of tracks that have passed the POI selection

  ClassDef(AliFlowEventSimple,1)
};

#endif


 AliFlowEventSimple.h:1
 AliFlowEventSimple.h:2
 AliFlowEventSimple.h:3
 AliFlowEventSimple.h:4
 AliFlowEventSimple.h:5
 AliFlowEventSimple.h:6
 AliFlowEventSimple.h:7
 AliFlowEventSimple.h:8
 AliFlowEventSimple.h:9
 AliFlowEventSimple.h:10
 AliFlowEventSimple.h:11
 AliFlowEventSimple.h:12
 AliFlowEventSimple.h:13
 AliFlowEventSimple.h:14
 AliFlowEventSimple.h:15
 AliFlowEventSimple.h:16
 AliFlowEventSimple.h:17
 AliFlowEventSimple.h:18
 AliFlowEventSimple.h:19
 AliFlowEventSimple.h:20
 AliFlowEventSimple.h:21
 AliFlowEventSimple.h:22
 AliFlowEventSimple.h:23
 AliFlowEventSimple.h:24
 AliFlowEventSimple.h:25
 AliFlowEventSimple.h:26
 AliFlowEventSimple.h:27
 AliFlowEventSimple.h:28
 AliFlowEventSimple.h:29
 AliFlowEventSimple.h:30
 AliFlowEventSimple.h:31
 AliFlowEventSimple.h:32
 AliFlowEventSimple.h:33
 AliFlowEventSimple.h:34
 AliFlowEventSimple.h:35
 AliFlowEventSimple.h:36
 AliFlowEventSimple.h:37
 AliFlowEventSimple.h:38
 AliFlowEventSimple.h:39
 AliFlowEventSimple.h:40
 AliFlowEventSimple.h:41
 AliFlowEventSimple.h:42
 AliFlowEventSimple.h:43
 AliFlowEventSimple.h:44
 AliFlowEventSimple.h:45
 AliFlowEventSimple.h:46
 AliFlowEventSimple.h:47
 AliFlowEventSimple.h:48
 AliFlowEventSimple.h:49
 AliFlowEventSimple.h:50
 AliFlowEventSimple.h:51
 AliFlowEventSimple.h:52
 AliFlowEventSimple.h:53
 AliFlowEventSimple.h:54
 AliFlowEventSimple.h:55
 AliFlowEventSimple.h:56
 AliFlowEventSimple.h:57
 AliFlowEventSimple.h:58
 AliFlowEventSimple.h:59
 AliFlowEventSimple.h:60
 AliFlowEventSimple.h:61
 AliFlowEventSimple.h:62
 AliFlowEventSimple.h:63
 AliFlowEventSimple.h:64
 AliFlowEventSimple.h:65
 AliFlowEventSimple.h:66
 AliFlowEventSimple.h:67
 AliFlowEventSimple.h:68
 AliFlowEventSimple.h:69
 AliFlowEventSimple.h:70
 AliFlowEventSimple.h:71
 AliFlowEventSimple.h:72
 AliFlowEventSimple.h:73
 AliFlowEventSimple.h:74
 AliFlowEventSimple.h:75
 AliFlowEventSimple.h:76
 AliFlowEventSimple.h:77
 AliFlowEventSimple.h:78
 AliFlowEventSimple.h:79
 AliFlowEventSimple.h:80
 AliFlowEventSimple.h:81
 AliFlowEventSimple.h:82
 AliFlowEventSimple.h:83
 AliFlowEventSimple.h:84
 AliFlowEventSimple.h:85
 AliFlowEventSimple.h:86
 AliFlowEventSimple.h:87
 AliFlowEventSimple.h:88
 AliFlowEventSimple.h:89
 AliFlowEventSimple.h:90
 AliFlowEventSimple.h:91
 AliFlowEventSimple.h:92
 AliFlowEventSimple.h:93
 AliFlowEventSimple.h:94
 AliFlowEventSimple.h:95
 AliFlowEventSimple.h:96
 AliFlowEventSimple.h:97
 AliFlowEventSimple.h:98
 AliFlowEventSimple.h:99
 AliFlowEventSimple.h:100
 AliFlowEventSimple.h:101
 AliFlowEventSimple.h:102
 AliFlowEventSimple.h:103
 AliFlowEventSimple.h:104
 AliFlowEventSimple.h:105
 AliFlowEventSimple.h:106
 AliFlowEventSimple.h:107
 AliFlowEventSimple.h:108
 AliFlowEventSimple.h:109
 AliFlowEventSimple.h:110
 AliFlowEventSimple.h:111
 AliFlowEventSimple.h:112
 AliFlowEventSimple.h:113
 AliFlowEventSimple.h:114
 AliFlowEventSimple.h:115
 AliFlowEventSimple.h:116
 AliFlowEventSimple.h:117
 AliFlowEventSimple.h:118
 AliFlowEventSimple.h:119
 AliFlowEventSimple.h:120
 AliFlowEventSimple.h:121
 AliFlowEventSimple.h:122
 AliFlowEventSimple.h:123
 AliFlowEventSimple.h:124
 AliFlowEventSimple.h:125
 AliFlowEventSimple.h:126
 AliFlowEventSimple.h:127
 AliFlowEventSimple.h:128
 AliFlowEventSimple.h:129
 AliFlowEventSimple.h:130
 AliFlowEventSimple.h:131
 AliFlowEventSimple.h:132
 AliFlowEventSimple.h:133
 AliFlowEventSimple.h:134
 AliFlowEventSimple.h:135
 AliFlowEventSimple.h:136
 AliFlowEventSimple.h:137
 AliFlowEventSimple.h:138
 AliFlowEventSimple.h:139
 AliFlowEventSimple.h:140
 AliFlowEventSimple.h:141
 AliFlowEventSimple.h:142
 AliFlowEventSimple.h:143
 AliFlowEventSimple.h:144
 AliFlowEventSimple.h:145
 AliFlowEventSimple.h:146
 AliFlowEventSimple.h:147
 AliFlowEventSimple.h:148
 AliFlowEventSimple.h:149
 AliFlowEventSimple.h:150
 AliFlowEventSimple.h:151
 AliFlowEventSimple.h:152
 AliFlowEventSimple.h:153
 AliFlowEventSimple.h:154
 AliFlowEventSimple.h:155
 AliFlowEventSimple.h:156
 AliFlowEventSimple.h:157
 AliFlowEventSimple.h:158
 AliFlowEventSimple.h:159
 AliFlowEventSimple.h:160
 AliFlowEventSimple.h:161
 AliFlowEventSimple.h:162
 AliFlowEventSimple.h:163
 AliFlowEventSimple.h:164
 AliFlowEventSimple.h:165
 AliFlowEventSimple.h:166