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 ---
72 #include "AliRunDigitizer.h"
73 #include "AliRunLoader.h"
74 #include "AliEMCALDigit.h"
76 #include "AliEMCALLoader.h"
77 #include "AliEMCALDigitizer.h"
78 #include "AliEMCALSDigitizer.h"
79 #include "AliEMCALGeometry.h"
80 #include "AliEMCALTick.h"
81 #include "AliEMCALJetMicroDst.h"
83 ClassImp(AliEMCALDigitizer)
86 //____________________________________________________________________________
87 AliEMCALDigitizer::AliEMCALDigitizer():AliDigitizer("",""),
94 fDefaultInit = kTRUE ;
95 fManager = 0 ; // We work in the standalong mode
96 fEventFolderName = "" ;
99 //____________________________________________________________________________
100 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName):
101 AliDigitizer("EMCAL"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
102 fInputFileNames(0), fEventNames(0), fEventFolderName(eventFolderName)
108 fDefaultInit = kFALSE ;
109 fManager = 0 ; // We work in the standalong mode
112 //____________________________________________________________________________
113 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d) : AliDigitizer(d)
117 SetName(d.GetName()) ;
118 SetTitle(d.GetTitle()) ;
119 fDigitThreshold = d.fDigitThreshold ;
120 fMeanPhotonElectron = d.fMeanPhotonElectron ;
121 fPedestal = d.fPedestal ;
123 fPinNoise = d.fPinNoise ;
124 fTimeResolution = d.fTimeResolution ;
125 fTimeThreshold = d.fTimeThreshold ;
126 fTimeSignalLength = d.fTimeSignalLength ;
127 fADCchannelEC = d.fADCchannelEC ;
128 fADCpedestalEC = d.fADCpedestalEC ;
129 fNADCEC = d.fNADCEC ;
130 fEventFolderName = d.fEventFolderName;
133 //____________________________________________________________________________
134 AliEMCALDigitizer::AliEMCALDigitizer(AliRunDigitizer * rd):
135 AliDigitizer(rd,"EMCAL"+AliConfig::Instance()->GetDigitizerTaskName()),
138 // ctor Init() is called by RunDigitizer
140 fEventFolderName = fManager->GetInputFolderName(0) ;
141 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
143 fDefaultInit = kFALSE ;
146 //____________________________________________________________________________
147 AliEMCALDigitizer::~AliEMCALDigitizer()
149 if (AliRunLoader::GetRunLoader()) {
150 AliLoader *emcalLoader=0;
151 if ((emcalLoader = AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL")))
152 emcalLoader->CleanDigitizer();
155 AliDebug(1," no runloader present");
156 delete [] fInputFileNames ;
157 delete [] fEventNames ;
160 //____________________________________________________________________________
161 void AliEMCALDigitizer::Digitize(Int_t event)
164 // Makes the digitization of the collected summable digits
165 // for this it first creates the array of all EMCAL modules
166 // filled with noise (different for EMC, CPV and PPSD) and
167 // after that adds contributions from SDigits. This design
168 // helps to avoid scanning over the list of digits to add
169 // contribution of any new SDigit.
170 static int isTrd1Geom = -1; // -1 - mean undefined
171 static int nEMC=0; //max number of digits possible
173 AliRunLoader *rl = AliRunLoader::GetRunLoader();
174 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
175 Int_t ReadEvent = event ;
176 // fManager is data member from AliDigitizer
178 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
179 AliDebug(1,Form("Adding event %d from input stream 0 %s %s",
180 ReadEvent, GetTitle(), fEventFolderName.Data())) ;
181 rl->GetEvent(ReadEvent);
183 TClonesArray * digits = emcalLoader->Digits() ;
186 const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
188 TString ng(geom->GetName());
190 if(ng.Contains("SHISH") && ng.Contains("TRD1")) isTrd1Geom = 1;
192 if(isTrd1Geom == 0) nEMC = geom->GetNPhi()*geom->GetNZ();
193 else nEMC = geom->GetNCells();
194 AliDebug(1,Form("nEMC %i (number cells in EMCAL) | %s | isTrd1Geom %i\n", nEMC, geom->GetName(), isTrd1Geom));
198 digits->Expand(nEMC) ;
200 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
201 if ( !emcalLoader->SDigitizer() )
202 emcalLoader->LoadSDigitizer();
203 AliEMCALSDigitizer * sDigitizer = dynamic_cast<AliEMCALSDigitizer*>(emcalLoader->SDigitizer());
206 Fatal("Digitize", "SDigitizer with name %s %s not found", fEventFolderName.Data(), GetTitle() ) ;
208 //take all the inputs to add together and load the SDigits
209 TObjArray * sdigArray = new TObjArray(fInput) ;
210 sdigArray->AddAt(emcalLoader->SDigits(), 0) ;
213 for(i = 1 ; i < fInput ; i++){
214 TString tempo(fEventNames[i]) ;
216 rl = AliRunLoader::Open(fInputFileNames[i], tempo) ;
218 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
219 Info("Digitize", "Adding event %d from input stream %d %s %s", ReadEvent, i, fInputFileNames[i].Data(), tempo.Data()) ;
221 rl->GetEvent(ReadEvent);
222 sdigArray->AddAt(emcalLoader->SDigits(), i) ;
225 //Find the first tower with signal
226 Int_t nextSig = nEMC + 1 ;
227 TClonesArray * sdigits ;
228 for(i = 0 ; i < fInput ; i++){
229 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
230 if ( !sdigits->GetEntriesFast() )
232 Int_t curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(0))->GetId() ;
233 if(curNext < nextSig)
235 AliDebug(1,Form("input %i : #sdigits %i \n",
236 i, sdigits->GetEntriesFast()));
238 AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
240 TArrayI index(fInput) ;
241 index.Reset() ; //Set all indexes to zero
243 AliEMCALDigit * digit ;
244 AliEMCALDigit * curSDigit ;
246 TClonesArray * ticks = new TClonesArray("AliEMCALTick",1000) ;
248 //Put Noise contribution
249 for(absID = 1; absID <= nEMC; absID++){
251 // amplitude set to zero, noise will be added later
252 new((*digits)[absID-1]) AliEMCALDigit( -1, -1, absID, 0, TimeOfNoise() ) ;
253 //look if we have to add signal?
254 digit = dynamic_cast<AliEMCALDigit *>(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++]) AliEMCALTick(digit->GetTime(),0, b);
264 //Mark the end of the signal
265 new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime()+fTimeSignalLength, -a, -b);
268 for(i = 0; i< fInput ; i++){ //loop over (possible) merge sources
269 if(dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
270 curSDigit = dynamic_cast<AliEMCALDigit*>(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) ;
281 curSDigit->ShiftPrimary(primaryoffset) ;
283 a = curSDigit->GetAmp() ;
284 b = a /fTimeSignalLength ;
285 new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime(),0, b);
286 new((*ticks)[contrib++]) AliEMCALTick(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<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
297 // add fluctuations for photo-electron creation
298 amp = sDigitizer->Calibrate(digit->GetAmp()) ; // GeV
299 amp *= static_cast<Float_t>(gRandom->Poisson(fMeanPhotonElectron)) / static_cast<Float_t>(fMeanPhotonElectron) ;
301 //calculate and set time
302 Float_t time = FrontEdgeTime(ticks) ;
303 digit->SetTime(time) ;
305 //Find next signal module
307 for(i = 0 ; i < fInput ; i++){
308 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
309 Int_t curNext = nextSig ;
310 if(sdigits->GetEntriesFast() > index[i] ){
311 curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))->GetId() ;
313 if(curNext < nextSig) nextSig = curNext ;
317 amp += TMath::Abs(gRandom->Gaus(0., fPinNoise)) ;
318 digit->SetAmp(sDigitizer->Digitize(amp)) ;
319 AliDebug(10,Form(" absID %5i amp %f nextSig %5i\n",
320 absID, amp, nextSig));
321 } // for(absID = 1; absID <= nEMC; absID++)
326 delete sdigArray ; //We should not delete its contents
328 //remove digits below thresholds
329 for(i = 0 ; i < nEMC ; i++){
330 digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
331 Float_t threshold = fDigitThreshold ;
332 if(sDigitizer->Calibrate( digit->GetAmp() ) < threshold)
333 digits->RemoveAt(i) ;
335 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
340 Int_t ndigits = digits->GetEntriesFast() ;
341 digits->Expand(ndigits) ;
343 //Set indexes in list of digits and fill hists.
344 sv::FillH1(fHists, 0, Double_t(ndigits));
345 Float_t energy=0., esum=0.;
346 for (i = 0 ; i < ndigits ; i++) {
347 digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ;
348 digit->SetIndexInList(i) ;
349 energy = sDigitizer->Calibrate(digit->GetAmp()) ;
351 digit->SetAmp(DigitizeEnergy(energy) ) ; // for what ??
352 sv::FillH1(fHists, 2, double(digit->GetAmp()));
353 sv::FillH1(fHists, 3, double(energy));
354 sv::FillH1(fHists, 4, double(digit->GetId()));
356 sv::FillH1(fHists, 1, esum);
359 //____________________________________________________________________________
361 Int_t AliEMCALDigitizer::DigitizeEnergy(Float_t energy)
363 Int_t channel = -999;
364 channel = static_cast<Int_t> (TMath::Ceil( (energy + fADCpedestalEC)/fADCchannelEC )) ;
365 if(channel > fNADCEC )
370 //____________________________________________________________________________
371 void AliEMCALDigitizer::Exec(Option_t *option)
373 // Steering method to process digitization for events
374 // in the range from fFirstEvent to fLastEvent.
375 // This range is optionally set by SetEventRange().
376 // if fLastEvent=-1, then process events until the end.
377 // by default fLastEvent = fFirstEvent (process only one event)
379 if (!fInit) { // to prevent overwrite existing file
380 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
384 if (strstr(option,"print")) {
389 if(strstr(option,"tim"))
390 gBenchmark->Start("EMCALDigitizer");
392 AliRunLoader *rl = AliRunLoader::GetRunLoader();
393 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
395 // Post Digitizer to the white board
396 emcalLoader->PostDigitizer(this) ;
398 if (fLastEvent == -1)
399 fLastEvent = rl->GetNumberOfEvents() - 1 ;
401 fLastEvent = fFirstEvent ; // what is this ??
403 Int_t nEvents = fLastEvent - fFirstEvent + 1;
406 rl->LoadSDigits("EMCAL");
407 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
408 //gime->Event(ievent,"S") ;
409 rl->GetEvent(ievent);
411 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
415 if(strstr(option,"deb"))
417 if(strstr(option,"table")) gObjectTable->Print();
419 //increment the total number of Digits per run
420 fDigitsInRun += emcalLoader->Digits()->GetEntriesFast() ;
422 // gime->WriteDigitizer("OVERWRITE");
424 emcalLoader->CleanDigitizer() ;
426 if(strstr(option,"tim")){
427 gBenchmark->Stop("EMCALDigitizer");
428 printf("Exec: took %f seconds for Digitizing %f seconds per event",
429 gBenchmark->GetCpuTime("EMCALDigitizer"), gBenchmark->GetCpuTime("EMCALDigitizer")/nEvents ) ;
433 //____________________________________________________________________________
434 Float_t AliEMCALDigitizer::FrontEdgeTime(TClonesArray * ticks)
436 // Returns the shortest time among all time ticks
438 ticks->Sort() ; //Sort in accordance with times of ticks
440 AliEMCALTick * ctick = (AliEMCALTick *) it.Next() ;
441 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
444 while((t=(AliEMCALTick*) it.Next())){
445 if(t->GetTime() < time) //This tick starts before crossing
450 time = ctick->CrossingTime(fTimeThreshold) ;
455 //____________________________________________________________________________
456 Bool_t AliEMCALDigitizer::Init()
458 // Makes all memory allocations
460 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
462 if ( emcalLoader == 0 ) {
463 Fatal("Init", "Could not obtain the AliEMCALLoader");
468 fLastEvent = fFirstEvent ;
471 fInput = fManager->GetNinputs() ;
475 fInputFileNames = new TString[fInput] ;
476 fEventNames = new TString[fInput] ;
477 fInputFileNames[0] = GetTitle() ;
478 fEventNames[0] = fEventFolderName.Data() ;
480 for (index = 1 ; index < fInput ; index++) {
481 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
482 TString tempo = fManager->GetInputFolderName(index) ;
483 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fManager
486 //to prevent cleaning of this object while GetEvent is called
487 emcalLoader->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
492 //____________________________________________________________________________
493 void AliEMCALDigitizer::InitParameters()
494 { // Tune parameters - 24-nov-04
496 fMeanPhotonElectron = 3300 ; // electrons per GeV
498 if (fPinNoise == 0. )
499 Warning("InitParameters", "No noise added\n") ;
500 fDigitThreshold = fPinNoise * 3; // 3 * sigma
501 fTimeResolution = 0.3e-9 ; // 300 psc
502 fTimeSignalLength = 1.0e-9 ;
504 fADCchannelEC = 0.00305; // 200./65536 - width of one ADC channel in GeV
505 fADCpedestalEC = 0.009 ; // GeV
506 fNADCEC = (Int_t) TMath::Power(2,16) ; // number of channels in Tower ADC - 65536
508 fTimeThreshold = 0.001*10000000 ; // Means 1 MeV in terms of SDigits amplitude ??
509 // hists. for control; no hists on default
514 //__________________________________________________________________
515 void AliEMCALDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
517 // Allows to produce digits by superimposing background and signal event.
518 // It is assumed, that headers file with SIGNAL events is opened in
520 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
521 // Thus we avoid writing (changing) huge and expensive
522 // backgound files: all output will be writen into SIGNAL, i.e.
523 // opened in constructor file.
525 // One can open as many files to mix with as one needs.
526 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
529 if( strcmp(GetName(), "") == 0 )
533 Error("MixWith", "Cannot use this method under AliRunDigitizer") ;
536 // looking for file which contains AliRun
537 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
538 Error("MixWith", "File %s does not exist!", alirunFileName.Data()) ;
541 // looking for the file which contains SDigits
542 //AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
543 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
544 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
545 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
546 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
547 if ( (gSystem->AccessPathName(fileName)) ) {
548 Error("MixWith", "The file %s does not exist!", fileName.Data()) ;
551 // need to increase the arrays
552 TString tempo = fInputFileNames[fInput-1] ;
553 delete [] fInputFileNames ;
554 fInputFileNames = new TString[fInput+1] ;
555 fInputFileNames[fInput-1] = tempo ;
557 tempo = fEventNames[fInput-1] ;
558 delete [] fEventNames ;
559 fEventNames = new TString[fInput+1] ;
560 fEventNames[fInput-1] = tempo ;
562 fInputFileNames[fInput] = alirunFileName ;
563 fEventNames[fInput] = eventFolderName ;
567 void AliEMCALDigitizer::Print1(Option_t * option)
568 { // 19-nov-04 - just for convinience
573 //__________________________________________________________________
574 void AliEMCALDigitizer::Print(Option_t*)const
576 // Print Digitizer's parameters
577 printf("Print: \n------------------- %s -------------", GetName() ) ;
578 if( strcmp(fEventFolderName.Data(), "") != 0 ){
579 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
583 nStreams = GetNInputStreams() ;
591 for (index = 0 ; index < nStreams ; index++) {
592 TString tempo(fEventNames[index]) ;
594 if ((rl = AliRunLoader::GetRunLoader(tempo)) == 0)
595 rl = AliRunLoader::Open(fInputFileNames[index], tempo) ;
596 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
597 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
598 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
599 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
600 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
603 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
605 printf("\nWriting digits to %s", emcalLoader->GetDigitsFileName().Data()) ;
607 printf("\nWith following parameters:\n") ;
609 printf(" Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise) ;
610 printf(" Threshold in EMC (fDigitThreshold) = %f\n", fDigitThreshold) ;
611 printf("---------------------------------------------------\n") ;
614 printf("Print: AliEMCALDigitizer not initialized") ;
617 //__________________________________________________________________
618 void AliEMCALDigitizer::PrintDigits(Option_t * option){
619 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
620 TClonesArray * digits = emcalLoader->Digits() ;
621 TClonesArray * sdigits = emcalLoader->SDigits() ;
623 printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ;
624 printf("\n event %d", emcalLoader->GetRunLoader()->GetEventNumber());
626 if(strstr(option,"all")){
628 AliEMCALDigit * digit;
629 printf("\nEMC digits (with primaries):\n") ;
630 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
632 for (index = 0 ; index < digits->GetEntries() ; index++) {
633 digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
634 printf("\n%6d %8d %6.5e %4d %2d : ",
635 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
637 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
638 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
645 //__________________________________________________________________
646 Float_t AliEMCALDigitizer::TimeOfNoise(void)
648 // Calculates the time signal generated by noise
649 //PH Info("TimeOfNoise", "Change me") ;
650 return gRandom->Rndm() * 1.28E-5;
653 //__________________________________________________________________
654 void AliEMCALDigitizer::Unload()
656 // Unloads the SDigits and Digits
660 for(i = 1 ; i < fInput ; i++){
661 TString tempo(fEventNames[i]) ;
663 if ((rl = AliRunLoader::GetRunLoader(tempo)))
664 rl->GetDetectorLoader("EMCAL")->UnloadSDigits() ;
666 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
667 emcalLoader->UnloadDigits() ;
670 //_________________________________________________________________________________________
671 void AliEMCALDigitizer::WriteDigits()
674 // Makes TreeD in the output file.
675 // Check if branch already exists:
676 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
677 // already existing branches.
678 // else creates branch with Digits, named "EMCAL", title "...",
679 // and branch "AliEMCALDigitizer", with the same title to keep all the parameters
680 // and names of files, from which digits are made.
682 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
684 const TClonesArray * digits = emcalLoader->Digits() ;
685 TTree * treeD = emcalLoader->TreeD();
687 emcalLoader->MakeDigitsContainer();
688 treeD = emcalLoader->TreeD();
691 // -- create Digits branch
692 Int_t bufferSize = 32000 ;
693 TBranch * digitsBranch = 0;
694 if ((digitsBranch = treeD->GetBranch("EMCAL")))
695 digitsBranch->SetAddress(&digits);
697 treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
698 //digitsBranch->SetTitle(fEventFolderName);
701 emcalLoader->WriteDigits("OVERWRITE");
702 emcalLoader->WriteDigitizer("OVERWRITE");
708 void AliEMCALDigitizer::Browse(TBrowser* b)
710 if(fHists) b->Add(fHists);
714 TList *AliEMCALDigitizer::BookControlHists(int var)
716 Info("BookControlHists"," started ");
718 const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
720 new TH1F("hDigiN", "#EMCAL digits with fAmp > fDigitThreshold",
721 fNADCEC+1, -0.5, Double_t(fNADCEC));
722 new TH1F("HDigiSumEnergy","Sum.EMCAL energy from digi", 1000, 0.0, 200.);
723 new TH1F("hDigiAmp", "EMCAL digital amplitude", fNADCEC+1, -0.5, Double_t(fNADCEC));
724 new TH1F("hDigiEnergy","EMCAL cell energy", 2000, 0.0, 200.);
725 new TH1F("hDigiAbsId","EMCAL absId cells with fAmp > fDigitThreshold ",
726 geom->GetNCells(), 0.5, Double_t(geom->GetNCells())+0.5);
728 fHists = sv::MoveHistsToList("EmcalDigiControlHists", kFALSE);
733 void AliEMCALDigitizer::SaveHists(const char* name, Bool_t kSingleKey, const char* opt)
735 sv::SaveListOfHists(fHists, name, kSingleKey, opt);