]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliCaloRawAnalyzerFakeALTRO.cxx
correcting baryon mass subtraction for visible energy in MC
[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  
18  
19 Author: R. GUERNANE LPSC Grenoble CNRS/IN2P3
20 */
21
22 #include "AliCaloRawAnalyzerFakeALTRO.h"
23 #include "AliCaloBunchInfo.h"
24 #include "AliCaloFitResults.h"
25 #include "AliLog.h"
26 #include "TMath.h"
27 #include <stdexcept>
28 #include <iostream>
29 #include "TF1.h"
30 #include "TGraph.h"
31 #include "AliCaloConstants.h"
32
33 using namespace std;
34
35
36 #define BAD 4  //CRAP PTH
37  
38 ClassImp( AliCaloRawAnalyzerFakeALTRO )
39
40
41 AliCaloRawAnalyzerFakeALTRO::AliCaloRawAnalyzerFakeALTRO() : AliCaloRawAnalyzer("Chi Square Fit", "LMS"),
42                                                  fkEulerSquared(7.389056098930650227),
43                                                  fTf1(0),
44                                                  fTau(2.35),
45                                                  fFixTau(kTRUE)
46 {
47   //comment
48   for(int i=0; i < MAXSAMPLES; i++)
49     {
50       fXaxis[i] = i;
51     }
52   
53   fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])",  0, 30 ); 
54   if (fFixTau) {
55     fTf1->FixParameter(2, fTau);
56   }
57   else {
58     fTf1->ReleaseParameter(2); // allow par. to vary
59     fTf1->SetParameter(2, fTau);
60   }
61  
62 }
63
64
65 AliCaloRawAnalyzerFakeALTRO::~AliCaloRawAnalyzerFakeALTRO()
66 {
67   delete fTf1;
68 }
69
70
71
72 AliCaloFitResults
73 AliCaloRawAnalyzerFakeALTRO::Evaluate( const vector<AliCaloBunchInfo>  &bunchvector, const UInt_t altrocfg1,  const UInt_t altrocfg2 )
74 {
75   // Extracting signal parameters using fitting
76   short maxampindex; //index of maximum amplitude
77   short maxamp; //Maximum amplitude
78   int index = SelectBunch( bunchvector,  &maxampindex,  &maxamp );
79   
80   if( index >= 0)
81     {
82       Float_t  ped  = ReverseAndSubtractPed( &(bunchvector.at(index))  ,  altrocfg1, altrocfg2, fReversed  );
83       Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(),  fReversed );
84       short maxrev = maxampindex  -  bunchvector.at(index).GetStartBin();
85       // timebinOffset is timebin value at maximum (maxrev)
86       short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
87       if(  maxf < fAmpCut  ||  ( maxamp - ped) > fOverflowCut  ) // (maxamp - ped) > fOverflowCut = Close to saturation (use low gain then)
88         {
89           return  AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset);
90         }            
91       else if ( maxf >= fAmpCut )
92         {
93           int first = 0;
94           int last = 0;
95           SelectSubarray( fReversed,  bunchvector.at(index).GetLength(),  maxrev, &first, &last);
96           int nsamples =  last - first + 1;
97           
98           if( ( nsamples  )  >= fNsampleCut )
99             {
100               Float_t tmax = (maxrev - first); // local tmax estimate
101               TGraph *graph =  new TGraph(  nsamples, fXaxis,  &fReversed[first] );
102               fTf1->SetParameter(0, maxf*fkEulerSquared );
103               fTf1->SetParameter(1, tmax - fTau); 
104               // set rather loose parameter limits
105               fTf1->SetParLimits(0, 0.5*maxf*fkEulerSquared, 2*maxf*fkEulerSquared );
106               fTf1->SetParLimits(1, tmax - fTau - 4, tmax - fTau + 4); 
107
108               if (fFixTau) {
109                 fTf1->FixParameter(2, fTau);
110               }
111               else {
112                 fTf1->ReleaseParameter(2); // allow par. to vary
113                 fTf1->SetParameter(2, fTau);
114               }
115
116               Short_t tmpStatus = 0;
117               try {
118                 tmpStatus =  graph->Fit(fTf1, "Q0RW");
119               }
120               catch (const std::exception & e) {
121                 AliError( Form("TGraph Fit exception %s", e.what()) ); 
122                 return AliCaloFitResults( maxamp, ped, Ret::kNoFit, maxf, timebinOffset,
123                                           timebinOffset, Ret::kDummy, Ret::kDummy,  Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
124               }
125
126               if( fVerbose == true )
127                 {
128                   AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) ); 
129                   PrintFitResult( fTf1 ) ;
130                 }  
131               // global tmax
132               tmax = fTf1->GetParameter(1) + timebinOffset - (maxrev - first) // abs. t0
133                 + fTf1->GetParameter(2); // +tau, makes sum tmax
134               
135                 delete graph;
136                 return AliCaloFitResults( maxamp, ped , Ret::kFitPar,
137                                           fTf1->GetParameter(0)/fkEulerSquared, 
138                                           tmax,
139                                           timebinOffset,  
140                                           fTf1->GetChisquare(), 
141                                           fTf1->GetNDF(),
142                                           Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
143                                 
144                 //     delete graph;
145         
146             }
147           else
148             {
149               Float_t chi2 = CalculateChi2(maxf, maxrev, first, last);
150               Int_t ndf = last - first - 1; // nsamples - 2
151               return AliCaloFitResults( maxamp, ped, Ret::kCrude, maxf, timebinOffset,
152                                         timebinOffset, chi2, ndf, Ret::kDummy, AliCaloFitSubarray(index, maxrev, first, last) ); 
153             }
154         } // ampcut
155     }
156   return AliCaloFitResults(  Ret::kInvalid,  Ret::kInvalid );
157   
158 }
159
160
161 void 
162 AliCaloRawAnalyzerFakeALTRO::PrintFitResult(const TF1 *f) const
163 {
164   //comment
165   cout << endl;
166   cout << __FILE__ << __LINE__ << "Using this samplerange we get" << endl;
167   cout << __FILE__ << __LINE__ << "AMPLITUDE = " << f->GetParameter(0)/fkEulerSquared << ",.. !!!!" << endl;
168   cout << __FILE__ << __LINE__ << "TOF = " << f->GetParameter(1) << ",.. !!!!" << endl;
169   cout << __FILE__ << __LINE__ << "NDF = " << f->GetNDF() << ",.. !!!!" << endl;
170   //  cout << __FILE__ << __LINE__ << "STATUS = " << f->GetStatus()  << ",.. !!!!" << endl << endl;
171   cout << endl << endl;
172 }
173