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.
24 // The name of the TTask is also the title of the branch that will contain
25 // the created SDigits
26 // The title of the TTAsk is the name of the file that contains the hits from
27 // which the SDigits are created
29 // For each event two branches are created in TreeD:
30 // "PHOS" - list of digits
31 // "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
33 // Note, that one can set a title for new digits branch, and repeat digitization with
34 // another set of parameters.
37 // root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
38 // root[1] d->ExecuteTask()
39 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
40 // //Digitizes SDigitis in all events found in file galice.root
42 // root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;
43 // // Will read sdigits from galice1.root
44 // root[3] d1->MixWith("galice2.root")
45 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
46 // // Reads another set of sdigits from galice2.root
47 // root[3] d1->MixWith("galice3.root")
48 // // Reads another set of sdigits from galice3.root
49 // root[4] d->ExecuteTask("deb timing")
50 // // Reads SDigits from files galice1.root, galice2.root ....
51 // // mixes them and stores produced Digits in file galice1.root
52 // // deb - prints number of produced digits
53 // // deb all - prints list of produced digits
54 // // timing - prints time used for digitization
57 // --- ROOT system ---
63 #include "TObjString.h"
64 #include "TBenchmark.h"
66 // --- Standard library ---
69 // --- AliRoot header files ---
72 #include "AliRunDigitizer.h"
73 #include "AliPHOSDigit.h"
75 #include "AliPHOSGetter.h"
76 #include "AliPHOSDigitizer.h"
77 #include "AliPHOSSDigitizer.h"
78 #include "AliPHOSGeometry.h"
79 #include "AliPHOSTick.h"
81 ClassImp(AliPHOSDigitizer)
84 //____________________________________________________________________________
85 AliPHOSDigitizer::AliPHOSDigitizer()
90 fEMCDigitThreshold = 0.01 ;
92 fCPVDigitThreshold = 0.09 ;
93 fTimeResolution = 0.5e-9 ;
94 fTimeSignalLength = 1.0e-9 ;
96 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
97 fADCpedestalEmc = 0.005 ; //
98 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
100 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
101 fADCpedestalCpv = 0.012 ; //
102 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
104 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
105 fManager = 0 ; // We work in the standalong mode
108 //____________________________________________________________________________
109 AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile,const char * name)
113 SetTitle(headerFile) ;
114 fManager = 0 ; // We work in the standalong mode
119 //____________________________________________________________________________
120 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * ard):AliDigitizer(ard)
123 SetName(""); //Will call init in the digitizing
124 SetTitle("aliroot") ;
127 //____________________________________________________________________________
128 AliPHOSDigitizer::~AliPHOSDigitizer()
135 //____________________________________________________________________________
136 void AliPHOSDigitizer::Digitize(const Int_t event)
139 // Makes the digitization of the collected summable digits.
140 // It first creates the array of all PHOS modules
141 // filled with noise (different for EMC, CPV and PPSD) and
142 // then adds contributions from SDigits.
143 // This design avoids scanning over the list of digits to add
144 // contribution to new SDigits only.
146 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
147 TClonesArray * digits = gime->Digits(GetName()) ;
151 const AliPHOSGeometry *geom = gime->PHOSGeometry() ;
153 //Making digits with noise, first EMC
154 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
158 TString name = geom->GetName() ;
160 nCPV = nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
161 geom->GetNModules() ;
163 digits->Expand(nCPV) ;
165 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
166 const AliPHOSSDigitizer * sDigitizer = gime->SDigitizer(GetName());
168 cerr << "ERROR: AliPHOSDigitizer::Digitize -> SDigitizer with name " << GetName() << " not found " << endl ;
172 // loop through the sdigits posted to the White Board and add them to the noise
173 TCollection * folderslist = gime->SDigitsFolder()->GetListOfFolders() ;
174 TIter next(folderslist) ;
175 TFolder * folder = 0 ;
176 TClonesArray * sdigits = 0 ;
178 TObjArray * sdigArray = new TObjArray(2) ;
179 while ( (folder = (TFolder*)next()) )
180 if ( (sdigits = (TClonesArray*)folder->FindObject(GetName()) ) ) {
181 cout << "INFO: AliPHOSDigitizer::Digitize -> Adding SDigits "
182 << GetName() << " from " << folder->GetName() << endl ;
183 sdigArray->AddAt(sdigits, input) ;
187 //Find the first crystall with signal
188 Int_t nextSig = 200000 ;
190 for(i=0; i<input; i++){
191 sdigits = (TClonesArray *)sdigArray->At(i) ;
192 if ( !sdigits->GetEntriesFast() )
194 Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
195 if(curNext < nextSig)
199 TArrayI index(input) ;
200 index.Reset() ; //Set all indexes to zero
202 AliPHOSDigit * digit ;
203 AliPHOSDigit * curSDigit ;
205 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
207 //Put Noise contribution
208 for(absID = 1; absID <= nEMC; absID++){
209 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
210 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
211 //look if we have to add signal?
213 //Add SDigits from all inputs
214 digit = (AliPHOSDigit *) digits->At(absID-1) ;
218 Float_t a = digit->GetAmp() ;
219 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
220 //Mark the beginnign of the signal
221 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
222 //Mark the end of the ignal
223 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
226 for(i=0; i<input; i++){
227 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
228 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
231 //May be several digits will contribute from the same input
232 while(curSDigit && curSDigit->GetId() == absID){
233 //Shift primary to separate primaries belonging different inputs
234 Int_t primaryoffset ;
236 primaryoffset = fManager->GetMask(i) ;
239 curSDigit->ShiftPrimary(i) ;
241 a = curSDigit->GetAmp() ;
242 b = a /fTimeSignalLength ;
243 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
244 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
246 *digit = *digit + *curSDigit ; //add energies
249 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
250 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
256 //calculate and set time
257 Float_t time = FrontEdgeTime(ticks) ;
258 digit->SetTime(time) ;
260 //Find next signal module
262 for(i=0; i<input; i++){
263 sdigits = ((TClonesArray *)sdigArray->At(i)) ;
264 Int_t curNext = nextSig ;
265 if(sdigits->GetEntriesFast() > index[i] ){
266 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
269 if(curNext < nextSig) nextSig = curNext ;
278 //Now CPV digits (different noise and no timing)
279 for(absID = nEMC+1; absID <= nCPV; absID++){
280 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
281 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
282 //look if we have to add signal?
284 digit = (AliPHOSDigit *) digits->At(absID-1) ;
285 //Add SDigits from all inputs
286 for(i=0; i<input; i++){
287 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
288 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
292 //May be several digits will contribute from the same input
293 while(curSDigit && curSDigit->GetId() == absID){
294 //Shift primary to separate primaries belonging different inputs
295 Int_t primaryoffset ;
297 primaryoffset = fManager->GetMask(i) ;
299 primaryoffset = 10000000*i ;
300 curSDigit->ShiftPrimary(primaryoffset) ;
303 *digit = *digit + *curSDigit ;
305 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
306 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
312 //Find next signal module
314 for(i=0; i<input; i++){
315 sdigits = (TClonesArray *)sdigArray->At(i) ;
316 Int_t curNext = nextSig ;
317 if(sdigits->GetEntriesFast() > index[i] )
318 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
319 if(curNext < nextSig) nextSig = curNext ;
324 delete sdigArray ; //We should not delete its contents
328 //remove digits below thresholds
329 for(i = 0; i < nEMC ; i++){
330 digit = (AliPHOSDigit*) digits->At(i) ;
331 if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
332 digits->RemoveAt(i) ;
334 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
338 for(i = nEMC; i < nCPV ; i++)
339 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(i))->GetAmp()) < fCPVDigitThreshold)
340 digits->RemoveAt(i) ;
344 Int_t ndigits = digits->GetEntriesFast() ;
345 digits->Expand(ndigits) ;
347 //Set indexes in list of digits and make true digitization of the energy
348 for (i = 0 ; i < ndigits ; i++) {
349 AliPHOSDigit * digit = (AliPHOSDigit *) digits->At(i) ;
350 digit->SetIndexInList(i) ;
351 Float_t energy = sDigitizer->Calibrate(digit->GetAmp()) ;
352 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
356 //____________________________________________________________________________
357 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
360 if(absId <= fEmcCrystals){ //digitize as EMC
361 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
362 if(chanel > fNADCemc ) chanel = fNADCemc ;
364 else{ //Digitize as CPV
365 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
366 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
370 //____________________________________________________________________________
371 void AliPHOSDigitizer::Exec(Option_t *option)
375 if(strcmp(GetName(), "") == 0 )
378 if (strstr(option,"print")) {
383 if(strstr(option,"tim"))
384 gBenchmark->Start("PHOSDigitizer");
386 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
393 treeD = fManager->GetTreeD() ;
394 nevents = 1 ; // Will process only one event
397 gAlice->GetEvent(0) ;
398 nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
399 treeD=gAlice->TreeD() ;
402 cerr << " AliPHOSDigitizer :: Can not find TreeD " << endl ;
406 //Check, if this branch already exits
407 TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
409 TBranch * branch = 0 ;
410 Bool_t phosfound = kFALSE, digitizerfound = kFALSE ;
412 while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
413 if ( (strcmp(branch->GetName(), "PHOS")==0) &&
414 (strcmp(branch->GetTitle(), GetName())==0) )
417 else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) &&
418 (strcmp(branch->GetTitle(), GetName())==0) )
419 digitizerfound = kTRUE ;
423 cerr << "WARNING: AliPHOSDigitizer -> Digits branch with name " << GetName()
424 << " already exits" << endl ;
427 if ( digitizerfound ) {
428 cerr << "WARNING: AliPHOSDigitizer -> Digitizer branch with name " << GetName()
429 << " already exits" << endl ;
435 for(ievent = 0; ievent < nevents; ievent++){
439 for(input = 0 ; input < fManager->GetNinputs(); input ++){
440 TTree * treeS = fManager->GetInputTreeS(input) ;
442 cerr << "AliPHOSDigitizer -> No Input " << endl ;
445 gime->ReadTreeS(treeS,input) ;
449 gime->Event(ievent,"S") ;
451 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
453 WriteDigits(ievent) ;
455 if(strstr(option,"deb"))
458 //increment the total number of Digits per run
459 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
462 if(strstr(option,"tim")){
463 gBenchmark->Stop("PHOSDigitizer");
464 cout << "AliPHOSDigitizer:" << endl ;
465 cout << " took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for Digitizing "
466 << gBenchmark->GetCpuTime("PHOSDigitizer")/nevents << " seconds per event " << endl ;
472 //____________________________________________________________________________
473 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks)
475 ticks->Sort() ; //Sort in accordance with times of ticks
477 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
478 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
481 while((t=(AliPHOSTick*) it.Next())){
482 if(t->GetTime() < time) //This tick starts before crossing
487 time = ctick->CrossingTime(fTimeThreshold) ;
491 //____________________________________________________________________________
492 Bool_t AliPHOSDigitizer::Init()
495 fEMCDigitThreshold = 0.01 ;
497 fCPVDigitThreshold = 0.09 ;
498 fTimeResolution = 0.5e-9 ;
499 fTimeSignalLength = 1.0e-9 ;
501 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
502 fADCpedestalEmc = 0.005 ; //
503 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
505 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
506 fADCpedestalCpv = 0.012 ; //
507 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
509 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
512 SetTitle("aliroot") ;
514 // SetTitle("galice.root") ;
516 if( strcmp(GetName(), "") == 0 )
519 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
521 cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
525 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
526 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
528 // Post Digits to the white board
529 gime->PostDigits(GetName() ) ;
531 // Post Digitizer to the white board
532 gime->PostDigitizer(this) ;
534 //Mark that we will use current header file
536 gime->PostSDigits(GetName(),GetTitle()) ;
537 gime->PostSDigitizer(GetName(),GetTitle()) ;
542 //__________________________________________________________________
543 void AliPHOSDigitizer::MixWith(const char* headerFile)
545 // Allows to produce digits by superimposing background and signal event.
546 // It is assumed, that headers file with SIGNAL events is opened in
548 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
549 // Thus we avoid writing (changing) huge and expensive
550 // backgound files: all output will be writen into SIGNAL, i.e.
551 // opened in constructor file.
553 // One can open as many files to mix with as one needs.
554 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
557 if( strcmp(GetName(), "") == 0 )
561 cout << "Can not use this method under AliRunDigitizer " << endl ;
565 // check if the specified SDigits do not already exist on the White Board:
566 // //Folders/RunMC/Event/Data/PHOS/SDigits/headerFile/sdigitsname
568 TString path = "Folders/RunMC/Event/Data/PHOS/SDigits" ;
572 if ( gROOT->FindObjectAny(path.Data()) ) {
573 cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
577 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
578 gime->PostSDigits(GetName(),headerFile) ;
580 // check if the requested file is already open or exist and if SDigits Branch exist
581 TFile * file = (TFile*)gROOT->FindObject(headerFile);
583 file = new TFile(headerFile, "READ") ;
585 cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ;
592 //__________________________________________________________________
593 void AliPHOSDigitizer::Print(Option_t* option)const {
594 // Print Digitizer's parameters
595 if( strcmp(GetName(), "") != 0 ){
597 cout << "------------------- "<< GetName() << " -------------" << endl ;
598 cout << "Digitizing sDigits from file(s): " <<endl ;
600 TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("Folders/RunMC/Event/Data/PHOS/SDigits"))->GetListOfFolders() ;
601 TIter next(folderslist) ;
602 TFolder * folder = 0 ;
604 while ( (folder = (TFolder*)next()) ) {
605 if ( folder->FindObject(GetName()) )
606 cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ;
609 cout << "Writing digits to " << GetTitle() << endl ;
612 cout << "With following parameters: " << endl ;
613 cout << " Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
614 cout << " Threshold in EMC (fEMCDigitThreshold) = " << fEMCDigitThreshold << endl ; ;
615 cout << " Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ;
616 cout << " Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ;
617 cout << "---------------------------------------------------" << endl ;
620 cout << "AliPHOSDigitizer not initialized " << endl ;
624 //__________________________________________________________________
625 void AliPHOSDigitizer::PrintDigits(Option_t * option){
626 // Print a table of digits
628 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
629 TClonesArray * digits = gime->Digits() ;
631 cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
632 cout << " Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
634 if(strstr(option,"all")||strstr(option,"EMC")){
637 AliPHOSDigit * digit;
638 cout << "EMC digits (with primaries): " << endl ;
639 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
640 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
642 for (index = 0 ; (index < digits->GetEntriesFast()) &&
643 (((AliPHOSDigit * ) digits->At(index))->GetId() <= maxEmc) ; index++) {
644 digit = (AliPHOSDigit * ) digits->At(index) ;
645 if(digit->GetNprimary() == 0) continue;
646 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
647 << setw(6) << digit->GetIndexInList() << " "
648 << setw(5) << digit->GetNprimary() <<" ";
651 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
652 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
658 if(strstr(option,"all")||strstr(option,"CPV")){
660 //loop over CPV digits
661 AliPHOSDigit * digit;
662 cout << "CPV digits: " << endl ;
663 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
664 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
666 for (index = 0 ; index < digits->GetEntriesFast(); index++) {
667 digit = (AliPHOSDigit * ) digits->At(index) ;
668 if(digit->GetId() > maxEmc){
669 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
670 << setw(6) << digit->GetIndexInList() << " "
671 << setw(5) << digit->GetNprimary() <<" ";
674 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
675 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
683 //__________________________________________________________________
684 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
686 // we set title (comment) of the SDigits branch in the first! header file
687 if( strcmp(GetName(), "") == 0 )
690 AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ;
693 //__________________________________________________________________
694 Float_t AliPHOSDigitizer::TimeOfNoise(void)
695 { // Calculates the time signal generated by noise
696 //to be rewritten, now returns just big number
700 //____________________________________________________________________________
701 void AliPHOSDigitizer::Reset()
703 // sets current event number to the first simulated event
705 if( strcmp(GetName(), "") == 0 )
709 // for(inputs = 0; inputs < fNinputs ;inputs++)
710 // fIevent->AddAt(-1, inputs ) ;
714 //____________________________________________________________________________
715 void AliPHOSDigitizer::WriteDigits(Int_t event)
718 // Makes TreeD in the output file.
719 // Check if branch already exists:
720 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
721 // already existing branches.
722 // else creates branch with Digits, named "PHOS", title "...",
723 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
724 // and names of files, from which digits are made.
726 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
727 const TClonesArray * digits = gime->Digits(GetName()) ;
731 treeD = fManager->GetTreeD() ;
733 treeD = gAlice->TreeD();
735 // create new branches
736 // -- generate file name if necessary
738 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
739 file = new char[strlen(gAlice->GetBaseFile())+20] ;
740 sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
743 TDirectory *cwd = gDirectory;
744 // -- create Digits branch
745 Int_t bufferSize = 32000 ;
746 TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
747 digitsBranch->SetTitle(GetName());
749 digitsBranch->SetFile(file);
750 TIter next( digitsBranch->GetListOfBranches());
752 while ((sbr=(TBranch*)next())) {
758 // -- Create Digitizer branch
759 Int_t splitlevel = 0 ;
760 AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
761 TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel);
762 digitizerBranch->SetTitle(GetName());
764 digitizerBranch->SetFile(file);
765 TIter next( digitizerBranch->GetListOfBranches());
767 while ((sbr=(TBranch*)next())) {
773 digitsBranch->Fill() ;
774 digitizerBranch->Fill() ;
776 treeD->Write(0,kOverwrite) ;