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:
22 * Revision 1.48 2006/04/22 10:30:17 hristov
23 * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
25 * Revision 1.47 2005/05/28 14:19:04 schutz
26 * Compilation warnings fixed by T.P.
30 //_________________________________________________________________________
31 // This is a TTask that makes SDigits out of Hits
32 // The name of the TTask is also the title of the branch that will contain
33 // the created SDigits
34 // The title of the TTAsk is the name of the file that contains the hits from
35 // which the SDigits are created
36 // A Summable Digits is the sum of all hits originating
37 // from one primary in one active cell
38 // A threshold for assignment of the primary to SDigit is applied
39 // SDigits are written to TreeS, branch "PHOS"
40 // AliPHOSSDigitizer with all current parameters is written
41 // to TreeS branch "AliPHOSSDigitizer".
42 // Both branches have the same title. If necessary one can produce
43 // another set of SDigits with different parameters. Two versions
44 // can be distunguished using titles of the branches.
46 // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
47 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
48 // root [1] s->ExecuteTask()
49 // // Makes SDigitis for all events stored in galice.root
50 // root [2] s->SetPedestalParameter(0.001)
51 // // One can change parameters of digitization
52 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
53 // // and write them into the new branch
54 // root [4] s->ExecuteTask("deb all tim")
55 // // available parameters:
56 // deb - print # of produced SDigitis
57 // deb all - print # and list of produced SDigits
58 // tim - print benchmarking information
60 //*-- Author : Dmitri Peressounko (SUBATECH & KI)
61 //////////////////////////////////////////////////////////////////////////////
64 // --- ROOT system ---
65 #include "TBenchmark.h"
66 //#include "TObjectTable.h"
68 // --- Standard library ---
70 // --- AliRoot header files ---
72 #include "AliPHOSGeometry.h"
73 #include "AliPHOSDigit.h"
74 #include "AliPHOSGetter.h"
75 #include "AliPHOSHit.h"
76 #include "AliPHOSSDigitizer.h"
77 //#include "AliMemoryWatcher.h"
79 ClassImp(AliPHOSSDigitizer)
82 //____________________________________________________________________________
83 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","")
86 fFirstEvent = fLastEvent = 0 ;
87 fDefaultInit = kTRUE ;
90 //____________________________________________________________________________
91 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName,
92 const char * eventFolderName):
93 TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
94 fEventFolderName(eventFolderName)
98 fFirstEvent = fLastEvent = 0 ; // runs one event by defaut
101 fDefaultInit = kFALSE ;
104 //____________________________________________________________________________
105 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer & sd)
110 fFirstEvent = sd.fFirstEvent ;
111 fLastEvent = sd.fLastEvent ;
114 fPrimThreshold = sd.fPrimThreshold ;
115 fSDigitsInRun = sd.fSDigitsInRun ;
116 SetName(sd.GetName()) ;
117 SetTitle(sd.GetTitle()) ;
118 fEventFolderName = sd.fEventFolderName;
122 //____________________________________________________________________________
123 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
125 // AliPHOSGetter * gime =
126 // AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());
127 AliPHOSGetter * gime =
128 AliPHOSGetter::Instance();
129 gime->PhosLoader()->CleanSDigitizer();
131 //____________________________________________________________________________
132 void AliPHOSSDigitizer::Init()
134 // Uses the getter to access the required files
138 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
140 Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
144 TString opt("SDigits") ;
145 if(gime->VersionExists(opt) ) {
146 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
150 gime->PostSDigitizer(this);
151 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
155 //____________________________________________________________________________
156 void AliPHOSSDigitizer::InitParameters()
158 // initializes the parameters for digitization
161 fPrimThreshold = 0.01 ;
165 //____________________________________________________________________________
166 void AliPHOSSDigitizer::Exec(Option_t *option)
168 // Steering method to produce summable digits for events
169 // in the range from fFirstEvent to fLastEvent.
170 // This range is optionally set by SetEventRange().
171 // if fLastEvent=-1 (by default), then process events until the end.
173 // Summable digit is a sum of all hits in the same active
176 if (strstr(option, "print") ) {
181 if(strstr(option,"tim"))
182 gBenchmark->Start("PHOSSDigitizer");
184 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
186 //switch off reloading of this task while getting event
187 if (!fInit) { // to prevent overwrite existing file
188 AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
192 if (fLastEvent == -1)
193 fLastEvent = gime->MaxEvent() - 1 ;
195 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time
196 Int_t nEvents = fLastEvent - fFirstEvent + 1;
200 //AliMemoryWatcher memwatcher;
202 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
203 gime->Event(ievent,"H") ;
204 TTree * treeS = gime->TreeS();
205 TClonesArray * hits = gime->Hits() ;
206 TClonesArray * sdigits = gime->SDigits() ;
209 //Now make SDigits from hits, for PHOS it is the same, so just copy
210 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
211 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
212 // Assign primary number only if contribution is significant
214 if( hit->GetEnergy() > fPrimThreshold)
215 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
216 hit->GetEnergy() ,hit->GetTime()) ;
218 new((*sdigits)[nSdigits]) AliPHOSDigit(-1 ,hit->GetId(),
219 hit->GetEnergy() ,hit->GetTime()) ;
226 nSdigits = sdigits->GetEntriesFast() ;
228 fSDigitsInRun += nSdigits ;
229 sdigits->Expand(nSdigits) ;
231 for (i = 0 ; i < nSdigits ; i++) {
232 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
233 digit->SetIndexInList(i) ;
239 //First list of sdigits
241 Int_t bufferSize = 32000 ;
242 TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
244 sdigitsBranch->Fill() ;
246 gime->WriteSDigits("OVERWRITE");
250 gime->WriteSDigitizer("OVERWRITE");
251 //gObjectTable->Print() ;
253 if(strstr(option,"deb"))
254 PrintSDigits(option) ;
256 //memwatcher.Watch(ievent);
261 // gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
263 if(strstr(option,"tim")){
264 gBenchmark->Stop("PHOSSDigitizer");
265 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
266 gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
269 //TFile f("out.root","RECREATE");
270 //memwatcher.WriteToFile();
274 //__________________________________________________________________
275 void AliPHOSSDigitizer::Print(const Option_t *)const
277 // Prints parameters of SDigitizer
278 Info("Print", "\n------------------- %s -------------", GetName() ) ;
279 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
280 printf(" with digitization parameters A = %f\n", fA) ;
281 printf(" B = %f\n", fB) ;
282 printf(" Threshold for Primary assignment= %f\n", fPrimThreshold) ;
283 printf("---------------------------------------------------\n") ;
287 //__________________________________________________________________
288 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
291 // SDititizers are equal if their pedestal, slope and threshold are equal
293 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
299 //__________________________________________________________________
300 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
302 // Prints list of digits produced in the current pass of AliPHOSDigitizer
305 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
306 const TClonesArray * sdigits = gime->SDigits() ;
308 Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ;
310 if(strstr(option,"all")||strstr(option,"EMC")){
313 AliPHOSDigit * digit;
314 printf("\nEMC sdigits\n") ;
315 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
317 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
318 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
319 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
320 // if(digit->GetNprimary() == 0)
322 // printf("%6d %8d %6.5e %4d %2d :\n", // YVK
323 printf("%6d %.4f %6.5e %4d %2d :\n",
324 digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
326 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
327 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
332 if(strstr(option,"all")||strstr(option,"CPV")){
334 //loop over CPV digits
335 AliPHOSDigit * digit;
336 printf("\nCPV sdigits\n") ;
337 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
339 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
340 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
341 if(digit->GetId() > maxEmc){
342 printf("\n%6d %8d %4d %2d :",
343 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
345 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
346 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
353 //____________________________________________________________________________
354 void AliPHOSSDigitizer::Unload() const
356 // Unloads the objects from the folder
357 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
358 AliPHOSLoader * loader = gime->PhosLoader() ;
359 loader->UnloadHits() ;
360 loader->UnloadSDigits() ;