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, const char * eventFolderName):
79 TTask("PHOS"+AliConfig::fgkSDigitizerTaskName, alirunFileName),
80 fEventFolderName(eventFolderName)
84 fFirstEvent = fLastEvent = 0 ; // runs one event by defaut
87 fDefaultInit = kFALSE ;
90 //____________________________________________________________________________
91 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer & sd)
96 fFirstEvent = sd.fFirstEvent ;
97 fLastEvent = sd.fLastEvent ;
100 fPrimThreshold = sd.fPrimThreshold ;
101 fSDigitsInRun = sd.fSDigitsInRun ;
102 SetName(sd.GetName()) ;
103 SetTitle(sd.GetTitle()) ;
104 fEventFolderName = sd.fEventFolderName;
108 //____________________________________________________________________________
109 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
111 AliPHOSGetter * gime =
112 AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());
113 gime->PhosLoader()->CleanSDigitizer();
115 //____________________________________________________________________________
116 void AliPHOSSDigitizer::Init()
118 // Uses the getter to access the required files
122 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
124 Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
128 TString opt("SDigits") ;
129 if(gime->VersionExists(opt) ) {
130 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
134 gime->PostSDigitizer(this);
135 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
139 //____________________________________________________________________________
140 void AliPHOSSDigitizer::InitParameters()
142 // initializes the parameters for digitization
145 fPrimThreshold = 0.01 ;
149 //____________________________________________________________________________
150 void AliPHOSSDigitizer::Exec(Option_t *option)
152 // Steering method to produce summable digits for events
153 // in the range from fFirstEvent to fLastEvent.
154 // This range is optionally set by SetEventRange().
155 // if fLastEvent=-1 (by default), then process events until the end.
157 // Summable digit is a sum of all hits in the same active
160 if (strstr(option, "print") ) {
165 if(strstr(option,"tim"))
166 gBenchmark->Start("PHOSSDigitizer");
168 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
170 //switch off reloading of this task while getting event
171 if (!fInit) { // to prevent overwrite existing file
172 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
176 if (fLastEvent == -1)
177 fLastEvent = gime->MaxEvent() - 1 ;
179 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time
180 Int_t nEvents = fLastEvent - fFirstEvent + 1;
183 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
184 gime->Event(ievent,"H") ;
185 TTree * treeS = gime->TreeS();
186 TClonesArray * hits = gime->Hits() ;
187 TClonesArray * sdigits = gime->SDigits() ;
190 //Now make SDigits from hits, for PHOS it is the same, so just copy
191 Int_t nPrim = static_cast<Int_t>((gime->TreeH())->GetEntries()) ;
192 // Attention nPrim is the number of primaries tracked by Geant
193 // and this number could be different to the number of Primaries in TreeK;
196 for (iprim = 0 ; iprim < nPrim ; iprim ++) {
197 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
200 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
201 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
202 // Assign primary number only if contribution is significant
204 if( hit->GetEnergy() > fPrimThreshold)
205 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
206 Digitize(hit->GetEnergy()), hit->GetTime()) ;
208 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
209 Digitize(hit->GetEnergy()), hit->GetTime()) ;
218 nSdigits = sdigits->GetEntriesFast() ;
220 fSDigitsInRun += nSdigits ;
221 sdigits->Expand(nSdigits) ;
224 for (i = 0 ; i < nSdigits ; i++) {
225 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
226 digit->SetIndexInList(i) ;
232 //First list of sdigits
234 Int_t bufferSize = 32000 ;
235 TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
237 sdigitsBranch->Fill() ;
239 gime->WriteSDigits("OVERWRITE");
243 gime->WriteSDigitizer("OVERWRITE");
245 if(strstr(option,"deb"))
246 PrintSDigits(option) ;
251 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
253 if(strstr(option,"tim")){
254 gBenchmark->Stop("PHOSSDigitizer");
255 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
256 gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
260 //__________________________________________________________________
261 void AliPHOSSDigitizer::Print()const
263 // Prints parameters of SDigitizer
264 Info("Print", "\n------------------- %s -------------", GetName() ) ;
265 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
266 printf(" with digitization parameters A = %f\n", fA) ;
267 printf(" B = %f\n", fB) ;
268 printf(" Threshold for Primary assignment= %f\n", fPrimThreshold) ;
269 printf("---------------------------------------------------\n") ;
273 //__________________________________________________________________
274 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
277 // SDititizers are equal if their pedestal, slope and threshold are equal
279 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
285 //__________________________________________________________________
286 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
288 // Prints list of digits produced in the current pass of AliPHOSDigitizer
291 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
292 const TClonesArray * sdigits = gime->SDigits() ;
294 Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ;
296 if(strstr(option,"all")||strstr(option,"EMC")){
299 AliPHOSDigit * digit;
300 printf("\nEMC sdigits\n") ;
301 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
303 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
304 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
305 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
306 // if(digit->GetNprimary() == 0)
308 printf("%6d %8d %6.5e %4d %2d :\n",
309 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
311 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
312 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
317 if(strstr(option,"all")||strstr(option,"CPV")){
319 //loop over CPV digits
320 AliPHOSDigit * digit;
321 printf("\nCPV sdigits\n") ;
322 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
324 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
325 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
326 if(digit->GetId() > maxEmc){
327 printf("\n%6d %8d %4d %2d :",
328 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
330 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
331 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
338 //____________________________________________________________________________
339 void AliPHOSSDigitizer::Unload() const
341 // Unloads the objects from the folder
342 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
343 AliPHOSLoader * loader = gime->PhosLoader() ;
344 loader->UnloadHits() ;
345 loader->UnloadSDigits() ;