]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSDigitizer.cxx
Adding includes now needed by ROOT
[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 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 
23 // SDigits are written to TreeS, branch "EMCAL"
24 // AliEMCALSDigitizer with all current parameters is written 
25 // to TreeS branch "AliEMCALSDigitizer".
26 // Both branches have the same title. If necessary one can produce 
27 // another set of SDigits with different parameters. Two versions
28 // can be distunguished using titles of the branches.
29 // User case:
30 // root [0] AliEMCALSDigitizer * s = new AliEMCALSDigitizer("galice.root")
31 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
32 // root [1] s->ExecuteTask()
33 //             // Makes SDigitis for all events stored in galice.root
34 // root [2] s->SetPedestalParameter(0.001)
35 //             // One can change parameters of digitization
36 // root [3] s->SetSDigitsBranch("Redestal 0.001")
37 //             // and write them into the new branch
38 // root [4] s->ExeciteTask("deb all tim")
39 //             // available parameters:
40 //             deb - print # of produced SDigitis
41 //             deb all  - print # and list of produced SDigits
42 //             tim - print benchmarking information
43 //
44 //*-- Author : Sahal Yacoob (LBL)
45 // based on  : AliPHOSSDigitzer 
46 // Modif: 
47 //  August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
48 //                           of new  IO (à la PHOS)
49 //////////////////////////////////////////////////////////////////////////////
50
51 // --- ROOT system ---
52 #include <TBenchmark.h>
53 #include <TH1.h>
54 #include <TBrowser.h>
55 #include <Riostream.h>
56 #include <TMath.h>
57 #include <TROOT.h>
58
59 // --- Standard library ---
60 #include "stdlib.h"
61
62 // --- AliRoot header files ---
63 #include "AliLog.h"
64 #include "AliRunLoader.h"
65 #include "AliStack.h"
66 #include "AliEMCALDigit.h"
67 #include "AliEMCALLoader.h"
68 #include "AliEMCALHit.h"
69 #include "AliEMCALSDigitizer.h"
70 #include "AliEMCALGeometry.h"
71 #include "AliEMCALHistoUtilities.h"
72
73 ClassImp(AliEMCALSDigitizer)
74            
75 //____________________________________________________________________________ 
76 AliEMCALSDigitizer::AliEMCALSDigitizer()
77   : TTask("",""),
78     fA(0.),fB(0.),fECPrimThreshold(0.),
79     fDefaultInit(kTRUE),
80     fEventFolderName(0),
81     fInit(0),
82     fSDigitsInRun(0),
83     fFirstEvent(0),
84     fLastEvent(0),
85     fSampling(0.),
86     fControlHists(0),
87     fHists(0)
88 {
89   // ctor
90   InitParameters();
91 }
92
93 //____________________________________________________________________________ 
94 AliEMCALSDigitizer::AliEMCALSDigitizer(const char * alirunFileName, 
95                                        const char * eventFolderName)
96   : TTask("EMCAL"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
97     fA(0.),fB(0.),fECPrimThreshold(0.),
98     fDefaultInit(kFALSE),
99     fEventFolderName(eventFolderName),
100     fInit(0),
101     fSDigitsInRun(0),
102     fFirstEvent(0),
103     fLastEvent(0),
104     fSampling(0.),
105     fControlHists(1),
106     fHists(0)
107 {
108   // ctor
109   Init();
110   InitParameters() ; 
111   if(fControlHists) BookControlHists(1);
112 }
113
114
115 //____________________________________________________________________________ 
116 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd) 
117   : TTask(sd.GetName(),sd.GetTitle()),
118     fA(sd.fA),
119     fB(sd.fB),
120     fECPrimThreshold(sd.fECPrimThreshold),
121     fDefaultInit(sd.fDefaultInit),
122     fEventFolderName(sd.fEventFolderName),
123     fInit(sd.fInit),
124     fSDigitsInRun(sd.fSDigitsInRun),
125     fFirstEvent(sd.fFirstEvent),
126     fLastEvent(sd.fLastEvent),
127     fSampling(sd.fSampling),
128     fControlHists(sd.fControlHists),
129     fHists(sd.fHists)
130 {
131   //cpy ctor 
132 }
133
134
135 //____________________________________________________________________________ 
136 AliEMCALSDigitizer::~AliEMCALSDigitizer() {
137   //dtor
138   AliLoader *emcalLoader=0;
139   if ((emcalLoader = AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL")))
140       emcalLoader->CleanSDigitizer();
141 }
142
143 //____________________________________________________________________________ 
144 void AliEMCALSDigitizer::Init(){
145   // Initialization: open root-file, allocate arrays for hits and sdigits,
146   // attach task SDigitizer to the list of EMCAL tasks
147   // 
148   // Initialization can not be done in the default constructor
149   //============================================================= YS
150   //  The initialisation is now done by the getter
151
152   fInit = kTRUE; 
153    
154   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
155
156   if ( emcalLoader == 0 ) {
157     Fatal("Init", "Could not obtain the AliEMCALLoader");
158     return ;
159   } 
160   
161   emcalLoader->PostSDigitizer(this);
162   emcalLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
163   
164 }
165
166 //____________________________________________________________________________ 
167 void AliEMCALSDigitizer::InitParameters()
168 {
169   //initialize parameters for sdigitization
170
171   const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
172   if (geom->GetSampling() == 0.) {
173     Fatal("InitParameters", "Sampling factor not set !") ; 
174   }
175   // Initializes parameters
176   fA         = 0;
177   fB         = 1.e+7;
178   fSampling  = geom->GetSampling();
179
180  // threshold for deposit energy of hit
181   fECPrimThreshold  = 0.; // 24-nov-04 - was 1.e-6;
182
183   AliDebug(2,Form("Print: \n------------------- %s -------------\n",GetName()));
184   AliDebug(2,Form("   fInit                                 %i\n", int(fInit)));
185   AliDebug(2,Form("   fFirstEvent                           %i\n", fFirstEvent));
186   AliDebug(2,Form("   fLastEvent                            %i\n", fLastEvent));
187   AliDebug(2,Form("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()));
188   AliDebug(2,Form("   with digitization parameters       A = %f\n", fA));
189   AliDebug(2,Form("                                      B = %f\n", fB));
190   AliDebug(2,Form("   Threshold for EC Primary assignment  = %f\n", fECPrimThreshold));
191   AliDebug(2,Form("   Sampling                             = %f\n", fSampling));
192   AliDebug(2,Form("---------------------------------------------------\n"));
193
194   //  Print("");
195
196 }
197
198 //____________________________________________________________________________
199 void AliEMCALSDigitizer::Exec(Option_t *option) 
200
201   // Collects all hit of the same tower into digits
202   TString o(option); o.ToUpper();
203   if (strstr(option, "print") ) {
204
205     AliDebug(2,Form("Print: \n------------------- %s -------------\n",GetName()));
206     AliDebug(2,Form("   fInit                                 %i\n", int(fInit)));
207     AliDebug(2,Form("   fFirstEvent                           %i\n", fFirstEvent));
208     AliDebug(2,Form("   fLastEvent                            %i\n", fLastEvent));
209     AliDebug(2,Form("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()));
210     AliDebug(2,Form("   with digitization parameters       A = %f\n", fA));
211     AliDebug(2,Form("                                      B = %f\n", fB));
212     AliDebug(2,Form("   Threshold for EC Primary assignment  = %f\n", fECPrimThreshold));
213     AliDebug(2,Form("   Sampling                             = %f\n", fSampling));
214     AliDebug(2,Form("---------------------------------------------------\n"));
215
216     //    Print();
217     return ; 
218   }
219   
220   if(strstr(option,"tim"))
221     gBenchmark->Start("EMCALSDigitizer");
222
223   AliRunLoader *rl = AliRunLoader::GetRunLoader();
224   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
225
226   //switch off reloading of this task while getting event
227   if (!fInit) { // to prevent overwrite existing file
228     AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
229     return ;
230     }
231
232   if (fLastEvent == -1) 
233     fLastEvent = rl->GetNumberOfEvents() - 1 ;
234   else {
235     fLastEvent = TMath::Min(fLastEvent, rl->GetNumberOfEvents()-1);
236   }
237   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
238
239   Int_t ievent;
240   Float_t energy=0.; // de * fSampling - 23-nov-04
241   rl->LoadKinematics();
242   rl->LoadHits("EMCAL");
243
244   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
245     rl->GetEvent(ievent);
246     TTree * treeS = emcalLoader->TreeS();
247     if ( !treeS ) { 
248       emcalLoader->MakeSDigitsContainer();
249       treeS = emcalLoader->TreeS();
250     }
251     TClonesArray * hits = emcalLoader->Hits() ; 
252     TClonesArray * sdigits = emcalLoader->SDigits() ;
253     sdigits->Clear();
254
255     Int_t nSdigits = 0 ;
256     Int_t i;
257     AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance(); 
258     for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
259       AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
260       AliEMCALDigit * curSDigit = 0 ;
261       AliEMCALDigit * sdigit = 0 ;
262       Bool_t newsdigit = kTRUE; 
263       
264       // hit->GetId() - Absolute Id number EMCAL segment
265       if(geom->CheckAbsCellId(hit->GetId())) { // was IsInECA(hit->GetId())
266         energy = hit->GetEnergy() * fSampling; // 23-nov-04
267         if(energy >  fECPrimThreshold )
268           // Assign primary number only if deposited energy is significant
269           curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
270                                           hit->GetIparent(), hit->GetId(), 
271                                           Digitize(energy), hit->GetTime() ) ;
272           else
273             curSDigit =  new AliEMCALDigit( -1               , 
274                                             -1               ,
275                                             hit->GetId(), 
276                                             Digitize(energy), hit->GetTime() ) ;
277       } else {
278         Warning("Exec"," abs id %i is bad \n", hit->GetId());
279         newsdigit = kFALSE;
280         curSDigit = 0;
281       }
282       
283       if(curSDigit != 0){
284         for(Int_t check= 0; check < nSdigits ; check++) {
285           sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
286           
287           if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL tower ?              
288             *sdigit = *sdigit + *curSDigit;
289             newsdigit = kFALSE;
290           }
291         }
292       }
293       if (newsdigit) {
294         new((*sdigits)[nSdigits])  AliEMCALDigit(*curSDigit);
295         nSdigits++ ;  
296       }
297       delete curSDigit ; 
298     }  // loop over all hits (hit = deposited energy/entering particle)
299     sdigits->Sort() ;
300     
301     nSdigits = sdigits->GetEntriesFast() ;
302     fSDigitsInRun += nSdigits ;  
303     
304     Double_t e=0.,esum=0.;
305     AliEMCALHistoUtilities::FillH1(fHists, 0, double(sdigits->GetEntriesFast()));
306     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { 
307       AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
308       sdigit->SetIndexInList(i) ;
309       
310       AliEMCALHistoUtilities::FillH1(fHists, 2, double(sdigit->GetAmp()));
311       e = double(Calibrate(sdigit->GetAmp()));
312       esum += e;
313       AliEMCALHistoUtilities::FillH1(fHists, 3, e);
314       AliEMCALHistoUtilities::FillH1(fHists, 4, double(sdigit->GetId()));
315     }
316     if(esum>0.) AliEMCALHistoUtilities::FillH1(fHists, 1, esum);
317     
318     // Now write SDigits    
319     
320     Int_t bufferSize = 32000 ;    
321     TBranch * sdigitsBranch = treeS->GetBranch("EMCAL");
322     if (sdigitsBranch)
323       sdigitsBranch->SetAddress(&sdigits);
324     else
325       treeS->Branch("EMCAL",&sdigits,bufferSize);
326     
327     treeS->Fill();
328     
329     emcalLoader->WriteSDigits("OVERWRITE");
330     
331     //NEXT - SDigitizer
332     emcalLoader->WriteSDigitizer("OVERWRITE");  // why in event cycle ?
333     
334     if(strstr(option,"deb"))
335       PrintSDigits(option) ;  
336   }
337   
338   Unload();
339   
340   emcalLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
341   
342   if(strstr(option,"tim")){
343     gBenchmark->Stop("EMCALSDigitizer"); 
344     printf("\n Exec: took %f seconds for SDigitizing %f seconds per event\n", 
345            gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer")/nEvents ) ; 
346   }
347 }
348
349 //__________________________________________________________________
350
351 Int_t AliEMCALSDigitizer::Digitize(Float_t energy)const {
352   // Digitize the energy
353     Double_t aSignal = fA + energy*fB;
354     if (TMath::Abs(aSignal)>2147483647.0) { 
355       //PH 2147483647 is the max. integer
356       //PH This apparently is a problem which needs investigation
357       AliWarning(Form("Too big or too small energy %f",aSignal));
358       aSignal = TMath::Sign((Double_t)2147483647,aSignal);
359     }
360     return (Int_t ) aSignal;
361   }
362  
363
364 //__________________________________________________________________
365
366 void AliEMCALSDigitizer::Print1(Option_t * option)
367 {
368   Print(); 
369   PrintSDigits(option);
370 }
371
372 //__________________________________________________________________
373 void AliEMCALSDigitizer::Print(Option_t *option) const
374
375   // Prints parameters of SDigitizer
376   printf("Print: \n------------------- %s ------------- option %s\n", GetName() , option) ; 
377   printf("   fInit                                 %i\n", int(fInit));
378   printf("   fFirstEvent                           %i\n", fFirstEvent);
379   printf("   fLastEvent                            %i\n", fLastEvent);
380   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
381   printf("   with digitization parameters       A = %f\n", fA) ; 
382   printf("                                      B = %f\n", fB) ;
383   printf("   Threshold for EC Primary assignment  = %f\n", fECPrimThreshold)  ;
384   printf("   Sampling                             = %f\n", fSampling);
385   printf("---------------------------------------------------\n") ;
386 }
387
388 //__________________________________________________________________
389 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
390 {
391   // Equal operator.
392   // SDititizers are equal if their pedestal, slope and threshold are equal
393   if( (fA==sd.fA)&&(fB==sd.fB)&&
394       (fECPrimThreshold==sd.fECPrimThreshold))
395     return kTRUE ;
396   else
397     return kFALSE ;
398 }
399
400 //__________________________________________________________________
401 void AliEMCALSDigitizer::PrintSDigits(Option_t * option)
402 {
403   //Prints list of digits produced at the current pass of AliEMCALDigitizer
404   
405   
406   // AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
407   AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
408   const TClonesArray * sdigits = rl->SDigits() ; 
409   
410   printf("\n") ;  
411   printf("event %i", rl->GetRunLoader()->GetEventNumber());
412   printf(" Number of entries in SDigits list %i", sdigits->GetEntriesFast()); 
413   if(strstr(option,"all")||strstr(option,"EMC")){
414     
415     //loop over digits
416     AliEMCALDigit * digit;
417     printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
418     Int_t index, isum=0;
419     char * tempo = new char[8192]; 
420     for (index = 0 ; index < sdigits->GetEntries() ; index++) {
421       digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
422       sprintf(tempo, "\n%6d  %8d    %6.5e %4d      %2d :",
423               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
424       printf(tempo);
425       isum += digit->GetAmp();
426       
427       Int_t iprimary;
428       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
429         sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
430         printf(tempo); 
431       }          
432     }
433     delete tempo ;
434     printf("\n** Sum %i : %10.3f GeV/c **\n ", isum, double(isum)*1.e-6);
435   } else printf("\n");
436 }
437
438 //____________________________________________________________________________ 
439 void AliEMCALSDigitizer::Unload() const
440 {
441   // Unload Hits and SDigits from the folder
442   AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
443   rl->UnloadHits() ; 
444   rl->UnloadSDigits() ; 
445 }
446
447 void AliEMCALSDigitizer::Browse(TBrowser* b)
448 {
449   if(fHists) b->Add(fHists);
450   TTask::Browse(b);
451 }
452
453 TList *AliEMCALSDigitizer::BookControlHists(int var)
454
455   //book histograms for monitoring sdigitization
456   // 22-nov-04
457   gROOT->cd();
458   const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance() ;
459   if(var>=1){
460     AliDebug(1, " BookControlHists() in action ");
461     new TH1F("HSDigiN",  "#EMCAL  sdigits ", 1001, -0.5, 1000.5);
462     new TH1F("HSDigiSumEnergy","Sum.EMCAL energy", 1000, 0.0, 100.);
463     new TH1F("HSDigiAmp",  "EMCAL sdigits amplitude", 1000, 0., 2.e+9);
464     new TH1F("HSDigiEnergy","EMCAL cell energy", 1000, 0.0, 100.);
465     new TH1F("HSDigiAbsId","EMCAL absID for sdigits",
466     geom->GetNCells(), 0.5, Double_t(geom->GetNCells())+0.5);
467   }
468
469   fHists = AliEMCALHistoUtilities::MoveHistsToList("EmcalSDigiControlHists", kFALSE);
470   //  fHists = 0; ??
471
472   return fHists;
473 }
474
475 void AliEMCALSDigitizer::SaveHists(const char* name, Bool_t kSingleKey, const char* opt)
476 {
477   AliEMCALHistoUtilities::SaveListOfHists(fHists, name, kSingleKey, opt); 
478 }