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