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 AliFMDRecPoint objects from of Digits
24 // This class reads either digits from a TClonesArray or raw data from
25 // a DDL file (or similar), and stores the read ADC counts in an
26 // internal cache (fAdcs). The rec-points are made via the naiive
29 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
30 // Latest changes by Christian Holm Christensen <cholm@nbi.dk>
33 //____________________________________________________________________
35 // #include <AliLog.h> // ALILOG_H
36 // #include <AliRun.h> // ALIRUN_H
37 #include "AliFMDDebug.h"
38 // to be removed as soon as we remove it from the base class
39 #include "AliRunLoader.h"
40 #include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
41 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
42 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
43 #include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
44 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
45 #include "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
46 #include "AliESDEvent.h" // ALIESDEVENT_H
47 #include "AliESDVertex.h" // ALIESDVERTEX_H
48 #include <AliESDFMD.h> // ALIESDFMD_H
55 //____________________________________________________________________
56 ClassImp(AliFMDReconstructor)
58 ; // This is here to keep Emacs for indenting the next line
61 //____________________________________________________________________
62 AliFMDReconstructor::AliFMDReconstructor()
71 fVertexType(kNoVertex),
80 // Make a new FMD reconstructor object - default CTOR.
86 //____________________________________________________________________
87 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other)
92 fCurrentVertex(other.fCurrentVertex),
93 fESDObj(other.fESDObj),
94 fNoiseFactor(other.fNoiseFactor),
95 fAngleCorrect(other.fAngleCorrect),
96 fVertexType(other.fVertexType),
98 fDiagnostics(other.fDiagnostics),
99 fDiagStep1(other.fDiagStep1),
100 fDiagStep2(other.fDiagStep2),
101 fDiagStep3(other.fDiagStep3),
102 fDiagStep4(other.fDiagStep4),
103 fDiagAll(other.fDiagAll)
109 //____________________________________________________________________
111 AliFMDReconstructor::operator=(const AliFMDReconstructor& other)
113 // Assignment operator
115 fNMult = other.fNMult;
116 fTreeR = other.fTreeR;
117 fCurrentVertex = other.fCurrentVertex;
118 fESDObj = other.fESDObj;
119 fNoiseFactor = other.fNoiseFactor;
120 fAngleCorrect = other.fAngleCorrect;
121 fVertexType = other.fVertexType;
123 fDiagnostics = other.fDiagnostics;
124 fDiagStep1 = other.fDiagStep1;
125 fDiagStep2 = other.fDiagStep2;
126 fDiagStep3 = other.fDiagStep3;
127 fDiagStep4 = other.fDiagStep4;
128 fDiagAll = other.fDiagAll;
132 //____________________________________________________________________
133 AliFMDReconstructor::~AliFMDReconstructor()
136 if (fMult) fMult->Delete();
137 if (fMult) delete fMult;
138 if (fESDObj) delete fESDObj;
141 //____________________________________________________________________
143 AliFMDReconstructor::Init(AliRunLoader* /*runLoader*/)
145 // Initialize the reconstructor
147 // Initialize the geometry
148 AliFMDGeometry* geom = AliFMDGeometry::Instance();
150 geom->InitTransformations();
152 // Initialize the parameters
153 AliFMDParameters* param = AliFMDParameters::Instance();
156 // Current vertex position
158 // Create array of reconstructed strip multiplicities
159 fMult = new TClonesArray("AliFMDRecPoint", 51200);
160 // Create ESD output object
161 fESDObj = new AliESDFMD;
163 // Check if we need diagnostics histograms
164 if (!fDiagnostics) return;
165 fDiagStep1 = new TH2I("diagStep1", "Read ADC vs. Noise surpressed ADC",
166 1024, -.5, 1023.5, 1024, -.5, 1023.5);
167 fDiagStep1->SetDirectory(0);
168 fDiagStep1->GetXaxis()->SetTitle("ADC (read)");
169 fDiagStep1->GetYaxis()->SetTitle(Form("ADC (noise surpressed %4.f)",
171 fDiagStep2 = new TH2F("diagStep2", "ADC vs Edep deduced",
172 1024, -.5, 1023.5, 100, 0, 2);
173 fDiagStep2->SetDirectory(0);
174 fDiagStep2->GetXaxis()->SetTitle("ADC (noise surpressed)");
175 fDiagStep2->GetYaxis()->SetTitle("#Delta E [GeV]");
176 fDiagStep3 = new TH2F("diagStep3", "Edep vs Edep path corrected",
177 100, 0., 2., 100, 0., 2.);
178 fDiagStep3->SetDirectory(0);
179 fDiagStep3->GetXaxis()->SetTitle("#Delta E [GeV]");
180 fDiagStep3->GetYaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
181 fDiagStep4 = new TH2F("diagStep4", "Edep vs Multiplicity deduced",
182 100, 0., 2., 100, -.1, 19.9);
183 fDiagStep4->SetDirectory(0);
184 fDiagStep4->GetXaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
185 fDiagStep4->GetYaxis()->SetTitle("Multiplicity");
186 fDiagAll = new TH2F("diagAll", "Read ADC vs Multiplicity deduced",
187 1024, -.5, 1023.5, 100, -.1, 19.9);
188 fDiagAll->SetDirectory(0);
189 fDiagAll->GetXaxis()->SetTitle("ADC (read)");
190 fDiagAll->GetYaxis()->SetTitle("Multiplicity");
193 //____________________________________________________________________
195 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
196 TTree* digitsTree) const
198 // Convert Raw digits to AliFMDDigit's in a tree
199 AliFMDDebug(2, ("Reading raw data into digits tree"));
200 AliFMDRawReader rawRead(reader, digitsTree);
201 // rawRead.SetSampleRate(fFMD->GetSampleRate());
205 //____________________________________________________________________
207 AliFMDReconstructor::GetVertex() const
209 // Return the vertex to use.
210 // This is obtained from the ESD object.
211 // If not found, a warning is issued.
212 fVertexType = kNoVertex;
215 const AliESDVertex* vertex = fESD->GetVertex();
217 AliFMDDebug(2, ("Got vertex from ESD: %f", vertex->GetZv()));
218 fCurrentVertex = vertex->GetZv();
219 fVertexType = kESDVertex;
223 AliWarning("Didn't get any vertex from ESD or generator");
227 //____________________________________________________________________
229 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
230 TTree* clusterTree) const
232 // Reconstruct event from digits in tree
233 // Get the FMD branch holding the digits.
234 // FIXME: The vertex may not be known yet, so we may have to move
235 // some of this to FillESD.
236 AliFMDDebug(2, ("Reconstructing from digits in a tree"));
239 TBranch *digitBranch = digitsTree->GetBranch("FMD");
241 Error("Exec", "No digit branch for the FMD found");
244 TClonesArray* digits = new TClonesArray("AliFMDDigit");
245 digitBranch->SetAddress(&digits);
247 if (fMult) fMult->Clear();
248 if (fESDObj) fESDObj->Clear();
251 fTreeR = clusterTree;
252 fTreeR->Branch("FMD", &fMult);
254 AliFMDDebug(5, ("Getting entry 0 from digit branch"));
255 digitBranch->GetEntry(0);
257 AliFMDDebug(5, ("Processing digits"));
258 ProcessDigits(digits);
260 Int_t written = clusterTree->Fill();
261 AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
267 //____________________________________________________________________
269 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
271 // For each digit, find the pseudo rapdity, azimuthal angle, and
272 // number of corrected ADC counts, and pass it on to the algorithms
274 Int_t nDigits = digits->GetEntries();
275 AliFMDDebug(1, ("Got %d digits", nDigits));
276 fESDObj->SetNoiseFactor(fNoiseFactor);
277 fESDObj->SetAngleCorrected(fAngleCorrect);
278 for (Int_t i = 0; i < nDigits; i++) {
279 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
280 AliFMDParameters* param = AliFMDParameters::Instance();
281 // Check that the strip is not marked as dead
282 if (param->IsDead(digit->Detector(), digit->Ring(),
283 digit->Sector(), digit->Strip())) {
284 AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", digit->Detector(),
285 digit->Ring(), digit->Sector(), digit->Strip()));
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 AliFMDDebug(10, ("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]) AliFMDRecPoint(digit->Detector(),
316 (void)m; // Suppress warnings about unused variables.
319 fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(),
320 digit->Sector(), digit->Strip(), mult);
321 fESDObj->SetEta(digit->Detector(), digit->Ring(),
322 digit->Sector(), digit->Strip(), eta);
324 if (fDiagAll) fDiagAll->Fill(digit->Counts(), mult);
328 //____________________________________________________________________
330 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
332 // Member function to subtract the pedestal from a digit
333 // This implementation does nothing, but a derived class could over
334 // load this to subtract a pedestal that was given in a database or
335 // something like that.
339 AliFMDParameters* param = AliFMDParameters::Instance();
340 Float_t ped = param->GetPedestal(digit->Detector(),
344 Float_t noise = param->GetPedestalWidth(digit->Detector(),
348 AliFMDDebug(15, ("Subtracting pedestal %f from signal %d",
349 ped, digit->Counts()));
350 if (digit->Count3() > 0) adc = digit->Count3();
351 else if (digit->Count2() > 0) adc = digit->Count2();
352 else adc = digit->Count1();
353 counts = TMath::Max(Int_t(adc - ped), 0);
354 if (counts < noise * fNoiseFactor) counts = 0;
355 if (counts > 0) AliFMDDebug(15, ("Got a hit strip"));
356 if (fDiagStep1) fDiagStep1->Fill(adc, counts);
358 return UShort_t(counts);
361 //____________________________________________________________________
363 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit,
365 UShort_t count) const
367 // Converts number of ADC counts to energy deposited.
368 // Note, that this member function can be overloaded by derived
369 // classes to do strip-specific look-ups in databases or the like,
370 // to find the proper gain for a strip.
372 // In this simple version, we calculate the energy deposited as
374 // EnergyDeposited = cos(theta) * gain * count
379 // gain = ----------------- * Energy_deposited_per_MIP
382 // is constant and the same for all strips.
383 if (count <= 0) return 0;
384 AliFMDParameters* param = AliFMDParameters::Instance();
385 Float_t gain = param->GetPulseGain(digit->Detector(),
389 AliFMDDebug(15, ("Converting counts %d to energy via factor %f",
392 Double_t edep = count * gain;
393 if (fDiagStep2) fDiagStep2->Fill(count, edep);
395 Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
396 Double_t corr = TMath::Abs(TMath::Cos(theta));
397 Double_t cedep = corr * edep;
398 AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
399 edep, corr, cedep, eta, theta));
400 if (fDiagStep3) fDiagStep3->Fill(edep, cedep);
406 //____________________________________________________________________
408 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */,
411 // Converts an energy signal to number of particles.
412 // Note, that this member function can be overloaded by derived
413 // classes to do strip-specific look-ups in databases or the like,
414 // to find the proper gain for a strip.
416 // In this simple version, we calculate the multiplicity as
418 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
422 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
424 // is constant and the same for all strips
425 AliFMDParameters* param = AliFMDParameters::Instance();
426 Double_t edepMIP = param->GetEdepMip();
427 Float_t mult = edep / edepMIP;
429 AliFMDDebug(15, ("Translating energy %f to multiplicity via "
430 "divider %f->%f", edep, edepMIP, mult));
431 if (fDiagStep4) fDiagStep4->Fill(edep, mult);
435 //____________________________________________________________________
437 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit,
441 // Get the eta and phi of a digit
444 AliFMDGeometry* geom = AliFMDGeometry::Instance();
445 Double_t x, y, z, r, theta;
446 geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(),
447 digit->Strip(), x, y, z);
448 // Correct for vertex offset.
450 phi = TMath::ATan2(y, x);
451 r = TMath::Sqrt(y * y + x * x);
452 theta = TMath::ATan2(r, z);
453 eta = -TMath::Log(TMath::Tan(theta / 2));
458 //____________________________________________________________________
460 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
461 TTree* /* clusterTree */,
462 AliESDEvent* esd) const
464 // nothing to be done
465 // FIXME: The vertex may not be known when Reconstruct is executed,
466 // so we may have to move some of that member function here.
467 AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
471 AliFMDDebug(2, ("Writing FMD data to ESD tree"));
472 esd->SetFMDData(fESDObj);
475 if (!fDiagnostics || !esd) return;
476 static bool first = true;
477 Int_t evno = esd->GetEventNumberInFile(); // This is most likely NOT the event number you'd like to use. It has nothing to do with the 'real' event number.
478 AliFMDDebug(1, ("Writing diagnostics histograms to FMD.Diag.root/%03d",
480 TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
483 TDirectory* d = f.mkdir(Form("%03d", evno),
484 Form("Diagnostics histograms for event # %d", evno));
486 if (fDiagStep1) fDiagStep1->Write();
487 if (fDiagStep2) fDiagStep2->Write();
488 if (fDiagStep3) fDiagStep3->Write();
489 if (fDiagStep4) fDiagStep4->Write();
490 if (fDiagAll) fDiagAll->Write();
495 if (fDiagStep1) fDiagStep1->Reset();
496 if (fDiagStep2) fDiagStep2->Reset();
497 if (fDiagStep3) fDiagStep3->Reset();
498 if (fDiagStep4) fDiagStep4->Reset();
499 if (fDiagAll) fDiagAll->Reset();
503 //____________________________________________________________________
505 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const
507 // Cannot be used. See member function with same name but with 2
508 // TTree arguments. Make sure you do local reconstrucion
509 AliFMDDebug(2, ("Calling FillESD with loader and tree"));
510 AliError("MayNotUse");
512 //____________________________________________________________________
514 AliFMDReconstructor::Reconstruct(AliRunLoader*) const
516 // Cannot be used. See member function with same name but with 2
517 // TTree arguments. Make sure you do local reconstrucion
518 AliFMDDebug(2, ("Calling FillESD with loader"));
519 AliError("MayNotUse");
521 //____________________________________________________________________
523 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const
525 // Cannot be used. See member function with same name but with 2
526 // TTree arguments. Make sure you do local reconstrucion
527 AliFMDDebug(2, ("Calling FillESD with loader and raw reader"));
528 AliError("MayNotUse");
530 //____________________________________________________________________
532 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESDEvent*) const
534 // Cannot be used. See member function with same name but with 2
535 // TTree arguments. Make sure you do local reconstrucion
536 AliFMDDebug(2, ("Calling FillESD with raw reader, tree, and ESD"));
537 AliError("MayNotUse");
539 //____________________________________________________________________
541 AliFMDReconstructor::FillESD(AliRunLoader*,AliESDEvent*) const
543 // Cannot be used. See member function with same name but with 2
544 // TTree arguments. Make sure you do local reconstrucion
545 AliFMDDebug(2, ("Calling FillESD with loader and ESD"));
546 AliError("MayNotUse");
548 //____________________________________________________________________
550 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESDEvent*) const
552 // Cannot be used. See member function with same name but with 2
553 // TTree arguments. Make sure you do local reconstrucion
554 AliFMDDebug(2, ("Calling FillESD with loader, raw reader, and ESD"));
555 AliError("MayNotUse");
558 //____________________________________________________________________