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("","") {
83 fDefaultInit = kTRUE ;
86 //____________________________________________________________________________
87 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
92 fDefaultInit = kFALSE ;
95 //____________________________________________________________________________
96 AliPHOSSDigitizer::~AliPHOSSDigitizer()
99 // fDefaultInit = kTRUE if SDigitizer created by default ctor (to get just the parameters)
103 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
105 // remove the task from the folder list
106 gime->RemoveTask("S",GetName()) ;
108 TString name(GetName()) ;
109 if (! name.IsNull() )
110 name.Remove(name.Index(":")) ;
112 // remove the Hits from the folder list
113 gime->RemoveObjects("H",name) ;
115 // remove the SDigits from the folder list
116 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::InitParameters()
158 fPrimThreshold = 0.01 ;
163 //____________________________________________________________________________
164 void AliPHOSSDigitizer::Exec(Option_t *option)
166 // Collects all hits in the same active volume into digit
167 if( strcmp(GetName(), "") == 0 )
170 if (strstr(option, "print") ) {
175 if(strstr(option,"tim"))
176 gBenchmark->Start("PHOSSDigitizer");
178 //Check, if this branch already exits
179 gAlice->GetEvent(0) ;
181 TString sdname(GetName()) ;
182 sdname.Remove(sdname.Index(GetTitle())-1) ;
184 if(gAlice->TreeS() ) {
185 TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
187 TBranch * branch = 0 ;
188 Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ;
190 while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
191 TString thisName( GetName() ) ;
192 TString branchName( branch->GetTitle() ) ;
193 branchName.Append(":") ;
194 if ( (strcmp(branch->GetName(), "PHOS")==0) && thisName.BeginsWith(branchName) )
196 else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && thisName.BeginsWith(branchName) )
197 sdigitizerfound = kTRUE ;
200 if ( phosfound || sdigitizerfound ) {
201 cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName()
202 << " already exits" << endl ;
207 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
209 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
212 for(ievent = 0; ievent < nevents; ievent++){
213 gime->Event(ievent,"H") ;
214 const TClonesArray * hits = gime->Hits() ;
215 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
220 //Now make SDigits from hits, for PHOS it is the same, so just copy
222 Int_t Nprim = (Int_t) (gAlice->TreeH())->GetEntries();
223 // Attention Nprim is the number of primaries tracked by Geant and this number could be different to the number of Primaries in TreeK;
225 for (iprim=0; iprim<Nprim; iprim++) {
226 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
229 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
230 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
231 // Assign primary number only if contribution is significant
233 if( hit->GetEnergy() > fPrimThreshold)
234 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
235 Digitize(hit->GetEnergy()), hit->GetTime()) ;
237 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
238 Digitize(hit->GetEnergy()), hit->GetTime()) ;
247 nSdigits = sdigits->GetEntriesFast() ;
248 fSDigitsInRun += nSdigits ;
249 sdigits->Expand(nSdigits) ;
252 for (i = 0 ; i < nSdigits ; i++) {
253 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
254 digit->SetIndexInList(i) ;
257 if(gAlice->TreeS() == 0)
258 gAlice->MakeTree("S", fSplitFile);
260 //First list of sdigits
261 Int_t bufferSize = 32000 ;
262 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
263 sdigitsBranch->SetTitle(sdname);
266 Int_t splitlevel = 0 ;
267 AliPHOSSDigitizer * sd = this ;
268 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
269 &sd,bufferSize,splitlevel);
270 sdigitizerBranch->SetTitle(sdname);
272 sdigitsBranch->Fill() ;
273 sdigitizerBranch->Fill() ;
274 gAlice->TreeS()->AutoSave() ;
276 if(strstr(option,"deb"))
277 PrintSDigits(option) ;
281 if(strstr(option,"tim")){
282 gBenchmark->Stop("PHOSSDigitizer");
283 cout << "AliPHOSSDigitizer:" << endl ;
284 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
285 << gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
292 //__________________________________________________________________
293 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
295 // Setting title to branch SDigits
297 TString stitle(title) ;
299 // check if branch with title already exists
300 TBranch * sdigitsBranch =
301 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ;
302 TBranch * sdigitizerBranch =
303 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
304 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
305 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
306 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
307 cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
311 cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
315 // Post to the WhiteBoard
316 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
317 gime->PostSDigits( title, GetTitle()) ;
320 //__________________________________________________________________
321 void AliPHOSSDigitizer::SetSplitFile(const TString splitFileName)
323 // Diverts the SDigits in a file separate from the hits file
325 TDirectory * cwd = gDirectory ;
327 if ( !(gAlice->GetTreeSFileName() == splitFileName) ) {
328 if (gAlice->GetTreeSFile() )
329 gAlice->GetTreeSFile()->Close() ;
332 fSplitFile = gAlice->InitTreeFile("S",splitFileName.Data());
334 gAlice->Write(0, TObject::kOverwrite);
336 TTree *treeE = gAlice->TreeE();
338 cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> No TreeE found "<<endl;
343 AliHeader *header = new AliHeader();
344 treeE->SetBranchAddress("Header", &header);
345 treeE->SetBranchStatus("*",1);
346 TTree *treeENew = treeE->CloneTree();
347 treeENew->Write(0, TObject::kOverwrite);
351 TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
353 cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> AliceGeom was not found in the input file "<<endl;
356 AliceGeom->Write(0, TObject::kOverwrite);
358 gAlice->MakeTree("S",fSplitFile);
360 cout << "INFO: AliPHOSSDigitizer::SetSPlitMode -> SDigits will be stored in " << splitFileName.Data() << endl ;
364 //__________________________________________________________________
365 void AliPHOSSDigitizer::Print(Option_t* option)const
367 // Prints parameters of SDigitizer
368 cout << "------------------- "<< GetName() << " -------------" << endl ;
369 cout << " Writing SDigits to branch with title " << GetName() << endl ;
370 cout << " with digitization parameters A = " << fA << endl ;
371 cout << " B = " << fB << endl ;
372 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
373 cout << "---------------------------------------------------"<<endl ;
376 //__________________________________________________________________
377 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
380 // SDititizers are equal if their pedestal, slope and threshold are equal
382 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
387 //__________________________________________________________________
388 //__________________________________________________________________
389 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
391 // Prints list of digits produced in the current pass of AliPHOSDigitizer
394 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
395 TString sdname(GetName()) ;
396 sdname.Remove(sdname.Index(GetTitle())-1) ;
397 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
399 cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
400 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
402 if(strstr(option,"all")||strstr(option,"EMC")){
405 AliPHOSDigit * digit;
406 cout << "EMC sdigits " << endl ;
407 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
408 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
410 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
411 (((AliPHOSDigit * ) sdigits->At(index))->GetId() <= maxEmc) ; index++) {
412 digit = (AliPHOSDigit * ) sdigits->At(index) ;
413 if(digit->GetNprimary() == 0) continue;
414 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
415 << setw(6) << digit->GetIndexInList() << " "
416 << setw(5) << digit->GetNprimary() <<" ";
419 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
420 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
426 if(strstr(option,"all")||strstr(option,"CPV")){
428 //loop over CPV digits
429 AliPHOSDigit * digit;
430 cout << "CPV sdigits " << endl ;
431 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
432 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
434 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
435 digit = (AliPHOSDigit * ) sdigits->At(index) ;
436 if(digit->GetId() > maxEmc){
437 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
438 << setw(6) << digit->GetIndexInList() << " "
439 << setw(5) << digit->GetNprimary() <<" ";
442 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
443 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
451 //____________________________________________________________________________
452 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)