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"
35 #define BAD 4 //CRAP PTH
37 ClassImp( AliCaloRawAnalyzerFakeALTRO )
40 AliCaloRawAnalyzerFakeALTRO::AliCaloRawAnalyzerFakeALTRO() : AliCaloRawAnalyzer("Chi Square Fit", "LMS"),
41 fkEulerSquared(7.389056098930650227),
47 for(int i=0; i < MAXSAMPLES; i++)
52 fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])", 0, 30 );
54 fTf1->FixParameter(2, fTau);
57 fTf1->ReleaseParameter(2); // allow par. to vary
58 fTf1->SetParameter(2, fTau);
64 AliCaloRawAnalyzerFakeALTRO::~AliCaloRawAnalyzerFakeALTRO()
72 AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
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 );
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)
88 return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kCrude, maxf, timebinOffset);
90 else if ( maxf >= fAmpCut )
94 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last);
95 int nsamples = last - first + 1;
97 if( ( nsamples ) >= fNsampleCut )
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);
108 fTf1->FixParameter(2, fTau);
111 fTf1->ReleaseParameter(2); // allow par. to vary
112 fTf1->SetParameter(2, fTau);
115 Short_t tmpStatus = 0;
117 tmpStatus = graph->Fit(fTf1, "Q0RW");
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) );
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 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) );
155 return AliCaloFitResults( AliCaloFitResults::kInvalid, AliCaloFitResults::kInvalid );
161 AliCaloRawAnalyzerFakeALTRO::PrintFitResult(const TF1 *f) const
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;