]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliCaloRawAnalyzerFakeALTRO.cxx
Move AliEMCALGeometry to libEMCALUtils class with the other geometry classes for...
[u/mrichter/AliRoot.git] / EMCAL / AliCaloRawAnalyzerFakeALTRO.cxx
CommitLineData
de39a0ff 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/*
17
18
19Author: R. GUERNANE LPSC Grenoble CNRS/IN2P3
20*/
21
22#include "AliCaloRawAnalyzerFakeALTRO.h"
23#include "AliCaloBunchInfo.h"
24#include "AliCaloFitResults.h"
25#include "AliLog.h"
26#include "TMath.h"
27#include <stdexcept>
28#include <iostream>
29#include "TF1.h"
30#include "TGraph.h"
31
32using namespace std;
33
34
35#define BAD 4 //CRAP PTH
36
37ClassImp( AliCaloRawAnalyzerFakeALTRO )
38
39
40AliCaloRawAnalyzerFakeALTRO::AliCaloRawAnalyzerFakeALTRO() : AliCaloRawAnalyzer("Chi Square Fit", "LMS"),
41 fkEulerSquared(7.389056098930650227),
42 fTf1(0),
43 fTau(2.35),
44 fFixTau(kTRUE)
45{
46 //comment
47 for(int i=0; i < MAXSAMPLES; i++)
48 {
49 fXaxis[i] = i;
50 }
51
52 fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
53 if (fFixTau) {
54 fTf1->FixParameter(2, fTau);
55 }
56 else {
57 fTf1->ReleaseParameter(2); // allow par. to vary
58 fTf1->SetParameter(2, fTau);
59 }
60
61}
62
63
64AliCaloRawAnalyzerFakeALTRO::~AliCaloRawAnalyzerFakeALTRO()
65{
66 delete fTf1;
67}
68
69
70
71AliCaloFitResults
72AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
73{
74 // Extracting signal parameters using fitting
75 short maxampindex; //index of maximum amplitude
76 short maxamp; //Maximum amplitude
77 int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
78
79 if( index >= 0)
80 {
81 Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
82 Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
83 short maxrev = maxampindex - bunchvector.at(index).GetStartBin();
84 // timebinOffset is timebin value at maximum (maxrev)
85 short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
86 if( maxf < fAmpCut || ( maxamp - ped) > fOverflowCut ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
87 {
88 return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kCrude, maxf, timebinOffset);
89 }
90 else if ( maxf >= fAmpCut )
91 {
92 int first = 0;
93 int last = 0;
94 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last);
95 int nsamples = last - first + 1;
96
97 if( ( nsamples ) >= fNsampleCut )
98 {
99 Float_t tmax = (maxrev - first); // local tmax estimate
100 TGraph *graph = new TGraph( nsamples, fXaxis, &fReversed[first] );
101 fTf1->SetParameter(0, maxf*fkEulerSquared );
102 fTf1->SetParameter(1, tmax - fTau);
103 // set rather loose parameter limits
104 fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
105 fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4);
106
107 if (fFixTau) {
108 fTf1->FixParameter(2, fTau);
109 }
110 else {
111 fTf1->ReleaseParameter(2); // allow par. to vary
112 fTf1->SetParameter(2, fTau);
113 }
114
115 Short_t tmpStatus = 0;
116 try {
117 tmpStatus = graph->Fit(fTf1, "Q0RW");
118 }
119 catch (const std::exception & e) {
120 AliError( Form("TGraph Fit exception %s", e.what()) );
121 return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kNoFit, maxf, timebinOffset,
122 timebinOffset, AliCaloFitResults::kDummy, AliCaloFitResults::kDummy, AliCaloFitResults::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
123 }
124
125 if( fVerbose == true )
126 {
127 AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
128 PrintFitResult( fTf1 ) ;
129 }
130 // global tmax
131 tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
132 + fTf1->GetParameter(2); // +tau, makes sum tmax
133
134 delete graph;
135 return AliCaloFitResults( maxamp, ped , AliCaloFitResults::kFitPar,
136 fTf1->GetParameter(0)/fkEulerSquared,
137 tmax,
138 timebinOffset,
139 fTf1->GetChisquare(),
140 fTf1->GetNDF(),
141 AliCaloFitResults::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
142
143 // delete graph;
144
145 }
146 else
147 {
148 Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
149 Int_t ndf = last - first - 1; // nsamples - 2
150 return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kCrude, maxf, timebinOffset,
151 timebinOffset, chi2, ndf, AliCaloFitResults::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
152 }
153 } // ampcut
154 }
155 return AliCaloFitResults( AliCaloFitResults::kInvalid, AliCaloFitResults::kInvalid );
156
157}
158
159
160void
161AliCaloRawAnalyzerFakeALTRO::PrintFitResult(const TF1 *f) const
162{
163 //comment
164 cout << endl;
165 cout << __FILE__ << __LINE__ << "Using this samplerange we get" << endl;
166 cout << __FILE__ << __LINE__ << "AMPLITUDE = " << f->GetParameter(0)/fkEulerSquared << ",.. !!!!" << endl;
167 cout << __FILE__ << __LINE__ << "TOF = " << f->GetParameter(1) << ",.. !!!!" << endl;
168 cout << __FILE__ << __LINE__ << "NDF = " << f->GetNDF() << ",.. !!!!" << endl;
169 // cout << __FILE__ << __LINE__ << "STATUS = " << f->GetStatus() << ",.. !!!!" << endl << endl;
170 cout << endl << endl;
171}
172