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 **************************************************************************/
18 //____________________________________________________________________
20 // This is a class that constructs AliFMDMult (reconstructed
21 // multiplicity) from of Digits
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).
27 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
28 // Latest changes by Christian Holm Christensen <cholm@nbi.dk>
31 //____________________________________________________________________
33 #include <AliLog.h> // ALILOG_H
34 #include <AliRun.h> // ALIRUN_H
35 #include <AliRunLoader.h> // ALIRUNLOADER_H
36 #include <AliLoader.h> // ALILOADER_H
37 #include <AliHeader.h> // ALIHEADER_H
38 #include <AliRawReader.h> // ALIRAWREADER_H
39 #include <AliGenEventHeader.h> // ALIGENEVENTHEADER_H
40 #include "AliFMD.h" // ALIFMD_H
41 #include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
42 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
43 #include "AliFMDDetector.h" // ALIFMDDETECTOR_H
44 #include "AliFMDRing.h" // ALIFMDRING_H
45 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
46 #include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
47 #include "AliFMDRawStream.h" // ALIFMDRAWSTREAM_H
48 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
49 #include "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
50 #include "AliESD.h" // ALIESD_H
51 #include <AliESDFMD.h> // ALIESDFMD_H
54 //____________________________________________________________________
55 ClassImp(AliFMDReconstructor)
57 ; // This is here to keep Emacs for indenting the next line
60 //____________________________________________________________________
61 AliFMDReconstructor::AliFMDReconstructor()
70 // Make a new FMD reconstructor object - default CTOR.
74 //____________________________________________________________________
75 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other)
80 fCurrentVertex(other.fCurrentVertex),
81 fESDObj(other.fESDObj),
88 //____________________________________________________________________
90 AliFMDReconstructor::operator=(const AliFMDReconstructor& other)
92 // Assignment operator
94 fNMult = other.fNMult;
95 fTreeR = other.fTreeR;
96 fCurrentVertex = other.fCurrentVertex;
97 fESDObj = other.fESDObj;
102 //____________________________________________________________________
103 AliFMDReconstructor::~AliFMDReconstructor()
106 if (fMult) fMult->Delete();
107 if (fMult) delete fMult;
108 if (fESDObj) delete fESDObj;
111 //____________________________________________________________________
113 AliFMDReconstructor::Init(AliRunLoader* runLoader)
115 // Initialize the reconstructor
116 AliDebug(1, Form("Init called with runloader 0x%x", runLoader));
117 // Initialize the geometry
118 AliFMDGeometry* geom = AliFMDGeometry::Instance();
120 geom->InitTransformations();
122 // Initialize the parameters
123 AliFMDParameters* param = AliFMDParameters::Instance();
126 // Current vertex position
128 // Create array of reconstructed strip multiplicities
129 fMult = new TClonesArray("AliFMDRecPoint", 51200);
130 // Create ESD output object
131 fESDObj = new AliESDFMD;
133 // Check that we have a run loader
135 Warning("Init", "No run loader");
139 // Check if we can get the header tree
140 AliHeader* header = runLoader->GetHeader();
142 Warning("Init", "No header");
146 // Check if we can get a simulation header
147 AliGenEventHeader* eventHeader = header->GenEventHeader();
150 eventHeader->PrimaryVertex(vtx);
151 fCurrentVertex = vtx[2];
152 AliDebug(1, Form("Primary vertex Z coordinate for event # %d/%d is %f",
153 header->GetRun(), header->GetEvent(), fCurrentVertex));
154 Warning("Init", "no generator event header");
157 Warning("Init", "No generator event header - "
158 "perhaps we get the vertex from ESD?");
162 //____________________________________________________________________
164 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
165 TTree* digitsTree) const
167 // Convert Raw digits to AliFMDDigit's in a tree
168 AliDebug(1, "Reading raw data into digits tree");
169 AliFMDRawReader rawRead(reader, digitsTree);
170 // rawRead.SetSampleRate(fFMD->GetSampleRate());
174 //____________________________________________________________________
176 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
177 TTree* clusterTree) const
179 // Reconstruct event from digits in tree
180 // Get the FMD branch holding the digits.
181 // FIXME: The vertex may not be known yet, so we may have to move
182 // some of this to FillESD.
183 AliDebug(1, "Reconstructing from digits in a tree");
186 const AliESDVertex* vertex = fESD->GetVertex();
188 AliDebug(1, Form("Got vertex from ESD: %f", vertex->GetZv()));
189 fCurrentVertex = vertex->GetZv();
193 TBranch *digitBranch = digitsTree->GetBranch("FMD");
195 Error("Exec", "No digit branch for the FMD found");
198 TClonesArray* digits = new TClonesArray("AliFMDDigit");
199 digitBranch->SetAddress(&digits);
201 if (fMult) fMult->Clear();
202 if (fESDObj) fESDObj->Clear();
205 fTreeR = clusterTree;
206 fTreeR->Branch("FMD", &fMult);
208 AliDebug(5, "Getting entry 0 from digit branch");
209 digitBranch->GetEntry(0);
211 AliDebug(5, "Processing digits");
212 ProcessDigits(digits);
214 Int_t written = clusterTree->Fill();
215 AliDebug(10, Form("Filled %d bytes into cluster tree", written));
221 //____________________________________________________________________
223 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
225 // For each digit, find the pseudo rapdity, azimuthal angle, and
226 // number of corrected ADC counts, and pass it on to the algorithms
228 Int_t nDigits = digits->GetEntries();
229 AliDebug(1, Form("Got %d digits", nDigits));
230 for (Int_t i = 0; i < nDigits; i++) {
231 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
232 AliFMDParameters* param = AliFMDParameters::Instance();
233 // Check that the strip is not marked as dead
234 if (param->IsDead(digit->Detector(), digit->Ring(),
235 digit->Sector(), digit->Strip())) continue;
240 PhysicalCoordinates(digit, eta, phi);
242 // Substract pedestal.
243 UShort_t counts = SubtractPedestal(digit);
245 // Gain match digits.
246 Double_t edep = Adc2Energy(digit, eta, counts);
248 // Make rough multiplicity
249 Double_t mult = Energy2Multiplicity(digit, edep);
251 AliDebug(10, Form("FMD%d%c[%2d,%3d]: "
252 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
253 digit->Detector(), digit->Ring(), digit->Sector(),
254 digit->Strip(), digit->Counts(), counts, edep, mult));
256 // Create a `RecPoint' on the output branch.
258 new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(),
264 (void)m; // Suppress warnings about unused variables.
267 fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(),
268 digit->Sector(), digit->Strip(), mult);
269 fESDObj->SetEta(digit->Detector(), digit->Ring(),
270 digit->Sector(), digit->Strip(), eta);
274 //____________________________________________________________________
276 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
278 // Member function to subtract the pedestal from a digit
279 // This implementation does nothing, but a derived class could over
280 // load this to subtract a pedestal that was given in a database or
281 // something like that.
284 AliFMDParameters* param = AliFMDParameters::Instance();
285 Float_t pedM = param->GetPedestal(digit->Detector(),
289 AliDebug(10, Form("Subtracting pedestal %f from signal %d",
290 pedM, digit->Counts()));
291 if (digit->Count3() > 0) counts = digit->Count3();
292 else if (digit->Count2() > 0) counts = digit->Count2();
293 else counts = digit->Count1();
294 counts = TMath::Max(Int_t(counts - pedM), 0);
295 if (counts > 0) AliDebug(10, "Got a hit strip");
297 return UShort_t(counts);
300 //____________________________________________________________________
302 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit,
304 UShort_t count) const
306 // Converts number of ADC counts to energy deposited.
307 // Note, that this member function can be overloaded by derived
308 // classes to do strip-specific look-ups in databases or the like,
309 // to find the proper gain for a strip.
311 // In this simple version, we calculate the energy deposited as
313 // EnergyDeposited = cos(theta) * gain * count
318 // gain = ----------------- * Energy_deposited_per_MIP
321 // is constant and the same for all strips.
323 // Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
324 // Double_t edep = TMath::Abs(TMath::Cos(theta)) * fGain * count;
325 AliFMDParameters* param = AliFMDParameters::Instance();
326 Float_t gain = param->GetPulseGain(digit->Detector(),
330 Double_t edep = count * gain;
331 AliDebug(10, Form("Converting counts %d to energy via factor %f",
336 //____________________________________________________________________
338 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */,
341 // Converts an energy signal to number of particles.
342 // Note, that this member function can be overloaded by derived
343 // classes to do strip-specific look-ups in databases or the like,
344 // to find the proper gain for a strip.
346 // In this simple version, we calculate the multiplicity as
348 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
352 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
354 // is constant and the same for all strips
355 AliFMDParameters* param = AliFMDParameters::Instance();
356 Double_t edepMIP = param->GetEdepMip();
357 Float_t mult = edep / edepMIP;
359 AliDebug(10, Form("Translating energy %f to multiplicity via "
360 "divider %f->%f", edep, edepMIP, mult));
364 //____________________________________________________________________
366 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit,
370 // Get the eta and phi of a digit
373 AliFMDGeometry* geom = AliFMDGeometry::Instance();
374 Double_t x, y, z, r, theta;
375 geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(),
376 digit->Strip(), x, y, z);
377 // Correct for vertex offset.
379 phi = TMath::ATan2(y, x);
380 r = TMath::Sqrt(y * y + x * x);
381 theta = TMath::ATan2(r, z);
382 eta = -TMath::Log(TMath::Tan(theta / 2));
387 //____________________________________________________________________
389 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
390 TTree* /* clusterTree */,
393 // nothing to be done
394 // FIXME: The vertex may not be known when Reconstruct is executed,
395 // so we may have to move some of that member function here.
396 AliDebug(1, Form("Calling FillESD with two trees and one ESD"));
400 AliDebug(1, Form("Writing FMD data to ESD tree"));
401 esd->SetFMDData(fESDObj);
406 //____________________________________________________________________
408 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const
410 // Cannot be used. See member function with same name but with 2
411 // TTree arguments. Make sure you do local reconstrucion
412 AliDebug(1, Form("Calling FillESD with loader and tree"));
413 AliError("MayNotUse");
415 //____________________________________________________________________
417 AliFMDReconstructor::Reconstruct(AliRunLoader*) const
419 // Cannot be used. See member function with same name but with 2
420 // TTree arguments. Make sure you do local reconstrucion
421 AliDebug(1, Form("Calling FillESD with loader"));
422 AliError("MayNotUse");
424 //____________________________________________________________________
426 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const
428 // Cannot be used. See member function with same name but with 2
429 // TTree arguments. Make sure you do local reconstrucion
430 AliDebug(1, Form("Calling FillESD with loader and raw reader"));
431 AliError("MayNotUse");
433 //____________________________________________________________________
435 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESD*) const
437 // Cannot be used. See member function with same name but with 2
438 // TTree arguments. Make sure you do local reconstrucion
439 AliDebug(1, Form("Calling FillESD with raw reader, tree, and ESD"));
440 AliError("MayNotUse");
442 //____________________________________________________________________
444 AliFMDReconstructor::FillESD(AliRunLoader*,AliESD*) const
446 // Cannot be used. See member function with same name but with 2
447 // TTree arguments. Make sure you do local reconstrucion
448 AliDebug(1, Form("Calling FillESD with loader and ESD"));
449 AliError("MayNotUse");
451 //____________________________________________________________________
453 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESD*) const
455 // Cannot be used. See member function with same name but with 2
456 // TTree arguments. Make sure you do local reconstrucion
457 AliDebug(1, Form("Calling FillESD with loader, raw reader, and ESD"));
458 AliError("MayNotUse");
461 //____________________________________________________________________