1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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.
25 // For each event two branches are created in TreeD:
26 // "PHOS" - list of digits
27 // "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
29 // Note, that one can set a title for new digits branch, and repeat digitization with
30 // another set of parameters.
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
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
53 // --- ROOT system ---
59 #include "TObjString.h"
60 #include "TBenchmark.h"
61 // --- Standard library ---
64 // --- AliRoot header files ---
67 #include "AliPHOSDigit.h"
69 #include "AliPHOSDigitizer.h"
70 #include "AliPHOSSDigitizer.h"
71 #include "AliPHOSGeometry.h"
73 ClassImp(AliPHOSDigitizer)
76 //____________________________________________________________________________
77 AliPHOSDigitizer::AliPHOSDigitizer():TTask("AliPHOSDigitizer","")
84 fEMCDigitThreshold = 0.01 ;
86 fCPVDigitThreshold = 0.09 ;
87 fPPSDNoise = 0.0000001;
88 fPPSDDigitThreshold = 0.0000002 ;
89 fInitialized = kFALSE ;
97 //____________________________________________________________________________
98 void AliPHOSDigitizer::Init()
100 // Makes all memory allocations
104 fHeaderFiles = new TClonesArray("TObjString",1) ;
105 new((*fHeaderFiles)[0]) TObjString("galice.root") ;
107 //Test, if this file already open
109 TFile *file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(0))->GetString() ) ;
112 file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString(),"update") ;
113 gAlice = (AliRun *) file->Get("gAlice") ;
116 file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString()) ;
120 fSDigitsTitles = new TClonesArray("TObjString",1);
121 new((*fSDigitsTitles)[0]) TObjString("") ;
123 fSDigits = new TClonesArray("TClonesArray",1) ;
124 new((*fSDigits)[0]) TClonesArray("AliPHOSDigit",1000) ;
130 fDigits = new TClonesArray("AliPHOSDigit",200000) ;
132 fIevent = new TArrayI(1) ;
133 fIevent->AddAt(-1,0 ) ;
134 fIeventMax = new TArrayI(1) ;
136 fIeventMax->AddAt((Int_t) gAlice->TreeE()->GetEntries(), 0 );
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) ;
143 fInitialized = kTRUE ;
147 //____________________________________________________________________________
148 AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile,const char *sDigitsTitle):
149 TTask("AliPHOSDigitizer","")
152 fHeaderFiles = new TClonesArray("TObjString",1) ;
153 new((*fHeaderFiles)[0]) TObjString(headerFile) ;
155 // Header file, where result will be stored
156 TFile * file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(0))->GetString() ) ;
158 if(((TObjString *) fHeaderFiles->At(0))->GetString().Contains("rfio"))
159 file = TFile::Open(((TObjString *) fHeaderFiles->At(0))->GetString(),"update") ;
161 file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString(),"update") ;
162 gAlice = (AliRun *) file->Get("gAlice") ; //If not read yet
167 fSDigitsTitles = new TClonesArray("TObjString",1); // Title name of the SDigits branch
168 new((*fSDigitsTitles)[0]) TObjString(sDigitsTitle) ;
170 fSDigits = new TClonesArray("TClonesArray",1) ; // here list of SDigits wil be stored
171 new((*fSDigits)[0]) TClonesArray("AliPHOSDigit",1000) ;
173 fDigits = new TClonesArray("AliPHOSDigit",200000) ;
180 fIevent = new TArrayI(1) ;
181 fIevent->AddAt(-1,0 ) ;
182 fIeventMax = new TArrayI(1) ;
184 // Get number of events to process
185 fIeventMax->AddAt((Int_t) gAlice->TreeE()->GetEntries(), 0 );
190 fEMCDigitThreshold = 0.01 ;
192 fCPVDigitThreshold = 0.09 ;
193 fPPSDNoise = 0.0000001;
194 fPPSDDigitThreshold = 0.0000002 ;
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 ;
204 //____________________________________________________________________________
205 AliPHOSDigitizer::~AliPHOSDigitizer()
209 if(fHeaderFiles) delete fHeaderFiles ;
210 if(fSDigitsTitles) delete fSDigitsTitles ;
211 if(fSDigits) delete fSDigits ;
212 if(fDigits) delete fDigits ;
214 //____________________________________________________________________________
215 void AliPHOSDigitizer::Reset()
217 // sets current event number to the first simulated event
223 for(inputs = 0; inputs < fNinputs ;inputs++)
224 fIevent->AddAt(-1, inputs ) ;
227 //____________________________________________________________________________
228 Bool_t AliPHOSDigitizer::Combinator()
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...
239 Bool_t endNotReached = kTRUE ;
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 ) ;
246 endNotReached = kFALSE ;
247 else //for inputs other than base one start from the beginning
248 fIevent->AddAt(0, inputs ) ;
251 return endNotReached ;
255 //____________________________________________________________________________
256 void AliPHOSDigitizer::Digitize(Option_t *option)
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.
271 AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;
272 AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance( phos->GetGeometry()->GetName(), phos->GetGeometry()->GetTitle() );
274 //Making digits with noise, first EMC
275 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
280 TString name = geom->GetName() ;
282 if ( name == "IHEP" || name == "MIXT" )
283 nCPV =nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
284 geom->GetNCPVModules()*geom->GetNumberOfCPVLayers() ;
288 if ( name == "GPS2" || name == "MIXT" )
289 nPPSD =nCPV+2*geom->GetNPPSDModules()*geom->GetNumberOfModulesPhi()*geom->GetNumberOfModulesZ()*
290 geom->GetNumberOfPadsPhi()*geom->GetNumberOfPadsZ() ;
295 fDigits->Expand(nPPSD) ;
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) ) ;
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) ) ;
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) ) ;
314 // Now look throught (unsorted) list of SDigits and add corresponding digits
315 AliPHOSDigit *curSDigit ;
316 AliPHOSDigit *digit ;
319 for(inputs = 0; inputs< fNinputs ; inputs++){ //loop over (possible) merge sources
321 TClonesArray * sdigits= (TClonesArray *)fSDigits->At(inputs) ;
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 ;
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) ;
346 fDigits->Compress() ;
348 Int_t ndigits = fDigits->GetEntriesFast() ;
350 fDigits->Expand(ndigits) ;
353 //Set indexes in list of digits
355 for (i = 0 ; i < ndigits ; i++) {
356 AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ;
357 digit->SetIndexInList(i) ;
360 //____________________________________________________________________________
361 void AliPHOSDigitizer::WriteDigits()
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.
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
375 if(gAlice->TreeD()==0)
376 gAlice->MakeTree("D") ;
378 //Check, if this branch already exits?
379 TBranch * digitsBranch = 0;
380 TBranch * digitizerBranch = 0;
382 TObjArray * branches = gAlice->TreeD()->GetListOfBranches() ;
384 Bool_t phosNotFound = kTRUE ;
385 Bool_t digitizerNotFound = kTRUE ;
387 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
390 digitsBranch=(TBranch *) branches->At(ibranch) ;
391 if( (strcmp("PHOS",digitsBranch->GetName())==0 ) &&
392 (fDigitsTitle.CompareTo(digitsBranch->GetTitle()) == 0) )
393 phosNotFound = kFALSE ;
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 ;
404 if(!(digitizerNotFound && phosNotFound)){
405 cout << "AliPHOSDigitizer error: " << endl ;
406 cout << " can not update/overwrite existing branches "<< endl ;
407 cout << " do not write " << endl ;
411 // create new branches
413 //First generate file name
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()) ;
420 TDirectory *cwd = gDirectory;
422 //First create list of sdigits
423 Int_t bufferSize = 32000 ;
424 digitsBranch = gAlice->TreeD()->Branch("PHOS",&fDigits,bufferSize);
425 digitsBranch->SetTitle(fDigitsTitle.Data());
427 digitsBranch->SetFile(file);
428 TIter next( digitsBranch->GetListOfBranches());
430 while ((sbr=(TBranch*)next())) {
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());
443 digitizerBranch->SetFile(file);
444 TIter next( digitizerBranch->GetListOfBranches());
446 while ((sbr=(TBranch*)next())) {
452 digitsBranch->Fill() ;
453 digitizerBranch->Fill() ;
455 // gAlice->TreeD()->Fill() ; // YK 28.05.2001
456 gAlice->TreeD()->Write(0,kOverwrite) ;
458 //remove fSDigitizer before new event.
467 //____________________________________________________________________________
468 void AliPHOSDigitizer::Exec(Option_t *option)
472 if(!fInitialized) Init() ;
474 if(strstr(option,"tim"))
475 gBenchmark->Start("PHOSDigitizer");
477 //reset events numbers to start from the beginnig
482 if(!ReadSDigits()) //read sdigits event(s) evaluated by Combinator() from file(s)
485 Digitize(option) ; //Add prepared SDigits to digits and add the noise
488 if(strstr(option,"deb"))
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 ;
503 //__________________________________________________________________
504 Bool_t AliPHOSDigitizer::ReadSDigits()
506 // Reads summable digits from the opened files for the particular set of events given by fIevent
508 if(!fInitialized) Init() ;
512 for(inputs = fNinputs-1; inputs >= 0; inputs --){
514 Int_t event = fIevent->At(inputs) ;
516 TFile * file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(inputs))->GetString() ) ;
519 // Get SDigits Tree header from file
521 sprintf(treeName,"TreeS%d",event);
522 TTree * treeS = (TTree*)file->Get(treeName);
525 cout << "Error at AliPHOSDigitizer: no "<<treeName << " in file " << file->GetName() << endl ;
526 cout << "Do nothing " << endl ;
530 TBranch * sdigitsBranch = 0;
531 TBranch * sdigitizerBranch = 0;
533 TObjArray * branches = treeS->GetListOfBranches() ;
535 Bool_t phosNotFound = kTRUE ;
536 Bool_t sdigitizerNotFound = kTRUE ;
538 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
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 ;
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 ;
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 ;
562 cout << ((TObjString*)fSDigitsTitles->At(inputs))->GetString().Data() << endl ;
563 cout << "Do nothing" <<endl ;
567 TClonesArray * sdigits = (TClonesArray*) fSDigits->At(inputs) ;
568 sdigitsBranch->SetAddress(&sdigits) ;
570 AliPHOSSDigitizer *sDigitizer = new AliPHOSSDigitizer();
571 sdigitizerBranch->SetAddress(&sDigitizer) ;
573 sdigitsBranch->GetEntry(0) ;
574 sdigitizerBranch->GetEntry(0) ;
577 fSDigitizer = sDigitizer ;
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 ;
590 fPedestal = fSDigitizer->GetPedestalParameter() ;
591 fSlope = fSDigitizer->GetCalibrationParameter() ;
596 //__________________________________________________________________
597 void AliPHOSDigitizer::MixWith(char* HeaderFile, char* sDigitsTitle)
599 // Alows to produce digits by superimposing background and signal event.
600 // It is assumed, that headers file with SIGNAL events is opened in
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.
607 // One can open as many files to mix with as one needs.
615 cout << "Specify at least header file to merge"<< endl ;
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 ;
629 if(((TObjString*)fSDigitsTitles->At(inputs))->GetString().CompareTo(sDigitsTitle)){
630 cout << "Entry already exists, do not add" << endl ;
636 fHeaderFiles->Expand(fNinputs+1) ;
637 new((*fHeaderFiles)[fNinputs]) TObjString(HeaderFile) ;
640 TFile * file = new TFile(((TObjString *) fHeaderFiles->At(fNinputs))->GetString()) ;
644 fSDigitsTitles->Expand(fNinputs+1) ;
645 new((*fSDigitsTitles)[fNinputs]) TObjString(sDigitsTitle) ;
647 fSDigits->Expand(fNinputs+1) ;
648 new((*fSDigits)[fNinputs]) TClonesArray("AliPHOSDigit",1000) ;
650 fIevent->Set(fNinputs+1) ;
651 fIevent->AddAt(-1, fNinputs) ;
653 fIeventMax->Set(fNinputs+1) ;
655 TTree * te = (TTree *) file->Get("TE") ;
656 fIeventMax->AddAt((Int_t) te->GetEntries(), fNinputs );
661 //__________________________________________________________________
662 void AliPHOSDigitizer::Print(Option_t* option)const {
663 // Print Digitizer's parameters
666 cout << "------------------- "<< GetName() << " -------------" << endl ;
667 cout << "Digitizing sDigits from file(s): " <<endl ;
669 for(input = 0; input < fNinputs ; input++) {
670 cout << " " << ((TObjString *) fHeaderFiles->At(input))->GetString() <<
671 " Branch title:" << ((TObjString *) fSDigitsTitles->At(input))->GetString() << endl ;
674 cout << "Writing digits to " << ((TObjString *) fHeaderFiles->At(0))->GetString() << 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 ;
687 cout << "AliPHOSDigitizer not initialized " << endl ;
690 //__________________________________________________________________
691 void AliPHOSDigitizer::PrintDigits(Option_t * option){
692 // Print a table of digits
694 cout << "AliPHOSDigitiser:"<< endl ;
695 cout << " Number of entries in Digits list " << fDigits->GetEntriesFast() << endl ;
697 if(strstr(option,"all")){
700 AliPHOSDigit * digit;
701 cout << "Digit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
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() <<" ";
710 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
711 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
717 //__________________________________________________________________
718 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
720 // we set title (comment) of the SDigits branch in the first! header file
721 if(!fInitialized) Init() ;
723 ((TObjString*) fSDigitsTitles->At(0) )->SetString((char*)title) ;
726 //__________________________________________________________________
727 void AliPHOSDigitizer::SetDigitsBranch(const char* title)
729 //Sets the title (comment) of the branch to which Digits branch
730 if(!fInitialized) Init() ;
732 fDigitsTitle = title ;
735 //__________________________________________________________________