]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDAnaRing.cxx
Adding a class to create dN/deta from FMD analysis output
[u/mrichter/AliRoot.git] / FMD / AliFMDAnaRing.cxx
CommitLineData
9b98d361 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15//____________________________________________________________________
16//
17// Utility class for analysing ESD data.
18// This class does sharing and background correction
6ce810fc 19// It can form a base class for other things too.
20// This line is here to make the silly code checker happy!
21// This line is here to make the silly code checker happy!
9b98d361 22//
23#include "AliFMDAnaRing.h"
24#include <TH2.h>
25#include <TMath.h>
26#include <TBrowser.h>
27// #include <AliLog.h>
6ce810fc 28#define ne 80
29#define emin -3.7
30#define emax 5.3
31#define nm 100
32#define mmin -.5
33#define mmax 9.5
34// namespace {
35// Int_t ne = 80;
36// Float_t emin = -3.7;
37// Float_t emax = 5.3;
38// Int_t nm = 100;
39// Int_t mmin = -.5;
40// Int_t mmax = 9.5;
41// }
9b98d361 42
43//____________________________________________________________________
44AliFMDAnaRing::AliFMDAnaRing()
45 : fDet(0),
46 fRing('\0'),
47 fBg(0),
48 fCut0(0),
49 fCut1(0),
50 fUseBgCor(kFALSE),
6ce810fc 51 fNSeq(0),
52 fBareMult(),
53 fMergedMult(),
54 fRemovedMult(),
55 fMult(),
56 fStep1(),
57 fStep2(),
58 fStep3(),
59 fNEvents(0)
9b98d361 60{
61 // Default constructor.
62 // Do not use.
63 // Used by ROOT I/O
64}
65
66//____________________________________________________________________
67AliFMDAnaRing::AliFMDAnaRing(UShort_t det, Char_t ring,
68 TH2* bg, Float_t c0, Float_t c1)
69 : fDet(det),
70 fRing(ring),
71 fBg(bg),
72 fCut0(c0),
73 fCut1(c1),
74 fUseBgCor(bg ? kTRUE : kFALSE),
75 fNSeq(ring == 'I' || ring == 'i' ? 20 : 40),
76 fBareMult(Form("bareMult%d%c", fDet, fRing), "Bare multiplicity",
77 ne, emin, emax, fNSeq, 0, 2*TMath::Pi()),
78 fMergedMult(Form("mergedMult%d%c", fDet, fRing), "Merged multiplicity",
79 ne, emin, emax, fNSeq, 0, 2*TMath::Pi()),
80 fRemovedMult(Form("removedMult%d%c", fDet, fRing), "Removed multiplicity",
81 ne, emin, emax, fNSeq, 0, 2*TMath::Pi()),
82 fMult(Form("mult%d%c", fDet, fRing), "Multiplicity",
83 ne, emin, emax, fNSeq, 0, 2*TMath::Pi()),
84 fStep1(Form("step1_%d%c",fDet,fRing), "Step 1", nm,mmin,mmax,nm,mmin,mmax),
85 fStep2(Form("step2_%d%c",fDet,fRing), "Step 2", nm,mmin,mmax,nm,mmin,mmax),
6ce810fc 86 fStep3(Form("step3_%d%c",fDet,fRing), "Step 3", nm,mmin,mmax,nm,mmin,mmax),
87 fNEvents(0)
9b98d361 88{
89 // Constructor
6ce810fc 90 //
9b98d361 91 // Parameters
92 // det Detector
93 // ring Ring
94 // bg Background
95 // c0 Lower cut
96 // c1 higher cut */
97 fName[0] = 'F'; fName[1] = 'M'; fName[2] = 'D';
98 fName[3] = Char_t(det+48); fName[4] = ring;
99 fName[5] = '\0';
100 fBareMult.SetDirectory(0);
101 fBareMult.SetXTitle("#eta");
102 fBareMult.SetYTitle("#varphi");
103 fBareMult.SetZTitle("d^2M_{ch}'/d#etad#varphi");
104 fMergedMult.SetDirectory(0);
105 fMergedMult.SetXTitle("#eta");
106 fMergedMult.SetYTitle("#varphi");
107 fMergedMult.SetZTitle("d^2M_{ch}''/d#etad#varphi");
108 fRemovedMult.SetDirectory(0);
109 fRemovedMult.SetXTitle("#eta");
110 fRemovedMult.SetYTitle("#varphi");
111 fRemovedMult.SetZTitle("d^2M_{ch}'''/d#etad#varphi");
112 fMult.SetDirectory(0);
113 fMult.SetXTitle("#eta");
114 fMult.SetYTitle("#varphi");
115 fMult.SetZTitle("d^2M_{ch}/d#etad#varphi");
116 fStep1.SetDirectory(0);
117 fStep1.SetXTitle("Bare multiplicity");
118 fStep1.SetYTitle("Merged multiplicity");
119 fStep2.SetDirectory(0);
120 fStep2.SetXTitle("Merged multiplicity");
121 fStep2.SetYTitle("Cut multiplicity");
122 fStep3.SetDirectory(0);
123 fStep3.SetXTitle("Cut multiplicity");
124 fStep3.SetYTitle("Multiplicity");
125}
126
127//____________________________________________________________________
128AliFMDAnaRing::AliFMDAnaRing(const AliFMDAnaRing& o)
6ce810fc 129 : TObject(o),
130 fDet(o.fDet),
9b98d361 131 fRing(o.fRing),
132 fBg(o.fBg),
133 fCut0(o.fCut0),
134 fCut1(o.fCut1),
135 fUseBgCor(o.fUseBgCor),
136 fNSeq(o.fNSeq),
137 fBareMult(o.fBareMult),
138 fMergedMult(o.fMergedMult),
139 fRemovedMult(o.fRemovedMult),
6ce810fc 140 fMult(o.fMult),
141 fStep1(o.fStep1),
142 fStep2(o.fStep2),
143 fStep3(o.fStep3),
144 fNEvents(o.fNEvents)
9b98d361 145{
146 // Copy constructor.
147 // Do not use.
148 // Used by ROOT I/O
149 fBareMult.SetDirectory(0);
150 fBareMult.SetXTitle("#eta");
151 fBareMult.SetYTitle("#varphi");
152 fBareMult.SetZTitle("dM_{ch}'/d#eta");
153 fMergedMult.SetDirectory(0);
154 fMergedMult.SetXTitle("#eta");
155 fMergedMult.SetYTitle("#varphi");
156 fMergedMult.SetZTitle("dM_{ch}''/d#eta");
157 fRemovedMult.SetDirectory(0);
158 fRemovedMult.SetXTitle("#eta");
159 fRemovedMult.SetYTitle("#varphi");
160 fRemovedMult.SetZTitle("dM_{ch}'''/d#eta");
161 fMult.SetDirectory(0);
162 fMult.SetXTitle("#eta");
163 fMult.SetYTitle("#varphi");
164 fMult.SetZTitle("dM_{ch}/d#eta");
165}
166
167//____________________________________________________________________
168AliFMDAnaRing&
169AliFMDAnaRing::operator=(const AliFMDAnaRing& o)
170{
171 // Assignment operator
6ce810fc 172 //
173 // Parameters:
174 // o Object to assign from
175 // Returns:
176 // Reference to this object.
9b98d361 177 this->fDet = o.fDet;
178 this->fRing = o.fRing;
179 this->fBg = o.fBg;
180 this->fCut0 = o.fCut0;
181 this->fCut1 = o.fCut1;
182 this->fUseBgCor = o.fUseBgCor;
183 this->fNSeq = o.fNSeq;
184 fBareMult.Reset(); fBareMult.Add(&o.fBareMult);
185 fMergedMult.Reset(); fMergedMult.Add(&o.fMergedMult);
186 fRemovedMult.Reset(); fRemovedMult.Add(&o.fRemovedMult);
187 fMult.Reset(); fMult.Add(&o.fMult);
188 return *this;
189}
190
191
192//____________________________________________________________________
193Bool_t
194AliFMDAnaRing::ProcessESD(Float_t phi, Float_t eta, Float_t& m1, Float_t m2)
195{
196 // Process ESD
197 // Parameters
198 // phi Azimuthal angle @f$ \varphi @f$ of the hit
199 // eta Psuedo-rapidity @f$ \eta@f$ of hit
200 // m1 Multiplicity of this strip. C ontains the corrected
201 // value on return
202 // m2 Multiplicity of neighbor strip
203 // Return
204 // true if hits are merged */
205
206 // Merge shared hits
207 Bool_t merged = false;
208 Float_t sig = m1;
209 fBareMult.Fill(eta, phi, m1);
210 if ((m1 < fCut0 || m2 < fCut0) && m2 > 0.001) {
211 sig = m1 + m2;
212 merged = true;
213 fRemovedMult.Fill(eta, phi, TMath::Min(m1,m2));
214 }
215 fMergedMult.Fill(eta, phi, sig);
216 fStep1.Fill(m1, sig);
217
218 // Real multiplicity
219 Double_t cmult = 1;
220 if (sig < fCut0 / 2) cmult = 0; // return merged; // cmult = 0;
221 if (sig > fCut1) cmult = 2;
222 fStep2.Fill(sig, cmult);
223
224 // Background correction
225 Float_t bmult = cmult;
226 if (fUseBgCor && fBg) {
227 // Float_t bgPhi = phi + (phi > TMath::Pi() ? -2*TMath::Pi() : 0);
228 Float_t bgPhi = phi - TMath::Pi();
229 Int_t bgBin = fBg->FindBin(eta, bgPhi);
230 Float_t bgCor = (bgBin > 0 ? fBg->GetBinContent(bgBin) : 0);
231 bmult = (bgCor < 0.01 ? 0 : cmult / bgCor);
232 }
233 fStep3.Fill(cmult, bmult);
234 fMult.Fill(eta, phi, bmult);
235
236 // Call user code
237 Fill(phi, eta, bmult);
238 m1 = bmult;
239
240 return merged;
241}
242//____________________________________________________________________
243void
244AliFMDAnaRing::Finish()
245{
6ce810fc 246 // Finish this task
247 //
248 // Parameters:
249 // none
250 // Returns:
251 // nothing
9b98d361 252 Double_t de = (fMult.GetXaxis()->GetXmax() -
253 fMult.GetXaxis()->GetXmin())/fMult.GetNbinsX();
254 Double_t dp = 2*TMath::Pi() / fNSeq;
255 Double_t scale = ((fNEvents == 0 ? 1. : 1. / fNEvents) * 1 / de * 1 / dp);
256 fBareMult.Scale(scale);
257 fMergedMult.Scale(scale);
258 fRemovedMult.Scale(scale);
259 fMult.Scale(scale);
260}
261//____________________________________________________________________
262Int_t
263AliFMDAnaRing::Color() const
264{
265 // Get the ring specific color
6ce810fc 266 //
267 // Returns:
268 // The color of the current ring
9b98d361 269 switch (fDet) {
270 case 1: return kGreen + 2;
271 case 2: return kRed + (fRing == 'I' || fRing == 'i' ? 2 : -7);
272 case 3: return kBlue + (fRing == 'I' || fRing == 'i' ? 2 : -7);
273 }
274 return 0;
275}
276//____________________________________________________________________
277void
278AliFMDAnaRing::Browse(TBrowser* b)
279{
6ce810fc 280 // Browse this object
281 //
282 // Parameters:
283 // b Browser to use
284 // Returns:
285 // nothing
9b98d361 286 b->Add(&fBareMult);
287 b->Add(&fMergedMult);
288 b->Add(&fRemovedMult);
289 b->Add(&fMult);
290 b->Add(&fStep1);
291 b->Add(&fStep2);
292 b->Add(&fStep3);
293 if (fBg) b->Add(fBg);
294}
295
296//____________________________________________________________________
297//
298// EOF
299//
300