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("","")
83 fDefaultInit = kTRUE ;
86 //____________________________________________________________________________
87 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle, const Bool_t toSplit):
88 TTask(sDigitsTitle, headerFile)
94 fDefaultInit = kFALSE ;
97 //____________________________________________________________________________
98 AliPHOSSDigitizer::~AliPHOSSDigitizer()
105 //____________________________________________________________________________
106 void AliPHOSSDigitizer::Init()
108 // Initialization: open root-file, allocate arrays for hits and sdigits,
109 // attach task SDigitizer to the list of PHOS tasks
111 // Initialization can not be done in the default constructor
112 //============================================================= YS
113 // The initialisation is now done by AliPHOSGetter
115 if( strcmp(GetTitle(), "") == 0 )
116 SetTitle("galice.root") ;
118 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName(),fToSplit) ;
120 cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !"
125 gime->PostSDigits( GetName(), GetTitle() ) ;
129 // construct the name of the file as /path/PHOS.SDigits.root
130 // First - extract full path if necessary
131 TString sDigitsFileName(GetTitle()) ;
132 Ssiz_t islash = sDigitsFileName.Last('/') ;
133 if(islash<sDigitsFileName.Length())
134 sDigitsFileName.Remove(islash+1,sDigitsFileName.Length()) ;
137 // Next - append the file name
138 sDigitsFileName+="PHOS.SDigits." ;
139 if((strcmp(GetName(),"Default")!=0)&&(strcmp(GetName(),"")!=0)){
140 sDigitsFileName+=GetName() ;
141 sDigitsFileName+="." ;
143 sDigitsFileName+="root" ;
144 // Finally - check if the file already opened or open the file
145 fSplitFile = static_cast<TFile*>(gROOT->GetFile(sDigitsFileName.Data()));
147 fSplitFile = TFile::Open(sDigitsFileName.Data(),"update") ;
150 TString sdname(GetName() );
152 sdname.Append(GetTitle() ) ;
154 gime->PostSDigitizer(this) ;
157 //____________________________________________________________________________
158 void AliPHOSSDigitizer::InitParameters()
162 fPrimThreshold = 0.01 ;
168 //____________________________________________________________________________
169 void AliPHOSSDigitizer::Exec(Option_t *option)
171 // Collects all hits in the same active volume into digit
173 if( strcmp(GetName(), "") == 0 )
176 if (strstr(option, "print") ) {
181 if(strstr(option,"tim"))
182 gBenchmark->Start("PHOSSDigitizer");
184 //Check, if this branch already exits
185 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
186 if(gime->BranchExists("SDigits") )
189 TString sdname(GetName()) ;
190 sdname.Remove(sdname.Index(GetTitle())-1) ;
192 Int_t nevents = gime->MaxEvent() ;
194 for(ievent = 0; ievent < nevents; ievent++){
195 gime->Event(ievent,"H") ;
196 const TClonesArray * hits = gime->Hits() ;
197 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
201 //Now make SDigits from hits, for PHOS it is the same, so just copy
202 Int_t nPrim = static_cast<Int_t>((gAlice->TreeH())->GetEntries()) ;
203 // Attention nPrim is the number of primaries tracked by Geant
204 // and this number could be different to the number of Primaries in TreeK;
206 for (iprim = 0 ; iprim < nPrim ; iprim ++) {
207 //=========== Get the PHOS branch from Hits Tree for the Primary iprim
210 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
211 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
212 // Assign primary number only if contribution is significant
214 if( hit->GetEnergy() > fPrimThreshold)
215 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
216 Digitize(hit->GetEnergy()), hit->GetTime()) ;
218 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
219 Digitize(hit->GetEnergy()), hit->GetTime()) ;
228 nSdigits = sdigits->GetEntriesFast() ;
229 fSDigitsInRun += nSdigits ;
230 sdigits->Expand(nSdigits) ;
233 for (i = 0 ; i < nSdigits ; i++) {
234 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
235 digit->SetIndexInList(i) ;
240 if((gAlice->TreeS() == 0)|| (fSplitFile))
241 gAlice->MakeTree("S", fSplitFile);
246 //First list of sdigits
247 Int_t bufferSize = 32000 ;
248 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
249 sdigitsBranch->SetTitle(sdname);
252 Int_t splitlevel = 0 ;
253 AliPHOSSDigitizer * sd = this ;
254 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
255 &sd,bufferSize,splitlevel);
256 sdigitizerBranch->SetTitle(sdname);
258 sdigitsBranch->Fill() ;
259 sdigitizerBranch->Fill() ;
261 gAlice->TreeS()->AutoSave() ;
263 if(strstr(option,"deb"))
264 PrintSDigits(option) ;
267 if(strstr(option,"tim")){
268 gBenchmark->Stop("PHOSSDigitizer");
269 cout << "AliPHOSSDigitizer:" << endl ;
270 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
271 << gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
276 //__________________________________________________________________
277 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
279 // Setting title to branch SDigits
281 TString stitle(title) ;
283 // check if branch with title already exists
284 TBranch * sdigitsBranch =
285 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ;
286 TBranch * sdigitizerBranch =
287 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
288 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
289 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
290 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
291 cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
295 cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
299 // Post to the WhiteBoard
300 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
301 gime->PostSDigits( title, GetTitle()) ;
305 //__________________________________________________________________
306 void AliPHOSSDigitizer::Print(Option_t* option)const
308 // Prints parameters of SDigitizer
309 cout << "------------------- "<< GetName() << " -------------" << endl ;
310 cout << " Writing SDigits to branch with title " << GetName() << endl ;
311 cout << " with digitization parameters A = " << fA << endl ;
312 cout << " B = " << fB << endl ;
313 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
314 cout << "---------------------------------------------------"<<endl ;
318 //__________________________________________________________________
319 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
322 // SDititizers are equal if their pedestal, slope and threshold are equal
324 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
330 //__________________________________________________________________
331 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
333 // Prints list of digits produced in the current pass of AliPHOSDigitizer
336 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
337 TString sdname(GetName()) ;
338 sdname.Remove(sdname.Index(GetTitle())-1) ;
339 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
341 cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
342 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
344 if(strstr(option,"all")||strstr(option,"EMC")){
347 AliPHOSDigit * digit;
348 cout << "EMC sdigits " << endl ;
349 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
350 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
352 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
353 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
354 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
355 if(digit->GetNprimary() == 0) continue;
356 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
357 << setw(6) << digit->GetIndexInList() << " "
358 << setw(5) << digit->GetNprimary() <<" ";
361 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
362 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
368 if(strstr(option,"all")||strstr(option,"CPV")){
370 //loop over CPV digits
371 AliPHOSDigit * digit;
372 cout << "CPV sdigits " << endl ;
373 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
374 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
376 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
377 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
378 if(digit->GetId() > maxEmc){
379 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
380 << setw(6) << digit->GetIndexInList() << " "
381 << setw(5) << digit->GetNprimary() <<" ";
384 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
385 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
393 //____________________________________________________________________________
394 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)