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 "TGeometry.h"
65 #include "TBenchmark.h"
67 // --- Standard library ---
69 // --- AliRoot header files ---
71 #include "AliHeader.h"
72 #include "AliStream.h"
73 #include "AliRunDigitizer.h"
74 #include "AliPHOSDigit.h"
76 #include "AliPHOSGetter.h"
77 #include "AliPHOSDigitizer.h"
78 #include "AliPHOSSDigitizer.h"
79 #include "AliPHOSGeometry.h"
80 #include "AliPHOSTick.h"
82 ClassImp(AliPHOSDigitizer)
85 //____________________________________________________________________________
86 AliPHOSDigitizer::AliPHOSDigitizer()
90 fDefaultInit = kTRUE ;
91 fManager = 0 ; // We work in the standalong mode
95 //____________________________________________________________________________
96 AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile, const char * name, const Bool_t toSplit)
100 SetTitle(headerFile) ;
102 fManager = 0 ; // We work in the standalong mode
107 fDefaultInit = kFALSE ;
110 //____________________________________________________________________________
111 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * ard):AliDigitizer(ard)
114 SetTitle(ard->GetInputFileName(0,0)) ;
116 fDefaultInit = kFALSE ;
119 if (ard->GetOutputFile()) {
120 SetName(ard->GetOutputFile().Data());
128 //____________________________________________________________________________
129 AliPHOSDigitizer::~AliPHOSDigitizer()
136 //____________________________________________________________________________
137 void AliPHOSDigitizer::Digitize(const Int_t event)
140 // Makes the digitization of the collected summable digits.
141 // It first creates the array of all PHOS modules
142 // filled with noise (different for EMC, CPV and PPSD) and
143 // then adds contributions from SDigits.
144 // This design avoids scanning over the list of digits to add
145 // contribution to new SDigits only.
147 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
148 TClonesArray * digits = gime->Digits(GetName()) ;
152 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 Error("Digitize", "SDigitizer with name %s not found", GetName() ) ;
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 TString fileName(folder->GetName()) ;
182 fileName.ReplaceAll("_","/") ;
183 // Info("Digitize", "Adding SDigits %s from %s", GetName(), fileName() ) ;
184 sdigArray->AddAt(sdigits, input) ;
189 //Find the first crystall with signal
190 Int_t nextSig = 200000 ;
192 for(i=0; i<input; i++){
193 sdigits = (TClonesArray *)sdigArray->At(i) ;
194 if ( !sdigits->GetEntriesFast() )
196 Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
197 if(curNext < nextSig)
201 TArrayI index(input) ;
202 index.Reset() ; //Set all indexes to zero
204 AliPHOSDigit * digit ;
205 AliPHOSDigit * curSDigit ;
207 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
209 //Put Noise contribution
210 for(absID = 1; absID <= nEMC; absID++){
211 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
212 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
213 //look if we have to add signal?
214 digit = (AliPHOSDigit *) digits->At(absID-1) ;
217 //Add SDigits from all inputs
220 Float_t a = digit->GetAmp() ;
221 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
222 //Mark the beginnign of the signal
223 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
224 //Mark the end of the ignal
225 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
228 for(i=0; i<input; i++){
229 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
230 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
233 //May be several digits will contribute from the same input
234 while(curSDigit && curSDigit->GetId() == absID){
235 //Shift primary to separate primaries belonging different inputs
236 Int_t primaryoffset ;
238 primaryoffset = fManager->GetMask(i) ;
240 primaryoffset = 10000000*i ;
241 curSDigit->ShiftPrimary(primaryoffset) ;
243 a = curSDigit->GetAmp() ;
244 b = a /fTimeSignalLength ;
245 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
246 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
248 *digit = *digit + *curSDigit ; //add energies
251 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
252 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
258 //calculate and set time
259 Float_t time = FrontEdgeTime(ticks) ;
260 digit->SetTime(time) ;
262 //Find next signal module
264 for(i=0; i<input; i++){
265 sdigits = ((TClonesArray *)sdigArray->At(i)) ;
266 Int_t curNext = nextSig ;
267 if(sdigits->GetEntriesFast() > index[i] ){
268 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
271 if(curNext < nextSig) nextSig = curNext ;
280 //Now CPV digits (different noise and no timing)
281 for(absID = nEMC+1; absID <= nCPV; absID++){
282 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
283 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
284 //look if we have to add signal?
286 digit = (AliPHOSDigit *) digits->At(absID-1) ;
287 //Add SDigits from all inputs
288 for(i=0; i<input; i++){
289 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
290 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
294 //May be several digits will contribute from the same input
295 while(curSDigit && curSDigit->GetId() == absID){
296 //Shift primary to separate primaries belonging different inputs
297 Int_t primaryoffset ;
299 primaryoffset = fManager->GetMask(i) ;
301 primaryoffset = 10000000*i ;
302 curSDigit->ShiftPrimary(primaryoffset) ;
305 *digit = *digit + *curSDigit ;
307 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
308 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
314 //Find next signal module
316 for(i=0; i<input; i++){
317 sdigits = (TClonesArray *)sdigArray->At(i) ;
318 Int_t curNext = nextSig ;
319 if(sdigits->GetEntriesFast() > index[i] )
320 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
321 if(curNext < nextSig) nextSig = curNext ;
326 delete sdigArray ; //We should not delete its contents
330 //remove digits below thresholds
331 for(i = 0; i < nEMC ; i++){
332 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
333 if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
334 digits->RemoveAt(i) ;
336 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
340 for(i = nEMC; i < nCPV ; i++)
341 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(i))->GetAmp()) < fCPVDigitThreshold)
342 digits->RemoveAt(i) ;
346 Int_t ndigits = digits->GetEntriesFast() ;
347 digits->Expand(ndigits) ;
349 //Set indexes in list of digits and make true digitization of the energy
350 for (i = 0 ; i < ndigits ; i++) {
351 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
352 digit->SetIndexInList(i) ;
353 Float_t energy = sDigitizer->Calibrate(digit->GetAmp()) ;
354 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
359 //____________________________________________________________________________
360 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
363 if(absId <= fEmcCrystals){ //digitize as EMC
364 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
365 if(chanel > fNADCemc ) chanel = fNADCemc ;
367 else{ //Digitize as CPV
368 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
369 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
374 //____________________________________________________________________________
375 void AliPHOSDigitizer::Exec(Option_t *option)
379 if(strcmp(GetName(), "") == 0 )
381 if (strstr(option,"print")) {
386 if(strstr(option,"tim"))
387 gBenchmark->Start("PHOSDigitizer");
389 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
396 treeD = fManager->GetTreeD() ;
397 nevents = 1 ; // Will process only one event
399 //Check, if this branch already exits
401 TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
403 TBranch * branch = 0 ;
404 Bool_t phosfound = kFALSE, digitizerfound = kFALSE ;
406 while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
407 if ( (strcmp(branch->GetName(), "PHOS")==0) &&
408 (strcmp(branch->GetTitle(), GetName())==0) )
411 else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) &&
412 (strcmp(branch->GetTitle(), GetName())==0) )
413 digitizerfound = kTRUE ;
417 Error("Exec", "Digits branch with name %s already exists", GetName() ) ;
420 if ( digitizerfound ) {
421 Error("Exec", "Digitizer branch with name %s already exists", GetName() );
426 else { //PHOS standalone
427 if(gime->BranchExists("Digits") )
429 nevents=gime->MaxEvent() ;
434 for(ievent = 0; ievent < nevents; ievent++){
439 for(input = 0 ; input < fManager->GetNinputs(); input ++){
440 TTree * treeS = fManager->GetInputTreeS(input) ;
442 Error("Exec", "No Input") ;
445 gime->ReadTreeS(treeS,input) ;
450 gime->Event(ievent,"S") ;
452 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
454 WriteDigits(ievent) ;
456 if(strstr(option,"deb"))
459 //increment the total number of Digits per run
460 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
463 if(strstr(option,"tim")){
464 gBenchmark->Stop("PHOSDigitizer");
466 message = " took %f seconds for Digitizing %f seconds per event\n" ;
467 Info("Exec", message.Data(),
468 gBenchmark->GetCpuTime("PHOSDigitizer"),
469 gBenchmark->GetCpuTime("PHOSDigitizer")/nevents );
473 //____________________________________________________________________________
474 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks)
476 ticks->Sort() ; //Sort in accordance with times of ticks
478 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
479 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
482 while((t=(AliPHOSTick*) it.Next())){
483 if(t->GetTime() < time) //This tick starts before crossing
488 time = ctick->CrossingTime(fTimeThreshold) ;
493 //____________________________________________________________________________
494 Bool_t AliPHOSDigitizer::Init()
496 // Makes all memory allocations
498 if( strcmp(GetTitle(), "") == 0 )
499 SetTitle("galice.root") ;
501 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName(), fToSplit) ;
503 Error("Init", "Could not obtain the Getter object !") ;
507 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
509 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
511 // Post Digits to the white board
512 gime->PostDigits(GetName() ) ;
514 // Post Digitizer to the white board
515 gime->PostDigitizer(this) ;
519 // construct the name of the file as /path/PHOS.SDigits.root
520 //First - extract full path if necessary
521 TString digitsFileName(GetTitle()) ;
522 Ssiz_t islash = digitsFileName.Last('/') ;
523 if(islash<digitsFileName.Length())
524 digitsFileName.Remove(islash+1,digitsFileName.Length()) ;
527 // Next - append the file name
528 digitsFileName+="PHOS.Digits." ;
529 if((strcmp(GetName(),"Default")!=0)&&(strcmp(GetName(),"")!=0)){
530 digitsFileName+=GetName() ;
531 digitsFileName+="." ;
533 digitsFileName+="root" ;
534 // Finally - check if the file already opened or open the file
535 fSplitFile = static_cast<TFile*>(gROOT->GetFile(digitsFileName.Data()));
537 fSplitFile = TFile::Open(digitsFileName.Data(),"update") ;
540 //Mark that we will use current header file
542 gime->PostSDigits(GetName(),GetTitle()) ;
543 gime->PostSDigitizer(GetName(),GetTitle()) ;
548 //____________________________________________________________________________
549 void AliPHOSDigitizer::InitParameters()
552 fEMCDigitThreshold = 0.012 ;
554 fCPVDigitThreshold = 0.09 ;
555 fTimeResolution = 0.5e-9 ;
556 fTimeSignalLength = 1.0e-9 ;
558 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
559 fADCpedestalEmc = 0.005 ; //
560 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
562 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
563 fADCpedestalCpv = 0.012 ; //
564 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
566 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
570 //__________________________________________________________________
571 void AliPHOSDigitizer::MixWith(const char* headerFile)
573 // Allows to produce digits by superimposing background and signal event.
574 // It is assumed, that headers file with SIGNAL events is opened in
576 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
577 // Thus we avoid writing (changing) huge and expensive
578 // backgound files: all output will be writen into SIGNAL, i.e.
579 // opened in constructor file.
581 // One can open as many files to mix with as one needs.
582 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
585 if( strcmp(GetName(), "") == 0 )
589 Warning("MixWith", "Can not use this method under AliRunDigitizer\n" ) ;
593 // check if the specified SDigits do not already exist on the White Board:
594 // //Folders/RunMC/Event/Data/PHOS/SDigits/headerFile/sdigitsname
596 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
597 TString path = gime->SDigitsFolder()->GetName() ;
599 // before it was ???? "Folders/RunMC/Event/Data/PHOS/SDigits" ;
603 if ( gROOT->FindObjectAny(path.Data()) ) {
604 Warning("MixWith", "Entry already exists, do not add\n" ) ;
608 gime->PostSDigits(GetName(),headerFile) ;
610 // check if the requested file is already open or exist and if SDigits Branch exist
611 TFile * file = (TFile*)gROOT->FindObject(headerFile);
613 file = new TFile(headerFile, "READ") ;
615 Error("MixWith", "File %s does not exist\n", headerFile );
621 //__________________________________________________________________
622 void AliPHOSDigitizer::Print(Option_t* option)const
624 // Print Digitizer's parameters
625 if( strcmp(GetName(), "") != 0 ){
627 TString message("\n-------------------") ;
628 message += GetName() ;
629 message += "-------------\n" ;
631 const Int_t nStreams = GetNInputStreams() ;
634 for (index = 0 ; index < nStreams ; index++) {
635 message += "Adding SDigits " ;
636 message += GetName() ;
638 message += fManager->GetInputFileName(index, 0) ;
640 message += "\nWriting digits to " ;
641 message += fManager->GetInputFileName(0, 0) ;
643 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
644 // gime->Folder("sdigits") ;
645 // TCollection * folderslist = gime->Folder("sdigits")->GetListOfFolders() ;
646 // TIter next(folderslist) ;
647 // TFolder * folder = 0 ;
649 // while ( (folder = (TFolder*)next()) ) {
650 // if ( folder->FindObject(GetName()) )
652 message += "\nWriting digits to " ;
653 message += GetTitle() ;
655 message += "\nWith following parameters:\n" ;
656 message += " Electronics noise in EMC (fPinNoise) = %f\n" ;
657 message += " Threshold in EMC (fEMCDigitThreshold) = %f\n" ;
658 message += " Noise in CPV (fCPVNoise) = %f\n" ;
659 message += " Threshold in CPV (fCPVDigitThreshold) = %f\n" ;
660 message += "---------------------------------------------------\n" ;
661 Info("Print", message.Data(),
665 fCPVDigitThreshold ) ;
668 Info("Print", "AliPHOSDigitizer not initialized\n" ) ;
672 //__________________________________________________________________
673 void AliPHOSDigitizer::PrintDigits(Option_t * option)
675 // Print a table of digits
677 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
678 TClonesArray * digits = gime->Digits() ;
681 message = "event %d\n" ;
682 message += " Number of entries in Digits list %d\n" ;
683 Info("Print", message.Data(), gAlice->GetEvNumber(), digits->GetEntriesFast() ) ;
684 if(strstr(option,"all")||strstr(option,"EMC")){
687 AliPHOSDigit * digit;
688 message = "EMC digits (with primaries):\n" ;
689 message += "Digit Id Amplitude Index Nprim Primaries list\n" ;
690 Info("Print", message.Data()) ;
691 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
693 for (index = 0 ; (index < digits->GetEntriesFast()) &&
694 (((AliPHOSDigit * ) digits->At(index))->GetId() <= maxEmc) ; index++) {
695 digit = (AliPHOSDigit * ) digits->At(index) ;
696 if(digit->GetNprimary() == 0)
698 message = " %d %d %d %d\n" ;
699 Info("Print", message.Data(),
702 digit->GetIndexInList(),
703 digit->GetNprimary() ) ;
706 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
708 Info("Print", message.Data(), digit->GetPrimary(iprimary+1) ) ;
713 if(strstr(option,"all")||strstr(option,"CPV")){
715 //loop over CPV digits
716 AliPHOSDigit * digit;
717 message = "CPV digits:\n" ;
718 message += "Digit Id Amplitude Index Nprim Primaries list\n" ;
719 Info("Print", message.Data()) ;
720 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
722 for (index = 0 ; index < digits->GetEntriesFast(); index++) {
723 digit = (AliPHOSDigit * ) digits->At(index) ;
724 if(digit->GetId() > maxEmc){
725 message = " %d %d %d %d\n" ;
726 Info("Print", message.Data(),
729 digit->GetIndexInList(),
730 digit->GetNprimary() ) ;
733 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
735 Info("Print", message.Data(), digit->GetPrimary(iprimary+1) ) ;
742 //__________________________________________________________________
743 Float_t AliPHOSDigitizer::TimeOfNoise(void)
744 { // Calculates the time signal generated by noise
745 //to be rewritten, now returns just big number
750 //____________________________________________________________________________
751 void AliPHOSDigitizer::WriteDigits(Int_t event)
754 // Makes TreeD in the output file.
755 // Check if branch already exists:
756 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
757 // already existing branches.
758 // else creates branch with Digits, named "PHOS", title "...",
759 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
760 // and names of files, from which digits are made.
762 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
763 const TClonesArray * digits = gime->Digits(GetName()) ;
767 treeD = fManager->GetTreeD() ;
769 if((gAlice->TreeD() == 0) || (fSplitFile)) // we should not create TreeD if it is already here
770 gAlice->MakeTree("D", fSplitFile); // We overwrite TreeD in split file in the case of second reconstruction
773 treeD = gAlice->TreeD();
776 // -- create Digits branch
777 Int_t bufferSize = 32000 ;
778 TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
779 digitsBranch->SetTitle(GetName());
781 // -- Create Digitizer branch
782 Int_t splitlevel = 0 ;
783 const AliPHOSDigitizer * d = dynamic_cast<const AliPHOSDigitizer *>(gime->Digitizer()) ;
784 TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel);
785 digitizerBranch->SetTitle(GetName());
787 digitsBranch->Fill() ;
788 digitizerBranch->Fill() ;