Using default task names. Cleaning the task in the destructor if it was posted
[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 //_________________________________________________________________________
20 // This is a TTask that makes SDigits out of Hits
21 // The name of the TTask is also the title of the branch that will contain 
22 // the created SDigits
23 // The title of the TTAsk is the name of the file that contains the hits from
24 // which the SDigits are created
25 // A Summable Digits is the sum of all hits originating 
26 // from one primary in one active cell
27 // A threshold for assignment of the primary to SDigit is applied 
28 // SDigits are written to TreeS, branch "PHOS"
29 // AliPHOSSDigitizer with all current parameters is written 
30 // to TreeS branch "AliPHOSSDigitizer".
31 // Both branches have the same title. If necessary one can produce 
32 // another set of SDigits with different parameters. Two versions
33 // can be distunguished using titles of the branches.
34 // User case:
35 //  root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
36 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
37 //  root [1] s->ExecuteTask()
38 //             // Makes SDigitis for all events stored in galice.root
39 //  root [2] s->SetPedestalParameter(0.001)
40 //             // One can change parameters of digitization
41 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
42 //             // and write them into the new branch
43 // root [4] s->ExecuteTask("deb all tim")
44 //             // available parameters:
45 //             deb - print # of produced SDigitis
46 //             deb all  - print # and list of produced SDigits
47 //             tim - print benchmarking information
48 //
49 //*-- Author :  Dmitri Peressounko (SUBATECH & KI) 
50 //////////////////////////////////////////////////////////////////////////////
51
52
53 // --- ROOT system ---
54 #include "TBenchmark.h"
55 #include "TRandom.h"
56
57 // --- Standard library ---
58
59 // --- AliRoot header files ---
60 #include "AliPHOSGeometry.h" 
61 #include "AliPHOSDigit.h"
62 #include "AliPHOSGetter.h"
63 #include "AliPHOSHit.h"
64 #include "AliPHOSSDigitizer.h"
65
66 ClassImp(AliPHOSSDigitizer)
67
68            
69 //____________________________________________________________________________ 
70   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","")
71 {
72   // ctor
73   InitParameters() ;
74   fDefaultInit = kTRUE ; 
75 }
76
77 //____________________________________________________________________________ 
78 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName, const char * eventFolderName):
79   TTask("PHOS"+AliConfig::fgkSDigitizerTaskName, alirunFileName),
80   fEventFolderName(eventFolderName)
81 {
82
83   // ctor
84   InitParameters() ; 
85   Init();
86   fDefaultInit = kFALSE ; 
87 }
88
89 //____________________________________________________________________________ 
90 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer & sd)
91   : TTask(sd)
92 {
93   //cpy ctor 
94
95   fA             = sd.fA ;
96   fB             = sd.fB ;
97   fPrimThreshold = sd.fPrimThreshold ;
98   fSDigitsInRun  = sd.fSDigitsInRun ;
99   SetName(sd.GetName()) ; 
100   SetTitle(sd.GetTitle()) ; 
101   fEventFolderName = sd.fEventFolderName;
102 }
103
104
105 //____________________________________________________________________________ 
106 AliPHOSSDigitizer::~AliPHOSSDigitizer() {
107   //dtor
108   AliPHOSGetter * gime =
109     AliPHOSGetter::Instance(GetTitle(),fEventFolderName.Data());  
110   gime->PhosLoader()->CleanSDigitizer();
111 }
112 //____________________________________________________________________________ 
113 void AliPHOSSDigitizer::Init()
114 {
115   // Uses the getter to access the required files
116   
117   fInit = kTRUE ; 
118   
119   AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());  
120   if ( gime == 0 ) {
121     Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;  
122     return ;
123   } 
124   
125   TString opt("SDigits") ; 
126   if(gime->VersionExists(opt) ) { 
127     Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
128     fInit = kFALSE ; 
129   }
130
131   gime->PostSDigitizer(this);
132   gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
133   
134 }
135
136 //____________________________________________________________________________ 
137 void AliPHOSSDigitizer::InitParameters()
138
139   // initializes the parameters for digitization
140   fA             = 0;
141   fB             = 10000000.;
142   fPrimThreshold = 0.01 ;
143   fSDigitsInRun  = 0 ;
144 }
145
146 //____________________________________________________________________________
147 void AliPHOSSDigitizer::Exec(Option_t *option) 
148
149   // Collects all hits in the same active volume into digit
150   
151   if (strstr(option, "print") ) {
152     Print() ; 
153     return ; 
154   }
155
156   if(strstr(option,"tim"))
157     gBenchmark->Start("PHOSSDigitizer");
158   
159   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
160
161   //switch off reloading of this task while getting event
162   if (!fInit) { // to prevent overwrite existing file
163     Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
164     return ;
165   }
166
167
168   Int_t nevents = gime->MaxEvent() ; 
169   Int_t ievent ;
170   for(ievent = 0; ievent < nevents; ievent++){
171
172     gime->Event(ievent,"H") ;
173
174     TTree * treeS = gime->TreeS(); 
175     TClonesArray * hits = gime->Hits() ;
176     TClonesArray * sdigits = gime->SDigits() ;
177     sdigits->Clear();
178     Int_t nSdigits = 0 ;
179     //Now make SDigits from hits, for PHOS it is the same, so just copy    
180     Int_t nPrim =  static_cast<Int_t>((gime->TreeH())->GetEntries()) ; 
181     // Attention nPrim is the number of primaries tracked by Geant 
182     // and this number could be different to the number of Primaries in TreeK;
183     Int_t iprim ;
184
185     for (iprim = 0 ; iprim < nPrim ; iprim ++) { 
186       //=========== Get the PHOS branch from Hits Tree for the Primary iprim
187       gime->Track(iprim) ;
188       Int_t i;
189       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
190         AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
191         // Assign primary number only if contribution is significant
192         
193         if( hit->GetEnergy() > fPrimThreshold)
194           new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
195                                                  Digitize(hit->GetEnergy()), hit->GetTime()) ;
196         else
197           new((*sdigits)[nSdigits]) AliPHOSDigit( -1              , hit->GetId(), 
198                                                   Digitize(hit->GetEnergy()), hit->GetTime()) ;
199         nSdigits++ ;    
200         
201       }
202  
203     } // loop over iprim
204
205     sdigits->Sort() ;
206
207     nSdigits = sdigits->GetEntriesFast() ;
208
209     fSDigitsInRun += nSdigits ;  
210     sdigits->Expand(nSdigits) ;
211
212     Int_t i ;
213     for (i = 0 ; i < nSdigits ; i++) { 
214       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
215       digit->SetIndexInList(i) ;     
216     }
217
218     //Now write SDigits
219
220     
221     //First list of sdigits
222
223     Int_t bufferSize = 32000 ;
224     TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
225
226     sdigitsBranch->Fill() ;
227
228     gime->WriteSDigits("OVERWRITE");
229
230     //Next - SDigitizer
231
232     gime->WriteSDigitizer("OVERWRITE");
233
234     if(strstr(option,"deb"))
235       PrintSDigits(option) ;
236   }
237   
238   Unload();
239
240   gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
241   
242   if(strstr(option,"tim")){
243     gBenchmark->Stop("PHOSSDigitizer");
244     Info("Exec","   took %f seconds for SDigitizing  %f seconds per event",
245          gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents) ;
246   }
247 }
248
249 //__________________________________________________________________
250 void AliPHOSSDigitizer::Print()const
251 {
252   // Prints parameters of SDigitizer
253   Info("Print", "\n------------------- %s -------------", GetName() ) ; 
254   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
255   printf("   with digitization parameters  A = %f\n", fA) ; 
256   printf("                                 B = %f\n", fB) ;
257   printf("   Threshold for Primary assignment= %f\n", fPrimThreshold)  ; 
258   printf("---------------------------------------------------\n") ;
259   
260 }
261
262 //__________________________________________________________________
263 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
264 {
265   // Equal operator.
266   // SDititizers are equal if their pedestal, slope and threshold are equal
267
268   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
269     return kTRUE ;
270   else
271     return kFALSE ;
272 }
273
274 //__________________________________________________________________
275 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
276 {
277   // Prints list of digits produced in the current pass of AliPHOSDigitizer
278
279
280   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
281   const TClonesArray * sdigits = gime->SDigits() ;
282   
283   Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ; 
284
285   if(strstr(option,"all")||strstr(option,"EMC")){
286     
287     //loop over digits
288     AliPHOSDigit * digit;
289     printf("\nEMC sdigits\n") ; 
290     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
291     Int_t index ;
292     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
293          ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
294       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
295       //  if(digit->GetNprimary() == 0) 
296       //        continue;
297       printf("%6d  %8d    %6.5e %4d      %2d :",
298               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
299       Int_t iprimary;
300       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
301         printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
302       }  
303     }    
304   }
305
306   if(strstr(option,"all")||strstr(option,"CPV")){
307     
308     //loop over CPV digits
309     AliPHOSDigit * digit;
310     printf("\nCPV sdigits\n") ; 
311     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
312     Int_t index ;
313     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
314       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
315       if(digit->GetId() > maxEmc){
316         printf("\n%6d  %8d    %4d      %2d :",
317                 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;  
318         Int_t iprimary;
319         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
320           printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
321         }
322       }    
323     }
324   }
325 }
326
327 //____________________________________________________________________________ 
328 void AliPHOSSDigitizer::Unload() const
329 {
330   // Unloads the objects from the folder
331   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
332   AliPHOSLoader * loader = gime->PhosLoader() ; 
333   loader->UnloadHits() ; 
334   loader->UnloadSDigits() ; 
335 }