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

Last change on this file since 3614 was 3612, checked in by gaug, 21 years ago
*** empty log message ***
File size: 10.6 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;
58const Float_t MHCalibrationChargePix::fgPickupLimit = 5.;
59// --------------------------------------------------------------------------
60//
61// Default Constructor.
62//
63// Sets:
64// - the default number for fChargeNbins (fgChargeNbins)
65// - the default number for fChargeFirst (fgChargeFirst)
66// - the default number for fChargeLast (fgChargeLast)
67// - the default number for fAbsTimeNbins (fgAbstTimeNbins)
68// - the default number for fAbsTimeFirst (fgAbsTimeFirst)
69// - the default number for fAbsTimeLast (fgAbsTimeLast)
70// - the default number for fPickupLimit (fgPickupLimit)
71//
72// - the default name of the fHGausHist ("HCalibrationCharge")
73// - the default title of the fHGausHist ("Distribution of Summed FADC slices Pixel ")
74// - the default x-axis title for fHGausHist ("Sum FADC Slices")
75// - the default y-axis title for fHGausHist ("Nr. of events")
76//
77// - the default name of the fHAbsTime ("HAbsTimePixel")
78// - the default title of the fHAbsTime ("Distribution of Absolute Arrival Times Pixel ")
79// - the default x-axis title for fHAbsTime ("Absolute Arrival Time [FADC slice nr]")
80// - the default y-axis title for fHAbsTime ("Nr. of events");
81// - the default directory of the fHAbsTime (NULL)
82//
83// Initializes:
84// - fHAbsTime()
85// - fPixId to -1
86// - all variables to 0.
87// - all flags to kFALSE
88//
89MHCalibrationChargePix::MHCalibrationChargePix(const char *name, const char *title)
90 : fHAbsTime(), fFlags(0)
91{
92
93 fName = name ? name : "MHCalibrationChargePix";
94 fTitle = title ? title : "Statistics of the FADC sums of calibration events";
95
96 SetChargeNbins();
97 SetChargeFirst();
98 SetChargeLast();
99
100 SetAbsTimeNbins();
101 SetAbsTimeFirst();
102 SetAbsTimeLast();
103
104 SetPickupLimit();
105
106 fHGausHist.SetName("HCalibrationCharge");
107 fHGausHist.SetTitle("Distribution of Summed FADC slices Pixel");
108 fHGausHist.SetXTitle("Sum FADC Slices");
109 fHGausHist.SetYTitle("Nr. of events");
110
111 fHAbsTime.SetName("HAbsTimePixel");
112 fHAbsTime.SetTitle("Distribution of Absolute Arrival Times Pixel ");
113 fHAbsTime.SetXTitle("Absolute Arrival Time [FADC slice nr]");
114 fHAbsTime.SetYTitle("Nr. of events");
115
116 fHAbsTime.UseCurrentStyle();
117 fHAbsTime.SetDirectory(NULL);
118
119 Clear();
120
121}
122
123// --------------------------------------------------------------------------
124//
125// Sets:
126// - fHGausHist.SetBins(fChargeNbins,fChargeFirst,fChargeLast);
127// - fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
128//
129void MHCalibrationChargePix::Init()
130{
131
132 fHGausHist.SetBins(fChargeNbins,fChargeFirst,fChargeLast);
133 fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
134}
135
136// --------------------------------------------------------------------------
137//
138// Sets:
139// - all variables to 0., except fPixId to -1
140// - all flags to kFALSE
141//
142// - executes MHGausEvents::Clear()
143//
144void MHCalibrationChargePix::Clear(Option_t *o)
145{
146
147 fPixId = -1;
148 fSaturated = 0.;
149 fPickup = 0.;
150
151 SetExcluded ( kFALSE );
152
153 MHGausEvents::Clear();
154 return;
155}
156
157// --------------------------------------------------------------------------
158//
159// Empty function to overload MHGausEvents::Reset()
160//
161void MHCalibrationChargePix::Reset()
162{
163}
164
165// --------------------------------------------------------------------------
166//
167// - Set fPixId to id
168// - Add id to name and title of this
169// - Add id to names of fHGausHist and fHAbsTime
170// - Add id to the titles of fHGausHist and fHAbsTime
171//
172void MHCalibrationChargePix::ChangeHistId(Int_t id)
173{
174
175 fPixId = id;
176
177 fHGausHist.SetName(Form("%s%d", fHGausHist.GetName(), id));
178 fHGausHist.SetTitle(Form("%s%d", fHGausHist.GetTitle(), id));
179
180 fHAbsTime.SetName(Form("%s%d", fHAbsTime.GetName(), id));
181 fHAbsTime.SetTitle(Form("%s%d", fHAbsTime.GetTitle(), id));
182
183 fName = Form("%s%d", fName.Data(), id);
184 fTitle = Form("%s%d", fTitle.Data(), id);
185}
186
187
188void MHCalibrationChargePix::SetExcluded(const Bool_t b)
189{
190 b ? SETBIT(fFlags,kExcluded) : CLRBIT(fFlags,kExcluded);
191}
192
193// --------------------------------------------------------------------------
194//
195// returns fHGausHist.Integral("width")
196//
197const Float_t MHCalibrationChargePix::GetIntegral() const
198{
199 return fHGausHist.Integral("width");
200}
201
202const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
203{
204 return fHAbsTime.GetMean();
205}
206
207const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
208{
209 return fHAbsTime.GetRMS();
210}
211
212const Bool_t MHCalibrationChargePix::IsExcluded() const
213{
214 return TESTBIT(fFlags,kExcluded);
215}
216
217// --------------------------------------------------------------------------
218//
219// Fills fHAbsTime with t
220// Returns kFALSE, if overflow or underflow occurred, else kTRUE
221//
222Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
223{
224 return fHAbsTime.Fill(t) > -1;
225}
226
227// -----------------------------------------------------------------------------
228//
229// Default draw:
230//
231// The following options can be chosen:
232//
233// "": displays the fHGausHist and the fHAbsTime
234// "all": executes additionally MHGausEvents::Draw(), with options
235//
236void MHCalibrationChargePix::Draw(const Option_t *opt)
237{
238
239 TString option(opt);
240 option.ToLower();
241
242 Int_t win = 1;
243
244 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
245 TVirtualPad *pad = NULL;
246
247 if (option.Contains("all"))
248 {
249 option.ReplaceAll("all","");
250 oldpad->Divide(2,1);
251 win = 2;
252 oldpad->cd(1);
253 TVirtualPad *newpad = gPad;
254 pad = newpad;
255 pad->Divide(1,2);
256 pad->cd(1);
257 }
258 else
259 {
260 pad = oldpad;
261 pad->Divide(1,2);
262 pad->cd(1);
263 }
264
265 if (!IsEmpty())
266 gPad->SetLogy();
267
268 gPad->SetTicks();
269
270 fHGausHist.GetXaxis()->SetLabelSize(0.06);
271 fHGausHist.GetYaxis()->SetLabelSize(0.07);
272 fHGausHist.GetXaxis()->SetLabelOffset(0.01);
273 fHGausHist.GetYaxis()->SetLabelOffset(0.01);
274 fHGausHist.GetXaxis()->SetTitleSize(0.065);
275 fHGausHist.GetYaxis()->SetTitleSize(0.07);
276 fHGausHist.GetXaxis()->SetTitleOffset(0.6);
277 fHGausHist.GetYaxis()->SetTitleOffset(0.6);
278 fHGausHist.Draw(opt);
279 if (fFGausFit)
280 {
281 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
282 fFGausFit->Draw("same");
283 }
284
285 pad->cd(2);
286 gPad->SetTicks();
287
288 fHAbsTime.GetXaxis()->SetLabelSize(0.06);
289 fHAbsTime.GetYaxis()->SetLabelSize(0.07);
290 fHAbsTime.GetXaxis()->SetLabelOffset(0.01);
291 fHAbsTime.GetYaxis()->SetLabelOffset(0.01);
292 fHAbsTime.GetXaxis()->SetTitleSize(0.065);
293 fHAbsTime.GetYaxis()->SetTitleSize(0.07);
294 fHAbsTime.GetXaxis()->SetTitleOffset(0.6);
295 fHAbsTime.GetYaxis()->SetTitleOffset(0.6);
296 fHAbsTime.Draw(opt);
297
298 if (win < 2)
299 return;
300
301 oldpad->cd(2);
302 MHGausEvents::Draw("fourierevents");
303
304}
305
306// -----------------------------------------------------------------------------
307//
308// Repeats the Gauss fit in a smaller range, defined by:
309//
310// min = GetMean() - fPickupLimit * GetSigma();
311// max = GetMean() + fPickupLimit * GetSigma();
312//
313Bool_t MHCalibrationChargePix::RepeatFit(const Option_t *option)
314{
315
316 //
317 // Get new fitting ranges
318 //
319 Axis_t rmin = GetMean() - fPickupLimit * GetSigma();
320 Axis_t rmax = GetMean() + fPickupLimit * GetSigma();
321
322 GetFGausFit()->SetRange(rmin,rmax);
323
324 GetHGausHist()->Fit(GetFGausFit(),option);
325
326 SetMean ( GetFGausFit()->GetParameter(1) );
327 SetMeanErr ( GetFGausFit()->GetParameter(2) );
328 SetSigma ( GetFGausFit()->GetParError(1) );
329 SetSigmaErr ( GetFGausFit()->GetParError(2) );
330 SetProb ( GetFGausFit()->GetProb() );
331
332 //
333 // The fit result is accepted under condition:
334 // 1) The results are not nan's
335 // 2) The NDF is not smaller than fNDFLimit (5)
336 // 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
337 //
338 if ( TMath::IsNaN ( GetMean() )
339 || TMath::IsNaN ( GetMeanErr() )
340 || TMath::IsNaN ( GetProb() )
341 || TMath::IsNaN ( GetSigma() )
342 || TMath::IsNaN ( GetSigmaErr() )
343 || GetFGausFit()->GetNDF() < fNDFLimit
344 || GetProb() < fProbLimit )
345 return kFALSE;
346
347 SetGausFitOK(kTRUE);
348 return kTRUE;
349
350}
351
352
353// -----------------------------------------------------------------------------
354//
355// Bypasses the Gauss fit by taking mean and RMS from the histogram
356//
357// Errors are determined in the following way:
358// MeanErr = RMS / Sqrt(entries)
359// SigmaErr = RMS / (2.*Sqrt(entries) )
360//
361void MHCalibrationChargePix::BypassFit()
362{
363
364 //
365 // In case, we do not obtain reasonable values
366 // with the fit, we take the histogram values
367 //
368 SetMean ( fHGausHist.GetMean() );
369 SetMeanErr ( fHGausHist.GetRMS() / TMath::Sqrt(fHGausHist.GetEntries()) );
370 SetSigma ( fHGausHist.GetRMS() );
371 SetSigmaErr ( fHGausHist.GetRMS() / TMath::Sqrt(fHGausHist.GetEntries()) / 2. );
372}
373
374void MHCalibrationChargePix::CountPickup()
375{
376 fPickup = (Int_t)GetHGausHist()->Integral(GetHGausHist()->GetXaxis()->FindBin(GetMean()+fPickupLimit*GetSigma()),
377 GetHGausHist()->GetXaxis()->GetLast(),
378 "width");
379}
380
381
382
383
384
385
Note: See TracBrowser for help on using the repository browser.