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

Last change on this file since 4787 was 4636, checked in by gaug, 20 years ago
*** empty log message ***
File size: 8.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 the charge calibration.
29// Stores and fits the charges and stores the location of the maximum FADC
30// slice. Charges are taken from MExtractedSignalPix.
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// --------------------------------------------------------------------------
59//
60// Default Constructor.
61//
62// Sets:
63// - the default number for fNbins (fgChargeNbins)
64// - the default number for fFirst (fgChargeFirst)
65// - the default number for fLast (fgChargeLast)
66// - the default number for fAbsTimeNbins (fgAbsTimeNbins)
67// - the default number for fAbsTimeFirst (fgAbsTimeFirst)
68// - the default number for fAbsTimeLast (fgAbsTimeLast)
69//
70// - the default name of the fHGausHist ("HCalibrationCharge")
71// - the default title of the fHGausHist ("Distribution of Summed FADC slices Pixel ")
72// - the default x-axis title for fHGausHist ("Sum FADC Slices")
73// - the default y-axis title for fHGausHist ("Nr. of events")
74//
75// - the default name of the fHAbsTime ("HAbsTimePixel")
76// - the default title of the fHAbsTime ("Distribution of Absolute Arrival Times Pixel ")
77// - the default x-axis title for fHAbsTime ("Absolute Arrival Time [FADC slice nr]")
78// - the default y-axis title for fHAbsTime ("Nr. of events");
79// - the default directory of the fHAbsTime (NULL)
80// - the current style for fHAbsTime
81//
82// Initializes:
83// - fHAbsTime()
84//
85// Calls:
86// - Clear();
87//
88MHCalibrationChargePix::MHCalibrationChargePix(const char *name, const char *title)
89 : fHAbsTime()
90{
91
92 fName = name ? name : "MHCalibrationChargePix";
93 fTitle = title ? title : "Statistics of the FADC sums of calibration events";
94
95 SetNbins ( fgChargeNbins );
96 SetFirst ( fgChargeFirst );
97 SetLast ( fgChargeLast );
98
99 SetAbsTimeNbins();
100 SetAbsTimeFirst();
101 SetAbsTimeLast();
102
103 fHGausHist.SetName("HCalibrationCharge");
104 fHGausHist.SetTitle("Distribution of Summed FADC slices Pixel");
105 fHGausHist.SetXTitle("Sum FADC Slices");
106 fHGausHist.SetYTitle("Nr. of events");
107
108 fHAbsTime.SetName("HAbsTimePixel");
109 fHAbsTime.SetTitle("Distribution of Absolute Arrival Times Pixel ");
110 fHAbsTime.SetXTitle("Absolute Arrival Time [FADC slice nr]");
111 fHAbsTime.SetYTitle("Nr. of events");
112
113 fHAbsTime.UseCurrentStyle();
114 fHAbsTime.SetDirectory(NULL);
115
116 Clear();
117
118}
119
120// --------------------------------------------------------------------------
121//
122// Sets:
123// - fHGausHist.SetBins(fNbins,fFirst,fLast);
124// - fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
125//
126void MHCalibrationChargePix::InitBins()
127{
128 MHGausEvents::InitBins();
129 fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
130}
131
132// --------------------------------------------------------------------------
133//
134// Empty function to overload MHGausEvents::Reset()
135//
136void MHCalibrationChargePix::Reset()
137{
138}
139
140// --------------------------------------------------------------------------
141//
142// Calls MHGausEvents::ChangeHistId()
143//
144// Add id to names and titles of:
145// - fHAbsTime
146//
147void MHCalibrationChargePix::ChangeHistId(Int_t id)
148{
149
150 MHGausEvents::ChangeHistId(id);
151
152 fHAbsTime.SetName (Form("%s%d", fHAbsTime.GetName(), id));
153 fHAbsTime.SetTitle(Form("%s%d", fHAbsTime.GetTitle(), id));
154
155}
156
157// --------------------------------------------------------------------------
158//
159// returns fHGausHist.Integral("width")
160//
161const Float_t MHCalibrationChargePix::GetIntegral() const
162{
163 return fHGausHist.Integral("width");
164}
165
166// --------------------------------------------------------------------------
167//
168// returns fHAbsTime.GetMean()
169//
170const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
171{
172 return fHAbsTime.GetMean();
173}
174
175// --------------------------------------------------------------------------
176//
177// returns fHAbsTime.GetRMS()
178//
179const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
180{
181 return fHAbsTime.GetRMS();
182}
183
184// --------------------------------------------------------------------------
185//
186// Fills fHAbsTime with t
187// Returns kFALSE, if overflow or underflow occurred, else kTRUE
188//
189Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
190{
191 return fHAbsTime.Fill(t) > -1;
192}
193
194// -----------------------------------------------------------------------------
195//
196// Default draw:
197//
198// The following options can be chosen:
199//
200// "": displays the fHGausHist and the fHAbsTime
201// "all": executes additionally MHGausEvents::Draw(), with options
202//
203// The following picture shows a typical outcome of call to Draw("all"):
204// One the left side:
205// - The distribution of the values with the Gauss fit are shown (points connected
206// with each other). The line is green, thus the Gauss fit has a probability higher
207// than 0.5%.
208// - The distribution of the positions of the maximum (abs. arrival times)
209// is displayed. Most of the events have an arrival time of slice 7 (==hardware:8)
210//
211// On the right side:
212// - The first plot shows the distribution of the values with the Gauss fit
213// with error bars
214// - The second plot shows the TGraph with the events vs. time
215// - The third plot shows the fourier transform and a peak at 100 Hz.
216// - The fourth plot shows the projection of the fourier components and an exponential
217// fit, with the result that the observed deviation is not statistical, but signficant with a
218// probability smaller than 0.5%.
219//
220//Begin_Html
221/*
222<img src="images/MHCalibrationChargePixDraw.gif">
223*/
224//End_Html
225//
226void MHCalibrationChargePix::Draw(const Option_t *opt)
227{
228
229 TString option(opt);
230 option.ToLower();
231
232 Int_t win = 1;
233
234 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
235 TVirtualPad *pad = NULL;
236
237 if (option.Contains("all"))
238 {
239 option.ReplaceAll("all","");
240 oldpad->Divide(2,1);
241 win = 2;
242 oldpad->cd(1);
243 TVirtualPad *newpad = gPad;
244 pad = newpad;
245 pad->Divide(1,2);
246 pad->cd(1);
247 }
248 else if (option.Contains("datacheck"))
249 {
250 MHGausEvents::Draw("events");
251 return;
252 }
253 else
254 {
255 pad = oldpad;
256 pad->Divide(1,2);
257 pad->cd(1);
258 }
259 /*
260 else
261 {
262 option.ReplaceAll("time","");
263 pad = oldpad;
264 pad->Divide(1,2);
265 pad->cd(1);
266 }
267 */
268 if (!IsEmpty() && !IsOnlyOverflow())
269 gPad->SetLogy();
270
271 gPad->SetTicks();
272
273 fHGausHist.GetXaxis()->SetLabelSize(0.06);
274 fHGausHist.GetYaxis()->SetLabelSize(0.07);
275 fHGausHist.GetXaxis()->SetLabelOffset(0.01);
276 fHGausHist.GetYaxis()->SetLabelOffset(0.01);
277 fHGausHist.GetXaxis()->SetTitleSize(0.065);
278 fHGausHist.GetYaxis()->SetTitleSize(0.07);
279 fHGausHist.GetXaxis()->SetTitleOffset(0.6);
280 fHGausHist.GetYaxis()->SetTitleOffset(0.6);
281 fHGausHist.Draw();
282 if (fFGausFit)
283 {
284 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
285 fFGausFit->Draw("same");
286 }
287
288 pad->cd(2);
289 gPad->SetTicks();
290
291 fHAbsTime.GetXaxis()->SetLabelSize(0.06);
292 fHAbsTime.GetYaxis()->SetLabelSize(0.07);
293 fHAbsTime.GetXaxis()->SetLabelOffset(0.01);
294 fHAbsTime.GetYaxis()->SetLabelOffset(0.01);
295 fHAbsTime.GetXaxis()->SetTitleSize(0.065);
296 fHAbsTime.GetYaxis()->SetTitleSize(0.07);
297 fHAbsTime.GetXaxis()->SetTitleOffset(0.6);
298 fHAbsTime.GetYaxis()->SetTitleOffset(0.6);
299 fHAbsTime.Draw();
300
301 if (win < 2)
302 return;
303
304 oldpad->cd(2);
305 MHGausEvents::Draw("fourierevents");
306
307}
Note: See TracBrowser for help on using the repository browser.