]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDReconstructor.cxx
Added code to do flow analysis.
[u/mrichter/AliRoot.git] / FMD / AliFMDReconstructor.cxx
CommitLineData
37c55dc0 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
0d0e6995 15/* $Id$ */
c2fc1258 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
20*/
4347b38f 21//____________________________________________________________________
42403906 22//
02a27b50 23// This is a class that constructs AliFMDRecPoint objects from of Digits
6169f936 24// This class reads either digits from a TClonesArray or raw data from
4347b38f 25// a DDL file (or similar), and stores the read ADC counts in an
6169f936 26// internal cache (fAdcs). The rec-points are made via the naiive
27// method.
4347b38f 28//
37c55dc0 29//-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
4347b38f 30// Latest changes by Christian Holm Christensen <cholm@nbi.dk>
31//
32//
33//____________________________________________________________________
37c55dc0 34
f95a63c4 35// #include <AliLog.h> // ALILOG_H
6169f936 36// #include <AliRun.h> // ALIRUN_H
f95a63c4 37#include "AliFMDDebug.h"
f011bd38 38// to be removed as soon as we remove it from the base class
39#include "AliRunLoader.h"
1a1fdef7 40#include "AliFMDGeometry.h" // ALIFMDGEOMETRY_H
41#include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
56b1929b 42#include "AliFMDDigit.h" // ALIFMDDIGIT_H
43#include "AliFMDReconstructor.h" // ALIFMDRECONSTRUCTOR_H
56b1929b 44#include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
bf000c32 45#include "AliFMDRecPoint.h" // ALIFMDMULTNAIIVE_H
af885e0f 46#include "AliESDEvent.h" // ALIESDEVENT_H
aad72f45 47#include "AliESDVertex.h" // ALIESDVERTEX_H
8f6ee336 48#include <AliESDFMD.h> // ALIESDFMD_H
aad72f45 49#include <TMath.h>
9684be2f 50#include <TH1.h>
51#include <TH2.h>
52#include <TFile.h>
02a27b50 53class AliRawReader;
4347b38f 54
55//____________________________________________________________________
925e6570 56ClassImp(AliFMDReconstructor)
1a1fdef7 57#if 0
58 ; // This is here to keep Emacs for indenting the next line
59#endif
37c55dc0 60
4347b38f 61//____________________________________________________________________
62AliFMDReconstructor::AliFMDReconstructor()
63 : AliReconstructor(),
e9fd1e20 64 fMult(0x0),
65 fNMult(0),
66 fTreeR(0x0),
67 fCurrentVertex(0),
68 fESDObj(0x0),
9684be2f 69 fNoiseFactor(0),
70 fAngleCorrect(kTRUE),
71 fVertexType(kNoVertex),
9684be2f 72 fESD(0x0),
73 fDiagnostics(kFALSE),
74 fDiagStep1(0),
75 fDiagStep2(0),
76 fDiagStep3(0),
77 fDiagStep4(0),
78 fDiagAll(0)
4347b38f 79{
8f6ee336 80 // Make a new FMD reconstructor object - default CTOR.
a9579262 81 SetNoiseFactor();
82 SetAngleCorrect();
4347b38f 83}
84
42403906 85
86//____________________________________________________________________
87AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other)
8f6ee336 88 : AliReconstructor(),
89 fMult(other.fMult),
e9fd1e20 90 fNMult(other.fNMult),
91 fTreeR(other.fTreeR),
92 fCurrentVertex(other.fCurrentVertex),
93 fESDObj(other.fESDObj),
a9579262 94 fNoiseFactor(other.fNoiseFactor),
9684be2f 95 fAngleCorrect(other.fAngleCorrect),
96 fVertexType(other.fVertexType),
9684be2f 97 fESD(other.fESD),
98 fDiagnostics(other.fDiagnostics),
99 fDiagStep1(other.fDiagStep1),
100 fDiagStep2(other.fDiagStep2),
101 fDiagStep3(other.fDiagStep3),
102 fDiagStep4(other.fDiagStep4),
103 fDiagAll(other.fDiagAll)
42403906 104{
56b1929b 105 // Copy constructor
42403906 106}
107
108
109//____________________________________________________________________
110AliFMDReconstructor&
111AliFMDReconstructor::operator=(const AliFMDReconstructor& other)
112{
56b1929b 113 // Assignment operator
a9579262 114 fMult = other.fMult;
115 fNMult = other.fNMult;
116 fTreeR = other.fTreeR;
e9fd1e20 117 fCurrentVertex = other.fCurrentVertex;
a9579262 118 fESDObj = other.fESDObj;
a9579262 119 fNoiseFactor = other.fNoiseFactor;
120 fAngleCorrect = other.fAngleCorrect;
9684be2f 121 fVertexType = other.fVertexType;
9684be2f 122 fESD = other.fESD;
123 fDiagnostics = other.fDiagnostics;
124 fDiagStep1 = other.fDiagStep1;
125 fDiagStep2 = other.fDiagStep2;
126 fDiagStep3 = other.fDiagStep3;
127 fDiagStep4 = other.fDiagStep4;
128 fDiagAll = other.fDiagAll;
42403906 129 return *this;
130}
56b1929b 131
132//____________________________________________________________________
133AliFMDReconstructor::~AliFMDReconstructor()
134{
135 // Destructor
8f6ee336 136 if (fMult) fMult->Delete();
137 if (fMult) delete fMult;
138 if (fESDObj) delete fESDObj;
56b1929b 139}
4347b38f 140
141//____________________________________________________________________
142void
f011bd38 143AliFMDReconstructor::Init(AliRunLoader* /*runLoader*/)
1a1fdef7 144{
145 // Initialize the reconstructor
f011bd38 146
bf000c32 147 // Initialize the geometry
9b48326f 148 AliFMDGeometry* geom = AliFMDGeometry::Instance();
149 geom->Init();
150 geom->InitTransformations();
151
152 // Initialize the parameters
153 AliFMDParameters* param = AliFMDParameters::Instance();
154 param->Init();
bf000c32 155
8f6ee336 156 // Current vertex position
1a1fdef7 157 fCurrentVertex = 0;
8f6ee336 158 // Create array of reconstructed strip multiplicities
bf000c32 159 fMult = new TClonesArray("AliFMDRecPoint", 51200);
8f6ee336 160 // Create ESD output object
161 fESDObj = new AliESDFMD;
162
9684be2f 163 // Check if we need diagnostics histograms
164 if (!fDiagnostics) return;
165 fDiagStep1 = new TH2I("diagStep1", "Read ADC vs. Noise surpressed ADC",
166 1024, -.5, 1023.5, 1024, -.5, 1023.5);
167 fDiagStep1->SetDirectory(0);
168 fDiagStep1->GetXaxis()->SetTitle("ADC (read)");
169 fDiagStep1->GetYaxis()->SetTitle(Form("ADC (noise surpressed %4.f)",
170 fNoiseFactor));
171 fDiagStep2 = new TH2F("diagStep2", "ADC vs Edep deduced",
172 1024, -.5, 1023.5, 100, 0, 2);
173 fDiagStep2->SetDirectory(0);
174 fDiagStep2->GetXaxis()->SetTitle("ADC (noise surpressed)");
175 fDiagStep2->GetYaxis()->SetTitle("#Delta E [GeV]");
176 fDiagStep3 = new TH2F("diagStep3", "Edep vs Edep path corrected",
177 100, 0., 2., 100, 0., 2.);
178 fDiagStep3->SetDirectory(0);
179 fDiagStep3->GetXaxis()->SetTitle("#Delta E [GeV]");
180 fDiagStep3->GetYaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
181 fDiagStep4 = new TH2F("diagStep4", "Edep vs Multiplicity deduced",
182 100, 0., 2., 100, -.1, 19.9);
183 fDiagStep4->SetDirectory(0);
184 fDiagStep4->GetXaxis()->SetTitle("#Delta E/#Delta x #times #delta x [GeV]");
185 fDiagStep4->GetYaxis()->SetTitle("Multiplicity");
186 fDiagAll = new TH2F("diagAll", "Read ADC vs Multiplicity deduced",
187 1024, -.5, 1023.5, 100, -.1, 19.9);
188 fDiagAll->SetDirectory(0);
189 fDiagAll->GetXaxis()->SetTitle("ADC (read)");
190 fDiagAll->GetYaxis()->SetTitle("Multiplicity");
4347b38f 191}
dc8af42e 192
4347b38f 193//____________________________________________________________________
194void
1a1fdef7 195AliFMDReconstructor::ConvertDigits(AliRawReader* reader,
196 TTree* digitsTree) const
197{
198 // Convert Raw digits to AliFMDDigit's in a tree
f95a63c4 199 AliFMDDebug(2, ("Reading raw data into digits tree"));
1a1fdef7 200 AliFMDRawReader rawRead(reader, digitsTree);
201 // rawRead.SetSampleRate(fFMD->GetSampleRate());
202 rawRead.Exec();
4347b38f 203}
204
4347b38f 205//____________________________________________________________________
206void
9684be2f 207AliFMDReconstructor::GetVertex() const
4347b38f 208{
9684be2f 209 fVertexType = kNoVertex;
210 fCurrentVertex = 0;
a3537838 211 if (fESD) {
8f6ee336 212 const AliESDVertex* vertex = fESD->GetVertex();
213 if (vertex) {
f95a63c4 214 AliFMDDebug(2, ("Got vertex from ESD: %f", vertex->GetZv()));
8f6ee336 215 fCurrentVertex = vertex->GetZv();
9684be2f 216 fVertexType = kESDVertex;
217 return;
8f6ee336 218 }
a3537838 219 }
9684be2f 220 AliWarning("Didn't get any vertex from ESD or generator");
221}
222
223
224//____________________________________________________________________
225void
226AliFMDReconstructor::Reconstruct(TTree* digitsTree,
227 TTree* clusterTree) const
228{
229 // Reconstruct event from digits in tree
230 // Get the FMD branch holding the digits.
231 // FIXME: The vertex may not be known yet, so we may have to move
232 // some of this to FillESD.
f95a63c4 233 AliFMDDebug(2, ("Reconstructing from digits in a tree"));
9684be2f 234 GetVertex();
235
1a1fdef7 236 TBranch *digitBranch = digitsTree->GetBranch("FMD");
4347b38f 237 if (!digitBranch) {
1a1fdef7 238 Error("Exec", "No digit branch for the FMD found");
239 return;
4347b38f 240 }
0abe7182 241 TClonesArray* digits = new TClonesArray("AliFMDDigit");
1a1fdef7 242 digitBranch->SetAddress(&digits);
4347b38f 243
8f6ee336 244 if (fMult) fMult->Clear();
245 if (fESDObj) fESDObj->Clear();
246
247 fNMult = 0;
248 fTreeR = clusterTree;
249 fTreeR->Branch("FMD", &fMult);
250
f95a63c4 251 AliFMDDebug(5, ("Getting entry 0 from digit branch"));
1a1fdef7 252 digitBranch->GetEntry(0);
253
f95a63c4 254 AliFMDDebug(5, ("Processing digits"));
e802be3e 255 ProcessDigits(digits);
256
8f6ee336 257 Int_t written = clusterTree->Fill();
f95a63c4 258 AliFMDDebug(10, ("Filled %d bytes into cluster tree", written));
0abe7182 259 digits->Delete();
260 delete digits;
4347b38f 261}
1a1fdef7 262
8f6ee336 263
4347b38f 264//____________________________________________________________________
e802be3e 265void
266AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
4347b38f 267{
69b696b9 268 // For each digit, find the pseudo rapdity, azimuthal angle, and
269 // number of corrected ADC counts, and pass it on to the algorithms
270 // used.
e802be3e 271 Int_t nDigits = digits->GetEntries();
f95a63c4 272 AliFMDDebug(1, ("Got %d digits", nDigits));
a9579262 273 fESDObj->SetNoiseFactor(fNoiseFactor);
274 fESDObj->SetAngleCorrected(fAngleCorrect);
e802be3e 275 for (Int_t i = 0; i < nDigits; i++) {
276 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
8f6ee336 277 AliFMDParameters* param = AliFMDParameters::Instance();
278 // Check that the strip is not marked as dead
279 if (param->IsDead(digit->Detector(), digit->Ring(),
15b17c89 280 digit->Sector(), digit->Strip())) {
f95a63c4 281 AliFMDDebug(10, ("FMD%d%c[%2d,%3d] is dead", digit->Detector(),
15b17c89 282 digit->Ring(), digit->Sector(), digit->Strip()));
283 continue;
284 }
8f6ee336 285
286 // digit->Print();
287 // Get eta and phi
288 Float_t eta, phi;
289 PhysicalCoordinates(digit, eta, phi);
4347b38f 290
8f6ee336 291 // Substract pedestal.
e802be3e 292 UShort_t counts = SubtractPedestal(digit);
4347b38f 293
8f6ee336 294 // Gain match digits.
295 Double_t edep = Adc2Energy(digit, eta, counts);
296
297 // Make rough multiplicity
298 Double_t mult = Energy2Multiplicity(digit, edep);
299
f95a63c4 300 AliFMDDebug(10, ("FMD%d%c[%2d,%3d]: "
8f6ee336 301 "ADC: %d, Counts: %d, Energy: %f, Mult: %f",
302 digit->Detector(), digit->Ring(), digit->Sector(),
303 digit->Strip(), digit->Counts(), counts, edep, mult));
304
305 // Create a `RecPoint' on the output branch.
bf000c32 306 AliFMDRecPoint* m =
307 new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(),
308 digit->Ring(),
309 digit->Sector(),
310 digit->Strip(),
311 eta, phi,
312 edep, mult);
8f6ee336 313 (void)m; // Suppress warnings about unused variables.
314 fNMult++;
315
316 fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(),
317 digit->Sector(), digit->Strip(), mult);
318 fESDObj->SetEta(digit->Detector(), digit->Ring(),
319 digit->Sector(), digit->Strip(), eta);
9684be2f 320
321 if (fDiagAll) fDiagAll->Fill(digit->Counts(), mult);
4347b38f 322 }
4347b38f 323}
8f6ee336 324
1a1fdef7 325//____________________________________________________________________
326UShort_t
327AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
328{
329 // Member function to subtract the pedestal from a digit
330 // This implementation does nothing, but a derived class could over
331 // load this to subtract a pedestal that was given in a database or
332 // something like that.
333
8f6ee336 334 Int_t counts = 0;
9684be2f 335 Int_t adc = 0;
8f6ee336 336 AliFMDParameters* param = AliFMDParameters::Instance();
a9579262 337 Float_t ped = param->GetPedestal(digit->Detector(),
8f6ee336 338 digit->Ring(),
339 digit->Sector(),
340 digit->Strip());
a9579262 341 Float_t noise = param->GetPedestalWidth(digit->Detector(),
342 digit->Ring(),
343 digit->Sector(),
344 digit->Strip());
f95a63c4 345 AliFMDDebug(15, ("Subtracting pedestal %f from signal %d",
a9579262 346 ped, digit->Counts()));
9684be2f 347 if (digit->Count3() > 0) adc = digit->Count3();
348 else if (digit->Count2() > 0) adc = digit->Count2();
349 else adc = digit->Count1();
350 counts = TMath::Max(Int_t(adc - ped), 0);
a9579262 351 if (counts < noise * fNoiseFactor) counts = 0;
f95a63c4 352 if (counts > 0) AliFMDDebug(15, ("Got a hit strip"));
9684be2f 353 if (fDiagStep1) fDiagStep1->Fill(adc, counts);
8f6ee336 354
1a1fdef7 355 return UShort_t(counts);
356}
357
8f6ee336 358//____________________________________________________________________
359Float_t
360AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit,
a9579262 361 Float_t eta,
8f6ee336 362 UShort_t count) const
363{
364 // Converts number of ADC counts to energy deposited.
365 // Note, that this member function can be overloaded by derived
366 // classes to do strip-specific look-ups in databases or the like,
367 // to find the proper gain for a strip.
368 //
369 // In this simple version, we calculate the energy deposited as
370 //
371 // EnergyDeposited = cos(theta) * gain * count
372 //
373 // where
374 //
375 // Pre_amp_MIP_Range
376 // gain = ----------------- * Energy_deposited_per_MIP
377 // ADC_channel_size
378 //
379 // is constant and the same for all strips.
a9579262 380 if (count <= 0) return 0;
8f6ee336 381 AliFMDParameters* param = AliFMDParameters::Instance();
382 Float_t gain = param->GetPulseGain(digit->Detector(),
383 digit->Ring(),
384 digit->Sector(),
385 digit->Strip());
f95a63c4 386 AliFMDDebug(15, ("Converting counts %d to energy via factor %f",
8f6ee336 387 count, gain));
a9579262 388
389 Double_t edep = count * gain;
9684be2f 390 if (fDiagStep2) fDiagStep2->Fill(count, edep);
a9579262 391 if (fAngleCorrect) {
9684be2f 392 Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
393 Double_t corr = TMath::Abs(TMath::Cos(theta));
394 Double_t cedep = corr * edep;
f95a63c4 395 AliFMDDebug(10, ("correcting for path %f * %f = %f (eta=%f, theta=%f)",
9684be2f 396 edep, corr, cedep, eta, theta));
397 if (fDiagStep3) fDiagStep3->Fill(edep, cedep);
398 edep = cedep;
a9579262 399 }
8f6ee336 400 return edep;
401}
402
403//____________________________________________________________________
404Float_t
405AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */,
406 Float_t edep) const
407{
408 // Converts an energy signal to number of particles.
409 // Note, that this member function can be overloaded by derived
410 // classes to do strip-specific look-ups in databases or the like,
411 // to find the proper gain for a strip.
412 //
413 // In this simple version, we calculate the multiplicity as
414 //
415 // multiplicity = Energy_deposited / Energy_deposited_per_MIP
416 //
417 // where
418 //
419 // Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness
420 //
421 // is constant and the same for all strips
422 AliFMDParameters* param = AliFMDParameters::Instance();
423 Double_t edepMIP = param->GetEdepMip();
424 Float_t mult = edep / edepMIP;
425 if (edep > 0)
f95a63c4 426 AliFMDDebug(15, ("Translating energy %f to multiplicity via "
8f6ee336 427 "divider %f->%f", edep, edepMIP, mult));
9684be2f 428 if (fDiagStep4) fDiagStep4->Fill(edep, mult);
8f6ee336 429 return mult;
430}
431
432//____________________________________________________________________
433void
434AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit,
435 Float_t& eta,
436 Float_t& phi) const
437{
438 // Get the eta and phi of a digit
439 //
440 // Get geometry.
9b48326f 441 AliFMDGeometry* geom = AliFMDGeometry::Instance();
bf000c32 442 Double_t x, y, z, r, theta;
9b48326f 443 geom->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(),
bf000c32 444 digit->Strip(), x, y, z);
445 // Correct for vertex offset.
446 z += fCurrentVertex;
447 phi = TMath::ATan2(y, x);
448 r = TMath::Sqrt(y * y + x * x);
449 theta = TMath::ATan2(r, z);
450 eta = -TMath::Log(TMath::Tan(theta / 2));
8f6ee336 451}
452
453
454
4347b38f 455//____________________________________________________________________
456void
1a1fdef7 457AliFMDReconstructor::FillESD(TTree* /* digitsTree */,
458 TTree* /* clusterTree */,
af885e0f 459 AliESDEvent* esd) const
121a60bd 460{
42403906 461 // nothing to be done
69b696b9 462 // FIXME: The vertex may not be known when Reconstruct is executed,
463 // so we may have to move some of that member function here.
f95a63c4 464 AliFMDDebug(2, ("Calling FillESD with two trees and one ESD"));
8f6ee336 465 // fESDObj->Print();
466
467 if (esd) {
f95a63c4 468 AliFMDDebug(2, ("Writing FMD data to ESD tree"));
8f6ee336 469 esd->SetFMDData(fESDObj);
a3537838 470 }
9684be2f 471
472 if (!fDiagnostics || !esd) return;
473 static bool first = true;
31fd97b2 474 Int_t evno = esd->GetEventNumberInFile(); // This is most likely NOT the event number you'd like to use. It has nothing to do with the 'real' event number.
f95a63c4 475 AliFMDDebug(1, ("Writing diagnostics histograms to FMD.Diag.root/%03d",
9684be2f 476 evno));
477 TFile f("FMD.Diag.root", (first ? "RECREATE" : "UPDATE"));
478 first = false;
479 f.cd();
480 TDirectory* d = f.mkdir(Form("%03d", evno),
481 Form("Diagnostics histograms for event # %d", evno));
482 d->cd();
483 if (fDiagStep1) fDiagStep1->Write();
484 if (fDiagStep2) fDiagStep2->Write();
485 if (fDiagStep3) fDiagStep3->Write();
486 if (fDiagStep4) fDiagStep4->Write();
487 if (fDiagAll) fDiagAll->Write();
488 d->Write();
489 f.Write();
490 f.Close();
491
492 if (fDiagStep1) fDiagStep1->Reset();
493 if (fDiagStep2) fDiagStep2->Reset();
494 if (fDiagStep3) fDiagStep3->Reset();
495 if (fDiagStep4) fDiagStep4->Reset();
496 if (fDiagAll) fDiagAll->Reset();
121a60bd 497}
498
54e415a8 499
500//____________________________________________________________________
501void
502AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const
503{
504 // Cannot be used. See member function with same name but with 2
505 // TTree arguments. Make sure you do local reconstrucion
f95a63c4 506 AliFMDDebug(2, ("Calling FillESD with loader and tree"));
54e415a8 507 AliError("MayNotUse");
508}
509//____________________________________________________________________
510void
511AliFMDReconstructor::Reconstruct(AliRunLoader*) const
512{
513 // Cannot be used. See member function with same name but with 2
514 // TTree arguments. Make sure you do local reconstrucion
f95a63c4 515 AliFMDDebug(2, ("Calling FillESD with loader"));
54e415a8 516 AliError("MayNotUse");
517}
518//____________________________________________________________________
519void
520AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const
521{
522 // Cannot be used. See member function with same name but with 2
523 // TTree arguments. Make sure you do local reconstrucion
f95a63c4 524 AliFMDDebug(2, ("Calling FillESD with loader and raw reader"));
54e415a8 525 AliError("MayNotUse");
526}
527//____________________________________________________________________
528void
af885e0f 529AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESDEvent*) const
54e415a8 530{
531 // Cannot be used. See member function with same name but with 2
532 // TTree arguments. Make sure you do local reconstrucion
f95a63c4 533 AliFMDDebug(2, ("Calling FillESD with raw reader, tree, and ESD"));
54e415a8 534 AliError("MayNotUse");
535}
536//____________________________________________________________________
537void
af885e0f 538AliFMDReconstructor::FillESD(AliRunLoader*,AliESDEvent*) const
54e415a8 539{
540 // Cannot be used. See member function with same name but with 2
541 // TTree arguments. Make sure you do local reconstrucion
f95a63c4 542 AliFMDDebug(2, ("Calling FillESD with loader and ESD"));
54e415a8 543 AliError("MayNotUse");
544}
545//____________________________________________________________________
546void
af885e0f 547AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESDEvent*) const
54e415a8 548{
549 // Cannot be used. See member function with same name but with 2
550 // TTree arguments. Make sure you do local reconstrucion
f95a63c4 551 AliFMDDebug(2, ("Calling FillESD with loader, raw reader, and ESD"));
54e415a8 552 AliError("MayNotUse");
553}
554
42403906 555//____________________________________________________________________
556//
557// EOF
558//