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
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.
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
44 //*-- Author : Sahal Yacoob (LBL)
45 // based on : AliPHOSSDigitzer
47 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
48 // of new IO (à la PHOS)
49 //////////////////////////////////////////////////////////////////////////////
51 // --- ROOT system ---
54 #include "TBenchmark.h"
56 // --- Standard library ---
58 // --- AliRoot header files ---
60 #include "AliEMCALDigit.h"
61 #include "AliEMCALGetter.h"
62 #include "AliEMCALHit.h"
63 #include "AliEMCALSDigitizer.h"
64 #include "AliEMCALGeometry.h"
66 ClassImp(AliEMCALSDigitizer)
68 //____________________________________________________________________________
69 AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("","")
73 fDefaultInit = kTRUE ;
76 //____________________________________________________________________________
77 AliEMCALSDigitizer::AliEMCALSDigitizer(const char * alirunFileName, const char * eventFolderName):
78 TTask("EMCAL"+AliConfig::fgkSDigitizerTaskName, alirunFileName),
79 fEventFolderName(eventFolderName)
84 fDefaultInit = kFALSE ;
88 //____________________________________________________________________________
89 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd) : TTask(sd) {
94 fPREPrimThreshold = sd.fPREPrimThreshold ;
95 fECPrimThreshold = sd.fECPrimThreshold ;
96 fHCPrimThreshold = sd.fHCPrimThreshold ;
97 fSDigitsInRun = sd.fSDigitsInRun ;
98 SetName(sd.GetName()) ;
99 SetTitle(sd.GetTitle()) ;
100 fEventFolderName = sd.fEventFolderName;
104 //____________________________________________________________________________
105 void AliEMCALSDigitizer::Init(){
106 // Initialization: open root-file, allocate arrays for hits and sdigits,
107 // attach task SDigitizer to the list of EMCAL tasks
109 // Initialization can not be done in the default constructor
110 //============================================================= YS
111 // The initialisation is now done by the getter
115 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());
117 Fatal("Init", "Could not obtain the Getter objectfor file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
121 TString opt("SDigits") ;
122 if(gime->VersionExists(opt) ) {
123 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
127 gime->PostSDigitizer(this);
128 gime->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
132 //____________________________________________________________________________
133 void AliEMCALSDigitizer::InitParameters()
138 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
139 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
140 if (geom->GetSampling() == 0.) {
141 Error("InitParameters", "Sampling factor not set !") ;
145 Info("InitParameters", "Sampling factor set to %f\n", geom->GetSampling()) ;
147 // this threshold corresponds approximately to 100 MeV
148 fECPrimThreshold = 100E-3 / ( geom->GetSampling() * ( geom->GetNPRLayers() + geom->GetNECLayers()) ) * geom->GetNECLayers() ;
149 fPREPrimThreshold = 100E-3 / ( geom->GetSampling() * ( geom->GetNPRLayers() + geom->GetNECLayers()) ) * geom->GetNPRLayers() ;
150 fHCPrimThreshold = fECPrimThreshold/5. ; // 5 is totally arbitrary
154 //____________________________________________________________________________
155 void AliEMCALSDigitizer::Exec(Option_t *option)
157 // Collects all hits in the section (PRE/ECAL/HCAL) of the same tower into digit
159 if (strstr(option, "print") ) {
164 if(strstr(option,"tim"))
165 gBenchmark->Start("EMCALSDigitizer");
167 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
169 //switch off reloading of this task while getting event
170 if (!fInit) { // to prevent overwrite existing file
171 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
175 Int_t nevents = gime->MaxEvent() ;
177 for(ievent = 0; ievent < nevents; ievent++){
179 gime->Event(ievent,"H") ;
181 TTree * treeS = gime->TreeS();
182 TClonesArray * hits = gime->Hits() ;
183 TClonesArray * sdigits = gime->SDigits() ;
186 //Now make SDigits from hits, for EMCAL it is the same, so just copy
187 Int_t nPrim = static_cast<Int_t>((gime->TreeH())->GetEntries()) ;
188 // Attention nPrim is the number of primaries tracked by Geant
189 // and this number could be different to the number of Primaries in TreeK;
192 for ( iprim = 0 ; iprim < nPrim ; iprim++ ) {
193 //=========== Get the EMCAL branch from Hits Tree for the Primary iprim
196 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
197 AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
198 AliEMCALDigit * curSDigit = 0 ;
199 AliEMCALDigit * sdigit = 0 ;
200 Bool_t newsdigit = kTRUE;
202 // Assign primary number only if deposited energy is significant
204 AliEMCALGeometry * geom = gime->EMCALGeometry() ;
206 if( geom->IsInPRE(hit->GetId()) )
207 if( hit->GetEnergy() > fPREPrimThreshold )
208 curSDigit = new AliEMCALDigit( hit->GetPrimary(),
209 hit->GetIparent(), hit->GetId(),
210 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
212 curSDigit = new AliEMCALDigit( -1 ,
215 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
216 else if( geom->IsInECA(hit->GetId()) )
217 if( hit->GetEnergy() > fECPrimThreshold )
218 curSDigit = new AliEMCALDigit( hit->GetPrimary(),
219 hit->GetIparent(), hit->GetId(),
220 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
222 curSDigit = new AliEMCALDigit( -1 ,
225 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
226 else if( geom->IsInHCA(hit->GetId()) )
227 if( hit->GetEnergy() > fHCPrimThreshold )
229 curSDigit = new AliEMCALDigit( hit->GetPrimary(),
230 hit->GetIparent(), hit->GetId(),
231 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
233 curSDigit = new AliEMCALDigit( -1 ,
236 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
239 for(check= 0; check < nSdigits ; check++) {
240 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
241 if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL/HCAL/preshower tower ?
242 *sdigit = *sdigit + *curSDigit;
247 new((*sdigits)[nSdigits]) AliEMCALDigit(*curSDigit);
251 } // loop over all hits (hit = deposited energy/layer/entering particle)
256 nSdigits = sdigits->GetEntriesFast() ;
257 fSDigitsInRun += nSdigits ;
258 sdigits->Expand(nSdigits) ;
260 Int_t NPrimarymax = -1 ;
262 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
263 AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
264 sdigit->SetIndexInList(i) ;
267 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
268 if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > NPrimarymax)
269 NPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
274 //First list of sdigits
276 Int_t bufferSize = 32000 ;
277 TBranch * sdigitsBranch = treeS->Branch("EMCAL",&sdigits,bufferSize);
279 sdigitsBranch->Fill() ;
281 gime->WriteSDigits("OVERWRITE");
285 gime->WriteSDigitizer("OVERWRITE");
287 if(strstr(option,"deb"))
288 PrintSDigits(option) ;
293 gime->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
295 if(strstr(option,"tim")){
296 gBenchmark->Stop("EMCALSDigitizer");
297 Info("Exec", "took %f seconds for SDigitizing %f seconds per event",
298 gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer") ) ;
303 //__________________________________________________________________
304 void AliEMCALSDigitizer::Print(Option_t*)const
306 // Prints parameters of SDigitizer
307 Info("Print", "\n------------------- %s -------------", GetName() ) ;
308 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
309 printf(" with digitization parameters A = %f\n", fA) ;
310 printf(" B = %f\n", fB) ;
311 printf(" Threshold for PRE Primary assignment= %f\n", fPREPrimThreshold) ;
312 printf(" Threshold for EC Primary assignment= %f\n", fECPrimThreshold) ;
313 printf(" Threshold for HC Primary assignment= %f\n", fHCPrimThreshold) ;
314 printf("---------------------------------------------------\n") ;
318 //__________________________________________________________________
319 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
322 // SDititizers are equal if their pedestal, slope and threshold are equal
324 if( (fA==sd.fA)&&(fB==sd.fB)&&
325 (fECPrimThreshold==sd.fECPrimThreshold) &&
326 (fHCPrimThreshold==sd.fHCPrimThreshold) &&
327 (fPREPrimThreshold==sd.fPREPrimThreshold))
333 //__________________________________________________________________
334 void AliEMCALSDigitizer::PrintSDigits(Option_t * option){
335 //Prints list of digits produced at the current pass of AliEMCALDigitizer
337 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
338 const TClonesArray * sdigits = gime->SDigits() ;
340 TString message("\n") ;
341 message += "event " ;
342 message += gAlice->GetEvNumber() ;
343 message += "\n Number of entries in SDigits list " ;
344 message += sdigits->GetEntriesFast() ;
345 if(strstr(option,"all")||strstr(option,"EMC")){
348 AliEMCALDigit * digit;
349 message += "\n Id Amplitude Time Index Nprim: Primaries list \n" ;
351 char * tempo = new char[8192];
352 for (index = 0 ; index < sdigits->GetEntries() ; index++) {
353 digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
354 sprintf(tempo, "\n%6d %8d %6.5e %4d %2d :",
355 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
359 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
360 sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ;
366 Info("PrintSDigits", message.Data() ) ;
369 //____________________________________________________________________________
370 void AliEMCALSDigitizer::Unload() const
372 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
373 AliEMCALLoader * loader = gime->EmcalLoader() ;
374 loader->UnloadHits() ;
375 loader->UnloadSDigits() ;