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") ;
129 //____________________________________________________________________________
130 AliPHOSDigitizer::~AliPHOSDigitizer()
137 //____________________________________________________________________________
138 void AliPHOSDigitizer::Digitize(const Int_t event)
141 // Makes the digitization of the collected summable digits.
142 // It first creates the array of all PHOS modules
143 // filled with noise (different for EMC, CPV and PPSD) and
144 // then adds contributions from SDigits.
145 // This design avoids scanning over the list of digits to add
146 // contribution to new SDigits only.
148 if( !AliPHOSGetter::GetInstance())
151 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
152 TClonesArray * digits = gime->Digits(GetName()) ;
153 AliPHOSSDigitizer * sdiz = gime->SDigitizer(GetName()) ;
157 const AliPHOSGeometry *geom = gime->PHOSGeometry() ;
159 //Making digits with noise, first EMC
160 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
164 TString name = geom->GetName() ;
166 nCPV = nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
167 geom->GetNModules() ;
169 digits->Expand(nCPV) ;
172 // sdigitize random gaussian noise and add it to all cells (EMCA+CPV+PPSD)
173 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
174 const AliPHOSSDigitizer * sDigitizer = gime->SDigitizer(GetName());
176 cerr << "ERROR: AliPHOSDigitizer::Digitize -> SDigitizer with name " << GetName() << " not found " << endl ;
180 // loop through the sdigits posted to the White Board and add them to the noise
181 TCollection * folderslist = gime->SDigitsFolder()->GetListOfFolders() ;
182 TIter next(folderslist) ;
183 TFolder * folder = 0 ;
184 TClonesArray * sdigits = 0 ;
186 TObjArray * sdigArray = new TObjArray(2) ;
187 while ( (folder = (TFolder*)next()) )
188 if ( (sdigits = (TClonesArray*)folder->FindObject(GetName()) ) ) {
189 cout << "INFO: AliPHOSDigitizer::Digitize -> Adding SDigits "
190 << GetName() << " from " << folder->GetName() << endl ;
191 sdigArray->AddAt(sdigits, input) ;
195 //Find the first crystall with signal
196 Int_t nextSig = 200000 ;
198 for(i=0; i<input; i++){
199 sdigits = (TClonesArray *)sdigArray->At(i) ;
200 if ( !sdigits->GetEntries() )
202 Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
203 if(curNext < nextSig)
207 TArrayI index(input) ;
208 index.Reset() ; //Set all indexes to zero
210 AliPHOSDigit * digit ;
211 AliPHOSDigit * curSDigit ;
213 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
215 //Put Noise contribution
216 for(absID = 1; absID <= nEMC; absID++){
217 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
218 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
219 //look if we have to add signal?
221 //Add SDigits from all inputs
222 digit = (AliPHOSDigit *) digits->At(absID-1) ;
226 Float_t a = digit->GetAmp() ;
227 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
228 //Mark the beginnign of the signal
229 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
230 //Mark the end of the ignal
231 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
234 for(i=0; i<input; i++){
235 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
236 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
239 //May be several digits will contribute from the same input
240 while(curSDigit && curSDigit->GetId() == absID){
241 //Shift primary to separate primaries belonging different inputs
242 Int_t primaryoffset ;
244 primaryoffset = fARD->GetMask(i) ;
246 primaryoffset = 10000000*i ;
247 curSDigit->ShiftPrimary(primaryoffset) ;
249 a = curSDigit->GetAmp() ;
250 b = a /fTimeSignalLength ;
251 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
252 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
254 *digit = *digit + *curSDigit ; //add energies
257 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
258 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
264 //calculate and set time
265 Float_t time = FrontEdgeTime(ticks) ;
266 digit->SetTime(time) ;
268 //Find next signal module
270 for(i=0; i<input; i++){
271 sdigits = ((TClonesArray *)sdigArray->At(i)) ;
272 Int_t curNext = nextSig ;
273 if(sdigits->GetEntriesFast() > index[i] ){
274 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
277 if(curNext < nextSig) nextSig = curNext ;
286 //Now CPV digits (different noise and no timing)
287 for(absID = nEMC+1; absID <= nCPV; absID++){
288 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
289 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
290 //look if we have to add signal?
292 digit = (AliPHOSDigit *) digits->At(absID-1) ;
293 //Add SDigits from all inputs
294 for(i=0; i<input; i++){
295 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
296 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
299 //May be several digits will contribute from the same input
300 while(curSDigit && curSDigit->GetId() == absID){
301 //Shift primary to separate primaries belonging different inputs
302 Int_t primaryoffset ;
304 primaryoffset = fARD->GetMask(i) ;
306 primaryoffset = 10000000*i ;
307 curSDigit->ShiftPrimary(primaryoffset) ;
310 *digit = *digit + *curSDigit ;
312 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
313 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
319 //Find next signal module
320 for(i=0; i<input; i++){
321 sdigits = (TClonesArray *)sdigArray->At(i) ;
322 Int_t curNext = nextSig ;
323 if(sdigits->GetEntriesFast() > index[i] )
324 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
325 if(curNext < nextSig) nextSig = curNext ;
330 delete sdigArray ; //We should not delete its contents
334 //remove digits below thresholds
335 for(absID = 0; absID < nEMC ; absID++){
336 digit = (AliPHOSDigit*) digits->At(absID) ;
337 if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
338 digits->RemoveAt(absID) ;
340 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
343 for(absID = nEMC; absID < nCPV ; absID++)
344 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(absID))->GetAmp()) < fCPVDigitThreshold)
345 digits->RemoveAt(absID) ;
349 Int_t ndigits = digits->GetEntriesFast() ;
350 digits->Expand(ndigits) ;
352 //Set indexes in list of digits and make true digitization of the energy
353 for (i = 0 ; i < ndigits ; i++) {
354 AliPHOSDigit * digit = (AliPHOSDigit *) digits->At(i) ;
355 digit->SetIndexInList(i) ;
356 Float_t energy = sdiz->Calibrate(digit->GetAmp()) ;
357 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
361 //____________________________________________________________________________
362 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
365 if(absId <= fEmcCrystals){ //digitize as EMC
366 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
367 if(chanel > fNADCemc ) chanel = fNADCemc ;
369 else{ //Digitize as CPV
370 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
371 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
375 //____________________________________________________________________________
376 void AliPHOSDigitizer::Exec(Option_t *option)
380 if( strcmp(GetName(), "") == 0 )
383 if (strstr(option,"print")) {
388 if(strstr(option,"tim"))
389 gBenchmark->Start("PHOSDigitizer");
391 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
398 treeD = fARD->GetTreeD() ;
399 nevents = 1 ; // Will process only one event
402 gAlice->GetEvent(0) ;
403 nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
404 treeD=gAlice->TreeD() ;
407 cerr << " AliPHOSDigitizer :: Can not find TreeD " << endl ;
411 //Check, if this branch already exits
412 TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
414 TBranch * branch = 0 ;
415 Bool_t phosfound = kFALSE, digitizerfound = kFALSE ;
417 while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
418 if ( (strcmp(branch->GetName(), "PHOS")==0) &&
419 (strcmp(branch->GetTitle(), GetName())==0) )
422 else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) &&
423 (strcmp(branch->GetTitle(), GetName())==0) )
424 digitizerfound = kTRUE ;
428 cerr << "WARNING: AliPHOSDigitizer -> Digits branch with name " << GetName()
429 << " already exits" << endl ;
432 if ( digitizerfound ) {
433 cerr << "WARNING: AliPHOSDigitizer -> Digitizer branch with name " << GetName()
434 << " already exits" << endl ;
440 for(ievent = 0; ievent < nevents; ievent++){
444 for(input = 0 ; input < fARD->GetNinputs(); input ++){
445 TTree * treeS = fARD->GetInputTreeS(input) ;
447 cerr << "AliPHOSDigitizer -> No Input " << endl ;
450 gime->ReadTreeS(treeS,input) ;
454 gime->Event(ievent,"S") ;
456 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
458 WriteDigits(ievent) ;
460 if(strstr(option,"deb"))
463 //increment the total number of Digits per run
464 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
467 if(strstr(option,"tim")){
468 gBenchmark->Stop("PHOSDigitizer");
469 cout << "AliPHOSDigitizer:" << endl ;
470 cout << " took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for Digitizing "
471 << gBenchmark->GetCpuTime("PHOSDigitizer")/nevents << " seconds per event " << endl ;
477 //____________________________________________________________________________
478 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks)
480 ticks->Sort() ; //Sort in accordance with times of ticks
482 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
483 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
486 while((t=(AliPHOSTick*) it.Next())){
487 if(t->GetTime() < time) //This tick starts before crossing
492 time = ctick->CrossingTime(fTimeThreshold) ;
496 //____________________________________________________________________________
497 Bool_t AliPHOSDigitizer::Init()
500 fEMCDigitThreshold = 0.01 ;
502 fCPVDigitThreshold = 0.09 ;
503 fTimeResolution = 0.5e-9 ;
504 fTimeSignalLength = 1.0e-9 ;
506 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
507 fADCpedestalEmc = 0.005 ; //
508 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
510 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
511 fADCpedestalCpv = 0.012 ; //
512 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
514 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
516 // Makes all memory allocations
517 // Adds Digitizer task to the folder of PHOS tasks
519 //============================================================= YS
520 // The initialisation is now done by AliPHOSGetter
522 if( strcmp(GetTitle(), "") == 0 )
523 SetTitle("galice.root") ;
525 // the SDigits name is stored by AliPHOSGetter as the name of the TClones Array
526 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/branchname and has branchTitle as title.
528 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
530 cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
534 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
535 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
537 // create a folder on the white board //YSAlice/WhiteBoard/Digits/PHOS/headerFile/digitsTitle
538 gime->PostDigits(GetName() ) ;
540 //add Task to //YSAlice/tasks/Digitizer/PHOS
541 gime->PostDigitizer(this) ;
543 //Mark that we will use current header file
545 gime->PostSDigits(GetName(),GetTitle()) ;
546 gime->PostSDigitizer(GetName(),GetTitle()) ;
551 //__________________________________________________________________
552 void AliPHOSDigitizer::MixWith(const char* headerFile)
554 // Allows to produce digits by superimposing background and signal event.
555 // It is assumed, that headers file with SIGNAL events is opened in
557 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
558 // Thus we avoid writing (changing) huge and expensive
559 // backgound files: all output will be writen into SIGNAL, i.e.
560 // opened in constructor file.
562 // One can open as many files to mix with as one needs.
563 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
566 if( strcmp(GetName(), "") == 0 )
570 cout << "Can not use this method under AliRunDigitizer " << endl ;
574 // check if the specified SDigits do not already exist on the White Board:
575 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/sDigitsTitle
577 TString path = "YSAlice/WhiteBoard/SDigits/PHOS/" ;
581 if ( gROOT->FindObjectAny(path.Data()) ) {
582 cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
586 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
587 gime->PostSDigits(GetName(),headerFile) ;
589 // check if the requested file is already open or exist and if SDigits Branch exist
590 TFile * file = (TFile*)gROOT->FindObject(headerFile);
592 file = new TFile(headerFile, "READ") ;
594 cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ;
601 //__________________________________________________________________
602 void AliPHOSDigitizer::Print(Option_t* option)const {
603 // Print Digitizer's parameters
604 if( strcmp(GetName(), "") != 0 ){
606 cout << "------------------- "<< GetName() << " -------------" << endl ;
607 cout << "Digitizing sDigits from file(s): " <<endl ;
609 TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("YSAlice/WhiteBoard/SDigits/PHOS"))->GetListOfFolders() ;
610 TIter next(folderslist) ;
611 TFolder * folder = 0 ;
613 while ( (folder = (TFolder*)next()) ) {
614 if ( folder->FindObject(GetName()) )
615 cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ;
618 cout << "Writing digits to " << GetTitle() << endl ;
621 cout << "With following parameters: " << endl ;
622 cout << " Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
623 cout << " Threshold in EMC (fEMCDigitThreshold) = " << fEMCDigitThreshold << endl ; ;
624 cout << " Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ;
625 cout << " Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ;
626 cout << "---------------------------------------------------" << endl ;
629 cout << "AliPHOSDigitizer not initialized " << endl ;
633 //__________________________________________________________________
634 void AliPHOSDigitizer::PrintDigits(Option_t * option){
635 // Print a table of digits
637 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
638 TClonesArray * digits = gime->Digits() ;
640 cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
641 cout << " Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
643 if(strstr(option,"all")||strstr(option,"EMC")){
646 AliPHOSDigit * digit;
647 cout << "Digit Id Amplitude Index " << " Nprim " << " Primaries list " << endl;
648 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
650 for (index = 0 ; (index < digits->GetEntries()) &&
651 (((AliPHOSDigit * ) digits->At(index))->GetId() <= maxEmc) ; index++) {
652 digit = (AliPHOSDigit * ) digits->At(index) ;
653 if(digit->GetNprimary() == 0) continue;
654 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
655 << setw(6) << digit->GetIndexInList() << " "
656 << setw(5) << digit->GetNprimary() <<" ";
659 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
660 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
666 if(strstr(option,"all")||strstr(option,"CPV")){
668 //loop over CPV digits
669 AliPHOSDigit * digit;
670 cout << "Digit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
671 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
673 for (index = 0 ; index < digits->GetEntries(); index++) {
674 digit = (AliPHOSDigit * ) digits->At(index) ;
675 if(digit->GetId() > maxEmc){
676 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
677 << setw(6) << digit->GetIndexInList() << " "
678 << setw(5) << digit->GetNprimary() <<" ";
681 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
682 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
690 //__________________________________________________________________
691 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
693 // we set title (comment) of the SDigits branch in the first! header file
694 if( strcmp(GetName(), "") == 0 )
697 AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ;
700 //__________________________________________________________________
701 Float_t AliPHOSDigitizer::TimeOfNoise(void)
702 { // Calculates the time signal generated by noise
703 //to be rewritten, now returns just big number
707 //____________________________________________________________________________
708 void AliPHOSDigitizer::Reset()
710 // sets current event number to the first simulated event
712 if( strcmp(GetName(), "") == 0 )
716 // for(inputs = 0; inputs < fNinputs ;inputs++)
717 // fIevent->AddAt(-1, inputs ) ;
721 //____________________________________________________________________________
722 void AliPHOSDigitizer::WriteDigits(Int_t event)
725 // Makes TreeD in the output file.
726 // Check if branch already exists:
727 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
728 // already existing branches.
729 // else creates branch with Digits, named "PHOS", title "...",
730 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
731 // and names of files, from which digits are made.
733 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
734 TClonesArray * digits = gime->Digits(GetName()) ;
739 treeD = fARD->GetTreeD() ;
741 treeD = gAlice->TreeD();
743 // create new branches
744 // -- generate file name if necessary
746 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
747 file = new char[strlen(gAlice->GetBaseFile())+20] ;
748 sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
751 TDirectory *cwd = gDirectory;
752 // -- create Digits branch
753 Int_t bufferSize = 32000 ;
754 TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
755 digitsBranch->SetTitle(GetName());
757 digitsBranch->SetFile(file);
758 TIter next( digitsBranch->GetListOfBranches());
760 while ((sbr=(TBranch*)next())) {
766 // -- Create Digitizer branch
767 Int_t splitlevel = 0 ;
768 AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
769 TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel);
770 digitizerBranch->SetTitle(GetName());
772 digitizerBranch->SetFile(file);
773 TIter next( digitizerBranch->GetListOfBranches());
775 while ((sbr=(TBranch*)next())) {
781 digitsBranch->Fill() ;
782 digitizerBranch->Fill() ;
784 treeD->Write(0,kOverwrite) ;