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 from simulated data
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>
67 #include <TObjectTable.h>
72 // --- AliRoot header files ---
75 #include "AliRunDigitizer.h"
76 #include "AliRunLoader.h"
77 #include "AliCDBManager.h"
78 #include "AliCDBEntry.h"
79 #include "AliEMCALDigit.h"
81 #include "AliEMCALLoader.h"
82 #include "AliEMCALDigitizer.h"
83 #include "AliEMCALSDigitizer.h"
84 #include "AliEMCALGeometry.h"
85 #include "AliEMCALTick.h"
86 #include "AliEMCALCalibData.h"
87 #include "AliEMCALSimParam.h"
88 #include "AliEMCALRawDigit.h"
92 Double_t HeavisideTheta(Double_t x)
96 if (x > 0.) signal = 1.;
101 Double_t AnalogFastORFunction(Double_t *x, Double_t *par)
103 Double_t v0 = par[0];
104 Double_t t0 = par[1];
105 Double_t tr = par[2];
108 Double_t C1 = 33e-12;
110 Double_t C2 = 22e-12;
114 return (((0.8*(-((TMath::Power(C1,2)*C2*TMath::Power(TMath::E(),(-t + t0)/(C1*R1))*
115 TMath::Power(R1,2)*R2)/(C1*R1 - C2*R2)) +
116 C1*C2*R1*R2*(1 - (C2*TMath::Power(TMath::E(),(-t + t0)/(C2*R2))*R2)/(-(C1*R1) + C2*R2)))*v0*
117 HeavisideTheta(t - t0))/tr
118 - (0.8*(C1*C2*R1*R2 -
119 (TMath::Power(C1,2)*C2*TMath::Power(TMath::E(),(-1.*t + t0 + 1.25*tr)/(C1*R1))*
120 TMath::Power(R1,2)*R2)/(C1*R1 - C2*R2) +
121 (C1*TMath::Power(C2,2)*TMath::Power(TMath::E(),(-1.*t + t0 + 1.25*tr)/(C2*R2))*
122 R1*TMath::Power(R2,2))/(C1*R1 - C2*R2))*v0*
123 HeavisideTheta(t - t0 - 1.25*tr))/tr)/(C2*R1));
127 ClassImp(AliEMCALDigitizer)
130 //____________________________________________________________________________
131 AliEMCALDigitizer::AliEMCALDigitizer()
132 : AliDigitizer("",""),
137 fInputFileNames(0x0),
140 fMeanPhotonElectron(0),
141 // fPedestal(0), //Not used, remove?
142 // fSlope(0), //Not used, remove?
145 // fTimeThreshold(0), //Not used, remove?
146 // fTimeSignalLength(0), //Not used, remove?
150 fEventFolderName(""),
157 fManager = 0 ; // We work in the standalone mode
160 //____________________________________________________________________________
161 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName)
162 : AliDigitizer("EMCAL"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
163 fDefaultInit(kFALSE),
170 fMeanPhotonElectron(0),
171 // fPedestal(0),//Not used, remove?
172 // fSlope(0), //Not used, remove?
175 // fTimeThreshold(0), //Not used, remove?
176 // fTimeSignalLength(0), //Not used, remove?
180 fEventFolderName(eventFolderName),
188 fManager = 0 ; // We work in the standalone mode
191 //____________________________________________________________________________
192 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d)
193 : AliDigitizer(d.GetName(),d.GetTitle()),
194 fDefaultInit(d.fDefaultInit),
195 fDigitsInRun(d.fDigitsInRun),
198 fInputFileNames(d.fInputFileNames),
199 fEventNames(d.fEventNames),
200 fDigitThreshold(d.fDigitThreshold),
201 fMeanPhotonElectron(d.fMeanPhotonElectron),
202 // fPedestal(d.fPedestal), //Not used, remove?
203 // fSlope(d.fSlope), //Not used, remove?
204 fPinNoise(d.fPinNoise),
205 fTimeResolution(d.fTimeResolution),
206 // fTimeThreshold(d.fTimeThreshold), //Not used, remove?
207 // fTimeSignalLength(d.fTimeSignalLength), //Not used, remove?
208 fADCchannelEC(d.fADCchannelEC),
209 fADCpedestalEC(d.fADCpedestalEC),
211 fEventFolderName(d.fEventFolderName),
212 fFirstEvent(d.fFirstEvent),
213 fLastEvent(d.fLastEvent),
214 fCalibData(d.fCalibData)
219 //____________________________________________________________________________
220 AliEMCALDigitizer::AliEMCALDigitizer(AliRunDigitizer * rd)
221 : AliDigitizer(rd,"EMCAL"+AliConfig::Instance()->GetDigitizerTaskName()),
222 fDefaultInit(kFALSE),
229 fMeanPhotonElectron(0),
230 // fPedestal(0), //Not used, remove?
231 // fSlope(0.), //Not used, remove?
234 // fTimeThreshold(0), //Not used, remove?
235 // fTimeSignalLength(0), //Not used, remove?
244 // ctor Init() is called by RunDigitizer
246 fEventFolderName = fManager->GetInputFolderName(0) ;
247 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
251 //____________________________________________________________________________
252 AliEMCALDigitizer::~AliEMCALDigitizer()
255 if (AliRunLoader::Instance()) {
256 AliLoader *emcalLoader=0;
257 if ((emcalLoader = AliRunLoader::Instance()->GetDetectorLoader("EMCAL")))
258 emcalLoader->CleanDigitizer();
261 AliDebug(1," no runloader present");
262 delete [] fInputFileNames ;
263 delete [] fEventNames ;
267 //____________________________________________________________________________
268 void AliEMCALDigitizer::Digitize(Int_t event)
271 // Makes the digitization of the collected summable digits
272 // for this it first creates the array of all EMCAL modules
273 // filled with noise and after that adds contributions from
274 // SDigits. This design helps to avoid scanning over the
275 // list of digits to add contribution of any new SDigit.
278 // Note that SDigit energy info is stored as an amplitude, so we
279 // must call the Calibrate() method of the SDigitizer to convert it
280 // back to an energy in GeV before adding it to the Digit
282 static int nEMC=0; //max number of digits possible
284 AliRunLoader *rl = AliRunLoader::Instance();
285 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
286 Int_t readEvent = event ;
287 // fManager is data member from AliDigitizer
289 readEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
290 AliDebug(1,Form("Adding event %d from input stream 0 %s %s",
291 readEvent, GetTitle(), fEventFolderName.Data())) ;
292 rl->GetEvent(readEvent);
294 TClonesArray * digits = emcalLoader->Digits() ;
295 digits->Delete() ; //correct way to clear array when memory is
296 //allocated by objects stored in array
299 AliEMCALGeometry *geom = 0;
300 if (rl->GetAliRun()) {
301 AliEMCAL * emcal = (AliEMCAL*)rl->GetAliRun()->GetDetector("EMCAL");
302 geom = emcal->GetGeometry();
305 AliFatal("Could not get AliRun from runLoader");
307 nEMC = geom->GetNCells();
308 AliDebug(1,Form("nEMC %i (number cells in EMCAL) | %s \n", nEMC, geom->GetName()));
312 digits->Expand(nEMC) ;
314 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
315 if ( !emcalLoader->SDigitizer() )
316 emcalLoader->LoadSDigitizer();
317 AliEMCALSDigitizer * sDigitizer = dynamic_cast<AliEMCALSDigitizer*>(emcalLoader->SDigitizer());
320 Fatal("Digitize", "SDigitizer with name %s %s not found", fEventFolderName.Data(), GetTitle() ) ;
322 //take all the inputs to add together and load the SDigits
323 TObjArray * sdigArray = new TObjArray(fInput) ;
324 sdigArray->AddAt(emcalLoader->SDigits(), 0) ;
327 for(i = 1 ; i < fInput ; i++){
328 TString tempo(fEventNames[i]) ;
331 AliRunLoader * rl2 = AliRunLoader::GetRunLoader(tempo) ;
334 rl2 = AliRunLoader::Open(fInputFileNames[i], tempo) ;
337 readEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
338 Info("Digitize", "Adding event %d from input stream %d %s %s", readEvent, i, fInputFileNames[i].Data(), tempo.Data()) ;
340 rl2->GetEvent(readEvent);
341 AliEMCALLoader *emcalLoader2 = dynamic_cast<AliEMCALLoader*>(rl2->GetDetectorLoader("EMCAL"));
342 sdigArray->AddAt(emcalLoader2->SDigits(), i) ;
345 //Find the first tower with signal
346 Int_t nextSig = nEMC + 1 ;
347 TClonesArray * sdigits ;
348 for(i = 0 ; i < fInput ; i++){
349 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
350 if ( !sdigits->GetEntriesFast() )
352 Int_t curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(0))->GetId() ;
353 if(curNext < nextSig)
355 AliDebug(1,Form("input %i : #sdigits %i \n",
356 i, sdigits->GetEntriesFast()));
358 AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
360 TArrayI index(fInput) ;
361 index.Reset() ; //Set all indexes to zero
363 AliEMCALDigit * digit ;
364 AliEMCALDigit * curSDigit ;
366 // TClonesArray * ticks = new TClonesArray("AliEMCALTick",1000) ;
368 //Put Noise contribution
369 for(absID = 0; absID < nEMC; absID++){ // Nov 30, 2006 by PAI; was from 1 to nEMC
371 // amplitude set to zero, noise will be added later
372 new((*digits)[absID]) AliEMCALDigit( -1, -1, absID, 0, TimeOfNoise() ); // absID-1->absID
373 //look if we have to add signal?
374 digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID)); // absID-1->absID
377 //Add SDigits from all inputs
379 //Int_t contrib = 0 ;
381 //Follow PHOS and comment out this timing model til a better one
382 //can be developed - JLK 28-Apr-2008
384 //Float_t a = digit->GetAmp() ;
385 //Float_t b = TMath::Abs( a /fTimeSignalLength) ;
386 //Mark the beginning of the signal
387 //new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime(),0, b);
388 //Mark the end of the signal
389 //new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime()+fTimeSignalLength, -a, -b);
391 // Calculate time as time of the largest digit
392 Float_t time = digit->GetTime() ;
393 Float_t aTime= digit->GetAmp() ;
396 for(i = 0; i< fInput ; i++){ //loop over (possible) merge sources
397 if(dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
398 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
401 //May be several digits will contribute from the same input
402 while(curSDigit && (curSDigit->GetId() == absID)){
403 //Shift primary to separate primaries belonging different inputs
404 Int_t primaryoffset ;
406 primaryoffset = fManager->GetMask(i) ;
409 curSDigit->ShiftPrimary(primaryoffset) ;
411 //Remove old timing model - JLK 28-April-2008
412 //a = curSDigit->GetAmp() ;
413 //b = a /fTimeSignalLength ;
414 //new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime(),0, b);
415 //new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
416 if(curSDigit->GetAmp()>aTime) {
417 aTime = curSDigit->GetAmp();
418 time = curSDigit->GetTime();
421 *digit = *digit + *curSDigit ; //adds amplitudes of each digit
424 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
425 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
430 //Here we convert the summed amplitude to an energy in GeV
431 energy = sDigitizer->Calibrate(digit->GetAmp()) ; // GeV
432 // add fluctuations for photo-electron creation
433 energy *= static_cast<Float_t>(gRandom->Poisson(fMeanPhotonElectron)) / static_cast<Float_t>(fMeanPhotonElectron) ;
435 //calculate and set time
436 //New timing model needed - JLK 28-April-2008
437 //Float_t time = FrontEdgeTime(ticks) ;
438 digit->SetTime(time) ;
440 //Find next signal module
442 for(i = 0 ; i < fInput ; i++){
443 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
444 Int_t curNext = nextSig ;
445 if(sdigits->GetEntriesFast() > index[i] ){
446 curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))->GetId() ;
448 if(curNext < nextSig) nextSig = curNext ;
452 energy += TMath::Abs(gRandom->Gaus(0., fPinNoise)) ;
454 //Now digitize the energy according to the sDigitizer method,
455 //which merely converts the energy to an integer. Later we will
456 //check that the stored value matches our allowed dynamic ranges
457 digit->SetAmp(sDigitizer->Digitize(energy)) ;
458 AliDebug(10,Form(" absID %5i energy %f nextSig %5i\n",
459 absID, energy, nextSig));
460 } // for(absID = 0; absID < nEMC; absID++)
465 //JLK is it better to call Clear() here?
466 delete sdigArray ; //We should not delete its contents
468 //remove digits below thresholds
469 // until 10-02-2010 remove digits with energy smaller than fDigitThreshold 3*fPinNoise
470 // now, remove digits with Digitized ADC smaller than fDigitThreshold = 3
472 for(i = 0 ; i < nEMC ; i++){
473 digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
474 //First get the energy in GeV units.
475 energy = sDigitizer->Calibrate(digit->GetAmp()) ;
476 //Then digitize using the calibration constants of the ocdb
477 Int_t ampADC = DigitizeEnergy(energy, digit->GetId()) ;
478 //if(ampADC>2)printf("Digit energy %f, amp %d, amp cal %d, threshold %d\n",energy,digit->GetAmp(),ampADC,fDigitThreshold);
479 if(ampADC < fDigitThreshold)
480 digits->RemoveAt(i) ;
482 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
487 Int_t ndigits = digits->GetEntriesFast() ;
490 //After we have done the summing and digitizing to create the
491 //digits, now we want to calibrate the resulting amplitude to match
492 //the dynamic range of our real data.
493 for (i = 0 ; i < ndigits ; i++) {
494 digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ;
495 digit->SetIndexInList(i) ;
496 energy = sDigitizer->Calibrate(digit->GetAmp()) ;
497 digit->SetAmp(DigitizeEnergy(energy, digit->GetId()) ) ;
502 // //_____________________________________________________________________
503 Int_t AliEMCALDigitizer::DigitizeEnergy(Float_t energy, Int_t AbsId)
506 // Returns digitized value of the energy in a cell absId
507 // using the calibration constants stored in the OCDB
508 // or default values if no CalibData object is found.
509 // This effectively converts everything to match the dynamic range
510 // of the real data we will collect
513 const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
516 AliFatal("Did not get geometry from EMCALLoader");
524 Int_t channel = -999;
526 Bool_t bCell = geom->GetCellIndex(AbsId, iSupMod, nModule, nIphi, nIeta) ;
528 Error("DigitizeEnergy","Wrong cell id number : AbsId %i ", AbsId) ;
529 geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
532 fADCpedestalEC = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
533 fADCchannelEC = fCalibData->GetADCchannel(iSupMod,ieta,iphi);
536 channel = static_cast<Int_t> (TMath::Ceil( (energy + fADCpedestalEC)/fADCchannelEC )) ;
538 if(channel > fNADCEC )
544 //____________________________________________________________________________
545 void AliEMCALDigitizer::Exec(Option_t *option)
547 // Steering method to process digitization for events
548 // in the range from fFirstEvent to fLastEvent.
549 // This range is optionally set by SetEventRange().
550 // if fLastEvent=-1, then process events until the end.
551 // by default fLastEvent = fFirstEvent (process only one event)
553 if (!fInit) { // to prevent overwrite existing file
554 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
558 if (strstr(option,"print")) {
564 if(strstr(option,"tim"))
565 gBenchmark->Start("EMCALDigitizer");
567 AliRunLoader *rl = AliRunLoader::Instance();
568 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
570 // Post Digitizer to the white board
571 emcalLoader->PostDigitizer(this) ;
573 if (fLastEvent == -1) {
574 fLastEvent = rl->GetNumberOfEvents() - 1 ;
577 fLastEvent = fFirstEvent ; // what is this ??
579 Int_t nEvents = fLastEvent - fFirstEvent + 1;
582 TClonesArray* digitsTRG = new TClonesArray("AliEMCALRawDigit", 32 * 96);
583 TClonesArray* digitsTMP = new TClonesArray("AliEMCALDigit", 32 * 96);
584 rl->LoadSDigits("EMCAL");
585 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
587 rl->GetEvent(ievent);
589 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
594 //-------------------------------------
595 Digits2FastOR(digitsTMP, digitsTRG);
597 WriteDigits(digitsTRG);
599 emcalLoader->WriteDigits( "OVERWRITE");
600 emcalLoader->WriteDigitizer("OVERWRITE");
606 //-------------------------------------
608 if(strstr(option,"deb"))
610 if(strstr(option,"table")) gObjectTable->Print();
612 //increment the total number of Digits per run
613 fDigitsInRun += emcalLoader->Digits()->GetEntriesFast() ;
616 emcalLoader->CleanDigitizer() ;
618 if(strstr(option,"tim")){
619 gBenchmark->Stop("EMCALDigitizer");
620 AliInfo(Form("Exec: took %f seconds for Digitizing %f seconds per event",
621 gBenchmark->GetCpuTime("EMCALDigitizer"), gBenchmark->GetCpuTime("EMCALDigitizer")/nEvents )) ;
625 //____________________________________________________________________________
626 void AliEMCALDigitizer::Digits2FastOR(TClonesArray* digitsTMP, TClonesArray* digitsTRG)
628 // FEE digits afterburner to produce TRG digits
629 // we are only interested in the FEE digit deposited energy
630 // to be converted later into a voltage value
632 // push the FEE digit to its associated FastOR (numbered from 0:95)
633 // TRU is in charge of summing module digits
635 AliRunLoader *runLoader = AliRunLoader::Instance();
637 AliRun* run = runLoader->GetAliRun();
639 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
641 AliEMCALGeometry* geom = dynamic_cast<AliEMCAL*>(run->GetDetector("EMCAL"))->GetGeometry();
643 // build FOR from simulated digits
644 // and xfer to the corresponding TRU input (mapping)
646 TClonesArray* digits = emcalLoader->Digits();
648 TIter NextDigit(digits);
649 while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
651 Int_t id = digit->GetId();
653 Int_t iSupMod, nModule, nIphi, nIeta, iphi, ieta, iphim, ietam;
655 geom->GetCellIndex( id, iSupMod, nModule, nIphi, nIeta );
656 geom->GetModulePhiEtaIndexInSModule( iSupMod, nModule, iphim, ietam );
657 geom->GetCellPhiEtaIndexInSModule( iSupMod, nModule, nIphi, nIeta, iphi, ieta);
659 // identify to which TRU this FEE digit belong
660 //Int_t itru = (iSupMod < 11) ? iphim / 4 + 3 * iSupMod : 31;
663 itru = (iSupMod % 2) ? (2 - int(iphim / 4)) + 3 * iSupMod : iphim / 4 + 3 * iSupMod;
669 // FIXME: bad numbering solution to deal w/ the last 2 SM which have only 1 TRU each
670 // using the AliEMCALGeometry official numbering
671 // only 1 TRU/SM in SM 10 & SM 11
674 if ((itru == 31 && iphim < 2) || (itru == 30 && iphim > 5)) continue;
676 // to be compliant with %4 per TRU
677 if (itru == 31) iphim -= 2;
680 Bool_t isOK = geom->GetAbsFastORIndexFromPositionInTRU(itru, ietam, iphim % 4, trgid);
682 AliDebug(2,Form("trigger digit id: %d itru: %d isOK: %d\n",trgid,itru,isOK));
686 AliEMCALDigit* d = static_cast<AliEMCALDigit*>(digitsTMP->At(trgid));
690 new((*digitsTMP)[trgid]) AliEMCALDigit(*digit);
691 d = (AliEMCALDigit*)digitsTMP->At(trgid);
702 Int_t *timeSamples = new Int_t[nSamples];
704 NextDigit = TIter(digitsTMP);
705 while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
709 Int_t id = digit->GetId();
710 Float_t time = digit->GetTime();
712 Double_t depositedEnergy = 0.;
713 for (Int_t j = 1; j <= digit->GetNprimary(); j++) depositedEnergy += digit->GetDEPrimary(j);
715 // FIXME: Check digit time!
718 DigitalFastOR(time, depositedEnergy, timeSamples, nSamples);
720 for (Int_t j=0;j<nSamples;j++)
722 timeSamples[j] = ((j << 12) & 0xFF000) | (timeSamples[j] & 0xFFF);
725 new((*digitsTRG)[digitsTRG->GetEntriesFast()]) AliEMCALRawDigit(id, timeSamples, nSamples);
731 //____________________________________________________________________________
732 void AliEMCALDigitizer::DigitalFastOR( Double_t time, Double_t dE, Int_t timeSamples[], Int_t nSamples )
736 const Int_t reso = 11; // 11-bit resolution ADC
737 const Double_t vFSR = 1; // Full scale input voltage range
738 const Double_t Ne = 125; // signal of the APD per MeV of energy deposit in a tower: 125 photo-e-/MeV @ M=30
739 const Double_t vA = .136e-6; // CSP output range: 0.136uV/e-
740 const Double_t rise = 40e-9; // rise time (10-90%) of the FastOR signal before shaping
742 const Double_t kTimeBinWidth = 25E-9; // sampling frequency (40MHz)
744 Double_t vV = 1000. * dE * Ne * vA; // GeV 2 MeV
746 TF1 signalF("signal", AnalogFastORFunction, 0, nSamples * kTimeBinWidth, 3);
747 signalF.SetParameter( 0, vV );
748 signalF.SetParameter( 1, time ); // FIXME: when does the signal arrive? Might account for cable lengths
749 signalF.SetParameter( 2, rise );
751 for (Int_t iTime=0; iTime<nSamples; iTime++)
753 // FIXME: add noise (probably not simply Gaussian) according to DA measurements
754 // probably plan an access to OCDB
756 timeSamples[iTime] = int((TMath::Power(2, reso) / vFSR) * signalF.Eval(iTime * kTimeBinWidth) + 0.5);
761 //____________________________________________________________________________
762 //Float_t AliEMCALDigitizer::FrontEdgeTime(TClonesArray * ticks)
764 // // Returns the shortest time among all time ticks
766 // ticks->Sort() ; //Sort in accordance with times of ticks
768 // AliEMCALTick * ctick = (AliEMCALTick *) it.Next() ;
769 // Float_t time = ctick->CrossingTime(fTimeThreshold) ;
771 // AliEMCALTick * t ;
772 // while((t=(AliEMCALTick*) it.Next())){
773 // if(t->GetTime() < time) //This tick starts before crossing
778 // time = ctick->CrossingTime(fTimeThreshold) ;
784 //____________________________________________________________________________
785 Bool_t AliEMCALDigitizer::Init()
787 // Makes all memory allocations
789 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
791 if ( emcalLoader == 0 ) {
792 Fatal("Init", "Could not obtain the AliEMCALLoader");
797 fLastEvent = fFirstEvent ;
800 fInput = fManager->GetNinputs() ;
804 fInputFileNames = new TString[fInput] ;
805 fEventNames = new TString[fInput] ;
806 fInputFileNames[0] = GetTitle() ;
807 fEventNames[0] = fEventFolderName.Data() ;
809 for (index = 1 ; index < fInput ; index++) {
810 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
811 TString tempo = fManager->GetInputFolderName(index) ;
812 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fManager
815 //to prevent cleaning of this object while GetEvent is called
816 emcalLoader->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
818 //Calibration instance
819 fCalibData = emcalLoader->CalibData();
823 //____________________________________________________________________________
824 void AliEMCALDigitizer::InitParameters()
826 // Parameter initialization for digitizer
828 // Get the parameters from the OCDB via the loader
829 AliRunLoader *rl = AliRunLoader::Instance();
830 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
831 AliEMCALSimParam * simParam = 0x0;
832 if(emcalLoader) simParam = emcalLoader->SimulationParameters();
835 simParam = AliEMCALSimParam::GetInstance();
836 AliWarning("Simulation Parameters not available in OCDB?");
839 fMeanPhotonElectron = simParam->GetMeanPhotonElectron();//4400; // electrons per GeV
840 fPinNoise = simParam->GetPinNoise();//0.012; // pin noise in GeV from analysis test beam data
841 if (fPinNoise < 0.0001 )
842 Warning("InitParameters", "No noise added\n") ;
843 fDigitThreshold = simParam->GetDigitThreshold(); //fPinNoise * 3; // 3 * sigma
844 fTimeResolution = simParam->GetTimeResolution(); //0.6e-9 ; // 600 psc
846 // These defaults are normally not used.
847 // Values are read from calibration database instead
848 fADCchannelEC = 0.0153; // Update 24 Apr 2007: 250./16/1024 - width of one ADC channel in GeV
849 fADCpedestalEC = 0.0 ; // GeV
851 fNADCEC = simParam->GetNADCEC();//(Int_t) TMath::Power(2,16) ; // number of channels in Tower ADC - 65536
853 AliDebug(2,Form("Mean Photon Electron %d, Noise %f, E Threshold %f,Time Resolution %g,NADCEC %d",
854 fMeanPhotonElectron,fPinNoise,fDigitThreshold,fTimeResolution,fNADCEC));
856 // Not used anymore, remove?
857 // fTimeSignalLength = 1.0e-9 ;
858 // fTimeThreshold = 0.001*10000000 ; // Means 1 MeV in terms of SDigits amplitude ??
862 //__________________________________________________________________
863 void AliEMCALDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
865 // Allows to produce digits by superimposing background and signal event.
866 // It is assumed, that headers file with SIGNAL events is opened in
868 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
869 // Thus we avoid writing (changing) huge and expensive
870 // backgound files: all output will be writen into SIGNAL, i.e.
871 // opened in constructor file.
873 // One can open as many files to mix with as one needs.
874 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
877 if( strcmp(GetName(), "") == 0 )
881 Error("MixWith", "Cannot use this method under AliRunDigitizer") ;
884 // looking for file which contains AliRun
885 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
886 Error("MixWith", "File %s does not exist!", alirunFileName.Data()) ;
889 // looking for the file which contains SDigits
890 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
891 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
892 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
893 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
894 if ( (gSystem->AccessPathName(fileName)) ) {
895 Error("MixWith", "The file %s does not exist!", fileName.Data()) ;
898 // need to increase the arrays
899 // MvL: This code only works when fInput == 1, I think.
900 TString tempo = fInputFileNames[fInput-1] ;
901 delete [] fInputFileNames ;
902 fInputFileNames = new TString[fInput+1] ;
903 fInputFileNames[fInput-1] = tempo ;
905 tempo = fEventNames[fInput-1] ;
906 delete [] fEventNames ;
907 fEventNames = new TString[fInput+1] ;
908 fEventNames[fInput-1] = tempo ;
910 fInputFileNames[fInput] = alirunFileName ;
911 fEventNames[fInput] = eventFolderName ;
915 //__________________________________________________________________
916 void AliEMCALDigitizer::Print1(Option_t * option)
917 { // 19-nov-04 - just for convinience
922 //__________________________________________________________________
923 void AliEMCALDigitizer::Print(Option_t*)const
925 // Print Digitizer's parameters
926 printf("Print: \n------------------- %s -------------", GetName() ) ;
927 if( strcmp(fEventFolderName.Data(), "") != 0 ){
928 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
932 nStreams = GetNInputStreams() ;
940 for (index = 0 ; index < nStreams ; index++) {
941 TString tempo(fEventNames[index]) ;
943 if ((rl = AliRunLoader::GetRunLoader(tempo)) == 0)
944 rl = AliRunLoader::Open(fInputFileNames[index], tempo) ;
945 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
946 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
947 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
948 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
949 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
952 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
954 printf("\nWriting digits to %s", emcalLoader->GetDigitsFileName().Data()) ;
956 printf("\nWith following parameters:\n") ;
958 printf(" Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise) ;
959 printf(" Threshold in Tower (fDigitThreshold) = %d\n", fDigitThreshold) ;
960 printf("---------------------------------------------------\n") ;
963 printf("Print: AliEMCALDigitizer not initialized") ;
966 //__________________________________________________________________
967 void AliEMCALDigitizer::PrintDigits(Option_t * option)
969 //utility method for printing digit information
971 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
972 TClonesArray * digits = emcalLoader->Digits() ;
973 TClonesArray * sdigits = emcalLoader->SDigits() ;
975 printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ;
976 printf("\n event %d", emcalLoader->GetRunLoader()->GetEventNumber());
978 if(strstr(option,"all")){
980 AliEMCALDigit * digit;
981 printf("\nEMC digits (with primaries):\n") ;
982 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
984 for (index = 0 ; index < digits->GetEntries() ; index++) {
985 digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
986 printf("\n%6d %8d %6.5e %4d %2d : ",
987 digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
989 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
990 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
997 //__________________________________________________________________
998 Float_t AliEMCALDigitizer::TimeOfNoise(void)
1000 // Calculates the time signal generated by noise
1001 //PH Info("TimeOfNoise", "Change me") ;
1002 return gRandom->Rndm() * 1.28E-5;
1005 //__________________________________________________________________
1006 void AliEMCALDigitizer::Unload()
1008 // Unloads the SDigits and Digits
1012 for(i = 1 ; i < fInput ; i++){
1013 TString tempo(fEventNames[i]) ;
1015 if ((rl = AliRunLoader::GetRunLoader(tempo)))
1016 rl->GetDetectorLoader("EMCAL")->UnloadSDigits() ;
1018 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1019 emcalLoader->UnloadDigits() ;
1022 //_________________________________________________________________________________________
1023 void AliEMCALDigitizer::WriteDigits()
1026 // Makes TreeD in the output file.
1027 // Check if branch already exists:
1028 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
1029 // already existing branches.
1030 // else creates branch with Digits, named "EMCAL", title "...",
1031 // and branch "AliEMCALDigitizer", with the same title to keep all the parameters
1032 // and names of files, from which digits are made.
1034 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1036 const TClonesArray * digits = emcalLoader->Digits() ;
1037 TTree * treeD = emcalLoader->TreeD();
1039 emcalLoader->MakeDigitsContainer();
1040 treeD = emcalLoader->TreeD();
1043 // -- create Digits branch
1044 Int_t bufferSize = 32000 ;
1045 TBranch * digitsBranch = 0;
1046 if ((digitsBranch = treeD->GetBranch("EMCAL"))) {
1047 digitsBranch->SetAddress(&digits);
1048 AliWarning("Digits Branch already exists. Not all digits will be visible");
1051 treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
1052 //digitsBranch->SetTitle(fEventFolderName);
1055 emcalLoader->WriteDigits("OVERWRITE");
1056 emcalLoader->WriteDigitizer("OVERWRITE");
1062 //__________________________________________________________________
1063 void AliEMCALDigitizer::WriteDigits(TClonesArray* digits, const char* branchName)
1066 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1068 TTree* treeD = emcalLoader->TreeD();
1071 emcalLoader->MakeDigitsContainer();
1072 treeD = emcalLoader->TreeD();
1075 // -- create Digits branch
1076 Int_t bufferSize = 32000;
1078 if (TBranch* triggerBranch = treeD->GetBranch(branchName))
1080 triggerBranch->SetAddress(&digits);
1084 treeD->Branch(branchName,"TClonesArray",&digits,bufferSize);
1090 //__________________________________________________________________
1091 void AliEMCALDigitizer::Browse(TBrowser* b)