source: trunk/MagicSoft/Mars/mhflux/MMcSpectrumWeight.cc@ 8785

Last change on this file since 8785 was 8746, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 19.0 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): Thomas Bretz 5/2005 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Marcos Lopez 10/2003 <mailto:marcos@gae.ucm.es>
20!
21! Copyright: MAGIC Software Development, 2000-2007
22!
23!
24\* ======================================================================== */
25
26//////////////////////////////////////////////////////////////////////////////
27//
28// MMcSpectrumWeight
29//
30// Change the spectrum of the MC showers simulated with Corsika (a power law)
31// to a new one, which can be either, again a power law but with a different
32// spectral index, or a generalizeed spectrum. The new spectrum can be
33// pass to this class in different ways:
34//
35// 1. If the new spectrum will be a power law, just introduce the slope
36// of this power law.
37// 2. If the new spectrum will have a general shape:
38// The new spectrum is passed as a char* (SetFormula())
39//
40// Method:
41// -------
42//
43// - Corsika spectrun: dN/dE = A * E^(a)
44// with a = fOldSlope, and A = N/integral{E*de} from ELowLim to EUppLim
45//
46// - New spectrum: dN/dE = B * g(E)
47// where B = N/integral{g*dE} from ELowLim to EUppLim, and N=NumEvents
48//
49// For converting the spectrum simulated with Corsika to the new one, we
50// apply a weight to each event, given by:
51//
52// W(E) = B/A * g(E)/E^(a)
53//
54// In the case the new spectrum is simply a power law: dN/dE = B * E^(b), we
55// have:
56//
57// W(E) = B/A * E^(b-a)
58//
59// (The factor B/A is used in order both the original and new spectrum have
60// the same area (i.e. in order they represent the same number of showers))
61//
62//
63// If using SetFormula you can specify formulas accepted by TF1, eg:
64// pow(X, -2.6)
65// (Rem: all capital (!) 'X' are replaced by the corresponding %s.fEnergy
66// automatically)
67//
68// For more details of the setup see MMcSpectrumWeight::ReadEnv
69//
70//
71// Input Containers:
72// MMcEvt
73// MMcCorsikaRunHeader
74// [MPointingPos]
75// [MHillas]
76//
77// Output Container:
78// MWeight [MParameterD]
79//
80//////////////////////////////////////////////////////////////////////////////
81#include "MMcSpectrumWeight.h"
82
83#include <TF1.h>
84#include <TH1.h>
85#include <TH2.h>
86#include <TH3.h>
87#include <TSpline.h>
88
89#include "MLog.h"
90#include "MLogManip.h"
91
92#include "MParList.h"
93#include "MParameters.h"
94
95#include "MHillas.h"
96#include "MPointingPos.h"
97
98#include "MMcEvt.hxx"
99#include "MMcCorsikaRunHeader.h"
100
101ClassImp(MMcSpectrumWeight);
102
103using namespace std;
104
105void MMcSpectrumWeight::Init(const char *name, const char *title)
106{
107 fName = name ? name : "MMcSpectrumWeight";
108 fTitle = title ? title : "Task to calculate weights to change the energy spectrum";
109
110 AddToBranchList("MMcEvt.fEnergy");
111
112 fNameWeight = "MWeight";
113 fNameMcEvt = "MMcEvt";
114
115 fNewSlope = -99;
116 fOldSlope = -99;
117
118 fEnergyMin = -1;
119 fEnergyMax = -2;
120
121 fNorm = 1;
122 fNormEnergy = -1;
123
124 fAllowChange = kFALSE;
125
126 fFunc = NULL;
127 fMcEvt = NULL;
128 fHillas = NULL;
129 fWeight = NULL;
130 fWeightsZd = NULL;
131 fWeightsSize = NULL;
132 fPointing = NULL;
133}
134
135// ---------------------------------------------------------------------------
136//
137// Default Constructor.
138//
139MMcSpectrumWeight::MMcSpectrumWeight(const char *name, const char *title)
140{
141 Init(name,title);
142}
143
144// ---------------------------------------------------------------------------
145//
146// Destructor. If necessary delete fFunc
147//
148MMcSpectrumWeight::~MMcSpectrumWeight()
149{
150 if (fFunc)
151 delete fFunc;
152// if (fWeightsSize)
153// delete fWeightsSize;
154}
155
156// ---------------------------------------------------------------------------
157//
158// Search for
159// - fNameMcEvt [MMcEvtBasic]
160//
161// Find/Create:
162// - fNameWeight [MWeight]
163//
164Int_t MMcSpectrumWeight::PreProcess(MParList *pList)
165{
166 fMcEvt = (MMcEvt*)pList->FindObject(fNameMcEvt, "MMcEvtBasic");
167 if (!fMcEvt)
168 {
169 *fLog << err << fNameMcEvt << " [MMcEvtBasic] not found... abort." << endl;
170 return kFALSE;
171 }
172
173 fWeight = (MParameterD*)pList->FindCreateObj("MParameterD", fNameWeight);
174 if (!fWeight)
175 return kFALSE;
176
177 if (fWeightsZd)
178 {
179 fPointing = (MPointingPos*)pList->FindObject("MPointingPos");
180 if (!fPointing)
181 {
182 *fLog << err << "MPointingPos not found... abort." << endl;
183 return kFALSE;
184 }
185 }
186
187 if (fWeightsSize)
188 {
189 fHillas = (MHillas*)pList->FindObject("MHillas");
190 if (!fHillas)
191 {
192 *fLog << err << "MHillas not found... abort." << endl;
193 return kFALSE;
194 }
195 }
196
197 return kTRUE;
198}
199
200// ---------------------------------------------------------------------------
201//
202// Replace {fNameMcEvt}.fEnergy by "(x)" and return the result.
203//
204TString MMcSpectrumWeight::ReplaceX(TString str) const
205{
206 return str.ReplaceAll(Form("%s.fEnergy", fNameMcEvt.Data()), "(x)");
207}
208
209// ---------------------------------------------------------------------------
210//
211// Return the function corresponding to the mc spectrum with
212// slope fOldSlope: pow({fNameMcEvt}.fEnergy, fOldSlope)
213//
214// The slope is returned as %.3f
215//
216TString MMcSpectrumWeight::GetFormulaSpecOld(const char *name) const
217{
218 return Form("pow(%s.fEnergy, %.3f)", name, fOldSlope);
219}
220
221// ---------------------------------------------------------------------------
222//
223// Return the function corresponding to the new spectrum with
224// slope fNewSlope: pow({fNameMcEvt}.fEnergy, fNewSlope)
225//
226// The slope is returned as %.3f
227//
228// If a different formula is set (SetFormula()) this formula is returned
229// unchanged.
230//
231TString MMcSpectrumWeight::GetFormulaSpecNew(const char *name) const
232{
233 TString str = fFormula.IsNull() ? Form("pow(%s.fEnergy, %.3f)", name, fNewSlope) : fFormula.Data();
234 if (!fFormula.IsNull())
235 str.ReplaceAll("X", Form("(%s.fEnergy)", name));
236
237 return str;
238}
239
240// ---------------------------------------------------------------------------
241//
242// Return the formula to calculate weights.
243// Is is compiled by
244// o1 = integral(fEnergyMin, fEnergyMax, GetFormulaSpecOldX());
245// n1 = integral(fEnergyMin, fEnergyMax, GetFormulaSpecNewX());
246// o2 = CalcSpecOld(fNormEnergy);
247// n2 = CalcSpecNew(fNormEnergy);
248//
249// result (fNormEnergy<0):
250// fNorm*o1/n1*GetFormulaNewSpec()/GetFormulaOldSpec()
251//
252// result (fNormEnergy>=0):
253// fNorm*o2/n2*GetFormulaNewSpec()/GetFormulaOldSpec()
254//
255// fNorm is 1 by default but can be overwritten using SetNorm()
256//
257// If the formulas GetFormulaSpecOldX() and GetFormulaSpecNewX()
258// are equal only fNorm is returned.
259//
260// The normalization constant is returned as %.16e
261//
262// Example: 0.3712780019*(pow(MMcEvt.fEnergy,-2.270))/(pow(MMcEvt.fEnergy,-2.600))
263//
264TString MMcSpectrumWeight::GetFormulaWeights(const char *name) const
265{
266 if (GetFormulaSpecOld()==GetFormulaSpecNew())
267 return Form("%.16e", fNorm);
268
269 const Double_t iold = fNormEnergy<0 ? GetSpecOldIntegral() : CalcSpecOld(fNormEnergy);
270 const Double_t inew = fNormEnergy<0 ? GetSpecNewIntegral() : CalcSpecNew(fNormEnergy);
271
272 const Double_t norm = fNorm*iold/inew;
273
274 return Form("%.16e*(%s)/(%s)", norm, GetFormulaSpecNew(name).Data(), GetFormulaSpecOld(name).Data());
275}
276
277// ---------------------------------------------------------------------------
278//
279// Returns the integral between fEnergyMin and fEnergyMax of
280// GetFormulaSpecNewX() describing the destination spectrum
281//
282Double_t MMcSpectrumWeight::GetSpecNewIntegral(Double_t emin, Double_t emax) const
283{
284 TF1 funcnew("Dummy", GetFormulaSpecNewX());
285 return funcnew.Integral(emin, emax);
286}
287
288// ---------------------------------------------------------------------------
289//
290// Returns the integral between fEnergyMin and fEnergyMax of
291// GetFormulaSpecOldX() describing the simulated spectrum
292//
293Double_t MMcSpectrumWeight::GetSpecOldIntegral(Double_t emin, Double_t emax) const
294{
295 TF1 funcold("Dummy", GetFormulaSpecOldX());
296 return funcold.Integral(emin, emax);
297}
298
299// ---------------------------------------------------------------------------
300//
301// Returns the value of GetFormulaSpecNewX() at the energy e describing
302// the destination spectrum
303//
304Double_t MMcSpectrumWeight::CalcSpecNew(Double_t e) const
305{
306 TF1 funcnew("Dummy", GetFormulaSpecNewX());
307 return funcnew.Eval(e);
308}
309
310// ---------------------------------------------------------------------------
311//
312// Returns the value of GetFormulaSpecOldX() at the energy e describing
313// the simulated spectrum
314//
315Double_t MMcSpectrumWeight::CalcSpecOld(Double_t e) const
316{
317 TF1 funcnew("Dummy", GetFormulaSpecOldX());
318 return funcnew.Eval(e);
319}
320
321void MMcSpectrumWeight::SetWeightsSize(TH1D *h)
322{
323 fWeightsSize=h;
324 /*
325 if (h==0)
326 {
327 fWeightsSize=0;
328 return;
329 }
330
331 if (fWeightsSize)
332 delete fWeightsSize;
333
334 const Double_t xmin = TMath::Log10(h->GetXaxis()->GetXmin());
335 const Double_t xmax = TMath::Log10(h->GetXaxis()->GetXmax());
336 const Double_t xnum = h->GetNbinsX()+1;
337
338 fWeightsSize = new TSpline3("WeightsSize", xmin, xmax,
339 h->GetArray()+1, xnum);*/
340}
341
342// ---------------------------------------------------------------------------
343//
344// Initialize fEnergyMin, fEnergymax and fOldSlope from MMcCorsikaRunHeader
345// by GetELowLim(), GetEUppLim() and GetSlopeSpec().
346//
347// If fEnergyMax>fEnergyMin (means: the values have already been
348// initialized) and !fAllowChange the consistency of the new values
349// with the present values is checked with a numerical precision of 1e-10.
350// If one doesn't match kFALSE is returned.
351//
352// If the mc slope is -1 kFALSE is returned.
353//
354// If the new slope for the spectrum is -1 it is set to the original MC
355// slope.
356//
357// fFunc is set to the formula returned by GetFormulaWeightsX()
358//
359Bool_t MMcSpectrumWeight::Set(const MMcCorsikaRunHeader &rh)
360{
361 if (fEnergyMax>fEnergyMin && !fAllowChange)
362 {
363 if (TMath::Abs(fEnergyMax-rh.GetEUppLim())>1e-10)
364 {
365 *fLog << err;
366 *fLog << "ERROR - The maximum simulated Monte Carlo energy is not allowed to change ";
367 *fLog << "(" << fEnergyMax << " --> " << rh.GetEUppLim() << ")... abort." << endl;
368 return kFALSE;
369 }
370
371 if (TMath::Abs(fOldSlope-rh.GetSlopeSpec())>1e-10)
372 {
373 *fLog << err;
374 *fLog << "ERROR - The slope of the Monte Carlo is not allowed to change ";
375 *fLog << "(" << fOldSlope << " --> " << rh.GetSlopeSpec() << ")... abort." << endl;
376 return kFALSE;
377 }
378
379 // No change happened
380 if (TMath::Abs(fEnergyMin-rh.GetELowLim())<=1e-10)
381 return kTRUE;
382
383 // The lower energy limit has changed
384 *fLog << warn;
385 *fLog << "The minimum simulated Monte Carlo energy has changed from ";
386 *fLog << fEnergyMin << "GeV to " << rh.GetELowLim() << "GeV." << endl;
387 fEnergyMin = rh.GetELowLim();
388 }
389
390 fOldSlope = rh.GetSlopeSpec();
391 fEnergyMin = rh.GetELowLim();
392 fEnergyMax = rh.GetEUppLim();
393
394 if (fNewSlope==-99 && fFormula.IsNull())
395 {
396 *fLog << inf << "A new slope for the power law has not yet been defined... using " << fOldSlope << "." << endl;
397 fNewSlope = fOldSlope;
398 }
399
400 if (fFunc)
401 delete fFunc;
402
403 if (GetFormulaSpecOld()==GetFormulaSpecNew())
404 *fLog << inf << "No spectral change requested..." << endl;
405 else
406 {
407 *fLog << inf << "Weighting from slope " << fOldSlope << " to ";
408 if (fFormula.IsNull())
409 *fLog << "slope " << fNewSlope << "." << endl;
410 else
411 *fLog << GetFormulaSpecNewX() << endl;
412 }
413
414 fFunc = new TF1("", GetFormulaWeightsX());
415 gROOT->GetListOfFunctions()->Remove(fFunc);
416 fFunc->SetName("SpectralWeighs");
417
418 return kTRUE;
419}
420
421// ---------------------------------------------------------------------------
422//
423// completes a simulated spectrum starting at an energy fEnergyMin down to
424// an energy emin.
425//
426// It is assumed that the contents of MMcSpectrumWeight for the new spectrum
427// correctly describe the spectrum within the histogram, and fEnergyMin
428// and fEnergyMax correctly describe the range.
429//
430// In the 1D case it is assumed that the x-axis is a zenith angle binning.
431// In the 2D case the x-axis is assumed to be zenith angle, the y-axis
432// to be energy.
433//
434void MMcSpectrumWeight::CompleteEnergySpectrum(TH1 &h, Double_t emin) const
435{
436 if (h.InheritsFrom(TH3::Class()))
437 {
438 *fLog << "ERROR - MMcSpctrumWeight::CompleteEnergySpectrum doesn't support TH3." << endl;
439 return;
440 }
441
442 if (fEnergyMin <= emin)
443 return;
444
445 const Double_t norm = GetSpecNewIntegral();
446
447 if (!h.InheritsFrom(TH2::Class()))
448 {
449 h.Scale(GetSpecNewIntegral(emin, fEnergyMax)/norm);
450 return;
451 }
452
453 const TAxis &axey = *h.GetYaxis();
454
455 const Int_t first = axey.FindFixBin(emin);
456 const Int_t last = axey.FindFixBin(fEnergyMin);
457
458 for (int x=1; x<=h.GetNbinsX(); x++)
459 {
460 const Double_t f = h.Integral(x, x, -1, 9999)/norm;
461
462 for (int y=first; y<=last; y++)
463 {
464 const Double_t lo = axey.GetBinLowEdge(y) <emin ? emin : axey.GetBinLowEdge(y);
465 const Double_t hi = axey.GetBinLowEdge(y+1)>fEnergyMin ? fEnergyMin : axey.GetBinLowEdge(y+1);
466
467 h.AddBinContent(h.GetBin(x, y), f*GetSpecNewIntegral(lo, hi));
468 }
469 }
470}
471
472// ---------------------------------------------------------------------------
473//
474// The current contants are printed
475//
476void MMcSpectrumWeight::Print(Option_t *o) const
477{
478 const TString opt(o);
479
480 const Bool_t hasnew = opt.Contains("new") || opt.IsNull();
481 const Bool_t hasold = opt.Contains("old") || opt.IsNull();
482
483 *fLog << all << GetDescriptor() << endl;
484
485 if (hasold)
486 {
487 *fLog << " Simulated energy range: " << fEnergyMin << "GeV - " << fEnergyMax << "GeV" << endl;
488 *fLog << " Simulated spectral slope: ";
489 if (fOldSlope==-99)
490 *fLog << "undefined" << endl;
491 else
492 *fLog << fOldSlope << endl;
493 }
494 if (hasnew)
495 {
496 *fLog << " New spectral slope: ";
497 if (fNewSlope==-99)
498 *fLog << "undefined/no change" << endl;
499 else
500 *fLog << fNewSlope << endl;
501 }
502 *fLog << " Additional user norm.: " << fNorm << endl;
503 *fLog << " Spectra are normalized: " << (fNormEnergy<0?"by integral":Form("at %.1fGeV", fNormEnergy)) << endl;
504 if (hasold)
505 {
506 *fLog << " Old Spectrum: ";
507 if (fNewSlope==-99)
508 *fLog << "undefined";
509 else
510 *fLog << GetFormulaSpecOldX();
511 if (fEnergyMin>=0 && fEnergyMax>0)
512 *fLog << " (I=" << GetSpecOldIntegral() << ")";
513 *fLog << endl;
514 }
515 if (hasnew)
516 {
517 *fLog << " New Spectrum: ";
518 if (fNewSlope==-99 && fFormula.IsNull())
519 *fLog << "undefined/no change";
520 else
521 *fLog << GetFormulaSpecNewX();
522 if (fEnergyMin>=0 && fEnergyMax>0)
523 *fLog << " (I=" << GetSpecNewIntegral() << ")";
524 *fLog << endl;
525 }
526 if (fFunc)
527 *fLog << " Weight func: " << fFunc->GetTitle() << endl;
528}
529
530// ----------------------------------------------------------------------------
531//
532// Executed each time a new root file is loaded
533// We will need fOldSlope and fE{Upp,Low}Lim to calculate the weights
534//
535Bool_t MMcSpectrumWeight::ReInit(MParList *plist)
536{
537 MMcCorsikaRunHeader *rh = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
538 if (!rh)
539 {
540 *fLog << err << "MMcCorsikaRunHeader not found... abort." << endl;
541 return kFALSE;
542 }
543
544 return Set(*rh);
545}
546
547/*
548 * This could be used to improve the Zd-weighting within a bin.
549 * Another option is to use more bins, or both.
550 * Note that it seems unnecessary, because the shape within the
551 * theta-bins should be similar in data and Monte Carlo... hopefully.
552 *
553void MMcSpectrumWeight::InitZdWeights()
554{
555 TH2D w(*fWeightsZd);
556
557 for (int i=1; i<=w.GetNbinsX(); i++)
558 {
559 const Double_t tmin = w.GetBinLowEdge(i) *TMath::DegToRad();
560 const Double_t tmax = w.GetBinLowEdge(i+1)*TMath::DegToRad();
561
562 const Double_t wdth = tmax-tmin;
563 const Double_t integ = cos(tmin)-cos(tmax);
564
565 w.SetBinContent(i, w.GetBinContent(i)*wdth/integ);
566 }
567
568 // const Int_t i = fWeightsZd->GetXaxis()->FindFixBin(fPointing->GetZd());
569 // const Double_t theta = fPointing->GetZd()*TMath::DegToRad();
570 // w = sin(theta)*w.GetBinContent(i);
571}
572*/
573
574// ----------------------------------------------------------------------------
575//
576// Fill the result of the evaluation of fFunc at fEvEvt->GetEnergy
577// into the weights container.
578//
579Int_t MMcSpectrumWeight::Process()
580{
581 Double_t w = 1;
582
583 if (fWeightsZd)
584 {
585 const Int_t i = fWeightsZd->GetXaxis()->FindFixBin(fPointing->GetZd());
586 w = fWeightsZd->GetBinContent(i);
587 }
588 if (fWeightsSize)
589 {
590 const Int_t i = fWeightsSize->GetXaxis()->FindFixBin(fHillas->GetSize());
591 w *= fWeightsSize->GetBinContent(i);
592 // w *= fWeightsSize->Eval(TMath::Log10(fHillas->GetSize()));
593 }
594
595 const Double_t e = fMcEvt->GetEnergy();
596
597 fWeight->SetVal(fFunc->Eval(e)*w);
598
599 return kTRUE;
600}
601
602// --------------------------------------------------------------------------
603//
604// Read the setup from a TEnv, eg:
605//
606// MMcSpectrumWeight.NewSlope: -2.6
607// The new slope of the spectrum
608//
609// MMcSpectrumWeight.Norm: 1.0
610// An additional artificial scale factor
611//
612// MMcSpectrumWeight.NormEnergy: 200
613// To normalize at a given energy instead of the integral
614//
615// MMcSpectrumWeight.Formula: pow(X, -2.6)
616// A formula to which the spectrum is weighted (use a capital X for
617// the energy)
618//
619Int_t MMcSpectrumWeight::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
620{
621 Bool_t rc = kFALSE;
622 if (IsEnvDefined(env, prefix, "NewSlope", print))
623 {
624 rc = kTRUE;
625 SetNewSlope(GetEnvValue(env, prefix, "NewSlope", fNewSlope));
626 }
627 if (IsEnvDefined(env, prefix, "Norm", print))
628 {
629 rc = kTRUE;
630 SetNorm(GetEnvValue(env, prefix, "Norm", fNorm));
631 }
632 if (IsEnvDefined(env, prefix, "NormEnergy", print))
633 {
634 rc = kTRUE;
635 SetNormEnergy(GetEnvValue(env, prefix, "NormEnergy", fNormEnergy));
636 }
637 if (IsEnvDefined(env, prefix, "Formula", print))
638 {
639 rc = kTRUE;
640 SetFormula(GetEnvValue(env, prefix, "Formula", fFormula));
641 }
642
643 return rc;
644}
Note: See TracBrowser for help on using the repository browser.