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 **************************************************************************/
19 //_________________________________________________________________________
20 // This is a TTask that makes SDigits out of Hits
21 // The name of the TTask is also the title of the branch that will contain
22 // the created SDigits
23 // The title of the TTAsk is the name of the file that contains the hits from
24 // which the SDigits are created
25 // A Summable Digits is the sum of all hits originating
26 // from one primary in one active cell
27 // A threshold for assignment of the primary to SDigit is applied
28 // SDigits are written to TreeS, branch "PHOS"
29 // AliPHOSSDigitizer with all current parameters is written
30 // to TreeS branch "AliPHOSSDigitizer".
31 // Both branches have the same title. If necessary one can produce
32 // another set of SDigits with different parameters. Two versions
33 // can be distunguished using titles of the branches.
35 // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
36 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
37 // root [1] s->ExecuteTask()
38 // // Makes SDigitis for all events stored in galice.root
39 // root [2] s->SetPedestalParameter(0.001)
40 // // One can change parameters of digitization
41 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
42 // // and write them into the new branch
43 // root [4] s->ExecuteTask("deb all tim")
44 // // available parameters:
45 // deb - print # of produced SDigitis
46 // deb all - print # and list of produced SDigits
47 // tim - print benchmarking information
49 //*-- Author : Dmitri Peressounko (SUBATECH & KI)
50 //////////////////////////////////////////////////////////////////////////////
53 // --- ROOT system ---
54 #include "TBenchmark.h"
57 // --- Standard library ---
59 // --- AliRoot header files ---
60 #include "AliPHOSGeometry.h"
61 #include "AliPHOSDigit.h"
62 #include "AliPHOSGetter.h"
63 #include "AliPHOSHit.h"
64 #include "AliPHOSSDigitizer.h"
66 ClassImp(AliPHOSSDigitizer)
69 //____________________________________________________________________________
70 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","")
73 fFirstEvent = fLastEvent = 0 ;
74 fDefaultInit = kTRUE ;
77 //____________________________________________________________________________
78 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName,
79 const char * eventFolderName):
80 TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
81 fEventFolderName(eventFolderName)
85 fFirstEvent = fLastEvent = 0 ; // runs one event by defaut
88 fDefaultInit = kFALSE ;
91 //____________________________________________________________________________
92 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer & sd)
97 fFirstEvent = sd.fFirstEvent ;
98 fLastEvent = sd.fLastEvent ;
101 fPrimThreshold = sd.fPrimThreshold ;
102 fSDigitsInRun = sd.fSDigitsInRun ;
103 SetName(sd.GetName()) ;
104 SetTitle(sd.GetTitle()) ;
105 fEventFolderName = sd.fEventFolderName;
109 //____________________________________________________________________________
110 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
112 AliPHOSGetter * gime =
113 AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());
114 gime->PhosLoader()->CleanSDigitizer();
116 //____________________________________________________________________________
117 void AliPHOSSDigitizer::Init()
119 // Uses the getter to access the required files
123 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
125 Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
129 TString opt("SDigits") ;
130 if(gime->VersionExists(opt) ) {
131 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
135 gime->PostSDigitizer(this);
136 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
140 //____________________________________________________________________________
141 void AliPHOSSDigitizer::InitParameters()
143 // initializes the parameters for digitization
146 fPrimThreshold = 0.01 ;
150 //____________________________________________________________________________
151 void AliPHOSSDigitizer::Exec(Option_t *option)
153 // Steering method to produce summable digits for events
154 // in the range from fFirstEvent to fLastEvent.
155 // This range is optionally set by SetEventRange().
156 // if fLastEvent=-1 (by default), then process events until the end.
158 // Summable digit is a sum of all hits in the same active
161 if (strstr(option, "print") ) {
166 if(strstr(option,"tim"))
167 gBenchmark->Start("PHOSSDigitizer");
169 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
171 //switch off reloading of this task while getting event
172 if (!fInit) { // to prevent overwrite existing file
173 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
177 if (fLastEvent == -1)
178 fLastEvent = gime->MaxEvent() - 1 ;
180 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time
181 Int_t nEvents = fLastEvent - fFirstEvent + 1;
184 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
185 gime->Event(ievent,"H") ;
186 TTree * treeS = gime->TreeS();
187 TClonesArray * hits = gime->Hits() ;
188 TClonesArray * sdigits = gime->SDigits() ;
191 //Now make SDigits from hits, for PHOS it is the same, so just copy
192 Int_t nPrim = static_cast<Int_t>((gime->TreeH())->GetEntries()) ;
193 // Attention nPrim is the number of primaries tracked by Geant
194 // and this number could be different to the number of Primaries in TreeK;
197 for (iprim = 0 ; iprim < nPrim ; iprim ++) {
198 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
201 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
202 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
203 // Assign primary number only if contribution is significant
205 if( hit->GetEnergy() > fPrimThreshold)
206 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
207 Digitize(hit->GetEnergy()), hit->GetTime()) ;
209 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
210 Digitize(hit->GetEnergy()), hit->GetTime()) ;
219 nSdigits = sdigits->GetEntriesFast() ;
221 fSDigitsInRun += nSdigits ;
222 sdigits->Expand(nSdigits) ;
225 for (i = 0 ; i < nSdigits ; i++) {
226 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
227 digit->SetIndexInList(i) ;
233 //First list of sdigits
235 Int_t bufferSize = 32000 ;
236 TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
238 sdigitsBranch->Fill() ;
240 gime->WriteSDigits("OVERWRITE");
244 gime->WriteSDigitizer("OVERWRITE");
246 if(strstr(option,"deb"))
247 PrintSDigits(option) ;
252 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
254 if(strstr(option,"tim")){
255 gBenchmark->Stop("PHOSSDigitizer");
256 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
257 gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
261 //__________________________________________________________________
262 void AliPHOSSDigitizer::Print()const
264 // Prints parameters of SDigitizer
265 Info("Print", "\n------------------- %s -------------", GetName() ) ;
266 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
267 printf(" with digitization parameters A = %f\n", fA) ;
268 printf(" B = %f\n", fB) ;
269 printf(" Threshold for Primary assignment= %f\n", fPrimThreshold) ;
270 printf("---------------------------------------------------\n") ;
274 //__________________________________________________________________
275 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
278 // SDititizers are equal if their pedestal, slope and threshold are equal
280 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
286 //__________________________________________________________________
287 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
289 // Prints list of digits produced in the current pass of AliPHOSDigitizer
292 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
293 const TClonesArray * sdigits = gime->SDigits() ;
295 Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ;
297 if(strstr(option,"all")||strstr(option,"EMC")){
300 AliPHOSDigit * digit;
301 printf("\nEMC sdigits\n") ;
302 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
304 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
305 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
306 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
307 // if(digit->GetNprimary() == 0)
309 printf("%6d %8d %6.5e %4d %2d :\n",
310 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
312 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
313 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
318 if(strstr(option,"all")||strstr(option,"CPV")){
320 //loop over CPV digits
321 AliPHOSDigit * digit;
322 printf("\nCPV sdigits\n") ;
323 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
325 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
326 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
327 if(digit->GetId() > maxEmc){
328 printf("\n%6d %8d %4d %2d :",
329 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
331 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
332 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
339 //____________________________________________________________________________
340 void AliPHOSSDigitizer::Unload() const
342 // Unloads the objects from the folder
343 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
344 AliPHOSLoader * loader = gime->PhosLoader() ;
345 loader->UnloadHits() ;
346 loader->UnloadSDigits() ;