]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliCaloRawAnalyzerFakeALTRO.cxx
Add gain fluctuations in simulation (Evi)
[u/mrichter/AliRoot.git] / EMCAL / AliCaloRawAnalyzerFakeALTRO.cxx
CommitLineData
de39a0ff 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/*
652d1723 17 Author: R. GUERNANE LPSC Grenoble CNRS/IN2P3
de39a0ff 18*/
19
652d1723 20
de39a0ff 21#include "AliCaloRawAnalyzerFakeALTRO.h"
22#include "AliCaloBunchInfo.h"
23#include "AliCaloFitResults.h"
24#include "AliLog.h"
25#include "TMath.h"
26#include <stdexcept>
27#include <iostream>
28#include "TF1.h"
29#include "TGraph.h"
168c7b3c 30#include "AliCaloConstants.h"
de39a0ff 31
32using namespace std;
33
de39a0ff 34ClassImp( AliCaloRawAnalyzerFakeALTRO )
35
36
396baaf6 37AliCaloRawAnalyzerFakeALTRO::AliCaloRawAnalyzerFakeALTRO() : AliCaloRawAnalyzerFitter("Chi Square Fit", "LMS")
de39a0ff 38{
41da8832 39 fAlgo= Algo::kFakeAltro;
de39a0ff 40}
41
42
43AliCaloRawAnalyzerFakeALTRO::~AliCaloRawAnalyzerFakeALTRO()
44{
396baaf6 45 //delete fTf1;
de39a0ff 46}
47
48
de39a0ff 49AliCaloFitResults
50AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1, const UInt_t altrocfg2 )
51{
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 );
56
57 if( index >= 0)
58 {
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)
65 {
168c7b3c 66 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset);
de39a0ff 67 }
68 else if ( maxf >= fAmpCut )
69 {
70 int first = 0;
71 int last = 0;
92d9f317 72 SelectSubarray( fReversed, bunchvector.at(index).GetLength(), maxrev, &first, &last, fFitArrayCut );
de39a0ff 73 int nsamples = last - first + 1;
74
75 if( ( nsamples ) >= fNsampleCut )
76 {
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);
84
85 if (fFixTau) {
86 fTf1->FixParameter(2, fTau);
87 }
88 else {
89 fTf1->ReleaseParameter(2); // allow par. to vary
90 fTf1->SetParameter(2, fTau);
91 }
92
93 Short_t tmpStatus = 0;
94 try {
95 tmpStatus = graph->Fit(fTf1, "Q0RW");
96 }
97 catch (const std::exception & e) {
98 AliError( Form("TGraph Fit exception %s", e.what()) );
168c7b3c 99 return AliCaloFitResults( maxamp, ped, Ret::kNoFit, maxf, timebinOffset,
100 timebinOffset, Ret::kDummy, Ret::kDummy, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
de39a0ff 101 }
102
103 if( fVerbose == true )
104 {
105 AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
106 PrintFitResult( fTf1 ) ;
107 }
108 // global tmax
109 tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
110 + fTf1->GetParameter(2); // +tau, makes sum tmax
111
112 delete graph;
168c7b3c 113 return AliCaloFitResults( maxamp, ped , Ret::kFitPar,
de39a0ff 114 fTf1->GetParameter(0)/fkEulerSquared,
115 tmax,
116 timebinOffset,
117 fTf1->GetChisquare(),
118 fTf1->GetNDF(),
168c7b3c 119 Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
de39a0ff 120
121 // delete graph;
122
123 }
124 else
125 {
126 Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
127 Int_t ndf = last - first - 1; // nsamples - 2
168c7b3c 128 return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset,
129 timebinOffset, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
de39a0ff 130 }
131 } // ampcut
132 }
168c7b3c 133 return AliCaloFitResults( Ret::kInvalid, Ret::kInvalid );
de39a0ff 134}
135