AliRoot Core  3dc7879 (3dc7879)
AliCaloRawAnalyzerFastFit.cxx
Go to the documentation of this file.
1 // -*- mode: c++ -*-
2 /**************************************************************************
3  * This file is property of and copyright by the Experimental Nuclear *
4  * Physics Group, Dep. of Physics *
5  * University of Oslo, Norway, 2007 *
6  * *
7  * Author: Per Thomas Hille <perthi@fys.uio.no> for the ALICE HLT Project.*
8  * Contributors are mentioned in the code where appropriate. *
9  * Please report bugs to perthi@fys.uio.no *
10  * *
11  * Permission to use, copy, modify and distribute this software and its *
12  * documentation strictly for non-commercial purposes is hereby granted *
13  * without fee, provided that the above copyright notice appears in all *
14  * copies and that both the copyright notice and this permission notice *
15  * appear in the supporting documentation. The authors make no claims *
16  * about the suitability of this software for any purpose. It is *
17  * provided "as is" without express or implied warranty. *
18  **************************************************************************/
19 
20 // EMCal system
22 #include "AliCaloFastAltroFitv0.h"
23 #include "AliCaloFitResults.h"
24 #include "AliCaloBunchInfo.h"
25 #include "AliCaloConstants.h"
26 
27 // Root
28 #include "TMath.h"
29 
30 // Standard libraries
31 #include <iostream>
32 using namespace std;
33 
34 
36 ClassImp( AliCaloRawAnalyzerFastFit ) ;
38 
41 //_______________________________________________________________________
43 {
45 }
46 
49 //_______________________________________________________________________
51 AliCaloRawAnalyzerFastFit::Evaluate( const vector<AliCaloBunchInfo> &bunchvector,
52  UInt_t altrocfg1, UInt_t altrocfg2 )
53 {
54  short maxampindex; //index of maximum amplitude
55  short maxamp; //Maximum amplitude
56  int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
57 
58  if( index >= 0)
59  {
60  Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
61  Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
62  short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
63  Float_t time = (timebinOffset*TIMEBINWITH)-fL1Phase;
64 
65  if( maxf < fAmpCut || maxamp > fOverflowCut ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
66  {
67  return AliCaloFitResults( maxamp, ped, Algo::kCrude, maxf, time, (int)time, 0, 0, Ret::kDummy); //Time scale 19/08/2014 (Antônio)
68  }
69  else if ( maxf >= fAmpCut ) // no if statement needed really; keep for readability
70  {
71  int first = 0;
72  int last = 0;
73  int maxrev = maxampindex - bunchvector.at(index).GetStartBin();
74 
75  SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev , &first, &last, fFitArrayCut);
76 
77  int nsamples = last - first + 1;
78 
79  if( ( nsamples ) >= fNsampleCut )
80  {
81  Double_t ordered[1008];
82 
83  for(int i=0; i < nsamples ; i++ )
84  {
85  ordered[i] = fReversed[first + i];
86  }
87 
88  Double_t eSignal = 1; // nominal 1 ADC error
89  Double_t dAmp = maxf;
90  Double_t eAmp = 0;
91  Double_t dTime0 = 0;
92  Double_t eTime = 0;
93  Double_t chi2 = 0;
94  Double_t dTau = 2.35; // time-bin units
95 
96  AliCaloFastAltroFitv0::FastFit(fXaxis, ordered , nsamples,
97  eSignal, dTau, dAmp, eAmp, dTime0, eTime, chi2);
98  Double_t dTimeMax = dTime0 + timebinOffset - (maxrev - first) // abs. t0
99  + dTau; // +tau, makes sum tmax
100  Float_t timemax = (dTimeMax*TIMEBINWITH)-fL1Phase;
101  return AliCaloFitResults(maxamp,ped,Ret::kFitPar,dAmp,timemax,(int)timemax,chi2,Ret::kDummy,Ret::kDummy,AliCaloFitSubarray(index,maxrev,first,last)); //Time scale 19/08/2014 (Antônio)
102  } // samplecut
103  else
104  {
105 
106  Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
107  Int_t ndf = last - first - 1; // nsamples - 2
108  return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, time, (int)time, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) ); //Time scale 19/08/2014 (Antônio)
109  }
110  } // ampcut
111  } // bunch index
112 
114 }
int fFitArrayCut
Cut on ADC value (after ped. subtraction) for signals used for fit.
Raw data fitting: special fast fit.
Float_t fAmpCut
Max ADC - pedestal must be higher than this befor attemting to extract the amplitude.
int fOverflowCut
Value when ADC starts to saturate.
int fNsampleCut
Minimum number of sample require before attemting to extract signal parameters.
virtual void FastFit(Int_t *t, Int_t *y, Int_t nPoints, Double_t sig, Double_t tau, Double_t n, Double_t ped, Double_t tMax)
Fast fit.
Container class to hold info from bunches/samples.
int SelectBunch(const std::vector< AliCaloBunchInfo > &bunchvector, short *maxampbin, short *maxamplitude)
We select the bunch with the highest amplitude unless any time constraints is set.
const double TIMEBINWITH
each sample is 100 ns
Container class to hold results from fitting.
Double_t chi2
Definition: AnalyzeLaser.C:7
Raw data fitters base class.
double fXaxis[ALTROMAXSAMPLES]
Axis if time bins, ( used by TGraph )
Algo::fitAlgorithm fAlgo
Which algorithm to use.
void SelectSubarray(const Double_t *date, int length, short maxindex, int *first, int *last, int cut) const
Float_t ReverseAndSubtractPed(const AliCaloBunchInfo *bunch, UInt_t altrocfg1, UInt_t altrocfg2, double *outarray) const
Double_t fReversed[ALTROMAXSAMPLES]
Reversed sequence of samples (pedestalsubtracted)
Double_t CalculateChi2(const Double_t amp, const Double_t time, const Int_t first, const Int_t last, const Double_t adcErr=1, const Double_t tau=2.35) const
virtual AliCaloFitResults Evaluate(const std::vector< AliCaloBunchInfo > &bunchvector, UInt_t altrocfg1, UInt_t altrocfg2)
Execute algorithm.
Double_t fL1Phase
Phase of the ADC sampling clock relative to the LHC clock.