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