]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliCaloRawAnalyzerComparison.cxx
merge
[u/mrichter/AliRoot.git] / EMCAL / AliCaloRawAnalyzerComparison.cxx
1 /**************************************************************************
2  * This file is property of and copyright by the Experimental Nuclear     *
3  * Physics Group, Dep. of Physics                                         *
4  * University of Oslo, Norway, 2007                                       *
5  *                                                                        *
6  * Author: Per Thomas Hille <perthi@fys.uio.no> for the ALICE HLT Project.*
7  * Contributors are mentioned in the code where appropriate.              *
8  * Please report bugs to perthi@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 #include "AliCaloRawAnalyzerComparison.h"
19
20 #include "AliCaloRawAnalyzerCrude.h"
21 #include "AliCaloRawAnalyzerLMS.h"
22 #include "AliCaloRawAnalyzerFastFit.h"
23 #include "AliCaloRawAnalyzerNN.h"
24 #include "AliCaloRawAnalyzerPeakFinder.h"
25
26 #include "AliCaloFitResults.h"
27 #include "TH2D.h"
28 #include "TMath.h"
29 #include "TFile.h"
30
31 #include <iostream>
32
33 using namespace std;
34
35 AliCaloRawAnalyzerComparison::AliCaloRawAnalyzerComparison() : fMod(0),
36                                                                fMonCol1(1),
37                                                                fMonRow1(15),
38                                                                fMonCol2(1),
39                                                                fMonRow2(16)
40 { // ctor
41   fReferenceAnalyzer = new  AliCaloRawAnalyzerLMS();
42
43   AliCaloRawAnalyzerCrude *crude = new AliCaloRawAnalyzerCrude();
44   AliCaloRawAnalyzerLMS *lms = new AliCaloRawAnalyzerLMS();
45   AliCaloRawAnalyzerFastFit *fastfit = new AliCaloRawAnalyzerFastFit();
46   AliCaloRawAnalyzerNN *neuralnet = new AliCaloRawAnalyzerNN();
47   AliCaloRawAnalyzerPeakFinder *peakfinder = new AliCaloRawAnalyzerPeakFinder();
48
49   fReferenceAnalyzer->SetIsZeroSuppressed();
50   peakfinder->SetIsZeroSuppressed();
51   lms->SetIsZeroSuppressed();
52   fastfit->SetIsZeroSuppressed();
53   neuralnet->SetIsZeroSuppressed();
54   crude->SetIsZeroSuppressed(); 
55   
56   fRawAnalyzers.push_back( (AliCaloRawAnalyzer*)crude);
57   fRawAnalyzers.push_back( (AliCaloRawAnalyzer*)lms);
58   fRawAnalyzers.push_back( (AliCaloRawAnalyzer*)fastfit);
59   fRawAnalyzers.push_back( (AliCaloRawAnalyzer*)neuralnet);
60   
61   fRawAnalyzers.push_back( (AliCaloRawAnalyzer*)peakfinder);
62   InitHistograms( fRawAnalyzers, fReferenceAnalyzer );
63 }
64
65
66 AliCaloRawAnalyzerComparison::~AliCaloRawAnalyzerComparison()
67 { // dtor
68
69 }
70  
71
72 void 
73 AliCaloRawAnalyzerComparison::Evaluate( const vector<AliCaloBunchInfo> &bunchvector, const UInt_t altrocfg1,  
74                                         const UInt_t altrocfg2, const int event, const int col, const int row )
75 { // evaluate methods
76   //  cout << __FILE__ << __LINE__ << endl;
77
78   AliCaloFitResults ref = fReferenceAnalyzer->Evaluate( bunchvector, altrocfg1, altrocfg2 );
79   //  AliCaloFitResults ref;
80
81   for(int i=0; i < fRawAnalyzers.size(); i++ )
82     {
83       AliCaloFitResults an = fRawAnalyzers.at(i)->Evaluate( bunchvector, altrocfg1,  altrocfg2 ); 
84       
85       //     cout << __FILE__ << __LINE__ << ":" << fRawAnalyzers.at(i)->GetAlgoAbbr() << "  col = " << col << " row " << row << "  amp= " <<  an.GetAmp() << endl ;
86
87       //   fAmplitudeVsEvent[i]->Fill( event, an.GetAmp() );
88       //   fTofVsEvent[i]->Fill( event, an.GetTof() );
89     
90       fRefAmpVsAnalyzers[i]->Fill(ref.GetAmp(), an.GetAmp() );
91       fRefTofVsAnalyzers[i]->Fill(ref.GetTof(), an.GetTof() ); 
92       fAmpDiff[i]->Fill(ref.GetAmp() - an.GetAmp() ); 
93       fTofDiff[i]->Fill(ref.GetTof() - an.GetTof() );  
94      
95       if(col  == fMonCol1 && fMonRow1 == row )
96         {
97           fAmplitudeVsEvent[i]->Fill( event, an.GetAmp() );
98           fTofVsEvent[i]->Fill( event, an.GetTof() );
99           fMon1[i] = an;
100         }
101       
102       if(col  == fMonCol2 &&  row == fMonRow2  )
103         {
104           fMon2[i] = an;
105         }
106  
107       if(row >= 0 && col >= 0 )
108         {
109           fAmpHistograms[i][col][row]->Fill( an.GetAmp() );
110         }
111     }
112 }  
113
114
115
116 void 
117 AliCaloRawAnalyzerComparison::EventChanged()
118 { // new event
119   for(int i=0; i < fRawAnalyzers.size(); i++ )
120     {
121       /*
122       if( ( fMon1[i].GetAmp() > 50 &&   fMon2[i].GetAmp() > 50 ) &&  (  fMon1[i].GetAmp() < 1023  &&   fMon2[i].GetAmp() < 1023 )  )
123         {
124           cout << __FILE__ << __LINE__ << fRawAnalyzers.at(i)->GetAlgoAbbr() << "\tamp1=" <<  fMon1[i].GetAmp()  << "\tamp2 = " <<  fMon2[i].GetAmp() << endl;
125         }
126       */
127
128       if( fMon1[i].GetAmp()  > 50  &&  fMon2[i].GetAmp() > 50  )
129         {
130           //  if(  fMon1[i].GetTof() != fMon2[i].GetTof())
131           {
132             //    fTofResDifferential[i]->Fill( ( fMon1[i].GetTof() - fMon2[i].GetTof())/TMath::Sqrt(2) );
133             fTofResDifferential[i]->Fill( ( fMon1[i].GetTof() - fMon2[i].GetTof()) );
134             fTofResAbsolute[i]->Fill(  fMon1[i].GetTof() ); 
135           }
136         }
137     }
138 }
139
140
141 void 
142 AliCaloRawAnalyzerComparison::WriteHistograms()
143 { // write histograms 
144   TFile *f = new TFile("comparison2.root", "recreate");
145   
146   /*
147   for(int col=0; col < NZCOLSSMOD; col ++ )
148     {
149       for(int row=0; row < NXROWSSMOD; row ++ )
150         {
151           fAmpHistograms[col][row]->Write();
152         }
153     }
154   */  
155
156   for(int i=0; i < fRawAnalyzers.size(); i++ )
157     {
158       for(int col=0; col < NZCOLSSMOD; col ++ )
159         {
160           for(int row=0; row < NXROWSSMOD; row ++ )
161             {
162               fAmpHistograms[i][col][row]->Write();
163             }
164         }
165       
166       fAmplitudeVsEvent[i]->Write();
167       fTofVsEvent[i]->Write();
168       fRefAmpVsAnalyzers[i]->Write();
169       fRefTofVsAnalyzers[i]->Write(); 
170       fAmpDiff[i]->Write(); 
171       fTofDiff[i]->Write();  
172       fTofResDifferential[i]->Write();
173       fTofResAbsolute[i]->Write(); 
174
175    }
176  
177   f->Close();
178 }
179
180
181 void
182 AliCaloRawAnalyzerComparison::InitHistograms( vector <AliCaloRawAnalyzer*> analyzers, AliCaloRawAnalyzer* ref )
183 { // init histograms 
184   char tmpname[256];
185  
186   /* 
187   for(int col=0; col < NZCOLSSMOD; col ++ )
188     {
189       for(int row=0; row < NXROWSSMOD; row ++ )
190         {
191           sprintf(tmpname, "z(col)%d_x(row)%d_amplitude;Counts;Amplitude/ADC counts", col, row);
192           fAmpHistograms[col][row] = new TH1D(tmpname, tmpname, 1024, 0, 1023 );
193         }
194     }
195   */
196
197   // TH1D *fAmpHistograms[NZCOLSSMOD][NXROWSSMOD];
198
199   for(int i=0; i < analyzers.size(); i++ )
200     {
201       for(int col=0; col < NZCOLSSMOD; col ++ )
202         {
203           for(int row=0; row < NXROWSSMOD; row ++ )
204             {
205               sprintf(tmpname, "z(col)%d_x(row)%d_amplitude_%s;Amplitude/ADC counts;Counts", col, row,  analyzers.at(i)->GetAlgoAbbr() );
206               fAmpHistograms[i][col][row] = new TH1D(tmpname, tmpname, 1024, 0, 1023 );
207             }
208         }
209       
210       sprintf(tmpname, "%s_amplitude_vs_event_row%d_col%d;Event;Amplitude/ADC counts", analyzers.at(i)->GetAlgoAbbr(), fMonRow1, fMonCol1 );
211       fAmplitudeVsEvent[i]=new TH2D(tmpname, tmpname, 8000, 0, 7999, 1024, 0, 1023 ); 
212      
213       sprintf(tmpname, "%s_tof_vs_event_row%d_col%d;Event;Amplitude/ADC counts", analyzers.at(i)->GetAlgoAbbr(), fMonRow1, fMonCol1 );  
214       
215       //     fTofVsEvent[i]=new TH2D(tmpname, tmpname, 8000, 0, 7999, 2000, 3000, 5000); 
216       fTofVsEvent[i]=new TH2D(tmpname, tmpname, 8000, 0, 7999, 2000, 1000, 4000);  
217
218
219       sprintf(tmpname, "%s_vs_%s_ampltude; Amplitude_{%s}/ADC counts; Amplitude_{%s}/ADC counts", ref->GetAlgoAbbr(), 
220               analyzers.at(i)->GetAlgoAbbr(), ref->GetAlgoAbbr(), analyzers.at(i)->GetAlgoAbbr() );
221       fRefAmpVsAnalyzers[i] = new TH2D(tmpname, tmpname, 1024, 0, 1023, 1024, 0, 1023); 
222       sprintf(tmpname, "%s_vs_%s_tof; tof_{%s}/ns; tof_{%s}/ns", ref->GetAlgoAbbr(), analyzers.at(i)->GetAlgoAbbr(), ref->GetAlgoAbbr(), analyzers.at(i)->GetAlgoAbbr() );
223       fRefTofVsAnalyzers[i] = new TH2D(tmpname, tmpname, 500, 2000, 4999, 500, 2000, 4999 ); 
224       sprintf( tmpname, "%s-%s amplitude;counts;A_{%s} - A_{%s}",   ref->GetAlgoAbbr(), 
225                analyzers.at(i)->GetAlgoAbbr(), ref->GetAlgoAbbr(), analyzers.at(i)->GetAlgoAbbr());
226       fAmpDiff[i] = new TH1D(tmpname, tmpname, 100, -10, 10 );
227       sprintf( tmpname, "%s-%s tof;counts;A_{%s} - A_{%s}",   ref->GetAlgoAbbr(), 
228                analyzers.at(i)->GetAlgoAbbr(), ref->GetAlgoAbbr(), analyzers.at(i)->GetAlgoAbbr());
229       fTofDiff[i] = new TH1D(tmpname, tmpname, 1000, -5000, 5000 );
230       sprintf( tmpname, "%s Differential tof resolution (%d, %d) vs (%d, %d);#sigma_{tof}^{%s}/ns;Counts",  analyzers.at(i)->GetAlgoAbbr(), fMonCol1, fMonRow1, fMonCol2, fMonRow2, 
231                analyzers.at(i)->GetAlgoAbbr() );
232       fTofResDifferential[i] = new TH1D(tmpname, tmpname, 1000, -250, 250 );
233       sprintf( tmpname, "%s Absolute tof distribution (%d, %d);#sigma_{tof}^{%s}/ns;Counts",  analyzers.at(i)->GetAlgoAbbr(), fMonCol1, fMonRow1, analyzers.at(i)->GetAlgoAbbr() );
234       
235       fTofResAbsolute[i] = new TH1D(tmpname, tmpname, 2000 , -1000, 7000 );
236
237
238     }
239 }
240