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 if (fBad(d, r, s, t)) {
316 AliDebug(5, Form("Marking FMD%d%c[%2d,%3d] as bad", d, r, s, t));
317 esd->SetMultiplicity(d, r, s, t, AliESDFMD::kInvalidMult);
319 if (param->IsDead(d, r, s, t)) {
320 AliDebug(5, Form("Marking FMD%d%c[%2d,%3d] as dead", d, r, s, t));
321 esd->SetMultiplicity(d, r, s, t, AliESDFMD::kInvalidMult);
322 // esd->SetEta(d, r, s, t, AliESDFMD::kInvalidEta);
324 else if (esd->Multiplicity(d, r, s, t) == AliESDFMD::kInvalidMult) {
325 AliDebug(10, Form("Setting null signal in FMD%d%c[%2d,%3d]",
327 esd->SetMultiplicity(d, r, s, t, 0);
335 //____________________________________________________________________
337 AliFMDReconstructor::Reconstruct(AliFMDRawReader& rawReader) const
339 AliFMDDebug(1, ("Reconstructing from FMD raw reader"));
341 AliWarning("I'm a zombie - cannot do anything");
345 UShort_t det, sec, str, fac;
346 Short_t adc, oldDet = -1;
351 while (rawReader.NextSignal(det, rng, sec, str, adc, zs, fac)) {
354 fZSFactor[det-1] = fac;
357 ProcessSignal(det, rng, sec, str, adc);
359 UseRecoParam(kFALSE);
363 //____________________________________________________________________
365 AliFMDReconstructor::Reconstruct(AliRawReader* reader, TTree*) const
367 // Reconstruct directly from raw data (no intermediate output on
368 // digit tree or rec point tree).
371 // reader Raw event reader
372 // ctree Not used - 'cluster tree' to store rec-points in.
373 AliFMDDebug(1, ("Reconstructing from raw reader"));
375 AliWarning("I'm a zombie - cannot do anything");
378 AliFMDRawReader rawReader(reader, 0);
379 Reconstruct(rawReader);
382 //____________________________________________________________________
384 AliFMDReconstructor::Digitize(AliRawReader* reader, TClonesArray* sdigits) const
386 // Reconstruct directly from raw data (no intermediate output on
387 // digit tree or rec point tree).
390 // reader Raw event reader
393 AliWarning("I'm a zombie - cannot do anything");
396 AliFMDRawReader rawReader(reader, 0);
398 UShort_t det, sec, str, sam, rat, fac;
399 Short_t adc, oldDet = -1;
404 while (rawReader.NextSample(det, rng, sec, str, sam, rat, adc, zs, fac)) {
405 if (!rawReader.SelectSample(sam, rat)) continue;
408 fZSFactor[det-1] = fac;
411 DigitizeSignal(sdigits, det, rng, sec, str, sam, adc);
413 UseRecoParam(kFALSE);
416 //____________________________________________________________________
418 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
419 TTree* clusterTree) const
421 // Reconstruct event from digits in tree
422 // Get the FMD branch holding the digits.
423 // FIXME: The vertex may not be known yet, so we may have to move
424 // some of this to FillESD.
427 // digitsTree Pointer to a tree containing digits
428 // clusterTree Pointer to output tree
431 AliWarning("I'm a zombie - cannot do anything");
434 if (!fMult) fMult = new TClonesArray("AliFMDRecPoint");
436 AliFMDDebug(1, ("Reconstructing from digits in a tree"));
441 static TClonesArray* digits = new TClonesArray("AliFMDDigit");
442 TBranch* digitBranch = digitsTree->GetBranch("FMD");
444 Error("Exec", "No digit branch for the FMD found");
447 digitBranch->SetAddress(&digits);
449 if (digits) digits->Clear();
450 if (fMult) fMult->Clear();
451 if (fESDObj) fESDObj->Clear();
454 fTreeR = clusterTree;
455 fTreeR->Branch("FMD", &fMult);
457 AliFMDDebug(5, ("Getting entry 0 from digit branch"));
458 digitBranch->GetEntry(0);
460 AliFMDDebug(5, ("Processing digits"));
462 ProcessDigits(digits);
463 UseRecoParam(kFALSE);
465 Int_t written = clusterTree->Fill();
466 AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
472 //____________________________________________________________________
474 AliFMDReconstructor::ProcessDigits(TClonesArray* digits,
475 const AliFMDRawReader& rawRead) const
477 // For each digit, find the pseudo rapdity, azimuthal angle, and
478 // number of corrected ADC counts, and pass it on to the algorithms
482 // digits Array of digits
485 AliWarning("I'm a zombie - cannot do anything");
488 AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
489 for (size_t i = 1; i <= 3; i++) {
490 fZS[i-1] = rawRead.IsZeroSuppressed(map->Detector2DDL(i));
491 fZSFactor[i-1] = rawRead.NoiseFactor(map->Detector2DDL(i));
494 ProcessDigits(digits);
495 UseRecoParam(kFALSE);
498 //____________________________________________________________________
500 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
502 // For each digit, find the pseudo rapdity, azimuthal angle, and
503 // number of corrected ADC counts, and pass it on to the algorithms
507 // digits Array of digits
509 Int_t nDigits = digits->GetEntries();
510 AliFMDDebug(2, ("Got %d digits", nDigits));
511 fESDObj->SetNoiseFactor(fNoiseFactor);
512 fESDObj->SetAngleCorrected(fAngleCorrect);
514 for (Int_t i = 0; i < nDigits; i++) {
515 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
516 if (!digit) continue;
521 //____________________________________________________________________
523 AliFMDReconstructor::ProcessDigit(AliFMDDigit* digit) const
526 // Process a single digit
529 // digit Digiti to process
531 UShort_t det = digit->Detector();
532 Char_t rng = digit->Ring();
533 UShort_t sec = digit->Sector();
534 UShort_t str = digit->Strip();
535 Short_t adc = digit->Counts();
537 ProcessSignal(det, rng, sec, str, adc);
540 //____________________________________________________________________
542 AliFMDReconstructor::ProcessSignal(UShort_t det,
548 // Process the signal from a single strip
557 if (adc >= AliFMDRawReader::kBadSignal) {
558 AliFMDDebug(3, ("FMD%d%c[%2d,%3d] is marked bad", det, rng, sec, str));
559 fBad(det,rng,sec,str) = true;
563 // Check that the strip is not marked as dead
564 AliFMDParameters* param = AliFMDParameters::Instance();
565 if (param->IsDead(det, rng, sec, str)) {
566 AliFMDDebug(3, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
567 fBad(det,rng,sec,str) = true;
574 PhysicalCoordinates(det, rng, sec, str, eta, phi);
576 // Substract pedestal.
577 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
578 if(counts == USHRT_MAX) return;
580 // Gain match digits.
581 Double_t edep = Adc2Energy(det, rng, sec, str, eta, counts);
582 // Get rid of nonsense energy
585 // Make rough multiplicity
586 Double_t mult = Energy2Multiplicity(det, rng, sec, str, edep);
587 // Get rid of nonsense mult
589 // AliWarning(Form("The mutliplicity in FMD%d%c[%2d,%3d]=%f > 20 "
590 // "(ADC: %d, Energy: %f)", det, rng, sec, str, mult,
593 if (mult < 0) return;
594 AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
595 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
596 det, rng, sec, str, adc, counts, edep, mult));
598 // Create a `RecPoint' on the output branch.
601 new ((*fMult)[fNMult]) AliFMDRecPoint(det, rng, sec, str,
602 eta, phi, edep, mult);
603 (void)m; // Suppress warnings about unused variables.
607 fESDObj->SetMultiplicity(det, rng, sec, str, mult);
608 fESDObj->SetEta(det, rng, sec, str, eta);
610 if (fDiagAll) fDiagAll->Fill(adc, mult);
614 //____________________________________________________________________
616 AliFMDReconstructor::DigitizeSignal(TClonesArray* sdigits,
624 // Process the signal from a single strip
633 AliFMDParameters* param = AliFMDParameters::Instance();
634 // Check that the strip is not marked as dead
635 if (param->IsDead(det, rng, sec, str)) {
636 AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
640 // Substract pedestal.
641 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
642 if(counts == USHRT_MAX || counts == 0) return;
644 // Gain match digits.
645 Double_t edep = Adc2Energy(det, rng, sec, str, counts);
646 // Get rid of nonsense energy
649 Int_t n = sdigits->GetEntriesFast();
650 // AliFMDSDigit* sdigit =
652 AliFMDSDigit(det, rng, sec, str, edep, counts, counts, counts, counts);
653 // sdigit->SetCount(sam, counts);
656 //____________________________________________________________________
658 AliFMDReconstructor::SubtractPedestal(UShort_t det,
665 UShort_t zsNoiseFactor) const
668 // Subtract the pedestal off the ADC counts.
671 // det Detector number
672 // rng Ring identifier
676 // noiseFactor If pedestal substracted pedestal is less then
677 // this times the noise, then consider this to be 0.
678 // zsEnabled Whether zero-suppression is on.
679 // zsNoiseFactor Noise factor used in on-line pedestal
683 // The pedestal subtracted ADC counts (possibly 0), or @c
684 // USHRT_MAX in case of problems.
686 AliFMDParameters* param = AliFMDParameters::Instance();
687 Float_t ped = (zsEnabled ? 0 :
688 param->GetPedestal(det, rng, sec, str));
689 Float_t noise = param->GetPedestalWidth(det, rng, sec, str);
690 if(ped < 0 || noise < 0) {
691 AliWarningClass(Form("Invalid pedestal (%f) or noise (%f) "
692 "for FMD%d%c[%02d,%03d]",
693 ped, noise, det, rng, sec, str));
696 AliDebugClass(15, Form("Subtracting pedestal for FMD%d%c[%2d,%3d]=%4d "
697 "(%s w/factor %d, noise factor %f, "
698 "pedestal %8.2f+/-%8.2f)",
699 det, rng, sec, str, adc,
700 (zsEnabled ? "zs'ed" : "straight"),
701 zsNoiseFactor, noiseFactor, ped, noise));
703 Int_t counts = adc + Int_t(zsEnabled ? zsNoiseFactor * noise : - ped);
704 counts = TMath::Max(Int_t(counts), 0);
705 // Calculate the noise factor for suppressing remenants of the noise
706 // peak. If we have done on-line zero suppression, we only check
707 // for noise signals that are larger than the suppressed noise. If
708 // the noise factor used on line is larger than the factor used
709 // here, we do not do this check at all.
712 // Online factor | Read factor | Result
713 // ---------------+--------------+-------------------------------
714 // 2 | 3 | Check if signal > 1 * noise
715 // 3 | 3 | Check if signal > 0
716 // 3 | 2 | Check if signal > 0
718 // In this way, we make sure that we do not suppress away too much
719 // data, and that the read-factor is the most stringent cut.
720 Float_t nf = TMath::Max(0.F, noiseFactor - (zsEnabled ? zsNoiseFactor : 0));
721 if (counts < noise * nf) counts = 0;
722 if (counts > 0) AliDebugClass(15, "Got a hit strip");
724 UShort_t ret = counts < 0 ? 0 : counts;
729 //____________________________________________________________________
731 AliFMDReconstructor::SubtractPedestal(UShort_t det,
737 // Member function to subtract the pedestal from a digit
744 // adc # of ADC counts
746 // Pedestal subtracted signal or USHRT_MAX in case of problems
748 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc,
749 fNoiseFactor, fZS[det-1],
751 if (fDiagStep1) fDiagStep1->Fill(adc, counts);
756 //____________________________________________________________________
758 AliFMDReconstructor::Adc2Energy(UShort_t det,
762 UShort_t count) const
764 // Converts number of ADC counts to energy deposited.
765 // Note, that this member function can be overloaded by derived
766 // classes to do strip-specific look-ups in databases or the like,
767 // to find the proper gain for a strip.
769 // In the first simple version, we calculate the energy deposited as
771 // EnergyDeposited = cos(theta) * gain * count
776 // gain = ----------------- * Energy_deposited_per_MIP
779 // is constant and the same for all strips.
781 // For the production we use the conversion measured in the NBI lab.
782 // The total conversion is then:
787 // => energy = ----------------
795 // counts Number of ADC counts over pedestal
797 // The energy deposited in a single strip, or -1 in case of problems
799 if (count <= 0) return 0;
800 AliFMDParameters* param = AliFMDParameters::Instance();
801 Float_t gain = param->GetPulseGain(det, rng, sec, str);
802 // 'Tagging' bad gains as bad energy
804 AliWarning(Form("Invalid gain (%f) for FMD%d%c[%02d,%03d]",
805 gain, det, rng, sec, str));
808 AliFMDDebug(5, ("Converting counts %d to energy (factor=%f, DAC2MIP=%f)",
809 count, gain,param->GetDACPerMIP()));
811 Double_t edep = ((count * param->GetEdepMip())
812 / (gain * param->GetDACPerMIP()));
816 //____________________________________________________________________
818 AliFMDReconstructor::Adc2Energy(UShort_t det,
823 UShort_t count) const
825 // Converts number of ADC counts to energy deposited.
826 // Note, that this member function can be overloaded by derived
827 // classes to do strip-specific look-ups in databases or the like,
828 // to find the proper gain for a strip.
830 // In the first simple version, we calculate the energy deposited as
832 // EnergyDeposited = cos(theta) * gain * count
837 // gain = ----------------- * Energy_deposited_per_MIP
840 // is constant and the same for all strips.
842 // For the production we use the conversion measured in the NBI lab.
843 // The total conversion is then:
848 // => energy = ----------------
856 // eta Psuedo-rapidity
857 // counts Number of ADC counts over pedestal
859 // The energy deposited in a single strip, or -1 in case of problems
861 Double_t edep = Adc2Energy(det, rng, sec, str, count);
863 if (fDiagStep2) fDiagStep2->Fill(count, edep);
865 Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
866 Double_t corr = TMath::Abs(TMath::Cos(theta));
867 Double_t cedep = corr * edep;
868 AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
869 edep, corr, cedep, eta, theta));
870 if (fDiagStep3) fDiagStep3->Fill(edep, cedep);
876 //____________________________________________________________________
878 AliFMDReconstructor::Energy2Multiplicity(UShort_t /*det*/,
884 // Converts an energy signal to number of particles.
885 // Note, that this member function can be overloaded by derived
886 // classes to do strip-specific look-ups in databases or the like,
887 // to find the proper gain for a strip.
889 // In this simple version, we calculate the multiplicity as
891 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
895 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
897 // is constant and the same for all strips
904 // edep Energy deposited in a single strip
906 // The "bare" multiplicity corresponding to the energy deposited
907 AliFMDParameters* param = AliFMDParameters::Instance();
908 Double_t edepMIP = param->GetEdepMip();
909 Float_t mult = edep / edepMIP;
912 AliFMDDebug(15, ("Translating energy %f to multiplicity via "
913 "divider %f->%f", edep, edepMIP, mult));
915 if (fDiagStep4) fDiagStep4->Fill(edep, mult);
919 //____________________________________________________________________
921 AliFMDReconstructor::PhysicalCoordinates(UShort_t det,
928 // Get the eta and phi of a digit
935 // eta On return, contains the psuedo-rapidity of the strip
936 // phi On return, contains the azimuthal angle of the strip
938 AliFMDGeometry* geom = AliFMDGeometry::Instance();
939 Double_t x, y, z, r, theta, deta, dphi;
940 geom->Detector2XYZ(det, rng, sec, str, x, y, z);
942 // Correct for vertex offset.
944 AliFMDGeometry::XYZ2REtaPhiTheta(x, y, z, r, deta, dphi, theta);
950 class ESDPrinter : public AliESDFMD::ForOne
954 Bool_t operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t,
955 Float_t m, Float_t e)
957 if (m > 0 && m != AliESDFMD::kInvalidMult)
958 printf(" FMD%d%c[%2d,%3d] = %6.3f / %6.3f\n", d, r, s, t, m, e);
965 //____________________________________________________________________
967 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
968 TTree* /* clusterTree */,
969 AliESDEvent* esd) const
971 // nothing to be done
972 // FIXME: The vertex may not be known when Reconstruct is executed,
973 // so we may have to move some of that member function here.
974 AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
976 AliWarning("I'm a zombie - cannot do anything");
981 // Fix up ESD so that only truely dead channels get the kInvalidMult flag.
982 MarkDeadChannels(fESDObj);
984 Double_t oldVz = fCurrentVertex;
986 if (fVertexType != kNoVertex) {
987 AliFMDDebug(2, ("Revertexing the ESD data to vz=%f (was %f)",
988 fCurrentVertex, oldVz));
989 AliFMDESDRevertexer revertexer;
990 revertexer.Revertex(fESDObj, fCurrentVertex);
993 if (AliDebugLevel() > 10) {
999 AliFMDDebug(2, ("Writing FMD data to ESD tree"));
1000 esd->SetFMDData(fESDObj);
1003 if (!fDiagnostics || !esd) return;
1004 static bool first = true;
1005 // This is most likely NOT the event number you'd like to use. It
1006 // has nothing to do with the 'real' event number.
1007 // - That's OK. We just use it for the name of the directory -
1008 // nothing else. Christian
1009 Int_t evno = esd->GetEventNumberInFile();
1010 AliFMDDebug(3, ("Writing diagnostics histograms to FMD.Diag.root/%03d",evno));
1011 TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
1014 TDirectory* d = f.mkdir(Form("%03d", evno),
1015 Form("Diagnostics histograms for event # %d", evno));
1017 if (fDiagStep1) fDiagStep1->Write();
1018 if (fDiagStep2) fDiagStep2->Write();
1019 if (fDiagStep3) fDiagStep3->Write();
1020 if (fDiagStep4) fDiagStep4->Write();
1021 if (fDiagAll) fDiagAll->Write();
1026 if (fDiagStep1) fDiagStep1->Reset();
1027 if (fDiagStep2) fDiagStep2->Reset();
1028 if (fDiagStep3) fDiagStep3->Reset();
1029 if (fDiagStep4) fDiagStep4->Reset();
1030 if (fDiagAll) fDiagAll->Reset();
1033 //____________________________________________________________________
1035 AliFMDReconstructor::FillESD(AliRawReader*, TTree* clusterTree,
1036 AliESDEvent* esd) const
1039 // Forwards to above member function
1042 AliWarning("I'm a zombie - cannot do anything");
1046 FillESD(dummy, clusterTree, esd);
1048 //____________________________________________________________________