]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/GammaConv/AliGammaConversionHistograms.cxx
Fix coding convention violations.
[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.1                                                            *\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   fBackgroundContainer(NULL),\r
48   fDebugContainer(NULL),\r
49   fResolutionContainer(NULL),\r
50   fMatchContainer(NULL),\r
51   fESDContainer(NULL),\r
52   fMCContainer(NULL),\r
53   fOtherContainer(NULL)\r
54 {\r
55   // see header file for documenation\r
56 }\r
57 \r
58 \r
59 AliGammaConversionHistograms::AliGammaConversionHistograms(const AliGammaConversionHistograms & original) :\r
60   fHistogramMap(original.fHistogramMap),\r
61   fNPhiIndex(original.fNPhiIndex),\r
62   fNRIndex(original.fNRIndex),\r
63   fMinRadius(original.fMinRadius),\r
64   fMaxRadius(original.fMaxRadius),\r
65   fDeltaR(original.fDeltaR),\r
66   fMinPhi(original.fMinPhi),\r
67   fMaxPhi(original.fMaxPhi),\r
68   fDeltaPhi(original.fDeltaPhi),\r
69   fMappingContainer(original.fMappingContainer),\r
70   fBackgroundContainer(original.fBackgroundContainer),\r
71   fDebugContainer(original.fDebugContainer),\r
72   fResolutionContainer(original.fResolutionContainer),\r
73   fMatchContainer(original.fMatchContainer),\r
74   fESDContainer(original.fESDContainer),\r
75   fMCContainer(original.fMCContainer),\r
76   fOtherContainer(original.fOtherContainer)\r
77 {    \r
78   //see header file for documentation\r
79 }\r
80 \r
81 \r
82 AliGammaConversionHistograms & AliGammaConversionHistograms::operator = (const AliGammaConversionHistograms & /*original*/)\r
83 {\r
84   // assignment operator\r
85   return *this;\r
86 }\r
87 \r
88 \r
89 AliGammaConversionHistograms::~AliGammaConversionHistograms() {\r
90   //destructor\r
91         \r
92         \r
93 }\r
94 \r
95 void AliGammaConversionHistograms::AddHistogram(TString histogramName, TString histogramTitle, Int_t nXBins, Double_t firstX,Double_t lastX,TString xAxisTitle, TString yAxisTitle){\r
96   // see header file for documentation\r
97   TH1F *tmp = new TH1F(histogramName, histogramTitle,nXBins,firstX,lastX);\r
98   tmp->GetXaxis()->SetTitle(xAxisTitle);\r
99   tmp->GetYaxis()->SetTitle(yAxisTitle);\r
100   TObjString* tobjstring = new TObjString(histogramName.Data());\r
101   fHistogramMap->Add((TObject*)tobjstring,(TObject*)tmp);\r
102 }\r
103 \r
104 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
105   // see header file for documentation\r
106   TH2F *tmp = new TH2F(histogramName, histogramTitle,nXBins,firstX,lastX,nYBins,firstY,lastY);\r
107   tmp->GetXaxis()->SetTitle(xAxisTitle);\r
108   tmp->GetYaxis()->SetTitle(yAxisTitle);\r
109   TObjString *tobjstring = new TObjString(histogramName.Data());\r
110   fHistogramMap->Add((TObject*)tobjstring,(TObject*)tmp);\r
111 }\r
112 \r
113 void AliGammaConversionHistograms::FillHistogram(TString histogramName, Double_t xValue) const{\r
114   //see header file for documentation\r
115   TH1 *tmp = (TH1*)fHistogramMap->GetValue(histogramName.Data());\r
116   if(tmp){\r
117     tmp->Fill(xValue);\r
118   }\r
119 }\r
120 \r
121 void AliGammaConversionHistograms::FillHistogram(TString histogramName, Double_t xValue, Double_t yValue) const{\r
122   //see header file for documentation\r
123   TH1 *tmp = (TH1*)fHistogramMap->GetValue(histogramName.Data());\r
124   if(tmp){\r
125     tmp->Fill(xValue, yValue);\r
126   }\r
127 }\r
128 \r
129 void AliGammaConversionHistograms::GetOutputContainer(TList *fOutputContainer){\r
130   //checking if the container is alrerady created\r
131         \r
132   if(fOutputContainer == NULL){\r
133     //print warning\r
134     return;\r
135   }\r
136         \r
137   if(fHistogramMap != NULL){\r
138     TIter iter(fHistogramMap);\r
139     TObjString *histogramName;\r
140     while ((histogramName = (TObjString*) iter.Next())) {\r
141       TString histogramString = histogramName->GetString();\r
142       if(histogramString.Contains("Mapping")){// means it should be put in the mapping folder\r
143         if(fMappingContainer == NULL){\r
144           fMappingContainer = new TList();\r
145           fMappingContainer->SetName("Mapping histograms");\r
146         }\r
147         if(fMappingContainer != NULL){\r
148           fMappingContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
149         }\r
150       }\r
151       else if(histogramString.Contains("Background")){// means it should be put in the background folder\r
152         if(fBackgroundContainer == NULL){\r
153           fBackgroundContainer = new TList();\r
154           fBackgroundContainer->SetName("Background histograms");\r
155         }\r
156         if(fBackgroundContainer != NULL){\r
157           fBackgroundContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
158         }\r
159       }\r
160       else if(histogramString.Contains("Debug")){// means it should be put in the debug folder\r
161         if(fDebugContainer == NULL){\r
162           fDebugContainer = new TList();\r
163           fDebugContainer->SetName("Debug histograms");\r
164         }\r
165         if(fDebugContainer != NULL){\r
166           fDebugContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
167         }\r
168       }\r
169       else if(histogramString.Contains("Resolution")){// means it should be put in the resolution folder\r
170         if(fResolutionContainer == NULL){\r
171           fResolutionContainer = new TList();\r
172           fResolutionContainer->SetName("Resolution histograms");\r
173         }\r
174         if(fResolutionContainer != NULL){\r
175           fResolutionContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
176         }\r
177       }\r
178       else if(histogramString.Contains("Match")){// means it should be put in the mapping folder\r
179         if(fMatchContainer == NULL){\r
180           fMatchContainer = new TList();\r
181           fMatchContainer->SetName("Match histograms");\r
182         }\r
183         if(fMatchContainer != NULL){\r
184           fMatchContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
185         }\r
186       }\r
187       else if(histogramString.Contains("ESD")){// means it should be put in the ESD folder\r
188         if(fESDContainer == NULL){\r
189           fESDContainer = new TList();\r
190           fESDContainer->SetName("ESD histograms");\r
191         }\r
192         if(fESDContainer != NULL){\r
193           fESDContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
194         }\r
195       }\r
196       else if(histogramString.Contains("MC")){// means it should be put in the MC folder\r
197         if(fMCContainer == NULL){\r
198           fMCContainer = new TList();\r
199           fMCContainer->SetName("MC histograms");\r
200         }\r
201         if(fMCContainer != NULL){\r
202           fMCContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
203         }\r
204       }\r
205       else{\r
206         if(fOtherContainer == NULL){\r
207           fOtherContainer = new TList();\r
208           fOtherContainer->SetName("Other histograms");\r
209         }\r
210         if(fOtherContainer != NULL){\r
211           fOtherContainer->Add((TH1*)fHistogramMap->GetValue(histogramString.Data()));\r
212         }\r
213       }\r
214       histogramName = NULL;\r
215     } // end while\r
216     if(fMappingContainer != NULL){\r
217       fOutputContainer->Add(fMappingContainer);\r
218     }\r
219     if(fBackgroundContainer != NULL){\r
220       fOutputContainer->Add(fBackgroundContainer);\r
221     }\r
222     if(fDebugContainer != NULL){\r
223       fOutputContainer->Add(fDebugContainer);\r
224     }\r
225     if(fResolutionContainer != NULL){\r
226       fOutputContainer->Add(fResolutionContainer);\r
227     }\r
228     if(fMatchContainer != NULL){\r
229       fOutputContainer->Add(fMatchContainer);\r
230     }\r
231     if(fESDContainer != NULL){\r
232       fOutputContainer->Add(fESDContainer);\r
233     }\r
234     if(fMCContainer != NULL){\r
235       fOutputContainer->Add(fMCContainer);\r
236     }\r
237     if(fOtherContainer != NULL){\r
238       fOutputContainer->Add(fMCContainer);\r
239     }\r
240   }\r
241 }\r
242 \r
243 Int_t AliGammaConversionHistograms::GetRBin(Double_t radius) const{\r
244   // see header file for documentation\r
245   Int_t iResult=0;\r
246   if(fDeltaR>0){\r
247     iResult = (Int_t)((radius - fMinRadius)/fDeltaR);\r
248   }\r
249   return iResult;\r
250 }\r
251 \r
252 Int_t AliGammaConversionHistograms::GetPhiBin(Double_t phi) const{\r
253   // see header file for documentation\r
254   Int_t iResult=0;\r
255   if(fDeltaPhi>0){\r
256     if(phi>TMath::Pi()){\r
257       phi-=2*TMath::Pi();\r
258     }\r
259     iResult = (Int_t)((phi - fMinPhi)/fDeltaPhi);\r
260   }\r
261   return iResult;\r
262 }\r
263 \r
264 \r
265 \r
266 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
267   // Initializing the valuse for the mapping\r
268         \r
269   fNPhiIndex = nPhiIndex;\r
270   fNRIndex   = nRIndex;\r
271   fMinRadius      = minRadius;\r
272   fMaxRadius      = maxRadius;\r
273   if(nBinsR>0 && nRIndex!=0){\r
274     fDeltaR       = (fMaxRadius - fMinRadius)/nRIndex;\r
275   }\r
276   fMinPhi         = minPhi;\r
277   fMaxPhi         = maxPhi;\r
278   if(nBinsPhi>0 && nPhiIndex!=0){\r
279     fDeltaPhi     = (fMaxPhi-fMinPhi)/nPhiIndex;\r
280   }\r
281 }\r
282 \r
283 \r
284 //mapping\r
285 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
286   // see header file for documentation\r
287         \r
288   for(Int_t phi =0; phi<=fNPhiIndex;phi++){\r
289                 \r
290     for(Int_t r =0; r<fNRIndex;r++){\r
291                         \r
292       // setting axis to "" changes below\r
293       xAxisTitle="";\r
294       yAxisTitle="";\r
295       //Creating the axis titles\r
296       if(xAxisTitle.Length() == 0){\r
297         xAxisTitle.Form("Phi %02d",phi);\r
298       }\r
299                         \r
300       if(yAxisTitle.Length() == 0){\r
301         yAxisTitle.Form("R %02d",phi);\r
302       }\r
303                         \r
304       //MC\r
305       TString nameMC="";\r
306       nameMC.Form("MC_Conversion_Mapping-Phi%02d-R%02d",phi,r);\r
307       TString titleMC="";\r
308       titleMC.Form("Electron-Positron MC Mapping-Phi%02d-R%02d",phi,r);\r
309                         \r
310       AddHistogram(nameMC, titleMC, nXBins, firstX, lastX, nYBins, firstY, lastY, xAxisTitle, yAxisTitle);\r
311                         \r
312       //ESD\r
313       TString nameESD="";\r
314       nameESD.Form("ESD_Conversion_Mapping-Phi%02d-R%02d",phi,r);\r
315       TString titleESD="";\r
316       titleESD.Form("Electron-Positron ESD Mapping-Phi%02d-R%02d",phi,r);\r
317                         \r
318       AddHistogram(nameESD, titleESD, nXBins, firstX, lastX, nYBins, firstY, lastY, xAxisTitle, yAxisTitle);\r
319     }\r
320   }\r
321         \r
322         \r
323   for(Int_t phi =0; phi<=nPhiIndex;phi++){ \r
324                 \r
325     // setting axis to "" changes below\r
326     xAxisTitle="";\r
327     yAxisTitle="";\r
328     //Creating the axis titles\r
329     if(xAxisTitle.Length() == 0){\r
330       xAxisTitle.Form("Phi %02d",phi);\r
331     }\r
332     if(yAxisTitle.Length() == 0){\r
333       yAxisTitle = "Counts";\r
334     }\r
335                 \r
336     //MC\r
337     TString nameMC="";\r
338     nameMC.Form("MC_Conversion_Mapping-Phi%02d",phi);\r
339     TString titleMC="";\r
340     titleMC.Form("Electron-Positron MC Mapping-Phi%02d",phi);\r
341                 \r
342     AddHistogram(nameMC, titleMC, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
343                 \r
344     //MC\r
345     TString nameESD="";\r
346     nameESD.Form("ESD_Conversion_Mapping-Phi%02d",phi);\r
347     TString titleESD="";\r
348     titleESD.Form("Electron-Positron ESD Mapping-Phi%02d",phi);\r
349                 \r
350     AddHistogram(nameESD, titleESD, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
351   }\r
352         \r
353         \r
354   for(Int_t r =0; r<=nRIndex;r++){\r
355                 \r
356     // setting axis to "" changes below\r
357     xAxisTitle="";\r
358     yAxisTitle="";\r
359     //Creating the axis titles\r
360     if(xAxisTitle.Length() == 0){\r
361       xAxisTitle.Form("R %02d",r);\r
362     }\r
363     if(yAxisTitle.Length() == 0){\r
364       yAxisTitle = "Counts";\r
365     }\r
366                 \r
367     //MC\r
368     TString nameMC="";\r
369     nameMC.Form("MC_Conversion_Mapping-R%02d",r);\r
370     TString titleMC="";\r
371     titleMC.Form("Electron-Positron MC Mapping-R%02d",r);\r
372                 \r
373     AddHistogram(nameMC, titleMC, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
374                 \r
375     //ESD\r
376     TString nameESD="";\r
377     nameESD.Form("ESD_Conversion_Mapping-R%02d",r);\r
378     TString titleESD="";\r
379     titleESD.Form("Electron-Positron ESD Mapping-R%02d",r);\r
380                 \r
381     AddHistogram(nameESD, titleESD, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
382                 \r
383     //Mapping Phi in R\r
384     TString nameMCPhiInR="";\r
385     nameMCPhiInR.Form("MC_Conversion_Mapping_Phi_R-%02d",r);\r
386     TString titleMCPhiInR="";\r
387     titleMCPhiInR.Form("MC Mapping of Phi in R%02d",r);\r
388     AddHistogram(nameMCPhiInR, titleMCPhiInR, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);\r
389                 \r
390     //Mapping Phi in R\r
391     TString nameESDPhiInR="";\r
392     nameESDPhiInR.Form("ESD_Conversion_Mapping_Phi_R-%02d",r);\r
393     TString titleESDPhiInR="";\r
394     titleESDPhiInR.Form("ESD Mapping of Phi in R%02d",r);\r
395     AddHistogram(nameESDPhiInR, titleESDPhiInR, nXBins, firstX, lastX, xAxisTitle, yAxisTitle);    \r
396   }\r
397 }\r