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