AliRoot Core  3dc7879 (3dc7879)
AliCaloRawAnalyzerFakeALTRO.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
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  **************************************************************************/
15 
16 // AliRoot/EMCal system
18 #include "AliCaloBunchInfo.h"
19 #include "AliCaloFitResults.h"
20 #include "AliCaloConstants.h"
21 #include "AliLog.h"
22 
23 // ROOT system
24 #include "TMath.h"
25 #include "TF1.h"
26 #include "TGraph.h"
27 
28 // Standard libraries
29 #include <stdexcept>
30 #include <iostream>
31 
32 using namespace std;
33 
35 ClassImp( AliCaloRawAnalyzerFakeALTRO ) ;
37 
40 //_______________________________________________________________________
42 {
44 }
45 
48 //_______________________________________________________________________
50 {
51  //delete fTf1;
52 }
53 
56 //_______________________________________________________________________
58 AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo> &bunchvector,
59  UInt_t altrocfg1, UInt_t altrocfg2 )
60 {
61  short maxampindex; //index of maximum amplitude
62  short maxamp; //Maximum amplitude
63  int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
64 
65  if( index >= 0)
66  {
67  Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
68  Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
69  short maxrev = maxampindex - bunchvector.at(index).GetStartBin();
70  // timebinOffset is timebin value at maximum (maxrev)
71  short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
72  Float_t time = (timebinOffset*TIMEBINWITH)-fL1Phase;
73 
74  if( maxf < fAmpCut || maxamp > fOverflowCut ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
75  {
76  return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, time, (int)time, 0, 0, Ret::kDummy);
77  }
78  else if ( maxf >= fAmpCut )
79  {
80  int first = 0;
81  int last = 0;
82  SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last, fFitArrayCut );
83  int nsamples = last - first + 1;
84 
85  if( ( nsamples ) >= fNsampleCut )
86  {
87  Float_t tmax = (maxrev - first); // local tmax estimate
88  TGraph *graph = new TGraph( nsamples, fXaxis, &fReversed[first] );
89  fTf1->SetParameter(0, maxf*fkEulerSquared );
90  fTf1->SetParameter(1, tmax - fTau);
91  // set rather loose parameter limits
92  fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
93  fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4);
94 
95  if (fFixTau)
96  {
97  fTf1->FixParameter(2, fTau);
98  }
99  else
100  {
101  fTf1->ReleaseParameter(2); // allow par. to vary
102  fTf1->SetParameter(2, fTau);
103  }
104 
105  Short_t tmpStatus = 0;
106  try
107  {
108  tmpStatus = graph->Fit(fTf1, "Q0RW");
109  }
110  catch (const std::exception & e)
111  {
112  AliError( Form("TGraph Fit exception %s, fit status %d", e.what(),tmpStatus) );
113  return AliCaloFitResults( maxamp, ped, Ret::kNoFit, maxf, time,
114  (int)time, Ret::kDummy, Ret::kDummy, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
115  }
116 
117  if( fVerbose == true )
118  {
119  AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
120  PrintFitResult( fTf1 ) ;
121  }
122 
123  // global tmax
124  tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
125  + fTf1->GetParameter(2); // +tau, makes sum tmax
126  Float_t timemax = (tmax*TIMEBINWITH)-fL1Phase;
127 
128  delete graph;
129 
130  return AliCaloFitResults( maxamp, ped , Ret::kFitPar,
131  fTf1->GetParameter(0)/fkEulerSquared,
132  timemax,
133  (int)timemax,
134  fTf1->GetChisquare(),
135  fTf1->GetNDF(),
136  Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
137 
138  // delete graph;
139 
140  }
141  else
142  {
143  Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
144  Int_t ndf = last - first - 1; // nsamples - 2
145 
146  return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, time,
147  (int)time, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
148  }
149  } // ampcut
150  }
151 
153 }
154 
int fFitArrayCut
Cut on ADC value (after ped. subtraction) for signals used for fit.
void PrintFitResult(const TF1 *f) const
Print fit results.
Float_t fAmpCut
Max ADC - pedestal must be higher than this befor attemting to extract the amplitude.
const double fkEulerSquared
e^2 = 7.389056098930650227
virtual ~AliCaloRawAnalyzerFakeALTRO()
Destructor.
int fOverflowCut
Value when ADC starts to saturate.
int fNsampleCut
Minimum number of sample require before attemting to extract signal parameters.
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.
Float_t fTau
Rise time of the signal (peak position = t0 +tau), by defauly it is 235 ns.
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
static void PrintBunch(const AliCaloBunchInfo &bunch)
Print bunch information.
TF1 * fTf1
Analytical formula of the Semi Gaussian to be fitted.
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
bool fVerbose
Print debug information to std out if set to true.
Double_t fL1Phase
Phase of the ADC sampling clock relative to the LHC clock.
virtual AliCaloFitResults Evaluate(const std::vector< AliCaloBunchInfo > &bunchvector, UInt_t altrocfg1, UInt_t altrocfg2)
Extract signal parameters using fitting.
#define AliError(message)
Definition: AliLog.h:591
Bool_t fFixTau
Fixed fit parameter or not, used in AliCaloRawAnalyzerFitter.