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