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

Last change on this file since 6938 was 6932, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 13.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): 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 lg(E)}{lg(E_{mc})} vs E_{est} and E_{mc}");
81 fHResolution.SetXTitle("E_{est} [GeV]");
82 fHResolution.SetYTitle("E_{mc} [GeV]");
83 fHResolution.SetZTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(E_{mc})}");
84
85 fHImpact.SetDirectory(NULL);
86 fHImpact.SetName("Impact");
87 fHImpact.SetTitle("\\frac{\\Delta lg(E)}{lg(E)} vs Impact parameter");
88 fHImpact.SetXTitle("Impact parameter [m]");
89 fHImpact.SetYTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(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(50, -1.3, 1.3);
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 = (MParameterD*)plist->FindObject("MEnergyEst", "MParameterD");
123 if (!fEnergy)
124 {
125 *fLog << err << "MEnergyEst [MParameterD] 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->GetVal();
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 = (log10(eest)-log10(etru));///log10(etru);
167 const Double_t resE = (eest-etru)/etru;
168
169 fHEnergy.Fill(eest, etru, theta, w);
170 fHResolution.Fill(eest, etru, resE, w);
171 fHImpact.Fill(imp, resE, w);
172
173 fChisq += res*res;
174 fBias += res;
175
176 return kTRUE;
177}
178
179Bool_t MHEnergyEst::Finalize()
180{
181 fChisq /= GetNumExecutions();
182 fBias /= GetNumExecutions();
183
184 Double_t sigma = TMath::Sqrt(fChisq);//+TMath::Abs(fBias);
185 fResult->SetVal(sigma);
186
187 Print();
188
189 return kTRUE;
190}
191
192void MHEnergyEst::Print(Option_t *o) const
193{
194 *fLog << all;
195 *fLog << "Mean log10(Energy) Resoltion: " << Form("%5.2f%%", TMath::Sqrt(fChisq-fBias*fBias)) << endl;
196 *fLog << "Mean log10(Energy) Bias: " << Form("%5.2f%%", fBias) << endl;
197}
198
199void MHEnergyEst::Paint(Option_t *opt)
200{
201 TVirtualPad *pad = gPad;
202
203 pad->cd(1);
204
205 TH1D *hx=0;
206 TH1D *hy=0;
207
208 if (pad->GetPad(1))
209 {
210 pad->GetPad(1)->cd(1);
211
212 if ((hx=(TH1D*)gPad->FindObject("EnergyEst_ex")))
213 {
214 TH1D *h2 = (TH1D*)fHEnergy.Project3D("dum_ex");
215 hx->Reset();
216 hx->Add(h2);
217 delete h2;
218 }
219
220 if ((hy=(TH1D*)gPad->FindObject("EnergyEst_ey")))
221 {
222 TH1D *h2 = (TH1D*)fHEnergy.Project3D("dum_ey");
223 hy->Reset();
224 hy->Add(h2);
225 delete h2;
226 }
227
228 if (hx && hy)
229 {
230 hy->SetMaximum();
231 hy->SetMaximum(TMath::Max(hx->GetMaximum(), hy->GetMaximum())*1.2);
232 if (hy->GetMaximum()>0)
233 gPad->SetLogy();
234 }
235
236 if (pad->GetPad(1)->GetPad(2))
237 {
238 pad->GetPad(1)->GetPad(2)->cd(1);
239 if ((hx=(TH1D*)gPad->FindObject("EnergyEst_ez")))
240 {
241 TH1D *h2 = (TH1D*)fHEnergy.Project3D("dum_ez");
242 hx->Reset();
243 hx->Add(h2);
244 delete h2;
245 }
246
247 pad->GetPad(1)->GetPad(2)->cd(2);
248 /*=*/fHResolution.ProjectionZ("Resolution");
249 ///*h =*/ fHImpact.ProjectionX("Impact", -1, 9999, "e");
250 }
251 }
252
253 if (pad->GetPad(2))
254 {
255 pad->GetPad(2)->cd(1);
256 UpdatePlot(fHEnergy, "yx", kTRUE);
257
258 TLine *l = (TLine*)gPad->FindObject("TLine");
259 if (l)
260 {
261 const Float_t min = TMath::Max(fHEnergy.GetXaxis()->GetXmin(), fHEnergy.GetYaxis()->GetXmin());
262 const Float_t max = TMath::Min(fHEnergy.GetXaxis()->GetXmax(), fHEnergy.GetYaxis()->GetXmax());
263
264 l->SetX1(min);
265 l->SetX2(max);
266 l->SetY1(min);
267 l->SetY2(max);
268 }
269
270 pad->GetPad(2)->cd(2);
271 UpdatePlot(fHResolution, "zy");
272
273 pad->GetPad(2)->cd(3);
274 UpdatePlot(fHResolution, "zx");
275 }
276}
277
278void MHEnergyEst::UpdatePlot(TH3 &h, const char *how, Bool_t logy)
279{
280 TH2D *hyx=0;
281 if (!(hyx=(TH2D*)gPad->FindObject(MString::Form("%s_%s", h.GetName(), how))))
282 return;
283
284 TH2D *h2 = (TH2D*)h.Project3D(MString::Form("dum_%s", how));
285 hyx->Reset();
286 hyx->Add(h2);
287 delete h2;
288
289 TH1D *hx = 0;
290 if ((hx=(TH1D*)gPad->FindObject(Form("Prof%s", h.GetName()))))
291 {
292 hx = hyx->ProfileX(Form("Prof%s", h.GetName()), -1, 9999, "s");
293
294 if (logy && hx->GetMaximum()>0)
295 gPad->SetLogy();
296 }
297}
298
299TH1 *MHEnergyEst::MakePlot(TH3 &h, const char *how)
300{
301 gPad->SetBorderMode(0);
302 gPad->SetLogx();
303 gPad->SetGridx();
304 gPad->SetGridy();
305
306 // Results in crashes....
307 //gROOT->GetListOfCleanups()->Add(gPad); // WHY?
308
309 TH2D *h2 = (TH2D*)h.Project3D(how);
310 h2->SetDirectory(NULL);
311 h2->SetBit(kCanDelete);
312 h2->SetFillColor(kBlue);
313 h2->SetLineColor(kRed);
314
315 TH1D *h1 = h2->ProfileX(Form("Prof%s", h.GetName()), -1, 9999, "s");
316 h1->SetDirectory(NULL);
317 h1->SetBit(kCanDelete);
318 h1->SetLineWidth(2);
319 h1->SetLineColor(kRed);
320 h1->SetFillStyle(4000);
321 h1->SetStats(kFALSE);
322
323 h2->Draw("");
324 h1->Draw("E3 hist C same");
325// h1->Draw("Chistsame");
326
327 return h1;
328}
329
330
331// --------------------------------------------------------------------------
332//
333// Draw the histogram
334//
335void MHEnergyEst::Draw(Option_t *opt)
336{
337 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
338
339 // Do the projection before painting the histograms into
340 // the individual pads
341 AppendPad("");
342
343 pad->SetBorderMode(0);
344 pad->Divide(2, 1, 1e-10, 1e-10);
345
346 TH1 *h;
347
348 pad->cd(1);
349 gPad->SetBorderMode(0);
350
351 gPad->Divide(1, 2, 1e-10, 1e-10);
352
353 TVirtualPad *pad2 = gPad;
354
355 pad2->cd(1);
356 gPad->SetBorderMode(0);
357
358 gPad->SetGridx();
359 gPad->SetGridy();
360 gPad->SetLogx();
361 h = (TH1D*)fHEnergy.Project3D("ey");
362 h->SetBit(TH1::kNoStats);
363 h->SetTitle("Energy disribution: Monte Carlo E_{mc} (black), Estimated E_{est} (green)");
364 h->SetXTitle("E [GeV]"); // E_mc
365 h->SetYTitle("Counts");
366 h->SetBit(kCanDelete);
367 h->SetDirectory(NULL);
368 h->SetMarkerStyle(kFullDotMedium);
369 h->Draw();
370
371 h = (TH1D*)fHEnergy.Project3D("ex");
372 h->SetBit(TH1::kNoTitle|TH1::kNoStats);
373 h->SetXTitle("E [GeV]"); // E_est
374 h->SetYTitle("Counts");
375 h->SetBit(kCanDelete);
376 h->SetDirectory(NULL);
377 h->SetMarkerStyle(kFullDotMedium);
378 h->SetLineColor(kGreen);
379 h->SetMarkerColor(kGreen);
380 h->Draw("same");
381
382 // FIXME: LEGEND
383
384 pad2->cd(2);
385 gPad->SetBorderMode(0);
386
387 TVirtualPad *pad3 = gPad;
388 pad3->Divide(2, 1, 1e-10, 1e-10);
389 pad3->cd(1);
390 gPad->SetBorderMode(0);
391 gPad->SetGridx();
392 gPad->SetGridy();
393 h = fHEnergy.Project3D("ez");
394 h->SetTitle("Zenith Angle Distribution");
395 h->SetBit(TH1::kNoStats);
396 h->SetDirectory(NULL);
397 h->SetXTitle("\\Theta [\\circ]");
398 h->SetBit(kCanDelete);
399 h->Draw();
400
401 pad3->cd(2);
402 gPad->SetBorderMode(0);
403 gPad->SetGridx();
404 gPad->SetGridy();
405 /*
406 h = fHImpact.ProjectionX("Impact", -1, 9999, "e");
407 h->SetBit(TH1::kNoStats);
408 h->SetTitle("Distribution of Impact");
409 h->SetDirectory(NULL);
410 h->SetXTitle("Impact [m]");
411 h->SetBit(kCanDelete);
412 h->Draw();*/
413 // ----------------------------------------
414 h = fHResolution.ProjectionZ("Resolution");
415 //h->SetXTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(E_{mc})}");
416 h->SetYTitle("Counts");
417 h->SetTitle("Distribution of \\Delta lg(E) / lg(E_{mc})");
418 h->SetDirectory(NULL);
419 h->SetBit(kCanDelete);
420 h->GetXaxis()->SetRangeUser(-1.3, 1.3);
421 h->Draw("");
422 // ----------------------------------------
423
424 pad->cd(2);
425 gPad->SetBorderMode(0);
426
427 gPad->Divide(1, 3, 1e-10, 1e-10);
428 pad2 = gPad;
429
430 pad2->cd(1);
431 h = MakePlot(fHEnergy, "yx");
432 h->SetXTitle("E_{mc} [GeV]");
433 h->SetYTitle("E_{est} [GeV]");
434 h->SetTitle("Estimated energy E_{est} vs Monte Carlo energy E_{mc}");
435
436 TLine line;
437 line.DrawLine(0,0,1,1);
438
439 pad2->cd(2);
440 h = MakePlot(fHResolution, "zy");
441 h->SetXTitle("E_{mc} [GeV]");
442 h->SetYTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(E_{mc})}");
443 h->SetTitle("Energy resolution \\Delta lg(E) / lg(E) vs Monte Carlo energy E_{mc}");
444 h->SetMinimum(-1.3);
445 h->SetMaximum(1.3);
446
447 pad2->cd(3);
448 h = MakePlot(fHResolution, "zx");
449 h->SetXTitle("E_{est} [GeV]");
450 h->SetYTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(E_{mc})}");
451 h->SetTitle("Energy Resolution \\Delta lg(E) / lg(E) vs estimated Energy E_{est}");
452 h->SetMinimum(-1.3);
453 h->SetMaximum(1.3);
454}
455
456// --------------------------------------------------------------------------
457//
458// Returns the mapped value from the Matrix
459//
460Double_t MHEnergyEst::GetVal(Int_t i) const
461{
462 return (*fMatrix)[fMap[i]];
463}
464
465// --------------------------------------------------------------------------
466//
467// You can use this function if you want to use a MHMatrix instead of the
468// given containers. This function adds all necessary columns to the
469// given matrix. Afterward you should fill the matrix with the corresponding
470// data (eg from a file by using MHMatrix::Fill). If you now loop
471// through the matrix (eg using MMatrixLoop) MEnergyEstParam2::Process
472// will take the values from the matrix instead of the containers.
473//
474void MHEnergyEst::InitMapping(MHMatrix *mat)
475{
476 if (fMatrix)
477 return;
478
479 fMatrix = mat;
480
481 fMap[0] = fMatrix->AddColumn("MMcEvt.fEnergy");
482 fMap[1] = fMatrix->AddColumn("MMcEvt.fImpact/100");
483 fMap[2] = fMatrix->AddColumn("MMcEvt.fTelescopeTheta*kRad2Deg");
484}
485
486void MHEnergyEst::StopMapping()
487{
488 fMatrix = NULL;
489}
490
Note: See TracBrowser for help on using the repository browser.