1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 Author: R. GUERNANE LPSC Grenoble CNRS/IN2P3
22 #include "AliCaloRawAnalyzerFakeALTRO.h"
23 #include "AliCaloBunchInfo.h"
24 #include "AliCaloFitResults.h"
31 #include "AliCaloConstants.h"
36 //#define BAD 4 //CRAP PTH
38 ClassImp( AliCaloRawAnalyzerFakeALTRO )
41 AliCaloRawAnalyzerFakeALTRO::AliCaloRawAnalyzerFakeALTRO() : AliCaloRawAnalyzerFitter("Chi Square Fit", "LMS")
42 // fkEulerSquared(7.389056098930650227),
50 for(int i=0; i < ALTROMAXSAMPLES; i++)
55 fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
57 fTf1->FixParameter(2, fTau);
60 fTf1->ReleaseParameter(2); // allow par. to vary
61 fTf1->SetParameter(2, fTau);
68 AliCaloRawAnalyzerFakeALTRO::~AliCaloRawAnalyzerFakeALTRO()
76 AliCaloRawAnalyzerFakeALTRO::InitFormula( TF1* f)
78 f = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
83 AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
85 // Extracting signal parameters using fitting
86 short maxampindex; //index of maximum amplitude
87 short maxamp; //Maximum amplitude
88 int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
92 Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
93 Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
94 short maxrev = maxampindex - bunchvector.at(index).GetStartBin();
95 // timebinOffset is timebin value at maximum (maxrev)
96 short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
97 if( maxf < fAmpCut || ( maxamp - ped) > fOverflowCut ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
99 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset);
101 else if ( maxf >= fAmpCut )
105 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last, fFitArrayCut );
106 int nsamples = last - first + 1;
108 if( ( nsamples ) >= fNsampleCut )
110 Float_t tmax = (maxrev - first); // local tmax estimate
111 TGraph *graph = new TGraph( nsamples, fXaxis, &fReversed[first] );
112 fTf1->SetParameter(0, maxf*fkEulerSquared );
113 fTf1->SetParameter(1, tmax - fTau);
114 // set rather loose parameter limits
115 fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
116 fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4);
119 fTf1->FixParameter(2, fTau);
122 fTf1->ReleaseParameter(2); // allow par. to vary
123 fTf1->SetParameter(2, fTau);
126 Short_t tmpStatus = 0;
128 tmpStatus = graph->Fit(fTf1, "Q0RW");
130 catch (const std::exception & e) {
131 AliError( Form("TGraph Fit exception %s", e.what()) );
132 return AliCaloFitResults( maxamp, ped, Ret::kNoFit, maxf, timebinOffset,
133 timebinOffset, Ret::kDummy, Ret::kDummy, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
136 if( fVerbose == true )
138 AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
139 PrintFitResult( fTf1 ) ;
142 tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
143 + fTf1->GetParameter(2); // +tau, makes sum tmax
146 return AliCaloFitResults( maxamp, ped , Ret::kFitPar,
147 fTf1->GetParameter(0)/fkEulerSquared,
150 fTf1->GetChisquare(),
152 Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
159 Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
160 Int_t ndf = last - first - 1; // nsamples - 2
161 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset,
162 timebinOffset, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
166 return AliCaloFitResults( Ret::kInvalid, Ret::kInvalid );
174 AliCaloRawAnalyzerFakeALTRO::PrintFitResult(const TF1 *f) const
178 cout << __FILE__ << __LINE__ << "Using this samplerange we get" << endl;
179 cout << __FILE__ << __LINE__ << "AMPLITUDE = " << f->GetParameter(0)/fkEulerSquared << ",.. !!!!" << endl;
180 cout << __FILE__ << __LINE__ << "TOF = " << f->GetParameter(1) << ",.. !!!!" << endl;
181 cout << __FILE__ << __LINE__ << "NDF = " << f->GetNDF() << ",.. !!!!" << endl;
182 // cout << __FILE__ << __LINE__ << "STATUS = " << f->GetStatus() << ",.. !!!!" << endl << endl;
183 cout << endl << endl;