]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliCaloRawAnalyzerFakeALTRO.cxx
5d21ed809d3863f77dc81b59d47a25cdbfda5f25
[u/mrichter/AliRoot.git] / EMCAL / AliCaloRawAnalyzerFakeALTRO.cxx
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 /*
17   Author: R. GUERNANE LPSC Grenoble CNRS/IN2P3
18 */
19
20
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"
30 #include "AliCaloConstants.h"
31
32 using namespace std;
33
34 ClassImp( AliCaloRawAnalyzerFakeALTRO )
35
36
37 AliCaloRawAnalyzerFakeALTRO::AliCaloRawAnalyzerFakeALTRO() : AliCaloRawAnalyzerFitter("Chi Square Fit", "FakeAltro")
38 {
39   // constructor
40   
41   fAlgo= Algo::kFakeAltro;
42 }
43
44 AliCaloRawAnalyzerFakeALTRO::~AliCaloRawAnalyzerFakeALTRO()
45 {
46   // destructor
47   
48   //delete fTf1;
49 }
50
51 AliCaloFitResults
52 AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo>  &bunchvector, const UInt_t altrocfg1,  const UInt_t altrocfg2 )
53 {
54   // Extracting signal parameters using fitting
55   short maxampindex; //index of maximum amplitude
56   short maxamp; //Maximum amplitude
57   int index = SelectBunch( bunchvector,  &maxampindex,  &maxamp );
58   
59   if( index >= 0)
60   {
61     Float_t  ped  = ReverseAndSubtractPed( &(bunchvector.at(index))  ,  altrocfg1, altrocfg2, fReversed  );
62     Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(),  fReversed );
63     short maxrev = maxampindex  -  bunchvector.at(index).GetStartBin();
64     // timebinOffset is timebin value at maximum (maxrev)
65     short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
66     if(  maxf < fAmpCut  ||  ( maxamp - ped) > fOverflowCut  ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
67     {
68       return  AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset);
69     }
70     else if ( maxf >= fAmpCut )
71     {
72       int first = 0;
73       int last = 0;
74       SelectSubarray( fReversed,  bunchvector.at(index).GetLength(),  maxrev, &first, &last, fFitArrayCut );
75       int nsamples =  last - first + 1;
76       
77       if( ( nsamples  )  >= fNsampleCut )
78             {
79               Float_t tmax = (maxrev - first); // local tmax estimate
80               TGraph *graph =  new TGraph(  nsamples, fXaxis,  &fReversed[first] );
81               fTf1->SetParameter(0, maxf*fkEulerSquared );
82               fTf1->SetParameter(1, tmax - fTau);
83               // set rather loose parameter limits
84               fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
85               fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4);
86         
87               if (fFixTau) {
88           fTf1->FixParameter(2, fTau);
89               }
90               else {
91           fTf1->ReleaseParameter(2); // allow par. to vary
92           fTf1->SetParameter(2, fTau);
93               }
94         
95               Short_t tmpStatus = 0;
96               try {
97           tmpStatus =  graph->Fit(fTf1, "Q0RW");
98               }
99               catch (const std::exception & e) {
100           AliError( Form("TGraph Fit exception %s, fit status %d", e.what(),tmpStatus) );
101           return AliCaloFitResults( maxamp, ped, Ret::kNoFit, maxf, timebinOffset,
102                                    timebinOffset, Ret::kDummy, Ret::kDummy,  Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
103               }
104         
105               if( fVerbose == true )
106         {
107           AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) );
108           PrintFitResult( fTf1 ) ;
109         }
110               // global tmax
111               tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
112         + fTf1->GetParameter(2); // +tau, makes sum tmax
113               
114         delete graph;
115         return AliCaloFitResults( maxamp, ped , Ret::kFitPar,
116                                  fTf1->GetParameter(0)/fkEulerSquared,
117                                  tmax,
118                                  timebinOffset,
119                                  fTf1->GetChisquare(),
120                                  fTf1->GetNDF(),
121                                  Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
122                                 
123         //     delete graph;
124         
125             }
126       else
127             {
128               Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
129               Int_t ndf = last - first - 1; // nsamples - 2
130               return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset,
131                                  timebinOffset, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) ); 
132             }
133     } // ampcut
134   }
135   return AliCaloFitResults(  Ret::kInvalid,  Ret::kInvalid );
136 }
137