]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
The standalone QA data maker is called from AliSimulation and AliReconstruction outsi...
[u/mrichter/AliRoot.git] / PHOS / AliPHOSSDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16
17 /* $Id$ */
18
19 /* History of cvs commits:
20  *
21  * $Log$
22  * Revision 1.55  2007/10/14 21:08:10  schutz
23  * Introduced the checking of QA results from previous step before entering the event loop
24  *
25  * Revision 1.54  2007/10/10 09:05:10  schutz
26  * Changing name QualAss to QA
27  *
28  * Revision 1.53  2007/09/30 17:08:20  schutz
29  * Introducing the notion of QA data acquisition cycle (needed by online)
30  *
31  * Revision 1.52  2007/09/26 14:22:18  cvetan
32  * 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.
33  *
34  * Revision 1.51  2007/08/07 14:12:03  kharlov
35  * Quality assurance added (Yves Schutz)
36  *
37  * Revision 1.50  2006/08/28 10:01:56  kharlov
38  * Effective C++ warnings fixed (Timur Pocheptsov)
39  *
40  * Revision 1.49  2006/05/10 06:42:53  kharlov
41  * Remove redundant loop over primaries
42  *
43  * Revision 1.48  2006/04/22 10:30:17  hristov
44  * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
45  *
46  * Revision 1.47  2005/05/28 14:19:04  schutz
47  * Compilation warnings fixed by T.P.
48  *
49  */
50
51 //_________________________________________________________________________
52 // This is a TTask that makes SDigits out of Hits
53 // The name of the TTask is also the title of the branch that will contain 
54 // the created SDigits
55 // The title of the TTAsk is the name of the file that contains the hits from
56 // which the SDigits are created
57 // A Summable Digits is the sum of all hits originating 
58 // from one primary in one active cell
59 // A threshold for assignment of the primary to SDigit is applied 
60 // SDigits are written to TreeS, branch "PHOS"
61 // AliPHOSSDigitizer with all current parameters is written 
62 // to TreeS branch "AliPHOSSDigitizer".
63 // Both branches have the same title. If necessary one can produce 
64 // another set of SDigits with different parameters. Two versions
65 // can be distunguished using titles of the branches.
66 // User case:
67 //  root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
68 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
69 //  root [1] s->ExecuteTask()
70 //             // Makes SDigitis for all events stored in galice.root
71 //  root [2] s->SetPedestalParameter(0.001)
72 //             // One can change parameters of digitization
73 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
74 //             // and write them into the new branch
75 // root [4] s->ExecuteTask("deb all tim")
76 //             // available parameters:
77 //             deb - print # of produced SDigitis
78 //             deb all  - print # and list of produced SDigits
79 //             tim - print benchmarking information
80 //
81 //*-- Author :  Dmitri Peressounko (SUBATECH & KI) 
82 //////////////////////////////////////////////////////////////////////////////
83
84
85 // --- ROOT system ---
86 #include "TBenchmark.h"
87                 //#include "TObjectTable.h"
88
89 // --- Standard library ---
90
91 // --- AliRoot header files ---
92 #include "AliLog.h"
93 #include "AliPHOSGeometry.h" 
94 #include "AliPHOSDigit.h"
95 #include "AliPHOSGetter.h"
96 #include "AliPHOSHit.h"
97 #include "AliPHOSSDigitizer.h"
98 #include "AliPHOSQADataMaker.h" 
99
100 ClassImp(AliPHOSSDigitizer)
101
102            
103 //____________________________________________________________________________ 
104 AliPHOSSDigitizer::AliPHOSSDigitizer() : 
105   TTask("",""),
106   fA(0.f), fB(0.f),
107   fPrimThreshold(0.f),
108   fDefaultInit(kTRUE),
109   fEventFolderName(""),
110   fInit(kFALSE),
111   fSDigitsInRun(0),
112   fFirstEvent(0),
113   fLastEvent(0)
114 //  , fQADM (0x0)
115 {
116   // ctor
117   // Intialize the quality assurance data maker         
118 }
119
120 //____________________________________________________________________________ 
121 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName, 
122                                      const char * eventFolderName):
123   TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
124   fA(0.f), fB(0.f),
125   fPrimThreshold(0.f),
126   fDefaultInit(kFALSE),
127   fEventFolderName(eventFolderName),
128   fInit(kFALSE),
129   fSDigitsInRun(0),
130   fFirstEvent(0),
131   fLastEvent(0)
132 //  , fQADM (0x0)
133 {
134   // ctor
135   InitParameters() ; 
136   Init();
137   fDefaultInit = kFALSE ; 
138 //  // Intialize the quality assurance data maker 
139 //  //FIXME: get the run number
140 //  Int_t run = 0 ;
141 //  //EMXIF     
142 //  GetQADataMaker()->Init(AliQA::kHITS, run, fgkCycles) ;
143 //  GetQADataMaker()->StartOfCycle(AliQA::kHITS) ;
144 //  GetQADataMaker()->Init(AliQA::kSDIGITS, run) ; 
145 //  GetQADataMaker()->StartOfCycle(AliQA::kSDIGITS, kTRUE) ;
146 }
147
148 //____________________________________________________________________________
149 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer& sd) :
150   TTask(sd.GetName(), sd.GetTitle()),
151   fA(sd.fA), fB(sd.fB),
152   fPrimThreshold(sd.fPrimThreshold),
153   fDefaultInit(kFALSE),
154   fEventFolderName(sd.fEventFolderName),
155   fInit(kFALSE),
156   fSDigitsInRun(sd.fSDigitsInRun),
157   fFirstEvent(sd.fFirstEvent),
158   fLastEvent(sd.fLastEvent)
159 //  , fQADM (sd.fQADM)
160
161   // cpy ctor
162 //  // Intialize the quality assurance data maker       
163 //  //FIXME: get the run number
164 //  Int_t run = 0 ;
165 //  //EMXIF     
166 //  GetQADataMaker()->Init(AliQA::kHITS, run, fgkCycles) ;
167 //  GetQADataMaker()->StartOfCycle(AliQA::kHITS) ;
168 //  GetQADataMaker()->Init(AliQA::kSDIGITS, run) ; 
169 //  GetQADataMaker()->StartOfCycle(AliQA::kSDIGITS, kTRUE) ;
170 }
171
172 //_____________________________________________________________________________
173 AliPHOSSDigitizer& AliPHOSSDigitizer::operator = (const AliPHOSSDigitizer& qa)
174 {
175 // assignment operator
176
177   this->~AliPHOSSDigitizer();
178   new(this) AliPHOSSDigitizer(qa);
179   return *this;
180 }
181
182 //____________________________________________________________________________ 
183 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
184   //dtor
185   //  AliPHOSGetter * gime =
186   //  AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());  
187   AliPHOSGetter * gime =
188     AliPHOSGetter::Instance();  
189   gime->PhosLoader()->CleanSDigitizer();
190 //  delete fQADM ; 
191 }
192
193 //____________________________________________________________________________ 
194 void AliPHOSSDigitizer::Init()
195 {
196   // Uses the getter to access the required files
197   
198   fInit = kTRUE ; 
199   
200   AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());  
201   if ( gime == 0 ) {
202     Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;  
203     return ;
204   } 
205   
206   TString opt("SDigits") ; 
207   if(gime->VersionExists(opt) ) { 
208     Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
209     fInit = kFALSE ; 
210   }
211
212   gime->PostSDigitizer(this);
213   gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
214  
215 //  fQADM = new AliPHOSQADataMaker() ;  
216
217 }
218
219 //____________________________________________________________________________ 
220 void AliPHOSSDigitizer::InitParameters()
221
222   // initializes the parameters for digitization
223   fA             = 0;
224   fB             = 10000000.;
225   fPrimThreshold = 0.01 ;
226   fSDigitsInRun  = 0 ;
227 }
228
229 //____________________________________________________________________________
230 void AliPHOSSDigitizer::Exec(Option_t *option) 
231
232   // Steering method to produce summable digits for events
233   // in the range from fFirstEvent to fLastEvent.
234   // This range is optionally set by SetEventRange().
235   // if fLastEvent=-1 (by default), then process events until the end.
236   //
237   // Summable digit is a sum of all hits in the same active
238   // volume into digit
239   
240   if (strstr(option, "print") ) {
241     Print() ; 
242     return ; 
243   }
244
245   if(strstr(option,"tim"))
246     gBenchmark->Start("PHOSSDigitizer");
247   
248   // check the QA result for RAWS
249   AliQA * qa = AliQA::Instance(AliQA::kPHOS) ; 
250   if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kFATAL)) {
251         AliFatal("QA status in RAW was Fatal") ;
252   } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kERROR)) {
253         AliError("QA status in RAW was Error") ;
254   } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kWARNING) ) {
255         AliWarning("QA status in RAW was Warning") ;
256   } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kRAW, AliQA::kINFO) ) {
257         AliInfo("QA status in RAW was Info") ;
258   }
259
260   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
261
262   //switch off reloading of this task while getting event
263   if (!fInit) { // to prevent overwrite existing file
264     AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
265     return ;
266   }
267
268   if (fLastEvent == -1) 
269     fLastEvent = gime->MaxEvent() - 1 ;
270   else 
271     fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time 
272   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
273   
274   Int_t ievent, i;
275
276   //AliMemoryWatcher memwatcher;
277
278   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {    
279     gime->Event(ievent,"H") ;   
280     TTree * treeS = gime->TreeS(); 
281     TClonesArray * hits = gime->Hits() ;
282     TClonesArray * sdigits = gime->SDigits() ;
283     sdigits->Clear();
284     Int_t nSdigits = 0 ;
285     //Now make SDigits from hits, for PHOS it is the same, so just copy    
286     for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
287       AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
288       // Assign primary number only if contribution is significant
289       
290       if( hit->GetEnergy() > fPrimThreshold)
291         new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
292                                                hit->GetEnergy() ,hit->GetTime()) ;
293       else
294         new((*sdigits)[nSdigits]) AliPHOSDigit(-1               ,hit->GetId(), 
295                                                hit->GetEnergy() ,hit->GetTime()) ;
296       nSdigits++ ;      
297       
298     }
299  
300     sdigits->Sort() ;
301
302     nSdigits = sdigits->GetEntriesFast() ;
303
304     fSDigitsInRun += nSdigits ;  
305     sdigits->Expand(nSdigits) ;
306
307     for (i = 0 ; i < nSdigits ; i++) { 
308       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
309       digit->SetIndexInList(i) ;     
310     }
311
312 //    // make Quality Assurance data
313 //
314 //    if (GetQADataMaker()->IsCycleDone() ) {
315 //      GetQADataMaker()->EndOfCycle(AliQA::kHITS) ; 
316 //        GetQADataMaker()->EndOfCycle(AliQA::kSDIGITS) ; 
317 //      GetQADataMaker()->StartOfCycle(AliQA::kHITS) ; 
318 //        GetQADataMaker()->StartOfCycle(AliQA::kSDIGITS, kTRUE) ; 
319 //   }
320 //    GetQADataMaker()->Exec(AliQA::kHITS, hits) ; 
321 //    GetQADataMaker()->Exec(AliQA::kSDIGITS, sdigits) ; 
322 //    GetQADataMaker()->Increment() ;
323         
324     //Now write SDigits
325
326     
327     //First list of sdigits
328
329     Int_t bufferSize = 32000 ;
330     TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
331
332     sdigitsBranch->Fill() ;
333
334     gime->WriteSDigits("OVERWRITE");
335
336     //Next - SDigitizer
337
338     gime->WriteSDigitizer("OVERWRITE");
339     //gObjectTable->Print() ; 
340   
341     if(strstr(option,"deb"))
342       PrintSDigits(option) ;
343       
344     //memwatcher.Watch(ievent); 
345   }// event loop
346   
347 //  //Write the quality assurance data 
348 //  GetQADataMaker()->EndOfCycle(AliQA::kHITS) ;    
349 //  GetQADataMaker()->EndOfCycle(AliQA::kSDIGITS) ;    
350 //  GetQADataMaker()->Finish() ;
351
352   Unload();
353
354   //  gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
355   
356   if(strstr(option,"tim")){
357     gBenchmark->Stop("PHOSSDigitizer");
358     Info("Exec","   took %f seconds for SDigitizing  %f seconds per event",
359          gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
360   }
361   
362   //TFile f("out.root","RECREATE");
363   //memwatcher.WriteToFile(); 
364   //f.Close();
365 }
366
367 //__________________________________________________________________
368 void AliPHOSSDigitizer::Print(const Option_t *)const
369 {
370   // Prints parameters of SDigitizer
371   Info("Print", "\n------------------- %s -------------", GetName() ) ; 
372   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
373   printf("   with digitization parameters  A = %f\n", fA) ; 
374   printf("                                 B = %f\n", fB) ;
375   printf("   Threshold for Primary assignment= %f\n", fPrimThreshold)  ; 
376   printf("---------------------------------------------------\n") ;
377   
378 }
379
380 //__________________________________________________________________
381 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
382 {
383   // Equal operator.
384   // SDititizers are equal if their pedestal, slope and threshold are equal
385
386   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
387     return kTRUE ;
388   else
389     return kFALSE ;
390 }
391
392 //__________________________________________________________________
393 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
394 {
395   // Prints list of digits produced in the current pass of AliPHOSDigitizer
396
397
398   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
399   const TClonesArray * sdigits = gime->SDigits() ;
400   
401   Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ; 
402
403   if(strstr(option,"all")||strstr(option,"EMC")){
404     
405     //loop over digits
406     AliPHOSDigit * digit;
407     printf("\nEMC sdigits\n") ; 
408     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
409     Int_t index ;
410     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
411          ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
412       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
413       //  if(digit->GetNprimary() == 0) 
414       //        continue;
415 //       printf("%6d  %8d    %6.5e %4d      %2d :\n", // YVK
416       printf("%6d  %.4f    %6.5e %4d      %2d :\n",
417               digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
418       Int_t iprimary;
419       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
420         printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
421       }  
422     }    
423   }
424
425   if(strstr(option,"all")||strstr(option,"CPV")){
426     
427     //loop over CPV digits
428     AliPHOSDigit * digit;
429     printf("\nCPV sdigits\n") ; 
430     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
431     Int_t index ;
432     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
433       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
434       if(digit->GetId() > maxEmc){
435         printf("\n%6d  %8d    %4d      %2d :",
436                 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;  
437         Int_t iprimary;
438         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
439           printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
440         }
441       }    
442     }
443   }
444 }
445
446 //____________________________________________________________________________ 
447 void AliPHOSSDigitizer::Unload() const
448 {
449   // Unloads the objects from the folder
450   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
451   AliPHOSLoader * loader = gime->PhosLoader() ; 
452   loader->UnloadHits() ; 
453   loader->UnloadSDigits() ; 
454 }