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.49 2006/05/10 06:42:53 kharlov
23 * Remove redundant loop over primaries
25 * Revision 1.48 2006/04/22 10:30:17 hristov
26 * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
28 * Revision 1.47 2005/05/28 14:19:04 schutz
29 * Compilation warnings fixed by T.P.
33 //_________________________________________________________________________
34 // This is a TTask that makes SDigits out of Hits
35 // The name of the TTask is also the title of the branch that will contain
36 // the created SDigits
37 // The title of the TTAsk is the name of the file that contains the hits from
38 // which the SDigits are created
39 // A Summable Digits is the sum of all hits originating
40 // from one primary in one active cell
41 // A threshold for assignment of the primary to SDigit is applied
42 // SDigits are written to TreeS, branch "PHOS"
43 // AliPHOSSDigitizer with all current parameters is written
44 // to TreeS branch "AliPHOSSDigitizer".
45 // Both branches have the same title. If necessary one can produce
46 // another set of SDigits with different parameters. Two versions
47 // can be distunguished using titles of the branches.
49 // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
50 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
51 // root [1] s->ExecuteTask()
52 // // Makes SDigitis for all events stored in galice.root
53 // root [2] s->SetPedestalParameter(0.001)
54 // // One can change parameters of digitization
55 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
56 // // and write them into the new branch
57 // root [4] s->ExecuteTask("deb all tim")
58 // // available parameters:
59 // deb - print # of produced SDigitis
60 // deb all - print # and list of produced SDigits
61 // tim - print benchmarking information
63 //*-- Author : Dmitri Peressounko (SUBATECH & KI)
64 //////////////////////////////////////////////////////////////////////////////
67 // --- ROOT system ---
68 #include "TBenchmark.h"
69 //#include "TObjectTable.h"
71 // --- Standard library ---
73 // --- AliRoot header files ---
75 #include "AliPHOSGeometry.h"
76 #include "AliPHOSDigit.h"
77 #include "AliPHOSGetter.h"
78 #include "AliPHOSHit.h"
79 #include "AliPHOSSDigitizer.h"
80 //#include "AliMemoryWatcher.h"
82 ClassImp(AliPHOSSDigitizer)
85 //____________________________________________________________________________
86 AliPHOSSDigitizer::AliPHOSSDigitizer() :
100 //____________________________________________________________________________
101 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName,
102 const char * eventFolderName):
103 TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
106 fDefaultInit(kFALSE),
107 fEventFolderName(eventFolderName),
116 fDefaultInit = kFALSE ;
119 //____________________________________________________________________________
120 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
122 // AliPHOSGetter * gime =
123 // AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());
124 AliPHOSGetter * gime =
125 AliPHOSGetter::Instance();
126 gime->PhosLoader()->CleanSDigitizer();
128 //____________________________________________________________________________
129 void AliPHOSSDigitizer::Init()
131 // Uses the getter to access the required files
135 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
137 Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
141 TString opt("SDigits") ;
142 if(gime->VersionExists(opt) ) {
143 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
147 gime->PostSDigitizer(this);
148 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
152 //____________________________________________________________________________
153 void AliPHOSSDigitizer::InitParameters()
155 // initializes the parameters for digitization
158 fPrimThreshold = 0.01 ;
162 //____________________________________________________________________________
163 void AliPHOSSDigitizer::Exec(Option_t *option)
165 // Steering method to produce summable digits for events
166 // in the range from fFirstEvent to fLastEvent.
167 // This range is optionally set by SetEventRange().
168 // if fLastEvent=-1 (by default), then process events until the end.
170 // Summable digit is a sum of all hits in the same active
173 if (strstr(option, "print") ) {
178 if(strstr(option,"tim"))
179 gBenchmark->Start("PHOSSDigitizer");
181 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
183 //switch off reloading of this task while getting event
184 if (!fInit) { // to prevent overwrite existing file
185 AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
189 if (fLastEvent == -1)
190 fLastEvent = gime->MaxEvent() - 1 ;
192 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time
193 Int_t nEvents = fLastEvent - fFirstEvent + 1;
197 //AliMemoryWatcher memwatcher;
199 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
200 gime->Event(ievent,"H") ;
201 TTree * treeS = gime->TreeS();
202 TClonesArray * hits = gime->Hits() ;
203 TClonesArray * sdigits = gime->SDigits() ;
206 //Now make SDigits from hits, for PHOS it is the same, so just copy
207 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
208 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
209 // Assign primary number only if contribution is significant
211 if( hit->GetEnergy() > fPrimThreshold)
212 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
213 hit->GetEnergy() ,hit->GetTime()) ;
215 new((*sdigits)[nSdigits]) AliPHOSDigit(-1 ,hit->GetId(),
216 hit->GetEnergy() ,hit->GetTime()) ;
223 nSdigits = sdigits->GetEntriesFast() ;
225 fSDigitsInRun += nSdigits ;
226 sdigits->Expand(nSdigits) ;
228 for (i = 0 ; i < nSdigits ; i++) {
229 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
230 digit->SetIndexInList(i) ;
236 //First list of sdigits
238 Int_t bufferSize = 32000 ;
239 TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
241 sdigitsBranch->Fill() ;
243 gime->WriteSDigits("OVERWRITE");
247 gime->WriteSDigitizer("OVERWRITE");
248 //gObjectTable->Print() ;
250 if(strstr(option,"deb"))
251 PrintSDigits(option) ;
253 //memwatcher.Watch(ievent);
258 // gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
260 if(strstr(option,"tim")){
261 gBenchmark->Stop("PHOSSDigitizer");
262 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
263 gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
266 //TFile f("out.root","RECREATE");
267 //memwatcher.WriteToFile();
271 //__________________________________________________________________
272 void AliPHOSSDigitizer::Print(const Option_t *)const
274 // Prints parameters of SDigitizer
275 Info("Print", "\n------------------- %s -------------", GetName() ) ;
276 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
277 printf(" with digitization parameters A = %f\n", fA) ;
278 printf(" B = %f\n", fB) ;
279 printf(" Threshold for Primary assignment= %f\n", fPrimThreshold) ;
280 printf("---------------------------------------------------\n") ;
284 //__________________________________________________________________
285 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
288 // SDititizers are equal if their pedestal, slope and threshold are equal
290 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
296 //__________________________________________________________________
297 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
299 // Prints list of digits produced in the current pass of AliPHOSDigitizer
302 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
303 const TClonesArray * sdigits = gime->SDigits() ;
305 Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ;
307 if(strstr(option,"all")||strstr(option,"EMC")){
310 AliPHOSDigit * digit;
311 printf("\nEMC sdigits\n") ;
312 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
314 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
315 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
316 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
317 // if(digit->GetNprimary() == 0)
319 // printf("%6d %8d %6.5e %4d %2d :\n", // YVK
320 printf("%6d %.4f %6.5e %4d %2d :\n",
321 digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
323 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
324 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
329 if(strstr(option,"all")||strstr(option,"CPV")){
331 //loop over CPV digits
332 AliPHOSDigit * digit;
333 printf("\nCPV sdigits\n") ;
334 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
336 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
337 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
338 if(digit->GetId() > maxEmc){
339 printf("\n%6d %8d %4d %2d :",
340 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
342 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
343 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
350 //____________________________________________________________________________
351 void AliPHOSSDigitizer::Unload() const
353 // Unloads the objects from the folder
354 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
355 AliPHOSLoader * loader = gime->PhosLoader() ;
356 loader->UnloadHits() ;
357 loader->UnloadSDigits() ;