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?
146 // fTimeThreshold(0), //Not used, remove?
147 // fTimeSignalLength(0), //Not used, remove?
151 fEventFolderName(""),
158 fManager = 0 ; // We work in the standalone mode
161 //____________________________________________________________________________
162 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName)
163 : AliDigitizer("EMCAL"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
164 fDefaultInit(kFALSE),
171 fMeanPhotonElectron(0),
172 // fPedestal(0),//Not used, remove?
173 // fSlope(0), //Not used, remove?
177 // fTimeThreshold(0), //Not used, remove?
178 // fTimeSignalLength(0), //Not used, remove?
182 fEventFolderName(eventFolderName),
190 fManager = 0 ; // We work in the standalone mode
193 //____________________________________________________________________________
194 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d)
195 : AliDigitizer(d.GetName(),d.GetTitle()),
196 fDefaultInit(d.fDefaultInit),
197 fDigitsInRun(d.fDigitsInRun),
200 fInputFileNames(d.fInputFileNames),
201 fEventNames(d.fEventNames),
202 fDigitThreshold(d.fDigitThreshold),
203 fMeanPhotonElectron(d.fMeanPhotonElectron),
204 // fPedestal(d.fPedestal), //Not used, remove?
205 // fSlope(d.fSlope), //Not used, remove?
206 fPinNoise(d.fPinNoise),
207 fTimeDelay(d.fTimeDelay),
208 fTimeResolution(d.fTimeResolution),
209 // fTimeThreshold(d.fTimeThreshold), //Not used, remove?
210 // fTimeSignalLength(d.fTimeSignalLength), //Not used, remove?
211 fADCchannelEC(d.fADCchannelEC),
212 fADCpedestalEC(d.fADCpedestalEC),
214 fEventFolderName(d.fEventFolderName),
215 fFirstEvent(d.fFirstEvent),
216 fLastEvent(d.fLastEvent),
217 fCalibData(d.fCalibData)
222 //____________________________________________________________________________
223 AliEMCALDigitizer::AliEMCALDigitizer(AliRunDigitizer * rd)
224 : AliDigitizer(rd,"EMCAL"+AliConfig::Instance()->GetDigitizerTaskName()),
225 fDefaultInit(kFALSE),
232 fMeanPhotonElectron(0),
233 // fPedestal(0), //Not used, remove?
234 // fSlope(0.), //Not used, remove?
238 // fTimeThreshold(0), //Not used, remove?
239 // fTimeSignalLength(0), //Not used, remove?
248 // ctor Init() is called by RunDigitizer
250 fEventFolderName = fManager->GetInputFolderName(0) ;
251 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
255 //____________________________________________________________________________
256 AliEMCALDigitizer::~AliEMCALDigitizer()
259 if (AliRunLoader::Instance()) {
260 AliLoader *emcalLoader=0;
261 if ((emcalLoader = AliRunLoader::Instance()->GetDetectorLoader("EMCAL")))
262 emcalLoader->CleanDigitizer();
265 AliDebug(1," no runloader present");
266 delete [] fInputFileNames ;
267 delete [] fEventNames ;
271 //____________________________________________________________________________
272 void AliEMCALDigitizer::Digitize(Int_t event)
275 // Makes the digitization of the collected summable digits
276 // for this it first creates the array of all EMCAL modules
277 // filled with noise and after that adds contributions from
278 // SDigits. This design helps to avoid scanning over the
279 // list of digits to add contribution of any new SDigit.
282 // Note that SDigit energy info is stored as an amplitude, so we
283 // must call the Calibrate() method of the SDigitizer to convert it
284 // back to an energy in GeV before adding it to the Digit
286 static int nEMC=0; //max number of digits possible
288 AliRunLoader *rl = AliRunLoader::Instance();
289 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
290 Int_t readEvent = event ;
291 // fManager is data member from AliDigitizer
293 readEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
294 AliDebug(1,Form("Adding event %d from input stream 0 %s %s",
295 readEvent, GetTitle(), fEventFolderName.Data())) ;
296 rl->GetEvent(readEvent);
298 TClonesArray * digits = emcalLoader->Digits() ;
299 digits->Delete() ; //correct way to clear array when memory is
300 //allocated by objects stored in array
303 AliEMCALGeometry *geom = 0;
304 if (rl->GetAliRun()) {
305 AliEMCAL * emcal = (AliEMCAL*)rl->GetAliRun()->GetDetector("EMCAL");
306 geom = emcal->GetGeometry();
309 AliFatal("Could not get AliRun from runLoader");
311 nEMC = geom->GetNCells();
312 AliDebug(1,Form("nEMC %i (number cells in EMCAL) | %s \n", nEMC, geom->GetName()));
316 digits->Expand(nEMC) ;
318 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
319 if ( !emcalLoader->SDigitizer() )
320 emcalLoader->LoadSDigitizer();
321 AliEMCALSDigitizer * sDigitizer = dynamic_cast<AliEMCALSDigitizer*>(emcalLoader->SDigitizer());
324 Fatal("Digitize", "SDigitizer with name %s %s not found", fEventFolderName.Data(), GetTitle() ) ;
326 //take all the inputs to add together and load the SDigits
327 TObjArray * sdigArray = new TObjArray(fInput) ;
328 sdigArray->AddAt(emcalLoader->SDigits(), 0) ;
331 for(i = 1 ; i < fInput ; i++){
332 TString tempo(fEventNames[i]) ;
335 AliRunLoader * rl2 = AliRunLoader::GetRunLoader(tempo) ;
338 rl2 = AliRunLoader::Open(fInputFileNames[i], tempo) ;
341 readEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
342 Info("Digitize", "Adding event %d from input stream %d %s %s", readEvent, i, fInputFileNames[i].Data(), tempo.Data()) ;
344 rl2->GetEvent(readEvent);
345 AliEMCALLoader *emcalLoader2 = dynamic_cast<AliEMCALLoader*>(rl2->GetDetectorLoader("EMCAL"));
346 sdigArray->AddAt(emcalLoader2->SDigits(), i) ;
349 //Find the first tower with signal
350 Int_t nextSig = nEMC + 1 ;
351 TClonesArray * sdigits ;
352 for(i = 0 ; i < fInput ; i++){
353 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
354 if ( !sdigits->GetEntriesFast() )
356 Int_t curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(0))->GetId() ;
357 if(curNext < nextSig)
359 AliDebug(1,Form("input %i : #sdigits %i \n",
360 i, sdigits->GetEntriesFast()));
362 AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
364 TArrayI index(fInput) ;
365 index.Reset() ; //Set all indexes to zero
367 AliEMCALDigit * digit ;
368 AliEMCALDigit * curSDigit ;
370 // TClonesArray * ticks = new TClonesArray("AliEMCALTick",1000) ;
372 //Put Noise contribution
373 for(absID = 0; absID < nEMC; absID++){ // Nov 30, 2006 by PAI; was from 1 to nEMC
375 // amplitude set to zero, noise will be added later
376 new((*digits)[absID]) AliEMCALDigit( -1, -1, absID, 0., TimeOfNoise(),kFALSE); // absID-1->absID
377 //look if we have to add signal?
378 digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID)); // absID-1->absID
381 //Add SDigits from all inputs
383 //Int_t contrib = 0 ;
385 //Follow PHOS and comment out this timing model til a better one
386 //can be developed - JLK 28-Apr-2008
388 //Float_t a = digit->GetAmplitude() ;
389 //Float_t b = TMath::Abs( a /fTimeSignalLength) ;
390 //Mark the beginning of the signal
391 //new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime(),0, b);
392 //Mark the end of the signal
393 //new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime()+fTimeSignalLength, -a, -b);
395 // Calculate time as time of the largest digit
396 Float_t time = digit->GetTime() ;
397 Float_t aTime= digit->GetAmplitude() ;
400 for(i = 0; i< fInput ; i++){ //loop over (possible) merge sources
401 if(dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
402 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
405 //May be several digits will contribute from the same input
406 while(curSDigit && (curSDigit->GetId() == absID)){
407 //Shift primary to separate primaries belonging different inputs
408 Int_t primaryoffset ;
410 primaryoffset = fManager->GetMask(i) ;
413 curSDigit->ShiftPrimary(primaryoffset) ;
415 //Remove old timing model - JLK 28-April-2008
416 //a = curSDigit->GetAmplitude() ;
417 //b = a /fTimeSignalLength ;
418 //new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime(),0, b);
419 //new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
420 if(curSDigit->GetAmplitude()>aTime) {
421 aTime = curSDigit->GetAmplitude();
422 time = curSDigit->GetTime();
425 *digit = *digit + *curSDigit ; //adds amplitudes of each digit
428 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
429 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
434 //Here we convert the summed amplitude to an energy in GeV
435 energy = sDigitizer->Calibrate(digit->GetAmplitude()) ; // GeV
436 // add fluctuations for photo-electron creation
437 energy *= static_cast<Float_t>(gRandom->Poisson(fMeanPhotonElectron)) / static_cast<Float_t>(fMeanPhotonElectron) ;
439 //calculate and set time
440 //New timing model needed - JLK 28-April-2008
441 //Float_t time = FrontEdgeTime(ticks) ;
442 digit->SetTime(time) ;
444 //Find next signal module
446 for(i = 0 ; i < fInput ; i++){
447 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
448 Int_t curNext = nextSig ;
449 if(sdigits->GetEntriesFast() > index[i] ){
450 curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))->GetId() ;
452 if(curNext < nextSig) nextSig = curNext ;
456 energy += TMath::Abs(gRandom->Gaus(0., fPinNoise)) ;
458 //Now digitize the energy according to the sDigitizer method,
459 //which merely converts the energy to an integer. Later we will
460 //check that the stored value matches our allowed dynamic ranges
461 digit->SetAmplitude(sDigitizer->Digitize(energy)) ;
462 AliDebug(10,Form(" absID %5i energy %f nextSig %5i\n",
463 absID, energy, nextSig));
464 } // for(absID = 0; absID < nEMC; absID++)
469 //JLK is it better to call Clear() here?
470 delete sdigArray ; //We should not delete its contents
472 //remove digits below thresholds
473 // until 10-02-2010 remove digits with energy smaller than fDigitThreshold 3*fPinNoise
474 // now, remove digits with Digitized ADC smaller than fDigitThreshold = 3
476 for(i = 0 ; i < nEMC ; i++){
477 digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
478 //First get the energy in GeV units.
479 energy = sDigitizer->Calibrate(digit->GetAmplitude()) ;
480 //Then digitize using the calibration constants of the ocdb
481 Float_t ampADC = DigitizeEnergy(energy, digit->GetId()) ;
482 //if(ampADC>2)printf("Digit energy %f, amp %d, amp cal %d, threshold %d\n",energy,digit->GetAmplitude(),ampADC,fDigitThreshold);
483 if(ampADC < fDigitThreshold)
484 digits->RemoveAt(i) ;
486 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
491 Int_t ndigits = digits->GetEntriesFast() ;
494 //After we have done the summing and digitizing to create the
495 //digits, now we want to calibrate the resulting amplitude to match
496 //the dynamic range of our real data.
497 for (i = 0 ; i < ndigits ; i++) {
498 digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ;
499 digit->SetIndexInList(i) ;
500 energy = sDigitizer->Calibrate(digit->GetAmplitude()) ;
501 digit->SetAmplitude(DigitizeEnergy(energy, digit->GetId()) ) ;
503 digit->SetTime(digit->GetTime()+fTimeDelay) ;
504 // printf("digit amplitude set at end: i %d, amp %f\n",i,digit->GetAmplitude());
509 // //_____________________________________________________________________
510 Float_t AliEMCALDigitizer::DigitizeEnergy(Float_t energy, Int_t AbsId)
513 // Returns digitized value of the energy in a cell absId
514 // using the calibration constants stored in the OCDB
515 // or default values if no CalibData object is found.
516 // This effectively converts everything to match the dynamic range
517 // of the real data we will collect
520 const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
523 AliFatal("Did not get geometry from EMCALLoader");
531 Float_t channel = -999;
533 Bool_t bCell = geom->GetCellIndex(AbsId, iSupMod, nModule, nIphi, nIeta) ;
535 Error("DigitizeEnergy","Wrong cell id number : AbsId %i ", AbsId) ;
536 geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
539 fADCpedestalEC = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
540 fADCchannelEC = fCalibData->GetADCchannel(iSupMod,ieta,iphi);
543 //channel = static_cast<Int_t> (TMath::Floor( (energy + fADCpedestalEC)/fADCchannelEC )) ;
544 channel = (energy + fADCpedestalEC)/fADCchannelEC ;
546 if(channel > fNADCEC )
552 //____________________________________________________________________________
553 void AliEMCALDigitizer::Exec(Option_t *option)
555 // Steering method to process digitization for events
556 // in the range from fFirstEvent to fLastEvent.
557 // This range is optionally set by SetEventRange().
558 // if fLastEvent=-1, then process events until the end.
559 // by default fLastEvent = fFirstEvent (process only one event)
561 if (!fInit) { // to prevent overwrite existing file
562 Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
566 if (strstr(option,"print")) {
572 if(strstr(option,"tim"))
573 gBenchmark->Start("EMCALDigitizer");
575 AliRunLoader *rl = AliRunLoader::Instance();
576 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
578 // Post Digitizer to the white board
579 emcalLoader->PostDigitizer(this) ;
581 if (fLastEvent == -1) {
582 fLastEvent = rl->GetNumberOfEvents() - 1 ;
585 fLastEvent = fFirstEvent ; // what is this ??
587 Int_t nEvents = fLastEvent - fFirstEvent + 1;
590 TClonesArray* digitsTRG = new TClonesArray("AliEMCALRawDigit", 32 * 96);
591 TClonesArray* digitsTMP = new TClonesArray("AliEMCALDigit", 32 * 96);
592 rl->LoadSDigits("EMCAL");
593 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
595 rl->GetEvent(ievent);
597 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
602 //-------------------------------------
603 Digits2FastOR(digitsTMP, digitsTRG);
605 WriteDigits(digitsTRG);
607 (emcalLoader->TreeD())->Fill();
609 emcalLoader->WriteDigits( "OVERWRITE");
610 emcalLoader->WriteDigitizer("OVERWRITE");
614 digitsTRG->Clear("C");
615 digitsTMP->Clear("C");
616 //-------------------------------------
618 if(strstr(option,"deb"))
620 if(strstr(option,"table")) gObjectTable->Print();
622 //increment the total number of Digits per run
623 fDigitsInRun += emcalLoader->Digits()->GetEntriesFast() ;
626 emcalLoader->CleanDigitizer() ;
628 if(strstr(option,"tim")){
629 gBenchmark->Stop("EMCALDigitizer");
630 AliInfo(Form("Exec: took %f seconds for Digitizing %f seconds per event",
631 gBenchmark->GetCpuTime("EMCALDigitizer"), gBenchmark->GetCpuTime("EMCALDigitizer")/nEvents )) ;
635 //____________________________________________________________________________
636 void AliEMCALDigitizer::Digits2FastOR(TClonesArray* digitsTMP, TClonesArray* digitsTRG)
638 // FEE digits afterburner to produce TRG digits
639 // we are only interested in the FEE digit deposited energy
640 // to be converted later into a voltage value
642 // push the FEE digit to its associated FastOR (numbered from 0:95)
643 // TRU is in charge of summing module digits
645 AliRunLoader *runLoader = AliRunLoader::Instance();
647 AliRun* run = runLoader->GetAliRun();
649 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
651 AliEMCALGeometry* geom = dynamic_cast<AliEMCAL*>(run->GetDetector("EMCAL"))->GetGeometry();
653 // build FOR from simulated digits
654 // and xfer to the corresponding TRU input (mapping)
656 TClonesArray* digits = emcalLoader->Digits();
658 TIter NextDigit(digits);
659 while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
661 Int_t id = digit->GetId();
663 Int_t iSupMod, nModule, nIphi, nIeta, iphi, ieta, iphim, ietam;
665 geom->GetCellIndex( id, iSupMod, nModule, nIphi, nIeta );
666 geom->GetModulePhiEtaIndexInSModule( iSupMod, nModule, iphim, ietam );
667 geom->GetCellPhiEtaIndexInSModule( iSupMod, nModule, nIphi, nIeta, iphi, ieta);
669 // identify to which TRU this FEE digit belong
670 //Int_t itru = (iSupMod < 11) ? iphim / 4 + 3 * iSupMod : 31;
673 itru = (iSupMod % 2) ? (2 - int(iphim / 4)) + 3 * iSupMod : iphim / 4 + 3 * iSupMod;
679 // FIXME: bad numbering solution to deal w/ the last 2 SM which have only 1 TRU each
680 // using the AliEMCALGeometry official numbering
681 // only 1 TRU/SM in SM 10 & SM 11
684 if ((itru == 31 && iphim < 2) || (itru == 30 && iphim > 5)) continue;
686 // to be compliant with %4 per TRU
687 if (itru == 31) iphim -= 2;
690 Bool_t isOK = geom->GetAbsFastORIndexFromPositionInTRU(itru, ietam, iphim % 4, trgid);
692 AliDebug(2,Form("trigger digit id: %d itru: %d isOK: %d\n",trgid,itru,isOK));
696 AliEMCALDigit* d = static_cast<AliEMCALDigit*>(digitsTMP->At(trgid));
700 new((*digitsTMP)[trgid]) AliEMCALDigit(*digit);
701 d = (AliEMCALDigit*)digitsTMP->At(trgid);
712 Int_t *timeSamples = new Int_t[nSamples];
714 NextDigit = TIter(digitsTMP);
715 while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
719 Int_t id = digit->GetId();
720 Float_t time = digit->GetTime();
722 Double_t depositedEnergy = 0.;
723 for (Int_t j = 1; j <= digit->GetNprimary(); j++) depositedEnergy += digit->GetDEPrimary(j);
725 // FIXME: Check digit time!
728 DigitalFastOR(time, depositedEnergy, timeSamples, nSamples);
730 for (Int_t j=0;j<nSamples;j++)
732 timeSamples[j] = ((j << 12) & 0xFF000) | (timeSamples[j] & 0xFFF);
735 new((*digitsTRG)[digitsTRG->GetEntriesFast()]) AliEMCALRawDigit(id, timeSamples, nSamples);
740 delete [] timeSamples;
743 //____________________________________________________________________________
744 void AliEMCALDigitizer::DigitalFastOR( Double_t time, Double_t dE, Int_t timeSamples[], Int_t nSamples )
748 const Int_t reso = 11; // 11-bit resolution ADC
749 const Double_t vFSR = 1; // Full scale input voltage range
750 const Double_t Ne = 125; // signal of the APD per MeV of energy deposit in a tower: 125 photo-e-/MeV @ M=30
751 const Double_t vA = .136e-6; // CSP output range: 0.136uV/e-
752 const Double_t rise = 40e-9; // rise time (10-90%) of the FastOR signal before shaping
754 const Double_t kTimeBinWidth = 25E-9; // sampling frequency (40MHz)
756 Double_t vV = 1000. * dE * Ne * vA; // GeV 2 MeV
758 TF1 signalF("signal", AnalogFastORFunction, 0, nSamples * kTimeBinWidth, 3);
759 signalF.SetParameter( 0, vV );
760 signalF.SetParameter( 1, time ); // FIXME: when does the signal arrive? Might account for cable lengths
761 signalF.SetParameter( 2, rise );
763 for (Int_t iTime=0; iTime<nSamples; iTime++)
765 // FIXME: add noise (probably not simply Gaussian) according to DA measurements
766 // probably plan an access to OCDB
768 timeSamples[iTime] = int((TMath::Power(2, reso) / vFSR) * signalF.Eval(iTime * kTimeBinWidth) + 0.5);
773 //____________________________________________________________________________
774 //Float_t AliEMCALDigitizer::FrontEdgeTime(TClonesArray * ticks)
776 // // Returns the shortest time among all time ticks
778 // ticks->Sort() ; //Sort in accordance with times of ticks
780 // AliEMCALTick * ctick = (AliEMCALTick *) it.Next() ;
781 // Float_t time = ctick->CrossingTime(fTimeThreshold) ;
783 // AliEMCALTick * t ;
784 // while((t=(AliEMCALTick*) it.Next())){
785 // if(t->GetTime() < time) //This tick starts before crossing
790 // time = ctick->CrossingTime(fTimeThreshold) ;
796 //____________________________________________________________________________
797 Bool_t AliEMCALDigitizer::Init()
799 // Makes all memory allocations
801 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
803 if ( emcalLoader == 0 ) {
804 Fatal("Init", "Could not obtain the AliEMCALLoader");
809 fLastEvent = fFirstEvent ;
812 fInput = fManager->GetNinputs() ;
816 fInputFileNames = new TString[fInput] ;
817 fEventNames = new TString[fInput] ;
818 fInputFileNames[0] = GetTitle() ;
819 fEventNames[0] = fEventFolderName.Data() ;
821 for (index = 1 ; index < fInput ; index++) {
822 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
823 TString tempo = fManager->GetInputFolderName(index) ;
824 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fManager
827 //to prevent cleaning of this object while GetEvent is called
828 emcalLoader->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
830 //Calibration instance
831 fCalibData = emcalLoader->CalibData();
835 //____________________________________________________________________________
836 void AliEMCALDigitizer::InitParameters()
838 // Parameter initialization for digitizer
840 // Get the parameters from the OCDB via the loader
841 AliRunLoader *rl = AliRunLoader::Instance();
842 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
843 AliEMCALSimParam * simParam = 0x0;
844 if(emcalLoader) simParam = emcalLoader->SimulationParameters();
847 simParam = AliEMCALSimParam::GetInstance();
848 AliWarning("Simulation Parameters not available in OCDB?");
851 fMeanPhotonElectron = simParam->GetMeanPhotonElectron();//4400; // electrons per GeV
852 fPinNoise = simParam->GetPinNoise();//0.012; // pin noise in GeV from analysis test beam data
853 if (fPinNoise < 0.0001 )
854 Warning("InitParameters", "No noise added\n") ;
855 fDigitThreshold = simParam->GetDigitThreshold(); //fPinNoise * 3; // 3 * sigma
856 fTimeResolution = simParam->GetTimeResolution(); //0.6e-9 ; // 600 pc
857 fTimeDelay = simParam->GetTimeDelay(); //600e-9 ; // 600 ns
859 // These defaults are normally not used.
860 // Values are read from calibration database instead
861 fADCchannelEC = 0.0153; // Update 24 Apr 2007: 250./16/1024 - width of one ADC channel in GeV
862 fADCpedestalEC = 0.0 ; // GeV
864 fNADCEC = simParam->GetNADCEC();//(Int_t) TMath::Power(2,16) ; // number of channels in Tower ADC - 65536
866 AliDebug(2,Form("Mean Photon Electron %d, Noise %f, Digit Threshold %d,Time Resolution %g,NADCEC %d",
867 fMeanPhotonElectron,fPinNoise,fDigitThreshold,fTimeResolution,fNADCEC));
869 // Not used anymore, remove?
870 // fTimeSignalLength = 1.0e-9 ;
871 // fTimeThreshold = 0.001*10000000 ; // Means 1 MeV in terms of SDigits amplitude ??
875 //__________________________________________________________________
876 void AliEMCALDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
878 // Allows to produce digits by superimposing background and signal event.
879 // It is assumed, that headers file with SIGNAL events is opened in
881 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
882 // Thus we avoid writing (changing) huge and expensive
883 // backgound files: all output will be writen into SIGNAL, i.e.
884 // opened in constructor file.
886 // One can open as many files to mix with as one needs.
887 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
890 if( strcmp(GetName(), "") == 0 )
894 Error("MixWith", "Cannot use this method under AliRunDigitizer") ;
897 // looking for file which contains AliRun
898 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
899 Error("MixWith", "File %s does not exist!", alirunFileName.Data()) ;
902 // looking for the file which contains SDigits
903 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
904 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
905 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
906 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
907 if ( (gSystem->AccessPathName(fileName)) ) {
908 Error("MixWith", "The file %s does not exist!", fileName.Data()) ;
911 // need to increase the arrays
912 // MvL: This code only works when fInput == 1, I think.
913 TString tempo = fInputFileNames[fInput-1] ;
914 delete [] fInputFileNames ;
915 fInputFileNames = new TString[fInput+1] ;
916 fInputFileNames[fInput-1] = tempo ;
918 tempo = fEventNames[fInput-1] ;
919 delete [] fEventNames ;
920 fEventNames = new TString[fInput+1] ;
921 fEventNames[fInput-1] = tempo ;
923 fInputFileNames[fInput] = alirunFileName ;
924 fEventNames[fInput] = eventFolderName ;
928 //__________________________________________________________________
929 void AliEMCALDigitizer::Print1(Option_t * option)
930 { // 19-nov-04 - just for convinience
935 //__________________________________________________________________
936 void AliEMCALDigitizer::Print(Option_t*)const
938 // Print Digitizer's parameters
939 printf("Print: \n------------------- %s -------------", GetName() ) ;
940 if( strcmp(fEventFolderName.Data(), "") != 0 ){
941 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
945 nStreams = GetNInputStreams() ;
952 for (index = 0 ; index < nStreams ; index++) {
953 TString tempo(fEventNames[index]) ;
955 if ((rl = AliRunLoader::GetRunLoader(tempo)) == 0)
956 rl = AliRunLoader::Open(fInputFileNames[index], tempo) ;
957 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
958 TString fileName( emcalLoader->GetSDigitsFileName() ) ;
959 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
960 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
961 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
964 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
966 printf("\nWriting digits to %s", emcalLoader->GetDigitsFileName().Data()) ;
968 printf("\nWith following parameters:\n") ;
970 printf(" Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise) ;
971 printf(" Threshold in Tower (fDigitThreshold) = %d\n", fDigitThreshold) ;
972 printf("---------------------------------------------------\n") ;
975 printf("Print: AliEMCALDigitizer not initialized") ;
978 //__________________________________________________________________
979 void AliEMCALDigitizer::PrintDigits(Option_t * option)
981 //utility method for printing digit information
983 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
984 TClonesArray * digits = emcalLoader->Digits() ;
985 TClonesArray * sdigits = emcalLoader->SDigits() ;
987 printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ;
988 printf("\n event %d", emcalLoader->GetRunLoader()->GetEventNumber());
990 if(strstr(option,"all")){
992 AliEMCALDigit * digit;
993 printf("\nEMC digits (with primaries):\n") ;
994 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
996 for (index = 0 ; index < digits->GetEntries() ; index++) {
997 digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
998 printf("\n%6d %8f %6.5e %4d %2d : ",
999 digit->GetId(), digit->GetAmplitude(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
1001 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
1002 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
1009 //__________________________________________________________________
1010 Float_t AliEMCALDigitizer::TimeOfNoise(void)
1012 // Calculates the time signal generated by noise
1013 //PH Info("TimeOfNoise", "Change me") ;
1014 return gRandom->Rndm() * 1.28E-5;
1017 //__________________________________________________________________
1018 void AliEMCALDigitizer::Unload()
1020 // Unloads the SDigits and Digits
1024 for(i = 1 ; i < fInput ; i++){
1025 TString tempo(fEventNames[i]) ;
1027 if ((rl = AliRunLoader::GetRunLoader(tempo)))
1028 rl->GetDetectorLoader("EMCAL")->UnloadSDigits() ;
1030 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1031 emcalLoader->UnloadDigits() ;
1034 //_________________________________________________________________________________________
1035 void AliEMCALDigitizer::WriteDigits()
1038 // Makes TreeD in the output file.
1039 // Check if branch already exists:
1040 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
1041 // already existing branches.
1042 // else creates branch with Digits, named "EMCAL", title "...",
1043 // and branch "AliEMCALDigitizer", with the same title to keep all the parameters
1044 // and names of files, from which digits are made.
1046 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1048 const TClonesArray * digits = emcalLoader->Digits() ;
1049 TTree * treeD = emcalLoader->TreeD();
1051 emcalLoader->MakeDigitsContainer();
1052 treeD = emcalLoader->TreeD();
1055 // -- create Digits branch
1056 Int_t bufferSize = 32000 ;
1057 TBranch * digitsBranch = 0;
1058 if ((digitsBranch = treeD->GetBranch("EMCAL"))) {
1059 digitsBranch->SetAddress(&digits);
1060 AliWarning("Digits Branch already exists. Not all digits will be visible");
1063 treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
1064 //digitsBranch->SetTitle(fEventFolderName);
1068 emcalLoader->WriteDigits("OVERWRITE");
1069 emcalLoader->WriteDigitizer("OVERWRITE");
1075 //__________________________________________________________________
1076 void AliEMCALDigitizer::WriteDigits(TClonesArray* digits, const char* branchName)
1079 AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1081 TTree* treeD = emcalLoader->TreeD();
1084 emcalLoader->MakeDigitsContainer();
1085 treeD = emcalLoader->TreeD();
1088 // -- create Digits branch
1089 Int_t bufferSize = 32000;
1091 if (TBranch* triggerBranch = treeD->GetBranch(branchName))
1093 triggerBranch->SetAddress(&digits);
1097 treeD->Branch(branchName,"TClonesArray",&digits,bufferSize);
1103 //__________________________________________________________________
1104 void AliEMCALDigitizer::Browse(TBrowser* b)