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