source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationRelTimeCam.cc@ 5134

Last change on this file since 5134 was 5098, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 21.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// MHCalibrationRelTimeCam
27//
28// Fills the extracted relative arrival times of MArrivalTimeCam into
29// the MHCalibrationPix-classes MHCalibrationPix for every:
30//
31// - Pixel, stored in the TObjArray's MHCalibrationCam::fHiGainArray
32// or MHCalibrationCam::fHiGainArray, respectively, depending if
33// MArrivalTimePix::IsLoGainUsed() is set.
34//
35// - Average pixel per AREA index (e.g. inner and outer for the MAGIC camera),
36// stored in the TObjArray's MHCalibrationCam::fAverageHiGainAreas and
37// MHCalibrationCam::fAverageHiGainAreas
38//
39// - Average pixel per camera SECTOR (e.g. sectors 1-6 for the MAGIC camera),
40// stored in the TObjArray's MHCalibrationCam::fAverageHiGainSectors
41// and MHCalibrationCam::fAverageHiGainSectors
42//
43// Every relative time is calculated as the difference between the individual
44// pixel arrival time and the one of pixel 1 (hardware number: 2).
45// The relative times are filled into a histogram and an array, in order to perform
46// a Fourier analysis (see MHGausEvents). The signals are moreover averaged on an
47// event-by-event basis and written into the corresponding average pixels.
48//
49// The histograms are fitted to a Gaussian, mean and sigma with its errors
50// and the fit probability are extracted. If none of these values are NaN's and
51// if the probability is bigger than MHGausEvents::fProbLimit (default: 0.5%),
52// the fit is declared valid.
53// Otherwise, the fit is repeated within ranges of the previous mean
54// +- MHCalibrationPix::fPickupLimit (default: 5) sigma (see MHCalibrationPix::RepeatFit())
55// In case this does not make the fit valid, the histogram means and RMS's are
56// taken directly (see MHCalibrationPix::BypassFit()) and the following flags are set:
57// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kRelTimeNotFitted ) and
58// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
59//
60// Outliers of more than MHCalibrationPix::fPickupLimit (default: 5) sigmas
61// from the mean are counted as Pickup events (stored in MHCalibrationPix::fPickup)
62//
63// The class also fills arrays with the signal vs. event number, creates a fourier
64// spectrum (see MHGausEvents::CreateFourierSpectrum()) and investigates if the
65// projected fourier components follow an exponential distribution.
66// In case that the probability of the exponential fit is less than
67// MHGausEvents::fProbLimit (default: 0.5%), the following flags are set:
68// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kRelTimeOscillating ) and
69// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
70//
71// This same procedure is performed for the average pixels.
72//
73// The following results are written into MCalibrationRelTimeCam:
74//
75// - MCalibrationPix::SetMean()
76// - MCalibrationPix::SetMeanErr()
77// - MCalibrationPix::SetSigma()
78// - MCalibrationPix::SetSigmaErr()
79// - MCalibrationPix::SetProb()
80// - MCalibrationPix::SetNumPickup()
81//
82// For all averaged areas, the fitted sigma is multiplied with the square root of
83// the number involved pixels in order to be able to compare it to the average of
84// sigmas in the camera.
85//
86/////////////////////////////////////////////////////////////////////////////
87#include "MHCalibrationRelTimeCam.h"
88#include "MHCalibrationPix.h"
89
90#include "MLog.h"
91#include "MLogManip.h"
92
93#include "MParList.h"
94
95#include "MCalibrationIntensityRelTimeCam.h"
96
97#include "MCalibrationRelTimeCam.h"
98#include "MCalibrationRelTimePix.h"
99#include "MCalibrationPix.h"
100
101#include "MArrivalTimeCam.h"
102#include "MArrivalTimePix.h"
103
104#include "MGeomCam.h"
105#include "MGeomPix.h"
106
107#include "MBadPixelsIntensityCam.h"
108#include "MBadPixelsCam.h"
109#include "MBadPixelsPix.h"
110
111#include <TPad.h>
112#include <TVirtualPad.h>
113#include <TCanvas.h>
114#include <TStyle.h>
115#include <TF1.h>
116#include <TLine.h>
117#include <TLatex.h>
118#include <TLegend.h>
119#include <TGraph.h>
120
121ClassImp(MHCalibrationRelTimeCam);
122
123using namespace std;
124
125const Float_t MHCalibrationRelTimeCam::fgNumHiGainSaturationLimit = 0.25;
126const UInt_t MHCalibrationRelTimeCam::fgReferencePixel = 1;
127const Int_t MHCalibrationRelTimeCam::fgNbins = 900;
128const Axis_t MHCalibrationRelTimeCam::fgFirst = -5.;
129const Axis_t MHCalibrationRelTimeCam::fgLast = 5.;
130const TString MHCalibrationRelTimeCam::gsHistName = "RelTime";
131const TString MHCalibrationRelTimeCam::gsHistTitle = "Rel. Arr. Times";
132const TString MHCalibrationRelTimeCam::gsHistXTitle = "Rel. Arr. Time [FADC slices]";
133const TString MHCalibrationRelTimeCam::gsHistYTitle = "Nr. events";
134// --------------------------------------------------------------------------
135//
136// Default Constructor.
137//
138// Sets:
139// - fReferencePixel to fgReferencePixel
140// - fNbins to fgNbins
141// - fFirst to fgFirst
142// - fLast to fgLast
143//
144// - fHistName to gsHistName
145// - fHistTitle to gsHistTitle
146// - fHistXTitle to gsHistXTitle
147// - fHistYTitle to gsHistYTitle
148//
149MHCalibrationRelTimeCam::MHCalibrationRelTimeCam(const char *name, const char *title)
150{
151
152 fName = name ? name : "MHCalibrationRelTimeCam";
153 fTitle = title ? title : "Histogram class for the relative time calibration of the camera";
154
155 SetNumHiGainSaturationLimit(fgNumHiGainSaturationLimit);
156
157 SetReferencePixel();
158
159 SetNbins(fgNbins);
160 SetFirst(fgFirst);
161 SetLast (fgLast );
162
163 SetHistName (gsHistName .Data());
164 SetHistTitle (gsHistTitle .Data());
165 SetHistXTitle(gsHistXTitle.Data());
166 SetHistYTitle(gsHistYTitle.Data());
167
168}
169
170// --------------------------------------------------------------------------
171//
172// Gets or creates the pointers to:
173// - MCalibrationRelTimeCam
174//
175// Searches pointer to:
176// - MArrivalTimeCam
177//
178// Calls:
179// - MHCalibrationCam::InitHiGainArrays()
180// - MHCalibrationCam::InitLoGainArrays()
181//
182// Sets:
183// - fSumareahi to nareas
184// - fSumarealo to nareas
185// - fSumsectorhi to nareas
186// - fSumsectorlo to nareas
187// - fNumareahi to nareas
188// - fNumarealo to nareas
189// - fNumsectorhi to nareas
190// - fNumsectorlo to nareas
191//
192Bool_t MHCalibrationRelTimeCam::ReInitHists(MParList *pList)
193{
194
195 fIntensCam = (MCalibrationIntensityCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityRelTimeCam"));
196 if (fIntensCam)
197 *fLog << inf << "Found MCalibrationIntensityRelTimeCam ... " << endl;
198 else
199 {
200 fCam = (MCalibrationCam*)pList->FindObject(AddSerialNumber("MCalibrationRelTimeCam"));
201 if (!fCam)
202 {
203 fCam = (MCalibrationCam*)pList->FindCreateObj(AddSerialNumber("MCalibrationRelTimeCam"));
204 if (!fCam)
205 {
206 *fLog << err << "Cannot find nor create MCalibrationRelTimeCam ... abort." << endl;
207 return kFALSE;
208 }
209 fCam->Init(*fGeom);
210 }
211 }
212
213 MArrivalTimeCam *signal = (MArrivalTimeCam*)pList->FindObject("MArrivalTimeCam");
214 if (!signal)
215 {
216 *fLog << err << "MArrivalTimeCam not found... abort." << endl;
217 return kFALSE;
218 }
219
220 const Int_t npixels = fGeom->GetNumPixels();
221 const Int_t nsectors = fGeom->GetNumSectors();
222 const Int_t nareas = fGeom->GetNumAreas();
223
224 InitHiGainArrays(npixels,nareas,nsectors);
225 InitLoGainArrays(npixels,nareas,nsectors);
226
227 fSumareahi .Set(nareas);
228 fSumarealo .Set(nareas);
229 fSumsectorhi.Set(nsectors);
230 fSumsectorlo.Set(nsectors);
231 fNumareahi .Set(nareas);
232 fNumarealo .Set(nareas);
233 fNumsectorhi.Set(nsectors);
234 fNumsectorlo.Set(nsectors);
235
236 return kTRUE;
237}
238
239
240// -------------------------------------------------------------------------------
241//
242// Retrieves pointer to MArrivalTimeCam:
243//
244// Retrieves from MGeomCam:
245// - number of pixels
246// - number of pixel areas
247// - number of sectors
248//
249// Fills HiGain or LoGain histograms (MHGausEvents::FillHistAndArray()), respectively
250// depending on MArrivalTimePix::IsLoGainUsed(), with:
251// - MArrivalTimePix::GetArrivalTime(pixid) - MArrivalTimePix::GetArrivalTime(1);
252// (i.e. the time difference between pixel i and pixel 1 (hardware number: 2) )
253//
254Bool_t MHCalibrationRelTimeCam::FillHists(const MParContainer *par, const Stat_t w)
255{
256
257 MArrivalTimeCam *arrtime = (MArrivalTimeCam*)par;
258 if (!arrtime)
259 {
260 gLog << err << "No argument in MArrivalTime::Fill... abort." << endl;
261 return kFALSE;
262 }
263
264 const Int_t npixels = fGeom->GetNumPixels();
265 const Int_t nareas = fGeom->GetNumAreas();
266 const Int_t nsectors = fGeom->GetNumSectors();
267
268 fSumareahi .Reset();
269 fSumarealo .Reset();
270 fSumsectorhi.Reset();
271 fSumsectorlo.Reset();
272 fNumareahi .Reset();
273 fNumarealo .Reset();
274 fNumsectorhi.Reset();
275 fNumsectorlo.Reset();
276
277 const MArrivalTimePix &refpix = (*arrtime)[fReferencePixel];
278 const Float_t reftime = refpix.IsLoGainUsed()
279 ? refpix.GetArrivalTimeLoGain() : refpix.GetArrivalTimeHiGain();
280
281 for (Int_t i=0; i<npixels; i++)
282 {
283
284 MHCalibrationPix &histhi = (*this)[i];
285 MHCalibrationPix &histlo = (*this)(i);
286
287 if (histhi.IsExcluded())
288 continue;
289
290 const MArrivalTimePix &pix = (*arrtime)[i];
291 const Int_t aidx = (*fGeom)[i].GetAidx();
292 const Int_t sector = (*fGeom)[i].GetSector();
293
294 if (pix.IsLoGainUsed() && IsLoGain())
295 {
296 const Float_t reltime = pix.GetArrivalTimeLoGain() - reftime;
297 histhi.AddSaturated(1);
298 histlo.FillHistAndArray(reltime);
299 fSumarealo [aidx] += reltime;
300 fNumarealo [aidx] ++;
301 fSumsectorlo[sector] += reltime;
302 fNumsectorlo[sector] ++;
303 }
304 else
305 {
306 const Float_t reltime = pix.GetArrivalTimeHiGain() - reftime;
307
308 histhi.FillHistAndArray(reltime) ;
309 fSumareahi [aidx] += reltime;
310 fNumareahi [aidx] ++;
311 fSumsectorhi[sector] += reltime;
312 fNumsectorhi[sector] ++;
313 }
314 }
315
316 for (Int_t j=0; j<nareas; j++)
317 {
318 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
319 histhi.FillHistAndArray(fNumareahi[j] == 0 ? 0. : fSumareahi[j]/fNumareahi[j]);
320
321 if (IsLoGain())
322 {
323 MHCalibrationPix &histlo = GetAverageLoGainArea(j);
324 histlo.FillHistAndArray(fNumarealo[j] == 0 ? 0. : fSumarealo[j]/fNumarealo[j]);
325 }
326 }
327
328 for (Int_t j=0; j<nsectors; j++)
329 {
330 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
331 histhi.FillHistAndArray(fNumsectorhi[j] == 0 ? 0. : fSumsectorhi[j]/fNumsectorhi[j]);
332
333 if (IsLoGain())
334 {
335 MHCalibrationPix &histlo = GetAverageLoGainSector(j);
336 histlo.FillHistAndArray(fNumsectorlo[j] == 0 ? 0. : fSumsectorlo[j]/fNumsectorlo[j]);
337 }
338 }
339
340 return kTRUE;
341}
342
343// --------------------------------------------------------------------------
344//
345// Calls:
346// - MHCalibrationCam::FitHiGainArrays() with flags:
347// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
348// - MHCalibrationCam::FitLoGainArrays() with flags:
349// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
350//
351Bool_t MHCalibrationRelTimeCam::FinalizeHists()
352{
353 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
354 {
355
356 MHCalibrationPix &histhi = (*this)[i];
357
358 if (histhi.IsExcluded())
359 continue;
360
361 MCalibrationRelTimePix &pix = fIntensCam
362 ? (MCalibrationRelTimePix&)(*fIntensCam)[i]
363 : (MCalibrationRelTimePix&)(*fCam)[i];
364
365 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
366 {
367 pix.SetHiGainSaturation();
368 histhi.SetExcluded();
369 }
370 else
371 if (IsLoGain())
372 (*this)(i).SetExcluded();
373
374 Stat_t overflow = histhi.GetHGausHist()->GetBinContent(histhi.GetHGausHist()->GetNbinsX()+1);
375 if (overflow > 0.1)
376 {
377 *fLog << warn << GetDescriptor()
378 << ": HiGain Histogram Overflow occurred " << overflow
379 << " times in pixel: " << i << " (without saturation!) " << endl;
380 // bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
381 }
382
383 overflow = histhi.GetHGausHist()->GetBinContent(0);
384 if (overflow > 0.1)
385 {
386 *fLog << warn << GetDescriptor()
387 << ": HiGain Histogram Underflow occurred " << overflow
388 << " times in pixel: " << i << " (without saturation!) " << endl;
389 // bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
390 }
391 }
392
393 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
394 {
395
396 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
397
398 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
399 {
400 MCalibrationRelTimePix &pix = fIntensCam
401 ? (MCalibrationRelTimePix&)fIntensCam->GetAverageArea(j)
402 : (MCalibrationRelTimePix&)fCam->GetAverageArea(j);
403 pix.SetHiGainSaturation();
404 histhi.SetExcluded();
405 }
406 else
407 if (IsLoGain())
408 GetAverageLoGainArea(j).SetExcluded();
409
410 }
411
412 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
413 {
414
415 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
416 MCalibrationRelTimePix &pix = fIntensCam
417 ? (MCalibrationRelTimePix&)fIntensCam->GetAverageSector(j)
418 : (MCalibrationRelTimePix&)fCam->GetAverageSector(j);
419
420 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
421 {
422 pix.SetHiGainSaturation();
423 histhi.SetExcluded();
424 }
425 else
426 if (IsLoGain())
427 GetAverageLoGainSector(j).SetExcluded();
428 }
429
430 FitHiGainArrays(fIntensCam ? (MCalibrationCam&)(*fIntensCam->GetCam()) : (MCalibrationCam&)(*fCam),
431 fIntensBad ? (*fIntensBad->GetCam()) : *fBadPixels,
432 MBadPixelsPix::kRelTimeNotFitted,
433 MBadPixelsPix::kRelTimeOscillating);
434
435 if (IsLoGain())
436 FitLoGainArrays(fIntensCam ? (MCalibrationCam&)(*fIntensCam->GetCam()) : (MCalibrationCam&)(*fCam),
437 fIntensBad ? (*fIntensBad->GetCam()) : *fBadPixels,
438 MBadPixelsPix::kRelTimeNotFitted,
439 MBadPixelsPix::kRelTimeOscillating);
440
441 return kTRUE;
442}
443
444// --------------------------------------------------------------------------
445//
446// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
447// - MBadPixelsPix::kRelTimeNotFitted
448// - MBadPixelsPix::kRelTimeOscillating
449//
450void MHCalibrationRelTimeCam::FinalizeBadPixels()
451{
452
453 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
454
455 for (Int_t i=0; i<badcam->GetSize(); i++)
456 {
457 MBadPixelsPix &bad = (*badcam)[i];
458
459 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeNotFitted ))
460 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
461
462 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeOscillating))
463 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
464
465 }
466}
467
468// --------------------------------------------------------------------------
469//
470// The types are as follows:
471//
472// Fitted values:
473// ==============
474//
475// 0: Fitted Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMean()
476// 1: Error Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMeanErr()
477// 2: Sigma fitted Relative Arrival Time in FADC slices (MHGausEvents::GetSigma()
478// 3: Error Sigma Relative Arrival Time in FADC slices (MHGausEvents::GetSigmaErr()
479//
480// Useful variables derived from the fit results:
481// =============================================
482//
483// 4: Returned probability of Gauss fit (calls: MHGausEvents::GetProb())
484//
485// Localized defects:
486// ==================
487//
488// 5: Gaus fit not OK (calls: MHGausEvents::IsGausFitOK())
489// 6: Fourier spectrum not OK (calls: MHGausEvents::IsFourierSpectrumOK())
490//
491Bool_t MHCalibrationRelTimeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
492{
493
494 if (fHiGainArray->GetSize() <= idx)
495 return kFALSE;
496
497 const MHCalibrationPix &pix = (*this)[idx];
498
499 switch (type)
500 {
501 case 0:
502 val = pix.GetMean();
503 break;
504 case 1:
505 val = pix.GetMeanErr();
506 break;
507 case 2:
508 val = pix.GetSigma();
509 break;
510 case 3:
511 val = pix.GetSigmaErr();
512 break;
513 case 4:
514 val = pix.GetProb();
515 break;
516 case 5:
517 if (!pix.IsGausFitOK())
518 val = 1.;
519 break;
520 case 6:
521 if (!pix.IsFourierSpectrumOK())
522 val = 1.;
523 break;
524 default:
525 return kFALSE;
526 }
527 return kTRUE;
528}
529
530// --------------------------------------------------------------------------
531//
532// Calls MHCalibrationPix::DrawClone() for pixel idx
533//
534void MHCalibrationRelTimeCam::DrawPixelContent(Int_t idx) const
535{
536 (*this)[idx].DrawClone();
537}
538
539// -----------------------------------------------------------------------------
540//
541// Default draw:
542//
543// Displays the averaged areas, both High Gain and Low Gain
544//
545// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
546//
547void MHCalibrationRelTimeCam::Draw(const Option_t *opt)
548{
549
550 const Int_t nareas = fAverageHiGainAreas->GetEntries();
551 if (nareas == 0)
552 return;
553
554 TString option(opt);
555 option.ToLower();
556
557 if (!option.Contains("datacheck"))
558 {
559 MHCalibrationCam::Draw(opt);
560 return;
561 }
562
563 //
564 // From here on , the datacheck - Draw
565 //
566 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
567 pad->SetBorderMode(0);
568 pad->Divide(1,nareas);
569
570 //
571 // Loop over inner and outer pixels
572 //
573 for (Int_t i=0; i<nareas;i++)
574 {
575
576 pad->cd(i+1);
577
578 MHCalibrationPix &hipix = GetAverageHiGainArea(i);
579 //
580 // Ask for Hi-Gain saturation
581 //
582 if (hipix.GetSaturated() > fNumHiGainSaturationLimit*hipix.GetHGausHist()->GetEntries() && IsLoGain())
583 {
584 MHCalibrationPix &lopix = GetAverageLoGainArea(i);
585 DrawDataCheckPixel(lopix,0.);
586 }
587 else
588 DrawDataCheckPixel(hipix,0.);
589 }
590}
591
592// -----------------------------------------------------------------------------
593//
594// Draw the average pixel for the datacheck:
595//
596// Displays the averaged areas, both High Gain and Low Gain
597//
598// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
599//
600void MHCalibrationRelTimeCam::DrawDataCheckPixel(MHCalibrationPix &pix, const Float_t refline)
601{
602
603 TVirtualPad *newpad = gPad;
604 newpad->Divide(1,2);
605 newpad->cd(1);
606
607 gPad->SetTicks();
608 if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
609 gPad->SetLogy();
610
611 TH1F *hist = pix.GetHGausHist();
612
613 TH1F *null = new TH1F("Null",hist->GetTitle(),100,pix.GetFirst() < -3.0 ? -3.0 : pix.GetFirst(),
614 pix.GetLast() > 3.0 ? 3.0 : pix.GetLast());
615
616 null->SetMaximum(1.1*hist->GetMaximum());
617 null->SetDirectory(NULL);
618 null->SetBit(kCanDelete);
619 null->SetStats(kFALSE);
620 //
621 // set the labels bigger
622 //
623 TAxis *xaxe = null->GetXaxis();
624 TAxis *yaxe = null->GetYaxis();
625 xaxe->CenterTitle();
626 yaxe->CenterTitle();
627 xaxe->SetTitleSize(0.07);
628 yaxe->SetTitleSize(0.07);
629 xaxe->SetTitleOffset(0.7);
630 yaxe->SetTitleOffset(0.55);
631 xaxe->SetLabelSize(0.06);
632 yaxe->SetLabelSize(0.06);
633
634 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
635 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
636
637 null->Draw();
638 hist->Draw("same");
639
640 gStyle->SetOptFit();
641
642 TF1 *fit = pix.GetFGausFit();
643
644 if (fit)
645 {
646 switch ( fColor )
647 {
648 case MCalibrationCam::kGREEN:
649 fit->SetLineColor(kGreen);
650 break;
651 case MCalibrationCam::kBLUE:
652 fit->SetLineColor(kBlue);
653 break;
654 case MCalibrationCam::kUV:
655 fit->SetLineColor(106);
656 break;
657 case MCalibrationCam::kCT1:
658 fit->SetLineColor(006);
659 break;
660 default:
661 fit->SetLineColor(kRed);
662 }
663 fit->Draw("same");
664 }
665
666 // DisplayRefLines(null,refline);
667
668 gPad->Modified();
669 gPad->Update();
670
671 newpad->cd(2);
672 gPad->SetTicks();
673
674 TH1F *null2 = new TH1F("Null2",hist->GetTitle(),100,0.,pix.GetEvents()->GetSize()/pix.GetEventFrequency());
675
676 null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
677 null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
678 null2->SetDirectory(NULL);
679 null2->SetBit(kCanDelete);
680 null2->SetStats(kFALSE);
681 //
682 // set the labels bigger
683 //
684 TAxis *xaxe2 = null2->GetXaxis();
685 TAxis *yaxe2 = null2->GetYaxis();
686 xaxe2->CenterTitle();
687 yaxe2->CenterTitle();
688 xaxe2->SetTitleSize(0.07);
689 yaxe2->SetTitleSize(0.07);
690 xaxe2->SetTitleOffset(0.7);
691 yaxe2->SetTitleOffset(0.55);
692 xaxe2->SetLabelSize(0.06);
693 yaxe2->SetLabelSize(0.06);
694
695 pix.CreateGraphEvents();
696 TGraph *gr = pix.GetGraphEvents();
697
698 xaxe2->SetTitle(gr->GetXaxis()->GetTitle());
699 yaxe2->SetTitle(gr->GetYaxis()->GetTitle());
700
701 null2->Draw();
702
703 pix.DrawEvents("same");
704
705 return;
706}
707
Note: See TracBrowser for help on using the repository browser.