]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDigitizer.cxx
coding convention corrections
[u/mrichter/AliRoot.git] / PHOS / AliPHOSDigitizer.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 /* $Id$ */
17
18 //_________________________________________________________________________
19 //*-- Author :  Dmitri Peressounko (SUBATECH & Kurchatov Institute) 
20 //////////////////////////////////////////////////////////////////////////////
21 // This TTask performs digitization of Summable digits (in the PHOS case it is just
22 // the sum of contributions from all primary particles into a given cell). 
23 // In addition it performs mixing of summable digits from different events.
24 // The name of the TTask is also the title of the branch that will contain 
25 // the created SDigits
26 // The title of the TTAsk is the name of the file that contains the hits from
27 // which the SDigits are created
28 //
29 // For each event two branches are created in TreeD:
30 //   "PHOS" - list of digits
31 //   "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
32 //
33 // Note, that one can set a title for new digits branch, and repeat digitization with
34 // another set of parameters.
35 //
36 // Use case:
37 // root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
38 // root[1] d->ExecuteTask()             
39 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
40 //                       //Digitizes SDigitis in all events found in file galice.root 
41 //
42 // root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;  
43 //                       // Will read sdigits from galice1.root
44 // root[3] d1->MixWith("galice2.root")       
45 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
46 //                       // Reads another set of sdigits from galice2.root
47 // root[3] d1->MixWith("galice3.root")       
48 //                       // Reads another set of sdigits from galice3.root
49 // root[4] d->ExecuteTask("deb timing")    
50 //                       // Reads SDigits from files galice1.root, galice2.root ....
51 //                       // mixes them and stores produced Digits in file galice1.root          
52 //                       // deb - prints number of produced digits
53 //                       // deb all - prints list of produced digits
54 //                       // timing  - prints time used for digitization
55 //
56
57 // --- ROOT system ---
58 #include "TFile.h"
59 #include "TTree.h"
60 #include "TSystem.h"
61 #include "TROOT.h"
62 #include "TFolder.h"
63 #include "TObjString.h"
64 #include "TGeometry.h"
65 #include "TBenchmark.h"
66
67 // --- Standard library ---
68
69 // --- AliRoot header files ---
70 #include "AliRun.h"
71 #include "AliHeader.h"
72 #include "AliStream.h"
73 #include "AliRunDigitizer.h"
74 #include "AliPHOSDigit.h"
75 #include "AliPHOS.h"
76 #include "AliPHOSGetter.h"
77 #include "AliPHOSDigitizer.h"
78 #include "AliPHOSSDigitizer.h"
79 #include "AliPHOSGeometry.h"
80 #include "AliPHOSTick.h"
81
82 ClassImp(AliPHOSDigitizer)
83
84
85 //____________________________________________________________________________ 
86   AliPHOSDigitizer::AliPHOSDigitizer() 
87 {
88   // ctor
89   InitParameters() ; 
90   fDefaultInit = kTRUE ;
91   fManager = 0 ;                     // We work in the standalong mode
92
93  }
94
95 //____________________________________________________________________________ 
96 AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile, const char * name, const Bool_t toSplit)
97 {
98   // ctor
99
100   SetTitle(headerFile) ;
101   SetName(name) ;
102   fManager = 0 ;                     // We work in the standalong mode
103   fSplitFile= 0 ; 
104   InitParameters() ; 
105   fToSplit = toSplit ;
106   Init() ;
107   fDefaultInit = kFALSE ; 
108 }
109
110 //____________________________________________________________________________ 
111 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * ard):AliDigitizer(ard)
112 {
113   // ctor
114   SetTitle(ard->GetInputFileName(0,0)) ;
115   InitParameters() ; 
116   fDefaultInit = kFALSE ; 
117   fSplitFile   = 0 ; 
118
119   if (ard->GetOutputFile()) {
120     SetName(ard->GetOutputFile().Data());
121     fToSplit = kTRUE ;
122   } else {
123     SetName("Default") ;
124     fToSplit = kFALSE ;
125   }
126 }
127
128 //____________________________________________________________________________ 
129   AliPHOSDigitizer::~AliPHOSDigitizer()
130 {
131   // dtor
132      
133    fSplitFile = 0 ; 
134 }
135
136 //____________________________________________________________________________
137 void AliPHOSDigitizer::Digitize(const Int_t event) 
138
139   
140   // Makes the digitization of the collected summable digits.
141   //  It first creates the array of all PHOS modules
142   //  filled with noise (different for EMC, CPV and PPSD) and
143   //  then adds contributions from SDigits. 
144   // This design avoids scanning over the list of digits to add 
145   // contribution to new SDigits only.
146
147   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
148   TClonesArray * digits = gime->Digits(GetName()) ; 
149
150   digits->Clear() ;
151
152   const AliPHOSGeometry *geom = gime->PHOSGeometry() ; 
153   //Making digits with noise, first EMC
154   Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
155   
156   Int_t nCPV ;
157   Int_t absID ;
158   TString name      =  geom->GetName() ;
159   
160   nCPV = nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
161     geom->GetNModules() ;
162
163   digits->Expand(nCPV) ;
164
165   // get first the sdigitizer from the tasks list (must have same name as the digitizer)
166   const AliPHOSSDigitizer * sDigitizer = gime->SDigitizer(GetName()); 
167   if ( !sDigitizer) {
168     Fatal("Digitize", "SDigitizer with name %s not found", GetName() ) ; 
169   }
170     
171   // loop through the sdigits posted to the White Board and add them to the noise
172   TCollection * folderslist = gime->SDigitsFolder()->GetListOfFolders() ; 
173   TIter next(folderslist) ; 
174   TFolder * folder = 0 ; 
175   TClonesArray * sdigits = 0 ;
176   Int_t input = 0 ;
177   TObjArray * sdigArray = new TObjArray(2) ;
178   while ( (folder = (TFolder*)next()) ) { 
179     if ( (sdigits = (TClonesArray*)folder->FindObject(GetName()) ) ) {
180       TString fileName(folder->GetName()) ;
181       fileName.ReplaceAll("_","/") ;
182 //       Info("Digitize", "Adding SDigits %s from %s", GetName(), fileName() ) ;  
183       sdigArray->AddAt(sdigits, input) ;
184       input++ ;
185     }
186   }
187
188   //Find the first crystall with signal
189   Int_t nextSig = 200000 ; 
190   Int_t i;
191   for(i=0; i<input; i++){
192     sdigits = (TClonesArray *)sdigArray->At(i) ;
193     if ( !sdigits->GetEntriesFast() )
194       continue ; 
195     Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
196      if(curNext < nextSig) 
197        nextSig = curNext ;
198   }
199
200   TArrayI index(input) ;
201   index.Reset() ;  //Set all indexes to zero
202
203   AliPHOSDigit * digit ;
204   AliPHOSDigit * curSDigit ;
205
206   TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
207
208   //Put Noise contribution
209   for(absID = 1; absID <= nEMC; absID++){
210     Float_t noise = gRandom->Gaus(0., fPinNoise) ; 
211     new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
212     //look if we have to add signal?
213     digit = (AliPHOSDigit *) digits->At(absID-1) ;
214  
215     if(absID==nextSig){
216       //Add SDigits from all inputs 
217       ticks->Clear() ;
218       Int_t contrib = 0 ;
219       Float_t a = digit->GetAmp() ;
220       Float_t b = TMath::Abs( a /fTimeSignalLength) ;
221       //Mark the beginnign of the signal
222       new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);  
223       //Mark the end of the ignal     
224       new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b); 
225
226       //loop over inputs
227       for(i=0; i<input; i++){
228         if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
229           curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;         
230         else
231           curSDigit = 0 ;
232         //May be several digits will contribute from the same input
233         while(curSDigit && curSDigit->GetId() == absID){           
234           //Shift primary to separate primaries belonging different inputs
235           Int_t primaryoffset ;
236           if(fManager)
237             primaryoffset = fManager->GetMask(i) ; 
238           else
239             primaryoffset = 10000000*i ;
240           curSDigit->ShiftPrimary(primaryoffset) ;
241           
242           a = curSDigit->GetAmp() ;
243           b = a /fTimeSignalLength ;
244           new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);  
245           new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b); 
246
247           *digit = *digit + *curSDigit ;  //add energies
248
249           index[i]++ ;
250           if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
251             curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;       
252           else
253             curSDigit = 0 ;
254         }
255       }
256
257       //calculate and set time
258       Float_t time = FrontEdgeTime(ticks) ;
259       digit->SetTime(time) ;
260
261       //Find next signal module
262       nextSig = 200000 ;
263       for(i=0; i<input; i++){
264         sdigits = ((TClonesArray *)sdigArray->At(i)) ;
265         Int_t curNext = nextSig ;
266         if(sdigits->GetEntriesFast() > index[i] ){
267           curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
268           
269         }
270         if(curNext < nextSig) nextSig = curNext ;
271       }
272     }
273   }
274   
275   ticks->Delete() ;
276   delete ticks ;
277
278
279   //Now CPV digits (different noise and no timing)
280   for(absID = nEMC+1; absID <= nCPV; absID++){
281     Float_t noise = gRandom->Gaus(0., fCPVNoise) ; 
282     new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
283     //look if we have to add signal?
284     if(absID==nextSig){
285       digit = (AliPHOSDigit *) digits->At(absID-1) ;
286       //Add SDigits from all inputs
287       for(i=0; i<input; i++){
288         if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
289           curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;         
290         else
291           curSDigit = 0 ;
292
293         //May be several digits will contribute from the same input
294         while(curSDigit && curSDigit->GetId() == absID){           
295           //Shift primary to separate primaries belonging different inputs
296           Int_t primaryoffset ;
297           if(fManager)
298             primaryoffset = fManager->GetMask(i) ; 
299           else
300             primaryoffset = 10000000*i ;
301           curSDigit->ShiftPrimary(primaryoffset) ;
302
303           //add energies
304           *digit = *digit + *curSDigit ;  
305           index[i]++ ;
306           if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
307             curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;       
308           else
309             curSDigit = 0 ;
310         }
311       }
312
313       //Find next signal module
314       nextSig = 200000 ;
315       for(i=0; i<input; i++){
316         sdigits = (TClonesArray *)sdigArray->At(i) ;
317         Int_t curNext = nextSig ;
318         if(sdigits->GetEntriesFast() > index[i] )
319           curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
320         if(curNext < nextSig) nextSig = curNext ;
321       }
322       
323     }
324   }
325   delete sdigArray ; //We should not delete its contents
326   
327   
328   
329   //remove digits below thresholds
330   for(i = 0; i < nEMC ; i++){
331     digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
332     if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
333       digits->RemoveAt(i) ;
334     else
335       digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
336   }
337
338
339   for(i = nEMC; i < nCPV ; i++)
340     if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(i))->GetAmp()) < fCPVDigitThreshold)
341       digits->RemoveAt(i) ;
342     
343   digits->Compress() ;  
344   
345   Int_t ndigits = digits->GetEntriesFast() ;
346   digits->Expand(ndigits) ;
347
348   //Set indexes in list of digits and make true digitization of the energy
349   for (i = 0 ; i < ndigits ; i++) { 
350     digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ; 
351     digit->SetIndexInList(i) ;     
352     Float_t energy = sDigitizer->Calibrate(digit->GetAmp()) ;
353     digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
354   }
355
356 }
357
358 //____________________________________________________________________________
359 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
360 {
361   Int_t chanel ;
362   if(absId <= fEmcCrystals){ //digitize as EMC 
363     chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;       
364     if(chanel > fNADCemc ) chanel =  fNADCemc ;
365   }
366   else{ //Digitize as CPV
367     chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;       
368     if(chanel > fNADCcpv ) chanel =  fNADCcpv ;
369   }
370   return chanel ;
371 }
372
373 //____________________________________________________________________________
374 void AliPHOSDigitizer::Exec(Option_t *option) 
375
376   // Managing method
377
378   if(strcmp(GetName(), "") == 0 )   
379     Init() ;
380   if (strstr(option,"print")) {
381     Print("");
382     return ; 
383   }
384   
385   if(strstr(option,"tim"))
386     gBenchmark->Start("PHOSDigitizer");
387
388   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
389   
390   Int_t nevents ;
391   
392   TTree * treeD ;
393   
394   if(fManager){
395     treeD = fManager->GetTreeD() ;
396     nevents = 1 ;    // Will process only one event
397     
398     //Check, if this branch already exits
399     if (treeD) { 
400       TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
401       TIter next(lob) ; 
402       TBranch * branch = 0 ;  
403       Bool_t phosfound = kFALSE, digitizerfound = kFALSE ; 
404       
405       while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
406         if ( (strcmp(branch->GetName(), "PHOS")==0) && 
407              (strcmp(branch->GetTitle(), GetName())==0) ) 
408           phosfound = kTRUE ;
409         
410         else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) && 
411                   (strcmp(branch->GetTitle(), GetName())==0) ) 
412           digitizerfound = kTRUE ; 
413       }
414       
415       if ( phosfound ) {
416         Error("Exec", "Digits branch with name %s already exists", GetName() ) ;
417         return ; 
418       }   
419       if ( digitizerfound ) {
420         Error("Exec", "Digitizer branch with name %s already exists", GetName() );
421         return ; 
422       }
423     }   
424   }
425   else { //PHOS standalone
426     if(gime->BranchExists("Digits") ) 
427       return ;
428     nevents=gime->MaxEvent() ;
429   }
430   
431   Int_t ievent ;
432   
433   for(ievent = 0; ievent < nevents; ievent++){
434   
435     if(fManager){
436
437       Int_t input ;
438       for(input = 0 ; input < fManager->GetNinputs(); input ++){
439         TTree * treeS = fManager->GetInputTreeS(input) ;
440         if(!treeS){
441           Error("Exec", "No Input") ;  
442           return ;
443         }
444         gime->ReadTreeS(treeS,input) ;
445       }
446
447     }
448     else
449       gime->Event(ievent,"S") ;
450     
451     Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
452     
453     WriteDigits(ievent) ;
454     
455     if(strstr(option,"deb"))
456       PrintDigits(option);
457     
458     //increment the total number of Digits per run 
459     fDigitsInRun += gime->Digits()->GetEntriesFast() ;  
460   }
461    
462   if(strstr(option,"tim")){
463     gBenchmark->Stop("PHOSDigitizer");
464     TString message ; 
465     message = "  took %f seconds for Digitizing %f seconds per event\n" ; 
466     Info("Exec", message.Data(), 
467          gBenchmark->GetCpuTime("PHOSDigitizer"), 
468          gBenchmark->GetCpuTime("PHOSDigitizer")/nevents ); 
469   } 
470 }
471
472 //____________________________________________________________________________ 
473 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks) 
474 { // 
475   ticks->Sort() ; //Sort in accordance with times of ticks
476   TIter it(ticks) ;
477   AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
478   Float_t time = ctick->CrossingTime(fTimeThreshold) ;    
479
480   AliPHOSTick * t ;  
481   while((t=(AliPHOSTick*) it.Next())){
482     if(t->GetTime() < time)  //This tick starts before crossing
483       *ctick+=*t ;
484     else
485       return time ;
486
487     time = ctick->CrossingTime(fTimeThreshold) ;    
488   }
489   return time ;
490 }
491
492 //____________________________________________________________________________ 
493 Bool_t AliPHOSDigitizer::Init()
494 {
495   // Makes all memory allocations
496
497   if( strcmp(GetTitle(), "") == 0 )
498     SetTitle("galice.root") ;
499   
500   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName(), fToSplit) ; 
501   if ( gime == 0 ) {
502     Error("Init", "Could not obtain the Getter object !") ; 
503     return kFALSE;
504   } 
505   
506   const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
507
508   fEmcCrystals = geom->GetNModules() *  geom->GetNCristalsInModule() ;
509   
510   // Post Digits to the white board
511   gime->PostDigits(GetName() ) ;   
512   
513   // Post Digitizer to the white board
514   gime->PostDigitizer(this) ;
515
516   fSplitFile = 0 ;
517   if(fToSplit){
518     // construct the name of the file as /path/PHOS.SDigits.root
519     //First - extract full path if necessary
520     TString digitsFileName(GetTitle()) ;
521     Ssiz_t islash = digitsFileName.Last('/') ;
522     if(islash<digitsFileName.Length())
523       digitsFileName.Remove(islash+1,digitsFileName.Length()) ;
524     else
525       digitsFileName="" ;
526     // Next - append the file name 
527     digitsFileName+="PHOS.Digits." ;
528     if((strcmp(GetName(),"Default")!=0)&&(strcmp(GetName(),"")!=0)){
529       digitsFileName+=GetName() ;
530       digitsFileName+="." ;
531     }
532     digitsFileName+="root" ;
533     // Finally - check if the file already opened or open the file
534     fSplitFile = static_cast<TFile*>(gROOT->GetFile(digitsFileName.Data()));   
535     if(!fSplitFile)
536       fSplitFile =  TFile::Open(digitsFileName.Data(),"update") ;
537   }
538   
539   //Mark that we will use current header file
540   if(!fManager){
541     gime->PostSDigits(GetName(),GetTitle()) ;
542     gime->PostSDigitizer(GetName(),GetTitle()) ;
543   }
544   return kTRUE ;
545 }
546
547 //____________________________________________________________________________ 
548 void AliPHOSDigitizer::InitParameters()
549 {
550   fPinNoise           = 0.004 ;
551   fEMCDigitThreshold  = 0.012 ;
552   fCPVNoise           = 0.01;
553   fCPVDigitThreshold  = 0.09 ;
554   fTimeResolution     = 0.5e-9 ;
555   fTimeSignalLength   = 1.0e-9 ;
556   fDigitsInRun  = 0 ; 
557   fADCchanelEmc = 0.0015;        // width of one ADC channel in GeV
558   fADCpedestalEmc = 0.005 ;      //
559   fNADCemc = (Int_t) TMath::Power(2,16) ;  // number of channels in EMC ADC
560
561   fADCchanelCpv = 0.0012 ;          // width of one ADC channel in CPV 'popugais'
562   fADCpedestalCpv = 0.012 ;         // 
563   fNADCcpv = (Int_t) TMath::Power(2,12);      // number of channels in CPV ADC
564
565   fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
566     
567 }
568
569 //__________________________________________________________________
570 void AliPHOSDigitizer::MixWith(const char* headerFile,const char* brname)
571 {
572   // Allows to produce digits by superimposing background and signal event.
573   // It is assumed, that headers file with SIGNAL events is opened in 
574   // the constructor. 
575   // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed 
576   // Thus we avoid writing (changing) huge and expensive 
577   // backgound files: all output will be writen into SIGNAL, i.e. 
578   // opened in constructor file. 
579   //
580   // One can open as many files to mix with as one needs.
581   // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
582   // can be mixed.
583
584   if( strcmp(GetName(), "") == 0 )
585     Init() ;
586
587   if(fManager){
588     Warning("MixWith", "Can not use this method under AliRunDigitizer\n" ) ;
589     return ;
590   }
591   
592   // check if the specified SDigits do not already exist on the White Board:
593   // //Folders/RunMC/Event/Data/PHOS/SDigits/headerFile/sdigitsname
594
595   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
596   TString path = gime->SDigitsFolder()->GetName() ; 
597
598   // before it was ???? "Folders/RunMC/Event/Data/PHOS/SDigits" ; 
599   path += headerFile ; 
600   path += "/" ; 
601   if(brname)
602      path += brname ; 
603   else
604     path += GetName() ;
605   if ( gROOT->FindObjectAny(path.Data()) ) {
606     Warning("MixWith", "Entry already exists, do not add\n" ) ;
607     return;
608   }
609   gime->PostSDigits(GetName(),headerFile) ;
610   // check if the requested file is already open or exist and if SDigits Branch exist
611   TFile * file = (TFile*)gROOT->GetFile(headerFile); 
612   if ( !file ) { 
613     file = new TFile(headerFile, "READ") ; 
614     if (!file) { 
615       Error("MixWith", "File %s does not exist\n", headerFile ); 
616       return ; 
617     }
618   }
619 }
620
621 //__________________________________________________________________
622 void AliPHOSDigitizer::Print(Option_t* option)const 
623 {
624   // Print Digitizer's parameters
625   if( strcmp(GetName(), "") != 0 ){
626
627     TString message("\n-------------------") ; 
628     message += GetName() ; 
629     message += "-------------\n" ;
630
631     const Int_t nStreams = GetNInputStreams() ; 
632     if (nStreams) {
633       Int_t index = 0 ;  
634       for (index = 0 ; index < nStreams ; index++) {  
635         message += "Adding SDigits " ; 
636         message += GetName() ;
637         message += "from " ; 
638         message += fManager->GetInputFileName(index, 0) ; 
639       }
640       message += "\nWriting digits to " ; 
641       message += fManager->GetInputFileName(0, 0) ;   
642     } else { 
643       //       AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;  
644       //       gime->Folder("sdigits")  ;
645       //       TCollection * folderslist = gime->Folder("sdigits")->GetListOfFolders() ; 
646       //       TIter next(folderslist) ; 
647       //       TFolder * folder = 0 ; 
648       
649       //       while ( (folder = (TFolder*)next()) ) {
650       //        if ( folder->FindObject(GetName())  ) 
651       //      }
652       message += "\nWriting digits to " ;
653       message += GetTitle() ;
654     }    
655     message += "\nWith following parameters:\n" ;
656     message += "     Electronics noise in EMC (fPinNoise) = %f\n" ; 
657     message += "  Threshold  in EMC  (fEMCDigitThreshold) = %f\n" ; 
658     message += "                 Noise in CPV (fCPVNoise) = %f\n" ;
659     message += "    Threshold in CPV (fCPVDigitThreshold) = %f\n" ;
660     message += "---------------------------------------------------\n" ; 
661     Info("Print", message.Data(),  
662          fPinNoise, 
663          fEMCDigitThreshold, 
664          fCPVNoise, 
665          fCPVDigitThreshold ) ; 
666   }
667   else
668     Info("Print", "AliPHOSDigitizer not initialized\n" ) ;
669   
670 }
671  
672 //__________________________________________________________________
673  void AliPHOSDigitizer::PrintDigits(Option_t * option)
674 {
675   // Print a table of digits
676   
677   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
678   TClonesArray * digits = gime->Digits() ; 
679   
680   TString message ; 
681   message  = "\nevent " ;
682   message += gAlice->GetEvNumber() ;
683   message += "\n       Number of entries in Digits list " ;  
684   message += digits->GetEntriesFast()  ;  
685   char * tempo = new char[8192]; 
686  
687   if(strstr(option,"all")||strstr(option,"EMC")){  
688     //loop over digits
689     AliPHOSDigit * digit;
690     message += "\nEMC digits (with primaries):\n"  ;
691     message += "\n   Id  Amplitude    Time          Index Nprim: Primaries list \n" ;    
692     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
693     Int_t index ;
694     for (index = 0 ; (index < digits->GetEntriesFast()) && 
695            (((AliPHOSDigit * )  digits->At(index))->GetId() <= maxEmc) ; index++) {
696       digit = (AliPHOSDigit * )  digits->At(index) ;
697       if(digit->GetNprimary() == 0) 
698         continue;
699       sprintf(tempo, "\n%6d  %8d    %6.5e %4d      %2d :",
700               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
701       message += tempo ; 
702       Int_t iprimary;
703       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
704         sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
705         message += tempo ; 
706       }    
707     }
708   }
709   
710   if(strstr(option,"all")||strstr(option,"CPV")){
711     
712     //loop over CPV digits
713     AliPHOSDigit * digit;
714     message += "\nCPV digits:\n" ;
715     message += "\n   Id  Amplitude    Index Nprim: Primaries list \n" ;    
716
717     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
718     Int_t index ;
719     for (index = 0 ; index < digits->GetEntriesFast(); index++) {
720       digit = (AliPHOSDigit * )  digits->At(index) ;
721       if(digit->GetId() > maxEmc){
722         sprintf(tempo, "\n%6d  %8d    %4d      %2d :",
723                 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;  
724         message += tempo ; 
725         Int_t iprimary;
726         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
727           sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
728           message += tempo ; 
729         }    
730       }
731     }
732   }
733   delete tempo ; 
734   Info("Print", message.Data() ) ; 
735 }
736
737 //__________________________________________________________________
738 Float_t AliPHOSDigitizer::TimeOfNoise(void)
739 {  // Calculates the time signal generated by noise
740   //to be rewritten, now returns just big number
741   return 1. ;
742
743 }
744
745 //____________________________________________________________________________
746 void AliPHOSDigitizer::WriteDigits(Int_t event)
747 {
748
749   // Makes TreeD in the output file. 
750   // Check if branch already exists: 
751   //   if yes, exit without writing: ROOT TTree does not support overwriting/updating of 
752   //      already existing branches. 
753   //   else creates branch with Digits, named "PHOS", title "...",
754   //      and branch "AliPHOSDigitizer", with the same title to keep all the parameters
755   //      and names of files, from which digits are made.
756
757   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
758   const TClonesArray * digits = gime->Digits(GetName()) ; 
759   TTree * treeD ;
760   
761   if(fManager) 
762     treeD = fManager->GetTreeD() ;
763   else {
764     if((gAlice->TreeD() == 0) || (fSplitFile))  // we should not create TreeD if it is already here
765       gAlice->MakeTree("D", fSplitFile);       // We overwrite TreeD in split file in the case of second reconstruction
766     if(fSplitFile)
767       fSplitFile->cd() ;
768     treeD = gAlice->TreeD();
769   }
770   
771   // -- create Digits branch
772   Int_t bufferSize = 32000 ;    
773   TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
774   digitsBranch->SetTitle(GetName());
775   
776   // -- Create Digitizer branch
777   Int_t splitlevel = 0 ;
778   const AliPHOSDigitizer * d = dynamic_cast<const AliPHOSDigitizer *>(gime->Digitizer()) ;
779   TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel); 
780   digitizerBranch->SetTitle(GetName());
781
782   digitsBranch->Fill() ;
783   digitizerBranch->Fill() ; 
784   treeD->AutoSave() ; 
785 }