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"
75 ClassImp(AliPHOSSDigitizer)
78 //____________________________________________________________________________
79 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") {
82 fDefaultInit = kTRUE ;
85 //____________________________________________________________________________
86 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle, const Bool_t toSplit):
87 TTask(sDigitsTitle, headerFile)
93 fDefaultInit = kFALSE ;
96 //____________________________________________________________________________
97 AliPHOSSDigitizer::~AliPHOSSDigitizer()
104 //____________________________________________________________________________
105 void AliPHOSSDigitizer::Init()
107 // Initialization: open root-file, allocate arrays for hits and sdigits,
108 // attach task SDigitizer to the list of PHOS tasks
110 // Initialization can not be done in the default constructor
111 //============================================================= YS
112 // The initialisation is now done by AliPHOSGetter
114 if( strcmp(GetTitle(), "") == 0 )
115 SetTitle("galice.root") ;
117 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName(),fToSplit) ;
119 cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !"
124 gime->PostSDigits( GetName(), GetTitle() ) ;
128 // construct the name of the file as /path/PHOS.SDigits.root
129 // First - extract full path if necessary
130 TString sDigitsFileName(GetTitle()) ;
131 Ssiz_t islash = sDigitsFileName.Last('/') ;
132 if(islash<sDigitsFileName.Length())
133 sDigitsFileName.Remove(islash+1,sDigitsFileName.Length()) ;
136 // Next - append the file name
137 sDigitsFileName+="PHOS.SDigits." ;
138 if((strcmp(GetName(),"Default")!=0)&&(strcmp(GetName(),"")!=0)){
139 sDigitsFileName+=GetName() ;
140 sDigitsFileName+="." ;
142 sDigitsFileName+="root" ;
143 // Finally - check if the file already opened or open the file
144 fSplitFile = static_cast<TFile*>(gROOT->GetFile(sDigitsFileName.Data()));
146 fSplitFile = TFile::Open(sDigitsFileName.Data(),"update") ;
149 TString sdname(GetName() );
151 sdname.Append(GetTitle() ) ;
153 gime->PostSDigitizer(this) ;
156 //____________________________________________________________________________
157 void AliPHOSSDigitizer::InitParameters()
161 fPrimThreshold = 0.01 ;
167 //____________________________________________________________________________
168 void AliPHOSSDigitizer::Exec(Option_t *option)
170 // Collects all hits in the same active volume into digit
171 if( strcmp(GetName(), "") == 0 )
174 if (strstr(option, "print") ) {
179 if(strstr(option,"tim"))
180 gBenchmark->Start("PHOSSDigitizer");
182 //Check, if this branch already exits
183 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
184 if(gime->BranchExists("SDigits") )
187 TString sdname(GetName()) ;
188 sdname.Remove(sdname.Index(GetTitle())-1) ;
190 Int_t nevents = gime->MaxEvent() ;
192 for(ievent = 0; ievent < nevents; ievent++){
193 gime->Event(ievent,"H") ;
194 const TClonesArray * hits = gime->Hits() ;
195 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
199 //Now make SDigits from hits, for PHOS it is the same, so just copy
200 Int_t nPrim = static_cast<Int_t>((gAlice->TreeH())->GetEntries()) ;
201 // Attention nPrim is the number of primaries tracked by Geant
202 // and this number could be different to the number of Primaries in TreeK;
204 for (iprim = 0 ; iprim < nPrim ; iprim ++) {
205 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
208 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
209 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
210 // Assign primary number only if contribution is significant
212 if( hit->GetEnergy() > fPrimThreshold)
213 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
214 Digitize(hit->GetEnergy()), hit->GetTime()) ;
216 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
217 Digitize(hit->GetEnergy()), hit->GetTime()) ;
226 nSdigits = sdigits->GetEntriesFast() ;
227 fSDigitsInRun += nSdigits ;
228 sdigits->Expand(nSdigits) ;
231 for (i = 0 ; i < nSdigits ; i++) {
232 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
233 digit->SetIndexInList(i) ;
238 if((gAlice->TreeS() == 0)|| (fSplitFile))
239 gAlice->MakeTree("S", fSplitFile);
244 //First list of sdigits
245 Int_t bufferSize = 32000 ;
246 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
247 sdigitsBranch->SetTitle(sdname);
250 Int_t splitlevel = 0 ;
251 AliPHOSSDigitizer * sd = this ;
252 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
253 &sd,bufferSize,splitlevel);
254 sdigitizerBranch->SetTitle(sdname);
256 sdigitsBranch->Fill() ;
257 sdigitizerBranch->Fill() ;
259 gAlice->TreeS()->AutoSave() ;
261 if(strstr(option,"deb"))
262 PrintSDigits(option) ;
265 if(strstr(option,"tim")){
266 gBenchmark->Stop("PHOSSDigitizer");
267 cout << "AliPHOSSDigitizer:" << endl ;
268 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
269 << gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
274 //__________________________________________________________________
275 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
277 // Setting title to branch SDigits
279 TString stitle(title) ;
281 // check if branch with title already exists
282 TBranch * sdigitsBranch =
283 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ;
284 TBranch * sdigitizerBranch =
285 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
286 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
287 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
288 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
289 cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
293 cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
297 // Post to the WhiteBoard
298 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
299 gime->PostSDigits( title, GetTitle()) ;
302 // //__________________________________________________________________
303 // void AliPHOSSDigitizer::SetSplitFile(const TString splitFileName)
305 // // Diverts the SDigits in a file separate from the hits file
307 // TDirectory * cwd = gDirectory ;
309 // if ( !(gAlice->GetTreeSFileName() == splitFileName) ) {
310 // if (gAlice->GetTreeSFile() )
311 // gAlice->GetTreeSFile()->Close() ;
314 // fSplitFile = gAlice->InitTreeFile("S",splitFileName.Data());
315 // fSplitFile->cd() ;
316 // gAlice->Write(0, TObject::kOverwrite);
318 // TTree *treeE = gAlice->TreeE();
320 // cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> No TreeE found "<<endl;
325 // AliHeader *header = new AliHeader();
326 // treeE->SetBranchAddress("Header", &header);
327 // treeE->SetBranchStatus("*",1);
328 // TTree *treeENew = treeE->CloneTree();
329 // treeENew->Write(0, TObject::kOverwrite);
333 // TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
335 // cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> AliceGeom was not found in the input file "<<endl;
338 // AliceGeom->Write(0, TObject::kOverwrite);
340 // gAlice->MakeTree("S",fSplitFile);
342 // cout << "INFO: AliPHOSSDigitizer::SetSPlitMode -> SDigits will be stored in " << splitFileName.Data() << endl ;
346 //__________________________________________________________________
347 void AliPHOSSDigitizer::Print(Option_t* option)const
349 // Prints parameters of SDigitizer
350 cout << "------------------- "<< GetName() << " -------------" << endl ;
351 cout << " Writing SDigits to branch with title " << GetName() << endl ;
352 cout << " with digitization parameters A = " << fA << endl ;
353 cout << " B = " << fB << endl ;
354 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
355 cout << "---------------------------------------------------"<<endl ;
358 //__________________________________________________________________
359 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
362 // SDititizers are equal if their pedestal, slope and threshold are equal
364 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
369 //__________________________________________________________________
370 //__________________________________________________________________
371 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
373 // Prints list of digits produced in the current pass of AliPHOSDigitizer
376 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
377 TString sdname(GetName()) ;
378 sdname.Remove(sdname.Index(GetTitle())-1) ;
379 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
381 cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
382 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
384 if(strstr(option,"all")||strstr(option,"EMC")){
387 AliPHOSDigit * digit;
388 cout << "EMC sdigits " << endl ;
389 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
390 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
392 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
393 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
394 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
395 if(digit->GetNprimary() == 0) continue;
396 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
397 << setw(6) << digit->GetIndexInList() << " "
398 << setw(5) << digit->GetNprimary() <<" ";
401 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
402 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
408 if(strstr(option,"all")||strstr(option,"CPV")){
410 //loop over CPV digits
411 AliPHOSDigit * digit;
412 cout << "CPV sdigits " << endl ;
413 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
414 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
416 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
417 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
418 if(digit->GetId() > maxEmc){
419 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
420 << setw(6) << digit->GetIndexInList() << " "
421 << setw(5) << digit->GetNprimary() <<" ";
424 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
425 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
433 //____________________________________________________________________________
434 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)