source: trunk/MagicSoft/Mars/mhflux/MHEnergyEst.cc@ 6917

Last change on this file since 6917 was 6917, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 12.6 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): Wolfgang Wittek 4/2002 <mailto:wittek@mppmu.mpg.de>
19! Author(s): Abelardo Moralejo 5/2003 <mailto:moralejo@pd.infn.it>
20! Author(s): Thomas Bretz 1/2005 <mailto:tbretz@astro.uni-wuerzburg.de>
21!
22! Copyright: MAGIC Software Development, 2000-2005
23!
24!
25\* ======================================================================== */
26
27//////////////////////////////////////////////////////////////////////////////
28//
29// MHEnergyEst
30//
31// calculates the migration matrix E-est vs. E-true
32// for different bins in Theta
33//
34//////////////////////////////////////////////////////////////////////////////
35#include "MHEnergyEst.h"
36
37#include <TLine.h>
38#include <TCanvas.h>
39#include <TProfile.h>
40
41#include "MString.h"
42
43#include "MMcEvt.hxx"
44
45#include "MEnergyEst.h"
46#include "MBinning.h"
47#include "MParList.h"
48#include "MParameters.h"
49#include "MHMatrix.h"
50
51#include "MLog.h"
52#include "MLogManip.h"
53
54ClassImp(MHEnergyEst);
55
56using namespace std;
57
58// --------------------------------------------------------------------------
59//
60// Default Constructor. It sets name and title of the histogram.
61//
62MHEnergyEst::MHEnergyEst(const char *name, const char *title)
63 : fMcEvt(0), fEnergy(0), fResult(0), fMatrix(0)
64{
65 //
66 // set the name and title of this object
67 //
68 fName = name ? name : "MHEnergyEst";
69 fTitle = title ? title : "3-D histogram E-true E-est Theta";
70
71 fHEnergy.SetDirectory(NULL);
72 fHEnergy.SetName("EnergyEst");
73 fHEnergy.SetTitle("Histogram in E_{est}, E_{mc} and \\Theta");
74 fHEnergy.SetXTitle("E_{est} [GeV]");
75 fHEnergy.SetYTitle("E_{mc} [GeV]");
76 fHEnergy.SetZTitle("\\Theta [\\circ]");
77
78 fHResolution.SetDirectory(NULL);
79 fHResolution.SetName("EnergyRes");
80 fHResolution.SetTitle("Histogram in \\frac{\\Delta E}{E_{mc}} vs E_{est} and E_{mc}");
81 fHResolution.SetXTitle("E_{est} [GeV]");
82 fHResolution.SetYTitle("E_{mc} [GeV]");
83 fHResolution.SetZTitle("\\frac{E_{est} - E_{mc}}{E_{mc}}");
84
85 fHImpact.SetDirectory(NULL);
86 fHImpact.SetName("Impact");
87 fHImpact.SetTitle("\\frac{\\Delta E}{E} vs Impact parameter");
88 fHImpact.SetXTitle("Impact parameter [m]");
89 fHImpact.SetYTitle("\\frac{E_{est} - E_{mc}}{E_{mc}}");
90
91 fHEnergy.Sumw2();
92 fHImpact.Sumw2();
93 fHResolution.Sumw2();
94
95 MBinning binsi, binse, binst, binsr;
96 binse.SetEdgesLog(15, 10, 1000000);
97 binst.SetEdgesCos(50, 0, 60);
98 binsi.SetEdges(10, 0, 400);
99 binsr.SetEdges(10, 0, 1);
100
101 SetBinning(&fHEnergy, &binse, &binse, &binst);
102 SetBinning(&fHImpact, &binsi, &binsr);
103 SetBinning(&fHResolution, &binse, &binse, &binsr);
104}
105
106// --------------------------------------------------------------------------
107//
108// Set the binnings and prepare the filling of the histograms
109//
110Bool_t MHEnergyEst::SetupFill(const MParList *plist)
111{
112 if (!fMatrix)
113 {
114 fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
115 if (!fMcEvt)
116 {
117 *fLog << err << "MMcEvt not found... aborting." << endl;
118 return kFALSE;
119 }
120 }
121
122 fEnergy = (MEnergyEst*)plist->FindObject("MEnergyEst");
123 if (!fEnergy)
124 {
125 *fLog << err << "MEnergyEst not found... aborting." << endl;
126 return kFALSE;
127 }
128
129 fResult = (MParameterD*)const_cast<MParList*>(plist)->FindCreateObj("MParameterD", "MinimizationValue");
130 if (!fResult)
131 return kFALSE;
132
133 MBinning binst, binse, binsi, binsr;
134 binse.SetEdges(fHEnergy, 'x');
135 binst.SetEdges(fHEnergy, 'z');
136 binsi.SetEdges(fHImpact, 'x');
137 binsr.SetEdges(fHImpact, 'y');
138
139 binst.SetEdges(*plist, "BinningTheta");
140 binse.SetEdges(*plist, "BinningEnergyEst");
141 binsi.SetEdges(*plist, "BinningImpact");
142 binsr.SetEdges(*plist, "BinningEnergyRes");
143
144 SetBinning(&fHEnergy, &binse, &binse, &binst);
145 SetBinning(&fHImpact, &binsi, &binsr);
146 SetBinning(&fHResolution, &binse, &binse, &binsr);
147
148 fChisq = 0;
149 fHEnergy.Reset();
150 fHImpact.Reset();
151 fHResolution.Reset();
152
153 return kTRUE;
154}
155
156// --------------------------------------------------------------------------
157//
158// Fill the histogram
159//
160Bool_t MHEnergyEst::Fill(const MParContainer *par, const Stat_t w)
161{
162 const Double_t eest = fEnergy->GetEnergy();
163 const Double_t etru = fMatrix ? GetVal(0) : fMcEvt->GetEnergy();
164 const Double_t imp = fMatrix ? GetVal(1) : fMcEvt->GetImpact()/100;
165 const Double_t theta = fMatrix ? GetVal(2) : fMcEvt->GetTelescopeTheta()*TMath::RadToDeg();
166 const Double_t res = (eest-etru)/etru;
167
168 fHEnergy.Fill(eest, etru, theta, w);
169 fHResolution.Fill(eest, etru, TMath::Abs(res), w);
170 fHImpact.Fill(imp, TMath::Abs(res), w);
171
172 fChisq += TMath::Abs(res);//*res;
173 fBias += res;
174
175 return kTRUE;
176}
177
178Bool_t MHEnergyEst::Finalize()
179{
180 fChisq /= GetNumExecutions();
181 fBias /= GetNumExecutions();
182
183 Double_t res = fChisq; //TMath::Sqrt(fChisq - fBias*fBias);
184
185 fResult->SetVal(TMath::IsNaN(res)?0:res);/// GetNumExecutions());
186
187 *fLog << all << "Mean Energy Resoltuion: " << Form("%.1f%%", fResult->GetVal()*100) << endl;
188 *fLog << all << "Energy Bias at: " << Form("%.1f%%", fBias*100) << endl;
189
190 return kTRUE;
191}
192
193void MHEnergyEst::Paint(Option_t *opt)
194{
195 TVirtualPad *pad = gPad;
196
197 pad->cd(1);
198
199 TH1D *hx=0;
200 TH1D *hy=0;
201
202 if (pad->GetPad(1))
203 {
204 pad->GetPad(1)->cd(1);
205
206 if ((hx=(TH1D*)gPad->FindObject("EnergyEst_ex")))
207 {
208 TH1D *h2 = (TH1D*)fHEnergy.Project3D("dum_ex");
209 hx->Reset();
210 hx->Add(h2);
211 delete h2;
212 }
213
214 if ((hy=(TH1D*)gPad->FindObject("EnergyEst_ey")))
215 {
216 TH1D *h2 = (TH1D*)fHEnergy.Project3D("dum_ey");
217 hy->Reset();
218 hy->Add(h2);
219 delete h2;
220 }
221
222 if (hx && hy)
223 {
224 hy->SetMaximum();
225 hy->SetMaximum(TMath::Max(hx->GetMaximum(), hy->GetMaximum())*1.2);
226 if (hy->GetMaximum()>0)
227 gPad->SetLogy();
228 }
229
230 if (pad->GetPad(1)->GetPad(2))
231 {
232 pad->GetPad(1)->GetPad(2)->cd(1);
233 if ((hx=(TH1D*)gPad->FindObject("EnergyEst_ez")))
234 {
235 TH1D *h2 = (TH1D*)fHEnergy.Project3D("dum_ez");
236 hx->Reset();
237 hx->Add(h2);
238 delete h2;
239 }
240
241 //pad->GetPad(1)->GetPad(2)->cd(2);
242 ///*h =*/ fHImpact.ProjectionX("Impact", -1, 9999, "e");
243 }
244 }
245
246 if (pad->GetPad(2))
247 {
248 pad->GetPad(2)->cd(1);
249 UpdatePlot(fHEnergy, "yx", kTRUE);
250
251 TLine *l = (TLine*)gPad->FindObject("TLine");
252 if (l)
253 {
254 const Float_t min = TMath::Max(fHEnergy.GetXaxis()->GetXmin(), fHEnergy.GetYaxis()->GetXmin());
255 const Float_t max = TMath::Min(fHEnergy.GetXaxis()->GetXmax(), fHEnergy.GetYaxis()->GetXmax());
256
257 l->SetX1(min);
258 l->SetX2(max);
259 l->SetY1(min);
260 l->SetY2(max);
261 }
262
263 pad->GetPad(2)->cd(2);
264 UpdatePlot(fHResolution, "zy");
265
266 pad->GetPad(2)->cd(3);
267 UpdatePlot(fHResolution, "zx");
268 }
269}
270
271void MHEnergyEst::UpdatePlot(TH3 &h, const char *how, Bool_t logy)
272{
273 TH2D *hyx=0;
274 if (!(hyx=(TH2D*)gPad->FindObject(MString::Form("%s_%s", h.GetName(), how))))
275 return;
276
277 TH2D *h2 = (TH2D*)h.Project3D(MString::Form("dum_%s", how));
278 hyx->Reset();
279 hyx->Add(h2);
280 delete h2;
281
282 TH1D *hx = 0;
283 if ((hx=(TH1D*)gPad->FindObject(Form("Prof%s", h.GetName()))))
284 {
285 hx = hyx->ProfileX(Form("Prof%s", h.GetName()), -1, 9999, "s");
286
287 if (logy && hx->GetMaximum()>0)
288 gPad->SetLogy();
289 }
290}
291
292TH1 *MHEnergyEst::MakePlot(TH3 &h, const char *how)
293{
294 gPad->SetBorderMode(0);
295 gPad->SetLogx();
296
297 //gROOT->GetListOfCleanups()->Add(gPad); // WHY?
298
299 TH2D *h2 = (TH2D*)h.Project3D(how);
300 h2->SetDirectory(NULL);
301 h2->SetBit(kCanDelete);
302 h2->SetFillColor(kBlue);
303
304 TH1D *h1 = h2->ProfileX(Form("Prof%s", h.GetName()), -1, 9999, "s");
305 h1->SetDirectory(NULL);
306 h1->SetBit(kCanDelete);
307 h1->SetLineWidth(2);
308 h1->SetLineColor(kRed);
309 h1->SetFillStyle(4000);
310 h1->SetStats(kFALSE);
311
312
313 //h1->Draw("E3");
314 h2->Draw();
315 h1->Draw("Chistsame");
316
317 return h1;
318}
319
320
321// --------------------------------------------------------------------------
322//
323// Draw the histogram
324//
325void MHEnergyEst::Draw(Option_t *opt)
326{
327 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
328
329 // Do the projection before painting the histograms into
330 // the individual pads
331 AppendPad("");
332
333 pad->SetBorderMode(0);
334 pad->Divide(2, 1, 1e-10, 1e-10);
335
336 TH1 *h;
337
338 pad->cd(1);
339 gPad->SetBorderMode(0);
340
341 gPad->Divide(1, 2, 1e-10, 1e-10);
342
343 TVirtualPad *pad2 = gPad;
344
345 pad2->cd(1);
346 gPad->SetBorderMode(0);
347
348 gPad->SetLogx();
349 h = (TH1D*)fHEnergy.Project3D("ey");
350 h->SetBit(TH1::kNoStats);
351 h->SetTitle("Energy disribution: Monte Carlo E_{mc} (black), Estimated E_{est} (green)");
352 h->SetXTitle("E [GeV]"); // E_mc
353 h->SetYTitle("Counts");
354 h->SetBit(kCanDelete);
355 h->SetDirectory(NULL);
356 h->SetMarkerStyle(kFullDotMedium);
357 h->Draw();
358
359 h = (TH1D*)fHEnergy.Project3D("ex");
360 h->SetBit(TH1::kNoTitle|TH1::kNoStats);
361 h->SetXTitle("E [GeV]"); // E_est
362 h->SetYTitle("Counts");
363 h->SetBit(kCanDelete);
364 h->SetDirectory(NULL);
365 h->SetMarkerStyle(kFullDotMedium);
366 h->SetLineColor(kGreen);
367 h->SetMarkerColor(kGreen);
368 h->Draw("same");
369
370 // FIXME: LEGEND
371
372 pad2->cd(2);
373 gPad->SetBorderMode(0);
374
375 TVirtualPad *pad3 = gPad;
376 pad3->Divide(2, 1, 1e-10, 1e-10);
377 pad3->cd(1);
378 gPad->SetBorderMode(0);/*
379 h = fHImpact.ProjectionX("Impact", -1, 9999, "e");
380 h->SetBit(TH1::kNoStats);
381 h->SetTitle("Distribution of Impact");
382 h->SetDirectory(NULL);
383 h->SetXTitle("Impact [m]");
384 h->SetBit(kCanDelete);
385 h->Draw();*/
386 h = fHEnergy.Project3D("ez");
387 h->SetTitle("Distribution of Theta");
388 h->SetBit(TH1::kNoStats);
389 h->SetDirectory(NULL);
390 h->SetXTitle("\\Theta [\\circ]");
391 h->SetBit(kCanDelete);
392 h->Draw();
393
394 pad3->cd(2);
395 gPad->SetBorderMode(0);
396
397 pad->cd(2);
398 gPad->SetBorderMode(0);
399
400 gPad->Divide(1, 3, 1e-10, 1e-10);
401 pad2 = gPad;
402
403 pad2->cd(1);
404 h = MakePlot(fHEnergy, "yx");
405 h->SetXTitle("E_{mc} [GeV]");
406 h->SetYTitle("E_{est} [GeV]");
407 h->SetTitle("Estimated energy E_{est} vs Monte Carlo energy E_{mc}");
408
409 TLine line;
410 line.DrawLine(0,0,1,1);
411
412 pad2->cd(2);
413 h = MakePlot(fHResolution, "zy");
414 h->SetXTitle("E_{mc} [GeV]");
415 h->SetYTitle("\\frac{E_{est} - E_{mc}}{E_{mc}}");
416 h->SetTitle("Energy resolution \\Delta E / E vs Monte Carlo energy E_{mc}");
417 h->SetMinimum(0);
418 h->SetMaximum(1);
419
420 pad2->cd(3);
421 h = MakePlot(fHResolution, "zx");
422 h->SetXTitle("E_{est} [GeV]");
423 h->SetYTitle("\\frac{E_{est} - E_{mc}}{E_{mc}}");
424 h->SetTitle("Energy Resolution \\Delta E / E vs estimated Energy E_{est}");
425 h->SetMinimum(0);
426 h->SetMaximum(1);
427}
428
429// --------------------------------------------------------------------------
430//
431// Returns the mapped value from the Matrix
432//
433Double_t MHEnergyEst::GetVal(Int_t i) const
434{
435 return (*fMatrix)[fMap[i]];
436}
437
438// --------------------------------------------------------------------------
439//
440// You can use this function if you want to use a MHMatrix instead of the
441// given containers. This function adds all necessary columns to the
442// given matrix. Afterward you should fill the matrix with the corresponding
443// data (eg from a file by using MHMatrix::Fill). If you now loop
444// through the matrix (eg using MMatrixLoop) MEnergyEstParam2::Process
445// will take the values from the matrix instead of the containers.
446//
447void MHEnergyEst::InitMapping(MHMatrix *mat)
448{
449 if (fMatrix)
450 return;
451
452 fMatrix = mat;
453
454 fMap[0] = fMatrix->AddColumn("MMcEvt.fEnergy");
455 fMap[1] = fMatrix->AddColumn("MMcEvt.fImpact/100");
456 fMap[2] = fMatrix->AddColumn("MMcEvt.fTelescopeTheta*kRad2Deg");
457}
458
459void MHEnergyEst::StopMapping()
460{
461 fMatrix = NULL;
462}
463
Note: See TracBrowser for help on using the repository browser.