]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrDep/AliAnaPi0.cxx
add histograms for pi0 analysis, both clusters in same module
[u/mrichter/AliRoot.git] / PWG4 / PartCorrDep / AliAnaPi0.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15 /* $Id: $ */
16
17 //_________________________________________________________________________
18 // Class to collect two-photon invariant mass distributions for
19 // extractin raw pi0 yield.
20 //
21 //-- Author: Dmitri Peressounko (RRC "KI") 
22 //-- Adapted to PartCorr frame by Lamia Benhabib (SUBATECH)
23 //-- and Gustavo Conesa (INFN-Frascati)
24 //_________________________________________________________________________
25
26
27 // --- ROOT system ---
28 #include "TH3.h"
29 //#include "Riostream.h"
30 #include "TCanvas.h"
31 #include "TPad.h"
32 #include "TROOT.h"
33 #include "TClonesArray.h"
34 #include "TObjString.h"
35
36 //---- AliRoot system ----
37 #include "AliAnaPi0.h"
38 #include "AliCaloTrackReader.h"
39 #include "AliCaloPID.h"
40 #include "AliStack.h"
41 #include "AliFiducialCut.h"
42 #include "TParticle.h"
43 #include "AliAODCaloCluster.h"
44 #include "AliVEvent.h"
45 #include "AliPHOSGeoUtils.h"
46 #include "AliEMCALGeoUtils.h"
47 #include "AliESDCaloCluster.h"
48 #include "AliESDEvent.h"
49 #include "AliAODEvent.h"
50
51 ClassImp(AliAnaPi0)
52
53 //________________________________________________________________________________________________________________________________________________  
54 AliAnaPi0::AliAnaPi0() : AliAnaPartCorrBaseClass(),
55 fNCentrBin(0),fNZvertBin(0),fNrpBin(0),
56 fNPID(0),fNmaxMixEv(0), fZvtxCut(0.),fCalorimeter(""),
57 fEMCALGeoName("EMCAL_COMPLETE"), fNModules(12), fEventsList(0x0), fhEtalon(0x0), fhReMod(0x0),
58 fhRe1(0x0),fhMi1(0x0),fhRe2(0x0),fhMi2(0x0),fhRe3(0x0),fhMi3(0x0),fhEvents(0x0),
59 fhPrimPt(0x0), fhPrimAccPt(0x0), fhPrimY(0x0), fhPrimAccY(0x0), fhPrimPhi(0x0), fhPrimAccPhi(0x0), 
60 fPHOSGeo(0x0),fEMCALGeo(0x0)
61 {
62 //Default Ctor
63  InitParameters();
64  
65 }
66
67 //________________________________________________________________________________________________________________________________________________
68 AliAnaPi0::AliAnaPi0(const AliAnaPi0 & ex) : AliAnaPartCorrBaseClass(ex),  
69 fNCentrBin(ex.fNCentrBin),fNZvertBin(ex.fNZvertBin),fNrpBin(ex.fNrpBin),
70 fNPID(ex.fNPID),fNmaxMixEv(ex.fNmaxMixEv),fZvtxCut(ex.fZvtxCut), fCalorimeter(ex.fCalorimeter),
71 fEMCALGeoName(ex.fEMCALGeoName), fNModules(ex.fNModules), fEventsList(ex.fEventsList), fhEtalon(ex.fhEtalon), fhReMod(ex.fhReMod),
72 fhRe1(ex.fhRe1),fhMi1(ex.fhMi1),fhRe2(ex.fhRe2),fhMi2(ex.fhMi2),fhRe3(ex.fhRe3),fhMi3(ex.fhMi3),fhEvents(ex.fhEvents),
73 fhPrimPt(ex.fhPrimPt), fhPrimAccPt(ex.fhPrimAccPt), fhPrimY(ex.fhPrimY), 
74 fhPrimAccY(ex.fhPrimAccY), fhPrimPhi(ex.fhPrimPhi), fhPrimAccPhi(ex.fhPrimAccPhi),
75 fPHOSGeo(ex.fPHOSGeo),fEMCALGeo(ex.fEMCALGeo)
76 {
77   // cpy ctor
78   //Do not need it
79 }
80
81 //________________________________________________________________________________________________________________________________________________
82 AliAnaPi0 & AliAnaPi0::operator = (const AliAnaPi0 & ex)
83 {
84   // assignment operator
85   
86   if(this == &ex)return *this;
87   ((AliAnaPartCorrBaseClass *)this)->operator=(ex);
88   
89   fNCentrBin = ex.fNCentrBin  ; fNZvertBin = ex.fNZvertBin  ; fNrpBin = ex.fNrpBin  ; 
90   fNPID = ex.fNPID  ; fNmaxMixEv = ex.fNmaxMixEv  ; fZvtxCut = ex.fZvtxCut  ;  fCalorimeter = ex.fCalorimeter  ;  
91   fEMCALGeoName = ex.fEMCALGeoName ; fNModules = ex.fNModules; fEventsList = ex.fEventsList  ;  fhEtalon = ex.fhEtalon  ; 
92   fhRe1 = ex.fhRe1  ; fhMi1 = ex.fhMi1  ; fhRe2 = ex.fhRe2  ; fhMi2 = ex.fhMi2  ; fhReMod = ex.fhReMod; 
93   fhRe3 = ex.fhRe3  ; fhMi3 = ex.fhMi3  ; fhEvents = ex.fhEvents  ; 
94   fhPrimPt = ex.fhPrimPt  ;  fhPrimAccPt = ex.fhPrimAccPt  ;  fhPrimY = ex.fhPrimY  ;  
95   fhPrimAccY = ex.fhPrimAccY  ;  fhPrimPhi = ex.fhPrimPhi  ;  fhPrimAccPhi = ex.fhPrimAccPhi ;
96   
97   return *this;
98   
99 }
100
101 //________________________________________________________________________________________________________________________________________________
102 AliAnaPi0::~AliAnaPi0() {
103   // Remove event containers
104   if(fEventsList){
105     for(Int_t ic=0; ic<fNCentrBin; ic++){
106       for(Int_t iz=0; iz<fNZvertBin; iz++){
107         for(Int_t irp=0; irp<fNrpBin; irp++){
108           fEventsList[ic*fNZvertBin*fNrpBin+iz*fNrpBin+irp]->Delete() ;
109           delete fEventsList[ic*fNZvertBin*fNrpBin+iz*fNrpBin+irp] ;
110         }
111       }
112     }
113     delete[] fEventsList; 
114     fEventsList=0 ;
115   }
116   
117   if(fPHOSGeo)  delete fPHOSGeo  ;
118   if(fEMCALGeo) delete fEMCALGeo ;
119         
120 }
121
122 //________________________________________________________________________________________________________________________________________________
123 void AliAnaPi0::InitParameters()
124 {
125 //Init parameters when first called the analysis
126 //Set default parameters
127   SetInputAODName("PWG4Particle");
128   
129   AddToHistogramsName("AnaPi0_");
130   fNModules = 12; // set maximum to maximum number of EMCAL modules
131   fNCentrBin = 1;
132   fNZvertBin = 1;
133   fNrpBin    = 1;
134   fNPID      = 9;
135   fNmaxMixEv = 10;
136   fZvtxCut   = 40;
137   fCalorimeter  = "PHOS";
138   fEMCALGeoName = "EMCAL_COMPLETE";
139 }
140 //________________________________________________________________________________________________________________________________________________
141 void AliAnaPi0::Init()
142 {  
143   //Init some data members needed in analysis
144   
145   //Histograms binning and range
146   if(!fhEtalon){                                                   //  p_T      alpha   d m_gg    
147     fhEtalon = new TH3D("hEtalon","Histo with binning parameters",50,0.,25.,10,0.,1.,200,0.,1.) ; 
148     fhEtalon->SetXTitle("P_{T} (GeV)") ;
149     fhEtalon->SetYTitle("#alpha") ;
150     fhEtalon->SetZTitle("m_{#gamma#gamma} (GeV)") ;
151   }
152   
153 }
154
155 //________________________________________________________________________________________________________________________________________________
156 TList * AliAnaPi0::GetCreateOutputObjects()
157 {  
158   // Create histograms to be saved in output file and 
159   // store them in fOutputContainer
160   
161   //create event containers
162   fEventsList = new TList*[fNCentrBin*fNZvertBin*fNrpBin] ;
163         
164   for(Int_t ic=0; ic<fNCentrBin; ic++){
165     for(Int_t iz=0; iz<fNZvertBin; iz++){
166       for(Int_t irp=0; irp<fNrpBin; irp++){
167         fEventsList[ic*fNZvertBin*fNrpBin+iz*fNrpBin+irp] = new TList() ;
168       }
169     }
170   }
171   
172   //If Geometry library loaded, do geometry selection during analysis.
173   printf("AliAnaPi0::GetCreateOutputObjects() - PHOS geometry initialized!\n");
174   fPHOSGeo = new AliPHOSGeoUtils("PHOSgeo") ;
175   
176   printf("AliAnaPi0::GetCreateOutputObjects() - EMCAL geometry initialized!\n");
177   fEMCALGeo = new AliEMCALGeoUtils(fEMCALGeoName) ;
178
179   TList * outputContainer = new TList() ; 
180   outputContainer->SetName(GetName()); 
181         
182   fhReMod = new TH3D*[fNModules] ;
183   fhRe1 = new TH3D*[fNCentrBin*fNPID] ;
184   fhRe2 = new TH3D*[fNCentrBin*fNPID] ;
185   fhRe3 = new TH3D*[fNCentrBin*fNPID] ;
186   fhMi1 = new TH3D*[fNCentrBin*fNPID] ;
187   fhMi2 = new TH3D*[fNCentrBin*fNPID] ;
188   fhMi3 = new TH3D*[fNCentrBin*fNPID] ;
189   
190   char key[255] ;
191   char title[255] ;
192         
193   for(Int_t ic=0; ic<fNCentrBin; ic++){
194     for(Int_t ipid=0; ipid<fNPID; ipid++){
195                 
196       //Distance to bad module 1
197       sprintf(key,"hRe_cen%d_pid%d_dist1",ic,ipid) ;
198       sprintf(title,"Real m_{#gamma#gamma} distr. for centrality=%d and PID=%d",ic,ipid) ;
199       
200       fhEtalon->Clone(key);
201       fhRe1[ic*fNPID+ipid]=(TH3D*)fhEtalon->Clone(key) ;
202       fhRe1[ic*fNPID+ipid]->SetName(key) ;
203       fhRe1[ic*fNPID+ipid]->SetTitle(title) ;
204       outputContainer->Add(fhRe1[ic*fNPID+ipid]) ;
205       
206       sprintf(key,"hMi_cen%d_pid%d_dist1",ic,ipid) ;
207       sprintf(title,"Mixed m_{#gamma#gamma} distr. for centrality=%d and PID=%d",ic,ipid) ;
208       fhMi1[ic*fNPID+ipid]=(TH3D*)fhEtalon->Clone(key) ;
209       fhMi1[ic*fNPID+ipid]->SetName(key) ;
210       fhMi1[ic*fNPID+ipid]->SetTitle(title) ;
211       outputContainer->Add(fhMi1[ic*fNPID+ipid]) ;
212       
213       //Distance to bad module 2
214       sprintf(key,"hRe_cen%d_pid%d_dist2",ic,ipid) ;
215       sprintf(title,"Real m_{#gamma#gamma} distr. for centrality=%d and PID=%d",ic,ipid) ;
216       fhRe2[ic*fNPID+ipid]=(TH3D*)fhEtalon->Clone(key) ;
217       fhRe2[ic*fNPID+ipid]->SetName(key) ;
218       fhRe2[ic*fNPID+ipid]->SetTitle(title) ;
219       outputContainer->Add(fhRe2[ic*fNPID+ipid]) ;
220       
221       sprintf(key,"hMi_cen%d_pid%d_dist2",ic,ipid) ;
222       sprintf(title,"Mixed m_{#gamma#gamma} distr. for centrality=%d and PID=%d",ic,ipid) ;
223       fhMi2[ic*fNPID+ipid]=(TH3D*)fhEtalon->Clone(key) ;
224       fhMi2[ic*fNPID+ipid]->SetName(key) ;
225       fhMi2[ic*fNPID+ipid]->SetTitle(title) ;
226       outputContainer->Add(fhMi2[ic*fNPID+ipid]) ;
227       
228       //Distance to bad module 3
229       sprintf(key,"hRe_cen%d_pid%d_dist3",ic,ipid) ;
230       sprintf(title,"Real m_{#gamma#gamma} distr. for centrality=%d and PID=%d",ic,ipid) ;
231       fhRe3[ic*fNPID+ipid]=(TH3D*)fhEtalon->Clone(key) ;
232       fhRe3[ic*fNPID+ipid]->SetName(key) ; 
233       fhRe3[ic*fNPID+ipid]->SetTitle(title) ;
234       outputContainer->Add(fhRe3[ic*fNPID+ipid]) ;
235       
236       sprintf(key,"hMi_cen%d_pid%d_dist3",ic,ipid) ;
237       sprintf(title,"Mixed m_{#gamma#gamma} distr. for centrality=%d and PID=%d",ic,ipid) ;
238       fhMi3[ic*fNPID+ipid]=(TH3D*)fhEtalon->Clone(key) ;
239       fhMi3[ic*fNPID+ipid]->SetName(key) ;
240       fhMi3[ic*fNPID+ipid]->SetTitle(title) ;
241       outputContainer->Add(fhMi3[ic*fNPID+ipid]) ;
242     }
243   }
244   
245   
246   fhEvents=new TH3D("hEvents","Number of events",fNCentrBin,0.,1.*fNCentrBin,
247                     fNZvertBin,0.,1.*fNZvertBin,fNrpBin,0.,1.*fNrpBin) ;
248   outputContainer->Add(fhEvents) ;
249   
250   //Histograms filled only if MC data is requested      
251   if(IsDataMC() || (GetReader()->GetDataType() == AliCaloTrackReader::kMC) ){
252     if(fhEtalon->GetXaxis()->GetXbins() && fhEtalon->GetXaxis()->GetXbins()->GetSize()){ //Variable bin size
253       fhPrimPt = new TH1D("hPrimPt","Primary pi0 pt",fhEtalon->GetXaxis()->GetNbins(),fhEtalon->GetXaxis()->GetXbins()->GetArray()) ;
254       fhPrimAccPt = new TH1D("hPrimAccPt","Primary pi0 pt with both photons in acceptance",fhEtalon->GetXaxis()->GetNbins(),
255                              fhEtalon->GetXaxis()->GetXbins()->GetArray()) ;
256     }
257     else{
258       fhPrimPt = new TH1D("hPrimPt","Primary pi0 pt",fhEtalon->GetXaxis()->GetNbins(),fhEtalon->GetXaxis()->GetXmin(),fhEtalon->GetXaxis()->GetXmax()) ;
259       fhPrimAccPt = new TH1D("hPrimAccPt","Primary pi0 pt with both photons in acceptance",
260                              fhEtalon->GetXaxis()->GetNbins(),fhEtalon->GetXaxis()->GetXmin(),fhEtalon->GetXaxis()->GetXmax()) ;
261     }
262     outputContainer->Add(fhPrimPt) ;
263     outputContainer->Add(fhPrimAccPt) ;
264     
265     fhPrimY = new TH1D("hPrimaryRapidity","Rapidity of primary pi0",100,-5.,5.) ; 
266     outputContainer->Add(fhPrimY) ;
267     
268     fhPrimAccY = new TH1D("hPrimAccRapidity","Rapidity of primary pi0",100,-5.,5.) ; 
269     outputContainer->Add(fhPrimAccY) ;
270     
271     fhPrimPhi = new TH1D("hPrimaryPhi","Azimithal of primary pi0",180,0.,360.) ; 
272     outputContainer->Add(fhPrimPhi) ;
273     
274     fhPrimAccPhi = new TH1D("hPrimAccPhi","Azimithal of primary pi0 with accepted daughters",180,-0.,360.) ; 
275     outputContainer->Add(fhPrimAccPhi) ;
276   }
277
278   for(Int_t imod=0; imod<fNModules; imod++){
279         //Module dependent invariant mass
280         sprintf(key,"hReMod_%d",imod) ;
281         sprintf(title,"Real m_{#gamma#gamma} distr. for Module %d",imod) ;
282         fhEtalon->Clone(key);
283         fhReMod[imod]=(TH3D*)fhEtalon->Clone(key) ;
284         fhReMod[imod]->SetName(key) ;
285         fhReMod[imod]->SetTitle(title) ;
286         outputContainer->Add(fhReMod[imod]) ;
287   }
288         
289         
290   //Save parameters used for analysis
291   TString parList ; //this will be list of parameters used for this analysis.
292   char onePar[255] ;
293   sprintf(onePar,"--- AliAnaPi0 ---\n") ;
294   parList+=onePar ;     
295   sprintf(onePar,"Number of bins in Centrality:  %d \n",fNCentrBin) ;
296   parList+=onePar ;
297   sprintf(onePar,"Number of bins in Z vert. pos: %d \n",fNZvertBin) ;
298   parList+=onePar ;
299   sprintf(onePar,"Number of bins in Reac. Plain: %d \n",fNrpBin) ;
300   parList+=onePar ;
301   sprintf(onePar,"Depth of event buffer: %d \n",fNmaxMixEv) ;
302   parList+=onePar ;
303   sprintf(onePar,"Number of different PID used:  %d \n",fNPID) ;
304   parList+=onePar ;
305   sprintf(onePar,"Cuts: \n") ;
306   parList+=onePar ;
307   sprintf(onePar,"Z vertex position: -%f < z < %f \n",fZvtxCut,fZvtxCut) ;
308   parList+=onePar ;
309   sprintf(onePar,"Calorimeter: %s \n",fCalorimeter.Data()) ;
310   parList+=onePar ;
311   sprintf(onePar,"Number of modules: %d \n",fNModules) ;
312   parList+=onePar ;
313   
314   TObjString *oString= new TObjString(parList) ;
315   outputContainer->Add(oString);
316   
317   return outputContainer;
318 }
319
320 //_________________________________________________________________________________________________________________________________________________
321 void AliAnaPi0::Print(const Option_t * /*opt*/) const
322 {
323   //Print some relevant parameters set for the analysis
324   printf("**** Print %s %s ****\n", GetName(), GetTitle() ) ;
325   AliAnaPartCorrBaseClass::Print(" ");
326
327   printf("Number of bins in Centrality:  %d \n",fNCentrBin) ;
328   printf("Number of bins in Z vert. pos: %d \n",fNZvertBin) ;
329   printf("Number of bins in Reac. Plain: %d \n",fNrpBin) ;
330   printf("Depth of event buffer: %d \n",fNmaxMixEv) ;
331   printf("Number of different PID used:  %d \n",fNPID) ;
332   printf("Cuts: \n") ;
333   printf("Z vertex position: -%2.3f < z < %2.3f \n",fZvtxCut,fZvtxCut) ;
334   printf("Number of modules: %d \n",fNModules) ;
335   printf("------------------------------------------------------\n") ;
336
337
338 //____________________________________________________________________________________________________________________________________________________
339 Int_t AliAnaPi0::GetModuleNumber(AliAODPWG4Particle * particle) 
340 {
341         //Get the EMCAL/PHOS module number that corresponds to this particle
342         
343         Int_t absId = -1;
344         if(fCalorimeter=="EMCAL"){
345                 fEMCALGeo->GetAbsCellIdFromEtaPhi(particle->Eta(),particle->Phi(), absId);
346                 if(GetDebug() > 2) 
347                         printf("EMCAL: cluster eta %f, phi %f, absid %d, SuperModule %d\n",
348                                    particle->Eta(), particle->Phi()*TMath::RadToDeg(),absId, fEMCALGeo->GetSuperModuleNumber(absId));
349                 return fEMCALGeo->GetSuperModuleNumber(absId) ;
350         }//EMCAL
351         else{//PHOS
352                 Int_t    relId[4];
353                 if(!strcmp((GetReader()->GetInputEvent())->GetName(),"AliESDEvent"))   {
354                         AliESDCaloCluster *cluster = ((AliESDEvent*)GetReader()->GetInputEvent())->GetCaloCluster(particle->GetCaloLabel(0));
355                         if ( cluster->GetNCells() > 0) {
356                                 absId = cluster->GetCellAbsId(0);
357                                 if(GetDebug() > 2) 
358                                         printf("PHOS: cluster eta %f, phi %f, e %f, e cluster %f, absId %d\n",
359                                                    particle->Eta(), particle->Phi()*TMath::RadToDeg(), particle->E(), cluster->E(), absId);
360                         }
361                         else return -1;
362                 }//ESDs
363                 else{
364                         AliAODCaloCluster *cluster = ((AliAODEvent*)GetReader()->GetInputEvent())->GetCaloCluster(particle->GetCaloLabel(0));
365                         if ( cluster->GetNCells() > 0) {
366                                 absId = cluster->GetCellAbsId(0);
367                                 if(GetDebug() > 2) 
368                                         printf("PHOS: cluster eta %f, phi %f, e %f, e cluster %f, absId %d\n",
369                                                    particle->Eta(), particle->Phi()*TMath::RadToDeg(), particle->E(), cluster->E(), absId);
370                         }
371                         else return -1;
372                 }//AODs
373
374                 if ( absId >= 0) {
375                         fPHOSGeo->AbsToRelNumbering(absId,relId);
376                         if(GetDebug() > 2) 
377                                 printf("PHOS: Module %d\n",relId[0]-1);
378                         return relId[0]-1;
379                 }
380                 else return -1;
381         }//PHOS
382         
383         return -1;
384 }
385
386 //____________________________________________________________________________________________________________________________________________________
387 void AliAnaPi0::MakeAnalysisFillHistograms() 
388 {
389   //Process one event and extract photons from AOD branch 
390   // filled with AliAnaPhoton and fill histos with invariant mass
391   
392   //Apply some cuts on event: vertex position and centrality range  
393   Int_t iRun=(GetReader()->GetInputEvent())->GetRunNumber() ;
394   if(IsBadRun(iRun)) return ;   
395   
396   Double_t vert[]={0,0,0} ; //vertex ;
397   GetReader()->GetVertex(vert);
398   if(vert[2]<-fZvtxCut || vert[2]> fZvtxCut) return ; //Event can not be used (vertex, centrality,... cuts not fulfilled)
399   
400   //Get Centrality and calculate centrality bin
401   //Does not exist in ESD yet???????
402   Int_t curCentrBin=0 ;
403   
404   //Get Reaction Plain position and calculate RP bin
405   //does not exist in ESD yet????
406   Int_t curRPBin=0 ;    
407   
408   Int_t curZvertBin=(Int_t)(0.5*fNZvertBin*(vert[2]+fZvtxCut)/fZvtxCut) ;
409   
410   fhEvents->Fill(curCentrBin+0.5,curZvertBin+0.5,curRPBin+0.5) ;
411   
412   Int_t nPhot = GetInputAODBranch()->GetEntriesFast() ;
413   if(GetDebug() > 1) printf("AliAnaPi0::MakeAnalysisFillHistograms() - Photon entries %d\n", nPhot);
414   Int_t module1 = -1;
415   Int_t module2 = -1;
416   for(Int_t i1=0; i1<nPhot-1; i1++){
417     AliAODPWG4Particle * p1 = (AliAODPWG4Particle*) (GetInputAODBranch()->At(i1)) ;
418     TLorentzVector photon1(p1->Px(),p1->Py(),p1->Pz(),p1->E());
419         //Get Module number
420         module1 = GetModuleNumber(p1);
421     for(Int_t i2=i1+1; i2<nPhot; i2++){
422       AliAODPWG4Particle * p2 = (AliAODPWG4Particle*) (GetInputAODBranch()->At(i2)) ;
423       TLorentzVector photon2(p2->Px(),p2->Py(),p2->Pz(),p2->E());
424           //Get module number
425           module2 = GetModuleNumber(p2);
426       Double_t m  = (photon1 + photon2).M() ;
427       Double_t pt = (photon1 + photon2).Pt();
428       Double_t a  = TMath::Abs(p1->E()-p2->E())/(p1->E()+p2->E()) ;
429       if(GetDebug() > 2)
430         printf("AliAnaPi0::MakeAnalysisFillHistograms() - Current Event: pT: photon1 %2.2f, photon2 %2.2f; Pair: pT %2.2f, mass %2.3f, a %f2.3\n",
431                p1->Pt(), p2->Pt(), pt,m,a);
432                 //Fill module dependent histograms
433                 //if(module1==module2) printf("mod1 %d\n",module1);
434                 if(module1==module2 && module1 >=0 && module1<fNModules)
435                         fhReMod[module1]->Fill(pt,a,m) ;
436                 
437       for(Int_t ipid=0; ipid<fNPID; ipid++)
438         {
439           if((p1->IsPIDOK(ipid,AliCaloPID::kPhoton)) && (p2->IsPIDOK(ipid,AliCaloPID::kPhoton))){ 
440             fhRe1[curCentrBin*fNPID+ipid]->Fill(pt,a,m) ;
441             if(p1->DistToBad()>0 && p2->DistToBad()>0){
442               fhRe2[curCentrBin*fNPID+ipid]->Fill(pt,a,m) ;
443               if(p1->DistToBad()>1 && p2->DistToBad()>1){
444                 fhRe3[curCentrBin*fNPID+ipid]->Fill(pt,a,m) ;
445               }
446             }
447           }
448         } 
449     }
450   }
451   
452   //Fill mixed
453   
454   TList * evMixList=fEventsList[curCentrBin*fNZvertBin*fNrpBin+curZvertBin*fNrpBin+curRPBin] ;
455   Int_t nMixed = evMixList->GetSize() ;
456   for(Int_t ii=0; ii<nMixed; ii++){  
457     TClonesArray* ev2= (TClonesArray*) (evMixList->At(ii));
458     Int_t nPhot2=ev2->GetEntriesFast() ;
459     Double_t m = -999;
460     if(GetDebug() > 1) printf("AliAnaPi0::MakeAnalysisFillHistograms() - Mixed event %d photon entries %d\n", ii, nPhot);
461     
462     for(Int_t i1=0; i1<nPhot; i1++){
463       AliAODPWG4Particle * p1 = (AliAODPWG4Particle*) (GetInputAODBranch()->At(i1)) ;
464       TLorentzVector photon1(p1->Px(),p1->Py(),p1->Pz(),p1->E());
465       for(Int_t i2=0; i2<nPhot2; i2++){
466         AliAODPWG4Particle * p2 = (AliAODPWG4Particle*) (ev2->At(i2)) ;
467         
468         TLorentzVector photon2(p2->Px(),p2->Py(),p2->Pz(),p2->E());
469         m =           (photon1+photon2).M() ; 
470         Double_t pt = (photon1 + photon2).Pt();
471         Double_t a  = TMath::Abs(p1->E()-p2->E())/(p1->E()+p2->E()) ;
472         if(GetDebug() > 2)
473           printf("AliAnaPi0::MakeAnalysisFillHistograms() - Mixed Event: pT: photon1 %2.2f, photon2 %2.2f; Pair: pT %2.2f, mass %2.3f, a %f2.3\n",
474                  p1->Pt(), p2->Pt(), pt,m,a);                   
475         for(Int_t ipid=0; ipid<fNPID; ipid++){ 
476           if((p1->IsPIDOK(ipid,AliCaloPID::kPhoton)) && (p2->IsPIDOK(ipid,AliCaloPID::kPhoton))){ 
477             fhMi1[curCentrBin*fNPID+ipid]->Fill(pt,a,m) ;
478             if(p1->DistToBad()>0 && p2->DistToBad()>0){
479               fhMi2[curCentrBin*fNPID+ipid]->Fill(pt,a,m) ;
480               if(p1->DistToBad()>1 && p2->DistToBad()>1){
481                 fhMi3[curCentrBin*fNPID+ipid]->Fill(pt,a,m) ;
482               }
483               
484             }
485           }
486         }
487       }
488     }
489   }
490   
491   TClonesArray *currentEvent = new TClonesArray(*GetInputAODBranch());
492   //Add current event to buffer and Remove redandant events 
493   if(currentEvent->GetEntriesFast()>0){
494     evMixList->AddFirst(currentEvent) ;
495     currentEvent=0 ; //Now list of particles belongs to buffer and it will be deleted with buffer
496     if(evMixList->GetSize()>=fNmaxMixEv)
497       {
498         TClonesArray * tmp = (TClonesArray*) (evMixList->Last()) ;
499         evMixList->RemoveLast() ;
500         delete tmp ;
501       }
502   } 
503   else{ //empty event
504     delete currentEvent ;
505     currentEvent=0 ; 
506   }
507   
508   //Acceptance
509   if(GetReader()->ReadStack()){ 
510           AliStack * stack = GetMCStack();
511           if(stack && (IsDataMC() || (GetReader()->GetDataType() == AliCaloTrackReader::kMC)) ){
512              for(Int_t i=0 ; i<stack->GetNprimary(); i++){
513                  TParticle * prim = stack->Particle(i) ;
514                          if(prim->GetPdgCode() == 111){
515                      Double_t pi0Pt = prim->Pt() ;
516                      //printf("pi0, pt %2.2f\n",pi0Pt);
517                     if(prim->Energy() == TMath::Abs(prim->Pz()))  continue ; //Protection against floating point exception        
518                     Double_t pi0Y  = 0.5*TMath::Log((prim->Energy()-prim->Pz())/(prim->Energy()+prim->Pz())) ;
519                     Double_t phi   = TMath::RadToDeg()*prim->Phi() ;
520                 if(TMath::Abs(pi0Y) < 0.5){
521                    fhPrimPt->Fill(pi0Pt) ;
522                 }
523                 fhPrimY  ->Fill(pi0Y) ;
524                     fhPrimPhi->Fill(phi) ;
525         
526                 //Check if both photons hit Calorimeter
527                 Int_t iphot1=prim->GetFirstDaughter() ;
528                 Int_t iphot2=prim->GetLastDaughter() ;
529                 if(iphot1>-1 && iphot1<stack->GetNtrack() && iphot2>-1 && iphot2<stack->GetNtrack()){
530                    TParticle * phot1 = stack->Particle(iphot1) ;
531                    TParticle * phot2 = stack->Particle(iphot2) ;
532                    if(phot1 && phot2 && phot1->GetPdgCode()==22 && phot2->GetPdgCode()==22){
533                    //printf("2 photons: photon 1: pt %2.2f, phi %3.2f, eta %1.2f; photon 2: pt %2.2f, phi %3.2f, eta %1.2f\n",
534                            //   phot1->Pt(), phot1->Phi()*180./3.1415, phot1->Eta(), phot2->Pt(), phot2->Phi()*180./3.1415, phot2->Eta());
535                    Bool_t inacceptance = kFALSE;
536
537                            if(fCalorimeter == "PHOS"){
538                                    if(fPHOSGeo){
539                                            Int_t mod ;
540                                            Double_t x,z ;
541                                            if(fPHOSGeo->ImpactOnEmc(phot1,mod,z,x) && fPHOSGeo->ImpactOnEmc(phot2,mod,z,x)) 
542                                                    inacceptance = kTRUE;
543                                            if(GetDebug() > 2) printf("In %s Real acceptance? %d\n",fCalorimeter.Data(),inacceptance);
544                                    }
545                                    else{
546                                            TLorentzVector lv1, lv2;
547                                            phot1->Momentum(lv1);
548                                            phot2->Momentum(lv2);
549                                            if(GetFiducialCut()->IsInFiducialCut(lv1,fCalorimeter) && GetFiducialCut()->IsInFiducialCut(lv2,fCalorimeter)) 
550                                            inacceptance = kTRUE ;
551                                            if(GetDebug() > 2) printf("In %s fiducial cut acceptance? %d\n",fCalorimeter.Data(),inacceptance);
552                                    }
553                                    
554                            }       
555                            else if(fCalorimeter == "EMCAL"){
556                                    if(fEMCALGeo){
557                                            if(fEMCALGeo->Impact(phot1) && fEMCALGeo->Impact(phot2)) 
558                                                    inacceptance = kTRUE;
559                                            if(GetDebug() > 2) printf("In %s Real acceptance? %d\n",fCalorimeter.Data(),inacceptance);
560                                         }
561                                    else{
562                                            TLorentzVector lv1, lv2;
563                                            phot1->Momentum(lv1);
564                                            phot2->Momentum(lv2);
565                                            if(GetFiducialCut()->IsInFiducialCut(lv1,fCalorimeter) && GetFiducialCut()->IsInFiducialCut(lv2,fCalorimeter)) 
566                                                    inacceptance = kTRUE ;
567                                            if(GetDebug() > 2) printf("In %s fiducial cut acceptance? %d\n",fCalorimeter.Data(),inacceptance);
568                                    }
569                    }      
570                 
571                    if(inacceptance){
572                       fhPrimAccPt->Fill(pi0Pt) ;
573                       fhPrimAccPhi->Fill(phi) ;
574                       fhPrimAccY->Fill(pi0Y) ;
575                    }//Accepted
576                }// 2 photons      
577              }//Check daughters exist
578        }// Primary pi0
579      }//loop on primaries       
580    }//stack exists and data is MC
581   }//read stack
582   else if(GetReader()->ReadAODMCParticles()){
583                 if(GetDebug() >= 0)  printf("AliAnaPi0::MakeAnalysisFillHistograms() - Acceptance calculation with MCParticles not implemented yet\n");
584   }     
585   
586 }       
587
588 //________________________________________________________________________
589 void AliAnaPi0::ReadHistograms(TList* outputList)
590 {
591         // Needed when Terminate is executed in distributed environment
592         // Refill analysis histograms of this class with corresponding histograms in output list. 
593         
594         // Histograms of this analsys are kept in the same list as other analysis, recover the position of
595         // the first one and then add the next.
596         Int_t index = outputList->IndexOf(outputList->FindObject(GetAddedHistogramsStringToName()+"hRe_cen0_pid0_dist1"));
597         
598         if(!fhRe1) fhRe1 = new TH3D*[fNCentrBin*fNPID] ;
599         if(!fhRe2) fhRe2 = new TH3D*[fNCentrBin*fNPID] ;
600         if(!fhRe3) fhRe3 = new TH3D*[fNCentrBin*fNPID] ;
601         if(!fhMi1) fhMi1 = new TH3D*[fNCentrBin*fNPID] ;
602         if(!fhMi2) fhMi2 = new TH3D*[fNCentrBin*fNPID] ;
603         if(!fhMi3) fhMi3 = new TH3D*[fNCentrBin*fNPID] ;        
604         if(!fhReMod) fhReMod = new TH3D*[fNModules] ;   
605         
606     for(Int_t ic=0; ic<fNCentrBin; ic++){
607         for(Int_t ipid=0; ipid<fNPID; ipid++){
608                         fhRe1[ic*fNPID+ipid] = (TH3D*) outputList->At(index++);
609                         fhMi1[ic*fNPID+ipid] = (TH3D*) outputList->At(index++);
610                         fhRe2[ic*fNPID+ipid] = (TH3D*) outputList->At(index++);
611                         fhMi2[ic*fNPID+ipid] = (TH3D*) outputList->At(index++);
612                         fhRe3[ic*fNPID+ipid] = (TH3D*) outputList->At(index++);
613                         fhMi3[ic*fNPID+ipid] = (TH3D*) outputList->At(index++);
614                 }
615         }
616         
617         fhEvents = (TH3D *) outputList->At(index++); 
618         
619         //Histograms filled only if MC data is requested        
620         if(IsDataMC() || (GetReader()->GetDataType() == AliCaloTrackReader::kMC) ){
621                  fhPrimPt     = (TH1D*)  outputList->At(index++);
622                  fhPrimAccPt  = (TH1D*)  outputList->At(index++);
623                  fhPrimY      = (TH1D*)  outputList->At(index++);
624                  fhPrimAccY   = (TH1D*)  outputList->At(index++);
625                  fhPrimPhi    = (TH1D*)  outputList->At(index++);
626                  fhPrimAccPhi = (TH1D*)  outputList->At(index++);
627         }
628         
629         for(Int_t imod=0; imod < fNModules; imod++)
630                         fhReMod[imod] = (TH3D*) outputList->At(index++);
631                         
632 }
633
634
635 //____________________________________________________________________________________________________________________________________________________
636 void AliAnaPi0::Terminate(TList* outputList) 
637 {
638   //Do some calculations and plots from the final histograms.
639   
640   printf(" *** %s Terminate:\n", GetName()) ; 
641  
642   //Recover histograms from output histograms list, needed for distributed analysis.    
643   ReadHistograms(outputList);
644         
645   if (!fhRe1) {
646      printf("AliAnaPi0::Terminate() - Error: Remote output histograms not imported in AliAnaPi0 object");
647      return;
648   }
649         
650   printf("AliAnaPi0::Terminate()         Mgg Real        : %5.3f , RMS : %5.3f \n", fhRe1[0]->GetMean(),   fhRe1[0]->GetRMS() ) ;
651  
652   char nameIM[128];
653   sprintf(nameIM,"AliAnaPi0_%s_cPt",fCalorimeter.Data());
654   TCanvas  * cIM = new TCanvas(nameIM, "", 400, 10, 600, 700) ;
655   cIM->Divide(2, 2);
656
657   cIM->cd(1) ; 
658   //gPad->SetLogy();
659   TH1D * hIMAllPt = (TH1D*) fhRe1[0]->ProjectionZ();
660   hIMAllPt->SetLineColor(2);
661   hIMAllPt->SetTitle("No cut on  p_{T, #gamma#gamma} ");
662   hIMAllPt->Draw();
663
664   cIM->cd(2) ; 
665   TH3F * hRe1Pt5 = (TH3F*)fhRe1[0]->Clone("IMPt5");
666   hRe1Pt5->GetXaxis()->SetRangeUser(0,5);
667   TH1D * hIMPt5 = (TH1D*) hRe1Pt5->Project3D("z");
668   hIMPt5->SetLineColor(2);  
669   hIMPt5->SetTitle("0 < p_{T, #gamma#gamma} < 5 GeV/c");
670   hIMPt5->Draw();
671   
672   cIM->cd(3) ; 
673   TH3F * hRe1Pt10 =  (TH3F*)fhRe1[0]->Clone("IMPt10");
674   hRe1Pt10->GetXaxis()->SetRangeUser(5,10);
675   TH1D * hIMPt10 = (TH1D*) hRe1Pt10->Project3D("z");
676   hIMPt10->SetLineColor(2);  
677   hIMPt10->SetTitle("5 < p_{T, #gamma#gamma} < 10 GeV/c");
678   hIMPt10->Draw();
679   
680   cIM->cd(4) ; 
681   TH3F * hRe1Pt20 =  (TH3F*)fhRe1[0]->Clone("IMPt20");
682   hRe1Pt20->GetXaxis()->SetRangeUser(10,20);
683   TH1D * hIMPt20 = (TH1D*) hRe1Pt20->Project3D("z");
684   hIMPt20->SetLineColor(2);  
685   hIMPt20->SetTitle("10 < p_{T, #gamma#gamma} < 20 GeV/c");
686   hIMPt20->Draw();
687    
688   char nameIMF[128];
689   sprintf(nameIMF,"AliAnaPi0_%s_Mgg.eps",fCalorimeter.Data());
690   cIM->Print(nameIMF);
691
692   char namePt[128];
693   sprintf(namePt,"AliAnaPi0_%s_cPt",fCalorimeter.Data());
694   TCanvas  * cPt = new TCanvas(namePt, "", 400, 10, 600, 700) ;
695   cPt->Divide(2, 2);
696
697   cPt->cd(1) ; 
698   //gPad->SetLogy();
699   TH1D * hPt = (TH1D*) fhRe1[0]->Project3D("x");
700   hPt->SetLineColor(2);
701   hPt->SetTitle("No cut on  M_{#gamma#gamma} ");
702   hPt->Draw();
703
704   cPt->cd(2) ; 
705   TH3F * hRe1IM1 = (TH3F*)fhRe1[0]->Clone("PtIM1");
706   hRe1IM1->GetZaxis()->SetRangeUser(0.05,0.21);
707   TH1D * hPtIM1 = (TH1D*) hRe1IM1->Project3D("x");
708   hPtIM1->SetLineColor(2);  
709   hPtIM1->SetTitle("0.05 < M_{#gamma#gamma} < 0.21 GeV/c^{2}");
710   hPtIM1->Draw();
711   
712   cPt->cd(3) ; 
713   TH3F * hRe1IM2 = (TH3F*)fhRe1[0]->Clone("PtIM2");
714   hRe1IM2->GetZaxis()->SetRangeUser(0.09,0.17);
715   TH1D * hPtIM2 = (TH1D*) hRe1IM2->Project3D("x");
716   hPtIM2->SetLineColor(2);  
717   hPtIM2->SetTitle("0.09 < M_{#gamma#gamma} < 0.17 GeV/c^{2}");
718   hPtIM2->Draw();
719
720   cPt->cd(4) ; 
721   TH3F * hRe1IM3 = (TH3F*)fhRe1[0]->Clone("PtIM3");
722   hRe1IM3->GetZaxis()->SetRangeUser(0.11,0.15);
723   TH1D * hPtIM3 = (TH1D*) hRe1IM1->Project3D("x");
724   hPtIM3->SetLineColor(2);  
725   hPtIM3->SetTitle("0.11 < M_{#gamma#gamma} < 0.15 GeV/c^{2}");
726   hPtIM3->Draw();
727    
728   char namePtF[128];
729   sprintf(namePtF,"AliAnaPi0_%s_Pt.eps",fCalorimeter.Data());
730   cPt->Print(namePtF);
731
732  
733   char line[1024] ; 
734   sprintf(line, ".!tar -zcf %s_%s.tar.gz *.eps", GetName(),fCalorimeter.Data()) ; 
735   gROOT->ProcessLine(line);
736   sprintf(line, ".!rm -fR AliAnaPi0_%s*.eps",fCalorimeter.Data()); 
737   gROOT->ProcessLine(line);
738  
739   printf(" AliAnaPi0::Terminate() - !! All the eps files are in %s_%s.tar.gz !!!\n", GetName(), fCalorimeter.Data());
740
741 }
742
743
744
745
746