]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliCaloRawAnalyzerLMS.cxx
skip unused TGraph + add access to TF1
[u/mrichter/AliRoot.git] / EMCAL / AliCaloRawAnalyzerLMS.cxx
1 /**************************************************************************
2  * This file is property of and copyright by                              *
3  * the Relativistic Heavy Ion Group (RHIG), Yale University, US, 2009     *
4  *                                                                        *
5  * Primary Author: Per Thomas Hille <p.t.hille@fys.uio.no>                *
6  *                                                                        *
7  * Contributors are mentioned in the code where appropriate.              *
8  * Please report bugs to p.t.hille@fys.uio.no                             *
9  *                                                                        *
10  * Permission to use, copy, modify and distribute this software and its   *
11  * documentation strictly for non-commercial purposes is hereby granted   *
12  * without fee, provided that the above copyright notice appears in all   *
13  * copies and that both the copyright notice and this permission notice   *
14  * appear in the supporting documentation. The authors make no claims     *
15  * about the suitability of this software for any purpose. It is          *
16  * provided "as is" without express or implied warranty.                  *
17  **************************************************************************/
18
19
20 // Extraction of amplitude and peak position
21 // FRom CALO raw data using
22 // least square fit for the
23 // Moment assuming identical and 
24 // independent errors (equivalent with chi square)
25 // 
26
27 #include "AliCaloRawAnalyzerLMS.h"
28 #include "AliCaloBunchInfo.h"
29 #include "AliCaloFitResults.h"
30 #include "AliLog.h"
31 #include "TMath.h"
32 #include <iostream>
33 #include "TF1.h"
34 #include "TGraph.h"
35
36 using namespace std;
37
38
39 #define BAD 4  //CRAP PTH
40  
41 ClassImp( AliCaloRawAnalyzerLMS )
42
43
44 AliCaloRawAnalyzerLMS::AliCaloRawAnalyzerLMS() : AliCaloRawAnalyzer("Chi Square Fit", "LMS"),
45                                                  fkEulerSquared(7.389056098930650227),
46                                                  fTf1(0),
47                                                  fTau(2.35),
48                                                  fFixTau(kTRUE)
49 {
50   //comment
51   for(int i=0; i < MAXSAMPLES; i++)
52     {
53       fXaxis[i] = i;
54     }
55   
56   fTf1 = new TF1( "myformula", "[0]*((x - [1])/[2])^2*exp(-2*(x -[1])/[2])",  0, 30 ); 
57   if (fFixTau) {
58     fTf1->FixParameter(2, fTau);
59   }
60   else {
61     fTf1->ReleaseParameter(2); // allow par. to vary
62     fTf1->SetParameter(2, fTau);
63   }
64  
65 }
66
67
68 AliCaloRawAnalyzerLMS::~AliCaloRawAnalyzerLMS()
69 {
70   delete fTf1;
71 }
72
73
74
75 AliCaloFitResults
76 AliCaloRawAnalyzerLMS::Evaluate( const vector<AliCaloBunchInfo>  &bunchvector, const UInt_t altrocfg1,  const UInt_t altrocfg2 )
77 {
78   // Extracting signal parameters using fitting
79   short maxampindex; //index of maximum amplitude
80   short maxamp; //Maximum amplitude
81   int index = SelectBunch( bunchvector,  &maxampindex,  &maxamp );
82   
83   if( index >= 0)
84     {
85       Float_t  ped  = ReverseAndSubtractPed( &(bunchvector.at(index))  ,  altrocfg1, altrocfg2, fReversed  );
86       int first;
87       int last;
88       Float_t maxf = TMath::MaxElement( bunchvector.at(index).GetLength(),  fReversed );
89       short maxrev = maxampindex  -  bunchvector.at(index).GetStartBin();
90       short timebinOffset = maxampindex - (bunchvector.at(index).GetLength()-1);
91
92       if ( maxf > fAmpCut )
93         {
94           SelectSubarray( fReversed,  bunchvector.at(index).GetLength(),  maxrev, &first, &last);
95           int nsamples =  last - first + 1;
96           
97           if( ( nsamples  )  >= fNsampleCut )
98             {
99               
100               TGraph *graph =  new TGraph(  nsamples, fXaxis,  &fReversed[first] );
101               fTf1->SetParameter(0, maxf*fkEulerSquared );
102               fTf1->SetParameter(1, 0.2);
103
104               if (fFixTau) {
105                 fTf1->FixParameter(2, fTau);
106               }
107               else {
108                 fTf1->ReleaseParameter(2); // allow par. to vary
109                 fTf1->SetParameter(2, fTau);
110               }
111
112               Short_t tmpStatus =  graph->Fit(fTf1, "Q0RW");
113              
114               if( fVerbose == true )
115                 {
116                   AliCaloRawAnalyzer::PrintBunch( bunchvector.at(index) ); 
117                   PrintFitResult( fTf1 ) ;
118                 }  
119               
120                 delete graph;
121                 return AliCaloFitResults( maxamp, ped ,    tmpStatus,  
122                                           fTf1->GetParameter(0)/fkEulerSquared, 
123                                           fTf1->GetParameter(1) + timebinOffset,  
124                                           fTf1->GetChisquare(), 
125                                           fTf1->GetNDF(),
126                                           AliCaloFitResults::kDummy, AliCaloFitSubarray(index, maxrev, first, last) );
127                                 
128                 //     delete graph;
129         
130             }
131           else
132             {
133               return AliCaloFitResults( maxamp, ped, AliCaloFitResults::kNoFit, maxf, maxrev+timebinOffset, AliCaloFitResults::kNoFit, AliCaloFitResults::kNoFit,
134                                         AliCaloFitResults::kNoFit, AliCaloFitSubarray(index, maxrev, first, last) ); 
135             }
136         }
137       else
138         {
139           return AliCaloFitResults( maxamp , ped, AliCaloFitResults::kNoFit, maxf, maxrev+timebinOffset, AliCaloFitResults::kNoFit, AliCaloFitResults::kNoFit);
140         }       
141     }
142
143   return AliCaloFitResults( AliCaloFitResults::kInvalid, AliCaloFitResults::kInvalid );
144   
145 }
146
147
148 void 
149 AliCaloRawAnalyzerLMS::PrintFitResult(const TF1 *f) const
150 {
151   //comment
152   cout << endl;
153   cout << __FILE__ << __LINE__ << "Using this samplerange we get" << endl;
154   cout << __FILE__ << __LINE__ << "AMPLITUDE = " << f->GetParameter(0)/fkEulerSquared << ",.. !!!!" << endl;
155   cout << __FILE__ << __LINE__ << "TOF = " << f->GetParameter(1) << ",.. !!!!" << endl;
156   cout << __FILE__ << __LINE__ << "NDF = " << f->GetNDF() << ",.. !!!!" << endl;
157   //  cout << __FILE__ << __LINE__ << "STATUS = " << f->GetStatus()  << ",.. !!!!" << endl << endl;
158   cout << endl << endl;
159 }
160