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 **************************************************************************/
18 //_________________________________________________________________________
19 // This is a TTask that makes SDigits out of Hits
20 // The name of the TTask is also the title of the branch that will contain
21 // the created SDigits
22 // The title of the TTAsk is the name of the file that contains the hits from
23 // which the SDigits are created
24 // A Summable Digits is the sum of all hits originating
25 // from one primary in one active cell
26 // A threshold for assignment of the primary to SDigit is applied
27 // SDigits are written to TreeS, branch "PHOS"
28 // AliPHOSSDigitizer with all current parameters is written
29 // to TreeS branch "AliPHOSSDigitizer".
30 // Both branches have the same title. If necessary one can produce
31 // another set of SDigits with different parameters. Two versions
32 // can be distunguished using titles of the branches.
34 // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
35 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
36 // root [1] s->ExecuteTask()
37 // // Makes SDigitis for all events stored in galice.root
38 // root [2] s->SetPedestalParameter(0.001)
39 // // One can change parameters of digitization
40 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
41 // // and write them into the new branch
42 // root [4] s->ExecuteTask("deb all tim")
43 // // available parameters:
44 // deb - print # of produced SDigitis
45 // deb all - print # and list of produced SDigits
46 // tim - print benchmarking information
48 //*-- Author : Dmitri Peressounko (SUBATECH & KI)
49 //////////////////////////////////////////////////////////////////////////////
52 // --- ROOT system ---
59 #include "TBenchmark.h"
60 // --- Standard library ---
63 // --- AliRoot header files ---
65 #include "AliPHOSDigit.h"
66 #include "AliPHOSGeometry.h"
67 #include "AliPHOSGetter.h"
68 #include "AliPHOSHit.h"
69 #include "AliPHOSSDigitizer.h"
72 ClassImp(AliPHOSSDigitizer)
75 //____________________________________________________________________________
76 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","")
81 fPrimThreshold = 0.01 ;
85 //____________________________________________________________________________
86 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
91 fPrimThreshold = 0.01 ;
97 //____________________________________________________________________________
98 void AliPHOSSDigitizer::Init()
100 // Initialization: open root-file, allocate arrays for hits and sdigits,
101 // attach task SDigitizer to the list of PHOS tasks
103 // Initialization can not be done in the default constructor
104 //============================================================= YS
105 // The initialisation is now done by AliPHOSGetter
107 if( strcmp(GetTitle(), "") == 0 )
108 SetTitle("galice.root") ;
110 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
112 cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
116 gime->PostSDigits( GetName(), GetTitle() ) ;
118 TString sdname(GetName() );
120 sdname.Append(GetTitle() ) ;
122 gime->PostSDigitizer(this) ;
124 // create a folder on the white board //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/sdigitsTitle
128 //____________________________________________________________________________
129 void AliPHOSSDigitizer::Exec(Option_t *option)
131 // Collects all hits in the same active volume into digit
132 if( strcmp(GetName(), "") == 0 )
135 if (strstr(option, "print") ) {
140 if(strstr(option,"tim"))
141 gBenchmark->Start("PHOSSDigitizer");
143 //Check, if this branch already exits
144 gAlice->GetEvent(0) ;
145 if(gAlice->TreeS() ) {
146 TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
148 TBranch * branch = 0 ;
149 Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ;
151 while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
152 if ( (strcmp(branch->GetName(), "PHOS")==0) && (strcmp(branch->GetTitle(), GetName())==0) )
155 else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) )
156 sdigitizerfound = kTRUE ;
159 if ( phosfound || sdigitizerfound ) {
160 cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName()
161 << " already exits" << endl ;
166 TString sdname(GetName()) ;
167 sdname.Remove(sdname.Index(GetTitle())-1) ;
169 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
171 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
174 for(ievent = 0; ievent < nevents; ievent++){
175 gime->Event(ievent,"H") ;
177 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
182 //Now make SDigits from hits, for PHOS it is the same, so just copy
184 for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
186 //=========== Get the PHOS branch from Hits Tree for the Primary track itrack
187 gime->Track(itrack) ;
188 TClonesArray * hits = gime->Hits() ;
190 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
191 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
192 // Assign primary number only if contribution is significant
194 if( hit->GetEnergy() > fPrimThreshold)
195 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
196 Digitize(hit->GetEnergy()), hit->GetTime()) ;
198 new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(),
199 Digitize(hit->GetEnergy()), hit->GetTime()) ;
203 } // loop over tracks
207 nSdigits = sdigits->GetEntriesFast() ;
208 fSDigitsInRun += nSdigits ;
209 sdigits->Expand(nSdigits) ;
211 for (i = 0 ; i < nSdigits ; i++) {
212 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
213 digit->SetIndexInList(i) ;
216 if(gAlice->TreeS() == 0)
217 gAlice->MakeTree("S") ;
219 //Make (if necessary) branches
221 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
222 file = new char[strlen(gAlice->GetBaseFile())+20] ;
223 sprintf(file,"%s/PHOS.SDigits.root",gAlice->GetBaseFile()) ;
226 TDirectory *cwd = gDirectory;
228 //First list of sdigits
229 Int_t bufferSize = 32000 ;
230 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
231 sdigitsBranch->SetTitle(sdname);
233 sdigitsBranch->SetFile(file);
234 TIter next( sdigitsBranch->GetListOfBranches());
236 while ((subbr=static_cast<TBranch*>(next()))) {
237 subbr->SetFile(file);
243 Int_t splitlevel = 0 ;
244 AliPHOSSDigitizer * sd = this ;
245 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
246 &sd,bufferSize,splitlevel);
247 sdigitizerBranch->SetTitle(sdname);
249 sdigitizerBranch->SetFile(file);
250 TIter next( sdigitizerBranch->GetListOfBranches());
252 while ((subbr=static_cast<TBranch*>(next()))) {
253 subbr->SetFile(file);
259 sdigitsBranch->Fill() ;
260 sdigitizerBranch->Fill() ;
261 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
263 if(strstr(option,"deb"))
264 PrintSDigits(option) ;
268 if(strstr(option,"tim")){
269 gBenchmark->Stop("PHOSSDigitizer");
270 cout << "AliPHOSSDigitizer:" << endl ;
271 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
272 << gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
278 //__________________________________________________________________
279 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
281 // Setting title to branch SDigits
283 TString stitle(title) ;
285 // check if branch with title already exists
286 TBranch * sdigitsBranch =
287 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ;
288 TBranch * sdigitizerBranch =
289 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
290 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
291 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
292 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
293 cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
297 cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
301 // Post to the WhiteBoard
302 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
303 gime->PostSDigits( title, GetTitle()) ;
306 //__________________________________________________________________
307 void AliPHOSSDigitizer::Print(Option_t* option)const
309 // Prints parameters of SDigitizer
310 cout << "------------------- "<< GetName() << " -------------" << endl ;
311 cout << " Writing SDigits to branch with title " << GetName() << endl ;
312 cout << " with digitization parameters A = " << fA << endl ;
313 cout << " B = " << fB << endl ;
314 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
315 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))
329 //__________________________________________________________________
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 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 (((AliPHOSDigit * ) sdigits->At(index))->GetId() <= maxEmc) ; index++) {
354 digit = (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 = (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)