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 /* History of cvs commits:
21 * Revision 1.85 2005/05/28 14:19:04 schutz
22 * Compilation warnings fixed by T.P.
26 //_________________________________________________________________________
27 //*-- Author : Dmitri Peressounko (SUBATECH & Kurchatov Institute)
28 //////////////////////////////////////////////////////////////////////////////
29 // This TTask performs digitization of Summable digits (in the PHOS case it is just
30 // the sum of contributions from all primary particles into a given cell).
31 // In addition it performs mixing of summable digits from different events.
32 // The name of the TTask is also the title of the branch that will contain
33 // the created SDigits
34 // The title of the TTAsk is the name of the file that contains the hits from
35 // which the SDigits are created
37 // For each event two branches are created in TreeD:
38 // "PHOS" - list of digits
39 // "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
41 // Note, that one can set a title for new digits branch, and repeat digitization with
42 // another set of parameters.
45 // root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
46 // root[1] d->ExecuteTask()
47 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
48 // //Digitizes SDigitis in all events found in file galice.root
50 // root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;
51 // // Will read sdigits from galice1.root
52 // root[3] d1->MixWith("galice2.root")
53 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
54 // // Reads another set of sdigits from galice2.root
55 // root[3] d1->MixWith("galice3.root")
56 // // Reads another set of sdigits from galice3.root
57 // root[4] d->ExecuteTask("deb timing")
58 // // Reads SDigits from files galice1.root, galice2.root ....
59 // // mixes them and stores produced Digits in file galice1.root
60 // // deb - prints number of produced digits
61 // // deb all - prints list of produced digits
62 // // timing - prints time used for digitization
65 // --- ROOT system ---
68 #include "TBenchmark.h"
71 // --- Standard library ---
73 // --- AliRoot header files ---
75 #include "AliRunDigitizer.h"
76 #include "AliPHOSDigit.h"
77 #include "AliPHOSGetter.h"
78 #include "AliPHOSDigitizer.h"
79 #include "AliPHOSSDigitizer.h"
80 #include "AliPHOSGeometry.h"
81 #include "AliPHOSTick.h"
83 ClassImp(AliPHOSDigitizer)
86 //____________________________________________________________________________
87 AliPHOSDigitizer::AliPHOSDigitizer():AliDigitizer("",""),
94 fDefaultInit = kTRUE ;
95 fManager = 0 ; // We work in the standalong mode
96 fEventFolderName = "" ;
99 //____________________________________________________________________________
100 AliPHOSDigitizer::AliPHOSDigitizer(TString alirunFileName,
101 TString eventFolderName):
102 AliDigitizer("PHOS"+AliConfig::Instance()->GetDigitizerTaskName(),
104 fInputFileNames(0), fEventNames(0), fEventFolderName(eventFolderName)
109 fDefaultInit = kFALSE ;
110 fManager = 0 ; // We work in the standalong mode
113 //____________________________________________________________________________
114 AliPHOSDigitizer::AliPHOSDigitizer(const AliPHOSDigitizer & d)
119 SetName(d.GetName()) ;
120 SetTitle(d.GetTitle()) ;
121 fPinNoise = d.fPinNoise ;
122 fEMCDigitThreshold = d.fEMCDigitThreshold ;
123 fCPVNoise = d.fCPVNoise ;
124 fCPVDigitThreshold = d.fCPVDigitThreshold ;
125 fTimeResolution = d.fTimeResolution ;
126 fTimeThreshold = d.fTimeThreshold ;
127 fTimeSignalLength = d.fTimeSignalLength ;
128 fADCchanelEmc = d.fADCchanelEmc ;
129 fADCpedestalEmc = d.fADCpedestalEmc ;
130 fNADCemc = d.fNADCemc ;
131 fADCchanelCpv = d.fADCchanelCpv ;
132 fADCpedestalCpv = d.fADCpedestalCpv ;
133 fNADCcpv = d.fNADCcpv ;
134 fEventFolderName = d.fEventFolderName;
137 //____________________________________________________________________________
138 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * rd):
139 AliDigitizer(rd,"PHOS"+AliConfig::Instance()->GetDigitizerTaskName()),
142 // ctor Init() is called by RunDigitizer
144 fEventFolderName = fManager->GetInputFolderName(0) ;
145 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
147 fDefaultInit = kFALSE ;
150 //____________________________________________________________________________
151 AliPHOSDigitizer::~AliPHOSDigitizer()
153 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
155 // Clean Digitizer from the white board
156 gime->PhosLoader()->CleanDigitizer() ;
158 delete [] fInputFileNames ;
159 delete [] fEventNames ;
163 //____________________________________________________________________________
164 void AliPHOSDigitizer::Digitize(Int_t event)
167 // Makes the digitization of the collected summable digits.
168 // It first creates the array of all PHOS modules
169 // filled with noise (different for EMC, and CPV) and
170 // then adds contributions from SDigits.
171 // This design avoids scanning over the list of digits to add
172 // contribution to new SDigits only.
174 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
175 Int_t ReadEvent = event ;
177 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
178 AliInfo(Form("Adding event %d from input stream 0 %s %s",
179 ReadEvent, GetTitle(), fEventFolderName.Data())) ;
180 gime->Event(ReadEvent, "S") ;
181 TClonesArray * digits = gime->Digits() ;
184 const AliPHOSGeometry *geom = gime->PHOSGeometry() ;
185 //Making digits with noise, first EMC
186 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
191 nCPV = nEMC + geom->GetNumberOfCPVPadsZ() * geom->GetNumberOfCPVPadsPhi() * geom->GetNModules() ;
193 digits->Expand(nCPV) ;
195 // get first the sdigitizer from the tasks list
196 if ( !gime->SDigitizer() )
197 gime->LoadSDigitizer();
198 AliPHOSSDigitizer * sDigitizer = gime->SDigitizer();
201 AliFatal(Form("SDigitizer with name %s %s not found",
202 GetTitle(), fEventFolderName.Data() )) ;
204 //take all the inputs to add together and load the SDigits
205 TObjArray * sdigArray = new TObjArray(fInput) ;
206 sdigArray->AddAt(gime->SDigits(), 0) ;
208 for(i = 1 ; i < fInput ; i++){
209 TString tempo(fEventNames[i]) ;
211 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[i], tempo) ;
213 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
214 AliInfo(Form("Adding event %d from input stream %d %s %s",
215 ReadEvent, i, fInputFileNames[i].Data(), tempo.Data())) ;
216 gime->Event(ReadEvent,"S");
217 sdigArray->AddAt(gime->SDigits(), i) ;
220 //Find the first crystall with signal
221 Int_t nextSig = 200000 ;
222 TClonesArray * sdigits ;
223 for(i = 0 ; i < fInput ; i++){
224 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
225 if ( !sdigits->GetEntriesFast() )
227 Int_t curNext = dynamic_cast<AliPHOSDigit *>(sdigits->At(0))->GetId() ;
228 if(curNext < nextSig)
232 TArrayI index(fInput) ;
233 index.Reset() ; //Set all indexes to zero
235 AliPHOSDigit * digit ;
236 AliPHOSDigit * curSDigit ;
238 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
240 //Put Noise contribution
241 for(absID = 1 ; absID <= nEMC ; absID++){
242 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
243 new((*digits)[absID-1]) AliPHOSDigit( -1, absID, sDigitizer->Digitize(noise), TimeOfNoise() ) ;
244 //look if we have to add signal?
245 digit = dynamic_cast<AliPHOSDigit *>(digits->At(absID-1)) ;
248 //Add SDigits from all inputs
251 Float_t a = digit->GetAmp() ;
252 Float_t b = TMath::Abs( a / fTimeSignalLength) ;
253 //Mark the beginning of the signal
254 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
255 //Mark the end of the ignal
256 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
259 for(i = 0 ; i < fInput ; i++){
260 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
261 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
264 //May be several digits will contribute from the same input
265 while(curSDigit && curSDigit->GetId() == absID){
266 //Shift primary to separate primaries belonging different inputs
267 Int_t primaryoffset ;
269 primaryoffset = fManager->GetMask(i) ;
271 primaryoffset = 10000000*i ;
272 curSDigit->ShiftPrimary(primaryoffset) ;
274 a = curSDigit->GetAmp() ;
275 b = a /fTimeSignalLength ;
276 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
277 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
279 *digit = *digit + *curSDigit ; //add energies
282 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
283 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
289 //calculate and set time
290 Float_t time = FrontEdgeTime(ticks) ;
291 digit->SetTime(time) ;
293 //Find next signal module
295 for(i = 0 ; i < fInput ; i++){
296 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
297 Int_t curNext = nextSig ;
298 if(sdigits->GetEntriesFast() > index[i] ){
299 curNext = dynamic_cast<AliPHOSDigit *>(sdigits->At(index[i]))->GetId() ;
301 if(curNext < nextSig) nextSig = curNext ;
309 //Now CPV digits (different noise and no timing)
310 for(absID = nEMC+1; absID <= nCPV; absID++){
311 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
312 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
313 //look if we have to add signal?
315 digit = dynamic_cast<AliPHOSDigit *>(digits->At(absID-1)) ;
316 //Add SDigits from all inputs
317 for(i = 0 ; i < fInput ; i++){
318 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
319 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
323 //May be several digits will contribute from the same input
324 while(curSDigit && curSDigit->GetId() == absID){
325 //Shift primary to separate primaries belonging different inputs
326 Int_t primaryoffset ;
328 primaryoffset = fManager->GetMask(i) ;
330 primaryoffset = 10000000*i ;
331 curSDigit->ShiftPrimary(primaryoffset) ;
334 *digit = *digit + *curSDigit ;
336 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
337 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i]) ) ;
343 //Find next signal module
345 for(i = 0 ; i < fInput ; i++){
346 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
347 Int_t curNext = nextSig ;
348 if(sdigits->GetEntriesFast() > index[i] )
349 curNext = dynamic_cast<AliPHOSDigit *>( sdigits->At(index[i]) )->GetId() ;
350 if(curNext < nextSig) nextSig = curNext ;
356 delete sdigArray ; //We should not delete its contents
358 //remove digits below thresholds
359 for(i = 0 ; i < nEMC ; i++){
360 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
361 if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
362 digits->RemoveAt(i) ;
364 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
368 for(i = nEMC; i < nCPV ; i++)
369 if( sDigitizer->Calibrate( dynamic_cast<AliPHOSDigit*>(digits->At(i))->GetAmp() ) < fCPVDigitThreshold )
370 digits->RemoveAt(i) ;
374 Int_t ndigits = digits->GetEntriesFast() ;
375 digits->Expand(ndigits) ;
377 //Set indexes in list of digits and make true digitization of the energy
378 for (i = 0 ; i < ndigits ; i++) {
379 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
380 digit->SetIndexInList(i) ;
381 Float_t energy = sDigitizer->Calibrate(digit->GetAmp()) ;
382 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
386 //____________________________________________________________________________
387 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
389 // Returns digitized value of the energy in a cell absId
392 if(absId <= fEmcCrystals){ //digitize as EMC
393 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
394 if(chanel > fNADCemc ) chanel = fNADCemc ;
396 else{ //Digitize as CPV
397 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
398 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
403 //____________________________________________________________________________
404 void AliPHOSDigitizer::Exec(Option_t *option)
406 // Steering method to process digitization for events
407 // in the range from fFirstEvent to fLastEvent.
408 // This range is optionally set by SetEventRange().
409 // if fLastEvent=-1, then process events until the end.
410 // by default fLastEvent = fFirstEvent (process only one event)
412 if (!fInit) { // to prevent overwrite existing file
413 AliError(Form("Give a version name different from %s",
414 fEventFolderName.Data() )) ;
418 if (strstr(option,"print")) {
423 if(strstr(option,"tim"))
424 gBenchmark->Start("PHOSDigitizer");
426 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
428 // Post Digitizer to the white board
429 gime->PostDigitizer(this) ;
431 if (fLastEvent == -1)
432 fLastEvent = gime->MaxEvent() - 1 ;
434 fLastEvent = fFirstEvent ;
436 Int_t nEvents = fLastEvent - fFirstEvent + 1;
440 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
442 gime->Event(ievent,"S") ;
444 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
448 if(strstr(option,"deb"))
451 //increment the total number of Digits per run
452 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
455 gime->PhosLoader()->CleanDigitizer();
457 if(strstr(option,"tim")){
458 gBenchmark->Stop("PHOSDigitizer");
460 message = " took %f seconds for Digitizing %f seconds per event\n" ;
461 AliInfo(Form( message.Data(),
462 gBenchmark->GetCpuTime("PHOSDigitizer"),
463 gBenchmark->GetCpuTime("PHOSDigitizer")/nEvents ));
467 //____________________________________________________________________________
468 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks) const
470 // Returns the shortest time among all time ticks
472 ticks->Sort() ; //Sort in accordance with times of ticks
474 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
475 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
478 while((t=(AliPHOSTick*) it.Next())){
479 if(t->GetTime() < time) //This tick starts before crossing
484 time = ctick->CrossingTime(fTimeThreshold) ;
489 //____________________________________________________________________________
490 Bool_t AliPHOSDigitizer::Init()
492 // Makes all memory allocations
494 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
496 AliFatal(Form("Could not obtain the Getter object for file %s and event %s !",
497 GetTitle(), fEventFolderName.Data()));
501 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
503 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
505 TString opt("Digits") ;
506 if(gime->VersionExists(opt) ) {
507 AliError(Form("Give a version name different from %s",
508 fEventFolderName.Data() )) ;
513 fLastEvent = fFirstEvent ;
515 fInput = fManager->GetNinputs() ;
519 fInputFileNames = new TString[fInput] ;
520 fEventNames = new TString[fInput] ;
521 fInputFileNames[0] = GetTitle() ;
522 fEventNames[0] = fEventFolderName.Data() ;
524 for (index = 1 ; index < fInput ; index++) {
525 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
526 TString tempo = fManager->GetInputFolderName(index) ;
527 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added by fManager
530 //to prevent cleaning of this object while GetEvent is called
531 gime->PhosLoader()->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
536 //____________________________________________________________________________
537 void AliPHOSDigitizer::InitParameters()
539 // Set initial parameters Digitizer
542 fEMCDigitThreshold = 0.012 ;
544 fCPVDigitThreshold = 0.09 ;
545 fTimeResolution = 0.5e-9 ;
546 fTimeSignalLength = 1.0e-9 ;
548 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
549 fADCpedestalEmc = 0.005 ; //
550 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
552 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
553 fADCpedestalCpv = 0.012 ; //
554 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
556 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
557 SetEventRange(0,-1) ;
561 //__________________________________________________________________
562 void AliPHOSDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
564 // Allows to produce digits by superimposing background and signal event.
565 // It is assumed, that headers file with SIGNAL events is opened in
567 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
568 // Thus we avoid writing (changing) huge and expensive
569 // backgound files: all output will be writen into SIGNAL, i.e.
570 // opened in constructor file.
572 // One can open as many files to mix with as one needs.
573 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
576 if( strcmp(fEventFolderName, "") == 0 )
580 Warning("MixWith", "Cannot use this method with AliRunDigitizer\n" ) ;
583 // looking for file which contains AliRun
584 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
585 AliError(Form("File %s does not exist!", alirunFileName.Data())) ;
588 // looking for the file which contains SDigits
589 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
590 TString fileName( gime->GetSDigitsFileName() ) ;
591 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
592 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
593 if ( (gSystem->AccessPathName(fileName)) ) {
594 AliError(Form("The file %s does not exist!", fileName.Data())) ;
597 // need to increase the arrays
598 TString tempo = fInputFileNames[fInput-1] ;
599 delete [] fInputFileNames ;
600 fInputFileNames = new TString[fInput+1] ;
601 fInputFileNames[fInput-1] = tempo ;
603 tempo = fEventNames[fInput-1] ;
604 delete [] fEventNames ;
605 fEventNames = new TString[fInput+1] ;
606 fEventNames[fInput-1] = tempo ;
608 fInputFileNames[fInput] = alirunFileName ;
609 fEventNames[fInput] = eventFolderName ;
613 //__________________________________________________________________
614 void AliPHOSDigitizer::Print(const Option_t *)const
616 // Print Digitizer's parameters
617 AliInfo(Form("\n------------------- %s -------------", GetName() )) ;
618 if( strcmp(fEventFolderName.Data(), "") != 0 ){
619 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
623 nStreams = GetNInputStreams() ;
628 for (index = 0 ; index < nStreams ; index++) {
629 TString tempo(fEventNames[index]) ;
631 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[index], tempo) ;
632 TString fileName( gime->GetSDigitsFileName() ) ;
633 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
634 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
635 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
637 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
638 printf("\nWriting digits to %s", gime->GetDigitsFileName().Data()) ;
640 printf("\nWith following parameters:\n") ;
641 printf(" Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise ) ;
642 printf(" Threshold in EMC (fEMCDigitThreshold) = %f\n", fEMCDigitThreshold ) ;
643 printf(" Noise in CPV (fCPVNoise) = %f\n", fCPVNoise ) ;
644 printf(" Threshold in CPV (fCPVDigitThreshold) = %f\n",fCPVDigitThreshold ) ;
645 printf(" ---------------------------------------------------\n") ;
648 AliInfo(Form("AliPHOSDigitizer not initialized" )) ;
652 //__________________________________________________________________
653 void AliPHOSDigitizer::PrintDigits(Option_t * option)
655 // Print a table of digits
657 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
658 TClonesArray * digits = gime->Digits() ;
660 AliInfo(Form("%d", digits->GetEntriesFast())) ;
661 printf("\nevent %d", gAlice->GetEvNumber()) ;
662 printf("\n Number of entries in Digits list %d", digits->GetEntriesFast() ) ;
665 if(strstr(option,"all")||strstr(option,"EMC")){
667 AliPHOSDigit * digit;
668 printf("\nEMC digits (with primaries):\n") ;
669 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
670 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
672 for (index = 0 ; (index < digits->GetEntriesFast()) &&
673 (dynamic_cast<AliPHOSDigit *>(digits->At(index))->GetId() <= maxEmc) ; index++) {
674 digit = (AliPHOSDigit * ) digits->At(index) ;
675 if(digit->GetNprimary() == 0)
677 printf("%6d %8d %6.5e %4d %2d :",
678 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
680 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
681 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
687 if(strstr(option,"all")||strstr(option,"CPV")){
689 //loop over CPV digits
690 AliPHOSDigit * digit;
691 printf("\nCPV digits (with primaries):\n") ;
692 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
693 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
695 for (index = 0 ; index < digits->GetEntriesFast(); index++) {
696 digit = (AliPHOSDigit * ) digits->At(index) ;
697 if(digit->GetId() > maxEmc){
698 printf("%6d %8d %4d %2d :",
699 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
701 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
702 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
711 //__________________________________________________________________
712 Float_t AliPHOSDigitizer::TimeOfNoise(void) const
713 { // Calculates the time signal generated by noise
714 //PH Info("TimeOfNoise", "Change me") ;
715 return gRandom->Rndm() * 1.28E-5;
718 //__________________________________________________________________
719 void AliPHOSDigitizer::Unload()
723 for(i = 1 ; i < fInput ; i++){
724 TString tempo(fEventNames[i]) ;
726 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[i], tempo) ;
727 gime->PhosLoader()->UnloadSDigits() ;
730 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
731 gime->PhosLoader()->UnloadDigits() ;
734 //____________________________________________________________________________
735 void AliPHOSDigitizer::WriteDigits()
738 // Makes TreeD in the output file.
739 // Check if branch already exists:
740 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
741 // already existing branches.
742 // else creates branch with Digits, named "PHOS", title "...",
743 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
744 // and names of files, from which digits are made.
746 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
747 const TClonesArray * digits = gime->Digits() ;
748 TTree * treeD = gime->TreeD();
750 // -- create Digits branch
751 Int_t bufferSize = 32000 ;
752 TBranch * digitsBranch = treeD->Branch("PHOS","TClonesArray",&digits,bufferSize);
753 digitsBranch->SetTitle(fEventFolderName);
754 digitsBranch->Fill() ;
756 gime->WriteDigits("OVERWRITE");
757 gime->WriteDigitizer("OVERWRITE");