source: trunk/Mars/mhist/MHFlux.cc@ 18158

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