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"
40 #define BAD 4 //CRAP PTH
42 ClassImp( AliCaloRawAnalyzerLMS )
45 AliCaloRawAnalyzerLMS::AliCaloRawAnalyzerLMS() : AliCaloRawAnalyzer("Chi Square Fit", "LMS"),
46 fkEulerSquared(7.389056098930650227),
52 for(int i=0; i < MAXSAMPLES; i++)
57 fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
59 fTf1->FixParameter(2, fTau);
62 fTf1->ReleaseParameter(2); // allow par. to vary
63 fTf1->SetParameter(2, fTau);
69 AliCaloRawAnalyzerLMS::~AliCaloRawAnalyzerLMS()
77 AliCaloRawAnalyzerLMS::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
79 // Extracting signal parameters using fitting
80 short maxampindex; //index of maximum amplitude
81 short maxamp; //Maximum amplitude
82 int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
86 Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
89 Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
90 short maxrev = maxampindex - bunchvector.at(index).GetStartBin();
91 // timebinOffset is timebin value at maximum (maxrev)
92 short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
93 if ( maxf >= fAmpCut )
95 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last);
96 int nsamples = last - first + 1;
98 if( ( nsamples ) >= fNsampleCut )
100 Float_t tmax = (maxrev - first); // local tmax estimate
101 TGraph *graph = new TGraph( nsamples, fXaxis, &fReversed[first] );
102 fTf1->SetParameter(0, maxf*fkEulerSquared );
103 fTf1->SetParameter(1, tmax - fTau);
104 // set rather loose parameter limits
105 fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
106 fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4);
109 fTf1->FixParameter(2, fTau);
112 fTf1->ReleaseParameter(2); // allow par. to vary
113 fTf1->SetParameter(2, fTau);
116 Short_t tmpStatus = 0;
118 tmpStatus = graph->Fit(fTf1, "Q0RW");
120 catch (const std::exception & e) {
121 AliError( Form("TGraph Fit exception %s", e.what()) );
122 return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kNoFit, maxf, timebinOffset);
125 if( fVerbose == true )
127 AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
128 PrintFitResult( fTf1 ) ;
131 tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
132 + fTf1->GetParameter(2); // +tau, makes sum tmax
135 return AliCaloFitResults( maxamp, ped , AliCaloFitResults::kFitPar,
136 fTf1->GetParameter(0)/fkEulerSquared,
139 fTf1->GetChisquare(),
141 AliCaloFitResults::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
148 return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kCrude, maxf, timebinOffset);
153 return AliCaloFitResults( maxamp , ped, AliCaloFitResults::kCrude, maxf, timebinOffset);
157 return AliCaloFitResults( AliCaloFitResults::kInvalid, AliCaloFitResults::kInvalid );
163 AliCaloRawAnalyzerLMS::PrintFitResult(const TF1 *f) const
167 cout << __FILE__ << __LINE__ << "Using this samplerange we get" << endl;
168 cout << __FILE__ << __LINE__ << "AMPLITUDE = " << f->GetParameter(0)/fkEulerSquared << ",.. !!!!" << endl;
169 cout << __FILE__ << __LINE__ << "TOF = " << f->GetParameter(1) << ",.. !!!!" << endl;
170 cout << __FILE__ << __LINE__ << "NDF = " << f->GetNDF() << ",.. !!!!" << endl;
171 // cout << __FILE__ << __LINE__ << "STATUS = " << f->GetStatus() << ",.. !!!!" << endl << endl;
172 cout << endl << endl;