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