]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDigitizer.cxx
0ac1e94b6e1093041d959799730f26da8f136cb1
[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 }
109
110 //____________________________________________________________________________ 
111 AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile,const char * name)
112 {
113   // ctor
114   SetName(name) ;
115   SetTitle(headerFile) ;
116   fManager = 0 ;                     // We work in the standalong mode
117   Init() ;
118   
119 }
120
121 //____________________________________________________________________________ 
122 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * ard):AliDigitizer(ard)
123 {
124   // ctor
125   SetName("");     //Will call init in the digitizing
126   SetTitle("aliroot") ;  
127 }
128
129 //____________________________________________________________________________ 
130   AliPHOSDigitizer::~AliPHOSDigitizer()
131 {
132   // dtor
133
134
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
155   //Making digits with noise, first EMC
156   Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
157   
158   Int_t nCPV ;
159   Int_t absID ;
160   TString name      =  geom->GetName() ;
161   
162   nCPV = nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
163     geom->GetNModules() ;
164
165   digits->Expand(nCPV) ;
166
167   // get first the sdigitizer from the tasks list (must have same name as the digitizer)
168   const AliPHOSSDigitizer * sDigitizer = gime->SDigitizer(GetName()); 
169   if ( !sDigitizer) {
170     cerr << "ERROR: AliPHOSDigitizer::Digitize -> SDigitizer with name " << GetName() << " not found " << endl ; 
171     abort() ; 
172   }
173     
174   // loop through the sdigits posted to the White Board and add them to the noise
175   TCollection * folderslist = gime->SDigitsFolder()->GetListOfFolders() ; 
176   TIter next(folderslist) ; 
177   TFolder * folder = 0 ; 
178   TClonesArray * sdigits = 0 ;
179   Int_t input = 0 ;
180   TObjArray * sdigArray = new TObjArray(2) ;
181   while ( (folder = (TFolder*)next()) ) 
182     if ( (sdigits = (TClonesArray*)folder->FindObject(GetName()) ) ) {
183       TString fileName(folder->GetName()) ;
184       fileName.ReplaceAll("_","/") ;
185       cout << "INFO: AliPHOSDigitizer::Digitize -> Adding SDigits " 
186            << GetName() << " from " << fileName << endl ; 
187       sdigArray->AddAt(sdigits, input) ;
188       input++ ;
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 = (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     AliPHOSDigit * digit = (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 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
362 {
363   Int_t chanel ;
364   if(absId <= fEmcCrystals){ //digitize as EMC 
365     chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;       
366     if(chanel > fNADCemc ) chanel =  fNADCemc ;
367   }
368   else{ //Digitize as CPV
369     chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;       
370     if(chanel > fNADCcpv ) chanel =  fNADCcpv ;
371   }
372   return chanel ;
373 }
374 //____________________________________________________________________________
375 void AliPHOSDigitizer::Exec(Option_t *option) 
376
377   // Managing method
378
379   if(strcmp(GetName(), "") == 0 )   
380     Init() ;
381   
382   if (strstr(option,"print")) {
383     Print("");
384     return ; 
385   }
386   
387   if(strstr(option,"tim"))
388     gBenchmark->Start("PHOSDigitizer");
389
390   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
391   
392   Int_t nevents ;
393   
394   TTree * treeD ;
395   
396   if(fManager){
397     treeD = fManager->GetTreeD() ;
398     nevents = 1 ;    // Will process only one event
399   }
400   else {
401     gAlice->GetEvent(0) ;
402     nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
403     treeD=gAlice->TreeD() ;
404   }
405
406
407   //Check, if this branch already exits
408   if (treeD) { 
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
436   Int_t ievent ;
437   
438   for(ievent = 0; ievent < nevents; ievent++){
439   
440     if(fManager){
441       Int_t input ;
442       for(input = 0 ; input < fManager->GetNinputs(); input ++){
443         TTree * treeS = fManager->GetInputTreeS(input) ;
444         if(!treeS){
445           cerr << "AliPHOSDigitizer -> No Input " << endl ;
446           return ;
447         }
448         gime->ReadTreeS(treeS,input) ;
449       }
450     }
451     else
452       gime->Event(ievent,"S") ;
453     
454     Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
455     
456     WriteDigits(ievent) ;
457     
458     if(strstr(option,"deb"))
459       PrintDigits(option);
460     
461     //increment the total number of Digits per run 
462     fDigitsInRun += gime->Digits()->GetEntriesFast() ;  
463   }
464   
465   if(strstr(option,"tim")){
466     gBenchmark->Stop("PHOSDigitizer");
467     cout << "AliPHOSDigitizer:" << endl ;
468     cout << "  took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for Digitizing " 
469          <<  gBenchmark->GetCpuTime("PHOSDigitizer")/nevents << " seconds per event " << endl ;
470     cout << endl ;
471   }
472   
473 }
474
475 //____________________________________________________________________________ 
476 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks) 
477 { // 
478   ticks->Sort() ; //Sort in accordance with times of ticks
479   TIter it(ticks) ;
480   AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
481   Float_t time = ctick->CrossingTime(fTimeThreshold) ;    
482
483   AliPHOSTick * t ;  
484   while((t=(AliPHOSTick*) it.Next())){
485     if(t->GetTime() < time)  //This tick starts before crossing
486       *ctick+=*t ;
487     else
488       return time ;
489
490     time = ctick->CrossingTime(fTimeThreshold) ;    
491   }
492   return time ;
493 }
494 //____________________________________________________________________________ 
495 Bool_t AliPHOSDigitizer::Init()
496 {
497   fPinNoise           = 0.004 ;
498   fEMCDigitThreshold  = 0.012 ;
499   fCPVNoise           = 0.01;
500   fCPVDigitThreshold  = 0.09 ;
501   fTimeResolution     = 0.5e-9 ;
502   fTimeSignalLength   = 1.0e-9 ;
503   fDigitsInRun  = 0 ; 
504   fADCchanelEmc = 0.0015;        // width of one ADC channel in GeV
505   fADCpedestalEmc = 0.005 ;      //
506   fNADCemc = (Int_t) TMath::Power(2,16) ;  // number of channels in EMC ADC
507
508   fADCchanelCpv = 0.0012 ;          // width of one ADC channel in CPV 'popugais'
509   fADCpedestalCpv = 0.012 ;         // 
510   fNADCcpv = (Int_t) TMath::Power(2,12);      // number of channels in CPV ADC
511
512   fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
513     
514   if(fManager)
515     SetTitle("aliroot") ;
516   else if (strcmp(GetTitle(),"")==0) 
517    SetTitle("galice.root") ;
518
519   if( strcmp(GetName(), "") == 0 )
520     SetName("Default") ;
521   
522   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ; 
523   if ( gime == 0 ) {
524     cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ; 
525     return kFALSE;
526   } 
527   
528   const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
529   fEmcCrystals = geom->GetNModules() *  geom->GetNCristalsInModule() ;
530   
531   // Post Digits to the white board
532   gime->PostDigits(GetName() ) ;   
533   
534   // Post Digitizer to the white board
535   gime->PostDigitizer(this) ;
536   
537   //Mark that we will use current header file
538   if(!fManager){
539     gime->PostSDigits(GetName(),GetTitle()) ;
540     gime->PostSDigitizer(GetName(),GetTitle()) ;
541   }
542   return kTRUE ;
543 }
544
545 //__________________________________________________________________
546 void AliPHOSDigitizer::MixWith(const char* headerFile)
547 {
548   // Allows to produce digits by superimposing background and signal event.
549   // It is assumed, that headers file with SIGNAL events is opened in 
550   // the constructor. 
551   // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed 
552   // Thus we avoid writing (changing) huge and expensive 
553   // backgound files: all output will be writen into SIGNAL, i.e. 
554   // opened in constructor file. 
555   //
556   // One can open as many files to mix with as one needs.
557   // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
558   // can be mixed.
559
560   if( strcmp(GetName(), "") == 0 )
561     Init() ;
562
563   if(fManager){
564     cout << "Can not use this method under AliRunDigitizer " << endl ;
565     return ;
566   }
567   
568   // check if the specified SDigits do not already exist on the White Board:
569   // //Folders/RunMC/Event/Data/PHOS/SDigits/headerFile/sdigitsname
570
571   TString path = "Folders/RunMC/Event/Data/PHOS/SDigits" ; 
572   path += headerFile ; 
573   path += "/" ; 
574   path += GetName() ;
575   if ( gROOT->FindObjectAny(path.Data()) ) {
576     cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
577     return;
578   }
579
580   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
581   gime->PostSDigits(GetName(),headerFile) ;
582   
583   // check if the requested file is already open or exist and if SDigits Branch exist
584   TFile * file = (TFile*)gROOT->FindObject(headerFile); 
585   if ( !file ) { 
586     file = new TFile(headerFile, "READ") ; 
587     if (!file) { 
588       cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ; 
589       return ; 
590     }
591   }
592   
593 }
594
595 //__________________________________________________________________
596 void AliPHOSDigitizer::SetSplitFile(const TString splitFileName) const
597 {
598   // Diverts the Digits in a file separate from the hits file
599   
600   // I guess it is not going to work if we do merging
601   if (fManager) {
602     cerr << "ERROR: AliPHOSDigitizer::SetSplitFile -> Not yet available in case of merging activated " << endl ;  
603     return ; 
604   }
605
606   TDirectory * cwd = gDirectory ;
607   TFile * splitFile = gAlice->InitTreeFile("D",splitFileName.Data());
608   splitFile->cd() ; 
609   if ( !splitFile->Get("gAlice") ) 
610     gAlice->Write();
611   
612   TTree *treeE  = gAlice->TreeE();
613   if (!treeE) {
614     cerr << "ERROR: AliPHOSDigitizer::SetSplitFile -> No TreeE found "<<endl;
615     abort() ;
616   }      
617   
618   // copy TreeE
619   if ( !splitFile->Get("TreeE") ) {
620     AliHeader *header = new AliHeader();
621     treeE->SetBranchAddress("Header", &header);
622     treeE->SetBranchStatus("*",1);
623     TTree *treeENew =  treeE->CloneTree();
624     treeENew->Write();
625   }
626
627   // copy AliceGeom
628   if ( !splitFile->Get("AliceGeom") ) {
629     TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
630     if (!AliceGeom) {
631       cerr << "ERROR: AliPHOSDigitizer::SetSplitFile -> AliceGeom was not found in the input file "<<endl;
632       abort() ;
633     }
634     AliceGeom->Write();
635   }
636   
637   cwd->cd() ; 
638   gAlice->MakeTree("D",splitFile);
639   cout << "INFO: AliPHOSDigitizer::SetSPlitMode -> Digits will be stored in " << splitFileName.Data() << endl ; 
640 }
641
642 //__________________________________________________________________
643 void AliPHOSDigitizer::Print(Option_t* option)const {
644   // Print Digitizer's parameters
645   if( strcmp(GetName(), "") != 0 ){
646     
647     cout << "------------------- "<< GetName() << " -------------" << endl ;
648     cout << "Digitizing sDigits from file(s): " <<endl ;
649     
650      TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("Folders/RunMC/Event/Data/PHOS/SDigits"))->GetListOfFolders() ; 
651     TIter next(folderslist) ; 
652     TFolder * folder = 0 ; 
653     
654     while ( (folder = (TFolder*)next()) ) {
655       if ( folder->FindObject(GetName())  ) 
656         cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ; 
657     }
658     cout << endl ;
659     cout << "Writing digits to " << GetTitle() << endl ;
660     
661     cout << endl ;
662     cout << "With following parameters: " << endl ;
663     cout << "     Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
664     cout << "  Threshold  in EMC  (fEMCDigitThreshold) = " << fEMCDigitThreshold  << endl ; ;
665     cout << "                 Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ; 
666     cout << "    Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ; 
667     cout << "---------------------------------------------------" << endl ;
668   }
669   else
670     cout << "AliPHOSDigitizer not initialized " << endl ;
671   
672 }
673
674 //__________________________________________________________________
675 void AliPHOSDigitizer::PrintDigits(Option_t * option){
676   // Print a table of digits
677
678   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
679   TClonesArray * digits = gime->Digits() ; 
680
681   cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
682   cout << "       Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
683   cout << endl ;
684   if(strstr(option,"all")||strstr(option,"EMC")){
685     
686     //loop over digits
687     AliPHOSDigit * digit;
688     cout << "EMC digits (with primaries): " << endl ;
689     cout << "Digit Id     Amplitude   Index       Nprim  Primaries list " <<  endl;      
690     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
691     Int_t index ;
692     for (index = 0 ; (index < digits->GetEntriesFast()) && 
693          (((AliPHOSDigit * )  digits->At(index))->GetId() <= maxEmc) ; index++) {
694       digit = (AliPHOSDigit * )  digits->At(index) ;
695       if(digit->GetNprimary() == 0) continue;
696       cout << setw(6)  <<  digit->GetId() << "   "  <<  setw(10)  <<  digit->GetAmp() <<   "    "  
697            << setw(6)  <<  digit->GetIndexInList() << "    "   
698            << setw(5)  <<  digit->GetNprimary() <<"    ";
699       
700       Int_t iprimary;
701       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
702         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
703       cout << endl;      
704     }    
705     cout << endl;
706   }
707
708   if(strstr(option,"all")||strstr(option,"CPV")){
709     
710     //loop over CPV digits
711     AliPHOSDigit * digit;
712     cout << "CPV digits: " << 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(); index++) {
717       digit = (AliPHOSDigit * )  digits->At(index) ;
718       if(digit->GetId() > maxEmc){
719         cout << setw(6)  <<  digit->GetId() << "   "  <<        setw(10)  <<  digit->GetAmp() <<   "    "  
720              << setw(6)  <<  digit->GetIndexInList() << "    "   
721              << setw(5)  <<  digit->GetNprimary() <<"    ";
722         
723         Int_t iprimary;
724         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
725           cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
726         cout << endl;    
727       }    
728     }
729   }
730
731 }
732
733 //__________________________________________________________________
734 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
735 {
736   // we set title (comment) of the SDigits branch in the first! header file
737   if( strcmp(GetName(), "") == 0 )
738     Init() ;
739
740   AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ; 
741  
742 }
743 //__________________________________________________________________
744 Float_t AliPHOSDigitizer::TimeOfNoise(void)
745 {  // Calculates the time signal generated by noise
746   //to be rewritten, now returns just big number
747   return 1. ;
748
749 }
750 //____________________________________________________________________________
751 void AliPHOSDigitizer::Reset() 
752
753   // sets current event number to the first simulated event
754
755   if( strcmp(GetName(), "") == 0 )
756     Init() ;
757
758  //  Int_t inputs ;
759 //   for(inputs = 0; inputs < fNinputs ;inputs++)
760 //       fIevent->AddAt(-1, inputs ) ;
761   
762 }
763
764 //____________________________________________________________________________
765 void AliPHOSDigitizer::WriteDigits(Int_t event)
766 {
767
768   // Makes TreeD in the output file. 
769   // Check if branch already exists: 
770   //   if yes, exit without writing: ROOT TTree does not support overwriting/updating of 
771   //      already existing branches. 
772   //   else creates branch with Digits, named "PHOS", title "...",
773   //      and branch "AliPHOSDigitizer", with the same title to keep all the parameters
774   //      and names of files, from which digits are made.
775
776   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
777   const TClonesArray * digits = gime->Digits(GetName()) ; 
778  TTree * treeD ;
779
780  
781  if(fManager) 
782    treeD = fManager->GetTreeD() ;
783  else {
784    if (!gAlice->TreeD() ) 
785      gAlice->MakeTree("D");
786    treeD = gAlice->TreeD();
787  }
788
789
790   // -- create Digits branch
791   Int_t bufferSize = 32000 ;    
792   TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
793   digitsBranch->SetTitle(GetName());
794     
795   // -- Create Digitizer branch
796   Int_t splitlevel = 0 ;
797   AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
798   TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel); 
799   digitizerBranch->SetTitle(GetName());
800
801   digitsBranch->Fill() ;
802   treeD->AutoSave() ; 
803  
804 }