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

Last change on this file since 6924 was 6924, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 13.4 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 = (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 = (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 *fLog << all << "Mean log10(Energy) Resoltion: " << Form("%.1f%%", TMath::Sqrt(fChisq-fBias*fBias)*100) << endl;
188 *fLog << all << "Mean log10(Energy) Bias: " << 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 /*=*/fHResolution.ProjectionZ("Resolution");
243 ///*h =*/ fHImpact.ProjectionX("Impact", -1, 9999, "e");
244 }
245 }
246
247 if (pad->GetPad(2))
248 {
249 pad->GetPad(2)->cd(1);
250 UpdatePlot(fHEnergy, "yx", kTRUE);
251
252 TLine *l = (TLine*)gPad->FindObject("TLine");
253 if (l)
254 {
255 const Float_t min = TMath::Max(fHEnergy.GetXaxis()->GetXmin(), fHEnergy.GetYaxis()->GetXmin());
256 const Float_t max = TMath::Min(fHEnergy.GetXaxis()->GetXmax(), fHEnergy.GetYaxis()->GetXmax());
257
258 l->SetX1(min);
259 l->SetX2(max);
260 l->SetY1(min);
261 l->SetY2(max);
262 }
263
264 pad->GetPad(2)->cd(2);
265 UpdatePlot(fHResolution, "zy");
266
267 pad->GetPad(2)->cd(3);
268 UpdatePlot(fHResolution, "zx");
269 }
270}
271
272void MHEnergyEst::UpdatePlot(TH3 &h, const char *how, Bool_t logy)
273{
274 TH2D *hyx=0;
275 if (!(hyx=(TH2D*)gPad->FindObject(MString::Form("%s_%s", h.GetName(), how))))
276 return;
277
278 TH2D *h2 = (TH2D*)h.Project3D(MString::Form("dum_%s", how));
279 hyx->Reset();
280 hyx->Add(h2);
281 delete h2;
282
283 TH1D *hx = 0;
284 if ((hx=(TH1D*)gPad->FindObject(Form("Prof%s", h.GetName()))))
285 {
286 hx = hyx->ProfileX(Form("Prof%s", h.GetName()), -1, 9999, "s");
287
288 if (logy && hx->GetMaximum()>0)
289 gPad->SetLogy();
290 }
291}
292
293TH1 *MHEnergyEst::MakePlot(TH3 &h, const char *how)
294{
295 gPad->SetBorderMode(0);
296 gPad->SetLogx();
297 gPad->SetGridx();
298 gPad->SetGridy();
299
300 //gROOT->GetListOfCleanups()->Add(gPad); // WHY?
301
302 TH2D *h2 = (TH2D*)h.Project3D(how);
303 h2->SetDirectory(NULL);
304 h2->SetBit(kCanDelete);
305 h2->SetFillColor(kBlue);
306 h2->SetLineColor(kRed);
307
308 TH1D *h1 = h2->ProfileX(Form("Prof%s", h.GetName()), -1, 9999, "s");
309 h1->SetDirectory(NULL);
310 h1->SetBit(kCanDelete);
311 h1->SetLineWidth(2);
312 h1->SetLineColor(kRed);
313 h1->SetFillStyle(4000);
314 h1->SetStats(kFALSE);
315
316 h2->Draw("");
317 h1->Draw("E3same");
318 h1->Draw("Chistsame");
319
320 return h1;
321}
322
323
324// --------------------------------------------------------------------------
325//
326// Draw the histogram
327//
328void MHEnergyEst::Draw(Option_t *opt)
329{
330 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
331
332 // Do the projection before painting the histograms into
333 // the individual pads
334 AppendPad("");
335
336 pad->SetBorderMode(0);
337 pad->Divide(2, 1, 1e-10, 1e-10);
338
339 TH1 *h;
340
341 pad->cd(1);
342 gPad->SetBorderMode(0);
343
344 gPad->Divide(1, 2, 1e-10, 1e-10);
345
346 TVirtualPad *pad2 = gPad;
347
348 pad2->cd(1);
349 gPad->SetBorderMode(0);
350
351 gPad->SetGridx();
352 gPad->SetGridy();
353 gPad->SetLogx();
354 h = (TH1D*)fHEnergy.Project3D("ey");
355 h->SetBit(TH1::kNoStats);
356 h->SetTitle("Energy disribution: Monte Carlo E_{mc} (black), Estimated E_{est} (green)");
357 h->SetXTitle("E [GeV]"); // E_mc
358 h->SetYTitle("Counts");
359 h->SetBit(kCanDelete);
360 h->SetDirectory(NULL);
361 h->SetMarkerStyle(kFullDotMedium);
362 h->Draw();
363
364 h = (TH1D*)fHEnergy.Project3D("ex");
365 h->SetBit(TH1::kNoTitle|TH1::kNoStats);
366 h->SetXTitle("E [GeV]"); // E_est
367 h->SetYTitle("Counts");
368 h->SetBit(kCanDelete);
369 h->SetDirectory(NULL);
370 h->SetMarkerStyle(kFullDotMedium);
371 h->SetLineColor(kGreen);
372 h->SetMarkerColor(kGreen);
373 h->Draw("same");
374
375 // FIXME: LEGEND
376
377 pad2->cd(2);
378 gPad->SetBorderMode(0);
379
380 TVirtualPad *pad3 = gPad;
381 pad3->Divide(2, 1, 1e-10, 1e-10);
382 pad3->cd(1);
383 gPad->SetBorderMode(0);
384 gPad->SetGridx();
385 gPad->SetGridy();
386 h = fHEnergy.Project3D("ez");
387 h->SetTitle("Zenith Angle Distribution");
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 gPad->SetGridx();
397 gPad->SetGridy();
398 /*
399 h = fHImpact.ProjectionX("Impact", -1, 9999, "e");
400 h->SetBit(TH1::kNoStats);
401 h->SetTitle("Distribution of Impact");
402 h->SetDirectory(NULL);
403 h->SetXTitle("Impact [m]");
404 h->SetBit(kCanDelete);
405 h->Draw();*/
406 // ----------------------------------------
407 h = fHResolution.ProjectionZ("Resolution");
408 //h->SetXTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(E_{mc})}");
409 h->SetYTitle("Counts");
410 h->SetTitle("Distribution of \\Delta lg(E) / lg(E_{mc})");
411 h->SetDirectory(NULL);
412 h->SetBit(kCanDelete);
413 h->GetXaxis()->SetRangeUser(-1.3, 1.3);
414 h->Draw("");
415 // ----------------------------------------
416
417 pad->cd(2);
418 gPad->SetBorderMode(0);
419
420 gPad->Divide(1, 3, 1e-10, 1e-10);
421 pad2 = gPad;
422
423 pad2->cd(1);
424 h = MakePlot(fHEnergy, "yx");
425 h->SetXTitle("E_{mc} [GeV]");
426 h->SetYTitle("E_{est} [GeV]");
427 h->SetTitle("Estimated energy E_{est} vs Monte Carlo energy E_{mc}");
428
429 TLine line;
430 line.DrawLine(0,0,1,1);
431
432 pad2->cd(2);
433 h = MakePlot(fHResolution, "zy");
434 h->SetXTitle("E_{mc} [GeV]");
435 h->SetYTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(E_{mc})}");
436 h->SetTitle("Energy resolution \\Delta lg(E) / lg(E) vs Monte Carlo energy E_{mc}");
437 h->SetMinimum(-1.3);
438 h->SetMaximum(1.3);
439
440 pad2->cd(3);
441 h = MakePlot(fHResolution, "zx");
442 h->SetXTitle("E_{est} [GeV]");
443 h->SetYTitle("\\frac{lg(E_{est}) - lg(E_{mc})}{lg(E_{mc})}");
444 h->SetTitle("Energy Resolution \\Delta lg(E) / lg(E) vs estimated Energy E_{est}");
445 h->SetMinimum(-1.3);
446 h->SetMaximum(1.3);
447}
448
449// --------------------------------------------------------------------------
450//
451// Returns the mapped value from the Matrix
452//
453Double_t MHEnergyEst::GetVal(Int_t i) const
454{
455 return (*fMatrix)[fMap[i]];
456}
457
458// --------------------------------------------------------------------------
459//
460// You can use this function if you want to use a MHMatrix instead of the
461// given containers. This function adds all necessary columns to the
462// given matrix. Afterward you should fill the matrix with the corresponding
463// data (eg from a file by using MHMatrix::Fill). If you now loop
464// through the matrix (eg using MMatrixLoop) MEnergyEstParam2::Process
465// will take the values from the matrix instead of the containers.
466//
467void MHEnergyEst::InitMapping(MHMatrix *mat)
468{
469 if (fMatrix)
470 return;
471
472 fMatrix = mat;
473
474 fMap[0] = fMatrix->AddColumn("MMcEvt.fEnergy");
475 fMap[1] = fMatrix->AddColumn("MMcEvt.fImpact/100");
476 fMap[2] = fMatrix->AddColumn("MMcEvt.fTelescopeTheta*kRad2Deg");
477}
478
479void MHEnergyEst::StopMapping()
480{
481 fMatrix = NULL;
482}
483
Note: See TracBrowser for help on using the repository browser.