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 // From the cached values it then calculates the number of particles
28 // that hit a region of the FMDs, as specified by the user.
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.
35 // +---------------------+ +---------------------+
36 // | AliFMDReconstructor |<>-----| AliFMDMultAlgorithm |
37 // +---------------------+ +---------------------+
40 // +-----------+---------+
42 // +-------------------+ +------------------+
43 // | AliFMDMultPoisson | | AliFMDMultNaiive |
44 // +-------------------+ +------------------+
46 // AliFMDReconstructor acts as a manager class. It contains a list of
47 // AliFMDMultAlgorithm objects. The call graph looks something like
50 // +----------------------+ +----------------------+
51 // | :AliFMDReconstructor | | :AliFMDMultAlgorithm |
52 // +----------------------+ +----------------------+
55 // ------------>| | PreRun +-+
56 // | |------------------------------->| |
58 // | |-----+ (for each event) |
59 // | | | *ProcessEvent |
61 // || |<---+ PreEvent +-+
62 // || |------------------------------>| |
65 // || | | ProcessDigits |
68 // ||| | *ProcessDigit(digit) +-+
69 // ||| |----------------------------->| |
73 // || |------------------------------>| |
77 // | |------------------------------->| |
84 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
85 // Latest changes by Christian Holm Christensen <cholm@nbi.dk>
88 //____________________________________________________________________
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
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
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 "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
107 #include "AliESD.h" // ALIESD_H
108 #include <AliESDFMD.h> // ALIESDFMD_H
111 //____________________________________________________________________
112 ClassImp(AliFMDReconstructor)
114 ; // This is here to keep Emacs for indenting the next line
117 //____________________________________________________________________
118 AliFMDReconstructor::AliFMDReconstructor()
119 : AliReconstructor(),
123 // Make a new FMD reconstructor object - default CTOR.
127 //____________________________________________________________________
128 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other)
129 : AliReconstructor(),
131 fESDObj(other.fESDObj)
137 //____________________________________________________________________
139 AliFMDReconstructor::operator=(const AliFMDReconstructor& other)
141 // Assignment operator
143 fESDObj = other.fESDObj;
147 //____________________________________________________________________
148 AliFMDReconstructor::~AliFMDReconstructor()
151 if (fMult) fMult->Delete();
152 if (fMult) delete fMult;
153 if (fESDObj) delete fESDObj;
156 //____________________________________________________________________
158 AliFMDReconstructor::Init(AliRunLoader* runLoader)
160 // Initialize the reconstructor
161 AliDebug(1, Form("Init called with runloader 0x%x", runLoader));
162 // Initialize the geometry
163 AliFMDGeometry* fmd = AliFMDGeometry::Instance();
165 fmd->InitTransformations();
167 // Current vertex position
169 // Create array of reconstructed strip multiplicities
170 fMult = new TClonesArray("AliFMDRecPoint", 51200);
171 // Create ESD output object
172 fESDObj = new AliESDFMD;
174 // Check that we have a run loader
176 Warning("Init", "No run loader");
180 // Check if we can get the header tree
181 AliHeader* header = runLoader->GetHeader();
183 Warning("Init", "No header");
187 // Check if we can get a simulation header
188 AliGenEventHeader* eventHeader = header->GenEventHeader();
191 eventHeader->PrimaryVertex(vtx);
192 fCurrentVertex = vtx[2];
193 AliDebug(1, Form("Primary vertex Z coordinate for event # %d/%d is %f",
194 header->GetRun(), header->GetEvent(), fCurrentVertex));
195 Warning("Init", "no generator event header");
198 Warning("Init", "No generator event header - "
199 "perhaps we get the vertex from ESD?");
203 //____________________________________________________________________
205 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
206 TTree* digitsTree) const
208 // Convert Raw digits to AliFMDDigit's in a tree
209 AliDebug(1, "Reading raw data into digits tree");
210 AliFMDRawReader rawRead(reader, digitsTree);
211 // rawRead.SetSampleRate(fFMD->GetSampleRate());
215 //____________________________________________________________________
217 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
218 TTree* clusterTree) const
220 // Reconstruct event from digits in tree
221 // Get the FMD branch holding the digits.
222 // FIXME: The vertex may not be known yet, so we may have to move
223 // some of this to FillESD.
224 AliDebug(1, "Reconstructing from digits in a tree");
227 const AliESDVertex* vertex = fESD->GetVertex();
229 AliDebug(1, Form("Got vertex from ESD: %f", vertex->GetZv()));
230 fCurrentVertex = vertex->GetZv();
234 TBranch *digitBranch = digitsTree->GetBranch("FMD");
236 Error("Exec", "No digit branch for the FMD found");
239 TClonesArray* digits = new TClonesArray("AliFMDDigit");
240 digitBranch->SetAddress(&digits);
242 // TIter next(&fAlgorithms);
243 // AliFMDMultAlgorithm* algorithm = 0;
244 // while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next()))) {
245 // AliDebug(10, Form("PreEvent called for algorithm %s",
246 // algorithm->GetName()));
247 // algorithm->PreEvent(clusterTree, fCurrentVertex);
249 if (fMult) fMult->Clear();
250 if (fESDObj) fESDObj->Clear();
253 fTreeR = clusterTree;
254 fTreeR->Branch("FMD", &fMult);
256 AliDebug(5, "Getting entry 0 from digit branch");
257 digitBranch->GetEntry(0);
259 AliDebug(5, "Processing digits");
260 ProcessDigits(digits);
264 // while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next()))) {
265 // AliDebug(10, Form("PostEvent called for algorithm %s",
266 // algorithm->GetName()));
267 // algorithm->PostEvent();
269 Int_t written = clusterTree->Fill();
270 AliDebug(10, Form("Filled %d bytes into cluster tree", written));
276 //____________________________________________________________________
278 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
280 // For each digit, find the pseudo rapdity, azimuthal angle, and
281 // number of corrected ADC counts, and pass it on to the algorithms
283 Int_t nDigits = digits->GetEntries();
284 AliDebug(1, Form("Got %d digits", nDigits));
285 for (Int_t i = 0; i < nDigits; i++) {
286 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
287 AliFMDParameters* param = AliFMDParameters::Instance();
288 // Check that the strip is not marked as dead
289 if (param->IsDead(digit->Detector(), digit->Ring(),
290 digit->Sector(), digit->Strip())) continue;
295 PhysicalCoordinates(digit, eta, phi);
297 // Substract pedestal.
298 UShort_t counts = SubtractPedestal(digit);
300 // Gain match digits.
301 Double_t edep = Adc2Energy(digit, eta, counts);
303 // Make rough multiplicity
304 Double_t mult = Energy2Multiplicity(digit, edep);
306 AliDebug(10, Form("FMD%d%c[%2d,%3d]: "
307 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
308 digit->Detector(), digit->Ring(), digit->Sector(),
309 digit->Strip(), digit->Counts(), counts, edep, mult));
311 // Create a `RecPoint' on the output branch.
313 new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(),
319 (void)m; // Suppress warnings about unused variables.
322 fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(),
323 digit->Sector(), digit->Strip(), mult);
324 fESDObj->SetEta(digit->Detector(), digit->Ring(),
325 digit->Sector(), digit->Strip(), eta);
329 //____________________________________________________________________
331 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
333 // Member function to subtract the pedestal from a digit
334 // This implementation does nothing, but a derived class could over
335 // load this to subtract a pedestal that was given in a database or
336 // something like that.
339 AliFMDParameters* param = AliFMDParameters::Instance();
340 Float_t pedM = param->GetPedestal(digit->Detector(),
344 AliDebug(10, Form("Subtracting pedestal %f from signal %d",
345 pedM, digit->Counts()));
346 if (digit->Count3() > 0) counts = digit->Count3();
347 else if (digit->Count2() > 0) counts = digit->Count2();
348 else counts = digit->Count1();
349 counts = TMath::Max(Int_t(counts - pedM), 0);
350 if (counts > 0) AliDebug(10, "Got a hit strip");
352 return UShort_t(counts);
355 //____________________________________________________________________
357 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit,
359 UShort_t count) const
361 // Converts number of ADC counts to energy deposited.
362 // Note, that this member function can be overloaded by derived
363 // classes to do strip-specific look-ups in databases or the like,
364 // to find the proper gain for a strip.
366 // In this simple version, we calculate the energy deposited as
368 // EnergyDeposited = cos(theta) * gain * count
373 // gain = ----------------- * Energy_deposited_per_MIP
376 // is constant and the same for all strips.
378 // Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
379 // Double_t edep = TMath::Abs(TMath::Cos(theta)) * fGain * count;
380 AliFMDParameters* param = AliFMDParameters::Instance();
381 Float_t gain = param->GetPulseGain(digit->Detector(),
385 Double_t edep = count * gain;
386 AliDebug(10, Form("Converting counts %d to energy via factor %f",
391 //____________________________________________________________________
393 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */,
396 // Converts an energy signal to number of particles.
397 // Note, that this member function can be overloaded by derived
398 // classes to do strip-specific look-ups in databases or the like,
399 // to find the proper gain for a strip.
401 // In this simple version, we calculate the multiplicity as
403 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
407 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
409 // is constant and the same for all strips
410 AliFMDParameters* param = AliFMDParameters::Instance();
411 Double_t edepMIP = param->GetEdepMip();
412 Float_t mult = edep / edepMIP;
414 AliDebug(10, Form("Translating energy %f to multiplicity via "
415 "divider %f->%f", edep, edepMIP, mult));
419 //____________________________________________________________________
421 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit,
425 // Get the eta and phi of a digit
428 AliFMDGeometry* fmd = AliFMDGeometry::Instance();
430 AliFMDDetector* subDetector = fmd->GetDetector(digit->Detector());
432 Warning("ProcessDigits", "Unknown detector: FMD%d" , digit->Detector());
436 // Get the ring - we should really use geometry->Detector2XYZ(...)
438 AliFMDRing* ring = subDetector->GetRing(digit->Ring());
439 Float_t ringZ = subDetector->GetRingZ(digit->Ring());
441 Warning("ProcessDigits", "Unknown ring: FMD%d%c", digit->Detector(),
445 Float_t realZ = fCurrentVertex + ringZ;
446 Float_t stripR = ((ring->GetHighR() - ring->GetLowR())
447 / ring->GetNStrips() * (digit->Strip() + .5)
449 Float_t theta = TMath::ATan2(stripR, realZ);
451 Double_t x, y, z, r, theta;
452 fmd->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(),
453 digit->Strip(), x, y, z);
454 // Correct for vertex offset.
456 phi = TMath::ATan2(y, x);
457 r = TMath::Sqrt(y * y + x * x);
458 theta = TMath::ATan2(r, z);
459 eta = -TMath::Log(TMath::Tan(theta / 2));
464 //____________________________________________________________________
466 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
467 TTree* /* clusterTree */,
470 // nothing to be done
471 // FIXME: The vertex may not be known when Reconstruct is executed,
472 // so we may have to move some of that member function here.
473 AliDebug(1, Form("Calling FillESD with two trees and one ESD"));
477 AliDebug(1, Form("Writing FMD data to ESD tree"));
478 esd->SetFMDData(fESDObj);
479 // Let's check the data in the ESD
481 AliESDFMD* fromEsd = esd->GetFMDData();
483 AliWarning("No FMD object in ESD!");
486 for (UShort_t det = 1; det <= fESDObj->MaxDetectors(); det++) {
487 for (UShort_t ir = 0; ir < fESDObj->MaxRings(); ir++) {
488 Char_t ring = (ir == 0 ? 'I' : 'O');
489 for (UShort_t sec = 0; sec < fESDObj->MaxSectors(); sec++) {
490 for (UShort_t str = 0; str < fESDObj->MaxStrips(); str++) {
491 if (fESDObj->Multiplicity(det, ring, sec, str) !=
492 fromEsd->Multiplicity(det, ring, sec, str))
493 AliWarning(Form("Mult for FMD%d%c[%2d,%3d]",det,ring,sec,str));
494 if (fESDObj->Eta(det, ring, sec, str) !=
495 fromEsd->Eta(det, ring, sec, str))
496 AliWarning(Form("Eta for FMD%d%c[%2d,%3d]", det,ring,sec,str));
497 if (fESDObj->Multiplicity(det, ring, sec, str) > 0 &&
498 fESDObj->Multiplicity(det, ring, sec, str)
499 != AliESDFMD::kInvalidMult)
500 AliInfo(Form("Mult in FMD%d%c[%2d,%3d] is %f", det,ring,sec,str,
501 fESDObj->Multiplicity(det, ring, sec, str)));
510 static Int_t evNo = -1;
512 if (esd) evNo = esd->GetEventNumber();
513 TString fname(Form("FMD.ESD.%03d.root", evNo));
514 TFile* file = TFile::Open(fname.Data(), "RECREATE");
516 AliError(Form("Failed to open file %s", fname.Data()));
524 TClonesArray* multStrips = 0;
525 TClonesArray* multRegions = 0;
526 TTree* treeR = fmdLoader->TreeR();
527 TBranch* branchRegions = treeR->GetBranch("FMDPoisson");
528 TBranch* branchStrips = treeR->GetBranch("FMDNaiive");
529 branchRegions->SetAddress(&multRegions);
530 branchStrips->SetAddress(&multStrips);
533 Int_t nEntries = clusterTree->GetEntries();
534 for (Int_t entry = 0; entry < nEntries; entry++) {
535 AliDebug(5, Form("Entry # %d in cluster tree", entry));
536 treeR->GetEntry(entry);
539 Int_t nMults = multRegions->GetLast();
540 for (Int_t i = 0; i <= nMults; i++) {
541 AliFMDMultRegion* multR =
542 static_cast<AliFMDMultRegion*>(multRegions->UncheckedAt(i));
543 Int_t nParticles=multR->Particles();
544 if (i>=0 && i<=13) hEtaPoissonI1->AddBinContent(i+1,nParticles);
545 if (i>=14 && i<=27 ) hEtaPoissonI2->AddBinContent(i-13,nParticles);
546 if (i>=28 && i<=33 );
547 if (i>=34 && i<=47 ) hEtaPoissonI3->AddBinContent(48-i,nParticles);
548 if (i>=48 && i<=53) hEtaPoissonO3->AddBinContent(54-i,nParticles);
555 //____________________________________________________________________
557 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const
559 // Cannot be used. See member function with same name but with 2
560 // TTree arguments. Make sure you do local reconstrucion
561 AliDebug(1, Form("Calling FillESD with loader and tree"));
562 AliError("MayNotUse");
564 //____________________________________________________________________
566 AliFMDReconstructor::Reconstruct(AliRunLoader*) const
568 // Cannot be used. See member function with same name but with 2
569 // TTree arguments. Make sure you do local reconstrucion
570 AliDebug(1, Form("Calling FillESD with loader"));
571 AliError("MayNotUse");
573 //____________________________________________________________________
575 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const
577 // Cannot be used. See member function with same name but with 2
578 // TTree arguments. Make sure you do local reconstrucion
579 AliDebug(1, Form("Calling FillESD with loader and raw reader"));
580 AliError("MayNotUse");
582 //____________________________________________________________________
584 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESD*) const
586 // Cannot be used. See member function with same name but with 2
587 // TTree arguments. Make sure you do local reconstrucion
588 AliDebug(1, Form("Calling FillESD with raw reader, tree, and ESD"));
589 AliError("MayNotUse");
591 //____________________________________________________________________
593 AliFMDReconstructor::FillESD(AliRunLoader*,AliESD*) const
595 // Cannot be used. See member function with same name but with 2
596 // TTree arguments. Make sure you do local reconstrucion
597 AliDebug(1, Form("Calling FillESD with loader and ESD"));
598 AliError("MayNotUse");
600 //____________________________________________________________________
602 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESD*) const
604 // Cannot be used. See member function with same name but with 2
605 // TTree arguments. Make sure you do local reconstrucion
606 AliDebug(1, Form("Calling FillESD with loader, raw reader, and ESD"));
607 AliError("MayNotUse");
610 //____________________________________________________________________