]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSDigitizer.cxx
from M Weber: L1Gamma sigma check value relaxed
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALSDigitizer.cxx
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$ */
17
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 
23 //
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
29 //
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.
36 // User case:
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
50 //
51 //*-- Author : Sahal Yacoob (LBL)
52 // based on  : AliPHOSSDigitzer 
53 //////////////////////////////////////////////////////////////////////////////
54
55 // --- ROOT system ---
56 #include <TBenchmark.h>
57 #include <TBrowser.h>
58 //#include <Riostream.h>
59 #include <TMath.h>
60 #include <TROOT.h>
61
62 // --- Standard library ---
63 #include "stdlib.h"
64
65 // --- AliRoot header files ---
66 #include "AliLog.h"
67 #include "AliRunLoader.h"
68 #include "AliStack.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"
75
76 ClassImp(AliEMCALSDigitizer)
77            
78 //____________________________________________________________________________ 
79 AliEMCALSDigitizer::AliEMCALSDigitizer()
80   : TNamed("",""),
81     fA(0.),fB(0.),fECPrimThreshold(0.),
82     fDefaultInit(kTRUE),
83     fEventFolderName(0),
84     fInit(0),
85     fSDigitsInRun(0),
86     fFirstEvent(0),
87     fLastEvent(0),
88     fSampling(0.),
89         fHits(0)
90 {
91   // ctor
92   InitParameters();
93 }
94
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),
102     fInit(0),
103     fSDigitsInRun(0),
104     fFirstEvent(0),
105     fLastEvent(0),
106     fSampling(0.),
107     fHits(0)
108 {
109   // ctor
110   Init();
111   InitParameters() ; 
112
113 }
114
115
116 //____________________________________________________________________________ 
117 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd) 
118   : TNamed(sd.GetName(),sd.GetTitle()),
119     fA(sd.fA),
120     fB(sd.fB),
121     fECPrimThreshold(sd.fECPrimThreshold),
122     fDefaultInit(sd.fDefaultInit),
123     fEventFolderName(sd.fEventFolderName),
124     fInit(sd.fInit),
125     fSDigitsInRun(sd.fSDigitsInRun),
126     fFirstEvent(sd.fFirstEvent),
127     fLastEvent(sd.fLastEvent),
128     fSampling(sd.fSampling),
129     fHits(0)
130 {
131   //cpy ctor 
132 }
133
134 //_____________________________________________________________________
135 AliEMCALSDigitizer& AliEMCALSDigitizer::operator = (const AliEMCALSDigitizer& source)
136 { // assignment operator; use copy ctor
137   if (&source == this) return *this;
138
139   new (this) AliEMCALSDigitizer(source);
140   return *this;
141 }
142
143 //____________________________________________________________________________ 
144 AliEMCALSDigitizer::~AliEMCALSDigitizer() {
145   //dtor
146         
147   if(fHits){
148           fHits->Clear();
149           delete fHits;
150   }
151 }
152
153 //____________________________________________________________________________ 
154 void AliEMCALSDigitizer::Init(){
155   // Initialization: open root-file, allocate arrays for hits and sdigits
156   // 
157   // Initialization can not be done in the default constructor
158   //============================================================= YS
159   //  The initialisation is now done by the getter
160
161   fInit = kTRUE; 
162    
163   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
164
165   if ( emcalLoader == 0 ) {
166     Fatal("Init", "Could not obtain the AliEMCALLoader");
167     return ;
168   } 
169   
170 }
171
172 //____________________________________________________________________________ 
173 void AliEMCALSDigitizer::InitParameters()
174 {
175   //initialize parameters for sdigitization
176
177   const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
178   if (geom->GetSampling() == 0.) {
179     Fatal("InitParameters", "Sampling factor not set !") ; 
180   }
181
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();
187         
188   if(!simParam){
189         simParam = AliEMCALSimParam::GetInstance();
190         AliWarning("Simulation Parameters not available in OCDB?");
191   }
192         
193   //
194   //JLK 26-Jun-2008 THIS SHOULD HAVE BEEN EXPLAINED AGES AGO:
195   //
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
201   //
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();
206
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;
209   
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"));
220
221 }
222
223 //____________________________________________________________________________
224 void AliEMCALSDigitizer::Digitize(Option_t *option) 
225
226         // Collects all hit of the same tower into digits
227         TString o(option); o.ToUpper();
228         if (strstr(option, "print") ) {
229                 
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"));
240                 
241                 return ; 
242         }
243         
244         
245         if(strstr(option,"tim"))
246                 gBenchmark->Start("EMCALSDigitizer");
247         
248         AliRunLoader *rl = AliRunLoader::Instance();
249         AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
250         
251         if (!fInit) { // to prevent overwrite existing file
252           AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
253           return ;
254         }
255         
256         if (fLastEvent == -1) 
257                 fLastEvent = rl->GetNumberOfEvents() - 1 ;
258         else {
259                 fLastEvent = TMath::Min(fLastEvent, rl->GetNumberOfEvents()-1);
260         }
261         Int_t nEvents   = fLastEvent - fFirstEvent + 1;
262         
263         Float_t energy=0.; // de * fSampling - 23-nov-04
264         rl->LoadKinematics();
265         rl->LoadHits("EMCAL");
266   
267   Int_t ievent;
268         for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
269                 rl->GetEvent(ievent);
270                 TTree * treeS = emcalLoader->TreeS();
271                 if ( !treeS ) { 
272                         emcalLoader->MakeSDigitsContainer();
273                         treeS = emcalLoader->TreeS();
274                 }
275     
276                 TClonesArray * sdigits = emcalLoader->SDigits() ;
277                 sdigits->Clear("C");
278                 
279                 Int_t nSdigits = 0 ;
280                 Int_t iHit, iTrack, iSDigit;
281                 
282                 AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance(); 
283                 
284                 TTree *treeH = emcalLoader->TreeH();    
285                 if (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); 
292         
293         if(fHits){
294           
295           Int_t nHit = fHits->GetEntriesFast();
296           for(iHit = 0; iHit< nHit;iHit++){
297             
298             AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(fHits->At(iHit)) ;
299             AliEMCALDigit * curSDigit = 0 ;
300             AliEMCALDigit * sdigit = 0 ;
301             Bool_t newsdigit = kTRUE; 
302             
303             // hit->GetId() - Absolute Id number EMCAL segment
304             if(hit){
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,
312                                                  -1, 0,0,energy ) ;
313                 else
314                   curSDigit =  new AliEMCALDigit(-1, 
315                                                  -1,
316                                                  hit->GetId(), 
317                                                  Digitize(energy), hit->GetTime(),kFALSE,
318                                                  -1, 0,0,energy ) ;
319               } else {
320                 Warning("Digitize"," abs id %i is bad \n", hit->GetId());
321                 newsdigit = kFALSE;
322                 curSDigit = 0;
323               }
324               
325               if(curSDigit != 0){
326                 for(Int_t check= 0; check < nSdigits ; check++) {
327                   sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
328                   if(sdigit){
329                     if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL tower ?              
330                       *sdigit = *sdigit + *curSDigit;
331                       newsdigit = kFALSE;
332                     }
333                   }// sdigit exists
334                   else {
335                     AliWarning("Sdigit do not exist");
336                     newsdigit = kFALSE;  
337                   }// sdigit does not exist
338                 }//sdigit loop
339               }// currsdigit exists
340               
341               if (newsdigit) {
342                 new((*sdigits)[nSdigits])  AliEMCALDigit(*curSDigit);
343                 nSdigits++ ;  
344               }
345               delete curSDigit ;
346               
347             }// hit exists
348             else AliFatal("Hit is NULL!");
349             
350           }  // loop over all hits (hit = deposited energy/entering particle)
351           
352         }//fHits is not NULL
353         else AliFatal("fHit is NULL!");
354         
355                                 sdigits->Sort() ;
356                                 
357                                 nSdigits = sdigits->GetEntriesFast() ;
358                                 fSDigitsInRun += nSdigits ;  
359                                 
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!");
364                                 }       
365                                 if(fHits)fHits->Clear();
366                         }//track loop
367                 }// tree exists
368                 
369                 // Now write SDigits    
370                 
371                 Int_t bufferSize = 32000 ;    
372                 TBranch * sdigitsBranch = treeS->GetBranch("EMCAL");
373                 if (sdigitsBranch)
374                         sdigitsBranch->SetAddress(&sdigits);
375                 else
376                         treeS->Branch("EMCAL",&sdigits,bufferSize);
377                 
378                 treeS->Fill();
379                 
380                 emcalLoader->WriteSDigits("OVERWRITE");
381                 
382                 //NEXT - SDigitizer
383                 //emcalLoader->WriteSDigitizer("OVERWRITE");  // why in event cycle ?
384                 
385                 if(strstr(option,"deb"))
386                         PrintSDigits(option) ;  
387         }
388         
389         Unload();
390                 
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 ) ; 
395         }
396 }
397
398 //__________________________________________________________________
399 Float_t AliEMCALSDigitizer::Digitize(Float_t energy)const {
400   // Digitize the energy
401   //
402   //JLK 26-Jun-2008 EXPLANATION LONG OVERDUE:
403   //
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 
409   //
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);
416   }
417
418   return (Float_t ) aSignal;
419 }
420
421 //__________________________________________________________________
422 Float_t AliEMCALSDigitizer::Calibrate(Float_t amp)const {
423   //
424   // Convert the amplitude back to energy in GeV
425   //
426   //JLK 26-Jun-2008 EXPLANATION LONG OVERDUE:
427   //
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
434   //
435   return (Float_t)(amp - fA)/fB;
436
437 }
438  
439
440 //__________________________________________________________________
441 void AliEMCALSDigitizer::Print1(Option_t * option)
442 {
443   Print(); 
444   PrintSDigits(option);
445 }
446
447 //__________________________________________________________________
448 void AliEMCALSDigitizer::Print(Option_t *option) const
449
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") ;
461 }
462
463 //__________________________________________________________________
464 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
465 {
466   // Equal operator.
467   // SDititizers are equal if their pedestal, slope and threshold are equal
468   if( (fA==sd.fA)&&(fB==sd.fB)&&
469       (fECPrimThreshold==sd.fECPrimThreshold))
470     return kTRUE ;
471   else
472     return kFALSE ;
473 }
474 //__________________________________________________________________
475 void AliEMCALSDigitizer::PrintSDigits(Option_t * option)
476 {
477   //Prints list of digits produced at the current pass of AliEMCALDigitizer
478   
479   AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
480   if(rl){
481     const TClonesArray * sdigits = rl->SDigits() ; 
482     
483     printf("\n") ;  
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")){
487       
488       //loop over digits
489       AliEMCALDigit * digit;
490       printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
491       Int_t   index = 0;
492       Float_t isum  = 0.;
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) ) ;
497         if(digit){
498           snprintf(tempo, bufferSize,"\n%6d  %8f    %6.5e %4d      %2d :",
499                   digit->GetId(), digit->GetAmplitude(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
500           printf("%s",tempo);
501           isum += digit->GetAmplitude();
502           
503           Int_t iprimary;
504           for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
505             snprintf(tempo,bufferSize, "%d ",digit->GetPrimary(iprimary+1) ) ; 
506             printf("%s",tempo); 
507           }  
508         } //sdigit exists
509         else AliFatal("SDigit is NULL!");
510       }//loop
511       delete [] tempo ;
512       printf("\n** Sum %2.3f : %10.3f GeV/c **\n ", isum, Calibrate(isum));
513     } else printf("\n");
514   }
515   else AliFatal("EMCALLoader is NULL!");
516 }
517
518 //____________________________________________________________________________ 
519 void AliEMCALSDigitizer::Unload() const
520 {
521   // Unload Hits and SDigits from the folder
522   AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
523   if(rl){
524     rl->UnloadHits() ; 
525     rl->UnloadSDigits() ;
526   }
527   else AliFatal("EMCALLoader is NULL!");
528 }
529
530 //____________________________________________________________________________ 
531 void AliEMCALSDigitizer::Browse(TBrowser* b)
532 {
533   TNamed::Browse(b);
534 }