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"
39 #define BAD 4 //CRAP PTH
41 ClassImp( AliCaloRawAnalyzerLMS )
44 AliCaloRawAnalyzerLMS::AliCaloRawAnalyzerLMS() : AliCaloRawAnalyzer("Chi Square Fit", "LMS"),
45 fkEulerSquared(7.389056098930650227),
51 for(int i=0; i < MAXSAMPLES; i++)
56 fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
58 fTf1->FixParameter(2, fTau);
61 fTf1->ReleaseParameter(2); // allow par. to vary
62 fTf1->SetParameter(2, fTau);
68 AliCaloRawAnalyzerLMS::~AliCaloRawAnalyzerLMS()
76 AliCaloRawAnalyzerLMS::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
78 // Extracting signal parameters using fitting
79 short maxampindex; //index of maximum amplitude
80 short maxamp; //Maximum amplitude
81 int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
85 Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
88 Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
89 short maxrev = maxampindex - bunchvector.at(index).GetStartBin();
90 short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
94 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last);
95 int nsamples = last - first + 1;
97 if( ( nsamples ) >= fNsampleCut )
100 TGraph *graph = new TGraph( nsamples, fXaxis, &fReversed[first] );
101 fTf1->SetParameter(0, maxf*fkEulerSquared );
102 fTf1->SetParameter(1, 0.2);
105 fTf1->FixParameter(2, fTau);
108 fTf1->ReleaseParameter(2); // allow par. to vary
109 fTf1->SetParameter(2, fTau);
112 Short_t tmpStatus = graph->Fit(fTf1, "Q0RW");
114 if( fVerbose == true )
116 AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
117 PrintFitResult( fTf1 ) ;
121 return AliCaloFitResults( maxamp, ped , tmpStatus,
122 fTf1->GetParameter(0)/fkEulerSquared,
123 fTf1->GetParameter(1) + timebinOffset,
124 fTf1->GetChisquare(),
126 AliCaloFitResults::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
133 return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kNoFit, maxf, maxrev+timebinOffset, AliCaloFitResults::kNoFit, AliCaloFitResults::kNoFit,
134 AliCaloFitResults::kNoFit, AliCaloFitSubarray(index, maxrev, first, last) );
139 return AliCaloFitResults( maxamp , ped, AliCaloFitResults::kNoFit, maxf, maxrev+timebinOffset, AliCaloFitResults::kNoFit, AliCaloFitResults::kNoFit);
143 return AliCaloFitResults( AliCaloFitResults::kInvalid, AliCaloFitResults::kInvalid );
149 AliCaloRawAnalyzerLMS::PrintFitResult(const TF1 *f) const
153 cout << __FILE__ << __LINE__ << "Using this samplerange we get" << endl;
154 cout << __FILE__ << __LINE__ << "AMPLITUDE = " << f->GetParameter(0)/fkEulerSquared << ",.. !!!!" << endl;
155 cout << __FILE__ << __LINE__ << "TOF = " << f->GetParameter(1) << ",.. !!!!" << endl;
156 cout << __FILE__ << __LINE__ << "NDF = " << f->GetNDF() << ",.. !!!!" << endl;
157 // cout << __FILE__ << __LINE__ << "STATUS = " << f->GetStatus() << ",.. !!!!" << endl << endl;
158 cout << endl << endl;