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 /* History of cvs commits:
22 * Revision 1.56 2007/10/19 18:04:29 schutz
23 * The standalone QA data maker is called from AliSimulation and AliReconstruction outside the event loop; i.e. re-reading the data. The QA data making in the event loop has been commented out.
25 * Revision 1.55 2007/10/14 21:08:10 schutz
26 * Introduced the checking of QA results from previous step before entering the event loop
28 * Revision 1.54 2007/10/10 09:05:10 schutz
29 * Changing name QualAss to QA
31 * Revision 1.53 2007/09/30 17:08:20 schutz
32 * Introducing the notion of QA data acquisition cycle (needed by online)
34 * Revision 1.52 2007/09/26 14:22:18 cvetan
35 * Important changes to the reconstructor classes. Complete elimination of the run-loaders, which are now steered only from AliReconstruction. Removal of the corresponding Reconstruct() and FillESD() methods.
37 * Revision 1.51 2007/08/07 14:12:03 kharlov
38 * Quality assurance added (Yves Schutz)
40 * Revision 1.50 2006/08/28 10:01:56 kharlov
41 * Effective C++ warnings fixed (Timur Pocheptsov)
43 * Revision 1.49 2006/05/10 06:42:53 kharlov
44 * Remove redundant loop over primaries
46 * Revision 1.48 2006/04/22 10:30:17 hristov
47 * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
49 * Revision 1.47 2005/05/28 14:19:04 schutz
50 * Compilation warnings fixed by T.P.
54 //_________________________________________________________________________
55 // This is a TTask that makes SDigits out of Hits
56 // The name of the TTask is also the title of the branch that will contain
57 // the created SDigits
58 // The title of the TTAsk is the name of the file that contains the hits from
59 // which the SDigits are created
60 // A Summable Digits is the sum of all hits originating
61 // from one primary in one active cell
62 // A threshold for assignment of the primary to SDigit is applied
63 // SDigits are written to TreeS, branch "PHOS"
64 // AliPHOSSDigitizer with all current parameters is written
65 // to TreeS branch "AliPHOSSDigitizer".
66 // Both branches have the same title. If necessary one can produce
67 // another set of SDigits with different parameters. Two versions
68 // can be distunguished using titles of the branches.
70 // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
71 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
72 // root [1] s->ExecuteTask()
73 // // Makes SDigitis for all events stored in galice.root
74 // root [2] s->SetPedestalParameter(0.001)
75 // // One can change parameters of digitization
76 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
77 // // and write them into the new branch
78 // root [4] s->ExecuteTask("deb all tim")
79 // // available parameters:
80 // deb - print # of produced SDigitis
81 // deb all - print # and list of produced SDigits
82 // tim - print benchmarking information
84 //*-- Author : Dmitri Peressounko (SUBATECH & KI)
85 //////////////////////////////////////////////////////////////////////////////
88 // --- ROOT system ---
89 #include "TBenchmark.h"
90 //#include "TObjectTable.h"
92 // --- Standard library ---
94 // --- AliRoot header files ---
96 #include "AliPHOSGeometry.h"
97 #include "AliPHOSDigit.h"
98 #include "AliPHOSGetter.h"
99 #include "AliPHOSHit.h"
100 #include "AliPHOSSDigitizer.h"
102 ClassImp(AliPHOSSDigitizer)
105 //____________________________________________________________________________
106 AliPHOSSDigitizer::AliPHOSSDigitizer() :
110 fEventFolderName(""),
117 // Intialize the quality assurance data maker
120 //____________________________________________________________________________
121 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName,
122 const char * eventFolderName):
123 TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
125 fDefaultInit(kFALSE),
126 fEventFolderName(eventFolderName),
135 fDefaultInit = kFALSE ;
138 //____________________________________________________________________________
139 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer& sd) :
140 TTask(sd.GetName(), sd.GetTitle()),
141 fPrimThreshold(sd.fPrimThreshold),
142 fDefaultInit(kFALSE),
143 fEventFolderName(sd.fEventFolderName),
145 fSDigitsInRun(sd.fSDigitsInRun),
146 fFirstEvent(sd.fFirstEvent),
147 fLastEvent(sd.fLastEvent)
152 //_____________________________________________________________________________
153 AliPHOSSDigitizer& AliPHOSSDigitizer::operator = (const AliPHOSSDigitizer& qa)
155 // assignment operator
157 this->~AliPHOSSDigitizer();
158 new(this) AliPHOSSDigitizer(qa);
162 //____________________________________________________________________________
163 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
165 // AliPHOSGetter * gime =
166 // AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());
167 AliPHOSGetter * gime =
168 AliPHOSGetter::Instance();
169 gime->PhosLoader()->CleanSDigitizer();
173 //____________________________________________________________________________
174 void AliPHOSSDigitizer::Init()
176 // Uses the getter to access the required files
180 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
182 Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
186 TString opt("SDigits") ;
187 if(gime->VersionExists(opt) ) {
188 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
192 gime->PostSDigitizer(this);
193 gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
195 // fQADM = new AliPHOSQADataMaker() ;
199 //____________________________________________________________________________
200 void AliPHOSSDigitizer::InitParameters()
202 // initializes the parameters for digitization
203 fPrimThreshold = 0.01 ;
207 //____________________________________________________________________________
208 void AliPHOSSDigitizer::Exec(Option_t *option)
210 // Steering method to produce summable digits for events
211 // in the range from fFirstEvent to fLastEvent.
212 // This range is optionally set by SetEventRange().
213 // if fLastEvent=-1 (by default), then process events until the end.
215 // Summable digit is a sum of all hits in the same active
218 if (strstr(option, "print") ) {
223 if(strstr(option,"tim"))
224 gBenchmark->Start("PHOSSDigitizer");
227 // check the QA result for RAWS
228 AliQA * qa = AliQA::Instance(AliQA::kPHOS) ;
229 if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kFATAL)) {
230 AliFatal("QA status in RAW was Fatal") ;
231 } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kERROR)) {
232 AliError("QA status in RAW was Error") ;
233 } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kWARNING) ) {
234 AliWarning("QA status in RAW was Warning") ;
235 } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kINFO) ) {
236 AliInfo("QA status in RAW was Info") ;
239 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
241 //switch off reloading of this task while getting event
242 if (!fInit) { // to prevent overwrite existing file
243 AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
247 if (fLastEvent == -1)
248 fLastEvent = gime->MaxEvent() - 1 ;
250 fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time
251 Int_t nEvents = fLastEvent - fFirstEvent + 1;
255 //AliMemoryWatcher memwatcher;
257 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
258 gime->Event(ievent,"H") ;
259 TTree * treeS = gime->TreeS();
260 TClonesArray * hits = gime->Hits() ;
261 TClonesArray * sdigits = gime->SDigits() ;
264 //Now make SDigits from hits, for PHOS it is the same, so just copy
265 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
266 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
267 // Assign primary number only if contribution is significant
269 if( hit->GetEnergy() > fPrimThreshold)
270 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
271 hit->GetEnergy() ,hit->GetTime()) ;
273 new((*sdigits)[nSdigits]) AliPHOSDigit(-1 ,hit->GetId(),
274 hit->GetEnergy() ,hit->GetTime()) ;
281 nSdigits = sdigits->GetEntriesFast() ;
283 fSDigitsInRun += nSdigits ;
284 sdigits->Expand(nSdigits) ;
286 for (i = 0 ; i < nSdigits ; i++) {
287 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
288 digit->SetIndexInList(i) ;
291 // // make Quality Assurance data
293 // if (GetQADataMaker()->IsCycleDone() ) {
294 // GetQADataMaker()->EndOfCycle(AliQA::kHITS) ;
295 // GetQADataMaker()->EndOfCycle(AliQA::kSDIGITS) ;
296 // GetQADataMaker()->StartOfCycle(AliQA::kHITS) ;
297 // GetQADataMaker()->StartOfCycle(AliQA::kSDIGITS, kTRUE) ;
299 // GetQADataMaker()->Exec(AliQA::kHITS, hits) ;
300 // GetQADataMaker()->Exec(AliQA::kSDIGITS, sdigits) ;
301 // GetQADataMaker()->Increment() ;
306 //First list of sdigits
308 Int_t bufferSize = 32000 ;
309 TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
311 sdigitsBranch->Fill() ;
313 gime->WriteSDigits("OVERWRITE");
317 gime->WriteSDigitizer("OVERWRITE");
318 //gObjectTable->Print() ;
320 if(strstr(option,"deb"))
321 PrintSDigits(option) ;
323 //memwatcher.Watch(ievent);
326 // //Write the quality assurance data
327 // GetQADataMaker()->EndOfCycle(AliQA::kHITS) ;
328 // GetQADataMaker()->EndOfCycle(AliQA::kSDIGITS) ;
329 // GetQADataMaker()->Finish() ;
333 // gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
335 if(strstr(option,"tim")){
336 gBenchmark->Stop("PHOSSDigitizer");
337 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
338 gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
341 //TFile f("out.root","RECREATE");
342 //memwatcher.WriteToFile();
346 //__________________________________________________________________
347 void AliPHOSSDigitizer::Print(const Option_t *)const
349 // Prints parameters of SDigitizer
350 Info("Print", "\n------------------- %s -------------", GetName() ) ;
351 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
352 printf(" Threshold for Primary assignment= %f\n", fPrimThreshold) ;
353 printf("---------------------------------------------------\n") ;
357 //__________________________________________________________________
358 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
361 // SDititizers are equal if their pedestal, slope and threshold are equal
363 if(fPrimThreshold==sd.fPrimThreshold)
369 //__________________________________________________________________
370 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
372 // Prints list of digits produced in the current pass of AliPHOSDigitizer
375 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
376 const TClonesArray * sdigits = gime->SDigits() ;
378 Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ;
380 if(strstr(option,"all")||strstr(option,"EMC")){
383 AliPHOSDigit * digit;
384 printf("\nEMC sdigits\n") ;
385 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
387 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
388 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
389 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
390 // if(digit->GetNprimary() == 0)
392 // printf("%6d %8d %6.5e %4d %2d :\n", // YVK
393 printf("%6d %.4f %6.5e %4d %2d :\n",
394 digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
396 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
397 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
402 if(strstr(option,"all")||strstr(option,"CPV")){
404 //loop over CPV digits
405 AliPHOSDigit * digit;
406 printf("\nCPV sdigits\n") ;
407 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
409 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
410 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
411 if(digit->GetId() > maxEmc){
412 printf("\n%6d %8d %4d %2d :",
413 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
415 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
416 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
423 //____________________________________________________________________________
424 void AliPHOSSDigitizer::Unload() const
426 // Unloads the objects from the folder
427 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
428 AliPHOSLoader * loader = gime->PhosLoader() ;
429 loader->UnloadHits() ;
430 loader->UnloadSDigits() ;