]> git.uio.no Git - u/mrichter/AliRoot.git/blob - AD/AliADDigitizer.cxx
Add siscone libs (does not work without)
[u/mrichter/AliRoot.git] / AD / 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
57 ClassImp(AliADDigitizer)
58
59 //____________________________________________________________________________
60  AliADDigitizer::AliADDigitizer()
61                    :AliDigitizer(),
62                     fCalibData(GetCalibData()),
63                     fNdigits(0),
64                     fDigits(0),
65                     fSignalShape(NULL),
66                     fPMResponse(NULL),
67                     fSinglePhESpectrum(NULL),
68                     fEvenOrOdd(kFALSE),
69                     fTask(kHits2Digits),
70                     fAD(NULL)
71 {
72   // default constructor
73   // Initialize OCDB and containers used in the digitization
74
75   Init();
76 }
77
78 //____________________________________________________________________________ 
79   AliADDigitizer::AliADDigitizer(AliAD *AD, DigiTask_t task)
80                     :AliDigitizer(),
81                      fCalibData(GetCalibData()),
82                      fNdigits(0),
83                      fDigits(0),
84                      fSignalShape(NULL),
85                      fPMResponse(NULL),
86                      fSinglePhESpectrum(NULL),
87                      fEvenOrOdd(kFALSE),
88                      fTask(task),
89                      fAD(AD)
90 {
91   // constructor
92   // Initialize OCDB and containers used in the digitization
93
94   Init();
95 }
96            
97 //____________________________________________________________________________ 
98   AliADDigitizer::AliADDigitizer(AliDigitizationInput* digInput)
99                     :AliDigitizer(digInput),
100                      fCalibData(GetCalibData()),
101                      fNdigits(0),
102                      fDigits(0),
103                      fSignalShape(NULL),
104                      fPMResponse(NULL),
105                      fSinglePhESpectrum(NULL),
106                      fEvenOrOdd(kFALSE),
107                      fTask(kHits2Digits),
108                      fAD(NULL)
109 {
110   // constructor
111   // Initialize OCDB and containers used in the digitization
112
113   Init();
114 }
115            
116 //____________________________________________________________________________ 
117   AliADDigitizer::~AliADDigitizer()
118 {
119   // destructor
120   
121   if (fDigits) {
122     fDigits->Delete();
123     delete fDigits;
124     fDigits=0; 
125   }
126
127   if (fSignalShape) {
128     delete fSignalShape;
129     fSignalShape = NULL;
130   }
131   if (fPMResponse) {
132     delete fPMResponse;
133     fPMResponse = NULL;
134   }
135   if (fSinglePhESpectrum) {
136     delete fSinglePhESpectrum;
137     fSinglePhESpectrum = NULL;
138   }
139
140   for(Int_t i = 0 ; i < 16; ++i) {
141     if (fTime[i]) delete [] fTime[i];
142   }
143 }
144
145 //____________________________________________________________________________ 
146 Bool_t AliADDigitizer::Init()
147 {
148   // Initialises the digitizer
149   // Initialize OCDB and containers used in the digitization
150
151   // check if the digitizer was already initialized
152   if (fSignalShape) return kTRUE;
153
154   fSignalShape = new TF1("ADSignalShape",this,&AliADDigitizer::SignalShape,0,200,6,"AliADDigitizer","SignalShape");
155   //  fSignalShape->SetParameters(0,1.57345e1,-4.25603e-1,2.9,6.40982,3.69339e-01);
156   //  fSignalShape->SetParameters(1.34330e+00,1.13007e+02,-4.95705e-01,
157   //                          3.68911e+00,1.01040e+00, 3.94675e-01);
158   fSignalShape->SetParameters(-1.07335e+00,2.16002e+01,-1.26133e-01,
159                               1.41619e+00,5.50334e-01,3.86111e-01);
160   fPMResponse = new TF1("ADPMResponse",this,&AliADDigitizer::PMResponse,-kPMRespTime,2.*kPMRespTime,0,"AliADDigitizer","PMResponse");
161   fSinglePhESpectrum = new TF1("ADSinglePhESpectrum",this,&AliADDigitizer::SinglePhESpectrum,0,20,0,"AliADDigitizer","SinglePhESpectrum");
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("VZERO/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     fHptdcOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
202                         (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
203                        fCalibData->GetTimeOffset(i)-
204                        l1Delay-
205                        phase->GetMeanPhase()+
206                        delays->GetBinContent(i+1)+
207                        kV0Offset);
208     fClockOffset[i] = (((Float_t)fCalibData->GetRollOver(board)-
209                         (Float_t)fCalibData->GetTriggerCountOffset(board))*25.0+
210                        fCalibData->GetTimeOffset(i)-
211                        l1Delay+
212                        kV0Offset);
213
214     fTime[i] = new Float_t[fNBins[i]];
215     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
216     
217    // AliWarning(Form("PMT %d,PM gain %f, fNBins %d, TimeBinSize %f,",i, fPmGain[i], fNBins[i],fBinSize[i]));
218     
219   }
220
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            //if (pmt < 32) t += kV0CDelayCables;
314            t += fHptdcOffset[pmt];
315            Int_t nPhE;
316            Float_t prob = fCalibData->GetLightYields(pmt)*kPhotoCathodeEfficiency; // Optical losses included!
317            if (nPhot > 100)
318              nPhE = (Int_t)gRandom->Gaus(prob*Float_t(nPhot)+0.5,
319                                          sqrt(Float_t(nPhot)*prob*(1.-prob)));
320            else
321              nPhE = gRandom->Binomial(nPhot,prob);
322            Float_t charge = TMath::Qe()*fPmGain[pmt]*fBinSize[pmt]/integral;
323            
324            
325            for (Int_t iPhE = 0; iPhE < nPhE; ++iPhE) {       
326              Float_t tPhE = t + fSignalShape->GetRandom(0,fBinSize[pmt]*Float_t(fNBins[pmt]));
327              Float_t gainVar = fSinglePhESpectrum->GetRandom(0,20)/meansPhE;
328              Int_t firstBin = TMath::Max(0,(Int_t)((tPhE-kPMRespTime)/fBinSize[pmt]));
329              Int_t lastBin = TMath::Min(fNBins[pmt]-1,(Int_t)((tPhE+2.*kPMRespTime)/fBinSize[pmt]));
330              for(Int_t iBin = firstBin; iBin <= lastBin; ++iBin) {
331                Float_t tempT = fBinSize[pmt]*(0.5+iBin)-tPhE;
332                fTime[pmt][iBin] += gainVar*charge*fPMResponse->Eval(tempT);
333              }
334            }         // ph.e. loop
335          }           // hit loop
336      }               // track loop
337      loader->UnloadHits();
338 }
339
340 //____________________________________________________________________________ 
341 void AliADDigitizer::DigitizeSDigits()
342 {
343   // Digitize the fTime arrays (SDigits) to the level of
344   // Digits (fAdc arrays)
345   for(Int_t i = 0 ; i < 16; ++i) {
346     for(Int_t j = 0; j < kNClocks; ++j) fAdc[i][j] = 0;
347     fLeadingTime[i] = fTimeWidth[i] = 0;
348   }
349
350   Float_t maximum = 0.9*fSignalShape->GetMaximum(0,200); // Not exact, one needs to do this on the convoluted
351   Float_t integral2 = fSignalShape->Integral(0,200); // function. Anyway the effect is small <10% on the 2.5 ADC thr
352   for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
353     Float_t thr = fCalibData->GetCalibDiscriThr(ipmt,kFALSE)*kChargePerADC*maximum*fBinSize[ipmt]/integral2;
354     //Float_t thr = 1e-25;
355        
356     Bool_t ltFound = kFALSE, ttFound = kFALSE;
357     for (Int_t iBin = 0; iBin < fNBins[ipmt]; ++iBin) {
358       Float_t t = fBinSize[ipmt]*Float_t(iBin);
359       
360       if (fTime[ipmt][iBin] > thr) {
361         if (!ltFound && (iBin < fNBinsLT[ipmt])) {
362           ltFound = kTRUE;
363           fLeadingTime[ipmt] = t;
364         }
365       }
366       else {
367         if (ltFound) {
368           if (!ttFound) {
369             ttFound = kTRUE;
370             fTimeWidth[ipmt] = t - fLeadingTime[ipmt];
371           }
372         }
373       }
374       Float_t tadc = t - fClockOffset[ipmt];
375       Int_t clock = kNClocks/2 - Int_t(tadc/25.0);
376       if (clock >= 0 && clock < kNClocks)
377         fAdc[ipmt][clock] += fTime[ipmt][iBin]/kChargePerADC;
378     }
379     AliDebug(1,Form("Channel %d Offset %f Time %f",ipmt,fClockOffset[ipmt],fLeadingTime[ipmt]));
380     Int_t board = AliADCalibData::GetBoardNumber(ipmt);
381     if (ltFound && ttFound) {
382       fTimeWidth[ipmt] = fCalibData->GetWidthResolution(board)*
383         Float_t(Int_t(fTimeWidth[ipmt]/fCalibData->GetWidthResolution(board)));
384       if (fTimeWidth[ipmt] < Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board))
385         fTimeWidth[ipmt] = Float_t(kMinTDCWidth)*fCalibData->GetWidthResolution(board);
386       if (fTimeWidth[ipmt] > Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board))
387         fTimeWidth[ipmt] = Float_t(kMaxTDCWidth)*fCalibData->GetWidthResolution(board);
388     }
389   }
390
391   fEvenOrOdd = gRandom->Integer(2);
392   for (Int_t j=0; j<16; ++j){
393     for (Int_t iClock = 0; iClock < kNClocks; ++iClock) {
394       Int_t integrator = (iClock + fEvenOrOdd) % 2;
395       AliDebug(1,Form("ADC %d %d %f",j,iClock,fAdc[j][iClock]));
396       fAdc[j][iClock]  += gRandom->Gaus(fAdcPedestal[j][integrator], fAdcSigma[j][integrator]);
397     }
398   }
399         
400 }
401
402 //____________________________________________________________________________ 
403 void AliADDigitizer::ReadSDigits()
404 {
405   // Read SDigits which are then to precessed
406   // in the following method
407   for(Int_t i = 0 ; i < 16; ++i) {
408     memset(fTime[i],0,fNBins[i]*sizeof(Float_t));
409     fLabels[i][0] = fLabels[i][1] = fLabels[i][2] = -1;
410   }
411
412   // Loop over input files
413   Int_t nFiles= fDigInput->GetNinputs();
414   for (Int_t inputFile = 0; inputFile < nFiles; inputFile++) {
415     // Get the current loader 
416     AliRunLoader* currentLoader = 
417       AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
418
419     AliLoader *loader = currentLoader->GetLoader("ADLoader");
420     loader->LoadSDigits("READ");
421   
422     // Get the tree of summable digits
423     TTree* sdigitsTree = loader->TreeS();
424     if (!sdigitsTree)  {
425       AliError("No sdigit tree from digInput");
426       continue;
427     }
428
429     // Get the branch 
430     TBranch* sdigitsBranch = sdigitsTree->GetBranch("ADSDigit");
431     if (!sdigitsBranch) {
432       AliError("Failed to get sdigit branch");
433       return;
434     }
435
436     // Set the branch address
437     TClonesArray *sdigitsArray = NULL;
438     sdigitsBranch->SetAddress(&sdigitsArray);
439
440     // Sum contributions from the sdigits
441     // Get number of entries in the tree 
442     Int_t nentries  = Int_t(sdigitsBranch->GetEntries());
443     for (Int_t entry = 0; entry < nentries; ++entry)  {
444       sdigitsBranch->GetEntry(entry);
445       // Get the number of sdigits 
446       Int_t nsdigits = sdigitsArray->GetEntries();
447       
448       for (Int_t sdigit = 0; sdigit < nsdigits; sdigit++) {
449         AliADSDigit* sDigit = static_cast<AliADSDigit*>(sdigitsArray->UncheckedAt(sdigit));
450         Int_t pmNumber = sDigit->PMNumber();
451         Int_t nbins = sDigit->GetNBins();
452         if (nbins != fNBins[pmNumber]) {
453           AliError(Form("Incompatible number of bins between digitizer (%d) and sdigit (%d) for PM %d! Skipping sdigit!",
454                         fNBins[pmNumber],nbins,pmNumber));
455           continue;
456         }
457         // Sum the charges
458         Float_t *charges = sDigit->GetCharges();
459         for(Int_t iBin = 0; iBin < nbins; ++iBin) fTime[pmNumber][iBin] += charges[iBin];
460         // and the labels
461         Int_t *labels = sDigit->GetTracks();
462         Int_t j = 0;
463         for(Int_t i = 0; i < 3; ++i) {
464           if (fLabels[pmNumber][i] < 0) {
465             if (labels[j] < 0) break;
466             fLabels[pmNumber][i] = labels[j];
467             j++;
468           }
469         }
470       }
471     }
472     loader->UnloadSDigits();
473   }
474 }
475
476
477 //____________________________________________________________________________
478 void AliADDigitizer::WriteDigits(AliLoader *loader)
479 {
480   // Take fAdc arrays filled by the previous
481   // method and produce and add digits to the digit Tree
482
483   loader->LoadDigits("UPDATE");
484
485   if (!loader->TreeD()) loader->MakeTree("D");
486   loader->MakeDigitsContainer();
487   TTree* treeD  = loader->TreeD();
488   DigitsArray();
489   treeD->Branch("ADDigit", &fDigits); 
490   //fAD->MakeBranchInTree(treeD,"AD",&fDigits,1000,"");
491   
492   Short_t *chargeADC = new Short_t[kNClocks];
493   for (Int_t i=0; i<16; i++) {      
494     for (Int_t j = 0; j < kNClocks; ++j) {
495       Int_t tempadc = Int_t(fAdc[i][j]);
496       if (tempadc > 1023) tempadc = 1023;
497       chargeADC[j] = tempadc;
498     }
499     AddDigit(i, fLeadingTime[i], fTimeWidth[i], Bool_t((10+fEvenOrOdd)%2), chargeADC, fLabels[i]);
500   }
501   delete [] chargeADC;
502
503   treeD->Fill();
504   loader->WriteDigits("OVERWRITE");  
505   loader->UnloadDigits();     
506   ResetDigits();
507 }
508
509 //____________________________________________________________________________
510 void AliADDigitizer::WriteSDigits(AliLoader *loader)
511 {
512   // Take fTime arrays filled by the previous
513   // method and produce and add sdigits to the sdigit Tree
514
515   loader->LoadSDigits("UPDATE");
516
517   if (!loader->TreeS()) loader->MakeTree("S");
518   loader->MakeSDigitsContainer();
519   TTree* treeS  = loader->TreeS();
520   SDigitsArray();
521   treeS->Branch("ADSDigit", &fDigits); 
522   //fAD->MakeBranchInTree(treeS,"AD",&fDigits,8000,"");
523   
524   for (Int_t ipmt = 0; ipmt < 16; ++ipmt) {
525     AddSDigit(ipmt,fNBins[ipmt],fTime[ipmt],fLabels[ipmt]);
526   }
527
528   treeS->Fill();
529   loader->WriteSDigits("OVERWRITE");  
530   loader->UnloadSDigits();     
531   ResetDigits();
532 }
533
534
535
536 //____________________________________________________________________________
537 void AliADDigitizer::AddDigit(Int_t pmnumber, Float_t time, Float_t width, Bool_t integrator, Short_t *chargeADC, Int_t *labels) 
538  { 
539  
540 // Adds Digit 
541  
542   TClonesArray &ldigits = *fDigits;  
543          
544   new(ldigits[fNdigits++]) AliADdigit(pmnumber,time,width,integrator,chargeADC,labels);
545          
546 }
547 //____________________________________________________________________________
548 void AliADDigitizer::AddSDigit(Int_t pmnumber, Int_t nbins, Float_t *charges, Int_t *labels) 
549  { 
550  
551 // Adds SDigit 
552  
553   TClonesArray &ldigits = *fDigits;  
554          
555   new(ldigits[fNdigits++]) AliADSDigit(pmnumber,nbins,charges,labels);
556          
557 }
558 //____________________________________________________________________________
559 void AliADDigitizer::ResetDigits()
560 {
561
562 // Clears Digits
563
564   fNdigits = 0;
565   if (fDigits) fDigits->Clear();
566 }
567
568 //____________________________________________________________________________
569 AliADCalibData* AliADDigitizer::GetCalibData() const
570
571 {
572 /*/
573   AliCDBManager *man = AliCDBManager::Instance();
574
575   AliCDBEntry *entry=0;
576
577   entry = man->Get("AD/Calib/Data");
578
579   AliADCalibData *calibdata = 0;
580
581   if (entry) calibdata = (AliADCalibData*) entry->GetObject();
582   if (!calibdata)  AliFatal("No calibration data from calibration database !");
583 /*/
584   AliADCalibData *calibdata = new AliADCalibData();
585   
586   return calibdata;
587
588 }
589
590 //____________________________________________________________________________
591 double AliADDigitizer::SignalShape(double *x, double *par)
592 {
593   // this function simulates the time
594   // of arrival of the photons at the
595   // photocathode
596   Double_t xx = x[0];
597   if (xx <= par[0]) return 0;
598   Double_t a = 1./TMath::Power((xx-par[0])/par[1],1./par[2]);
599   if (xx <= par[3]) return a;
600   Double_t b = 1./TMath::Power((xx-par[3])/par[4],1./par[5]);
601   Double_t f = a*b/(a+b);
602   AliDebug(100,Form("x=%f func=%f",xx,f));
603   return f;
604 }
605
606 //____________________________________________________________________________
607 double AliADDigitizer::PMResponse(double *x, double * /* par */)
608 {
609   // this function describes the
610   // PM time response to a single
611   // photoelectron
612   Double_t xx = x[0]+kPMRespTime;
613   return xx*xx*TMath::Exp(-xx*xx/(kPMRespTime*kPMRespTime));
614 }
615
616 //____________________________________________________________________________
617 double AliADDigitizer::SinglePhESpectrum(double *x, double * /* par */)
618 {
619   // this function describes the
620   // PM amplitude response to a single
621   // photoelectron
622   Double_t xx = x[0];
623   if (xx < 0) return 0;
624   return (TMath::Poisson(xx,kPMNbOfSecElec)+kPMTransparency*TMath::Poisson(xx,1.0));
625 }
626 //____________________________________________________________________
627 TClonesArray* AliADDigitizer::DigitsArray() 
628 {
629   // Initialize digit array if not already and
630   // return pointer to it. 
631   if (!fDigits) { 
632     fDigits = new TClonesArray("AliADdigit", 16);
633     fNdigits = 0;
634   }
635   return fDigits;
636 }
637
638 //____________________________________________________________________
639 TClonesArray* AliADDigitizer::SDigitsArray() 
640 {
641   // Initialize sdigit array if not already and
642   // return pointer to it. 
643   if (!fDigits) { 
644     fDigits = new TClonesArray("AliADSDigit", 16);
645     fNdigits = 0;
646   }
647   return fDigits;
648 }
649