]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDReconstructor.cxx
Fixing copy/pase pub (K.Shileev)
[u/mrichter/AliRoot.git] / FMD / AliFMDReconstructor.cxx
CommitLineData
37c55dc0 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
0d0e6995 16/* $Id$ */
17
4347b38f 18//____________________________________________________________________
42403906 19//
7684b53c 20// This is a class that constructs AliFMDMult (reconstructed
21// multiplicity) from of Digits
4347b38f 22//
23// This class reads either digits from a TClonesArray or raw data from
24// a DDL file (or similar), and stores the read ADC counts in an
25// internal cache (fAdcs).
26//
27// From the cached values it then calculates the number of particles
28// that hit a region of the FMDs, as specified by the user.
29//
30// The reconstruction can be done in two ways: Either via counting the
31// number of empty strips (Poisson method), or by converting the ADC
32// signal to an energy deposition, and then dividing by the typical
33// energy loss of a particle.
34//
7684b53c 35// +---------------------+ +---------------------+
36// | AliFMDReconstructor |<>-----| AliFMDMultAlgorithm |
37// +---------------------+ +---------------------+
38// ^
39// |
40// +-----------+---------+
41// | |
42// +-------------------+ +------------------+
43// | AliFMDMultPoisson | | AliFMDMultNaiive |
44// +-------------------+ +------------------+
45//
46// AliFMDReconstructor acts as a manager class. It contains a list of
47// AliFMDMultAlgorithm objects. The call graph looks something like
48//
49//
50// +----------------------+ +----------------------+
51// | :AliFMDReconstructor | | :AliFMDMultAlgorithm |
52// +----------------------+ +----------------------+
53// | |
54// Reconstruct +-+ |
55// ------------>| | PreRun +-+
56// | |------------------------------->| |
57// | | +-+
58// | |-----+ (for each event) |
59// | | | *ProcessEvent |
60// |+-+ | |
61// || |<---+ PreEvent +-+
62// || |------------------------------>| |
63// || | +-+
64// || |-----+ |
65// || | | ProcessDigits |
66// ||+-+ | |
67// ||| |<---+ |
68// ||| | *ProcessDigit(digit) +-+
69// ||| |----------------------------->| |
70// ||| | +-+
71// ||+-+ |
72// || | PostEvent +-+
73// || |------------------------------>| |
74// || | +-+
75// |+-+ |
76// | | PostRun +-+
77// | |------------------------------->| |
78// | | +-+
79// +-+ |
80// | |
81//
4347b38f 82//
37c55dc0 83//
84//-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
4347b38f 85// Latest changes by Christian Holm Christensen <cholm@nbi.dk>
86//
87//
88//____________________________________________________________________
37c55dc0 89
56b1929b 90#include <AliLog.h> // ALILOG_H
91#include <AliRun.h> // ALIRUN_H
92#include <AliRunLoader.h> // ALIRUNLOADER_H
93#include <AliLoader.h> // ALILOADER_H
94#include <AliHeader.h> // ALIHEADER_H
95#include <AliRawReader.h> // ALIRAWREADER_H
96#include <AliGenEventHeader.h> // ALIGENEVENTHEADER_H
1a1fdef7 97#include "AliFMD.h" // ALIFMD_H
98#include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
99#include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
100#include "AliFMDDetector.h" // ALIFMDDETECTOR_H
101#include "AliFMDRing.h" // ALIFMDRING_H
56b1929b 102#include "AliFMDDigit.h" // ALIFMDDIGIT_H
103#include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
104#include "AliFMDRawStream.h" // ALIFMDRAWSTREAM_H
105#include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
106#include "AliFMDMultAlgorithm.h" // ALIFMDMULTALGORITHM_H
107#include "AliFMDMultPoisson.h" // ALIFMDMULTPOISSON_H
108#include "AliFMDMultNaiive.h" // ALIFMDMULTNAIIVE_H
4347b38f 109
110//____________________________________________________________________
925e6570 111ClassImp(AliFMDReconstructor)
1a1fdef7 112#if 0
113 ; // This is here to keep Emacs for indenting the next line
114#endif
37c55dc0 115
4347b38f 116//____________________________________________________________________
117AliFMDReconstructor::AliFMDReconstructor()
118 : AliReconstructor(),
4347b38f 119 fPedestal(0),
e802be3e 120 fPedestalWidth(0),
121 fPedestalFactor(0)
4347b38f 122{
42403906 123 // Make a new FMD reconstructor object - default CTOR.
1a1fdef7 124 AliFMDParameters* pars = AliFMDParameters::Instance();
125 fPedestal = pars->GetPedestal();
126 fPedestalWidth = pars->GetPedestalWidth();
127 fPedestalFactor = pars->GetPedestalFactor();
128
56b1929b 129 fAlgorithms.Add(new AliFMDMultNaiive);
130 fAlgorithms.Add(new AliFMDMultPoisson);
4347b38f 131}
132
42403906 133
134//____________________________________________________________________
135AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other)
136 : AliReconstructor(),
42403906 137 fPedestal(0),
e802be3e 138 fPedestalWidth(0),
139 fPedestalFactor(0)
42403906 140{
56b1929b 141 // Copy constructor
1a1fdef7 142 fPedestal = other.fPedestal;
143 fPedestalWidth = other.fPedestalWidth;
144 fPedestalFactor = other.fPedestalFactor;
42403906 145
56b1929b 146
147 fAlgorithms.Delete();
148 TIter next(&(other.fAlgorithms));
149 AliFMDMultAlgorithm* algorithm = 0;
150 while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next())))
151 fAlgorithms.Add(algorithm);
152 fAlgorithms.SetOwner(kFALSE);
42403906 153}
154
155
156//____________________________________________________________________
157AliFMDReconstructor&
158AliFMDReconstructor::operator=(const AliFMDReconstructor& other)
159{
56b1929b 160 // Assignment operator
1a1fdef7 161 fPedestal = other.fPedestal;
162 fPedestalWidth = other.fPedestalWidth;
163 fPedestalFactor = other.fPedestalFactor;
42403906 164
56b1929b 165 fAlgorithms.Delete();
166 TIter next(&(other.fAlgorithms));
167 AliFMDMultAlgorithm* algorithm = 0;
168 while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next())))
169 fAlgorithms.Add(algorithm);
170 fAlgorithms.SetOwner(kFALSE);
171
42403906 172 return *this;
173}
56b1929b 174
175//____________________________________________________________________
176AliFMDReconstructor::~AliFMDReconstructor()
177{
178 // Destructor
179 fAlgorithms.Delete();
180}
4347b38f 181
182//____________________________________________________________________
183void
1a1fdef7 184AliFMDReconstructor::Init(AliRunLoader* runLoader)
185{
186 // Initialize the reconstructor
187 AliDebug(1, Form("Init called with runloader 0x%x", runLoader));
188 fCurrentVertex = 0;
189 if (!runLoader) {
190 Warning("Init", "No run loader");
4347b38f 191 return;
192 }
1a1fdef7 193 AliHeader* header = runLoader->GetHeader();
194 if (!header) {
195 Warning("Init", "No header");
4347b38f 196 return;
197 }
1a1fdef7 198 AliGenEventHeader* eventHeader = header->GenEventHeader();
199 if (!eventHeader) {
200 Warning("Init", "no event header");
201 return;
4347b38f 202 }
1a1fdef7 203 TArrayF vtx;
204 eventHeader->PrimaryVertex(vtx);
205 fCurrentVertex = vtx[2];
206 AliDebug(1, Form("Primary vertex Z coordinate for event # %d/%d is %f",
207 header->GetRun(), header->GetEvent(), fCurrentVertex));
4347b38f 208}
dc8af42e 209
4347b38f 210//____________________________________________________________________
211void
1a1fdef7 212AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
213 TTree* digitsTree) const
214{
215 // Convert Raw digits to AliFMDDigit's in a tree
216 AliFMDRawReader rawRead(reader, digitsTree);
217 // rawRead.SetSampleRate(fFMD->GetSampleRate());
218 rawRead.Exec();
4347b38f 219}
220
4347b38f 221//____________________________________________________________________
222void
1a1fdef7 223AliFMDReconstructor::Reconstruct(TTree* digitsTree,
224 TTree* clusterTree) const
4347b38f 225{
1a1fdef7 226 // Reconstruct event from digits in tree
4347b38f 227 // Get the FMD branch holding the digits.
1a1fdef7 228 AliDebug(1, "Reconstructing from digits in a tree");
229
230 TBranch *digitBranch = digitsTree->GetBranch("FMD");
231 TClonesArray* digits = new TClonesArray("AliFMDDigit");
4347b38f 232 if (!digitBranch) {
1a1fdef7 233 Error("Exec", "No digit branch for the FMD found");
234 return;
4347b38f 235 }
1a1fdef7 236 digitBranch->SetAddress(&digits);
4347b38f 237
e802be3e 238 TIter next(&fAlgorithms);
56b1929b 239 AliFMDMultAlgorithm* algorithm = 0;
240 while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next())))
1a1fdef7 241 algorithm->PreEvent(clusterTree, fCurrentVertex);
242 digitBranch->GetEntry(0);
243
e802be3e 244 ProcessDigits(digits);
245
246 next.Reset();
247 algorithm = 0;
56b1929b 248 while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next())))
e802be3e 249 algorithm->PostEvent();
1a1fdef7 250 clusterTree->Fill();
4347b38f 251}
1a1fdef7 252
4347b38f 253//____________________________________________________________________
e802be3e 254void
255AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
4347b38f 256{
e802be3e 257 Int_t nDigits = digits->GetEntries();
1a1fdef7 258 AliDebug(1, Form("Got %d digits", nDigits));
e802be3e 259 for (Int_t i = 0; i < nDigits; i++) {
260 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
1a1fdef7 261 AliFMDGeometry* fmd = AliFMDGeometry::Instance();
262 AliFMDDetector* subDetector = fmd->GetDetector(digit->Detector());
e802be3e 263 if (!subDetector) {
264 Warning("ProcessDigits", "Unknown detector: FMD%d" , digit->Detector());
265 continue;
266 }
4347b38f 267
1a1fdef7 268 AliFMDRing* ring = subDetector->GetRing(digit->Ring());
269 Float_t ringZ = subDetector->GetRingZ(digit->Ring());
e802be3e 270 if (!ring) {
271 Warning("ProcessDigits", "Unknown ring: FMD%d%c", digit->Detector(),
272 digit->Ring());
273 break;
4347b38f 274 }
275
e802be3e 276 Float_t realZ = fCurrentVertex + ringZ;
56b1929b 277 Float_t stripR = ((ring->GetHighR() - ring->GetLowR())
278 / ring->GetNStrips() * (digit->Strip() + .5)
279 + ring->GetLowR());
e802be3e 280 Float_t theta = TMath::ATan2(stripR, realZ);
281 Float_t phi = (2 * TMath::Pi() / ring->GetNSectors()
282 * (digit->Sector() + .5));
283 Float_t eta = -TMath::Log(TMath::Tan(theta / 2));
284 UShort_t counts = SubtractPedestal(digit);
4347b38f 285
e802be3e 286 TIter next(&fAlgorithms);
56b1929b 287 AliFMDMultAlgorithm* algorithm = 0;
288 while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next())))
e802be3e 289 algorithm->ProcessDigit(digit, eta, phi, counts);
4347b38f 290 }
4347b38f 291}
4347b38f 292
1a1fdef7 293//____________________________________________________________________
294UShort_t
295AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
296{
297 // Member function to subtract the pedestal from a digit
298 // This implementation does nothing, but a derived class could over
299 // load this to subtract a pedestal that was given in a database or
300 // something like that.
301
302 Int_t counts = 0;
303 Float_t ped = fPedestal + fPedestalFactor * fPedestalWidth;
304 if (digit->Count3() > 0) counts = digit->Count3();
305 else if (digit->Count2() > 0) counts = digit->Count2();
306 else counts = digit->Count1();
307 counts = TMath::Max(Int_t(counts - ped), 0);
308 return UShort_t(counts);
309}
310
4347b38f 311//____________________________________________________________________
312void
1a1fdef7 313AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
314 TTree* /* clusterTree */,
315 AliESD* /* esd*/) const
121a60bd 316{
42403906 317 // nothing to be done
121a60bd 318
319}
320
42403906 321//____________________________________________________________________
322//
323// EOF
324//