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