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