]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDReconstructor.cxx
0. General code clean-up, including messages, and the like.
[u/mrichter/AliRoot.git] / FMD / AliFMDReconstructor.cxx
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  **************************************************************************/
15
16 /* $Id$ */
17
18 //____________________________________________________________________
19 //
20 // This is a class that constructs AliFMDMult (reconstructed
21 // multiplicity) from of Digits
22 //
23 // This class reads either digits from a TClonesArray or raw data from
24 // a DDL file (or similar), and stores the read ADC counts in an
25 // internal cache (fAdcs). 
26 //
27 // From the cached values it then calculates the number of particles
28 // that hit a region of the FMDs, as specified by the user. 
29 //
30 // The reconstruction can be done in two ways: Either via counting the
31 // number of empty strips (Poisson method), or by converting the ADC
32 // signal to an energy deposition, and then dividing by the typical
33 // energy loss of a particle.
34 // 
35 //      +---------------------+       +---------------------+
36 //      | AliFMDReconstructor |<>-----| AliFMDMultAlgorithm |
37 //      +---------------------+       +---------------------+
38 //                                               ^
39 //                                               |
40 //                                   +-----------+---------+
41 //                                   |                     |
42 //                         +-------------------+   +------------------+
43 //                         | AliFMDMultPoisson |   | AliFMDMultNaiive |
44 //                         +-------------------+   +------------------+
45 //
46 // AliFMDReconstructor acts as a manager class.  It contains a list of
47 // AliFMDMultAlgorithm objects.  The call graph looks something like 
48 //
49 //
50 //       +----------------------+            +----------------------+
51 //       | :AliFMDReconstructor |            | :AliFMDMultAlgorithm |
52 //       +----------------------+            +----------------------+
53 //                  |                                  |
54 //    Reconstruct  +-+                                 |
55 //    ------------>| |                         PreRun +-+
56 //                 | |------------------------------->| |   
57 //                 | |                                +-+
58 //                 | |-----+ (for each event)          |
59 //                 | |     | *ProcessEvent             |
60 //                 |+-+    |                           |
61 //                 || |<---+                 PreEvent +-+
62 //                 || |------------------------------>| |      
63 //                 || |                               +-+
64 //                 || |-----+                          |
65 //                 || |     | ProcessDigits            |
66 //                 ||+-+    |                          |
67 //                 ||| |<---+                          |
68 //                 ||| |         *ProcessDigit(digit) +-+
69 //                 ||| |----------------------------->| |
70 //                 ||| |                              +-+
71 //                 ||+-+                               |
72 //                 || |                     PostEvent +-+
73 //                 || |------------------------------>| |
74 //                 || |                               +-+
75 //                 |+-+                                |
76 //                 | |                        PostRun +-+
77 //                 | |------------------------------->| |
78 //                 | |                                +-+
79 //                 +-+                                 |
80 //                  |                                  |
81 //
82 //
83 // 
84 //-- Authors: Evgeny Karpechev(INR) and Alla Maevsksia
85 //  Latest changes by Christian Holm Christensen <cholm@nbi.dk>
86 //
87 //
88 //____________________________________________________________________
89
90 #include <AliLog.h>                        // ALILOG_H
91 #include <AliRun.h>                        // ALIRUN_H
92 #include <AliRunLoader.h>                  // ALIRUNLOADER_H
93 #include <AliLoader.h>                     // ALILOADER_H
94 #include <AliHeader.h>                     // ALIHEADER_H
95 #include <AliRawReader.h>                  // ALIRAWREADER_H
96 #include <AliGenEventHeader.h>             // ALIGENEVENTHEADER_H
97 #include "AliFMD.h"                        // ALIFMD_H
98 #include "AliFMDGeometry.h"                // ALIFMDGEOMETRY_H
99 #include "AliFMDParameters.h"              // ALIFMDPARAMETERS_H
100 #include "AliFMDDetector.h"                // ALIFMDDETECTOR_H
101 #include "AliFMDRing.h"                    // ALIFMDRING_H
102 #include "AliFMDDigit.h"                   // ALIFMDDIGIT_H
103 #include "AliFMDReconstructor.h"           // ALIFMDRECONSTRUCTOR_H
104 #include "AliFMDRawStream.h"               // ALIFMDRAWSTREAM_H
105 #include "AliFMDRawReader.h"               // ALIFMDRAWREADER_H
106 #include "AliFMDRecPoint.h"                // ALIFMDMULTNAIIVE_H
107 #include "AliESD.h"                        // ALIESD_H
108 #include <AliESDFMD.h>                     // ALIESDFMD_H
109 #include <TFile.h>
110
111 //____________________________________________________________________
112 ClassImp(AliFMDReconstructor)
113 #if 0
114   ; // This is here to keep Emacs for indenting the next line
115 #endif
116
117 //____________________________________________________________________
118 AliFMDReconstructor::AliFMDReconstructor() 
119   : AliReconstructor(),
120     fMult(0), 
121     fESDObj(0)
122 {
123   // Make a new FMD reconstructor object - default CTOR.  
124 }
125   
126
127 //____________________________________________________________________
128 AliFMDReconstructor::AliFMDReconstructor(const AliFMDReconstructor& other) 
129   : AliReconstructor(), 
130     fMult(other.fMult),
131     fESDObj(other.fESDObj)
132 {
133   // Copy constructor 
134 }
135   
136
137 //____________________________________________________________________
138 AliFMDReconstructor&
139 AliFMDReconstructor::operator=(const AliFMDReconstructor& other) 
140 {
141   // Assignment operator
142   fMult   = other.fMult;
143   fESDObj = other.fESDObj;
144   return *this;
145 }
146
147 //____________________________________________________________________
148 AliFMDReconstructor::~AliFMDReconstructor() 
149 {
150   // Destructor 
151   if (fMult)   fMult->Delete();
152   if (fMult)   delete fMult;
153   if (fESDObj) delete fESDObj;
154 }
155
156 //____________________________________________________________________
157 void 
158 AliFMDReconstructor::Init(AliRunLoader* runLoader) 
159 {
160   // Initialize the reconstructor 
161   AliDebug(1, Form("Init called with runloader 0x%x", runLoader));
162   // Initialize the geometry 
163   AliFMDGeometry* fmd = AliFMDGeometry::Instance();
164   fmd->Init();
165   fmd->InitTransformations();
166   
167   // Current vertex position
168   fCurrentVertex = 0;
169   // Create array of reconstructed strip multiplicities 
170   fMult = new TClonesArray("AliFMDRecPoint", 51200);
171   // Create ESD output object 
172   fESDObj = new AliESDFMD;
173   
174   // Check that we have a run loader
175   if (!runLoader) { 
176     Warning("Init", "No run loader");
177     return;
178   }
179
180   // Check if we can get the header tree 
181   AliHeader* header = runLoader->GetHeader();
182   if (!header) {
183     Warning("Init", "No header");
184     return;
185   }
186
187   // Check if we can get a simulation header 
188   AliGenEventHeader* eventHeader = header->GenEventHeader();
189   if (eventHeader) {
190     TArrayF vtx;
191     eventHeader->PrimaryVertex(vtx);
192     fCurrentVertex = vtx[2];
193     AliDebug(1, Form("Primary vertex Z coordinate for event # %d/%d is %f", 
194                      header->GetRun(), header->GetEvent(), fCurrentVertex));
195     Warning("Init", "no generator event header");
196   }
197   else {
198     Warning("Init", "No generator event header - "
199             "perhaps we get the vertex from ESD?");
200   }
201   // Get the ESD tree 
202   SetESD(new AliESD);
203 }
204
205 //____________________________________________________________________
206 void 
207 AliFMDReconstructor::ConvertDigits(AliRawReader* reader, 
208                                    TTree* digitsTree) const
209 {
210   // Convert Raw digits to AliFMDDigit's in a tree 
211   AliDebug(1, "Reading raw data into digits tree");
212   AliFMDRawReader rawRead(reader, digitsTree);
213   // rawRead.SetSampleRate(fFMD->GetSampleRate());
214   rawRead.Exec();
215 }
216
217 //____________________________________________________________________
218 void 
219 AliFMDReconstructor::Reconstruct(TTree* digitsTree, 
220                                  TTree* clusterTree) const 
221 {
222   // Reconstruct event from digits in tree 
223   // Get the FMD branch holding the digits. 
224   // FIXME: The vertex may not be known yet, so we may have to move
225   // some of this to FillESD. 
226   AliDebug(1, "Reconstructing from digits in a tree");
227 #if 0
228   if (fESD) {
229     const AliESDVertex* vertex = fESD->GetVertex();
230     if (vertex) {
231       AliDebug(1, Form("Got vertex from ESD: %f", vertex->GetZv()));
232       fCurrentVertex = vertex->GetZv();
233     }
234   }
235 #endif  
236   TBranch *digitBranch = digitsTree->GetBranch("FMD");
237   if (!digitBranch) {
238     Error("Exec", "No digit branch for the FMD found");
239     return;
240   }
241   TClonesArray* digits = new TClonesArray("AliFMDDigit");
242   digitBranch->SetAddress(&digits);
243
244   // TIter next(&fAlgorithms);
245   // AliFMDMultAlgorithm* algorithm = 0;
246   // while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next()))) {
247   //   AliDebug(10, Form("PreEvent called for algorithm %s", 
248   //                     algorithm->GetName()));    
249   //   algorithm->PreEvent(clusterTree, fCurrentVertex);
250   // }
251   if (fMult)   fMult->Clear();
252   if (fESDObj) fESDObj->Clear();
253   
254   fNMult = 0;
255   fTreeR = clusterTree;
256   fTreeR->Branch("FMD", &fMult);
257   
258   AliDebug(5, "Getting entry 0 from digit branch");
259   digitBranch->GetEntry(0);
260   
261   AliDebug(5, "Processing digits");
262   ProcessDigits(digits);
263
264   // next.Reset();
265   // algorithm = 0;
266   // while ((algorithm = static_cast<AliFMDMultAlgorithm*>(next()))) {
267   //   AliDebug(10, Form("PostEvent called for algorithm %s", 
268   //                     algorithm->GetName()));
269   // algorithm->PostEvent();
270   // }
271   Int_t written = clusterTree->Fill();
272   AliDebug(10, Form("Filled %d bytes into cluster tree", written));
273   digits->Delete();
274   delete digits;
275 }
276  
277
278 //____________________________________________________________________
279 void
280 AliFMDReconstructor::ProcessDigits(TClonesArray* digits) const
281 {
282   // For each digit, find the pseudo rapdity, azimuthal angle, and
283   // number of corrected ADC counts, and pass it on to the algorithms
284   // used. 
285   Int_t nDigits = digits->GetEntries();
286   AliDebug(1, Form("Got %d digits", nDigits));
287   for (Int_t i = 0; i < nDigits; i++) {
288     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
289     AliFMDParameters* param  = AliFMDParameters::Instance();
290     // Check that the strip is not marked as dead 
291     if (param->IsDead(digit->Detector(), digit->Ring(), 
292                       digit->Sector(), digit->Strip())) continue;
293
294     // digit->Print();
295     // Get eta and phi 
296     Float_t eta, phi;
297     PhysicalCoordinates(digit, eta, phi);
298     
299     // Substract pedestal. 
300     UShort_t counts   = SubtractPedestal(digit);
301     
302     // Gain match digits. 
303     Double_t edep     = Adc2Energy(digit, eta, counts);
304     
305     // Make rough multiplicity 
306     Double_t mult     = Energy2Multiplicity(digit, edep);
307     
308     AliDebug(10, Form("FMD%d%c[%2d,%3d]: "
309                       "ADC: %d, Counts: %d, Energy: %f, Mult: %f", 
310                       digit->Detector(), digit->Ring(), digit->Sector(),
311                       digit->Strip(), digit->Counts(), counts, edep, mult));
312     
313     // Create a `RecPoint' on the output branch. 
314     AliFMDRecPoint* m = 
315       new ((*fMult)[fNMult]) AliFMDRecPoint(digit->Detector(), 
316                                             digit->Ring(), 
317                                             digit->Sector(),
318                                             digit->Strip(),
319                                             eta, phi, 
320                                             edep, mult);
321     (void)m; // Suppress warnings about unused variables. 
322     fNMult++;
323
324     fESDObj->SetMultiplicity(digit->Detector(), digit->Ring(), 
325                              digit->Sector(),  digit->Strip(), mult);
326     fESDObj->SetEta(digit->Detector(), digit->Ring(), 
327                     digit->Sector(),  digit->Strip(), eta);
328   }
329 }
330
331 //____________________________________________________________________
332 UShort_t
333 AliFMDReconstructor::SubtractPedestal(AliFMDDigit* digit) const
334 {
335   // Member function to subtract the pedestal from a digit
336   // This implementation does nothing, but a derived class could over
337   // load this to subtract a pedestal that was given in a database or
338   // something like that. 
339
340   Int_t             counts = 0;
341   AliFMDParameters* param  = AliFMDParameters::Instance();
342   Float_t           pedM   = param->GetPedestal(digit->Detector(), 
343                                                 digit->Ring(), 
344                                                 digit->Sector(), 
345                                                 digit->Strip());
346   AliDebug(10, Form("Subtracting pedestal %f from signal %d", 
347                    pedM, digit->Counts()));
348   if (digit->Count3() > 0)      counts = digit->Count3();
349   else if (digit->Count2() > 0) counts = digit->Count2();
350   else                          counts = digit->Count1();
351   counts = TMath::Max(Int_t(counts - pedM), 0);
352   if (counts > 0) AliDebug(10, "Got a hit strip");
353   
354   return  UShort_t(counts);
355 }
356
357 //____________________________________________________________________
358 Float_t
359 AliFMDReconstructor::Adc2Energy(AliFMDDigit* digit, 
360                                 Float_t      /* eta */, 
361                                 UShort_t     count) const
362 {
363   // Converts number of ADC counts to energy deposited. 
364   // Note, that this member function can be overloaded by derived
365   // classes to do strip-specific look-ups in databases or the like,
366   // to find the proper gain for a strip. 
367   // 
368   // In this simple version, we calculate the energy deposited as 
369   // 
370   //    EnergyDeposited = cos(theta) * gain * count
371   // 
372   // where 
373   // 
374   //           Pre_amp_MIP_Range
375   //    gain = ----------------- * Energy_deposited_per_MIP
376   //           ADC_channel_size    
377   // 
378   // is constant and the same for all strips.
379
380   // Double_t theta = 2 * TMath::ATan(TMath::Exp(-eta));
381   // Double_t edep  = TMath::Abs(TMath::Cos(theta)) * fGain * count;
382   AliFMDParameters* param = AliFMDParameters::Instance();
383   Float_t           gain  = param->GetPulseGain(digit->Detector(), 
384                                                 digit->Ring(), 
385                                                 digit->Sector(), 
386                                                 digit->Strip());
387   Double_t          edep  = count * gain;
388   AliDebug(10, Form("Converting counts %d to energy via factor %f", 
389                     count, gain));
390   return edep;
391 }
392
393 //____________________________________________________________________
394 Float_t
395 AliFMDReconstructor::Energy2Multiplicity(AliFMDDigit* /* digit */, 
396                                          Float_t      edep) const
397 {
398   // Converts an energy signal to number of particles. 
399   // Note, that this member function can be overloaded by derived
400   // classes to do strip-specific look-ups in databases or the like,
401   // to find the proper gain for a strip. 
402   // 
403   // In this simple version, we calculate the multiplicity as 
404   // 
405   //   multiplicity = Energy_deposited / Energy_deposited_per_MIP
406   // 
407   // where 
408   //
409   //   Energy_deposited_per_MIP = 1.664 * SI_density * SI_thickness 
410   // 
411   // is constant and the same for all strips 
412   AliFMDParameters* param   = AliFMDParameters::Instance();
413   Double_t          edepMIP = param->GetEdepMip();
414   Float_t           mult    = edep / edepMIP;
415   if (edep > 0) 
416     AliDebug(10, Form("Translating energy %f to multiplicity via "
417                      "divider %f->%f", edep, edepMIP, mult));
418   return mult;
419 }
420
421 //____________________________________________________________________
422 void
423 AliFMDReconstructor::PhysicalCoordinates(AliFMDDigit* digit, 
424                                          Float_t& eta, 
425                                          Float_t& phi) const
426 {
427   // Get the eta and phi of a digit 
428   // 
429   // Get geometry. 
430   AliFMDGeometry* fmd = AliFMDGeometry::Instance();
431 #if 0
432   AliFMDDetector* subDetector = fmd->GetDetector(digit->Detector());
433   if (!subDetector) { 
434     Warning("ProcessDigits", "Unknown detector: FMD%d" , digit->Detector());
435     return;
436   }
437     
438   // Get the ring - we should really use geometry->Detector2XYZ(...)
439   // here. 
440   AliFMDRing* ring  = subDetector->GetRing(digit->Ring());
441   Float_t     ringZ = subDetector->GetRingZ(digit->Ring());
442   if (!ring) {
443     Warning("ProcessDigits", "Unknown ring: FMD%d%c", digit->Detector(), 
444             digit->Ring());
445     return;
446   }
447   Float_t  realZ    = fCurrentVertex + ringZ;
448   Float_t  stripR   = ((ring->GetHighR() - ring->GetLowR()) 
449                        / ring->GetNStrips() * (digit->Strip() + .5) 
450                        + ring->GetLowR());
451   Float_t  theta    = TMath::ATan2(stripR, realZ);
452 #endif    
453   Double_t x, y, z, r, theta;
454   fmd->Detector2XYZ(digit->Detector(), digit->Ring(), digit->Sector(), 
455                     digit->Strip(), x, y, z);
456   // Correct for vertex offset. 
457   z     += fCurrentVertex;
458   phi   =  TMath::ATan2(y, x);
459   r     =  TMath::Sqrt(y * y + x * x);
460   theta =  TMath::ATan2(r, z);
461   eta   = -TMath::Log(TMath::Tan(theta / 2));
462 }
463
464       
465
466 //____________________________________________________________________
467 void 
468 AliFMDReconstructor::FillESD(TTree*  /* digitsTree */, 
469                              TTree*  /* clusterTree */,
470                              AliESD* esd) const
471 {
472   // nothing to be done
473   // FIXME: The vertex may not be known when Reconstruct is executed,
474   // so we may have to move some of that member function here. 
475   AliDebug(1, Form("Calling FillESD with two trees and one ESD"));
476   // fESDObj->Print();
477
478   if (esd) { 
479     AliDebug(1, Form("Writing FMD data to ESD tree"));
480     esd->SetFMDData(fESDObj);
481     // Let's check the data in the ESD
482 #if 0
483     AliESDFMD* fromEsd = esd->GetFMDData();
484     if (!fromEsd) {
485       AliWarning("No FMD object in ESD!");
486       return;
487     }
488     for (UShort_t det = 1; det <= fESDObj->MaxDetectors(); det++) {
489       for (UShort_t ir = 0; ir < fESDObj->MaxRings(); ir++) {
490         Char_t ring = (ir == 0 ? 'I' : 'O');
491         for (UShort_t sec = 0; sec < fESDObj->MaxSectors(); sec++) {
492           for (UShort_t str = 0; str < fESDObj->MaxStrips(); str++) {
493             if (fESDObj->Multiplicity(det, ring, sec, str) != 
494                 fromEsd->Multiplicity(det, ring, sec, str))
495               AliWarning(Form("Mult for FMD%d%c[%2d,%3d]",det,ring,sec,str));
496             if (fESDObj->Eta(det, ring, sec, str) != 
497                 fromEsd->Eta(det, ring, sec, str))
498               AliWarning(Form("Eta for FMD%d%c[%2d,%3d]", det,ring,sec,str));
499             if (fESDObj->Multiplicity(det, ring, sec, str) > 0 && 
500                 fESDObj->Multiplicity(det, ring, sec, str) 
501                 != AliESDFMD::kInvalidMult) 
502               AliInfo(Form("Mult in FMD%d%c[%2d,%3d] is %f", det,ring,sec,str,
503                            fESDObj->Multiplicity(det, ring, sec, str)));
504           }
505         }
506       }
507     }
508 #endif
509   }
510
511 #if 0  
512   static Int_t evNo = -1;
513   evNo++;
514   if (esd) evNo = esd->GetEventNumber();
515   TString fname(Form("FMD.ESD.%03d.root", evNo));
516   TFile* file = TFile::Open(fname.Data(), "RECREATE");
517   if (!file) {
518     AliError(Form("Failed to open file %s", fname.Data()));
519     return;
520   }
521   fESDObj->Write();
522   file->Close();
523 #endif
524     
525 #if 0
526   TClonesArray* multStrips  = 0;
527   TClonesArray* multRegions = 0;
528   TTree*        treeR  = fmdLoader->TreeR();
529   TBranch*      branchRegions = treeR->GetBranch("FMDPoisson");
530   TBranch*      branchStrips  = treeR->GetBranch("FMDNaiive");
531   branchRegions->SetAddress(&multRegions);
532   branchStrips->SetAddress(&multStrips);
533   
534   Int_t total = 0;
535   Int_t nEntries  = clusterTree->GetEntries();
536   for (Int_t entry = 0; entry < nEntries; entry++) {
537     AliDebug(5, Form("Entry # %d in cluster tree", entry));
538     treeR->GetEntry(entry);
539     
540     
541     Int_t nMults = multRegions->GetLast();
542     for (Int_t i = 0; i <= nMults; i++) {
543       AliFMDMultRegion* multR =
544         static_cast<AliFMDMultRegion*>(multRegions->UncheckedAt(i));
545       Int_t nParticles=multR->Particles();
546       if (i>=0 && i<=13)   hEtaPoissonI1->AddBinContent(i+1,nParticles);
547       if (i>=14 && i<=27 ) hEtaPoissonI2->AddBinContent(i-13,nParticles);
548       if (i>=28 && i<=33 );
549       if (i>=34 && i<=47 ) hEtaPoissonI3->AddBinContent(48-i,nParticles);
550       if (i>=48 && i<=53)  hEtaPoissonO3->AddBinContent(54-i,nParticles);
551     }
552   }
553 #endif   
554 }
555
556
557 //____________________________________________________________________
558 void 
559 AliFMDReconstructor::Reconstruct(AliRawReader*,TTree*) const 
560 {
561   // Cannot be used.  See member function with same name but with 2
562   // TTree arguments.   Make sure you do local reconstrucion 
563   AliDebug(1, Form("Calling FillESD with loader and tree"));
564   AliError("MayNotUse");
565 }
566 //____________________________________________________________________
567 void 
568 AliFMDReconstructor::Reconstruct(AliRunLoader*) const 
569 {
570   // Cannot be used.  See member function with same name but with 2
571   // TTree arguments.   Make sure you do local reconstrucion 
572   AliDebug(1, Form("Calling FillESD with loader"));
573   AliError("MayNotUse");
574 }
575 //____________________________________________________________________
576 void 
577 AliFMDReconstructor::Reconstruct(AliRunLoader*, AliRawReader*) const 
578 {
579   // Cannot be used.  See member function with same name but with 2
580   // TTree arguments.   Make sure you do local reconstrucion 
581   AliDebug(1, Form("Calling FillESD with loader and raw reader"));
582   AliError("MayNotUse");
583 }
584 //____________________________________________________________________
585 void 
586 AliFMDReconstructor::FillESD(AliRawReader*,TTree*,AliESD*) const 
587 {
588   // Cannot be used.  See member function with same name but with 2
589   // TTree arguments.   Make sure you do local reconstrucion 
590   AliDebug(1, Form("Calling FillESD with raw reader, tree, and ESD"));
591   AliError("MayNotUse");
592 }
593 //____________________________________________________________________
594 void 
595 AliFMDReconstructor::FillESD(AliRunLoader*,AliESD*) const
596 {
597   // Cannot be used.  See member function with same name but with 2
598   // TTree arguments.   Make sure you do local reconstrucion 
599   AliDebug(1, Form("Calling FillESD with loader and ESD"));
600   AliError("MayNotUse");
601 }
602 //____________________________________________________________________
603 void 
604 AliFMDReconstructor::FillESD(AliRunLoader*,AliRawReader*,AliESD*) const 
605 {
606   // Cannot be used.  See member function with same name but with 2
607   // TTree arguments.   Make sure you do local reconstrucion 
608   AliDebug(1, Form("Calling FillESD with loader, raw reader, and ESD"));
609   AliError("MayNotUse");
610 }
611
612 //____________________________________________________________________
613 //
614 // EOF
615 //