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

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