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 /* History of cvs commits:
24 //_________________________________________________________________________
25 // This is a TTask that makes SDigits out of Hits
26 // The name of the TTask is also the title of the branch that will contain
27 // the created SDigits
28 // The title of the TTAsk is the name of the file that contains the hits from
29 // which the SDigits are created
30 // A Summable Digits is the sum of all hits originating
31 // from one primary in one active cell
32 // A threshold for assignment of the primary to SDigit is applied
33 // SDigits are written to TreeS, branch "PHOS"
34 // AliPHOSSDigitizer with all current parameters is written
35 // to TreeS branch "AliPHOSSDigitizer".
36 // Both branches have the same title. If necessary one can produce
37 // another set of SDigits with different parameters. Two versions
38 // can be distunguished using titles of the branches.
40 // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
41 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 // root [1] s->ExecuteTask()
43 // // Makes SDigitis for all events stored in galice.root
44 // root [2] s->SetPedestalParameter(0.001)
45 // // One can change parameters of digitization
46 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
47 // // and write them into the new branch
48 // root [4] s->ExecuteTask("deb all tim")
49 // // available parameters:
50 // deb - print # of produced SDigitis
51 // deb all - print # and list of produced SDigits
52 // tim - print benchmarking information
54 //*-- Author : Dmitri Peressounko (SUBATECH & KI)
55 //////////////////////////////////////////////////////////////////////////////
58 // --- ROOT system ---
59 #include "TBenchmark.h"
60 //#include "TObjectTable.h"
62 // --- Standard library ---
64 // --- AliRoot header files ---
66 #include "AliPHOSGeometry.h"
67 #include "AliPHOSDigit.h"
68 #include "AliPHOSGetter.h"
69 #include "AliPHOSHit.h"
70 #include "AliPHOSSDigitizer.h"
71 //#include "AliMemoryWatcher.h"
73 ClassImp(AliPHOSSDigitizer)
76 //____________________________________________________________________________
77 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","")
80 fFirstEvent = fLastEvent = 0 ;
81 fDefaultInit = kTRUE ;
84 //____________________________________________________________________________
85 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName,
86 const char * eventFolderName):
87 TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
88 fEventFolderName(eventFolderName)
92 fFirstEvent = fLastEvent = 0 ; // runs one event by defaut
95 fDefaultInit = kFALSE ;
98 //____________________________________________________________________________
99 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer & sd)
104 fFirstEvent = sd.fFirstEvent ;
105 fLastEvent = sd.fLastEvent ;
108 fPrimThreshold = sd.fPrimThreshold ;
109 fSDigitsInRun = sd.fSDigitsInRun ;
110 SetName(sd.GetName()) ;
111 SetTitle(sd.GetTitle()) ;
112 fEventFolderName = sd.fEventFolderName;
116 //____________________________________________________________________________
117 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
119 // AliPHOSGetter * gime =
120 // AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());
121 AliPHOSGetter * gime =
122 AliPHOSGetter::Instance();
123 gime->PhosLoader()->CleanSDigitizer();
125 //____________________________________________________________________________
126 void AliPHOSSDigitizer::Init()
128 // Uses the getter to access the required files
132 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
134 Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
138 TString opt("SDigits") ;
139 if(gime->VersionExists(opt) ) {
140 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
144 gime->PostSDigitizer(this);
145 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
149 //____________________________________________________________________________
150 void AliPHOSSDigitizer::InitParameters()
152 // initializes the parameters for digitization
155 fPrimThreshold = 0.01 ;
159 //____________________________________________________________________________
160 void AliPHOSSDigitizer::Exec(Option_t *option)
162 // Steering method to produce summable digits for events
163 // in the range from fFirstEvent to fLastEvent.
164 // This range is optionally set by SetEventRange().
165 // if fLastEvent=-1 (by default), then process events until the end.
167 // Summable digit is a sum of all hits in the same active
170 if (strstr(option, "print") ) {
175 if(strstr(option,"tim"))
176 gBenchmark->Start("PHOSSDigitizer");
178 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
180 //switch off reloading of this task while getting event
181 if (!fInit) { // to prevent overwrite existing file
182 AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
186 if (fLastEvent == -1)
187 fLastEvent = gime->MaxEvent() - 1 ;
189 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time
190 Int_t nEvents = fLastEvent - fFirstEvent + 1;
194 //AliMemoryWatcher memwatcher;
196 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
197 gime->Event(ievent,"H") ;
198 TTree * treeS = gime->TreeS();
199 TClonesArray * hits = gime->Hits() ;
200 TClonesArray * sdigits = gime->SDigits() ;
203 //Now make SDigits from hits, for PHOS it is the same, so just copy
204 Int_t nPrim = static_cast<Int_t>((gime->TreeH())->GetEntries()) ;
205 // Attention nPrim is the number of primaries tracked by Geant
206 // and this number could be different to the number of Primaries in TreeK;
209 for (iprim = 0 ; iprim < nPrim ; iprim ++) {
210 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
213 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
214 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
215 // Assign primary number only if contribution is significant
217 if( hit->GetEnergy() > fPrimThreshold)
218 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
219 Digitize(hit->GetEnergy()), hit->GetTime()) ;
221 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
222 Digitize(hit->GetEnergy()), hit->GetTime()) ;
231 nSdigits = sdigits->GetEntriesFast() ;
233 fSDigitsInRun += nSdigits ;
234 sdigits->Expand(nSdigits) ;
237 for (i = 0 ; i < nSdigits ; i++) {
238 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
239 digit->SetIndexInList(i) ;
245 //First list of sdigits
247 Int_t bufferSize = 32000 ;
248 TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
250 sdigitsBranch->Fill() ;
252 gime->WriteSDigits("OVERWRITE");
256 gime->WriteSDigitizer("OVERWRITE");
257 //gObjectTable->Print() ;
259 if(strstr(option,"deb"))
260 PrintSDigits(option) ;
262 //memwatcher.Watch(ievent);
267 // gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
269 if(strstr(option,"tim")){
270 gBenchmark->Stop("PHOSSDigitizer");
271 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
272 gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
275 //TFile f("out.root","RECREATE");
276 //memwatcher.WriteToFile();
280 //__________________________________________________________________
281 void AliPHOSSDigitizer::Print(const Option_t *)const
283 // Prints parameters of SDigitizer
284 Info("Print", "\n------------------- %s -------------", GetName() ) ;
285 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
286 printf(" with digitization parameters A = %f\n", fA) ;
287 printf(" B = %f\n", fB) ;
288 printf(" Threshold for Primary assignment= %f\n", fPrimThreshold) ;
289 printf("---------------------------------------------------\n") ;
293 //__________________________________________________________________
294 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
297 // SDititizers are equal if their pedestal, slope and threshold are equal
299 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
305 //__________________________________________________________________
306 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
308 // Prints list of digits produced in the current pass of AliPHOSDigitizer
311 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
312 const TClonesArray * sdigits = gime->SDigits() ;
314 Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ;
316 if(strstr(option,"all")||strstr(option,"EMC")){
319 AliPHOSDigit * digit;
320 printf("\nEMC sdigits\n") ;
321 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
323 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
324 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
325 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
326 // if(digit->GetNprimary() == 0)
328 printf("%6d %8d %6.5e %4d %2d :\n",
329 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
331 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
332 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
337 if(strstr(option,"all")||strstr(option,"CPV")){
339 //loop over CPV digits
340 AliPHOSDigit * digit;
341 printf("\nCPV sdigits\n") ;
342 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
344 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
345 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
346 if(digit->GetId() > maxEmc){
347 printf("\n%6d %8d %4d %2d :",
348 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
350 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
351 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
358 //____________________________________________________________________________
359 void AliPHOSSDigitizer::Unload() const
361 // Unloads the objects from the folder
362 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
363 AliPHOSLoader * loader = gime->PhosLoader() ;
364 loader->UnloadHits() ;
365 loader->UnloadSDigits() ;