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