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