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)
41 //_________________________________________________________________________
42 // This is a TTask that makes SDigits out of Hits
43 // The name of the TTask is also the title of the branch that will contain
44 // the created SDigits
45 // The title of the TTAsk is the name of the file that contains the hits from
46 // which the SDigits are created
47 // A Summable Digits is the sum of all hits originating
48 // from one primary in one active cell
49 // A threshold for assignment of the primary to SDigit is applied
50 // SDigits are written to TreeS, branch "PHOS"
51 // AliPHOSSDigitizer with all current parameters is written
52 // to TreeS branch "AliPHOSSDigitizer".
53 // Both branches have the same title. If necessary one can produce
54 // another set of SDigits with different parameters. Two versions
55 // can be distunguished using titles of the branches.
57 // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
58 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
59 // root [1] s->ExecuteTask()
60 // // Makes SDigitis for all events stored in galice.root
61 // root [2] s->SetPedestalParameter(0.001)
62 // // One can change parameters of digitization
63 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
64 // // and write them into the new branch
65 // root [4] s->ExecuteTask("deb all tim")
66 // // available parameters:
67 // deb - print # of produced SDigitis
68 // deb all - print # and list of produced SDigits
69 // tim - print benchmarking information
71 //-- Author : Dmitri Peressounko (SUBATECH & KI)
72 //////////////////////////////////////////////////////////////////////////////
75 // --- ROOT system ---
76 #include "TBenchmark.h"
78 // --- Standard library ---
80 // --- AliRoot header files ---
82 #include "AliPHOSGeometry.h"
83 #include "AliPHOSDigit.h"
84 #include "AliPHOSLoader.h"
85 #include "AliPHOSHit.h"
86 #include "AliPHOSSDigitizer.h"
88 ClassImp(AliPHOSSDigitizer)
91 //____________________________________________________________________________
92 AliPHOSSDigitizer::AliPHOSSDigitizer() :
103 // Intialize the quality assurance data maker
106 //____________________________________________________________________________
107 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName,
108 const char * eventFolderName):
109 TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
111 fDefaultInit(kFALSE),
112 fEventFolderName(eventFolderName),
121 fDefaultInit = kFALSE ;
124 //____________________________________________________________________________
125 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer& sd) :
126 TTask(sd.GetName(), sd.GetTitle()),
127 fPrimThreshold(sd.fPrimThreshold),
128 fDefaultInit(kFALSE),
129 fEventFolderName(sd.fEventFolderName),
131 fSDigitsInRun(sd.fSDigitsInRun),
132 fFirstEvent(sd.fFirstEvent),
133 fLastEvent(sd.fLastEvent)
138 //_____________________________________________________________________________
139 AliPHOSSDigitizer& AliPHOSSDigitizer::operator = (const AliPHOSSDigitizer& qa)
141 // assignment operator
143 this->~AliPHOSSDigitizer();
144 new(this) AliPHOSSDigitizer(qa);
148 //____________________________________________________________________________
149 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
151 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
153 AliPHOSLoader * phosLoader =
154 dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
156 phosLoader->CleanSDigitizer() ;
160 //____________________________________________________________________________
161 void AliPHOSSDigitizer::Init()
163 // Uses the Loader to access the required files
167 //to prevent cleaning of this object while GetEvent is called
168 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
170 rl = AliRunLoader::Open(GetTitle(), fEventFolderName) ;
172 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
173 phosLoader->PostSDigitizer(this);
174 phosLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
177 //____________________________________________________________________________
178 void AliPHOSSDigitizer::InitParameters()
180 // initializes the parameters for digitization
181 fPrimThreshold = 0.01 ;
185 //____________________________________________________________________________
186 void AliPHOSSDigitizer::Exec(Option_t *option)
188 // Steering method to produce summable digits for events
189 // in the range from fFirstEvent to fLastEvent.
190 // This range is optionally set by SetEventRange().
191 // if fLastEvent=-1 (by default), then process events until the end.
193 // Summable digit is a sum of all hits in the same active
196 if (strstr(option, "print") ) {
201 if(strstr(option,"tim"))
202 gBenchmark->Start("PHOSSDigitizer");
205 // check the QA result for RAWS
206 AliQAv1 * qa = AliQAv1::Instance(AliQAv1::kPHOS) ;
207 if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kFATAL)) {
208 AliFatal("QA status in RAW was Fatal") ;
209 } else if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kERROR)) {
210 AliError("QA status in RAW was Error") ;
211 } else if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kWARNING) ) {
212 AliWarning("QA status in RAW was Warning") ;
213 } else if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kINFO) ) {
214 AliInfo("QA status in RAW was Info") ;
217 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
218 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
220 //switch off reloading of this task while getting event
221 if (!fInit) { // to prevent overwrite existing file
222 AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
226 if (fLastEvent == -1)
227 fLastEvent = rl->GetNumberOfEvents() - 1 ;
229 fLastEvent = TMath::Min(fFirstEvent, rl->GetNumberOfEvents()); // only one event at the time
230 Int_t nEvents = fLastEvent - fFirstEvent + 1;
234 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
235 rl->GetEvent(ievent) ;
236 TTree * treeS = phosLoader->TreeS();
238 phosLoader->MakeTree("S");
239 treeS = phosLoader->TreeS();
242 phosLoader->CleanHits() ;
243 phosLoader->LoadHits("READ") ;
245 TClonesArray * hits = phosLoader->Hits() ;
246 TClonesArray * sdigits = phosLoader->SDigits() ;
248 phosLoader->MakeSDigitsArray() ;
249 sdigits = phosLoader->SDigits() ;
254 //Now make SDigits from hits, for PHOS it is the same, so just copy
255 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
256 AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
257 // Assign primary number only if contribution is significant
259 if( hit->GetEnergy() > fPrimThreshold)
260 new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
261 hit->GetEnergy() ,hit->GetTime()) ;
263 new((*sdigits)[nSdigits]) AliPHOSDigit(-1 ,hit->GetId(),
264 hit->GetEnergy() ,hit->GetTime()) ;
271 nSdigits = sdigits->GetEntriesFast() ;
273 fSDigitsInRun += nSdigits ;
274 sdigits->Expand(nSdigits) ;
276 for (i = 0 ; i < nSdigits ; i++) {
277 AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ;
278 digit->SetIndexInList(i) ;
281 // // make Quality Assurance data
283 // if (GetQADataMaker()->IsCycleDone() ) {
284 // GetQADataMaker()->EndOfCycle(AliQAv1::kHITS) ;
285 // GetQADataMaker()->EndOfCycle(AliQAv1::kSDIGITS) ;
286 // GetQADataMaker()->StartOfCycle(AliQAv1::kHITS) ;
287 // GetQADataMaker()->StartOfCycle(AliQAv1::kSDIGITS, kTRUE) ;
289 // GetQADataMaker()->Exec(AliQAv1::kHITS, hits) ;
290 // GetQADataMaker()->Exec(AliQAv1::kSDIGITS, sdigits) ;
291 // GetQADataMaker()->Increment() ;
296 //First list of sdigits
298 Int_t bufferSize = 32000 ;
299 TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
300 sdigitsBranch->Fill() ;
302 phosLoader->WriteSDigits("OVERWRITE");
303 phosLoader->WriteSDigitizer("OVERWRITE");
305 if(strstr(option,"deb"))
306 PrintSDigits(option) ;
310 // //Write the quality assurance data
311 // GetQADataMaker()->EndOfCycle(AliQAv1::kHITS) ;
312 // GetQADataMaker()->EndOfCycle(AliQAv1::kSDIGITS) ;
313 // GetQADataMaker()->Finish() ;
317 if(strstr(option,"tim")){
318 gBenchmark->Stop("PHOSSDigitizer");
319 Info("Exec"," took %f seconds for SDigitizing %f seconds per event",
320 gBenchmark->GetCpuTime("PHOSSDigitizer"),
321 gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
326 //__________________________________________________________________
327 void AliPHOSSDigitizer::Print(const Option_t *)const
329 // Prints parameters of SDigitizer
330 Info("Print", "\n------------------- %s -------------", GetName() ) ;
331 printf(" Writing SDigits to branch with title %s\n", fEventFolderName.Data()) ;
332 printf(" Threshold for Primary assignment= %f\n", fPrimThreshold) ;
333 printf("---------------------------------------------------\n") ;
337 //__________________________________________________________________
338 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
341 // SDititizers are equal if their pedestal, slope and threshold are equal
343 if(fPrimThreshold==sd.fPrimThreshold)
349 //__________________________________________________________________
350 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
352 // Prints list of digits produced in the current pass of AliPHOSDigitizer
354 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
355 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
357 // Get PHOS Geometry object
358 AliPHOSGeometry *geom;
359 if (!(geom = AliPHOSGeometry::GetInstance()))
360 geom = AliPHOSGeometry::GetInstance("IHEP","");
362 const TClonesArray * sdigits = phosLoader->SDigits() ;
364 Info( "\nPrintSDigits", "event # %d %d sdigits",
365 gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ;
367 if(strstr(option,"all")||strstr(option,"EMC")){
370 AliPHOSDigit * digit;
371 printf("\nEMC sdigits\n") ;
372 Int_t maxEmc = geom->GetNModules() * geom->GetNCristalsInModule() ;
374 for (index = 0 ; (index < sdigits->GetEntriesFast()) &&
375 ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
376 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
377 // if(digit->GetNprimary() == 0)
379 // printf("%6d %8d %6.5e %4d %2d :\n", // YVK
380 printf("%6d %.4f %6.5e %4d %2d :\n",
381 digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
383 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
384 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
389 if(strstr(option,"all")||strstr(option,"CPV")){
391 //loop over CPV digits
392 AliPHOSDigit * digit;
393 printf("\nCPV sdigits\n") ;
394 Int_t maxEmc = geom->GetNModules() * geom->GetNCristalsInModule() ;
396 for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
397 digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
398 if(digit->GetId() > maxEmc){
399 printf("\n%6d %8d %4d %2d :",
400 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
402 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
403 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
410 //____________________________________________________________________________
411 void AliPHOSSDigitizer::Unload() const
413 // Unloads the objects from the folder
414 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
415 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
416 phosLoader->UnloadHits() ;
417 phosLoader->UnloadSDigits() ;