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 ---
60 #include "TBenchmark.h"
61 #include "TGeometry.h"
63 // --- Standard library ---
66 // --- AliRoot header files ---
68 #include "AliHeader.h"
69 #include "AliPHOSDigit.h"
70 #include "AliPHOSGeometry.h"
71 #include "AliPHOSGetter.h"
72 #include "AliPHOSHit.h"
73 #include "AliPHOSSDigitizer.h"
76 ClassImp(AliPHOSSDigitizer)
79 //____________________________________________________________________________
80 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") {
85 //____________________________________________________________________________
86 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
93 //____________________________________________________________________________
94 AliPHOSSDigitizer::~AliPHOSSDigitizer()
97 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
99 // remove the task from the folder list
100 gime->RemoveTask("S",GetName()) ;
102 TString name(GetName()) ;
103 if (! name.IsNull() )
104 name.Remove(name.Index(":")) ;
106 // remove the Hits from the folder list
107 gime->RemoveObjects("H",name) ;
109 // remove the SDigits from the folder list
110 gime->RemoveObjects("S", name) ;
119 //____________________________________________________________________________
120 void AliPHOSSDigitizer::Init()
122 // Initialization: open root-file, allocate arrays for hits and sdigits,
123 // attach task SDigitizer to the list of PHOS tasks
125 // Initialization can not be done in the default constructor
126 //============================================================= YS
127 // The initialisation is now done by AliPHOSGetter
129 if( strcmp(GetTitle(), "") == 0 )
130 SetTitle("galice.root") ;
132 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
134 cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
138 gime->PostSDigits( GetName(), GetTitle() ) ;
140 TString sdname(GetName() );
142 sdname.Append(GetTitle() ) ;
144 gime->PostSDigitizer(this) ;
147 //____________________________________________________________________________
148 void AliPHOSSDigitizer::InitParameters()
152 fPrimThreshold = 0.01 ;
157 //____________________________________________________________________________
158 void AliPHOSSDigitizer::Exec(Option_t *option)
160 // Collects all hits in the same active volume into digit
161 if( strcmp(GetName(), "") == 0 )
164 if (strstr(option, "print") ) {
169 if(strstr(option,"tim"))
170 gBenchmark->Start("PHOSSDigitizer");
172 //Check, if this branch already exits
173 gAlice->GetEvent(0) ;
175 TString sdname(GetName()) ;
176 sdname.Remove(sdname.Index(GetTitle())-1) ;
178 if(gAlice->TreeS() ) {
179 TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
181 TBranch * branch = 0 ;
182 Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ;
184 while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
185 TString thisName( GetName() ) ;
186 TString branchName( branch->GetTitle() ) ;
187 branchName.Append(":") ;
188 if ( (strcmp(branch->GetName(), "PHOS")==0) && thisName.BeginsWith(branchName) )
190 else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && thisName.BeginsWith(branchName) )
191 sdigitizerfound = kTRUE ;
194 if ( phosfound || sdigitizerfound ) {
195 cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName()
196 << " already exits" << endl ;
201 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
203 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
206 for(ievent = 0; ievent < nevents; ievent++){
207 gime->Event(ievent,"H") ;
208 const TClonesArray * hits = gime->Hits() ;
209 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
214 //Now make SDigits from hits, for PHOS it is the same, so just copy
216 Int_t Nprim = (Int_t) (gAlice->TreeH())->GetEntries();
217 // Attention Nprim is the number of primaries tracked by Geant and this number could be different to the number of Primaries in TreeK;
219 for (iprim=0; iprim<Nprim; iprim++) {
220 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
223 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
224 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
225 // Assign primary number only if contribution is significant
227 if( hit->GetEnergy() > fPrimThreshold)
228 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
229 Digitize(hit->GetEnergy()), hit->GetTime()) ;
231 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
232 Digitize(hit->GetEnergy()), hit->GetTime()) ;
241 nSdigits = sdigits->GetEntriesFast() ;
242 fSDigitsInRun += nSdigits ;
243 sdigits->Expand(nSdigits) ;
246 for (i = 0 ; i < nSdigits ; i++) {
247 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
248 digit->SetIndexInList(i) ;
251 if(gAlice->TreeS() == 0)
252 gAlice->MakeTree("S", fSplitFile);
254 //First list of sdigits
255 Int_t bufferSize = 32000 ;
256 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
257 sdigitsBranch->SetTitle(sdname);
260 Int_t splitlevel = 0 ;
261 AliPHOSSDigitizer * sd = this ;
262 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
263 &sd,bufferSize,splitlevel);
264 sdigitizerBranch->SetTitle(sdname);
266 sdigitsBranch->Fill() ;
267 sdigitizerBranch->Fill() ;
268 gAlice->TreeS()->AutoSave() ;
270 if(strstr(option,"deb"))
271 PrintSDigits(option) ;
275 if(strstr(option,"tim")){
276 gBenchmark->Stop("PHOSSDigitizer");
277 cout << "AliPHOSSDigitizer:" << endl ;
278 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
279 << gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
286 //__________________________________________________________________
287 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
289 // Setting title to branch SDigits
291 TString stitle(title) ;
293 // check if branch with title already exists
294 TBranch * sdigitsBranch =
295 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ;
296 TBranch * sdigitizerBranch =
297 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
298 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
299 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
300 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
301 cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
305 cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
309 // Post to the WhiteBoard
310 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
311 gime->PostSDigits( title, GetTitle()) ;
314 //__________________________________________________________________
315 void AliPHOSSDigitizer::SetSplitFile(const TString splitFileName)
317 // Diverts the SDigits in a file separate from the hits file
319 TDirectory * cwd = gDirectory ;
321 if ( !(gAlice->GetTreeSFileName() == splitFileName) ) {
322 if (gAlice->GetTreeSFile() )
323 gAlice->GetTreeSFile()->Close() ;
326 fSplitFile = gAlice->InitTreeFile("S",splitFileName.Data());
328 gAlice->Write(0, TObject::kOverwrite);
330 TTree *treeE = gAlice->TreeE();
332 cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> No TreeE found "<<endl;
337 AliHeader *header = new AliHeader();
338 treeE->SetBranchAddress("Header", &header);
339 treeE->SetBranchStatus("*",1);
340 TTree *treeENew = treeE->CloneTree();
341 treeENew->Write(0, TObject::kOverwrite);
345 TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
347 cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> AliceGeom was not found in the input file "<<endl;
350 AliceGeom->Write(0, TObject::kOverwrite);
352 gAlice->MakeTree("S",fSplitFile);
354 cout << "INFO: AliPHOSSDigitizer::SetSPlitMode -> SDigits will be stored in " << splitFileName.Data() << endl ;
358 //__________________________________________________________________
359 void AliPHOSSDigitizer::Print(Option_t* option)const
361 // Prints parameters of SDigitizer
362 cout << "------------------- "<< GetName() << " -------------" << endl ;
363 cout << " Writing SDigits to branch with title " << GetName() << endl ;
364 cout << " with digitization parameters A = " << fA << endl ;
365 cout << " B = " << fB << endl ;
366 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
367 cout << "---------------------------------------------------"<<endl ;
370 //__________________________________________________________________
371 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
374 // SDititizers are equal if their pedestal, slope and threshold are equal
376 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
381 //__________________________________________________________________
382 //__________________________________________________________________
383 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
385 // Prints list of digits produced in the current pass of AliPHOSDigitizer
388 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
389 TString sdname(GetName()) ;
390 sdname.Remove(sdname.Index(GetTitle())-1) ;
391 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
393 cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
394 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
396 if(strstr(option,"all")||strstr(option,"EMC")){
399 AliPHOSDigit * digit;
400 cout << "EMC sdigits " << endl ;
401 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
402 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
404 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
405 (((AliPHOSDigit * ) sdigits->At(index))->GetId() <= maxEmc) ; index++) {
406 digit = (AliPHOSDigit * ) sdigits->At(index) ;
407 if(digit->GetNprimary() == 0) continue;
408 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
409 << setw(6) << digit->GetIndexInList() << " "
410 << setw(5) << digit->GetNprimary() <<" ";
413 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
414 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
420 if(strstr(option,"all")||strstr(option,"CPV")){
422 //loop over CPV digits
423 AliPHOSDigit * digit;
424 cout << "CPV sdigits " << endl ;
425 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
426 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
428 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
429 digit = (AliPHOSDigit * ) sdigits->At(index) ;
430 if(digit->GetId() > maxEmc){
431 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
432 << setw(6) << digit->GetIndexInList() << " "
433 << setw(5) << digit->GetNprimary() <<" ";
436 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
437 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
445 //____________________________________________________________________________
446 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)