]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliV0ReaderV1.h
added new addtask + major modifications for pPb
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliV0ReaderV1.h
1 #ifndef ALIV0READERV1_H
2 #define ALIV0READERV1_H
3
4 #include "AliAnalysisTaskSE.h"
5 #include "AliAODv0.h"
6 #include "AliESDv0.h"
7 #include "AliConversionCuts.h"
8 #include "AliExternalTrackParam.h"
9 #include "TObject.h"
10 #include "AliMCEvent.h"   // for CF
11 #include "AliESDEvent.h"
12 #include "AliKFParticle.h"
13 #include "TParticle.h"
14 #include <vector>
15 #include "AliESDpid.h"
16 #include "TF1.h"
17 #include "TRandom3.h"
18 #include "AliAnalysisManager.h"
19
20 class AliConversionPhotonBase;
21 class TRandom3;
22 class AliStack;
23 class TList;
24 class AliKFConversionPhoton;
25 class TString;
26 class TClonesArray;
27 class TH1F;
28 class TH2F;
29 class AliAODConversionPhoton;
30
31 using namespace std;
32
33 class AliV0ReaderV1 : public AliAnalysisTaskSE {
34
35  public:
36
37   AliV0ReaderV1(const char *name="V0ReaderV1");
38   virtual ~AliV0ReaderV1();                            //virtual destructor
39   void UserCreateOutputObjects();
40   virtual Bool_t Notify();
41   virtual void UserExec(Option_t *option);
42   virtual void Terminate(Option_t *);
43   virtual void Init();
44
45   Bool_t ProcessEvent(AliVEvent *inputEvent,AliMCEvent *mcEvent=NULL);
46   Bool_t IsEventSelected(){return fEventIsSelected;}
47
48   // Return Reconstructed Gammas
49   TClonesArray *GetReconstructedGammas(){return fConversionGammas;}
50   Int_t GetNReconstructedGammas(){if(fConversionGammas){return fConversionGammas->GetEntriesFast();}else{return 0;}}
51   AliConversionCuts *GetConversionCuts(){return fConversionCuts;}
52   TList *GetCutHistograms(){if(fConversionCuts){return fConversionCuts->GetCutHistograms();}return NULL;}
53   // Set Options
54
55   void SetConversionCuts(const TString cut);
56   void SetConversionCuts(AliConversionCuts *cuts){fConversionCuts=cuts;}
57   void SetUseOwnXYZCalculation(Bool_t flag){fUseOwnXYZCalculation=flag;}
58   void SetUseConstructGamma(Bool_t flag){fUseConstructGamma=flag;}
59   void SetUseAODConversionPhoton(Bool_t b){if(b){cout<<"Setting Outputformat to AliAODConversionPhoton "<<endl;}else{cout<<"Setting Outputformat to AliKFConversionPhoton "<<endl;};kUseAODConversionPhoton=b;}
60   void SetCreateAODs(Bool_t k=kTRUE){fCreateAOD=k;}
61   void SetDeltaAODFilename(TString s){fDeltaAODFilename=s;}
62   void SetDeltaAODBranchName(TString string) { fDeltaAODBranchName = string;AliInfo(Form("Set DeltaAOD BranchName to: %s",fDeltaAODBranchName.Data()));}
63   TString GetPeriodName(){return fPeriodName;}
64   
65 protected:
66     // Reconstruct Gammas
67     Bool_t ProcessESDV0s();
68     AliKFConversionPhoton *ReconstructV0(AliESDv0* fCurrentV0,Int_t currentV0Index);
69     void FillAODOutput();
70     void FindDeltaAODBranchName();
71     Bool_t GetAODConversionGammas();
72
73     // Getter Functions
74
75     const AliExternalTrackParam *GetExternalTrackParam(AliESDv0 *fCurrentV0,Int_t &tracklabel,Int_t charge);
76     const AliExternalTrackParam *GetExternalTrackParamP(AliESDv0 *fCurrentV0,Int_t &tracklabel){return GetExternalTrackParam(fCurrentV0,tracklabel,1);};
77     const AliExternalTrackParam *GetExternalTrackParamN(AliESDv0 *fCurrentV0,Int_t &tracklabel){return GetExternalTrackParam(fCurrentV0,tracklabel,-1);};
78     AliKFParticle *GetPositiveKFParticle(AliAODv0 *fCurrentV0,Int_t fTrackLabel[2]);
79     AliKFParticle *GetNegativeKFParticle(AliAODv0 *fCurrentV0,Int_t fTrackLabel[2]);
80     AliKFParticle *GetPositiveKFParticle(AliESDv0 *fCurrentV0,Int_t fTrackLabel[2]);
81     AliKFParticle *GetNegativeKFParticle(AliESDv0 *fCurrentV0,Int_t fTrackLabel[2]);
82
83     Bool_t GetConversionPoint(const AliExternalTrackParam *pparam,const AliExternalTrackParam *nparam,Double_t convpos[3],Double_t dca[2]);
84     Bool_t GetHelixCenter(const AliExternalTrackParam *track,Double_t center[2]);
85     Double_t GetPsiPair(const AliESDv0* v0, const AliExternalTrackParam *positiveparam,const AliExternalTrackParam *negativeparam) const;
86
87     AliConversionCuts *fConversionCuts; // Pointer to the ConversionCut Selection
88     TClonesArray *fConversionGammas; // TClonesArray holding the reconstructed photons
89     Bool_t fUseImprovedVertex; // set flag to improve primary vertex estimation by adding photons
90     Bool_t fUseOwnXYZCalculation; //flag that determines if we use our own calculation of xyz (markus)
91     Bool_t fUseConstructGamma; //flag that determines if we use ConstructGamma method from AliKF
92     Bool_t kUseAODConversionPhoton; // set flag to use AOD instead of KF output format for photons
93     Bool_t fCreateAOD; // set flag for AOD creation
94     TString     fDeltaAODBranchName;// File where Gamma Conv AOD is located, if not in default AOD
95     TString fDeltaAODFilename; // set filename for delta/satellite aod
96     Bool_t fEventIsSelected;
97     TString fPeriodName;
98     
99 private:
100     AliV0ReaderV1(AliV0ReaderV1 &original);
101     AliV0ReaderV1 &operator=(const AliV0ReaderV1 &ref);
102
103
104     ClassDef(AliV0ReaderV1,2)
105 };
106
107 inline void AliV0ReaderV1::SetConversionCuts(const TString cut){
108     if(fConversionCuts != NULL){
109         delete fConversionCuts;
110         fConversionCuts=NULL;
111     }
112     if(fConversionCuts == NULL){
113         fConversionCuts=new AliConversionCuts("V0ReaderCuts","V0ReaderCuts");
114         fConversionCuts->InitializeCutsFromCutString(cut.Data());
115     }
116 }
117
118
119 #endif