1 //____________________________________________________________________
3 // This is a class that constructs AliFMDRecPoint objects from of Digits
4 // This class reads either digits from a TClonesArray or raw data from
5 // a DDL file (or similar), and stores the read ADC counts in an
6 // internal cache (fAdcs). The rec-points are made via the naiive
9 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
10 // Latest changes by Christian Holm Christensen <cholm@nbi.dk>
13 //____________________________________________________________________
14 /**************************************************************************
15 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
17 * Author: The ALICE Off-line Project. *
18 * Contributors are mentioned in the code where appropriate. *
20 * Permission to use, copy, modify and distribute this software and its *
21 * documentation strictly for non-commercial purposes is hereby granted *
22 * without fee, provided that the above copyright notice appears in all *
23 * copies and that both the copyright notice and this permission notice *
24 * appear in the supporting documentation. The authors make no claims *
25 * about the suitability of this software for any purpose. It is *
26 * provided "as is" without express or implied warranty. *
27 **************************************************************************/
30 * @file AliFMDReconstructor.cxx
31 * @author Christian Holm Christensen <cholm@nbi.dk>
32 * @date Mon Mar 27 12:47:09 2006
33 * @brief FMD reconstruction
36 // #include <AliLog.h> // ALILOG_H
37 // #include <AliRun.h> // ALIRUN_H
38 #include "AliFMDDebug.h"
39 #include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
40 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
41 #include "AliFMDAltroMapping.h" // ALIFMDALTROMAPPING_H
42 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
43 #include "AliFMDSDigit.h" // ALIFMDDIGIT_H
44 #include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
45 #include "AliFMDRecoParam.h" // ALIFMDRECOPARAM_H
46 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
47 #include "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
48 #include "AliESDEvent.h" // ALIESDEVENT_H
49 #include "AliESDVertex.h" // ALIESDVERTEX_H
50 #include "AliESDTZERO.h" // ALIESDVERTEX_H
51 #include <AliESDFMD.h> // ALIESDFMD_H
57 #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),
89 // Make a new FMD reconstructor object - default CTOR.
92 if (AliDebugLevel() > 0) fDiagnostics = kTRUE;
93 for(Int_t det = 1; det<=3; det++) {
99 //____________________________________________________________________
100 AliFMDReconstructor::~AliFMDReconstructor()
103 if (fMult) fMult->Delete();
104 if (fMult) delete fMult;
105 if (fESDObj) delete fESDObj;
108 //____________________________________________________________________
110 AliFMDReconstructor::Init()
112 // Initialize the reconstructor
114 // Initialize the geometry
115 AliFMDGeometry* geom = AliFMDGeometry::Instance();
117 geom->InitTransformations();
119 // Initialize the parameters
120 AliFMDParameters* param = AliFMDParameters::Instance();
121 if (param->Init() != 0) {
122 AliError("Failed to initialize parameters, making zombie");
128 // Current vertex position
130 // Create array of reconstructed strip multiplicities
131 // fMult = new TClonesArray("AliFMDRecPoint", 51200);
132 // Create ESD output object
133 fESDObj = new AliESDFMD;
135 // Check if we need diagnostics histograms
136 if (!fDiagnostics) return;
137 AliInfo("Making diagnostics histograms");
139 fDiagStep1 = new TH2I("diagStep1", "Read ADC vs. Noise surpressed ADC",
140 1024, -.5, 1023.5, 1024, -.5, 1023.5);
141 fDiagStep1->SetDirectory(0);
142 fDiagStep1->GetXaxis()->SetTitle("ADC (read)");
143 fDiagStep1->GetYaxis()->SetTitle(Form("ADC (noise surpressed %4.f)",
147 fDiagStep2 = new TH2F("diagStep2", "ADC vs Edep deduced",
148 1024, -.5, 1023.5, 100, 0, 2);
149 fDiagStep2->SetDirectory(0);
150 fDiagStep2->GetXaxis()->SetTitle("ADC (noise surpressed)");
151 fDiagStep2->GetYaxis()->SetTitle("#Delta E [GeV]");
154 fDiagStep3 = new TH2F("diagStep3", "Edep vs Edep path corrected",
155 100, 0., 2., 100, 0., 2.);
156 fDiagStep3->SetDirectory(0);
157 fDiagStep3->GetXaxis()->SetTitle("#Delta E [GeV]");
158 fDiagStep3->GetYaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
161 fDiagStep4 = new TH2F("diagStep4", "Edep vs Multiplicity deduced",
162 100, 0., 2., 100, -.1, 19.9);
163 fDiagStep4->SetDirectory(0);
164 fDiagStep4->GetXaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
165 fDiagStep4->GetYaxis()->SetTitle("Multiplicity");
166 fDiagAll = new TH2F("diagAll", "Read ADC vs Multiplicity deduced",
167 1024, -.5, 1023.5, 100, -.1, 19.9);
170 fDiagAll->SetDirectory(0);
171 fDiagAll->GetXaxis()->SetTitle("ADC (read)");
172 fDiagAll->GetYaxis()->SetTitle("Multiplicity");
176 //____________________________________________________________________
178 AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
179 TTree* digitsTree) const
181 // Convert Raw digits to AliFMDDigit's in a tree
183 AliWarning("I'm a zombie - cannot do anything");
186 AliFMDDebug(1, ("Reading raw data into digits tree"));
188 AliError("No digits tree passed");
191 static TClonesArray* array = new TClonesArray("AliFMDDigit");
192 digitsTree->Branch("FMD", &array);
195 AliFMDRawReader rawRead(reader, digitsTree);
196 // rawRead.SetSampleRate(fFMD->GetSampleRate());
198 rawRead.ReadAdcs(array);
200 Int_t nWrite = digitsTree->Fill();
201 AliFMDDebug(1, ("Got a grand total of %d digits, wrote %d bytes to tree",
202 array->GetEntriesFast(), nWrite));
205 AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
206 for (size_t i = 1; i <= 3; i++) {
207 fZS[i-1] = rawRead.IsZeroSuppressed(map->Detector2DDL(i));
208 fZSFactor[i-1] = rawRead.NoiseFactor(map->Detector2DDL(i));
212 //____________________________________________________________________
214 AliFMDReconstructor::GetVertex(AliESDEvent* esd) const
216 // Return the vertex to use.
217 // This is obtained from the ESD object.
218 // If not found, a warning is issued.
219 fVertexType = kNoVertex;
223 const AliESDVertex* vertex = esd->GetPrimaryVertex();
224 if (!vertex) vertex = esd->GetPrimaryVertexSPD();
225 if (!vertex) vertex = esd->GetPrimaryVertexTPC();
226 if (!vertex) vertex = esd->GetVertex();
229 AliFMDDebug(2, ("Got %s (%s) from ESD: %f",
230 vertex->GetName(), vertex->GetTitle(), vertex->GetZv()));
231 fCurrentVertex = vertex->GetZv();
232 fVertexType = kESDVertex;
235 else if (esd->GetESDTZERO()) {
236 AliFMDDebug(2, ("Got primary vertex from T0: %f", esd->GetT0zVertex()));
237 fCurrentVertex = esd->GetT0zVertex();
238 fVertexType = kESDVertex;
241 AliWarning("Didn't get any vertex from ESD or generator");
244 //____________________________________________________________________
246 AliFMDReconstructor::GetIdentifier() const
248 // Get the detector identifier.
249 // Note the actual value is cached so that we do not
250 // need to do many expensive string comparisons.
251 static Int_t idx = AliReconstruction::GetDetIndex(GetDetectorName());
255 //____________________________________________________________________
256 const AliFMDRecoParam*
257 AliFMDReconstructor::GetParameters() const
259 // Get the reconstruction parameters.
262 // Pointer to reconstruction parameters or null if not found or wrong type
263 Int_t iDet = GetIdentifier(); // Was 12 - but changed on Cvetans request
264 const AliDetectorRecoParam* params = AliReconstructor::GetRecoParam(iDet);
265 if (!params || params->IsA() != AliFMDRecoParam::Class()) return 0;
266 return static_cast<const AliFMDRecoParam*>(params);
269 //____________________________________________________________________
271 AliFMDReconstructor::UseRecoParam(Bool_t set) const
274 // Set-up reconstructor to use values from reconstruction
275 // parameters, if present, for this event. If the argument @a set
276 // is @c false, then restore preset values.
281 static Float_t savedNoiseFactor = fNoiseFactor;
282 static Bool_t savedAngleCorrect = fAngleCorrect;
284 const AliFMDRecoParam* params = GetParameters();
286 fNoiseFactor = params->NoiseFactor();
287 fAngleCorrect = params->AngleCorrect();
291 fNoiseFactor = savedNoiseFactor;
292 fAngleCorrect = savedAngleCorrect;
296 //____________________________________________________________________
298 AliFMDReconstructor::MarkDeadChannels(AliESDFMD* esd) const
300 // Loop over all entries of the ESD and mark
301 // those that are dead as such
302 // - otherwise put in the zero signal.
303 AliFMDParameters* param = AliFMDParameters::Instance();
305 for (UShort_t d = 1; d <= 3; d++) {
306 UShort_t nR = (d == 1 ? 1 : 2);
308 for (UShort_t q = 0; q < nR; q++) {
309 Char_t r = (q == 0 ? 'I' : 'O');
310 UShort_t nS = (q == 0 ? 20 : 40);
311 UShort_t nT = (q == 0 ? 512 : 256);
313 for (UShort_t s = 0; s < nS; s++) {
314 for (UShort_t t = 0; t < nT; t++) {
315 // A strip can be marked `bad' for two reasons:
317 // - The raw reader fails to read the value
318 // - The strip is marked in OCDB as bad (IsDead)
320 // Hence, a dead strip will always be marked kInvalid here,
321 // while a non-dead bad-read strip will be filled with 0.
322 if (fBad(d, r, s, t)) {
323 AliDebug(5, Form("Marking FMD%d%c[%2d,%3d] as bad", d, r, s, t));
324 esd->SetMultiplicity(d, r, s, t, AliESDFMD::kInvalidMult);
326 if (param->IsDead(d, r, s, t)) {
327 AliDebug(5, Form("Marking FMD%d%c[%2d,%3d] as dead", d, r, s, t));
328 esd->SetMultiplicity(d, r, s, t, AliESDFMD::kInvalidMult);
329 // esd->SetEta(d, r, s, t, AliESDFMD::kInvalidEta);
331 else if (esd->Multiplicity(d, r, s, t) == AliESDFMD::kInvalidMult) {
332 AliDebug(10, Form("Setting null signal in FMD%d%c[%2d,%3d]",
334 esd->SetMultiplicity(d, r, s, t, 0);
342 //____________________________________________________________________
344 AliFMDReconstructor::PreReconstruct() const
346 AliFMDDebug(1, ("Before reoconstructing"));
348 AliWarning("I'm a zombie - cannot do anything");
358 // Reset the output ESD
362 // Pre-set eta values
363 for (UShort_t d=1; d<=3; d++) {
364 UShort_t nQ = (d == 1 ? 1 : 2);
365 for (UShort_t q=0; q<nQ; q++) {
366 UShort_t nStr = (q == 0 ? 512 : 256);
367 Char_t r = (q == 0 ? 'I' : 'O');
369 for (UShort_t t = 0; t < nStr; t++) {
371 // Always use sector 0
372 PhysicalCoordinates(d, r, 0, t, eta, phi);
373 fESDObj->SetEta(d, r, 0, t, eta);
383 //____________________________________________________________________
385 AliFMDReconstructor::Reconstruct(AliFMDRawReader& rawReader) const
387 // Reconstruct directly from raw data (no intermediate output on
388 // digit tree or rec point tree).
391 // rawReader FMD Raw event reader
392 AliFMDDebug(1, ("Reconstructing from FMD raw reader"));
393 if (!PreReconstruct()) return;
395 UShort_t det, sec, str, fac;
396 Short_t adc, oldDet = -1;
401 while (rawReader.NextSignal(det, rng, sec, str, adc, zs, fac)) {
404 fZSFactor[det-1] = fac;
407 ProcessSignal(det, rng, sec, str, adc);
409 UseRecoParam(kFALSE);
413 //____________________________________________________________________
415 AliFMDReconstructor::Reconstruct(AliRawReader* reader, TTree*) const
417 // Reconstruct directly from raw data (no intermediate output on
418 // digit tree or rec point tree).
421 // reader Raw event reader
422 // ctree Not used - 'cluster tree' to store rec-points in.
423 AliFMDDebug(1, ("Reconstructing from raw reader"));
425 AliWarning("I'm a zombie - cannot do anything");
428 AliFMDRawReader rawReader(reader, 0);
429 Reconstruct(rawReader);
432 //____________________________________________________________________
434 AliFMDReconstructor::Digitize(AliRawReader* reader, TClonesArray* sdigits) const
436 // Reconstruct directly from raw data (no intermediate output on
437 // digit tree or rec point tree).
440 // reader Raw event reader
443 AliWarning("I'm a zombie - cannot do anything");
446 AliFMDRawReader rawReader(reader, 0);
448 UShort_t det, sec, str, sam, rat, fac;
449 Short_t adc, oldDet = -1;
454 while (rawReader.NextSample(det, rng, sec, str, sam, rat, adc, zs, fac)) {
455 if (!rawReader.SelectSample(sam, rat)) continue;
458 fZSFactor[det-1] = fac;
461 DigitizeSignal(sdigits, det, rng, sec, str, sam, adc);
463 UseRecoParam(kFALSE);
466 //____________________________________________________________________
468 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
469 TTree* clusterTree) const
471 // Reconstruct event from digits in tree
472 // Get the FMD branch holding the digits.
473 // FIXME: The vertex may not be known yet, so we may have to move
474 // some of this to FillESD.
477 // digitsTree Pointer to a tree containing digits
478 // clusterTree Pointer to output tree
480 if (!PreReconstruct()) return;
482 if (!fMult) fMult = new TClonesArray("AliFMDRecPoint");
484 AliFMDDebug(1, ("Reconstructing from digits in a tree"));
486 // Get the digitis array
487 static TClonesArray* digits = new TClonesArray("AliFMDDigit");
488 TBranch* digitBranch = digitsTree->GetBranch("FMD");
490 Error("Exec", "No digit branch for the FMD found");
493 digitBranch->SetAddress(&digits);
495 if (digits) digits->Clear();
496 if (fMult) fMult->Clear();
498 // Create rec-point output branch
500 fTreeR = clusterTree;
501 fTreeR->Branch("FMD", &fMult);
503 AliFMDDebug(5, ("Getting entry 0 from digit branch"));
504 digitBranch->GetEntry(0);
506 AliFMDDebug(5, ("Processing digits"));
508 ProcessDigits(digits);
509 UseRecoParam(kFALSE);
511 Int_t written = clusterTree->Fill();
512 AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
518 //____________________________________________________________________
520 AliFMDReconstructor::ProcessDigits(TClonesArray* digits,
521 const AliFMDRawReader& rawRead) const
523 // For each digit, find the pseudo rapdity, azimuthal angle, and
524 // number of corrected ADC counts, and pass it on to the algorithms
528 // digits Array of digits
531 AliWarning("I'm a zombie - cannot do anything");
534 AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
535 for (size_t i = 1; i <= 3; i++) {
536 fZS[i-1] = rawRead.IsZeroSuppressed(map->Detector2DDL(i));
537 fZSFactor[i-1] = rawRead.NoiseFactor(map->Detector2DDL(i));
540 ProcessDigits(digits);
541 UseRecoParam(kFALSE);
544 //____________________________________________________________________
546 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
548 // For each digit, find the pseudo rapdity, azimuthal angle, and
549 // number of corrected ADC counts, and pass it on to the algorithms
553 // digits Array of digits
555 Int_t nDigits = digits->GetEntries();
556 AliFMDDebug(2, ("Got %d digits", nDigits));
557 fESDObj->SetNoiseFactor(fNoiseFactor);
558 fESDObj->SetAngleCorrected(fAngleCorrect);
560 for (Int_t i = 0; i < nDigits; i++) {
561 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
562 if (!digit) continue;
567 //____________________________________________________________________
569 AliFMDReconstructor::ProcessDigit(AliFMDDigit* digit) const
572 // Process a single digit
575 // digit Digiti to process
577 UShort_t det = digit->Detector();
578 Char_t rng = digit->Ring();
579 UShort_t sec = digit->Sector();
580 UShort_t str = digit->Strip();
581 Short_t adc = digit->Counts();
583 ProcessSignal(det, rng, sec, str, adc);
586 //____________________________________________________________________
588 AliFMDReconstructor::ProcessSignal(UShort_t det,
594 // Process the signal from a single strip
603 if (adc >= AliFMDRawReader::kBadSignal) {
604 AliFMDDebug(3, ("FMD%d%c[%2d,%3d] is marked bad", det, rng, sec, str));
605 fBad(det,rng,sec,str) = true;
609 // Check that the strip is not marked as dead
610 AliFMDParameters* param = AliFMDParameters::Instance();
611 if (param->IsDead(det, rng, sec, str)) {
612 AliFMDDebug(3, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
613 fBad(det,rng,sec,str) = true;
619 Float_t eta = fESDObj->Eta(det, rng, 0, str);
621 // Substract pedestal.
622 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
623 if(counts == USHRT_MAX) return;
625 // Gain match digits.
626 Double_t edep = Adc2Energy(det, rng, sec, str, eta, counts);
627 // Get rid of nonsense energy
630 // Make rough multiplicity
631 Double_t mult = Energy2Multiplicity(det, rng, sec, str, edep);
632 // Get rid of nonsense mult
634 // AliWarning(Form("The mutliplicity in FMD%d%c[%2d,%3d]=%f > 20 "
635 // "(ADC: %d, Energy: %f)", det, rng, sec, str, mult,
638 if (mult < 0) return;
639 AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
640 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
641 det, rng, sec, str, adc, counts, edep, mult));
643 // Create a `RecPoint' on the output branch.
646 PhysicalCoordinates(det, rng, sec, str, eta, phi);
649 new ((*fMult)[fNMult]) AliFMDRecPoint(det, rng, sec, str,
650 eta, phi, edep, mult);
651 (void)m; // Suppress warnings about unused variables.
655 fESDObj->SetMultiplicity(det, rng, sec, str, mult);
656 // fESDObj->SetEta(det, rng, sec, str, eta);
658 if (fDiagAll) fDiagAll->Fill(adc, mult);
662 //____________________________________________________________________
664 AliFMDReconstructor::DigitizeSignal(TClonesArray* sdigits,
672 // Process the signal from a single strip
681 AliFMDParameters* param = AliFMDParameters::Instance();
682 // Check that the strip is not marked as dead
683 if (param->IsDead(det, rng, sec, str)) {
684 AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
688 // Substract pedestal.
689 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
690 if(counts == USHRT_MAX || counts == 0) return;
692 // Gain match digits.
693 Double_t edep = Adc2Energy(det, rng, sec, str, counts);
694 // Get rid of nonsense energy
697 Int_t n = sdigits->GetEntriesFast();
698 // AliFMDSDigit* sdigit =
700 AliFMDSDigit(det, rng, sec, str, edep, counts, counts, counts, counts);
701 // sdigit->SetCount(sam, counts);
704 //____________________________________________________________________
706 AliFMDReconstructor::SubtractPedestal(UShort_t det,
713 UShort_t zsNoiseFactor) const
716 // Subtract the pedestal off the ADC counts.
719 // det Detector number
720 // rng Ring identifier
724 // noiseFactor If pedestal substracted pedestal is less then
725 // this times the noise, then consider this to be 0.
726 // zsEnabled Whether zero-suppression is on.
727 // zsNoiseFactor Noise factor used in on-line pedestal
731 // The pedestal subtracted ADC counts (possibly 0), or @c
732 // USHRT_MAX in case of problems.
734 AliFMDParameters* param = AliFMDParameters::Instance();
735 Float_t ped = (zsEnabled ? 0 :
736 param->GetPedestal(det, rng, sec, str));
737 Float_t noise = param->GetPedestalWidth(det, rng, sec, str);
738 if(ped < 0 || noise < 0) {
739 AliWarningClass(Form("Invalid pedestal (%f) or noise (%f) "
740 "for FMD%d%c[%02d,%03d]",
741 ped, noise, det, rng, sec, str));
744 AliDebugClass(15, Form("Subtracting pedestal for FMD%d%c[%2d,%3d]=%4d "
745 "(%s w/factor %d, noise factor %f, "
746 "pedestal %8.2f+/-%8.2f)",
747 det, rng, sec, str, adc,
748 (zsEnabled ? "zs'ed" : "straight"),
749 zsNoiseFactor, noiseFactor, ped, noise));
751 Int_t counts = adc + Int_t(zsEnabled ? zsNoiseFactor * noise : - ped);
752 counts = TMath::Max(Int_t(counts), 0);
753 // Calculate the noise factor for suppressing remenants of the noise
754 // peak. If we have done on-line zero suppression, we only check
755 // for noise signals that are larger than the suppressed noise. If
756 // the noise factor used on line is larger than the factor used
757 // here, we do not do this check at all.
760 // Online factor | Read factor | Result
761 // ---------------+--------------+-------------------------------
762 // 2 | 3 | Check if signal > 1 * noise
763 // 3 | 3 | Check if signal > 0
764 // 3 | 2 | Check if signal > 0
766 // In this way, we make sure that we do not suppress away too much
767 // data, and that the read-factor is the most stringent cut.
768 Float_t nf = TMath::Max(0.F, noiseFactor - (zsEnabled ? zsNoiseFactor : 0));
769 if (counts < noise * nf) counts = 0;
770 if (counts > 0) AliDebugClass(15, "Got a hit strip");
772 UShort_t ret = counts < 0 ? 0 : counts;
777 //____________________________________________________________________
779 AliFMDReconstructor::SubtractPedestal(UShort_t det,
785 // Member function to subtract the pedestal from a digit
792 // adc # of ADC counts
794 // Pedestal subtracted signal or USHRT_MAX in case of problems
796 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc,
797 fNoiseFactor, fZS[det-1],
799 if (fDiagStep1) fDiagStep1->Fill(adc, counts);
804 //____________________________________________________________________
806 AliFMDReconstructor::Adc2Energy(UShort_t det,
810 UShort_t count) const
812 // Converts number of ADC counts to energy deposited.
813 // Note, that this member function can be overloaded by derived
814 // classes to do strip-specific look-ups in databases or the like,
815 // to find the proper gain for a strip.
817 // In the first simple version, we calculate the energy deposited as
819 // EnergyDeposited = cos(theta) * gain * count
824 // gain = ----------------- * Energy_deposited_per_MIP
827 // is constant and the same for all strips.
829 // For the production we use the conversion measured in the NBI lab.
830 // The total conversion is then:
835 // => energy = ----------------
843 // counts Number of ADC counts over pedestal
845 // The energy deposited in a single strip, or -1 in case of problems
847 if (count <= 0) return 0;
848 AliFMDParameters* param = AliFMDParameters::Instance();
849 Float_t gain = param->GetPulseGain(det, rng, sec, str);
850 // 'Tagging' bad gains as bad energy
852 AliWarning(Form("Invalid gain (%f) for FMD%d%c[%02d,%03d]",
853 gain, det, rng, sec, str));
856 AliFMDDebug(5, ("Converting counts %d to energy (factor=%f, DAC2MIP=%f)",
857 count, gain,param->GetDACPerMIP()));
859 Double_t edep = ((count * param->GetEdepMip())
860 / (gain * param->GetDACPerMIP()));
864 //____________________________________________________________________
866 AliFMDReconstructor::Adc2Energy(UShort_t det,
871 UShort_t count) const
873 // Converts number of ADC counts to energy deposited.
874 // Note, that this member function can be overloaded by derived
875 // classes to do strip-specific look-ups in databases or the like,
876 // to find the proper gain for a strip.
878 // In the first simple version, we calculate the energy deposited as
880 // EnergyDeposited = cos(theta) * gain * count
885 // gain = ----------------- * Energy_deposited_per_MIP
888 // is constant and the same for all strips.
890 // For the production we use the conversion measured in the NBI lab.
891 // The total conversion is then:
896 // => energy = ----------------
904 // eta Psuedo-rapidity
905 // counts Number of ADC counts over pedestal
907 // The energy deposited in a single strip, or -1 in case of problems
909 Double_t edep = Adc2Energy(det, rng, sec, str, count);
911 if (fDiagStep2) fDiagStep2->Fill(count, edep);
913 Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
914 Double_t corr = TMath::Abs(TMath::Cos(theta));
915 Double_t cedep = corr * edep;
916 AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
917 edep, corr, cedep, eta, theta));
918 if (fDiagStep3) fDiagStep3->Fill(edep, cedep);
924 //____________________________________________________________________
926 AliFMDReconstructor::Energy2Multiplicity(UShort_t /*det*/,
932 // Converts an energy signal to number of particles.
933 // Note, that this member function can be overloaded by derived
934 // classes to do strip-specific look-ups in databases or the like,
935 // to find the proper gain for a strip.
937 // In this simple version, we calculate the multiplicity as
939 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
943 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
945 // is constant and the same for all strips
952 // edep Energy deposited in a single strip
954 // The "bare" multiplicity corresponding to the energy deposited
955 AliFMDParameters* param = AliFMDParameters::Instance();
956 Double_t edepMIP = param->GetEdepMip();
957 Float_t mult = edep / edepMIP;
960 AliFMDDebug(15, ("Translating energy %f to multiplicity via "
961 "divider %f->%f", edep, edepMIP, mult));
963 if (fDiagStep4) fDiagStep4->Fill(edep, mult);
967 //____________________________________________________________________
969 AliFMDReconstructor::PhysicalCoordinates(UShort_t det,
976 // Get the eta and phi of a digit
983 // eta On return, contains the psuedo-rapidity of the strip
984 // phi On return, contains the azimuthal angle of the strip
986 AliFMDGeometry* geom = AliFMDGeometry::Instance();
987 Double_t x, y, z, r, theta, deta, dphi;
988 geom->Detector2XYZ(det, rng, sec, str, x, y, z);
990 // Correct for vertex offset.
992 AliFMDGeometry::XYZ2REtaPhiTheta(x, y, z, r, deta, dphi, theta);
998 class ESDPrinter : public AliESDFMD::ForOne
1002 Bool_t operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t,
1003 Float_t m, Float_t e)
1005 if (m > 0 && m != AliESDFMD::kInvalidMult)
1006 printf(" FMD%d%c[%2d,%3d] = %6.3f / %6.3f\n", d, r, s, t, m, e);
1013 //____________________________________________________________________
1015 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
1016 TTree* /* clusterTree */,
1017 AliESDEvent* esd) const
1019 // nothing to be done
1020 // FIXME: The vertex may not be known when Reconstruct is executed,
1021 // so we may have to move some of that member function here.
1022 AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
1024 AliWarning("I'm a zombie - cannot do anything");
1027 // fESDObj->Print();
1029 // Fix up ESD so that only truely dead channels get the kInvalidMult flag.
1030 MarkDeadChannels(fESDObj);
1032 Double_t oldVz = fCurrentVertex;
1034 if (fVertexType != kNoVertex) {
1035 AliFMDDebug(2, ("Revertexing the ESD data to vz=%f (was %f)",
1036 fCurrentVertex, oldVz));
1037 AliFMDESDRevertexer revertexer;
1038 revertexer.Revertex(fESDObj, fCurrentVertex);
1041 if (AliDebugLevel() > 10) {
1043 fESDObj->ForEach(p);
1047 AliFMDDebug(2, ("Writing FMD data to ESD tree"));
1048 esd->SetFMDData(fESDObj);
1051 if (!fDiagnostics || !esd) return;
1052 static bool first = true;
1053 // This is most likely NOT the event number you'd like to use. It
1054 // has nothing to do with the 'real' event number.
1055 // - That's OK. We just use it for the name of the directory -
1056 // nothing else. Christian
1057 Int_t evno = esd->GetEventNumberInFile();
1058 AliFMDDebug(3, ("Writing diagnostics histograms to FMD.Diag.root/%03d",evno));
1059 TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
1062 TDirectory* d = f.mkdir(Form("%03d", evno),
1063 Form("Diagnostics histograms for event # %d", evno));
1065 if (fDiagStep1) fDiagStep1->Write();
1066 if (fDiagStep2) fDiagStep2->Write();
1067 if (fDiagStep3) fDiagStep3->Write();
1068 if (fDiagStep4) fDiagStep4->Write();
1069 if (fDiagAll) fDiagAll->Write();
1074 if (fDiagStep1) fDiagStep1->Reset();
1075 if (fDiagStep2) fDiagStep2->Reset();
1076 if (fDiagStep3) fDiagStep3->Reset();
1077 if (fDiagStep4) fDiagStep4->Reset();
1078 if (fDiagAll) fDiagAll->Reset();
1081 //____________________________________________________________________
1083 AliFMDReconstructor::FillESD(AliRawReader*, TTree* clusterTree,
1084 AliESDEvent* esd) const
1087 // Forwards to above member function
1090 AliWarning("I'm a zombie - cannot do anything");
1094 FillESD(dummy, clusterTree, esd);
1096 //____________________________________________________________________