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 "AliFMDSDigit.h" // ALIFMDDIGIT_H
43 #include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
44 #include "AliFMDRecoParam.h" // ALIFMDRECOPARAM_H
45 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
46 #include "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
47 #include "AliESDEvent.h" // ALIESDEVENT_H
48 #include "AliESDVertex.h" // ALIESDVERTEX_H
49 #include "AliESDTZERO.h" // ALIESDVERTEX_H
50 #include <AliESDFMD.h> // ALIESDFMD_H
56 #include "AliFMDESDRevertexer.h"
61 //____________________________________________________________________
62 ClassImp(AliFMDReconstructor)
64 ; // This is here to keep Emacs for indenting the next line
67 //____________________________________________________________________
68 AliFMDReconstructor::AliFMDReconstructor()
77 fVertexType(kNoVertex),
86 // Make a new FMD reconstructor object - default CTOR.
89 if (AliDebugLevel() > 0) fDiagnostics = kTRUE;
90 for(Int_t det = 1; det<=3; det++) {
96 //____________________________________________________________________
97 AliFMDReconstructor::~AliFMDReconstructor()
100 if (fMult) fMult->Delete();
101 if (fMult) delete fMult;
102 if (fESDObj) delete fESDObj;
105 //____________________________________________________________________
107 AliFMDReconstructor::Init()
109 // Initialize the reconstructor
111 // Initialize the geometry
112 AliFMDGeometry* geom = AliFMDGeometry::Instance();
114 geom->InitTransformations();
116 // Initialize the parameters
117 AliFMDParameters* param = AliFMDParameters::Instance();
120 // Current vertex position
122 // Create array of reconstructed strip multiplicities
123 fMult = new TClonesArray("AliFMDRecPoint", 51200);
124 // Create ESD output object
125 fESDObj = new AliESDFMD;
127 // Check if we need diagnostics histograms
128 if (!fDiagnostics) return;
129 AliInfo("Making diagnostics histograms");
130 fDiagStep1 = new TH2I("diagStep1", "Read ADC vs. Noise surpressed ADC",
131 1024, -.5, 1023.5, 1024, -.5, 1023.5);
132 fDiagStep1->SetDirectory(0);
133 fDiagStep1->GetXaxis()->SetTitle("ADC (read)");
134 fDiagStep1->GetYaxis()->SetTitle(Form("ADC (noise surpressed %4.f)",
136 fDiagStep2 = new TH2F("diagStep2", "ADC vs Edep deduced",
137 1024, -.5, 1023.5, 100, 0, 2);
138 fDiagStep2->SetDirectory(0);
139 fDiagStep2->GetXaxis()->SetTitle("ADC (noise surpressed)");
140 fDiagStep2->GetYaxis()->SetTitle("#Delta E [GeV]");
141 fDiagStep3 = new TH2F("diagStep3", "Edep vs Edep path corrected",
142 100, 0., 2., 100, 0., 2.);
143 fDiagStep3->SetDirectory(0);
144 fDiagStep3->GetXaxis()->SetTitle("#Delta E [GeV]");
145 fDiagStep3->GetYaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
146 fDiagStep4 = new TH2F("diagStep4", "Edep vs Multiplicity deduced",
147 100, 0., 2., 100, -.1, 19.9);
148 fDiagStep4->SetDirectory(0);
149 fDiagStep4->GetXaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
150 fDiagStep4->GetYaxis()->SetTitle("Multiplicity");
151 fDiagAll = new TH2F("diagAll", "Read ADC vs Multiplicity deduced",
152 1024, -.5, 1023.5, 100, -.1, 19.9);
153 fDiagAll->SetDirectory(0);
154 fDiagAll->GetXaxis()->SetTitle("ADC (read)");
155 fDiagAll->GetYaxis()->SetTitle("Multiplicity");
158 //____________________________________________________________________
160 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
161 TTree* digitsTree) const
163 // Convert Raw digits to AliFMDDigit's in a tree
164 AliFMDDebug(1, ("Reading raw data into digits tree"));
165 AliFMDRawReader rawRead(reader, digitsTree);
166 // rawRead.SetSampleRate(fFMD->GetSampleRate());
168 AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
169 for (size_t i = 1; i <= 3; i++) {
170 fZS[i] = rawRead.IsZeroSuppressed(map->Detector2DDL(i));
171 fZSFactor[i] = rawRead.NoiseFactor(map->Detector2DDL(i));
175 //____________________________________________________________________
177 AliFMDReconstructor::GetVertex(AliESDEvent* esd) const
179 // Return the vertex to use.
180 // This is obtained from the ESD object.
181 // If not found, a warning is issued.
182 fVertexType = kNoVertex;
186 const AliESDVertex* vertex = esd->GetPrimaryVertex();
187 if (!vertex) vertex = esd->GetPrimaryVertexSPD();
188 if (!vertex) vertex = esd->GetPrimaryVertexTPC();
189 if (!vertex) vertex = esd->GetVertex();
192 AliFMDDebug(2, ("Got %s (%s) from ESD: %f",
193 vertex->GetName(), vertex->GetTitle(), vertex->GetZv()));
194 fCurrentVertex = vertex->GetZv();
195 fVertexType = kESDVertex;
198 else if (esd->GetESDTZERO()) {
199 AliFMDDebug(2, ("Got primary vertex from T0: %f", esd->GetT0zVertex()));
200 fCurrentVertex = esd->GetT0zVertex();
201 fVertexType = kESDVertex;
204 AliWarning("Didn't get any vertex from ESD or generator");
207 //____________________________________________________________________
209 AliFMDReconstructor::GetIdentifier() const
211 return AliReconstruction::GetDetIndex(GetDetectorName());
214 //____________________________________________________________________
215 const AliFMDRecoParam*
216 AliFMDReconstructor::GetParameters() const
218 Int_t iDet = 12; // GetIdentifier();
219 const AliDetectorRecoParam* params = AliReconstructor::GetRecoParam(iDet);
220 if (!params || params->IsA() != AliFMDRecoParam::Class()) return 0;
221 return static_cast<const AliFMDRecoParam*>(params);
224 //____________________________________________________________________
226 AliFMDReconstructor::UseRecoParam(Bool_t set) const
228 static Float_t savedNoiseFactor = fNoiseFactor;
229 static Bool_t savedAngleCorrect = fAngleCorrect;
231 const AliFMDRecoParam* params = GetParameters();
233 fNoiseFactor = params->NoiseFactor();
234 fAngleCorrect = params->AngleCorrect();
238 fNoiseFactor = savedNoiseFactor;
239 fAngleCorrect = savedAngleCorrect;
244 //____________________________________________________________________
246 AliFMDReconstructor::Reconstruct(AliRawReader* reader, TTree*) const
248 // Reconstruct directly from raw data (no intermediate output on
249 // digit tree or rec point tree).
252 // reader Raw event reader
254 AliFMDDebug(1, ("Reconstructing from raw reader"));
255 AliFMDRawReader rawReader(reader, 0);
257 UShort_t det, sec, str, fac;
258 Short_t adc, oldDet = -1;
263 while (rawReader.NextSignal(det, rng, sec, str, adc, zs, fac)) {
266 fZSFactor[det-1] = fac;
269 ProcessSignal(det, rng, sec, str, adc);
271 UseRecoParam(kFALSE);
275 //____________________________________________________________________
277 AliFMDReconstructor::Digitize(AliRawReader* reader, TClonesArray* sdigits) const
279 // Reconstruct directly from raw data (no intermediate output on
280 // digit tree or rec point tree).
283 // reader Raw event reader
285 AliFMDRawReader rawReader(reader, 0);
287 UShort_t det, sec, str, sam, rat, fac;
288 Short_t adc, oldDet = -1;
293 while (rawReader.NextSample(det, rng, sec, str, sam, rat, adc, zs, fac)) {
294 if (!rawReader.SelectSample(sam, rat)) continue;
297 fZSFactor[det-1] = fac;
300 DigitizeSignal(sdigits, det, rng, sec, str, sam, adc);
302 UseRecoParam(kFALSE);
305 //____________________________________________________________________
307 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
308 TTree* clusterTree) const
310 // Reconstruct event from digits in tree
311 // Get the FMD branch holding the digits.
312 // FIXME: The vertex may not be known yet, so we may have to move
313 // some of this to FillESD.
316 // digitsTree Pointer to a tree containing digits
317 // clusterTree Pointer to output tree
319 AliFMDDebug(1, ("Reconstructing from digits in a tree"));
322 static TClonesArray* digits = new TClonesArray("AliFMDDigit");
323 TBranch *digitBranch = digitsTree->GetBranch("FMD");
325 Error("Exec", "No digit branch for the FMD found");
328 // TClonesArray* digits = new TClonesArray("AliFMDDigit");
329 digitBranch->SetAddress(&digits);
331 if (fMult) fMult->Clear();
332 if (fESDObj) fESDObj->Clear();
335 fTreeR = clusterTree;
336 fTreeR->Branch("FMD", &fMult);
338 AliFMDDebug(5, ("Getting entry 0 from digit branch"));
339 digitBranch->GetEntry(0);
341 AliFMDDebug(5, ("Processing digits"));
343 ProcessDigits(digits);
344 UseRecoParam(kFALSE);
346 Int_t written = clusterTree->Fill();
347 AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
353 //____________________________________________________________________
355 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
357 // For each digit, find the pseudo rapdity, azimuthal angle, and
358 // number of corrected ADC counts, and pass it on to the algorithms
362 // digits Array of digits
364 Int_t nDigits = digits->GetEntries();
365 AliFMDDebug(2, ("Got %d digits", nDigits));
366 fESDObj->SetNoiseFactor(fNoiseFactor);
367 fESDObj->SetAngleCorrected(fAngleCorrect);
368 for (Int_t i = 0; i < nDigits; i++) {
369 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
370 if (!digit) continue;
375 //____________________________________________________________________
377 AliFMDReconstructor::ProcessDigit(AliFMDDigit* digit) const
379 UShort_t det = digit->Detector();
380 Char_t rng = digit->Ring();
381 UShort_t sec = digit->Sector();
382 UShort_t str = digit->Strip();
383 Short_t adc = digit->Counts();
385 ProcessSignal(det, rng, sec, str, adc);
388 //____________________________________________________________________
390 AliFMDReconstructor::ProcessSignal(UShort_t det,
396 // Process the signal from a single strip
405 AliFMDParameters* param = AliFMDParameters::Instance();
406 // Check that the strip is not marked as dead
407 if (param->IsDead(det, rng, sec, str)) {
408 AliFMDDebug(1, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
415 PhysicalCoordinates(det, rng, sec, str, eta, phi);
417 // Substract pedestal.
418 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
419 if(counts == USHRT_MAX) return;
421 // Gain match digits.
422 Double_t edep = Adc2Energy(det, rng, sec, str, eta, counts);
423 // Get rid of nonsense energy
426 // Make rough multiplicity
427 Double_t mult = Energy2Multiplicity(det, rng, sec, str, edep);
428 // Get rid of nonsense mult
430 // AliWarning(Form("The mutliplicity in FMD%d%c[%2d,%3d]=%f > 20 "
431 // "(ADC: %d, Energy: %f)", det, rng, sec, str, mult,
434 if (mult < 0) return;
435 AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
436 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
437 det, rng, sec, str, adc, counts, edep, mult));
439 // Create a `RecPoint' on the output branch.
442 new ((*fMult)[fNMult]) AliFMDRecPoint(det, rng, sec, str,
443 eta, phi, edep, mult);
444 (void)m; // Suppress warnings about unused variables.
448 fESDObj->SetMultiplicity(det, rng, sec, str, mult);
449 fESDObj->SetEta(det, rng, sec, str, eta);
451 if (fDiagAll) fDiagAll->Fill(adc, mult);
455 //____________________________________________________________________
457 AliFMDReconstructor::DigitizeSignal(TClonesArray* sdigits,
465 // Process the signal from a single strip
474 AliFMDParameters* param = AliFMDParameters::Instance();
475 // Check that the strip is not marked as dead
476 if (param->IsDead(det, rng, sec, str)) {
477 AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
481 // Substract pedestal.
482 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
483 if(counts == USHRT_MAX || counts == 0) return;
485 // Gain match digits.
486 Double_t edep = Adc2Energy(det, rng, sec, str, counts);
487 // Get rid of nonsense energy
490 Int_t n = sdigits->GetEntriesFast();
491 // AliFMDSDigit* sdigit =
493 AliFMDSDigit(det, rng, sec, str, edep, counts, counts, counts, counts);
494 // sdigit->SetCount(sam, counts);
497 //____________________________________________________________________
499 AliFMDReconstructor::SubtractPedestal(UShort_t det,
506 UShort_t zsNoiseFactor) const
508 AliFMDParameters* param = AliFMDParameters::Instance();
509 Float_t ped = (zsEnabled ? 0 :
510 param->GetPedestal(det, rng, sec, str));
511 Float_t noise = param->GetPedestalWidth(det, rng, sec, str);
512 if(ped < 0 || noise < 0) {
513 AliWarningClass(Form("Invalid pedestal (%f) or noise (%f) "
514 "for FMD%d%c[%02d,%03d]",
515 ped, noise, det, rng, sec, str));
518 AliDebugClass(15, Form("Subtracting pedestal for FMD%d%c[%2d,%3d]=%4d "
519 "(%s w/factor %d, noise factor %f, "
520 "pedestal %8.2f+/-%8.2f)",
521 det, rng, sec, str, adc,
522 (zsEnabled ? "zs'ed" : "straight"),
523 zsNoiseFactor, noiseFactor, ped, noise));
525 Int_t counts = adc + Int_t(zsEnabled ? zsNoiseFactor * noise : - ped);
526 counts = TMath::Max(Int_t(counts), 0);
527 // Calculate the noise factor for suppressing remenants of the noise
528 // peak. If we have done on-line zero suppression, we only check
529 // for noise signals that are larger than the suppressed noise. If
530 // the noise factor used on line is larger than the factor used
531 // here, we do not do this check at all.
534 // Online factor | Read factor | Result
535 // ---------------+--------------+-------------------------------
536 // 2 | 3 | Check if signal > 1 * noise
537 // 3 | 3 | Check if signal > 0
538 // 3 | 2 | Check if signal > 0
540 // In this way, we make sure that we do not suppress away too much
541 // data, and that the read-factor is the most stringent cut.
542 Float_t nf = TMath::Max(0.F, noiseFactor - (zsEnabled ? zsNoiseFactor : 0));
543 if (counts < noise * nf) counts = 0;
544 if (counts > 0) AliDebugClass(15, "Got a hit strip");
546 UShort_t ret = counts < 0 ? 0 : counts;
551 //____________________________________________________________________
553 AliFMDReconstructor::SubtractPedestal(UShort_t det,
559 // Member function to subtract the pedestal from a digit
566 // adc # of ADC counts
568 // Pedestal subtracted signal or USHRT_MAX in case of problems
570 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc,
571 fNoiseFactor, fZS[det-1],
573 if (fDiagStep1) fDiagStep1->Fill(adc, counts);
578 //____________________________________________________________________
580 AliFMDReconstructor::Adc2Energy(UShort_t det,
584 UShort_t count) const
586 // Converts number of ADC counts to energy deposited.
587 // Note, that this member function can be overloaded by derived
588 // classes to do strip-specific look-ups in databases or the like,
589 // to find the proper gain for a strip.
591 // In the first simple version, we calculate the energy deposited as
593 // EnergyDeposited = cos(theta) * gain * count
598 // gain = ----------------- * Energy_deposited_per_MIP
601 // is constant and the same for all strips.
603 // For the production we use the conversion measured in the NBI lab.
604 // The total conversion is then:
609 // => energy = ----------------
617 // counts Number of ADC counts over pedestal
619 // The energy deposited in a single strip, or -1 in case of problems
621 if (count <= 0) return 0;
622 AliFMDParameters* param = AliFMDParameters::Instance();
623 Float_t gain = param->GetPulseGain(det, rng, sec, str);
624 // 'Tagging' bad gains as bad energy
626 AliWarning(Form("Invalid gain (%f) for FMD%d%c[%02d,%03d]",
627 gain, det, rng, sec, str));
630 AliFMDDebug(5, ("Converting counts %d to energy (factor=%f, DAC2MIP=%f)",
631 count, gain,param->GetDACPerMIP()));
633 Double_t edep = ((count * param->GetEdepMip())
634 / (gain * param->GetDACPerMIP()));
638 //____________________________________________________________________
640 AliFMDReconstructor::Adc2Energy(UShort_t det,
645 UShort_t count) const
647 // Converts number of ADC counts to energy deposited.
648 // Note, that this member function can be overloaded by derived
649 // classes to do strip-specific look-ups in databases or the like,
650 // to find the proper gain for a strip.
652 // In the first simple version, we calculate the energy deposited as
654 // EnergyDeposited = cos(theta) * gain * count
659 // gain = ----------------- * Energy_deposited_per_MIP
662 // is constant and the same for all strips.
664 // For the production we use the conversion measured in the NBI lab.
665 // The total conversion is then:
670 // => energy = ----------------
678 // eta Psuedo-rapidity
679 // counts Number of ADC counts over pedestal
681 // The energy deposited in a single strip, or -1 in case of problems
683 Double_t edep = Adc2Energy(det, rng, sec, str, count);
685 if (fDiagStep2) fDiagStep2->Fill(count, edep);
687 Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
688 Double_t corr = TMath::Abs(TMath::Cos(theta));
689 Double_t cedep = corr * edep;
690 AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
691 edep, corr, cedep, eta, theta));
692 if (fDiagStep3) fDiagStep3->Fill(edep, cedep);
698 //____________________________________________________________________
700 AliFMDReconstructor::Energy2Multiplicity(UShort_t /*det*/,
706 // Converts an energy signal to number of particles.
707 // Note, that this member function can be overloaded by derived
708 // classes to do strip-specific look-ups in databases or the like,
709 // to find the proper gain for a strip.
711 // In this simple version, we calculate the multiplicity as
713 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
717 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
719 // is constant and the same for all strips
726 // edep Energy deposited in a single strip
728 // The "bare" multiplicity corresponding to the energy deposited
729 AliFMDParameters* param = AliFMDParameters::Instance();
730 Double_t edepMIP = param->GetEdepMip();
731 Float_t mult = edep / edepMIP;
734 AliFMDDebug(15, ("Translating energy %f to multiplicity via "
735 "divider %f->%f", edep, edepMIP, mult));
737 if (fDiagStep4) fDiagStep4->Fill(edep, mult);
741 //____________________________________________________________________
743 AliFMDReconstructor::PhysicalCoordinates(UShort_t det,
750 // Get the eta and phi of a digit
757 // eta On return, contains the psuedo-rapidity of the strip
758 // phi On return, contains the azimuthal angle of the strip
760 AliFMDGeometry* geom = AliFMDGeometry::Instance();
761 Double_t x, y, z, r, theta;
762 geom->Detector2XYZ(det, rng, sec, str, x, y, z);
763 // Correct for vertex offset.
765 phi = TMath::ATan2(y, x);
766 r = TMath::Sqrt(y * y + x * x);
767 theta = TMath::ATan2(r, z);
768 eta = -TMath::Log(TMath::Tan(theta / 2));
773 //____________________________________________________________________
775 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
776 TTree* /* clusterTree */,
777 AliESDEvent* esd) const
779 // nothing to be done
780 // FIXME: The vertex may not be known when Reconstruct is executed,
781 // so we may have to move some of that member function here.
782 AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
785 Double_t oldVz = fCurrentVertex;
787 if (fVertexType != kNoVertex) {
788 AliFMDDebug(2, ("Revertexing the ESD data to vz=%f (was %f)",
789 fCurrentVertex, oldVz));
790 AliFMDESDRevertexer revertexer;
791 revertexer.Revertex(fESDObj, fCurrentVertex);
795 AliFMDDebug(2, ("Writing FMD data to ESD tree"));
796 esd->SetFMDData(fESDObj);
799 if (!fDiagnostics || !esd) return;
800 static bool first = true;
801 // This is most likely NOT the event number you'd like to use. It
802 // has nothing to do with the 'real' event number.
803 // - That's OK. We just use it for the name of the directory -
804 // nothing else. Christian
805 Int_t evno = esd->GetEventNumberInFile();
806 AliFMDDebug(3, ("Writing diagnostics histograms to FMD.Diag.root/%03d",evno));
807 TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
810 TDirectory* d = f.mkdir(Form("%03d", evno),
811 Form("Diagnostics histograms for event # %d", evno));
813 if (fDiagStep1) fDiagStep1->Write();
814 if (fDiagStep2) fDiagStep2->Write();
815 if (fDiagStep3) fDiagStep3->Write();
816 if (fDiagStep4) fDiagStep4->Write();
817 if (fDiagAll) fDiagAll->Write();
822 if (fDiagStep1) fDiagStep1->Reset();
823 if (fDiagStep2) fDiagStep2->Reset();
824 if (fDiagStep3) fDiagStep3->Reset();
825 if (fDiagStep4) fDiagStep4->Reset();
826 if (fDiagAll) fDiagAll->Reset();
829 //____________________________________________________________________
831 AliFMDReconstructor::FillESD(AliRawReader*, TTree* clusterTree,
832 AliESDEvent* esd) const
835 FillESD(dummy, clusterTree, esd);
837 //____________________________________________________________________