]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero.h
remaining histograms with hardcoded binning parameters fixed
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero.h
1 #ifndef ALIANALYSISTASKNEUTRALMESONTOPIPLPIMIPIZERO_H
2 #define ALIANALYSISTASKNEUTRALMESONTOPIPLPIMIPIZERO_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7
8 #include "AliAnalysisTaskSE.h"
9 #include "AliV0ReaderV1.h"
10 #include "AliKFConversionPhoton.h"
11 #include "AliPrimaryPionSelector.h"
12 #include "AliConversionMesonCuts.h"
13 #include "AliConvEventCuts.h"
14 #include "AliCaloPhotonCuts.h"
15 #include "AliGammaConversionAODBGHandler.h"
16 #include "TProfile2D.h"
17
18 class AliESDInputHandler;
19 class AliMCEventHandler;
20 class AliESDEvent;
21 class AliESDtrack;
22 class AliESDtrackCuts;
23 class AliESDpidCuts;
24 class AliTriggerAnalysis;
25
26 class AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero: public AliAnalysisTaskSE
27 {
28         public:
29
30                 AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero();
31                 AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero( const char* name );
32                 virtual ~AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero();
33
34                 virtual void UserExec(Option_t *);
35                 virtual void UserCreateOutputObjects();
36                 virtual Bool_t Notify();
37                 virtual void Terminate(const Option_t *);
38
39                         
40                 void SetMoveParticleAccordingToVertex(Bool_t flag){fMoveParticleAccordingToVertex = flag;}
41                 void SetIsHeavyIon(Int_t flag){
42                         if (flag == 1 || flag ==2 ){
43                                 fIsHeavyIon = 1;    
44                         } else {
45                                 fIsHeavyIon = 0;    
46                         }
47                 }
48                 
49                 void SetIsMC(Bool_t isMC){fIsMC=isMC;}
50                 void SetEventCutList(Int_t nCuts, TList *CutArray){ 
51                         fnCuts= nCuts;  
52                         fEventCutArray = CutArray;
53                 }
54                 void SetConversionCutList(TList *CutArray){ fGammaCutArray = CutArray;}
55                 void SetClusterCutList(TList *CutArray){ fClusterCutArray = CutArray;}
56                 void SetPionCutList(TList *CutArray){ fPionCutArray = CutArray;}
57                 void SetNeutralPionCutList(TList *CutArray){ fNeutralPionMesonCutArray = CutArray; }
58                 void SetMesonCutList(TList *CutArray){ fMesonCutArray = CutArray; }
59                 void SetDoMesonQA(Bool_t flag){ fDoMesonQA = flag; }
60                 void SetNeutralPionMode(Int_t mode){fNeutralPionMode = mode; }
61         
62
63         private:
64
65                 void InitBack();
66                 
67                 // routines for photon selection from conversions
68                 void ProcessConversionPhotonCandidates();
69                 void ProcessTrueConversionPhotonCandidates(AliAODConversionPhoton*);
70                 
71                 // routines for photon selection from clusters
72                 void ProcessCaloPhotonCandidates();
73                 void ProcessTrueCaloPhotonCandidates(AliAODConversionPhoton *TruePhotonCandidate);
74                 
75                 void ProcessTrueMesonCandidates(AliAODConversionMother *Pi0Candidate, AliAODConversionMother *TrueNeutralPionCandidate, AliAODConversionPhoton *TrueVirtualGammaCandidate);
76                 void MoveParticleAccordingToVertex(AliAODConversionMother* particle,const AliGammaConversionAODBGHandler::GammaConversionVertex *vertex);
77         
78                 // routines for neutral pion candidates from pure conversion
79                 void ProcessNeutralPionCandidatesPureConversions();     
80                 void ProcessTrueNeutralPionCandidatesPureConversions(AliAODConversionMother *Pi0Candidate, AliAODConversionPhoton *TrueGammaCandidate0, AliAODConversionPhoton *TrueGammaCandidate1);
81                 void ProcessTrueNeutralPionCandidatesPureConversionsAOD(AliAODConversionMother *Pi0Candidate, AliAODConversionPhoton *TrueGammaCandidate0, AliAODConversionPhoton *TrueGammaCandidate1);
82                 
83                 // routines for neutral pion candidates from pure calo
84                 void ProcessNeutralPionCandidatesPureCalo();
85                 void ProcessTrueNeutralPionCandidatesPureCalo(AliAODConversionMother *Pi0Candidate, AliAODConversionPhoton *TrueGammaCandidate0, AliAODConversionPhoton *TrueGammaCandidate1);
86
87                 // routines for neutral pion candidates from mixed conv + calo
88                 void ProcessNeutralPionCandidatesMixedConvCalo();
89                 void ProcessTrueNeutralPionCandidatesMixedConvCalo( AliAODConversionMother *Pi0Candidate, AliAODConversionPhoton *TrueGammaCandidate0, AliAODConversionPhoton *TrueGammaCandidate1);
90                 
91                 void ProcessPionCandidates();
92                 void ProcessMCParticles();
93                 void CalculateMesonCandidates();
94         void CalculateBackground();
95                 void UpdateEventByEventData();
96         
97                 Bool_t IsPiPlPiMiPiZeroDecay(TParticle *fMCMother) const;
98                 Bool_t IsEtaPiPlPiMiPiZeroDaughter( Int_t label ) const;
99                 Bool_t IsOmegaPiPlPiMiPiZeroDaughter( Int_t label ) const;
100                 Bool_t GammaIsNeutralMesonPiPlPiMiPiZeroDaughter( Int_t label ) const;
101
102                 AliV0ReaderV1                                   *fV0Reader;                                                                     // V0Reader for basic conversion photon selection
103                 AliPrimaryPionSelector                  *fPionSelector;                                                         // primary charged pion selector, basic selection of pi+,pi-
104                 AliGammaConversionAODBGHandler  **fBGHandler;                                                           // BG handler
105                 AliESDEvent                                     *fESDEvent;                                                                     // current event
106                 AliMCEvent                                              *fMCEvent;                                                                      // current MC event
107                 AliStack                                                *fMCStack;                                                                      // current MC stack
108                 TList                                                   **fCutFolder;                                                           // list of output folders with main cut name 
109                 TList                                                   **fESDList;                                                                     // list with main output histograms for data
110                 TList                                                   **fBackList;                                                            // list with THnSparseF for BG 
111                 TList                                                   **fMotherList;                                                          // list with THnSparseF for FG 
112                 TList                                                   **fTrueList;                                                            // list with validated reconstructed MC histograms
113                 TList                                                   **fMCList;                                                                      // list with pure MC histograms
114                 TList                                                   *fOutputContainer;                                                      // output container
115                 TClonesArray                                    *fReaderGammas;                                                         // array with photon from fV0Reader
116                 vector<Int_t>                                   fSelectorNegPionIndex;                                          // array with pion indices for negative pions from fPionSelector
117                 vector<Int_t>                                   fSelectorPosPionIndex;                                          // array with pion indices for positive pions from fPionSelector
118                 TList                                                   *fGoodConvGammas;                                                       // good conv gammas after selection
119                 TList                                                   *fClusterCandidates;                                            //! good calo gammas after selection 
120                 TList                                                   *fNeutralPionCandidates;                                        // good neutral pion candidates
121                 TList                                                   *fGoodVirtualParticles;                                         // combination of pi+pi- candidates
122                 TList                                                   *fEventCutArray;                                                        // array with event cuts
123                 TList                                                   *fGammaCutArray;                                                        // array with Conversion Cuts
124                 TList                                                   *fClusterCutArray;                                                      // array with Cluster Cuts
125                 TList                                                   *fPionCutArray;                                                         // array with charged pion cuts
126                 TList                                                   *fNeutralPionMesonCutArray;                                     // array with neutral pion cuts
127                 TList                                                   *fMesonCutArray;                                                        // array with neutral meson cuts
128                 AliConvEventCuts                                *fEventCuts;                                                            // current event cuts
129                 AliConversionPhotonCuts                 *fConversionCuts;                                                       // current conversion cuts
130                 AliCaloPhotonCuts                               *fClusterCuts;                                                          // current cluster cuts
131                 
132                 // reconstructed particles
133                 TH1F                                                    **fHistoConvGammaPt;                                            // array of histos of conversion photon, pt
134                 TH1F                                                    **fHistoConvGammaEta;                                           // array of histos of conversion photon, eta
135                 TH1F                                                    **fHistoClusterGammaPt;                                         // array of histos of Cluster photon, pt
136                 TH1F                                                    **fHistoClusterGammaEta;                                        // array of histos of Cluster photon, eta
137                 TH1F                                                    **fHistoNegPionPt;                                                      // array of histos of negative pion, pt
138                 TH1F                                                    **fHistoPosPionPt;                                                      // array of histos of positive pion, pt
139                 TH1F                                                    **fHistoNegPionPhi;                                                     // array of histos of negative pion, phi
140                 TH1F                                                    **fHistoPosPionPhi;                                                     // array of histos of positive pion, phi
141                 TH1F                                                    **fHistoNegPionEta;                                                     // array of histos of negative pion, eta
142                 TH1F                                                    **fHistoPosPionEta;                                                     // array of histos of positive pion, eta
143                 TH2F                                                    **fHistoNegPionClsTPC;                                          // array of histos of negative pion, findable tpc cluster, pT
144                 TH2F                                                    **fHistoPosPionClsTPC;                                          // array of histos of positive pion, findable tpc cluster, pT
145                 TH2F                                                    **fHistoPionDCAxy;                                                      // array of histos of pion, dca_xy, pT
146                 TH2F                                                    **fHistoPionDCAz;                                                       // array of histos of pion, dca_z, pT
147                 TH2F                                                    **fHistoPionTPCdEdxNSigma;                                      // array of histos of pion, p, TPC nSigma dEdx pion
148                 TH2F                                                    **fHistoPionTPCdEdx;                                            // array of histos of pion, p, TPC dEdx
149                 TH2F                                                    **fHistoPionPionInvMassPt;                                      // array of histos of pion pion, invMass, pT_{pi+pi-}
150                 TH2F                                                    **fHistoGammaGammaInvMassPt;                            // array of histos of gamma-gamma, invMass, pT_{gamma gamma}
151                 TH2F                                                    **fHistoMotherInvMassPt;                                        // array of histos of pi+pi-pi0 same event, invMass, pT_{pi+pi-pi0}
152                 THnSparseF                                              **fTHnSparseMotherInvMassPtZM;                          // array of THnSparseF of pi+pi-pi0 same event, invMass, pT_{pi+pi-pi0}, Z, M
153                 TH2F                                                    **fHistoMotherBackInvMassPt;                            // array of histos of pi+pi-pi0 mixed event, invMass, pT_{pi+pi-pi0}
154                 THnSparseF                                              **fTHnSparseMotherBackInvMassPtZM;                      // array of THnSparseF of pi+pi-pi0 mixed event, invMass, pT_{pi+pi-pi0}, Z, M
155                 
156                 // pure MC properties
157                 TH1F                                                    **fHistoMCAllGammaPt;                                           // array of histos of all produced gammas in the specified y range
158                 TH1F                                                    **fHistoMCConvGammaPt;                                          // array of histos of all converted gammas in the specified y range 
159                 TH1F                                                    **fHistoMCAllPosPionsPt;                                        // array of histos with all produced primary positive pions in the specified y range
160                 TH1F                                                    **fHistoMCAllNegPionsPt;                                        // array of histos with all produced primary negative pions in the specified y range
161                 TH1F                                                    **fHistoMCGammaFromNeutralMesonPt;                      // array of histos of all produced gammas from omega or eta via pi+pi-pi0 in the specified y range/
162                 TH1F                                                    **fHistoMCPosPionsFromNeutralMesonPt;           // array of histos of all produced positive pions from omega or eta via pi+pi-pi0 in the specified y range/
163                 TH1F                                                    **fHistoMCNegPionsFromNeutralMesonPt;           // array of histos of all produced negative pions from omega or eta via pi+pi-pi0 in the specified y range/
164                 TH1F                                                    **fHistoMCEtaPiPlPiMiPiZeroPt;                          // array of histos of produced etas via pi+pi-pi0 in the specified y range
165                 TH1F                                                    **fHistoMCEtaPiPlPiMiPiZeroInAccPt;                     // array of histos of produced etas via pi+pi-pi0 in the specified y range, with decay products in respective y, eta ranges 
166                 TH1F                                                    **fHistoMCOmegaPiPlPiMiPiZeroPt;                        // array of histos of produced omegas via pi+pi-pi0 in the specified y range
167                 TH1F                                                    **fHistoMCOmegaPiPlPiMiPiZeroInAccPt;           // array of histos of produced omegas via pi+pi-pi0 in the specified y range, with decay products in respective y, eta ranges 
168
169                 // reconstructed particles MC validated
170                 TH2F                                                    **fHistoTrueMotherPiPlPiMiPiZeroInvMassPt;      // histos with reconstructed validated eta or omega, inv mass, pT
171                 TH2F                                                    **fHistoTrueMotherGammaGammaInvMassPt;          // histos with reconstructed validated pi0, inv mass, pT
172                 TH2F                                                    **fHistoTrueMotherGammaGammaFromEtaInvMassPt;   // histos with reconstructed validated pi0, inv mass, pT
173                 TH2F                                                    **fHistoTrueMotherGammaGammaFromOmegaInvMassPt; // histos with reconstructed validated pi0, inv mass, pT
174                 TH1F                                                    **fHistoTrueConvGammaPt;                                        // histos with reconstructed validated conv gamma, pT
175                 TH1F                                                    **fHistoTrueConvGammaFromNeutralMesonPt;        // histos with reconstructed validated conv gamma from eta or omega via pi0, pT
176                 TH1F                                                    **fHistoTrueClusterGammaPt;                                     // histos with reconstructed validated cluster gamma, pT
177                 TH1F                                                    **fHistoTrueClusterGammaFromNeutralMesonPt;     // histos with reconstructed validated cluster gamma from eta or omega via pi0, pT
178                 TH1F                                                    **fHistoTruePosPionPt;                                          // histos with reconstructed validated positive pion, pT
179                 TH1F                                                    **fHistoTruePosPionFromNeutralMesonPt;          // histos with reconstructed validated positive pion from eta or omega, pT
180                 TH1F                                                    **fHistoTrueNegPionPt;                                          // histos with reconstructed validated negative pion, pT
181                 TH1F                                                    **fHistoTrueNegPionFromNeutralMesonPt;          // histos with reconstructed validated negative pion from eta or omega, pT
182                 TH2F                                                    **fHistoTruePionPionInvMassPt;                          // histos with reconstructed validated two pion, invariant mass, pT
183                 TH2F                                                    **fHistoTruePionPionFromSameMotherInvMassPt;// histos with reconstructed validated two pion from same mother, invariant mass, pT
184                 TH2F                                                    **fHistoTruePionPionFromEtaInvMassPt;           // histos with reconstructed validated two pion from eta , invariant mass, pT
185                 TH2F                                                    **fHistoTruePionPionFromOmegaInvMassPt;         // histos with reconstructed validated two pion from omega, invariant mass, pT
186                 // Event properties
187                 TH1I                                                    **fHistoNEvents;                                                        // histo for event counting
188                 TH1I                                                    **fHistoNGoodESDTracks;                                         // histo number of reconstructed primary tracks
189                 TProfile                                                **fProfileEtaShift;                                                     // profile for eta shift bookkeeping
190                         
191                 TRandom3                                                fRandom;                                                                        // random number
192                 Int_t                                                   fnCuts;                                                                         // number of cuts to be run in parallel
193                 Int_t                                                   fiCut;                                                                          // current cut
194                 Int_t                                                   fNumberOfESDTracks;                                                     // integer with number of primary tracks in this event
195                 Bool_t                                                  fMoveParticleAccordingToVertex;                         // Flag to move parice to the vertex
196                 Int_t                                                   fIsHeavyIon;                                                            // Flag for collision system 0: pp, 1: PbPb, 2: pPb
197                 Bool_t                                                  fDoMesonAnalysis;                                                       // Flag for switching on meson analysis
198                 Bool_t                                                  fDoMesonQA;                                                                     // Flag for switching on small meson QA
199                 Bool_t                                                  fIsFromMBHeader;                                                        // Flag for particle whether it belongs to accepted header
200                 Bool_t                                                  fIsMC;                                                                          // Flag for MC  
201                 Int_t                                                   fNeutralPionMode;                                                       // Flag how neutral pion is reconstructed 0=PCM-PCM, 1=PCM-Calo, 2=Calo-Calo
202
203         private:
204                 AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero( const AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero& ); // Not implemented
205                 AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero& operator=( const AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero& ); // Not implemented
206
207                 ClassDef( AliAnalysisTaskNeutralMesonToPiPlPiMiPiZero, 2 );
208 };
209
210 #endif // ALIANALYSISTASKNEUTRALMESONTOPIPLPIMIPIZERO_H
211