]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
Fixes for the creation of PAR files on Mac (Chiara, Gerri)
[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->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
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   TTask("",""),
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   TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), 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   TTask(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   //dtor
151   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
152   if(rl){
153     AliPHOSLoader * phosLoader = 
154       static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
155     if(phosLoader)
156       phosLoader->CleanSDigitizer() ;
157   }
158 }
159
160 //____________________________________________________________________________ 
161 void AliPHOSSDigitizer::Init()
162 {
163   // Uses the Loader to access the required files
164   
165   fInit = kTRUE ; 
166   
167   //to prevent cleaning of this object while GetEvent is called
168   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
169   if (!rl)
170     rl = AliRunLoader::Open(GetTitle(), fEventFolderName) ; 
171
172   AliPHOSLoader * phosLoader = static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
173   phosLoader->PostSDigitizer(this);
174   phosLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
175 }
176
177 //____________________________________________________________________________ 
178 void AliPHOSSDigitizer::InitParameters()
179
180   // initializes the parameters for digitization
181   fPrimThreshold = 0.01 ;
182   fSDigitsInRun  = 0 ;
183 }
184
185 //____________________________________________________________________________
186 void AliPHOSSDigitizer::Exec(Option_t *option) 
187
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.
192   //
193   // Summable digit is a sum of all hits in the same active
194   // volume into digit
195   
196   if (strstr(option, "print") ) {
197     Print() ; 
198     return ; 
199   }
200
201   if(strstr(option,"tim"))
202     gBenchmark->Start("PHOSSDigitizer");
203   
204 /*
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") ;
215   }
216 */
217   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
218   AliPHOSLoader * phosLoader = static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
219   
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()) ) ;
223     return ;
224   }
225   
226   if (fLastEvent == -1) 
227     fLastEvent = rl->GetNumberOfEvents() - 1 ;
228   else 
229     fLastEvent = TMath::Min(fFirstEvent, rl->GetNumberOfEvents()); // only one event at the time 
230   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
231   
232   Int_t ievent, i;
233
234   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {    
235     rl->GetEvent(ievent) ;
236     TTree * treeS = phosLoader->TreeS(); 
237     if(!treeS){
238       phosLoader->MakeTree("S");
239       treeS = phosLoader->TreeS();
240     }
241
242     phosLoader->CleanHits() ; 
243     phosLoader->LoadHits("READ") ;
244
245     TClonesArray * hits    = phosLoader->Hits() ;
246     TClonesArray * sdigits = phosLoader->SDigits() ;
247     if( !sdigits ) {
248       phosLoader->MakeSDigitsArray() ; 
249       sdigits = phosLoader->SDigits() ;
250     }
251     sdigits->Clear();
252     Int_t nSdigits = 0 ;
253
254     //Now make SDigits from hits, for PHOS it is the same, so just copy    
255     for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
256       
257       AliPHOSHit * hit = static_cast<AliPHOSHit *>(hits->At(i)) ;
258       // Assign primary number only if contribution is significant
259       
260       if( hit->GetEnergy() > fPrimThreshold)
261         new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
262                                                hit->GetEnergy() ,hit->GetTime()) ;
263       else
264         new((*sdigits)[nSdigits]) AliPHOSDigit(-1               ,hit->GetId(), 
265                                                hit->GetEnergy() ,hit->GetTime()) ;
266       nSdigits++ ;      
267       
268     }
269  
270     sdigits->Sort() ;
271
272     nSdigits = sdigits->GetEntriesFast() ;
273
274     fSDigitsInRun += nSdigits ;  
275     sdigits->Expand(nSdigits) ;
276     
277     for (i = 0 ; i < nSdigits ; i++) {       
278       AliPHOSDigit * digit = static_cast<AliPHOSDigit *>(sdigits->At(i)) ;
279       digit->SetIndexInList(i) ;     
280     }
281     
282 //    // make Quality Assurance data
283 //
284 //    if (GetQADataMaker()->IsCycleDone() ) {
285 //      GetQADataMaker()->EndOfCycle(AliQAv1::kHITS) ; 
286 //        GetQADataMaker()->EndOfCycle(AliQAv1::kSDIGITS) ; 
287 //      GetQADataMaker()->StartOfCycle(AliQAv1::kHITS) ; 
288 //        GetQADataMaker()->StartOfCycle(AliQAv1::kSDIGITS, kTRUE) ; 
289 //   }
290 //    GetQADataMaker()->Exec(AliQAv1::kHITS, hits) ; 
291 //    GetQADataMaker()->Exec(AliQAv1::kSDIGITS, sdigits) ; 
292 //    GetQADataMaker()->Increment() ;
293         
294     //Now write SDigits
295
296     
297     //First list of sdigits
298
299     Int_t bufferSize = 32000 ;
300     TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
301     sdigitsBranch->Fill() ;
302
303     phosLoader->WriteSDigits("OVERWRITE");
304     phosLoader->WriteSDigitizer("OVERWRITE");
305
306     if(strstr(option,"deb"))
307       PrintSDigits(option) ;
308       
309   }// event loop
310   
311 //  //Write the quality assurance data 
312 //  GetQADataMaker()->EndOfCycle(AliQAv1::kHITS) ;    
313 //  GetQADataMaker()->EndOfCycle(AliQAv1::kSDIGITS) ;    
314 //  GetQADataMaker()->Finish() ;
315
316   Unload();
317
318   if(strstr(option,"tim")){
319     gBenchmark->Stop("PHOSSDigitizer");
320     Info("Exec","   took %f seconds for SDigitizing  %f seconds per event",
321          gBenchmark->GetCpuTime("PHOSSDigitizer"), 
322          gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
323   }
324   
325 }
326
327 //__________________________________________________________________
328 void AliPHOSSDigitizer::Print(const Option_t *)const
329 {
330   // Prints parameters of SDigitizer
331   Info("Print", "\n------------------- %s -------------", GetName() ) ; 
332   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
333   printf("   Threshold for Primary assignment= %f\n", fPrimThreshold)  ; 
334   printf("---------------------------------------------------\n") ;
335   
336 }
337
338 //__________________________________________________________________
339 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
340 {
341   // Equal operator.
342   // SDititizers are equal if their pedestal, slope and threshold are equal
343
344   if(fPrimThreshold==sd.fPrimThreshold)
345     return kTRUE ;
346   else
347     return kFALSE ;
348 }
349
350 //__________________________________________________________________
351 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
352 {
353   // Prints list of digits produced in the current pass of AliPHOSDigitizer
354
355   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
356   AliPHOSLoader * phosLoader = static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
357
358   // Get PHOS Geometry object
359   AliPHOSGeometry *geom;
360   if (!(geom = AliPHOSGeometry::GetInstance())) 
361         geom = AliPHOSGeometry::GetInstance("IHEP","");
362
363   const TClonesArray * sdigits = phosLoader->SDigits() ;
364   
365   Info( "\nPrintSDigits", "event # %d %d sdigits", 
366         gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ; 
367
368   if(strstr(option,"all")||strstr(option,"EMC")){
369     
370     //loop over digits
371     AliPHOSDigit * digit;
372     printf("\nEMC sdigits\n") ; 
373     Int_t maxEmc = geom->GetNModules() * geom->GetNCristalsInModule() ;
374     Int_t index ;
375     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
376          ((static_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
377       digit = static_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
378       //  if(digit->GetNprimary() == 0) 
379       //        continue;
380 //       printf("%6d  %8d    %6.5e %4d      %2d :\n", // YVK
381       printf("%6d  %.4f    %6.5e %4d      %2d :\n",
382               digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
383       Int_t iprimary;
384       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
385         printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
386       }  
387     }    
388   }
389
390   if(strstr(option,"all")||strstr(option,"CPV")){
391     
392     //loop over CPV digits
393     AliPHOSDigit * digit;
394     printf("\nCPV sdigits\n") ; 
395     Int_t maxEmc = geom->GetNModules() * geom->GetNCristalsInModule() ;
396     Int_t index ;
397     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
398       digit = static_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
399       if(digit->GetId() > maxEmc){
400         printf("\n%6d  %8d    %4d      %2d :",
401                 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;  
402         Int_t iprimary;
403         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
404           printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
405         }
406       }    
407     }
408   }
409 }
410
411 //____________________________________________________________________________ 
412 void AliPHOSSDigitizer::Unload() const
413 {
414   // Unloads the objects from the folder
415   AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
416   AliPHOSLoader * phosLoader = static_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
417   
418   phosLoader->UnloadHits() ; 
419   phosLoader->UnloadSDigits() ; 
420 }