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 = 1.0e-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 fARD = 0 ; // We work in the standalong mode
108 //____________________________________________________________________________
109 AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile,const char * name)
113 SetTitle(headerFile) ;
114 fARD = 0 ; // We work in the standalong mode
119 //____________________________________________________________________________
120 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * ard)
125 SetTitle("aliroot") ;
130 //____________________________________________________________________________
131 AliPHOSDigitizer::~AliPHOSDigitizer()
138 //____________________________________________________________________________
139 void AliPHOSDigitizer::Digitize(const Int_t event)
142 // Makes the digitization of the collected summable digits.
143 // It first creates the array of all PHOS modules
144 // filled with noise (different for EMC, CPV and PPSD) and
145 // then adds contributions from SDigits.
146 // This design avoids scanning over the list of digits to add
147 // contribution to new SDigits only.
149 if( strcmp(GetName(), "") == 0 )
152 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
153 TClonesArray * digits = gime->Digits(GetName()) ;
154 AliPHOSSDigitizer * sdiz = gime->SDigitizer(GetName()) ;
158 const AliPHOSGeometry *geom = gime->PHOSGeometry() ;
160 //Making digits with noise, first EMC
161 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
165 TString name = geom->GetName() ;
167 nCPV = nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
168 geom->GetNModules() ;
170 digits->Expand(nCPV) ;
173 // sdigitize random gaussian noise and add it to all cells (EMCA+CPV+PPSD)
174 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
175 const AliPHOSSDigitizer * sDigitizer = gime->SDigitizer(GetName());
177 cerr << "ERROR: AliPHOSDigitizer::Digitize -> SDigitizer with name " << GetName() << " not found " << endl ;
181 // loop through the sdigits posted to the White Board and add them to the noise
182 TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("YSAlice/WhiteBoard/SDigits/PHOS"))->GetListOfFolders() ;
183 TIter next(folderslist) ;
184 TFolder * folder = 0 ;
185 TClonesArray * sdigits = 0 ;
187 TObjArray * sdigArray = new TObjArray(2) ;
188 while ( (folder = (TFolder*)next()) )
189 if ( (sdigits = (TClonesArray*)folder->FindObject(GetName()) ) ) {
190 cout << "INFO: AliPHOSDigitizer::Digitize -> Adding SDigits "
191 << GetName() << " from " << folder->GetName() << endl ;
192 sdigArray->AddAt(sdigits, input) ;
196 //Find the first crystall with sygnal
197 Int_t nextSig = 200000 ;
199 for(i=0; i<input; i++){
200 sdigits = (TClonesArray *)sdigArray->At(i) ;
201 Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
202 if(curNext < nextSig) nextSig = curNext ;
205 TArrayI index(input) ;
206 index.Reset() ; //Set all indexes to zero
208 AliPHOSDigit * digit ;
209 AliPHOSDigit * curSDigit ;
211 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
213 //Put Noise contribution
214 for(absID = 1; absID <= nEMC; absID++){
215 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
216 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
217 //look if we have to add signal?
219 //Add SDigits from all inputs
220 digit = (AliPHOSDigit *) digits->At(absID-1) ;
224 Float_t a = digit->GetAmp() ;
225 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
226 //Mark the beginnign of the signal
227 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
228 //Mark the end of the ignal
229 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
232 for(i=0; i<input; i++){
233 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
234 //May be several digits will contribute from the same input
235 while(curSDigit && curSDigit->GetId() == absID){
236 //Shift primary to separate primaries belonging different inputs
237 Int_t primaryoffset ;
239 primaryoffset = fARD->GetMask(i) ;
241 primaryoffset = 10000000*i ;
242 curSDigit->ShiftPrimary(primaryoffset) ;
244 a = curSDigit->GetAmp() ;
245 b = a /fTimeSignalLength ;
246 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
247 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
249 *digit = *digit + *curSDigit ; //add energies
252 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 = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
265 if(curNext < nextSig) nextSig = curNext ;
274 //Now CPV digits (different noise and no timing)
275 for(absID = nEMC+1; absID <= nCPV; absID++){
276 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
277 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
278 //look if we have to add signal?
280 digit = (AliPHOSDigit *) digits->At(absID-1) ;
281 //Add SDigits from all inputs
282 for(i=0; i<input; i++){
283 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
284 //May be several digits will contribute from the same input
285 while(curSDigit && curSDigit->GetId() == absID){
286 //Shift primary to separate primaries belonging different inputs
287 Int_t primaryoffset ;
289 primaryoffset = fARD->GetMask(i) ;
291 primaryoffset = 10000000*i ;
292 curSDigit->ShiftPrimary(primaryoffset) ;
295 *digit = *digit + *curSDigit ;
297 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
301 //Find next signal module
302 for(i=0; i<input; i++){
303 sdigits = (TClonesArray *)sdigArray->At(i) ;
304 Int_t curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
305 if(curNext < nextSig) nextSig = curNext ;
310 delete sdigArray ; //We should not delete its contents
314 //remove digits below thresholds
315 for(absID = 0; absID < nEMC ; absID++)
316 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(absID))->GetAmp()) < fEMCDigitThreshold)
317 digits->RemoveAt(absID) ;
319 for(absID = nEMC; absID < nCPV ; absID++)
320 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(absID))->GetAmp()) < fCPVDigitThreshold)
321 digits->RemoveAt(absID) ;
325 Int_t ndigits = digits->GetEntriesFast() ;
326 digits->Expand(ndigits) ;
328 //Set indexes in list of digits and make true digitization of the energy
329 for (i = 0 ; i < ndigits ; i++) {
330 AliPHOSDigit * digit = (AliPHOSDigit *) digits->At(i) ;
331 digit->SetIndexInList(i) ;
332 Float_t energy = sdiz->Calibrate(digit->GetAmp()) ;
333 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
337 //____________________________________________________________________________
338 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
341 if(absId <= fEmcCrystals){ //digitize as EMC
342 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
343 if(chanel > fNADCemc ) chanel = fNADCemc ;
345 else{ //Digitize as CPV
346 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
347 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
351 //____________________________________________________________________________
352 void AliPHOSDigitizer::Exec(Option_t *option)
356 if( strcmp(GetName(), "") == 0 )
359 if (strstr(option,"print")) {
364 if(strstr(option,"tim"))
365 gBenchmark->Start("PHOSDigitizer");
367 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
374 treeD = fARD->GetTreeD() ;
375 nevents = 1 ; // Will process only one event
378 gAlice->GetEvent(0) ;
379 nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
380 treeD=gAlice->TreeD() ;
383 cerr << " AliPHOSDigitizer :: Can not find TreeD " << endl ;
387 //Check, if this branch already exits
388 TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
390 TBranch * branch = 0 ;
391 Bool_t phosfound = kFALSE, digitizerfound = kFALSE ;
393 while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
394 if ( (strcmp(branch->GetName(), "PHOS")==0) &&
395 (strcmp(branch->GetTitle(), GetName())==0) )
398 else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) &&
399 (strcmp(branch->GetTitle(), GetName())==0) )
400 digitizerfound = kTRUE ;
404 cerr << "WARNING: AliPHOSDigitizer -> Digits branch with name " << GetName()
405 << " already exits" << endl ;
408 if ( digitizerfound ) {
409 cerr << "WARNING: AliPHOSDigitizer -> Digitizer branch with name " << GetName()
410 << " already exits" << endl ;
416 for(ievent = 0; ievent < nevents; ievent++){
420 for(input = 0 ; input < fARD->GetNinputs(); input ++){
421 TTree * treeS = fARD->GetInputTreeS(input) ;
423 cerr << "AliPHOSDigitizer -> No Input " << endl ;
426 gime->ReadTreeS(treeS,input) ;
430 gime->Event(ievent,"S") ;
432 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
434 WriteDigits(ievent) ;
436 if(strstr(option,"deb"))
440 if(strstr(option,"tim")){
441 gBenchmark->Stop("PHOSDigitizer");
442 cout << "AliPHOSDigitizer:" << endl ;
443 cout << " took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for Digitizing "
444 << gBenchmark->GetCpuTime("PHOSDigitizer")/nevents << " seconds per event " << endl ;
450 //____________________________________________________________________________
451 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks)
453 ticks->Sort() ; //Sort in accordance with times of ticks
455 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
456 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
459 while((t=(AliPHOSTick*) it.Next())){
460 if(t->GetTime() < time) //This tick starts before crossing
465 time = ctick->CrossingTime(fTimeThreshold) ;
469 //____________________________________________________________________________
470 Bool_t AliPHOSDigitizer::Init()
473 fEMCDigitThreshold = 0.01 ;
475 fCPVDigitThreshold = 0.09 ;
476 fTimeResolution = 1.0e-9 ;
477 fTimeSignalLength = 1.0e-9 ;
479 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
480 fADCpedestalEmc = 0.005 ; //
481 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
483 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
484 fADCpedestalCpv = 0.012 ; //
485 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
487 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
489 // Makes all memory allocations
490 // Adds Digitizer task to the folder of PHOS tasks
492 //============================================================= YS
493 // The initialisation is now done by AliPHOSGetter
495 if( strcmp(GetTitle(), "") == 0 )
496 SetTitle("galice.root") ;
498 // the SDigits name is stored by AliPHOSGetter as the name of the TClones Array
499 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/branchname and has branchTitle as title.
501 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
503 cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
507 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
508 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
510 // create a folder on the white board //YSAlice/WhiteBoard/Digits/PHOS/headerFile/digitsTitle
511 gime->PostDigits(GetName() ) ;
513 //add Task to //YSAlice/tasks/Digitizer/PHOS
514 gime->PostDigitizer(this) ;
516 //Mark that we will use current header file
518 gime->PostSDigits(GetName(),GetTitle()) ;
519 gime->PostSDigitizer(GetName(),GetTitle()) ;
524 //__________________________________________________________________
525 void AliPHOSDigitizer::MixWith(const char* headerFile)
527 // Allows to produce digits by superimposing background and signal event.
528 // It is assumed, that headers file with SIGNAL events is opened in
530 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
531 // Thus we avoid writing (changing) huge and expensive
532 // backgound files: all output will be writen into SIGNAL, i.e.
533 // opened in constructor file.
535 // One can open as many files to mix with as one needs.
536 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
539 if( strcmp(GetName(), "") == 0 )
543 cout << "Can not use this method under AliRunDigitizer " << endl ;
547 // check if the specified SDigits do not already exist on the White Board:
548 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/sDigitsTitle
550 TString path = "YSAlice/WhiteBoard/SDigits/PHOS/" ;
554 if ( gROOT->FindObjectAny(path.Data()) ) {
555 cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
559 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
560 gime->PostSDigits(GetName(),headerFile) ;
562 // check if the requested file is already open or exist and if SDigits Branch exist
563 TFile * file = (TFile*)gROOT->FindObject(headerFile);
565 file = new TFile(headerFile, "READ") ;
567 cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ;
574 //__________________________________________________________________
575 void AliPHOSDigitizer::Print(Option_t* option)const {
576 // Print Digitizer's parameters
577 if( strcmp(GetName(), "") != 0 ){
579 cout << "------------------- "<< GetName() << " -------------" << endl ;
580 cout << "Digitizing sDigits from file(s): " <<endl ;
582 TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("YSAlice/WhiteBoard/SDigits/PHOS"))->GetListOfFolders() ;
583 TIter next(folderslist) ;
584 TFolder * folder = 0 ;
586 while ( (folder = (TFolder*)next()) ) {
587 if ( folder->FindObject(GetName()) )
588 cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ;
591 cout << "Writing digits to " << GetTitle() << endl ;
594 cout << "With following parameters: " << endl ;
595 cout << " Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
596 cout << " Threshold in EMC (fEMCDigitThreshold) = " << fEMCDigitThreshold << endl ; ;
597 cout << " Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ;
598 cout << " Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ;
599 cout << "---------------------------------------------------" << endl ;
602 cout << "AliPHOSDigitizer not initialized " << endl ;
606 //__________________________________________________________________
607 void AliPHOSDigitizer::PrintDigits(Option_t * option){
608 // Print a table of digits
610 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
611 TClonesArray * digits = gime->Digits() ;
613 cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
614 cout << " Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
616 if(strstr(option,"all")||strstr(option,"EMC")){
619 AliPHOSDigit * digit;
620 cout << "Digit Id Amplitude Index " << " Nprim " << " Primaries list " << endl;
621 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
623 for (index = 0 ; (index < digits->GetEntries()) &&
624 (((AliPHOSDigit * ) digits->At(index))->GetId() <= maxEmc) ; index++) {
625 digit = (AliPHOSDigit * ) digits->At(index) ;
626 if(digit->GetNprimary() == 0) continue;
627 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
628 << setw(6) << digit->GetIndexInList() << " "
629 << setw(5) << digit->GetNprimary() <<" ";
632 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
633 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
639 if(strstr(option,"all")||strstr(option,"CPV")){
641 //loop over CPV digits
642 AliPHOSDigit * digit;
643 cout << "Digit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
644 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
646 for (index = 0 ; index < digits->GetEntries(); index++) {
647 digit = (AliPHOSDigit * ) digits->At(index) ;
648 if(digit->GetId() > maxEmc){
649 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
650 << setw(6) << digit->GetIndexInList() << " "
651 << setw(5) << digit->GetNprimary() <<" ";
654 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
655 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
663 //__________________________________________________________________
664 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
666 // we set title (comment) of the SDigits branch in the first! header file
667 if( strcmp(GetName(), "") == 0 )
670 AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ;
673 //__________________________________________________________________
674 Float_t AliPHOSDigitizer::TimeOfNoise(void)
675 { // Calculates the time signal generated by noise
676 //to be rewritten, now returns just big number
680 //____________________________________________________________________________
681 void AliPHOSDigitizer::Reset()
683 // sets current event number to the first simulated event
685 if( strcmp(GetName(), "") == 0 )
689 // for(inputs = 0; inputs < fNinputs ;inputs++)
690 // fIevent->AddAt(-1, inputs ) ;
694 //____________________________________________________________________________
695 void AliPHOSDigitizer::WriteDigits(Int_t event)
698 // Makes TreeD in the output file.
699 // Check if branch already exists:
700 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
701 // already existing branches.
702 // else creates branch with Digits, named "PHOS", title "...",
703 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
704 // and names of files, from which digits are made.
706 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
707 TClonesArray * digits = gime->Digits(GetName()) ;
712 treeD = fARD->GetTreeD() ;
714 treeD = gAlice->TreeD();
716 // create new branches
717 // -- generate file name if necessary
719 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
720 file = new char[strlen(gAlice->GetBaseFile())+20] ;
721 sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
724 TDirectory *cwd = gDirectory;
725 // -- create Digits branch
726 Int_t bufferSize = 32000 ;
727 TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
728 digitsBranch->SetTitle(GetName());
730 digitsBranch->SetFile(file);
731 TIter next( digitsBranch->GetListOfBranches());
733 while ((sbr=(TBranch*)next())) {
739 // -- Create Digitizer branch
740 Int_t splitlevel = 0 ;
741 AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
742 TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel);
743 digitizerBranch->SetTitle(GetName());
745 digitizerBranch->SetFile(file);
746 TIter next( digitizerBranch->GetListOfBranches());
748 while ((sbr=(TBranch*)next())) {
754 digitsBranch->Fill() ;
755 digitizerBranch->Fill() ;
757 treeD->Write(0,kOverwrite) ;