Changes for Root6 (Mikolaj)
[u/mrichter/AliRoot.git] / AD / ADsim / AliADDigitizer.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: AliADDigitizer.cxx  $ */
17
18 ///_________________________________________________________________________
19 ///
20 /// This class constructs Digits out of Hits
21 ///
22 ///
23
24 // --- Standard library ---
25
26 // --- ROOT system ---
27 #include <TMath.h>
28 #include <TTree.h>
29 #include <TMap.h>
30 #include <TGeoManager.h>
31 #include <TGeoPhysicalNode.h>
32 #include <AliGeomManager.h>
33 #include <TRandom.h>
34 #include <TF1.h>
35 #include <TH1F.h>
36
37 // --- AliRoot header files ---
38 #include "AliRun.h"
39 #include "AliDetector.h"
40 #include "AliAD.h"
41 #include "AliADhit.h"
42 #include "AliADConst.h"
43 #include "AliRunLoader.h"
44 #include "AliLoader.h"
45 #include "AliGRPObject.h"
46 #include "AliDigitizationInput.h"
47 #include "AliCDBManager.h"
48 #include "AliCDBStorage.h"
49 #include "AliCDBEntry.h"
50 #include "AliADCalibData.h"
51 #include "AliCTPTimeParams.h"
52 #include "AliLHCClockPhase.h"
53 #include "AliADdigit.h"
54 #include "AliADDigitizer.h"
55 #include "AliADSDigit.h"
56 #include "AliADTriggerSimulator.h"
57
58 ClassImp(AliADDigitizer)
59
60 //____________________________________________________________________________
61  AliADDigitizer::AliADDigitizer()
62                    :AliDigitizer(),
63                     fCalibData(GetCalibData()),
64                     fNdigits(0),
65                     fDigits(0),
66                     fSignalShape(NULL),
67                     fPMResponse(NULL),
68                     fSinglePhESpectrum(NULL),
69                     fEvenOrOdd(kFALSE),
70                     fTask(kHits2Digits),
71                     fAD(NULL)
72 {
73   // default constructor
74   // Initialize OCDB and containers used in the digitization
75
76   Init();
77 }
78
79 //____________________________________________________________________________ 
80   AliADDigitizer::AliADDigitizer(AliAD *AD, DigiTask_t task)
81                     :AliDigitizer(),
82                      fCalibData(GetCalibData()),
83                      fNdigits(0),
84                      fDigits(0),
85                      fSignalShape(NULL),
86                      fPMResponse(NULL),
87                      fSinglePhESpectrum(NULL),
88                      fEvenOrOdd(kFALSE),
89                      fTask(task),
90                      fAD(AD)
91 {
92   // constructor
93   // Initialize OCDB and containers used in the digitization
94
95   Init();
96 }
97            
98 //____________________________________________________________________________ 
99   AliADDigitizer::AliADDigitizer(AliDigitizationInput* digInput)
100                     :AliDigitizer(digInput),
101                      fCalibData(GetCalibData()),
102                      fNdigits(0),
103                      fDigits(0),
104                      fSignalShape(NULL),
105                      fPMResponse(NULL),
106                      fSinglePhESpectrum(NULL),
107                      fEvenOrOdd(kFALSE),
108                      fTask(kHits2Digits),
109                      fAD(NULL)
110 {
111   // constructor
112   // Initialize OCDB and containers used in the digitization
113
114   Init();
115 }
116            
117 //____________________________________________________________________________ 
118   AliADDigitizer::~AliADDigitizer()
119 {
120   // destructor
121   
122   if (fDigits) {
123     fDigits->Delete();
124     delete fDigits;
125     fDigits=0; 
126   }
127
128   if (fSignalShape) {
129     delete fSignalShape;
130     fSignalShape = NULL;
131   }
132   if (fPMResponse) {
133     delete fPMResponse;
134     fPMResponse = NULL;
135   }
136   if (fSinglePhESpectrum) {
137     delete fSinglePhESpectrum;
138     fSinglePhESpectrum = NULL;
139   }
140
141   for(Int_t i = 0 ; i < 16; ++i) {
142     if (fTime[i]) delete [] fTime[i];
143   }
144 }
145
146 //____________________________________________________________________________ 
147 Bool_t AliADDigitizer::Init()
148 {
149   // Initialises the digitizer
150   // Initialize OCDB and containers used in the digitization
151
152   // check if the digitizer was already initialized
153   if (fSignalShape) return kTRUE;
154   fSignalShape = new TF1("ADSignalShape",this,&AliADDigitizer::SignalShape,0,200,6,"AliADDigitizer","SignalShape");
155   fPMResponse = new TF1("ADPMResponse",this,&AliADDigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliADDigitizer","PMResponse");
156   fSinglePhESpectrum = new TF1("ADSinglePhESpectrum",this,&AliADDigitizer::SinglePhESpectrum,0,20,0,"AliADDigitizer","SinglePhESpectrum");
157   //  fSignalShape->SetParameters(0,1.57345e1,-4.25603e-1,2.9,6.40982,3.69339e-01);
158   //  fSignalShape->SetParameters(1.34330e+00,1.13007e+02,-4.95705e-01,
159   //                          3.68911e+00,1.01040e+00, 3.94675e-01);
160   fSignalShape->SetParameters(-1.07335e+00,2.16002e+01,-1.26133e-01,
161                               1.41619e+00,5.50334e-01,3.86111e-01);
162   
163   // Now get the CTP L0->L1 delay
164   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
165   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
166   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
167   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
168
169   AliCDBEntry *entry1 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
170   if (!entry1) AliFatal("CTP time-alignment is not found in OCDB !");
171   AliCTPTimeParams *ctpTimeAlign = (AliCTPTimeParams*)entry1->GetObject();
172   l1Delay += ((Float_t)ctpTimeAlign->GetDelayL1L0()*25.0);
173
174   AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("AD/Calib/TimeDelays");
175   if (!entry2) AliFatal("AD time delays are not found in OCDB !");
176   TH1F *delays = (TH1F*)entry2->GetObject();
177
178   AliCDBEntry *entry3 = AliCDBManager::Instance()->Get("GRP/Calib/LHCClockPhase");
179   if (!entry3) AliFatal("LHC clock-phase shift is not found in OCDB !");
180   AliLHCClockPhase *phase = (AliLHCClockPhase*)entry3->GetObject();
181
182   for(Int_t i = 0 ; i < 16; ++i) {
183
184     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
185     fLeadingTime[i] = fTimeWidth[i] = 0;
186
187     fPmGain[i] = fCalibData->GetGain(i);
188
189     fAdcPedestal[i][0] = fCalibData->GetPedestal(i);
190     fAdcSigma[i][0]    = fCalibData->GetSigma(i); 
191     fAdcPedestal[i][1] = fCalibData->GetPedestal(i+16);
192     fAdcSigma[i][1]    = fCalibData->GetSigma(i+16); 
193
194     Int_t board = AliADCalibData::GetBoardNumber(i);
195     fNBins[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0+
196                              (Float_t)kMaxTDCWidth*fCalibData->GetWidthResolution(board))/
197                             fCalibData->GetTimeResolution(board));
198     fNBinsLT[i] = TMath::Nint(((Float_t)(fCalibData->GetMatchWindow(board)+1)*25.0)/
199                               fCalibData->GetTimeResolution(board));
200     fBinSize[i] = fCalibData->GetTimeResolution(board);
201     
202     fHptdcOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
203                         (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0
204                        //+fCalibData->GetTimeOffset(i)
205                        +delays->GetBinContent(16-i) - 250);
206                        //-l1Delay-
207                        //-phase->GetMeanPhase()
208                        //+kADOffset);
209                        
210     fClockOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
211                         (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0
212                        +fCalibData->GetTimeOffset(i) - 250);
213                        //-l1Delay
214                        //+kADOffset);
215
216     fTime[i] = new Float_t[fNBins[i]];
217     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
218       
219   }
220   //std::cout<<"AD: "<<" fNBins = "<<fNBins[0]<<" fNBinsLT = "<<fNBinsLT[0]<<" fHptdcOffset = "<<fHptdcOffset[0]<<" fClockOffset = "<<fClockOffset[0]<<std::endl;
221   return kTRUE;
222
223 }
224
225 //____________________________________________________________________________ 
226 void AliADDigitizer::Digitize(Option_t* /*option*/) 
227 {   
228    // Creates digits from hits
229   fNdigits = 0;  
230
231   if (fAD && !fDigInput) {
232     AliLoader *loader = fAD->GetLoader();
233     if (!loader) {
234       AliError("Can not get AD Loader via AliAD object!");
235       return;
236     }
237     AliRunLoader* runLoader = AliRunLoader::Instance();
238     for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); ++iEvent) {
239       runLoader->GetEvent(iEvent);
240       if (fTask == kHits2Digits) {
241         DigitizeHits();
242         DigitizeSDigits();
243         WriteDigits(loader);
244       }
245       else {
246         DigitizeHits();
247         WriteSDigits(loader);
248       }
249     }
250   }
251   else if (fDigInput) {
252       ReadSDigits();
253       DigitizeSDigits();
254       AliRunLoader *currentLoader = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
255       AliLoader *loader = currentLoader->GetLoader("ADLoader");
256       if (!loader) { 
257         AliError("Cannot get AD Loader via RunDigitizer!");
258         return;
259       }
260       WriteDigits(loader);
261   }
262   else {
263     AliFatal("Invalid digitization task! Exiting!");
264   }
265 }
266
267 //____________________________________________________________________________ 
268 void AliADDigitizer::DigitizeHits()
269 {
270   // Digitize the hits to the level of
271   // SDigits (fTime arrays)
272
273   for(Int_t i = 0 ; i < 16; ++i) {
274     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
275     fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
276   }
277   Float_t integral = fPMResponse->Integral(-kPMRespTime,2.*kPMRespTime);
278   Float_t meansPhE = fSinglePhESpectrum->Mean(0,20);
279   
280      AliLoader* loader = fAD->GetLoader();
281      if (!loader) {
282        AliError("Can not get AD Loader!");
283        return;
284      }
285      loader->LoadHits();
286      TTree* treeH = loader->TreeH();
287      if (!treeH) {
288        AliError("Cannot get TreeH!");
289        return;
290      }
291      TClonesArray* hits = fAD->Hits();
292
293 //  Now makes Digits from hits
294      Int_t nTracks = (Int_t) treeH->GetEntries();
295      for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
296          fAD->ResetHits();
297          treeH->GetEvent(iTrack);
298          Int_t nHits = hits->GetEntriesFast();
299          for (Int_t iHit = 0; iHit < nHits; iHit++) {
300            AliADhit* hit = (AliADhit *)hits->UncheckedAt(iHit);
301            Int_t nPhot = hit->GetNphot();
302            Int_t pmt  = hit->GetCell();//One PM per cell in AD                          
303            if (pmt < 0) continue;
304            Int_t trackLabel = hit->GetTrack();
305            for(Int_t l = 0; l < 3; ++l) {
306              if (fLabels[pmt][l] < 0) {
307                fLabels[pmt][l] = trackLabel;
308                break;
309              }
310            }
311            Float_t dt_scintillator = gRandom->Gaus(0,kIntTimeRes);
312            Float_t t = dt_scintillator + hit->GetTof();
313            t += fHptdcOffset[pmt];
314            
315            Float_t charge = nPhot*fPmGain[pmt]*fBinSize[pmt]/integral;
316              
317              Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
318              
319              Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
320              Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
321              //std::cout<<"Bins: "<<firstBin*fBinSize[pmt]<<" - "<<lastBin*fBinSize[pmt]<<std::endl;
322              //std::cout<<"Bins: "<<firstBin<<" - "<<lastBin<<std::endl;
323              
324              for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
325                Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
326                fTime[pmt][iBin] += charge*fPMResponse->Eval(tempT);
327              }
328          }           // hit loop
329      }               // track loop
330      loader->UnloadHits();
331 }
332
333 //____________________________________________________________________________ 
334 void AliADDigitizer::DigitizeSDigits()
335 {
336   // Digitize the fTime arrays (SDigits) to the level of
337   // Digits (fAdc arrays)
338   for(Int_t i = 0 ; i < 16; ++i) {
339     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
340     fLeadingTime[i] = fTimeWidth[i] = 0;
341   }
342
343   Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200); // Not exact, one needs to do this on the convoluted
344   Float_t integral2 = fSignalShape->Integral(0,200); // function. Anyway the effect is small <10% on the 2.5 ADC thr
345   for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
346     Float_t thr = fCalibData->GetCalibDiscriThr(ipmt,kFALSE)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
347        
348     Bool_t ltFound = kFALSE, ttFound = kFALSE;
349     for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
350       Float_t t = fBinSize[ipmt]*Float_t(iBin);
351       if (fTime[ipmt][iBin] > thr) {
352         if (!ltFound && (iBin < fNBinsLT[ipmt])) {
353           ltFound = kTRUE;
354           fLeadingTime[ipmt] = t;
355         }
356       }
357       else {
358         if (ltFound) {
359           if (!ttFound) {
360             ttFound = kTRUE;
361             fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
362           }
363         }
364       }
365       Float_t tadc = t - fClockOffset[ipmt];
366       Int_t clock = kNClocks/2 - Int_t(tadc/25.0);
367       if (clock >= 0 && clock < kNClocks)
368         fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
369     }
370     AliDebug(1,Form("Channel %d Offset %f Time %f",ipmt,fClockOffset[ipmt],fLeadingTime[ipmt]));
371     Int_t board = AliADCalibData::GetBoardNumber(ipmt);
372     if (ltFound && ttFound) {
373       fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
374         Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
375       if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
376         fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
377       if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
378         fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
379     }
380   }
381
382   fEvenOrOdd = gRandom->Integer(2);
383   for (Int_t j=0; j<16; ++j){
384     for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
385       Int_t integrator = (iClock + fEvenOrOdd) % 2;
386       AliDebug(1,Form("ADC %d %d %f",j,iClock,fAdc[j][iClock]));
387       fAdc[j][iClock]  += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
388     }
389   }
390   //Fill BB and BG flags in trigger simulator
391   AliADTriggerSimulator * triggerSimulator = new AliADTriggerSimulator();
392   triggerSimulator->FillFlags(fBBFlag,fBGFlag,fLeadingTime);
393         
394 }
395
396 //____________________________________________________________________________ 
397 void AliADDigitizer::ReadSDigits()
398 {
399   // Read SDigits which are then to precessed
400   // in the following method
401   for(Int_t i = 0 ; i < 16; ++i) {
402     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
403     fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
404   }
405
406   // Loop over input files
407   Int_t nFiles= fDigInput->GetNinputs();
408   for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
409     // Get the current loader 
410     AliRunLoader* currentLoader = 
411       AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
412
413     AliLoader *loader = currentLoader->GetLoader("ADLoader");
414     loader->LoadSDigits("READ");
415   
416     // Get the tree of summable digits
417     TTree* sdigitsTree = loader->TreeS();
418     if (!sdigitsTree)  {
419       AliError("No sdigit tree from digInput");
420       continue;
421     }
422
423     // Get the branch 
424     TBranch* sdigitsBranch = sdigitsTree->GetBranch("ADSDigit");
425     if (!sdigitsBranch) {
426       AliError("Failed to get sdigit branch");
427       return;
428     }
429
430     // Set the branch address
431     TClonesArray *sdigitsArray = NULL;
432     sdigitsBranch->SetAddress(&sdigitsArray);
433
434     // Sum contributions from the sdigits
435     // Get number of entries in the tree 
436     Int_t nentries  = Int_t(sdigitsBranch->GetEntries());
437     for (Int_t entry = 0; entry < nentries; ++entry)  {
438       sdigitsBranch->GetEntry(entry);
439       // Get the number of sdigits 
440       Int_t nsdigits = sdigitsArray->GetEntries();
441       
442       for (Int_t sdigit = 0; sdigit < nsdigits; sdigit++) {
443         AliADSDigit* sDigit = static_cast<AliADSDigit*>(sdigitsArray->UncheckedAt(sdigit));
444         Int_t pmNumber = sDigit->PMNumber();
445         Int_t nbins = sDigit->GetNBins();
446         if (nbins != fNBins[pmNumber]) {
447           AliError(Form("Incompatible number of bins between digitizer (%d) and sdigit (%d) for PM %d! Skipping sdigit!",
448                         fNBins[pmNumber],nbins,pmNumber));
449           continue;
450         }
451         // Sum the charges
452         Float_t *charges = sDigit->GetCharges();
453         for(Int_t iBin = 0; iBin < nbins; ++iBin) fTime[pmNumber][iBin] += charges[iBin];
454         // and the labels
455         Int_t *labels = sDigit->GetTracks();
456         Int_t j = 0;
457         for(Int_t i = 0; i < 3; ++i) {
458           if (fLabels[pmNumber][i] < 0) {
459             if (labels[j] < 0) break;
460             fLabels[pmNumber][i] = labels[j];
461             j++;
462           }
463         }
464       }
465     }
466     loader->UnloadSDigits();
467   }
468 }
469
470
471 //____________________________________________________________________________
472 void AliADDigitizer::WriteDigits(AliLoader *loader)
473 {
474   // Take fAdc arrays filled by the previous
475   // method and produce and add digits to the digit Tree
476
477   loader->LoadDigits("UPDATE");
478
479   if (!loader->TreeD()) loader->MakeTree("D");
480   loader->MakeDigitsContainer();
481   TTree* treeD  = loader->TreeD();
482   DigitsArray();
483   treeD->Branch("ADDigit", &fDigits); 
484   
485   Short_t *chargeADC = new Short_t[kNClocks];
486   for (Int_t i=0; i<16; i++) {      
487     for (Int_t j = 0; j < kNClocks; ++j) {
488       Int_t tempadc = Int_t(fAdc[i][j]);
489       if (tempadc > 1023) tempadc = 1023;
490       chargeADC[j] = tempadc;
491     }
492     AddDigit(i, fLeadingTime[i], fTimeWidth[i], Bool_t((10+fEvenOrOdd)%2), chargeADC, fBBFlag[i], fBGFlag[i], fLabels[i]);
493   }
494   delete [] chargeADC;
495
496   treeD->Fill();
497   loader->WriteDigits("OVERWRITE");  
498   loader->UnloadDigits();     
499   ResetDigits();
500 }
501
502 //____________________________________________________________________________
503 void AliADDigitizer::WriteSDigits(AliLoader *loader)
504 {
505   // Take fTime arrays filled by the previous
506   // method and produce and add sdigits to the sdigit Tree
507
508   loader->LoadSDigits("UPDATE");
509
510   if (!loader->TreeS()) loader->MakeTree("S");
511   loader->MakeSDigitsContainer();
512   TTree* treeS  = loader->TreeS();
513   SDigitsArray();
514   treeS->Branch("ADSDigit", &fDigits); 
515   //fAD->MakeBranchInTree(treeS,"AD",&fDigits,8000,"");
516   
517   for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
518     AddSDigit(ipmt,fNBins[ipmt],fTime[ipmt],fLabels[ipmt]);
519   }
520
521   treeS->Fill();
522   loader->WriteSDigits("OVERWRITE");  
523   loader->UnloadSDigits();     
524   ResetDigits();
525 }
526
527
528
529 //____________________________________________________________________________
530 void AliADDigitizer::AddDigit(Int_t pmnumber, Float_t time, Float_t width, Bool_t integrator, Short_t *chargeADC, Bool_t bbFlag, Bool_t bgFlag, Int_t *labels) 
531  { 
532  
533 // Adds Digit 
534  
535   TClonesArray &ldigits = *fDigits;  
536          
537   new(ldigits[fNdigits++]) AliADdigit(pmnumber,time,width,integrator,chargeADC,bbFlag,bgFlag,labels);
538          
539 }
540 //____________________________________________________________________________
541 void AliADDigitizer::AddSDigit(Int_t pmnumber, Int_t nbins, Float_t *charges, Int_t *labels) 
542  { 
543  
544 // Adds SDigit 
545  
546   TClonesArray &ldigits = *fDigits;  
547          
548   new(ldigits[fNdigits++]) AliADSDigit(pmnumber,nbins,charges,labels);
549          
550 }
551 //____________________________________________________________________________
552 void AliADDigitizer::ResetDigits()
553 {
554
555 // Clears Digits
556
557   fNdigits = 0;
558   if (fDigits) fDigits->Clear();
559 }
560
561 //____________________________________________________________________________
562 AliADCalibData* AliADDigitizer::GetCalibData() const
563
564 {
565 AliCDBManager *man = AliCDBManager::Instance();
566
567   AliCDBEntry *entry=0;
568
569   entry = man->Get("AD/Calib/Data");
570   if(!entry){
571     AliWarning("Load of calibration data from default storage failed!");
572     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
573         
574     man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
575     man->SetRun(1);
576     entry = man->Get("AD/Calib/Data");
577   }
578   // Retrieval of data in directory AD/Calib/Data:
579
580   AliADCalibData *calibdata = 0;
581
582   if (entry) calibdata = (AliADCalibData*) entry->GetObject();
583   if (!calibdata)  AliFatal("No calibration data from calibration database !");
584
585   return calibdata;
586
587 }
588
589 //____________________________________________________________________________
590 double AliADDigitizer::SignalShape(double *x, double *par)
591 {
592   // this function simulates the time
593   // of arrival of the photons at the
594   // photocathode
595   Double_t xx = x[0];
596   if (xx <= par[0]) return 0;
597   Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
598   if (xx <= par[3]) return a;
599   Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
600   Double_t f = a*b/(a+b);
601   AliDebug(100,Form("x=%f func=%f",xx,f));
602   return f;
603 }
604
605 //____________________________________________________________________________
606 double AliADDigitizer::PMResponse(double *x, double * /* par */)
607 {
608   // this function describes the
609   // PM time response to a single
610   // photoelectron
611   Double_t xx = x[0]+kPMRespTime;
612   return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
613 }
614
615 //____________________________________________________________________________
616 double AliADDigitizer::SinglePhESpectrum(double *x, double * /* par */)
617 {
618   // this function describes the
619   // PM amplitude response to a single
620   // photoelectron
621   Double_t xx = x[0];
622   if (xx < 0) return 0;
623   return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
624 }
625 //____________________________________________________________________
626 TClonesArray* AliADDigitizer::DigitsArray() 
627 {
628   // Initialize digit array if not already and
629   // return pointer to it. 
630   if (!fDigits) { 
631     fDigits = new TClonesArray("AliADdigit", 16);
632     fNdigits = 0;
633   }
634   return fDigits;
635 }
636
637 //____________________________________________________________________
638 TClonesArray* AliADDigitizer::SDigitsArray() 
639 {
640   // Initialize sdigit array if not already and
641   // return pointer to it. 
642   if (!fDigits) { 
643     fDigits = new TClonesArray("AliADSDigit", 16);
644     fNdigits = 0;
645   }
646   return fDigits;
647 }
648