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"));
166 AliError("No digits tree passed");
169 static TClonesArray* array = new TClonesArray("AliFMDDigit");
170 digitsTree->Branch("FMD", &array);
172 AliFMDRawReader rawRead(reader, digitsTree);
173 // rawRead.SetSampleRate(fFMD->GetSampleRate());
175 rawRead.ReadAdcs(array);
177 Int_t nWrite = digitsTree->Fill();
178 AliFMDDebug(1, ("Got a grand total of %d digits, wrote %d bytes to tree",
179 array->GetEntriesFast(), nWrite));
182 AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
183 for (size_t i = 1; i <= 3; i++) {
184 fZS[i-1] = rawRead.IsZeroSuppressed(map->Detector2DDL(i));
185 fZSFactor[i-1] = rawRead.NoiseFactor(map->Detector2DDL(i));
189 //____________________________________________________________________
191 AliFMDReconstructor::GetVertex(AliESDEvent* esd) const
193 // Return the vertex to use.
194 // This is obtained from the ESD object.
195 // If not found, a warning is issued.
196 fVertexType = kNoVertex;
200 const AliESDVertex* vertex = esd->GetPrimaryVertex();
201 if (!vertex) vertex = esd->GetPrimaryVertexSPD();
202 if (!vertex) vertex = esd->GetPrimaryVertexTPC();
203 if (!vertex) vertex = esd->GetVertex();
206 AliFMDDebug(2, ("Got %s (%s) from ESD: %f",
207 vertex->GetName(), vertex->GetTitle(), vertex->GetZv()));
208 fCurrentVertex = vertex->GetZv();
209 fVertexType = kESDVertex;
212 else if (esd->GetESDTZERO()) {
213 AliFMDDebug(2, ("Got primary vertex from T0: %f", esd->GetT0zVertex()));
214 fCurrentVertex = esd->GetT0zVertex();
215 fVertexType = kESDVertex;
218 AliWarning("Didn't get any vertex from ESD or generator");
221 //____________________________________________________________________
223 AliFMDReconstructor::GetIdentifier() const
225 return AliReconstruction::GetDetIndex(GetDetectorName());
228 //____________________________________________________________________
229 const AliFMDRecoParam*
230 AliFMDReconstructor::GetParameters() const
232 Int_t iDet = 12; // GetIdentifier();
233 const AliDetectorRecoParam* params = AliReconstructor::GetRecoParam(iDet);
234 if (!params || params->IsA() != AliFMDRecoParam::Class()) return 0;
235 return static_cast<const AliFMDRecoParam*>(params);
238 //____________________________________________________________________
240 AliFMDReconstructor::UseRecoParam(Bool_t set) const
242 static Float_t savedNoiseFactor = fNoiseFactor;
243 static Bool_t savedAngleCorrect = fAngleCorrect;
245 const AliFMDRecoParam* params = GetParameters();
247 fNoiseFactor = params->NoiseFactor();
248 fAngleCorrect = params->AngleCorrect();
252 fNoiseFactor = savedNoiseFactor;
253 fAngleCorrect = savedAngleCorrect;
257 //____________________________________________________________________
259 AliFMDReconstructor::MarkDeadChannels(AliESDFMD* esd) const
261 // Loop over all entries of the ESD and mark
262 // those that are dead as such
263 // - otherwise put in the zero signal.
264 AliFMDParameters* param = AliFMDParameters::Instance();
266 for (UShort_t d = 1; d <= 3; d++) {
267 UShort_t nR = (d == 1 ? 1 : 2);
269 for (UShort_t q = 0; q < nR; q++) {
270 Char_t r = (q == 0 ? 'I' : 'O');
271 UShort_t nS = (q == 0 ? 20 : 40);
272 UShort_t nT = (q == 0 ? 512 : 256);
274 for (UShort_t s = 0; s < nS; s++) {
275 for (UShort_t t = 0; t < nT; t++) {
276 if (param->IsDead(d, r, s, t))
277 esd->SetMultiplicity(d, r, s, t, AliESDFMD::kInvalidMult);
278 else if (esd->Multiplicity(d, r, s, t) == AliESDFMD::kInvalidMult)
279 esd->SetMultiplicity(d, r, s, t, 0);
286 //____________________________________________________________________
288 AliFMDReconstructor::Reconstruct(AliRawReader* reader, TTree*) const
290 // Reconstruct directly from raw data (no intermediate output on
291 // digit tree or rec point tree).
294 // reader Raw event reader
295 // ctree Not used - 'cluster tree' to store rec-points in.
296 AliFMDDebug(1, ("Reconstructing from raw reader"));
297 AliFMDRawReader rawReader(reader, 0);
299 UShort_t det, sec, str, fac;
300 Short_t adc, oldDet = -1;
305 while (rawReader.NextSignal(det, rng, sec, str, adc, zs, fac)) {
308 fZSFactor[det-1] = fac;
311 ProcessSignal(det, rng, sec, str, adc);
313 UseRecoParam(kFALSE);
317 //____________________________________________________________________
319 AliFMDReconstructor::Digitize(AliRawReader* reader, TClonesArray* sdigits) const
321 // Reconstruct directly from raw data (no intermediate output on
322 // digit tree or rec point tree).
325 // reader Raw event reader
327 AliFMDRawReader rawReader(reader, 0);
329 UShort_t det, sec, str, sam, rat, fac;
330 Short_t adc, oldDet = -1;
335 while (rawReader.NextSample(det, rng, sec, str, sam, rat, adc, zs, fac)) {
336 if (!rawReader.SelectSample(sam, rat)) continue;
339 fZSFactor[det-1] = fac;
342 DigitizeSignal(sdigits, det, rng, sec, str, sam, adc);
344 UseRecoParam(kFALSE);
347 //____________________________________________________________________
349 AliFMDReconstructor::Reconstruct(TTree* digitsTree,
350 TTree* clusterTree) const
352 // Reconstruct event from digits in tree
353 // Get the FMD branch holding the digits.
354 // FIXME: The vertex may not be known yet, so we may have to move
355 // some of this to FillESD.
358 // digitsTree Pointer to a tree containing digits
359 // clusterTree Pointer to output tree
361 if (!fMult) fMult = new TClonesArray("AliFMDRecPoint");
363 AliFMDDebug(1, ("Reconstructing from digits in a tree"));
368 static TClonesArray* digits = new TClonesArray("AliFMDDigit");
369 TBranch* digitBranch = digitsTree->GetBranch("FMD");
371 Error("Exec", "No digit branch for the FMD found");
374 digitBranch->SetAddress(&digits);
376 if (digits) digits->Clear();
377 if (fMult) fMult->Clear();
378 if (fESDObj) fESDObj->Clear();
381 fTreeR = clusterTree;
382 fTreeR->Branch("FMD", &fMult);
384 AliFMDDebug(5, ("Getting entry 0 from digit branch"));
385 digitBranch->GetEntry(0);
387 AliFMDDebug(5, ("Processing digits"));
389 ProcessDigits(digits);
390 UseRecoParam(kFALSE);
392 Int_t written = clusterTree->Fill();
393 AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
399 //____________________________________________________________________
401 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
403 // For each digit, find the pseudo rapdity, azimuthal angle, and
404 // number of corrected ADC counts, and pass it on to the algorithms
408 // digits Array of digits
410 Int_t nDigits = digits->GetEntries();
411 AliFMDDebug(2, ("Got %d digits", nDigits));
412 fESDObj->SetNoiseFactor(fNoiseFactor);
413 fESDObj->SetAngleCorrected(fAngleCorrect);
414 for (Int_t i = 0; i < nDigits; i++) {
415 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
416 if (!digit) continue;
421 //____________________________________________________________________
423 AliFMDReconstructor::ProcessDigit(AliFMDDigit* digit) const
425 UShort_t det = digit->Detector();
426 Char_t rng = digit->Ring();
427 UShort_t sec = digit->Sector();
428 UShort_t str = digit->Strip();
429 Short_t adc = digit->Counts();
431 ProcessSignal(det, rng, sec, str, adc);
434 //____________________________________________________________________
436 AliFMDReconstructor::ProcessSignal(UShort_t det,
442 // Process the signal from a single strip
451 AliFMDParameters* param = AliFMDParameters::Instance();
452 // Check that the strip is not marked as dead
453 if (param->IsDead(det, rng, sec, str)) {
454 AliFMDDebug(1, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
461 PhysicalCoordinates(det, rng, sec, str, eta, phi);
463 // Substract pedestal.
464 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
465 if(counts == USHRT_MAX) return;
467 // Gain match digits.
468 Double_t edep = Adc2Energy(det, rng, sec, str, eta, counts);
469 // Get rid of nonsense energy
472 // Make rough multiplicity
473 Double_t mult = Energy2Multiplicity(det, rng, sec, str, edep);
474 // Get rid of nonsense mult
476 // AliWarning(Form("The mutliplicity in FMD%d%c[%2d,%3d]=%f > 20 "
477 // "(ADC: %d, Energy: %f)", det, rng, sec, str, mult,
480 if (mult < 0) return;
481 AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
482 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
483 det, rng, sec, str, adc, counts, edep, mult));
485 // Create a `RecPoint' on the output branch.
488 new ((*fMult)[fNMult]) AliFMDRecPoint(det, rng, sec, str,
489 eta, phi, edep, mult);
490 (void)m; // Suppress warnings about unused variables.
494 fESDObj->SetMultiplicity(det, rng, sec, str, mult);
495 fESDObj->SetEta(det, rng, sec, str, eta);
497 if (fDiagAll) fDiagAll->Fill(adc, mult);
501 //____________________________________________________________________
503 AliFMDReconstructor::DigitizeSignal(TClonesArray* sdigits,
511 // Process the signal from a single strip
520 AliFMDParameters* param = AliFMDParameters::Instance();
521 // Check that the strip is not marked as dead
522 if (param->IsDead(det, rng, sec, str)) {
523 AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", det, rng, sec, str));
527 // Substract pedestal.
528 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc);
529 if(counts == USHRT_MAX || counts == 0) return;
531 // Gain match digits.
532 Double_t edep = Adc2Energy(det, rng, sec, str, counts);
533 // Get rid of nonsense energy
536 Int_t n = sdigits->GetEntriesFast();
537 // AliFMDSDigit* sdigit =
539 AliFMDSDigit(det, rng, sec, str, edep, counts, counts, counts, counts);
540 // sdigit->SetCount(sam, counts);
543 //____________________________________________________________________
545 AliFMDReconstructor::SubtractPedestal(UShort_t det,
552 UShort_t zsNoiseFactor) const
554 AliFMDParameters* param = AliFMDParameters::Instance();
555 Float_t ped = (zsEnabled ? 0 :
556 param->GetPedestal(det, rng, sec, str));
557 Float_t noise = param->GetPedestalWidth(det, rng, sec, str);
558 if(ped < 0 || noise < 0) {
559 AliWarningClass(Form("Invalid pedestal (%f) or noise (%f) "
560 "for FMD%d%c[%02d,%03d]",
561 ped, noise, det, rng, sec, str));
564 AliDebugClass(15, Form("Subtracting pedestal for FMD%d%c[%2d,%3d]=%4d "
565 "(%s w/factor %d, noise factor %f, "
566 "pedestal %8.2f+/-%8.2f)",
567 det, rng, sec, str, adc,
568 (zsEnabled ? "zs'ed" : "straight"),
569 zsNoiseFactor, noiseFactor, ped, noise));
571 Int_t counts = adc + Int_t(zsEnabled ? zsNoiseFactor * noise : - ped);
572 counts = TMath::Max(Int_t(counts), 0);
573 // Calculate the noise factor for suppressing remenants of the noise
574 // peak. If we have done on-line zero suppression, we only check
575 // for noise signals that are larger than the suppressed noise. If
576 // the noise factor used on line is larger than the factor used
577 // here, we do not do this check at all.
580 // Online factor | Read factor | Result
581 // ---------------+--------------+-------------------------------
582 // 2 | 3 | Check if signal > 1 * noise
583 // 3 | 3 | Check if signal > 0
584 // 3 | 2 | Check if signal > 0
586 // In this way, we make sure that we do not suppress away too much
587 // data, and that the read-factor is the most stringent cut.
588 Float_t nf = TMath::Max(0.F, noiseFactor - (zsEnabled ? zsNoiseFactor : 0));
589 if (counts < noise * nf) counts = 0;
590 if (counts > 0) AliDebugClass(15, "Got a hit strip");
592 UShort_t ret = counts < 0 ? 0 : counts;
597 //____________________________________________________________________
599 AliFMDReconstructor::SubtractPedestal(UShort_t det,
605 // Member function to subtract the pedestal from a digit
612 // adc # of ADC counts
614 // Pedestal subtracted signal or USHRT_MAX in case of problems
616 UShort_t counts = SubtractPedestal(det, rng, sec, str, adc,
617 fNoiseFactor, fZS[det-1],
619 if (fDiagStep1) fDiagStep1->Fill(adc, counts);
624 //____________________________________________________________________
626 AliFMDReconstructor::Adc2Energy(UShort_t det,
630 UShort_t count) const
632 // Converts number of ADC counts to energy deposited.
633 // Note, that this member function can be overloaded by derived
634 // classes to do strip-specific look-ups in databases or the like,
635 // to find the proper gain for a strip.
637 // In the first simple version, we calculate the energy deposited as
639 // EnergyDeposited = cos(theta) * gain * count
644 // gain = ----------------- * Energy_deposited_per_MIP
647 // is constant and the same for all strips.
649 // For the production we use the conversion measured in the NBI lab.
650 // The total conversion is then:
655 // => energy = ----------------
663 // counts Number of ADC counts over pedestal
665 // The energy deposited in a single strip, or -1 in case of problems
667 if (count <= 0) return 0;
668 AliFMDParameters* param = AliFMDParameters::Instance();
669 Float_t gain = param->GetPulseGain(det, rng, sec, str);
670 // 'Tagging' bad gains as bad energy
672 AliWarning(Form("Invalid gain (%f) for FMD%d%c[%02d,%03d]",
673 gain, det, rng, sec, str));
676 AliFMDDebug(5, ("Converting counts %d to energy (factor=%f, DAC2MIP=%f)",
677 count, gain,param->GetDACPerMIP()));
679 Double_t edep = ((count * param->GetEdepMip())
680 / (gain * param->GetDACPerMIP()));
684 //____________________________________________________________________
686 AliFMDReconstructor::Adc2Energy(UShort_t det,
691 UShort_t count) const
693 // Converts number of ADC counts to energy deposited.
694 // Note, that this member function can be overloaded by derived
695 // classes to do strip-specific look-ups in databases or the like,
696 // to find the proper gain for a strip.
698 // In the first simple version, we calculate the energy deposited as
700 // EnergyDeposited = cos(theta) * gain * count
705 // gain = ----------------- * Energy_deposited_per_MIP
708 // is constant and the same for all strips.
710 // For the production we use the conversion measured in the NBI lab.
711 // The total conversion is then:
716 // => energy = ----------------
724 // eta Psuedo-rapidity
725 // counts Number of ADC counts over pedestal
727 // The energy deposited in a single strip, or -1 in case of problems
729 Double_t edep = Adc2Energy(det, rng, sec, str, count);
731 if (fDiagStep2) fDiagStep2->Fill(count, edep);
733 Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
734 Double_t corr = TMath::Abs(TMath::Cos(theta));
735 Double_t cedep = corr * edep;
736 AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
737 edep, corr, cedep, eta, theta));
738 if (fDiagStep3) fDiagStep3->Fill(edep, cedep);
744 //____________________________________________________________________
746 AliFMDReconstructor::Energy2Multiplicity(UShort_t /*det*/,
752 // Converts an energy signal to number of particles.
753 // Note, that this member function can be overloaded by derived
754 // classes to do strip-specific look-ups in databases or the like,
755 // to find the proper gain for a strip.
757 // In this simple version, we calculate the multiplicity as
759 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
763 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
765 // is constant and the same for all strips
772 // edep Energy deposited in a single strip
774 // The "bare" multiplicity corresponding to the energy deposited
775 AliFMDParameters* param = AliFMDParameters::Instance();
776 Double_t edepMIP = param->GetEdepMip();
777 Float_t mult = edep / edepMIP;
780 AliFMDDebug(15, ("Translating energy %f to multiplicity via "
781 "divider %f->%f", edep, edepMIP, mult));
783 if (fDiagStep4) fDiagStep4->Fill(edep, mult);
787 //____________________________________________________________________
789 AliFMDReconstructor::PhysicalCoordinates(UShort_t det,
796 // Get the eta and phi of a digit
803 // eta On return, contains the psuedo-rapidity of the strip
804 // phi On return, contains the azimuthal angle of the strip
806 AliFMDGeometry* geom = AliFMDGeometry::Instance();
807 Double_t x, y, z, r, theta, deta, dphi;
808 geom->Detector2XYZ(det, rng, sec, str, x, y, z);
810 // Correct for vertex offset.
812 AliFMDGeometry::XYZ2REtaPhiTheta(x, y, z, r, deta, dphi, theta);
819 //____________________________________________________________________
821 AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
822 TTree* /* clusterTree */,
823 AliESDEvent* esd) const
825 // nothing to be done
826 // FIXME: The vertex may not be known when Reconstruct is executed,
827 // so we may have to move some of that member function here.
828 AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
831 // Fix up ESD so that only truely dead channels get the kInvalidMult flag.
832 MarkDeadChannels(fESDObj);
834 Double_t oldVz = fCurrentVertex;
836 if (fVertexType != kNoVertex) {
837 AliFMDDebug(2, ("Revertexing the ESD data to vz=%f (was %f)",
838 fCurrentVertex, oldVz));
839 AliFMDESDRevertexer revertexer;
840 revertexer.Revertex(fESDObj, fCurrentVertex);
844 AliFMDDebug(2, ("Writing FMD data to ESD tree"));
845 esd->SetFMDData(fESDObj);
848 if (!fDiagnostics || !esd) return;
849 static bool first = true;
850 // This is most likely NOT the event number you'd like to use. It
851 // has nothing to do with the 'real' event number.
852 // - That's OK. We just use it for the name of the directory -
853 // nothing else. Christian
854 Int_t evno = esd->GetEventNumberInFile();
855 AliFMDDebug(3, ("Writing diagnostics histograms to FMD.Diag.root/%03d",evno));
856 TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
859 TDirectory* d = f.mkdir(Form("%03d", evno),
860 Form("Diagnostics histograms for event # %d", evno));
862 if (fDiagStep1) fDiagStep1->Write();
863 if (fDiagStep2) fDiagStep2->Write();
864 if (fDiagStep3) fDiagStep3->Write();
865 if (fDiagStep4) fDiagStep4->Write();
866 if (fDiagAll) fDiagAll->Write();
871 if (fDiagStep1) fDiagStep1->Reset();
872 if (fDiagStep2) fDiagStep2->Reset();
873 if (fDiagStep3) fDiagStep3->Reset();
874 if (fDiagStep4) fDiagStep4->Reset();
875 if (fDiagAll) fDiagAll->Reset();
878 //____________________________________________________________________
880 AliFMDReconstructor::FillESD(AliRawReader*, TTree* clusterTree,
881 AliESDEvent* esd) const
884 FillESD(dummy, clusterTree, esd);
886 //____________________________________________________________________