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 //_________________________________________________________________________
19 // Utility Class for handling Raw data
20 // Does all transitions from Digits to Raw and vice versa,
21 // for simu and reconstruction
23 // Note: the current version is still simplified. Only
24 // one raw signal per digit is generated; either high-gain or low-gain
25 // Need to add concurrent high and low-gain info in the future
26 // No pedestal is added to the raw signal.
27 //*-- Author: Marco van Leeuwen (LBL)
29 #include "AliEMCALRawUtils.h"
38 #include "AliRunLoader.h"
39 class AliCaloAltroMapping;
40 #include "AliAltroBuffer.h"
41 #include "AliRawReader.h"
42 #include "AliCaloRawStreamV3.h"
45 #include "AliEMCALRecParam.h"
46 #include "AliEMCALLoader.h"
47 #include "AliEMCALGeometry.h"
48 class AliEMCALDigitizer;
49 #include "AliEMCALDigit.h"
50 #include "AliEMCALRawDigit.h"
52 #include "AliCaloCalibPedestal.h"
53 #include "AliCaloFastAltroFitv0.h"
54 #include "AliCaloNeuralFit.h"
55 #include "AliCaloBunchInfo.h"
56 #include "AliCaloFitResults.h"
57 #include "AliCaloRawAnalyzerFastFit.h"
58 #include "AliCaloRawAnalyzerNN.h"
59 #include "AliCaloRawAnalyzerLMS.h"
60 #include "AliCaloRawAnalyzerPeakFinder.h"
61 #include "AliCaloRawAnalyzerCrude.h"
63 ClassImp(AliEMCALRawUtils)
65 // Signal shape parameters
66 Int_t AliEMCALRawUtils::fgTimeBins = 256; // number of sampling bins of the raw RO signal (we typically use 15-50; theoretical max is 1k+)
67 Double_t AliEMCALRawUtils::fgTimeBinWidth = 100E-9 ; // each sample is 100 ns
68 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ; // 15 time bins ~ 1.5 musec
70 // some digitization constants
71 Int_t AliEMCALRawUtils::fgThreshold = 1;
72 Int_t AliEMCALRawUtils::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
73 Int_t AliEMCALRawUtils::fgPedestalValue = 0; // pedestal value for digits2raw, default generate ZS data
74 Double_t AliEMCALRawUtils::fgFEENoise = 3.; // 3 ADC channels of noise (sampled)
76 AliEMCALRawUtils::AliEMCALRawUtils(fitAlgorithm fitAlgo)
77 : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
78 fNPedSamples(0), fGeom(0), fOption(""),
79 fRemoveBadChannels(kTRUE),fFittingAlgorithm(0),fUseFALTRO(kFALSE),fRawAnalyzer(0)
82 //These are default parameters.
83 //Can be re-set from without with setter functions
84 //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
85 fHighLowGainFactor = 16. ; // Adjusted for a low gain range of 82 GeV (10 bits)
86 fOrder = 2; // Order of gamma fn
87 fTau = 2.35; // in units of timebin, from CERN 2007 testbeam
88 fNoiseThreshold = 3; // 3 ADC counts is approx. noise level
89 fNPedSamples = 4; // Less than this value => likely pedestal samples
90 fRemoveBadChannels = kFALSE; // Do not remove bad channels before fitting
91 fUseFALTRO = kTRUE; // Get the trigger FALTRO information and pass it to digits.
92 SetFittingAlgorithm(fitAlgo);
94 //Get Mapping RCU files from the AliEMCALRecParam
95 const TObjArray* maps = AliEMCALRecParam::GetMappings();
96 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
98 for(Int_t i = 0; i < 4; i++) {
99 fMapping[i] = (AliAltroMapping*)maps->At(i);
102 //To make sure we match with the geometry in a simulation file,
103 //let's try to get it first. If not, take the default geometry
104 AliRunLoader *rl = AliRunLoader::Instance();
105 if (rl && rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
106 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
108 AliInfo(Form("Using default geometry in raw reco"));
109 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
112 if(!fGeom) AliFatal(Form("Could not get geometry!"));
116 //____________________________________________________________________________
117 AliEMCALRawUtils::AliEMCALRawUtils(AliEMCALGeometry *pGeometry, fitAlgorithm fitAlgo)
118 : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
119 fNPedSamples(0), fGeom(pGeometry), fOption(""),
120 fRemoveBadChannels(kTRUE),fFittingAlgorithm(0),fUseFALTRO(kFALSE),fRawAnalyzer()
123 // Initialize with the given geometry - constructor required by HLT
124 // HLT does not use/support AliRunLoader(s) instances
125 // This is a minimum intervention solution
126 // Comment by MPloskon@lbl.gov
129 //These are default parameters.
130 //Can be re-set from without with setter functions
131 //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
132 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
133 fOrder = 2; // order of gamma fn
134 fTau = 2.35; // in units of timebin, from CERN 2007 testbeam
135 fNoiseThreshold = 3; // 3 ADC counts is approx. noise level
136 fNPedSamples = 4; // Less than this value => likely pedestal samples
137 fRemoveBadChannels = kFALSE; // Do not remove bad channels before fitting
138 fUseFALTRO = kTRUE; // Get the trigger FALTRO information and pass it to digits.
139 SetFittingAlgorithm(fitAlgo);
141 //Get Mapping RCU files from the AliEMCALRecParam
142 const TObjArray* maps = AliEMCALRecParam::GetMappings();
143 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
145 for(Int_t i = 0; i < 4; i++) {
146 fMapping[i] = (AliAltroMapping*)maps->At(i);
149 if(!fGeom) AliFatal(Form("Could not get geometry!"));
153 //____________________________________________________________________________
154 AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
156 fHighLowGainFactor(rawU.fHighLowGainFactor),
159 fNoiseThreshold(rawU.fNoiseThreshold),
160 fNPedSamples(rawU.fNPedSamples),
162 fOption(rawU.fOption),
163 fRemoveBadChannels(rawU.fRemoveBadChannels),
164 fFittingAlgorithm(rawU.fFittingAlgorithm),
165 fUseFALTRO(rawU.fUseFALTRO),
166 fRawAnalyzer(rawU.fRawAnalyzer)
169 fMapping[0] = rawU.fMapping[0];
170 fMapping[1] = rawU.fMapping[1];
171 fMapping[2] = rawU.fMapping[2];
172 fMapping[3] = rawU.fMapping[3];
175 //____________________________________________________________________________
176 AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
178 //assignment operator
181 fHighLowGainFactor = rawU.fHighLowGainFactor;
182 fOrder = rawU.fOrder;
184 fNoiseThreshold = rawU.fNoiseThreshold;
185 fNPedSamples = rawU.fNPedSamples;
187 fOption = rawU.fOption;
188 fRemoveBadChannels = rawU.fRemoveBadChannels;
189 fFittingAlgorithm = rawU.fFittingAlgorithm;
190 fUseFALTRO = rawU.fUseFALTRO;
191 fRawAnalyzer = rawU.fRawAnalyzer;
192 fMapping[0] = rawU.fMapping[0];
193 fMapping[1] = rawU.fMapping[1];
194 fMapping[2] = rawU.fMapping[2];
195 fMapping[3] = rawU.fMapping[3];
202 //____________________________________________________________________________
203 AliEMCALRawUtils::~AliEMCALRawUtils() {
208 //____________________________________________________________________________
209 void AliEMCALRawUtils::Digits2Raw()
211 // convert digits of the current event to raw data
213 AliRunLoader *rl = AliRunLoader::Instance();
214 AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
217 loader->LoadDigits("EMCAL");
219 TClonesArray* digits = loader->Digits() ;
222 Warning("Digits2Raw", "no digits found !");
226 static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
227 AliAltroBuffer* buffers[nDDL];
228 for (Int_t i=0; i < nDDL; i++)
231 TArrayI adcValuesLow(fgTimeBins);
232 TArrayI adcValuesHigh(fgTimeBins);
234 // loop over digits (assume ordered digits)
235 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
236 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
237 if (digit->GetAmp() < fgThreshold)
247 fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
248 fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
250 //Check which is the RCU, 0 or 1, of the cell.
253 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
254 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
257 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
259 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
261 if (nSM%2==1) iRCU = 1 - iRCU; // swap for odd=C side, to allow us to cable both sides the same
264 Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
267 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
269 Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
271 if (buffers[iDDL] == 0) {
272 // open new file and write dummy header
273 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
274 //Select mapping file RCU0A, RCU0C, RCU1A, RCU1C
275 Int_t iRCUside=iRCU+(nSM%2)*2;
276 //iRCU=0 and even (0) SM -> RCU0A.data 0
277 //iRCU=1 and even (0) SM -> RCU1A.data 1
278 //iRCU=0 and odd (1) SM -> RCU0C.data 2
279 //iRCU=1 and odd (1) SM -> RCU1C.data 3
280 //cout<<" nSM "<<nSM<<"; iRCU "<<iRCU<<"; iRCUside "<<iRCUside<<endl;
281 buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCUside]);
282 buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE); //Dummy;
285 // out of time range signal (?)
286 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
287 AliInfo("Signal is out of time range.\n");
288 buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
289 buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() ); // time bin
290 buffers[iDDL]->FillBuffer(3); // bunch length
291 buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM); // trailer
292 // calculate the time response function
294 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh.GetArray(), adcValuesLow.GetArray()) ;
296 buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow.GetArray(), fgThreshold);
298 buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh.GetArray(), fgThreshold);
302 // write headers and close files
303 for (Int_t i=0; i < nDDL; i++) {
306 buffers[i]->WriteDataHeader(kFALSE, kFALSE);
311 loader->UnloadDigits();
314 //____________________________________________________________________________
315 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr, const AliCaloCalibPedestal* pedbadmap, TClonesArray *digitsTRG)
317 // convert raw data of the current event to digits
322 Error("Raw2Digits", "no digits found !");
326 Error("Raw2Digits", "no raw reader found !");
330 AliCaloRawStreamV3 in(reader,"EMCAL",fMapping);
331 // Select EMCAL DDL's;
332 reader->Select("EMCAL",0,43); // 43 = AliEMCALGeoParams::fgkLastAltroDDL
334 // fRawAnalyzer setup
335 fRawAnalyzer->SetNsampleCut(5); // requirement for fits to be done
336 fRawAnalyzer->SetAmpCut(fNoiseThreshold);
337 fRawAnalyzer->SetFitArrayCut(fNoiseThreshold);
338 fRawAnalyzer->SetIsZeroSuppressed(true); // TMP - should use stream->IsZeroSuppressed(), or altro cfg registers later
340 // channel info parameters
342 Int_t caloFlag = 0; // low, high gain, or TRU, or LED ref.
344 // start loop over input stream
345 while (in.NextDDL()) {
347 // if ( in.GetDDLNumber() != 0 && in.GetDDLNumber() != 2 ) continue;
349 while (in.NextChannel()) {
352 Int_t hhwAdd = in.GetHWAddress();
353 UShort_t iiBranch = ( hhwAdd >> 11 ) & 0x1; // 0/1
354 UShort_t iiFEC = ( hhwAdd >> 7 ) & 0xF;
355 UShort_t iiChip = ( hhwAdd >> 4 ) & 0x7;
356 UShort_t iiChannel = hhwAdd & 0xF;
358 if ( !( iiBranch == 0 && iiFEC == 1 && iiChip == 3 && ( iiChannel >= 8 && iiChannel <= 15 ) ) && !( iiBranch == 1 && iiFEC == 0 && in.GetColumn() == 0 ) ) continue;
361 //Check if the signal is high or low gain and then do the fit,
362 //if it is from TRU or LEDMon do not fit
363 caloFlag = in.GetCaloFlag();
364 // if (caloFlag != 0 && caloFlag != 1) continue;
365 if (caloFlag > 2) continue; // Work with ALTRO and FALTRO
367 //Do not fit bad channels of ALTRO
368 if(caloFlag < 2 && fRemoveBadChannels && pedbadmap->IsBadChannel(in.GetModule(),in.GetColumn(),in.GetRow())) {
369 //printf("Tower from SM %d, column %d, row %d is BAD!!! Skip \n", in.GetModule(),in.GetColumn(),in.GetRow());
373 vector<AliCaloBunchInfo> bunchlist;
374 while (in.NextBunch()) {
375 bunchlist.push_back( AliCaloBunchInfo(in.GetStartTimeBin(), in.GetBunchLength(), in.GetSignals() ) );
376 } // loop over bunches
379 if ( caloFlag < 2 ){ // ALTRO
384 if ( fFittingAlgorithm == kFastFit || fFittingAlgorithm == kNeuralNet || fFittingAlgorithm == kLMS || fFittingAlgorithm == kPeakFinder || fFittingAlgorithm == kCrude) {
385 // all functionality to determine amp and time etc is encapsulated inside the Evaluate call for these methods
386 AliCaloFitResults fitResults = fRawAnalyzer->Evaluate( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2());
388 amp = fitResults.GetAmp();
389 time = fitResults.GetTof();
391 else { // for the other methods we for now use the functionality of
392 // AliCaloRawAnalyzer as well, to select samples and prepare for fits,
393 // if it looks like there is something to fit
396 Float_t ampEstimate = 0;
398 short timeEstimate = 0;
399 Float_t pedEstimate = 0;
402 Int_t bunchIndex = 0;
404 // The PreFitEvaluateSamples + later call to FitRaw will hopefully
405 // be replaced by a single Evaluate call or so soon, like for the other
406 // methods, but this should be good enough for evaluation of
407 // the methods for now (Jan. 2010)
409 int nsamples = fRawAnalyzer->PreFitEvaluateSamples( bunchlist, in.GetAltroCFG1(), in.GetAltroCFG2(), bunchIndex, ampEstimate, maxADC, timeEstimate, pedEstimate, first, last);
411 if (ampEstimate > fNoiseThreshold) { // something worth looking at
413 time = timeEstimate; // maxrev in AliCaloRawAnalyzer speak; comes with an offset w.r.t. real timebin
414 Int_t timebinOffset = bunchlist.at(bunchIndex).GetStartBin() - (bunchlist.at(bunchIndex).GetLength()-1);
418 if ( nsamples > 1 ) { // possibly something to fit
419 FitRaw(first, last, amp, time);
420 time += timebinOffset;
423 if ( amp>0 && time>0 ) { // brief sanity check of fit results
425 // check fit results: should be consistent with initial estimates
426 // more magic numbers, but very loose cuts, for now..
427 // We have checked that amp and ampEstimate values are positive so division for assymmetry
428 // calculation should be OK/safe
429 Float_t ampAsymm = (amp - ampEstimate)/(amp + ampEstimate);
430 if ( (TMath::Abs(ampAsymm) > 0.1) ) {
431 AliDebug(2,Form("Fit results amp %f time %f not consistent with expectations ped %f max-ped %f time %d",
432 amp, time, pedEstimate, ampEstimate, timeEstimate));
434 // what should do we do then? skip this channel or assign the simple estimate?
435 // for now just overwrite the fit results with the simple estimate
439 } // amp & time check
440 } // ampEstimate check
441 } // method selection
443 if (amp > fNoiseThreshold && amp<fgkRawSignalOverflow) { // something to be stored
444 Int_t id = fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
445 lowGain = in.IsLowGain();
447 // go from time-bin units to physical time fgtimetrigger
448 time = time * GetRawFormatTimeBinWidth(); // skip subtraction of fgTimeTrigger?
450 AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
451 // printf("Added tower: SM %d, row %d, column %d, amp %3.2f\n",in.GetModule(), in.GetRow(), in.GetColumn(),amp);
452 // round off amplitude value to nearest integer
453 AddDigit(digitsArr, id, lowGain, TMath::Nint(amp), time);
459 // if (maxTimeBin && gSig->GetN() > maxTimeBin + 10) gSig->Set(maxTimeBin + 10); // set actual max size of TGraph
460 Int_t hwAdd = in.GetHWAddress();
461 UShort_t iRCU = in.GetDDLNumber() % 2; // 0/1
462 UShort_t iBranch = ( hwAdd >> 11 ) & 0x1; // 0/1
464 // Now find TRU number
465 Int_t itru = 3 * in.GetModule() + ( (iRCU << 1) | iBranch ) - 1;
467 AliDebug(1,Form("Found TRG digit in TRU: %2d ADC: %2d",itru,in.GetColumn()));
471 Bool_t isOK = fGeom->GetAbsFastORIndexFromTRU(itru, in.GetColumn(), idtrg);
473 Int_t timeSamples[256]; for (Int_t j=0;j<256;j++) timeSamples[j] = 0;
476 for (std::vector<AliCaloBunchInfo>::iterator itVectorData = bunchlist.begin(); itVectorData != bunchlist.end(); itVectorData++)
478 AliCaloBunchInfo bunch = *(itVectorData);
480 const UShort_t* sig = bunch.GetData();
481 Int_t startBin = bunch.GetStartBin();
483 for (Int_t iS = 0; iS < bunch.GetLength(); iS++)
485 Int_t time = startBin--;
488 if ( amp ) timeSamples[nSamples++] = ( ( time << 12 ) & 0xFF000 ) | ( amp & 0xFFF );
492 if (nSamples && isOK) AddDigit(digitsTRG, idtrg, timeSamples, nSamples);
494 } // end while over channel
495 } //end while over DDL's, of input stream
500 //____________________________________________________________________________
501 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t timeSamples[], Int_t nSamples)
503 new((*digitsArr)[digitsArr->GetEntriesFast()]) AliEMCALRawDigit(id, timeSamples, nSamples);
505 // Int_t idx = digitsArr->GetEntriesFast()-1;
506 // AliEMCALRawDigit* d = (AliEMCALRawDigit*)digitsArr->At(idx);
509 //____________________________________________________________________________
510 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
513 // This routine checks whether a digit exists already for this tower
514 // and then decides whether to use the high or low gain info
516 // Called by Raw2Digits
518 AliEMCALDigit *digit = 0, *tmpdigit = 0;
519 TIter nextdigit(digitsArr);
520 while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
521 if (tmpdigit->GetId() == id)
525 if (!digit) { // no digit existed for this tower; create one
526 if (lowGain && amp > fgkOverflowCut)
527 amp = Int_t(fHighLowGainFactor * amp);
528 Int_t idigit = digitsArr->GetEntries();
529 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;
531 else { // a digit already exists, check range
532 // (use high gain if signal < cut value, otherwise low gain)
533 if (lowGain) { // new digit is low gain
534 if (digit->GetAmp() > fgkOverflowCut) { // use if stored digit is out of range
535 digit->SetAmp(Int_t(fHighLowGainFactor * amp));
536 digit->SetTime(time);
539 else if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
541 digit->SetTime(time);
546 //____________________________________________________________________________
547 void AliEMCALRawUtils::FitRaw(const Int_t firstTimeBin, const Int_t lastTimeBin, Float_t & amp, Float_t & time) const
548 { // Fits the raw signal time distribution
550 //--------------------------------------------------
551 //Do the fit, different fitting algorithms available
552 //--------------------------------------------------
553 int nsamples = lastTimeBin - firstTimeBin + 1;
555 switch(fFittingAlgorithm) {
558 if (nsamples < 3) { return; } // nothing much to fit
559 //printf("Standard fitter \n");
561 // Create Graph to hold data we will fit
562 TGraph *gSig = new TGraph( nsamples);
563 for (int i=0; i<nsamples; i++) {
564 Int_t timebin = firstTimeBin + i;
565 gSig->SetPoint(i, timebin, fRawAnalyzer->GetReversed(timebin));
568 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
569 signalF->SetParameters(10.,5.,fTau,fOrder,0.); //set all defaults once, just to be safe
570 signalF->SetParNames("amp","t0","tau","N","ped");
571 signalF->FixParameter(2,fTau); // tau in units of time bin
572 signalF->FixParameter(3,fOrder); // order
573 signalF->FixParameter(4, 0); // pedestal should be subtracted when we get here
574 signalF->SetParameter(1, time);
575 signalF->SetParameter(0, amp);
577 gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
579 // assign fit results
580 amp = signalF->GetParameter(0);
581 time = signalF->GetParameter(1);
585 // cross-check with ParabolaFit to see if the results make sense
586 FitParabola(gSig, amp); // amp is possibly updated
588 //printf("Std : Amp %f, time %g\n",amp, time);
589 delete gSig; // delete TGraph
593 //----------------------------
596 if (nsamples < 3) { return; } // nothing much to fit
597 //printf("LogFit \n");
599 // Create Graph to hold data we will fit
600 TGraph *gSigLog = new TGraph( nsamples);
601 for (int i=0; i<nsamples; i++) {
602 Int_t timebin = firstTimeBin + i;
603 gSigLog->SetPoint(timebin, timebin, TMath::Log(fRawAnalyzer->GetReversed(timebin) ) );
606 TF1 * signalFLog = new TF1("signalLog", RawResponseFunctionLog, 0, GetRawFormatTimeBins(), 5);
607 signalFLog->SetParameters(2.3, 5.,fTau,fOrder,0.); //set all defaults once, just to be safe
608 signalFLog->SetParNames("amplog","t0","tau","N","ped");
609 signalFLog->FixParameter(2,fTau); // tau in units of time bin
610 signalFLog->FixParameter(3,fOrder); // order
611 signalFLog->FixParameter(4, 0); // pedestal should be subtracted when we get here
612 signalFLog->SetParameter(1, time);
614 signalFLog->SetParameter(0, TMath::Log(amp));
617 gSigLog->Fit(signalFLog, "QROW"); // Note option 'W': equal errors on all points
619 // assign fit results
620 Double_t amplog = signalFLog->GetParameter(0); //Not Amp, but Log of Amp
621 amp = TMath::Exp(amplog);
622 time = signalFLog->GetParameter(1);
625 //printf("LogFit: Amp %f, time %g\n",amp, time);
629 //----------------------------
631 //----------------------------
632 }//switch fitting algorithms
637 //__________________________________________________________________
638 void AliEMCALRawUtils::FitParabola(const TGraph *gSig, Float_t & amp) const
640 //BEG YS alternative methods to calculate the amplitude
641 Double_t * ymx = gSig->GetX() ;
642 Double_t * ymy = gSig->GetY() ;
644 Double_t ymMaxX[kN] = {0., 0., 0.} ;
645 Double_t ymMaxY[kN] = {0., 0., 0.} ;
647 // find the maximum amplitude
649 for (Int_t ymi = 0; ymi < gSig->GetN(); ymi++) {
650 if (ymy[ymi] > ymMaxY[0] ) {
651 ymMaxY[0] = ymy[ymi] ; //<========== This is the maximum amplitude
652 ymMaxX[0] = ymx[ymi] ;
656 // find the maximum by fitting a parabola through the max and the two adjacent samples
657 if ( ymiMax < gSig->GetN()-1 && ymiMax > 0) {
658 ymMaxY[1] = ymy[ymiMax+1] ;
659 ymMaxY[2] = ymy[ymiMax-1] ;
660 ymMaxX[1] = ymx[ymiMax+1] ;
661 ymMaxX[2] = ymx[ymiMax-1] ;
662 if (ymMaxY[0]*ymMaxY[1]*ymMaxY[2] > 0) {
663 //fit a parabola through the 3 points y= a+bx+x*x*x
671 for (Int_t i = 0; i < kN ; i++) {
674 sx2 += ymMaxX[i]*ymMaxX[i] ;
675 sx3 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ;
676 sx4 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ;
677 sxy += ymMaxX[i]*ymMaxY[i] ;
678 sx2y += ymMaxX[i]*ymMaxX[i]*ymMaxY[i] ;
680 Double_t cN = (sx2y*kN-sy*sx2)*(sx3*sx-sx2*sx2)-(sx2y*sx-sxy*sx2)*(sx3*kN-sx*sx2);
681 Double_t cD = (sx4*kN-sx2*sx2)*(sx3*sx-sx2*sx2)-(sx4*sx-sx3*sx2)*(sx3*kN-sx*sx2) ;
682 Double_t c = cN / cD ;
683 Double_t b = ((sx2y*kN-sy*sx2)-c*(sx4*kN-sx2*sx2))/(sx3*kN-sx*sx2) ;
684 Double_t a = (sy-b*sx-c*sx2)/kN ;
685 Double_t xmax = -b/(2*c) ;
686 ymax = a + b*xmax + c*xmax*xmax ;//<========== This is the maximum amplitude
691 Double_t diff = TMath::Abs(1-ymMaxY[0]/amp) ;
694 //printf("Yves : Amp %f, time %g\n",amp, time);
699 //__________________________________________________________________
700 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
702 // Matches version used in 2007 beam test
704 // Shape of the electronics raw reponse:
705 // It is a semi-gaussian, 2nd order Gamma function of the general form
707 // xx = (t - t0 + tau) / tau [xx is just a convenient help variable]
708 // F = A * (xx**N * exp( N * ( 1 - xx) ) for xx >= 0
712 // A: par[0] // Amplitude = peak value
719 Double_t tau =par[2];
721 Double_t ped = par[4];
722 Double_t xx = ( x[0] - par[1] + tau ) / tau ;
727 signal = ped + par[0] * TMath::Power(xx , n) * TMath::Exp(n * (1 - xx )) ;
732 //__________________________________________________________________
733 Double_t AliEMCALRawUtils::RawResponseFunctionLog(Double_t *x, Double_t *par)
735 // Matches version used in 2007 beam test
737 // Shape of the electronics raw reponse:
738 // It is a semi-gaussian, 2nd order Gamma function of the general form
740 // xx = (t - t0 + tau) / tau [xx is just a convenient help variable]
741 // F = A * (xx**N * exp( N * ( 1 - xx) ) for xx >= 0
745 // Log[A]: par[0] // Amplitude = peak value
752 Double_t tau =par[2];
754 //Double_t ped = par[4]; // not used
755 Double_t xx = ( x[0] - par[1] + tau ) / tau ;
758 signal = par[0] - n*TMath::Log(TMath::Abs(xx)) + n * (1 - xx ) ;
760 signal = par[0] + n*TMath::Log(xx) + n * (1 - xx ) ;
765 //__________________________________________________________________
766 Bool_t AliEMCALRawUtils::RawSampledResponse(const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL, const Int_t keyErr) const
768 // for a start time dtime and an amplitude damp given by digit,
769 // calculates the raw sampled response AliEMCAL::RawResponseFunction
771 Bool_t lowGain = kFALSE ;
773 // A: par[0] // Amplitude = peak value
779 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
780 signalF.SetParameter(0, damp) ;
781 signalF.SetParameter(1, (dtime + fgTimeTrigger)/fgTimeBinWidth) ;
782 signalF.SetParameter(2, fTau) ;
783 signalF.SetParameter(3, fOrder);
784 signalF.SetParameter(4, fgPedestalValue);
786 Double_t signal=0.0, noise=0.0;
787 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
788 signal = signalF.Eval(iTime) ;
790 // Next lines commeted for the moment but in principle it is not necessary to add
791 // extra noise since noise already added at the digits level.
793 //According to Terry Awes, 13-Apr-2008
794 //add gaussian noise in quadrature to each sample
795 //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
796 //signal = sqrt(signal*signal + noise*noise);
798 // March 17,09 for fast fit simulations by Alexei Pavlinov.
799 // Get from PHOS analysis. In some sense it is open questions.
801 noise = gRandom->Gaus(0.,fgFEENoise);
805 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
806 if ( adcH[iTime] > fgkRawSignalOverflow ){ // larger than 10 bits
807 adcH[iTime] = fgkRawSignalOverflow ;
811 signal /= fHighLowGainFactor;
813 adcL[iTime] = static_cast<Int_t>(signal + 0.5) ;
814 if ( adcL[iTime] > fgkRawSignalOverflow) // larger than 10 bits
815 adcL[iTime] = fgkRawSignalOverflow ;
820 //__________________________________________________________________
821 void AliEMCALRawUtils::SetFittingAlgorithm(Int_t fitAlgo)
823 //Set fitting algorithm and initialize it if this same algorithm was not set before.
824 //printf("**** Set Algorithm , number %d ****\n",fitAlgo);
826 if(fitAlgo == fFittingAlgorithm && fRawAnalyzer) {
827 //Do nothing, this same algorithm already set before.
828 //printf("**** Algorithm already set before, number %d, %s ****\n",fitAlgo, fRawAnalyzer->GetName());
831 //Initialize the requested algorithm
832 if(fitAlgo != fFittingAlgorithm || !fRawAnalyzer) {
833 //printf("**** Init Algorithm , number %d ****\n",fitAlgo);
835 fFittingAlgorithm = fitAlgo;
836 if (fRawAnalyzer) delete fRawAnalyzer; // delete prev. analyzer if existed.
838 if (fitAlgo == kFastFit) {
839 fRawAnalyzer = new AliCaloRawAnalyzerFastFit();
841 else if (fitAlgo == kNeuralNet) {
842 fRawAnalyzer = new AliCaloRawAnalyzerNN();
844 else if (fitAlgo == kLMS) {
845 fRawAnalyzer = new AliCaloRawAnalyzerLMS();
847 else if (fitAlgo == kPeakFinder) {
848 fRawAnalyzer = new AliCaloRawAnalyzerPeakFinder();
850 else if (fitAlgo == kCrude) {
851 fRawAnalyzer = new AliCaloRawAnalyzerCrude();
854 fRawAnalyzer = new AliCaloRawAnalyzer();