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