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

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