]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FORWARD/analysis2/AliFMDCorrector.cxx
Moved reading corrections into base class AliForwardMultiplicityBase.
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliFMDCorrector.cxx
CommitLineData
7984e5f7 1//
2// This class calculates the exclusive charged particle density
3// in each for the 5 FMD rings.
4//
72cc12cd 5#include "AliFMDCorrector.h"
7e4038b5 6#include <AliESDFMD.h>
7#include <TAxis.h>
8#include <TList.h>
9#include <TMath.h>
0bd4b00f 10#include "AliForwardCorrectionManager.h"
7e4038b5 11#include "AliLog.h"
12#include <TH2D.h>
0bd4b00f 13#include <TROOT.h>
14#include <iostream>
15#include <iomanip>
7e4038b5 16
72cc12cd 17ClassImp(AliFMDCorrector)
7e4038b5 18#if 0
19; // For Emacs
20#endif
21
22//____________________________________________________________________
72cc12cd 23AliFMDCorrector::AliFMDCorrector()
7e4038b5 24 : TNamed(),
25 fRingHistos(),
7ec4d843 26 fUseSecondaryMap(true),
27 fUseVertexBias(true),
28 fUseAcceptance(true),
81eda625 29 fUseMergingEfficiency(true),
ea3e5d95 30 fDebug(0)
7984e5f7 31{
32 // Constructor
33}
7e4038b5 34
35//____________________________________________________________________
72cc12cd 36AliFMDCorrector::AliFMDCorrector(const char* title)
cc83fca2 37 : TNamed("fmdCorrector", title),
7e4038b5 38 fRingHistos(),
7ec4d843 39 fUseSecondaryMap(true),
40 fUseVertexBias(true),
41 fUseAcceptance(true),
81eda625 42 fUseMergingEfficiency(true),
ea3e5d95 43 fDebug(0)
7e4038b5 44{
7984e5f7 45 // Constructor
46 //
47 // Parameters:
48 // title Title
7e4038b5 49 fRingHistos.SetName(GetName());
50 fRingHistos.Add(new RingHistos(1, 'I'));
51 fRingHistos.Add(new RingHistos(2, 'I'));
52 fRingHistos.Add(new RingHistos(2, 'O'));
53 fRingHistos.Add(new RingHistos(3, 'I'));
54 fRingHistos.Add(new RingHistos(3, 'O'));
55}
56
57//____________________________________________________________________
72cc12cd 58AliFMDCorrector::AliFMDCorrector(const AliFMDCorrector& o)
7e4038b5 59 : TNamed(o),
60 fRingHistos(),
7ec4d843 61 fUseSecondaryMap(o.fUseSecondaryMap),
62 fUseVertexBias(o.fUseVertexBias),
63 fUseAcceptance(o.fUseAcceptance),
81eda625 64 fUseMergingEfficiency(o.fUseMergingEfficiency),
ea3e5d95 65 fDebug(o.fDebug)
7e4038b5 66{
7984e5f7 67 // Copy constructor
68 //
69 // Parameters:
70 // o Object to copy from
7e4038b5 71 TIter next(&o.fRingHistos);
72 TObject* obj = 0;
73 while ((obj = next())) fRingHistos.Add(obj);
74}
75
76//____________________________________________________________________
72cc12cd 77AliFMDCorrector::~AliFMDCorrector()
7e4038b5 78{
7984e5f7 79 // Destructor
80 //
81 //
7e4038b5 82 fRingHistos.Delete();
83}
84
85//____________________________________________________________________
72cc12cd 86AliFMDCorrector&
87AliFMDCorrector::operator=(const AliFMDCorrector& o)
7e4038b5 88{
7984e5f7 89 // Assignment operator
90 //
91 // Parameters:
92 // o Object to assign from
ea3e5d95 93 TNamed::operator=(o);
7e4038b5 94
ea3e5d95 95 fDebug = o.fDebug;
7e4038b5 96 fRingHistos.Delete();
7ec4d843 97 fUseSecondaryMap = o.fUseSecondaryMap;
98 fUseVertexBias = o.fUseVertexBias;
99 fUseAcceptance = o.fUseAcceptance;
81eda625 100 fUseMergingEfficiency = o.fUseMergingEfficiency;
7e4038b5 101 TIter next(&o.fRingHistos);
102 TObject* obj = 0;
103 while ((obj = next())) fRingHistos.Add(obj);
104
105 return *this;
106}
107
7ec4d843 108//____________________________________________________________________
109void
110AliFMDCorrector::Init(const TAxis&)
111{
112 //
113 // Initialize this object
114 //
115 // Parameters:
116 // etaAxis Eta axis to use
117 //
118 if (!fUseSecondaryMap)
119 AliWarning("Secondary maps not used - BE CAREFUL");
120 if (!fUseVertexBias)
121 AliWarning("Vertex bias not used");
122 if (!fUseAcceptance)
123 AliWarning("Acceptance from dead-channels not used");
124}
125
7e4038b5 126//____________________________________________________________________
72cc12cd 127AliFMDCorrector::RingHistos*
128AliFMDCorrector::GetRingHistos(UShort_t d, Char_t r) const
7e4038b5 129{
7984e5f7 130 //
131 // Get the ring histogram container
132 // Parameters:
133 // d Detector
134 // r Ring
135 //
136 // Return:
137 // Ring histogram container
138 //
7e4038b5 139 Int_t idx = -1;
140 switch (d) {
141 case 1: idx = 0; break;
142 case 2: idx = 1 + (r == 'I' || r == 'i' ? 0 : 1); break;
143 case 3: idx = 3 + (r == 'I' || r == 'i' ? 0 : 1); break;
144 }
145 if (idx < 0 || idx >= fRingHistos.GetEntries()) return 0;
146
147 return static_cast<RingHistos*>(fRingHistos.At(idx));
148}
149
150//____________________________________________________________________
151Bool_t
72cc12cd 152AliFMDCorrector::Correct(AliForwardUtil::Histos& hists,
7ec4d843 153 UShort_t vtxbin)
7e4038b5 154{
7984e5f7 155 //
156 // Do the calculations
157 // Parameters:
158 // hists Cache of histograms
159 // vtxBin Vertex bin
160 //
161 // Return:
162 // true on successs
163 //
0bd4b00f 164 AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
7e4038b5 165
0bd4b00f 166 UShort_t uvb = vtxbin;
7e4038b5 167 for (UShort_t d=1; d<=3; d++) {
168 UShort_t nr = (d == 1 ? 1 : 2);
169 for (UShort_t q=0; q<nr; q++) {
72cc12cd 170 Char_t r = (q == 0 ? 'I' : 'O');
171 TH2D* h = hists.Get(d,r);
172 RingHistos* rh = GetRingHistos(d,r);
7ec4d843 173
174 if (fUseSecondaryMap) {
175 TH2D* bg = fcm.GetSecondaryMap()->GetCorrection(d, r, uvb);
176 if (!bg) {
177 AliWarning(Form("No secondary correction for FMDM%d%c in vertex bin %d",
178 d, r, uvb));
179 continue;
180 }
181 // Divide by primary/total ratio
182 h->Divide(bg);
7e4038b5 183 }
7ec4d843 184 if (fUseVertexBias) {
185 TH2D* ef = fcm.GetVertexBias()->GetCorrection(r, uvb);
186 if (!ef) {
187 AliWarning(Form("No event %s vertex bias correction in vertex bin %d",
188 (r == 'I' || r == 'i' ? "inner" : "outer"), uvb));
189 continue;
190 }
191 // Divide by the event selection efficiency
192 h->Divide(ef);
72cc12cd 193 }
7ec4d843 194 if (fUseAcceptance) {
195 TH2D* ac = fcm.GetAcceptance()->GetCorrection(d, r, uvb);
196 if (!ac) {
197 AliWarning(Form("No acceptance correction for FMD%d%c in vertex bin %d",
72cc12cd 198 d, r, uvb));
7ec4d843 199 continue;
200 }
201 // Divide by the acceptance correction
202 h->Divide(ac);
7e4038b5 203 }
204
7e4038b5 205
7ec4d843 206
207
72cc12cd 208
81eda625 209 if (fUseMergingEfficiency) {
210 if (!fcm.GetMergingEfficiency()) {
211 AliWarning("No merging efficiencies");
212 continue;
213 }
214 TH1D* sf = fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb);
215 if (!fcm.GetMergingEfficiency()->GetCorrection(d,r,uvb)) {
216 AliWarning(Form("No merging efficiency for FMD%d%c at vertex bin %d",
217 d, r, uvb));
218 continue;
219 }
7e4038b5 220
81eda625 221
222 for (Int_t ieta = 1; ieta <= h->GetNbinsX(); ieta++) {
223 Float_t c = sf->GetBinContent(ieta);
224 Float_t ec = sf->GetBinError(ieta);
225
226 if (c == 0) continue;
7e4038b5 227
81eda625 228 for (Int_t iphi = 1; iphi <= h->GetNbinsY(); iphi++) {
229 Double_t m = h->GetBinContent(ieta, iphi) / c;
230 Double_t em = h->GetBinError(ieta, iphi);
7e4038b5 231
81eda625 232 Double_t e = TMath::Sqrt(em * em + (m * ec) * (m * ec));
233
234 h->SetBinContent(ieta,iphi,m);
235 h->SetBinError(ieta,iphi,e);
236 }
7e4038b5 237 }
238 }
7e4038b5 239 rh->fDensity->Add(h);
240 }
241 }
242
243 return kTRUE;
244}
245
246//____________________________________________________________________
247void
72cc12cd 248AliFMDCorrector::ScaleHistograms(TList* dir, Int_t nEvents)
7e4038b5 249{
7984e5f7 250 //
251 // Scale the histograms to the total number of events
252 // Parameters:
253 // dir Where the output is stored
254 // nEvents Number of events
255 //
7e4038b5 256 if (nEvents <= 0) return;
9d99b0dd 257 TList* d = static_cast<TList*>(dir->FindObject(GetName()));
258 if (!d) return;
7e4038b5 259
260 TIter next(&fRingHistos);
261 RingHistos* o = 0;
9d99b0dd 262 while ((o = static_cast<RingHistos*>(next())))
263 o->ScaleHistograms(d, nEvents);
7e4038b5 264}
7e4038b5 265//____________________________________________________________________
266void
72cc12cd 267AliFMDCorrector::DefineOutput(TList* dir)
7e4038b5 268{
269 TList* d = new TList;
270 d->SetName(GetName());
271 dir->Add(d);
272 TIter next(&fRingHistos);
273 RingHistos* o = 0;
274 while ((o = static_cast<RingHistos*>(next()))) {
275 o->Output(d);
276 }
277}
278
0bd4b00f 279//____________________________________________________________________
280void
72cc12cd 281AliFMDCorrector::Print(Option_t* /* option */) const
0bd4b00f 282{
7984e5f7 283 //
284 // Print information
285 // Parameters:
286 // option Not used
287 //
0bd4b00f 288 char ind[gROOT->GetDirLevel()+1];
289 for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
290 ind[gROOT->GetDirLevel()] = '\0';
7ec4d843 291 std::cout << ind << "AliFMDCorrector: " << GetName() << "\n"
292 << std::boolalpha
293 << ind << " Use secondary maps: " << fUseSecondaryMap << "\n"
294 << ind << " Use vertex bias: " << fUseVertexBias << "\n"
295 << ind << " Use acceptance: " << fUseAcceptance << "\n"
296 << ind << " Use merging efficiency: " << fUseMergingEfficiency
297 << std::endl;
0bd4b00f 298}
299
7e4038b5 300//====================================================================
72cc12cd 301AliFMDCorrector::RingHistos::RingHistos()
9d99b0dd 302 : AliForwardUtil::RingHistos(),
7e4038b5 303 fDensity(0)
7984e5f7 304{
305 // Constructor
306 //
307 //
308}
7e4038b5 309
310//____________________________________________________________________
72cc12cd 311AliFMDCorrector::RingHistos::RingHistos(UShort_t d, Char_t r)
9d99b0dd 312 : AliForwardUtil::RingHistos(d,r),
7e4038b5 313 fDensity(0)
314{
7984e5f7 315 //
316 // Constructor
317 // Parameters:
318 // d detector
319 // r ring
320 //
7e4038b5 321 fDensity = new TH2D(Form("FMD%d%c_Primary_Density", d, r),
322 Form("in FMD%d%c", d, r),
323 200, -4, 6, (r == 'I' || r == 'i' ? 20 : 40),
324 0, 2*TMath::Pi());
325 fDensity->SetDirectory(0);
326 fDensity->SetXTitle("#eta");
327 fDensity->SetYTitle("#phi [radians]");
328 fDensity->SetZTitle("Primary N_{ch} density");
329}
330//____________________________________________________________________
72cc12cd 331AliFMDCorrector::RingHistos::RingHistos(const RingHistos& o)
9d99b0dd 332 : AliForwardUtil::RingHistos(o),
7e4038b5 333 fDensity(o.fDensity)
7984e5f7 334{
335 //
336 // Copy constructor
337 // Parameters:
338 // o Object to copy from
339 //
340}
7e4038b5 341
342//____________________________________________________________________
72cc12cd 343AliFMDCorrector::RingHistos&
344AliFMDCorrector::RingHistos::operator=(const RingHistos& o)
7e4038b5 345{
7984e5f7 346 //
347 // Assignment operator
348 // Parameters:
349 // o Object to assign from
350 //
351 // Return:
352 // Reference to this
353 //
9d99b0dd 354 AliForwardUtil::RingHistos::operator=(o);
7e4038b5 355
356 if (fDensity) delete fDensity;
357
358 fDensity = static_cast<TH2D*>(o.fDensity->Clone());
359
360 return *this;
361}
362//____________________________________________________________________
72cc12cd 363AliFMDCorrector::RingHistos::~RingHistos()
7e4038b5 364{
7984e5f7 365 //
366 // Destructor
367 //
7e4038b5 368 if (fDensity) delete fDensity;
369}
370
371//____________________________________________________________________
372void
72cc12cd 373AliFMDCorrector::RingHistos::Output(TList* dir)
7e4038b5 374{
7984e5f7 375 //
376 // Make output
377 // Parameters:
378 // dir Where to put it
379 //
9d99b0dd 380 TList* d = DefineOutputList(dir);
7e4038b5 381 d->Add(fDensity);
9d99b0dd 382}
383
384//____________________________________________________________________
385void
72cc12cd 386AliFMDCorrector::RingHistos::ScaleHistograms(TList* dir, Int_t nEvents)
9d99b0dd 387{
7984e5f7 388 //
389 // Scale the histograms to the total number of events
390 // Parameters:
391 // dir where the output is stored
392 // nEvents Number of events
393 //
9d99b0dd 394 TList* l = GetOutputList(dir);
395 if (!l) return;
396
397 TH1* density = GetOutputHist(l,Form("%s_Primary_Density", fName.Data()));
398 if (density) density->Scale(1./nEvents);
7e4038b5 399}
400
401//____________________________________________________________________
402//
403// EOF
404//
405
406
407