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