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 #include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
39 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
40 #include "AliFMDAltroMapping.h" // ALIFMDALTROMAPPING_H
41 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
42 #include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
43 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
44 #include "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
45 #include "AliESDEvent.h" // ALIESDEVENT_H
46 #include "AliESDVertex.h" // ALIESDVERTEX_H
47 #include "AliESDTZERO.h" // ALIESDVERTEX_H
48 #include <AliESDFMD.h> // ALIESDFMD_H
54 // Import revertexer into a private namespace (to prevent conflicts)
56 # include "AliFMDESDRevertexer.h"
62 //____________________________________________________________________
63 ClassImp(AliFMDReconstructor)
65 ; // This is here to keep Emacs for indenting the next line
68 //____________________________________________________________________
69 AliFMDReconstructor::AliFMDReconstructor()
78 fVertexType(kNoVertex),
87 // Make a new FMD reconstructor object - default CTOR.
90 if (AliDebugLevel() > 0) fDiagnostics = kTRUE;
94 //____________________________________________________________________
95 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other)
100 fCurrentVertex(other.fCurrentVertex),
101 fESDObj(other.fESDObj),
102 fNoiseFactor(other.fNoiseFactor),
103 fAngleCorrect(other.fAngleCorrect),
104 fVertexType(other.fVertexType),
106 fDiagnostics(other.fDiagnostics),
107 fDiagStep1(other.fDiagStep1),
108 fDiagStep2(other.fDiagStep2),
109 fDiagStep3(other.fDiagStep3),
110 fDiagStep4(other.fDiagStep4),
111 fDiagAll(other.fDiagAll)
117 //____________________________________________________________________
119 AliFMDReconstructor::operator=(const AliFMDReconstructor& other)
121 // Assignment operator
123 fNMult = other.fNMult;
124 fTreeR = other.fTreeR;
125 fCurrentVertex = other.fCurrentVertex;
126 fESDObj = other.fESDObj;
127 fNoiseFactor = other.fNoiseFactor;
128 fAngleCorrect = other.fAngleCorrect;
129 fVertexType = other.fVertexType;
131 fDiagnostics = other.fDiagnostics;
132 fDiagStep1 = other.fDiagStep1;
133 fDiagStep2 = other.fDiagStep2;
134 fDiagStep3 = other.fDiagStep3;
135 fDiagStep4 = other.fDiagStep4;
136 fDiagAll = other.fDiagAll;
140 //____________________________________________________________________
141 AliFMDReconstructor::~AliFMDReconstructor()
144 if (fMult) fMult->Delete();
145 if (fMult) delete fMult;
146 if (fESDObj) delete fESDObj;
149 //____________________________________________________________________
151 AliFMDReconstructor::Init()
153 // Initialize the reconstructor
155 // Initialize the geometry
156 AliFMDGeometry* geom = AliFMDGeometry::Instance();
158 geom->InitTransformations();
160 // Initialize the parameters
161 AliFMDParameters* param = AliFMDParameters::Instance();
164 // Current vertex position
166 // Create array of reconstructed strip multiplicities
167 fMult = new TClonesArray("AliFMDRecPoint", 51200);
168 // Create ESD output object
169 fESDObj = new AliESDFMD;
171 // Check if we need diagnostics histograms
172 if (!fDiagnostics) return;
173 AliInfo("Making diagnostics histograms");
174 fDiagStep1 = new TH2I("diagStep1", "Read ADC vs. Noise surpressed ADC",
175 1024, -.5, 1023.5, 1024, -.5, 1023.5);
176 fDiagStep1->SetDirectory(0);
177 fDiagStep1->GetXaxis()->SetTitle("ADC (read)");
178 fDiagStep1->GetYaxis()->SetTitle(Form("ADC (noise surpressed %4.f)",
180 fDiagStep2 = new TH2F("diagStep2", "ADC vs Edep deduced",
181 1024, -.5, 1023.5, 100, 0, 2);
182 fDiagStep2->SetDirectory(0);
183 fDiagStep2->GetXaxis()->SetTitle("ADC (noise surpressed)");
184 fDiagStep2->GetYaxis()->SetTitle("#Delta E [GeV]");
185 fDiagStep3 = new TH2F("diagStep3", "Edep vs Edep path corrected",
186 100, 0., 2., 100, 0., 2.);
187 fDiagStep3->SetDirectory(0);
188 fDiagStep3->GetXaxis()->SetTitle("#Delta E [GeV]");
189 fDiagStep3->GetYaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
190 fDiagStep4 = new TH2F("diagStep4", "Edep vs Multiplicity deduced",
191 100, 0., 2., 100, -.1, 19.9);
192 fDiagStep4->SetDirectory(0);
193 fDiagStep4->GetXaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
194 fDiagStep4->GetYaxis()->SetTitle("Multiplicity");
195 fDiagAll = new TH2F("diagAll", "Read ADC vs Multiplicity deduced",
196 1024, -.5, 1023.5, 100, -.1, 19.9);
197 fDiagAll->SetDirectory(0);
198 fDiagAll->GetXaxis()->SetTitle("ADC (read)");
199 fDiagAll->GetYaxis()->SetTitle("Multiplicity");
202 //____________________________________________________________________
204 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
205 TTree* digitsTree) const
207 // Convert Raw digits to AliFMDDigit's in a tree
208 AliFMDDebug(2, ("Reading raw data into digits tree"));
209 AliFMDRawReader rawRead(reader, digitsTree);
210 // rawRead.SetSampleRate(fFMD->GetSampleRate());
212 AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
213 for (size_t i = 1; i <= 3; i++) {
214 fZS[i] = rawRead.IsZeroSuppressed(map->Detector2DDL(i));
215 fZSFactor[i] = rawRead.NoiseFactor(map->Detector2DDL(i));
219 //____________________________________________________________________
221 AliFMDReconstructor::GetVertex() const
223 // Return the vertex to use.
224 // This is obtained from the ESD object.
225 // If not found, a warning is issued.
226 fVertexType = kNoVertex;
229 const AliESDVertex* vertex = fESD->GetPrimaryVertex();
230 if (!vertex) vertex = fESD->GetPrimaryVertexSPD();
231 if (!vertex) vertex = fESD->GetPrimaryVertexTPC();
232 if (!vertex) vertex = fESD->GetVertex();
235 AliFMDDebug(2, ("Got %s (%s) from ESD: %f",
236 vertex->GetName(), vertex->GetTitle(), vertex->GetZv()));
237 fCurrentVertex = vertex->GetZv();
238 fVertexType = kESDVertex;
241 else if (fESD->GetESDTZERO()) {
242 AliFMDDebug(2, ("Got primary vertex from T0: %f", fESD->GetT0zVertex()));
243 fCurrentVertex = fESD->GetT0zVertex();
244 fVertexType = kESDVertex;
248 AliWarning("Didn't get any vertex from ESD or generator");
252 //____________________________________________________________________
254 AliFMDReconstructor::Reconstruct(AliRawReader* reader, TTree*) const
256 // Reconstruct directly from raw data (no intermediate output on
257 // digit tree or rec point tree).
260 // reader Raw event reader
262 AliFMDRawReader rawReader(reader, 0);
264 UShort_t det, sec, str, fac;
265 Short_t adc, oldDet = -1;
269 while (rawReader.NextSignal(det, rng, sec, str, adc, zs, fac)) {
272 fZSFactor[det-1] = fac;
275 ProcessSignal(det, rng, sec, str, adc);
279 //____________________________________________________________________
281 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
282 TTree* clusterTree) const
284 // Reconstruct event from digits in tree
285 // Get the FMD branch holding the digits.
286 // FIXME: The vertex may not be known yet, so we may have to move
287 // some of this to FillESD.
290 // digitsTree Pointer to a tree containing digits
291 // clusterTree Pointer to output tree
293 AliFMDDebug(2, ("Reconstructing from digits in a tree"));
296 TBranch *digitBranch = digitsTree->GetBranch("FMD");
298 Error("Exec", "No digit branch for the FMD found");
301 TClonesArray* digits = new TClonesArray("AliFMDDigit");
302 digitBranch->SetAddress(&digits);
304 if (fMult) fMult->Clear();
305 if (fESDObj) fESDObj->Clear();
308 fTreeR = clusterTree;
309 fTreeR->Branch("FMD", &fMult);
311 AliFMDDebug(5, ("Getting entry 0 from digit branch"));
312 digitBranch->GetEntry(0);
314 AliFMDDebug(5, ("Processing digits"));
315 ProcessDigits(digits);
317 Int_t written = clusterTree->Fill();
318 AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
324 //____________________________________________________________________
326 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
328 // For each digit, find the pseudo rapdity, azimuthal angle, and
329 // number of corrected ADC counts, and pass it on to the algorithms
333 // digits Array of digits
335 Int_t nDigits = digits->GetEntries();
336 AliFMDDebug(1, ("Got %d digits", nDigits));
337 fESDObj->SetNoiseFactor(fNoiseFactor);
338 fESDObj->SetAngleCorrected(fAngleCorrect);
339 for (Int_t i = 0; i < nDigits; i++) {
340 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
341 if (!digit) continue;
346 //____________________________________________________________________
348 AliFMDReconstructor::ProcessDigit(AliFMDDigit* digit) const
350 UShort_t det = digit->Detector();
351 Char_t rng = digit->Ring();
352 UShort_t sec = digit->Sector();
353 UShort_t str = digit->Strip();
354 Short_t adc = digit->Counts();
356 ProcessSignal(det, rng, sec, str, adc);
359 //____________________________________________________________________
361 AliFMDReconstructor::ProcessSignal(UShort_t det,
367 // Process the signal from a single strip
376 AliFMDParameters* param = AliFMDParameters::Instance();
377 // Check that the strip is not marked as dead
378 if (param->IsDead(det, rng, sec, str)) {
379 AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
386 PhysicalCoordinates(det, rng, sec, str, eta, phi);
388 // Substract pedestal.
389 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
390 if(counts == USHRT_MAX) return;
392 // Gain match digits.
393 Double_t edep = Adc2Energy(det, rng, sec, str, eta, counts);
394 // Get rid of nonsense energy
397 // Make rough multiplicity
398 Double_t mult = Energy2Multiplicity(det, rng, sec, str, edep);
399 // Get rid of nonsense mult
400 if (mult < 0) return;
401 AliFMDDebug(5, ("FMD%d%c[%2d,%3d]: "
402 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
403 det, rng, sec, str, adc, counts, edep, mult));
405 // Create a `RecPoint' on the output branch.
408 new ((*fMult)[fNMult]) AliFMDRecPoint(det, rng, sec, str,
409 eta, phi, edep, mult);
410 (void)m; // Suppress warnings about unused variables.
414 fESDObj->SetMultiplicity(det, rng, sec, str, mult);
415 fESDObj->SetEta(det, rng, sec, str, eta);
417 if (fDiagAll) fDiagAll->Fill(adc, mult);
423 //____________________________________________________________________
425 AliFMDReconstructor::SubtractPedestal(UShort_t det,
431 // Member function to subtract the pedestal from a digit
438 // adc # of ADC counts
440 // Pedestal subtracted signal or USHRT_MAX in case of problems
442 AliFMDParameters* param = AliFMDParameters::Instance();
443 Bool_t zs = fZS[det-1];
444 UShort_t fac = fZSFactor[det-1];
445 Float_t ped = (zs ? 0 :
446 param->GetPedestal(det, rng, sec, str));
447 Float_t noise = param->GetPedestalWidth(det, rng, sec, str);
448 if(ped < 0 || noise < 0) {
449 AliWarning(Form("Invalid pedestal (%f) or noise (%f) "
450 "for FMD%d%c[%02d,%03d]", ped, noise, det, rng, sec, str));
454 AliFMDDebug(5, ("Subtracting pedestal %f from signal %d", ped, adc));
455 // if (digit->Count3() > 0) adc = digit->Count3();
456 // else if (digit->Count2() > 0) adc = digit->Count2();
457 // else adc = digit->Count1();
458 Int_t counts = adc + Int_t(zs ? fac * noise : - ped);
459 counts = TMath::Max(Int_t(counts), 0);
460 if (counts < noise * fNoiseFactor) counts = 0;
461 if (counts > 0) AliFMDDebug(15, ("Got a hit strip"));
462 if (fDiagStep1) fDiagStep1->Fill(adc, counts);
464 return UShort_t(counts);
467 //____________________________________________________________________
469 AliFMDReconstructor::Adc2Energy(UShort_t det,
474 UShort_t count) const
476 // Converts number of ADC counts to energy deposited.
477 // Note, that this member function can be overloaded by derived
478 // classes to do strip-specific look-ups in databases or the like,
479 // to find the proper gain for a strip.
481 // In the first simple version, we calculate the energy deposited as
483 // EnergyDeposited = cos(theta) * gain * count
488 // gain = ----------------- * Energy_deposited_per_MIP
491 // is constant and the same for all strips.
493 // For the production we use the conversion measured in the NBI lab.
494 // The total conversion is then:
499 // => energy = ----------------
507 // eta Psuedo-rapidity
508 // counts Number of ADC counts over pedestal
510 // The energy deposited in a single strip, or -1 in case of problems
512 if (count <= 0) return 0;
513 AliFMDParameters* param = AliFMDParameters::Instance();
514 Float_t gain = param->GetPulseGain(det, rng, sec, str);
515 // 'Tagging' bad gains as bad energy
517 AliWarning(Form("Invalid gain (%f) for FMD%d%c[%02d,%03d]",
518 gain, det, rng, sec, str));
521 AliFMDDebug(5, ("Converting counts %d to energy (factor=%f, DAC2MIP=%f)",
522 count, gain,param->GetDACPerMIP()));
524 Double_t edep = ((count * param->GetEdepMip())
525 / (gain * param->GetDACPerMIP()));
526 if (fDiagStep2) fDiagStep2->Fill(count, edep);
528 Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
529 Double_t corr = TMath::Abs(TMath::Cos(theta));
530 Double_t cedep = corr * edep;
531 AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
532 edep, corr, cedep, eta, theta));
533 if (fDiagStep3) fDiagStep3->Fill(edep, cedep);
539 //____________________________________________________________________
541 AliFMDReconstructor::Energy2Multiplicity(UShort_t /*det*/,
547 // Converts an energy signal to number of particles.
548 // Note, that this member function can be overloaded by derived
549 // classes to do strip-specific look-ups in databases or the like,
550 // to find the proper gain for a strip.
552 // In this simple version, we calculate the multiplicity as
554 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
558 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
560 // is constant and the same for all strips
567 // edep Energy deposited in a single strip
569 // The "bare" multiplicity corresponding to the energy deposited
570 AliFMDParameters* param = AliFMDParameters::Instance();
571 Double_t edepMIP = param->GetEdepMip();
572 Float_t mult = edep / edepMIP;
575 AliFMDDebug(15, ("Translating energy %f to multiplicity via "
576 "divider %f->%f", edep, edepMIP, mult));
578 if (fDiagStep4) fDiagStep4->Fill(edep, mult);
582 //____________________________________________________________________
584 AliFMDReconstructor::PhysicalCoordinates(UShort_t det,
591 // Get the eta and phi of a digit
598 // eta On return, contains the psuedo-rapidity of the strip
599 // phi On return, contains the azimuthal angle of the strip
601 AliFMDGeometry* geom = AliFMDGeometry::Instance();
602 Double_t x, y, z, r, theta;
603 geom->Detector2XYZ(det, rng, sec, str, x, y, z);
604 // Correct for vertex offset.
606 phi = TMath::ATan2(y, x);
607 r = TMath::Sqrt(y * y + x * x);
608 theta = TMath::ATan2(r, z);
609 eta = -TMath::Log(TMath::Tan(theta / 2));
614 //____________________________________________________________________
616 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
617 TTree* /* clusterTree */,
618 AliESDEvent* esd) const
620 // nothing to be done
621 // FIXME: The vertex may not be known when Reconstruct is executed,
622 // so we may have to move some of that member function here.
623 AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
626 Double_t oldVz = fCurrentVertex;
628 if (fVertexType != kNoVertex) {
629 AliFMDDebug(2, ("Revertexing the ESD data to vz=%f (was %f)",
630 fCurrentVertex, oldVz));
631 AliFMDESDRevertexer revertexer;
632 revertexer.Revertex(fESDObj, fCurrentVertex);
636 AliFMDDebug(2, ("Writing FMD data to ESD tree"));
637 esd->SetFMDData(fESDObj);
640 if (!fDiagnostics || !esd) return;
641 static bool first = true;
642 // This is most likely NOT the event number you'd like to use. It
643 // has nothing to do with the 'real' event number.
644 // - That's OK. We just use it for the name of the directory -
645 // nothing else. Christian
646 Int_t evno = esd->GetEventNumberInFile();
647 AliFMDDebug(1, ("Writing diagnostics histograms to FMD.Diag.root/%03d",evno));
648 TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
651 TDirectory* d = f.mkdir(Form("%03d", evno),
652 Form("Diagnostics histograms for event # %d", evno));
654 if (fDiagStep1) fDiagStep1->Write();
655 if (fDiagStep2) fDiagStep2->Write();
656 if (fDiagStep3) fDiagStep3->Write();
657 if (fDiagStep4) fDiagStep4->Write();
658 if (fDiagAll) fDiagAll->Write();
663 if (fDiagStep1) fDiagStep1->Reset();
664 if (fDiagStep2) fDiagStep2->Reset();
665 if (fDiagStep3) fDiagStep3->Reset();
666 if (fDiagStep4) fDiagStep4->Reset();
667 if (fDiagAll) fDiagAll->Reset();
670 //____________________________________________________________________
672 AliFMDReconstructor::FillESD(AliRawReader*, TTree* clusterTree,
673 AliESDEvent* esd) const
676 FillESD(dummy, clusterTree, esd);
679 //____________________________________________________________________