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