]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrDep/AliAnaOmegaToPi0Gamma.cxx
Updates from Renzhuo to enter the mini train
[u/mrichter/AliRoot.git] / PWG4 / PartCorrDep / AliAnaOmegaToPi0Gamma.cxx
1 /**************************************************************************\r
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: The ALICE Off-line Project.                                    *\r
5  * Contributors are mentioned in the code where appropriate.              *\r
6  *                                                                        *\r
7  * Permission to use, copy, modify and distribute this software and its   *\r
8  * documentation strictly for non-commercial purposes is hereby granted   *\r
9  * without fee, provided that the above copyright notice appears in all   *\r
10  * copies and that both the copyright notice and this permission notice   *\r
11  * appear in the supporting documentation. The authors make no claims     *\r
12  * about the suitability of this software for any purpose. It is          *\r
13  * provided "as is" without express or implied warranty.                  *\r
14  **************************************************************************/\r
15 /* $Id: $ */\r
16 //_________________________________________________________________________\r
17 // class to extract omega(782)->pi0+gamma->3gamma\r
18 //\r
19 //-- Author: Renzhuo Wan (IOPP-Wuhan, China)\r
20 //_________________________________________________________________________\r
21 \r
22 // --- ROOT system\r
23 class TROOT;\r
24 \r
25 // --- AliRoot system\r
26 //class AliVEvent;\r
27 // --- ROOT system ---\r
28 #include "TH2F.h"\r
29 #include "TLorentzVector.h"\r
30 #include "TParticle.h"\r
31 #include "TCanvas.h"\r
32 #include "TFile.h"\r
33 //---- AliRoot system ----\r
34 #include "AliAnaOmegaToPi0Gamma.h"\r
35 #include "AliCaloTrackReader.h"\r
36 #include "AliCaloPID.h"\r
37 #include "AliStack.h"\r
38 #include "AliVEvent.h"\r
39 #include "AliAODEvent.h"\r
40 #include "AliAODMCParticle.h"\r
41 ClassImp(AliAnaOmegaToPi0Gamma)\r
42 \r
43 //______________________________________________________________________________\r
44 AliAnaOmegaToPi0Gamma::AliAnaOmegaToPi0Gamma() : AliAnaPartCorrBaseClass(),\r
45 fInputAODGamma(0),fInputAODPi0(0), fInputAODGammaName(""),\r
46 fEventsList(0),fNVtxZBin(1), fNCentBin(1), fNRpBin(1), fNBadChDistBin(1), fNpid(3),\r
47 fNmaxMixEv(4), fVtxZCut(0), fCent(0), fRp(0), fBadChDist(0),\r
48 fPi0Mass(0.13498),fPi0MassWindow(0.015),fPi0OverOmegaPtCut(0.8),\r
49 fGammaOverOmegaPtCut(0.3),\r
50 fhEtalon(0),\r
51 fRealOmega0(0), fMixAOmega0(0),\r
52 fMixBOmega0(0), fMixCOmega0(0),\r
53 fRealOmega1(0), fMixAOmega1(0),\r
54 fMixBOmega1(0), fMixCOmega1(0),\r
55 fRealOmega2(0), fMixAOmega2(0),\r
56 fMixBOmega2(0), fMixCOmega2(0),\r
57 fhOmegaPriPt(0)\r
58 {\r
59  //Default Ctor\r
60 }\r
61 \r
62 //______________________________________________________________________________\r
63 AliAnaOmegaToPi0Gamma::AliAnaOmegaToPi0Gamma(const AliAnaOmegaToPi0Gamma & ex) : AliAnaPartCorrBaseClass(ex),\r
64 fInputAODGamma(new TClonesArray (*ex.fInputAODGamma)), \r
65 fInputAODPi0(new TClonesArray (*ex.fInputAODPi0)),\r
66 fInputAODGammaName(ex.fInputAODGammaName),\r
67 fEventsList(ex.fEventsList), \r
68 fNVtxZBin(ex.fNVtxZBin), fNCentBin(ex.fNCentBin), fNRpBin(ex.fNRpBin),\r
69 fNBadChDistBin(ex.fNBadChDistBin),fNpid(ex.fNpid),\r
70 fNmaxMixEv(ex.fNmaxMixEv),\r
71 fVtxZCut(ex.fVtxZCut), fCent(ex.fCent), fRp(ex.fRp), fBadChDist(ex.fBadChDist),\r
72 fPi0Mass(ex.fPi0Mass),\r
73 fPi0MassWindow(ex.fPi0MassWindow),\r
74 fPi0OverOmegaPtCut(ex.fPi0OverOmegaPtCut),\r
75 fGammaOverOmegaPtCut(ex.fGammaOverOmegaPtCut),\r
76 fhEtalon(ex.fhEtalon),\r
77 fRealOmega0(ex.fRealOmega0), fMixAOmega0(ex.fMixAOmega0),\r
78 fMixBOmega0(ex.fMixBOmega0), fMixCOmega0(ex.fMixCOmega0),\r
79 fRealOmega1(ex.fRealOmega1), fMixAOmega1(ex.fMixAOmega1),\r
80 fMixBOmega1(ex.fMixBOmega1), fMixCOmega1(ex.fMixCOmega1),\r
81 fRealOmega2(ex.fRealOmega2), fMixAOmega2(ex.fMixAOmega2),\r
82 fMixBOmega2(ex.fMixBOmega2), fMixCOmega2(ex.fMixCOmega2),\r
83 fhOmegaPriPt(ex.fhOmegaPriPt)\r
84 {\r
85  // cpy ctor\r
86  //Do not need it\r
87 }\r
88 \r
89 //______________________________________________________________________________\r
90 AliAnaOmegaToPi0Gamma & AliAnaOmegaToPi0Gamma::operator = (const AliAnaOmegaToPi0Gamma & ex)\r
91 {\r
92  // assignment operator\r
93 \r
94  if(this == &ex)return *this;\r
95    ((AliAnaPartCorrBaseClass *)this)->operator=(ex);\r
96    fInputAODGamma = new TClonesArray(*ex.fInputAODGamma);\r
97    fInputAODPi0 = new TClonesArray(*ex.fInputAODPi0);\r
98 //   fInputAODGamma=ex.fInputAODGamma;\r
99 //   fInputAODPi0=ex.fInputAODPi0;\r
100    fInputAODGammaName = ex.fInputAODGammaName;\r
101    fEventsList = ex.fEventsList;\r
102 \r
103    fNVtxZBin=ex.fNVtxZBin;\r
104    fNCentBin=ex.fNCentBin;\r
105    fNRpBin=ex.fNRpBin;\r
106    fNBadChDistBin=ex.fNBadChDistBin;\r
107    fNpid=ex.fNpid;\r
108    fNmaxMixEv =ex.fNmaxMixEv;\r
109 \r
110    fVtxZCut=ex.fVtxZCut;\r
111    fCent=ex.fCent;\r
112    fRp=ex.fRp;\r
113    fBadChDist=ex.fBadChDist;\r
114 \r
115    fPi0Mass=ex.fPi0Mass;\r
116    fPi0MassWindow=ex.fPi0MassWindow;\r
117    fPi0OverOmegaPtCut=ex.fPi0OverOmegaPtCut;\r
118    fGammaOverOmegaPtCut=ex.fGammaOverOmegaPtCut;\r
119 \r
120    fhEtalon=ex.fhEtalon;\r
121    fRealOmega0=ex.fRealOmega0;\r
122    fMixAOmega0=ex.fMixAOmega0;\r
123    fMixBOmega0=ex.fMixBOmega0;\r
124    fMixCOmega0=ex.fMixCOmega0;\r
125    fRealOmega1=ex.fRealOmega1;\r
126    fMixAOmega1=ex.fMixAOmega1;\r
127    fMixBOmega1=ex.fMixBOmega1;\r
128    fMixCOmega1=ex.fMixCOmega1;\r
129    fRealOmega2=ex.fRealOmega2;\r
130    fMixAOmega2=ex.fMixAOmega2;\r
131    fMixBOmega2=ex.fMixBOmega2;\r
132    fMixCOmega2=ex.fMixCOmega2;\r
133    fhOmegaPriPt=ex.fhOmegaPriPt;\r
134   return *this;\r
135         \r
136 }\r
137 \r
138 //______________________________________________________________________________\r
139 AliAnaOmegaToPi0Gamma::~AliAnaOmegaToPi0Gamma() {\r
140 \r
141   //dtor\r
142   if(fInputAODGamma){\r
143     fInputAODGamma->Clear();\r
144     delete fInputAODGamma;\r
145   }\r
146   \r
147   if(fInputAODPi0){\r
148     fInputAODPi0->Clear();\r
149     delete fInputAODPi0;\r
150   }  \r
151 \r
152   if(fEventsList){\r
153      for(Int_t i=0;i<fNVtxZBin;i++){\r
154         for(Int_t j=0;j<fNCentBin;j++){\r
155            for(Int_t k=0;k<fNRpBin;k++){\r
156                fEventsList[i*fNCentBin*fNRpBin+j*fNRpBin+k]->Clear();\r
157                delete fEventsList[i*fNCentBin*fNRpBin+j*fNRpBin+k];\r
158            }\r
159         }\r
160      }\r
161   }\r
162   delete [] fEventsList;\r
163   fEventsList=0;\r
164 \r
165  delete [] fVtxZCut;\r
166  delete [] fCent;\r
167  delete [] fRp;\r
168  delete [] fBadChDist;\r
169 \r
170 }\r
171 \r
172 //______________________________________________________________________________\r
173 void AliAnaOmegaToPi0Gamma::InitParameters()\r
174 {\r
175 //Init parameters when first called the analysis\r
176 //Set default parameters\r
177   fVtxZCut = new Double_t [fNVtxZBin];\r
178   for(Int_t i=0;i<fNVtxZBin;i++) fVtxZCut[i]=10*(i+1);\r
179 \r
180   fCent=new Double_t[fNCentBin];\r
181   fCent[0]=0;\r
182 \r
183   fRp=new Double_t[fNRpBin];\r
184   fRp[0]=0;\r
185 \r
186   fBadChDist=new Int_t [fNBadChDistBin];\r
187   for(Int_t j=0;j<fNBadChDistBin;j++)   fBadChDist[j] = j+2;\r
188 \r
189 }\r
190 \r
191 \r
192 //______________________________________________________________________________\r
193 TList * AliAnaOmegaToPi0Gamma::GetCreateOutputObjects()\r
194 {\r
195   InitParameters();\r
196 \r
197   Int_t nptbins   = GetHistoPtBins();\r
198   Float_t ptmax   = GetHistoPtMax();\r
199   Float_t ptmin  = GetHistoPtMin();\r
200 \r
201   Int_t nmassbins = GetHistoMassBins();\r
202   Float_t massmin = GetHistoMassMin();\r
203   Float_t massmax = GetHistoMassMax();\r
204 \r
205   fhEtalon = new TH2F("hEtalon","Histo with binning parameters", nptbins,ptmin,ptmax,nmassbins,massmin,massmax) ;\r
206   fhEtalon->SetXTitle("P_{T} (GeV)") ;\r
207   fhEtalon->SetYTitle("m_{inv} (GeV)") ;\r
208  \r
209  // store them in fOutputContainer\r
210  fEventsList = new TList*[fNVtxZBin*fNCentBin*fNRpBin];\r
211  for(Int_t i=0;i<fNVtxZBin;i++){\r
212      for(Int_t j=0;j<fNCentBin;j++){\r
213          for(Int_t k=0;k<fNRpBin;k++){\r
214              fEventsList[i*fNCentBin*fNRpBin+j*fNRpBin+k]=new TList();\r
215          }\r
216      }\r
217  }\r
218 \r
219         \r
220  TList * outputContainer = new TList() ; \r
221  outputContainer->SetName(GetName());\r
222  char key[255] ;\r
223  char title[255] ;\r
224  const char * detector= fInputAODGammaName.Data();\r
225  Int_t ndim=fNVtxZBin*fNCentBin*fNRpBin*fNBadChDistBin*fNpid;\r
226 \r
227  fRealOmega0 =new TH2F*[ndim];\r
228  fMixAOmega0 =new TH2F*[ndim];\r
229  fMixBOmega0 =new TH2F*[ndim];\r
230  fMixCOmega0 =new TH2F*[ndim];\r
231  \r
232  fRealOmega1 =new TH2F*[ndim];\r
233  fMixAOmega1 =new TH2F*[ndim];\r
234  fMixBOmega1 =new TH2F*[ndim];\r
235  fMixCOmega1 =new TH2F*[ndim];\r
236 \r
237  fRealOmega2 =new TH2F*[ndim];\r
238  fMixAOmega2 =new TH2F*[ndim];\r
239  fMixBOmega2 =new TH2F*[ndim];\r
240  fMixCOmega2 =new TH2F*[ndim];\r
241 \r
242  for(Int_t i=0;i<fNVtxZBin;i++){\r
243      for(Int_t j=0;j<fNCentBin;j++){\r
244          for(Int_t k=0;k<fNRpBin;k++){ //at event level\r
245              Int_t idim=i*fNCentBin*fNRpBin+j*fNRpBin+k;\r
246              for(Int_t ipid=0;ipid<fNpid;ipid++){\r
247                  for(Int_t idist=0;idist<fNBadChDistBin;idist++){ //at particle level\r
248 \r
249                      Int_t index=idim*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;\r
250 \r
251                      sprintf(key,"RealToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
252                      sprintf(title, "%s Real Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
253                      fhEtalon->Clone(key);\r
254                      fRealOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;\r
255                      fRealOmega0[index]->SetName(key) ;\r
256                      fRealOmega0[index]->SetTitle(title);\r
257                      outputContainer->Add(fRealOmega0[index]);\r
258 \r
259                      sprintf(key,"MixAToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
260                      sprintf(title, "%s MixA Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
261                      fhEtalon->Clone(key);\r
262                      fMixAOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;\r
263                      fMixAOmega0[index]->SetName(key) ;\r
264                      fMixAOmega0[index]->SetTitle(title);\r
265                      outputContainer->Add(fMixAOmega0[index]);\r
266 \r
267                      sprintf(key,"MixBToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
268                      sprintf(title, "%s MixB Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
269                      fhEtalon->Clone(key);\r
270                      fMixBOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;\r
271                      fMixBOmega0[index]->SetName(key) ;\r
272                      fMixBOmega0[index]->SetTitle(title);\r
273                      outputContainer->Add(fMixBOmega0[index]);\r
274 \r
275                      sprintf(key,"MixCToPi0Gamma_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
276                      sprintf(title, "%s MixC Pi0GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
277                      fhEtalon->Clone(key);\r
278                      fMixCOmega0[index]=(TH2F*)fhEtalon->Clone(key) ;\r
279                      fMixCOmega0[index]->SetName(key) ;\r
280                      fMixCOmega0[index]->SetTitle(title);\r
281                      outputContainer->Add(fMixCOmega0[index]);\r
282 \r
283                      sprintf(key,"RealToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
284                      sprintf(title, "%s Real Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
285                      fhEtalon->Clone(key);\r
286                      fRealOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;\r
287                      fRealOmega1[index]->SetName(key) ;\r
288                      fRealOmega1[index]->SetTitle(title);\r
289                      outputContainer->Add(fRealOmega1[index]);\r
290 \r
291                      sprintf(key,"MixAToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
292                      sprintf(title, "%s MixA Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
293                      fhEtalon->Clone(key);\r
294                      fMixAOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;\r
295                      fMixAOmega1[index]->SetName(key) ;\r
296                      fMixAOmega1[index]->SetTitle(title);\r
297                      outputContainer->Add(fMixAOmega1[index]);\r
298 \r
299                      sprintf(key,"MixBToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
300                      sprintf(title, "%s MixB Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
301                      fhEtalon->Clone(key);\r
302                      fMixBOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;\r
303                      fMixBOmega1[index]->SetName(key) ;\r
304                      fMixBOmega1[index]->SetTitle(title);\r
305                      outputContainer->Add(fMixBOmega1[index]);\r
306 \r
307                      sprintf(key,"MixCToPi0Gamma1_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
308                      sprintf(title, "%s MixC Pi0(A<0.7)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
309                      fhEtalon->Clone(key);\r
310                      fMixCOmega1[index]=(TH2F*)fhEtalon->Clone(key) ;\r
311                      fMixCOmega1[index]->SetName(key) ;\r
312                      fMixCOmega1[index]->SetTitle(title);\r
313                      outputContainer->Add(fMixCOmega1[index]);\r
314 \r
315                      sprintf(key,"RealToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
316                      sprintf(title, "%s Real Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
317                      fhEtalon->Clone(key);\r
318                      fRealOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;\r
319                      fRealOmega2[index]->SetName(key) ;\r
320                      fRealOmega2[index]->SetTitle(title);\r
321                      outputContainer->Add(fRealOmega2[index]);\r
322 \r
323                      sprintf(key,"MixAToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
324                      sprintf(title, "%s MixA Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
325                      fhEtalon->Clone(key);\r
326                      fMixAOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;\r
327                      fMixAOmega2[index]->SetName(key) ;\r
328                      fMixAOmega2[index]->SetTitle(title);\r
329                      outputContainer->Add(fMixAOmega2[index]);\r
330 \r
331                      sprintf(key,"MixBToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
332                      sprintf(title, "%s MixB Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
333                      fhEtalon->Clone(key);\r
334                      fMixBOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;\r
335                      fMixBOmega2[index]->SetName(key) ;\r
336                      fMixBOmega2[index]->SetTitle(title);\r
337                      outputContainer->Add(fMixBOmega2[index]);\r
338 \r
339                      sprintf(key,"MixCToPi0Gamma2_Vz%dC%dRp%dPid%dDist%d",i,j,k,ipid,idist);\r
340                      sprintf(title, "%s MixC Pi0(A<0.8)GammaIVM vz_%2.1f_ct_%2.1f_Rp_%2.1f_pid_%d_dist_%dcm",detector,fVtxZCut[i],fCent[j],fRp[k],ipid,fBadChDist[idist]);\r
341                      fhEtalon->Clone(key);\r
342                      fMixCOmega2[index]=(TH2F*)fhEtalon->Clone(key) ;\r
343                      fMixCOmega2[index]->SetName(key) ;\r
344                      fMixCOmega2[index]->SetTitle(title);\r
345                      outputContainer->Add(fMixCOmega2[index]);\r
346                  }\r
347              }\r
348          }\r
349      }  \r
350  }\r
351 \r
352  if(IsDataMC()){\r
353     sprintf(key, "%sOmegaPri",detector);\r
354     sprintf(title,"primary #omega in %s",detector);\r
355     fhOmegaPriPt=new TH1F(key, title,nptbins,ptmin,ptmax);\r
356     fhOmegaPriPt->GetXaxis()->SetTitle("P_{T}");\r
357     fhOmegaPriPt->GetYaxis()->SetTitle("dN/P_{T}");\r
358     outputContainer->Add(fhOmegaPriPt);\r
359  }\r
360 \r
361  delete fhEtalon;\r
362  return outputContainer;\r
363 }\r
364 \r
365 //______________________________________________________________________________\r
366 void AliAnaOmegaToPi0Gamma::Print(const Option_t * /*opt*/) const\r
367 {\r
368   //Print some relevant parameters set in the analysis\r
369   printf("**** Print %s %s ****\n", GetName(), GetTitle() ) ;\r
370   AliAnaPartCorrBaseClass::Print(" ");\r
371   printf("Omega->pi0+gamma->3gamma\n");\r
372   printf("Cuts at event level:            \n");\r
373   printf("Bins of vertex Z:                     %d \n", fNVtxZBin);\r
374   printf("Bins of centrality:                   %d \n",fNCentBin);\r
375   printf("Bins of Reaction plane:               %d\n",fNRpBin);\r
376   printf("Cuts at AOD particle level:\n");\r
377   printf("Number of PID:                        %d \n", fNpid);\r
378   printf("Number of DistToBadChannel cuts:      %d\n", fNBadChDistBin);\r
379   printf("number of events buffer to be mixed:  %d\n",fNmaxMixEv);\r
380\r
381 \r
382 //______________________________________________________________________________\r
383 void AliAnaOmegaToPi0Gamma::MakeAnalysisFillHistograms() \r
384 {\r
385  //fill the MC AOD if needed first\r
386  //-----------\r
387  //need to be further implemented\r
388  AliStack * stack = 0x0;\r
389 // TParticle * primary = 0x0;\r
390  TClonesArray * mcparticles0 = 0x0;\r
391  TClonesArray * mcparticles1 = 0x0;\r
392  AliAODMCParticle * aodprimary = 0x0;\r
393  Int_t pdg=0;\r
394  Double_t pt=0;\r
395  Double_t eta=0;\r
396 \r
397  if(IsDataMC()){\r
398      if(GetReader()->ReadStack()){\r
399         stack =  GetMCStack() ;\r
400         if(!stack)\r
401         printf("AliAnaAcceptance::MakeAnalysisFillHistograms() - There is no stack!\n");\r
402         for(Int_t i=0 ; i<stack->GetNtrack(); i++){\r
403             TParticle * prim = stack->Particle(i) ;\r
404             pdg = prim->GetPdgCode() ;\r
405             eta=prim->Eta();\r
406             pt=prim->Pt();\r
407             if(TMath::Abs(eta)<0.5) {\r
408                if(pdg==223) fhOmegaPriPt->Fill(pt);\r
409             }\r
410         }\r
411       }\r
412       else if(GetReader()->ReadAODMCParticles()){\r
413           //Get the list of MC particles\r
414           mcparticles0 = GetReader()->GetAODMCParticles(0);\r
415           if(!mcparticles0 && GetDebug() > 0)     {\r
416              printf("AliAnaAcceptance::MakeAnalysisFillHistograms() -  Standard MCParticles not available!\n");\r
417            }\r
418            if(GetReader()->GetSecondInputAODTree()){\r
419                mcparticles1 = GetReader()->GetAODMCParticles(1);\r
420                if(!mcparticles1 && GetDebug() > 0)     {\r
421                    printf("AliAnaAcceptance::MakeAnalysisFillHistograms() -  Second input MCParticles not available!\n");\r
422                 }\r
423            }\r
424            for(Int_t i=0;i<mcparticles0->GetEntries();i++){\r
425                aodprimary =(AliAODMCParticle*)mcparticles0->At(i);\r
426                pdg = aodprimary->GetPdgCode() ;\r
427                eta=aodprimary->Eta();\r
428                pt=aodprimary->Pt();\r
429                if(TMath::Abs(eta)<0.5) {\r
430                   if(pdg==223) fhOmegaPriPt->Fill(pt);\r
431                 }\r
432 \r
433            }\r
434      }\r
435  }// is data and MC\r
436   \r
437 \r
438  //process event from AOD brach \r
439  //extract pi0, eta and omega analysis\r
440  Int_t iRun=(GetReader()->GetInputEvent())->GetRunNumber() ;\r
441  if(IsBadRun(iRun)) return ;    \r
442 \r
443  //vertex z\r
444  Double_t vert[]={0,0,0} ;\r
445  GetReader()->GetVertex(vert);\r
446  Int_t curEventBin =0;\r
447 \r
448  Int_t ivtxzbin=(Int_t)TMath::Abs(vert[2])/10;\r
449  if(ivtxzbin>=fNVtxZBin)return;\r
450 \r
451  //centrality\r
452  Int_t icentbin=0;\r
453 \r
454  //reaction plane\r
455  Int_t irpbin=0;\r
456 \r
457  if(ivtxzbin==-1) return; \r
458  curEventBin = ivtxzbin*fNCentBin*fNRpBin + icentbin*fNRpBin + irpbin;\r
459 \r
460   fInputAODGamma = GetAODBranch(fInputAODGammaName); //photon array\r
461  //fInputAODGamma = (TClonesArray *) GetReader()->GetOutputEvent()->FindListObject(fInputAODGammaName); //photon array\r
462  Int_t nphotons = fInputAODGamma->GetEntries();\r
463 \r
464  fInputAODPi0 = (TClonesArray*)GetInputAODBranch();  //pi0 array\r
465  Int_t npi0s = fInputAODPi0 ->GetEntries();\r
466  \r
467  if(nphotons<3 || npi0s<1)return; //for pi0, eta and omega->pi0+gamma->3gamma reconstruction\r
468 \r
469  //reconstruction of omega(782)->pi0+gamma->3gamma\r
470  //loop for pi0 and photon\r
471  if(GetDebug() > 0) printf("omega->pi0+gamma->3gamma invariant mass analysis ! This event have %d photons and %d pi0 \n", nphotons, npi0s);\r
472  for(Int_t i=0;i<npi0s;i++){\r
473      AliAODPWG4Particle * pi0 = (AliAODPWG4Particle*) (fInputAODPi0->At(i)) ; //pi0\r
474      TLorentzVector vpi0(pi0->Px(),pi0->Py(),pi0->Pz(),pi0->E());\r
475      Int_t lab1=pi0->GetCaloLabel(0);  // photon1 from pi0 decay\r
476      Int_t lab2=pi0->GetCaloLabel(1);  // photon2 from pi0 decay\r
477      //for omega->pi0+gamma, it needs at least three photons per event\r
478      //Get the two decay photons from pi0\r
479      AliAODPWG4Particle * photon1 =0;\r
480      AliAODPWG4Particle * photon2 =0;\r
481      for(Int_t d1=0;d1<nphotons;d1++){\r
482          for(Int_t d2=0;d2<nphotons;d2++){\r
483              AliAODPWG4Particle * dp1 = (AliAODPWG4Particle*) (fInputAODGamma->At(d1));\r
484              AliAODPWG4Particle * dp2 = (AliAODPWG4Particle*) (fInputAODGamma->At(d2));\r
485              Int_t dlab1=dp1->GetCaloLabel(0);\r
486              Int_t dlab2=dp2->GetCaloLabel(0);\r
487              if(dlab1==lab1 && dlab2==lab2){\r
488                 photon1=dp1;\r
489                 photon2=dp2;\r
490              }\r
491              else continue;\r
492          }\r
493     }\r
494 \r
495     //caculate the asy and dist of the two photon from pi0 decay\r
496     TLorentzVector dph1(photon1->Px(),photon1->Py(),photon1->Pz(),photon1->E());\r
497     TLorentzVector dph2(photon2->Px(),photon2->Py(),photon2->Pz(),photon2->E());\r
498 \r
499     Double_t pi0asy= TMath::Abs(dph1.E()-dph2.E())/(dph1.E()+dph2.E());\r
500 //    Double_t phi1=dph1.Phi();\r
501 //    Double_t phi2=dph2.Phi();\r
502 //    Double_t eta1=dph1.Eta();\r
503 //    Double_t eta2=dph2.Eta();\r
504 //    Double_t pi0dist=TMath::Sqrt((phi1-phi2)*(phi1-phi2)+(eta1-eta2)*(eta1-eta2));\r
505 \r
506     if(pi0->GetPdg()==111  && nphotons>2 && npi0s\r
507            && TMath::Abs(vpi0.M()-fPi0Mass)<fPi0MassWindow) { //pi0 candidates\r
508 \r
509          //avoid the double counting\r
510          Int_t * dc1= new Int_t[nphotons];\r
511          Int_t * dc2= new Int_t[nphotons];\r
512          Int_t index1=0;\r
513          Int_t index2=0;\r
514          for(Int_t k=0;k<i;k++){\r
515              AliAODPWG4Particle * p3=(AliAODPWG4Particle*)(fInputAODPi0->At(k));\r
516              Int_t lab4=p3->GetCaloLabel(0);\r
517              Int_t lab5=p3->GetCaloLabel(1);\r
518              if(lab1==lab4){ dc1[index1]=lab5;  index1++;  }\r
519              if(lab2==lab5){ dc2[index2]=lab4;  index2++;  }\r
520          }\r
521 \r
522 \r
523          //loop the pi0 with third gamma\r
524          for(Int_t j=0;j<nphotons;j++){\r
525              AliAODPWG4Particle *photon3 = (AliAODPWG4Particle*) (fInputAODGamma->At(j));\r
526              TLorentzVector dph3(photon3->Px(),photon3->Py(),photon3->Pz(),photon3->E());\r
527              Int_t lab3=photon3->GetCaloLabel(0);\r
528              Double_t pi0gammapt=(vpi0+dph3).Pt();\r
529              Double_t pi0gammamass=(vpi0+dph3).M();\r
530              Double_t pi0OverOmegaPtRatio =vpi0.Pt()/pi0gammapt; \r
531              Double_t gammaOverOmegaPtRatio= dph3.Pt()/pi0gammapt;\r
532 \r
533              //pi0, gamma pt cut             \r
534              if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut || \r
535                 gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;\r
536 \r
537              for(Int_t l=0;l<index1;l++) if(lab3==dc1[l]) lab3=-1;\r
538              for(Int_t l=0;l<index2;l++) if(lab3==dc2[l]) lab3=-1;\r
539     \r
540              if(lab3>0 && lab3!=lab1 && lab3!=lab2){\r
541                 for(Int_t ipid=0;ipid<fNpid;ipid++){\r
542                     for(Int_t idist=0;idist<fNBadChDistBin;idist++){\r
543                         Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;\r
544                         if(photon1->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
545                            photon2->IsPIDOK(ipid,AliCaloPID::kPhoton) && \r
546                            photon3->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
547                            photon1->DistToBad()>=fBadChDist[idist] &&\r
548                            photon2->DistToBad()>=fBadChDist[idist] &&\r
549                            photon3->DistToBad()>=fBadChDist[idist] ){\r
550                            //fill the histograms\r
551                            if(GetDebug() > 2) printf("Real: index  %d  pt  %2.3f  mass   %2.3f \n", index, pi0gammapt, pi0gammamass);\r
552                            fRealOmega0[index]->Fill(pi0gammapt,pi0gammamass);\r
553                            if(pi0asy<0.7) fRealOmega1[index]->Fill(pi0gammapt,pi0gammamass);\r
554                            if(pi0asy<0.8) fRealOmega2[index]->Fill(pi0gammapt,pi0gammamass);\r
555                         }\r
556                      }\r
557                 }\r
558              }\r
559          }      \r
560          delete []dc1;\r
561          delete []dc2;\r
562  if(GetDebug() > 0) printf("MixA: (r1_event1+r2_event1)+r3_event2 \n");\r
563          //-------------------------\r
564          //background analysis\r
565          //three background\r
566          // --A   (r1_event1+r2_event1)+r3_event2\r
567          Int_t nMixed = fEventsList[curEventBin]->GetSize();\r
568          for(Int_t im=0;im<nMixed;im++){\r
569              TClonesArray* ev2= (TClonesArray*) (fEventsList[curEventBin]->At(im));\r
570              for(Int_t mix1=0;mix1<ev2->GetEntries();mix1++){\r
571                 AliAODPWG4Particle *mix1ph = (AliAODPWG4Particle*) (ev2->At(mix1));     \r
572                 TLorentzVector vmixph(mix1ph->Px(),mix1ph->Py(),mix1ph->Pz(),mix1ph->E());\r
573                 Double_t pi0gammapt=(vpi0+vmixph).Pt();\r
574                 Double_t pi0gammamass=(vpi0+vmixph).M();\r
575                 Double_t pi0OverOmegaPtRatio =vpi0.Pt()/pi0gammapt;\r
576                 Double_t gammaOverOmegaPtRatio= vmixph.Pt()/pi0gammapt;\r
577 \r
578               //pi0, gamma pt cut             \r
579                if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut || \r
580                   gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;\r
581 \r
582                 for(Int_t ipid=0;ipid<fNpid;ipid++){\r
583                    for(Int_t idist=0;idist<fNBadChDistBin;idist++){\r
584                        Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;\r
585                        if(photon1->IsPIDOK(ipid,AliCaloPID::kPhoton)&&\r
586                           photon2->IsPIDOK(ipid,AliCaloPID::kPhoton)&&\r
587                           mix1ph->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
588                           photon1->DistToBad()>=fBadChDist[idist] &&\r
589                           photon2->DistToBad()>=fBadChDist[idist] &&\r
590                           mix1ph->DistToBad()>=fBadChDist[idist] ){\r
591                           if(GetDebug() > 2) printf("MixA: index  %d   pt  %2.3f  mass   %2.3f \n",index, pi0gammapt, pi0gammamass);\r
592                           //fill the histograms\r
593                           fMixAOmega0[index]->Fill(pi0gammapt,pi0gammamass);\r
594                           if(pi0asy<0.7)fMixAOmega1[index]->Fill(pi0gammapt,pi0gammamass);\r
595                           if(pi0asy<0.8)fMixAOmega2[index]->Fill(pi0gammapt,pi0gammamass);\r
596 //printf("mix A  %d  %2.2f \n", index, pi0gammamass);\r
597                           \r
598                        }\r
599                    }\r
600                 }\r
601              }\r
602          }\r
603      }\r
604  }\r
605 \r
606  //\r
607  // --B   (r1_event1+r2_event2)+r3_event2\r
608  //\r
609  if(GetDebug() >0)printf("MixB:  (r1_event1+r2_event2)+r3_event2 \n");\r
610  for(Int_t i=0;i<nphotons;i++){\r
611      AliAODPWG4Particle *ph1 = (AliAODPWG4Particle*) (fInputAODGamma->At(i)); \r
612      TLorentzVector vph1(ph1->Px(),ph1->Py(),ph1->Pz(),ph1->E());\r
613 \r
614      Int_t nMixed = fEventsList[curEventBin]->GetSize();\r
615      for(Int_t ie=0;ie<nMixed;ie++){\r
616          TClonesArray* ev2= (TClonesArray*) (fEventsList[curEventBin]->At(ie));\r
617          for(Int_t mix1=0;mix1<ev2->GetEntries();mix1++){\r
618              AliAODPWG4Particle *ph2 = (AliAODPWG4Particle*) (ev2->At(mix1));\r
619              TLorentzVector vph2(ph2->Px(),ph2->Py(),ph2->Pz(),ph2->E());\r
620              Double_t pi0asy = TMath::Abs(vph1.E()-vph2.E())/(vph1.E()+vph2.E());            \r
621              Double_t pi0mass=(vph1+vph2).M();\r
622                 \r
623              if(TMath::Abs(pi0mass-fPi0Mass)<fPi0MassWindow){//for pi0 selection\r
624                 for(Int_t mix2=(mix1+1);mix2<ev2->GetEntries();mix2++){\r
625                     AliAODPWG4Particle *ph3 = (AliAODPWG4Particle*) (ev2->At(mix2));\r
626                     TLorentzVector vph3(ph3->Px(),ph3->Py(),ph3->Pz(),ph3->E());\r
627                     \r
628                     Double_t pi0gammapt=(vph1+vph2+vph3).Pt();\r
629                     Double_t pi0gammamass=(vph1+vph2+vph3).M(); \r
630                     Double_t pi0OverOmegaPtRatio =(vph1+vph2).Pt()/pi0gammapt;\r
631                     Double_t gammaOverOmegaPtRatio= vph3.Pt()/pi0gammapt;\r
632                    //pi0, gamma pt cut             \r
633                    if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut ||\r
634                        gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;\r
635 \r
636                     for(Int_t ipid=0;ipid<fNpid;ipid++){\r
637                          for(Int_t idist=0;idist<fNBadChDistBin;idist++){\r
638                             Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;\r
639                             if(ph1->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
640                                ph2->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
641                                ph3->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
642                                ph1->DistToBad()>=fBadChDist[idist] &&\r
643                                ph2->DistToBad()>=fBadChDist[idist] &&\r
644                                ph3->DistToBad()>=fBadChDist[idist] ){\r
645                                if(GetDebug() > 2) printf("MixB: index  %d   pt  %2.3f  mass   %2.3f \n", index, pi0gammapt, pi0gammamass);\r
646                                //fill histograms\r
647                                fMixBOmega0[index]->Fill(pi0gammapt,pi0gammamass);\r
648                                if(pi0asy<0.7) fMixBOmega1[index]->Fill(pi0gammapt,pi0gammamass);\r
649                                if(pi0asy<0.8) fMixBOmega2[index]->Fill(pi0gammapt,pi0gammamass);\r
650 //printf("mix B  %d  %2.2f \n", index, pi0gammamass);\r
651                             }\r
652                         }                   \r
653                     }\r
654                 }\r
655         \r
656                 //\r
657                 // --C   (r1_event1+r2_event2)+r3_event3\r
658                 //\r
659                 if(GetDebug() >0)printf("MixC: (r1_event1+r2_event2)+r3_event3\n");\r
660                 for(Int_t je=(ie+1);je<nMixed;je++){\r
661                     TClonesArray* ev3= (TClonesArray*) (fEventsList[curEventBin]->At(je));\r
662                     for(Int_t mix3=0;mix3<ev3->GetEntries();mix3++){\r
663                         AliAODPWG4Particle *ph3 = (AliAODPWG4Particle*) (ev3->At(mix3));\r
664                         TLorentzVector vph3(ph3->Px(),ph3->Py(),ph3->Pz(),ph3->E());\r
665                         \r
666                         Double_t pi0gammapt=(vph1+vph2+vph3).Pt();\r
667                         Double_t pi0gammamass=(vph1+vph2+vph3).M();\r
668                         Double_t pi0OverOmegaPtRatio =(vph1+vph2).Pt()/pi0gammapt;\r
669                         Double_t gammaOverOmegaPtRatio= vph3.Pt()/pi0gammapt;\r
670                         //pi0, gamma pt cut             \r
671                         if(pi0OverOmegaPtRatio>fPi0OverOmegaPtCut ||\r
672                           gammaOverOmegaPtRatio<fGammaOverOmegaPtCut) continue;\r
673         \r
674                         for(Int_t ipid=0;ipid<fNpid;ipid++){\r
675                             for(Int_t idist=0;idist<fNBadChDistBin;idist++){\r
676                                 Int_t index=curEventBin*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;\r
677                                 if(ph1->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
678                                    ph2->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
679                                    ph3->IsPIDOK(ipid,AliCaloPID::kPhoton) &&\r
680                                    ph1->DistToBad()>=fBadChDist[idist] &&\r
681                                    ph2->DistToBad()>=fBadChDist[idist] &&\r
682                                    ph3->DistToBad()>=fBadChDist[idist] ){\r
683                                    if(GetDebug() > 2) printf("MixC: index  %d  pt  %2.3f  mass   %2.3f \n", index, pi0gammapt, pi0gammamass);\r
684                                    //fill histograms\r
685                                   fMixCOmega0[index]->Fill(pi0gammapt,pi0gammamass);\r
686                                   if(pi0asy<0.7) fMixCOmega1[index]->Fill(pi0gammapt,pi0gammamass);\r
687                                   if(pi0asy<0.8) fMixCOmega2[index]->Fill(pi0gammapt,pi0gammamass);\r
688 //printf("mix C  %d  %2.2f \n", index, pi0gammamass);\r
689                                }\r
690                            }\r
691                         }\r
692                     }\r
693                 }\r
694              } //for pi0 selecton               \r
695          }\r
696      }\r
697  }\r
698 \r
699 \r
700  //event buffer \r
701  TClonesArray *currentEvent = new TClonesArray(*fInputAODGamma);\r
702  if(currentEvent->GetEntriesFast()>0){\r
703         fEventsList[curEventBin]->AddFirst(currentEvent) ;\r
704         currentEvent=0 ; \r
705         if(fEventsList[curEventBin]->GetSize()>=fNmaxMixEv) {\r
706             TClonesArray * tmp = (TClonesArray*) (fEventsList[curEventBin]->Last()) ;\r
707             fEventsList[curEventBin]->RemoveLast() ;\r
708             delete tmp ;\r
709         }\r
710  }\r
711  else{ \r
712        delete currentEvent ;\r
713        currentEvent=0 ;\r
714  }\r
715 }\r
716 \r
717 //______________________________________________________________________________\r
718 void AliAnaOmegaToPi0Gamma::ReadHistograms(TList * outputList)\r
719 {\r
720  //read the histograms \r
721  //for the finalization of the terminate analysis\r
722 \r
723  Int_t index = outputList->IndexOf(outputList->FindObject(GetAddedHistogramsStringToName()+"RealToPi0Gamma_Vz0C0Rp0Pid0Dist0"));\r
724 \r
725   Int_t ndim=fNVtxZBin*fNCentBin*fNRpBin*fNBadChDistBin*fNpid;\r
726 \r
727  if(!fRealOmega0) fRealOmega0 =new TH2F*[ndim];\r
728  if(!fMixAOmega0) fMixAOmega0 =new TH2F*[ndim];\r
729  if(!fMixBOmega0) fMixBOmega0 =new TH2F*[ndim];\r
730  if(!fMixCOmega0) fMixCOmega0 =new TH2F*[ndim];\r
731 \r
732  if(!fRealOmega1) fRealOmega1 =new TH2F*[ndim];\r
733  if(!fMixAOmega1) fMixAOmega1 =new TH2F*[ndim];\r
734  if(!fMixBOmega1) fMixBOmega1 =new TH2F*[ndim];\r
735  if(!fMixCOmega1) fMixCOmega1 =new TH2F*[ndim];\r
736 \r
737  if(!fRealOmega2) fRealOmega2 =new TH2F*[ndim];\r
738  if(!fMixAOmega2) fMixAOmega2 =new TH2F*[ndim];\r
739  if(!fMixBOmega2) fMixBOmega2 =new TH2F*[ndim];\r
740  if(!fMixCOmega2) fMixCOmega2 =new TH2F*[ndim];\r
741 \r
742   for(Int_t i=0;i<fNVtxZBin;i++){\r
743      for(Int_t j=0;j<fNCentBin;j++){\r
744          for(Int_t k=0;k<fNRpBin;k++){ //at event level\r
745              Int_t idim=i*fNCentBin*fNRpBin+j*fNRpBin+k;\r
746              for(Int_t ipid=0;ipid<fNpid;ipid++){ \r
747                 for(Int_t idist=0;idist<fNBadChDistBin;idist++){ //at particle\r
748                     Int_t ind=idim*fNpid*fNBadChDistBin+ipid*fNBadChDistBin+idist;\r
749                     fRealOmega0[ind]= (TH2F*) outputList->At(index++);\r
750                     fMixAOmega0[ind]= (TH2F*) outputList->At(index++);\r
751                     fMixBOmega0[ind]= (TH2F*) outputList->At(index++);\r
752                     fMixCOmega0[ind]= (TH2F*) outputList->At(index++);\r
753 \r
754                     fRealOmega1[ind]= (TH2F*) outputList->At(index++);\r
755                     fMixAOmega1[ind]= (TH2F*) outputList->At(index++);\r
756                     fMixBOmega1[ind]= (TH2F*) outputList->At(index++);\r
757                     fMixCOmega1[ind]= (TH2F*) outputList->At(index++);\r
758 \r
759                     fRealOmega2[ind]= (TH2F*) outputList->At(index++);\r
760                     fMixAOmega2[ind]= (TH2F*) outputList->At(index++);\r
761                     fMixBOmega2[ind]= (TH2F*) outputList->At(index++);\r
762                     fMixCOmega2[ind]= (TH2F*) outputList->At(index++);\r
763                     \r
764                  \r
765                 }\r
766               }\r
767           }\r
768       }\r
769   }\r
770   \r
771   if(IsDataMC()){\r
772      fhOmegaPriPt  = (TH1F*)  outputList->At(index++);\r
773   }\r
774 \r
775 }\r
776 \r
777 //______________________________________________________________________________\r
778 void AliAnaOmegaToPi0Gamma::Terminate(TList * outputList) \r
779 {\r
780 // //Do some calculations and plots from the final histograms.\r
781   if(GetDebug() >= 0) printf("AliAnaOmegaToPi0Gamma::Terminate() \n");\r
782   ReadHistograms(outputList);\r
783   char cvs1[128];  \r
784   sprintf(cvs1, "Neutral_%s_IVM",fInputAODGammaName.Data());\r
785 \r
786   TCanvas * cvsIVM = new TCanvas(cvs1, cvs1, 400, 10, 600, 700) ;\r
787   cvsIVM->Divide(2, 2);\r
788 \r
789   cvsIVM->cd(1);\r
790   char dec[128];\r
791   sprintf(dec,"h2Real_%s",fInputAODGammaName.Data());\r
792   TH2F * h2Real= (TH2F*)fRealOmega0[0]->Clone(dec);\r
793   h2Real->GetXaxis()->SetRangeUser(4,6);\r
794   TH1F * hRealOmega = (TH1F*) h2Real->ProjectionY();\r
795   hRealOmega->SetTitle("RealPi0Gamma 4<pt<6");\r
796   hRealOmega->SetLineColor(2);\r
797   hRealOmega->Draw();\r
798 \r
799   cvsIVM->cd(2);\r
800   sprintf(dec,"hMixA_%s",fInputAODGammaName.Data());\r
801   TH2F *h2MixA= (TH2F*)fMixAOmega0[0]->Clone(dec);\r
802   h2MixA->GetXaxis()->SetRangeUser(4,6);\r
803   TH1F * hMixAOmega = (TH1F*) h2MixA->ProjectionY();\r
804   hMixAOmega->SetTitle("MixA 4<pt<6");\r
805   hMixAOmega->SetLineColor(2);\r
806   hMixAOmega->Draw();\r
807 \r
808   cvsIVM->cd(3);\r
809   sprintf(dec,"hMixB_%s",fInputAODGammaName.Data());\r
810   TH2F * h2MixB= (TH2F*)fMixBOmega0[0]->Clone(dec);\r
811   h2MixB->GetXaxis()->SetRangeUser(4,6);\r
812   TH1F * hMixBOmega = (TH1F*) h2MixB->ProjectionY();\r
813   hMixBOmega->SetTitle("MixB 4<pt<6");\r
814   hMixBOmega->SetLineColor(2);\r
815   hMixBOmega->Draw();\r
816 \r
817   cvsIVM->cd(4);\r
818   sprintf(dec,"hMixC_%s",fInputAODGammaName.Data());\r
819   TH2F *h2MixC= (TH2F*)fMixCOmega0[0]->Clone(dec);\r
820   h2MixC->GetXaxis()->SetRangeUser(4,6);\r
821   TH1F * hMixCOmega = (TH1F*) h2MixC->ProjectionY();\r
822   hMixCOmega->SetTitle("MixC 4<pt<6");\r
823   hMixCOmega->SetLineColor(2);\r
824   hMixCOmega->Draw();\r
825 \r
826   char eps[128 ];\r
827   sprintf(eps,"CVS_%s_IVM.eps",fInputAODGammaName.Data());\r
828   cvsIVM->Print(eps);\r
829   cvsIVM->Modified();\r
830  \r
831 }\r