1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 //_________________________________________________________________________
17 // class to extract omega(782)->pi0+gamma->3gamma
18 // Mar. 22, 2011: Additional method, espeically for EMCAL. A high E cluster is assumpted as pi0 (two photons are overlapped) without unfolding
20 //-- Author: Renzhuo Wan (IOPP-Wuhan, China)
21 //_________________________________________________________________________
28 // --- ROOT system ---
30 #include "TLorentzVector.h"
31 #include "TParticle.h"
34 //---- AliRoot system ----
35 #include "AliAnaOmegaToPi0Gamma.h"
36 #include "AliCaloTrackReader.h"
37 #include "AliCaloPID.h"
39 #include "AliVEvent.h"
40 #include "AliAODEvent.h"
41 #include "AliAODMCParticle.h"
42 ClassImp(AliAnaOmegaToPi0Gamma)
44 //______________________________________________________________________________
45 AliAnaOmegaToPi0Gamma::AliAnaOmegaToPi0Gamma() : AliAnaPartCorrBaseClass(),
46 fInputAODPi0(0), fInputAODGammaName(""),
47 fEventsList(0x0),fNVtxZBin(0), fNCentBin(0), fNRpBin(0), fNBadChDistBin(0), fNpid(0),
48 fVtxZCut(0), fCent(0), fRp(0),
49 fPi0Mass(0),fPi0MassWindow(0),fPi0OverOmegaPtCut(0),
50 fGammaOverOmegaPtCut(0), fEOverlapCluster(0),
52 fRealOmega0(0), fMixAOmega0(0),
53 fMixBOmega0(0), fMixCOmega0(0),
54 fRealOmega1(0), fMixAOmega1(0),
55 fMixBOmega1(0), fMixCOmega1(0),
56 fRealOmega2(0), fMixAOmega2(0),
57 fMixBOmega2(0), fMixCOmega2(0),
65 //______________________________________________________________________________
66 AliAnaOmegaToPi0Gamma::~AliAnaOmegaToPi0Gamma() {
71 // fInputAODPi0->Clear();
72 // delete fInputAODPi0;
76 for(Int_t i=0;i<fNVtxZBin;i++){
77 for(Int_t j=0;j<fNCentBin;j++){
78 for(Int_t k=0;k<fNRpBin;k++){
79 fEventsList[i*fNCentBin*fNRpBin+j*fNRpBin+k]->Clear();
80 delete fEventsList[i*fNCentBin*fNRpBin+j*fNRpBin+k];
85 delete [] fEventsList;
94 //______________________________________________________________________________
95 void AliAnaOmegaToPi0Gamma::InitParameters()
97 //Init parameters when first called the analysis
98 //Set default parameters
99 fInputAODGammaName="PhotonsDetector";
107 fPi0MassWindow=0.015;
108 fPi0OverOmegaPtCut=0.8;
109 fGammaOverOmegaPtCut=0.2;
115 //______________________________________________________________________________
116 TList * AliAnaOmegaToPi0Gamma::GetCreateOutputObjects()
119 fVtxZCut = new Double_t [fNVtxZBin];
120 for(Int_t i=0;i<fNVtxZBin;i++) fVtxZCut[i]=10*(i+1);
122 fCent=new Double_t[fNCentBin];
123 for(int i = 0;i<fNCentBin;i++)fCent[i]=0;
125 fRp=new Double_t[fNRpBin];
126 for(int i = 0;i<fNRpBin;i++)fRp[i]=0;
128 Int_t nptbins = GetHistoPtBins();
129 Float_t ptmax = GetHistoPtMax();
130 Float_t ptmin = GetHistoPtMin();
132 Int_t nmassbins = GetHistoMassBins();
133 Float_t massmin = GetHistoMassMin();
134 Float_t massmax = GetHistoMassMax();
136 fhEtalon = new TH2F("hEtalon","Histo with binning parameters", nptbins,ptmin,ptmax,nmassbins,massmin,massmax) ;
137 fhEtalon->SetXTitle("P_{T} (GeV)") ;
138 fhEtalon->SetYTitle("m_{inv} (GeV)") ;
140 // store them in fOutputContainer
141 fEventsList = new TList*[fNVtxZBin*fNCentBin*fNRpBin];
142 for(Int_t i=0;i<fNVtxZBin;i++){
143 for(Int_t j=0;j<fNCentBin;j++){
144 for(Int_t k=0;k<fNRpBin;k++){
145 fEventsList[i*fNCentBin*fNRpBin+j*fNRpBin+k] = new TList();
146 fEventsList[i*fNCentBin*fNRpBin+j*fNRpBin+k]->SetOwner(kFALSE);
151 TList * outputContainer = new TList() ;
152 outputContainer->SetName(GetName());
153 const Int_t buffersize = 255;
154 char key[buffersize] ;
155 char title[buffersize] ;
156 const char * detector= fInputAODGammaName.Data();
157 Int_t ndim=fNVtxZBin*fNCentBin*fNRpBin*fNBadChDistBin*fNpid;
159 fRealOmega0 =new TH2F*[ndim];
160 fMixAOmega0 =new TH2F*[ndim];
161 fMixBOmega0 =new TH2F*[ndim];
162 fMixCOmega0 =new TH2F*[ndim];
164 fRealOmega1 =new TH2F*[ndim];
165 fMixAOmega1 =new TH2F*[ndim];
166 fMixBOmega1 =new TH2F*[ndim];
167 fMixCOmega1 =new TH2F*[ndim];
169 fRealOmega2 =new TH2F*[ndim];
170 fMixAOmega2 =new TH2F*[ndim];
171 fMixBOmega2 =new TH2F*[ndim];
172 fMixCOmega2 =new TH2F*[ndim];
174 fhFakeOmega = new TH2F*[fNCentBin];
176 for(Int_t i=0;i<fNVtxZBin;i++){
177 for(Int_t j=0;j<fNCentBin;j++){
178 for(Int_t k=0;k<fNRpBin;k++){ //at event level
179 Int_t idim=i*fNCentBin*fNRpBin+j*fNRpBin+k;
180 for(Int_t ipid=0;ipid<fNpid;ipid++){
181 for(Int_t idist=0;idist<fNBadChDistBin;idist++){ //at particle level
183 Int_t index=idim*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;
185 snprintf(key,buffersize,"RealToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
186 snprintf(title,buffersize, "%s Real Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
187 fRealOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;
188 fRealOmega0[index]->SetName(key) ;
189 fRealOmega0[index]->SetTitle(title);
190 outputContainer->Add(fRealOmega0[index]);
192 snprintf(key,buffersize,"MixAToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
193 snprintf(title,buffersize, "%s MixA Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
194 fMixAOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;
195 fMixAOmega0[index]->SetName(key) ;
196 fMixAOmega0[index]->SetTitle(title);
197 outputContainer->Add(fMixAOmega0[index]);
199 snprintf(key,buffersize,"MixBToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
200 snprintf(title,buffersize, "%s MixB Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
201 fMixBOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;
202 fMixBOmega0[index]->SetName(key) ;
203 fMixBOmega0[index]->SetTitle(title);
204 outputContainer->Add(fMixBOmega0[index]);
206 snprintf(key,buffersize,"MixCToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
207 snprintf(title,buffersize, "%s MixC Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
208 fMixCOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;
209 fMixCOmega0[index]->SetName(key) ;
210 fMixCOmega0[index]->SetTitle(title);
211 outputContainer->Add(fMixCOmega0[index]);
213 snprintf(key,buffersize,"RealToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
214 snprintf(title,buffersize, "%s Real Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
215 fRealOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;
216 fRealOmega1[index]->SetName(key) ;
217 fRealOmega1[index]->SetTitle(title);
218 outputContainer->Add(fRealOmega1[index]);
220 snprintf(key,buffersize,"MixAToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
221 snprintf(title,buffersize, "%s MixA Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
222 fMixAOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;
223 fMixAOmega1[index]->SetName(key) ;
224 fMixAOmega1[index]->SetTitle(title);
225 outputContainer->Add(fMixAOmega1[index]);
227 snprintf(key,buffersize,"MixBToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
228 snprintf(title,buffersize, "%s MixB Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
229 fMixBOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;
230 fMixBOmega1[index]->SetName(key) ;
231 fMixBOmega1[index]->SetTitle(title);
232 outputContainer->Add(fMixBOmega1[index]);
234 snprintf(key,buffersize,"MixCToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
235 snprintf(title,buffersize, "%s MixC Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
236 fMixCOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;
237 fMixCOmega1[index]->SetName(key) ;
238 fMixCOmega1[index]->SetTitle(title);
239 outputContainer->Add(fMixCOmega1[index]);
241 snprintf(key,buffersize,"RealToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
242 snprintf(title,buffersize, "%s Real Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
243 fRealOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;
244 fRealOmega2[index]->SetName(key) ;
245 fRealOmega2[index]->SetTitle(title);
246 outputContainer->Add(fRealOmega2[index]);
248 snprintf(key,buffersize,"MixAToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
249 snprintf(title,buffersize, "%s MixA Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
250 fMixAOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;
251 fMixAOmega2[index]->SetName(key) ;
252 fMixAOmega2[index]->SetTitle(title);
253 outputContainer->Add(fMixAOmega2[index]);
255 snprintf(key,buffersize,"MixBToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
256 snprintf(title,buffersize, "%s MixB Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
257 fMixBOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;
258 fMixBOmega2[index]->SetName(key) ;
259 fMixBOmega2[index]->SetTitle(title);
260 outputContainer->Add(fMixBOmega2[index]);
262 snprintf(key,buffersize,"MixCToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);
263 snprintf(title,buffersize, "%s MixC Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%d",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,idist);
264 fMixCOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;
265 fMixCOmega2[index]->SetName(key) ;
266 fMixCOmega2[index]->SetTitle(title);
267 outputContainer->Add(fMixCOmega2[index]);
274 for(Int_t i=0;i<fNCentBin;i++){
275 snprintf(key,buffersize, "fhFakeOmega%d",i);
276 snprintf(title,buffersize,"FakePi0(high pt cluster)+Gamma with Centrality%d",i);
277 fhFakeOmega[i] = (TH2F*)fhEtalon->Clone(key);
278 fhFakeOmega[i]->SetTitle(title);
279 outputContainer->Add(fhFakeOmega[i]);
282 snprintf(key,buffersize, "%sOmegaPri",detector);
283 snprintf(title,buffersize,"primary #omega in %s",detector);
284 fhOmegaPriPt=new TH1F(key, title,nptbins,ptmin,ptmax);
285 fhOmegaPriPt->GetXaxis()->SetTitle("P_{T}");
286 fhOmegaPriPt->GetYaxis()->SetTitle("dN/P_{T}");
287 outputContainer->Add(fhOmegaPriPt);
291 return outputContainer;
294 //______________________________________________________________________________
295 void AliAnaOmegaToPi0Gamma::Print(const Option_t * /*opt*/) const
297 //Print some relevant parameters set in the analysis
298 printf("**** Print %s %s ****\n", GetName(), GetTitle() ) ;
299 AliAnaPartCorrBaseClass::Print(" ");
300 printf("Omega->pi0+gamma->3gamma\n");
301 printf("Cuts at event level: \n");
302 printf("Bins of vertex Z: %d \n", fNVtxZBin);
303 printf("Bins of centrality: %d \n",fNCentBin);
304 printf("Bins of Reaction plane: %d\n",fNRpBin);
305 printf("Cuts at AOD particle level:\n");
306 printf("Number of PID: %d \n", fNpid);
307 printf("Number of DistToBadChannel cuts: %d\n", fNBadChDistBin);
310 //______________________________________________________________________________
311 void AliAnaOmegaToPi0Gamma::MakeAnalysisFillHistograms()
313 //fill the MC AOD if needed first
315 //need to be further implemented
316 AliStack * stack = 0x0;
317 // TParticle * primary = 0x0;
318 TClonesArray * mcparticles0 = 0x0;
319 //TClonesArray * mcparticles1 = 0x0;
320 AliAODMCParticle * aodprimary = 0x0;
326 if(GetReader()->ReadStack()){
327 stack = GetMCStack() ;
329 printf("AliAnaAcceptance::MakeAnalysisFillHistograms() - There is no stack!\n");
332 for(Int_t i=0 ; i<stack->GetNtrack(); i++){
333 TParticle * prim = stack->Particle(i) ;
334 pdg = prim->GetPdgCode() ;
337 if(TMath::Abs(eta)<0.5) {
338 if(pdg==223) fhOmegaPriPt->Fill(pt);
343 else if(GetReader()->ReadAODMCParticles()){
344 //Get the list of MC particles
345 mcparticles0 = GetReader()->GetAODMCParticles(0);
347 if(GetDebug() > 0) printf("AliAnaAcceptance::MakeAnalysisFillHistograms() - Standard MCParticles not available!\n");
349 // if(GetReader()->GetSecondInputAODTree()){
350 // mcparticles1 = GetReader()->GetAODMCParticles(1);
351 // if(!mcparticles1 && GetDebug() > 0) {
352 // printf("AliAnaAcceptance::MakeAnalysisFillHistograms() - Second input MCParticles not available!\n");
356 for(Int_t i=0;i<mcparticles0->GetEntries();i++){
357 aodprimary =(AliAODMCParticle*)mcparticles0->At(i);
358 pdg = aodprimary->GetPdgCode() ;
359 eta=aodprimary->Eta();
361 if(TMath::Abs(eta)<0.5) {
362 if(pdg==223) fhOmegaPriPt->Fill(pt);
366 }//mcparticles0 exists
371 //process event from AOD brach
372 //extract pi0, eta and omega analysis
373 Int_t iRun=(GetReader()->GetInputEvent())->GetRunNumber() ;
374 if(IsBadRun(iRun)) return ;
377 Double_t vert[]={0,0,0} ;
379 Int_t curEventBin =0;
381 Int_t ivtxzbin=(Int_t)TMath::Abs(vert[2])/10;
382 if(ivtxzbin>=fNVtxZBin)return;
385 Int_t currentCentrality = GetEventCentrality();
386 if(currentCentrality == -1) return;
387 Int_t optCent = GetReader()->GetCentralityOpt();
388 Int_t icentbin=currentCentrality/(optCent/fNCentBin) ; //GetEventCentrality();
390 printf("-------------- %d %d %d ",currentCentrality, optCent, icentbin);
394 if(ivtxzbin==-1) return;
395 curEventBin = ivtxzbin*fNCentBin*fNRpBin + icentbin*fNRpBin + irpbin;
396 TClonesArray *aodGamma = (TClonesArray*) GetAODBranch(fInputAODGammaName); //photon array
397 // TClonesArray *aodGamma = (TClonesArray *) GetReader()->GetOutputEvent()->FindListObject(fInputAODGammaName); //photon array
399 if(aodGamma) nphotons= aodGamma->GetEntries();
402 fInputAODPi0 = (TClonesArray*)GetInputAODBranch(); //pi0 array
404 if(fInputAODPi0) npi0s= fInputAODPi0 ->GetEntries();
407 if(nphotons<3 || npi0s<1)return; //for pi0, eta and omega->pi0+gamma->3gamma reconstruction
409 //reconstruction of omega(782)->pi0+gamma->3gamma
410 //loop for pi0 and photon
411 if(GetDebug() > 0) printf("omega->pi0+gamma->3gamma invariant mass analysis ! This event have %d photons and %d pi0 \n", nphotons, npi0s);
412 for(Int_t i=0;i<npi0s;i++){
413 AliAODPWG4Particle * pi0 = (AliAODPWG4Particle*) (fInputAODPi0->At(i)) ; //pi0
414 TLorentzVector vpi0(pi0->Px(),pi0->Py(),pi0->Pz(),pi0->E());
415 Int_t lab1=pi0->GetCaloLabel(0); // photon1 from pi0 decay
416 Int_t lab2=pi0->GetCaloLabel(1); // photon2 from pi0 decay
417 //for omega->pi0+gamma, it needs at least three photons per event
418 //Get the two decay photons from pi0
419 AliAODPWG4Particle * photon1 =0;
420 AliAODPWG4Particle * photon2 =0;
421 for(Int_t d1=0;d1<nphotons;d1++){
422 for(Int_t d2=0;d2<nphotons;d2++){
423 AliAODPWG4Particle * dp1 = (AliAODPWG4Particle*) (aodGamma->At(d1));
424 AliAODPWG4Particle * dp2 = (AliAODPWG4Particle*) (aodGamma->At(d2));
425 Int_t dlab1=dp1->GetCaloLabel(0);
426 Int_t dlab2=dp2->GetCaloLabel(0);
427 if(dlab1==lab1 && dlab2==lab2){
434 //caculate the asy and dist of the two photon from pi0 decay
435 TLorentzVector dph1(photon1->Px(),photon1->Py(),photon1->Pz(),photon1->E());
436 TLorentzVector dph2(photon2->Px(),photon2->Py(),photon2->Pz(),photon2->E());
438 Double_t pi0asy= TMath::Abs(dph1.E()-dph2.E())/(dph1.E()+dph2.E());
439 // Double_t phi1=dph1.Phi();
440 // Double_t phi2=dph2.Phi();
441 // Double_t eta1=dph1.Eta();
442 // Double_t eta2=dph2.Eta();
443 // Double_t pi0dist=TMath::Sqrt((phi1-phi2)*(phi1-phi2)+(eta1-eta2)*(eta1-eta2));
445 if(pi0->GetIdentifiedParticleType()==111 && nphotons>2 && npi0s
446 && TMath::Abs(vpi0.M()-fPi0Mass)<fPi0MassWindow) { //pi0 candidates
448 //avoid the double counting
449 Int_t * dc1= new Int_t[nphotons];
450 Int_t * dc2= new Int_t[nphotons];
453 for(Int_t k=0;k<i;k++){
454 AliAODPWG4Particle * p3=(AliAODPWG4Particle*)(fInputAODPi0->At(k));
455 Int_t lab4=p3->GetCaloLabel(0);
456 Int_t lab5=p3->GetCaloLabel(1);
457 if(lab1==lab4){ dc1[index1]=lab5; index1++; }
458 if(lab2==lab5){ dc2[index2]=lab4; index2++; }
462 //loop the pi0 with third gamma
463 for(Int_t j=0;j<nphotons;j++){
464 AliAODPWG4Particle *photon3 = (AliAODPWG4Particle*) (aodGamma->At(j));
465 TLorentzVector dph3(photon3->Px(),photon3->Py(),photon3->Pz(),photon3->E());
466 Int_t lab3=photon3->GetCaloLabel(0);
467 Double_t pi0gammapt=(vpi0+dph3).Pt();
468 Double_t pi0gammamass=(vpi0+dph3).M();
469 Double_t pi0OverOmegaPtRatio =vpi0.Pt()/pi0gammapt;
470 Double_t gammaOverOmegaPtRatio= dph3.Pt()/pi0gammapt;
473 if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut ||
474 gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;
476 for(Int_t l=0;l<index1;l++) if(lab3==dc1[l]) lab3=-1;
477 for(Int_t l=0;l<index2;l++) if(lab3==dc2[l]) lab3=-1;
479 if(lab3>0 && lab3!=lab1 && lab3!=lab2){
480 for(Int_t ipid=0;ipid<fNpid;ipid++){
481 for(Int_t idist=0;idist<fNBadChDistBin;idist++){
482 Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;
483 if(photon1->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
484 photon2->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
485 photon3->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
486 photon1->DistToBad()>=idist &&
487 photon2->DistToBad()>=idist &&
488 photon3->DistToBad()>=idist ){
489 //fill the histograms
490 if(GetDebug() > 2) printf("Real: index %d pt %2.3f mass %2.3f \n", index, pi0gammapt, pi0gammamass);
491 fRealOmega0[index]->Fill(pi0gammapt,pi0gammamass);
492 if(pi0asy<0.7) fRealOmega1[index]->Fill(pi0gammapt,pi0gammamass);
493 if(pi0asy<0.8) fRealOmega2[index]->Fill(pi0gammapt,pi0gammamass);
501 if(GetDebug() > 0) printf("MixA: (r1_event1+r2_event1)+r3_event2 \n");
502 //-------------------------
503 //background analysis
505 // --A (r1_event1+r2_event1)+r3_event2
506 Int_t nMixed = fEventsList[curEventBin]->GetSize();
507 for(Int_t im=0;im<nMixed;im++){
508 TClonesArray* ev2= (TClonesArray*) (fEventsList[curEventBin]->At(im));
509 for(Int_t mix1=0;mix1<ev2->GetEntries();mix1++){
510 AliAODPWG4Particle *mix1ph = (AliAODPWG4Particle*) (ev2->At(mix1));
511 TLorentzVector vmixph(mix1ph->Px(),mix1ph->Py(),mix1ph->Pz(),mix1ph->E());
512 Double_t pi0gammapt=(vpi0+vmixph).Pt();
513 Double_t pi0gammamass=(vpi0+vmixph).M();
514 Double_t pi0OverOmegaPtRatio =vpi0.Pt()/pi0gammapt;
515 Double_t gammaOverOmegaPtRatio= vmixph.Pt()/pi0gammapt;
518 if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut ||
519 gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;
521 for(Int_t ipid=0;ipid<fNpid;ipid++){
522 for(Int_t idist=0;idist<fNBadChDistBin;idist++){
523 Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;
524 if(photon1->IsPIDOK(ipid,AliCaloPID::kPhoton)&&
525 photon2->IsPIDOK(ipid,AliCaloPID::kPhoton)&&
526 mix1ph->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
527 photon1->DistToBad()>=idist &&
528 photon2->DistToBad()>=idist &&
529 mix1ph->DistToBad()>=idist ){
530 if(GetDebug() > 2) printf("MixA: index %d pt %2.3f mass %2.3f \n",index, pi0gammapt, pi0gammamass);
531 //fill the histograms
532 fMixAOmega0[index]->Fill(pi0gammapt,pi0gammamass);
533 if(pi0asy<0.7)fMixAOmega1[index]->Fill(pi0gammapt,pi0gammamass);
534 if(pi0asy<0.8)fMixAOmega2[index]->Fill(pi0gammapt,pi0gammamass);
535 //printf("mix A %d %2.2f \n", index, pi0gammamass);
546 // --B (r1_event1+r2_event2)+r3_event2
548 if(GetDebug() >0)printf("MixB: (r1_event1+r2_event2)+r3_event2 \n");
549 for(Int_t i=0;i<nphotons;i++){
550 AliAODPWG4Particle *ph1 = (AliAODPWG4Particle*) (aodGamma->At(i));
551 TLorentzVector vph1(ph1->Px(),ph1->Py(),ph1->Pz(),ph1->E());
552 //interrupt here...................
553 //especially for EMCAL
554 //we suppose the high pt clusters are overlapped pi0
556 for(Int_t j=i+1;j<nphotons;j++){
557 AliAODPWG4Particle *ph2 = (AliAODPWG4Particle*) (aodGamma->At(j));
558 TLorentzVector fakePi0, fakeOmega, vph;
560 if(ph1->E() > fEOverlapCluster && ph1->E() > ph2->E()) {
561 fakePi0.SetPxPyPzE(ph1->Px(),ph1->Py(),ph1->Pz(),TMath::Sqrt(ph1->Px()*ph1->Px()+ph1->Py()*ph1->Py()+ph1->Pz()*ph1->Pz()+0.135*0.135));
562 vph.SetPxPyPzE(ph2->Px(),ph2->Py(),ph2->Pz(),ph2->E());
564 else if(ph2->E() > fEOverlapCluster && ph2->E() > ph1->E()) {
565 fakePi0.SetPxPyPzE(ph2->Px(),ph2->Py(),ph2->Pz(),TMath::Sqrt(ph2->Px()*ph2->Px()+ph2->Py()*ph2->Py()+ph2->Pz()*ph2->Pz()+0.135*0.135));
566 vph.SetPxPyPzE(ph1->Px(), ph1->Py(),ph1->Pz(),ph1->E());
570 fakeOmega=fakePi0+vph;
571 for(Int_t ii=0;ii<fNCentBin;ii++){
572 fhFakeOmega[icentbin]->Fill(fakeOmega.Pt(), fakeOmega.M());
576 //continue ................
577 Int_t nMixed = fEventsList[curEventBin]->GetSize();
578 for(Int_t ie=0;ie<nMixed;ie++){
579 TClonesArray* ev2= (TClonesArray*) (fEventsList[curEventBin]->At(ie));
580 for(Int_t mix1=0;mix1<ev2->GetEntries();mix1++){
581 AliAODPWG4Particle *ph2 = (AliAODPWG4Particle*) (ev2->At(mix1));
582 TLorentzVector vph2(ph2->Px(),ph2->Py(),ph2->Pz(),ph2->E());
583 Double_t pi0asy = TMath::Abs(vph1.E()-vph2.E())/(vph1.E()+vph2.E());
584 Double_t pi0mass=(vph1+vph2).M();
586 if(TMath::Abs(pi0mass-fPi0Mass)<fPi0MassWindow){//for pi0 selection
587 for(Int_t mix2=(mix1+1);mix2<ev2->GetEntries();mix2++){
588 AliAODPWG4Particle *ph3 = (AliAODPWG4Particle*) (ev2->At(mix2));
589 TLorentzVector vph3(ph3->Px(),ph3->Py(),ph3->Pz(),ph3->E());
591 Double_t pi0gammapt=(vph1+vph2+vph3).Pt();
592 Double_t pi0gammamass=(vph1+vph2+vph3).M();
593 Double_t pi0OverOmegaPtRatio =(vph1+vph2).Pt()/pi0gammapt;
594 Double_t gammaOverOmegaPtRatio= vph3.Pt()/pi0gammapt;
596 if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut ||
597 gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;
599 for(Int_t ipid=0;ipid<fNpid;ipid++){
600 for(Int_t idist=0;idist<fNBadChDistBin;idist++){
601 Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;
602 if(ph1->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
603 ph2->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
604 ph3->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
605 ph1->DistToBad()>=idist &&
606 ph2->DistToBad()>=idist &&
607 ph3->DistToBad()>=idist ){
608 if(GetDebug() > 2) printf("MixB: index %d pt %2.3f mass %2.3f \n", index, pi0gammapt, pi0gammamass);
610 fMixBOmega0[index]->Fill(pi0gammapt,pi0gammamass);
611 if(pi0asy<0.7) fMixBOmega1[index]->Fill(pi0gammapt,pi0gammamass);
612 if(pi0asy<0.8) fMixBOmega2[index]->Fill(pi0gammapt,pi0gammamass);
613 //printf("mix B %d %2.2f \n", index, pi0gammamass);
620 // --C (r1_event1+r2_event2)+r3_event3
622 if(GetDebug() >0)printf("MixC: (r1_event1+r2_event2)+r3_event3\n");
623 for(Int_t je=(ie+1);je<nMixed;je++){
624 TClonesArray* ev3= (TClonesArray*) (fEventsList[curEventBin]->At(je));
625 for(Int_t mix3=0;mix3<ev3->GetEntries();mix3++){
626 AliAODPWG4Particle *ph3 = (AliAODPWG4Particle*) (ev3->At(mix3));
627 TLorentzVector vph3(ph3->Px(),ph3->Py(),ph3->Pz(),ph3->E());
629 Double_t pi0gammapt=(vph1+vph2+vph3).Pt();
630 Double_t pi0gammamass=(vph1+vph2+vph3).M();
631 Double_t pi0OverOmegaPtRatio =(vph1+vph2).Pt()/pi0gammapt;
632 Double_t gammaOverOmegaPtRatio= vph3.Pt()/pi0gammapt;
634 if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut ||
635 gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;
637 for(Int_t ipid=0;ipid<fNpid;ipid++){
638 for(Int_t idist=0;idist<fNBadChDistBin;idist++){
639 Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;
640 if(ph1->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
641 ph2->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
642 ph3->IsPIDOK(ipid,AliCaloPID::kPhoton) &&
643 ph1->DistToBad()>=idist &&
644 ph2->DistToBad()>=idist &&
645 ph3->DistToBad()>=idist ){
646 if(GetDebug() > 2) printf("MixC: index %d pt %2.3f mass %2.3f \n", index, pi0gammapt, pi0gammamass);
648 fMixCOmega0[index]->Fill(pi0gammapt,pi0gammamass);
649 if(pi0asy<0.7) fMixCOmega1[index]->Fill(pi0gammapt,pi0gammamass);
650 if(pi0asy<0.8) fMixCOmega2[index]->Fill(pi0gammapt,pi0gammamass);
651 //printf("mix C %d %2.2f \n", index, pi0gammamass);
664 TClonesArray *currentEvent = new TClonesArray(*aodGamma);
665 if(currentEvent->GetEntriesFast()>0){
666 fEventsList[curEventBin]->AddFirst(currentEvent) ;
668 if(fEventsList[curEventBin]->GetSize()>=GetNMaxEvMix()) {
669 TClonesArray * tmp = (TClonesArray*) (fEventsList[curEventBin]->Last()) ;
670 fEventsList[curEventBin]->RemoveLast() ;
675 delete currentEvent ;
681 //______________________________________________________________________________
682 void AliAnaOmegaToPi0Gamma::ReadHistograms(TList * outputList)
684 //read the histograms
685 //for the finalization of the terminate analysis
687 Int_t index = outputList->IndexOf(outputList->FindObject(GetAddedHistogramsStringToName()+"RealToPi0Gamma_Vz0C0Rp0Pid0Dist0"));
689 Int_t ndim=fNVtxZBin*fNCentBin*fNRpBin*fNBadChDistBin*fNpid;
691 if(!fRealOmega0) fRealOmega0 =new TH2F*[ndim];
692 if(!fMixAOmega0) fMixAOmega0 =new TH2F*[ndim];
693 if(!fMixBOmega0) fMixBOmega0 =new TH2F*[ndim];
694 if(!fMixCOmega0) fMixCOmega0 =new TH2F*[ndim];
696 if(!fRealOmega1) fRealOmega1 =new TH2F*[ndim];
697 if(!fMixAOmega1) fMixAOmega1 =new TH2F*[ndim];
698 if(!fMixBOmega1) fMixBOmega1 =new TH2F*[ndim];
699 if(!fMixCOmega1) fMixCOmega1 =new TH2F*[ndim];
701 if(!fRealOmega2) fRealOmega2 =new TH2F*[ndim];
702 if(!fMixAOmega2) fMixAOmega2 =new TH2F*[ndim];
703 if(!fMixBOmega2) fMixBOmega2 =new TH2F*[ndim];
704 if(!fMixCOmega2) fMixCOmega2 =new TH2F*[ndim];
706 for(Int_t i=0;i<fNVtxZBin;i++){
707 for(Int_t j=0;j<fNCentBin;j++){
708 for(Int_t k=0;k<fNRpBin;k++){ //at event level
709 Int_t idim=i*fNCentBin*fNRpBin+j*fNRpBin+k;
710 for(Int_t ipid=0;ipid<fNpid;ipid++){
711 for(Int_t idist=0;idist<fNBadChDistBin;idist++){ //at particle
712 Int_t ind=idim*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;
713 fRealOmega0[ind]= (TH2F*) outputList->At(index++);
714 fMixAOmega0[ind]= (TH2F*) outputList->At(index++);
715 fMixBOmega0[ind]= (TH2F*) outputList->At(index++);
716 fMixCOmega0[ind]= (TH2F*) outputList->At(index++);
718 fRealOmega1[ind]= (TH2F*) outputList->At(index++);
719 fMixAOmega1[ind]= (TH2F*) outputList->At(index++);
720 fMixBOmega1[ind]= (TH2F*) outputList->At(index++);
721 fMixCOmega1[ind]= (TH2F*) outputList->At(index++);
723 fRealOmega2[ind]= (TH2F*) outputList->At(index++);
724 fMixAOmega2[ind]= (TH2F*) outputList->At(index++);
725 fMixBOmega2[ind]= (TH2F*) outputList->At(index++);
726 fMixCOmega2[ind]= (TH2F*) outputList->At(index++);
736 fhOmegaPriPt = (TH1F*) outputList->At(index++);
741 //______________________________________________________________________________
742 void AliAnaOmegaToPi0Gamma::Terminate(TList * outputList)
744 // //Do some calculations and plots from the final histograms.
745 if(GetDebug() >= 0) printf("AliAnaOmegaToPi0Gamma::Terminate() \n");
746 ReadHistograms(outputList);
747 const Int_t buffersize = 255;
748 char cvs1[buffersize];
749 snprintf(cvs1, buffersize, "Neutral_%s_IVM",fInputAODGammaName.Data());
751 TCanvas * cvsIVM = new TCanvas(cvs1, cvs1, 400, 10, 600, 700) ;
752 cvsIVM->Divide(2, 2);
755 char dec[buffersize];
756 snprintf(dec,buffersize,"h2Real_%s",fInputAODGammaName.Data());
757 TH2F * h2Real= (TH2F*)fRealOmega0[0]->Clone(dec);
758 h2Real->GetXaxis()->SetRangeUser(4,6);
759 TH1F * hRealOmega = (TH1F*) h2Real->ProjectionY();
760 hRealOmega->SetTitle("RealPi0Gamma 4<pt<6");
761 hRealOmega->SetLineColor(2);
765 snprintf(dec,buffersize,"hMixA_%s",fInputAODGammaName.Data());
766 TH2F *h2MixA= (TH2F*)fMixAOmega0[0]->Clone(dec);
767 h2MixA->GetXaxis()->SetRangeUser(4,6);
768 TH1F * hMixAOmega = (TH1F*) h2MixA->ProjectionY();
769 hMixAOmega->SetTitle("MixA 4<pt<6");
770 hMixAOmega->SetLineColor(2);
774 snprintf(dec,buffersize,"hMixB_%s",fInputAODGammaName.Data());
775 TH2F * h2MixB= (TH2F*)fMixBOmega0[0]->Clone(dec);
776 h2MixB->GetXaxis()->SetRangeUser(4,6);
777 TH1F * hMixBOmega = (TH1F*) h2MixB->ProjectionY();
778 hMixBOmega->SetTitle("MixB 4<pt<6");
779 hMixBOmega->SetLineColor(2);
783 snprintf(dec,buffersize,"hMixC_%s",fInputAODGammaName.Data());
784 TH2F *h2MixC= (TH2F*)fMixCOmega0[0]->Clone(dec);
785 h2MixC->GetXaxis()->SetRangeUser(4,6);
786 TH1F * hMixCOmega = (TH1F*) h2MixC->ProjectionY();
787 hMixCOmega->SetTitle("MixC 4<pt<6");
788 hMixCOmega->SetLineColor(2);
791 char eps[buffersize];
792 snprintf(eps,buffersize,"CVS_%s_IVM.eps",fInputAODGammaName.Data());