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 //_________________________________________________________________________
20 //////////////////////////////////////////////////////////////////////////////
21 // Class performs digitization of Summable digits
23 // In addition it performs mixing of summable digits from different events.
25 // For each event two branches are created in TreeD:
26 // "EMCAL" - list of digits
27 // "AliEMCALDigitizer" - AliEMCALDigitizer with all parameters used in digitization
29 // Note, that one cset title for new digits branch, and repeat digitization with
30 // another set of parameters.
33 // root[0] AliEMCALDigitizer * d = new AliEMCALDigitizer() ;
34 // root[1] d->ExecuteTask()
35 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
36 // //Digitizes SDigitis in all events found in file galice.root
38 // root[2] AliEMCALDigitizer * d1 = new AliEMCALDigitizer("galice1.root") ;
39 // // Will read sdigits from galice1.root
40 // root[3] d1->MixWith("galice2.root")
41 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 // // Reads another portion of sdigits from galice2.root
43 // root[3] d1->MixWith("galice3.root")
44 // // Reads another portion of sdigits from galice3.root
45 // root[4] d->ExecuteTask("deb timing")
46 // // Reads SDigits from files galice1.root, galice2.root ....
47 // // mixes them and stores produced Digits in file galice1.root
48 // // deb - prints number of produced digits
49 // // deb all - prints list of produced digits
50 // // timing - prints time used for digitization
51 ////////////////////////////////////////////////////////////////////////////////////
53 //*-- Author: Sahal Yacoob (LBL)
54 // based on : AliEMCALDigitizer
56 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
57 // of new IO (à la PHOS)
58 // November 2003 Aleksei Pavlinov : adopted for Shish-Kebab geometry
59 //_________________________________________________________________________________
61 // --- ROOT system ---
65 #include "TBenchmark.h"
69 #include "TObjectTable.h"
71 // --- AliRoot header files ---
73 #include "AliRunDigitizer.h"
74 #include "AliRunLoader.h"
75 #include "AliEMCALDigit.h"
77 #include "AliEMCALLoader.h"
78 #include "AliEMCALDigitizer.h"
79 #include "AliEMCALSDigitizer.h"
80 #include "AliEMCALGeometry.h"
81 #include "AliEMCALTick.h"
82 #include "AliEMCALHistoUtilities.h"
84 ClassImp(AliEMCALDigitizer)
87 //____________________________________________________________________________
88 AliEMCALDigitizer::AliEMCALDigitizer():AliDigitizer("",""),
95 fDefaultInit = kTRUE ;
96 fManager = 0 ; // We work in the standalong mode
97 fEventFolderName = "" ;
100 //____________________________________________________________________________
101 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName):
102 AliDigitizer("EMCAL"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
103 fInputFileNames(0), fEventNames(0), fEventFolderName(eventFolderName)
109 fDefaultInit = kFALSE ;
110 fManager = 0 ; // We work in the standalong mode
113 //____________________________________________________________________________
114 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d) : AliDigitizer(d)
118 SetName(d.GetName()) ;
119 SetTitle(d.GetTitle()) ;
120 fDigitThreshold = d.fDigitThreshold ;
121 fMeanPhotonElectron = d.fMeanPhotonElectron ;
122 fPedestal = d.fPedestal ;
124 fPinNoise = d.fPinNoise ;
125 fTimeResolution = d.fTimeResolution ;
126 fTimeThreshold = d.fTimeThreshold ;
127 fTimeSignalLength = d.fTimeSignalLength ;
128 fADCchannelEC = d.fADCchannelEC ;
129 fADCpedestalEC = d.fADCpedestalEC ;
130 fNADCEC = d.fNADCEC ;
131 fEventFolderName = d.fEventFolderName;
134 //____________________________________________________________________________
135 AliEMCALDigitizer::AliEMCALDigitizer(AliRunDigitizer * rd):
136 AliDigitizer(rd,"EMCAL"+AliConfig::Instance()->GetDigitizerTaskName()),
139 // ctor Init() is called by RunDigitizer
141 fEventFolderName = fManager->GetInputFolderName(0) ;
142 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
144 fDefaultInit = kFALSE ;
147 //____________________________________________________________________________
148 AliEMCALDigitizer::~AliEMCALDigitizer()
150 if (AliRunLoader::GetRunLoader()) {
151 AliLoader *emcalLoader=0;
152 if ((emcalLoader = AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL")))
153 emcalLoader->CleanDigitizer();
156 AliDebug(1," no runloader present");
157 delete [] fInputFileNames ;
158 delete [] fEventNames ;
160 if(fHists) delete fHists;
163 //____________________________________________________________________________
164 void AliEMCALDigitizer::Digitize(Int_t event)
167 // Makes the digitization of the collected summable digits
168 // for this it first creates the array of all EMCAL modules
169 // filled with noise (different for EMC, CPV and PPSD) and
170 // after that adds contributions from SDigits. This design
171 // helps to avoid scanning over the list of digits to add
172 // contribution of any new SDigit.
173 static int isTrd1Geom = -1; // -1 - mean undefined
174 static int nEMC=0; //max number of digits possible
176 AliRunLoader *rl = AliRunLoader::GetRunLoader();
177 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
178 Int_t ReadEvent = event ;
179 // fManager is data member from AliDigitizer
181 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
182 AliDebug(1,Form("Adding event %d from input stream 0 %s %s",
183 ReadEvent, GetTitle(), fEventFolderName.Data())) ;
184 rl->GetEvent(ReadEvent);
186 TClonesArray * digits = emcalLoader->Digits() ;
190 // const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
192 AliRun * gAlice = rl->GetAliRun();
193 AliEMCAL * emcal = (AliEMCAL*)gAlice->GetDetector("EMCAL");
194 AliEMCALGeometry * geom = emcal->GetGeometry();
197 TString ng(geom->GetName());
199 if(ng.Contains("SHISH") && ng.Contains("TRD1")) isTrd1Geom = 1;
201 if(isTrd1Geom == 0) nEMC = geom->GetNPhi()*geom->GetNZ();
202 else nEMC = geom->GetNCells();
203 AliDebug(1,Form("nEMC %i (number cells in EMCAL) | %s | isTrd1Geom %i\n", nEMC, geom->GetName(), isTrd1Geom));
207 digits->Expand(nEMC) ;
209 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
210 if ( !emcalLoader->SDigitizer() )
211 emcalLoader->LoadSDigitizer();
212 AliEMCALSDigitizer * sDigitizer = dynamic_cast<AliEMCALSDigitizer*>(emcalLoader->SDigitizer());
215 Fatal("Digitize", "SDigitizer with name %s %s not found", fEventFolderName.Data(), GetTitle() ) ;
217 //take all the inputs to add together and load the SDigits
218 TObjArray * sdigArray = new TObjArray(fInput) ;
219 sdigArray->AddAt(emcalLoader->SDigits(), 0) ;
222 for(i = 1 ; i < fInput ; i++){
223 TString tempo(fEventNames[i]) ;
226 AliRunLoader * rl2 = AliRunLoader::GetRunLoader(tempo) ;
229 rl2 = AliRunLoader::Open(fInputFileNames[i], tempo) ;
232 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
233 Info("Digitize", "Adding event %d from input stream %d %s %s", ReadEvent, i, fInputFileNames[i].Data(), tempo.Data()) ;
235 rl2->GetEvent(ReadEvent);
236 AliEMCALLoader *emcalLoader2 = dynamic_cast<AliEMCALLoader*>(rl2->GetDetectorLoader("EMCAL"));
237 sdigArray->AddAt(emcalLoader2->SDigits(), i) ;
240 //Find the first tower with signal
241 Int_t nextSig = nEMC + 1 ;
242 TClonesArray * sdigits ;
243 for(i = 0 ; i < fInput ; i++){
244 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
245 if ( !sdigits->GetEntriesFast() )
247 Int_t curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(0))->GetId() ;
248 if(curNext < nextSig)
250 AliDebug(1,Form("input %i : #sdigits %i \n",
251 i, sdigits->GetEntriesFast()));
253 AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
255 TArrayI index(fInput) ;
256 index.Reset() ; //Set all indexes to zero
258 AliEMCALDigit * digit ;
259 AliEMCALDigit * curSDigit ;
261 TClonesArray * ticks = new TClonesArray("AliEMCALTick",1000) ;
263 //Put Noise contribution
264 for(absID = 1; absID <= nEMC; absID++){
266 // amplitude set to zero, noise will be added later
267 new((*digits)[absID-1]) AliEMCALDigit( -1, -1, absID, 0, TimeOfNoise() ) ;
268 //look if we have to add signal?
269 digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID-1)) ;
272 //Add SDigits from all inputs
275 Float_t a = digit->GetAmp() ;
276 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
277 //Mark the beginning of the signal
278 new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime(),0, b);
279 //Mark the end of the signal
280 new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime()+fTimeSignalLength, -a, -b);
283 for(i = 0; i< fInput ; i++){ //loop over (possible) merge sources
284 if(dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
285 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
288 //May be several digits will contribute from the same input
289 while(curSDigit && (curSDigit->GetId() == absID)){
290 //Shift primary to separate primaries belonging different inputs
291 Int_t primaryoffset ;
293 primaryoffset = fManager->GetMask(i) ;
296 curSDigit->ShiftPrimary(primaryoffset) ;
298 a = curSDigit->GetAmp() ;
299 b = a /fTimeSignalLength ;
300 new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime(),0, b);
301 new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
303 *digit = *digit + *curSDigit ; //add energies
306 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
307 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
312 // add fluctuations for photo-electron creation
313 amp = sDigitizer->Calibrate(digit->GetAmp()) ; // GeV
314 amp *= static_cast<Float_t>(gRandom->Poisson(fMeanPhotonElectron)) / static_cast<Float_t>(fMeanPhotonElectron) ;
316 //calculate and set time
317 Float_t time = FrontEdgeTime(ticks) ;
318 digit->SetTime(time) ;
320 //Find next signal module
322 for(i = 0 ; i < fInput ; i++){
323 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
324 Int_t curNext = nextSig ;
325 if(sdigits->GetEntriesFast() > index[i] ){
326 curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))->GetId() ;
328 if(curNext < nextSig) nextSig = curNext ;
332 amp += TMath::Abs(gRandom->Gaus(0., fPinNoise)) ;
333 digit->SetAmp(sDigitizer->Digitize(amp)) ;
334 AliDebug(10,Form(" absID %5i amp %f nextSig %5i\n",
335 absID, amp, nextSig));
336 } // for(absID = 1; absID <= nEMC; absID++)
341 delete sdigArray ; //We should not delete its contents
343 //remove digits below thresholds
344 for(i = 0 ; i < nEMC ; i++){
345 digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
346 Float_t threshold = fDigitThreshold ;
347 if(sDigitizer->Calibrate( digit->GetAmp() ) < threshold)
348 digits->RemoveAt(i) ;
350 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
355 Int_t ndigits = digits->GetEntriesFast() ;
356 digits->Expand(ndigits) ;
358 //Set indexes in list of digits and fill hists.
359 AliEMCALHistoUtilities::FillH1(fHists, 0, Double_t(ndigits));
360 Float_t energy=0., esum=0.;
361 for (i = 0 ; i < ndigits ; i++) {
362 digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ;
363 digit->SetIndexInList(i) ;
364 energy = sDigitizer->Calibrate(digit->GetAmp()) ;
366 digit->SetAmp(DigitizeEnergy(energy, digit->GetId()) ) ; // for what ??
367 AliEMCALHistoUtilities::FillH1(fHists, 2, double(digit->GetAmp()));
368 AliEMCALHistoUtilities::FillH1(fHists, 3, double(energy));
369 AliEMCALHistoUtilities::FillH1(fHists, 4, double(digit->GetId()));
371 AliEMCALHistoUtilities::FillH1(fHists, 1, esum);
374 // //_____________________________________________________________________
375 Int_t AliEMCALDigitizer::DigitizeEnergy(Float_t energy, Int_t AbsId)
377 // Returns digitized value of the energy in a cell absId
379 AliRunLoader *rl = AliRunLoader::GetRunLoader();
380 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>
381 (rl->GetDetectorLoader("EMCAL"));
383 // Load EMCAL Geometry
385 AliRun * gAlice = rl->GetAliRun();
386 AliEMCAL * emcal = (AliEMCAL*)gAlice->GetDetector("EMCAL");
387 AliEMCALGeometry * geom = emcal->GetGeometry();
390 AliFatal("Did not get geometry from EMCALLoader") ;
398 Int_t channel = -999;
400 Bool_t bCell = geom->GetCellIndex(AbsId, iSupMod, nTower, nIphi, nIeta) ;
402 Error("DigitizeEnergy","Wrong cell id number") ;
403 geom->GetCellPhiEtaIndexInSModule(iSupMod,nTower,nIphi, nIeta,iphi,ieta);
405 if(emcalLoader->CalibData()) {
406 fADCpedestalEC = emcalLoader->CalibData()
407 ->GetADCpedestal(iSupMod,ieta,iphi);
408 fADCchannelEC = emcalLoader->CalibData()
409 ->GetADCchannel(iSupMod,ieta,iphi);
412 channel = static_cast<Int_t> (TMath::Ceil( (energy + fADCpedestalEC)/fADCchannelEC )) ;
414 if(channel > fNADCEC )
420 //____________________________________________________________________________
421 void AliEMCALDigitizer::Exec(Option_t *option)
423 // Steering method to process digitization for events
424 // in the range from fFirstEvent to fLastEvent.
425 // This range is optionally set by SetEventRange().
426 // if fLastEvent=-1, then process events until the end.
427 // by default fLastEvent = fFirstEvent (process only one event)
429 if (!fInit) { // to prevent overwrite existing file
430 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
434 if (strstr(option,"print")) {
439 if(strstr(option,"tim"))
440 gBenchmark->Start("EMCALDigitizer");
442 AliRunLoader *rl = AliRunLoader::GetRunLoader();
443 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
445 // Post Digitizer to the white board
446 emcalLoader->PostDigitizer(this) ;
448 if (fLastEvent == -1)
449 fLastEvent = rl->GetNumberOfEvents() - 1 ;
451 fLastEvent = fFirstEvent ; // what is this ??
453 Int_t nEvents = fLastEvent - fFirstEvent + 1;
456 rl->LoadSDigits("EMCAL");
457 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
459 rl->GetEvent(ievent);
461 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
465 if(strstr(option,"deb"))
467 if(strstr(option,"table")) gObjectTable->Print();
469 //increment the total number of Digits per run
470 fDigitsInRun += emcalLoader->Digits()->GetEntriesFast() ;
473 emcalLoader->CleanDigitizer() ;
475 if(strstr(option,"tim")){
476 gBenchmark->Stop("EMCALDigitizer");
477 printf("Exec: took %f seconds for Digitizing %f seconds per event",
478 gBenchmark->GetCpuTime("EMCALDigitizer"), gBenchmark->GetCpuTime("EMCALDigitizer")/nEvents ) ;
482 //____________________________________________________________________________
483 Float_t AliEMCALDigitizer::FrontEdgeTime(TClonesArray * ticks)
485 // Returns the shortest time among all time ticks
487 ticks->Sort() ; //Sort in accordance with times of ticks
489 AliEMCALTick * ctick = (AliEMCALTick *) it.Next() ;
490 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
493 while((t=(AliEMCALTick*) it.Next())){
494 if(t->GetTime() < time) //This tick starts before crossing
499 time = ctick->CrossingTime(fTimeThreshold) ;
504 //____________________________________________________________________________
505 Bool_t AliEMCALDigitizer::Init()
507 // Makes all memory allocations
509 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
511 if ( emcalLoader == 0 ) {
512 Fatal("Init", "Could not obtain the AliEMCALLoader");
517 fLastEvent = fFirstEvent ;
520 fInput = fManager->GetNinputs() ;
524 fInputFileNames = new TString[fInput] ;
525 fEventNames = new TString[fInput] ;
526 fInputFileNames[0] = GetTitle() ;
527 fEventNames[0] = fEventFolderName.Data() ;
529 for (index = 1 ; index < fInput ; index++) {
530 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
531 TString tempo = fManager->GetInputFolderName(index) ;
532 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fManager
535 //to prevent cleaning of this object while GetEvent is called
536 emcalLoader->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
543 //____________________________________________________________________________
544 void AliEMCALDigitizer::InitParameters()
545 { // Tune parameters - 24-nov-04
547 fMeanPhotonElectron = 3300 ; // electrons per GeV
549 if (fPinNoise == 0. )
550 Warning("InitParameters", "No noise added\n") ;
551 fDigitThreshold = fPinNoise * 3; // 3 * sigma
552 fTimeResolution = 0.3e-9 ; // 300 psc
553 fTimeSignalLength = 1.0e-9 ;
555 fADCchannelEC = 0.00305; // 200./65536 - width of one ADC channel in GeV
556 fADCpedestalEC = 0.009 ; // GeV
557 fNADCEC = (Int_t) TMath::Power(2,16) ; // number of channels in Tower ADC - 65536
559 fTimeThreshold = 0.001*10000000 ; // Means 1 MeV in terms of SDigits amplitude ??
560 // hists. for control; no hists on default
565 //__________________________________________________________________
566 void AliEMCALDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
568 // Allows to produce digits by superimposing background and signal event.
569 // It is assumed, that headers file with SIGNAL events is opened in
571 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
572 // Thus we avoid writing (changing) huge and expensive
573 // backgound files: all output will be writen into SIGNAL, i.e.
574 // opened in constructor file.
576 // One can open as many files to mix with as one needs.
577 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
580 if( strcmp(GetName(), "") == 0 )
584 Error("MixWith", "Cannot use this method under AliRunDigitizer") ;
587 // looking for file which contains AliRun
588 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
589 Error("MixWith", "File %s does not exist!", alirunFileName.Data()) ;
592 // looking for the file which contains SDigits
593 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
594 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
595 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
596 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
597 if ( (gSystem->AccessPathName(fileName)) ) {
598 Error("MixWith", "The file %s does not exist!", fileName.Data()) ;
601 // need to increase the arrays
602 TString tempo = fInputFileNames[fInput-1] ;
603 delete [] fInputFileNames ;
604 fInputFileNames = new TString[fInput+1] ;
605 fInputFileNames[fInput-1] = tempo ;
607 tempo = fEventNames[fInput-1] ;
608 delete [] fEventNames ;
609 fEventNames = new TString[fInput+1] ;
610 fEventNames[fInput-1] = tempo ;
612 fInputFileNames[fInput] = alirunFileName ;
613 fEventNames[fInput] = eventFolderName ;
617 void AliEMCALDigitizer::Print1(Option_t * option)
618 { // 19-nov-04 - just for convinience
623 //__________________________________________________________________
624 void AliEMCALDigitizer::Print(Option_t*)const
626 // Print Digitizer's parameters
627 printf("Print: \n------------------- %s -------------", GetName() ) ;
628 if( strcmp(fEventFolderName.Data(), "") != 0 ){
629 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
633 nStreams = GetNInputStreams() ;
641 for (index = 0 ; index < nStreams ; index++) {
642 TString tempo(fEventNames[index]) ;
644 if ((rl = AliRunLoader::GetRunLoader(tempo)) == 0)
645 rl = AliRunLoader::Open(fInputFileNames[index], tempo) ;
646 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
647 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
648 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
649 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
650 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
653 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
655 printf("\nWriting digits to %s", emcalLoader->GetDigitsFileName().Data()) ;
657 printf("\nWith following parameters:\n") ;
659 printf(" Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise) ;
660 printf(" Threshold in EMC (fDigitThreshold) = %f\n", fDigitThreshold) ;
661 printf("---------------------------------------------------\n") ;
664 printf("Print: AliEMCALDigitizer not initialized") ;
667 //__________________________________________________________________
668 void AliEMCALDigitizer::PrintDigits(Option_t * option){
669 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
670 TClonesArray * digits = emcalLoader->Digits() ;
671 TClonesArray * sdigits = emcalLoader->SDigits() ;
673 printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ;
674 printf("\n event %d", emcalLoader->GetRunLoader()->GetEventNumber());
676 if(strstr(option,"all")){
678 AliEMCALDigit * digit;
679 printf("\nEMC digits (with primaries):\n") ;
680 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
682 for (index = 0 ; index < digits->GetEntries() ; index++) {
683 digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
684 printf("\n%6d %8d %6.5e %4d %2d : ",
685 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
687 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
688 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
695 //__________________________________________________________________
696 Float_t AliEMCALDigitizer::TimeOfNoise(void)
698 // Calculates the time signal generated by noise
699 //PH Info("TimeOfNoise", "Change me") ;
700 return gRandom->Rndm() * 1.28E-5;
703 //__________________________________________________________________
704 void AliEMCALDigitizer::Unload()
706 // Unloads the SDigits and Digits
710 for(i = 1 ; i < fInput ; i++){
711 TString tempo(fEventNames[i]) ;
713 if ((rl = AliRunLoader::GetRunLoader(tempo)))
714 rl->GetDetectorLoader("EMCAL")->UnloadSDigits() ;
716 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
717 emcalLoader->UnloadDigits() ;
720 //_________________________________________________________________________________________
721 void AliEMCALDigitizer::WriteDigits()
724 // Makes TreeD in the output file.
725 // Check if branch already exists:
726 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
727 // already existing branches.
728 // else creates branch with Digits, named "EMCAL", title "...",
729 // and branch "AliEMCALDigitizer", with the same title to keep all the parameters
730 // and names of files, from which digits are made.
732 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
734 const TClonesArray * digits = emcalLoader->Digits() ;
735 TTree * treeD = emcalLoader->TreeD();
737 emcalLoader->MakeDigitsContainer();
738 treeD = emcalLoader->TreeD();
741 // -- create Digits branch
742 Int_t bufferSize = 32000 ;
743 TBranch * digitsBranch = 0;
744 if ((digitsBranch = treeD->GetBranch("EMCAL")))
745 digitsBranch->SetAddress(&digits);
747 treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
748 //digitsBranch->SetTitle(fEventFolderName);
751 emcalLoader->WriteDigits("OVERWRITE");
752 emcalLoader->WriteDigitizer("OVERWRITE");
758 void AliEMCALDigitizer::Browse(TBrowser* b)
760 if(fHists) b->Add(fHists);
764 TList *AliEMCALDigitizer::BookControlHists(int var)
766 Info("BookControlHists"," started ");
768 const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
770 new TH1F("hDigiN", "#EMCAL digits with fAmp > fDigitThreshold",
771 fNADCEC+1, -0.5, Double_t(fNADCEC));
772 new TH1F("HDigiSumEnergy","Sum.EMCAL energy from digi", 1000, 0.0, 200.);
773 new TH1F("hDigiAmp", "EMCAL digital amplitude", fNADCEC+1, -0.5, Double_t(fNADCEC));
774 new TH1F("hDigiEnergy","EMCAL cell energy", 2000, 0.0, 200.);
775 new TH1F("hDigiAbsId","EMCAL absId cells with fAmp > fDigitThreshold ",
776 geom->GetNCells(), 0.5, Double_t(geom->GetNCells())+0.5);
779 fHists = AliEMCALHistoUtilities::MoveHistsToList("EmcalDigiControlHists", kFALSE);
780 fHists = 0; //huh? JLK 03-Mar-2006
784 void AliEMCALDigitizer::SaveHists(const char* name, Bool_t kSingleKey, const char* opt)
786 AliEMCALHistoUtilities::SaveListOfHists(fHists, name, kSingleKey, opt);