]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliAnalysisTaskPhiCorrelations.h
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskPhiCorrelations.h
1 #ifndef AliAnalysisTaskPhiCorrelations_H
2 #define AliAnalysisTaskPhiCorrelations_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 //
9 // Analysis class for Underlying Event studies w.r.t. leading track
10 //
11 // Look for correlations on the tranverse regions w.r.t
12 // the leading track in the event
13 //
14 // This class needs input AODs.
15 // The output is a list of analysis-specific containers.
16 //
17 // The AOD can be either connected to the InputEventHandler  
18 // for a chain of AOD files 
19 // or 
20 // to the OutputEventHandler
21 // for a chain of ESD files,
22 // in this case the class should be in the train after the jet-finder
23 //
24 //    Authors:
25 //    Jan Fiete Grosse-Oetringhaus
26 // 
27 ////////////////////////////////////////////////////////////////////////
28
29 #include "AliAnalysisTask.h"
30 #include "AliUEHist.h"
31 #include "TString.h"
32 #include "AliVParticle.h"
33 #include "AliLog.h"
34 #include "THn.h" // in cxx file causes .../THn.h:257: error: conflicting declaration ‘typedef class THnT<float> THnF’
35
36 class AliAODEvent;
37 class AliAnalyseLeadingTrackUE;
38 class AliInputEventHandler;
39 class AliMCEvent;
40 class AliMCEventHandler;
41 class AliUEHistograms;
42 class AliVParticle;
43 class TH1;
44 class TObjArray;
45 class AliEventPoolManager;
46 class AliESDEvent;
47 class AliHelperPID;
48 class AliAnalysisUtils;
49 class TFormula;
50 class TMap;
51 class AliGenEventHeader;
52 class AliVEvent;
53
54
55 class  AliAnalysisTaskPhiCorrelations : public AliAnalysisTask
56   {
57   public:
58     AliAnalysisTaskPhiCorrelations(const char* name="AliAnalysisTaskPhiCorrelations");
59     virtual           ~AliAnalysisTaskPhiCorrelations();
60        
61       
62     // Implementation of interace methods
63     virtual     void   ConnectInputData(Option_t *);
64     virtual     void   CreateOutputObjects();
65     virtual     void   Exec(Option_t *option);
66
67     // Setters/Getters
68     // general configuration
69     virtual     void    SetDebugLevel( Int_t level )  { fDebug = level; }
70     virtual     void    SetMode(Int_t mode)           { fMode  = mode;  }
71     virtual     void    SetReduceMemoryFootprint(Bool_t flag) { fReduceMemoryFootprint = flag; }
72     virtual     void    SetEventMixing(Bool_t flag) { fFillMixed = flag; }
73     virtual     void    SetMixingTracks(Int_t tracks) { fMixingTracks = tracks; }
74     virtual     void    SetTwoTrackEfficiencyStudy(Bool_t flag) { fTwoTrackEfficiencyStudy = flag; }
75     virtual     void    SetTwoTrackEfficiencyCut(Float_t value = 0.02, Float_t min = 0.8) { fTwoTrackEfficiencyCut = value; fTwoTrackCutMinRadius = min; }
76     virtual     void    SetUseVtxAxis(Int_t flag) { fUseVtxAxis = flag; }
77     virtual     void    SetCourseCentralityBinning(Bool_t flag) { fCourseCentralityBinning = flag; }
78     virtual     void    SetSkipTrigger(Bool_t flag) { fSkipTrigger = flag; }
79     virtual     void    SetInjectedSignals(Bool_t flag) { fInjectedSignals = flag; }
80     void SetRandomizeReactionPlane(Bool_t flag) { fRandomizeReactionPlane = flag; }
81     
82     // histogram settings
83     void SetEfficiencyCorrectionTriggers(THnF* hist) { fEfficiencyCorrectionTriggers = hist; }
84     void SetEfficiencyCorrectionAssociated(THnF* hist) { fEfficiencyCorrectionAssociated = hist; }
85     void SetCentralityWeights(TH1* hist) { fCentralityWeights = hist; }
86
87     // for event QA
88     void   SetTracksInVertex( Int_t val ){ fnTracksVertex = val; }
89     void   SetZVertex( Double_t val )    { fZVertex = val; }
90     void   SetAcceptOnlyMuEvents( Bool_t val )    { fAcceptOnlyMuEvents = val; }
91     
92     // track cuts
93     void   SetTrackEtaCut( Double_t val )    { fTrackEtaCut = val; }
94     void   SetTrackEtaCutMin( Double_t val )    { fTrackEtaCutMin = val; }
95     void   SetOnlyOneEtaSide(Int_t flag)     { fOnlyOneEtaSide = flag; }
96     void   SetTrackPhiCutEvPlMin(Double_t val)  { fTrackPhiCutEvPlMin = val; }
97     void   SetTrackPhiCutEvPlMax(Double_t val)  { fTrackPhiCutEvPlMax = val; }
98     void   SetPtMin(Double_t val)            { fPtMin = val; }
99     void   SetFilterBit( UInt_t val )        { fFilterBit = val;  }
100     void   SetDCAXYCut(TFormula* value)      { fDCAXYCut = value; }
101     void   SetSharedClusterCut(Float_t value) { fSharedClusterCut = value; }
102     void   SetCrossedRowsCut(Int_t value)    { fCrossedRowsCut = value; }
103     void   SetFoundFractionCut(Double_t value) { fFoundFractionCut = value; }
104     void   SetTrackStatus(UInt_t status)     { fTrackStatus = status; }
105     void   SetCheckMotherPDG(Bool_t checkpdg) { fCheckMotherPDG = checkpdg; }
106    
107     // track cuts
108     void   SetTrackletDphiCut( Double_t val )    { fTrackletDphiCut = val; }
109
110     void   SetEventSelectionBit( UInt_t val )        { fSelectBit = val;  }
111     void   SetUseChargeHadrons( Bool_t val ) { fUseChargeHadrons = val; }
112     void   SetSelectParticleSpecies( Int_t trigger, Int_t associated ) { fParticleSpeciesTrigger = trigger; fParticleSpeciesAssociated = associated; }
113     void   SetSelectCharge(Int_t selectCharge) { fSelectCharge = selectCharge; }
114     void   SetSelectTriggerCharge(Int_t selectCharge) { fTriggerSelectCharge = selectCharge; }
115     void   SetSelectAssociatedCharge(Int_t selectCharge) { fAssociatedSelectCharge = selectCharge; }
116     void   SetTriggerRestrictEta(Float_t eta) { fTriggerRestrictEta = eta; }
117     void   SetEtaOrdering(Bool_t flag) { fEtaOrdering = flag; }
118     void   SetPairCuts(Bool_t conversions, Bool_t resonances) { fCutConversions = conversions; fCutResonances = resonances; }
119     void   SetRejectResonanceDaughters(Int_t value) { fRejectResonanceDaughters = value; }
120     void   SetCentralityMethod(const char* method) { fCentralityMethod = method; }
121     void   SetFillpT(Bool_t flag) { fFillpT = flag; }
122     void   SetStepsFillSkip(Bool_t step0, Bool_t step6) { fFillOnlyStep0 = step0; fSkipStep6 = step6; }
123     void   SetRejectCentralityOutliers(Bool_t flag = kTRUE) { fRejectCentralityOutliers = flag; }
124     void   SetRejectZeroTrackEvents(Bool_t flag)    { fRejectZeroTrackEvents = flag; }
125     void   SetRemoveWeakDecays(Bool_t flag = kTRUE) { fRemoveWeakDecays = flag; }
126     void   SetRemoveDuplicates(Bool_t flag = kTRUE) { fRemoveDuplicates = flag; }
127     void   SetSkipFastCluster(Bool_t flag = kTRUE)  { fSkipFastCluster = flag; }
128     void   SetWeightPerEvent(Bool_t flag = kTRUE)   { fWeightPerEvent = flag; }
129     void   SetCustomBinning(const char* binningStr) { fCustomBinning = binningStr; }
130     void   SetPtOrder(Bool_t flag) { fPtOrder = flag; }
131     void   SetTriggersFromDetector(Int_t flag) { fTriggersFromDetector = flag; }
132     void   SetAssociatedFromDetector(Int_t flag) { fAssociatedFromDetector = flag; }
133     void   SetMCUseUncheckedCentrality(Bool_t flag) { fMCUseUncheckedCentrality = flag; }
134     
135     AliHelperPID* GetHelperPID() { return fHelperPID; }
136     void   SetHelperPID(AliHelperPID* pid){ fHelperPID = pid; }
137    
138     AliAnalysisUtils* GetAnalysisUtils() { return fAnalysisUtils; }
139     void   SetAnalysisUtils(AliAnalysisUtils* utils){ fAnalysisUtils = utils; }
140    
141     TMap* GetMap() { return fMap; }
142     void   SetMap(TMap* map){ fMap = map; }
143     
144   private:
145     AliAnalysisTaskPhiCorrelations(const  AliAnalysisTaskPhiCorrelations &det);
146     AliAnalysisTaskPhiCorrelations&   operator=(const  AliAnalysisTaskPhiCorrelations &det);
147     void            AddSettingsTree();                                  // add list of settings to output list
148     // Analysis methods
149     void            AnalyseCorrectionMode();                            // main algorithm to get correction maps
150     void            AnalyseDataMode();                                  // main algorithm to get raw distributions
151     void            Initialize();                                       // initialize some common pointer
152     TObjArray* CloneAndReduceTrackList(TObjArray* tracks);
153     void RemoveDuplicates(TObjArray* tracks);
154     void CleanUp(TObjArray* tracks, TObject* mcObj, Int_t maxLabel);
155     void SelectCharge(TObjArray* tracks);
156     AliGenEventHeader* GetFirstHeader();
157     Bool_t AcceptEventCentralityWeight(Double_t centrality);
158     void ShiftTracks(TObjArray* tracks, Double_t angle);
159     TObjArray* GetParticlesFromDetector(AliVEvent* inputEvent, Int_t idet);
160     Bool_t IsMuEvent();
161
162     // General configuration
163     Int_t               fDebug;           //  Debug flag
164     Int_t               fMode;            //  fMode = 0: data-like analysis 
165                                           //  fMode = 1: corrections analysis   
166     Bool_t              fReduceMemoryFootprint; // reduce memory consumption by writing less debug histograms
167     Bool_t              fFillMixed;             // enable event mixing (default: ON)
168     Int_t               fMixingTracks;          // size of track buffer for event mixing
169     Bool_t              fTwoTrackEfficiencyStudy; // two-track efficiency study on
170     Float_t             fTwoTrackEfficiencyCut;   // enable two-track efficiency cut
171     Float_t             fTwoTrackCutMinRadius;    // minimum radius for two-track efficiency cut
172     Int_t               fUseVtxAxis;              // use z vtx as axis (needs 7-10 times more memory!)
173     Bool_t              fCourseCentralityBinning; // less centrality bins
174     Bool_t              fSkipTrigger;             // skip trigger selection
175     Bool_t              fInjectedSignals;         // check header to skip injected signals in MC
176     Bool_t              fRandomizeReactionPlane;  // change the orientation of the RP by a random value by shifting all tracks
177     
178     AliHelperPID*     fHelperPID;      // points to class for PID
179     AliAnalysisUtils*     fAnalysisUtils;      // points to class with common analysis utilities
180     TMap*     fMap;                   // points to TMap class containing scaling factors for VZERO A signal
181  
182    // Pointers to external UE classes
183     AliAnalyseLeadingTrackUE*     fAnalyseUE;      //! points to class containing common analysis algorithms
184     AliUEHistograms*  fHistos;       //! points to class to handle histograms/containers  
185     AliUEHistograms*  fHistosMixed;       //! points to class to handle mixed histograms/containers  
186     
187     THnF* fEfficiencyCorrectionTriggers;     // if non-0 this efficiency correction is applied on the fly to the filling for trigger particles. The factor is multiplicative, i.e. should contain 1/efficiency. Axes: eta, pT, centrality, z-vtx
188     THnF* fEfficiencyCorrectionAssociated;   // if non-0 this efficiency correction is applied on the fly to the filling for associated particles. The factor is multiplicative, i.e. should contain 1/efficiency. Axes: eta, pT, centrality, z-vtx
189     TH1* fCentralityWeights;                 // for centrality flattening
190     
191     // Handlers and events
192     AliAODEvent*             fAOD;             //! AOD Event 
193     AliESDEvent*             fESD;             //! ESD Event 
194     TClonesArray*            fArrayMC;         //! Array of MC particles 
195     AliInputEventHandler*    fInputHandler;    //! Generic InputEventHandler 
196     AliMCEvent*              fMcEvent;         //! MC event
197     AliInputEventHandler*    fMcHandler;       //! MCEventHandler 
198     AliEventPoolManager*     fPoolMgr;         //! event pool manager
199     
200     // Histogram settings
201     TList*              fListOfHistos;    //  Output list of containers 
202     
203     // Event QA cuts
204     Int_t               fnTracksVertex;          // QA tracks pointing to principal vertex
205     Double_t            fZVertex;                 // Position of Vertex in Z direction
206     Bool_t              fAcceptOnlyMuEvents;   // Only Events with at least one muon are accepted
207     TString             fCentralityMethod;      // Method to determine centrality
208     
209     // Track cuts
210     Double_t            fTrackEtaCut;          // Maximum Eta cut on particles
211     Double_t            fTrackEtaCutMin;          // Minimum Eta cut on particles
212     Double_t            fTrackPhiCutEvPlMin;   // Minimum Phi cut on particles with respect to the Event Plane (values between 0 and Pi/2)
213     Double_t            fTrackPhiCutEvPlMax;   // Maximum Phi cut on particles with respect to the Event Plane (values between 0 and Pi/2), if = 0, then no cut is performed
214     Int_t               fOnlyOneEtaSide;       // decides that only trigger particle from one eta side are considered (0 = all; -1 = negative, 1 = positive)
215     Double_t            fPtMin;                // Min pT to start correlations
216     TFormula*           fDCAXYCut;             // additional pt dependent cut on DCA XY (only for AOD)
217     Double_t            fSharedClusterCut;  // cut on shared clusters (only for AOD)
218     Int_t               fCrossedRowsCut;   // cut on crossed rows (only for AOD)
219     Double_t            fFoundFractionCut;     // cut on crossed rows/findable clusters (only for AOD)
220     UInt_t              fFilterBit;            // Select tracks from an specific track cut 
221     UInt_t              fTrackStatus;          // if non-0, the bits set in this variable are required for each track
222     UInt_t              fSelectBit;            // Select events according to AliAnalysisTaskJetServices bit maps 
223     Bool_t              fUseChargeHadrons;     // Only use charge hadrons
224     Int_t               fParticleSpeciesTrigger; // Select which particle to use for the trigger [ -1 (all, default) 0 (pions) 1 (kaons) 2 (protons) 3 (others) particles ]
225     Int_t               fParticleSpeciesAssociated; // Select which particle to use for the associated [ -1 (all, default) 0 (pions) 1 (kaons) 2 (protons) 3 (others) particles ]
226     Bool_t             fCheckMotherPDG;     // Check the PDG code of mother for secondaries 
227
228     // Tracklets cuts
229     Double_t            fTrackletDphiCut;    //maximum Dphi cut on tracklets
230     
231     Int_t fSelectCharge;           // (un)like sign selection when building correlations: 0: no selection; 1: unlike sign; 2: like sign
232     Int_t fTriggerSelectCharge;    // select charge of trigger particle: 1: positive; -1 negative
233     Int_t fAssociatedSelectCharge; // select charge of associated particle: 1: positive; -1 negative
234     Float_t fTriggerRestrictEta;   // restrict eta range for trigger particle (default: -1 [off])
235     Bool_t fEtaOrdering;           // eta ordering, see AliUEHistograms.h for documentation
236     Bool_t fCutConversions;        // cut on conversions (inv mass)
237     Bool_t fCutResonances;         // cut on resonances (inv mass)
238     Int_t fRejectResonanceDaughters; // reject all daughters of all resonance candidates (1: test method (cut at m_inv=0.9); 2: k0; 3: lambda)
239     Bool_t fFillOnlyStep0;         // fill only step 0
240     Bool_t fSkipStep6;             // skip step 6 when filling
241     Bool_t fRejectCentralityOutliers;  // enable rejection of outliers in centrality vs no track correlation. Interferes with the event plane dependence code
242     Bool_t fRejectZeroTrackEvents;  // reject events which have no tracks (using the eta, pT cuts defined)
243     Bool_t fRemoveWeakDecays;      // remove secondaries from weak decays from tracks and particles
244     Bool_t fRemoveDuplicates;      // remove particles with the same label (double reconstruction)
245     Bool_t fSkipFastCluster;       // skip kFastOnly flagged events (only for data)
246     Bool_t fWeightPerEvent;        // weight with the number of trigger particles per event
247     TString fCustomBinning;        // supersedes default binning if set, see AliUEHist::GetBinning or AliUEHistograms::AliUEHistograms for syntax and examples
248     Bool_t fPtOrder;               // apply pT,a < pt,t condition; default: kTRUE
249     Int_t fTriggersFromDetector;   // 0 = tracks (default); 1 = VZERO_A; 2 = VZERO_C; 3 = SPD tracklets; 4 = forward muons
250     Int_t fAssociatedFromDetector;   // 0 = tracks (default); 1 = VZERO_A; 2 = VZERO_C; 3 = SPD tracklets; 4 = forward muons
251     Bool_t fMCUseUncheckedCentrality; // use unchecked centrality (only applies to MC); default: kFALSE
252     
253     Bool_t fFillpT;                // fill sum pT instead of number density
254     
255     ClassDef(AliAnalysisTaskPhiCorrelations, 47); // Analysis task for delta phi correlations
256   };
257
258 class AliDPhiBasicParticle : public AliVParticle
259 {
260   public:
261     AliDPhiBasicParticle(Float_t eta, Float_t phi, Float_t pt, Short_t charge)
262       : fEta(eta), fPhi(phi), fpT(pt), fCharge(charge)
263     {
264     }
265     ~AliDPhiBasicParticle() {}
266     
267     // kinematics
268     virtual Double_t Px() const { AliFatal("Not implemented"); return 0; }
269     virtual Double_t Py() const { AliFatal("Not implemented"); return 0; }
270     virtual Double_t Pz() const { AliFatal("Not implemented"); return 0; }
271     virtual Double_t Pt() const { return fpT; }
272     virtual Double_t P() const { AliFatal("Not implemented"); return 0; }
273     virtual Bool_t   PxPyPz(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
274
275     virtual Double_t Xv() const { AliFatal("Not implemented"); return 0; }
276     virtual Double_t Yv() const { AliFatal("Not implemented"); return 0; }
277     virtual Double_t Zv() const { AliFatal("Not implemented"); return 0; }
278     virtual Bool_t   XvYvZv(Double_t[3]) const { AliFatal("Not implemented"); return 0; }
279
280     virtual Double_t OneOverPt()  const { AliFatal("Not implemented"); return 0; }
281     virtual Double_t Phi()        const { return fPhi; }
282     virtual Double_t Theta()      const { AliFatal("Not implemented"); return 0; }
283
284
285     virtual Double_t E()          const { AliFatal("Not implemented"); return 0; }
286     virtual Double_t M()          const { AliFatal("Not implemented"); return 0; }
287     
288     virtual Double_t Eta()        const { return fEta; }
289     virtual Double_t Y()          const { AliFatal("Not implemented"); return 0; }
290     
291     virtual Short_t Charge()      const { return fCharge; }
292     virtual Int_t   GetLabel()    const { AliFatal("Not implemented"); return 0; }
293     // PID
294     virtual Int_t   PdgCode()     const { AliFatal("Not implemented"); return 0; }      
295     virtual const Double_t *PID() const { AliFatal("Not implemented"); return 0; }
296     
297     virtual Bool_t IsEqual(const TObject* obj) const { return (obj->GetUniqueID() == GetUniqueID()); }
298     
299     virtual void SetPhi(Double_t phi) { fPhi = phi; }
300     
301   private:
302     Float_t fEta;      // eta
303     Float_t fPhi;      // phi
304     Float_t fpT;       // pT
305     Short_t fCharge;   // charge
306     
307     ClassDef( AliDPhiBasicParticle, 1); // class which contains only quantities requires for this analysis to reduce memory consumption for event mixing
308 };
309   
310 #endif
311
312