source: trunk/MagicSoft/Mars/mhist/MHFlux.cc@ 1994

Last change on this file since 1994 was 1992, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 17.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 5/2002 <mailto:wittek@mppmu.mpg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26// //
27// MHFlux //
28// //
29// calculates absolute photon fluxes //
30// from the distributions of the estimated energy //
31// for the different bins in some variable 'Var' //
32// (Var = Theta or time) //
33// //
34//////////////////////////////////////////////////////////////////////////////
35
36#include "MHFlux.h"
37
38#include <TStyle.h>
39
40#include <TF1.h>
41#include <TH2.h>
42#include <TProfile.h>
43
44#include <TCanvas.h>
45
46#include "MTime.h"
47
48#include "MBinning.h"
49#include "MParList.h"
50
51#include "MLog.h"
52#include "MLogManip.h"
53
54#include "MHThetabarTheta.h"
55#include "MHEffOnTime.h"
56#include "MHGamma.h"
57
58ClassImp(MHFlux);
59
60MHFlux::MHFlux(const MHGamma &hist, const TString varname, const TString unit)
61 : fHOrig(), fHUnfold(), fHFlux()
62{
63 const TH2D &h2d = *hist.GetProject();
64
65 if (varname.IsNull() || unit.IsNull())
66 *fLog << warn << dbginf << "varname or unit not defined" << endl;
67
68 fVarname = varname;
69 fUnit = unit;
70
71 // char txt[100];
72
73 // original distribution of E-est for different bins
74 // of the variable (Theta or time)
75 // sprintf(txt, "gammas vs. E-est and %s",varname);
76
77 ((TH2D&)h2d).Copy(fHOrig);
78
79 fHOrig.SetName("E_est");
80 fHOrig.SetTitle(TString("No.of gammas vs. E-est and ")+fVarname);
81
82 fHOrig.SetDirectory(NULL);
83 fHOrig.SetXTitle("E_{est} [GeV]");
84 fHOrig.SetYTitle(fVarname+fUnit);
85 //fHOrig.Sumw2();
86
87 SetBinning((TH2*)&fHOrig, (TH2*)&h2d);
88
89 fHOrig.Copy(fHUnfold);
90
91 // unfolded distribution of E-unfold for different bins
92 // of the variable (Theta or time)
93 // sprintf(txt, "gammas vs. E-unfold and %s",varname);
94 fHUnfold.SetName("E-unfolded");
95 fHUnfold.SetTitle(TString("No.of gammas vs. E-unfold and ")+fVarname);
96
97 fHUnfold.SetDirectory(NULL);
98 fHUnfold.SetXTitle("E_{unfold} [GeV]");
99 fHUnfold.SetYTitle(fVarname+fUnit);
100 //fHUnfold.Sumw2();
101
102 SetBinning((TH2*)&fHUnfold, (TH2*)&fHOrig);
103
104
105 // absolute photon flux vs. E-unfold
106 // for different bins of the variable (Theta or time)
107 //
108 // sprintf(txt, "gamma flux [1/(s m2 GeV) vs. E-unfold and %s",varname);
109 fHFlux.SetName("photon flux");
110 fHFlux.SetTitle(TString("Gamma flux [1/(s m^2 GeV) vs. E-unfold and ")+fVarname);
111
112 fHFlux.SetDirectory(NULL);
113 fHFlux.SetXTitle("E_{unfold} [GeV]");
114 fHFlux.SetYTitle(fVarname+fUnit);
115 fHFlux.Sumw2();
116
117 SetBinning((TH2*)&fHFlux, (TH2*)&fHUnfold);
118}
119
120// --------------------------------------------------------------------------
121//
122// Default Constructor. It sets the variable name (Theta or time)
123// and the units for the variable
124//
125MHFlux::MHFlux(const TH2D &h2d, const TString varname, const TString unit)
126 : fHOrig(), fHUnfold(), fHFlux()
127{
128 if (varname.IsNull() || unit.IsNull())
129 *fLog << warn << dbginf << "varname or unit not defined" << endl;
130
131 fVarname = varname;
132 fUnit = unit;
133
134 // char txt[100];
135
136 // original distribution of E-est for different bins
137 // of the variable (Theta or time)
138 // sprintf(txt, "gammas vs. E-est and %s",varname);
139
140 ((TH2D&)h2d).Copy(fHOrig);
141
142 fHOrig.SetName("E_est");
143 fHOrig.SetTitle(TString("No.of gammas vs. E-est and ")+fVarname);
144
145 fHOrig.SetDirectory(NULL);
146 fHOrig.SetXTitle("E_{est} [GeV]");
147 fHOrig.SetYTitle(fVarname+fUnit);
148 //fHOrig.Sumw2();
149
150 // copy fHOrig into fHUnfold in case no unfolding is done
151 fHOrig.Copy(fHUnfold);
152
153 SetBinning((TH2*)&fHOrig, (TH2*)&h2d);
154
155
156 // unfolded distribution of E-unfold for different bins
157 // of the variable (Theta or time)
158 // sprintf(txt, "gammas vs. E-unfold and %s",varname);
159 fHUnfold.SetName("E-unfolded");
160 fHUnfold.SetTitle(TString("No.of gammas vs. E-unfold and ")+fVarname);
161
162 fHUnfold.SetDirectory(NULL);
163 fHUnfold.SetXTitle("E_{unfold} [GeV]");
164 fHUnfold.SetYTitle(fVarname+fUnit);
165 //fHUnfold.Sumw2();
166
167 SetBinning((TH2*)&fHUnfold, (TH2*)&fHOrig);
168
169
170 // absolute photon flux vs. E-unfold
171 // for different bins of the variable (Theta or time)
172 //
173 // sprintf(txt, "gamma flux [1/(s m2 GeV) vs. E-unfold and %s",varname);
174 fHFlux.SetName("photon flux");
175 fHFlux.SetTitle(TString("Gamma flux [1/(s m^{2} GeV)] vs. E-unfold and ")+fVarname);
176
177 fHFlux.SetDirectory(NULL);
178 fHFlux.SetXTitle("E_{unfold} [GeV]");
179 fHFlux.SetYTitle(fVarname+fUnit);
180 fHFlux.Sumw2();
181
182 SetBinning((TH2*)&fHFlux, (TH2*)&fHUnfold);
183}
184
185// -------------------------------------------------------------------------
186//
187// Unfold the distribution in E-est
188//
189void MHFlux::Unfold()
190{
191}
192
193void MHFlux::CalcFlux(const MHEffOnTime &teff, const MHThetabarTheta &thetabar,
194 const TH2D *aeff)
195{
196 CalcFlux(teff.GetHist(), thetabar.GetHist(), aeff);
197}
198
199Double_t MHFlux::ParabInterpolLog(const TAxis &axe, Int_t j,
200 Double_t y[], Double_t Ebar) const
201{
202 const double t1 = log10(axe.GetBinLowEdge(j-1)) + log10(axe.GetBinUpEdge(j-1));
203 const double t2 = log10(axe.GetBinLowEdge(j)) + log10(axe.GetBinUpEdge(j));
204 const double t3 = log10(axe.GetBinLowEdge(j+1)) + log10(axe.GetBinUpEdge(j+1));
205
206 const Double_t lebar = log10(Ebar);
207
208 return Parab(t1/2, t2/2, t3/2, y[j-2], y[j-1], y[j], lebar);
209}
210
211// --------------------------------------------------------------------
212//
213// determine bins for interpolation (k3 is the middle one) in bar.
214// k0 denotes the bin from which the error is copied
215//
216void MHFlux::FindBins(const TAxis &axe, const Double_t bar, Int_t &k3, Int_t &k0) const
217{
218 const Int_t n = axe.GetNbins();
219
220 k3 = axe.FindFixBin(bar);
221 k0 = k3;
222
223 if (k3<2)
224 {
225 k3 = 2;
226 if (bar<axe.GetBinLowEdge(2))
227 k0 = 1;
228 }
229
230 if (k3>n-1)
231 {
232 k3 = n-1;
233 if (bar>axe.GetBinLowEdge(n))
234 k0 = n;
235 }
236
237 if (bar>=axe.GetBinLowEdge(1) && bar<=axe.GetBinUpEdge(n))
238 return;
239
240 *fLog << dbginf << "extrapolation: bar = " << bar;
241 *fLog << ", min =" << axe.GetBinLowEdge(1);
242 *fLog << ", max =" << axe.GetBinUpEdge(n) << endl;
243}
244
245Double_t MHFlux::ParabInterpolCos(const TAxis &axe, const TH2D *aeff, Int_t j, Int_t k3, Double_t val) const
246{
247 const double t1 = cos( axe.GetBinCenter (k3-1) );
248 const double t2 = cos( axe.GetBinCenter (k3) );
249 const double t3 = cos( axe.GetBinCenter (k3+1) );
250
251 const double a1 = aeff->GetBinContent(j, k3-1);
252 const double a2 = aeff->GetBinContent(j, k3);
253 const double a3 = aeff->GetBinContent(j, k3+1);
254
255 return Parab(t1, t2, t3, a1, a2, a3, val);
256}
257
258// -------------------------------------------------------------------------
259//
260// Calculate photon flux by dividing the distribution in Eunf (fHUnfold) by
261// the width of the energy interval (deltaE)
262// the effective ontime (*teff)
263// and the effective collection area (*aeff)
264//
265void MHFlux::CalcFlux(const TH1D *teff, const TProfile *thetabar,
266 const TH2D *aeff)
267{
268 //
269 // Note that fHUnfold has bins in Eunf and Var
270 // *teff has bins in Var (the same bins in Var as fHUnfold)
271 // *thetabar has bins in Var (the same bins in Var as fHUnfold)
272 // *aeff has bins in Etru and Theta
273 // (where in general the binning in Etru is different
274 // from the binning in Eunf)
275 // The variable Var may be 'time' or 'Theta'
276
277 const TAxis &axex = *((TH2*)aeff)->GetXaxis();
278 const TAxis &axey = *((TH2*)aeff)->GetYaxis();
279
280 if (axex.GetNbins()<3)
281 {
282 *fLog << err << "ERROR - Number of Energy bins <3 not implemented!" << endl;
283 return;
284 }
285
286 if (axey.GetNbins()<3)
287 *fLog << warn << "WARNING - Less than 3 theta-bins not supported very well!" << endl;
288
289 //
290 // calculate effective collection area
291 // for the Eunf and Var bins of the histogram fHUnfold
292 // from the histogram *aeff, which has bins in Etru and Theta
293 // the result is the histogram fHAeff
294 //
295 TH2D fHAeff;
296 SetBinning((TH2*)&fHAeff, (TH2*)&fHUnfold);
297 fHAeff.Sumw2();
298
299 //
300 // ------ start loops ------
301 //
302 const Int_t nEtru = aeff->GetNbinsX();
303
304 Double_t *aeffbar = new Double_t[nEtru];
305 Double_t *daeffbar = new Double_t[nEtru];
306
307 const Int_t nVar = fHFlux.GetNbinsY();
308 for (int n=1; n<=nVar; n++)
309 {
310 const Double_t tbar = thetabar->GetBinContent(n);
311 const Double_t costbar = cos(tbar/kRad2Deg);
312
313 // determine bins for interpolation (k3, k0)
314 Int_t kv, ke;
315 FindBins(axey, tbar, kv, ke);
316
317 //
318 // calculate effective collection area at Theta = Thetabar
319 // by quadratic interpolation in cos(Theta);
320 // do this for each bin of Etru
321 //
322 for (int j=1; j<=nEtru; j++)
323 {
324 if (axey.GetNbins()<3)
325 {
326 // FIXME: Other interpolation?
327 aeffbar[j-1] = aeff->GetBinContent(j, n);
328 daeffbar[j-1] = aeff->GetBinError(j, n);
329 }
330 else
331 {
332 aeffbar[j-1] = ParabInterpolCos(axey, aeff, j, kv, costbar);
333 daeffbar[j-1] = aeff->GetBinError(j, ke);
334 }
335 }
336
337 //
338 // calculate effective collection area at (E = Ebar, Theta = Thetabar)
339 // by quadratic interpolation in log10(Etru)
340 // do this for each bin of Eunf
341 //
342 CalcEffCol(axex, fHAeff, n, aeffbar, daeffbar);
343 }
344
345 delete aeffbar;
346 delete daeffbar;
347
348 //
349 // now calculate the absolute gamma flux
350 //
351 CalcAbsGammaFlux(*teff, fHAeff);
352}
353
354// --------------------------------------------------------------------
355//
356// calculate effective collection area at (E = Ebar, Theta = Thetabar)
357// by quadratic interpolation in log10(Etru)
358// do this for each bin of Eunf
359//
360void MHFlux::CalcEffCol(const TAxis &axex, TH2D &fHAeff, Int_t n, Double_t aeffbar[], Double_t daeffbar[])
361{
362 const Int_t nEunf = fHFlux.GetNbinsX();
363
364 const TAxis &unfx = *fHUnfold.GetXaxis();
365
366 for (int m=1; m<=nEunf; m++)
367 {
368 const Double_t Ebar = GetBinCenterLog(unfx, m);
369
370 Int_t j0, j3;
371 FindBins(axex, Ebar, j3, j0);
372
373 const Double_t v = ParabInterpolLog(axex, j3, aeffbar, Ebar);
374
375 fHAeff.SetBinContent(m,n, v);
376 fHAeff.SetBinError(m,n, daeffbar[j0-1]);
377 }
378}
379
380// --------------------------------------------------------------------
381//
382// calculate the absolute gamma flux
383//
384void MHFlux::CalcAbsGammaFlux(const TH1D &teff, const TH2D &fHAeff)
385{
386 const Int_t nEunf = fHFlux.GetNbinsX();
387 const Int_t nVar = fHFlux.GetNbinsY();
388
389 for (int m=1; m<=nEunf; m++)
390 {
391 const Double_t DeltaE = fHFlux.GetXaxis()->GetBinWidth(m);
392
393 for (int n=1; n<=nVar; n++)
394 {
395 const Double_t Ngam = fHUnfold.GetBinContent(m,n);
396 const Double_t Aeff = fHAeff.GetBinContent(m,n);
397 const Double_t Effon = teff.GetBinContent(n);
398
399 const Double_t c1 = fHUnfold.GetBinError(m,n)/Ngam;
400 const Double_t c2 = teff.GetBinError(n) /Effon;
401 const Double_t c3 = fHAeff.GetBinError(m,n) /Aeff;
402
403 const Double_t cont = Ngam / (DeltaE * Effon * Aeff);
404 const Double_t dcont = sqrt(c1*c1 + c2*c2 + c3*c3);
405
406 //
407 // Out of Range
408 //
409 const Bool_t oor = Ngam<=0 || DeltaE<=0 || Effon<=0 || Aeff<=0;
410
411 if (oor)
412 *fLog << warn << "MHFlux::CalcAbsGammaFlux(" << m << "," << n << ") ";
413
414 if (Ngam<=0)
415 *fLog << " Ngam=0";
416 if (DeltaE<=0)
417 *fLog << " DeltaE=0";
418 if (Effon<=0)
419 *fLog << " Effon=0";
420 if (Aeff<=0)
421 *fLog << " Aeff=0";
422
423 if (oor)
424 *fLog << endl;
425
426 fHFlux.SetBinContent(m,n, oor ? 1e-20 : cont);
427 fHFlux.SetBinError(m,n, oor ? 1e-20 : dcont*cont);
428 }
429 }
430}
431
432// --------------------------------------------------------------------
433//
434// draw the differential photon flux vs. E-unf
435// for the individual bins of the variable Var
436//
437void MHFlux::DrawFluxProjectionX(Option_t *opt) const
438{
439 const Int_t nVar = fHFlux.GetNbinsY();
440
441 for (int n=1; n<=nVar; n++)
442 {
443 TString strg0("Flux-");
444
445 TH1D &h = *((TH2D)fHFlux).ProjectionX(strg0+fVarname, n, n, "E");
446
447 TString strg1 = "Photon flux vs. E_{unfold} for ";
448 TString strg2 = fVarname+"-bin #";
449 strg2 += n;
450
451 new TCanvas(strg2, strg1+strg2);
452 gPad->SetLogx();
453 gPad->SetLogy();
454
455 TString name = fVarname+"bin_";
456 name += n;
457
458 h.SetName(name);
459 h.SetTitle(strg1+strg2);
460 h.SetXTitle("E_{unfold} [GeV]");
461 h.SetYTitle("photons / (s m^{2} GeV)");
462 h.GetXaxis()->SetLabelOffset(-0.025);
463 h.GetXaxis()->SetTitleOffset(1.1);
464 h.GetXaxis()->SetNdivisions(1);
465 h.GetYaxis()->SetTitleOffset(1.25);
466 h.DrawCopy();
467 }
468}
469
470void MHFlux::DrawOrigProjectionX(Option_t *opt) const
471{
472 const Int_t nVar = fHOrig.GetNbinsY();
473
474 for (int n=1; n<=nVar; n++)
475 {
476 TString strg0 = "Orig-";
477 strg0 += fVarname;
478 strg0 += "_";
479 strg0 += n;
480
481 TH1D &h = *((TH2D)fHOrig).ProjectionX(strg0, n, n, "E");
482
483 TString strg1("No.of photons vs. E-est for ");
484 strg1 += fVarname+"-bin ";
485 strg1 += n;
486
487 new TCanvas(strg0, strg1);
488
489 gPad->SetLogx();
490 gPad->SetLogy();
491
492 h.SetName(strg0);
493 h.SetTitle(strg1);
494 h.SetXTitle("E_{est} [GeV]");
495 h.GetXaxis()->SetLabelOffset(-0.025);
496 h.GetXaxis()->SetTitleOffset(1.1);
497 h.GetXaxis()->SetNdivisions(1);
498 h.SetYTitle("No.of photons");
499 h.DrawCopy();
500 }
501}
502
503// -------------------------------------------------------------------------
504//
505// Draw copies of the histograms
506//
507TObject *MHFlux::DrawClone(Option_t *opt) const
508{
509 TCanvas &c = *MakeDefCanvas("flux", "Orig - Unfold - Flux plots");
510 c.Divide(2, 2);
511
512 gROOT->SetSelectedPad(NULL);
513
514 c.cd(1);
515 ((TH2*)&fHOrig)->DrawCopy("");
516 gPad->SetLogx();
517
518 c.cd(2);
519 ((TH2*)&fHUnfold)->DrawCopy("");
520 gPad->SetLogx();
521
522 c.cd(3);
523 ((TH2*)&fHFlux)->DrawCopy("");
524 gPad->SetLogx();
525
526 c.Modified();
527 c.Update();
528
529 return &c;
530}
531
532// -------------------------------------------------------------------------
533//
534// Draw the histograms
535//
536void MHFlux::Draw(Option_t *opt)
537{
538 if (!gPad)
539 MakeDefCanvas("flux", "orig-unfold-flux plots");
540
541 gPad->Divide(2,2);
542
543 gPad->cd(1);
544 fHOrig.Draw(opt);
545
546 gPad->cd(2);
547 fHUnfold.Draw(opt);
548
549 gPad->cd(3);
550 fHFlux.Draw(opt);
551
552 gPad->Modified();
553 gPad->Update();
554}
555
556Double_t MHFlux::Parab(Double_t x1, Double_t x2, Double_t x3,
557 Double_t y1, Double_t y2, Double_t y3,
558 Double_t val)
559{
560 Double_t c0, c1, c2;
561 Parab(x1, x2, x3, y1, y2, y3, &c0, &c1, &c2);
562 return c0 + c1*val + c2*val*val;
563}
564
565// -------------------------------------------------------------------------
566//
567// Quadratic interpolation
568//
569// *** calculate the parameters of a parabula
570// y = a + b*x + c*x**2 = F(x)
571// such that yi = F(xi) for (i=1,3)
572//
573Bool_t MHFlux::Parab(Double_t x1, Double_t x2, Double_t x3,
574 Double_t y1, Double_t y2, Double_t y3,
575 Double_t *a, Double_t *b, Double_t *c)
576{
577 const double det =
578 + x2*x3*x3 + x1*x2*x2 + x3*x1*x1
579 - x2*x1*x1 - x3*x2*x2 - x1*x3*x3;
580
581 if (det==0)
582 {
583 *a = 0;
584 *b = 0;
585 *c = 0;
586 return kFALSE;
587 }
588
589 const double det1 = 1.0/det;
590
591 const double ai11 = x2*x3*x3 - x3*x2*x2;
592 const double ai12 = x3*x1*x1 - x1*x3*x3;
593 const double ai13 = x1*x2*x2 - x2*x1*x1;
594
595 const double ai21 = x2*x2 - x3*x3;
596 const double ai22 = x3*x3 - x1*x1;
597 const double ai23 = x1*x1 - x2*x2;
598
599 const double ai31 = x3 - x2;
600 const double ai32 = x1 - x3;
601 const double ai33 = x2 - x1;
602
603 *a = (ai11*y1 + ai12*y2 + ai13*y3) * det1;
604 *b = (ai21*y1 + ai22*y2 + ai23*y3) * det1;
605 *c = (ai31*y1 + ai32*y2 + ai33*y3) * det1;
606
607 return kTRUE;
608}
Note: See TracBrowser for help on using the repository browser.