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