]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDInput.cxx
Fixes for reading zero-suppressed data. These should be propagated to
[u/mrichter/AliRoot.git] / FMD / AliFMDInput.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, 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 /* $Id$ */
16 /** @file    AliFMDInput.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:42:40 2006
19     @brief   FMD utility classes for reading FMD data
20 */
21 //___________________________________________________________________
22 //
23 // The classes defined here, are utility classes for reading in data
24 // for the FMD.  They are  put in a seperate library to not polute the
25 // normal libraries.  The classes are intended to be used as base
26 // classes for customized class that do some sort of analysis on the
27 // various types of data produced by the FMD. 
28 //
29 // Latest changes by Christian Holm Christensen
30 //
31 #include "AliFMDInput.h"        // ALIFMDHIT_H
32 #include "AliFMDDebug.h"        // ALIFMDDEBUG_H ALILOG_H
33 #include "AliLoader.h"          // ALILOADER_H
34 #include "AliRunLoader.h"       // ALIRUNLOADER_H
35 #include "AliRun.h"             // ALIRUN_H
36 #include "AliStack.h"           // ALISTACK_H
37 #include "AliRawReaderFile.h"   // ALIRAWREADERFILE_H
38 #include "AliRawReaderRoot.h"   // ALIRAWREADERROOT_H
39 #include "AliRawReaderDate.h"   // ALIRAWREADERDATE_H
40 #include "AliFMD.h"             // ALIFMD_H
41 #include "AliFMDHit.h"          // ALIFMDHIT_H
42 #include "AliFMDDigit.h"        // ALIFMDDigit_H
43 #include "AliFMDSDigit.h"       // ALIFMDDigit_H
44 #include "AliFMDRecPoint.h"     // ALIFMDRECPOINT_H
45 #include "AliFMDRawReader.h"    // ALIFMDRAWREADER_H
46 #include <AliESD.h>
47 #include <AliESDFMD.h>
48 #include <AliESDEvent.h>
49 #include <AliCDBManager.h>
50 #include <AliCDBEntry.h>
51 #include <AliAlignObjParams.h>
52 #include <TTree.h>              // ROOT_TTree
53 #include <TChain.h>             // ROOT_TChain
54 #include <TParticle.h>          // ROOT_TParticle
55 #include <TString.h>            // ROOT_TString
56 #include <TDatabasePDG.h>       // ROOT_TDatabasePDG
57 #include <TMath.h>              // ROOT_TMath
58 #include <TGeoManager.h>        // ROOT_TGeoManager 
59 #include <TSystemDirectory.h>   // ROOT_TSystemDirectory
60 #include <Riostream.h>          // ROOT_Riostream
61 #include <TFile.h>              // ROOT_TFile
62 #include <TStreamerInfo.h>
63 #include <TArrayF.h>
64
65 //____________________________________________________________________
66 ClassImp(AliFMDInput)
67 #if 0
68   ; // This is here to keep Emacs for indenting the next line
69 #endif
70
71
72 //____________________________________________________________________
73 AliFMDInput::AliFMDInput()
74   : TNamed("AliFMDInput", "Input handler for various FMD data"), 
75     fGAliceFile(""), 
76     fLoader(0),
77     fRun(0), 
78     fStack(0),
79     fFMDLoader(0), 
80     fReader(0),
81     fFMDReader(0),
82     fFMD(0),
83     fESD(0),
84     fESDEvent(0),
85     fTreeE(0),
86     fTreeH(0),
87     fTreeD(0),
88     fTreeS(0),
89     fTreeR(0), 
90     fTreeA(0), 
91     fChainE(0),
92     fArrayE(0),
93     fArrayH(0),
94     fArrayD(0),
95     fArrayS(0), 
96     fArrayR(0), 
97     fArrayA(0), 
98     fHeader(0),
99     fGeoManager(0),
100     fTreeMask(0), 
101     fRawFile(""),
102     fIsInit(kFALSE),
103     fEventCount(0)
104 {
105
106   // Constructor of an FMD input object.  Specify what data to read in
107   // using the AddLoad member function.   Sub-classes should at a
108   // minimum overload the member function Event.   A full job can be
109   // executed using the member function Run. 
110 }
111
112   
113
114 //____________________________________________________________________
115 AliFMDInput::AliFMDInput(const char* gAliceFile)
116   : TNamed("AliFMDInput", "Input handler for various FMD data"), 
117     fGAliceFile(gAliceFile),
118     fLoader(0),
119     fRun(0), 
120     fStack(0),
121     fFMDLoader(0), 
122     fReader(0),
123     fFMDReader(0),
124     fFMD(0),
125     fESD(0),
126     fESDEvent(0),
127     fTreeE(0),
128     fTreeH(0),
129     fTreeD(0),
130     fTreeS(0),
131     fTreeR(0), 
132     fTreeA(0), 
133     fChainE(0),
134     fArrayE(0),
135     fArrayH(0),
136     fArrayD(0),
137     fArrayS(0), 
138     fArrayR(0), 
139     fArrayA(0), 
140     fHeader(0),
141     fGeoManager(0),
142     fTreeMask(0), 
143     fRawFile(""),
144     fIsInit(kFALSE),
145     fEventCount(0)
146 {
147   
148   // Constructor of an FMD input object.  Specify what data to read in
149   // using the AddLoad member function.   Sub-classes should at a
150   // minimum overload the member function Event.   A full job can be
151   // executed using the member function Run. 
152 }
153
154 //____________________________________________________________________
155 Int_t
156 AliFMDInput::NEvents() const 
157 {
158   // Get number of events
159   if (TESTBIT(fTreeMask, kRaw)) return -1;
160   if (fTreeE) return fTreeE->GetEntries();
161   return -1;
162 }
163
164 //____________________________________________________________________
165 Bool_t
166 AliFMDInput::Init()
167 {
168   // Initialize the object.  Get the needed loaders, and such. 
169
170   // Check if we have been initialized
171   if (fIsInit) { 
172     AliWarning("Already initialized");
173     return fIsInit;
174   }
175   // Get the run 
176   if (!TESTBIT(fTreeMask, kRaw)) { 
177     if (fGAliceFile.IsNull()) fGAliceFile = "galice.root";
178     // Get the loader
179     fLoader = AliRunLoader::Open(fGAliceFile.Data(), "Alice", "read");
180     if (!fLoader) {
181       AliError(Form("Coulnd't read the file %s", fGAliceFile.Data()));
182       return kFALSE;
183     }
184   
185     if  (fLoader->LoadgAlice()) return kFALSE;
186     fRun = fLoader->GetAliRun();
187   
188     // Get the FMD 
189     fFMD = static_cast<AliFMD*>(fRun->GetDetector("FMD"));
190     if (!fFMD) {
191       AliError("Failed to get detector FMD from loader");
192       return kFALSE;
193     }
194   
195     // Get the FMD loader
196     fFMDLoader = fLoader->GetLoader("FMDLoader");
197     if (!fFMDLoader) {
198       AliError("Failed to get detector FMD loader from loader");
199       return kFALSE;
200     }
201     if (fLoader->LoadHeader()) { 
202       AliError("Failed to get event header information from loader");
203       return kFALSE;
204     }
205     fTreeE = fLoader->TreeE();
206   }
207
208   // Optionally, get the ESD files
209   if (TESTBIT(fTreeMask, kESD)) {
210     fChainE = new TChain("esdTree");
211     TSystemDirectory dir(".",".");
212     TList*           files = dir.GetListOfFiles();
213     TSystemFile*     file = 0;
214     if (!files) {
215       AliError("No files");
216       return kFALSE;
217     }
218     files->Sort();
219     TIter            next(files);
220     while ((file = static_cast<TSystemFile*>(next()))) {
221       TString fname(file->GetName());
222       if (fname.Contains("AliESDs")) fChainE->AddFile(fname.Data());
223     }
224     fESDEvent = new AliESDEvent();
225     fESDEvent->ReadFromTree(fChainE);
226     //    fChainE->SetBranchAddress("ESD", &fMainESD);
227     
228   }
229     
230   if (TESTBIT(fTreeMask, kRaw)) {
231     AliInfo("Getting FMD raw data digits");
232     fArrayA = new TClonesArray("AliFMDDigit");
233 #if 0
234     if (!fRawFile.IsNull() && fRawFile.EndsWith(".root"))
235       fReader = new AliRawReaderRoot(fRawFile.Data());
236     else if (!fRawFile.IsNull() && fRawFile.EndsWith(".raw"))
237       fReader = new AliRawReaderDate(fRawFile.Data());
238     else
239       fReader = new AliRawReaderFile(-1);
240 #else
241     if(!fRawFile.IsNull()) 
242       fReader = AliRawReader::Create(fRawFile.Data());
243     else 
244       fReader = new AliRawReaderFile(-1);
245 #endif
246     fFMDReader = new AliFMDRawReader(fReader, 0);
247   }
248   
249   // Optionally, get the geometry 
250   if (TESTBIT(fTreeMask, kGeometry)) {
251     if (fRun) {
252       TString fname(fRun->GetGeometryFileName());
253       if (fname.IsNull()) {
254         Warning("Init", "No file name for the geometry from AliRun");
255         fname = gSystem->DirName(fGAliceFile);
256         fname.Append("/geometry.root");
257       }
258       fGeoManager = TGeoManager::Import(fname.Data());
259       if (!fGeoManager) {
260         Fatal("Init", "No geometry manager found");
261         return kFALSE;
262       }
263     }
264     else { 
265       AliGeomManager::LoadGeometry();
266     }
267     AliCDBManager* cdb   = AliCDBManager::Instance();
268     AliCDBEntry*   align = cdb->Get("FMD/Align/Data");
269     if (align) {
270       AliInfo("Got alignment data from CDB");
271       TClonesArray* array = dynamic_cast<TClonesArray*>(align->GetObject());
272       if (!array) {
273         AliWarning("Invalid align data from CDB");
274       }
275       else {
276         Int_t nAlign = array->GetEntries();
277         for (Int_t i = 0; i < nAlign; i++) {
278           AliAlignObjParams* a = static_cast<AliAlignObjParams*>(array->At(i));
279           if (!a->ApplyToGeometry()) {
280             AliWarning(Form("Failed to apply alignment to %s", 
281                             a->GetSymName()));
282           }
283         }
284       }
285     }
286   }
287
288   fEventCount = 0;
289   fIsInit = kTRUE;
290   return fIsInit;
291 }
292
293 //____________________________________________________________________
294 Bool_t
295 AliFMDInput::Begin(Int_t event)
296 {
297   // Called at the begining of each event.  Per default, it gets the
298   // data trees and gets pointers to the output arrays.   Users can
299   // overload this, but should call this member function in the
300   // overloaded member function of the derived class. 
301
302   // Check if we have been initialized
303   if (!fIsInit) { 
304     AliError("Not initialized");
305     return fIsInit;
306   }
307
308   // Get the event 
309   if (fLoader && fLoader->GetEvent(event)) return kFALSE;
310   AliInfo(Form("Now in event %8d/%8d", event, NEvents()));
311
312   // Possibly load global kinematics information 
313   if (TESTBIT(fTreeMask, kKinematics) || TESTBIT(fTreeMask, kTracks)) {
314     // AliInfo("Getting kinematics");
315     if (fLoader->LoadKinematics("READ")) return kFALSE;
316     fStack = fLoader->Stack();
317   }
318
319   // Possibly load FMD Hit information 
320   if (TESTBIT(fTreeMask, kHits) || TESTBIT(fTreeMask, kTracks)) {
321     // AliInfo("Getting FMD hits");
322     if (!fFMDLoader || fFMDLoader->LoadHits("READ")) return kFALSE;
323     fTreeH = fFMDLoader->TreeH();
324     if (!fArrayH) fArrayH = fFMD->Hits(); 
325   }
326
327   // Possibly load heaedr information 
328   if (TESTBIT(fTreeMask, kHeader)) {
329     // AliInfo("Getting FMD hits");
330     if (!fLoader /* || fLoader->LoadHeader()*/) return kFALSE;
331     fHeader = fLoader->GetHeader();
332   }
333
334   // Possibly load FMD Digit information 
335   if (TESTBIT(fTreeMask, kDigits)) {
336     // AliInfo("Getting FMD digits");
337     if (!fFMDLoader || fFMDLoader->LoadDigits("READ")) return kFALSE;
338     fTreeD = fFMDLoader->TreeD();
339     if (fTreeD) {
340       if (!fArrayD) fArrayD = fFMD->Digits();
341     }
342     else {
343       fArrayD = 0;
344       AliWarning(Form("Failed to load FMD Digits"));
345     } 
346   }
347
348   // Possibly load FMD Sdigit information 
349   if (TESTBIT(fTreeMask, kSDigits)) {
350     // AliInfo("Getting FMD summable digits");
351     if (!fFMDLoader || fFMDLoader->LoadSDigits("READ")) return kFALSE;
352     fTreeS = fFMDLoader->TreeS();
353     if (!fArrayS) fArrayS = fFMD->SDigits();
354   }
355
356   // Possibly load FMD RecPoints information 
357   if (TESTBIT(fTreeMask, kRecPoints)) {
358     // AliInfo("Getting FMD reconstructed points");
359     if (!fFMDLoader || fFMDLoader->LoadRecPoints("READ")) return kFALSE;
360     fTreeR = fFMDLoader->TreeR();
361     if (!fArrayR) fArrayR = new TClonesArray("AliFMDRecPoint");
362     fTreeR->SetBranchAddress("FMD",  &fArrayR);
363   }  
364
365   // Possibly load FMD ESD information 
366   if (TESTBIT(fTreeMask, kESD)) {
367     // AliInfo("Getting FMD event summary data");
368     Int_t read = fChainE->GetEntry(event);
369     if (read <= 0) return kFALSE;
370     fESD = fESDEvent->GetFMDData();
371     if (!fESD) return kFALSE;
372 #if 0
373     TFile* f = fChainE->GetFile();
374     if (f) {
375       TObject* o = f->GetStreamerInfoList()->FindObject("AliFMDMap");
376       if (o) {
377         TStreamerInfo* info = static_cast<TStreamerInfo*>(o);
378         std::cout << "AliFMDMap class version read is " 
379                   <<  info->GetClassVersion() << std::endl;
380       }
381     }
382     // fESD->CheckNeedUShort(fChainE->GetFile());
383 #endif
384   }
385
386   // Possibly load FMD Digit information 
387   if (TESTBIT(fTreeMask, kRaw)) {
388     // AliInfo("Getting FMD raw data digits");
389     if (!fReader->NextEvent()) return kFALSE;
390     // AliFMDRawReader r(fReader, 0);
391     fArrayA->Clear();
392     fFMDReader->ReadAdcs(fArrayA);
393     AliFMDDebug(1, ("Got a total of %d digits", fArrayA->GetEntriesFast()));
394   }
395   fEventCount++;
396   return kTRUE;
397 }
398
399
400 //____________________________________________________________________
401 Bool_t 
402 AliFMDInput::Event()
403 {
404   // Process one event.  The default implementation one or more of 
405   //
406   //   -  ProcessHits       if the hits are loaded. 
407   //   -  ProcessDigits     if the digits are loaded. 
408   //   -  ProcessSDigits    if the sumbable digits are loaded. 
409   //   -  ProcessRecPoints  if the reconstructed points are loaded. 
410   //   -  ProcessESD        if the event summary data is loaded
411   // 
412   if (TESTBIT(fTreeMask, kHits)) 
413     if (!ProcessHits()) return kFALSE; 
414   if (TESTBIT(fTreeMask, kTracks)) 
415     if (!ProcessTracks()) return kFALSE; 
416   if (TESTBIT(fTreeMask, kDigits)) 
417     if (!ProcessDigits()) return kFALSE;
418   if (TESTBIT(fTreeMask, kSDigits)) 
419     if (!ProcessSDigits()) return kFALSE;
420   if (TESTBIT(fTreeMask, kRaw)) 
421     if (!ProcessRawDigits()) return kFALSE;
422   if (TESTBIT(fTreeMask, kRecPoints)) 
423     if (!ProcessRecPoints()) return kFALSE;
424   if (TESTBIT(fTreeMask, kESD))
425     if (!ProcessESDs()) return kFALSE;
426   
427   return kTRUE;
428 }
429
430 //____________________________________________________________________
431 Bool_t 
432 AliFMDInput::ProcessHits()
433 {
434   // Read the hit tree, and pass each hit to the member function
435   // ProcessHit.
436   if (!fTreeH) {
437     AliError("No hit tree defined");
438     return kFALSE;
439   }
440   if (!fArrayH) {
441     AliError("No hit array defined");
442     return kFALSE;
443   }
444
445   Int_t nTracks = fTreeH->GetEntries();
446   for (Int_t i = 0; i < nTracks; i++) {
447     Int_t hitRead  = fTreeH->GetEntry(i);
448     if (hitRead <= 0) continue;
449
450     Int_t nHit = fArrayH->GetEntries();
451     if (nHit <= 0) continue;
452   
453     for (Int_t j = 0; j < nHit; j++) {
454       AliFMDHit* hit = static_cast<AliFMDHit*>(fArrayH->At(j));
455       if (!hit) continue;
456
457       TParticle* track = 0;
458       if (TESTBIT(fTreeMask, kKinematics) && fStack) {
459         Int_t trackno = hit->Track();
460         track = fStack->Particle(trackno);
461       }
462       if (!ProcessHit(hit, track)) return kFALSE;
463     }    
464   }
465   return kTRUE;
466 }
467
468 //____________________________________________________________________
469 Bool_t 
470 AliFMDInput::ProcessTracks()
471 {
472   // Read the hit tree, and pass each hit to the member function
473   // ProcessTrack.
474   if (!fStack) {
475     AliError("No track tree defined");
476     return kFALSE;
477   }
478   if (!fTreeH) {
479     AliError("No hit tree defined");
480     return kFALSE;
481   }
482   if (!fArrayH) {
483     AliError("No hit array defined");
484     return kFALSE;
485   }
486
487   // Int_t nTracks = fStack->GetNtrack();
488   Int_t nTracks = fTreeH->GetEntries();
489   for (Int_t i = 0; i < nTracks; i++) {
490     Int_t      trackno = nTracks - i - 1;
491     TParticle* track   = fStack->Particle(trackno);
492     if (!track) continue;
493
494     // Get the hits for this track. 
495     Int_t hitRead  = fTreeH->GetEntry(i);
496     Int_t nHit     = fArrayH->GetEntries();
497     if (nHit == 0 || hitRead <= 0) { 
498       // Let user code see the track, even if there's no hits. 
499       if (!ProcessTrack(trackno, track, 0)) return kFALSE;
500       continue;
501     }
502
503     // Loop over the hits corresponding to this track.
504     for (Int_t j = 0; j < nHit; j++) {
505       AliFMDHit* hit = static_cast<AliFMDHit*>(fArrayH->At(j));
506       if (!ProcessTrack(trackno, track, hit)) return kFALSE;
507     }   
508   }
509   return kTRUE;
510 }
511
512 //____________________________________________________________________
513 Bool_t 
514 AliFMDInput::ProcessDigits()
515 {
516   // Read the digit tree, and pass each digit to the member function
517   // ProcessDigit.
518   if (!fTreeD) {
519     AliError("No digit tree defined");
520     return kFALSE;
521   }
522   if (!fArrayD) {
523     AliError("No digit array defined");
524     return kFALSE;
525   }
526
527   Int_t nEv = fTreeD->GetEntries();
528   for (Int_t i = 0; i < nEv; i++) {
529     Int_t digitRead  = fTreeD->GetEntry(i);
530     if (digitRead <= 0) continue;
531     Int_t nDigit = fArrayD->GetEntries();
532     AliFMDDebug(0, ("Got %5d digits for this event", nDigit));
533     if (nDigit <= 0) continue;
534     for (Int_t j = 0; j < nDigit; j++) {
535       AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayD->At(j));
536       if (!digit) continue;
537       if (!ProcessDigit(digit)) return kFALSE;
538     }    
539   }
540   return kTRUE;
541 }
542
543 //____________________________________________________________________
544 Bool_t 
545 AliFMDInput::ProcessSDigits()
546 {
547   // Read the summable digit tree, and pass each sumable digit to the
548   // member function ProcessSdigit.
549   if (!fTreeS) {
550     AliWarning("No sdigit tree defined");
551     return kTRUE; // Empty SDigits is fine
552   }
553   if (!fArrayS) {
554     AliWarning("No sdigit array defined");
555     return kTRUE; // Empty SDigits is fine
556   }
557
558   Int_t nEv = fTreeS->GetEntries();
559   for (Int_t i = 0; i < nEv; i++) {
560     Int_t sdigitRead  = fTreeS->GetEntry(i);
561     if (sdigitRead <= 0) continue;
562     Int_t nSdigit = fArrayS->GetEntries();
563     AliFMDDebug(0, ("Got %5d digits for this event", nSdigit));
564     if (nSdigit <= 0) continue;
565     for (Int_t j = 0; j < nSdigit; j++) {
566       AliFMDSDigit* sdigit = static_cast<AliFMDSDigit*>(fArrayS->At(j));
567       if (!sdigit) continue;
568       if (!ProcessSDigit(sdigit)) return kFALSE;
569     }    
570   }
571   return kTRUE;
572 }
573
574 //____________________________________________________________________
575 Bool_t 
576 AliFMDInput::ProcessRawDigits()
577 {
578   // Read the digit tree, and pass each digit to the member function
579   // ProcessDigit.
580   if (!fArrayA) {
581     AliError("No raw digit array defined");
582     return kFALSE;
583   }
584
585   Int_t nDigit = fArrayA->GetEntries();
586   if (nDigit <= 0) return kTRUE;
587   for (Int_t j = 0; j < nDigit; j++) {
588     AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayA->At(j));
589     if (!digit) continue;
590     if (AliLog::GetDebugLevel("FMD","") >= 40 && j < 30) 
591       digit->Print();
592     if (!ProcessRawDigit(digit)) return kFALSE;
593   }    
594   return kTRUE;
595 }
596
597 //____________________________________________________________________
598 Bool_t 
599 AliFMDInput::ProcessRecPoints()
600 {
601   // Read the reconstrcted points tree, and pass each reconstruction
602   // object (AliFMDRecPoint) to either ProcessRecPoint.
603   if (!fTreeR) {
604     AliError("No recpoint tree defined");
605     return kFALSE;
606   }
607   if (!fArrayR) {
608     AliError("No recpoints array defined");
609     return kFALSE;
610   }
611
612   Int_t nEv = fTreeR->GetEntries();
613   for (Int_t i = 0; i < nEv; i++) {
614     Int_t recRead  = fTreeR->GetEntry(i);
615     if (recRead <= 0) continue;
616     Int_t nRecPoint = fArrayR->GetEntries();
617     for (Int_t j = 0; j < nRecPoint; j++) {
618       AliFMDRecPoint* recPoint = static_cast<AliFMDRecPoint*>(fArrayR->At(j));
619       if (!recPoint) continue;
620       if (!ProcessRecPoint(recPoint)) return kFALSE;
621     }    
622   }
623   return kTRUE;
624 }
625
626 //____________________________________________________________________
627 Bool_t 
628 AliFMDInput::ProcessESDs()
629 {
630   // Process event summary data
631   if (!fESD) return kFALSE;
632   for (UShort_t det = 1; det <= 3; det++) {
633     Char_t rings[] = { 'I', (det == 1 ? '\0' : 'O'), '\0' };
634     for (Char_t* rng = rings; *rng != '\0'; rng++) {
635       UShort_t nsec = (*rng == 'I' ?  20 :  40);
636       UShort_t nstr = (*rng == 'I' ? 512 : 256);
637       for (UShort_t sec = 0; sec < nsec; sec++) {
638         for (UShort_t str = 0; str < nstr; str++) {
639           Float_t eta  = fESD->Eta(det,*rng,sec,str);
640           Float_t mult = fESD->Multiplicity(det,*rng,sec,str);
641           if (!fESD->IsAngleCorrected()) 
642             mult *= TMath::Abs(TMath::Cos(2.*TMath::ATan(TMath::Exp(-eta))));
643           if (!ProcessESD(det, *rng, sec, str, eta, mult)) continue;
644         }
645       }
646     }
647   }
648   return kTRUE;
649 }
650
651 //____________________________________________________________________
652 Bool_t
653 AliFMDInput::End()
654 {
655   // Called at the end of each event.  Per default, it unloads the
656   // data trees and resets the pointers to the output arrays.   Users
657   // can overload this, but should call this member function in the
658   // overloaded member function of the derived class. 
659
660   // Check if we have been initialized
661   if (!fIsInit) { 
662     AliError("Not initialized");
663     return fIsInit;
664   }
665   // Possibly unload global kinematics information 
666   if (TESTBIT(fTreeMask, kKinematics) || TESTBIT(fTreeMask, kTracks)) {
667     fLoader->UnloadKinematics();
668     // fTreeK = 0;
669     fStack = 0;
670   }
671   // Possibly unload FMD Hit information 
672   if (TESTBIT(fTreeMask, kHits) || TESTBIT(fTreeMask, kTracks)) {
673     fFMDLoader->UnloadHits();
674     fTreeH = 0;
675   }
676   // Possibly unload FMD Digit information 
677   if (TESTBIT(fTreeMask, kDigits)) {
678     fFMDLoader->UnloadDigits();
679     fTreeD = 0;
680   }
681   // Possibly unload FMD Sdigit information 
682   if (TESTBIT(fTreeMask, kSDigits)) {
683     fFMDLoader->UnloadSDigits();
684     fTreeS = 0;
685   }
686   // Possibly unload FMD RecPoints information 
687   if (TESTBIT(fTreeMask, kRecPoints)) {
688     fFMDLoader->UnloadRecPoints();
689     fTreeR = 0;
690   }
691   // AliInfo("Now out event");
692   return kTRUE;
693 }
694
695 //____________________________________________________________________
696 Bool_t
697 AliFMDInput::Run()
698 {
699   // Run over all events and files references in galice.root 
700
701   Bool_t retval;
702   if (!(retval = Init())) return retval;
703
704   Int_t nEvents = NEvents();
705   for (Int_t event = 0; nEvents < 0 || event < nEvents; event++) {
706     if (!(retval = Begin(event))) break;
707     if (!(retval = Event())) break;
708     if (!(retval = End())) break;
709   }
710   if (!retval) return retval;
711   retval = Finish();
712   return retval;
713 }
714
715 //__________________________________________________________________
716 TArrayF 
717 AliFMDInput::MakeLogScale(Int_t n, Double_t min, Double_t max) 
718 {
719   // Service function to define a logarithmic axis. 
720   // Parameters: 
721   //   n    Number of bins 
722   //   min  Minimum of axis 
723   //   max  Maximum of axis 
724   TArrayF bins(n+1);
725   bins[0]      = min;
726   if (n <= 20) {
727     for (Int_t i = 1; i < n+1; i++) bins[i] = bins[i-1] + (max-min)/n;
728     return bins;
729   }
730   Float_t dp   = n / TMath::Log10(max / min);
731   Float_t pmin = TMath::Log10(min);
732   for (Int_t i = 1; i < n+1; i++) {
733     Float_t p = pmin + i / dp;
734     bins[i]   = TMath::Power(10, p);
735   }
736   return bins;
737 }
738
739
740
741 //____________________________________________________________________
742 //
743 // EOF
744 //