]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSGetter.cxx
Moved MakePrimary to EMCRecPoint to rpduce correct order of primaries
[u/mrichter/AliRoot.git] / PHOS / AliPHOSGetter.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 AliPHOSAnalyser):
28 //  AliPHOSGetter * gime = AliPHOSGetter::GetInstance("galice.root","test") ;
29 //  for(Int_t irecp = 0; irecp < gime->NRecParticles() ; irecp++)
30 //     AliPHOSRecParticle * 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 AliPHOSIndexToObject 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 <TParticle.h>
47 #include <TF1.h>
48 #include <TGraph.h>
49 #include <TCanvas.h>
50 #include <TStyle.h>
51 //#include <TFrame.h>
52
53 // --- Standard library ---
54
55 // --- AliRoot header files ---
56 #include "AliESD.h"
57 #include "AliHeader.h"  
58 #include "AliMC.h"
59 #include "AliPHOS.h"
60 #include "AliPHOSBeamTestEvent.h"
61 #include "AliPHOSGetter.h"
62 #include "AliPHOSLoader.h"
63 #include "AliPHOSPulseGenerator.h"
64 #include "AliRunLoader.h"
65 #include "AliStack.h"  
66 #include "AliPHOSRawDecoder.h"
67 #include "AliRawReaderFile.h"
68 #include "AliLog.h"
69 #include "AliCDBLocal.h"
70 #include "AliCDBStorage.h"
71 #include "AliCDBManager.h"
72 #include "AliPHOSRawDigiProducer.h"
73 #include "AliPHOSReconstructor.h"
74 #include "AliPHOSRecoParam.h"
75
76 ClassImp(AliPHOSGetter)
77   
78 AliPHOSGetter   * AliPHOSGetter::fgObjGetter  = 0 ; 
79 AliPHOSLoader   * AliPHOSGetter::fgPhosLoader = 0;
80 AliPHOSCalibData* AliPHOSGetter::fgCalibData  = 0;
81 Int_t AliPHOSGetter::fgDebug = 0;
82
83 //  TFile * AliPHOSGetter::fgFile = 0 ; 
84
85 AliPHOSGetter::AliPHOSGetter() :
86   fBTE(0),
87   fLoadingStatus(),
88   fNPrimaries(0),
89   fPrimaries(0),
90   fESDFile(0),
91   fESDFileName(),
92   fESD(0),
93   fESDTree(0),
94   fRawDigits(kFALSE),
95   fcdb(0)
96 {
97   // ctor: this is a singleton, the ctor should never be called but cint needs it as public
98   Fatal("ctor", "AliPHOSGetter is a singleton default ctor not callable") ;
99
100
101
102 //____________________________________________________________________________ 
103 AliPHOSGetter::AliPHOSGetter(const char* headerFile, const char* version, Option_t * openingOption) :
104   fBTE(0),
105   fLoadingStatus(),
106   fNPrimaries(0),
107   fPrimaries(0),
108   fESDFile(0),
109   fESDFileName(),
110   fESD(0),
111   fESDTree(0),
112   fRawDigits(kFALSE),
113   fcdb(0)
114 {
115   // ctor only called by Instance()
116
117   AliRunLoader* rl = AliRunLoader::GetRunLoader(version) ; 
118   if (!rl) {
119     rl = AliRunLoader::Open(headerFile, version, openingOption);
120     if (!rl) {
121       Fatal("AliPHOSGetter", "Could not find the Run Loader for %s - %s",headerFile, version) ; 
122       return ;
123     } 
124     if (rl->GetAliRun() == 0x0) {
125       rl->LoadgAlice();
126       gAlice = rl->GetAliRun(); // should be removed
127       rl->LoadHeader();
128     }
129   }
130   fgPhosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
131   if ( !fgPhosLoader ) 
132     Error("AliPHOSGetter", "Could not find PHOSLoader") ; 
133   else 
134     fgPhosLoader->SetTitle(version);
135   
136   // initialize data members
137   SetDebug(0) ; 
138   fBTE = 0 ; 
139   fPrimaries = 0 ; 
140   fLoadingStatus = "" ; 
141  
142   fESDFileName = rl->GetFileName()  ; // this should be the galice.root file
143   fESDFileName.ReplaceAll("galice.root", "AliESDs.root") ;  
144   fESDFile = 0 ; 
145   fESD = 0 ; 
146   fESDTree = 0 ; 
147   fRawDigits = kFALSE ;
148
149 }
150
151 AliPHOSGetter::AliPHOSGetter(const AliPHOSGetter & obj) :
152   TObject(obj),
153   fBTE(0),
154   fLoadingStatus(),
155   fNPrimaries(0),
156   fPrimaries(0),
157   fESDFile(0),
158   fESDFileName(),
159   fESD(0),
160   fESDTree(0),
161   fRawDigits(kFALSE),
162   fcdb(0)
163 {
164   // cpy ctor requested by Coding Convention 
165   Fatal("cpy ctor", "not implemented") ;
166
167
168 //____________________________________________________________________________ 
169 AliPHOSGetter::AliPHOSGetter(Int_t /*i*/) :
170   fBTE(0),
171   fLoadingStatus(),
172   fNPrimaries(0),
173   fPrimaries(0),
174   fESDFile(0),
175   fESDFileName(),
176   fESD(0),
177   fESDTree(0),
178   fRawDigits(kFALSE),
179   fcdb(0)
180 {
181   // special constructor for onflight 
182
183
184
185 //____________________________________________________________________________ 
186 AliPHOSGetter::~AliPHOSGetter()
187 {
188   // dtor
189   if(fgPhosLoader){
190     delete fgPhosLoader ;
191     fgPhosLoader = 0 ;
192   }
193   if(fBTE){
194     delete fBTE ; 
195     fBTE = 0 ;
196   } 
197   if(fPrimaries){
198     fPrimaries->Delete() ; 
199     delete fPrimaries ;
200   } 
201   if (fESD) 
202     delete fESD ; 
203   if (fESDTree) 
204     delete fESDTree ;
205  
206   fgObjGetter = 0;
207 }
208
209 //____________________________________________________________________________ 
210 void AliPHOSGetter::Reset()
211 {
212   // resets things in case the getter is called consecutively with different files
213   // the PHOS Loader is already deleted by the Run Loader
214
215   if (fPrimaries) { 
216     fPrimaries->Delete() ; 
217     delete fPrimaries ;
218   } 
219   fgPhosLoader = 0; 
220   fgObjGetter = 0; 
221 }
222
223 //____________________________________________________________________________ 
224 AliPHOSClusterizer * AliPHOSGetter::Clusterizer()
225
226   // Returns pointer to the Clusterizer task 
227   AliPHOSClusterizer * rv ; 
228   rv =  dynamic_cast<AliPHOSClusterizer *>(PhosLoader()->Reconstructioner()) ;
229   if (!rv) {
230     Event(0, "R") ; 
231     rv =  dynamic_cast<AliPHOSClusterizer*>(PhosLoader()->Reconstructioner()) ;
232   }
233   return rv ; 
234 }
235
236 //____________________________________________________________________________ 
237 TObjArray * AliPHOSGetter::CpvRecPoints() const
238 {
239   // asks the Loader to return the CPV RecPoints container 
240
241   TObjArray * rv = 0 ; 
242   
243   rv = PhosLoader()->CpvRecPoints() ; 
244   if (!rv) {
245     PhosLoader()->MakeRecPointsArray() ;
246     rv = PhosLoader()->CpvRecPoints() ; 
247   }
248   return rv ; 
249 }
250
251 //____________________________________________________________________________ 
252 TClonesArray * AliPHOSGetter::Digits() const
253 {
254   // asks the Loader to return the Digits container 
255
256   TClonesArray * rv = 0 ; 
257   rv = PhosLoader()->Digits() ; 
258
259   if( !rv ) {
260     PhosLoader()->MakeDigitsArray() ; 
261     rv = PhosLoader()->Digits() ;
262   }
263   return rv ; 
264 }
265
266 //____________________________________________________________________________ 
267 AliPHOSDigitizer * AliPHOSGetter::Digitizer() 
268
269   // Returns pointer to the Digitizer task 
270   AliPHOSDigitizer * rv ; 
271   rv =  dynamic_cast<AliPHOSDigitizer *>(PhosLoader()->Digitizer()) ;
272   if (!rv) {
273     Event(0, "D") ; 
274     rv =  dynamic_cast<AliPHOSDigitizer *>(PhosLoader()->Digitizer()) ;
275   }
276   return rv ; 
277 }
278
279
280 //____________________________________________________________________________ 
281 TObjArray * AliPHOSGetter::EmcRecPoints() const
282 {
283   // asks the Loader to return the EMC RecPoints container 
284
285   TObjArray * rv = 0 ; 
286   
287   rv = PhosLoader()->EmcRecPoints() ; 
288   if (!rv) {
289     PhosLoader()->MakeRecPointsArray() ;
290     rv = PhosLoader()->EmcRecPoints() ; 
291   }
292   return rv ; 
293 }
294
295 //____________________________________________________________________________ 
296 TClonesArray * AliPHOSGetter::TrackSegments() const
297 {
298   // asks the Loader to return the TrackSegments container 
299
300   TClonesArray * rv = 0 ; 
301   
302   rv = PhosLoader()->TrackSegments() ; 
303   if (!rv) {
304     PhosLoader()->MakeTrackSegmentsArray() ;
305     rv = PhosLoader()->TrackSegments() ; 
306   }
307   return rv ; 
308 }
309
310 //____________________________________________________________________________ 
311 AliPHOSTrackSegmentMaker * AliPHOSGetter::TrackSegmentMaker()
312
313   // Returns pointer to the TrackSegmentMaker task 
314   AliPHOSTrackSegmentMaker * rv ; 
315   rv =  dynamic_cast<AliPHOSTrackSegmentMaker *>(PhosLoader()->TrackSegmentMaker()) ;
316   if (!rv) {
317     Event(0, "T") ; 
318     rv =  dynamic_cast<AliPHOSTrackSegmentMaker *>(PhosLoader()->TrackSegmentMaker()) ;
319   }
320   return rv ; 
321 }
322
323 //____________________________________________________________________________ 
324 TClonesArray * AliPHOSGetter::RecParticles() const
325 {
326   // asks the Loader to return the TrackSegments container 
327
328   TClonesArray * rv = 0 ; 
329   
330   rv = PhosLoader()->RecParticles() ; 
331   if (!rv) {
332     PhosLoader()->MakeRecParticlesArray() ;
333     rv = PhosLoader()->RecParticles() ; 
334   }
335   return rv ; 
336 }
337 //____________________________________________________________________________ 
338 void AliPHOSGetter::Event(Int_t event, const char* opt) 
339 {
340   // Reads the content of all Tree's S, D and R
341
342 //   if ( event >= MaxEvent() ) {
343 //     Error("Event", "%d not found in TreeE !", event) ; 
344 //     return ; 
345 //   }
346
347   AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
348
349 //   // checks if we are dealing with test-beam data
350 //   TBranch * btb = rl->TreeE()->GetBranch("AliPHOSBeamTestEvent") ;
351 //   if(btb){
352 //     if(!fBTE)
353 //       fBTE = new AliPHOSBeamTestEvent() ;
354 //     btb->SetAddress(&fBTE) ;
355 //     btb->GetEntry(event) ;
356 //   }
357 //   else{
358 //     if(fBTE){
359 //       delete fBTE ;
360 //       fBTE = 0 ;
361 //     }
362 //   }
363
364   // Loads the type of object(s) requested
365   
366   rl->GetEvent(event) ;
367
368   if(strstr(opt,"X") || (strcmp(opt,"")==0)){
369     ReadPrimaries() ;
370   }
371   
372   if(strstr(opt,"H")  ){
373     ReadTreeH();
374   }
375   
376   if(strstr(opt,"S")  ){
377     ReadTreeS() ;
378   }
379   
380   if(strstr(opt,"D") ){
381     ReadTreeD() ;
382   }
383   
384   if(strstr(opt,"R") ){
385     ReadTreeR() ;
386   }
387
388   if( strstr(opt,"T") ){
389     ReadTreeT() ;
390   }    
391
392   if( strstr(opt,"P") ){
393     ReadTreeP() ;
394   }    
395
396   if( strstr(opt,"E") ){
397     ReadTreeE(event) ;
398   }
399
400 }
401
402
403 //____________________________________________________________________________ 
404 void AliPHOSGetter::Event(AliRawReader *rawReader, const char* opt, Bool_t isOldRCUFormat) 
405 {
406   // Reads the raw event from rawReader
407   // isOldRCUFormat defines whenever to assume
408   // the old RCU format or not
409   
410   if( strstr(opt,"W")  ){
411     rawReader->NextEvent();
412     ReadRaw(rawReader,isOldRCUFormat) ;
413   }    
414  
415 }
416
417
418 //____________________________________________________________________________ 
419 Int_t AliPHOSGetter::EventNumber() const
420   {
421   // return the current event number
422   AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
423   return static_cast<Int_t>(rl->GetEventNumber()) ;   
424 }
425
426 //____________________________________________________________________________ 
427   TClonesArray * AliPHOSGetter::Hits() const
428 {
429   // asks the loader to return  the Hits container 
430   
431   TClonesArray * rv = 0 ; 
432   
433   rv = PhosLoader()->Hits() ; 
434   if ( !rv ) {
435     PhosLoader()->LoadHits("read"); 
436     rv = PhosLoader()->Hits() ; 
437   }
438   return rv ; 
439 }
440
441 //____________________________________________________________________________ 
442 AliPHOSGetter * AliPHOSGetter::Instance(const char* alirunFileName, const char* version, Option_t * openingOption) 
443 {
444   // Creates and returns the pointer of the unique instance
445   // Must be called only when the environment has changed
446   
447   if(!fgObjGetter){ // first time the getter is called 
448     fgObjGetter = new AliPHOSGetter(alirunFileName, version, openingOption) ;
449   }
450   else { // the getter has been called previously
451     AliRunLoader * rl = AliRunLoader::GetRunLoader(fgPhosLoader->GetTitle());
452     if ( rl->GetFileName() == alirunFileName ) {// the alirunFile has the same name
453       // check if the file is already open
454       TFile * galiceFile = dynamic_cast<TFile *>(gROOT->FindObject(rl->GetFileName()) ) ; 
455       
456       if ( !galiceFile ) 
457         fgObjGetter = new AliPHOSGetter(alirunFileName, version, openingOption) ;
458       
459       else {  // the file is already open check the version name
460         TString currentVersionName = rl->GetEventFolder()->GetName() ; 
461         TString newVersionName(version) ; 
462         if (currentVersionName == newVersionName) 
463           if(fgDebug)
464             ::Warning( "Instance", "Files with version %s already open", currentVersionName.Data() ) ;  
465         else {
466           fgObjGetter = new AliPHOSGetter(alirunFileName, version, openingOption) ;      
467         }
468       }
469     }
470     else {
471       rl = AliRunLoader::GetRunLoader(fgPhosLoader->GetTitle()) ; 
472       if ( strstr(version, AliConfig::GetDefaultEventFolderName()) ) // false in case of merging
473         delete rl ; 
474       fgObjGetter = new AliPHOSGetter(alirunFileName, version, openingOption) ;      
475     }
476   }
477   if (!fgObjGetter) 
478     ::Error("AliPHOSGetter::Instance", "Failed to create the PHOS Getter object") ;
479   else 
480     if (fgDebug)
481       Print() ;
482   
483   return fgObjGetter ;
484 }
485
486 //____________________________________________________________________________ 
487 AliPHOSGetter *  AliPHOSGetter::Instance()
488 {
489   // Returns the pointer of the unique instance already defined
490   
491   if(!fgObjGetter && fgDebug)
492      ::Warning("AliPHOSGetter::Instance", "Getter not initialized") ;
493
494    return fgObjGetter ;
495            
496 }
497
498 //____________________________________________________________________________ 
499 Int_t AliPHOSGetter::MaxEvent() const 
500 {
501   // returns the number of events in the run (from TE)
502
503   AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
504   return static_cast<Int_t>(rl->GetNumberOfEvents()) ; 
505 }
506
507 //____________________________________________________________________________ 
508 TParticle * AliPHOSGetter::Primary(Int_t index) const
509 {
510   AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
511   return rl->Stack()->Particle(index) ; 
512
513
514 //____________________________________________________________________________ 
515 AliPHOS * AliPHOSGetter:: PHOS() const  
516 {
517   // returns the PHOS object 
518   AliPHOSLoader *    loader = 0;
519   static AliPHOSLoader * oldloader = 0;
520   static AliPHOS * phos = 0;
521
522   loader = PhosLoader();
523
524   if ( loader != oldloader) {
525     phos = dynamic_cast<AliPHOS*>(loader->GetModulesFolder()->FindObject("PHOS")) ;
526     oldloader = loader;
527   }
528   if (!phos) 
529     if (fgDebug)
530       Warning("PHOS", "PHOS module not found in module folders: %s", PhosLoader()->GetModulesFolder()->GetName() ) ; 
531   return phos ; 
532 }  
533
534
535
536 //____________________________________________________________________________ 
537 AliPHOSPID * AliPHOSGetter::PID()
538
539   // Returns pointer to the PID task 
540   AliPHOSPID * rv ; 
541   rv =  dynamic_cast<AliPHOSPID *>(PhosLoader()->PIDTask()) ;
542   if (!rv) {
543     Event(0, "P") ; 
544     rv =  dynamic_cast<AliPHOSPID *>(PhosLoader()->PIDTask()) ;
545   }
546   return rv ; 
547 }
548
549 //____________________________________________________________________________ 
550 AliPHOSGeometry * AliPHOSGetter::PHOSGeometry() const 
551 {
552   // Returns PHOS geometry
553
554   AliPHOSGeometry * rv = 0 ; 
555   if (PHOS() )
556     rv =  PHOS()->GetGeometry() ;
557   else {
558     rv = AliPHOSGeometry::GetInstance();
559     if (!rv) {
560       AliError("Could not find PHOS geometry! Loading the default one !");
561       rv = AliPHOSGeometry::GetInstance("IHEP","");
562     }
563   }
564   return rv ; 
565
566
567 //____________________________________________________________________________ 
568 TClonesArray * AliPHOSGetter::Primaries()  
569 {
570   // creates the Primaries container if needed
571   if ( !fPrimaries ) {
572     if (fgDebug) 
573       Info("Primaries", "Creating a new TClonesArray for primaries") ; 
574     fPrimaries = new TClonesArray("TParticle", 1000) ;
575   } 
576   return fPrimaries ; 
577 }
578
579 //____________________________________________________________________________ 
580 void  AliPHOSGetter::Print() 
581 {
582   // Print usefull information about the getter
583     
584   AliRunLoader * rl = AliRunLoader::GetRunLoader(fgPhosLoader->GetTitle());
585   ::Info( "Print", "gAlice file is %s -- version name is %s", (rl->GetFileName()).Data(), rl->GetEventFolder()->GetName() ) ; 
586 }
587
588 //____________________________________________________________________________ 
589 void AliPHOSGetter::ReadPrimaries()  
590 {
591   // Read Primaries from Kinematics.root
592   
593   AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
594   
595   // gets kine tree from the root file (Kinematics.root)
596   if ( ! rl->TreeK() ) { // load treeK the first time
597     rl->LoadKinematics() ;
598   }
599   
600   fNPrimaries = (rl->GetHeader())->GetNtrack(); 
601   if (fgDebug) 
602     Info( "ReadTreeK", "Found %d particles in event # %d", fNPrimaries, EventNumber() ) ; 
603
604
605   // first time creates the container
606   if ( Primaries() ) 
607     fPrimaries->Clear() ; 
608   
609   Int_t index = 0 ; 
610   for (index = 0 ; index < fNPrimaries; index++) { 
611     new ((*fPrimaries)[index]) TParticle(*(Primary(index)));
612   }
613 }
614
615 //____________________________________________________________________________ 
616 Bool_t AliPHOSGetter::OpenESDFile() 
617 {
618   //Open the ESD file    
619   Bool_t rv = kTRUE ; 
620   if (!fESDFile) {
621     fESDFile = TFile::Open(fESDFileName) ;
622     if (!fESDFile ) 
623       return kFALSE ; 
624   }
625   else if (fESDFile->IsOpen()) {
626     fESDFile->Close() ; 
627     fESDFile = TFile::Open(fESDFileName) ;
628   }
629   if (!fESDFile->IsOpen())
630     rv = kFALSE ; 
631   return rv ; 
632 }
633
634 //____________________________________________________________________________ 
635 void AliPHOSGetter::FitRaw(Bool_t lowGainFlag, TGraph * gLowGain, TGraph * gHighGain, TF1* signalF, Double_t & energy, Double_t & time) const
636 {
637   // Fits the raw signal time distribution 
638
639   const Int_t kNoiseThreshold = 0 ;
640   Double_t timezero1 = 0., timezero2 = 0., timemax = 0. ;
641   Double_t signal = 0., signalmax = 0. ;       
642   time   = 0. ; 
643   energy = 0. ; 
644
645   // Create a shaper pulse object which contains all the shaper parameters
646   AliPHOSPulseGenerator pulse;
647
648   if (lowGainFlag) {
649     timezero1 = timezero2 = signalmax = timemax = 0. ;
650     signalF->FixParameter(0, pulse.GetRawFormatLowCharge()) ; 
651     signalF->FixParameter(1, pulse.GetRawFormatLowGain()) ; 
652     Int_t index ; 
653     for (index = 0; index < pulse.GetRawFormatTimeBins(); index++) {
654       gLowGain->GetPoint(index, time, signal) ; 
655       if (signal > kNoiseThreshold && timezero1 == 0.) 
656         timezero1 = time ;
657       if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
658         timezero2 = time ; 
659       if (signal > signalmax) {
660         signalmax = signal ; 
661         timemax   = time ; 
662       }
663     }
664     signalmax /= 
665       pulse.RawResponseFunctionMax(pulse.GetRawFormatLowCharge(), 
666                                    pulse.GetRawFormatLowGain()) ;
667     if ( timezero1 + pulse.GetRawFormatTimePeak() < pulse.GetRawFormatTimeMax() * 0.4 ) { // else its noise 
668       signalF->SetParameter(2, signalmax) ; 
669       signalF->SetParameter(3, timezero1) ;         
670       gLowGain->Fit(signalF, "QRO", "", 0., timezero2); //, "QRON") ; 
671       energy = signalF->GetParameter(2) ; 
672       time   = signalF->GetMaximumX() - pulse.GetRawFormatTimePeak() - pulse.GetRawFormatTimeTrigger() ;
673     }
674   } else {
675     timezero1 = timezero2 = signalmax = timemax = 0. ;
676     signalF->FixParameter(0, pulse.GetRawFormatHighCharge()) ; 
677     signalF->FixParameter(1, pulse.GetRawFormatHighGain()) ; 
678     Int_t index ; 
679     for (index = 0; index < pulse.GetRawFormatTimeBins(); index++) {
680       gHighGain->GetPoint(index, time, signal) ;               
681       if (signal > kNoiseThreshold && timezero1 == 0.) 
682         timezero1 = time ;
683       if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
684         timezero2 = time ; 
685       if (signal > signalmax) {
686         signalmax = signal ;   
687         timemax   = time ; 
688       }
689     }
690     signalmax /= pulse.RawResponseFunctionMax(pulse.GetRawFormatHighCharge(), 
691                                               pulse.GetRawFormatHighGain()) ;;
692     if ( timezero1 + pulse.GetRawFormatTimePeak() < pulse.GetRawFormatTimeMax() * 0.4 ) { // else its noise  
693       signalF->SetParameter(2, signalmax) ; 
694       signalF->SetParameter(3, timezero1) ;               
695       gHighGain->Fit(signalF, "QRO", "", 0., timezero2) ; 
696       energy = signalF->GetParameter(2) ; 
697       time   = signalF->GetMaximumX() - pulse.GetRawFormatTimePeak() - pulse.GetRawFormatTimeTrigger() ;
698     }
699   }
700   if (time == 0) energy = 0 ; 
701 }
702
703 //____________________________________________________________________________ 
704 Int_t AliPHOSGetter::CalibrateRaw(Double_t energy, Int_t *relId)
705 {
706   // Convert energy into digitized amplitude for a cell relId
707   // It is a user responsilibity to open CDB and set
708   // AliPHOSCalibData object by the following operators:
709   // 
710   // AliCDBLocal *loc = new AliCDBLocal("deCalibDB");
711   // AliPHOSCalibData* clb = (AliPHOSCalibData*)AliCDBStorage::Instance()
712   //    ->Get(path_to_calibdata,run_number);
713   // AliPHOSGetter* gime = AliPHOSGetter::Instance("galice.root");
714   // gime->SetCalibData(clb);
715
716   if (CalibData() == 0)
717     Warning("CalibrateRaw","Calibration DB is not initiated!");
718
719   Int_t   module = relId[0];
720   Int_t   column = relId[3];
721   Int_t   row    = relId[2];
722
723   Float_t gainFactor = 0.0015; // width of one Emc ADC channel in GeV
724   Float_t pedestal   = 0.005;  // Emc pedestals
725
726   if(CalibData()) {
727     gainFactor = CalibData()->GetADCchannelEmc (module,column,row);
728     pedestal   = CalibData()->GetADCpedestalEmc(module,column,row);
729   }
730   
731   Int_t   amp = static_cast<Int_t>( (energy - pedestal) / gainFactor + 0.5 ) ; 
732   return amp;
733 }
734 //____________________________________________________________________________ 
735 Int_t AliPHOSGetter::ReadRaw(AliRawReader *rawReader,Bool_t isOldRCUFormat)
736 {
737   // reads the raw format data, converts it into digits format and store digits in Digits()
738   // container.
739   // isOldRCUFormat = kTRUE in case of the old RCU
740   // format used in the raw data readout.
741   // Reimplemented by Boris Polichtchouk (Jul 2006)
742   // to make it working with the Jul-Aug 2006 beam test data.
743  
744   //Create raw decoder.
745
746   AliPHOSRawDecoder dc(rawReader);
747   dc.SetOldRCUFormat(isOldRCUFormat);
748   dc.SubtractPedestals(AliPHOSReconstructor::GetRecoParamEmc()->SubtractPedestals());
749
750   TClonesArray * digits = Digits() ;
751   AliPHOSRawDigiProducer pr;
752   pr.MakeDigits(digits,&dc);
753   
754   //ADC counts -> GeV
755   for(Int_t i=0; i<digits->GetEntries(); i++) {
756     AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(i);
757     digit->SetEnergy(digit->GetEnergy()/AliPHOSPulseGenerator::GeV2ADC());
758   }
759   
760   //!!!!for debug!!!
761   Int_t modMax=-111;
762   Int_t colMax=-111;
763   Int_t rowMax=-111;
764   Float_t eMax=-333;
765   //!!!for debug!!!
766
767   Int_t relId[4];
768   for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
769     AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
770     if(digit->GetEnergy()>eMax) {
771       PHOSGeometry()->AbsToRelNumbering(digit->GetId(),relId);
772       eMax=digit->GetEnergy();
773       modMax=relId[0];
774       rowMax=relId[2];
775       colMax=relId[3];
776     }
777   }
778
779   AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n\n",
780                   modMax,colMax,rowMax,eMax));
781
782   return digits->GetEntriesFast() ; 
783 }
784
785 //____________________________________________________________________________ 
786 Int_t AliPHOSGetter::ReadTreeD() const
787 {
788   // Read the Digits
789   
790   PhosLoader()->CleanDigits() ;    
791   PhosLoader()->LoadDigits("UPDATE") ;
792   PhosLoader()->LoadDigitizer("UPDATE") ;
793   return Digits()->GetEntries() ; 
794 }
795
796 //____________________________________________________________________________ 
797 Int_t AliPHOSGetter::ReadTreeH() const
798 {
799   // Read the Hits
800   PhosLoader()->CleanHits() ;
801   // gets TreeH from the root file (PHOS.Hit.root)
802   //if ( !IsLoaded("H") ) {
803     PhosLoader()->LoadHits("UPDATE") ;
804   //  SetLoaded("H") ; 
805   //}  
806   return Hits()->GetEntries() ; 
807 }
808
809 //____________________________________________________________________________ 
810 Int_t AliPHOSGetter::ReadTreeR() const
811 {
812   // Read the RecPoints
813   
814   PhosLoader()->CleanRecPoints() ;
815   // gets TreeR from the root file (PHOS.RecPoints.root)
816   //if ( !IsLoaded("R") ) {
817     PhosLoader()->LoadRecPoints("UPDATE") ;
818     //  SetLoaded("R") ; 
819     //}
820
821   return EmcRecPoints()->GetEntries() ; 
822 }
823
824 //____________________________________________________________________________ 
825 Int_t AliPHOSGetter::ReadTreeT() const
826 {
827   // Read the TrackSegments
828   
829   PhosLoader()->CleanTracks() ; 
830   // gets TreeT from the root file (PHOS.TrackSegments.root)
831   //if ( !IsLoaded("T") ) {
832     PhosLoader()->LoadTracks("UPDATE") ;
833     PhosLoader()->LoadTrackSegmentMaker("UPDATE") ;
834     //    SetLoaded("T") ; 
835     //}
836
837   return TrackSegments()->GetEntries() ; 
838 }
839 //____________________________________________________________________________ 
840 Int_t AliPHOSGetter::ReadTreeP() const
841 {
842   // Read the RecParticles
843   
844   PhosLoader()->CleanRecParticles() ; 
845
846   // gets TreeT from the root file (PHOS.TrackSegments.root)
847   //  if ( !IsLoaded("P") ) {
848     PhosLoader()->LoadRecParticles("UPDATE") ;
849     PhosLoader()->LoadPID("UPDATE") ;
850     //  SetLoaded("P") ; 
851     //}
852
853   return RecParticles()->GetEntries() ; 
854 }
855 //____________________________________________________________________________ 
856 Int_t AliPHOSGetter::ReadTreeS() const
857 {
858   // Read the SDigits
859   
860   PhosLoader()->CleanSDigits() ; 
861   // gets TreeS from the root file (PHOS.SDigits.root)
862   //if ( !IsLoaded("S") ) {
863     PhosLoader()->LoadSDigits("READ") ;
864     PhosLoader()->LoadSDigitizer("READ") ;
865     //  SetLoaded("S") ; 
866     //}
867
868   return SDigits()->GetEntries() ; 
869 }
870
871 //____________________________________________________________________________ 
872 Int_t AliPHOSGetter::ReadTreeE(Int_t event)
873 {
874   // Read the ESD
875   
876   // gets esdTree from the root file (AliESDs.root)
877   if (!fESDFile)
878     if ( !OpenESDFile() ) 
879       return -1 ; 
880
881   fESDTree = static_cast<TTree*>(fESDFile->Get("esdTree")) ; 
882   fESD = new AliESD;
883    if (!fESDTree) {
884
885      Error("ReadTreeE", "no ESD tree found");
886      return -1;
887    }
888    fESDTree->SetBranchAddress("ESD", &fESD);
889    fESDTree->GetEvent(event);
890
891    return event ; 
892 }
893
894 //____________________________________________________________________________ 
895 TClonesArray * AliPHOSGetter::SDigits() const
896 {
897   // asks the Loader to return the Digits container 
898
899   TClonesArray * rv = 0 ; 
900   
901   rv = PhosLoader()->SDigits() ; 
902   if (!rv) {
903     PhosLoader()->MakeSDigitsArray() ;
904     rv = PhosLoader()->SDigits() ; 
905   }
906   return rv ; 
907 }
908
909 //____________________________________________________________________________ 
910 AliPHOSSDigitizer * AliPHOSGetter::SDigitizer()
911
912   // Returns pointer to the SDigitizer task 
913   AliPHOSSDigitizer * rv ; 
914   rv =  dynamic_cast<AliPHOSSDigitizer *>(PhosLoader()->SDigitizer()) ;
915   if (!rv) {
916     Event(0, "S") ; 
917     rv =  dynamic_cast<AliPHOSSDigitizer *>(PhosLoader()->SDigitizer()) ;
918   }
919   return rv ; 
920 }
921
922 //____________________________________________________________________________ 
923 TParticle * AliPHOSGetter::Secondary(const TParticle* p, Int_t index) const
924 {
925   // Return first (index=1) or second (index=2) secondary particle of primary particle p 
926
927   if(index <= 0) 
928     return 0 ;
929   if(index > 2)
930     return 0 ;
931
932   if(p) {
933   Int_t daughterIndex = p->GetDaughter(index-1) ; 
934   AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
935   return  rl->GetAliRun()->GetMCApp()->Particle(daughterIndex) ; 
936   }
937   else
938     return 0 ;
939 }
940
941 //____________________________________________________________________________ 
942 void AliPHOSGetter::Track(Int_t itrack) 
943 {
944   // Read the first entry of PHOS branch in hit tree gAlice->TreeH()
945  
946  AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
947
948   if( !TreeH() ) // load treeH the first time
949     rl->LoadHits() ;
950
951   // first time create the container
952   TClonesArray * hits = Hits() ; 
953   if ( hits ) 
954     hits->Clear() ; 
955
956   TBranch * phosbranch = dynamic_cast<TBranch*>(TreeH()->GetBranch("PHOS")) ; 
957   phosbranch->SetAddress(&hits) ;
958   phosbranch->GetEntry(itrack) ;
959 }
960
961 //____________________________________________________________________________ 
962 TTree * AliPHOSGetter::TreeD() const 
963 {
964   // Returns pointer to the Digits Tree
965   TTree * rv = 0 ; 
966   rv = PhosLoader()->TreeD() ; 
967   if ( !rv ) {
968     PhosLoader()->MakeTree("D");
969     rv = PhosLoader()->TreeD() ;
970   } 
971   
972   return rv ; 
973 }
974
975 //____________________________________________________________________________ 
976 TTree * AliPHOSGetter::TreeH() const 
977 {
978   // Returns pointer to the Hits Tree
979   TTree * rv = 0 ; 
980   rv = PhosLoader()->TreeH() ; 
981   if ( !rv ) {
982     PhosLoader()->MakeTree("H");
983     rv = PhosLoader()->TreeH() ;
984   } 
985   
986   return rv ; 
987 }
988
989 //____________________________________________________________________________ 
990 TTree * AliPHOSGetter::TreeR() const 
991 {
992   // Returns pointer to the RecPoints Tree
993   TTree * rv = 0 ; 
994   rv = PhosLoader()->TreeR() ; 
995   if ( !rv ) {
996     PhosLoader()->MakeTree("R");
997     rv = PhosLoader()->TreeR() ;
998   } 
999   
1000   return rv ; 
1001 }
1002
1003 //____________________________________________________________________________ 
1004 TTree * AliPHOSGetter::TreeT() const 
1005 {
1006   // Returns pointer to the TrackSegments Tree
1007   TTree * rv = 0 ; 
1008   rv = PhosLoader()->TreeT() ; 
1009   if ( !rv ) {
1010     PhosLoader()->MakeTree("T");
1011     rv = PhosLoader()->TreeT() ;
1012   } 
1013   
1014   return rv ; 
1015 }
1016 //____________________________________________________________________________ 
1017 TTree * AliPHOSGetter::TreeP() const 
1018 {
1019   // Returns pointer to the RecParticles  Tree
1020   TTree * rv = 0 ; 
1021   rv = PhosLoader()->TreeP() ; 
1022   if ( !rv ) {
1023     PhosLoader()->MakeTree("P");
1024     rv = PhosLoader()->TreeP() ;
1025   } 
1026   
1027   return rv ; 
1028 }
1029
1030 //____________________________________________________________________________ 
1031 TTree * AliPHOSGetter::TreeS() const 
1032
1033  // Returns pointer to the SDigits Tree
1034   TTree * rv = 0 ; 
1035   rv = PhosLoader()->TreeS() ; 
1036   if ( !rv ) {
1037     PhosLoader()->MakeTree("S");
1038     rv = PhosLoader()->TreeS() ;
1039   } 
1040   
1041   return rv ; 
1042 }
1043
1044 //____________________________________________________________________________ 
1045 Bool_t AliPHOSGetter::VersionExists(TString & opt) const
1046 {
1047   // checks if the version with the present name already exists in the same directory
1048
1049   Bool_t rv = kFALSE ;
1050  
1051   AliRunLoader * rl = AliRunLoader::GetRunLoader(PhosLoader()->GetTitle());
1052   TString version( rl->GetEventFolder()->GetName() ) ; 
1053
1054   opt.ToLower() ; 
1055   
1056   if ( opt == "sdigits") {
1057     // add the version name to the root file name
1058     TString fileName( PhosLoader()->GetSDigitsFileName() ) ; 
1059     if (version != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
1060       fileName = fileName.ReplaceAll(".root", "") + "_" + version + ".root" ;
1061     if ( !(gSystem->AccessPathName(fileName)) ) { 
1062       Warning("VersionExists", "The file %s already exists", fileName.Data()) ;
1063       rv = kTRUE ; 
1064     }
1065     PhosLoader()->SetSDigitsFileName(fileName) ;
1066   }
1067
1068   if ( opt == "digits") {
1069     // add the version name to the root file name
1070     TString fileName( PhosLoader()->GetDigitsFileName() ) ; 
1071     if (version != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
1072       fileName = fileName.ReplaceAll(".root", "") + "_" + version + ".root" ;
1073     if ( !(gSystem->AccessPathName(fileName)) ) {
1074       Warning("VersionExists", "The file %s already exists", fileName.Data()) ;  
1075       rv = kTRUE ; 
1076     }
1077   }
1078
1079   return rv ;
1080
1081 }
1082
1083 //____________________________________________________________________________ 
1084 UShort_t AliPHOSGetter::EventPattern(void) const
1085 {
1086   // Return the pattern (trigger bit register) of the beam-test event
1087   if(fBTE)
1088     return fBTE->GetPattern() ;
1089   else
1090     return 0 ;
1091 }
1092 //____________________________________________________________________________ 
1093 Float_t AliPHOSGetter::BeamEnergy(void) const
1094 {
1095   // Return the beam energy of the beam-test event
1096   if(fBTE)
1097     return fBTE->GetBeamEnergy() ;
1098   else
1099     return 0 ;
1100 }
1101 //____________________________________________________________________________ 
1102
1103 AliPHOSCalibData* AliPHOSGetter::CalibData()
1104
1105   // Check if the instance of AliPHOSCalibData exists, and return it
1106
1107   return fgCalibData;
1108 }