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 fPrimThreshold = 0.01 ;
90 //____________________________________________________________________________
91 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
96 fPrimThreshold = 0.01 ;
102 //____________________________________________________________________________
103 AliPHOSSDigitizer::~AliPHOSSDigitizer()
106 if ( fSplitFile->IsOpen() )
107 fSplitFile->Close() ;
108 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
109 // Close the root file
112 // remove the task from the folder list
113 gime->RemoveTask("S",GetName()) ;
115 TString name(GetName()) ;
116 name.Remove(name.Index(":")) ;
118 // remove the Hits from the folder list
119 gime->RemoveObjects("H",name) ;
121 // remove the SDigits from the folder list
122 gime->RemoveObjects("S", name) ;
125 //____________________________________________________________________________
126 void AliPHOSSDigitizer::Init()
128 // Initialization: open root-file, allocate arrays for hits and sdigits,
129 // attach task SDigitizer to the list of PHOS tasks
131 // Initialization can not be done in the default constructor
132 //============================================================= YS
133 // The initialisation is now done by AliPHOSGetter
135 if( strcmp(GetTitle(), "") == 0 )
136 SetTitle("galice.root") ;
138 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
140 cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
144 gime->PostSDigits( GetName(), GetTitle() ) ;
146 TString sdname(GetName() );
148 sdname.Append(GetTitle() ) ;
150 gime->PostSDigitizer(this) ;
153 //____________________________________________________________________________
154 void AliPHOSSDigitizer::Exec(Option_t *option)
156 // Collects all hits in the same active volume into digit
157 if( strcmp(GetName(), "") == 0 )
160 if (strstr(option, "print") ) {
165 if(strstr(option,"tim"))
166 gBenchmark->Start("PHOSSDigitizer");
168 //Check, if this branch already exits
169 gAlice->GetEvent(0) ;
171 TString sdname(GetName()) ;
172 sdname.Remove(sdname.Index(GetTitle())-1) ;
174 if(gAlice->TreeS() ) {
175 TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
177 TBranch * branch = 0 ;
178 Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ;
180 while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
181 TString thisName( GetName() ) ;
182 TString branchName( branch->GetTitle() ) ;
183 branchName.Append(":") ;
184 if ( (strcmp(branch->GetName(), "PHOS")==0) && thisName.BeginsWith(branchName) )
186 else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && thisName.BeginsWith(branchName) )
187 sdigitizerfound = kTRUE ;
190 if ( phosfound || sdigitizerfound ) {
191 cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName()
192 << " already exits" << endl ;
197 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
199 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
202 for(ievent = 0; ievent < nevents; ievent++){
203 gime->Event(ievent,"H") ;
204 const TClonesArray * hits = gime->Hits() ;
205 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
210 //Now make SDigits from hits, for PHOS it is the same, so just copy
212 Int_t Nprim = (Int_t) (gAlice->TreeH())->GetEntries();
213 // Attention Nprim is the number of primaries tracked by Geant and this number could be different to the number of Primaries in TreeK;
215 for (iprim=0; iprim<Nprim; iprim++) {
216 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
219 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
220 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
221 // Assign primary number only if contribution is significant
223 if( hit->GetEnergy() > fPrimThreshold)
224 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
225 Digitize(hit->GetEnergy()), hit->GetTime()) ;
227 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
228 Digitize(hit->GetEnergy()), hit->GetTime()) ;
237 nSdigits = sdigits->GetEntriesFast() ;
238 fSDigitsInRun += nSdigits ;
239 sdigits->Expand(nSdigits) ;
242 for (i = 0 ; i < nSdigits ; i++) {
243 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
244 digit->SetIndexInList(i) ;
247 if(gAlice->TreeS() == 0)
248 gAlice->MakeTree("S", fSplitFile);
250 //First list of sdigits
251 Int_t bufferSize = 32000 ;
252 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
253 sdigitsBranch->SetTitle(sdname);
256 Int_t splitlevel = 0 ;
257 AliPHOSSDigitizer * sd = this ;
258 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
259 &sd,bufferSize,splitlevel);
260 sdigitizerBranch->SetTitle(sdname);
262 sdigitsBranch->Fill() ;
263 sdigitizerBranch->Fill() ;
264 gAlice->TreeS()->AutoSave() ;
266 if(strstr(option,"deb"))
267 PrintSDigits(option) ;
272 if ( fSplitFile->IsOpen() )
273 fSplitFile->Close() ;
276 if(strstr(option,"tim")){
277 gBenchmark->Stop("PHOSSDigitizer");
278 cout << "AliPHOSSDigitizer:" << endl ;
279 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
280 << gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
287 //__________________________________________________________________
288 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
290 // Setting title to branch SDigits
292 TString stitle(title) ;
294 // check if branch with title already exists
295 TBranch * sdigitsBranch =
296 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ;
297 TBranch * sdigitizerBranch =
298 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
299 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
300 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
301 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
302 cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
306 cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
310 // Post to the WhiteBoard
311 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
312 gime->PostSDigits( title, GetTitle()) ;
315 //__________________________________________________________________
316 void AliPHOSSDigitizer::SetSplitFile(const TString splitFileName)
318 // Diverts the SDigits in a file separate from the hits file
320 TDirectory * cwd = gDirectory ;
322 if ( !(gAlice->GetTreeSFileName() == splitFileName) ) {
323 if (gAlice->GetTreeSFile() )
324 gAlice->GetTreeSFile()->Close() ;
327 fSplitFile = gAlice->InitTreeFile("S",splitFileName.Data());
329 if ( !fSplitFile->Get("gAlice") )
332 TTree *treeE = gAlice->TreeE();
334 cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> No TreeE found "<<endl;
339 if ( !fSplitFile->Get("TreeE") ) {
340 AliHeader *header = new AliHeader();
341 treeE->SetBranchAddress("Header", &header);
342 treeE->SetBranchStatus("*",1);
343 TTree *treeENew = treeE->CloneTree();
348 if ( !fSplitFile->Get("AliceGeom") ) {
349 TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
351 cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> AliceGeom was not found in the input file "<<endl;
357 gAlice->MakeTree("S",fSplitFile);
359 cout << "INFO: AliPHOSSDigitizer::SetSPlitMode -> SDigits will be stored in " << splitFileName.Data() << endl ;
362 //__________________________________________________________________
363 void AliPHOSSDigitizer::Print(Option_t* option)const
365 // Prints parameters of SDigitizer
366 cout << "------------------- "<< GetName() << " -------------" << endl ;
367 cout << " Writing SDigits to branch with title " << GetName() << endl ;
368 cout << " with digitization parameters A = " << fA << endl ;
369 cout << " B = " << fB << endl ;
370 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
371 cout << "---------------------------------------------------"<<endl ;
374 //__________________________________________________________________
375 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
378 // SDititizers are equal if their pedestal, slope and threshold are equal
380 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
385 //__________________________________________________________________
386 //__________________________________________________________________
387 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
389 // Prints list of digits produced in the current pass of AliPHOSDigitizer
392 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
393 TString sdname(GetName()) ;
394 sdname.Remove(sdname.Index(GetTitle())-1) ;
395 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
397 cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
398 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
400 if(strstr(option,"all")||strstr(option,"EMC")){
403 AliPHOSDigit * digit;
404 cout << "EMC sdigits " << endl ;
405 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
406 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
408 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
409 (((AliPHOSDigit * ) sdigits->At(index))->GetId() <= maxEmc) ; index++) {
410 digit = (AliPHOSDigit * ) sdigits->At(index) ;
411 if(digit->GetNprimary() == 0) continue;
412 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
413 << setw(6) << digit->GetIndexInList() << " "
414 << setw(5) << digit->GetNprimary() <<" ";
417 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
418 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
424 if(strstr(option,"all")||strstr(option,"CPV")){
426 //loop over CPV digits
427 AliPHOSDigit * digit;
428 cout << "CPV sdigits " << endl ;
429 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
430 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
432 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
433 digit = (AliPHOSDigit * ) sdigits->At(index) ;
434 if(digit->GetId() > maxEmc){
435 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
436 << setw(6) << digit->GetIndexInList() << " "
437 << setw(5) << digit->GetNprimary() <<" ";
440 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
441 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
449 //____________________________________________________________________________
450 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)