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

Last change on this file since 3627 was 3625, checked in by gaug, 21 years ago
*** empty log message ***
File size: 10.3 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 and perform a Fourier analysis.
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:
143// - MHGausEvents::Clear()
144//
145void MHCalibrationChargePix::Clear(Option_t *o)
146{
147
148 fPixId = -1;
149 fSaturated = 0.;
150 fPickup = 0.;
151
152 SetExcluded ( kFALSE );
153
154 MHGausEvents::Clear();
155 return;
156}
157
158// --------------------------------------------------------------------------
159//
160// Empty function to overload MHGausEvents::Reset()
161//
162void MHCalibrationChargePix::Reset()
163{
164}
165
166// --------------------------------------------------------------------------
167//
168// Set fPixId to id
169//
170// Add id to names and titles of:
171// - this
172// - fHGausHist
173// - fHAbsTime
174//
175void MHCalibrationChargePix::ChangeHistId(Int_t id)
176{
177
178 fPixId = id;
179
180 fHGausHist.SetName(Form("%s%d", fHGausHist.GetName(), id));
181 fHGausHist.SetTitle(Form("%s%d", fHGausHist.GetTitle(), id));
182
183 fHAbsTime.SetName(Form("%s%d", fHAbsTime.GetName(), id));
184 fHAbsTime.SetTitle(Form("%s%d", fHAbsTime.GetTitle(), id));
185
186 fName = Form("%s%d", fName.Data(), id);
187 fTitle = Form("%s%d", fTitle.Data(), id);
188}
189
190// --------------------------------------------------------------------------
191//
192// Set Excluded bit from outside
193//
194void MHCalibrationChargePix::SetExcluded(const Bool_t b)
195{
196 b ? SETBIT(fFlags,kExcluded) : CLRBIT(fFlags,kExcluded);
197}
198
199const Bool_t MHCalibrationChargePix::IsExcluded() const
200{
201 return TESTBIT(fFlags,kExcluded);
202}
203
204// --------------------------------------------------------------------------
205//
206// returns fHGausHist.Integral("width")
207//
208const Float_t MHCalibrationChargePix::GetIntegral() const
209{
210 return fHGausHist.Integral("width");
211}
212
213// --------------------------------------------------------------------------
214//
215// returns fHAbsTime.GetMean()
216//
217const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
218{
219 return fHAbsTime.GetMean();
220}
221
222// --------------------------------------------------------------------------
223//
224// returns fHAbsTime.GetRMS()
225//
226const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
227{
228 return fHAbsTime.GetRMS();
229}
230
231// --------------------------------------------------------------------------
232//
233// Fills fHAbsTime with t
234// Returns kFALSE, if overflow or underflow occurred, else kTRUE
235//
236Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
237{
238 return fHAbsTime.Fill(t) > -1;
239}
240
241// -----------------------------------------------------------------------------
242//
243// Default draw:
244//
245// The following options can be chosen:
246//
247// "": displays the fHGausHist and the fHAbsTime
248// "all": executes additionally MHGausEvents::Draw(), with options
249//
250void MHCalibrationChargePix::Draw(const Option_t *opt)
251{
252
253 TString option(opt);
254 option.ToLower();
255
256 Int_t win = 1;
257
258 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
259 TVirtualPad *pad = NULL;
260
261 if (option.Contains("all"))
262 {
263 option.ReplaceAll("all","");
264 oldpad->Divide(2,1);
265 win = 2;
266 oldpad->cd(1);
267 TVirtualPad *newpad = gPad;
268 pad = newpad;
269 pad->Divide(1,2);
270 pad->cd(1);
271 }
272 else
273 {
274 pad = oldpad;
275 pad->Divide(1,2);
276 pad->cd(1);
277 }
278
279 if (!IsEmpty())
280 gPad->SetLogy();
281
282 gPad->SetTicks();
283
284 fHGausHist.GetXaxis()->SetLabelSize(0.06);
285 fHGausHist.GetYaxis()->SetLabelSize(0.07);
286 fHGausHist.GetXaxis()->SetLabelOffset(0.01);
287 fHGausHist.GetYaxis()->SetLabelOffset(0.01);
288 fHGausHist.GetXaxis()->SetTitleSize(0.065);
289 fHGausHist.GetYaxis()->SetTitleSize(0.07);
290 fHGausHist.GetXaxis()->SetTitleOffset(0.6);
291 fHGausHist.GetYaxis()->SetTitleOffset(0.6);
292 fHGausHist.Draw(opt);
293 if (fFGausFit)
294 {
295 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
296 fFGausFit->Draw("same");
297 }
298
299 pad->cd(2);
300 gPad->SetTicks();
301
302 fHAbsTime.GetXaxis()->SetLabelSize(0.06);
303 fHAbsTime.GetYaxis()->SetLabelSize(0.07);
304 fHAbsTime.GetXaxis()->SetLabelOffset(0.01);
305 fHAbsTime.GetYaxis()->SetLabelOffset(0.01);
306 fHAbsTime.GetXaxis()->SetTitleSize(0.065);
307 fHAbsTime.GetYaxis()->SetTitleSize(0.07);
308 fHAbsTime.GetXaxis()->SetTitleOffset(0.6);
309 fHAbsTime.GetYaxis()->SetTitleOffset(0.6);
310 fHAbsTime.Draw(opt);
311
312 if (win < 2)
313 return;
314
315 oldpad->cd(2);
316 MHGausEvents::Draw("fourierevents");
317
318}
319
320// -----------------------------------------------------------------------------
321//
322// Repeats the Gauss fit in a smaller range, defined by:
323//
324// min = GetMean() - fPickupLimit * GetSigma();
325// max = GetMean() + fPickupLimit * GetSigma();
326//
327Bool_t MHCalibrationChargePix::RepeatFit(const Option_t *option)
328{
329
330 //
331 // Get new fitting ranges
332 //
333 Axis_t rmin = GetMean() - fPickupLimit * GetSigma();
334 Axis_t rmax = GetMean() + fPickupLimit * GetSigma();
335
336 GetFGausFit()->SetRange(rmin,rmax);
337
338 GetHGausHist()->Fit(GetFGausFit(),option);
339
340 SetMean ( GetFGausFit()->GetParameter(1) );
341 SetMeanErr ( GetFGausFit()->GetParameter(2) );
342 SetSigma ( GetFGausFit()->GetParError(1) );
343 SetSigmaErr ( GetFGausFit()->GetParError(2) );
344 SetProb ( GetFGausFit()->GetProb() );
345
346 //
347 // The fit result is accepted under condition:
348 // 1) The results are not nan's
349 // 2) The NDF is not smaller than fNDFLimit (5)
350 // 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
351 //
352 if ( TMath::IsNaN ( GetMean() )
353 || TMath::IsNaN ( GetMeanErr() )
354 || TMath::IsNaN ( GetProb() )
355 || TMath::IsNaN ( GetSigma() )
356 || TMath::IsNaN ( GetSigmaErr() )
357 || GetFGausFit()->GetNDF() < fNDFLimit
358 || GetProb() < fProbLimit )
359 return kFALSE;
360
361 SetGausFitOK(kTRUE);
362 return kTRUE;
363
364}
365
366
367
368void MHCalibrationChargePix::CountPickup()
369{
370 fPickup = (Int_t)GetHGausHist()->Integral(GetHGausHist()->GetXaxis()->FindBin(GetMean()+fPickupLimit*GetSigma()),
371 GetHGausHist()->GetXaxis()->GetLast(),
372 "width");
373}
374
375
376
377
378
379
Note: See TracBrowser for help on using the repository browser.