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

Last change on this file since 3956 was 3765, 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// Sets:
135// - fSaturated to 0.
136//
137// Executes:
138// - MHGausEvents::Clear()
139//
140void MHCalibrationChargePix::Clear(Option_t *o)
141{
142
143 fSaturated = 0.;
144
145 MHGausEvents::Clear();
146 return;
147}
148
149// --------------------------------------------------------------------------
150//
151// Empty function to overload MHGausEvents::Reset()
152//
153void MHCalibrationChargePix::Reset()
154{
155}
156
157// --------------------------------------------------------------------------
158//
159// Calls MHGausEvents::ChangeHistId()
160//
161// Add id to names and titles of:
162// - fHAbsTime
163//
164void MHCalibrationChargePix::ChangeHistId(Int_t id)
165{
166
167 MHGausEvents::ChangeHistId(id);
168
169 fHAbsTime.SetName (Form("%s%d", fHAbsTime.GetName(), id));
170 fHAbsTime.SetTitle(Form("%s%d", fHAbsTime.GetTitle(), id));
171
172}
173
174// --------------------------------------------------------------------------
175//
176// returns fHGausHist.Integral("width")
177//
178const Float_t MHCalibrationChargePix::GetIntegral() const
179{
180 return fHGausHist.Integral("width");
181}
182
183// --------------------------------------------------------------------------
184//
185// returns fHAbsTime.GetMean()
186//
187const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
188{
189 return fHAbsTime.GetMean();
190}
191
192// --------------------------------------------------------------------------
193//
194// returns fHAbsTime.GetRMS()
195//
196const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
197{
198 return fHAbsTime.GetRMS();
199}
200
201// --------------------------------------------------------------------------
202//
203// Fills fHAbsTime with t
204// Returns kFALSE, if overflow or underflow occurred, else kTRUE
205//
206Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
207{
208 return fHAbsTime.Fill(t) > -1;
209}
210
211// -----------------------------------------------------------------------------
212//
213// Default draw:
214//
215// The following options can be chosen:
216//
217// "": displays the fHGausHist and the fHAbsTime
218// "all": executes additionally MHGausEvents::Draw(), with options
219//
220// The following picture shows a typical outcome of call to Draw("all"):
221// One the left side:
222// - The distribution of the values with the Gauss fit are shown (points connected
223// with each other). The line is green, thus the Gauss fit has a probability higher
224// than 0.5%.
225// - The distribution of the positions of the maximum (abs. arrival times)
226// is displayed. Most of the events have an arrival time of slice 7 (==hardware:8)
227//
228// On the right side:
229// - The first plot shows the distribution of the values with the Gauss fit
230// with error bars
231// - The second plot shows the TGraph with the events vs. time
232// - The third plot shows the fourier transform and a peak at 100 Hz.
233// - The fourth plot shows the projection of the fourier components and an exponential
234// fit, with the result that the observed deviation is not statistical, but signficant with a
235// probability smaller than 0.5%.
236//
237//Begin_Html
238/*
239<img src="images/MHCalibrationChargePixDraw.gif">
240*/
241//End_Html
242//
243void MHCalibrationChargePix::Draw(const Option_t *opt)
244{
245
246 TString option(opt);
247 option.ToLower();
248
249 Int_t win = 1;
250
251 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
252 TVirtualPad *pad = NULL;
253
254 if (option.Contains("all"))
255 {
256 option.ReplaceAll("all","");
257 oldpad->Divide(2,1);
258 win = 2;
259 oldpad->cd(1);
260 TVirtualPad *newpad = gPad;
261 pad = newpad;
262 pad->Divide(1,2);
263 pad->cd(1);
264 }
265 else
266 {
267 pad = oldpad;
268 pad->Divide(1,2);
269 pad->cd(1);
270 }
271
272 if (!IsEmpty())
273 gPad->SetLogy();
274
275 gPad->SetTicks();
276
277 fHGausHist.GetXaxis()->SetLabelSize(0.06);
278 fHGausHist.GetYaxis()->SetLabelSize(0.07);
279 fHGausHist.GetXaxis()->SetLabelOffset(0.01);
280 fHGausHist.GetYaxis()->SetLabelOffset(0.01);
281 fHGausHist.GetXaxis()->SetTitleSize(0.065);
282 fHGausHist.GetYaxis()->SetTitleSize(0.07);
283 fHGausHist.GetXaxis()->SetTitleOffset(0.6);
284 fHGausHist.GetYaxis()->SetTitleOffset(0.6);
285 fHGausHist.Draw(opt);
286 if (fFGausFit)
287 {
288 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
289 fFGausFit->Draw("same");
290 }
291
292 pad->cd(2);
293 gPad->SetTicks();
294
295 fHAbsTime.GetXaxis()->SetLabelSize(0.06);
296 fHAbsTime.GetYaxis()->SetLabelSize(0.07);
297 fHAbsTime.GetXaxis()->SetLabelOffset(0.01);
298 fHAbsTime.GetYaxis()->SetLabelOffset(0.01);
299 fHAbsTime.GetXaxis()->SetTitleSize(0.065);
300 fHAbsTime.GetYaxis()->SetTitleSize(0.07);
301 fHAbsTime.GetXaxis()->SetTitleOffset(0.6);
302 fHAbsTime.GetYaxis()->SetTitleOffset(0.6);
303 fHAbsTime.Draw(opt);
304
305 if (win < 2)
306 return;
307
308 oldpad->cd(2);
309 MHGausEvents::Draw("fourierevents");
310
311}
Note: See TracBrowser for help on using the repository browser.