]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALGetter.cxx
7543cda1bfdad762fcb9504f5ed878c1e63e78ba
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALGetter.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 //  A singleton. This class should be used in the analysis stage to get 
20 //  reconstructed objects: Digits, RecPoints, TrackSegments and RecParticles,
21 //  instead of directly reading them from galice.root file. This container 
22 //  ensures, that one reads Digits, made of these particular digits, RecPoints, 
23 //  made of these particular RecPoints, TrackSegments and RecParticles. 
24 //  This becomes non trivial if there are several identical branches, produced with
25 //  different set of parameters. 
26 //
27 //  An example of how to use (see also class AliEMCALAnalyser):
28 //  AliEMCALGetter * gime = AliEMCALGetter::GetInstance("galice.root","test") ;
29 //  for(Int_t irecp = 0; irecp < gime->NRecParticles() ; irecp++)
30 //     AliEMCALRecParticle * part = gime->RecParticle(1) ;
31 //     ................
32 //  gime->Event(event) ;    // reads new event from galice.root
33 //                  
34 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC KI & SUBATECH)
35 //*--         Completely redesigned by Dmitri Peressounko March 2001  
36 //
37 //*-- YS June 2001 : renamed the original AliEMCALIndexToObject and make
38 //*--         systematic usage of TFolders without changing the interface        
39 //////////////////////////////////////////////////////////////////////////////
40
41 // --- ROOT system ---
42
43 #include <TFile.h>
44 #include <TROOT.h>
45 #include <TSystem.h>
46 #include <TF1.h>
47 #include <TGraph.h>
48 //#include <TCanvas.h>
49 //#include <TFrame.h>
50
51 // --- Standard library ---
52
53 // --- AliRoot header files ---
54
55 #include "AliEMCAL.h"
56 #include "AliEMCALGetter.h"
57 #include "AliEMCALLoader.h"
58 #include "AliHeader.h"  
59 #include "AliMC.h"
60 #include "AliRunLoader.h"
61 #include "AliStack.h"  
62 #include "AliEMCALRawStream.h"
63 #include "AliRawReaderFile.h"
64
65 ClassImp(AliEMCALGetter)
66   
67 AliEMCALGetter * AliEMCALGetter::fgObjGetter = 0 ; 
68 AliEMCALLoader * AliEMCALGetter::fgEmcalLoader = 0;
69 Int_t AliEMCALGetter::fgDebug = 0;
70
71 //  TFile * AliEMCALGetter::fgFile = 0 ; 
72
73 //____________________________________________________________________________ 
74 AliEMCALGetter::AliEMCALGetter(const char* headerFile, const char* version, Option_t * openingOption) {
75   AliRunLoader* rl = AliRunLoader::GetRunLoader(version) ; 
76   if (!rl) {
77     rl = AliRunLoader::Open(headerFile, version, openingOption);
78     if (!rl) {
79       Fatal("AliEMCALGetter", "Could not find the Run Loader for %s - %s",headerFile, version) ; 
80       return ;
81     } 
82     if (rl->GetAliRun() == 0x0) {
83       rl->LoadgAlice();
84       gAlice = rl->GetAliRun(); // should be removed
85     }
86   } 
87   fgEmcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetLoader("EMCALLoader"));
88   if ( !fgEmcalLoader ) 
89     Error("AliEMCALGetter", "Could not find EMCALLoader") ; 
90   else 
91     fgEmcalLoader->SetTitle(version);
92
93    // initialize data members
94   SetDebug(0) ; 
95   fLoadingStatus = "" ; 
96
97 }
98
99 //____________________________________________________________________________ 
100 AliEMCALGetter::~AliEMCALGetter()
101 {
102   if (fgEmcalLoader) {
103     delete fgEmcalLoader;
104     fgEmcalLoader = 0 ;
105   }
106   fgObjGetter = 0;
107 }
108
109
110 //____________________________________________________________________________ 
111 void AliEMCALGetter::Reset()
112 {
113   // resets things in case the getter is called consecutively with different files
114   // the EMCAL Loader is already deleted by the Run Loader
115   fgEmcalLoader = 0;
116   fgObjGetter = 0;
117
118 }
119
120 //____________________________________________________________________________ 
121 AliEMCALClusterizer * AliEMCALGetter::Clusterizer()
122
123   // return pointer to Clusterizer Tree
124   AliEMCALClusterizer * rv ; 
125   rv =  dynamic_cast<AliEMCALClusterizer *>(EmcalLoader()->Reconstructioner()) ;
126   if (!rv) {
127     Event(0, "R") ; 
128     rv =  dynamic_cast<AliEMCALClusterizer*>(EmcalLoader()->Reconstructioner()) ;
129   }
130   return rv ; 
131 }
132
133
134 //____________________________________________________________________________ 
135 TClonesArray * AliEMCALGetter::Digits() const  
136 {
137   // asks the Loader to return the Digits container 
138
139   TClonesArray * rv = 0 ; 
140   rv = EmcalLoader()->Digits() ; 
141
142   if( !rv ) {
143     EmcalLoader()->MakeDigitsArray() ; 
144     rv = EmcalLoader()->Digits() ;
145   }
146   return rv ; 
147 }
148
149 //____________________________________________________________________________ 
150 AliEMCALDigitizer * AliEMCALGetter::Digitizer() 
151
152   // return pointer to Digitizer Tree
153   AliEMCALDigitizer * rv ; 
154   rv =  dynamic_cast<AliEMCALDigitizer *>(EmcalLoader()->Digitizer()) ;
155   if (!rv) {
156     Event(0, "D") ; 
157     rv =  dynamic_cast<AliEMCALDigitizer *>(EmcalLoader()->Digitizer()) ;
158   }
159   return rv ; 
160 }
161
162 //____________________________________________________________________________ 
163 TObjArray * AliEMCALGetter::ECARecPoints() const 
164 {
165   // asks the Loader to return the EMC RecPoints container 
166
167   TObjArray * rv = 0 ; 
168   
169   rv = EmcalLoader()->ECARecPoints() ; 
170   if (!rv) {
171     EmcalLoader()->MakeRecPointsArray() ;
172     rv = EmcalLoader()->ECARecPoints() ; 
173   }
174   return rv ; 
175 }
176
177 //____________________________________________________________________________ 
178 TClonesArray * AliEMCALGetter::TrackSegments() const 
179 {
180   // asks the Loader to return the TrackSegments container 
181
182   TClonesArray * rv = 0 ; 
183   
184   rv = EmcalLoader()->TrackSegments() ; 
185   if (!rv) {
186     EmcalLoader()->MakeTrackSegmentsArray() ;
187     rv = EmcalLoader()->TrackSegments() ; 
188   }
189   return rv ; 
190 }
191
192 //____________________________________________________________________________ 
193 AliEMCALTrackSegmentMaker * AliEMCALGetter::TrackSegmentMaker() 
194
195   // return pointer to TrackSegmentMaker Tree
196   AliEMCALTrackSegmentMaker * rv ; 
197   rv =  dynamic_cast<AliEMCALTrackSegmentMaker *>(EmcalLoader()->TrackSegmentMaker()) ;
198   if (!rv) {
199     Event(0, "T") ; 
200     rv =  dynamic_cast<AliEMCALTrackSegmentMaker *>(EmcalLoader()->TrackSegmentMaker()) ;
201   }
202   return rv ; 
203 }
204
205 //____________________________________________________________________________ 
206 TClonesArray * AliEMCALGetter::RecParticles() const  
207 {
208   // asks the Loader to return the TrackSegments container 
209
210   TClonesArray * rv = 0 ; 
211   
212   rv = EmcalLoader()->RecParticles() ; 
213   if (!rv) {
214     EmcalLoader()->MakeRecParticlesArray() ;
215     rv = EmcalLoader()->RecParticles() ; 
216   }
217   return rv ; 
218 }
219 //____________________________________________________________________________ 
220 void AliEMCALGetter::Event(Int_t event, const char* opt) 
221 {
222   // Reads the content of all Tree's S, D and R
223
224 //   if ( event >= MaxEvent() ) {
225 //     Error("Event", "%d not found in TreeE !", event) ; 
226 //     return ; 
227 //   }
228
229   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
230   // checks if we are dealing with test-beam data
231 //   TBranch * btb = rl->TreeE()->GetBranch("AliEMCALBeamTestEvent") ;
232 //   if(btb){
233 //     if(!fBTE)
234 //       fBTE = new AliEMCALBeamTestEvent() ;
235 //     btb->SetAddress(&fBTE) ;
236 //     btb->GetEntry(event) ;
237 //   }
238 //   else{
239 //     if(fBTE){
240 //       delete fBTE ;
241 //       fBTE = 0 ;
242 //     }
243 //   }
244
245   // Loads the type of object(s) requested
246   
247   rl->GetEvent(event) ;
248
249   if( strstr(opt,"X") || (strcmp(opt,"")==0) )
250     ReadPrimaries() ;
251
252   if(strstr(opt,"H") )
253     ReadTreeH();
254
255   if(strstr(opt,"S") )
256     ReadTreeS() ;
257
258   if( strstr(opt,"D") )
259     ReadTreeD() ;
260
261   if( strstr(opt,"R") )
262     ReadTreeR() ;
263
264   if( strstr(opt,"T") )
265     ReadTreeT() ;
266
267   if( strstr(opt,"P") )
268     ReadTreeP() ;
269
270   if( strstr(opt,"W") )
271     ReadRaw(event) ;
272   
273 }
274
275
276 //____________________________________________________________________________ 
277 Int_t AliEMCALGetter::EventNumber() const
278   {
279   // return the current event number
280   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
281   return static_cast<Int_t>(rl->GetEventNumber()) ;   
282 }
283
284 //____________________________________________________________________________ 
285   TClonesArray * AliEMCALGetter::Hits() const   
286 {
287   // asks the loader to return  the Hits container 
288   
289   TClonesArray * rv = 0 ; 
290   
291   rv = EmcalLoader()->Hits() ; 
292   if ( !rv ) {
293     EmcalLoader()->LoadHits("read"); 
294     rv = EmcalLoader()->Hits() ; 
295   }
296   return rv ; 
297 }
298
299 //____________________________________________________________________________ 
300 AliEMCALGetter * AliEMCALGetter::Instance(const char* alirunFileName, const char* version, Option_t * openingOption) 
301 {
302   // Creates and returns the pointer of the unique instance
303   // Must be called only when the environment has changed
304
305   if(!fgObjGetter){ // first time the getter is called 
306     fgObjGetter = new AliEMCALGetter(alirunFileName, version, openingOption) ;
307   }
308   else { // the getter has been called previously
309     AliRunLoader * rl = AliRunLoader::GetRunLoader(fgEmcalLoader->GetTitle());
310     if ( rl->GetFileName() == alirunFileName ) {// the alirunFile has the same name
311       // check if the file is already open
312       TFile * galiceFile = dynamic_cast<TFile *>(gROOT->FindObject(rl->GetFileName()) ) ; 
313       
314       if ( !galiceFile ) 
315         fgObjGetter = new AliEMCALGetter(alirunFileName, version, openingOption);
316       
317       else {  // the file is already open check the version name
318         TString currentVersionName = rl->GetEventFolder()->GetName() ; 
319         TString newVersionName(version) ; 
320         if (currentVersionName == newVersionName) 
321           if(fgDebug)
322             ::Warning( "Instance", "Files with version %s already open", currentVersionName.Data() ) ;  
323         else {    
324           fgObjGetter = new AliEMCALGetter(alirunFileName, version, openingOption) ; 
325         }
326       }
327     }
328     else { 
329       AliRunLoader * rl = AliRunLoader::GetRunLoader(fgEmcalLoader->GetTitle());
330       if ( strstr(version, AliConfig::GetDefaultEventFolderName()) ) // false in case of merging
331         delete rl ; 
332       fgObjGetter = new AliEMCALGetter(alirunFileName, version, openingOption) ;      
333     }
334   }
335   if (!fgObjGetter) 
336     ::Error("AliEMCALGetter::Instance", "Failed to create the EMCAL Getter object") ;
337   else 
338     if (fgDebug)
339       Print() ;
340
341   return fgObjGetter ;
342 }
343
344 //____________________________________________________________________________ 
345 AliEMCALGetter *  AliEMCALGetter::Instance()
346 {
347   // Returns the pointer of the unique instance already defined
348   
349   if(!fgObjGetter && fgDebug)
350      ::Warning("AliEMCALGetter::Instance", "Getter not initialized") ;
351
352    return fgObjGetter ;
353            
354 }
355
356 //____________________________________________________________________________ 
357 Int_t AliEMCALGetter::MaxEvent() const 
358 {
359   // returns the number of events in the run (from TE)
360
361   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
362   return static_cast<Int_t>(rl->GetNumberOfEvents()) ; 
363 }
364
365 //____________________________________________________________________________ 
366 TParticle * AliEMCALGetter::Primary(Int_t index) const
367 {
368   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
369   return rl->Stack()->Particle(index) ; 
370
371
372 //____________________________________________________________________________ 
373 Int_t AliEMCALGetter::NPrimaries() const
374 {
375   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
376   return (rl->GetHeader())->GetNtrack(); 
377
378
379 //____________________________________________________________________________ 
380 AliEMCAL * AliEMCALGetter:: EMCAL() const  
381 {
382   // returns the EMCAL object 
383   AliEMCALLoader *    loader = 0;
384   static AliEMCALLoader * oldloader = 0;
385   static AliEMCAL * emcal = 0;
386
387   loader = EmcalLoader();
388
389   if (loader != oldloader ) {
390     emcal = dynamic_cast<AliEMCAL*>(loader->GetModulesFolder()->FindObject("EMCAL")) ;  
391     oldloader = loader;
392   }
393   if (!emcal) 
394     if (fgDebug)
395       Warning("EMCAL", "EMCAL module not found in module folders: %s", EmcalLoader()->GetModulesFolder()->GetName() ) ; 
396   return emcal ; 
397 }  
398
399
400
401 //____________________________________________________________________________ 
402 AliEMCALPID * AliEMCALGetter::PID() 
403
404   // return pointer to PID Tree
405   AliEMCALPID * rv ; 
406   rv =  dynamic_cast<AliEMCALPID *>(EmcalLoader()->PIDTask()) ;
407   if (!rv) {
408     Event(0, "P") ; 
409     rv =  dynamic_cast<AliEMCALPID *>(EmcalLoader()->PIDTask()) ;
410   }
411   return rv ; 
412 }
413
414 //____________________________________________________________________________ 
415 AliEMCALGeometry * AliEMCALGetter::EMCALGeometry() const 
416 {
417   // Returns EMCAL geometry
418
419   AliEMCALGeometry * rv = 0 ; 
420   if (EMCAL() )
421     rv =  EMCAL()->GetGeometry() ;
422   return rv ; 
423
424
425 //____________________________________________________________________________ 
426 void  AliEMCALGetter::Print() 
427 {
428   // Print usefull information about the getter
429     
430   AliRunLoader * rl = AliRunLoader::GetRunLoader(fgEmcalLoader->GetTitle());
431   ::Info("Print", "gAlice file is %s -- version name is %s", (rl->GetFileName()).Data(), rl->GetEventFolder()->GetName() ) ; 
432 }
433
434 //____________________________________________________________________________ 
435 void AliEMCALGetter::ReadPrimaries()  
436 {
437   // Read Primaries from Kinematics.root
438   
439   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
440   
441   // gets kine tree from the root file (Kinematics.root)
442   if ( ! rl->TreeK() )  // load treeK the first time
443     rl->LoadKinematics() ;
444   
445   if (fgDebug) 
446     Info("ReadTreeK", "Found %d particles in event # %d", NPrimaries(), EventNumber() ) ; 
447 }
448
449 //____________________________________________________________________________ 
450 void AliEMCALGetter::FitRaw(Bool_t lowGainFlag, TGraph * gLowGain, TGraph * gHighGain, TF1* signalF, Int_t & amp, Double_t & time)
451 {
452   // Fits the raw signal time distribution 
453
454   const Int_t kNoiseThreshold = 0 ;
455   Double_t timezero1 = 0., timezero2 = 0., timemax = 0. ;
456   Double_t signal = 0., signalmax = 0. ;       
457   Double_t energy = time = 0. ; 
458
459   if (lowGainFlag) {
460     timezero1 = timezero2 = signalmax = timemax = 0. ;
461     signalF->FixParameter(0, EMCAL()->GetRawFormatLowCharge()) ; 
462     signalF->FixParameter(1, EMCAL()->GetRawFormatLowGain()) ; 
463     Int_t index ; 
464     for (index = 0; index < EMCAL()->GetRawFormatTimeBins(); index++) {
465       gLowGain->GetPoint(index, time, signal) ; 
466       if (signal > kNoiseThreshold && timezero1 == 0.) 
467         timezero1 = time ;
468       if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
469         timezero2 = time ; 
470       if (signal > signalmax) {
471         signalmax = signal ; 
472         timemax   = time ; 
473       }
474     }
475     signalmax /= EMCAL()->RawResponseFunctionMax(EMCAL()->GetRawFormatLowCharge(), 
476                                                 EMCAL()->GetRawFormatLowGain()) ;
477     if ( timezero1 + EMCAL()->GetRawFormatTimePeak() < EMCAL()->GetRawFormatTimeMax() * 0.4 ) { // else its noise 
478       signalF->SetParameter(2, signalmax) ; 
479       signalF->SetParameter(3, timezero1) ;         
480       gLowGain->Fit(signalF, "QRON", "", 0., timezero2); //, "QRON") ; 
481       energy = signalF->GetParameter(2) ; 
482       time   = signalF->GetMaximumX() - EMCAL()->GetRawFormatTimePeak() - EMCAL()->GetRawFormatTimeTrigger() ;
483     }
484   } else {
485     timezero1 = timezero2 = signalmax = timemax = 0. ;
486     signalF->FixParameter(0, EMCAL()->GetRawFormatHighCharge()) ; 
487     signalF->FixParameter(1, EMCAL()->GetRawFormatHighGain()) ; 
488     Int_t index ; 
489     for (index = 0; index < EMCAL()->GetRawFormatTimeBins(); index++) {
490       gHighGain->GetPoint(index, time, signal) ;               
491       if (signal > kNoiseThreshold && timezero1 == 0.) 
492         timezero1 = time ;
493       if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
494         timezero2 = time ; 
495       if (signal > signalmax) {
496         signalmax = signal ;   
497         timemax   = time ; 
498       }
499     }
500     signalmax /= EMCAL()->RawResponseFunctionMax(EMCAL()->GetRawFormatHighCharge(), 
501                                                 EMCAL()->GetRawFormatHighGain()) ;;
502     if ( timezero1 + EMCAL()->GetRawFormatTimePeak() < EMCAL()->GetRawFormatTimeMax() * 0.4 ) { // else its noise  
503       signalF->SetParameter(2, signalmax) ; 
504       signalF->SetParameter(3, timezero1) ;               
505       gHighGain->Fit(signalF, "QRON", "", 0., timezero2) ; 
506       energy = signalF->GetParameter(2) ; 
507       time   = signalF->GetMaximumX() - EMCAL()->GetRawFormatTimePeak() - EMCAL()->GetRawFormatTimeTrigger() ;
508     }
509   }
510   
511   if (time == 0. && energy == 0.) 
512     amp = 0 ; 
513   else {
514   AliEMCALDigitizer * digitizer = Digitizer() ; 
515   amp = static_cast<Int_t>( (energy - digitizer->GetECApedestal()) / digitizer->GetECAchannel() + 0.5 ) ; 
516   }
517   // dessin
518 //   TCanvas * c1 = new TCanvas("c1","A Simple Graph Example",200,10,700,500);
519 //   c1->SetFillColor(42);
520 //   c1->SetGrid();
521 //   gLowGain->SetLineColor(2);
522 //   gLowGain->SetLineWidth(4);
523 //   gLowGain->SetMarkerColor(4);
524 //   gLowGain->SetMarkerStyle(21);
525 //   gLowGain->SetTitle("Lowgain");
526 //   gLowGain->GetXaxis()->SetTitle("X title");
527 //   gLowGain->GetYaxis()->SetTitle("Y title");
528 //   gLowGain->Draw("ACP");
529   
530 //   c1->Update();
531 //   c1->GetFrame()->SetFillColor(21);
532 //   c1->GetFrame()->SetBorderSize(12);
533 //   c1->Modified();
534   
535 //   TCanvas * c2 = new TCanvas("c2","A Simple Graph Example",200,10,700,500);
536 //   c2->SetFillColor(42);
537 //   c2->SetGrid();
538 //   gHighGain->SetLineColor(2);
539 //   gHighGain->SetLineWidth(4);
540 //   gHighGain->SetMarkerColor(4);
541 //   gHighGain->SetMarkerStyle(21);
542 //   gHighGain->SetTitle("Highgain");
543 //   gHighGain->GetXaxis()->SetTitle("X title");
544 //   gHighGain->GetYaxis()->SetTitle("Y title");
545 //     gHighGain->Draw("ACP");
546   
547 //   c2->Update();
548 //   c2->GetFrame()->SetFillColor(21);
549 //   c2->GetFrame()->SetBorderSize(12);
550 //   c2->Modified();
551 }
552
553 //____________________________________________________________________________ 
554 Int_t AliEMCALGetter::ReadRaw(Int_t event)
555 {
556   // reads the raw format data, converts it into digits format and store digits in Digits()
557   // container.
558   
559   AliRawReaderFile rawReader(event) ; 
560   AliEMCALRawStream in(&rawReader);
561   
562   Bool_t first = kTRUE ;
563  
564   TF1 * signalF = new TF1("signal", AliEMCAL::RawResponseFunction, 0, EMCAL()->GetRawFormatTimeMax(), 4);
565   signalF->SetParNames("Charge", "Gain", "Amplitude", "TimeZero") ; 
566  
567   
568   Int_t id = -1;
569   Bool_t lowGainFlag = kFALSE ; 
570
571   TClonesArray * digits = Digits() ;
572   digits->Clear() ; 
573   Int_t idigit = 0 ; 
574   Int_t amp = 0 ; 
575   Double_t time = 0. ; 
576
577   TGraph * gLowGain = new TGraph(EMCAL()->GetRawFormatTimeBins()) ; 
578   TGraph * gHighGain= new TGraph(EMCAL()->GetRawFormatTimeBins()) ;  
579
580   while ( in.Next() ) { // EMCAL entries loop 
581     if ( in.IsNewId() ) {
582       if (!first) {
583         FitRaw(lowGainFlag, gLowGain, gHighGain, signalF, amp, time) ; 
584         if (amp > 0) {
585           new((*digits)[idigit]) AliEMCALDigit( -1, -1, id, amp, time) ;        
586           idigit++ ; 
587         }
588         Int_t index ; 
589         for (index = 0; index < EMCAL()->GetRawFormatTimeBins(); index++) {
590           gLowGain->SetPoint(index, index * EMCAL()->GetRawFormatTimeMax() / EMCAL()->GetRawFormatTimeBins(), 0) ;  
591           gHighGain->SetPoint(index, index * EMCAL()->GetRawFormatTimeMax() / EMCAL()->GetRawFormatTimeBins(), 0) ; 
592         } 
593       }  
594       first = kFALSE ; 
595       id = in.GetId() ; 
596       if (in.GetModule() == EMCAL()->GetRawFormatLowGainOffset() ) 
597         lowGainFlag = kTRUE ; 
598       else 
599         lowGainFlag = kFALSE ; 
600     }
601     if (lowGainFlag)
602       gLowGain->SetPoint(in.GetTime(), 
603                          in.GetTime()* EMCAL()->GetRawFormatTimeMax() / EMCAL()->GetRawFormatTimeBins(), 
604                          in.GetSignal()) ;
605     else 
606       gHighGain->SetPoint(in.GetTime(), 
607                           in.GetTime() * EMCAL()->GetRawFormatTimeMax() / EMCAL()->GetRawFormatTimeBins(), 
608                           in.GetSignal() ) ;
609     
610   } // EMCAL entries loop
611   digits->Sort() ; 
612
613   delete signalF ; 
614   delete gLowGain;
615   delete gHighGain ; 
616     
617   return Digits()->GetEntriesFast() ; 
618 }
619   
620   //____________________________________________________________________________ 
621 Int_t AliEMCALGetter::ReadTreeD()
622 {
623   // Read the Digits
624   
625    EmcalLoader()->CleanDigits() ; 
626    EmcalLoader()->LoadDigits("UPDATE") ;
627    EmcalLoader()->LoadDigitizer("UPDATE") ;
628    return Digits()->GetEntries() ; 
629 }
630
631 //____________________________________________________________________________ 
632 Int_t AliEMCALGetter::ReadTreeH()
633 {
634   // Read the Hits
635   EmcalLoader()->CleanHits() ; 
636   // gets TreeH from the root file (EMCAL.Hit.root)
637   //if ( !IsLoaded("H") ) {
638     EmcalLoader()->LoadHits("UPDATE") ;
639     //SetLoaded("H") ; 
640     //}  
641   return Hits()->GetEntries() ; 
642 }
643
644 //____________________________________________________________________________ 
645 Int_t AliEMCALGetter::ReadTreeR()
646 {
647   // Read the RecPoints
648
649    EmcalLoader()->CleanRecPoints() ; 
650   // gets TreeR from the root file (EMCAL.RecPoints.root)
651   //if ( !IsLoaded("R") ) {
652     EmcalLoader()->LoadRecPoints("UPDATE") ;
653     EmcalLoader()->LoadClusterizer("UPDATE") ;
654     //  SetLoaded("R") ; 
655     //}
656
657   return ECARecPoints()->GetEntries() ; 
658 }
659
660 //____________________________________________________________________________ 
661 Int_t AliEMCALGetter::ReadTreeT()
662 {
663   // Read the TrackSegments
664   
665   EmcalLoader()->CleanTracks() ; 
666   // gets TreeT from the root file (EMCAL.TrackSegments.root)
667   //if ( !IsLoaded("T") ) {
668     EmcalLoader()->LoadTracks("UPDATE") ;
669     EmcalLoader()->LoadTrackSegmentMaker("UPDATE") ;
670     //  SetLoaded("T") ; 
671     //}
672
673   return TrackSegments()->GetEntries() ; 
674 }
675 //____________________________________________________________________________ 
676 Int_t AliEMCALGetter::ReadTreeP()
677 {
678   // Read the RecParticles
679   
680   EmcalLoader()->CleanRecParticles() ; 
681   // gets TreeP from the root file (EMCAL.RecParticles.root)
682   //  if ( !IsLoaded("P") ) {
683     EmcalLoader()->LoadRecParticles("UPDATE") ;
684     EmcalLoader()->LoadPID("UPDATE") ;
685     //  SetLoaded("P") ; 
686     // }
687
688   return RecParticles()->GetEntries() ; 
689 }
690 //____________________________________________________________________________ 
691 Int_t AliEMCALGetter::ReadTreeS()
692 {
693   // Read the SDigits
694   
695   //  EmcalLoader()->CleanSDigits() ; 
696   // gets TreeS from the root file (EMCAL.SDigits.root)
697   // if ( !IsLoaded("S") ) {
698     EmcalLoader()->LoadSDigits("READ") ;
699     EmcalLoader()->LoadSDigitizer("READ") ;
700     //  SetLoaded("S") ; 
701     //}
702
703   return SDigits()->GetEntries() ; 
704 }
705
706 //____________________________________________________________________________ 
707 TClonesArray * AliEMCALGetter::SDigits() const  
708 {
709   // asks the Loader to return the Digits container 
710
711   TClonesArray * rv = 0 ; 
712   
713   rv = EmcalLoader()->SDigits() ; 
714   if (!rv) {
715     EmcalLoader()->MakeSDigitsArray() ;
716     rv = EmcalLoader()->SDigits() ; 
717   }
718   return rv ; 
719 }
720
721 //____________________________________________________________________________ 
722 AliEMCALSDigitizer * AliEMCALGetter::SDigitizer() 
723
724   // Return pointer to SDigitizer task
725   AliEMCALSDigitizer * rv ; 
726   rv =  dynamic_cast<AliEMCALSDigitizer *>(EmcalLoader()->SDigitizer()) ;
727   if (!rv) {
728     Event(0, "S") ; 
729     rv =  dynamic_cast<AliEMCALSDigitizer *>(EmcalLoader()->SDigitizer()) ;
730   }
731   return rv ; 
732 }
733
734 //____________________________________________________________________________ 
735 TParticle * AliEMCALGetter::Secondary(const TParticle* p, Int_t index) const
736 {
737   // Return first (index=1) or second (index=2) secondary particle of primary particle p 
738
739   if(index <= 0) 
740     return 0 ;
741   if(index > 2)
742     return 0 ;
743
744   if(p) {
745   Int_t daughterIndex = p->GetDaughter(index-1) ; 
746   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
747   return  rl->GetAliRun()->GetMCApp()->Particle(daughterIndex) ; 
748   }
749   else
750     return 0 ;
751 }
752
753 //____________________________________________________________________________ 
754 void AliEMCALGetter::Track(Int_t itrack) 
755 {
756   // Read the first entry of EMCAL branch in hit tree gAlice->TreeH()
757  
758  AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
759
760   if( !TreeH() ) // load treeH the first time
761     rl->LoadHits() ;
762
763   // first time create the container
764   TClonesArray * hits = Hits() ; 
765   if ( hits ) 
766     hits->Clear() ; 
767
768   TBranch * emcalbranch = dynamic_cast<TBranch*>(TreeH()->GetBranch("EMCAL")) ; 
769   emcalbranch->SetAddress(&hits) ;
770   emcalbranch->GetEntry(itrack) ;
771 }
772
773 //____________________________________________________________________________ 
774 TTree * AliEMCALGetter::TreeD() const 
775 {
776   // return pointer to Digits Tree
777   TTree * rv = 0 ; 
778   rv = EmcalLoader()->TreeD() ; 
779   if ( !rv ) {
780     EmcalLoader()->MakeTree("D");
781     rv = EmcalLoader()->TreeD() ;
782   } 
783   
784   return rv ; 
785 }
786
787 //____________________________________________________________________________ 
788 TTree * AliEMCALGetter::TreeH() const 
789 {
790   // return pointer to Hits Tree
791   TTree * rv = 0 ; 
792   rv = EmcalLoader()->TreeH() ; 
793   if ( !rv ) {
794     EmcalLoader()->MakeTree("H");
795     rv = EmcalLoader()->TreeH() ;
796   } 
797   
798   return rv ; 
799 }
800
801 //____________________________________________________________________________ 
802 TTree * AliEMCALGetter::TreeR() const 
803 {
804   // return pointer to RecPoints Tree
805
806   TTree * rv = 0 ; 
807   rv = EmcalLoader()->TreeR() ; 
808   if ( !rv ) {
809     EmcalLoader()->MakeTree("R");
810     rv = EmcalLoader()->TreeR() ;
811   } 
812   
813   return rv ; 
814 }
815
816 //____________________________________________________________________________ 
817 TTree * AliEMCALGetter::TreeT() const 
818 {  
819   // return pointer to TrackSegments Tree
820   TTree * rv = 0 ; 
821   rv = EmcalLoader()->TreeT() ; 
822   if ( !rv ) {
823     EmcalLoader()->MakeTree("T");
824     rv = EmcalLoader()->TreeT() ;
825   } 
826   
827   return rv ; 
828 }
829 //____________________________________________________________________________ 
830 TTree * AliEMCALGetter::TreeP() const 
831 {
832   // return pointer to RecParticles Tree
833   TTree * rv = 0 ; 
834   rv = EmcalLoader()->TreeP() ; 
835   if ( !rv ) {
836     EmcalLoader()->MakeTree("P");
837     rv = EmcalLoader()->TreeP() ;
838   } 
839   
840   return rv ; 
841 }
842
843 //____________________________________________________________________________ 
844 TTree * AliEMCALGetter::TreeS() const 
845 {
846   // return pointer to SDigits Tree
847   TTree * rv = 0 ; 
848   rv = EmcalLoader()->TreeS() ; 
849   if ( !rv ) {
850     EmcalLoader()->MakeTree("S");
851     rv = EmcalLoader()->TreeS() ;
852   } 
853   
854   return rv ; 
855 }
856
857 //____________________________________________________________________________ 
858 Bool_t AliEMCALGetter::VersionExists(TString & opt) const
859 {
860   // checks if the version with the present name already exists in the same directory
861
862   Bool_t rv = kFALSE ;
863  
864   AliRunLoader * rl = AliRunLoader::GetRunLoader(EmcalLoader()->GetTitle());
865   TString version( rl->GetEventFolder()->GetName() ) ; 
866
867   opt.ToLower() ; 
868   
869   if ( opt == "sdigits") {
870     // add the version name to the root file name
871     TString fileName( EmcalLoader()->GetSDigitsFileName() ) ; 
872     if (version != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
873       fileName = fileName.ReplaceAll(".root", "") + "_" + version + ".root" ;
874     if ( !(gSystem->AccessPathName(fileName)) ) { 
875       Warning("VersionExists", "The file %s already exists", fileName.Data()) ;
876       rv = kTRUE ; 
877     }
878     EmcalLoader()->SetSDigitsFileName(fileName) ;
879   }
880
881   if ( opt == "digits") {
882     // add the version name to the root file name
883     TString fileName( EmcalLoader()->GetDigitsFileName() ) ; 
884     if (version != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
885       fileName = fileName.ReplaceAll(".root", "") + "_" + version + ".root" ;
886     if ( !(gSystem->AccessPathName(fileName)) ) {
887       Warning("VersionExists", "The file %s already exists", fileName.Data()) ;  
888       rv = kTRUE ; 
889     }
890   }
891
892   return rv ;
893
894 }
895
896 //____________________________________________________________________________ 
897 UShort_t AliEMCALGetter::EventPattern(void) const
898 {
899   // Return the pattern (trigger bit register) of the beam-test event
900 //   if(fBTE)
901 //     return fBTE->GetPattern() ;
902 //   else
903     return 0 ;
904 }
905 //____________________________________________________________________________ 
906 Float_t AliEMCALGetter::BeamEnergy(void) const
907 {
908   // Return the beam energy of the beam-test event
909 //   if(fBTE)
910 //     return fBTE->GetBeamEnergy() ;
911 //   else
912     return 0 ;
913 }