]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/GammaConv/AliGammaConversionHistograms.cxx
Gamma conversion analysis moved to its own subdirectory from PartCorr
[u/mrichter/AliRoot.git] / PWG4 / GammaConv / AliGammaConversionHistograms.cxx
1 /**************************************************************************\r
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: Ana Marin, Kathrin Koch, Kenneth Aamodt                        *\r
5  * Version 1.0                                                            *\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 \r
16 ////////////////////////////////////////////////\r
17 //--------------------------------------------- \r
18 // Class used to do analysis on conversion pairs\r
19 //---------------------------------------------\r
20 ////////////////////////////////////////////////\r
21 \r
22 #include "AliGammaConversionHistograms.h"\r
23 #include "TMath.h"\r
24 #include "TObjString.h"\r
25 #include "TMap.h"\r
26 #include "TList.h"\r
27 #include "TH1F.h"\r
28 #include "TH2F.h"\r
29 \r
30 \r
31 using namespace std;\r
32 \r
33 ClassImp(AliGammaConversionHistograms)\r
34 \r
35 \r
36 AliGammaConversionHistograms::AliGammaConversionHistograms() :\r
37   fHistogramMap(new TMap()),\r
38   fNPhiIndex(0),\r
39   fNRIndex(0),\r
40   fMinRadius(0.),\r
41   fMaxRadius(0.),\r
42   fDeltaR(0.),\r
43   fMinPhi(0.),\r
44   fMaxPhi(0.),\r
45   fDeltaPhi(0.),\r
46   fMappingContainer(NULL)\r
47 {\r
48   // see header file for documenation\r
49 }\r
50 \r
51 \r
52 AliGammaConversionHistograms::AliGammaConversionHistograms(const AliGammaConversionHistograms & original) :\r
53   fHistogramMap(original.fHistogramMap),\r
54   fNPhiIndex(original.fNPhiIndex),\r
55   fNRIndex(original.fNRIndex),\r
56   fMinRadius(original.fMinRadius),\r
57   fMaxRadius(original.fMaxRadius),\r
58   fDeltaR(original.fDeltaR),\r
59   fMinPhi(original.fMinPhi),\r
60   fMaxPhi(original.fMaxPhi),\r
61   fDeltaPhi(original.fDeltaPhi),\r
62   fMappingContainer(original.fMappingContainer)\r
63 {    \r
64   //see header file for documentation\r
65 }\r
66 \r
67 \r
68 AliGammaConversionHistograms & AliGammaConversionHistograms::operator = (const AliGammaConversionHistograms & /*original*/)\r
69 {\r
70   // assignment operator\r
71   return *this;\r
72 }\r
73 \r
74 \r
75 AliGammaConversionHistograms::~AliGammaConversionHistograms() {\r
76   //destructor\r
77 \r
78 \r
79 }\r
80 \r
81 void AliGammaConversionHistograms::AddHistogram(TString histogramName, TString histogramTitle, Int_t nXBins, Double_t firstX,Double_t lastX,TString xAxisTitle, TString yAxisTitle){\r
82   // see header file for documentation\r
83   TH1F *tmp = new TH1F(histogramName, histogramTitle,nXBins,firstX,lastX);\r
84   tmp->GetXaxis()->SetTitle(xAxisTitle);\r
85   tmp->GetYaxis()->SetTitle(yAxisTitle);\r
86   TObjString* tobjstring = new TObjString(histogramName.Data());\r
87   fHistogramMap->Add((TObject*)tobjstring,(TObject*)tmp);\r
88 }\r
89 \r
90 void AliGammaConversionHistograms::AddHistogram(TString histogramName, TString histogramTitle, Int_t nXBins, Double_t firstX, Double_t lastX, Int_t nYBins, Double_t firstY, Double_t lastY, TString xAxisTitle, TString yAxisTitle){\r
91   // see header file for documentation\r
92   TH2F *tmp = new TH2F(histogramName, histogramTitle,nXBins,firstX,lastX,nYBins,firstY,lastY);\r
93   tmp->GetXaxis()->SetTitle(xAxisTitle);\r
94   tmp->GetYaxis()->SetTitle(yAxisTitle);\r
95   TObjString *tobjstring = new TObjString(histogramName.Data());\r
96   fHistogramMap->Add((TObject*)tobjstring,(TObject*)tmp);\r
97 }\r
98 \r
99 void AliGammaConversionHistograms::FillHistogram(TString histogramName, Double_t xValue) const{\r
100   //see header file for documentation\r
101   TH1 *tmp = (TH1*)fHistogramMap->GetValue(histogramName.Data());\r
102   if(tmp){\r
103       tmp->Fill(xValue);\r
104   }\r
105   else{\r
106     cout<<"Histogram does not exist: "<<histogramName.Data()<<endl;\r
107   }\r
108 }\r
109 \r
110 void AliGammaConversionHistograms::FillHistogram(TString histogramName, Double_t xValue, Double_t yValue) const{\r
111   //see header file for documentation\r
112   TH1 *tmp = (TH1*)fHistogramMap->GetValue(histogramName.Data());\r
113   if(tmp){\r
114     tmp->Fill(xValue, yValue);\r
115   }\r
116   else{\r
117     cout<<"Histogram does not exist: "<<histogramName.Data()<<endl;\r
118   }\r
119 }\r
120 \r
121 void AliGammaConversionHistograms::GetOutputContainer(TList *fOutputContainer){\r
122   //checking if the container is alrerady created\r
123 \r
124   if(fOutputContainer == NULL){\r
125     //print warning\r
126     return;\r
127   }\r
128 \r
129   if(fHistogramMap != NULL){\r
130     TIter iter(fHistogramMap);\r
131     TObjString *histogramName;\r
132     while ((histogramName = (TObjString*) iter.Next())) {\r
133       TString histogramString = histogramName->GetString();\r
134       if(histogramString.Contains("Mapping")){// means it should be put in the mapping folder\r
135         if(fMappingContainer == NULL){\r
136           fMappingContainer = new TList();\r
137           fMappingContainer->SetName("Mapping histograms");\r
138         }\r
139         if(fMappingContainer != NULL){\r
140           fMappingContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
141         }\r
142       }\r
143       fOutputContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
144       histogramName = NULL;\r
145     } // end while\r
146     //   fOutputContainer->Add(fMappingContainer);\r
147   }\r
148 }\r
149 \r
150 Int_t AliGammaConversionHistograms::GetRBin(Double_t radius) const{\r
151   // see header file for documentation\r
152   Int_t iResult=0;\r
153   if(fDeltaR>0){\r
154     iResult = (Int_t)((radius - fMinRadius)/fDeltaR);\r
155   }\r
156   return iResult;\r
157 }\r
158 \r
159 Int_t AliGammaConversionHistograms::GetPhiBin(Double_t phi) const{\r
160   // see header file for documentation\r
161   Int_t iResult=0;\r
162   if(fDeltaPhi>0){\r
163     if(phi>TMath::Pi()){\r
164       phi-=2*TMath::Pi();\r
165     }\r
166     iResult = (Int_t)((phi - fMinPhi)/fDeltaPhi);\r
167   }\r
168   return iResult;\r
169 }\r
170 \r
171 \r
172 \r
173 void AliGammaConversionHistograms::InitializeMappingValues(Int_t nPhiIndex, Int_t nRIndex, Int_t nBinsR, Double_t minRadius, Double_t maxRadius,Int_t nBinsPhi, Double_t minPhi, Double_t maxPhi){\r
174   // Initializing the valuse for the mapping\r
175 \r
176   fNPhiIndex = nPhiIndex;\r
177   fNRIndex   = nRIndex;\r
178   fMinRadius      = minRadius;\r
179   fMaxRadius      = maxRadius;\r
180   if(nBinsR>0 && nRIndex!=0){\r
181     fDeltaR       = (fMaxRadius - fMinRadius)/nRIndex;\r
182   }\r
183   fMinPhi         = minPhi;\r
184   fMaxPhi         = maxPhi;\r
185   if(nBinsPhi>0 && nPhiIndex!=0){\r
186     fDeltaPhi     = (fMaxPhi-fMinPhi)/nPhiIndex;\r
187   }\r
188 }\r
189 \r
190 \r
191 //mapping\r
192 void AliGammaConversionHistograms::AddMappingHistograms(Int_t nPhiIndex, Int_t nRIndex,Int_t nXBins, Double_t firstX, Double_t lastX, Int_t nYBins, Double_t firstY, Double_t lastY, TString xAxisTitle, TString yAxisTitle){\r
193   // see header file for documentation\r
194   \r
195   for(Int_t phi =0; phi<=fNPhiIndex;phi++){\r
196 \r
197     for(Int_t r =0; r<fNRIndex;r++){\r
198 \r
199       // setting axis to "" changes below\r
200       xAxisTitle="";\r
201       yAxisTitle="";\r
202       //Creating the axis titles\r
203       if(xAxisTitle.Length() == 0){\r
204         xAxisTitle.Form("Phi %02d",phi);\r
205       }\r
206       \r
207       if(yAxisTitle.Length() == 0){\r
208         yAxisTitle.Form("R %02d",phi);\r
209       }\r
210 \r
211       //MC\r
212       TString nameMC="";\r
213       nameMC.Form("MC_EP_Mapping-Phi%02d-R%02d",phi,r);\r
214       TString titleMC="";\r
215       titleMC.Form("Electron-Positron MC Mapping-Phi%02d-R%02d",phi,r);\r
216 \r
217       AddHistogram(nameMC, titleMC, nXBins, firstX, lastX, nYBins, firstY, lastY, xAxisTitle, yAxisTitle);\r
218 \r
219       //ESD\r
220       TString nameESD="";\r
221       nameESD.Form("ESD_EP_Mapping-Phi%02d-R%02d",phi,r);\r
222       TString titleESD="";\r
223       titleESD.Form("Electron-Positron ESD Mapping-Phi%02d-R%02d",phi,r);\r
224 \r
225       AddHistogram(nameESD, titleESD, nXBins, firstX, lastX, nYBins, firstY, lastY, xAxisTitle, yAxisTitle);\r
226     }\r
227   }\r
228 \r
229 \r
230   for(Int_t phi =0; phi<=nPhiIndex;phi++){ \r
231 \r
232     // setting axis to "" changes below\r
233     xAxisTitle="";\r
234     yAxisTitle="";\r
235     //Creating the axis titles\r
236     if(xAxisTitle.Length() == 0){\r
237       xAxisTitle.Form("Phi %02d",phi);\r
238     }\r
239     if(yAxisTitle.Length() == 0){\r
240       yAxisTitle = "Counts";\r
241     }\r
242     \r
243     //MC\r
244     TString nameMC="";\r
245     nameMC.Form("MC_EP_Mapping-Phi%02d",phi);\r
246     TString titleMC="";\r
247     titleMC.Form("Electron-Positron MC Mapping-Phi%02d",phi);\r
248     \r
249     AddHistogram(nameMC, titleMC, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
250 \r
251     //MC\r
252     TString nameESD="";\r
253     nameESD.Form("ESD_EP_Mapping-Phi%02d",phi);\r
254     TString titleESD="";\r
255     titleESD.Form("Electron-Positron ESD Mapping-Phi%02d",phi);\r
256     \r
257     AddHistogram(nameESD, titleESD, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
258   }\r
259 \r
260 \r
261   for(Int_t r =0; r<=nRIndex;r++){\r
262 \r
263     // setting axis to "" changes below\r
264     xAxisTitle="";\r
265     yAxisTitle="";\r
266     //Creating the axis titles\r
267     if(xAxisTitle.Length() == 0){\r
268       xAxisTitle.Form("R %02d",r);\r
269     }\r
270     if(yAxisTitle.Length() == 0){\r
271       yAxisTitle = "Counts";\r
272     }\r
273     \r
274     //MC\r
275     TString nameMC="";\r
276     nameMC.Form("MC_EP_Mapping-R%02d",r);\r
277     TString titleMC="";\r
278     titleMC.Form("Electron-Positron MC Mapping-R%02d",r);\r
279     \r
280     AddHistogram(nameMC, titleMC, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
281 \r
282     //ESD\r
283     TString nameESD="";\r
284     nameESD.Form("ESD_EP_Mapping-R%02d",r);\r
285     TString titleESD="";\r
286     titleESD.Form("Electron-Positron ESD Mapping-R%02d",r);\r
287     \r
288     AddHistogram(nameESD, titleESD, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
289 \r
290     //Mapping Phi in R\r
291     TString nameMCPhiInR="";\r
292     nameMCPhiInR.Form("MC_EP_Mapping_Phi_R-%02d",r);\r
293     TString titleMCPhiInR="";\r
294     titleMCPhiInR.Form("MC Mapping of Phi in R%02d",r);\r
295     AddHistogram(nameMCPhiInR, titleMCPhiInR, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
296 \r
297     //Mapping Phi in R\r
298     TString nameESDPhiInR="";\r
299     nameESDPhiInR.Form("ESD_EP_Mapping_Phi_R-%02d",r);\r
300     TString titleESDPhiInR="";\r
301     titleESDPhiInR.Form("ESD Mapping of Phi in R%02d",r);\r
302     AddHistogram(nameESDPhiInR, titleESDPhiInR, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);    \r
303   }\r
304 }\r