]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[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
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.
56 // User case:
57 //  root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
58 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
59 //  root [1] s->Digitize()
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
70 //
71 //-- Author :  Dmitri Peressounko (SUBATECH & KI) 
72 //////////////////////////////////////////////////////////////////////////////
73
74
75 // --- ROOT system ---
76 #include "TBenchmark.h"
77
78 // --- Standard library ---
79
80 // --- AliRoot header files ---
81 #include "AliLog.h"
82 #include "AliPHOSGeometry.h" 
83 #include "AliPHOSDigit.h"
84 #include "AliPHOSLoader.h"
85 #include "AliPHOSHit.h"
86 #include "AliPHOSSDigitizer.h"
87
88 ClassImp(AliPHOSSDigitizer)
89
90            
91 //____________________________________________________________________________ 
92 AliPHOSSDigitizer::AliPHOSSDigitizer() : 
93   TNamed("",""),
94   fPrimThreshold(0.f),
95   fDefaultInit(kTRUE),
96   fEventFolderName(""),
97   fInit(kFALSE),
98   fSDigitsInRun(0),
99   fFirstEvent(0),
100   fLastEvent(0)
101 {
102   // ctor
103   // Intialize the quality assurance data maker         
104 }
105
106 //____________________________________________________________________________ 
107 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName, 
108                                      const char * eventFolderName):
109   TNamed("PHOSSDigitizer", alirunFileName),
110   fPrimThreshold(0.f),
111   fDefaultInit(kFALSE),
112   fEventFolderName(eventFolderName),
113   fInit(kFALSE),
114   fSDigitsInRun(0),
115   fFirstEvent(0),
116   fLastEvent(0)
117 {
118   // ctor
119   InitParameters() ; 
120   Init();
121   fDefaultInit = kFALSE ; 
122 }
123
124 //____________________________________________________________________________
125 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer& sd) :
126   TNamed(sd.GetName(), sd.GetTitle()),
127   fPrimThreshold(sd.fPrimThreshold),
128   fDefaultInit(kFALSE),
129   fEventFolderName(sd.fEventFolderName),
130   fInit(kFALSE),
131   fSDigitsInRun(sd.fSDigitsInRun),
132   fFirstEvent(sd.fFirstEvent),
133   fLastEvent(sd.fLastEvent)
134
135   // cpy ctor
136 }
137
138 //_____________________________________________________________________________
139 AliPHOSSDigitizer& AliPHOSSDigitizer::operator = (const AliPHOSSDigitizer& qa)
140 {
141 // assignment operator
142
143   this->~AliPHOSSDigitizer();
144   new(this) AliPHOSSDigitizer(qa);
145   return *this;
146 }
147
148 //____________________________________________________________________________ 
149 AliPHOSSDigitizer::~AliPHOSSDigitizer() 
150 {
151   //dtor
152 }
153
154 //____________________________________________________________________________ 
155 void AliPHOSSDigitizer::Init()
156 {
157   // Uses the Loader to access the required files
158   
159   fInit = kTRUE ; 
160   
161   //to prevent cleaning of this object while GetEvent is called
162   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
163   if (!rl) AliRunLoader::Open(GetTitle(), fEventFolderName) ; 
164 }
165
166 //____________________________________________________________________________ 
167 void AliPHOSSDigitizer::InitParameters()
168
169   // initializes the parameters for digitization
170   fPrimThreshold = 0.01 ;
171   fSDigitsInRun  = 0 ;
172 }
173
174 //____________________________________________________________________________
175 void AliPHOSSDigitizer::Digitize(Option_t *option) 
176
177   // Steering method to produce summable digits for events
178   // in the range from fFirstEvent to fLastEvent.
179   // This range is optionally set by SetEventRange().
180   // if fLastEvent=-1 (by default), then process events until the end.
181   //
182   // Summable digit is a sum of all hits in the same active
183   // volume into digit
184   
185   if (strstr(option, "print") ) {
186     Print() ; 
187     return ; 
188   }
189
190   if(strstr(option,"tim"))
191     gBenchmark->Start("PHOSSDigitizer");
192   
193 /*
194         // check the QA result for RAWS
195   AliQAv1 * qa = AliQAv1::Instance(AliQAv1::kPHOS) ; 
196   if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kFATAL)) {
197         AliFatal("QA status in RAW was Fatal") ;
198   } else if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kERROR)) {
199         AliError("QA status in RAW was Error") ;
200   } else if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kWARNING) ) {
201         AliWarning("QA status in RAW was Warning") ;
202   } else if ( qa->IsSet(AliQAv1::kPHOS, AliQAv1::kRAW, AliQAv1::kINFO) ) {
203         AliInfo("QA status in RAW was Info") ;
204   }
205 */
206   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
207   AliPHOSLoader * phosLoader = static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
208   
209   //switch off reloading of this task while getting event
210   if (!fInit) { // to prevent overwrite existing file
211     AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
212     return ;
213   }
214   
215   if (fLastEvent == -1) 
216     fLastEvent = rl->GetNumberOfEvents() - 1 ;
217   else 
218     fLastEvent = TMath::Min(fFirstEvent, rl->GetNumberOfEvents()); // only one event at the time 
219   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
220   
221   Int_t ievent, i;
222
223   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {    
224     rl->GetEvent(ievent) ;
225     TTree * treeS = phosLoader->TreeS(); 
226     if(!treeS){
227       phosLoader->MakeTree("S");
228       treeS = phosLoader->TreeS();
229     }
230
231     phosLoader->CleanHits() ; 
232     phosLoader->LoadHits("READ") ;
233
234     TClonesArray * hits    = phosLoader->Hits() ;
235     TClonesArray * sdigits = phosLoader->SDigits() ;
236     if( !sdigits ) {
237       phosLoader->MakeSDigitsArray() ; 
238       sdigits = phosLoader->SDigits() ;
239     }
240     sdigits->Clear();
241     Int_t nSdigits = 0 ;
242
243     //Now make SDigits from hits, for PHOS it is the same, so just copy    
244     for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
245       
246       AliPHOSHit * hit = static_cast<AliPHOSHit *>(hits->At(i)) ;
247       // Assign primary number only if contribution is significant
248       
249       if( hit->GetEnergy() > fPrimThreshold)
250         new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
251                                                hit->GetEnergy() ,hit->GetTime()) ;
252       else
253         new((*sdigits)[nSdigits]) AliPHOSDigit(-1               ,hit->GetId(), 
254                                                hit->GetEnergy() ,hit->GetTime()) ;
255       nSdigits++ ;      
256       
257     }
258  
259     sdigits->Sort() ;
260
261     nSdigits = sdigits->GetEntriesFast() ;
262
263     fSDigitsInRun += nSdigits ;  
264     sdigits->Expand(nSdigits) ;
265     
266     for (i = 0 ; i < nSdigits ; i++) {       
267       AliPHOSDigit * digit = static_cast<AliPHOSDigit *>(sdigits->At(i)) ;
268       digit->SetIndexInList(i) ;     
269     }
270     
271 //    // make Quality Assurance data
272 //
273 //    if (GetQADataMaker()->IsCycleDone() ) {
274 //      GetQADataMaker()->EndOfCycle(AliQAv1::kHITS) ; 
275 //        GetQADataMaker()->EndOfCycle(AliQAv1::kSDIGITS) ; 
276 //      GetQADataMaker()->StartOfCycle(AliQAv1::kHITS) ; 
277 //        GetQADataMaker()->StartOfCycle(AliQAv1::kSDIGITS, kTRUE) ; 
278 //   }
279 //    GetQADataMaker()->Exec(AliQAv1::kHITS, hits) ; 
280 //    GetQADataMaker()->Exec(AliQAv1::kSDIGITS, sdigits) ; 
281 //    GetQADataMaker()->Increment() ;
282         
283     //Now write SDigits
284
285     
286     //First list of sdigits
287
288     Int_t bufferSize = 32000 ;
289     TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
290     sdigitsBranch->Fill() ;
291
292     phosLoader->WriteSDigits("OVERWRITE");
293
294     if(strstr(option,"deb"))
295       PrintSDigits(option) ;
296       
297   }// event loop
298   
299 //  //Write the quality assurance data 
300 //  GetQADataMaker()->EndOfCycle(AliQAv1::kHITS) ;    
301 //  GetQADataMaker()->EndOfCycle(AliQAv1::kSDIGITS) ;    
302 //  GetQADataMaker()->Finish() ;
303
304   Unload();
305
306   if(strstr(option,"tim")){
307     gBenchmark->Stop("PHOSSDigitizer");
308     Info("Exec","   took %f seconds for SDigitizing  %f seconds per event",
309          gBenchmark->GetCpuTime("PHOSSDigitizer"), 
310          gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
311   }
312   
313 }
314
315 //__________________________________________________________________
316 void AliPHOSSDigitizer::Print(const Option_t *)const
317 {
318   // Prints parameters of SDigitizer
319   Info("Print", "\n------------------- %s -------------", GetName() ) ; 
320   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
321   printf("   Threshold for Primary assignment= %f\n", fPrimThreshold)  ; 
322   printf("---------------------------------------------------\n") ;
323   
324 }
325
326 //__________________________________________________________________
327 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
328 {
329   // Equal operator.
330   // SDititizers are equal if their pedestal, slope and threshold are equal
331
332   if(fPrimThreshold==sd.fPrimThreshold)
333     return kTRUE ;
334   else
335     return kFALSE ;
336 }
337
338 //__________________________________________________________________
339 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
340 {
341   // Prints list of digits produced in the current pass of AliPHOSDigitizer
342
343   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
344   AliPHOSLoader * phosLoader = static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
345
346   // Get PHOS Geometry object
347   AliPHOSGeometry *geom;
348   if (!(geom = AliPHOSGeometry::GetInstance())) 
349         geom = AliPHOSGeometry::GetInstance("IHEP","");
350
351   const TClonesArray * sdigits = phosLoader->SDigits() ;
352   
353   Info( "\nPrintSDigits", "event # %d %d sdigits", 
354         gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ; 
355
356   if(strstr(option,"all")||strstr(option,"EMC")){
357     
358     //loop over digits
359     AliPHOSDigit * digit;
360     printf("\nEMC sdigits\n") ; 
361     Int_t maxEmc = geom->GetNModules() * geom->GetNCristalsInModule() ;
362     Int_t index ;
363     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
364          ((static_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
365       digit = static_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
366       //  if(digit->GetNprimary() == 0) 
367       //        continue;
368 //       printf("%6d  %8d    %6.5e %4d      %2d :\n", // YVK
369       printf("%6d  %.4f    %6.5e %4d      %2d :\n",
370               digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
371       Int_t iprimary;
372       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
373         printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
374       }  
375     }    
376   }
377
378   if(strstr(option,"all")||strstr(option,"CPV")){
379     
380     //loop over CPV digits
381     AliPHOSDigit * digit;
382     printf("\nCPV sdigits\n") ; 
383     Int_t maxEmc = geom->GetNModules() * geom->GetNCristalsInModule() ;
384     Int_t index ;
385     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
386       digit = static_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
387       if(digit->GetId() > maxEmc){
388         printf("\n%6d  %8d    %4d      %2d :",
389                 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;  
390         Int_t iprimary;
391         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
392           printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
393         }
394       }    
395     }
396   }
397 }
398
399 //____________________________________________________________________________ 
400 void AliPHOSSDigitizer::Unload() const
401 {
402   // Unloads the objects from the folder
403   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
404   AliPHOSLoader * phosLoader = static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
405   
406   phosLoader->UnloadHits() ; 
407   phosLoader->UnloadSDigits() ; 
408 }