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.88 2006/03/13 14:05:43 kharlov
22 * Calibration objects for EMC and CPV
24 * Revision 1.87 2005/08/24 15:33:49 kharlov
25 * Calibration data for raw digits
27 * Revision 1.86 2005/07/12 20:07:35 hristov
28 * Changes needed to run simulation and reconstrruction in the same AliRoot session
30 * Revision 1.85 2005/05/28 14:19:04 schutz
31 * Compilation warnings fixed by T.P.
35 //_________________________________________________________________________
36 //*-- Author : Dmitri Peressounko (SUBATECH & Kurchatov Institute)
37 //////////////////////////////////////////////////////////////////////////////
38 // This TTask performs digitization of Summable digits (in the PHOS case it is just
39 // the sum of contributions from all primary particles into a given cell).
40 // In addition it performs mixing of summable digits from different events.
41 // The name of the TTask is also the title of the branch that will contain
42 // the created SDigits
43 // The title of the TTAsk is the name of the file that contains the hits from
44 // which the SDigits are created
46 // For each event two branches are created in TreeD:
47 // "PHOS" - list of digits
48 // "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
50 // Note, that one can set a title for new digits branch, and repeat digitization with
51 // another set of parameters.
54 // root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
55 // root[1] d->ExecuteTask()
56 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
57 // //Digitizes SDigitis in all events found in file galice.root
59 // root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;
60 // // Will read sdigits from galice1.root
61 // root[3] d1->MixWith("galice2.root")
62 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
63 // // Reads another set of sdigits from galice2.root
64 // root[3] d1->MixWith("galice3.root")
65 // // Reads another set of sdigits from galice3.root
66 // root[4] d->ExecuteTask("deb timing")
67 // // Reads SDigits from files galice1.root, galice2.root ....
68 // // mixes them and stores produced Digits in file galice1.root
69 // // deb - prints number of produced digits
70 // // deb all - prints list of produced digits
71 // // timing - prints time used for digitization
74 // --- ROOT system ---
77 #include "TBenchmark.h"
80 // --- Standard library ---
82 // --- AliRoot header files ---
84 #include "AliRunDigitizer.h"
85 #include "AliPHOSDigit.h"
86 #include "AliPHOSGetter.h"
87 #include "AliPHOSDigitizer.h"
88 #include "AliPHOSSDigitizer.h"
89 #include "AliPHOSGeometry.h"
90 #include "AliPHOSTick.h"
92 ClassImp(AliPHOSDigitizer)
95 //____________________________________________________________________________
96 AliPHOSDigitizer::AliPHOSDigitizer():AliDigitizer("",""),
103 fDefaultInit = kTRUE ;
104 fManager = 0 ; // We work in the standalong mode
105 fEventFolderName = "" ;
108 //____________________________________________________________________________
109 AliPHOSDigitizer::AliPHOSDigitizer(TString alirunFileName,
110 TString eventFolderName):
111 AliDigitizer("PHOS"+AliConfig::Instance()->GetDigitizerTaskName(),
113 fInputFileNames(0), fEventNames(0), fEventFolderName(eventFolderName)
118 fDefaultInit = kFALSE ;
119 fManager = 0 ; // We work in the standalong mode
122 //____________________________________________________________________________
123 AliPHOSDigitizer::AliPHOSDigitizer(const AliPHOSDigitizer & d)
128 SetName(d.GetName()) ;
129 SetTitle(d.GetTitle()) ;
130 fPinNoise = d.fPinNoise ;
131 fEMCDigitThreshold = d.fEMCDigitThreshold ;
132 fCPVNoise = d.fCPVNoise ;
133 fCPVDigitThreshold = d.fCPVDigitThreshold ;
134 fTimeResolution = d.fTimeResolution ;
135 fTimeThreshold = d.fTimeThreshold ;
136 fTimeSignalLength = d.fTimeSignalLength ;
137 fADCchanelEmc = d.fADCchanelEmc ;
138 fADCpedestalEmc = d.fADCpedestalEmc ;
139 fNADCemc = d.fNADCemc ;
140 fADCchanelCpv = d.fADCchanelCpv ;
141 fADCpedestalCpv = d.fADCpedestalCpv ;
142 fNADCcpv = d.fNADCcpv ;
143 fEventFolderName = d.fEventFolderName;
146 //____________________________________________________________________________
147 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * rd):
148 AliDigitizer(rd,"PHOS"+AliConfig::Instance()->GetDigitizerTaskName()),
151 // ctor Init() is called by RunDigitizer
153 fEventFolderName = fManager->GetInputFolderName(0) ;
154 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
156 fDefaultInit = kFALSE ;
159 //____________________________________________________________________________
160 AliPHOSDigitizer::~AliPHOSDigitizer()
162 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
164 // Clean Digitizer from the white board
165 gime->PhosLoader()->CleanDigitizer() ;
167 delete [] fInputFileNames ;
168 delete [] fEventNames ;
172 //____________________________________________________________________________
173 void AliPHOSDigitizer::Digitize(Int_t event)
176 // Makes the digitization of the collected summable digits.
177 // It first creates the array of all PHOS modules
178 // filled with noise (different for EMC, and CPV) and
179 // then adds contributions from SDigits.
180 // This design avoids scanning over the list of digits to add
181 // contribution to new SDigits only.
183 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
184 Int_t ReadEvent = event ;
186 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
187 AliInfo(Form("Adding event %d from input stream 0 %s %s",
188 ReadEvent, GetTitle(), fEventFolderName.Data())) ;
189 gime->Event(ReadEvent, "S") ;
190 TClonesArray * digits = gime->Digits() ;
193 const AliPHOSGeometry *geom = gime->PHOSGeometry() ;
194 //Making digits with noise, first EMC
195 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
200 nCPV = nEMC + geom->GetNumberOfCPVPadsZ() * geom->GetNumberOfCPVPadsPhi() * geom->GetNModules() ;
202 digits->Expand(nCPV) ;
204 // get first the sdigitizer from the tasks list
205 if ( !gime->SDigitizer() )
206 gime->LoadSDigitizer();
207 AliPHOSSDigitizer * sDigitizer = gime->SDigitizer();
210 AliFatal(Form("SDigitizer with name %s %s not found",
211 GetTitle(), fEventFolderName.Data() )) ;
213 //take all the inputs to add together and load the SDigits
214 TObjArray * sdigArray = new TObjArray(fInput) ;
215 sdigArray->AddAt(gime->SDigits(), 0) ;
217 for(i = 1 ; i < fInput ; i++){
218 TString tempo(fEventNames[i]) ;
220 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[i], tempo) ;
222 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
223 AliInfo(Form("Adding event %d from input stream %d %s %s",
224 ReadEvent, i, fInputFileNames[i].Data(), tempo.Data())) ;
225 gime->Event(ReadEvent,"S");
226 sdigArray->AddAt(gime->SDigits(), i) ;
229 //Find the first crystall with signal
230 Int_t nextSig = 200000 ;
231 TClonesArray * sdigits ;
232 for(i = 0 ; i < fInput ; i++){
233 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
234 if ( !sdigits->GetEntriesFast() )
236 Int_t curNext = dynamic_cast<AliPHOSDigit *>(sdigits->At(0))->GetId() ;
237 if(curNext < nextSig)
241 TArrayI index(fInput) ;
242 index.Reset() ; //Set all indexes to zero
244 AliPHOSDigit * digit ;
245 AliPHOSDigit * curSDigit ;
247 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
249 //Put Noise contribution
250 for(absID = 1 ; absID <= nEMC ; absID++){
251 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
252 new((*digits)[absID-1]) AliPHOSDigit( -1, absID, sDigitizer->Digitize(noise), TimeOfNoise() ) ;
253 //look if we have to add signal?
254 digit = dynamic_cast<AliPHOSDigit *>(digits->At(absID-1)) ;
257 //Add SDigits from all inputs
260 Float_t a = digit->GetAmp() ;
261 Float_t b = TMath::Abs( a / fTimeSignalLength) ;
262 //Mark the beginning of the signal
263 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
264 //Mark the end of the ignal
265 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
268 for(i = 0 ; i < fInput ; i++){
269 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
270 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
273 //May be several digits will contribute from the same input
274 while(curSDigit && curSDigit->GetId() == absID){
275 //Shift primary to separate primaries belonging different inputs
276 Int_t primaryoffset ;
278 primaryoffset = fManager->GetMask(i) ;
280 primaryoffset = 10000000*i ;
281 curSDigit->ShiftPrimary(primaryoffset) ;
283 a = curSDigit->GetAmp() ;
284 b = a /fTimeSignalLength ;
285 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
286 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
288 *digit = *digit + *curSDigit ; //add energies
291 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
292 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
298 //calculate and set time
299 Float_t time = FrontEdgeTime(ticks) ;
300 digit->SetTime(time) ;
302 //Find next signal module
304 for(i = 0 ; i < fInput ; i++){
305 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
306 Int_t curNext = nextSig ;
307 if(sdigits->GetEntriesFast() > index[i] ){
308 curNext = dynamic_cast<AliPHOSDigit *>(sdigits->At(index[i]))->GetId() ;
310 if(curNext < nextSig) nextSig = curNext ;
318 //Now CPV digits (different noise and no timing)
319 for(absID = nEMC+1; absID <= nCPV; absID++){
320 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
321 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
322 //look if we have to add signal?
324 digit = dynamic_cast<AliPHOSDigit *>(digits->At(absID-1)) ;
325 //Add SDigits from all inputs
326 for(i = 0 ; i < fInput ; i++){
327 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
328 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
332 //May be several digits will contribute from the same input
333 while(curSDigit && curSDigit->GetId() == absID){
334 //Shift primary to separate primaries belonging different inputs
335 Int_t primaryoffset ;
337 primaryoffset = fManager->GetMask(i) ;
339 primaryoffset = 10000000*i ;
340 curSDigit->ShiftPrimary(primaryoffset) ;
343 *digit = *digit + *curSDigit ;
345 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
346 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i]) ) ;
352 //Find next signal module
354 for(i = 0 ; i < fInput ; i++){
355 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
356 Int_t curNext = nextSig ;
357 if(sdigits->GetEntriesFast() > index[i] )
358 curNext = dynamic_cast<AliPHOSDigit *>( sdigits->At(index[i]) )->GetId() ;
359 if(curNext < nextSig) nextSig = curNext ;
365 delete sdigArray ; //We should not delete its contents
367 //remove digits below thresholds
368 for(i = 0 ; i < nEMC ; i++){
369 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
370 if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
371 digits->RemoveAt(i) ;
373 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
377 for(i = nEMC; i < nCPV ; i++)
378 if( sDigitizer->Calibrate( dynamic_cast<AliPHOSDigit*>(digits->At(i))->GetAmp() ) < fCPVDigitThreshold )
379 digits->RemoveAt(i) ;
383 Int_t ndigits = digits->GetEntriesFast() ;
384 digits->Expand(ndigits) ;
386 //Set indexes in list of digits and make true digitization of the energy
387 for (i = 0 ; i < ndigits ; i++) {
388 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
389 digit->SetIndexInList(i) ;
390 Float_t energy = sDigitizer->Calibrate(digit->GetAmp()) ;
391 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
395 //____________________________________________________________________________
396 Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
398 // Returns digitized value of the energy in a cell absId
400 AliPHOSGetter* gime = AliPHOSGetter::Instance();
402 if(!gime->CalibData()) {
403 //AliPHOSCalibData* cdb = new AliPHOSCalibData(gAlice->GetRunNumber()); // original
404 AliPHOSCalibData* cdb = new AliPHOSCalibData(-1); // use AliCDBManager's run number
405 gime->SetCalibData(cdb);
408 //Determine rel.position of the cell absId
410 gime->PHOSGeometry()->AbsToRelNumbering(absId,relId);
411 Int_t module=relId[0];
413 Int_t column=relId[3];
417 if(absId <= fEmcCrystals){ //digitize as EMC
419 //reading calibration data for cell absId.
420 //If no calibration DB found, accept default values.
422 if(gime->CalibData()) {
423 fADCpedestalEmc = gime->CalibData()->GetADCpedestalEmc(module,column,raw);
424 fADCchanelEmc = gime->CalibData()->GetADCchannelEmc(module,column,raw);
426 // printf("\t\tabsId %d ==>>module=%d column=%d raw=%d ped=%.4f gain=%.4f\n",
427 // absId,module,column,raw,fADCpedestalEmc,fADCchanelEmc);
429 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
430 if(chanel > fNADCemc ) chanel = fNADCemc ;
432 else{ //Digitize as CPV
433 if(gime->CalibData()) {
434 fADCpedestalCpv = gime->CalibData()->GetADCpedestalCpv(module,column,raw);
435 fADCchanelCpv = gime->CalibData()->GetADCchannelCpv(module,column,raw);
438 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
439 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
444 //____________________________________________________________________________
445 void AliPHOSDigitizer::Exec(Option_t *option)
447 // Steering method to process digitization for events
448 // in the range from fFirstEvent to fLastEvent.
449 // This range is optionally set by SetEventRange().
450 // if fLastEvent=-1, then process events until the end.
451 // by default fLastEvent = fFirstEvent (process only one event)
453 if (!fInit) { // to prevent overwrite existing file
454 AliError(Form("Give a version name different from %s",
455 fEventFolderName.Data() )) ;
459 if (strstr(option,"print")) {
464 if(strstr(option,"tim"))
465 gBenchmark->Start("PHOSDigitizer");
467 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
469 // Post Digitizer to the white board
470 gime->PostDigitizer(this) ;
472 if (fLastEvent == -1)
473 fLastEvent = gime->MaxEvent() - 1 ;
475 fLastEvent = fFirstEvent ;
477 Int_t nEvents = fLastEvent - fFirstEvent + 1;
481 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
483 gime->Event(ievent,"S") ;
485 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
489 if(strstr(option,"deb"))
492 //increment the total number of Digits per run
493 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
496 gime->PhosLoader()->CleanDigitizer();
498 if(strstr(option,"tim")){
499 gBenchmark->Stop("PHOSDigitizer");
501 message = " took %f seconds for Digitizing %f seconds per event\n" ;
502 AliInfo(Form( message.Data(),
503 gBenchmark->GetCpuTime("PHOSDigitizer"),
504 gBenchmark->GetCpuTime("PHOSDigitizer")/nEvents ));
508 //____________________________________________________________________________
509 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks) const
511 // Returns the shortest time among all time ticks
513 ticks->Sort() ; //Sort in accordance with times of ticks
515 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
516 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
519 while((t=(AliPHOSTick*) it.Next())){
520 if(t->GetTime() < time) //This tick starts before crossing
525 time = ctick->CrossingTime(fTimeThreshold) ;
530 //____________________________________________________________________________
531 Bool_t AliPHOSDigitizer::Init()
533 // Makes all memory allocations
535 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
537 AliFatal(Form("Could not obtain the Getter object for file %s and event %s !",
538 GetTitle(), fEventFolderName.Data()));
542 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
544 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
546 TString opt("Digits") ;
547 if(gime->VersionExists(opt) ) {
548 AliError(Form("Give a version name different from %s",
549 fEventFolderName.Data() )) ;
554 fLastEvent = fFirstEvent ;
556 fInput = fManager->GetNinputs() ;
560 fInputFileNames = new TString[fInput] ;
561 fEventNames = new TString[fInput] ;
562 fInputFileNames[0] = GetTitle() ;
563 fEventNames[0] = fEventFolderName.Data() ;
565 for (index = 1 ; index < fInput ; index++) {
566 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
567 TString tempo = fManager->GetInputFolderName(index) ;
568 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added by fManager
571 //to prevent cleaning of this object while GetEvent is called
572 gime->PhosLoader()->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
577 //____________________________________________________________________________
578 void AliPHOSDigitizer::InitParameters()
580 // Set initial parameters Digitizer
583 fEMCDigitThreshold = 0.012 ;
585 fCPVDigitThreshold = 0.09 ;
586 fTimeResolution = 0.5e-9 ;
587 fTimeSignalLength = 1.0e-9 ;
589 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
590 fADCpedestalEmc = 0.005 ; //
591 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
593 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
594 fADCpedestalCpv = 0.012 ; //
595 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
597 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
598 SetEventRange(0,-1) ;
602 //__________________________________________________________________
603 void AliPHOSDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
605 // Allows to produce digits by superimposing background and signal event.
606 // It is assumed, that headers file with SIGNAL events is opened in
608 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
609 // Thus we avoid writing (changing) huge and expensive
610 // backgound files: all output will be writen into SIGNAL, i.e.
611 // opened in constructor file.
613 // One can open as many files to mix with as one needs.
614 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
617 if( strcmp(fEventFolderName, "") == 0 )
621 Warning("MixWith", "Cannot use this method with AliRunDigitizer\n" ) ;
624 // looking for file which contains AliRun
625 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
626 AliError(Form("File %s does not exist!", alirunFileName.Data())) ;
629 // looking for the file which contains SDigits
630 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
631 TString fileName( gime->GetSDigitsFileName() ) ;
632 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
633 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
634 if ( (gSystem->AccessPathName(fileName)) ) {
635 AliError(Form("The file %s does not exist!", fileName.Data())) ;
638 // need to increase the arrays
639 TString tempo = fInputFileNames[fInput-1] ;
640 delete [] fInputFileNames ;
641 fInputFileNames = new TString[fInput+1] ;
642 fInputFileNames[fInput-1] = tempo ;
644 tempo = fEventNames[fInput-1] ;
645 delete [] fEventNames ;
646 fEventNames = new TString[fInput+1] ;
647 fEventNames[fInput-1] = tempo ;
649 fInputFileNames[fInput] = alirunFileName ;
650 fEventNames[fInput] = eventFolderName ;
654 //__________________________________________________________________
655 void AliPHOSDigitizer::Print(const Option_t *)const
657 // Print Digitizer's parameters
658 AliInfo(Form("\n------------------- %s -------------", GetName() )) ;
659 if( strcmp(fEventFolderName.Data(), "") != 0 ){
660 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
664 nStreams = GetNInputStreams() ;
669 for (index = 0 ; index < nStreams ; index++) {
670 TString tempo(fEventNames[index]) ;
672 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[index], tempo) ;
673 TString fileName( gime->GetSDigitsFileName() ) ;
674 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
675 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
676 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
678 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
679 printf("\nWriting digits to %s", gime->GetDigitsFileName().Data()) ;
681 printf("\nWith following parameters:\n") ;
682 printf(" Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise ) ;
683 printf(" Threshold in EMC (fEMCDigitThreshold) = %f\n", fEMCDigitThreshold ) ;
684 printf(" Noise in CPV (fCPVNoise) = %f\n", fCPVNoise ) ;
685 printf(" Threshold in CPV (fCPVDigitThreshold) = %f\n",fCPVDigitThreshold ) ;
686 printf(" ---------------------------------------------------\n") ;
689 AliInfo(Form("AliPHOSDigitizer not initialized" )) ;
693 //__________________________________________________________________
694 void AliPHOSDigitizer::PrintDigits(Option_t * option)
696 // Print a table of digits
698 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
699 TClonesArray * digits = gime->Digits() ;
701 AliInfo(Form("%d", digits->GetEntriesFast())) ;
702 printf("\nevent %d", gAlice->GetEvNumber()) ;
703 printf("\n Number of entries in Digits list %d", digits->GetEntriesFast() ) ;
706 if(strstr(option,"all")||strstr(option,"EMC")){
708 AliPHOSDigit * digit;
709 printf("\nEMC digits (with primaries):\n") ;
710 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
711 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
713 for (index = 0 ; (index < digits->GetEntriesFast()) &&
714 (dynamic_cast<AliPHOSDigit *>(digits->At(index))->GetId() <= maxEmc) ; index++) {
715 digit = (AliPHOSDigit * ) digits->At(index) ;
716 if(digit->GetNprimary() == 0)
718 printf("%6d %8d %6.5e %4d %2d :",
719 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
721 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
722 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
728 if(strstr(option,"all")||strstr(option,"CPV")){
730 //loop over CPV digits
731 AliPHOSDigit * digit;
732 printf("\nCPV digits (with primaries):\n") ;
733 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
734 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
736 for (index = 0 ; index < digits->GetEntriesFast(); index++) {
737 digit = (AliPHOSDigit * ) digits->At(index) ;
738 if(digit->GetId() > maxEmc){
739 printf("%6d %8d %4d %2d :",
740 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
742 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
743 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
752 //__________________________________________________________________
753 Float_t AliPHOSDigitizer::TimeOfNoise(void) const
754 { // Calculates the time signal generated by noise
755 //PH Info("TimeOfNoise", "Change me") ;
756 return gRandom->Rndm() * 1.28E-5;
759 //__________________________________________________________________
760 void AliPHOSDigitizer::Unload()
764 for(i = 1 ; i < fInput ; i++){
765 TString tempo(fEventNames[i]) ;
767 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[i], tempo) ;
768 gime->PhosLoader()->UnloadSDigits() ;
771 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
772 gime->PhosLoader()->UnloadDigits() ;
775 //____________________________________________________________________________
776 void AliPHOSDigitizer::WriteDigits()
779 // Makes TreeD in the output file.
780 // Check if branch already exists:
781 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
782 // already existing branches.
783 // else creates branch with Digits, named "PHOS", title "...",
784 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
785 // and names of files, from which digits are made.
787 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
788 const TClonesArray * digits = gime->Digits() ;
789 TTree * treeD = gime->TreeD();
791 // -- create Digits branch
792 Int_t bufferSize = 32000 ;
793 TBranch * digitsBranch = treeD->Branch("PHOS","TClonesArray",&digits,bufferSize);
794 digitsBranch->SetTitle(fEventFolderName);
795 digitsBranch->Fill() ;
797 gime->WriteDigits("OVERWRITE");
798 gime->WriteDigitizer("OVERWRITE");