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