source: trunk/MagicSoft/Mars/mcalib/MHCalibrationChargePix.cc@ 3358

Last change on this file since 3358 was 3358, checked in by gaug, 21 years ago
*** empty log message ***
File size: 6.8 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Markus Gaug 02/2004 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24//////////////////////////////////////////////////////////////////////////////
25//
26// MHCalibrationChargePix
27//
28// Histogram class for charge calibration analysis. Holds the histogrammed
29// summed FADC slices and some rough absolute timing. Calculates the mean
30// sum of FADC slices and its sigma
31//
32//////////////////////////////////////////////////////////////////////////////
33#include "MHCalibrationChargePix.h"
34
35#include <TH1.h>
36#include <TF1.h>
37
38#include <TVirtualPad.h>
39#include <TCanvas.h>
40#include <TPad.h>
41#include <TGraph.h>
42
43#include "MH.h"
44
45#include "MLog.h"
46#include "MLogManip.h"
47
48ClassImp(MHCalibrationChargePix);
49
50using namespace std;
51
52const Int_t MHCalibrationChargePix::fgChargeNbins = 2000;
53const Axis_t MHCalibrationChargePix::fgChargeFirst = -0.5;
54const Axis_t MHCalibrationChargePix::fgChargeLast = 1999.5;
55const Int_t MHCalibrationChargePix::fgAbsTimeNbins = 15;
56const Axis_t MHCalibrationChargePix::fgAbsTimeFirst = -0.5;
57const Axis_t MHCalibrationChargePix::fgAbsTimeLast = 14.5;
58
59const Int_t MHCalibrationChargePix::fgPickupLimit = 5;
60const Int_t MHCalibrationChargePix::fgPulserFrequency = 500;
61// --------------------------------------------------------------------------
62//
63// Default Constructor.
64//
65MHCalibrationChargePix::MHCalibrationChargePix(const char *name, const char *title)
66 : fHAbsTime()
67{
68
69 fName = name ? name : "MHCalibrationChargePix";
70 fTitle = title ? title : "Fill the FADC sums of calibration events and perform the fits";
71
72 SetChargeNbins();
73 SetChargeFirst();
74 SetChargeLast();
75
76 SetAbsTimeNbins();
77 SetAbsTimeFirst();
78 SetAbsTimeLast();
79
80 SetPickupLimit();
81
82 fHGausHist.SetName("HCalibrationCharge");
83 fHGausHist.SetTitle("Distribution of Summed FADC slices Pixel");
84 fHGausHist.SetXTitle("Sum FADC Slices");
85 fHGausHist.SetYTitle("Nr. of events");
86
87 fHAbsTime.SetName("HAbsTimePixel");
88 fHAbsTime.SetTitle("Distribution of Absolute Arrival Times Pixel ");
89 fHAbsTime.SetXTitle("Absolute Arrival Time [FADC slice nr]");
90 fHAbsTime.SetYTitle("Nr. of events");
91
92 fHAbsTime.UseCurrentStyle();
93 fHAbsTime.SetDirectory(NULL);
94
95 Clear();
96
97 SetPulserFrequency();
98}
99
100void MHCalibrationChargePix::Init()
101{
102
103 fHGausHist.SetBins(fChargeNbins,fChargeFirst,fChargeLast);
104 fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
105}
106
107
108void MHCalibrationChargePix::Clear(Option_t *o)
109{
110
111 fPixId = -1;
112 fSaturated = 0;
113 fPickup = 0;
114
115 MHGausEvents::Clear();
116 return;
117}
118
119
120void MHCalibrationChargePix::Reset()
121{
122}
123
124
125void MHCalibrationChargePix::ChangeHistId(Int_t id)
126{
127
128 fPixId = id;
129
130 fHGausHist.SetName(Form("%s%d", fHGausHist.GetName(), id));
131 fHGausHist.SetTitle(Form("%s%d", fHGausHist.GetTitle(), id));
132
133 fHAbsTime.SetName(Form("%s%d", fHAbsTime.GetName(), id));
134 fHAbsTime.SetTitle(Form("%s%d", fHAbsTime.GetTitle(), id));
135
136 fName = Form("%s%d", fName.Data(), id);
137 fTitle = Form("%s%d", fTitle.Data(), id);
138}
139
140
141const Float_t MHCalibrationChargePix::GetIntegral() const
142{
143 return fHGausHist.Integral("width");
144}
145
146const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
147{
148 return fHAbsTime.GetMean();
149}
150
151const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
152{
153 return fHAbsTime.GetRMS();
154}
155
156void MHCalibrationChargePix::SetPulserFrequency(Float_t f)
157{
158 SetEventFrequency(f);
159}
160
161Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
162{
163 return fHAbsTime.Fill(t) > -1;
164}
165
166
167void MHCalibrationChargePix::Draw(const Option_t *opt)
168{
169
170 TString option(opt);
171 option.ToLower();
172
173 Int_t win = 1;
174
175 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
176 TVirtualPad *pad = NULL;
177
178 if (option.Contains("all"))
179 {
180 option.ReplaceAll("all","");
181 oldpad->Divide(2,1);
182 win = 2;
183 oldpad->cd(1);
184 TVirtualPad *newpad = gPad;
185 pad = newpad;
186 pad->Divide(1,2);
187 pad->cd(1);
188 }
189 else
190 {
191 pad = oldpad;
192 pad->Divide(1,2);
193 pad->cd(1);
194 }
195
196 if (!IsEmpty())
197 gPad->SetLogy();
198
199 gPad->SetTicks();
200
201 fHGausHist.Draw(opt);
202 if (fFGausFit)
203 {
204 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
205 fFGausFit->Draw("same");
206 }
207
208 pad->cd(2);
209 gPad->SetTicks();
210 fHAbsTime.Draw(opt);
211
212 if (win < 2)
213 return;
214
215 oldpad->cd(2);
216 MHGausEvents::Draw("fourierevents");
217
218}
219
220Bool_t MHCalibrationChargePix::RepeatFit(const Option_t *option)
221{
222
223 //
224 // Get new fitting ranges
225 //
226 Axis_t rmin = GetMean() - fPickupLimit * GetSigma();
227 Axis_t rmax = GetMean() + fPickupLimit * GetSigma();
228
229 GetFGausFit()->SetRange(rmin,rmax);
230
231 GetHGausHist()->Fit(GetFGausFit(),option);
232
233 SetMean ( GetFGausFit()->GetParameter(1) );
234 SetMeanErr ( GetFGausFit()->GetParameter(2) );
235 SetSigma ( GetFGausFit()->GetParError(1) );
236 SetSigmaErr ( GetFGausFit()->GetParError(2) );
237 SetProb ( GetFGausFit()->GetProb() );
238
239 //
240 // The fit result is accepted under condition:
241 // 1) The results are not nan's
242 // 2) The NDF is not smaller than fNDFLimit (5)
243 // 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
244 //
245 if ( TMath::IsNaN ( GetMean() )
246 || TMath::IsNaN ( GetMeanErr() )
247 || TMath::IsNaN ( GetProb() )
248 || TMath::IsNaN ( GetSigma() )
249 || TMath::IsNaN ( GetSigmaErr() )
250 || GetFGausFit()->GetNDF() < fNDFLimit
251 || GetProb() < fProbLimit )
252 return kFALSE;
253
254 SetGausFitOK(kTRUE);
255 return kTRUE;
256
257}
258
259void MHCalibrationChargePix::BypassFit()
260{
261
262 //
263 // In case, we do not obtain reasonable values
264 // with the fit, we take the histogram values
265 //
266 SetMean ( fHGausHist.GetMean() );
267 SetMeanErr ( fHGausHist.GetRMS() / fHGausHist.GetEntries() );
268 SetSigma ( fHGausHist.GetRMS() );
269 SetSigmaErr ( fHGausHist.GetRMS() / fHGausHist.GetEntries() / 2. );
270}
271
272void MHCalibrationChargePix::CountPickup()
273{
274 fPickup = (Int_t)GetHGausHist()->Integral(GetHGausHist()->GetXaxis()->FindBin(GetMean()+fPickupLimit*GetSigma()),
275 GetHGausHist()->GetXaxis()->GetLast(),
276 "width");
277}
278
279
280
281
282
283
Note: See TracBrowser for help on using the repository browser.