source: tags/Mars-V1.0/mhflux/MMcSpectrumWeight.cc

Last change on this file was 8047, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 15.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): 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-2006
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. Is the new spectrum will be a power law, just introduce the slope
36// of this power law.
37// 2. Is 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//
75// Output Container:
76// MWeight [MParameterD]
77//
78//////////////////////////////////////////////////////////////////////////////
79#include "MMcSpectrumWeight.h"
80
81#include <TF1.h>
82#include <TH1.h>
83#include <TSpline.h>
84
85#include "MLog.h"
86#include "MLogManip.h"
87
88#include "MParList.h"
89#include "MParameters.h"
90
91#include "MHillas.h"
92#include "MPointingPos.h"
93
94#include "MMcEvt.hxx"
95#include "MMcCorsikaRunHeader.h"
96
97ClassImp(MMcSpectrumWeight);
98
99using namespace std;
100
101void MMcSpectrumWeight::Init(const char *name, const char *title)
102{
103 fName = name ? name : "MMcSpectrumWeight";
104 fTitle = title ? title : "Task to calculate weights to change the energy spectrum";
105
106 AddToBranchList("MMcEvt.fEnergy");
107
108 fNameWeight = "MWeight";
109 fNameMcEvt = "MMcEvt";
110
111 fNewSlope = -9;
112 fOldSlope = -9;
113
114 fEnergyMin = -1;
115 fEnergyMax = -2;
116
117 fNorm = 1;
118 fNormEnergy = -1;
119
120 fAllowChange = kFALSE;
121
122 fFunc = NULL;
123 fMcEvt = NULL;
124 fHillas = NULL;
125 fWeight = NULL;
126 fWeightsZd = NULL;
127 fWeightsSize = NULL;
128 fPointing = NULL;
129}
130
131// ---------------------------------------------------------------------------
132//
133// Default Constructor.
134//
135MMcSpectrumWeight::MMcSpectrumWeight(const char *name, const char *title)
136{
137 Init(name,title);
138}
139
140// ---------------------------------------------------------------------------
141//
142// Destructor. If necessary delete fFunc
143//
144MMcSpectrumWeight::~MMcSpectrumWeight()
145{
146 if (fFunc)
147 delete fFunc;
148// if (fWeightsSize)
149// delete fWeightsSize;
150}
151
152// ---------------------------------------------------------------------------
153//
154// Search for
155// - fNameMcEvt [MMcEvtBasic]
156//
157// Find/Create:
158// - fNameWeight [MWeight]
159//
160Int_t MMcSpectrumWeight::PreProcess(MParList *pList)
161{
162 fMcEvt = (MMcEvt*)pList->FindObject(fNameMcEvt, "MMcEvtBasic");
163 if (!fMcEvt)
164 {
165 *fLog << err << fNameMcEvt << " [MMcEvtBasic] not found... abort." << endl;
166 return kFALSE;
167 }
168
169 fWeight = (MParameterD*)pList->FindCreateObj("MParameterD", fNameWeight);
170 if (!fWeight)
171 return kFALSE;
172
173 if (fWeightsZd)
174 {
175 fPointing = (MPointingPos*)pList->FindObject("MPointingPos");
176 if (!fPointing)
177 {
178 *fLog << err << "MPointingPos not found... abort." << endl;
179 return kFALSE;
180 }
181 }
182
183 if (fWeightsSize)
184 {
185 fHillas = (MHillas*)pList->FindObject("MHillas");
186 if (!fHillas)
187 {
188 *fLog << err << "MHillas not found... abort." << endl;
189 return kFALSE;
190 }
191 }
192
193 return kTRUE;
194}
195
196// ---------------------------------------------------------------------------
197//
198// Replace {fNameMcEvt}.fEnergy by "(x)" and return the result.
199//
200TString MMcSpectrumWeight::ReplaceX(TString str) const
201{
202 return str.ReplaceAll(Form("%s.fEnergy", fNameMcEvt.Data()), "(x)");
203}
204
205// ---------------------------------------------------------------------------
206//
207// Return the function corresponding to the mc spectrum with
208// slope fOldSlope: pow({fNameMcEvt}.fEnergy, fOldSlope)
209//
210// The slope is returned as %.3f
211//
212TString MMcSpectrumWeight::GetFormulaSpecOld() const
213{
214 return Form("pow(%s.fEnergy, %.3f)", fNameMcEvt.Data(), fOldSlope);
215}
216
217// ---------------------------------------------------------------------------
218//
219// Return the function corresponding to the new spectrum with
220// slope fNewSlope: pow({fNameMcEvt}.fEnergy, fNewSlope)
221//
222// The slope is returned as %.3f
223//
224// If a different formula is set (SetFormula()) this formula is returned
225// unchanged.
226//
227TString MMcSpectrumWeight::GetFormulaSpecNew() const
228{
229 TString str = fFormula.IsNull() ? Form("pow(%s.fEnergy, %.3f)", fNameMcEvt.Data(), fNewSlope) : fFormula.Data();
230 if (!fFormula.IsNull())
231 str.ReplaceAll("X", Form("(%s.fEnergy)", fNameMcEvt.Data()));
232
233 return str;
234}
235
236// ---------------------------------------------------------------------------
237//
238// Return the formula to calculate weights.
239// Is is compiled by
240// o1 = integral(fEnergyMin, fEnergyMax, GetFormulaSpecOldX());
241// n1 = integral(fEnergyMin, fEnergyMax, GetFormulaSpecNewX());
242// o2 = CalcSpecOld(fNormEnergy);
243// n2 = CalcSpecNew(fNormEnergy);
244//
245// result (fNormEnergy<0):
246// fNorm*o1/n1*GetFormulaNewSpec()/GetFormulaOldSpec()
247//
248// result (fNormEnergy>=0):
249// fNorm*o2/n2*GetFormulaNewSpec()/GetFormulaOldSpec()
250//
251// fNorm is 1 by default but can be overwritten using SetNorm()
252//
253// If the formulas GetFormulaSpecOldX() and GetFormulaSpecNewX()
254// are equal only fNorm is returned.
255//
256// The normalization constant is returned as %.16f
257//
258// Example: 0.3712780019*(pow(MMcEvt.fEnergy,-2.270))/(pow(MMcEvt.fEnergy,-2.600))
259//
260TString MMcSpectrumWeight::GetFormulaWeights() const
261{
262 if (GetFormulaSpecOld()==GetFormulaSpecNew())
263 return Form("%.16f", fNorm);
264
265 const Double_t iold = fNormEnergy<0 ? GetSpecOldIntegral() : CalcSpecOld(fNormEnergy);
266 const Double_t inew = fNormEnergy<0 ? GetSpecNewIntegral() : CalcSpecNew(fNormEnergy);
267
268 const Double_t norm = fNorm*iold/inew;
269
270 return Form("%.16f*(%s)/(%s)", norm, GetFormulaSpecNew().Data(), GetFormulaSpecOld().Data());
271}
272
273// ---------------------------------------------------------------------------
274//
275// Returns the integral between fEnergyMin and fEnergyMax of
276// GetFormulaSpecNewX() describing the destination spectrum
277//
278Double_t MMcSpectrumWeight::GetSpecNewIntegral() const
279{
280 TF1 funcnew("Dummy", GetFormulaSpecNewX());
281 return funcnew.Integral(fEnergyMin, fEnergyMax);
282}
283
284// ---------------------------------------------------------------------------
285//
286// Returns the integral between fEnergyMin and fEnergyMax of
287// GetFormulaSpecOldX() describing the simulated spectrum
288//
289Double_t MMcSpectrumWeight::GetSpecOldIntegral() const
290{
291 TF1 funcold("Dummy", GetFormulaSpecOldX());
292 return funcold.Integral(fEnergyMin, fEnergyMax);
293}
294
295// ---------------------------------------------------------------------------
296//
297// Returns the value of GetFormulaSpecNewX() at the energy e describing
298// the destination spectrum
299//
300Double_t MMcSpectrumWeight::CalcSpecNew(Double_t e) const
301{
302 TF1 funcnew("Dummy", GetFormulaSpecNewX());
303 return funcnew.Eval(e);
304}
305
306// ---------------------------------------------------------------------------
307//
308// Returns the value of GetFormulaSpecOldX() at the energy e describing
309// the simulated spectrum
310//
311Double_t MMcSpectrumWeight::CalcSpecOld(Double_t e) const
312{
313 TF1 funcnew("Dummy", GetFormulaSpecOldX());
314 return funcnew.Eval(e);
315}
316
317void MMcSpectrumWeight::SetWeightsSize(TH1D *h)
318{
319 fWeightsSize=h;
320 /*
321 if (h==0)
322 {
323 fWeightsSize=0;
324 return;
325 }
326
327 if (fWeightsSize)
328 delete fWeightsSize;
329
330 const Double_t xmin = TMath::Log10(h->GetXaxis()->GetXmin());
331 const Double_t xmax = TMath::Log10(h->GetXaxis()->GetXmax());
332 const Double_t xnum = h->GetNbinsX()+1;
333
334 fWeightsSize = new TSpline3("WeightsSize", xmin, xmax,
335 h->GetArray()+1, xnum);*/
336}
337
338// ---------------------------------------------------------------------------
339//
340// Initialize fEnergyMin, fEnergymax and fOldSlope from MMcCorsikaRunHeader
341// by GetELowLim(), GetEUppLim() and GetSlopeSpec().
342//
343// If fEnergyMax>fEnergyMin (means: the values have already been
344// initialized) and !fAllowChange the consistency of the new values
345// with the present values is checked with a numerical precision of 1e-10.
346// If one doesn't match kFALSE is returned.
347//
348// If the mc slope is -1 kFALSE is returned.
349//
350// If the new slope for the spectrum is -1 it is set to the original MC
351// slope.
352//
353// fFunc is set to the formula returned by GetFormulaWeightsX()
354//
355Bool_t MMcSpectrumWeight::Set(const MMcCorsikaRunHeader &rh)
356{
357 if (fEnergyMax>fEnergyMin && !fAllowChange)
358 {
359 if (TMath::Abs(fOldSlope-rh.GetSlopeSpec())>1e-10)
360 {
361 *fLog << err;
362 *fLog << "ERROR - The slope of the Monte Carlo is not allowed to change ";
363 *fLog << "(" << fOldSlope << " --> " << rh.GetSlopeSpec() << ")... abort." << endl;
364 return kFALSE;
365 }
366 if (TMath::Abs(fEnergyMin-rh.GetELowLim())>1e-10)
367 {
368 *fLog << err;
369 *fLog << "ERROR - The minimum simulated Monte Carlo energy is not allowed to change ";
370 *fLog << "(" << fEnergyMin << " --> " << rh.GetELowLim() << ")... abort." << endl;
371 return kFALSE;
372 }
373 if (TMath::Abs(fEnergyMax-rh.GetEUppLim())>1e-10)
374 {
375 *fLog << err;
376 *fLog << "ERROR - The maximum simulated Monte Carlo energy is not allowed to change (";
377 *fLog << "(" << fEnergyMax << " --> " << rh.GetEUppLim() << ")... abort." << endl;
378 return kFALSE;
379 }
380 return kTRUE;
381 }
382
383 fOldSlope = rh.GetSlopeSpec();
384 fEnergyMin = rh.GetELowLim();
385 fEnergyMax = rh.GetEUppLim();
386
387 if (fNewSlope==-9)
388 {
389 *fLog << inf << "The new slope of the power law is undefined (-9)... no weighting applied." << endl;
390 fNewSlope = fOldSlope;
391 }
392
393 if (fFunc)
394 delete fFunc;
395
396 fFunc = new TF1("", GetFormulaWeightsX());
397 gROOT->GetListOfFunctions()->Remove(fFunc);
398 fFunc->SetName("SpectralWeighs");
399
400 return kTRUE;
401}
402
403// ---------------------------------------------------------------------------
404//
405// The current contants are printed
406//
407void MMcSpectrumWeight::Print(Option_t *o) const
408{
409 *fLog << all << GetDescriptor() << endl;
410 *fLog << " Simulated energy range: " << fEnergyMin << "GeV - " << fEnergyMax << "GeV" << endl;
411 *fLog << " Simulated spectral slope: " << fOldSlope << endl;
412 *fLog << " New spectral slope: " << fNewSlope << endl;
413 *fLog << " Additional user norm.: " << fNorm << endl;
414 *fLog << " Spectra are normalized: " << (fNormEnergy<0?"by integral":Form("at %.1fGeV", fNormEnergy)) << endl;
415 *fLog << " Old Spectrum: " << GetFormulaSpecOldX() << " (I=" << GetSpecOldIntegral() << ")" << endl;
416 *fLog << " New Spectrum: " << GetFormulaSpecNewX() << " (I=" << GetSpecNewIntegral() << ")" << endl;
417 if (fFunc)
418 *fLog << " Weight function: " << fFunc->GetTitle() << endl;
419}
420
421// ----------------------------------------------------------------------------
422//
423// Executed each time a new root file is loaded
424// We will need fOldSlope and fE{Upp,Low}Lim to calculate the weights
425//
426Bool_t MMcSpectrumWeight::ReInit(MParList *plist)
427{
428 MMcCorsikaRunHeader *rh = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
429 if (!rh)
430 {
431 *fLog << err << "MMcCorsikaRunHeader not found... abort." << endl;
432 return kFALSE;
433 }
434
435 return Set(*rh);
436}
437
438// ----------------------------------------------------------------------------
439//
440// Fill the result of the evaluation of fFunc at fEvEvt->GetEnergy
441// into the weights container.
442//
443Int_t MMcSpectrumWeight::Process()
444{
445 Double_t w = 1;
446
447 if (fWeightsZd)
448 {
449 const Int_t i = fWeightsZd->GetXaxis()->FindFixBin(fPointing->GetZd());
450 w = fWeightsZd->GetBinContent(i);
451 }
452 if (fWeightsSize)
453 {
454 const Int_t i = fWeightsSize->GetXaxis()->FindFixBin(fHillas->GetSize());
455 w *= fWeightsSize->GetBinContent(i);
456 // w *= fWeightsSize->Eval(TMath::Log10(fHillas->GetSize()));
457 }
458
459 const Double_t e = fMcEvt->GetEnergy();
460 fWeight->SetVal(fFunc->Eval(e)*w);
461
462 return kTRUE;
463}
464
465// --------------------------------------------------------------------------
466//
467// Read the setup from a TEnv, eg:
468//
469// MMcSpectrumWeight.NewSlope: -2.6
470// The new slope of the spectrum
471//
472// MMcSpectrumWeight.Norm: 1.0
473// An additional artificial scale factor
474//
475// MMcSpectrumWeight.NormEnergy: 200
476// To normalize at a given energy instead of the integral
477//
478// MMcSpectrumWeight.Formula: pow(X, -2.6)
479// A formula to which the spectrum is weighted (use a capital X for
480// the energy)
481//
482Int_t MMcSpectrumWeight::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
483{
484 Bool_t rc = kFALSE;
485 if (IsEnvDefined(env, prefix, "NewSlope", print))
486 {
487 rc = kTRUE;
488 SetNewSlope(GetEnvValue(env, prefix, "NewSlope", fNewSlope));
489 }
490 if (IsEnvDefined(env, prefix, "Norm", print))
491 {
492 rc = kTRUE;
493 SetNorm(GetEnvValue(env, prefix, "Norm", fNorm));
494 }
495 if (IsEnvDefined(env, prefix, "NormEnergy", print))
496 {
497 rc = kTRUE;
498 SetNormEnergy(GetEnvValue(env, prefix, "NormEnergy", fNormEnergy));
499 }
500 if (IsEnvDefined(env, prefix, "Formula", print))
501 {
502 rc = kTRUE;
503 SetFormula(GetEnvValue(env, prefix, "Formula", fFormula));
504 }
505
506 return rc;
507}
Note: See TracBrowser for help on using the repository browser.