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