]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSDigitizer.cxx
fDigitizer set to 0 in ctr
[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 // Class performs digitization of Summable digits (in the PHOS case this is just
22 // sum of contributions of all primary particles into given cell). 
23 // In addition it performs mixing of summable digits from different events.
24 //
25 // For each event two branches are created in TreeD:
26 //   "PHOS" - list of digits
27 //   "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
28 //
29 // Note, that one cset title for new digits branch, and repeat digitization with
30 // another set of parameters.
31 //
32 // Examples of use:
33 // root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
34 // root[1] d->ExecuteTask()             
35 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
36 //                       //Digitizes SDigitis in all events found in file galice.root 
37 //
38 // root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;  
39 //                       // Will read sdigits from galice1.root
40 // root[3] d1->MixWith("galice2.root")       
41 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 //                       // Reads another portion of sdigits from galice2.root
43 // root[3] d1->MixWith("galice3.root")       
44 //                       // Reads another portion of sdigits from galice3.root
45 // root[4] d->ExecuteTask("deb timing")    
46 //                       // Reads SDigits from files galice1.root, galice2.root ....
47 //                       // mixes them and stores produced Digits in file galice1.root          
48 //                       // deb - prints number of produced digits
49 //                       // deb all - prints list of produced digits
50 //                       // timing  - prints time used for digitization
51 //
52
53 // --- ROOT system ---
54 #include "TFile.h"
55 #include "TTree.h"
56 #include "TSystem.h"
57 #include "TROOT.h"
58 #include "TFolder.h"
59 #include "TObjString.h"
60 #include "TBenchmark.h"
61 // --- Standard library ---
62 #include <iomanip.h>
63
64 // --- AliRoot header files ---
65
66 #include "AliRun.h"
67 #include "AliPHOSDigit.h"
68 #include "AliPHOSHit.h"
69 #include "AliPHOSv1.h"
70 #include "AliPHOSDigitizer.h"
71 #include "AliPHOSSDigitizer.h"
72 #include "AliPHOSGeometry.h"
73
74 ClassImp(AliPHOSDigitizer)
75
76
77 //____________________________________________________________________________ 
78   AliPHOSDigitizer::AliPHOSDigitizer():TTask("AliPHOSDigitizer","") 
79 {
80   // ctor
81
82   fSDigitizer = 0 ;
83   fNinputs = 1 ;
84   fPinNoise = 0.01 ;
85   fEMCDigitThreshold = 0.01 ;
86   fCPVNoise = 0.01;
87   fCPVDigitThreshold = 0.09 ;
88   fPPSDNoise = 0.0000001;
89   fPPSDDigitThreshold = 0.0000002 ;
90   fInitialized = kFALSE ;
91
92   fHeaderFiles = 0;
93   fSDigitsTitles = 0;
94   fSDigits  = 0 ;
95   fDigits = 0;
96
97 }
98 //____________________________________________________________________________ 
99 void AliPHOSDigitizer::Init(){
100 // Makes all memory allocations
101
102   if(!fInitialized){
103     
104     cout << "In Init" << endl ;
105
106     fHeaderFiles  = new TClonesArray("TObjString",1) ;
107     new((*fHeaderFiles)[0]) TObjString("galice.root") ;
108     
109     //Test, if this file already open
110     
111     TFile *file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(0))->GetString() ) ;
112     
113     if(file == 0){
114       file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString(),"update") ;
115       gAlice = (AliRun *) file->Get("gAlice") ;
116     }
117     else
118       file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString()) ;
119     
120     file->cd() ;
121     
122     fSDigitsTitles = new TClonesArray("TObjString",1);
123     new((*fSDigitsTitles)[0]) TObjString("") ;   
124     
125     fSDigits      = new TClonesArray("TClonesArray",1) ;
126     new((*fSDigits)[0]) TClonesArray("AliPHOSDigit",1000) ;
127
128     fSDigitizer = 0 ;
129     
130     fDigitsTitle = "" ;
131     
132     fDigits = new TClonesArray("AliPHOSDigit",200000) ;
133     
134     fIevent    = new TArrayI(1) ;
135     fIevent->AddAt(-1,0 ) ; 
136     fIeventMax = new TArrayI(1) ;
137     
138     fIeventMax->AddAt((Int_t) gAlice->TreeE()->GetEntries(), 0 );
139     
140     // add Task to //root/Tasks folder
141     TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
142     roottasks->Add(this) ; 
143     
144     fInitialized = kTRUE ;
145   }
146   
147 }
148
149 //____________________________________________________________________________ 
150 AliPHOSDigitizer::AliPHOSDigitizer(const char *HeaderFile,const char *sDigitsTitle):
151   TTask("AliPHOSDigitizer","")
152 {
153   // ctor
154   fHeaderFiles  = new TClonesArray("TObjString",1) ;          
155   new((*fHeaderFiles)[0]) TObjString(HeaderFile) ;
156   
157   // Header file, where result will be stored
158   TFile * file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(0))->GetString() ) ;
159   if(file==0){
160     file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString(),"update") ;      
161     gAlice = (AliRun *) file->Get("gAlice") ;  //If not read yet
162   }
163   
164   file->cd() ;
165   
166   fSDigitsTitles = new TClonesArray("TObjString",1);         // Title name of the SDigits branch
167   new((*fSDigitsTitles)[0]) TObjString(sDigitsTitle) ;  
168     
169   fSDigits      = new TClonesArray("TClonesArray",1) ;      // here list of SDigits wil be stored
170   new((*fSDigits)[0]) TClonesArray("AliPHOSDigit",1000) ;
171     
172   fDigits = new TClonesArray("AliPHOSDigit",200000) ;
173   
174   fDigitsTitle = "" ; 
175
176   
177   fSDigitizer = 0 ;
178   
179   fIevent    = new TArrayI(1) ;
180   fIevent->AddAt(-1,0 ) ; 
181   fIeventMax = new TArrayI(1) ;
182   
183   // Get number of events to process
184   fIeventMax->AddAt((Int_t) gAlice->TreeE()->GetEntries(), 0 );
185   
186   fNinputs = 1 ;
187   
188   fPinNoise = 0.01 ;
189   fEMCDigitThreshold = 0.01 ;
190   fCPVNoise = 0.01;
191   fCPVDigitThreshold = 0.09 ;
192   fPPSDNoise = 0.0000001;
193   fPPSDDigitThreshold = 0.0000002 ;  
194   
195   // add Task to //root/Tasks folder
196   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
197   roottasks->Add(this) ; 
198   fInitialized = kTRUE ;
199   
200 }
201
202 //____________________________________________________________________________ 
203   AliPHOSDigitizer::~AliPHOSDigitizer()
204 {
205   // dtor
206
207   if(fHeaderFiles)  delete fHeaderFiles ;
208   if(fSDigitsTitles) delete fSDigitsTitles ;
209   if(fSDigits)      delete fSDigits ;
210   if(fDigits)       delete fDigits ;
211 }
212 //____________________________________________________________________________
213 void AliPHOSDigitizer::Reset() { 
214   //sets current event number to the first simulated event
215
216   if(!fInitialized)
217     Init() ;
218
219   Int_t inputs ;
220   for(inputs = 0; inputs < fNinputs ;inputs++)
221       fIevent->AddAt(-1, inputs ) ;
222   
223 }
224 //____________________________________________________________________________
225 Bool_t AliPHOSDigitizer::Combinator() { 
226
227   //Makes all desirable combinations Signal+Background,
228   // returns kFALSE when all combinations are made
229   // May be useful to introduce some options like "One-to-One", "All-to-One" and "All-to-All" ?
230
231   //realizing "One-to-One" option...
232
233   if(!fInitialized)
234     Init() ;
235
236   Int_t inputs ;
237   Bool_t endNotReached = kTRUE ;
238
239   for(inputs = 0; (inputs < fNinputs) && endNotReached ;inputs++){
240     if(fIevent->At(inputs)+1 < fIeventMax->At(inputs))
241       fIevent->AddAt(fIevent->At(inputs)+1, inputs ) ;
242     else
243       if(inputs == 0)
244         endNotReached = kFALSE ;
245       else //for inputs other than base one start from the beginning
246         fIevent->AddAt(0, inputs ) ;
247     
248   }
249   return endNotReached ;
250   
251 }
252
253 //____________________________________________________________________________
254 void AliPHOSDigitizer::Digitize(Option_t *option) { 
255
256   // Makes the digitization of the collected summable digits
257   // for this it first creates the array of all PHOS modules
258   // filled with noise (different for EMC, CPV and PPSD) and
259   // after that adds contributions from SDigits. This design 
260   // helps to avoid scanning over the list of digits to add 
261   // contribution of any new SDigit.
262
263   if(!fInitialized)
264     Init() ;
265
266   fDigits->Clear() ;
267
268   AliPHOS * PHOS = (AliPHOS *) gAlice->GetDetector("PHOS") ;   
269   AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance( PHOS->GetGeometry()->GetName(), PHOS->GetGeometry()->GetTitle() );
270
271   //Making digits with noise, first EMC
272   Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
273   
274   Int_t nCPV ;
275   Int_t nPPSD ;
276   Int_t absID ;
277   TString name      =  geom->GetName() ;
278   
279   if ( name == "IHEP" || name == "MIXT" )    
280     nCPV =nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
281       geom->GetNCPVModules()*geom->GetNumberOfCPVLayers() ;
282   else
283     nCPV = nEMC; 
284   
285   if ( name == "GPS2" || name == "MIXT" )    
286     nPPSD =nCPV+2*geom->GetNPPSDModules()*geom->GetNumberOfModulesPhi()*geom->GetNumberOfModulesZ()*
287       geom->GetNumberOfPadsPhi()*geom->GetNumberOfPadsZ() ;
288   else
289     nPPSD = nCPV; 
290
291
292   fDigits->Expand(nPPSD) ;
293
294   
295   for(absID = 1; absID <= nEMC; absID++){
296     Float_t noise = gRandom->Gaus(0., fPinNoise) ; 
297     new((*fDigits)[absID-1]) AliPHOSDigit( -1,absID,fSDigitizer->Digitize(noise) ) ;
298   }
299   
300   for(absID = nEMC+1; absID <= nCPV; absID++){
301     Float_t noise = gRandom->Gaus(0., fCPVNoise) ; 
302     new((*fDigits)[absID-1]) AliPHOSDigit( -1,absID,fSDigitizer->Digitize(noise) ) ;
303   }
304   
305   for(absID = nCPV+1; absID <= nPPSD; absID++){
306     Float_t noise = gRandom->Gaus(0., fPPSDNoise) ; 
307     new((*fDigits)[absID-1]) AliPHOSDigit( -1,absID,fSDigitizer->Digitize(noise) ) ;
308   }
309   
310
311   // Now look throught (unsorted) list of SDigits and add corresponding digits  
312   AliPHOSDigit *curSDigit ;
313   AliPHOSDigit *digit ;
314     
315   Int_t inputs;
316   for(inputs = 0; inputs< fNinputs ; inputs++){  //loop over (possible) merge sources
317     
318     TClonesArray * sdigits= (TClonesArray *)fSDigits->At(inputs) ;
319     Int_t isdigit ;
320
321     Int_t nSDigits = sdigits->GetEntries() ;     
322     for(isdigit=0;isdigit< nSDigits; isdigit++){
323       curSDigit = (AliPHOSDigit *)sdigits->At(isdigit) ;
324       if(inputs)                                       //Shift primaries for non-background sdigits
325         curSDigit->ShiftPrimary(inputs) ;
326       digit = (AliPHOSDigit *)fDigits->At(curSDigit->GetId() - 1);
327       *digit = *digit + *curSDigit ;
328     }  
329   }
330
331
332   //remove digits below thresholds
333   for(absID = 0; absID < nEMC ; absID++)
334     if(fSDigitizer->Calibrate(((AliPHOSDigit*)fDigits->At(absID))->GetAmp()) < fEMCDigitThreshold)
335       fDigits->RemoveAt(absID) ;
336   for(absID = nEMC; absID < nCPV ; absID++)
337     if(fSDigitizer->Calibrate(((AliPHOSDigit*)fDigits->At(absID))->GetAmp()) < fCPVDigitThreshold)
338       fDigits->RemoveAt(absID) ;
339   for(absID = nCPV; absID < nPPSD ; absID++)
340     if(fSDigitizer->Calibrate(((AliPHOSDigit *)fDigits->At(absID))->GetAmp()) < fPPSDDigitThreshold)
341       fDigits->RemoveAt(absID) ;
342   
343   fDigits->Compress() ;  
344   
345   Int_t ndigits = fDigits->GetEntriesFast() ;
346
347   fDigits->Expand(ndigits) ;
348
349
350   //Set indexes in list of digits
351   Int_t i ;
352   for (i = 0 ; i < ndigits ; i++) { 
353     AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ; 
354     digit->SetIndexInList(i) ;     
355   }
356 }
357 //____________________________________________________________________________
358 void AliPHOSDigitizer::WriteDigits(){
359
360   // Made TreeD in the output file. Check if branch already exists: if yes, exits 
361   // without writing: ROOT TTree does not suppert overwriting/updating of the 
362   // already existing branches. Creates branch with Digits, named "PHOS", title "...",
363   // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
364   // and names of files, from which digits are made.
365
366   gAlice->GetEvent(fIevent->At(0)) ;  // Suitable only for One-To-One mixing
367   gAlice->SetEvent(fIevent->At(0)) ;  // for all-to-all will produce a lot of branches in TreeD
368
369   if(gAlice->TreeD()==0)
370     gAlice->MakeTree("D") ;
371
372   
373   //Check, if this branch already exits?
374   TBranch * digitsBranch = 0;
375   TBranch * digitizerBranch = 0;
376   
377   TObjArray * branches = gAlice->TreeD()->GetListOfBranches() ;
378   Int_t ibranch;
379   Bool_t phosNotFound = kTRUE ;
380   Bool_t digitizerNotFound = kTRUE ;
381   
382   for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
383     
384     if(phosNotFound){
385       digitsBranch=(TBranch *) branches->At(ibranch) ;
386       if( (strcmp("PHOS",digitsBranch->GetName())==0 ) &&
387           (fDigitsTitle.CompareTo(digitsBranch->GetTitle()) == 0) )
388         phosNotFound = kFALSE ;
389     }
390     if(digitizerNotFound){
391       digitizerBranch = (TBranch *) branches->At(ibranch) ;
392       if( (strcmp(digitizerBranch->GetName(),"AliPHOSDigitizer") == 0) &&
393           (fDigitsTitle.CompareTo(digitizerBranch->GetTitle()) == 0))
394         digitizerNotFound = kFALSE ;
395     }
396   }
397   
398   
399   if(!(digitizerNotFound && phosNotFound)){ 
400     cout << "AliPHOSDigitizer error: " << endl ;
401     cout << "       can not update/overwrite existing branches "<< endl ;
402     cout << "       do not write " << endl ;
403     return ;
404   }
405
406   // create new branches
407
408   //First generate file name
409   char * file =0;
410   if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
411     file = new char[strlen(gAlice->GetBaseFile())+20] ;
412     sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
413   }
414   
415   TDirectory *cwd = gDirectory;
416   
417   //First create list of sdigits
418   Int_t bufferSize = 32000 ;    
419   digitsBranch = gAlice->TreeD()->Branch("PHOS",&fDigits,bufferSize);
420   digitsBranch->SetTitle(fDigitsTitle.Data());
421   if (file) {
422     digitsBranch->SetFile(file);
423     TIter next( digitsBranch->GetListOfBranches());
424     while ((digitsBranch=(TBranch*)next())) {
425       digitsBranch->SetFile(file);
426     }   
427     cwd->cd();
428   } 
429     
430   //second - create Digitizer
431   Int_t splitlevel = 0 ;
432   AliPHOSDigitizer * d = this ;
433   digitizerBranch = gAlice->TreeD()->Branch("AliPHOSDigitizer","AliPHOSDigitizer",
434                                             &d,bufferSize,splitlevel); 
435   digitizerBranch->SetTitle(fDigitsTitle.Data());
436   if (file) {
437     digitizerBranch->SetFile(file);
438     TIter next( digitizerBranch->GetListOfBranches());
439     while ((digitizerBranch=(TBranch*)next())) {
440       digitizerBranch->SetFile(file);
441     }   
442     cwd->cd();
443   }
444
445   gAlice->TreeD()->Fill() ;
446   
447   gAlice->TreeD()->Write(0,kOverwrite) ;  
448
449   //remove fSDigitizer before new event.  
450   if(fSDigitizer){
451     delete fSDigitizer ;
452     fSDigitizer = 0 ;
453   }
454
455
456 }
457
458 //____________________________________________________________________________
459 void AliPHOSDigitizer::Exec(Option_t *option) { 
460   // Managing method
461
462   if(!fInitialized)    Init() ;
463
464   if(strstr(option,"tim"))
465     gBenchmark->Start("PHOSDigitizer");
466
467   //reset events numbers to start from the beginnig
468   Reset() ;
469   
470   while(Combinator()){  
471     
472     if(!ReadSDigits()) //read sdigits event(s) evaluated by Combinator() from file(s)
473       return ;    
474     
475     Digitize(option) ; //Add prepared SDigits to digits and add the noise
476     WriteDigits() ;
477     
478     if(strstr(option,"deb"))
479       PrintDigits(option);
480
481   }
482
483   if(strstr(option,"tim")){
484     gBenchmark->Stop("PHOSDigitizer");
485     cout << "AliPHOSDigitizer:" << endl ;
486     cout << "  took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for SDigitizing " 
487          <<  gBenchmark->GetCpuTime("PHOSDigitizer")/(fIeventMax->At(0)) << " seconds per event " << endl ;
488     cout << endl ;
489   }
490   
491 }
492
493 //__________________________________________________________________
494 Bool_t AliPHOSDigitizer::ReadSDigits(){
495 // Reads summable digits from the opened files for the particular set of events given by fIevent
496
497   if(!fInitialized)    Init() ;
498
499
500   Int_t inputs ;
501   for(inputs = fNinputs-1; inputs >= 0; inputs --){
502
503     Int_t event = fIevent->At(inputs) ;
504
505     TFile * file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(inputs))->GetString() ) ;
506     file->cd() ;
507
508     // Get SDigits Tree header from file
509     char treeName[20]; 
510     sprintf(treeName,"TreeS%d",event);
511     TTree * treeS = (TTree*)file->Get(treeName);
512    
513     if(treeS==0){
514       cout << "Error at AliPHOSDigitizer: no "<<treeName << "   in file " << file->GetName() << endl ;
515       cout << "Do nothing " << endl ;
516       return kFALSE ;
517     }
518
519     TBranch * sdigitsBranch = 0;
520     TBranch * sdigitizerBranch = 0;
521
522     TObjArray * branches = treeS->GetListOfBranches() ;
523     Int_t ibranch;
524     Bool_t phosNotFound = kTRUE ;
525     Bool_t sdigitizerNotFound = kTRUE ;
526   
527     for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
528             
529       if(phosNotFound){
530         sdigitsBranch=(TBranch *) branches->At(ibranch) ;
531         if(( strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
532            ((TObjString*) fSDigitsTitles->At(inputs))->GetString().CompareTo(sdigitsBranch->GetTitle())== 0 )
533               phosNotFound = kFALSE ;
534         
535       }
536       
537       if(sdigitizerNotFound){
538         sdigitizerBranch = (TBranch *) branches->At(ibranch) ;
539         if(( strcmp(sdigitizerBranch->GetName(),"AliPHOSSDigitizer") == 0) &&
540            ((TObjString*) fSDigitsTitles->At(inputs))->GetString().CompareTo(sdigitizerBranch->GetTitle())== 0 )
541               sdigitizerNotFound = kFALSE ;
542         
543       }
544     }
545     
546     if(sdigitizerNotFound || phosNotFound){
547       cout << "Can't find Branch with sdigits or SDigitizer in the file " ;
548       if( ((TObjString*)fSDigitsTitles->At(inputs))->GetString().IsNull() )
549         cout << file->GetName() << endl ;       
550       else
551         cout << ((TObjString*)fSDigitsTitles->At(inputs))->GetString().Data() << endl ;
552       cout << "Do nothing" <<endl  ;
553       return kFALSE ;
554     }
555     
556     TClonesArray * sdigits = (TClonesArray*) fSDigits->At(inputs) ;  
557     sdigitsBranch->SetAddress(&sdigits) ;
558     
559     AliPHOSSDigitizer *sDigitizer = new AliPHOSSDigitizer();
560     sdigitizerBranch->SetAddress(&sDigitizer) ;
561     treeS->GetEvent(0) ;
562     
563     if(fSDigitizer == 0)
564       fSDigitizer = sDigitizer ;
565     else
566       if(!((*fSDigitizer)==(*sDigitizer)) ){
567         cout << "AliPHOSDigitizer ERROR:" << endl ;
568         cout << "       you are using sdigits made with different SDigitizers" << endl ;
569         cout << "fSD " << fSDigitizer << "  SD" << sDigitizer << endl ;
570         fSDigitizer->Print("") ;
571         sDigitizer->Print("") ;
572         cout << "Do Nothing " << endl ;
573         return kFALSE ;
574       }
575     
576   }
577   fPedestal = fSDigitizer->GetPedestalParameter() ;
578   fSlope    = fSDigitizer->GetCalibrationParameter() ;
579   
580   return kTRUE ;
581
582 }
583 //__________________________________________________________________
584 void AliPHOSDigitizer::MixWith(char* HeaderFile, char* sDigitsTitle){
585   // Alows produce digits by superimposing background and signal event.
586   // It is assumed, that headers file with SIGNAL events is opened in 
587   // constructor, and now we set the BACKGROUND event, with which we 
588   // will mix. Thus we avoid writing (changing) huge and expencive 
589   // backgound files: all output will be writen into SIGNAL, i.e. 
590   // opened in constructor file. 
591   //
592   // One can open as many files to mix with as one wants.
593
594
595   if(!fInitialized)
596     Init() ;
597
598
599   if(HeaderFile == 0){
600     cout << "Specify at least header file to merge"<< endl ;
601     return ;
602   }
603   
604   Int_t inputs ;
605   for(inputs = 0; inputs < fNinputs ; inputs++){
606     if(strcmp(((TObjString *)fHeaderFiles->At(inputs))->GetString(),HeaderFile) == 0 ){
607       if(sDigitsTitle == 0){ 
608         if(((TObjString*)fSDigitsTitles->At(inputs))->GetString().CompareTo("")  == 0){
609           cout << "Entry already exists, do not add" << endl ;
610           return ;
611         }
612       }
613       else
614         if(((TObjString*)fSDigitsTitles->At(inputs))->GetString().CompareTo(sDigitsTitle)){
615           cout << "Entry already exists, do not add" << endl ;
616           return;
617         }
618     }   
619   }  
620   
621   fHeaderFiles->Expand(fNinputs+1) ;
622   new((*fHeaderFiles)[fNinputs]) TObjString(HeaderFile) ;
623   
624   
625   TFile * file = new TFile(((TObjString *) fHeaderFiles->At(fNinputs))->GetString()) ;  
626   
627   file->cd() ;
628   
629   fSDigitsTitles->Expand(fNinputs+1) ;
630   new((*fSDigitsTitles)[fNinputs]) TObjString(sDigitsTitle) ;
631   
632   fSDigits->Expand(fNinputs+1) ;
633   new((*fSDigits)[fNinputs]) TClonesArray("AliPHOSDigit",1000) ;
634   
635   fIevent->Set(fNinputs+1) ;
636   fIevent->AddAt(-1, fNinputs) ;
637   
638   fIeventMax->Set(fNinputs+1) ;  
639   
640   TTree * te = (TTree *) file->Get("TE") ;
641   fIeventMax->AddAt((Int_t) te->GetEntries(), fNinputs );
642   
643   fNinputs++ ;
644   
645 }
646 //__________________________________________________________________
647 void AliPHOSDigitizer::Print(Option_t* option)const {
648   
649   if(fInitialized){
650     
651     cout << "------------------- "<< GetName() << " -------------" << endl ;
652     cout << "Digitizing sDigits from file(s): " <<endl ;
653     Int_t input ;
654     for(input = 0; input < fNinputs ; input++) {
655       cout << "          " << ((TObjString *) fHeaderFiles->At(input))->GetString() << 
656         "   Branch title:" << ((TObjString *) fSDigitsTitles->At(input))->GetString() << endl ;
657     }
658     cout << endl ;
659     cout << "Writing digits to " << ((TObjString *) fHeaderFiles->At(0))->GetString() << 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 << "               Noise in PPSD (fPPSDNoise) = " << fPPSDNoise << endl ;
668     cout << "  Threshold in PPSD (fPPSDDigitThreshold) = " << fPPSDDigitThreshold << endl ;
669     cout << "---------------------------------------------------" << endl ;
670   }
671   else
672     cout << "AliPHOSDigitizer not initialized " << endl ;
673   
674 }
675 //__________________________________________________________________
676 void AliPHOSDigitizer::PrintDigits(Option_t * option){
677     
678   cout << "AliPHOSDigitiser:"<< endl ;
679   cout << "       Number of entries in Digits list " << fDigits->GetEntriesFast() << endl ;
680   cout << endl ;
681   if(strstr(option,"all")){
682     
683     //loop over digits
684     AliPHOSDigit * digit;
685     cout << "Digit Id " << " Amplitude " <<  " Index "  <<  " Nprim " << " Primaries list " <<  endl;      
686     Int_t index ;
687     for (index = 0 ; index < fDigits->GetEntries() ; index++) {
688       digit = (AliPHOSDigit * )  fDigits->At(index) ;
689       cout << setw(8)  <<  digit->GetId() << " "  <<    setw(3)  <<  digit->GetAmp() <<   "  "  
690            << setw(6)  <<  digit->GetIndexInList() << "  "   
691            << setw(5)  <<  digit->GetNprimary() <<"  ";
692       
693       Int_t iprimary;
694       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
695         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << " ";
696       cout << endl;      
697     }
698     
699   }
700 }
701 //__________________________________________________________________
702 void AliPHOSDigitizer::SetSDigitsBranch(const char* title){
703   // we set title (comment) of the SDigits branch in the first! header file
704   if(!fInitialized)    Init() ;
705
706   ((TObjString*) fSDigitsTitles->At(0) )->SetString((char*)title) ;
707
708 }
709 //__________________________________________________________________
710 void AliPHOSDigitizer::SetDigitsBranch(const char* title){
711   //Sets the title (comment) of the branch to which Digits branch
712   if(!fInitialized)    Init() ;
713   
714   fDigitsTitle = title ;
715
716 }