source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationChargePINDiode.cc@ 4962

Last change on this file since 4962 was 4950, checked in by gaug, 20 years ago
*** empty log message ***
File size: 14.7 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// MHCalibrationChargePINDiode
28//
29// Histogram class for the charge calibration of the PIN Diode.
30// Stores and fits the charges, the RMS of the charges and stores the
31// location of the maximum FADC slice. Charges are taken from MExtractedSignalPINDiode.
32//
33//////////////////////////////////////////////////////////////////////////////
34#include "MHCalibrationChargePINDiode.h"
35
36#include <TH1.h>
37#include <TF1.h>
38#include <TPad.h>
39#include <TVirtualPad.h>
40#include <TCanvas.h>
41
42#include "MLog.h"
43#include "MLogManip.h"
44
45#include "MParList.h"
46
47#include "MExtractedSignalPINDiode.h"
48#include "MCalibrationChargePINDiode.h"
49
50ClassImp(MHCalibrationChargePINDiode);
51
52using namespace std;
53
54const Axis_t MHCalibrationChargePINDiode::fgAbsTimeFirst = -0.5;
55const Axis_t MHCalibrationChargePINDiode::fgAbsTimeLast = 29.5;
56const Int_t MHCalibrationChargePINDiode::fgAbsTimeNbins = 30;
57const Axis_t MHCalibrationChargePINDiode::fgChargeFirst = -0.5;
58const Axis_t MHCalibrationChargePINDiode::fgChargeLast = 1999.5;
59const Int_t MHCalibrationChargePINDiode::fgChargeNbins = 2000;
60const Int_t MHCalibrationChargePINDiode::fgRmsChargeNbins = 200;
61const Axis_t MHCalibrationChargePINDiode::fgRmsChargeFirst = 0.;
62const Axis_t MHCalibrationChargePINDiode::fgRmsChargeLast = 200.;
63const Float_t MHCalibrationChargePINDiode::fgTimeLowerLimit = 3.;
64const Float_t MHCalibrationChargePINDiode::fgTimeUpperLimit = 4.;
65const TString MHCalibrationChargePINDiode::gsHistName = "Charge";
66const TString MHCalibrationChargePINDiode::gsHistTitle = "Signals";
67const TString MHCalibrationChargePINDiode::gsHistXTitle = "Signal [FADC counts]";
68const TString MHCalibrationChargePINDiode::gsHistYTitle = "Nr. events";
69const TString MHCalibrationChargePINDiode::gsAbsHistName = "AbsTime";
70const TString MHCalibrationChargePINDiode::gsAbsHistTitle = "Abs. Arr. Times";
71const TString MHCalibrationChargePINDiode::gsAbsHistXTitle = "Time [FADC slices]";
72const TString MHCalibrationChargePINDiode::gsAbsHistYTitle = "Nr. events";
73// --------------------------------------------------------------------------
74//
75// Default Constructor.
76//
77// Sets:
78// - the default number for fAbsTimeFirst (fgAbsTimeFirst)
79// - the default number for fAbsTimeLast (fgAbsTimeLast)
80// - the default number for fAbsTimeNbins (fgAbsTimeNbins)
81// - the default number for MHGausEvents::fNbins (fgChargeNbins)
82// - the default number for MHGausEvents::fFirst (fgChargeFirst)
83// - the default number for MHGausEvents::fLast (fgChargeLast)
84// - the default number for fRmsChargeNbins (fgRmsChargeNbins)
85// - the default number for fRmsChargeFirst (fgRmsChargeFirst)
86// - the default number for fRmsChargeLast (fgRmsChargeLast)
87// - the default number for fTimeLowerLimit (fgTimeLowerLimit)
88// - the default number for fTimeUpperLimit (fgTimeUpperLimit)
89//
90// - the default name of the fHGausHist ("HCalibrationChargePINDiode")
91// - the default title of the fHGausHist ("Distribution of Summed FADC slices PIN Diode")
92// - the default x-axis title for fHGausHist ("Sum FADC Slices")
93// - the default y-axis title for fHGausHist ("Nr. of events")
94// - the default name of the fHAbsTime ("HAbsTimePINDiode")
95// - the default title of the fHAbsTime ("Distribution of Absolute Arrival Times PIN Diode")
96// - the default x-axis title for fHAbsTime ("Absolute Arrival Time [FADC slice nr]")
97// - the default y-axis title for fHAbsTime ("Nr. of events")
98// - the default name of the fHRmsCharge ("HRmsChargePINDiode")
99// - the default title of the fHRmsCharge ("Distribution of Variances of summed FADC slices PIN Diode")
100// - the default x-axis title for fHRmsCharge ("RMS (sum) [FADC slices]")
101// - the default y-axis title for fHRmsCharge ("Nr. of events")
102// - the default directory of the fHRmsCharge (NULL)
103// - the current style for fHRmsCharge (NULL)
104//
105// - fHistName to gsHistName
106// - fHistTitle to gsHistTitle
107// - fHistXTitle to gsHistXTitle
108// - fHistYTitle to gsHistYTitle
109//
110// - fAbsHistName to gsAbsHistName
111// - fAbsHistTitle to gsAbsHistTitle
112// - fAbsHistXTitle to gsAbsHistXTitle
113// - fAbsHistYTitle to gsAbsHistYTitle
114//
115// Initializes:
116// - fHRmsCharge()
117// - all pointers to NULL
118//
119// Calls:
120// - Clear()
121//
122MHCalibrationChargePINDiode::MHCalibrationChargePINDiode(const char *name, const char *title)
123 : fPINDiode(NULL), fSigPIN(NULL), fHRmsCharge()
124{
125
126 fName = name ? name : "MHCalibrationChargePINDiode";
127 fTitle = title ? title : "Fill the FADC sums of the PINDiode events and perform the fits";
128
129 SetAbsTimeFirst();
130 SetAbsTimeLast();
131 SetAbsTimeNbins();
132
133 SetNbins( fgChargeNbins );
134 SetFirst( fgChargeFirst );
135 SetLast ( fgChargeLast );
136
137 SetRmsChargeNbins();
138 SetRmsChargeFirst();
139 SetRmsChargeLast();
140
141 SetTimeLowerLimit();
142 SetTimeUpperLimit();
143
144 SetHistName (gsHistName .Data());
145 SetHistTitle (gsHistTitle .Data());
146 SetHistXTitle(gsHistXTitle.Data());
147 SetHistYTitle(gsHistYTitle.Data());
148
149 SetAbsHistName (gsAbsHistName .Data());
150 SetAbsHistTitle (gsAbsHistTitle .Data());
151 SetAbsHistXTitle(gsAbsHistXTitle.Data());
152 SetAbsHistYTitle(gsAbsHistYTitle.Data());
153
154 fHRmsCharge.SetName("HRmsChargePINDiode");
155 fHRmsCharge.SetTitle("Distribution of Variances of summed FADC slices PIN Diode");
156 fHRmsCharge.SetXTitle("RMS (sum) [FADC slices]");
157 fHRmsCharge.SetYTitle("Nr. of events");
158 fHRmsCharge.UseCurrentStyle();
159 fHRmsCharge.SetDirectory(NULL);
160
161 Clear();
162}
163
164// --------------------------------------------------------------------------
165//
166// Initializes Binning of the following histograms:
167// - fHGausHist.SetBins(fNbins,fFirst,fLast);
168// - fHAbsTime.SetBins(fAbsTimeNbins,fAbsTimeFirst,fAbsTimeLast);
169// - fHRmsCharge.SetBins(fRmsChargeNbins,fRmsChargeFirst,fRmsChargeLast);
170//
171Bool_t MHCalibrationChargePINDiode::SetupFill(const MParList *pList)
172{
173
174 TH1F *h = GetHGausHist();
175
176 h->SetName (fHistName.Data());
177 h->SetTitle(fHistTitle.Data());
178 h->SetXTitle(fHistXTitle.Data());
179 h->SetYTitle(fHistYTitle.Data());
180
181 h = GetHAbsTime();
182
183 h->SetName (fAbsHistName.Data());
184 h->SetTitle(fAbsHistTitle.Data());
185 h->SetXTitle(fAbsHistXTitle.Data());
186 h->SetYTitle(fAbsHistYTitle.Data());
187
188 MHCalibrationPix::InitBins();
189
190 fHAbsTime. SetBins(fAbsTimeNbins, fAbsTimeFirst, fAbsTimeLast);
191 fHRmsCharge.SetBins(fRmsChargeNbins,fRmsChargeFirst,fRmsChargeLast);
192
193 return kTRUE;
194
195}
196
197// --------------------------------------------------------------------------
198//
199// Use the MHCalibrationChargePix::Clone function and clone additionally
200// the rest of the data members.
201//
202TObject *MHCalibrationChargePINDiode::Clone(const char *name) const
203{
204
205 MHCalibrationChargePINDiode &pix = (MHCalibrationChargePINDiode&)*MHCalibrationChargePix::Clone(name);
206 //
207 // Copy data members
208 //
209 pix.fHRmsCharge = fHRmsCharge;
210
211 pix.fRmsChargeFirst = fRmsChargeFirst;
212 pix.fRmsChargeLast = fRmsChargeLast;
213 pix.fRmsChargeNbins = fRmsChargeNbins;
214 pix.fRmsChargeMean = fRmsChargeMean;
215 pix.fRmsChargeMeanErr = fRmsChargeMeanErr;
216 pix.fRmsChargeSigma = fRmsChargeSigma;
217 pix.fRmsChargeSigmaErr = fRmsChargeSigmaErr;
218 pix.fTimeLowerLimit = fTimeLowerLimit;
219 pix.fTimeUpperLimit = fTimeUpperLimit;
220
221 return &pix;
222}
223
224
225// --------------------------------------------------------------------------
226//
227// Gets or creates the pointers to:
228// - MExtractedSignalPINDiode
229// - MCalibrationChargePINDiode
230//
231Bool_t MHCalibrationChargePINDiode::ReInit(MParList *pList)
232{
233
234 fSigPIN = (MExtractedSignalPINDiode*)pList->FindCreateObj("MExtractedSignalPINDiode");
235 if (!fSigPIN)
236 {
237 *fLog << err << "MExtractedSignalPINDiode not found... aborting " << endl;
238 return kFALSE;
239 }
240
241 fPINDiode = (MCalibrationChargePINDiode*)pList->FindCreateObj("MCalibrationChargePINDiode");
242 if (!fPINDiode)
243 {
244 *fLog << err << "MCalibrationChargePINDiode not found... aborting " << endl;
245 return kFALSE;
246 }
247
248 return kTRUE;
249}
250
251// --------------------------------------------------------------------------
252//
253// Retrieves from MExtractedSignalPINDiode:
254// - Number of used FADC samples via MExtractedSignalPINDiode::GetNumFADCSamples()
255// - Extracted signal via MExtractedSignalPINDiode::GetExtractedSignal()
256// - Signal Rms MExtractedSignalPINDiode::GetExtractedRms()
257// - Arrival Time MExtractedSignalPINDiode::GetExtractedTime()
258//
259// Fills the following histograms:
260// - MHGausEvents::FillHistAndArray(signal)
261// - MHCalibrationChargePix::FillAbsTime(time);
262// - FillRmsCharge(rms);
263//
264Bool_t MHCalibrationChargePINDiode::Fill(const MParContainer *par, const Stat_t w)
265{
266
267 MExtractedSignalPINDiode *extractor = (MExtractedSignalPINDiode*)par;
268
269 if (!extractor)
270 {
271 *fLog << err << "No argument in MExtractedSignalPINDiode::Fill... abort." << endl;
272 return kFALSE;
273 }
274
275 Float_t slices = (Float_t)extractor->GetNumFADCSamples();
276
277 if (slices == 0.)
278 {
279 *fLog << err << "Number of used signal slices in MExtractedSignalPINDiode is zero ... abort."
280 << endl;
281 return kFALSE;
282 }
283
284 const Float_t signal = (float)extractor->GetExtractedSignal();
285 const Float_t time = extractor->GetExtractedTime();
286 const Float_t rms = extractor->GetExtractedRms();
287
288 FillHistAndArray(signal);
289 FillAbsTime(time);
290 FillRmsCharge(rms);
291
292 return kTRUE;
293}
294
295// --------------------------------------------------------------------------
296//
297// Returns kTRUE, if empty
298//
299// Performs the following fits:
300// - MHGausEvents::FitGaus()
301// - FitRmsCharge()
302//
303// Creates the fourier spectrum (MHGausEvents::CreateFourierSpectrum()
304// and sets bit MCalibrationChargePINDiode::SetOscillating( MHGausEvents::IsFourierSpectrumOK() )
305// Retrieves the results of the following fits and stores them in MCalibrationChargePINDiode:
306// - Mean Charge and Error
307// - Sigma Charge and Error
308// - Fit Probability
309// - Abs Time Mean
310// - Abs Time Rms
311// - Rms Charge Mean and Error
312// - Rms Charge Sigma and Error
313//
314// Performs one consistency check on the arrival time:
315// The check returns kFALSE if:
316//
317// -The mean arrival time is in fTimeLowerLimit slices from the lower edge
318// and fUpperLimit slices from the upper edge
319//
320Bool_t MHCalibrationChargePINDiode::Finalize()
321{
322
323 if (IsGausFitOK() || IsEmpty())
324 return kTRUE;
325
326 FitGaus();
327 FitRmsCharge();
328
329 CreateFourierSpectrum();
330 fPINDiode->SetOscillating ( !IsFourierSpectrumOK() );
331
332 fPINDiode->SetMean ( fMean );
333 fPINDiode->SetMeanVar ( fMeanErr * fMeanErr );
334 fPINDiode->SetSigma ( fSigma );
335 fPINDiode->SetSigmaVar ( fSigmaErr * fMeanErr );
336 fPINDiode->SetProb ( fProb );
337
338 fPINDiode->SetAbsTimeMean( GetAbsTimeMean() );
339 fPINDiode->SetAbsTimeRms( GetAbsTimeRms() );
340
341 fPINDiode->SetRmsChargeMean( GetRmsChargeMean() );
342 fPINDiode->SetRmsChargeMeanErr( GetRmsChargeMeanErr() );
343 fPINDiode->SetRmsChargeSigma( GetRmsChargeSigma() );
344 fPINDiode->SetRmsChargeSigmaErr( GetRmsChargeSigmaErr() );
345
346 fPINDiode->SetValid(kTRUE);
347
348 const Byte_t loweredge = fSigPIN->GetFirstUsedSlice();
349 const Byte_t upperedge = fSigPIN->GetLastUsedSlice();
350 const Float_t lowerlimit = (Float_t)loweredge + fTimeLowerLimit;
351 const Float_t upperlimit = (Float_t)upperedge + fTimeUpperLimit;
352
353 if (GetAbsTimeMean() < lowerlimit)
354 {
355 *fLog << warn << GetDescriptor()
356 << Form("%s%3.1f%s%2.1f%s",": Mean ArrivalTime: ",GetAbsTimeMean()," smaller than ",
357 lowerlimit," FADC slices from lower edge in PIN Diode") << endl;
358 *fLog << warn << GetDescriptor() << ": No PIN Diode calibration!! " << endl;
359 fPINDiode->SetValid(kFALSE);
360 }
361
362 if ( GetAbsTimeMean() > upperlimit )
363 {
364 *fLog << warn << GetDescriptor()
365 << Form("%s%3.1f%s%2.1f%s",": Mean ArrivalTime: ",GetAbsTimeMean()," bigger than ",
366 upperlimit," FADC slices from upper edge in PIN Diode") << endl;
367 *fLog << warn << GetDescriptor() << ": No PIN Diode calibration!! " << endl;
368 fPINDiode->SetValid(kFALSE);
369 }
370
371 return kTRUE;
372}
373
374// --------------------------------------------------------------------------
375//
376// Fills fHRmsCharge with q
377// Returns kFALSE, if overflow or underflow occurred, else kTRUE
378//
379Bool_t MHCalibrationChargePINDiode::FillRmsCharge(const Float_t q)
380{
381 return fHRmsCharge.Fill(q) > -1;
382}
383
384// -----------------------------------------------------------
385//
386// Fits -- not yet implemented
387//
388Bool_t MHCalibrationChargePINDiode::FitRmsCharge(Option_t *option)
389{
390 return 1;
391}
392
393
394// -------------------------------------------------------------------------
395//
396// Draw the histogram
397//
398// The following options can be chosen:
399//
400// "": displays the fHGausHist with fits and fHRmsCharge
401// "all": executes additionally MHCalibrationPix::Draw(), with option "fourierevents"
402//
403void MHCalibrationChargePINDiode::Draw(const Option_t *opt)
404{
405
406 TString option(opt);
407 option.ToLower();
408
409 Int_t win = 1;
410
411 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,900, 600);
412 TVirtualPad *pad = NULL;
413
414 oldpad->SetBorderMode(0);
415
416 if (option.Contains("all"))
417 {
418 option.ReplaceAll("all","");
419 oldpad->Divide(2,1);
420 win = 2;
421 oldpad->cd(1);
422 TVirtualPad *newpad = gPad;
423 pad = newpad;
424 pad->Divide(1,2);
425 pad->cd(1);
426 }
427 else
428 {
429 pad = oldpad;
430 pad->Divide(1,2);
431 pad->cd(1);
432 }
433
434 if (IsEmpty())
435 return;
436
437 if (!IsOnlyOverflow() && !IsOnlyUnderflow())
438 gPad->SetLogy();
439
440 gPad->SetTicks();
441
442 fHGausHist.Draw(opt);
443 if (fFGausFit)
444 {
445 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
446 fFGausFit->Draw("same");
447 }
448
449 pad->cd(2);
450 fHRmsCharge.Draw(opt);
451
452 oldpad->cd(2);
453 MHCalibrationPix::Draw("fourierevents");
454}
455
456
457
458
Note: See TracBrowser for help on using the repository browser.