AliRoot Core  3dc7879 (3dc7879)
AliCaloRawAnalyzerCrude.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * This file is property of and copyright by the Experimental Nuclear *
3  * Physics Group, Dep. of Physics *
4  * University of Oslo, Norway, 2007 *
5  * *
6  * Author: Per Thomas Hille <perthi@fys.uio.no> for the ALICE HLT Project.*
7  * Contributors are mentioned in the code where appropriate. *
8  * Please report bugs to perthi@fys.uio.no *
9  * *
10  * Permission to use, copy, modify and distribute this software and its *
11  * documentation strictly for non-commercial purposes is hereby granted *
12  * without fee, provided that the above copyright notice appears in all *
13  * copies and that both the copyright notice and this permission notice *
14  * appear in the supporting documentation. The authors make no claims *
15  * about the suitability of this software for any purpose. It is *
16  * provided "as is" without express or implied warranty. *
17  **************************************************************************/
18 
20 #include "AliCaloFitResults.h"
21 #include "AliCaloBunchInfo.h"
22 #include "TMath.h"
23 using namespace std;
24 
26 ClassImp(AliCaloRawAnalyzerCrude) ;
28 
31 //_______________________________________________________________________
33 {
35 }
36 
39 //_______________________________________________________________________
41 AliCaloRawAnalyzerCrude::Evaluate(const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2)
42 {
43  short maxampindex; //index of maximum amplitude
44  short maxamp; //Maximum amplitude
45  int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
46 
47  if( index >= 0)
48  {
49  Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
50  Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
51  short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
52  Float_t time = (timebinOffset*TIMEBINWITH)-fL1Phase;
53 
54  if( maxf < fAmpCut || maxamp > fOverflowCut ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
55  //ped removed from the comparison (maybe temporarily)
56  {
57  return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, time, (int)time, 0, 0, Ret::kDummy);
58  }
59  else if ( maxf >= fAmpCut ) // no if statement needed really; keep for readability
60  {
61  int first = 0;
62  int last = 0;
63  int maxrev = maxampindex - bunchvector.at(index).GetStartBin();
64 
65  SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev , &first, &last, fFitArrayCut );
66 
67  Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
68  Int_t ndf = last - first - 1; // nsamples - 2
69 
70  return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, time,
71  (int)time, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
72  } // ampcut
73  } // bunch index
74 
76 
77 } //end Crude
78 
79 
int fFitArrayCut
Cut on ADC value (after ped. subtraction) for signals used for fit.
Raw data fitting: crude fit.
Base class for extraction of signal amplitude and peak position.
Float_t fAmpCut
Max ADC - pedestal must be higher than this befor attemting to extract the amplitude.
virtual AliCaloFitResults Evaluate(const std::vector< AliCaloBunchInfo > &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2)
Extract signal.
int fOverflowCut
Value when ADC starts to saturate.
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
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
Double_t fL1Phase
Phase of the ADC sampling clock relative to the LHC clock.