]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliAnalysisTaskDiJetCorrelations.h
2+1 analysis from Raghava included in PWGCFCorrelationsDPhi libraries
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskDiJetCorrelations.h
1 #ifndef ALIANALYSISTASKDIJETCORRELATIONS_H
2 #define ALIANALYSISTASKDIJETCORRELATIONS_H
3
4 #include <TROOT.h>
5 #include <TSystem.h>
6 #include <TNtuple.h>
7 #include <TH1F.h>
8 #include <TH2F.h>
9 #include <TH3F.h>
10 #include <THnSparse.h>
11 #include <TArrayD.h>
12 #include <TParticle.h>
13 #include <TClonesArray.h>
14 #include <TList.h>
15 #include <TObjArray.h>
16
17 #include "AliAnalysisTaskSE.h"
18 #include "AliEventPoolManager.h"
19 #include "AliAODMCParticle.h"
20 #include "AliCentrality.h"
21
22 class AliAODEvent;
23 class AliVParticle;
24 class TObjArray;
25 class AliEventPoolManager;
26 class AliESDEvent;
27 class AliESDtrackCuts;
28 class AliMultiplicity;
29 class AliAODTrack;
30 class AliAODHandler;
31 class AliAODInputHandler;
32 class AliMCParticle;
33 class TDatabasePDG;
34 class TParticlePDG;
35 class AliMCEvent;
36 class AliStack;
37 class AliInputEventHandler;
38 class AliAnalysisTaskSE;
39
40
41 class AliAnalysisTaskDiJetCorrelations : public AliAnalysisTaskSE
42 {
43
44 public:
45     
46     AliAnalysisTaskDiJetCorrelations();
47     AliAnalysisTaskDiJetCorrelations(const char *name);
48     virtual ~AliAnalysisTaskDiJetCorrelations();
49     
50     virtual void     UserCreateOutputObjects();
51     //virtual void   SetInputCorrection();
52     //virtual void   Init();
53     //virtual void   LocalInit() {Init();}
54     virtual void    UserExec(Option_t *option);
55     virtual void    Terminate(Option_t *);
56     
57     virtual void    SetSEorME(Bool_t flag) {fMixedEvent = flag; }
58     virtual void    SetMESettings(Int_t mevt, Int_t mtrk, Int_t nminMix){fMEMaxPoolEvent = mevt, fMEMinTracks = mtrk,  fMEMinEventToMix = nminMix;}
59     virtual void    SetSystem(Bool_t system) { fSetSystemValue = system; }
60     virtual void    SetTrigger1PTValue(Double_t pTmin1, Double_t pTmax1){fTrigger1pTLowThr = pTmin1, fTrigger1pTHighThr = pTmax1;}
61     virtual void    SetTrigger2PTValue(Double_t pTmin2, Double_t pTmax2){fTrigger2pTLowThr = pTmin2, fTrigger2pTHighThr = pTmax2;}
62     virtual void    SetCentralityRange(Double_t minCent, Double_t maxCent){fMinCentrality = minCent, fMaxCentrality = maxCent;}
63     virtual void    SetDataType(Bool_t DataOrPart){fRecoOrMontecarlo = DataOrPart;}
64     virtual void    SetFilterBit(Bool_t flag=kTRUE){fSetFilterBit=flag;}
65     virtual void    SetFilterType(Int_t bittype){fbit= bittype;}
66     virtual void    SetVarCentBin(Bool_t varCbin=kTRUE){fuseVarCentBins= varCbin;}
67     virtual void    SetVarPtBin(Bool_t varPtbin=kTRUE){fuseVarPtBins= varPtbin;}
68  
69 private:
70     
71     AliAnalysisTaskDiJetCorrelations(const AliAnalysisTaskDiJetCorrelations &source);
72     AliAnalysisTaskDiJetCorrelations& operator=(const AliAnalysisTaskDiJetCorrelations& source);
73   
74     void DefineHistoNames();
75     Bool_t ConversionResonanceCut(Double_t refmaxpT, Double_t phiMaxpT, Double_t etaMaxpT, Double_t Charge, AliAODTrack* AodTracks, TH2F*fControlConvResT, TH1F* fHistTCorrTrack);
76     Bool_t TwoTrackEfficiencyCut(Double_t refmaxpT, Double_t phiMaxpT, Double_t etaMaxpT, Double_t Charge, AliAODTrack* AodTracks, Float_t bSigntmp);
77  
78     //______________________________|  Track Efficiency
79     Double_t GetTrackbyTrackEffValue(AliAODTrack* track, Double_t CentrOrMult){
80       
81         Float_t effvalue = -999.;
82         Int_t binX = f3DEffCor->GetXaxis()->FindBin(CentrOrMult);
83         Int_t binY = f3DEffCor->GetYaxis()->FindBin(track->Eta());
84         Int_t binZ = f3DEffCor->GetZaxis()->FindBin(track->Pt());
85         effvalue = f3DEffCor->GetBinContent(binX, binY, binZ);
86         if(effvalue==0) effvalue=1.;
87         return effvalue;
88     }
89     
90     //______________________________|  DPhi Star
91     Float_t GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1, Float_t phi2, Float_t pt2, Float_t charge2, Float_t radius, Float_t bSign){
92         
93         //calculate dphistar for two track efficiency cut!
94         Float_t dphistar = phi1 - phi2 - charge1 * bSign * TMath::ASin(0.075 * radius / pt1) + charge2 * bSign * TMath::ASin(0.075 * radius / pt2);
95         static const Double_t kPi = TMath::Pi();
96         if (dphistar > kPi)
97             dphistar = kPi*2 - dphistar;
98         if (dphistar < -kPi)
99             dphistar = -kPi*2 - dphistar;
100         //if (dphistar > kPi)
101         //  dphistar = kPi * 2 - dphistar;
102         return dphistar;
103     }
104     
105     //______________________________|  Mass Squared
106     Float_t GetInvMassSquared(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2){        
107         //calculate inv ass squared
108         Float_t tantheta1 = 1e10;
109         if (eta1 < -1e-10 || eta1 > 1e-10){
110             Float_t expTmp = TMath::Exp(-eta1);
111             tantheta1 = 2.0 * expTmp/(1.0 -expTmp*expTmp);
112         }
113         
114         Float_t tantheta2 = 1e10;
115         if (eta2 < -1e-10 || eta2 > 1e-10){
116             Float_t expTmp = TMath::Exp(-eta2);
117             tantheta2 = 2* expTmp/(1.0 - expTmp*expTmp);
118         }
119         
120         Float_t e1squ = m0_1 * m0_1 + pt1 * pt1 * (1.0 + 1.0 / tantheta1 / tantheta1);
121         Float_t e2squ = m0_2 * m0_2 + pt2 * pt2 * (1.0 + 1.0 / tantheta2 / tantheta2);
122         Float_t mass2 = m0_1 * m0_1 + m0_2 * m0_2 + 2 * (TMath::Sqrt(e1squ * e2squ) - ( pt1 * pt2 * ( TMath::Cos(phi1 - phi2) + 1.0 / tantheta1 / tantheta2) ) );
123         return mass2;
124     }
125     
126     //______________________________|  Mass Cheap
127     Float_t GetInvMassSquaredCheap(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2){
128         
129         // calculate inv mass squared approximately
130         Float_t tantheta1 = 1e10;        
131         if (eta1 < -1e-10 || eta1 > 1e-10){
132             Float_t expTmp = 1.0-eta1+eta1*eta1/2-eta1*eta1*eta1/6+eta1*eta1*eta1*eta1/24;
133             tantheta1 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
134         }
135         
136         Float_t tantheta2 = 1e10;
137         if (eta2 < -1e-10 || eta2 > 1e-10){
138             Float_t expTmp = 1.0-eta2+eta2*eta2/2-eta2*eta2*eta2/6+eta2*eta2*eta2*eta2/24;
139             tantheta2 = 2.0 * expTmp / ( 1.0 - expTmp*expTmp);
140         }
141         
142         Float_t e1squ = m0_1 * m0_1 + pt1 * pt1 * (1.0 + 1.0 / tantheta1 / tantheta1);
143         Float_t e2squ = m0_2 * m0_2 + pt2 * pt2 * (1.0 + 1.0 / tantheta2 / tantheta2);
144         
145         // fold onto 0...pi
146         Float_t deltaPhi = TMath::Abs(phi1 - phi2);
147         while (deltaPhi > TMath::TwoPi())
148             deltaPhi -= TMath::TwoPi();
149         if (deltaPhi > TMath::Pi())
150             deltaPhi = TMath::TwoPi() - deltaPhi;
151         
152         Float_t cosDeltaPhi = 0;
153         if (deltaPhi < TMath::Pi()/3)
154             cosDeltaPhi = 1.0 - deltaPhi*deltaPhi/2 + deltaPhi*deltaPhi*deltaPhi*deltaPhi/24;
155         else if (deltaPhi < 2*TMath::Pi()/3)
156             cosDeltaPhi = -(deltaPhi - TMath::Pi()/2) + 1.0/6 * TMath::Power((deltaPhi - TMath::Pi()/2), 3);
157         else
158             cosDeltaPhi = -1.0 + 1.0/2.0*(deltaPhi - TMath::Pi())*(deltaPhi - TMath::Pi()) - 1.0/24.0 * TMath::Power(deltaPhi - TMath::Pi(), 4);
159                 
160         Float_t mass2 = m0_1 * m0_1 + m0_2 * m0_2 + 2 * ( TMath::Sqrt(e1squ * e2squ) - ( pt1 * pt2 * ( cosDeltaPhi + 1.0 / tantheta1 / tantheta2 ) ) );
161         //   Printf(Form("%f %f %f %f %f %f %f %f %f", pt1, eta1, phi1, pt2, eta2, phi2, m0_1, m0_2, mass2));
162         return mass2;
163     }
164     
165     //______________________________| Phi Settings
166     Double_t AssignCorrectPhiRange(Double_t phi){
167         if (phi > 1.5 * TMath::Pi())phi -= TMath::TwoPi();
168         if (phi < -0.5 * TMath::Pi())phi += TMath::TwoPi();
169         Double_t phiClone = 0.;
170         phiClone = phi;
171         return phiClone;
172     }
173     
174     //______________________________| Mixed Event Pool
175     Bool_t DefineMixedEventPool(){
176         Int_t  NofCentBins  = 4;
177         Double_t MBins[]={0, 10, 20, 40, 100.1};
178         Double_t * CentrORMultBins = MBins;
179         
180         Int_t NofZVrtxBins  = 3;
181         Double_t ZBins[]={-9.999, -2.5, 2.5, 10.1};
182         Double_t *ZVrtxBins = ZBins;
183         
184         fPoolMgr = new AliEventPoolManager(fMEMaxPoolEvent, fMEMinTracks, NofCentBins, CentrORMultBins, NofZVrtxBins, ZVrtxBins);
185         if(!fPoolMgr) return kFALSE;
186         return kTRUE;
187     }
188     
189     //______________________________| Mixed Event Pool Process
190     Bool_t ProcessMixedEventPool(){
191         if(!fMixedEvent) return kFALSE;
192         if(!fPool->IsReady()) return kFALSE;
193         if(fPool->GetCurrentNEvents()<fMEMinEventToMix) return kFALSE;
194         return kTRUE;
195     }
196     
197     
198     //______________________________| All Used Ojects
199     Bool_t    fSetSystemValue;
200     Bool_t    fRecoOrMontecarlo;
201     Bool_t    fReadMC;
202     Bool_t    fSetFilterBit; //
203     Int_t     fbit ;   //
204     TClonesArray *farrayMC;//!
205     
206     Double_t  fCentrOrMult; // Multiplicity of Event for D eff
207     Double_t  fMinCentrality; // Minimun Centrality Value
208     Double_t  fMaxCentrality; // Maximum Centrality Value
209     
210     Double_t  fTrigger1pTLowThr;
211     Double_t  fTrigger1pTHighThr;
212     Double_t  fTrigger2pTLowThr;
213     Double_t  fTrigger2pTHighThr;
214     
215     Bool_t    fCutResonances;
216     Bool_t    fCutConversions;
217     Bool_t    ftwoTrackEfficiencyCut;
218     Bool_t    fuseVarCentBins;
219     Bool_t    fuseVarPtBins;
220   
221     TH1F     *fHistNEvents; //!
222     TH1F     *fHistT1CorrTrack; //!
223     TH1F     *fHistT2CorrTrack; //!
224     TList    *fOutputQA; //! Output list
225     TList    *fOutputCorr; //! Output list
226     TH3F     *f3DEffCor; //!
227     
228     AliEventPool *fPool; //! Pool for event mixing
229     AliEventPoolManager  *fPoolMgr;         //! event pool manager
230     Bool_t   fMixedEvent;               // enable event mixing (default: ON)
231     Int_t  fMEMaxPoolEvent;
232     Int_t  fMEMinTracks;
233     Int_t  fMEMinEventToMix;
234
235     TH1F *fHistQA[9]; //!
236     TH1F *fHistTrigDPhi; //!
237     TH2F *fControlConvResT1; //!
238     TH2F *fControlConvResT2; //!
239     TH2F *fControlConvResMT1;//!
240     TH2F *fControlConvResMT2;//!
241
242
243     ClassDef(AliAnalysisTaskDiJetCorrelations, 1); // example of analysis
244 };
245
246 #endif