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 TTask 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->ExecuteTask()
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->ExeciteTask("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.),
94 //____________________________________________________________________________
95 AliEMCALSDigitizer::AliEMCALSDigitizer(const char * alirunFileName,
96 const char * eventFolderName)
97 : TTask("EMCAL"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
98 fA(0.),fB(0.),fECPrimThreshold(0.),
100 fEventFolderName(eventFolderName),
114 //____________________________________________________________________________
115 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd)
116 : TTask(sd.GetName(),sd.GetTitle()),
119 fECPrimThreshold(sd.fECPrimThreshold),
120 fDefaultInit(sd.fDefaultInit),
121 fEventFolderName(sd.fEventFolderName),
123 fSDigitsInRun(sd.fSDigitsInRun),
124 fFirstEvent(sd.fFirstEvent),
125 fLastEvent(sd.fLastEvent),
126 fSampling(sd.fSampling)
132 //____________________________________________________________________________
133 AliEMCALSDigitizer::~AliEMCALSDigitizer() {
135 AliLoader *emcalLoader=0;
136 if ((emcalLoader = AliRunLoader::Instance()->GetDetectorLoader("EMCAL")))
137 emcalLoader->CleanSDigitizer();
140 //____________________________________________________________________________
141 void AliEMCALSDigitizer::Init(){
142 // Initialization: open root-file, allocate arrays for hits and sdigits,
143 // attach task SDigitizer to the list of EMCAL tasks
145 // Initialization can not be done in the default constructor
146 //============================================================= YS
147 // The initialisation is now done by the getter
151 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
153 if ( emcalLoader == 0 ) {
154 Fatal("Init", "Could not obtain the AliEMCALLoader");
158 emcalLoader->PostSDigitizer(this);
159 emcalLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
163 //____________________________________________________________________________
164 void AliEMCALSDigitizer::InitParameters()
166 //initialize parameters for sdigitization
168 const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
169 if (geom->GetSampling() == 0.) {
170 Fatal("InitParameters", "Sampling factor not set !") ;
173 // Get the parameters from the OCDB via the loader
174 AliRunLoader *rl = AliRunLoader::Instance();
175 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
176 AliEMCALSimParam * simParam = 0x0;
177 if(emcalLoader) simParam = emcalLoader->SimulationParameters();
180 simParam = AliEMCALSimParam::GetInstance();
181 AliWarning("Simulation Parameters not available in OCDB?");
185 //JLK 26-Jun-2008 THIS SHOULD HAVE BEEN EXPLAINED AGES AGO:
187 //In order to be able to store SDigit Energy info into
188 //AliEMCALDigit, we need to convert it temporarily to an ADC amplitude
189 //and later when summing SDigits to form digits, convert it back to
190 //energy. These fA and fB parameters accomplish this through the
191 //Digitize() and Calibrate() methods
193 // Initializes parameters
194 fA = simParam->GetA(); //0;
195 fB = simParam->GetB(); //1.e+6; // Changed 24 Apr 2007. Dynamic range now 2 TeV
196 fSampling = geom->GetSampling();
198 // threshold for deposit energy of hit
199 fECPrimThreshold = simParam->GetECPrimaryThreshold();//0.05;// GeV // 22-may-07 was 0// 24-nov-04 - was 1.e-6;
201 AliDebug(2,Form("Print: \n------------------- %s -------------\n",GetName()));
202 AliDebug(2,Form(" fInit %i\n", int(fInit)));
203 AliDebug(2,Form(" fFirstEvent %i\n", fFirstEvent));
204 AliDebug(2,Form(" fLastEvent %i\n", fLastEvent));
205 AliDebug(2,Form(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()));
206 AliDebug(2,Form(" with digitization parameters A = %f\n", fA));
207 AliDebug(2,Form(" B = %f\n", fB));
208 AliDebug(2,Form(" Threshold for EC Primary assignment = %f\n", fECPrimThreshold));
209 AliDebug(2,Form(" Sampling = %f\n", fSampling));
210 AliDebug(2,Form("---------------------------------------------------\n"));
214 //____________________________________________________________________________
215 void AliEMCALSDigitizer::Exec(Option_t *option)
217 // Collects all hit of the same tower into digits
218 TString o(option); o.ToUpper();
219 if (strstr(option, "print") ) {
221 AliDebug(2,Form("Print: \n------------------- %s -------------\n",GetName()));
222 AliDebug(2,Form(" fInit %i\n", int(fInit)));
223 AliDebug(2,Form(" fFirstEvent %i\n", fFirstEvent));
224 AliDebug(2,Form(" fLastEvent %i\n", fLastEvent));
225 AliDebug(2,Form(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()));
226 AliDebug(2,Form(" with digitization parameters A = %f\n", fA));
227 AliDebug(2,Form(" B = %f\n", fB));
228 AliDebug(2,Form(" Threshold for EC Primary assignment = %f\n", fECPrimThreshold));
229 AliDebug(2,Form(" Sampling = %f\n", fSampling));
230 AliDebug(2,Form("---------------------------------------------------\n"));
236 if(strstr(option,"tim"))
237 gBenchmark->Start("EMCALSDigitizer");
239 AliRunLoader *rl = AliRunLoader::Instance();
240 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
242 //switch off reloading of this task while getting event
243 if (!fInit) { // to prevent overwrite existing file
244 AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
248 if (fLastEvent == -1)
249 fLastEvent = rl->GetNumberOfEvents() - 1 ;
251 fLastEvent = TMath::Min(fLastEvent, rl->GetNumberOfEvents()-1);
253 Int_t nEvents = fLastEvent - fFirstEvent + 1;
256 Float_t energy=0.; // de * fSampling - 23-nov-04
257 rl->LoadKinematics();
258 rl->LoadHits("EMCAL");
260 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
261 rl->GetEvent(ievent);
262 TTree * treeS = emcalLoader->TreeS();
264 emcalLoader->MakeSDigitsContainer();
265 treeS = emcalLoader->TreeS();
267 TClonesArray * hits = emcalLoader->Hits() ;
268 TClonesArray * sdigits = emcalLoader->SDigits() ;
274 AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
275 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
276 AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
277 AliEMCALDigit * curSDigit = 0 ;
278 AliEMCALDigit * sdigit = 0 ;
279 Bool_t newsdigit = kTRUE;
281 // hit->GetId() - Absolute Id number EMCAL segment
282 if(geom->CheckAbsCellId(hit->GetId())) { // was IsInECA(hit->GetId())
283 energy = hit->GetEnergy() * fSampling; // 23-nov-04
284 if(energy > fECPrimThreshold )
285 // Assign primary number only if deposited energy is significant
286 curSDigit = new AliEMCALDigit( hit->GetPrimary(),
287 hit->GetIparent(), hit->GetId(),
288 Digitize(energy), hit->GetTime(),kFALSE,
291 curSDigit = new AliEMCALDigit( -1,
294 Digitize(energy), hit->GetTime(),kFALSE,
297 Warning("Exec"," abs id %i is bad \n", hit->GetId());
303 for(Int_t check= 0; check < nSdigits ; check++) {
304 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
306 if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL tower ?
307 *sdigit = *sdigit + *curSDigit;
313 new((*sdigits)[nSdigits]) AliEMCALDigit(*curSDigit);
317 } // loop over all hits (hit = deposited energy/entering particle)
320 nSdigits = sdigits->GetEntriesFast() ;
321 fSDigitsInRun += nSdigits ;
323 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
324 AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
325 sdigit->SetIndexInList(i) ;
330 Int_t bufferSize = 32000 ;
331 TBranch * sdigitsBranch = treeS->GetBranch("EMCAL");
333 sdigitsBranch->SetAddress(&sdigits);
335 treeS->Branch("EMCAL",&sdigits,bufferSize);
339 emcalLoader->WriteSDigits("OVERWRITE");
342 emcalLoader->WriteSDigitizer("OVERWRITE"); // why in event cycle ?
344 if(strstr(option,"deb"))
345 PrintSDigits(option) ;
350 emcalLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
352 if(strstr(option,"tim")){
353 gBenchmark->Stop("EMCALSDigitizer");
354 printf("\n Exec: took %f seconds for SDigitizing %f seconds per event\n",
355 gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer")/nEvents ) ;
359 //__________________________________________________________________
360 Float_t AliEMCALSDigitizer::Digitize(Float_t energy)const {
361 // Digitize the energy
363 //JLK 26-Jun-2008 EXPLANATION LONG OVERDUE:
365 //We have to digitize the SDigit energy so that it can be stored in
366 //AliEMCALDigit, which has only an ADC amplitude member and
367 //(rightly) no energy member. This method converts the energy to an
368 //integer which can be re-converted back to an energy with the
369 //Calibrate(energy) method when it is time to create Digits from SDigits
371 Double_t aSignal = fA + energy*fB;
372 if (TMath::Abs(aSignal)>2147483647.0) {
373 //PH 2147483647 is the max. integer
374 //PH This apparently is a problem which needs investigation
375 AliWarning(Form("Too big or too small energy %f",aSignal));
376 aSignal = TMath::Sign((Double_t)2147483647,aSignal);
379 return (Float_t ) aSignal;
382 //__________________________________________________________________
383 Float_t AliEMCALSDigitizer::Calibrate(Float_t amp)const {
385 // Convert the amplitude back to energy in GeV
387 //JLK 26-Jun-2008 EXPLANATION LONG OVERDUE:
389 //We have to digitize the SDigit energy with the method Digitize()
390 //so that it can be stored in AliEMCALDigit, which has only an ADC
391 //amplitude member and (rightly) no energy member. This method is
392 //just the reverse of Digitize(): it converts the stored amplitude
393 //back to an energy value in GeV so that the SDigit energies can be
394 //summed before adding noise and creating digits out of them
396 return (Float_t)(amp - fA)/fB;
401 //__________________________________________________________________
402 void AliEMCALSDigitizer::Print1(Option_t * option)
405 PrintSDigits(option);
408 //__________________________________________________________________
409 void AliEMCALSDigitizer::Print(Option_t *option) const
411 // Prints parameters of SDigitizer
412 printf("Print: \n------------------- %s ------------- option %s\n", GetName() , option) ;
413 printf(" fInit %i\n", int(fInit));
414 printf(" fFirstEvent %i\n", fFirstEvent);
415 printf(" fLastEvent %i\n", fLastEvent);
416 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
417 printf(" with digitization parameters A = %f\n", fA) ;
418 printf(" B = %f\n", fB) ;
419 printf(" Threshold for EC Primary assignment = %f\n", fECPrimThreshold) ;
420 printf(" Sampling = %f\n", fSampling);
421 printf("---------------------------------------------------\n") ;
424 //__________________________________________________________________
425 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
428 // SDititizers are equal if their pedestal, slope and threshold are equal
429 if( (fA==sd.fA)&&(fB==sd.fB)&&
430 (fECPrimThreshold==sd.fECPrimThreshold))
436 //__________________________________________________________________
437 void AliEMCALSDigitizer::PrintSDigits(Option_t * option)
439 //Prints list of digits produced at the current pass of AliEMCALDigitizer
441 AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
442 const TClonesArray * sdigits = rl->SDigits() ;
445 printf("event %i", rl->GetRunLoader()->GetEventNumber());
446 printf(" Number of entries in SDigits list %i", sdigits->GetEntriesFast());
447 if(strstr(option,"all")||strstr(option,"EMC")){
450 AliEMCALDigit * digit;
451 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
454 char * tempo = new char[8192];
455 for (index = 0 ; index < sdigits->GetEntries() ; index++) {
456 digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
457 sprintf(tempo, "\n%6d %8f %6.5e %4d %2d :",
458 digit->GetId(), digit->GetAmplitude(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
460 isum += digit->GetAmplitude();
463 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
464 sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ;
469 printf("\n** Sum %2.3f : %10.3f GeV/c **\n ", isum, Calibrate(isum));
473 //____________________________________________________________________________
474 void AliEMCALSDigitizer::Unload() const
476 // Unload Hits and SDigits from the folder
477 AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
479 rl->UnloadSDigits() ;
482 //____________________________________________________________________________
483 void AliEMCALSDigitizer::Browse(TBrowser* b)