]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
The patch fixes
[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.47  2005/05/28 14:19:04  schutz
23  * Compilation warnings fixed by T.P.
24  *
25  */
26
27 //_________________________________________________________________________
28 // This is a TTask that makes SDigits out of Hits
29 // The name of the TTask is also the title of the branch that will contain 
30 // the created SDigits
31 // The title of the TTAsk is the name of the file that contains the hits from
32 // which the SDigits are created
33 // A Summable Digits is the sum of all hits originating 
34 // from one primary in one active cell
35 // A threshold for assignment of the primary to SDigit is applied 
36 // SDigits are written to TreeS, branch "PHOS"
37 // AliPHOSSDigitizer with all current parameters is written 
38 // to TreeS branch "AliPHOSSDigitizer".
39 // Both branches have the same title. If necessary one can produce 
40 // another set of SDigits with different parameters. Two versions
41 // can be distunguished using titles of the branches.
42 // User case:
43 //  root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
44 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
45 //  root [1] s->ExecuteTask()
46 //             // Makes SDigitis for all events stored in galice.root
47 //  root [2] s->SetPedestalParameter(0.001)
48 //             // One can change parameters of digitization
49 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
50 //             // and write them into the new branch
51 // root [4] s->ExecuteTask("deb all tim")
52 //             // available parameters:
53 //             deb - print # of produced SDigitis
54 //             deb all  - print # and list of produced SDigits
55 //             tim - print benchmarking information
56 //
57 //*-- Author :  Dmitri Peressounko (SUBATECH & KI) 
58 //////////////////////////////////////////////////////////////////////////////
59
60
61 // --- ROOT system ---
62 #include "TBenchmark.h"
63                 //#include "TObjectTable.h"
64
65 // --- Standard library ---
66
67 // --- AliRoot header files ---
68 #include "AliLog.h"
69 #include "AliPHOSGeometry.h" 
70 #include "AliPHOSDigit.h"
71 #include "AliPHOSGetter.h"
72 #include "AliPHOSHit.h"
73 #include "AliPHOSSDigitizer.h"
74                 //#include "AliMemoryWatcher.h"
75
76 ClassImp(AliPHOSSDigitizer)
77
78            
79 //____________________________________________________________________________ 
80   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","")
81 {
82   // ctor
83   fFirstEvent = fLastEvent  = 0 ;  
84   fDefaultInit = kTRUE ; 
85 }
86
87 //____________________________________________________________________________ 
88 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName, 
89                                      const char * eventFolderName):
90   TTask("PHOS"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
91   fEventFolderName(eventFolderName)
92 {
93
94   // ctor
95   fFirstEvent = fLastEvent  = 0 ; // runs one event by defaut  
96   InitParameters() ; 
97   Init();
98   fDefaultInit = kFALSE ; 
99 }
100
101 //____________________________________________________________________________ 
102 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer & sd)
103   : TTask(sd)
104 {
105   //cpy ctor 
106
107   fFirstEvent    = sd.fFirstEvent ; 
108   fLastEvent     = sd.fLastEvent ;
109   fA             = sd.fA ;
110   fB             = sd.fB ;
111   fPrimThreshold = sd.fPrimThreshold ;
112   fSDigitsInRun  = sd.fSDigitsInRun ;
113   SetName(sd.GetName()) ; 
114   SetTitle(sd.GetTitle()) ; 
115   fEventFolderName = sd.fEventFolderName;
116 }
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 ;
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     Int_t nPrim =  static_cast<Int_t>((gime->TreeH())->GetEntries()) ; 
208     // Attention nPrim is the number of primaries tracked by Geant 
209     // and this number could be different to the number of Primaries in TreeK;
210     Int_t iprim ;
211
212     for (iprim = 0 ; iprim < nPrim ; iprim ++) { 
213       //=========== Get the PHOS branch from Hits Tree for the Primary iprim
214       gime->Track(iprim) ;
215      Int_t i;
216        for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
217         AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
218         // Assign primary number only if contribution is significant
219         
220         if( hit->GetEnergy() > fPrimThreshold)
221           new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
222                                                  hit->GetEnergy() ,hit->GetTime()) ;
223         else
224           new((*sdigits)[nSdigits]) AliPHOSDigit(-1               ,hit->GetId(), 
225                                                  hit->GetEnergy() ,hit->GetTime()) ;
226         nSdigits++ ;    
227         
228        }
229  
230     } // loop over iprim
231
232     sdigits->Sort() ;
233
234     nSdigits = sdigits->GetEntriesFast() ;
235
236     fSDigitsInRun += nSdigits ;  
237     sdigits->Expand(nSdigits) ;
238
239     Int_t i ;
240     for (i = 0 ; i < nSdigits ; i++) { 
241       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
242       digit->SetIndexInList(i) ;     
243     }
244
245     //Now write SDigits
246
247     
248     //First list of sdigits
249
250     Int_t bufferSize = 32000 ;
251     TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
252
253     sdigitsBranch->Fill() ;
254
255     gime->WriteSDigits("OVERWRITE");
256
257     //Next - SDigitizer
258
259     gime->WriteSDigitizer("OVERWRITE");
260     //gObjectTable->Print() ; 
261   
262     if(strstr(option,"deb"))
263       PrintSDigits(option) ;
264       
265     //memwatcher.Watch(ievent); 
266   }// event loop
267   
268   Unload();
269
270   //  gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
271   
272   if(strstr(option,"tim")){
273     gBenchmark->Stop("PHOSSDigitizer");
274     Info("Exec","   took %f seconds for SDigitizing  %f seconds per event",
275          gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nEvents) ;
276   }
277   
278   //TFile f("out.root","RECREATE");
279   //memwatcher.WriteToFile(); 
280   //f.Close();
281 }
282
283 //__________________________________________________________________
284 void AliPHOSSDigitizer::Print(const Option_t *)const
285 {
286   // Prints parameters of SDigitizer
287   Info("Print", "\n------------------- %s -------------", GetName() ) ; 
288   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
289   printf("   with digitization parameters  A = %f\n", fA) ; 
290   printf("                                 B = %f\n", fB) ;
291   printf("   Threshold for Primary assignment= %f\n", fPrimThreshold)  ; 
292   printf("---------------------------------------------------\n") ;
293   
294 }
295
296 //__________________________________________________________________
297 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
298 {
299   // Equal operator.
300   // SDititizers are equal if their pedestal, slope and threshold are equal
301
302   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
303     return kTRUE ;
304   else
305     return kFALSE ;
306 }
307
308 //__________________________________________________________________
309 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
310 {
311   // Prints list of digits produced in the current pass of AliPHOSDigitizer
312
313
314   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
315   const TClonesArray * sdigits = gime->SDigits() ;
316   
317   Info( "\nPrintSDigits", "event # %d %d sdigits", gAlice->GetEvNumber(), sdigits->GetEntriesFast() ) ; 
318
319   if(strstr(option,"all")||strstr(option,"EMC")){
320     
321     //loop over digits
322     AliPHOSDigit * digit;
323     printf("\nEMC sdigits\n") ; 
324     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
325     Int_t index ;
326     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
327          ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
328       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
329       //  if(digit->GetNprimary() == 0) 
330       //        continue;
331 //       printf("%6d  %8d    %6.5e %4d      %2d :\n", // YVK
332       printf("%6d  %.4f    %6.5e %4d      %2d :\n",
333               digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
334       Int_t iprimary;
335       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
336         printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
337       }  
338     }    
339   }
340
341   if(strstr(option,"all")||strstr(option,"CPV")){
342     
343     //loop over CPV digits
344     AliPHOSDigit * digit;
345     printf("\nCPV sdigits\n") ; 
346     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
347     Int_t index ;
348     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
349       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
350       if(digit->GetId() > maxEmc){
351         printf("\n%6d  %8d    %4d      %2d :",
352                 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;  
353         Int_t iprimary;
354         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
355           printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
356         }
357       }    
358     }
359   }
360 }
361
362 //____________________________________________________________________________ 
363 void AliPHOSSDigitizer::Unload() const
364 {
365   // Unloads the objects from the folder
366   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
367   AliPHOSLoader * loader = gime->PhosLoader() ; 
368   loader->UnloadHits() ; 
369   loader->UnloadSDigits() ; 
370 }