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 //_________________________________________________________________________
19 // This is a Class that makes SDigits out of Hits
20 // A Summable Digits is the sum of all hits originating
21 // from one in one tower of the EMCAL
22 // A threshold for assignment of the primary to SDigit is applied
24 // JLK 26-Jun-2008 Added explanation:
25 // SDigits need to hold the energy sum of the hits, but AliEMCALDigit
26 // can (should) only store amplitude. Therefore, the SDigit energy is
27 // "digitized" before being stored and must be "calibrated" back to an
28 // energy before SDigits are summed to form true Digits
30 // SDigits are written to TreeS, branch "EMCAL"
31 // AliEMCALSDigitizer with all current parameters is written
32 // to TreeS branch "AliEMCALSDigitizer".
33 // Both branches have the same title. If necessary one can produce
34 // another set of SDigits with different parameters. Two versions
35 // can be distunguished using titles of the branches.
37 // root [0] AliEMCALSDigitizer * s = new AliEMCALSDigitizer("galice.root")
38 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
39 // root [1] s->Digitize()
40 // // Makes SDigitis for all events stored in galice.root
41 // root [2] s->SetPedestalParameter(0.001)
42 // // One can change parameters of digitization
43 // root [3] s->SetSDigitsBranch("Redestal 0.001")
44 // // and write them into the new branch
45 // root [4] s->Digitize("deb all tim")
46 // // available parameters:
47 // deb - print # of produced SDigitis
48 // deb all - print # and list of produced SDigits
49 // tim - print benchmarking information
51 //*-- Author : Sahal Yacoob (LBL)
52 // based on : AliPHOSSDigitzer
53 //////////////////////////////////////////////////////////////////////////////
55 // --- ROOT system ---
56 #include <TBenchmark.h>
58 //#include <Riostream.h>
62 // --- Standard library ---
65 // --- AliRoot header files ---
67 #include "AliRunLoader.h"
69 #include "AliEMCALDigit.h"
70 #include "AliEMCALLoader.h"
71 #include "AliEMCALHit.h"
72 #include "AliEMCALSDigitizer.h"
73 #include "AliEMCALGeometry.h"
74 #include "AliEMCALSimParam.h"
76 ClassImp(AliEMCALSDigitizer)
78 //____________________________________________________________________________
79 AliEMCALSDigitizer::AliEMCALSDigitizer()
81 fA(0.),fB(0.),fECPrimThreshold(0.),
95 //____________________________________________________________________________
96 AliEMCALSDigitizer::AliEMCALSDigitizer(const char * alirunFileName,
97 const char * eventFolderName)
98 : TNamed("EMCALSDigitizer", alirunFileName),
99 fA(0.),fB(0.),fECPrimThreshold(0.),
100 fDefaultInit(kFALSE),
101 fEventFolderName(eventFolderName),
116 //____________________________________________________________________________
117 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd)
118 : TNamed(sd.GetName(),sd.GetTitle()),
121 fECPrimThreshold(sd.fECPrimThreshold),
122 fDefaultInit(sd.fDefaultInit),
123 fEventFolderName(sd.fEventFolderName),
125 fSDigitsInRun(sd.fSDigitsInRun),
126 fFirstEvent(sd.fFirstEvent),
127 fLastEvent(sd.fLastEvent),
128 fSampling(sd.fSampling),
134 //_____________________________________________________________________
135 AliEMCALSDigitizer& AliEMCALSDigitizer::operator = (const AliEMCALSDigitizer& source)
136 { // assignment operator; use copy ctor
137 if (&source == this) return *this;
139 new (this) AliEMCALSDigitizer(source);
143 //____________________________________________________________________________
144 AliEMCALSDigitizer::~AliEMCALSDigitizer() {
153 //____________________________________________________________________________
154 void AliEMCALSDigitizer::Init(){
155 // Initialization: open root-file, allocate arrays for hits and sdigits
157 // Initialization can not be done in the default constructor
158 //============================================================= YS
159 // The initialisation is now done by the getter
163 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
165 if ( emcalLoader == 0 ) {
166 Fatal("Init", "Could not obtain the AliEMCALLoader");
172 //____________________________________________________________________________
173 void AliEMCALSDigitizer::InitParameters()
175 //initialize parameters for sdigitization
177 const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
178 if (geom->GetSampling() == 0.) {
179 Fatal("InitParameters", "Sampling factor not set !") ;
182 // Get the parameters from the OCDB via the loader
183 AliRunLoader *rl = AliRunLoader::Instance();
184 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
185 AliEMCALSimParam * simParam = 0x0;
186 if(emcalLoader) simParam = emcalLoader->SimulationParameters();
189 simParam = AliEMCALSimParam::GetInstance();
190 AliWarning("Simulation Parameters not available in OCDB?");
194 //JLK 26-Jun-2008 THIS SHOULD HAVE BEEN EXPLAINED AGES AGO:
196 //In order to be able to store SDigit Energy info into
197 //AliEMCALDigit, we need to convert it temporarily to an ADC amplitude
198 //and later when summing SDigits to form digits, convert it back to
199 //energy. These fA and fB parameters accomplish this through the
200 //Digitize() and Calibrate() methods
202 // Initializes parameters
203 fA = simParam->GetA(); //0;
204 fB = simParam->GetB(); //1.e+6; // Changed 24 Apr 2007. Dynamic range now 2 TeV
205 fSampling = geom->GetSampling();
207 // threshold for deposit energy of hit
208 fECPrimThreshold = simParam->GetECPrimaryThreshold();//0.05;// GeV // 22-may-07 was 0// 24-nov-04 - was 1.e-6;
210 AliDebug(2,Form("Print: \n------------------- %s -------------\n",GetName()));
211 AliDebug(2,Form(" fInit %i\n", int(fInit)));
212 AliDebug(2,Form(" fFirstEvent %i\n", fFirstEvent));
213 AliDebug(2,Form(" fLastEvent %i\n", fLastEvent));
214 AliDebug(2,Form(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()));
215 AliDebug(2,Form(" with digitization parameters A = %f\n", fA));
216 AliDebug(2,Form(" B = %f\n", fB));
217 AliDebug(2,Form(" Threshold for EC Primary assignment = %f\n", fECPrimThreshold));
218 AliDebug(2,Form(" Sampling = %f\n", fSampling));
219 AliDebug(2,Form("---------------------------------------------------\n"));
223 //____________________________________________________________________________
224 void AliEMCALSDigitizer::Digitize(Option_t *option)
226 // Collects all hit of the same tower into digits
227 TString o(option); o.ToUpper();
228 if (strstr(option, "print") ) {
230 AliDebug(2,Form("Print: \n------------------- %s -------------\n",GetName()));
231 AliDebug(2,Form(" fInit %i\n", int(fInit)));
232 AliDebug(2,Form(" fFirstEvent %i\n", fFirstEvent));
233 AliDebug(2,Form(" fLastEvent %i\n", fLastEvent));
234 AliDebug(2,Form(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()));
235 AliDebug(2,Form(" with digitization parameters A = %f\n", fA));
236 AliDebug(2,Form(" B = %f\n", fB));
237 AliDebug(2,Form(" Threshold for EC Primary assignment = %f\n", fECPrimThreshold));
238 AliDebug(2,Form(" Sampling = %f\n", fSampling));
239 AliDebug(2,Form("---------------------------------------------------\n"));
245 if(strstr(option,"tim"))
246 gBenchmark->Start("EMCALSDigitizer");
248 AliRunLoader *rl = AliRunLoader::Instance();
249 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
251 if (!fInit) { // to prevent overwrite existing file
252 AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
256 if (fLastEvent == -1)
257 fLastEvent = rl->GetNumberOfEvents() - 1 ;
259 fLastEvent = TMath::Min(fLastEvent, rl->GetNumberOfEvents()-1);
261 Int_t nEvents = fLastEvent - fFirstEvent + 1;
263 Float_t energy=0.; // de * fSampling - 23-nov-04
264 rl->LoadKinematics();
265 rl->LoadHits("EMCAL");
268 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
269 rl->GetEvent(ievent);
270 TTree * treeS = emcalLoader->TreeS();
272 emcalLoader->MakeSDigitsContainer();
273 treeS = emcalLoader->TreeS();
276 TClonesArray * sdigits = emcalLoader->SDigits() ;
280 Int_t iHit, iTrack, iSDigit;
282 AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
284 TTree *treeH = emcalLoader->TreeH();
286 Int_t nTrack = treeH->GetEntries(); // TreeH has array of hits for every primary
287 TBranch * branchH = treeH->GetBranch("EMCAL");
288 //if(fHits)fHits->Clear();
289 branchH->SetAddress(&fHits);
290 for (iTrack = 0; iTrack < nTrack; iTrack++) {
291 branchH->GetEntry(iTrack);
295 Int_t nHit = fHits->GetEntriesFast();
296 for(iHit = 0; iHit< nHit;iHit++){
298 AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(fHits->At(iHit)) ;
299 AliEMCALDigit * curSDigit = 0 ;
300 AliEMCALDigit * sdigit = 0 ;
301 Bool_t newsdigit = kTRUE;
303 // hit->GetId() - Absolute Id number EMCAL segment
305 if(geom->CheckAbsCellId(hit->GetId())) { // was IsInECA(hit->GetId())
306 energy = hit->GetEnergy() * fSampling; // 23-nov-04
307 if(energy > fECPrimThreshold )
308 // Assign primary number only if deposited energy is significant
309 curSDigit = new AliEMCALDigit(hit->GetPrimary(),
310 hit->GetIparent(), hit->GetId(),
311 Digitize(energy), hit->GetTime(),kFALSE,
314 curSDigit = new AliEMCALDigit(-1,
317 Digitize(energy), hit->GetTime(),kFALSE,
320 Warning("Digitize"," abs id %i is bad \n", hit->GetId());
326 for(Int_t check= 0; check < nSdigits ; check++) {
327 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
329 if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL tower ?
330 *sdigit = *sdigit + *curSDigit;
335 AliWarning("Sdigit do not exist");
337 }// sdigit does not exist
339 }// currsdigit exists
342 new((*sdigits)[nSdigits]) AliEMCALDigit(*curSDigit);
348 else AliFatal("Hit is NULL!");
350 } // loop over all hits (hit = deposited energy/entering particle)
353 else AliFatal("fHit is NULL!");
357 nSdigits = sdigits->GetEntriesFast() ;
358 fSDigitsInRun += nSdigits ;
360 for (iSDigit = 0 ; iSDigit < sdigits->GetEntriesFast() ; iSDigit++) {
361 AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(iSDigit)) ;
362 if(sdigit)sdigit->SetIndexInList(iSDigit) ;
363 else AliFatal("sdigit is NULL!");
365 if(fHits)fHits->Clear();
371 Int_t bufferSize = 32000 ;
372 TBranch * sdigitsBranch = treeS->GetBranch("EMCAL");
374 sdigitsBranch->SetAddress(&sdigits);
376 treeS->Branch("EMCAL",&sdigits,bufferSize);
380 emcalLoader->WriteSDigits("OVERWRITE");
383 //emcalLoader->WriteSDigitizer("OVERWRITE"); // why in event cycle ?
385 if(strstr(option,"deb"))
386 PrintSDigits(option) ;
391 if(strstr(option,"tim")){
392 gBenchmark->Stop("EMCALSDigitizer");
393 printf("\n Digitize: took %f seconds for SDigitizing %f seconds per event\n",
394 gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer")/nEvents ) ;
398 //__________________________________________________________________
399 Float_t AliEMCALSDigitizer::Digitize(Float_t energy)const {
400 // Digitize the energy
402 //JLK 26-Jun-2008 EXPLANATION LONG OVERDUE:
404 //We have to digitize the SDigit energy so that it can be stored in
405 //AliEMCALDigit, which has only an ADC amplitude member and
406 //(rightly) no energy member. This method converts the energy to an
407 //integer which can be re-converted back to an energy with the
408 //Calibrate(energy) method when it is time to create Digits from SDigits
410 Double_t aSignal = fA + energy*fB;
411 if (TMath::Abs(aSignal)>2147483647.0) {
412 //PH 2147483647 is the max. integer
413 //PH This apparently is a problem which needs investigation
414 AliWarning(Form("Too big or too small energy %f",aSignal));
415 aSignal = TMath::Sign((Double_t)2147483647,aSignal);
418 return (Float_t ) aSignal;
421 //__________________________________________________________________
422 Float_t AliEMCALSDigitizer::Calibrate(Float_t amp)const {
424 // Convert the amplitude back to energy in GeV
426 //JLK 26-Jun-2008 EXPLANATION LONG OVERDUE:
428 //We have to digitize the SDigit energy with the method Digitize()
429 //so that it can be stored in AliEMCALDigit, which has only an ADC
430 //amplitude member and (rightly) no energy member. This method is
431 //just the reverse of Digitize(): it converts the stored amplitude
432 //back to an energy value in GeV so that the SDigit energies can be
433 //summed before adding noise and creating digits out of them
435 return (Float_t)(amp - fA)/fB;
440 //__________________________________________________________________
441 void AliEMCALSDigitizer::Print1(Option_t * option)
444 PrintSDigits(option);
447 //__________________________________________________________________
448 void AliEMCALSDigitizer::Print(Option_t *option) const
450 // Prints parameters of SDigitizer
451 printf("Print: \n------------------- %s ------------- option %s\n", GetName() , option) ;
452 printf(" fInit %i\n", int(fInit));
453 printf(" fFirstEvent %i\n", fFirstEvent);
454 printf(" fLastEvent %i\n", fLastEvent);
455 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
456 printf(" with digitization parameters A = %f\n", fA) ;
457 printf(" B = %f\n", fB) ;
458 printf(" Threshold for EC Primary assignment = %f\n", fECPrimThreshold) ;
459 printf(" Sampling = %f\n", fSampling);
460 printf("---------------------------------------------------\n") ;
463 //__________________________________________________________________
464 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
467 // SDititizers are equal if their pedestal, slope and threshold are equal
468 if( (fA==sd.fA)&&(fB==sd.fB)&&
469 (fECPrimThreshold==sd.fECPrimThreshold))
474 //__________________________________________________________________
475 void AliEMCALSDigitizer::PrintSDigits(Option_t * option)
477 //Prints list of digits produced at the current pass of AliEMCALDigitizer
479 AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
481 const TClonesArray * sdigits = rl->SDigits() ;
484 printf("event %i", rl->GetRunLoader()->GetEventNumber());
485 printf(" Number of entries in SDigits list %i", sdigits->GetEntriesFast());
486 if(strstr(option,"all")||strstr(option,"EMC")){
489 AliEMCALDigit * digit;
490 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
493 const Int_t bufferSize = 8192;
494 char * tempo = new char[bufferSize];
495 for (index = 0 ; index < sdigits->GetEntries() ; index++) {
496 digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
498 snprintf(tempo, bufferSize,"\n%6d %8f %6.5e %4d %2d :",
499 digit->GetId(), digit->GetAmplitude(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
501 isum += digit->GetAmplitude();
504 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
505 snprintf(tempo,bufferSize, "%d ",digit->GetPrimary(iprimary+1) ) ;
509 else AliFatal("SDigit is NULL!");
512 printf("\n** Sum %2.3f : %10.3f GeV/c **\n ", isum, Calibrate(isum));
515 else AliFatal("EMCALLoader is NULL!");
518 //____________________________________________________________________________
519 void AliEMCALSDigitizer::Unload() const
521 // Unload Hits and SDigits from the folder
522 AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
525 rl->UnloadSDigits() ;
527 else AliFatal("EMCALLoader is NULL!");
530 //____________________________________________________________________________
531 void AliEMCALSDigitizer::Browse(TBrowser* b)