1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /** @file AliFMDReconstructor.cxx
17 @author Christian Holm Christensen <cholm@nbi.dk>
18 @date Mon Mar 27 12:47:09 2006
19 @brief FMD reconstruction
21 //____________________________________________________________________
23 // This is a class that constructs AliFMDMult (reconstructed
24 // multiplicity) from of Digits
26 // This class reads either digits from a TClonesArray or raw data from
27 // a DDL file (or similar), and stores the read ADC counts in an
28 // internal cache (fAdcs).
30 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
31 // Latest changes by Christian Holm Christensen <cholm@nbi.dk>
34 //____________________________________________________________________
36 #include <AliLog.h> // ALILOG_H
37 #include <AliRun.h> // ALIRUN_H
38 #include <AliRunLoader.h> // ALIRUNLOADER_H
39 #include <AliLoader.h> // ALILOADER_H
40 #include <AliHeader.h> // ALIHEADER_H
41 #include <AliRawReader.h> // ALIRAWREADER_H
42 #include <AliGenEventHeader.h> // ALIGENEVENTHEADER_H
43 #include "AliFMD.h" // ALIFMD_H
44 #include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
45 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
46 #include "AliFMDDetector.h" // ALIFMDDETECTOR_H
47 #include "AliFMDRing.h" // ALIFMDRING_H
48 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
49 #include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
50 #include "AliFMDRawStream.h" // ALIFMDRAWSTREAM_H
51 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
52 #include "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
53 #include "AliESD.h" // ALIESD_H
54 #include <AliESDFMD.h> // ALIESDFMD_H
57 //____________________________________________________________________
58 ClassImp(AliFMDReconstructor)
60 ; // This is here to keep Emacs for indenting the next line
63 //____________________________________________________________________
64 AliFMDReconstructor::AliFMDReconstructor()
73 // Make a new FMD reconstructor object - default CTOR.
77 //____________________________________________________________________
78 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other)
83 fCurrentVertex(other.fCurrentVertex),
84 fESDObj(other.fESDObj),
91 //____________________________________________________________________
93 AliFMDReconstructor::operator=(const AliFMDReconstructor& other)
95 // Assignment operator
97 fNMult = other.fNMult;
98 fTreeR = other.fTreeR;
99 fCurrentVertex = other.fCurrentVertex;
100 fESDObj = other.fESDObj;
105 //____________________________________________________________________
106 AliFMDReconstructor::~AliFMDReconstructor()
109 if (fMult) fMult->Delete();
110 if (fMult) delete fMult;
111 if (fESDObj) delete fESDObj;
114 //____________________________________________________________________
116 AliFMDReconstructor::Init(AliRunLoader* runLoader)
118 // Initialize the reconstructor
119 AliDebug(1, Form("Init called with runloader 0x%x", runLoader));
120 // Initialize the geometry
121 AliFMDGeometry* geom = AliFMDGeometry::Instance();
123 geom->InitTransformations();
125 // Initialize the parameters
126 AliFMDParameters* param = AliFMDParameters::Instance();
129 // Current vertex position
131 // Create array of reconstructed strip multiplicities
132 fMult = new TClonesArray("AliFMDRecPoint", 51200);
133 // Create ESD output object
134 fESDObj = new AliESDFMD;
136 // Check that we have a run loader
138 Warning("Init", "No run loader");
142 // Check if we can get the header tree
143 AliHeader* header = runLoader->GetHeader();
145 Warning("Init", "No header");
149 // Check if we can get a simulation header
150 AliGenEventHeader* eventHeader = header->GenEventHeader();
153 eventHeader->PrimaryVertex(vtx);
154 fCurrentVertex = vtx[2];
155 AliDebug(1, Form("Primary vertex Z coordinate for event # %d/%d is %f",
156 header->GetRun(), header->GetEvent(), fCurrentVertex));
157 Warning("Init", "no generator event header");
160 Warning("Init", "No generator event header - "
161 "perhaps we get the vertex from ESD?");
165 //____________________________________________________________________
167 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
168 TTree* digitsTree) const
170 // Convert Raw digits to AliFMDDigit's in a tree
171 AliDebug(1, "Reading raw data into digits tree");
172 AliFMDRawReader rawRead(reader, digitsTree);
173 // rawRead.SetSampleRate(fFMD->GetSampleRate());
177 //____________________________________________________________________
179 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
180 TTree* clusterTree) const
182 // Reconstruct event from digits in tree
183 // Get the FMD branch holding the digits.
184 // FIXME: The vertex may not be known yet, so we may have to move
185 // some of this to FillESD.
186 AliDebug(1, "Reconstructing from digits in a tree");
189 const AliESDVertex* vertex = fESD->GetVertex();
191 AliDebug(1, Form("Got vertex from ESD: %f", vertex->GetZv()));
192 fCurrentVertex = vertex->GetZv();
196 TBranch *digitBranch = digitsTree->GetBranch("FMD");
198 Error("Exec", "No digit branch for the FMD found");
201 TClonesArray* digits = new TClonesArray("AliFMDDigit");
202 digitBranch->SetAddress(&digits);
204 if (fMult) fMult->Clear();
205 if (fESDObj) fESDObj->Clear();
208 fTreeR = clusterTree;
209 fTreeR->Branch("FMD", &fMult);
211 AliDebug(5, "Getting entry 0 from digit branch");
212 digitBranch->GetEntry(0);
214 AliDebug(5, "Processing digits");
215 ProcessDigits(digits);
217 Int_t written = clusterTree->Fill();
218 AliDebug(10, Form("Filled %d bytes into cluster tree", written));
224 //____________________________________________________________________
226 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
228 // For each digit, find the pseudo rapdity, azimuthal angle, and
229 // number of corrected ADC counts, and pass it on to the algorithms
231 Int_t nDigits = digits->GetEntries();
232 AliDebug(1, Form("Got %d digits", nDigits));
233 for (Int_t i = 0; i < nDigits; i++) {
234 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
235 AliFMDParameters* param = AliFMDParameters::Instance();
236 // Check that the strip is not marked as dead
237 if (param->IsDead(digit->Detector(), digit->Ring(),
238 digit->Sector(), digit->Strip())) continue;
243 PhysicalCoordinates(digit, eta, phi);
245 // Substract pedestal.
246 UShort_t counts = SubtractPedestal(digit);
248 // Gain match digits.
249 Double_t edep = Adc2Energy(digit, eta, counts);
251 // Make rough multiplicity
252 Double_t mult = Energy2Multiplicity(digit, edep);
254 AliDebug(10, Form("FMD%d%c[%2d,%3d]: "
255 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
256 digit->Detector(), digit->Ring(), digit->Sector(),
257 digit->Strip(), digit->Counts(), counts, edep, mult));
259 // Create a `RecPoint' on the output branch.
261 new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(),
267 (void)m; // Suppress warnings about unused variables.
270 fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(),
271 digit->Sector(), digit->Strip(), mult);
272 fESDObj->SetEta(digit->Detector(), digit->Ring(),
273 digit->Sector(), digit->Strip(), eta);
277 //____________________________________________________________________
279 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
281 // Member function to subtract the pedestal from a digit
282 // This implementation does nothing, but a derived class could over
283 // load this to subtract a pedestal that was given in a database or
284 // something like that.
287 AliFMDParameters* param = AliFMDParameters::Instance();
288 Float_t pedM = param->GetPedestal(digit->Detector(),
292 AliDebug(10, Form("Subtracting pedestal %f from signal %d",
293 pedM, digit->Counts()));
294 if (digit->Count3() > 0) counts = digit->Count3();
295 else if (digit->Count2() > 0) counts = digit->Count2();
296 else counts = digit->Count1();
297 counts = TMath::Max(Int_t(counts - pedM), 0);
298 if (counts > 0) AliDebug(10, "Got a hit strip");
300 return UShort_t(counts);
303 //____________________________________________________________________
305 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit,
307 UShort_t count) const
309 // Converts number of ADC counts to energy deposited.
310 // Note, that this member function can be overloaded by derived
311 // classes to do strip-specific look-ups in databases or the like,
312 // to find the proper gain for a strip.
314 // In this simple version, we calculate the energy deposited as
316 // EnergyDeposited = cos(theta) * gain * count
321 // gain = ----------------- * Energy_deposited_per_MIP
324 // is constant and the same for all strips.
326 // Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
327 // Double_t edep = TMath::Abs(TMath::Cos(theta)) * fGain * count;
328 AliFMDParameters* param = AliFMDParameters::Instance();
329 Float_t gain = param->GetPulseGain(digit->Detector(),
333 Double_t edep = count * gain;
334 AliDebug(10, Form("Converting counts %d to energy via factor %f",
339 //____________________________________________________________________
341 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */,
344 // Converts an energy signal to number of particles.
345 // Note, that this member function can be overloaded by derived
346 // classes to do strip-specific look-ups in databases or the like,
347 // to find the proper gain for a strip.
349 // In this simple version, we calculate the multiplicity as
351 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
355 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
357 // is constant and the same for all strips
358 AliFMDParameters* param = AliFMDParameters::Instance();
359 Double_t edepMIP = param->GetEdepMip();
360 Float_t mult = edep / edepMIP;
362 AliDebug(10, Form("Translating energy %f to multiplicity via "
363 "divider %f->%f", edep, edepMIP, mult));
367 //____________________________________________________________________
369 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit,
373 // Get the eta and phi of a digit
376 AliFMDGeometry* geom = AliFMDGeometry::Instance();
377 Double_t x, y, z, r, theta;
378 geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(),
379 digit->Strip(), x, y, z);
380 // Correct for vertex offset.
382 phi = TMath::ATan2(y, x);
383 r = TMath::Sqrt(y * y + x * x);
384 theta = TMath::ATan2(r, z);
385 eta = -TMath::Log(TMath::Tan(theta / 2));
390 //____________________________________________________________________
392 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
393 TTree* /* clusterTree */,
396 // nothing to be done
397 // FIXME: The vertex may not be known when Reconstruct is executed,
398 // so we may have to move some of that member function here.
399 AliDebug(1, Form("Calling FillESD with two trees and one ESD"));
403 AliDebug(1, Form("Writing FMD data to ESD tree"));
404 esd->SetFMDData(fESDObj);
409 //____________________________________________________________________
411 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const
413 // Cannot be used. See member function with same name but with 2
414 // TTree arguments. Make sure you do local reconstrucion
415 AliDebug(1, Form("Calling FillESD with loader and tree"));
416 AliError("MayNotUse");
418 //____________________________________________________________________
420 AliFMDReconstructor::Reconstruct(AliRunLoader*) const
422 // Cannot be used. See member function with same name but with 2
423 // TTree arguments. Make sure you do local reconstrucion
424 AliDebug(1, Form("Calling FillESD with loader"));
425 AliError("MayNotUse");
427 //____________________________________________________________________
429 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const
431 // Cannot be used. See member function with same name but with 2
432 // TTree arguments. Make sure you do local reconstrucion
433 AliDebug(1, Form("Calling FillESD with loader and raw reader"));
434 AliError("MayNotUse");
436 //____________________________________________________________________
438 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESD*) const
440 // Cannot be used. See member function with same name but with 2
441 // TTree arguments. Make sure you do local reconstrucion
442 AliDebug(1, Form("Calling FillESD with raw reader, tree, and ESD"));
443 AliError("MayNotUse");
445 //____________________________________________________________________
447 AliFMDReconstructor::FillESD(AliRunLoader*,AliESD*) const
449 // Cannot be used. See member function with same name but with 2
450 // TTree arguments. Make sure you do local reconstrucion
451 AliDebug(1, Form("Calling FillESD with loader and ESD"));
452 AliError("MayNotUse");
454 //____________________________________________________________________
456 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESD*) const
458 // Cannot be used. See member function with same name but with 2
459 // TTree arguments. Make sure you do local reconstrucion
460 AliDebug(1, Form("Calling FillESD with loader, raw reader, and ESD"));
461 AliError("MayNotUse");
464 //____________________________________________________________________