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