source: trunk/Mars/mhcalib/MHCalibrationChargePix.cc@ 18811

Last change on this file since 18811 was 5137, checked in by gaug, 20 years ago
*** empty log message ***
File size: 6.9 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// returns fHGausHist.Integral("width")
108//
109const Float_t MHCalibrationChargePix::GetIntegral() const
110{
111 return fHGausHist.Integral("width");
112}
113
114// --------------------------------------------------------------------------
115//
116// returns fHAbsTime.GetMean()
117//
118const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
119{
120 return fHAbsTime.GetMean();
121}
122
123// --------------------------------------------------------------------------
124//
125// returns fHAbsTime.GetRMS()
126//
127const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
128{
129 return fHAbsTime.GetRMS();
130}
131
132// --------------------------------------------------------------------------
133//
134// Fills fHAbsTime with t
135// Returns kFALSE, if overflow or underflow occurred, else kTRUE
136//
137Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
138{
139 return fHAbsTime.Fill(t) > -1;
140}
141
142// -----------------------------------------------------------------------------
143//
144// Default draw:
145//
146// The following options can be chosen:
147//
148// "": displays the fHGausHist and the fHAbsTime
149// "all": executes additionally MHCalibrationPix::Draw(), with options
150//
151// The following picture shows a typical outcome of call to Draw("all"):
152// One the left side:
153// - The distribution of the values with the Gauss fit are shown (points connected
154// with each other). The line is green, thus the Gauss fit has a probability higher
155// than 0.5%.
156// - The distribution of the positions of the maximum (abs. arrival times)
157// is displayed. Most of the events have an arrival time of slice 7 (==hardware:8)
158//
159// On the right side:
160// - The first plot shows the distribution of the values with the Gauss fit
161// with error bars
162// - The second plot shows the TGraph with the events vs. time
163// - The third plot shows the fourier transform and a peak at 100 Hz.
164// - The fourth plot shows the projection of the fourier components and an exponential
165// fit, with the result that the observed deviation is not statistical, but signficant with a
166// probability smaller than 0.5%.
167//
168//Begin_Html
169/*
170<img src="images/MHCalibrationChargePixDraw.gif">
171*/
172//End_Html
173//
174void MHCalibrationChargePix::Draw(const Option_t *opt)
175{
176
177 TString option(opt);
178 option.ToLower();
179
180 Int_t win = 1;
181
182 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
183 TVirtualPad *pad = NULL;
184
185 if (option.Contains("all"))
186 {
187 option.ReplaceAll("all","");
188 oldpad->Divide(2,1);
189 win = 2;
190 oldpad->cd(1);
191 TVirtualPad *newpad = gPad;
192 pad = newpad;
193 pad->Divide(1,2);
194 pad->cd(1);
195 }
196 else if (option.Contains("datacheck"))
197 {
198 MHCalibrationPix::Draw("events");
199 return;
200 }
201 else
202 {
203 pad = oldpad;
204 pad->Divide(1,2);
205 pad->cd(1);
206 }
207 /*
208 else
209 {
210 option.ReplaceAll("time","");
211 pad = oldpad;
212 pad->Divide(1,2);
213 pad->cd(1);
214 }
215 */
216 if (!IsEmpty() && !IsOnlyOverflow() && !IsOnlyUnderflow())
217 gPad->SetLogy();
218
219 gPad->SetTicks();
220
221 fHGausHist.GetXaxis()->SetLabelSize(0.06);
222 fHGausHist.GetYaxis()->SetLabelSize(0.07);
223 fHGausHist.GetXaxis()->SetLabelOffset(0.01);
224 fHGausHist.GetYaxis()->SetLabelOffset(0.01);
225 fHGausHist.GetXaxis()->SetTitleSize(0.065);
226 fHGausHist.GetYaxis()->SetTitleSize(0.07);
227 fHGausHist.GetXaxis()->SetTitleOffset(0.6);
228 fHGausHist.GetYaxis()->SetTitleOffset(0.6);
229 fHGausHist.Draw();
230 if (fFGausFit)
231 {
232 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
233 fFGausFit->Draw("same");
234 }
235
236 pad->cd(2);
237 gPad->SetTicks();
238
239 fHAbsTime.GetXaxis()->SetLabelSize(0.06);
240 fHAbsTime.GetYaxis()->SetLabelSize(0.07);
241 fHAbsTime.GetXaxis()->SetLabelOffset(0.01);
242 fHAbsTime.GetYaxis()->SetLabelOffset(0.01);
243 fHAbsTime.GetXaxis()->SetTitleSize(0.065);
244 fHAbsTime.GetYaxis()->SetTitleSize(0.07);
245 fHAbsTime.GetXaxis()->SetTitleOffset(0.6);
246 fHAbsTime.GetYaxis()->SetTitleOffset(0.6);
247 fHAbsTime.Draw();
248
249 if (win < 2)
250 return;
251
252 oldpad->cd(2);
253 MHCalibrationPix::Draw("fourierevents");
254
255}
Note: See TracBrowser for help on using the repository browser.