source: trunk/Mars/mpedestal/MHPedestalCor.cc@ 10056

Last change on this file since 10056 was 9851, checked in by tbretz, 14 years ago
Changed MH::SetBinning and similar functions to take references instead of pointers as arguments. For convenience wrappers for the old style functions are provided.
File size: 7.2 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MHPedestalCor.cc,v 1.4 2008-11-11 11:49:49 tbretz Exp $
3! --------------------------------------------------------------------------
4!
5! *
6! * This file is part of MARS, the MAGIC Analysis and Reconstruction
7! * Software. It is distributed to you in the hope that it can be a useful
8! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
9! * It is distributed WITHOUT ANY WARRANTY.
10! *
11! * Permission to use, copy, modify and distribute this software and its
12! * documentation for any purpose is hereby granted without fee,
13! * provided that the above copyright notice appear in all copies and
14! * that both that copyright notice and this permission notice appear
15! * in supporting documentation. It is provided "as is" without express
16! * or implied warranty.
17! *
18!
19!
20! Author(s): Thomas Bretz, 10/2006 <mailto:tbretz@astro.uni-wuerzburg.de>
21!
22! Copyright: MAGIC Software Development, 2000-2006
23!
24!
25\* ======================================================================== */
26
27/////////////////////////////////////////////////////////////////////////////
28//
29// MHPedestalCor
30//
31// A histogram to get the pedestal autocorrelation matrix. Tests have shown
32// the the autocorrelation does not depend on the actual slice so a profile
33// with the autocorrelation depeding on the distance between two slices
34// is filled.
35//
36// Functions to provide the autocorrelation as a TH2D or a TMatrix are
37// provided.
38//
39// Input:
40// - MPedestalSubtractedEvt (from MFillH)
41//
42/////////////////////////////////////////////////////////////////////////////
43#include "MHPedestalCor.h"
44
45#include <TH2.h>
46#include <TCanvas.h>
47
48#include "MLog.h"
49#include "MLogManip.h"
50
51#include "MParList.h"
52#include "MCamEvent.h"
53
54#include "MGeomCam.h"
55
56#include "MBinning.h"
57#include "MPedestalSubtractedEvt.h"
58
59ClassImp(MHPedestalCor);
60
61using namespace std;
62
63const TString MHPedestalCor::gsDefName = "MHPedestalCor";
64const TString MHPedestalCor::gsDefTitle = "Histogram for autocorrelation";
65
66// --------------------------------------------------------------------------
67//
68// Initialize the name and title of the task. Set fType to 0
69//
70MHPedestalCor::MHPedestalCor(const char *name, const char *title)
71{
72 //
73 // set the name and title of this object
74 //
75 fName = name ? name : gsDefName.Data();
76 fTitle = title ? title : gsDefTitle.Data();
77 /*
78 fHist.SetNameTitle("AutoCor", "Autocorrelation of slices");
79 fHist.SetDirectory(0);
80
81 fHist2.SetNameTitle("AutoCorProf", "Autocorrelation from Profile");
82 fHist2.SetDirectory(0);
83
84 fHist3.SetNameTitle("Dev", "Deviation of Profile from plain correlation");
85 fHist3.SetDirectory(0);
86
87 MBinning binsx(61, -30.5, 30.5);
88 MBinning binsy(30, -0.5, 29.5);
89
90 MH::SetBinning(&fHist, &binsy, &binsy);
91 MH::SetBinning(&fHist2, &binsy, &binsy);
92 MH::SetBinning(&fHist3, &binsy, &binsy);
93 */
94 fProf.SetNameTitle("AutoCorP", "Profile of auto correlation");
95 fProf.SetDirectory(0);
96
97 const MBinning binsx(15, -0.5, 14.5);
98 MH::SetBinning(fProf, binsx);
99}
100
101// --------------------------------------------------------------------------
102//
103// Fill the histograms with data from a MCamEvent-Container.
104//
105Int_t MHPedestalCor::Fill(const MParContainer *par, const Stat_t w)
106{
107 const MPedestalSubtractedEvt *evt = dynamic_cast<const MPedestalSubtractedEvt*>(par);
108 if (!evt)
109 {
110 *fLog << err << dbginf << "Got no MCamEvent as argument of Fill()..." << endl;
111 return kERROR;
112 }
113
114 // Implement an extraction range
115 // Implement a check as in PedCalcFromLoGain
116
117 const Int_t np = evt->GetNumPixels();
118 const Int_t ns = evt->GetNumSamples();
119
120 Int_t fCheckWinFirst = 0;
121 Int_t fCheckWinLast = ns;
122
123 Int_t fExtractWinFirst = 17;
124 Int_t fExtractWinLast = ns;
125
126 Float_t fMaxSignalVar = 40;
127 Float_t fMaxSignalAbs = 250;
128
129 for (int k=0; k<np; k++)
130 {
131 // This is the fast workaround to put hi- and lo-gains together
132 USample_t *slices = evt->GetSamplesRaw(k);//pixel.GetSamples();
133
134 USample_t max = 0;
135 USample_t min = (USample_t)-1;
136
137 // Find the maximum and minimum signal per slice in the high gain window
138 for (USample_t *slice=slices+fCheckWinFirst; slice<slices+fCheckWinLast; slice++)
139 {
140 if (*slice > max)
141 max = *slice;
142 if (*slice < min)
143 min = *slice;
144 }
145
146 // If the maximum in the high gain window is smaller than
147 if (max-min>=fMaxSignalVar || max>=fMaxSignalAbs)
148 continue;
149
150 const Float_t *sig = evt->GetSamples(k);
151
152 for (int i=fExtractWinFirst; i<fExtractWinLast; i++)
153 {
154 const Double_t s2 = sig[i]*sig[i];
155 // fHist.Fill(i, i, s2);
156 fProf.Fill(0., s2);
157
158 for (int j=fExtractWinFirst; j<fExtractWinLast; j++)
159 {
160 const Double_t sij = sig[i]*sig[j];
161
162// fHist.Fill(i, j, sij);
163// fHist.Fill(j, i, sij);
164
165 fProf.Fill(i-j, sij);
166// fProf.Fill(j-i, sij);
167 }
168 }
169 }
170
171 return kTRUE;
172}
173
174// --------------------------------------------------------------------------
175//
176// Return the autocorrelation matrix as a TH2D
177//
178void MHPedestalCor::GetAutoCorrelation(TH2 &h) const
179{
180 const Int_t n = fProf.GetNbinsX()*2+1;
181
182 const Axis_t xmax = fProf.GetXaxis()->GetXmax();
183
184 const MBinning bins(n, -xmax, xmax);
185 MH::SetBinning(h, bins, bins);
186
187 for (int x=0; x<n; x++)
188 for (int y=0; y<n; y++)
189 h.SetBinContent(x+1, y+1, fProf.GetBinContent(TMath::Abs(x-y)));
190
191}
192
193// --------------------------------------------------------------------------
194//
195// Return the autocorrelation into the TMatrix
196//
197void MHPedestalCor::GetAutoCorrelation(TMatrix &m) const
198{
199 const Int_t n = fProf.GetNbinsX()*2+1;
200 m.ResizeTo(n, n);
201
202 for (int x=0; x<n; x++)
203 for (int y=0; y<n; y++)
204 m[x][y] = fProf.GetBinContent(TMath::Abs(x-y));
205}
206
207/*
208void MHPedestalCor::Paint(Option_t *)
209{
210 MH::SetPalette("pretty");
211
212 for (int x=0; x<fHist.GetNbinsX(); x++)
213 for (int y=0; y<fHist.GetNbinsX(); y++)
214 fHist2.SetBinContent(x+1, y+1, fProf.GetBinContent(fProf.GetXaxis()->FindFixBin(x-y)));
215
216 fHist.Copy(fHist3);
217 fHist3.Add(&fHist2, -(fHist.GetEntries()/30/30));
218 fHist3.Divide(&fHist);
219 fHist3.Scale(100);
220 fHist3.SetMaximum();
221}
222*/
223
224void MHPedestalCor::Draw(Option_t *)
225{
226 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
227 pad->SetBorderMode(0);
228 pad->SetFrameBorderMode(0);
229
230 //pad->Divide(2, 2, 0.001, 0.001);
231
232 AppendPad();
233/*
234 pad->cd(1);
235 gPad->SetBorderMode(0);
236 gPad->SetFrameBorderMode(0);
237 gPad->SetGridx();
238 gPad->SetGridy();
239 fHist.Draw("colz");
240
241 pad->cd(2);
242 gPad->SetBorderMode(0);
243 gPad->SetFrameBorderMode(0);
244 gPad->SetGridx();
245 gPad->SetGridy();
246 fHist2.Draw("colz");
247
248 pad->cd(3);
249 gPad->SetBorderMode(0);
250 gPad->SetFrameBorderMode(0);
251 gPad->SetGridx();
252 gPad->SetGridy();
253 fHist3.Draw("colz");
254
255 pad->cd(4);*/
256
257 gPad->SetBorderMode(0);
258 gPad->SetFrameBorderMode(0);
259 gPad->SetGridx();
260 gPad->SetGridy();
261 fProf.Draw();
262}
Note: See TracBrowser for help on using the repository browser.