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 "AliEMCALDigit.h"
75 #include "AliEMCALGetter.h"
76 #include "AliEMCALDigitizer.h"
77 #include "AliEMCALSDigitizer.h"
78 #include "AliEMCALGeometry.h"
79 #include "AliEMCALTick.h"
80 #include "AliEMCALJetMicroDst.h"
82 ClassImp(AliEMCALDigitizer)
85 //____________________________________________________________________________
86 AliEMCALDigitizer::AliEMCALDigitizer():AliDigitizer("",""),
93 fDefaultInit = kTRUE ;
94 fManager = 0 ; // We work in the standalong mode
95 fEventFolderName = "" ;
98 //____________________________________________________________________________
99 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName):
100 AliDigitizer("EMCAL"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
101 fInputFileNames(0), fEventNames(0), fEventFolderName(eventFolderName)
107 fDefaultInit = kFALSE ;
108 fManager = 0 ; // We work in the standalong mode
111 //____________________________________________________________________________
112 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d) : AliDigitizer(d)
116 SetName(d.GetName()) ;
117 SetTitle(d.GetTitle()) ;
118 fDigitThreshold = d.fDigitThreshold ;
119 fMeanPhotonElectron = d.fMeanPhotonElectron ;
120 fPedestal = d.fPedestal ;
122 fPinNoise = d.fPinNoise ;
123 fTimeResolution = d.fTimeResolution ;
124 fTimeThreshold = d.fTimeThreshold ;
125 fTimeSignalLength = d.fTimeSignalLength ;
126 fADCchannelEC = d.fADCchannelEC ;
127 fADCpedestalEC = d.fADCpedestalEC ;
128 fNADCEC = d.fNADCEC ;
129 fEventFolderName = d.fEventFolderName;
132 //____________________________________________________________________________
133 AliEMCALDigitizer::AliEMCALDigitizer(AliRunDigitizer * rd):
134 AliDigitizer(rd,"EMCAL"+AliConfig::Instance()->GetDigitizerTaskName()),
137 // ctor Init() is called by RunDigitizer
139 fEventFolderName = fManager->GetInputFolderName(0) ;
140 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
142 fDefaultInit = kFALSE ;
145 //____________________________________________________________________________
146 AliEMCALDigitizer::~AliEMCALDigitizer()
148 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ;
150 // Clean Digitizer from the white board
151 gime->EmcalLoader()->CleanDigitizer() ;
153 delete [] fInputFileNames ;
154 delete [] fEventNames ;
158 //____________________________________________________________________________
159 void AliEMCALDigitizer::Digitize(Int_t event)
162 // Makes the digitization of the collected summable digits
163 // for this it first creates the array of all EMCAL modules
164 // filled with noise (different for EMC, CPV and PPSD) and
165 // after that adds contributions from SDigits. This design
166 // helps to avoid scanning over the list of digits to add
167 // contribution of any new SDigit.
168 static int isTrd1Geom = -1; // -1 - mean undefined
169 static int nEMC=0; //max number of digits possible
171 //<<<<<<< AliEMCALDigitizer.cxx
173 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ;
175 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ;
177 Int_t ReadEvent = event ;
178 // fManager is data member from AliDigitizer
180 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
181 if(deb>0) Info("Digitize", "Adding event %d from input stream 0 %s %s",
182 ReadEvent, GetTitle(), fEventFolderName.Data()) ;
183 gime->Event(ReadEvent, "S") ;
184 TClonesArray * digits = gime->Digits() ;
187 const AliEMCALGeometry *geom = gime->EMCALGeometry() ;
189 TString ng(geom->GetName());
191 if(ng.Contains("SHISH") && ng.Contains("TRD1")) isTrd1Geom = 1;
193 if(isTrd1Geom == 0) nEMC = geom->GetNPhi()*geom->GetNZ();
194 else nEMC = geom->GetNCells();
195 printf(" nEMC %i (number cells in EMCAL) | %s | isTrd1Geom %i\n", nEMC, geom->GetName(), isTrd1Geom);
199 digits->Expand(nEMC) ;
201 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
202 if ( !gime->SDigitizer() )
203 gime->LoadSDigitizer();
204 AliEMCALSDigitizer * sDigitizer = gime->SDigitizer();
207 Fatal("Digitize", "SDigitizer with name %s %s not found", fEventFolderName.Data(), GetTitle() ) ;
209 //take all the inputs to add together and load the SDigits
210 TObjArray * sdigArray = new TObjArray(fInput) ;
211 sdigArray->AddAt(gime->SDigits(), 0) ;
213 for(i = 1 ; i < fInput ; i++){
214 TString tempo(fEventNames[i]) ;
216 AliEMCALGetter * gime = AliEMCALGetter::Instance(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()) ;
220 gime->Event(ReadEvent,"S");
221 sdigArray->AddAt(gime->SDigits(), i) ;
224 //Find the first tower with signal
225 Int_t nextSig = nEMC + 1 ;
226 TClonesArray * sdigits ;
227 for(i = 0 ; i < fInput ; i++){
228 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
229 if ( !sdigits->GetEntriesFast() )
231 Int_t curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(0))->GetId() ;
232 if(curNext < nextSig)
234 if(deb>0) printf("input %i : #sdigits %i \n", i, sdigits->GetEntriesFast());
236 if(deb>0) printf("FIRST tower with signal %i \n", nextSig);
238 TArrayI index(fInput) ;
239 index.Reset() ; //Set all indexes to zero
241 AliEMCALDigit * digit ;
242 AliEMCALDigit * curSDigit ;
244 TClonesArray * ticks = new TClonesArray("AliEMCALTick",1000) ;
246 //Put Noise contributiony
247 for(absID = 1; absID <= nEMC; absID++){
249 // amplitude set to zero, noise will be added later
250 new((*digits)[absID-1]) AliEMCALDigit( -1, -1, absID, 0, TimeOfNoise() ) ;
251 //look if we have to add signal?
252 digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID-1)) ;
255 //Add SDigits from all inputs
258 Float_t a = digit->GetAmp() ;
259 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
260 //Mark the beginning of the signal
261 new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime(),0, b);
262 //Mark the end of the signal
263 new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime()+fTimeSignalLength, -a, -b);
266 for(i = 0; i< fInput ; i++){ //loop over (possible) merge sources
267 if(dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
268 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
271 //May be several digits will contribute from the same input
272 while(curSDigit && (curSDigit->GetId() == absID)){
273 //Shift primary to separate primaries belonging different inputs
274 Int_t primaryoffset ;
276 primaryoffset = fManager->GetMask(i) ;
279 curSDigit->ShiftPrimary(primaryoffset) ;
281 a = curSDigit->GetAmp() ;
282 b = a /fTimeSignalLength ;
283 new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime(),0, b);
284 new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
286 *digit = *digit + *curSDigit ; //add energies
289 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
290 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
295 // add fluctuations for photo-electron creation
296 amp = sDigitizer->Calibrate(digit->GetAmp()) ; // GeV
297 amp *= static_cast<Float_t>(gRandom->Poisson(fMeanPhotonElectron)) / static_cast<Float_t>(fMeanPhotonElectron) ;
299 //calculate and set time
300 Float_t time = FrontEdgeTime(ticks) ;
301 digit->SetTime(time) ;
303 //Find next signal module
305 for(i = 0 ; i < fInput ; i++){
306 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
307 Int_t curNext = nextSig ;
308 if(sdigits->GetEntriesFast() > index[i] ){
309 curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))->GetId() ;
311 if(curNext < nextSig) nextSig = curNext ;
315 amp += TMath::Abs(gRandom->Gaus(0., fPinNoise)) ;
316 digit->SetAmp(sDigitizer->Digitize(amp)) ;
317 if(deb>=10) printf(" absID %5i amp %f nextSig %5i\n", absID, amp, nextSig);
318 } // for(absID = 1; absID <= nEMC; absID++)
323 delete sdigArray ; //We should not delete its contents
325 //remove digits below thresholds
326 for(i = 0 ; i < nEMC ; i++){
327 digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
328 Float_t threshold = fDigitThreshold ;
329 if(sDigitizer->Calibrate( digit->GetAmp() ) < threshold)
330 digits->RemoveAt(i) ;
332 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
337 Int_t ndigits = digits->GetEntriesFast() ;
338 digits->Expand(ndigits) ;
340 //Set indexes in list of digits and fill hists.
341 sv::FillH1(fHists, 0, Double_t(ndigits));
342 Float_t energy=0., esum=0.;
343 for (i = 0 ; i < ndigits ; i++) {
344 digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ;
345 digit->SetIndexInList(i) ;
346 energy = sDigitizer->Calibrate(digit->GetAmp()) ;
348 digit->SetAmp(DigitizeEnergy(energy) ) ; // for what ??
349 sv::FillH1(fHists, 2, double(digit->GetAmp()));
350 sv::FillH1(fHists, 3, double(energy));
351 sv::FillH1(fHists, 4, double(digit->GetId()));
353 sv::FillH1(fHists, 1, esum);
356 //____________________________________________________________________________
358 Int_t AliEMCALDigitizer::DigitizeEnergy(Float_t energy)
360 Int_t channel = -999;
361 channel = static_cast<Int_t> (TMath::Ceil( (energy + fADCpedestalEC)/fADCchannelEC )) ;
362 if(channel > fNADCEC )
367 //____________________________________________________________________________
368 void AliEMCALDigitizer::Exec(Option_t *option)
370 // Steering method to process digitization for events
371 // in the range from fFirstEvent to fLastEvent.
372 // This range is optionally set by SetEventRange().
373 // if fLastEvent=-1, then process events until the end.
374 // by default fLastEvent = fFirstEvent (process only one event)
376 if (!fInit) { // to prevent overwrite existing file
377 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
381 if (strstr(option,"print")) {
386 if(strstr(option,"tim"))
387 gBenchmark->Start("EMCALDigitizer");
389 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ;
391 // Post Digitizer to the white board
392 gime->PostDigitizer(this) ;
394 if (fLastEvent == -1)
395 fLastEvent = gime->MaxEvent() - 1 ;
397 fLastEvent = fFirstEvent ; // what is this ??
399 Int_t nEvents = fLastEvent - fFirstEvent + 1;
400 Int_t ievent, nfr=50;
402 fLastEvent = TMath::Min(fLastEvent, gime->MaxEvent());
403 printf("AliEMCALDigitizer::Exec : option: %s | %i -> %i events : Max events %i \n",
404 option, fFirstEvent, fLastEvent, gime->MaxEvent());
405 for (ievent = fFirstEvent; ievent < fLastEvent; ievent++) {
406 if(ievent%nfr==0 || ievent==fLastEvent-1);
407 printf(" processed event %i\n", ievent);
408 gime->Event(ievent,"S") ;
410 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
414 if(strstr(option,"deb"))
416 if(strstr(option,"table")) gObjectTable->Print();
418 //increment the total number of Digits per run
419 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
421 // gime->WriteDigitizer("OVERWRITE");
423 gime->EmcalLoader()->CleanDigitizer() ;
425 if(strstr(option,"tim")){
426 gBenchmark->Stop("EMCALDigitizer");
427 printf("Exec: took %f seconds for Digitizing %f seconds per event",
428 gBenchmark->GetCpuTime("EMCALDigitizer"), gBenchmark->GetCpuTime("EMCALDigitizer")/nEvents ) ;
432 //____________________________________________________________________________
433 Float_t AliEMCALDigitizer::FrontEdgeTime(TClonesArray * ticks)
435 // Returns the shortest time among all time ticks
437 ticks->Sort() ; //Sort in accordance with times of ticks
439 AliEMCALTick * ctick = (AliEMCALTick *) it.Next() ;
440 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
443 while((t=(AliEMCALTick*) it.Next())){
444 if(t->GetTime() < time) //This tick starts before crossing
449 time = ctick->CrossingTime(fTimeThreshold) ;
454 //____________________________________________________________________________
455 Bool_t AliEMCALDigitizer::Init()
457 // Makes all memory allocations
459 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ;
461 Error("Init", "Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;
464 TString opt("Digits") ;
465 if(gime->VersionExists(opt) ) {
466 Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
471 fLastEvent = fFirstEvent ;
474 fInput = fManager->GetNinputs() ;
478 fInputFileNames = new TString[fInput] ;
479 fEventNames = new TString[fInput] ;
480 fInputFileNames[0] = GetTitle() ;
481 fEventNames[0] = fEventFolderName.Data() ;
483 for (index = 1 ; index < fInput ; index++) {
484 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
485 TString tempo = fManager->GetInputFolderName(index) ;
486 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fManager
489 //to prevent cleaning of this object while GetEvent is called
490 gime->EmcalLoader()->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
495 //____________________________________________________________________________
496 void AliEMCALDigitizer::InitParameters()
497 { // Tune parameters - 24-nov-04
499 fMeanPhotonElectron = 3300 ; // electrons per GeV
501 if (fPinNoise == 0. )
502 Warning("InitParameters", "No noise added\n") ;
503 fDigitThreshold = fPinNoise * 3; // 3 * sigma
504 fTimeResolution = 0.3e-9 ; // 300 psc
505 fTimeSignalLength = 1.0e-9 ;
507 fADCchannelEC = 0.00305; // 200./65536 - width of one ADC channel in GeV
508 fADCpedestalEC = 0.009 ; // GeV
509 fNADCEC = (Int_t) TMath::Power(2,16) ; // number of channels in Tower ADC - 65536
511 fTimeThreshold = 0.001*10000000 ; // Means 1 MeV in terms of SDigits amplitude ??
512 // hists. for control; no hists on default
517 //__________________________________________________________________
518 void AliEMCALDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
520 // Allows to produce digits by superimposing background and signal event.
521 // It is assumed, that headers file with SIGNAL events is opened in
523 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
524 // Thus we avoid writing (changing) huge and expensive
525 // backgound files: all output will be writen into SIGNAL, i.e.
526 // opened in constructor file.
528 // One can open as many files to mix with as one needs.
529 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
532 if( strcmp(GetName(), "") == 0 )
536 Error("MixWith", "Cannot use this method under AliRunDigitizer") ;
539 // looking for file which contains AliRun
540 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
541 Error("MixWith", "File %s does not exist!", alirunFileName.Data()) ;
544 // looking for the file which contains SDigits
545 AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
546 TString fileName( gime->GetSDigitsFileName() ) ;
547 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
548 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
549 if ( (gSystem->AccessPathName(fileName)) ) {
550 Error("MixWith", "The file %s does not exist!", fileName.Data()) ;
553 // need to increase the arrays
554 TString tempo = fInputFileNames[fInput-1] ;
555 delete [] fInputFileNames ;
556 fInputFileNames = new TString[fInput+1] ;
557 fInputFileNames[fInput-1] = tempo ;
559 tempo = fEventNames[fInput-1] ;
560 delete [] fEventNames ;
561 fEventNames = new TString[fInput+1] ;
562 fEventNames[fInput-1] = tempo ;
564 fInputFileNames[fInput] = alirunFileName ;
565 fEventNames[fInput] = eventFolderName ;
569 void AliEMCALDigitizer::Print1(Option_t * option)
570 { // 19-nov-04 - just for convinience
575 //__________________________________________________________________
576 void AliEMCALDigitizer::Print()const
578 // Print Digitizer's parameters
579 printf("Print: \n------------------- %s -------------", GetName() ) ;
580 if( strcmp(fEventFolderName.Data(), "") != 0 ){
581 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
585 nStreams = GetNInputStreams() ;
590 for (index = 0 ; index < nStreams ; index++) {
591 TString tempo(fEventNames[index]) ;
593 AliEMCALGetter * gime = AliEMCALGetter::Instance(fInputFileNames[index], tempo) ;
594 TString fileName( gime->GetSDigitsFileName() ) ;
595 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
596 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
597 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
599 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ;
600 printf("\nWriting digits to %s", gime->GetDigitsFileName().Data()) ;
602 printf("\nWith following parameters:\n") ;
604 printf(" Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise) ;
605 printf(" Threshold in EMC (fDigitThreshold) = %f\n", fDigitThreshold) ;
606 printf("---------------------------------------------------\n") ;
609 printf("Print: AliEMCALDigitizer not initialized") ;
612 //__________________________________________________________________
613 void AliEMCALDigitizer::PrintDigits(Option_t * option){
615 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ;
616 TClonesArray * digits = gime->Digits() ;
617 TClonesArray * sdigits = gime->SDigits() ;
619 printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ;
620 printf("\n event %d", gAlice->GetEvNumber()) ;
622 if(strstr(option,"all")){
624 AliEMCALDigit * digit;
625 printf("\nEMC digits (with primaries):\n") ;
626 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
628 for (index = 0 ; index < digits->GetEntries() ; index++) {
629 digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
630 printf("\n%6d %8d %6.5e %4d %2d : ",
631 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
633 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
634 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
641 //__________________________________________________________________
642 Float_t AliEMCALDigitizer::TimeOfNoise(void)
644 // Calculates the time signal generated by noise
645 //PH Info("TimeOfNoise", "Change me") ;
646 return gRandom->Rndm() * 1.28E-5;
649 //__________________________________________________________________
650 void AliEMCALDigitizer::Unload()
652 // Unloads the SDigits and Digits
654 for(i = 1 ; i < fInput ; i++){
655 TString tempo(fEventNames[i]) ;
657 AliEMCALGetter * gime = AliEMCALGetter::Instance(fInputFileNames[i], tempo) ;
658 gime->EmcalLoader()->UnloadSDigits() ;
661 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ;
662 gime->EmcalLoader()->UnloadDigits() ;
665 //_________________________________________________________________________________________
666 void AliEMCALDigitizer::WriteDigits()
669 // Makes TreeD in the output file.
670 // Check if branch already exists:
671 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
672 // already existing branches.
673 // else creates branch with Digits, named "EMCAL", title "...",
674 // and branch "AliEMCALDigitizer", with the same title to keep all the parameters
675 // and names of files, from which digits are made.
677 //<<<<<<< AliEMCALDigitizer.cxx
678 static Int_t writeSdigitizer=0;
679 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ;
681 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ;
683 const TClonesArray * digits = gime->Digits() ;
684 TTree * treeD = gime->TreeD();
686 // -- create Digits branch
687 Int_t bufferSize = 32000 ;
688 TBranch * digitsBranch = treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
689 digitsBranch->SetTitle(fEventFolderName);
690 digitsBranch->Fill() ;
692 gime->WriteDigits("OVERWRITE");
693 if(writeSdigitizer==0) {
694 gime->WriteDigitizer("OVERWRITE");
702 void AliEMCALDigitizer::Browse(TBrowser* b)
704 if(fHists) b->Add(fHists);
708 TList *AliEMCALDigitizer::BookControlHists(const int var)
710 Info("BookControlHists"," started ");
712 AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName);
713 const AliEMCALGeometry *geom = gime->EMCALGeometry() ;
715 new TH1F("hDigiN", "#EMCAL digits with fAmp > fDigitThreshold",
716 fNADCEC+1, -0.5, Double_t(fNADCEC));
717 new TH1F("HDigiSumEnergy","Sum.EMCAL energy from digi", 1000, 0.0, 200.);
718 new TH1F("hDigiAmp", "EMCAL digital amplitude", fNADCEC+1, -0.5, Double_t(fNADCEC));
719 new TH1F("hDigiEnergy","EMCAL cell energy", 2000, 0.0, 200.);
720 new TH1F("hDigiAbsId","EMCAL absId cells with fAmp > fDigitThreshold ",
721 geom->GetNCells(), 0.5, Double_t(geom->GetNCells())+0.5);
723 fHists = sv::MoveHistsToList("EmcalDigiControlHists", kFALSE);
727 void AliEMCALDigitizer::SaveHists(const char* name, Bool_t kSingleKey, const char* opt)
729 sv::SaveListOfHists(fHists, name, kSingleKey, opt);