]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDigitizer.cxx
Do not write Tree's but do an AutoSave() to fit in the framework
[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     = 0.5e-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   fManager = 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   fManager = 0 ;                     // We work in the standalong mode
115   Init() ;
116   
117 }
118
119 //____________________________________________________________________________ 
120 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * ard):AliDigitizer(ard)
121 {
122   // ctor
123   SetName("");     //Will call init in the digitizing
124   SetTitle("aliroot") ;  
125 }
126
127 //____________________________________________________________________________ 
128   AliPHOSDigitizer::~AliPHOSDigitizer()
129 {
130   // dtor
131
132
133 }
134
135 //____________________________________________________________________________
136 void AliPHOSDigitizer::Digitize(const Int_t event) 
137
138   
139   // Makes the digitization of the collected summable digits.
140   //  It first creates the array of all PHOS modules
141   //  filled with noise (different for EMC, CPV and PPSD) and
142   //  then adds contributions from SDigits. 
143   // This design avoids scanning over the list of digits to add 
144   // contribution to new SDigits only.
145
146   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
147   TClonesArray * digits = gime->Digits(GetName()) ; 
148
149   digits->Clear() ;
150
151   const AliPHOSGeometry *geom = gime->PHOSGeometry() ; 
152
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     cerr << "ERROR: AliPHOSDigitizer::Digitize -> SDigitizer with name " << GetName() << " not found " << endl ; 
169     abort() ; 
170   }
171     
172   // loop through the sdigits posted to the White Board and add them to the noise
173   TCollection * folderslist = gime->SDigitsFolder()->GetListOfFolders() ; 
174   TIter next(folderslist) ; 
175   TFolder * folder = 0 ; 
176   TClonesArray * sdigits = 0 ;
177   Int_t input = 0 ;
178   TObjArray * sdigArray = new TObjArray(2) ;
179   while ( (folder = (TFolder*)next()) ) 
180     if ( (sdigits = (TClonesArray*)folder->FindObject(GetName()) ) ) {
181       TString fileName(folder->GetName()) ;
182       fileName.ReplaceAll("_","/") ;
183       cout << "INFO: AliPHOSDigitizer::Digitize -> Adding SDigits " 
184            << GetName() << " from " << fileName << endl ; 
185       sdigArray->AddAt(sdigits, input) ;
186       input++ ;
187     }
188
189   //Find the first crystall with signal
190   Int_t nextSig = 200000 ; 
191   Int_t i;
192   for(i=0; i<input; i++){
193     sdigits = (TClonesArray *)sdigArray->At(i) ;
194     if ( !sdigits->GetEntriesFast() )
195       continue ; 
196     Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
197      if(curNext < nextSig) 
198        nextSig = curNext ;
199   }
200
201   TArrayI index(input) ;
202   index.Reset() ;  //Set all indexes to zero
203
204   AliPHOSDigit * digit ;
205   AliPHOSDigit * curSDigit ;
206
207   TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
208
209   //Put Noise contribution
210   for(absID = 1; absID <= nEMC; absID++){
211     Float_t noise = gRandom->Gaus(0., fPinNoise) ; 
212     new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
213     //look if we have to add signal?
214     if(absID==nextSig){
215       //Add SDigits from all inputs 
216       digit = (AliPHOSDigit *) digits->At(absID-1) ;
217
218       ticks->Clear() ;
219       Int_t contrib = 0 ;
220       Float_t a = digit->GetAmp() ;
221       Float_t b = TMath::Abs( a /fTimeSignalLength) ;
222       //Mark the beginnign of the signal
223       new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);  
224       //Mark the end of the ignal     
225       new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b); 
226
227       //loop over inputs
228       for(i=0; i<input; i++){
229         if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
230           curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;         
231         else
232           curSDigit = 0 ;
233         //May be several digits will contribute from the same input
234         while(curSDigit && curSDigit->GetId() == absID){           
235           //Shift primary to separate primaries belonging different inputs
236           Int_t primaryoffset ;
237           if(fManager)
238             primaryoffset = fManager->GetMask(i) ; 
239           else
240             primaryoffset = 10000000*i ;
241           curSDigit->ShiftPrimary(primaryoffset) ;
242           
243           a = curSDigit->GetAmp() ;
244           b = a /fTimeSignalLength ;
245           new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);  
246           new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b); 
247
248           *digit = *digit + *curSDigit ;  //add energies
249
250           index[i]++ ;
251           if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
252             curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;       
253           else
254             curSDigit = 0 ;
255         }
256       }
257
258       //calculate and set time
259       Float_t time = FrontEdgeTime(ticks) ;
260       digit->SetTime(time) ;
261
262       //Find next signal module
263       nextSig = 200000 ;
264       for(i=0; i<input; i++){
265         sdigits = ((TClonesArray *)sdigArray->At(i)) ;
266         Int_t curNext = nextSig ;
267         if(sdigits->GetEntriesFast() > index[i] ){
268           curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
269           
270         }
271         if(curNext < nextSig) nextSig = curNext ;
272       }
273     }
274   }
275   
276   ticks->Delete() ;
277   delete ticks ;
278
279
280   //Now CPV digits (different noise and no timing)
281   for(absID = nEMC+1; absID <= nCPV; absID++){
282     Float_t noise = gRandom->Gaus(0., fCPVNoise) ; 
283     new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
284     //look if we have to add signal?
285     if(absID==nextSig){
286       digit = (AliPHOSDigit *) digits->At(absID-1) ;
287       //Add SDigits from all inputs
288       for(i=0; i<input; i++){
289         if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
290           curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;         
291         else
292           curSDigit = 0 ;
293
294         //May be several digits will contribute from the same input
295         while(curSDigit && curSDigit->GetId() == absID){           
296           //Shift primary to separate primaries belonging different inputs
297           Int_t primaryoffset ;
298           if(fManager)
299             primaryoffset = fManager->GetMask(i) ; 
300           else
301             primaryoffset = 10000000*i ;
302           curSDigit->ShiftPrimary(primaryoffset) ;
303
304           //add energies
305           *digit = *digit + *curSDigit ;  
306           index[i]++ ;
307           if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
308             curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;       
309           else
310             curSDigit = 0 ;
311         }
312       }
313
314       //Find next signal module
315       nextSig = 200000 ;
316       for(i=0; i<input; i++){
317         sdigits = (TClonesArray *)sdigArray->At(i) ;
318         Int_t curNext = nextSig ;
319         if(sdigits->GetEntriesFast() > index[i] )
320           curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
321         if(curNext < nextSig) nextSig = curNext ;
322       }
323       
324     }
325   }
326   delete sdigArray ; //We should not delete its contents
327   
328   
329   
330   //remove digits below thresholds
331   for(i = 0; i < nEMC ; i++){
332     digit = (AliPHOSDigit*) digits->At(i) ;
333     if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
334       digits->RemoveAt(i) ;
335     else
336       digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
337   }
338
339
340   for(i = nEMC; i < nCPV ; i++)
341     if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(i))->GetAmp()) < fCPVDigitThreshold)
342       digits->RemoveAt(i) ;
343     
344   digits->Compress() ;  
345   
346   Int_t ndigits = digits->GetEntriesFast() ;
347   digits->Expand(ndigits) ;
348
349   //Set indexes in list of digits and make true digitization of the energy
350   for (i = 0 ; i < ndigits ; i++) { 
351     AliPHOSDigit * digit = (AliPHOSDigit *) digits->At(i) ; 
352     digit->SetIndexInList(i) ;     
353     Float_t energy = sDigitizer->Calibrate(digit->GetAmp()) ;
354     digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
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 void AliPHOSDigitizer::Exec(Option_t *option) 
374
375   // Managing method
376
377   if(strcmp(GetName(), "") == 0 )   
378     Init() ;
379   
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   else {
399     gAlice->GetEvent(0) ;
400     nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
401     treeD=gAlice->TreeD() ;
402   }
403   if(treeD == 0 ){
404     cerr << " AliPHOSDigitizer :: Can not find TreeD " << endl ;
405     return ;
406   }
407
408   //Check, if this branch already exits
409   TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
410   TIter next(lob) ; 
411   TBranch * branch = 0 ;  
412   Bool_t phosfound = kFALSE, digitizerfound = kFALSE ; 
413   
414   while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
415     if ( (strcmp(branch->GetName(), "PHOS")==0) && 
416          (strcmp(branch->GetTitle(), GetName())==0) ) 
417       phosfound = kTRUE ;
418     
419     else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) && 
420               (strcmp(branch->GetTitle(), GetName())==0) ) 
421       digitizerfound = kTRUE ; 
422   }
423
424   if ( phosfound ) {
425     cerr << "WARNING: AliPHOSDigitizer -> Digits branch with name " << GetName() 
426          << " already exits" << endl ;
427     return ; 
428   }   
429   if ( digitizerfound ) {
430     cerr << "WARNING: AliPHOSDigitizer -> Digitizer branch with name " << GetName() 
431          << " already exits" << endl ;
432     return ; 
433   }   
434
435   Int_t ievent ;
436
437   for(ievent = 0; ievent < nevents; ievent++){
438     
439     if(fManager){
440       Int_t input ;
441       for(input = 0 ; input < fManager->GetNinputs(); input ++){
442         TTree * treeS = fManager->GetInputTreeS(input) ;
443         if(!treeS){
444           cerr << "AliPHOSDigitizer -> No Input " << endl ;
445           return ;
446         }
447         gime->ReadTreeS(treeS,input) ;
448       }
449     }
450     else
451       gime->Event(ievent,"S") ;
452     
453     Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
454     
455     WriteDigits(ievent) ;
456     
457     if(strstr(option,"deb"))
458       PrintDigits(option);
459     
460     //increment the total number of Digits per run 
461     fDigitsInRun += gime->Digits()->GetEntriesFast() ;  
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     = 0.5e-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   if(fManager)
514     SetTitle("aliroot") ;
515   else if (strcmp(GetTitle(),"")==0) 
516    SetTitle("galice.root") ;
517
518   if( strcmp(GetName(), "") == 0 )
519     SetName("Default") ;
520   
521   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ; 
522   if ( gime == 0 ) {
523     cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ; 
524     return kFALSE;
525   } 
526   
527   const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
528   fEmcCrystals = geom->GetNModules() *  geom->GetNCristalsInModule() ;
529   
530   // Post Digits to the white board
531   gime->PostDigits(GetName() ) ;   
532   
533   // Post Digitizer to the white board
534   gime->PostDigitizer(this) ;
535   
536   //Mark that we will use current header file
537   if(!fManager){
538     gime->PostSDigits(GetName(),GetTitle()) ;
539     gime->PostSDigitizer(GetName(),GetTitle()) ;
540   }
541   return kTRUE ;
542 }
543
544 //__________________________________________________________________
545 void AliPHOSDigitizer::MixWith(const char* headerFile)
546 {
547   // Allows to produce digits by superimposing background and signal event.
548   // It is assumed, that headers file with SIGNAL events is opened in 
549   // the constructor. 
550   // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed 
551   // Thus we avoid writing (changing) huge and expensive 
552   // backgound files: all output will be writen into SIGNAL, i.e. 
553   // opened in constructor file. 
554   //
555   // One can open as many files to mix with as one needs.
556   // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
557   // can be mixed.
558
559   if( strcmp(GetName(), "") == 0 )
560     Init() ;
561
562   if(fManager){
563     cout << "Can not use this method under AliRunDigitizer " << endl ;
564     return ;
565   }
566   
567   // check if the specified SDigits do not already exist on the White Board:
568   // //Folders/RunMC/Event/Data/PHOS/SDigits/headerFile/sdigitsname
569
570   TString path = "Folders/RunMC/Event/Data/PHOS/SDigits" ; 
571   path += headerFile ; 
572   path += "/" ; 
573   path += GetName() ;
574   if ( gROOT->FindObjectAny(path.Data()) ) {
575     cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
576     return;
577   }
578
579   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
580   gime->PostSDigits(GetName(),headerFile) ;
581   
582   // check if the requested file is already open or exist and if SDigits Branch exist
583   TFile * file = (TFile*)gROOT->FindObject(headerFile); 
584   if ( !file ) { 
585     file = new TFile(headerFile, "READ") ; 
586     if (!file) { 
587       cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ; 
588       return ; 
589     }
590   }
591   
592 }
593
594 //__________________________________________________________________
595 void AliPHOSDigitizer::Print(Option_t* option)const {
596   // Print Digitizer's parameters
597   if( strcmp(GetName(), "") != 0 ){
598     
599     cout << "------------------- "<< GetName() << " -------------" << endl ;
600     cout << "Digitizing sDigits from file(s): " <<endl ;
601     
602      TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("Folders/RunMC/Event/Data/PHOS/SDigits"))->GetListOfFolders() ; 
603     TIter next(folderslist) ; 
604     TFolder * folder = 0 ; 
605     
606     while ( (folder = (TFolder*)next()) ) {
607       if ( folder->FindObject(GetName())  ) 
608         cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ; 
609     }
610     cout << endl ;
611     cout << "Writing digits to " << GetTitle() << endl ;
612     
613     cout << endl ;
614     cout << "With following parameters: " << endl ;
615     cout << "     Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
616     cout << "  Threshold  in EMC  (fEMCDigitThreshold) = " << fEMCDigitThreshold  << endl ; ;
617     cout << "                 Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ; 
618     cout << "    Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ; 
619     cout << "---------------------------------------------------" << endl ;
620   }
621   else
622     cout << "AliPHOSDigitizer not initialized " << endl ;
623   
624 }
625
626 //__________________________________________________________________
627 void AliPHOSDigitizer::PrintDigits(Option_t * option){
628   // Print a table of digits
629
630   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
631   TClonesArray * digits = gime->Digits() ; 
632
633   cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
634   cout << "       Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
635   cout << endl ;
636   if(strstr(option,"all")||strstr(option,"EMC")){
637     
638     //loop over digits
639     AliPHOSDigit * digit;
640     cout << "EMC digits (with primaries): " << endl ;
641     cout << "Digit Id     Amplitude   Index       Nprim  Primaries list " <<  endl;      
642     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
643     Int_t index ;
644     for (index = 0 ; (index < digits->GetEntriesFast()) && 
645          (((AliPHOSDigit * )  digits->At(index))->GetId() <= maxEmc) ; index++) {
646       digit = (AliPHOSDigit * )  digits->At(index) ;
647       if(digit->GetNprimary() == 0) continue;
648       cout << setw(6)  <<  digit->GetId() << "   "  <<  setw(10)  <<  digit->GetAmp() <<   "    "  
649            << setw(6)  <<  digit->GetIndexInList() << "    "   
650            << setw(5)  <<  digit->GetNprimary() <<"    ";
651       
652       Int_t iprimary;
653       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
654         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
655       cout << endl;      
656     }    
657     cout << endl;
658   }
659
660   if(strstr(option,"all")||strstr(option,"CPV")){
661     
662     //loop over CPV digits
663     AliPHOSDigit * digit;
664     cout << "CPV digits: " << endl ;
665     cout << "Digit Id     Amplitude   Index       Nprim  Primaries list " <<  endl;      
666     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
667     Int_t index ;
668     for (index = 0 ; index < digits->GetEntriesFast(); index++) {
669       digit = (AliPHOSDigit * )  digits->At(index) ;
670       if(digit->GetId() > maxEmc){
671         cout << setw(6)  <<  digit->GetId() << "   "  <<        setw(10)  <<  digit->GetAmp() <<   "    "  
672              << setw(6)  <<  digit->GetIndexInList() << "    "   
673              << setw(5)  <<  digit->GetNprimary() <<"    ";
674         
675         Int_t iprimary;
676         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
677           cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
678         cout << endl;    
679       }    
680     }
681   }
682
683 }
684
685 //__________________________________________________________________
686 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
687 {
688   // we set title (comment) of the SDigits branch in the first! header file
689   if( strcmp(GetName(), "") == 0 )
690     Init() ;
691
692   AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ; 
693  
694 }
695 //__________________________________________________________________
696 Float_t AliPHOSDigitizer::TimeOfNoise(void)
697 {  // Calculates the time signal generated by noise
698   //to be rewritten, now returns just big number
699   return 1. ;
700
701 }
702 //____________________________________________________________________________
703 void AliPHOSDigitizer::Reset() 
704
705   // sets current event number to the first simulated event
706
707   if( strcmp(GetName(), "") == 0 )
708     Init() ;
709
710  //  Int_t inputs ;
711 //   for(inputs = 0; inputs < fNinputs ;inputs++)
712 //       fIevent->AddAt(-1, inputs ) ;
713   
714 }
715
716 //____________________________________________________________________________
717 void AliPHOSDigitizer::WriteDigits(Int_t event)
718 {
719
720   // Makes TreeD in the output file. 
721   // Check if branch already exists: 
722   //   if yes, exit without writing: ROOT TTree does not support overwriting/updating of 
723   //      already existing branches. 
724   //   else creates branch with Digits, named "PHOS", title "...",
725   //      and branch "AliPHOSDigitizer", with the same title to keep all the parameters
726   //      and names of files, from which digits are made.
727
728   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
729   const TClonesArray * digits = gime->Digits(GetName()) ; 
730  TTree * treeD ;
731
732   if(fManager)
733     treeD = fManager->GetTreeD() ;
734  else
735     treeD = gAlice->TreeD();
736   
737   // create new branches
738   // -- generate file name if necessary
739   char * file =0;
740   if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
741     file = new char[strlen(gAlice->GetBaseFile())+20] ;
742     sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
743   }
744
745   TDirectory *cwd = gDirectory;
746   // -- create Digits branch
747   Int_t bufferSize = 32000 ;    
748   TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
749   digitsBranch->SetTitle(GetName());
750   if (file) {
751     digitsBranch->SetFile(file);
752     TIter next( digitsBranch->GetListOfBranches());
753     TBranch * sbr ;
754     while ((sbr=(TBranch*)next())) {
755       sbr->SetFile(file);
756     }   
757     cwd->cd();
758   } 
759     
760   // -- Create Digitizer branch
761   Int_t splitlevel = 0 ;
762   AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
763   TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel); 
764   digitizerBranch->SetTitle(GetName());
765   if (file) {
766     digitizerBranch->SetFile(file);
767     TIter next( digitizerBranch->GetListOfBranches());
768     TBranch * sbr;
769     while ((sbr=(TBranch*)next())) {
770       sbr->SetFile(file);
771     }   
772     cwd->cd();
773   }
774   
775   digitsBranch->Fill() ;
776   treeD->AutoSave() ; // Write(0,kOverwrite) ;  
777  
778 }