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 = 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 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ATTENTION !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
202 // if there are no SDigits in this event ===> Dmitri should check
203 if ( !sdigits->GetEntries() )
205 Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
206 if(curNext < nextSig)
210 cout << " AliPHOSDigitizer::Digitize 2" << endl ;
211 TArrayI index(input) ;
212 index.Reset() ; //Set all indexes to zero
214 AliPHOSDigit * digit ;
215 AliPHOSDigit * curSDigit ;
217 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
219 //Put Noise contribution
220 for(absID = 1; absID <= nEMC; absID++){
221 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
222 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
223 //look if we have to add signal?
225 //Add SDigits from all inputs
226 digit = (AliPHOSDigit *) digits->At(absID-1) ;
230 Float_t a = digit->GetAmp() ;
231 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
232 //Mark the beginnign of the signal
233 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
234 //Mark the end of the ignal
235 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
238 for(i=0; i<input; i++){
239 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
240 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
243 //May be several digits will contribute from the same input
244 while(curSDigit && curSDigit->GetId() == absID){
245 //Shift primary to separate primaries belonging different inputs
246 Int_t primaryoffset ;
248 primaryoffset = fARD->GetMask(i) ;
250 primaryoffset = 10000000*i ;
251 curSDigit->ShiftPrimary(primaryoffset) ;
253 a = curSDigit->GetAmp() ;
254 b = a /fTimeSignalLength ;
255 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
256 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
258 *digit = *digit + *curSDigit ; //add energies
261 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
262 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
268 //calculate and set time
269 Float_t time = FrontEdgeTime(ticks) ;
270 digit->SetTime(time) ;
272 //Find next signal module
274 for(i=0; i<input; i++){
275 sdigits = ((TClonesArray *)sdigArray->At(i)) ;
276 Int_t curNext = nextSig ;
277 if(sdigits->GetEntriesFast() > index[i] ){
278 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
281 if(curNext < nextSig) nextSig = curNext ;
290 //Now CPV digits (different noise and no timing)
291 for(absID = nEMC+1; absID <= nCPV; absID++){
292 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
293 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
294 //look if we have to add signal?
296 digit = (AliPHOSDigit *) digits->At(absID-1) ;
297 //Add SDigits from all inputs
298 for(i=0; i<input; i++){
299 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
300 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
303 //May be several digits will contribute from the same input
304 while(curSDigit && curSDigit->GetId() == absID){
305 //Shift primary to separate primaries belonging different inputs
306 Int_t primaryoffset ;
308 primaryoffset = fARD->GetMask(i) ;
310 primaryoffset = 10000000*i ;
311 curSDigit->ShiftPrimary(primaryoffset) ;
314 *digit = *digit + *curSDigit ;
316 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
317 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
323 //Find next signal module
324 for(i=0; i<input; i++){
325 sdigits = (TClonesArray *)sdigArray->At(i) ;
326 Int_t curNext = nextSig ;
327 if(sdigits->GetEntriesFast() > index[i] )
328 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
329 if(curNext < nextSig) nextSig = curNext ;
334 delete sdigArray ; //We should not delete its contents
338 //remove digits below thresholds
339 for(absID = 0; absID < nEMC ; absID++)
340 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(absID))->GetAmp()) < fEMCDigitThreshold)
341 digits->RemoveAt(absID) ;
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"))
464 if(strstr(option,"tim")){
465 gBenchmark->Stop("PHOSDigitizer");
466 cout << "AliPHOSDigitizer:" << endl ;
467 cout << " took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for Digitizing "
468 << gBenchmark->GetCpuTime("PHOSDigitizer")/nevents << " seconds per event " << endl ;
474 //____________________________________________________________________________
475 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks)
477 ticks->Sort() ; //Sort in accordance with times of ticks
479 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
480 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
483 while((t=(AliPHOSTick*) it.Next())){
484 if(t->GetTime() < time) //This tick starts before crossing
489 time = ctick->CrossingTime(fTimeThreshold) ;
493 //____________________________________________________________________________
494 Bool_t AliPHOSDigitizer::Init()
497 fEMCDigitThreshold = 0.01 ;
499 fCPVDigitThreshold = 0.09 ;
500 fTimeResolution = 1.0e-9 ;
501 fTimeSignalLength = 1.0e-9 ;
503 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
504 fADCpedestalEmc = 0.005 ; //
505 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
507 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
508 fADCpedestalCpv = 0.012 ; //
509 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
511 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
513 // Makes all memory allocations
514 // Adds Digitizer task to the folder of PHOS tasks
516 //============================================================= YS
517 // The initialisation is now done by AliPHOSGetter
519 if( strcmp(GetTitle(), "") == 0 )
520 SetTitle("galice.root") ;
522 // the SDigits name is stored by AliPHOSGetter as the name of the TClones Array
523 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/branchname and has branchTitle as title.
525 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
527 cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
531 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
532 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
534 // create a folder on the white board //YSAlice/WhiteBoard/Digits/PHOS/headerFile/digitsTitle
535 gime->PostDigits(GetName() ) ;
537 //add Task to //YSAlice/tasks/Digitizer/PHOS
538 gime->PostDigitizer(this) ;
540 //Mark that we will use current header file
542 gime->PostSDigits(GetName(),GetTitle()) ;
543 gime->PostSDigitizer(GetName(),GetTitle()) ;
548 //__________________________________________________________________
549 void AliPHOSDigitizer::MixWith(const char* headerFile)
551 // Allows to produce digits by superimposing background and signal event.
552 // It is assumed, that headers file with SIGNAL events is opened in
554 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
555 // Thus we avoid writing (changing) huge and expensive
556 // backgound files: all output will be writen into SIGNAL, i.e.
557 // opened in constructor file.
559 // One can open as many files to mix with as one needs.
560 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
563 if( strcmp(GetName(), "") == 0 )
567 cout << "Can not use this method under AliRunDigitizer " << endl ;
571 // check if the specified SDigits do not already exist on the White Board:
572 // //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/sDigitsTitle
574 TString path = "YSAlice/WhiteBoard/SDigits/PHOS/" ;
578 if ( gROOT->FindObjectAny(path.Data()) ) {
579 cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
583 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
584 gime->PostSDigits(GetName(),headerFile) ;
586 // check if the requested file is already open or exist and if SDigits Branch exist
587 TFile * file = (TFile*)gROOT->FindObject(headerFile);
589 file = new TFile(headerFile, "READ") ;
591 cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ;
598 //__________________________________________________________________
599 void AliPHOSDigitizer::Print(Option_t* option)const {
600 // Print Digitizer's parameters
601 if( strcmp(GetName(), "") != 0 ){
603 cout << "------------------- "<< GetName() << " -------------" << endl ;
604 cout << "Digitizing sDigits from file(s): " <<endl ;
606 TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("YSAlice/WhiteBoard/SDigits/PHOS"))->GetListOfFolders() ;
607 TIter next(folderslist) ;
608 TFolder * folder = 0 ;
610 while ( (folder = (TFolder*)next()) ) {
611 if ( folder->FindObject(GetName()) )
612 cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ;
615 cout << "Writing digits to " << GetTitle() << endl ;
618 cout << "With following parameters: " << endl ;
619 cout << " Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
620 cout << " Threshold in EMC (fEMCDigitThreshold) = " << fEMCDigitThreshold << endl ; ;
621 cout << " Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ;
622 cout << " Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ;
623 cout << "---------------------------------------------------" << endl ;
626 cout << "AliPHOSDigitizer not initialized " << endl ;
630 //__________________________________________________________________
631 void AliPHOSDigitizer::PrintDigits(Option_t * option){
632 // Print a table of digits
634 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
635 TClonesArray * digits = gime->Digits() ;
637 cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
638 cout << " Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
640 if(strstr(option,"all")||strstr(option,"EMC")){
643 AliPHOSDigit * digit;
644 cout << "Digit Id Amplitude Index " << " Nprim " << " Primaries list " << endl;
645 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
647 for (index = 0 ; (index < digits->GetEntries()) &&
648 (((AliPHOSDigit * ) digits->At(index))->GetId() <= maxEmc) ; index++) {
649 digit = (AliPHOSDigit * ) digits->At(index) ;
650 if(digit->GetNprimary() == 0) continue;
651 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
652 << setw(6) << digit->GetIndexInList() << " "
653 << setw(5) << digit->GetNprimary() <<" ";
656 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
657 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
663 if(strstr(option,"all")||strstr(option,"CPV")){
665 //loop over CPV digits
666 AliPHOSDigit * digit;
667 cout << "Digit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
668 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
670 for (index = 0 ; index < digits->GetEntries(); index++) {
671 digit = (AliPHOSDigit * ) digits->At(index) ;
672 if(digit->GetId() > maxEmc){
673 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
674 << setw(6) << digit->GetIndexInList() << " "
675 << setw(5) << digit->GetNprimary() <<" ";
678 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
679 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
687 //__________________________________________________________________
688 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
690 // we set title (comment) of the SDigits branch in the first! header file
691 if( strcmp(GetName(), "") == 0 )
694 AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ;
697 //__________________________________________________________________
698 Float_t AliPHOSDigitizer::TimeOfNoise(void)
699 { // Calculates the time signal generated by noise
700 //to be rewritten, now returns just big number
704 //____________________________________________________________________________
705 void AliPHOSDigitizer::Reset()
707 // sets current event number to the first simulated event
709 if( strcmp(GetName(), "") == 0 )
713 // for(inputs = 0; inputs < fNinputs ;inputs++)
714 // fIevent->AddAt(-1, inputs ) ;
718 //____________________________________________________________________________
719 void AliPHOSDigitizer::WriteDigits(Int_t event)
722 // Makes TreeD in the output file.
723 // Check if branch already exists:
724 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
725 // already existing branches.
726 // else creates branch with Digits, named "PHOS", title "...",
727 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
728 // and names of files, from which digits are made.
730 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
731 TClonesArray * digits = gime->Digits(GetName()) ;
736 treeD = fARD->GetTreeD() ;
738 treeD = gAlice->TreeD();
740 // create new branches
741 // -- generate file name if necessary
743 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
744 file = new char[strlen(gAlice->GetBaseFile())+20] ;
745 sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
748 TDirectory *cwd = gDirectory;
749 // -- create Digits branch
750 Int_t bufferSize = 32000 ;
751 TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
752 digitsBranch->SetTitle(GetName());
754 digitsBranch->SetFile(file);
755 TIter next( digitsBranch->GetListOfBranches());
757 while ((sbr=(TBranch*)next())) {
763 // -- Create Digitizer branch
764 Int_t splitlevel = 0 ;
765 AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
766 TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel);
767 digitizerBranch->SetTitle(GetName());
769 digitizerBranch->SetFile(file);
770 TIter next( digitizerBranch->GetListOfBranches());
772 while ((sbr=(TBranch*)next())) {
778 digitsBranch->Fill() ;
779 digitizerBranch->Fill() ;
781 treeD->Write(0,kOverwrite) ;