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 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 = gime->SDigitsFolder()->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 signal
197 Int_t nextSig = 200000 ;
199 for(i=0; i<input; i++){
200 sdigits = (TClonesArray *)sdigArray->At(i) ;
201 if ( !sdigits->GetEntries() )
203 Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
204 if(curNext < nextSig)
208 TArrayI index(input) ;
209 index.Reset() ; //Set all indexes to zero
211 AliPHOSDigit * digit ;
212 AliPHOSDigit * curSDigit ;
214 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
216 //Put Noise contribution
217 for(absID = 1; absID <= nEMC; absID++){
218 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
219 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
220 //look if we have to add signal?
222 //Add SDigits from all inputs
223 digit = (AliPHOSDigit *) digits->At(absID-1) ;
227 Float_t a = digit->GetAmp() ;
228 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
229 //Mark the beginnign of the signal
230 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
231 //Mark the end of the ignal
232 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
235 for(i=0; i<input; i++){
236 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
237 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
240 //May be several digits will contribute from the same input
241 while(curSDigit && curSDigit->GetId() == absID){
242 //Shift primary to separate primaries belonging different inputs
243 Int_t primaryoffset ;
245 primaryoffset = fARD->GetMask(i) ;
247 primaryoffset = 10000000*i ;
248 curSDigit->ShiftPrimary(primaryoffset) ;
250 a = curSDigit->GetAmp() ;
251 b = a /fTimeSignalLength ;
252 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
253 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
255 *digit = *digit + *curSDigit ; //add energies
258 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
259 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
265 //calculate and set time
266 Float_t time = FrontEdgeTime(ticks) ;
267 digit->SetTime(time) ;
269 //Find next signal module
271 for(i=0; i<input; i++){
272 sdigits = ((TClonesArray *)sdigArray->At(i)) ;
273 Int_t curNext = nextSig ;
274 if(sdigits->GetEntriesFast() > index[i] ){
275 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
278 if(curNext < nextSig) nextSig = curNext ;
287 //Now CPV digits (different noise and no timing)
288 for(absID = nEMC+1; absID <= nCPV; absID++){
289 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
290 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
291 //look if we have to add signal?
293 digit = (AliPHOSDigit *) digits->At(absID-1) ;
294 //Add SDigits from all inputs
295 for(i=0; i<input; i++){
296 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
297 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
300 //May be several digits will contribute from the same input
301 while(curSDigit && curSDigit->GetId() == absID){
302 //Shift primary to separate primaries belonging different inputs
303 Int_t primaryoffset ;
305 primaryoffset = fARD->GetMask(i) ;
307 primaryoffset = 10000000*i ;
308 curSDigit->ShiftPrimary(primaryoffset) ;
311 *digit = *digit + *curSDigit ;
313 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
314 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
320 //Find next signal module
321 for(i=0; i<input; i++){
322 sdigits = (TClonesArray *)sdigArray->At(i) ;
323 Int_t curNext = nextSig ;
324 if(sdigits->GetEntriesFast() > index[i] )
325 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
326 if(curNext < nextSig) nextSig = curNext ;
331 delete sdigArray ; //We should not delete its contents
335 //remove digits below thresholds
336 for(absID = 0; absID < nEMC ; absID++){
337 digit = (AliPHOSDigit*) digits->At(absID) ;
338 if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
339 digits->RemoveAt(absID) ;
341 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
344 for(absID = nEMC; absID < nCPV ; absID++)
345 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(absID))->GetAmp()) < fCPVDigitThreshold)
346 digits->RemoveAt(absID) ;
350 Int_t ndigits = digits->GetEntriesFast() ;
351 digits->Expand(ndigits) ;
353 //Set indexes in list of digits and make true digitization of the energy
354 for (i = 0 ; i < ndigits ; i++) {
355 AliPHOSDigit * digit = (AliPHOSDigit *) digits->At(i) ;
356 digit->SetIndexInList(i) ;
357 Float_t energy = sdiz->Calibrate(digit->GetAmp()) ;
358 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
362 //____________________________________________________________________________
363 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
366 if(absId <= fEmcCrystals){ //digitize as EMC
367 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
368 if(chanel > fNADCemc ) chanel = fNADCemc ;
370 else{ //Digitize as CPV
371 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
372 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
376 //____________________________________________________________________________
377 void AliPHOSDigitizer::Exec(Option_t *option)
381 if( strcmp(GetName(), "") == 0 )
384 if (strstr(option,"print")) {
389 if(strstr(option,"tim"))
390 gBenchmark->Start("PHOSDigitizer");
392 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
399 treeD = fARD->GetTreeD() ;
400 nevents = 1 ; // Will process only one event
403 gAlice->GetEvent(0) ;
404 nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
405 treeD=gAlice->TreeD() ;
408 cerr << " AliPHOSDigitizer :: Can not find TreeD " << endl ;
412 //Check, if this branch already exits
413 TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
415 TBranch * branch = 0 ;
416 Bool_t phosfound = kFALSE, digitizerfound = kFALSE ;
418 while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
419 if ( (strcmp(branch->GetName(), "PHOS")==0) &&
420 (strcmp(branch->GetTitle(), GetName())==0) )
423 else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) &&
424 (strcmp(branch->GetTitle(), GetName())==0) )
425 digitizerfound = kTRUE ;
429 cerr << "WARNING: AliPHOSDigitizer -> Digits branch with name " << GetName()
430 << " already exits" << endl ;
433 if ( digitizerfound ) {
434 cerr << "WARNING: AliPHOSDigitizer -> Digitizer branch with name " << GetName()
435 << " already exits" << endl ;
441 for(ievent = 0; ievent < nevents; ievent++){
445 for(input = 0 ; input < fARD->GetNinputs(); input ++){
446 TTree * treeS = fARD->GetInputTreeS(input) ;
448 cerr << "AliPHOSDigitizer -> No Input " << endl ;
451 gime->ReadTreeS(treeS,input) ;
455 gime->Event(ievent,"S") ;
457 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
459 WriteDigits(ievent) ;
461 if(strstr(option,"deb"))
465 if(strstr(option,"tim")){
466 gBenchmark->Stop("PHOSDigitizer");
467 cout << "AliPHOSDigitizer:" << endl ;
468 cout << " took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for Digitizing "
469 << gBenchmark->GetCpuTime("PHOSDigitizer")/nevents << " seconds per event " << endl ;
475 //____________________________________________________________________________
476 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks)
478 ticks->Sort() ; //Sort in accordance with times of ticks
480 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
481 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
484 while((t=(AliPHOSTick*) it.Next())){
485 if(t->GetTime() < time) //This tick starts before crossing
490 time = ctick->CrossingTime(fTimeThreshold) ;
494 //____________________________________________________________________________
495 Bool_t AliPHOSDigitizer::Init()
498 fEMCDigitThreshold = 0.01 ;
500 fCPVDigitThreshold = 0.09 ;
501 fTimeResolution = 0.5e-9 ;
502 fTimeSignalLength = 1.0e-9 ;
504 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
505 fADCpedestalEmc = 0.005 ; //
506 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
508 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
509 fADCpedestalCpv = 0.012 ; //
510 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
512 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
514 // Makes all memory allocations
515 // Adds Digitizer task to the folder of PHOS tasks
517 //============================================================= YS
518 // The initialisation is now done by AliPHOSGetter
520 if( strcmp(GetTitle(), "") == 0 )
521 SetTitle("galice.root") ;
523 // the SDigits name is stored by AliPHOSGetter as the name of the TClones Array
524 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/branchname and has branchTitle as title.
526 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
528 cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
532 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
533 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
535 // create a folder on the white board //YSAlice/WhiteBoard/Digits/PHOS/headerFile/digitsTitle
536 gime->PostDigits(GetName() ) ;
538 //add Task to //YSAlice/tasks/Digitizer/PHOS
539 gime->PostDigitizer(this) ;
541 //Mark that we will use current header file
543 gime->PostSDigits(GetName(),GetTitle()) ;
544 gime->PostSDigitizer(GetName(),GetTitle()) ;
549 //__________________________________________________________________
550 void AliPHOSDigitizer::MixWith(const char* headerFile)
552 // Allows to produce digits by superimposing background and signal event.
553 // It is assumed, that headers file with SIGNAL events is opened in
555 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
556 // Thus we avoid writing (changing) huge and expensive
557 // backgound files: all output will be writen into SIGNAL, i.e.
558 // opened in constructor file.
560 // One can open as many files to mix with as one needs.
561 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
564 if( strcmp(GetName(), "") == 0 )
568 cout << "Can not use this method under AliRunDigitizer " << endl ;
572 // check if the specified SDigits do not already exist on the White Board:
573 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/sDigitsTitle
575 TString path = "YSAlice/WhiteBoard/SDigits/PHOS/" ;
579 if ( gROOT->FindObjectAny(path.Data()) ) {
580 cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
584 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
585 gime->PostSDigits(GetName(),headerFile) ;
587 // check if the requested file is already open or exist and if SDigits Branch exist
588 TFile * file = (TFile*)gROOT->FindObject(headerFile);
590 file = new TFile(headerFile, "READ") ;
592 cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ;
599 //__________________________________________________________________
600 void AliPHOSDigitizer::Print(Option_t* option)const {
601 // Print Digitizer's parameters
602 if( strcmp(GetName(), "") != 0 ){
604 cout << "------------------- "<< GetName() << " -------------" << endl ;
605 cout << "Digitizing sDigits from file(s): " <<endl ;
607 TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("YSAlice/WhiteBoard/SDigits/PHOS"))->GetListOfFolders() ;
608 TIter next(folderslist) ;
609 TFolder * folder = 0 ;
611 while ( (folder = (TFolder*)next()) ) {
612 if ( folder->FindObject(GetName()) )
613 cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ;
616 cout << "Writing digits to " << GetTitle() << endl ;
619 cout << "With following parameters: " << endl ;
620 cout << " Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
621 cout << " Threshold in EMC (fEMCDigitThreshold) = " << fEMCDigitThreshold << endl ; ;
622 cout << " Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ;
623 cout << " Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ;
624 cout << "---------------------------------------------------" << endl ;
627 cout << "AliPHOSDigitizer not initialized " << endl ;
631 //__________________________________________________________________
632 void AliPHOSDigitizer::PrintDigits(Option_t * option){
633 // Print a table of digits
635 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
636 TClonesArray * digits = gime->Digits() ;
638 cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
639 cout << " Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
641 if(strstr(option,"all")||strstr(option,"EMC")){
644 AliPHOSDigit * digit;
645 cout << "Digit Id Amplitude Index " << " Nprim " << " Primaries list " << endl;
646 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
648 for (index = 0 ; (index < digits->GetEntries()) &&
649 (((AliPHOSDigit * ) digits->At(index))->GetId() <= maxEmc) ; index++) {
650 digit = (AliPHOSDigit * ) digits->At(index) ;
651 if(digit->GetNprimary() == 0) continue;
652 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
653 << setw(6) << digit->GetIndexInList() << " "
654 << setw(5) << digit->GetNprimary() <<" ";
657 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
658 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
664 if(strstr(option,"all")||strstr(option,"CPV")){
666 //loop over CPV digits
667 AliPHOSDigit * digit;
668 cout << "Digit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
669 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
671 for (index = 0 ; index < digits->GetEntries(); index++) {
672 digit = (AliPHOSDigit * ) digits->At(index) ;
673 if(digit->GetId() > maxEmc){
674 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
675 << setw(6) << digit->GetIndexInList() << " "
676 << setw(5) << digit->GetNprimary() <<" ";
679 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
680 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
688 //__________________________________________________________________
689 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
691 // we set title (comment) of the SDigits branch in the first! header file
692 if( strcmp(GetName(), "") == 0 )
695 AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ;
698 //__________________________________________________________________
699 Float_t AliPHOSDigitizer::TimeOfNoise(void)
700 { // Calculates the time signal generated by noise
701 //to be rewritten, now returns just big number
705 //____________________________________________________________________________
706 void AliPHOSDigitizer::Reset()
708 // sets current event number to the first simulated event
710 if( strcmp(GetName(), "") == 0 )
714 // for(inputs = 0; inputs < fNinputs ;inputs++)
715 // fIevent->AddAt(-1, inputs ) ;
719 //____________________________________________________________________________
720 void AliPHOSDigitizer::WriteDigits(Int_t event)
723 // Makes TreeD in the output file.
724 // Check if branch already exists:
725 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
726 // already existing branches.
727 // else creates branch with Digits, named "PHOS", title "...",
728 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
729 // and names of files, from which digits are made.
731 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
732 TClonesArray * digits = gime->Digits(GetName()) ;
737 treeD = fARD->GetTreeD() ;
739 treeD = gAlice->TreeD();
741 // create new branches
742 // -- generate file name if necessary
744 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
745 file = new char[strlen(gAlice->GetBaseFile())+20] ;
746 sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
749 TDirectory *cwd = gDirectory;
750 // -- create Digits branch
751 Int_t bufferSize = 32000 ;
752 TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
753 digitsBranch->SetTitle(GetName());
755 digitsBranch->SetFile(file);
756 TIter next( digitsBranch->GetListOfBranches());
758 while ((sbr=(TBranch*)next())) {
764 // -- Create Digitizer branch
765 Int_t splitlevel = 0 ;
766 AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
767 TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel);
768 digitizerBranch->SetTitle(GetName());
770 digitizerBranch->SetFile(file);
771 TIter next( digitizerBranch->GetListOfBranches());
773 while ((sbr=(TBranch*)next())) {
779 digitsBranch->Fill() ;
780 digitizerBranch->Fill() ;
782 treeD->Write(0,kOverwrite) ;