1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
20 //////////////////////////////////////////////////////////////////////////////
21 // Class performs digitization of Summable digits from simulated data
23 // In addition it performs mixing/embedding of summable digits from different events.
25 // For each event 3 branches are created in TreeD:
26 // "EMCAL" - list of digits
27 // "EMCALTRG" - list of trigger digits
28 // "AliEMCALDigitizer" - AliEMCALDigitizer with all parameters used in digitization
31 ////////////////////////////////////////////////////////////////////////////////////
33 //*-- Author: Sahal Yacoob (LBL)
34 // based on : AliEMCALDigitizer
36 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
37 // of new IO (a la PHOS)
38 // November 2003 Aleksei Pavlinov : adopted for Shish-Kebab geometry
39 // July 2011 GCB: Digitizer modified to accomodate embedding.
40 // Time calibration added. Decalibration possibility of energy and time added
41 //_________________________________________________________________________________
43 // --- ROOT system ---
47 #include <TBenchmark.h>
49 #include <TObjectTable.h>
54 // --- AliRoot header files ---
57 #include "AliDigitizationInput.h"
58 #include "AliRunLoader.h"
59 #include "AliCDBManager.h"
60 #include "AliCDBEntry.h"
61 #include "AliEMCALDigit.h"
63 #include "AliEMCALLoader.h"
64 #include "AliEMCALDigitizer.h"
65 #include "AliEMCALSDigitizer.h"
66 #include "AliEMCALGeometry.h"
67 #include "AliEMCALCalibData.h"
68 #include "AliEMCALSimParam.h"
69 #include "AliEMCALRawDigit.h"
70 #include "AliCaloCalibPedestal.h"
74 Double_t HeavisideTheta(Double_t x)
78 if (x > 0.) signal = 1.;
83 Double_t AnalogFastORFunction(Double_t *x, Double_t *par)
96 return (((0.8*(-((TMath::Power(C1,2)*C2*TMath::Power(TMath::E(),(-t + t0)/(C1*R1))*
97 TMath::Power(R1,2)*R2)/(C1*R1 - C2*R2)) +
98 C1*C2*R1*R2*(1 - (C2*TMath::Power(TMath::E(),(-t + t0)/(C2*R2))*R2)/(-(C1*R1) + C2*R2)))*v0*
99 HeavisideTheta(t - t0))/tr
100 - (0.8*(C1*C2*R1*R2 -
101 (TMath::Power(C1,2)*C2*TMath::Power(TMath::E(),(-1.*t + t0 + 1.25*tr)/(C1*R1))*
102 TMath::Power(R1,2)*R2)/(C1*R1 - C2*R2) +
103 (C1*TMath::Power(C2,2)*TMath::Power(TMath::E(),(-1.*t + t0 + 1.25*tr)/(C2*R2))*
104 R1*TMath::Power(R2,2))/(C1*R1 - C2*R2))*v0*
105 HeavisideTheta(t - t0 - 1.25*tr))/tr)/(C2*R1));
109 ClassImp(AliEMCALDigitizer)
112 //____________________________________________________________________________
113 AliEMCALDigitizer::AliEMCALDigitizer()
114 : AliDigitizer("",""),
119 fInputFileNames(0x0),
122 fMeanPhotonElectron(0),
123 fGainFluctuations(0),
127 fTimeResolutionPar0(0),
128 fTimeResolutionPar1(0),
131 fADCchannelECDecal(0),
133 fTimeChannelDecal(0),
135 fEventFolderName(""),
143 fDigInput = 0 ; // We work in the standalone mode
146 //____________________________________________________________________________
147 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName)
148 : AliDigitizer("EMCALDigitizer", alirunFileName),
149 fDefaultInit(kFALSE),
156 fMeanPhotonElectron(0),
157 fGainFluctuations(0),
161 fTimeResolutionPar0(0),
162 fTimeResolutionPar1(0),
165 fADCchannelECDecal(0),
167 fTimeChannelDecal(0),
169 fEventFolderName(eventFolderName),
178 fDigInput = 0 ; // We work in the standalone mode
181 //____________________________________________________________________________
182 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d)
183 : AliDigitizer(d.GetName(),d.GetTitle()),
184 fDefaultInit(d.fDefaultInit),
185 fDigitsInRun(d.fDigitsInRun),
188 fInputFileNames(d.fInputFileNames),
189 fEventNames(d.fEventNames),
190 fDigitThreshold(d.fDigitThreshold),
191 fMeanPhotonElectron(d.fMeanPhotonElectron),
192 fGainFluctuations(d.fGainFluctuations),
193 fPinNoise(d.fPinNoise),
194 fTimeNoise(d.fTimeNoise),
195 fTimeDelay(d.fTimeDelay),
196 fTimeResolutionPar0(d.fTimeResolutionPar0),
197 fTimeResolutionPar1(d.fTimeResolutionPar1),
198 fADCchannelEC(d.fADCchannelEC),
199 fADCpedestalEC(d.fADCpedestalEC),
200 fADCchannelECDecal(d.fADCchannelECDecal),
201 fTimeChannel(d.fTimeChannel), fTimeChannelDecal(d.fTimeChannelDecal),
203 fEventFolderName(d.fEventFolderName),
204 fFirstEvent(d.fFirstEvent),
205 fLastEvent(d.fLastEvent),
206 fCalibData(d.fCalibData),
207 fSDigitizer(d.fSDigitizer ? new AliEMCALSDigitizer(*d.fSDigitizer) : 0)
212 //____________________________________________________________________________
213 AliEMCALDigitizer::AliEMCALDigitizer(AliDigitizationInput * rd)
214 : AliDigitizer(rd,"EMCALDigitizer"),
215 fDefaultInit(kFALSE),
222 fMeanPhotonElectron(0),
223 fGainFluctuations(0),
227 fTimeResolutionPar0(0.),
228 fTimeResolutionPar1(0.),
231 fADCchannelECDecal(0),
233 fTimeChannelDecal(0),
241 // ctor Init() is called by RunDigitizer
243 fEventFolderName = fDigInput->GetInputFolderName(0) ;
244 SetTitle(dynamic_cast<AliStream*>(fDigInput->GetInputStream(0))->GetFileName(0));
248 //____________________________________________________________________________
249 AliEMCALDigitizer::~AliEMCALDigitizer()
252 delete [] fInputFileNames ;
253 delete [] fEventNames ;
254 if (fSDigitizer) delete fSDigitizer;
257 //____________________________________________________________________________
258 void AliEMCALDigitizer::Digitize(Int_t event)
261 // Makes the digitization of the collected summable digits
262 // for this it first creates the array of all EMCAL modules
263 // filled with noise and after that adds contributions from
264 // SDigits. This design helps to avoid scanning over the
265 // list of digits to add contribution of any new SDigit.
268 // Note that SDigit energy info is stored as an amplitude, so we
269 // must call the Calibrate() method of the SDigitizer to convert it
270 // back to an energy in GeV before adding it to the Digit
272 static int nEMC=0; //max number of digits possible
273 AliRunLoader *rl = AliRunLoader::Instance();
274 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
277 Int_t readEvent = event ;
278 // fDigInput is data member from AliDigitizer
280 readEvent = dynamic_cast<AliStream*>(fDigInput->GetInputStream(0))->GetCurrentEventNumber() ;
281 AliDebug(1,Form("Adding event %d from input stream 0 %s %s",
282 readEvent, GetTitle(), fEventFolderName.Data())) ;
283 rl->GetEvent(readEvent);
285 TClonesArray * digits = emcalLoader->Digits() ;
286 digits->Delete() ; //correct way to clear array when memory is
287 //allocated by objects stored in array
290 AliEMCALGeometry *geom = 0;
291 if (!rl->GetAliRun()) {
292 AliFatal("Could not get AliRun from runLoader");
295 AliEMCAL * emcal = (AliEMCAL*)rl->GetAliRun()->GetDetector("EMCAL");
296 geom = emcal->GetGeometry();
297 nEMC = geom->GetNCells();
298 AliDebug(1,Form("nEMC %i (number cells in EMCAL) | %s \n", nEMC, geom->GetName()));
303 digits->Expand(nEMC) ;
305 // RS create a digitizer on the fly
306 if (!fSDigitizer) fSDigitizer = new AliEMCALSDigitizer(rl->GetFileName().Data());
307 fSDigitizer->SetEventRange(0, -1) ;
309 //take all the inputs to add together and load the SDigits
310 TObjArray * sdigArray = new TObjArray(fInput) ;
311 sdigArray->AddAt(emcalLoader->SDigits(), 0) ;
313 Int_t embed = kFALSE;
314 for(i = 1 ; i < fInput ; i++){
315 TString tempo(fEventNames[i]) ;
318 AliRunLoader * rl2 = AliRunLoader::GetRunLoader(tempo) ;
321 rl2 = AliRunLoader::Open(fInputFileNames[i], tempo) ;
324 readEvent = dynamic_cast<AliStream*>(fDigInput->GetInputStream(i))->GetCurrentEventNumber() ;
325 Info("Digitize", "Adding event %d from input stream %d %s %s", readEvent, i, fInputFileNames[i].Data(), tempo.Data()) ;
327 // rl2->LoadDigits();
328 rl2->GetEvent(readEvent);
329 if(rl2->GetDetectorLoader("EMCAL"))
331 AliEMCALLoader *emcalLoader2 = dynamic_cast<AliEMCALLoader*>(rl2->GetDetectorLoader("EMCAL"));
334 //Merge 2 simulated sdigits
335 if(emcalLoader2->SDigits()){
336 TClonesArray* sdigits2 = emcalLoader2->SDigits();
337 sdigArray->AddAt(sdigits2, i) ;
338 // Check if first sdigit is of embedded type, if so, handle the sdigits differently:
339 // do not smear energy of embedded, do not add noise to any sdigits
340 if(sdigits2->GetEntriesFast()>0){
341 //printf("Merged digit type: %d\n",dynamic_cast<AliEMCALDigit*> (sdigits2->At(0))->GetType());
342 AliEMCALDigit * digit2 = dynamic_cast<AliEMCALDigit*> (sdigits2->At(0));
343 if(digit2 && digit2->GetType()==AliEMCALDigit::kEmbedded){
350 else AliFatal("EMCAL Loader of second event not available!");
353 else Fatal("Digitize", "Loader of second input not found");
356 //Find the first tower with signal
357 Int_t nextSig = nEMC + 1 ;
358 TClonesArray * sdigits ;
359 for(i = 0 ; i < fInput ; i++){
360 if(i > 0 && embed) continue;
361 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
363 if (sdigits->GetEntriesFast() ){
364 AliEMCALDigit *sd = dynamic_cast<AliEMCALDigit *>(sdigits->At(0));
366 Int_t curNext = sd->GetId() ;
367 if(curNext < nextSig)
369 AliDebug(1,Form("input %i : #sdigits %i \n",
370 i, sdigits->GetEntriesFast()));
372 }//First entry is not NULL
374 Info("Digitize", "NULL sdigit pointer");
377 }//There is at least one entry
379 Info("Digitize", "NULL sdigits array");
382 }// SDigits array exists
384 Info("Digitizer","Null sdigit pointer");
388 AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
390 TArrayI index(fInput) ;
391 index.Reset() ; //Set all indexes to zero
393 AliEMCALDigit * digit ;
394 AliEMCALDigit * curSDigit ;
396 //Put Noise contribution, smear time and energy
397 Float_t timeResolution = 0;
398 for(absID = 0; absID < nEMC; absID++){ // Nov 30, 2006 by PAI; was from 1 to nEMC
400 if (IsDead(absID)) continue; // Don't instantiate dead digits
404 // amplitude set to zero, noise will be added later
405 Float_t noiseTime = 0.;
406 if(!embed) noiseTime = TimeOfNoise(); //No need for embedded events?
407 new((*digits)[absID]) AliEMCALDigit( -1, -1, absID, 0., noiseTime,kFALSE); // absID-1->absID
408 //look if we have to add signal?
409 digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID)); // absID-1->absID
415 // Calculate time as time of the largest digit
416 Float_t time = digit->GetTime() ;
417 Float_t aTime= digit->GetAmplitude() ;
420 Int_t nInputs = fInput;
421 if(embed) nInputs = 1; // In case of embedding, merge later real digits, do not smear energy and time of data
422 for(i = 0; i< nInputs ; i++){ //loop over (possible) merge sources
423 TClonesArray* sdtclarr = dynamic_cast<TClonesArray *>(sdigArray->At(i));
425 Int_t sDigitEntries = sdtclarr->GetEntriesFast();
427 if(sDigitEntries > index[i] )
428 curSDigit = dynamic_cast<AliEMCALDigit*>(sdtclarr->At(index[i])) ;
432 //May be several digits will contribute from the same input
433 while(curSDigit && (curSDigit->GetId() == absID)){
434 //Shift primary to separate primaries belonging different inputs
435 Int_t primaryoffset ;
437 primaryoffset = fDigInput->GetMask(i) ;
440 curSDigit->ShiftPrimary(primaryoffset) ;
442 if(curSDigit->GetAmplitude()>aTime) {
443 aTime = curSDigit->GetAmplitude();
444 time = curSDigit->GetTime();
447 *digit = *digit + *curSDigit ; //adds amplitudes of each digit
451 if( sDigitEntries > index[i] )
452 curSDigit = dynamic_cast<AliEMCALDigit*>(sdtclarr->At(index[i])) ;
457 }// loop over merging sources
460 //Here we convert the summed amplitude to an energy in GeV only for simulation or mixing of simulations
461 energy = fSDigitizer->Calibrate(digit->GetAmplitude()) ; // GeV
463 // add fluctuations for photo-electron creation
464 // corrected fluctuations after comparison with beam test, Paraskevi Ganoti (06/11/2011)
465 Float_t fluct = static_cast<Float_t>((energy*fMeanPhotonElectron)/fGainFluctuations);
466 energy *= static_cast<Float_t>(gRandom->Poisson(fluct)) / fluct ;
468 //calculate and set time
469 digit->SetTime(time) ;
471 //Find next signal module
473 for(i = 0 ; i < nInputs ; i++){
474 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
477 Int_t curNext = nextSig ;
478 if(sdigits->GetEntriesFast() > index[i])
480 AliEMCALDigit * tmpdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]));
483 curNext = tmpdigit->GetId() ;
486 if(curNext < nextSig) nextSig = curNext ;
492 // add the noise now, no need for embedded events with real data
494 energy += gRandom->Gaus(0., fPinNoise) ;
498 //Now digitize the energy according to the fSDigitizer method,
499 //which merely converts the energy to an integer. Later we will
500 //check that the stored value matches our allowed dynamic ranges
501 digit->SetAmplitude(fSDigitizer->Digitize(energy)) ;
503 //Set time resolution with final energy
504 timeResolution = GetTimeResolution(energy);
505 digit->SetTime(gRandom->Gaus(digit->GetTime(),timeResolution) ) ;
506 AliDebug(10,Form(" absID %5i energy %f nextSig %5i\n",
507 absID, energy, nextSig));
509 digit->SetTime(digit->GetTime()+fTimeDelay) ;
511 }// Digit pointer exists
513 Info("Digitizer","Digit pointer is null");
515 } // for(absID = 0; absID < nEMC; absID++)
520 //Embed simulated amplitude (and time?) to real data digits
522 for(Int_t input = 1; input<fInput; input++){
523 TClonesArray *realDigits = dynamic_cast<TClonesArray*> (sdigArray->At(input));
525 AliDebug(1,"No sdigits to merge\n");
528 for(Int_t i2 = 0 ; i2 < realDigits->GetEntriesFast() ; i2++){
529 AliEMCALDigit * digit2 = dynamic_cast<AliEMCALDigit*>( realDigits->At(i2) ) ;
531 digit = dynamic_cast<AliEMCALDigit*>( digits->At(digit2->GetId()) ) ;
534 // Put the embedded cell energy in same units as simulated sdigits -> transform to energy units
535 // and multiply to get a big int.
536 Float_t time2 = digit2->GetTime();
537 Float_t e2 = digit2->GetAmplitude();
538 CalibrateADCTime(e2,time2,digit2->GetId());
539 Float_t amp2 = fSDigitizer->Digitize(e2);
541 Float_t e0 = digit ->GetAmplitude();
543 AliDebug(1,Form("digit 1: Abs ID %d, amp %f, type %d, time %e; digit2: Abs ID %d, amp %f, type %d, time %e\n",
544 digit ->GetId(),digit ->GetAmplitude(), digit ->GetType(), digit->GetTime(),
545 digit2->GetId(),amp2, digit2->GetType(), time2 ));
547 // Sum amplitudes, change digit amplitude
548 digit->SetAmplitude( digit->GetAmplitude() + amp2);
549 // Rough assumption, assign time of the largest of the 2 digits,
550 // data or signal to the final digit.
551 if(amp2 > digit->GetAmplitude()) digit->SetTime(time2);
552 //Mark the new digit as embedded
553 digit->SetType(AliEMCALDigit::kEmbedded);
555 if(digit2->GetAmplitude()>0.01 && e0> 0.01 )
556 AliDebug(1,Form("Embedded digit: Abs ID %d, amp %f, type %d\n",
557 digit->GetId(), digit->GetAmplitude(), digit->GetType()));
564 //JLK is it better to call Clear() here?
565 delete sdigArray ; //We should not delete its contents
567 //remove digits below thresholds
568 // until 10-02-2010 remove digits with energy smaller than fDigitThreshold 3*fPinNoise
569 // now, remove digits with Digitized ADC smaller than fDigitThreshold = 3,
570 // before merge in the same loop real data digits if available
573 for(i = 0 ; i < nEMC ; i++){
574 digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
577 //Then get the energy in GeV units.
578 energy = fSDigitizer->Calibrate(digit->GetAmplitude()) ;
579 //Then digitize using the calibration constants of the ocdb
580 Float_t ampADC = energy;
581 DigitizeEnergyTime(ampADC, time, digit->GetId()) ;
582 if(ampADC < fDigitThreshold)
583 digits->RemoveAt(i) ;
590 Int_t ndigits = digits->GetEntriesFast() ;
593 //After we have done the summing and digitizing to create the
594 //digits, now we want to calibrate the resulting amplitude to match
595 //the dynamic range of our real data.
596 for (i = 0 ; i < ndigits ; i++) {
597 digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ;
599 digit->SetIndexInList(i) ;
600 energy = fSDigitizer->Calibrate(digit->GetAmplitude()) ;
601 time = digit->GetTime();
602 Float_t ampADC = energy;
603 DigitizeEnergyTime(ampADC, time, digit->GetId());
604 digit->SetAmplitude(ampADC) ;
605 digit->SetTime(time);
606 // printf("digit amplitude set at end: i %d, amp %f\n",i,digit->GetAmplitude());
611 else AliFatal("EMCALLoader is NULL!") ;
615 //_____________________________________________________________________
616 void AliEMCALDigitizer::DigitizeEnergyTime(Float_t & energy, Float_t & time, const Int_t absId)
619 // Returns digitized value of the energy in a cell absId
620 // using the calibration constants stored in the OCDB
621 // or default values if no CalibData object is found.
622 // This effectively converts everything to match the dynamic range
623 // of the real data we will collect
626 const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
629 AliFatal("Did not get geometry from EMCALLoader");
639 Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
641 Error("DigitizeEnergyTime","Wrong cell id number : absId %i ", absId) ;
642 geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
645 fADCpedestalEC = fCalibData->GetADCpedestal (iSupMod,ieta,iphi);
646 fADCchannelEC = fCalibData->GetADCchannel (iSupMod,ieta,iphi);
647 fADCchannelECDecal = fCalibData->GetADCchannelDecal (iSupMod,ieta,iphi);
648 fTimeChannel = fCalibData->GetTimeChannel (iSupMod,ieta,iphi,0); // Assign bunch crossing number equal to 0 (has simulation different bunch crossings?)
649 fTimeChannelDecal = fCalibData->GetTimeChannelDecal(iSupMod,ieta,iphi);
652 //Apply calibration to get ADC counts and partial decalibration as especified in OCDB
653 energy = (energy + fADCpedestalEC)/fADCchannelEC/fADCchannelECDecal ;
654 time += fTimeChannel-fTimeChannelDecal;
656 if(energy > fNADCEC )
660 //_____________________________________________________________________
661 void AliEMCALDigitizer::Decalibrate(AliEMCALDigit *digit)
664 const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
667 AliFatal("Did not get geometry from EMCALLoader");
671 Int_t absId = digit->GetId();
679 Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
681 if (!bCell) Error("Decalibrate","Wrong cell id number : absId %i ", absId) ;
682 geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
685 fADCchannelEC = fCalibData->GetADCchannel(iSupMod,ieta,iphi);
686 float factor = 1./(fADCchannelEC/0.0162);
688 *digit = *digit * factor;
692 //_____________________________________________________________________
693 void AliEMCALDigitizer::CalibrateADCTime(Float_t & adc, Float_t & time, const Int_t absId)
695 // Returns the energy in a cell absId with a given adc value
696 // using the calibration constants stored in the OCDB. Time also corrected from parameter in OCDB
697 // Used in case of embedding, transform ADC counts from real event into energy
698 // so that we can add the energy of the simulated sdigits which are in energy
700 // Same as in AliEMCALClusterizer::Calibrate() but here we do not reject channels being marked as hot
701 // or with time out of window
704 const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
707 AliFatal("Did not get geometry from EMCALLoader");
717 Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
719 Error("CalibrateADCTime","Wrong cell id number : absId %i ", absId) ;
720 geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
723 fADCpedestalEC = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
724 fADCchannelEC = fCalibData->GetADCchannel (iSupMod,ieta,iphi);
725 fTimeChannel = fCalibData->GetTimeChannel(iSupMod,ieta,iphi,0);// Assign bunch crossing number equal to 0 (has simulation different bunch crossings?)
728 adc = adc * fADCchannelEC - fADCpedestalEC;
729 time -= fTimeChannel;
735 //____________________________________________________________________________
736 void AliEMCALDigitizer::Digitize(Option_t *option)
738 // Steering method to process digitization for events
739 // in the range from fFirstEvent to fLastEvent.
740 // This range is optionally set by SetEventRange().
741 // if fLastEvent=-1, then process events until the end.
742 // by default fLastEvent = fFirstEvent (process only one event)
744 if (!fInit) { // to prevent overwrite existing file
745 Error( "Digitize", "Give a version name different from %s", fEventFolderName.Data() ) ;
749 if (strstr(option,"print")) {
755 if(strstr(option,"tim"))
756 gBenchmark->Start("EMCALDigitizer");
758 AliRunLoader *rl = AliRunLoader::Instance();
759 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
762 AliFatal("Did not get the Loader");
766 if (fLastEvent == -1) {
767 fLastEvent = rl->GetNumberOfEvents() - 1 ;
770 fLastEvent = fFirstEvent ; // what is this ??
772 nEvents = fLastEvent - fFirstEvent + 1;
775 TClonesArray* digitsTMP = new TClonesArray("AliEMCALDigit", 32*96);
776 TClonesArray* digitsTRG = new TClonesArray("AliEMCALRawDigit", 32*96);
778 rl->LoadSDigits("EMCAL");
779 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
781 rl->GetEvent(ievent);
783 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
788 //-------------------------------------
791 Digits2FastOR(digitsTMP, digitsTRG);
793 WriteDigits(digitsTRG);
795 (emcalLoader->TreeD())->Fill();
797 emcalLoader->WriteDigits( "OVERWRITE");
801 digitsTRG ->Delete();
802 digitsTMP ->Delete();
804 //-------------------------------------
806 if(strstr(option,"deb"))
808 if(strstr(option,"table")) gObjectTable->Print();
810 //increment the total number of Digits per run
811 fDigitsInRun += emcalLoader->Digits()->GetEntriesFast() ;
816 if(strstr(option,"tim")){
817 gBenchmark->Stop("EMCALDigitizer");
818 Float_t cputime = gBenchmark->GetCpuTime("EMCALDigitizer");
819 Float_t avcputime = cputime;
820 if(nEvents==0) avcputime = 0 ;
821 AliInfo(Form("Digitize: took %f seconds for Digitizing %f seconds per event", cputime, avcputime)) ;
825 //__________________________________________________________________
826 Float_t AliEMCALDigitizer::GetTimeResolution(Float_t energy) const
831 res = TMath::Sqrt(fTimeResolutionPar0 +
832 fTimeResolutionPar1/(energy*energy) );
834 // parametrization above is for ns. Convert to seconds:
838 //____________________________________________________________________________
839 void AliEMCALDigitizer::Digits2FastOR(TClonesArray* digitsTMP, TClonesArray* digitsTRG)
841 // FEE digits afterburner to produce TRG digits
842 // we are only interested in the FEE digit deposited energy
843 // to be converted later into a voltage value
845 // push the FEE digit to its associated FastOR (numbered from 0:95)
846 // TRU is in charge of summing module digits
848 AliRunLoader *runLoader = AliRunLoader::Instance();
850 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
851 if (!emcalLoader) AliFatal("Did not get the Loader");
853 const AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
855 // build FOR from simulated digits
856 // and xfer to the corresponding TRU input (mapping)
858 TClonesArray* sdigits = emcalLoader->SDigits();
860 TClonesArray *digits = (TClonesArray*)sdigits->Clone();
862 AliDebug(999,Form("=== %d SDigits to trigger digits ===",digits->GetEntriesFast()));
864 TIter NextDigit(digits);
865 while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
867 if (IsDead(digit)) continue;
871 Int_t id = digit->GetId();
874 if (geom && geom->GetFastORIndexFromCellIndex(id , trgid))
876 AliDebug(1,Form("trigger digit id: %d from cell id: %d\n",trgid,id));
878 AliEMCALDigit* d = static_cast<AliEMCALDigit*>(digitsTMP->At(trgid));
882 new((*digitsTMP)[trgid]) AliEMCALDigit(*digit);
883 d = (AliEMCALDigit*)digitsTMP->At(trgid);
893 if (AliDebugLevel()) printf("Number of TRG digits: %d\n",digitsTMP->GetEntriesFast());
896 Int_t *timeSamples = new Int_t[nSamples];
898 NextDigit = TIter(digitsTMP);
899 while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
903 Int_t id = digit->GetId();
904 Float_t time = 50.e-9;
906 Double_t depositedEnergy = 0.;
907 for (Int_t j = 1; j <= digit->GetNprimary(); j++) depositedEnergy += digit->GetDEPrimary(j);
909 if (AliDebugLevel()) printf("Deposited Energy: %f\n", depositedEnergy);
911 // FIXME: Check digit time!
912 if (depositedEnergy) {
913 depositedEnergy += gRandom->Gaus(0., .08);
914 DigitalFastOR(time, depositedEnergy, timeSamples, nSamples);
916 for (Int_t j=0;j<nSamples;j++) {
917 if (AliDebugLevel()) printf("timeSamples[%d]: %d\n",j,timeSamples[j]);
918 timeSamples[j] = ((j << 16) | (timeSamples[j] & 0xFFFF));
921 new((*digitsTRG)[digitsTRG->GetEntriesFast()]) AliEMCALRawDigit(id, timeSamples, nSamples);
923 if (AliDebugLevel()) ((AliEMCALRawDigit*)digitsTRG->At(digitsTRG->GetEntriesFast() - 1))->Print("");
928 delete [] timeSamples;
930 if (digits && digits->GetEntriesFast()) digits->Delete();
933 //____________________________________________________________________________
934 void AliEMCALDigitizer::DigitalFastOR( Double_t time, Double_t dE, Int_t timeSamples[], Int_t nSamples )
938 const Int_t reso = 12; // 11-bit resolution ADC
939 const Double_t vFSR = 2.; // Full scale input voltage range 2V (p-p)
940 // const Double_t dNe = 125; // signal of the APD per MeV of energy deposit in a tower: 125 photo-e-/MeV @ M=30
941 const Double_t dNe = 125/1.3; // F-ALTRO max V. FEE: factor 4
942 const Double_t vA = .136e-6; // CSP output range: 0.136uV/e-
943 const Double_t rise = 50e-9; // rise time (10-90%) of the FastOR signal before shaping
945 const Double_t kTimeBinWidth = 25E-9; // sampling frequency (40MHz)
947 Double_t vV = 1000. * dE * dNe * vA; // GeV 2 MeV
949 TF1 signalF("signal", AnalogFastORFunction, 0, nSamples * kTimeBinWidth, 3);
950 signalF.SetParameter( 0, vV );
951 signalF.SetParameter( 1, time ); // FIXME: when does the signal arrive? Might account for cable lengths
952 signalF.SetParameter( 2, rise );
954 for (Int_t iTime=0; iTime<nSamples; iTime++)
956 // FIXME: add noise (probably not simply Gaussian) according to DA measurements
957 // probably plan an access to OCDB
958 Double_t sig = signalF.Eval(iTime * kTimeBinWidth);
959 if (TMath::Abs(sig) > vFSR/2.) {
960 AliError("Signal overflow!");
961 timeSamples[iTime] = (1 << reso) - 1;
963 AliDebug(999,Form("iTime: %d sig: %f\n",iTime,sig));
964 timeSamples[iTime] = ((1 << reso) / vFSR) * sig + 0.5;
969 //____________________________________________________________________________
970 Bool_t AliEMCALDigitizer::Init()
972 // Makes all memory allocations
974 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
976 if ( emcalLoader == 0 ) {
977 Fatal("Init", "Could not obtain the AliEMCALLoader");
982 fLastEvent = fFirstEvent ;
985 fInput = fDigInput->GetNinputs() ;
989 fInputFileNames = new TString[fInput] ;
990 fEventNames = new TString[fInput] ;
991 fInputFileNames[0] = GetTitle() ;
992 fEventNames[0] = fEventFolderName.Data() ;
994 for (index = 1 ; index < fInput ; index++) {
995 fInputFileNames[index] = dynamic_cast<AliStream*>(fDigInput->GetInputStream(index))->GetFileName(0);
996 TString tempo = fDigInput->GetInputFolderName(index) ;
997 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fDigInput
1000 //Calibration instance
1001 fCalibData = emcalLoader->CalibData();
1005 //____________________________________________________________________________
1006 void AliEMCALDigitizer::InitParameters()
1008 // Parameter initialization for digitizer
1010 // Get the parameters from the OCDB via the loader
1011 AliRunLoader *rl = AliRunLoader::Instance();
1012 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
1013 AliEMCALSimParam * simParam = 0x0;
1014 if(emcalLoader) simParam = emcalLoader->SimulationParameters();
1017 simParam = AliEMCALSimParam::GetInstance();
1018 AliWarning("Simulation Parameters not available in OCDB?");
1021 fMeanPhotonElectron = simParam->GetMeanPhotonElectron() ; // 4400; // electrons per GeV
1022 fGainFluctuations = simParam->GetGainFluctuations() ; // 15.0;
1024 fPinNoise = simParam->GetPinNoise();//0.012; // pin noise in GeV from analysis test beam data
1025 if (fPinNoise < 0.0001 )
1026 Warning("InitParameters", "No noise added\n") ;
1027 fTimeNoise = simParam->GetTimeNoise(); // 1.28E-5 s
1028 fDigitThreshold = simParam->GetDigitThreshold(); //fPinNoise * 3; // 3 * sigma
1029 fTimeResolutionPar0 = simParam->GetTimeResolutionPar0();
1030 fTimeResolutionPar1 = simParam->GetTimeResolutionPar1();
1031 fTimeDelay = simParam->GetTimeDelay(); //600e-9 ; // 600 ns
1033 // These defaults are normally not used.
1034 // Values are read from calibration database instead
1035 fADCchannelEC = 0.0153; // Update 24 Apr 2007: 250./16/1024 - width of one ADC channel in GeV
1036 fADCpedestalEC = 0.0 ; // GeV
1037 fADCchannelECDecal = 1.0; // No decalibration by default
1038 fTimeChannel = 0.0; // No time calibration by default
1039 fTimeChannelDecal = 0.0; // No time decalibration by default
1041 fNADCEC = simParam->GetNADCEC();//(Int_t) TMath::Power(2,16) ; // number of channels in Tower ADC - 65536
1043 AliDebug(2,Form("Mean Photon Electron %d, Gain Fluct. %2.1f; Noise: APD %f, Time %f; Digit Threshold %d,Time Resolution Par0 %g Par1 %g,NADCEC %d",
1044 fMeanPhotonElectron, fGainFluctuations, fPinNoise,fTimeNoise, fDigitThreshold,fTimeResolutionPar0,fTimeResolutionPar1,fNADCEC));
1048 //__________________________________________________________________
1049 void AliEMCALDigitizer::Print1(Option_t * option)
1050 { // 19-nov-04 - just for convenience
1052 PrintDigits(option);
1055 //__________________________________________________________________
1056 void AliEMCALDigitizer::Print (Option_t * ) const
1058 // Print Digitizer's parameters
1059 printf("Print: \n------------------- %s -------------", GetName() ) ;
1060 if( strcmp(fEventFolderName.Data(), "") != 0 ){
1061 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
1065 nStreams = GetNInputStreams() ;
1072 for (index = 0 ; index < nStreams ; index++) {
1073 TString tempo(fEventNames[index]) ;
1075 if ((rl = AliRunLoader::GetRunLoader(tempo)) == 0)
1076 rl = AliRunLoader::Open(fInputFileNames[index], tempo) ;
1077 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
1079 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
1080 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
1081 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
1082 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
1086 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1088 if(emcalLoader) printf("\nWriting digits to %s", emcalLoader->GetDigitsFileName().Data()) ;
1089 else printf("\nNULL LOADER");
1091 printf("\nWith following parameters:\n") ;
1092 printf(" Electronics noise in EMC, APD (fPinNoise) = %f, Time = %f \n", fPinNoise, fTimeNoise) ;
1093 printf(" Threshold in Tower (fDigitThreshold) = %d\n", fDigitThreshold) ;
1094 printf("---------------------------------------------------\n") ;
1097 printf("Print: AliEMCALDigitizer not initialized") ;
1100 //__________________________________________________________________
1101 void AliEMCALDigitizer::PrintDigits(Option_t * option)
1103 //utility method for printing digit information
1105 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1107 TClonesArray * digits = emcalLoader->Digits() ;
1108 TClonesArray * sdigits = emcalLoader->SDigits() ;
1110 printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ;
1111 printf("\n event %d", emcalLoader->GetRunLoader()->GetEventNumber());
1113 if(strstr(option,"all")){
1115 AliEMCALDigit * digit;
1116 printf("\nEMC digits (with primaries):\n") ;
1117 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
1119 for (index = 0 ; index < digits->GetEntries() ; index++) {
1120 digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
1122 printf("\n%6d %8f %6.5e %4d %2d : ",
1123 digit->GetId(), digit->GetAmplitude(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
1125 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
1126 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
1133 else printf("NULL LOADER, cannot print\n");
1136 //__________________________________________________________________
1137 Float_t AliEMCALDigitizer::TimeOfNoise(void)
1139 // Calculates the time signal generated by noise
1140 //printf("Time noise %e\n",fTimeNoise);
1141 return gRandom->Rndm() * fTimeNoise;
1144 //__________________________________________________________________
1145 void AliEMCALDigitizer::Unload()
1147 // Unloads the SDigits and Digits
1151 for(i = 1 ; i < fInput ; i++){
1152 TString tempo(fEventNames[i]) ;
1154 if ((rl = AliRunLoader::GetRunLoader(tempo)))
1155 rl->GetDetectorLoader("EMCAL")->UnloadSDigits() ;
1157 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1158 if(emcalLoader)emcalLoader->UnloadDigits() ;
1159 else AliFatal("Did not get the loader");
1162 //_________________________________________________________________________________________
1163 void AliEMCALDigitizer::WriteDigits()
1164 { // Makes TreeD in the output file.
1165 // Check if branch already exists:
1166 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
1167 // already existing branches.
1168 // else creates branch with Digits, named "EMCAL", title "...",
1169 // and branch "AliEMCALDigitizer", with the same title to keep all the parameters
1170 // and names of files, from which digits are made.
1172 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1175 const TClonesArray * digits = emcalLoader->Digits() ;
1176 TTree * treeD = emcalLoader->TreeD();
1178 emcalLoader->MakeDigitsContainer();
1179 treeD = emcalLoader->TreeD();
1182 // -- create Digits branch
1183 Int_t bufferSize = 32000 ;
1184 TBranch * digitsBranch = 0;
1185 if ((digitsBranch = treeD->GetBranch("EMCAL"))) {
1186 digitsBranch->SetAddress(&digits);
1187 AliWarning("Digits Branch already exists. Not all digits will be visible");
1190 treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
1191 //digitsBranch->SetTitle(fEventFolderName);
1195 emcalLoader->WriteDigits("OVERWRITE");
1196 emcalLoader->WriteDigitizer("OVERWRITE");
1202 else AliFatal("Loader not available");
1205 //__________________________________________________________________
1206 void AliEMCALDigitizer::WriteDigits(TClonesArray* digits, const char* branchName)
1207 { // overloaded method
1208 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1211 TTree* treeD = emcalLoader->TreeD();
1214 emcalLoader->MakeDigitsContainer();
1215 treeD = emcalLoader->TreeD();
1218 // -- create Digits branch
1219 Int_t bufferSize = 32000;
1221 if (TBranch* triggerBranch = treeD->GetBranch(branchName))
1223 triggerBranch->SetAddress(&digits);
1227 treeD->Branch(branchName,"TClonesArray",&digits,bufferSize);
1232 else AliFatal("Loader not available");
1235 //__________________________________________________________________
1236 Bool_t AliEMCALDigitizer::IsDead(AliEMCALDigit *digit)
1238 AliRunLoader *runLoader = AliRunLoader::Instance();
1239 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
1240 if (!emcalLoader) AliFatal("Did not get the Loader");
1242 AliCaloCalibPedestal *caloPed = emcalLoader->PedestalData();
1244 AliWarning("Could not access pedestal data! No dead channel removal applied");
1249 const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
1250 if (!geom) AliFatal("Did not get geometry from EMCALLoader");
1252 Int_t absId = digit->GetId();
1260 Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
1262 if (!bCell) Error("IsDead","Wrong cell id number : absId %i ", absId) ;
1263 geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
1265 Int_t channelStatus = (Int_t)(caloPed->GetDeadMap(iSupMod))->GetBinContent(ieta,iphi);
1267 if (channelStatus == AliCaloCalibPedestal::kDead)
1274 //__________________________________________________________________
1275 Bool_t AliEMCALDigitizer::IsDead(Int_t absId)
1277 AliRunLoader *runLoader = AliRunLoader::Instance();
1278 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
1279 if (!emcalLoader) AliFatal("Did not get the Loader");
1281 AliCaloCalibPedestal *caloPed = emcalLoader->PedestalData();
1283 AliWarning("Could not access pedestal data! No dead channel removal applied");
1288 const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
1289 if (!geom) AliFatal("Did not get geometry from EMCALLoader");
1298 Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
1300 if (!bCell) Error("IsDead","Wrong cell id number : absId %i ", absId) ;
1301 geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
1303 Int_t channelStatus = (Int_t)(caloPed->GetDeadMap(iSupMod))->GetBinContent(ieta,iphi);
1305 if (channelStatus == AliCaloCalibPedestal::kDead)