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"
37 #include "AliRunLoader.h"
38 class AliCaloAltroMapping;
39 #include "AliAltroBuffer.h"
40 #include "AliRawReader.h"
41 #include "AliCaloRawStreamV3.h"
44 #include "AliEMCALRecParam.h"
45 #include "AliEMCALLoader.h"
46 #include "AliEMCALGeometry.h"
47 class AliEMCALDigitizer;
48 #include "AliEMCALDigit.h"
50 #include "AliCaloCalibPedestal.h"
51 #include "AliCaloFastAltroFitv0.h"
53 ClassImp(AliEMCALRawUtils)
55 // Signal shape parameters
56 Int_t AliEMCALRawUtils::fgTimeBins = 256; // number of sampling bins of the raw RO signal (we typically use 15-50; theoretical max is 1k+)
57 Double_t AliEMCALRawUtils::fgTimeBinWidth = 100E-9 ; // each sample is 100 ns
58 Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ; // 15 time bins ~ 1.5 musec
60 // some digitization constants
61 Int_t AliEMCALRawUtils::fgThreshold = 1;
62 Int_t AliEMCALRawUtils::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
63 Int_t AliEMCALRawUtils::fgPedestalValue = 32; // pedestal value for digits2raw
64 Double_t AliEMCALRawUtils::fgFEENoise = 3.; // 3 ADC channels of noise (sampled)
66 AliEMCALRawUtils::AliEMCALRawUtils()
67 : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
68 fNPedSamples(0), fGeom(0), fOption(""),
69 fRemoveBadChannels(kTRUE),fFittingAlgorithm(0)
72 //These are default parameters.
73 //Can be re-set from without with setter functions
74 //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
75 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
76 fOrder = 2; // order of gamma fn
77 fTau = 2.35; // in units of timebin, from CERN 2007 testbeam
78 fNoiseThreshold = 3; // 3 ADC counts is approx. noise level
79 fNPedSamples = 4; // less than this value => likely pedestal samples
80 fRemoveBadChannels = kTRUE; //Remove bad channels before fitting
81 fFittingAlgorithm = kFastFit;//kStandard; // Use default minuit fitter
83 //Get Mapping RCU files from the AliEMCALRecParam
84 const TObjArray* maps = AliEMCALRecParam::GetMappings();
85 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
87 for(Int_t i = 0; i < 4; i++) {
88 fMapping[i] = (AliAltroMapping*)maps->At(i);
91 //To make sure we match with the geometry in a simulation file,
92 //let's try to get it first. If not, take the default geometry
93 AliRunLoader *rl = AliRunLoader::Instance();
94 if(!rl) AliError("Cannot find RunLoader!");
95 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
96 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
98 AliInfo(Form("Using default geometry in raw reco"));
99 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
102 if(!fGeom) AliFatal(Form("Could not get geometry!"));
106 //____________________________________________________________________________
107 AliEMCALRawUtils::AliEMCALRawUtils(AliEMCALGeometry *pGeometry)
108 : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
109 fNPedSamples(0), fGeom(pGeometry), fOption(""),
110 fRemoveBadChannels(kTRUE),fFittingAlgorithm(0)
113 // Initialize with the given geometry - constructor required by HLT
114 // HLT does not use/support AliRunLoader(s) instances
115 // This is a minimum intervention solution
116 // Comment by MPloskon@lbl.gov
119 //These are default parameters.
120 //Can be re-set from without with setter functions
121 //Already set in the OCDB and passed via setter in the AliEMCALReconstructor
122 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
123 fOrder = 2; // order of gamma fn
124 fTau = 2.35; // in units of timebin, from CERN 2007 testbeam
125 fNoiseThreshold = 3; // 3 ADC counts is approx. noise level
126 fNPedSamples = 4; // less than this value => likely pedestal samples
127 fRemoveBadChannels = kTRUE; //Remove bad channels before fitting
128 fFittingAlgorithm = kStandard; // Use default minuit fitter
130 //Get Mapping RCU files from the AliEMCALRecParam
131 const TObjArray* maps = AliEMCALRecParam::GetMappings();
132 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
134 for(Int_t i = 0; i < 4; i++) {
135 fMapping[i] = (AliAltroMapping*)maps->At(i);
138 if(!fGeom) AliFatal(Form("Could not get geometry!"));
142 //____________________________________________________________________________
143 AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
145 fHighLowGainFactor(rawU.fHighLowGainFactor),
148 fNoiseThreshold(rawU.fNoiseThreshold),
149 fNPedSamples(rawU.fNPedSamples),
151 fOption(rawU.fOption),
152 fRemoveBadChannels(rawU.fRemoveBadChannels),
153 fFittingAlgorithm(rawU.fFittingAlgorithm)
156 fMapping[0] = rawU.fMapping[0];
157 fMapping[1] = rawU.fMapping[1];
158 fMapping[2] = rawU.fMapping[2];
159 fMapping[3] = rawU.fMapping[3];
162 //____________________________________________________________________________
163 AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
165 //assignment operator
168 fHighLowGainFactor = rawU.fHighLowGainFactor;
169 fOrder = rawU.fOrder;
171 fNoiseThreshold = rawU.fNoiseThreshold;
172 fNPedSamples = rawU.fNPedSamples;
174 fOption = rawU.fOption;
175 fRemoveBadChannels = rawU.fRemoveBadChannels;
176 fFittingAlgorithm = rawU.fFittingAlgorithm;
177 fMapping[0] = rawU.fMapping[0];
178 fMapping[1] = rawU.fMapping[1];
179 fMapping[2] = rawU.fMapping[2];
180 fMapping[3] = rawU.fMapping[3];
187 //____________________________________________________________________________
188 AliEMCALRawUtils::~AliEMCALRawUtils() {
193 //____________________________________________________________________________
194 void AliEMCALRawUtils::Digits2Raw()
196 // convert digits of the current event to raw data
198 AliRunLoader *rl = AliRunLoader::Instance();
199 AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
202 loader->LoadDigits("EMCAL");
204 TClonesArray* digits = loader->Digits() ;
207 Warning("Digits2Raw", "no digits found !");
211 static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
212 AliAltroBuffer* buffers[nDDL];
213 for (Int_t i=0; i < nDDL; i++)
216 TArrayI adcValuesLow(fgTimeBins);
217 TArrayI adcValuesHigh(fgTimeBins);
219 // loop over digits (assume ordered digits)
220 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
221 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
222 if (digit->GetAmp() < fgThreshold)
232 fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
233 fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
235 //Check which is the RCU, 0 or 1, of the cell.
238 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
239 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
242 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
244 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
246 if (nSM%2==1) iRCU = 1 - iRCU; // swap for odd=C side, to allow us to cable both sides the same
249 Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
252 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
254 Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
256 if (buffers[iDDL] == 0) {
257 // open new file and write dummy header
258 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
259 //Select mapping file RCU0A, RCU0C, RCU1A, RCU1C
260 Int_t iRCUside=iRCU+(nSM%2)*2;
261 //iRCU=0 and even (0) SM -> RCU0A.data 0
262 //iRCU=1 and even (0) SM -> RCU1A.data 1
263 //iRCU=0 and odd (1) SM -> RCU0C.data 2
264 //iRCU=1 and odd (1) SM -> RCU1C.data 3
265 //cout<<" nSM "<<nSM<<"; iRCU "<<iRCU<<"; iRCUside "<<iRCUside<<endl;
266 buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCUside]);
267 buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE); //Dummy;
270 // out of time range signal (?)
271 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
272 AliInfo("Signal is out of time range.\n");
273 buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
274 buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() ); // time bin
275 buffers[iDDL]->FillBuffer(3); // bunch length
276 buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM); // trailer
277 // calculate the time response function
279 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh.GetArray(), adcValuesLow.GetArray()) ;
281 buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow.GetArray(), fgThreshold);
283 buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh.GetArray(), fgThreshold);
287 // write headers and close files
288 for (Int_t i=0; i < nDDL; i++) {
291 buffers[i]->WriteDataHeader(kFALSE, kFALSE);
296 loader->UnloadDigits();
299 //____________________________________________________________________________
300 void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr, AliCaloCalibPedestal* pedbadmap)
302 // convert raw data of the current event to digits
307 Error("Raw2Digits", "no digits found !");
311 Error("Raw2Digits", "no raw reader found !");
315 AliCaloRawStreamV3 in(reader,"EMCAL",fMapping);
316 // Select EMCAL DDL's;
317 reader->Select("EMCAL",0,43); // 43 = AliEMCALGeoParams::fgkLastAltroDDL
319 //Updated fitting routine from 2007 beam test takes into account
320 //possibility of two peaks in data and selects first one for fitting
321 //Also sets some of the starting parameters based on the shape of the
322 //given raw signal being fit
324 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
325 signalF->SetParameters(10.,5.,fTau,fOrder,0.); //set all defaults once, just to be safe
326 signalF->SetParNames("amp","t0","tau","N","ped");
327 signalF->FixParameter(2,fTau); // tau in units of time bin
328 signalF->FixParameter(3,fOrder); // order
334 Float_t ampEstimate = 0;
335 Float_t timeEstimate = 0;
336 Float_t pedEstimate = 0;
340 //Graph to hold data we will fit (should be converted to an array
341 //later to speed up processing
342 TGraph * gSig = new TGraph(GetRawFormatTimeBins());
345 Int_t caloFlag = 0; // low, high gain, or TRU, or LED ref.
347 // start loop over input stream
348 while (in.NextDDL()) {
349 while (in.NextChannel()) {
351 //Check if the signal is high or low gain and then do the fit,
352 //if it is from TRU do not fit
353 caloFlag = in.GetCaloFlag();
354 if (caloFlag != 0 && caloFlag != 1) continue;
356 //Do not fit bad channels
357 if(fRemoveBadChannels && pedbadmap->IsBadChannel(in.GetModule(),in.GetColumn(),in.GetRow())) {
358 //printf("Tower from SM %d, column %d, row %d is BAD!!! Skip \n", in.GetModule(),in.GetColumn(),in.GetRow());
362 // There can be zero-suppression in the raw data,
363 // so set up the TGraph in advance
364 for (i=0; i < GetRawFormatTimeBins(); i++) {
365 gSig->SetPoint(i, i , -1); // init to out-of-range values
368 Int_t maxTimeBin = 0;
369 Int_t min = 0x3ff; // init to 10-bit max
370 Int_t max = 0; // init to 10-bit min
371 while (in.NextBunch()) {
373 const UShort_t *sig = in.GetSignals();
374 startBin = in.GetStartTimeBin();
375 if (maxTimeBin < startBin) {
376 maxTimeBin = startBin; // timebins come in reverse order
378 if (maxTimeBin < 0 || maxTimeBin >= GetRawFormatTimeBins()) {
379 AliWarning(Form("Invalid time bin %d",maxTimeBin));
380 maxTimeBin = GetRawFormatTimeBins();
383 for (i = 0; i < in.GetBunchLength(); i++) {
385 gSig->SetPoint((Int_t)time, time, (Double_t) sig[i]) ;
386 if (max < sig[i]) max = sig[i];
387 if (min > sig[i]) min = sig[i];
390 } // loop over bunches
392 gSig->Set(maxTimeBin+1); // set actual max size of TGraph
394 //Initialize the variables, do not keep previous values.
395 // not really necessary to reset all of them (only amp and time at the moment), but better safe than sorry
403 if ( (max - min) > fNoiseThreshold) {
404 FitRaw(gSig, signalF, maxTimeBin, amp, time, ped,
405 ampEstimate, timeEstimate, pedEstimate);
406 // switch(fFittingAlgorithm)
410 // //printf("Standard fitter \n");
411 // FitRaw(gSig, signalF, maxTimeBin, amp, time, ped,
412 // ampEstimate, timeEstimate, pedEstimate);
417 // //printf("FastFitter \n");
418 // Double_t eSignal = 0;
419 // Double_t dAmp = amp;
420 // Double_t dTimeEstimate = timeEstimate;
421 // Double_t eTimeEstimate = 0;
422 // Double_t eAmp = 0;
423 // Double_t chi2 = 0;
425 // AliCaloFastAltroFitv0::FastFit(gSig->GetX(), gSig->GetY(), gSig->GetN(),
427 // dAmp, eAmp, dTimeEstimate, eTimeEstimate, chi2);
429 // timeEstimate = dTimeEstimate;
430 // //printf("FastFitter: Amp %f, time %f, eAmp %f, eTimeEstimate %f, chi2 %f\n",amp, timeEstimate,eAmp,eTimeEstimate,chi2);
437 if ( amp>0 && amp<2000 && time>0 && time<(maxTimeBin*GetRawFormatTimeBinWidth()) ) { //check both high and low end of amplitude result, and time
438 //2000 is somewhat arbitrary - not nice with magic numbers in the code..
439 id = fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
440 lowGain = in.IsLowGain();
442 // check fit results: should be consistent with initial estimates
443 // more magic numbers, but very loose cuts, for now..
444 // We have checked that amp an time values are positive so division for assymmetry
445 // calculation should be OK/safe
446 Float_t ampAsymm = (amp - ampEstimate)/(amp + ampEstimate);
447 if ( (TMath::Abs(ampAsymm) > 0.1) ||
448 (TMath::Abs(time - timeEstimate) > 2*GetRawFormatTimeBinWidth()) ) {
449 AliDebug(2,Form("Fit results ped %f amp %f time %f not consistent with expectations ped %f max-ped %f time %d",
450 ped, amp, time, pedEstimate, ampEstimate, timeEstimate));
452 // what should do we do then? skip this channel or assign the simple estimate?
453 // for now just overwrite the fit results with the simple estimate
458 AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
459 // printf("Added tower: SM %d, row %d, column %d, amp %3.2f\n",in.GetModule(), in.GetRow(), in.GetColumn(),amp);
460 // round off amplitude value to nearest integer
461 AddDigit(digitsArr, id, lowGain, TMath::Nint(amp), time);
465 for (Int_t index = 0; index < gSig->GetN(); index++) {
466 gSig->SetPoint(index, index, -1) ;
468 // Reset starting parameters for fit function
469 signalF->SetParameters(10.,5.,fTau,fOrder,0.); //reset all defaults just to be safe
471 } // end while over channel
472 } //end while over DDL's, of input stream
480 //____________________________________________________________________________
481 void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
484 // This routine checks whether a digit exists already for this tower
485 // and then decides whether to use the high or low gain info
487 // Called by Raw2Digits
489 AliEMCALDigit *digit = 0, *tmpdigit = 0;
490 TIter nextdigit(digitsArr);
491 while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
492 if (tmpdigit->GetId() == id)
496 if (!digit) { // no digit existed for this tower; create one
497 if (lowGain && amp > fgkOverflowCut)
498 amp = Int_t(fHighLowGainFactor * amp);
499 Int_t idigit = digitsArr->GetEntries();
500 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;
502 else { // a digit already exists, check range
503 // (use high gain if signal < cut value, otherwise low gain)
504 if (lowGain) { // new digit is low gain
505 if (digit->GetAmp() > fgkOverflowCut) { // use if stored digit is out of range
506 digit->SetAmp(Int_t(fHighLowGainFactor * amp));
507 digit->SetTime(time);
510 else if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
512 digit->SetTime(time);
517 //____________________________________________________________________________
518 void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, const Int_t lastTimeBin, Float_t & amp, Float_t & time, Float_t & ped, Float_t & ampEstimate, Float_t & timeEstimate, Float_t & pedEstimate, const Float_t cut) const
520 // Fits the raw signal time distribution; from AliEMCALGetter
521 // last argument: Float_t cut = 0.0; // indicating how much of amplitude w.r.t. max value fit should be above noise and pedestal
523 // initialize return values
531 // 0th step: remove plateau / overflow candidates
532 // before trying to estimate amplitude, search for maxima etc.
534 Int_t nOrig = gSig->GetN(); // number of samples before we remove any overflows
535 // Values for readback from input graph
540 // start: tmp dump of all values
541 for (Int_t i=0; i<gSig->GetN(); i++) {
542 gSig->GetPoint(i, ttime, signal) ; // get values
543 printf("orig: i %d, time %f, signal %f\n",i, ttime, signal);
545 // end: tmp dump of all values
548 // start from back of TGraph since RemovePoint will downshift indices
549 for (Int_t i=nOrig-1; i>=0; i--) {
550 gSig->GetPoint(i, ttime, signal) ; // get values
551 if (signal >= (pedEstimate + fgkOverflowCut) ) {
552 gSig->RemovePoint(i);
556 // 1st step: we try to estimate the pedestal value
558 for (Int_t index = 0; index < gSig->GetN(); index++) {
559 gSig->GetPoint(index, ttime, signal) ;
560 // ttime < fNPedsamples used for pedestal estimate;
561 // ttime >= fNPedSamples used for signal checks
562 if (signal >= 0 && ttime<fNPedSamples) { // valid value
563 pedEstimate += signal;
571 //AliWarning("Could not determine pedestal");
572 AliDebug(1,"Could not determine pedestal");
573 pedEstimate = 0; // good estimate for ZeroSupp data (non ZS data should have no problem with pedestal estimate)
576 // 2nd step: we look through the rest of the time-bins/ADC values and
577 // see if we have something that looks like a signal.
578 // We look for a first local maxima, as well as for a global maxima
579 Int_t locMaxFound = 0;
580 Int_t locMaxId = 0; // time-bin index at first local max
581 Float_t locMaxSig = -1; // actual local max value
582 Int_t globMaxId = 0; // time-bin index at global max
583 Float_t globMaxSig = -1; // actual global max value
584 // We will also look for any values that look like they are in overflow region
585 for (Int_t i=0; i<gSig->GetN(); i++) {
586 gSig->GetPoint(i, ttime, signal) ; // get values
588 // ttime < fNPedsamples used for pedestal estimate;
589 // ttime >= fNPedSamples used for signal checks
590 if (ttime >= fNPedSamples) {
592 // look for first local maximum signal=ADC value
593 if (!locMaxFound && signal > locMaxSig) {
597 else if ( locMaxSig > (pedEstimate + fNoiseThreshold) ) {
598 // we enter this condition after signal<=max, but previous
599 // max value was large enough. I.e. at least a significant local
600 // maxima has been found (just before)
604 // also check for global maximum..
605 if (signal > globMaxSig) {
610 } // end for-loop over samples after pedestal
612 // OK, we have looked through the signal spectra, let's see if we should try to make the fit
613 ampEstimate = locMaxSig - pedEstimate; // estimate using first local maxima
614 if ( ampEstimate > fNoiseThreshold ) { // else it's just noise
616 //Check that the local maximum we will use is not at the end or beginning of time sample range
617 Double_t timeMax = -1;
618 Int_t iMax = locMaxId;
619 gSig->GetPoint(locMaxId, timeMax, signal) ;
620 if (timeMax < 2 || timeMax > lastTimeBin-1) { // lastTimeBin is the lowest kept time-sample; current (Dec 2009) case
621 // if (timeMax < 2 || timeMax > lastTimeBin-2) { // for when lastTimeBin is the lowest read-out time-sample, future (2010) case
622 AliDebug(1,Form("Skip fit, maximum of the sample close to the edges : timeMax %3.2f, ampEstimate %3.2f",timeMax, ampEstimate));
626 // Check if the local and global maximum disagree
627 if (locMaxId != globMaxId) {
628 AliDebug(1,Form("Warning, local first maximum %d does not agree with global maximum %d\n", locMaxId, globMaxId));
632 // Get the maximum and find the lowest timebin (tailmin) where the ADC value is not
633 // significantly different from the pedestal
634 // first lower times edge a.k.a. tailmin
636 Double_t tmptime = 0;
637 for (Int_t i=iMax-1; i > 0; i--) {
638 gSig->GetPoint(i, tmptime, signal) ;
639 if((signal-pedEstimate) < fNoiseThreshold){
644 // then same exercise for the higher times edge a.k.a. tailmax
645 Int_t tailMax = lastTimeBin;
646 for (Int_t i=iMax+1; i < gSig->GetN(); i++) {
647 gSig->GetPoint(i, tmptime, signal) ;
648 if ((signal-pedEstimate) <= (ampEstimate*cut + fNoiseThreshold)) { // stop fit at cut-fraction of amplitude above noise-threshold (cut>0 would mean avoid the pulse shape falling edge)
654 // remove all points which are not in the distribution around maximum
655 // i.e. up to tailmin, and from tailmax
656 if ( tailMax != (gSig->GetN()-1) ){ // else nothing to remove
657 nOrig = gSig->GetN(); // can't use GetN call in for loop below since gSig size changes..
658 for(int j = tailMax; j < nOrig; j++) gSig->RemovePoint(tailMax);
660 for(int j = 0; j<=tailMin; j++) gSig->RemovePoint(0);
662 if(gSig->GetN() < 3) {
663 AliDebug(2,Form("Skip fit, number of entries in sample smaller than number of fitting parameters: in sample %d, fitting param 3",
668 timeEstimate = timeMax * GetRawFormatTimeBinWidth();
670 //--------------------------------------------------
671 //Do the fit, different fitting algorithms available
672 //--------------------------------------------------
674 switch(fFittingAlgorithm) {
677 //printf("Standard fitter \n");
679 // determine what the valid fit range is
680 Double_t minFit = 9999;
682 for (Int_t i=0; i < gSig->GetN(); i++) {
683 gSig->GetPoint(i, ttime, signal);
684 if (minFit > ttime) minFit=ttime;
685 if (maxFit < ttime) maxFit=ttime;
686 //debug: printf("no tail: i %d, time %f, signal %f\n",i, ttime, signal);
688 signalF->SetRange(minFit, maxFit);
690 signalF->FixParameter(4, pedEstimate) ;
691 signalF->SetParameter(1, timeMax);
692 signalF->SetParameter(0, ampEstimate);
694 gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
696 // assign fit results
697 amp = signalF->GetParameter(0);
698 time = signalF->GetParameter(1) * GetRawFormatTimeBinWidth(); // skip subtraction of fgTimeTrigger?
699 ped = signalF->GetParameter(4);
700 //printf("Std : Amp %f, time %g\n",amp, time);
702 //BEG YS alternative methods to calculate the amplitude
703 Double_t * ymx = gSig->GetX() ;
704 Double_t * ymy = gSig->GetY() ;
706 Double_t ymMaxX[kN] = {0., 0., 0.} ;
707 Double_t ymMaxY[kN] = {0., 0., 0.} ;
709 // find the maximum amplitude
711 for (Int_t ymi = 0; ymi < gSig->GetN(); ymi++) {
712 if (ymy[ymi] > ymMaxY[0] ) {
713 ymMaxY[0] = ymy[ymi] ; //<========== This is the maximum amplitude
714 ymMaxX[0] = ymx[ymi] ;
718 // find the maximum by fitting a parabola through the max and the two adjacent samples
719 if ( ymiMax < gSig->GetN()-1 && ymiMax > 0) {
720 ymMaxY[1] = ymy[ymiMax+1] ;
721 ymMaxY[2] = ymy[ymiMax-1] ;
722 ymMaxX[1] = ymx[ymiMax+1] ;
723 ymMaxX[2] = ymx[ymiMax-1] ;
724 if (ymMaxY[0]*ymMaxY[1]*ymMaxY[2] > 0) {
725 //fit a parabola through the 3 points y= a+bx+x*x*x
733 for (Int_t i = 0; i < kN ; i++) {
736 sx2 += ymMaxX[i]*ymMaxX[i] ;
737 sx3 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ;
738 sx4 += ymMaxX[i]*ymMaxX[i]*ymMaxX[i]*ymMaxX[i] ;
739 sxy += ymMaxX[i]*ymMaxY[i] ;
740 sx2y += ymMaxX[i]*ymMaxX[i]*ymMaxY[i] ;
742 Double_t cN = (sx2y*kN-sy*sx2)*(sx3*sx-sx2*sx2)-(sx2y*sx-sxy*sx2)*(sx3*kN-sx*sx2);
743 Double_t cD = (sx4*kN-sx2*sx2)*(sx3*sx-sx2*sx2)-(sx4*sx-sx3*sx2)*(sx3*kN-sx*sx2) ;
744 Double_t c = cN / cD ;
745 Double_t b = ((sx2y*kN-sy*sx2)-c*(sx4*kN-sx2*sx2))/(sx3*kN-sx*sx2) ;
746 Double_t a = (sy-b*sx-c*sx2)/kN ;
747 Double_t xmax = -b/(2*c) ;
748 ymax = a + b*xmax + c*xmax*xmax ;//<========== This is the maximum amplitude
752 Double_t diff = TMath::Abs(1-ymMaxY[0]/amp) ;
755 //printf("Yves : Amp %f, time %g\n",amp, time);
759 //----------------------------
762 //printf("FastFitter \n");
763 Double_t eSignal = 0;
766 Double_t dTime = time;
770 AliCaloFastAltroFitv0::FastFit(gSig->GetX(), gSig->GetY(), gSig->GetN(),
772 dAmp, eAmp, dTime, eTime, chi2);
774 time = dTime * GetRawFormatTimeBinWidth();
775 //printf("FastFitter: Amp %f, time %g, eAmp %f, eTimeEstimate %g, chi2 %f\n",amp, time,eAmp,eTime,chi2);
778 //----------------------------
779 }//switch fitting algorithms
780 } // ampEstimate > fNoiseThreshold
783 //__________________________________________________________________
784 Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
786 // Matches version used in 2007 beam test
788 // Shape of the electronics raw reponse:
789 // It is a semi-gaussian, 2nd order Gamma function of the general form
791 // xx = (t - t0 + tau) / tau [xx is just a convenient help variable]
792 // F = A * (xx**N * exp( N * ( 1 - xx) ) for xx >= 0
796 // A: par[0] // Amplitude = peak value
803 Double_t tau =par[2];
805 Double_t ped = par[4];
806 Double_t xx = ( x[0] - par[1] + tau ) / tau ;
811 signal = ped + par[0] * TMath::Power(xx , n) * TMath::Exp(n * (1 - xx )) ;
816 //__________________________________________________________________
817 Bool_t AliEMCALRawUtils::RawSampledResponse(
818 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
820 // for a start time dtime and an amplitude damp given by digit,
821 // calculates the raw sampled response AliEMCAL::RawResponseFunction
823 Bool_t lowGain = kFALSE ;
825 // A: par[0] // Amplitude = peak value
831 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
832 signalF.SetParameter(0, damp) ;
833 signalF.SetParameter(1, (dtime + fgTimeTrigger)/fgTimeBinWidth) ;
834 signalF.SetParameter(2, fTau) ;
835 signalF.SetParameter(3, fOrder);
836 signalF.SetParameter(4, fgPedestalValue);
838 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
839 Double_t signal = signalF.Eval(iTime) ;
841 // Next lines commeted for the moment but in principle it is not necessary to add
842 // extra noise since noise already added at the digits level.
844 //According to Terry Awes, 13-Apr-2008
845 //add gaussian noise in quadrature to each sample
846 //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
847 //signal = sqrt(signal*signal + noise*noise);
849 // March 17,09 for fast fit simulations by Alexei Pavlinov.
850 // Get from PHOS analysis. In some sense it is open questions.
851 //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
854 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
855 if ( adcH[iTime] > fgkRawSignalOverflow ){ // larger than 10 bits
856 adcH[iTime] = fgkRawSignalOverflow ;
860 signal /= fHighLowGainFactor;
862 adcL[iTime] = static_cast<Int_t>(signal + 0.5) ;
863 if ( adcL[iTime] > fgkRawSignalOverflow) // larger than 10 bits
864 adcL[iTime] = fgkRawSignalOverflow ;