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