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 ---
65 // --- AliRoot header files ---
67 #include "AliHeader.h"
68 #include "AliPHOSDigit.h"
69 #include "AliPHOSGeometry.h"
70 #include "AliPHOSGetter.h"
71 #include "AliPHOSHit.h"
72 #include "AliPHOSSDigitizer.h"
74 ClassImp(AliPHOSSDigitizer)
77 //____________________________________________________________________________
78 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 Error("Init" ,"Could not obtain the Getter object !") ;
123 gime->PostSDigits( GetName(), GetTitle() ) ;
127 // construct the name of the file as /path/PHOS.SDigits.root
128 // First - extract full path if necessary
129 TString sDigitsFileName(GetTitle()) ;
130 Ssiz_t islash = sDigitsFileName.Last('/') ;
131 if(islash<sDigitsFileName.Length())
132 sDigitsFileName.Remove(islash+1,sDigitsFileName.Length()) ;
135 // Next - append the file name
136 sDigitsFileName+="PHOS.SDigits." ;
137 if((strcmp(GetName(),"Default")!=0)&&(strcmp(GetName(),"")!=0)){
138 sDigitsFileName+=GetName() ;
139 sDigitsFileName+="." ;
141 sDigitsFileName+="root" ;
142 // Finally - check if the file already opened or open the file
143 fSplitFile = static_cast<TFile*>(gROOT->GetFile(sDigitsFileName.Data()));
145 fSplitFile = TFile::Open(sDigitsFileName.Data(),"update") ;
148 TString sdname(GetName() );
150 sdname.Append(GetTitle() ) ;
152 gime->PostSDigitizer(this) ;
155 //____________________________________________________________________________
156 void AliPHOSSDigitizer::InitParameters()
160 fPrimThreshold = 0.01 ;
166 //____________________________________________________________________________
167 void AliPHOSSDigitizer::Exec(Option_t *option)
169 // 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 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
268 gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents) ;
272 //__________________________________________________________________
273 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
275 // Setting title to branch SDigits
277 TString stitle(title) ;
279 // check if branch with title already exists
280 TBranch * sdigitsBranch =
281 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ;
282 TBranch * sdigitizerBranch =
283 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
284 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
285 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
286 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
287 Error("SetSDigitsBranch", "Cannot overwrite existing branch with title %s", title) ;
291 Info("SetSDigitsBranch", "-> Changing SDigits file from %s to %s", GetName(), title) ;
295 // Post to the WhiteBoard
296 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
297 gime->PostSDigits( title, GetTitle()) ;
301 //__________________________________________________________________
302 void AliPHOSSDigitizer::Print(Option_t* option)const
304 // Prints parameters of SDigitizer
306 message = "\n------------------- %s -------------\n" ;
307 message += " Writing SDigits to branch with title %s\n" ;
308 message += " with digitization parameters A = %f\n" ;
309 message += " B = %f\n" ;
310 message += " Threshold for Primary assignment= %f\n" ;
311 message += "---------------------------------------------------\n" ;
312 Info("Print", message.Data(), GetName(), GetName(), fA, fB, fPrimThreshold ) ;
316 //__________________________________________________________________
317 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
320 // SDititizers are equal if their pedestal, slope and threshold are equal
322 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
328 //__________________________________________________________________
329 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
331 // Prints list of digits produced in the current pass of AliPHOSDigitizer
334 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
335 TString sdname(GetName()) ;
336 sdname.Remove(sdname.Index(GetTitle())-1) ;
337 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
340 message = "\nAliPHOSSDigitiser: event " ;
341 message += gAlice->GetEvNumber();
342 message += "\n Number of entries in SDigits list " ;
343 message += sdigits->GetEntriesFast() ;
344 char * tempo = new char[8192];
346 if(strstr(option,"all")||strstr(option,"EMC")){
349 AliPHOSDigit * digit;
350 message += "\nEMC sdigits\n" ;
351 message += "\n Id Amplitude Time Index Nprim: Primaries list \n" ;
352 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
354 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
355 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
356 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
357 if(digit->GetNprimary() == 0)
359 sprintf(tempo, "\n%6d %8d %6.5e %4d %2d :",
360 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
363 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
364 sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ;
370 if(strstr(option,"all")||strstr(option,"CPV")){
372 //loop over CPV digits
373 AliPHOSDigit * digit;
375 message += "CPV sdigits\n" ;
376 message += "\n Id Amplitude Index Nprim: Primaries list \n" ;
377 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
379 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
380 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
381 if(digit->GetId() > maxEmc){
382 sprintf(tempo, "\n%6d %8d %4d %2d :",
383 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
386 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
387 sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ;
394 Info("PrintSDigits", message.Data() ) ;
397 //____________________________________________________________________________
398 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)