1 /**************************************************************************
2 * This file is property of and copyright by *
3 * the Relativistic Heavy Ion Group (RHIG), Yale University, US, 2009 *
5 * Primary Author: Per Thomas Hille <p.t.hille@fys.uio.no> *
7 * Contributors are mentioned in the code where appropriate. *
8 * Please report bugs to p.t.hille@fys.uio.no *
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 **************************************************************************/
20 // Extraction of amplitude and peak position
21 // FRom CALO raw data using
22 // least square fit for the
23 // Moment assuming identical and
24 // independent errors (equivalent with chi square)
27 #include "AliCaloRawAnalyzerLMS.h"
28 #include "AliCaloBunchInfo.h"
29 #include "AliCaloFitResults.h"
41 //#define BAD 4 //CRAP PTH
43 ClassImp( AliCaloRawAnalyzerLMS )
46 AliCaloRawAnalyzerLMS::AliCaloRawAnalyzerLMS() : AliCaloRawAnalyzerFitter("Chi Square Fit", "LMS")
47 // fkEulerSquared(7.389056098930650227),
56 for(int i=0; i < ALTROMAXSAMPLES; i++)
62 // fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
68 fTf1->FixParameter(2, fTau);
72 fTf1->ReleaseParameter(2); // allow par. to vary
73 fTf1->SetParameter(2, fTau);
79 AliCaloRawAnalyzerLMS::~AliCaloRawAnalyzerLMS()
87 AliCaloRawAnalyzerLMS::InitFormula( TF1* f)
89 f = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
94 AliCaloRawAnalyzerLMS::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
96 // Extracting signal parameters using fitting
97 short maxampindex; //index of maximum amplitude
98 short maxamp; //Maximum amplitude
99 int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
103 Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
104 Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
105 short maxrev = maxampindex - bunchvector.at(index).GetStartBin();
106 // timebinOffset is timebin value at maximum (maxrev)
107 short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
108 if( maxf < fAmpCut || ( maxamp - ped) > fOverflowCut ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
110 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset);
112 else if ( maxf >= fAmpCut )
116 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last, fFitArrayCut);
117 int nsamples = last - first + 1;
119 if( ( nsamples ) >= fNsampleCut )
121 Float_t tmax = (maxrev - first); // local tmax estimate
122 TGraph *graph = new TGraph( nsamples, fXaxis, &fReversed[first] );
123 fTf1->SetParameter(0, maxf*fkEulerSquared );
124 fTf1->SetParameter(1, tmax - fTau);
125 // set rather loose parameter limits
126 fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
127 fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4);
130 fTf1->FixParameter(2, fTau);
133 fTf1->ReleaseParameter(2); // allow par. to vary
134 fTf1->SetParameter(2, fTau);
137 Short_t tmpStatus = 0;
139 tmpStatus = graph->Fit(fTf1, "Q0RW");
141 catch (const std::exception & e) {
142 AliError( Form("TGraph Fit exception %s", e.what()) );
143 return AliCaloFitResults( maxamp, ped, Ret::kNoFit, maxf, timebinOffset,
144 timebinOffset, Ret::kDummy, Ret::kDummy, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
147 if( fVerbose == true )
149 AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
150 PrintFitResult( fTf1 ) ;
153 tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
154 + fTf1->GetParameter(2); // +tau, makes sum tmax
157 return AliCaloFitResults( maxamp, ped , Ret::kFitPar,
158 fTf1->GetParameter(0)/fkEulerSquared,
161 fTf1->GetChisquare(),
163 Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
171 Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
172 Int_t ndf = last - first - 1; // nsamples - 2
173 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset,
174 timebinOffset, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
178 return AliCaloFitResults( Ret::kInvalid, Ret::kInvalid );
185 AliCaloRawAnalyzerLMS::PrintFitResult(const TF1 *f) const
189 cout << __FILE__ << __LINE__ << "Using this samplerange we get" << endl;
190 cout << __FILE__ << __LINE__ << "AMPLITUDE = " << f->GetParameter(0)/fkEulerSquared << ",.. !!!!" << endl;
191 cout << __FILE__ << __LINE__ << "TOF = " << f->GetParameter(1) << ",.. !!!!" << endl;
192 cout << __FILE__ << __LINE__ << "NDF = " << f->GetNDF() << ",.. !!!!" << endl;
193 // cout << __FILE__ << __LINE__ << "STATUS = " << f->GetStatus() << ",.. !!!!" << endl << endl;
194 cout << endl << endl;