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 **************************************************************************/
17 Author: R. GUERNANE LPSC Grenoble CNRS/IN2P3
21 #include "AliCaloRawAnalyzerFakeALTRO.h"
22 #include "AliCaloBunchInfo.h"
23 #include "AliCaloFitResults.h"
30 #include "AliCaloConstants.h"
34 ClassImp( AliCaloRawAnalyzerFakeALTRO )
37 AliCaloRawAnalyzerFakeALTRO::AliCaloRawAnalyzerFakeALTRO() : AliCaloRawAnalyzerFitter("Chi Square Fit", "LMS")
39 fAlgo= Algo::kFakeAltro;
43 AliCaloRawAnalyzerFakeALTRO::~AliCaloRawAnalyzerFakeALTRO()
50 AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
52 // Extracting signal parameters using fitting
53 short maxampindex; //index of maximum amplitude
54 short maxamp; //Maximum amplitude
55 int index = SelectBunch( bunchvector, &maxampindex, &maxamp );
59 Float_t ped = ReverseAndSubtractPed( &(bunchvector.at(index)) , altrocfg1, altrocfg2, fReversed );
60 Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(), fReversed );
61 short maxrev = maxampindex - bunchvector.at(index).GetStartBin();
62 // timebinOffset is timebin value at maximum (maxrev)
63 short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
64 if( maxf < fAmpCut || ( maxamp - ped) > fOverflowCut ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
66 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset);
68 else if ( maxf >= fAmpCut )
72 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last, fFitArrayCut );
73 int nsamples = last - first + 1;
75 if( ( nsamples ) >= fNsampleCut )
77 Float_t tmax = (maxrev - first); // local tmax estimate
78 TGraph *graph = new TGraph( nsamples, fXaxis, &fReversed[first] );
79 fTf1->SetParameter(0, maxf*fkEulerSquared );
80 fTf1->SetParameter(1, tmax - fTau);
81 // set rather loose parameter limits
82 fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
83 fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4);
86 fTf1->FixParameter(2, fTau);
89 fTf1->ReleaseParameter(2); // allow par. to vary
90 fTf1->SetParameter(2, fTau);
93 Short_t tmpStatus = 0;
95 tmpStatus = graph->Fit(fTf1, "Q0RW");
97 catch (const std::exception & e) {
98 AliError( Form("TGraph Fit exception %s", e.what()) );
99 return AliCaloFitResults( maxamp, ped, Ret::kNoFit, maxf, timebinOffset,
100 timebinOffset, Ret::kDummy, Ret::kDummy, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
103 if( fVerbose == true )
105 AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
106 PrintFitResult( fTf1 ) ;
109 tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
110 + fTf1->GetParameter(2); // +tau, makes sum tmax
113 return AliCaloFitResults( maxamp, ped , Ret::kFitPar,
114 fTf1->GetParameter(0)/fkEulerSquared,
117 fTf1->GetChisquare(),
119 Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
126 Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
127 Int_t ndf = last - first - 1; // nsamples - 2
128 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset,
129 timebinOffset, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
133 return AliCaloFitResults( Ret::kInvalid, Ret::kInvalid );