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

Last change on this file since 9483 was 9337, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 21.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-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 spectrum: 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 "MString.h"
93#include "MParList.h"
94#include "MParameters.h"
95
96#include "MHillas.h"
97#include "MPointingPos.h"
98
99#include "MMcEvt.hxx"
100#include "MMcCorsikaRunHeader.h"
101
102ClassImp(MMcSpectrumWeight);
103
104using namespace std;
105
106void MMcSpectrumWeight::Init(const char *name, const char *title)
107{
108 fName = name ? name : "MMcSpectrumWeight";
109 fTitle = title ? title : "Task to calculate weights to change the energy spectrum";
110
111 AddToBranchList("MMcEvt.fEnergy");
112
113 fNameWeight = "MWeight";
114 fNameMcEvt = "MMcEvt";
115
116 fNewSlope = -99;
117 fOldSlope = -99;
118
119 fEnergyMin = -1;
120 fEnergyMax = -2;
121
122 fNorm = 1;
123 fNormEnergy = -1;
124
125 fAllowChange = kFALSE;
126
127 fFunc = NULL;
128 fMcEvt = NULL;
129 fHillas = NULL;
130 fWeight = NULL;
131 fWeightsZd = NULL;
132 fWeightsSize = NULL;
133 fPointing = NULL;
134}
135
136// ---------------------------------------------------------------------------
137//
138// Default Constructor.
139//
140MMcSpectrumWeight::MMcSpectrumWeight(const char *name, const char *title)
141{
142 Init(name,title);
143}
144
145// ---------------------------------------------------------------------------
146//
147// Destructor. If necessary delete fFunc
148//
149MMcSpectrumWeight::~MMcSpectrumWeight()
150{
151 if (fFunc)
152 delete fFunc;
153// if (fWeightsSize)
154// delete fWeightsSize;
155}
156
157// ---------------------------------------------------------------------------
158//
159// Search for
160// - fNameMcEvt [MMcEvtBasic]
161//
162// Find/Create:
163// - fNameWeight [MWeight]
164//
165Int_t MMcSpectrumWeight::PreProcess(MParList *pList)
166{
167 fMcEvt = (MMcEvt*)pList->FindObject(fNameMcEvt, "MMcEvtBasic");
168 if (!fMcEvt)
169 {
170 *fLog << err << fNameMcEvt << " [MMcEvtBasic] not found... abort." << endl;
171 return kFALSE;
172 }
173
174 fWeight = (MParameterD*)pList->FindCreateObj("MParameterD", fNameWeight);
175 if (!fWeight)
176 return kFALSE;
177
178 if (fWeightsZd)
179 {
180 fPointing = (MPointingPos*)pList->FindObject("MPointingPos");
181 if (!fPointing)
182 {
183 *fLog << err << "MPointingPos not found... abort." << endl;
184 return kFALSE;
185 }
186 }
187
188 if (fWeightsSize)
189 {
190 fHillas = (MHillas*)pList->FindObject("MHillas");
191 if (!fHillas)
192 {
193 *fLog << err << "MHillas not found... abort." << endl;
194 return kFALSE;
195 }
196 }
197
198 return kTRUE;
199}
200
201// ---------------------------------------------------------------------------
202//
203// Replace {fNameMcEvt}.fEnergy by "(x)" and return the result.
204//
205TString MMcSpectrumWeight::ReplaceX(TString str) const
206{
207 return str.ReplaceAll(MString::Format("%s.fEnergy", fNameMcEvt.Data()), "(x)");
208}
209
210// ---------------------------------------------------------------------------
211//
212// Return the function corresponding to the mc spectrum with
213// slope fOldSlope: pow({fNameMcEvt}.fEnergy, fOldSlope)
214//
215// The slope is returned as %.3f
216//
217TString MMcSpectrumWeight::GetFormulaSpecOld(const char *name) const
218{
219 return MString::Format("pow(%s.fEnergy, %.3f)", name, fOldSlope);
220}
221
222// ---------------------------------------------------------------------------
223//
224// Return the function corresponding to the new spectrum with
225// slope fNewSlope: pow({fNameMcEvt}.fEnergy, fNewSlope)
226//
227// The slope is returned as %.3f
228//
229// If a different formula is set (SetFormula()) this formula is returned
230// unchanged.
231//
232TString MMcSpectrumWeight::GetFormulaSpecNew(const char *name) const
233{
234 TString str = fFormula.IsNull() ? MString::Format("pow(%s.fEnergy, %.3f)", name, fNewSlope) : fFormula;
235 if (!fFormula.IsNull())
236 str.ReplaceAll("X", MString::Format("(%s.fEnergy)", name));
237
238 return str;
239}
240
241// ---------------------------------------------------------------------------
242//
243// Return the formula to calculate weights.
244// Is is compiled by
245// o1 = integral(fEnergyMin, fEnergyMax, GetFormulaSpecOldX());
246// n1 = integral(fEnergyMin, fEnergyMax, GetFormulaSpecNewX());
247// o2 = CalcSpecOld(fNormEnergy);
248// n2 = CalcSpecNew(fNormEnergy);
249//
250// result (fNormEnergy<0):
251// fNorm*o1/n1*GetFormulaNewSpec()/GetFormulaOldSpec()
252//
253// result (fNormEnergy>=0):
254// fNorm*o2/n2*GetFormulaNewSpec()/GetFormulaOldSpec()
255//
256// fNorm is 1 by default but can be overwritten using SetNorm()
257//
258// If the formulas GetFormulaSpecOldX() and GetFormulaSpecNewX()
259// are equal only fNorm is returned.
260//
261// The normalization constant is returned as %.16e
262//
263// Example: 0.3712780019*(pow(MMcEvt.fEnergy,-2.270))/(pow(MMcEvt.fEnergy,-2.600))
264//
265TString MMcSpectrumWeight::GetFormulaWeights(const char *name) const
266{
267 if (GetFormulaSpecOld()==GetFormulaSpecNew())
268 return MString::Format("%.16e", fNorm);
269
270 const Double_t iold = fNormEnergy<0 ? GetSpecOldIntegral() : CalcSpecOld(fNormEnergy);
271 const Double_t inew = fNormEnergy<0 ? GetSpecNewIntegral() : CalcSpecNew(fNormEnergy);
272
273 const Double_t norm = fNorm*iold/inew;
274
275 return MString::Format("%.16e*(%s)/(%s)", norm, GetFormulaSpecNew(name).Data(), GetFormulaSpecOld(name).Data());
276}
277
278// ---------------------------------------------------------------------------
279//
280// Returns the integral between fEnergyMin and fEnergyMax of
281// GetFormulaSpecNewX() describing the destination spectrum
282//
283Double_t MMcSpectrumWeight::GetSpecNewIntegral(Double_t emin, Double_t emax) const
284{
285 TF1 funcnew("Dummy", GetFormulaSpecNewX().Data());
286 return funcnew.Integral(emin, emax);
287}
288
289// ---------------------------------------------------------------------------
290//
291// Returns the integral between fEnergyMin and fEnergyMax of
292// GetFormulaSpecOldX() describing the simulated spectrum
293//
294Double_t MMcSpectrumWeight::GetSpecOldIntegral(Double_t emin, Double_t emax) const
295{
296 TF1 funcold("Dummy", GetFormulaSpecOldX().Data());
297 return funcold.Integral(emin, emax);
298}
299
300// ---------------------------------------------------------------------------
301//
302// Returns the value of GetFormulaSpecNewX() at the energy e describing
303// the destination spectrum
304//
305Double_t MMcSpectrumWeight::CalcSpecNew(Double_t e) const
306{
307 TF1 funcnew("Dummy", GetFormulaSpecNewX().Data());
308 return funcnew.Eval(e);
309}
310
311// ---------------------------------------------------------------------------
312//
313// Returns the value of GetFormulaSpecOldX() at the energy e describing
314// the simulated spectrum
315//
316Double_t MMcSpectrumWeight::CalcSpecOld(Double_t e) const
317{
318 TF1 funcnew("Dummy", GetFormulaSpecOldX().Data());
319 return funcnew.Eval(e);
320}
321
322void MMcSpectrumWeight::SetWeightsSize(TH1D *h)
323{
324 fWeightsSize=h;
325 /*
326 if (h==0)
327 {
328 fWeightsSize=0;
329 return;
330 }
331
332 if (fWeightsSize)
333 delete fWeightsSize;
334
335 const Double_t xmin = TMath::Log10(h->GetXaxis()->GetXmin());
336 const Double_t xmax = TMath::Log10(h->GetXaxis()->GetXmax());
337 const Double_t xnum = h->GetNbinsX()+1;
338
339 fWeightsSize = new TSpline3("WeightsSize", xmin, xmax,
340 h->GetArray()+1, xnum);*/
341}
342
343// ---------------------------------------------------------------------------
344//
345// Initialize fEnergyMin, fEnergymax and fOldSlope from MMcCorsikaRunHeader
346// by GetELowLim(), GetEUppLim() and GetSlopeSpec().
347//
348// If fEnergyMax>fEnergyMin (means: the values have already been
349// initialized) and !fAllowChange the consistency of the new values
350// with the present values is checked with a numerical precision of 1e-10.
351// If one doesn't match kFALSE is returned.
352//
353// If the mc slope is -1 kFALSE is returned.
354//
355// If the new slope for the spectrum is -1 it is set to the original MC
356// slope.
357//
358// fFunc is set to the formula returned by GetFormulaWeightsX()
359//
360Bool_t MMcSpectrumWeight::Set(const MMcCorsikaRunHeader &rh)
361{
362 if (fEnergyMax>fEnergyMin && !fAllowChange)
363 {
364 if (TMath::Abs(fEnergyMax-rh.GetEUppLim())>1e-10)
365 {
366 *fLog << err;
367 *fLog << "ERROR - The maximum simulated Monte Carlo energy is not allowed to change ";
368 *fLog << "(" << fEnergyMax << " --> " << rh.GetEUppLim() << ")... abort." << endl;
369 return kFALSE;
370 }
371
372 if (TMath::Abs(fOldSlope-rh.GetSlopeSpec())>1e-10)
373 {
374 *fLog << err;
375 *fLog << "ERROR - The slope of the Monte Carlo is not allowed to change ";
376 *fLog << "(" << fOldSlope << " --> " << rh.GetSlopeSpec() << ")... abort." << endl;
377 return kFALSE;
378 }
379
380 // No change happened
381 if (TMath::Abs(fEnergyMin-rh.GetELowLim())<=1e-10)
382 return kTRUE;
383
384 // The lower energy limit has changed
385 *fLog << warn;
386 *fLog << "The minimum simulated Monte Carlo energy has changed from ";
387 *fLog << fEnergyMin << "GeV to " << rh.GetELowLim() << "GeV." << endl;
388 fEnergyMin = rh.GetELowLim();
389 }
390
391 if (fNormEnergy<0 && fEnergyMin>0 && TMath::Abs(fEnergyMin-rh.GetELowLim())>1e-10)
392 {
393 *fLog << err;
394 *fLog << "You try to use changing minimum simulated Monte Carlo energies" << endl;
395 *fLog << "together with a normalization calculated from the integral." << endl;
396 *fLog << "This is not yet supported. Please switch to a normalization" << endl;
397 *fLog << "at a dedicated energy by specifying the energy" << endl;
398 *fLog << " MMcSpectrumWeight.NormEnergy: 500" << endl;
399 *fLog << "in your sponde.rc." << endl;
400 return kFALSE;
401 }
402
403 fOldSlope = rh.GetSlopeSpec();
404 fEnergyMin = rh.GetELowLim();
405 fEnergyMax = rh.GetEUppLim();
406
407 if (fNewSlope==-99 && fFormula.IsNull())
408 {
409 *fLog << inf << "A new slope for the power law has not yet been defined... using " << fOldSlope << "." << endl;
410 fNewSlope = fOldSlope;
411 }
412
413 if (fFunc)
414 delete fFunc;
415
416 if (GetFormulaSpecOld()==GetFormulaSpecNew())
417 *fLog << inf << "No spectral change requested..." << endl;
418 else
419 {
420 *fLog << inf << "Weighting from slope " << fOldSlope << " to ";
421 if (fFormula.IsNull())
422 *fLog << "slope " << fNewSlope << "." << endl;
423 else
424 *fLog << GetFormulaSpecNewX() << endl;
425 }
426
427 fFunc = new TF1("", GetFormulaWeightsX().Data());
428 fFunc->SetName("SpectralWeighs");
429 gROOT->GetListOfFunctions()->Remove(fFunc);
430
431 return kTRUE;
432}
433
434// ---------------------------------------------------------------------------
435//
436// completes a simulated spectrum starting at an energy fEnergyMin down to
437// an energy emin.
438//
439// It is assumed that the contents of MMcSpectrumWeight for the new spectrum
440// correctly describe the spectrum within the histogram, and fEnergyMin
441// and fEnergyMax correctly describe the range.
442//
443// If scale is given the histogram statistics is further extended by the
444// new spectrum according to the scale factor (eg. 1.2: by 20%)
445//
446// In the 1D case it is assumed that the x-axis is a zenith angle binning.
447// In the 2D case the x-axis is assumed to be zenith angle, the y-axis
448// to be energy.
449//
450void MMcSpectrumWeight::CompleteEnergySpectrum(TH1 &h, Double_t emin, Double_t scale) const
451{
452 if (h.InheritsFrom(TH3::Class()))
453 {
454 return;
455 }
456
457 if (fEnergyMin < emin)
458 {
459 *fLog << err << "ERROR - MMcSpctrumWeight::CompleteEnergySpectrum: fEnergyMin (";
460 *fLog << fEnergyMin << ") smaller than emin (" << emin << ")." << endl;
461 return;
462 }
463
464 // Total number of events for the new spectrum in the same
465 // energy range as the current histogram is filled
466 const Double_t norm = GetSpecNewIntegral();
467
468 // Check if it is only a histogram in ZA
469 if (!h.InheritsFrom(TH2::Class()))
470 {
471 // Warning: Simply scaling the zenith angle distribution might
472 // increase fluctuations for low statistics.
473 const Double_t f = GetSpecNewIntegral(emin, fEnergyMax)/norm;
474 h.Scale(f*scale);
475 return;
476 }
477
478 const TAxis &axey = *h.GetYaxis();
479
480 // Find energy range between the minimum energy to be filled (emin)
481 // and the minimum energy corresponding to the data filled into
482 // this histogram (fEnergyMin)
483 const Int_t first = axey.FindFixBin(emin);
484 const Int_t last = axey.FindFixBin(fEnergyMin); // data range min energy
485 const Int_t max = axey.FindFixBin(fEnergyMax); // data range max energy
486
487 for (int x=1; x<=h.GetNbinsX(); x++)
488 {
489 // Ratio between the number of events in the zenith angle
490 // bin corresponding to x and the new spectrum.
491 const Double_t f = h.Integral(x, x, -1, 9999)/norm;
492
493 // Fill histogram with the "new spectrum" between
494 // emin and fEnergyMin.
495 if (emin<fEnergyMin)
496 for (int y=first; y<=last; y++)
497 {
498 // Check if the bin is only partly filled by the energy range
499 const Double_t lo = axey.GetBinLowEdge(y) <emin ? emin : axey.GetBinLowEdge(y);
500 const Double_t hi = axey.GetBinLowEdge(y+1)>fEnergyMin ? fEnergyMin : axey.GetBinLowEdge(y+1);
501
502 // Add the new spectrum extending the existing spectrum
503 h.AddBinContent(h.GetBin(x, y), f*GetSpecNewIntegral(lo, hi));
504 }
505
506 // If scale is >1 we also have to increse the statistics f the
507 // histogram according to scale.
508 if (scale>1)
509 for (int y=first; y<=max; y++)
510 {
511 // Check if the bin is only partly filled by the energy range
512 const Double_t lo = axey.GetBinLowEdge(y) <emin ? emin : axey.GetBinLowEdge(y);
513 const Double_t hi = axey.GetBinLowEdge(y+1)>fEnergyMax ? fEnergyMax : axey.GetBinLowEdge(y+1);
514
515 // Use the analytical solution to scale the histogram
516 h.AddBinContent(h.GetBin(x, y), f*GetSpecNewIntegral(lo, hi)*(scale-1));
517 }
518 }
519}
520
521// ---------------------------------------------------------------------------
522//
523// The current contants are printed
524//
525void MMcSpectrumWeight::Print(Option_t *o) const
526{
527 const TString opt(o);
528
529 const Bool_t hasnew = opt.Contains("new") || opt.IsNull();
530 const Bool_t hasold = opt.Contains("old") || opt.IsNull();
531
532 *fLog << all << GetDescriptor() << endl;
533
534 if (hasold)
535 {
536 *fLog << " Simulated energy range: " << fEnergyMin << "GeV - " << fEnergyMax << "GeV" << endl;
537 *fLog << " Simulated spectral slope: ";
538 if (fOldSlope==-99)
539 *fLog << "undefined" << endl;
540 else
541 *fLog << fOldSlope << endl;
542 }
543 if (hasnew)
544 {
545 *fLog << " New spectral slope: ";
546 if (fNewSlope==-99)
547 *fLog << "undefined/no change" << endl;
548 else
549 *fLog << fNewSlope << endl;
550 }
551 *fLog << " Additional user norm.: " << fNorm << endl;
552 *fLog << " Spectra are normalized: " << (fNormEnergy<0?"by integral":MString::Format("at %.1fGeV", fNormEnergy)) << endl;
553 if (hasold)
554 {
555 *fLog << " Old Spectrum: ";
556 if (fNewSlope==-99)
557 *fLog << "undefined";
558 else
559 *fLog << GetFormulaSpecOldX();
560 if (fEnergyMin>=0 && fEnergyMax>0)
561 *fLog << " (I=" << GetSpecOldIntegral() << ")";
562 *fLog << endl;
563 }
564 if (hasnew)
565 {
566 *fLog << " New Spectrum: ";
567 if (fNewSlope==-99 && fFormula.IsNull())
568 *fLog << "undefined/no change";
569 else
570 *fLog << GetFormulaSpecNewX();
571 if (fEnergyMin>=0 && fEnergyMax>0)
572 *fLog << " (I=" << GetSpecNewIntegral() << ")";
573 *fLog << endl;
574 }
575 if (fFunc)
576 *fLog << " Weight func: " << fFunc->GetTitle() << endl;
577}
578
579// ----------------------------------------------------------------------------
580//
581// Executed each time a new root file is loaded
582// We will need fOldSlope and fE{Upp,Low}Lim to calculate the weights
583//
584Bool_t MMcSpectrumWeight::ReInit(MParList *plist)
585{
586 MMcCorsikaRunHeader *rh = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
587 if (!rh)
588 {
589 *fLog << err << "MMcCorsikaRunHeader not found... abort." << endl;
590 return kFALSE;
591 }
592
593 return Set(*rh);
594}
595
596/*
597 * This could be used to improve the Zd-weighting within a bin.
598 * Another option is to use more bins, or both.
599 * Note that it seems unnecessary, because the shape within the
600 * theta-bins should be similar in data and Monte Carlo... hopefully.
601 *
602void MMcSpectrumWeight::InitZdWeights()
603{
604 TH2D w(*fWeightsZd);
605
606 for (int i=1; i<=w.GetNbinsX(); i++)
607 {
608 const Double_t tmin = w.GetBinLowEdge(i) *TMath::DegToRad();
609 const Double_t tmax = w.GetBinLowEdge(i+1)*TMath::DegToRad();
610
611 const Double_t wdth = tmax-tmin;
612 const Double_t integ = cos(tmin)-cos(tmax);
613
614 w.SetBinContent(i, w.GetBinContent(i)*wdth/integ);
615 }
616
617 // const Int_t i = fWeightsZd->GetXaxis()->FindFixBin(fPointing->GetZd());
618 // const Double_t theta = fPointing->GetZd()*TMath::DegToRad();
619 // w = sin(theta)*w.GetBinContent(i);
620}
621*/
622
623// ----------------------------------------------------------------------------
624//
625// Fill the result of the evaluation of fFunc at fEvEvt->GetEnergy
626// into the weights container.
627//
628Int_t MMcSpectrumWeight::Process()
629{
630 Double_t w = 1;
631
632 if (fWeightsZd)
633 {
634 const Int_t i = fWeightsZd->GetXaxis()->FindFixBin(fPointing->GetZd());
635 w = fWeightsZd->GetBinContent(i);
636 }
637 if (fWeightsSize)
638 {
639 const Int_t i = fWeightsSize->GetXaxis()->FindFixBin(fHillas->GetSize());
640 w *= fWeightsSize->GetBinContent(i);
641 // w *= fWeightsSize->Eval(TMath::Log10(fHillas->GetSize()));
642 }
643
644 const Double_t e = fMcEvt->GetEnergy();
645
646 fWeight->SetVal(fFunc->Eval(e)*w);
647
648 return kTRUE;
649}
650
651// --------------------------------------------------------------------------
652//
653// Read the setup from a TEnv, eg:
654//
655// MMcSpectrumWeight.NewSlope: -2.6
656// The new slope of the spectrum
657//
658// MMcSpectrumWeight.Norm: 1.0
659// An additional artificial scale factor
660//
661// MMcSpectrumWeight.NormEnergy: 200
662// To normalize at a given energy instead of the integral
663//
664// MMcSpectrumWeight.Formula: pow(X, -2.6)
665// A formula to which the spectrum is weighted (use a capital X for
666// the energy)
667//
668Int_t MMcSpectrumWeight::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
669{
670 Bool_t rc = kFALSE;
671 if (IsEnvDefined(env, prefix, "NewSlope", print))
672 {
673 rc = kTRUE;
674 SetNewSlope(GetEnvValue(env, prefix, "NewSlope", fNewSlope));
675 }
676 if (IsEnvDefined(env, prefix, "Norm", print))
677 {
678 rc = kTRUE;
679 SetNorm(GetEnvValue(env, prefix, "Norm", fNorm));
680 }
681 if (IsEnvDefined(env, prefix, "NormEnergy", print))
682 {
683 rc = kTRUE;
684 SetNormEnergy(GetEnvValue(env, prefix, "NormEnergy", fNormEnergy));
685 }
686 if (IsEnvDefined(env, prefix, "Formula", print))
687 {
688 rc = kTRUE;
689 SetFormula(GetEnvValue(env, prefix, "Formula", fFormula));
690 }
691
692 return rc;
693}
Note: See TracBrowser for help on using the repository browser.