source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationChargePix.cc@ 5102

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