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 "AliFMDMultStrip.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 // Current vertex position
164 // Create array of reconstructed strip multiplicities
165 fMult = new TClonesArray("AliFMDMultStrip", 51200);
166 // Create ESD output object
167 fESDObj = new AliESDFMD;
169 // Check that we have a run loader
171 Warning("Init", "No run loader");
175 // Check if we can get the header tree
176 AliHeader* header = runLoader->GetHeader();
178 Warning("Init", "No header");
182 // Check if we can get a simulation header
183 AliGenEventHeader* eventHeader = header->GenEventHeader();
186 eventHeader->PrimaryVertex(vtx);
187 fCurrentVertex = vtx[2];
188 AliDebug(1, Form("Primary vertex Z coordinate for event # %d/%d is %f",
189 header->GetRun(), header->GetEvent(), fCurrentVertex));
190 Warning("Init", "no generator event header");
193 Warning("Init", "No generator event header - "
194 "perhaps we get the vertex from ESD?");
200 //____________________________________________________________________
202 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
203 TTree* digitsTree) const
205 // Convert Raw digits to AliFMDDigit's in a tree
206 AliDebug(1, "Reading raw data into digits tree");
207 AliFMDRawReader rawRead(reader, digitsTree);
208 // rawRead.SetSampleRate(fFMD->GetSampleRate());
212 //____________________________________________________________________
214 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
215 TTree* clusterTree) const
217 // Reconstruct event from digits in tree
218 // Get the FMD branch holding the digits.
219 // FIXME: The vertex may not be known yet, so we may have to move
220 // some of this to FillESD.
221 AliDebug(1, "Reconstructing from digits in a tree");
224 const AliESDVertex* vertex = fESD->GetVertex();
226 AliDebug(1, Form("Got vertex from ESD: %f", vertex->GetZv()));
227 fCurrentVertex = vertex->GetZv();
231 TBranch *digitBranch = digitsTree->GetBranch("FMD");
233 Error("Exec", "No digit branch for the FMD found");
236 TClonesArray* digits = new TClonesArray("AliFMDDigit");
237 digitBranch->SetAddress(&digits);
239 // TIter next(&fAlgorithms);
240 // AliFMDMultAlgorithm* algorithm = 0;
241 // while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next()))) {
242 // AliDebug(10, Form("PreEvent called for algorithm %s",
243 // algorithm->GetName()));
244 // algorithm->PreEvent(clusterTree, fCurrentVertex);
246 if (fMult) fMult->Clear();
247 if (fESDObj) fESDObj->Clear();
250 fTreeR = clusterTree;
251 fTreeR->Branch("FMD", &fMult);
253 AliDebug(5, "Getting entry 0 from digit branch");
254 digitBranch->GetEntry(0);
256 AliDebug(5, "Processing digits");
257 ProcessDigits(digits);
261 // while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next()))) {
262 // AliDebug(10, Form("PostEvent called for algorithm %s",
263 // algorithm->GetName()));
264 // algorithm->PostEvent();
266 Int_t written = clusterTree->Fill();
267 AliDebug(10, Form("Filled %d bytes into cluster tree", written));
273 //____________________________________________________________________
275 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
277 // For each digit, find the pseudo rapdity, azimuthal angle, and
278 // number of corrected ADC counts, and pass it on to the algorithms
280 Int_t nDigits = digits->GetEntries();
281 AliDebug(1, Form("Got %d digits", nDigits));
282 for (Int_t i = 0; i < nDigits; i++) {
283 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
284 AliFMDParameters* param = AliFMDParameters::Instance();
285 // Check that the strip is not marked as dead
286 if (param->IsDead(digit->Detector(), digit->Ring(),
287 digit->Sector(), digit->Strip())) continue;
292 PhysicalCoordinates(digit, eta, phi);
294 // Substract pedestal.
295 UShort_t counts = SubtractPedestal(digit);
297 // Gain match digits.
298 Double_t edep = Adc2Energy(digit, eta, counts);
300 // Make rough multiplicity
301 Double_t mult = Energy2Multiplicity(digit, edep);
303 AliDebug(10, Form("FMD%d%c[%2d,%3d]: "
304 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
305 digit->Detector(), digit->Ring(), digit->Sector(),
306 digit->Strip(), digit->Counts(), counts, edep, mult));
308 // Create a `RecPoint' on the output branch.
310 new ((*fMult)[fNMult]) AliFMDMultStrip(digit->Detector(),
316 AliFMDMult::kNaiive);
317 (void)m; // Suppress warnings about unused variables.
320 fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(),
321 digit->Sector(), digit->Strip(), mult);
322 fESDObj->SetEta(digit->Detector(), digit->Ring(),
323 digit->Sector(), digit->Strip(), eta);
327 //____________________________________________________________________
329 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
331 // Member function to subtract the pedestal from a digit
332 // This implementation does nothing, but a derived class could over
333 // load this to subtract a pedestal that was given in a database or
334 // something like that.
337 AliFMDParameters* param = AliFMDParameters::Instance();
338 Float_t pedM = param->GetPedestal(digit->Detector(),
342 AliDebug(10, Form("Subtracting pedestal %f from signal %d",
343 pedM, digit->Counts()));
344 if (digit->Count3() > 0) counts = digit->Count3();
345 else if (digit->Count2() > 0) counts = digit->Count2();
346 else counts = digit->Count1();
347 counts = TMath::Max(Int_t(counts - pedM), 0);
348 if (counts > 0) AliDebug(10, "Got a hit strip");
350 return UShort_t(counts);
353 //____________________________________________________________________
355 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit,
357 UShort_t count) const
359 // Converts number of ADC counts to energy deposited.
360 // Note, that this member function can be overloaded by derived
361 // classes to do strip-specific look-ups in databases or the like,
362 // to find the proper gain for a strip.
364 // In this simple version, we calculate the energy deposited as
366 // EnergyDeposited = cos(theta) * gain * count
371 // gain = ----------------- * Energy_deposited_per_MIP
374 // is constant and the same for all strips.
376 // Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
377 // Double_t edep = TMath::Abs(TMath::Cos(theta)) * fGain * count;
378 AliFMDParameters* param = AliFMDParameters::Instance();
379 Float_t gain = param->GetPulseGain(digit->Detector(),
383 Double_t edep = count * gain;
384 AliDebug(10, Form("Converting counts %d to energy via factor %f",
389 //____________________________________________________________________
391 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */,
394 // Converts an energy signal to number of particles.
395 // Note, that this member function can be overloaded by derived
396 // classes to do strip-specific look-ups in databases or the like,
397 // to find the proper gain for a strip.
399 // In this simple version, we calculate the multiplicity as
401 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
405 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
407 // is constant and the same for all strips
408 AliFMDParameters* param = AliFMDParameters::Instance();
409 Double_t edepMIP = param->GetEdepMip();
410 Float_t mult = edep / edepMIP;
412 AliDebug(10, Form("Translating energy %f to multiplicity via "
413 "divider %f->%f", edep, edepMIP, mult));
417 //____________________________________________________________________
419 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit,
423 // Get the eta and phi of a digit
426 AliFMDGeometry* fmd = AliFMDGeometry::Instance();
427 AliFMDDetector* subDetector = fmd->GetDetector(digit->Detector());
429 Warning("ProcessDigits", "Unknown detector: FMD%d" , digit->Detector());
433 // Get the ring - we should really use geometry->Detector2XYZ(...)
435 AliFMDRing* ring = subDetector->GetRing(digit->Ring());
436 Float_t ringZ = subDetector->GetRingZ(digit->Ring());
438 Warning("ProcessDigits", "Unknown ring: FMD%d%c", digit->Detector(),
443 // Correct for vertex offset.
444 Float_t realZ = fCurrentVertex + ringZ;
445 Float_t stripR = ((ring->GetHighR() - ring->GetLowR())
446 / ring->GetNStrips() * (digit->Strip() + .5)
448 Float_t theta = TMath::ATan2(stripR, realZ);
449 phi = (2 * TMath::Pi() / ring->GetNSectors()
450 * (digit->Sector() + .5));
451 eta = -TMath::Log(TMath::Tan(theta / 2));
456 //____________________________________________________________________
458 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
459 TTree* /* clusterTree */,
462 // nothing to be done
463 // FIXME: The vertex may not be known when Reconstruct is executed,
464 // so we may have to move some of that member function here.
465 AliDebug(1, Form("Calling FillESD with two trees and one ESD"));
469 AliDebug(1, Form("Writing FMD data to ESD tree"));
470 esd->SetFMDData(fESDObj);
471 // Let's check the data in the ESD
473 AliESDFMD* fromEsd = esd->GetFMDData();
475 AliWarning("No FMD object in ESD!");
478 for (UShort_t det = 1; det <= fESDObj->MaxDetectors(); det++) {
479 for (UShort_t ir = 0; ir < fESDObj->MaxRings(); ir++) {
480 Char_t ring = (ir == 0 ? 'I' : 'O');
481 for (UShort_t sec = 0; sec < fESDObj->MaxSectors(); sec++) {
482 for (UShort_t str = 0; str < fESDObj->MaxStrips(); str++) {
483 if (fESDObj->Multiplicity(det, ring, sec, str) !=
484 fromEsd->Multiplicity(det, ring, sec, str))
485 AliWarning(Form("Mult for FMD%d%c[%2d,%3d]",det,ring,sec,str));
486 if (fESDObj->Eta(det, ring, sec, str) !=
487 fromEsd->Eta(det, ring, sec, str))
488 AliWarning(Form("Eta for FMD%d%c[%2d,%3d]", det,ring,sec,str));
489 if (fESDObj->Multiplicity(det, ring, sec, str) > 0 &&
490 fESDObj->Multiplicity(det, ring, sec, str)
491 != AliESDFMD::kInvalidMult)
492 AliInfo(Form("Mult in FMD%d%c[%2d,%3d] is %f", det,ring,sec,str,
493 fESDObj->Multiplicity(det, ring, sec, str)));
502 static Int_t evNo = -1;
504 if (esd) evNo = esd->GetEventNumber();
505 TString fname(Form("FMD.ESD.%03d.root", evNo));
506 TFile* file = TFile::Open(fname.Data(), "RECREATE");
508 AliError(Form("Failed to open file %s", fname.Data()));
516 TClonesArray* multStrips = 0;
517 TClonesArray* multRegions = 0;
518 TTree* treeR = fmdLoader->TreeR();
519 TBranch* branchRegions = treeR->GetBranch("FMDPoisson");
520 TBranch* branchStrips = treeR->GetBranch("FMDNaiive");
521 branchRegions->SetAddress(&multRegions);
522 branchStrips->SetAddress(&multStrips);
525 Int_t nEntries = clusterTree->GetEntries();
526 for (Int_t entry = 0; entry < nEntries; entry++) {
527 AliDebug(5, Form("Entry # %d in cluster tree", entry));
528 treeR->GetEntry(entry);
531 Int_t nMults = multRegions->GetLast();
532 for (Int_t i = 0; i <= nMults; i++) {
533 AliFMDMultRegion* multR =
534 static_cast<AliFMDMultRegion*>(multRegions->UncheckedAt(i));
535 Int_t nParticles=multR->Particles();
536 if (i>=0 && i<=13) hEtaPoissonI1->AddBinContent(i+1,nParticles);
537 if (i>=14 && i<=27 ) hEtaPoissonI2->AddBinContent(i-13,nParticles);
538 if (i>=28 && i<=33 );
539 if (i>=34 && i<=47 ) hEtaPoissonI3->AddBinContent(48-i,nParticles);
540 if (i>=48 && i<=53) hEtaPoissonO3->AddBinContent(54-i,nParticles);
547 //____________________________________________________________________
549 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const
551 // Cannot be used. See member function with same name but with 2
552 // TTree arguments. Make sure you do local reconstrucion
553 AliDebug(1, Form("Calling FillESD with loader and tree"));
554 AliError("MayNotUse");
556 //____________________________________________________________________
558 AliFMDReconstructor::Reconstruct(AliRunLoader*) const
560 // Cannot be used. See member function with same name but with 2
561 // TTree arguments. Make sure you do local reconstrucion
562 AliDebug(1, Form("Calling FillESD with loader"));
563 AliError("MayNotUse");
565 //____________________________________________________________________
567 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const
569 // Cannot be used. See member function with same name but with 2
570 // TTree arguments. Make sure you do local reconstrucion
571 AliDebug(1, Form("Calling FillESD with loader and raw reader"));
572 AliError("MayNotUse");
574 //____________________________________________________________________
576 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESD*) const
578 // Cannot be used. See member function with same name but with 2
579 // TTree arguments. Make sure you do local reconstrucion
580 AliDebug(1, Form("Calling FillESD with raw reader, tree, and ESD"));
581 AliError("MayNotUse");
583 //____________________________________________________________________
585 AliFMDReconstructor::FillESD(AliRunLoader*,AliESD*) const
587 // Cannot be used. See member function with same name but with 2
588 // TTree arguments. Make sure you do local reconstrucion
589 AliDebug(1, Form("Calling FillESD with loader and ESD"));
590 AliError("MayNotUse");
592 //____________________________________________________________________
594 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESD*) const
596 // Cannot be used. See member function with same name but with 2
597 // TTree arguments. Make sure you do local reconstrucion
598 AliDebug(1, Form("Calling FillESD with loader, raw reader, and ESD"));
599 AliError("MayNotUse");
602 //____________________________________________________________________