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

Last change on this file since 4994 was 4962, checked in by gaug, 20 years ago
*** empty log message ***
File size: 9.2 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#if 0
84// --------------------------------------------------------------------------
85//
86// ATTENTION: This nasty Clone function is necessary since the ROOT cloning
87// lead to crashes on SOME machines (unfortunately not mine...).
88// This function is a workaround in order to achieve the correct
89// DrawClone() behaviour.
90//
91TObject *MHCalibrationChargePix::Clone(const char *name) const
92{
93
94 MHCalibrationChargePix &pix =
95 *new MHCalibrationChargePix(name ? name : fName.Data(),fTitle.Data());
96 //
97 // Copy MHGausEvents data members
98 //
99 pix.fBinsAfterStripping = fBinsAfterStripping;
100 pix.fCurrentSize = fCurrentSize;
101 pix.fFlags = fFlags;
102 pix.fPowerProbabilityBins = fPowerProbabilityBins;
103
104 if (fHPowerProbability)
105 pix.fHPowerProbability=(TH1I*)fHPowerProbability->Clone();
106
107 if (fPowerSpectrum)
108 pix.fPowerSpectrum = new TArrayF(*fPowerSpectrum);
109
110 pix.fEvents = fEvents;
111
112 if (fFGausFit)
113 pix.fFGausFit=(TF1*)fFGausFit->Clone();
114 if (fFExpFit)
115 pix.fFExpFit=(TF1*)fFExpFit->Clone();
116
117 pix.fFirst = fFirst;
118
119 if (fGraphEvents)
120 pix.fGraphEvents=(TGraph*)fGraphEvents->Clone();
121 if (fGraphPowerSpectrum)
122 pix.fGraphPowerSpectrum=(TGraph*)fGraphPowerSpectrum->Clone();
123
124 pix.fHGausHist = fHGausHist;
125
126 pix.fLast = fLast;
127 pix.fMean = fMean;
128 pix.fMeanErr = fMeanErr;
129 pix.fNbins = fNbins;
130 pix.fNDFLimit = fNDFLimit;
131 pix.fSigma = fSigma;
132 pix.fSigmaErr = fSigmaErr;
133 pix.fProb = fProb;
134 pix.fProbLimit = fProbLimit;
135
136 //
137 // Copy MHCalibrationPix data members
138 //
139 pix.fEventFrequency = fEventFrequency;
140 pix.fBlackoutLimit = fBlackoutLimit;
141 pix.fSaturated = fSaturated;
142 pix.fPickupLimit = fPickupLimit;
143 pix.fPixId = fPixId;
144
145 //
146 // Copy MHCalibrationChargePix data members
147 //
148 pix.fHAbsTime = fHAbsTime;
149
150 return &pix;
151}
152#endif
153
154// --------------------------------------------------------------------------
155//
156// Sets:
157// - fHGausHist.SetBins(fNbins,fFirst,fLast);
158// - fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
159//
160void MHCalibrationChargePix::InitBins()
161{
162 MHGausEvents::InitBins();
163 fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
164}
165
166// --------------------------------------------------------------------------
167//
168// Empty function to overload MHGausEvents::Reset()
169//
170void MHCalibrationChargePix::Reset()
171{
172 MHGausEvents::Reset();
173 fHAbsTime.Reset();
174}
175
176// --------------------------------------------------------------------------
177//
178// Calls MHCalibrationPix::ChangeHistId()
179//
180// Add id to names and titles of:
181// - fHAbsTime
182//
183void MHCalibrationChargePix::ChangeHistId(Int_t id)
184{
185
186 MHCalibrationPix::ChangeHistId(id);
187
188 fHAbsTime.SetName (Form("%s%d", fHAbsTime.GetName(), id));
189 fHAbsTime.SetTitle(Form("%s%d", fHAbsTime.GetTitle(), id));
190
191}
192
193// --------------------------------------------------------------------------
194//
195// returns fHGausHist.Integral("width")
196//
197const Float_t MHCalibrationChargePix::GetIntegral() const
198{
199 return fHGausHist.Integral("width");
200}
201
202// --------------------------------------------------------------------------
203//
204// returns fHAbsTime.GetMean()
205//
206const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
207{
208 return fHAbsTime.GetMean();
209}
210
211// --------------------------------------------------------------------------
212//
213// returns fHAbsTime.GetRMS()
214//
215const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
216{
217 return fHAbsTime.GetRMS();
218}
219
220// --------------------------------------------------------------------------
221//
222// Fills fHAbsTime with t
223// Returns kFALSE, if overflow or underflow occurred, else kTRUE
224//
225Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
226{
227 return fHAbsTime.Fill(t) > -1;
228}
229
230// -----------------------------------------------------------------------------
231//
232// Default draw:
233//
234// The following options can be chosen:
235//
236// "": displays the fHGausHist and the fHAbsTime
237// "all": executes additionally MHCalibrationPix::Draw(), with options
238//
239// The following picture shows a typical outcome of call to Draw("all"):
240// One the left side:
241// - The distribution of the values with the Gauss fit are shown (points connected
242// with each other). The line is green, thus the Gauss fit has a probability higher
243// than 0.5%.
244// - The distribution of the positions of the maximum (abs. arrival times)
245// is displayed. Most of the events have an arrival time of slice 7 (==hardware:8)
246//
247// On the right side:
248// - The first plot shows the distribution of the values with the Gauss fit
249// with error bars
250// - The second plot shows the TGraph with the events vs. time
251// - The third plot shows the fourier transform and a peak at 100 Hz.
252// - The fourth plot shows the projection of the fourier components and an exponential
253// fit, with the result that the observed deviation is not statistical, but signficant with a
254// probability smaller than 0.5%.
255//
256//Begin_Html
257/*
258<img src="images/MHCalibrationChargePixDraw.gif">
259*/
260//End_Html
261//
262void MHCalibrationChargePix::Draw(const Option_t *opt)
263{
264
265 TString option(opt);
266 option.ToLower();
267
268 Int_t win = 1;
269
270 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
271 TVirtualPad *pad = NULL;
272
273 if (option.Contains("all"))
274 {
275 option.ReplaceAll("all","");
276 oldpad->Divide(2,1);
277 win = 2;
278 oldpad->cd(1);
279 TVirtualPad *newpad = gPad;
280 pad = newpad;
281 pad->Divide(1,2);
282 pad->cd(1);
283 }
284 else if (option.Contains("datacheck"))
285 {
286 MHCalibrationPix::Draw("events");
287 return;
288 }
289 else
290 {
291 pad = oldpad;
292 pad->Divide(1,2);
293 pad->cd(1);
294 }
295 /*
296 else
297 {
298 option.ReplaceAll("time","");
299 pad = oldpad;
300 pad->Divide(1,2);
301 pad->cd(1);
302 }
303 */
304 if (!IsEmpty() && !IsOnlyOverflow() && !IsOnlyUnderflow())
305 gPad->SetLogy();
306
307 gPad->SetTicks();
308
309 fHGausHist.GetXaxis()->SetLabelSize(0.06);
310 fHGausHist.GetYaxis()->SetLabelSize(0.07);
311 fHGausHist.GetXaxis()->SetLabelOffset(0.01);
312 fHGausHist.GetYaxis()->SetLabelOffset(0.01);
313 fHGausHist.GetXaxis()->SetTitleSize(0.065);
314 fHGausHist.GetYaxis()->SetTitleSize(0.07);
315 fHGausHist.GetXaxis()->SetTitleOffset(0.6);
316 fHGausHist.GetYaxis()->SetTitleOffset(0.6);
317 fHGausHist.Draw();
318 if (fFGausFit)
319 {
320 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
321 fFGausFit->Draw("same");
322 }
323
324 pad->cd(2);
325 gPad->SetTicks();
326
327 fHAbsTime.GetXaxis()->SetLabelSize(0.06);
328 fHAbsTime.GetYaxis()->SetLabelSize(0.07);
329 fHAbsTime.GetXaxis()->SetLabelOffset(0.01);
330 fHAbsTime.GetYaxis()->SetLabelOffset(0.01);
331 fHAbsTime.GetXaxis()->SetTitleSize(0.065);
332 fHAbsTime.GetYaxis()->SetTitleSize(0.07);
333 fHAbsTime.GetXaxis()->SetTitleOffset(0.6);
334 fHAbsTime.GetYaxis()->SetTitleOffset(0.6);
335 fHAbsTime.Draw();
336
337 if (win < 2)
338 return;
339
340 oldpad->cd(2);
341 MHCalibrationPix::Draw("fourierevents");
342
343}
Note: See TracBrowser for help on using the repository browser.