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

Last change on this file since 3211 was 3174, checked in by gaug, 21 years ago
*** empty log message ***
File size: 5.5 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//
27// MHCalibrationChargePix
28//
29// Histogram class for charge calibration analysis. Holds the histogrammed
30// summed FADC slices and some rough absolute timing. Calculates the mean
31// sum of FADC slices and its sigma
32//
33//////////////////////////////////////////////////////////////////////////////
34#include "MHCalibrationChargePix.h"
35
36#include <TH1.h>
37#include <TF1.h>
38
39#include <TVirtualPad.h>
40#include <TCanvas.h>
41#include <TPad.h>
42#include <TGraph.h>
43
44#include "MH.h"
45
46#include "MLog.h"
47#include "MLogManip.h"
48
49ClassImp(MHCalibrationChargePix);
50
51using namespace std;
52
53const Int_t MHCalibrationChargePix::fgChargeNbins = 2000;
54const Axis_t MHCalibrationChargePix::fgChargeFirst = -0.5;
55const Axis_t MHCalibrationChargePix::fgChargeLast = 1999.5;
56const Int_t MHCalibrationChargePix::fgAbsTimeNbins = 15;
57const Axis_t MHCalibrationChargePix::fgAbsTimeFirst = -0.5;
58const Axis_t MHCalibrationChargePix::fgAbsTimeLast = 14.5;
59
60const Int_t MHCalibrationChargePix::fgPulserFrequency = 200;
61// --------------------------------------------------------------------------
62//
63// Default Constructor.
64//
65MHCalibrationChargePix::MHCalibrationChargePix(const char *name, const char *title)
66 : fHAbsTime()
67{
68
69 fName = name ? name : "MHCalibrationChargePix";
70 fTitle = title ? title : "Fill the FADC sums of calibration events and perform the fits";
71
72 SetChargeNbins();
73 SetChargeFirst();
74 SetChargeLast();
75
76 SetAbsTimeNbins();
77 SetAbsTimeFirst();
78 SetAbsTimeLast();
79
80 fHAbsTime.UseCurrentStyle();
81 fHAbsTime.SetDirectory(NULL);
82
83 Clear();
84
85 SetPulserFrequency();
86}
87
88void MHCalibrationChargePix::Init()
89{
90
91 fHGausHist.SetName("HCalibrationCharge");
92 fHGausHist.SetTitle("Distribution of Summed FADC slices Pixel");
93 fHGausHist.SetXTitle("Sum FADC Slices");
94 fHGausHist.SetYTitle("Nr. of events");
95 fHGausHist.SetBins(fChargeNbins,fChargeFirst,fChargeLast);
96 fHGausHist.Sumw2();
97
98 fHAbsTime.SetName("HAbsTimePixel");
99 fHAbsTime.SetTitle("Distribution of Absolute Arrival Times Pixel ");
100 fHAbsTime.SetXTitle("Absolute Arrival Time [FADC slice nr]");
101 fHAbsTime.SetYTitle("Nr. of events");
102 fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
103
104}
105
106
107
108void MHCalibrationChargePix::Clear(Option_t *o)
109{
110
111 fPixId = -1;
112
113 SetMeanTimeInFirstBin( kFALSE );
114 SetMeanTimeInLastBin( kFALSE );
115
116 MHGausEvents::Clear();
117 return;
118}
119
120
121void MHCalibrationChargePix::Reset()
122{
123}
124
125
126void MHCalibrationChargePix::ChangeHistId(Int_t id)
127{
128
129 fPixId = id;
130
131 fHGausHist.SetName(Form("%s%d", fHGausHist.GetName(), id));
132 fHGausHist.SetTitle(Form("%s%d", fHGausHist.GetTitle(), id));
133
134 fHAbsTime.SetName(Form("%s%d", fHAbsTime.GetName(), id));
135 fHAbsTime.SetTitle(Form("%s%d", fHAbsTime.GetTitle(), id));
136}
137
138void MHCalibrationChargePix::SetMeanTimeInFirstBin(const Bool_t b)
139{
140 b ? SETBIT(fFlags,kMeanTimeInFirstBin) : CLRBIT(fFlags,kMeanTimeInFirstBin);
141}
142
143void MHCalibrationChargePix::SetMeanTimeInLastBin(const Bool_t b)
144{
145 b ? SETBIT(fFlags,kMeanTimeInLastBin) : CLRBIT(fFlags,kMeanTimeInLastBin);
146}
147
148const Bool_t MHCalibrationChargePix::IsMeanTimeInFirstBin() const
149{
150 return TESTBIT(fFlags,kMeanTimeInFirstBin);
151}
152
153const Bool_t MHCalibrationChargePix::IsMeanTimeInLastBin() const
154{
155 return TESTBIT(fFlags,kMeanTimeInLastBin);
156}
157
158
159const Float_t MHCalibrationChargePix::GetIntegral() const
160{
161 return fHGausHist.Integral("width");
162}
163
164const Float_t MHCalibrationChargePix::GetAbsTimeMean() const
165{
166 return fHAbsTime.GetMean();
167}
168
169const Float_t MHCalibrationChargePix::GetAbsTimeRms() const
170{
171 return fHAbsTime.GetRMS();
172}
173
174void MHCalibrationChargePix::SetPulserFrequency(Float_t f)
175{
176 SetEventFrequency(f);
177}
178
179Bool_t MHCalibrationChargePix::FillAbsTime(Float_t t)
180{
181 return fHAbsTime.Fill(t) > -1;
182}
183
184
185void MHCalibrationChargePix::Draw(const Option_t *opt)
186{
187
188 TString option(opt);
189 option.ToLower();
190
191 Int_t win = 1;
192
193 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this,600, 600);
194
195 pad->SetTicks();
196 pad->SetBorderMode(0);
197
198 if (option.Contains("time"))
199 win++;
200
201 pad->SetTicks();
202 pad->SetBorderMode(0);
203 pad->Divide(1,win);
204 pad->cd(1);
205
206 if (!IsEmpty())
207 pad->SetLogy();
208
209 MHGausEvents::Draw(opt);
210
211 if (win > 1)
212 {
213 pad->cd(2);
214 fHAbsTime.Draw(opt);
215 }
216}
217
218void MHCalibrationChargePix::BypassFit()
219{
220
221 //
222 // In case, we do not obtain reasonable values
223 // with the fit, we take the histogram values
224 //
225 SetMean(fHGausHist.GetMean());
226 SetMeanErr(fHGausHist.GetRMS()/fHGausHist.GetEntries());
227 SetSigma(fHGausHist.GetRMS());
228 SetSigmaErr(fHGausHist.GetRMS()/fHGausHist.GetEntries()/2.);
229}
230
Note: See TracBrowser for help on using the repository browser.