source: tags/Mars-V0.9.4/mhflux/MAlphaFitter.cc

Last change on this file was 7224, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 20.6 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, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MAlphaFitter
28//
29// Create a single Alpha-Plot. The alpha-plot is fitted online. You can
30// check the result when it is filles in the MStatusDisplay
31// For more information see MHFalseSource::FitSignificance
32//
33// For convinience (fit) the output significance is stored in a
34// container in the parlisrt
35//
36// PRELIMINARY!
37//
38//////////////////////////////////////////////////////////////////////////////
39#include "MAlphaFitter.h"
40
41#include <TF1.h>
42#include <TH1.h>
43#include <TH3.h>
44
45#include <TRandom.h>
46
47#include <TLine.h>
48#include <TLatex.h>
49#include <TVirtualPad.h>
50
51#include "MMath.h"
52
53#include "MLogManip.h"
54
55ClassImp(MAlphaFitter);
56
57using namespace std;
58
59void MAlphaFitter::Clear(Option_t *o)
60{
61 fSignificance=0;
62 fEventsExcess=0;
63 fEventsSignal=0;
64 fEventsBackground=0;
65
66 fChiSqSignal=0;
67 fChiSqBg=0;
68 fIntegralMax=0;
69 fScaleFactor=1;
70
71 fCoefficients.Reset();
72}
73
74// --------------------------------------------------------------------------
75//
76// This is a preliminary implementation of a alpha-fit procedure for
77// all possible source positions. It will be moved into its own
78// more powerfull class soon.
79//
80// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
81// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
82// or
83// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
84//
85// Parameter [1] is fixed to 0 while the alpha peak should be
86// symmetric around alpha=0.
87//
88// Parameter [4] is fixed to 0 because the first derivative at
89// alpha=0 should be 0, too.
90//
91// In a first step the background is fitted between bgmin and bgmax,
92// while the parameters [0]=0 and [2]=1 are fixed.
93//
94// In a second step the signal region (alpha<sigmax) is fittet using
95// the whole function with parameters [1], [3], [4] and [5] fixed.
96//
97// The number of excess and background events are calculated as
98// s = int(hist, 0, 1.25*sigint)
99// b = int(pol2(3), 0, 1.25*sigint)
100//
101// The Significance is calculated using the Significance() member
102// function.
103//
104Bool_t MAlphaFitter::Fit(TH1D &h, Bool_t paint)
105{
106 Clear();
107 if (h.GetEntries()==0)
108 return kFALSE;
109
110 Double_t sigmax=fSigMax;
111 Double_t bgmin =fBgMin;
112 Double_t bgmax =fBgMax;
113
114 //*fLog << inf << "Fit: " << sigmax << " " << fSigInt << " " << bgmin << " " << bgmax << endl;
115
116 //TF1 fFunc("", Form("gaus(0) + pol%d(3)", fPolynomOrder), 0, 90);
117
118 //fFunc->SetParameters(fCoefficients.GetArray());
119
120 fFunc->FixParameter(1, 0);
121 if (fPolynomOrder!=1)
122 fFunc->FixParameter(4, 0);
123 fFunc->SetParLimits(2, 0, 90);
124 fFunc->SetParLimits(3, -1, 1);
125
126 const Double_t alpha0 = h.GetBinContent(1);
127 const Double_t alphaw = h.GetXaxis()->GetBinWidth(1);
128
129 // Check for the regios which is not filled...
130 if (alpha0==0)
131 return kFALSE; //*fLog << warn << "Histogram empty." << endl;
132
133 // First fit a polynom in the off region
134 fFunc->FixParameter(0, 0);
135 fFunc->FixParameter(2, 1);
136 fFunc->ReleaseParameter(3);
137
138 for (int i=5; i<fFunc->GetNpar(); i++)
139 if (fFitBackground)
140 fFunc->ReleaseParameter(i);
141 else
142 fFunc->SetParameter(i, 0);
143
144 // options : N do not store the function, do not draw
145 // I use integral of function in bin rather than value at bin center
146 // R use the range specified in the function range
147 // Q quiet mode
148 // E Perform better Errors estimation using Minos technique
149 if (fFitBackground)
150 {
151 h.Fit(fFunc, "NQI", "", bgmin, bgmax);
152 fChiSqBg = fFunc->GetChisquare()/fFunc->GetNDF();
153 }
154
155
156 // ------------------------------------
157 if (paint && fFitBackground)
158 {
159 fFunc->SetRange(0, 90);
160 fFunc->SetLineColor(kRed);
161 fFunc->SetLineWidth(2);
162 fFunc->Paint("same");
163 }
164 // ------------------------------------
165
166 fFunc->ReleaseParameter(0);
167 //func.ReleaseParameter(1);
168 fFunc->ReleaseParameter(2);
169 fFunc->FixParameter(3, fFunc->GetParameter(3));
170 fFunc->FixParameter(4, fFunc->GetParameter(4));
171 for (int i=5; i<fFunc->GetNpar(); i++)
172 fFunc->FixParameter(i, fFunc->GetParameter(i));
173
174 // Do not allow signals smaller than the background
175 const Double_t A = alpha0-fFunc->GetParameter(3);
176 const Double_t dA = TMath::Abs(A);
177 fFunc->SetParLimits(0, -dA*4, dA*4);
178 fFunc->SetParLimits(2, 0, 90);
179
180 // Now fit a gaus in the on region on top of the polynom
181 fFunc->SetParameter(0, A);
182 fFunc->SetParameter(2, sigmax*0.75);
183
184 // options : N do not store the function, do not draw
185 // I use integral of function in bin rather than value at bin center
186 // R use the range specified in the function range
187 // Q quiet mode
188 // E Perform better Errors estimation using Minos technique
189 h.Fit(fFunc, "NQI", "", 0, sigmax);
190
191 fChiSqSignal = fFunc->GetChisquare()/fFunc->GetNDF();
192 fCoefficients.Set(fFunc->GetNpar(), fFunc->GetParameters());
193
194 //const Bool_t ok = NDF>0 && chi2<2.5*NDF;
195
196 // ------------------------------------
197 if (paint)
198 {
199 fFunc->SetLineColor(kGreen);
200 fFunc->SetLineWidth(2);
201 fFunc->Paint("same");
202 }
203 // ------------------------------------
204
205 //const Double_t s = fFunc->Integral(0, fSigInt)/alphaw;
206 fFunc->SetParameter(0, 0);
207 fFunc->SetParameter(2, 1);
208 //const Double_t b = fFunc->Integral(0, fSigInt)/alphaw;
209 //fSignificance = MMath::SignificanceLiMaSigned(s, b);
210
211 const Int_t bin = h.GetXaxis()->FindFixBin(fSigInt*0.999);
212
213 fIntegralMax = h.GetBinLowEdge(bin+1);
214 fEventsBackground = fFunc->Integral(0, fIntegralMax)/alphaw;
215 fEventsSignal = h.Integral(1, bin);
216 fEventsExcess = fEventsSignal-fEventsBackground;
217 fSignificance = MMath::SignificanceLiMaSigned(fEventsSignal, fEventsBackground);
218
219 if (TMath::IsNaN(fSignificance))
220 fSignificance=0;
221
222 if (fEventsExcess<0)
223 fEventsExcess=0;
224
225 return kTRUE;
226}
227
228Bool_t MAlphaFitter::Fit(const TH1D &hon, const TH1D &hof, Double_t alpha, Bool_t paint)
229{
230 TH1D h(hon);
231 h.Add(&hof, -1); // substracts also number of entries!
232 h.SetEntries(hon.GetEntries());
233
234 MAlphaFitter fit(*this);
235 fit.SetPolynomOrder(0);
236
237 if (alpha<=0 || !fit.Fit(h, paint))
238 return kFALSE;
239
240 fChiSqSignal = fit.GetChiSqSignal();
241 fChiSqBg = fit.GetChiSqBg();
242 fCoefficients = fit.GetCoefficients();
243
244 const Int_t bin = hon.GetXaxis()->FindFixBin(fSigInt*0.999);
245
246 fIntegralMax = hon.GetBinLowEdge(bin+1);
247 fEventsBackground = hof.Integral(1, bin);
248 fEventsSignal = hon.Integral(1, bin);
249 fEventsExcess = fEventsSignal-fEventsBackground;
250 fScaleFactor = alpha;
251 fSignificance = MMath::SignificanceLiMaSigned(fEventsSignal, fEventsBackground/alpha, alpha);
252
253 if (TMath::IsNaN(fSignificance))
254 fSignificance=0;
255 if (fEventsExcess<0)
256 fEventsExcess=0;
257
258 return kTRUE;
259}
260
261void MAlphaFitter::PaintResult(Float_t x, Float_t y, Float_t size) const
262{
263 const Double_t w = GetGausSigma();
264 const Double_t m = fIntegralMax;
265
266 const Int_t l1 = w<=0 ? 0 : (Int_t)TMath::Ceil(-TMath::Log10(w));
267 const Int_t l2 = m<=0 ? 0 : (Int_t)TMath::Ceil(-TMath::Log10(m));
268 const TString fmt = Form("\\sigma_{L/M}=%%.1f \\omega=%%.%df\\circ E=%%d B=%%d x<%%.%df \\tilde\\chi_{b}=%%.1f \\tilde\\chi_{s}=%%.1f c=%%.1f f=%%.2f",
269 l1<1?1:l1+1, l2<1?1:l2+1);
270
271 TLatex text(x, y, Form(fmt.Data(), fSignificance, w, (int)fEventsExcess,
272 (int)fEventsBackground, m, fChiSqBg, fChiSqSignal,
273 fCoefficients[3], fScaleFactor));
274
275 text.SetBit(TLatex::kTextNDC);
276 text.SetTextSize(size);
277 text.Paint();
278
279 TLine line;
280 line.SetLineColor(14);
281 line.PaintLine(m, gPad->GetUymin(), m, gPad->GetUymax());
282}
283
284void MAlphaFitter::Copy(TObject &o) const
285{
286 MAlphaFitter &f = static_cast<MAlphaFitter&>(o);
287
288 // Setup
289 f.fSigInt = fSigInt;
290 f.fSigMax = fSigMax;
291 f.fBgMin = fBgMin;
292 f.fBgMax = fBgMax;
293 f.fScaleMin = fScaleMin;
294 f.fScaleMax = fScaleMax;
295 f.fPolynomOrder = fPolynomOrder;
296 f.fFitBackground= fFitBackground;
297 f.fSignalFunc = fSignalFunc;
298 f.fScaleMode = fScaleMode;
299 f.fScaleUser = fScaleUser;
300 f.fStrategy = fStrategy;
301 f.fCoefficients.Set(fCoefficients.GetSize());
302 f.fCoefficients.Reset();
303
304 // Result
305 f.fSignificance = fSignificance;
306 f.fEventsExcess = fEventsExcess;
307 f.fEventsSignal = fEventsSignal;
308 f.fEventsBackground = fEventsBackground;
309 f.fChiSqSignal = fChiSqSignal;
310 f.fChiSqBg = fChiSqBg;
311 f.fIntegralMax = fIntegralMax;
312 f.fScaleFactor = fScaleFactor;
313
314 // Function
315 TF1 *fcn = f.fFunc;
316 f.fFunc = new TF1(*fFunc);
317 gROOT->GetListOfFunctions()->Remove(f.fFunc);
318 f.fFunc->SetName("Dummy");
319 delete fcn;
320}
321
322void MAlphaFitter::Print(Option_t *o) const
323{
324 *fLog << GetDescriptor() << ": Fitting..." << endl;
325 *fLog << " ...signal to " << fSigMax << " (integrate into bin at " << fSigInt << ")" << endl;
326 *fLog << " ...signal function: ";
327 switch (fSignalFunc)
328 {
329 case kGauss: *fLog << "gauss(x)"; break;
330 case kThetaSq: *fLog << "gauss(sqrt(x))"; break;
331 }
332 *fLog << endl;
333 if (!fFitBackground)
334 *fLog << " ...no background." << endl;
335 else
336 {
337 *fLog << " ...background from " << fBgMin << " to " << fBgMax << endl;
338 *fLog << " ...polynom order " << fPolynomOrder << endl;
339 *fLog << " ...scale mode: ";
340 switch (fScaleMode)
341 {
342 case kNone: *fLog << "none."; break;
343 case kEntries: *fLog << "entries."; break;
344 case kIntegral: *fLog << "integral."; break;
345 case kOffRegion: *fLog << "off region (integral between " << fScaleMin << " and " << fScaleMax << ")"; break;
346 case kBackground: *fLog << "background (integral between " << fBgMin << " and " << fBgMax << ")"; break;
347 case kLeastSquare: *fLog << "least square (N/A)"; break;
348 case kUserScale: *fLog << "user def (" << fScaleUser << ")"; break;
349 }
350 *fLog << endl;
351 }
352
353 if (TString(o).Contains("result"))
354 {
355 *fLog << "Result:" << endl;
356 *fLog << " - Significance (Li/Ma) " << fSignificance << endl;
357 *fLog << " - Excess Events " << fEventsExcess << endl;
358 *fLog << " - Signal Events " << fEventsSignal << endl;
359 *fLog << " - Background Events " << fEventsBackground << endl;
360 *fLog << " - Chi^2/ndf (Signal) " << fChiSqSignal << endl;
361 *fLog << " - Chi^2/ndf (Background) " << fChiSqBg << endl;
362 *fLog << " - Signal integrated up to " << fIntegralMax << "°" << endl;
363 *fLog << " - Scale Factor (Off) " << fScaleFactor << endl;
364 }
365}
366
367Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, UInt_t bin, Bool_t paint)
368{
369 const TString name(Form("TempAlphaEnergy%06d", gRandom->Integer(1000000)));
370 TH1D *h = hon.ProjectionZ(name, -1, 9999, bin, bin, "E");
371 h->SetDirectory(0);
372
373 const Bool_t rc = Fit(*h, paint);
374 delete h;
375 return rc;
376}
377
378Bool_t MAlphaFitter::FitTheta(const TH3D &hon, UInt_t bin, Bool_t paint)
379{
380 const TString name(Form("TempAlphaTheta%06d", gRandom->Integer(1000000)));
381 TH1D *h = hon.ProjectionZ(name, bin, bin, -1, 9999, "E");
382 h->SetDirectory(0);
383
384 const Bool_t rc = Fit(*h, paint);
385 delete h;
386 return rc;
387}
388/*
389Bool_t MAlphaFitter::FitTime(const TH3D &hon, UInt_t bin, Bool_t paint)
390{
391 const TString name(Form("TempAlphaTime%06d", gRandom->Integer(1000000)));
392
393 hon.GetZaxis()->SetRange(bin,bin);
394 TH1D *h = (TH1D*)hon.Project3D("ye");
395 hon.GetZaxis()->SetRange(-1,9999);
396
397 h->SetDirectory(0);
398
399 const Bool_t rc = Fit(*h, paint);
400 delete h;
401 return rc;
402}
403*/
404Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, Bool_t paint)
405{
406 const TString name(Form("TempAlpha%06d", gRandom->Integer(1000000)));
407 TH1D *h = hon.ProjectionZ(name, -1, 9999, -1, 9999, "E");
408 h->SetDirectory(0);
409
410 const Bool_t rc = Fit(*h, paint);
411 delete h;
412 return rc;
413}
414
415Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
416{
417 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
418 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
419
420 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, bin, bin, "E");
421 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, bin, bin, "E");
422 h1->SetDirectory(0);
423 h0->SetDirectory(0);
424
425 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
426
427 delete h0;
428 delete h1;
429
430 return rc;
431}
432
433Bool_t MAlphaFitter::FitTheta(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
434{
435 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
436 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
437
438 TH1D *h1 = hon.ProjectionZ(name1, bin, bin, -1, 9999, "E");
439 TH1D *h0 = hof.ProjectionZ(name0, bin, bin, -1, 9999, "E");
440 h1->SetDirectory(0);
441 h0->SetDirectory(0);
442
443 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
444
445 delete h0;
446 delete h1;
447
448 return rc;
449}
450/*
451Bool_t MAlphaFitter::FitTime(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
452{
453 const TString name1(Form("TempAlphaTime%06d_on", gRandom->Integer(1000000)));
454 const TString name0(Form("TempAlphaTime%06d_off", gRandom->Integer(1000000)));
455
456 hon.GetZaxis()->SetRange(bin,bin);
457 TH1D *h1 = (TH1D*)hon.Project3D("ye");
458 hon.GetZaxis()->SetRange(-1,9999);
459 h1->SetDirectory(0);
460
461 hof.GetZaxis()->SetRange(bin,bin);
462 TH1D *h0 = (TH1D*)hof.Project3D("ye");
463 hof.GetZaxis()->SetRange(-1,9999);
464 h0->SetDirectory(0);
465
466 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
467
468 delete h0;
469 delete h1;
470
471 return rc;
472}
473*/
474Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, const TH3D &hof, Bool_t paint)
475{
476 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
477 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
478
479 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, -1, 9999, "E");
480 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, -1, 9999, "E");
481 h1->SetDirectory(0);
482 h0->SetDirectory(0);
483
484 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
485
486 delete h0;
487 delete h1;
488
489 return rc;
490}
491
492Double_t MAlphaFitter::Scale(TH1D &of, const TH1D &on) const
493{
494 Float_t scaleon = 1;
495 Float_t scaleof = 1;
496 switch (fScaleMode)
497 {
498 case kNone:
499 return 1;
500
501 case kEntries:
502 scaleon = on.GetEntries();
503 scaleof = of.GetEntries();
504 break;
505
506 case kIntegral:
507 scaleon = on.Integral();
508 scaleof = of.Integral();
509 break;
510
511 case kOffRegion:
512 {
513 const Int_t min = on.GetXaxis()->FindFixBin(fScaleMin);
514 const Int_t max = on.GetXaxis()->FindFixBin(fScaleMax);
515 scaleon = on.Integral(min, max);
516 scaleof = of.Integral(min, max);
517 }
518 break;
519
520 case kBackground:
521 {
522 const Int_t min = on.GetXaxis()->FindFixBin(fBgMin);
523 const Int_t max = on.GetXaxis()->FindFixBin(fBgMax);
524 scaleon = on.Integral(min, max);
525 scaleof = of.Integral(min, max);
526 }
527 break;
528
529 case kUserScale:
530 scaleon = fScaleUser;
531 break;
532
533 // This is just to make some compiler happy
534 default:
535 return 1;
536 }
537
538 if (scaleof!=0)
539 {
540 of.Scale(scaleon/scaleof);
541 return scaleon/scaleof;
542 }
543 else
544 {
545 of.Reset();
546 return 0;
547 }
548}
549
550Double_t MAlphaFitter::GetMinimizationValue() const
551{
552 switch (fStrategy)
553 {
554 case kSignificance:
555 return -GetSignificance();
556 case kSignificanceChi2:
557 return -GetSignificance()/GetChiSqSignal();
558 case kSignificanceLogExcess:
559 if (GetEventsExcess()<1)
560 return 0;
561 return -GetSignificance()*TMath::Log10(GetEventsExcess());
562 case kSignificanceExcess:
563 return -GetSignificance()*GetEventsExcess();
564 case kExcess:
565 return -GetEventsExcess();
566 case kGaussSigma:
567 return GetGausSigma();
568 }
569 return 0;
570}
571
572Int_t MAlphaFitter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
573{
574 Bool_t rc = kFALSE;
575
576 //void SetScaleUser(Float_t scale) { fScaleUser = scale; fScaleMode=kUserScale; }
577 //void SetScaleMode(ScaleMode_t mode) { fScaleMode = mode; }
578
579 if (IsEnvDefined(env, prefix, "SignalIntegralMax", print))
580 {
581 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalIntegralMax", fSigInt));
582 rc = kTRUE;
583 }
584 if (IsEnvDefined(env, prefix, "SignalFitMax", print))
585 {
586 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalFitMax", fSigMax));
587 rc = kTRUE;
588 }
589 if (IsEnvDefined(env, prefix, "BackgroundFitMax", print))
590 {
591 SetBackgroundFitMax(GetEnvValue(env, prefix, "BackgroundFitMax", fBgMax));
592 rc = kTRUE;
593 }
594 if (IsEnvDefined(env, prefix, "BackgroundFitMin", print))
595 {
596 SetBackgroundFitMin(GetEnvValue(env, prefix, "BackgroundFitMin", fBgMin));
597 rc = kTRUE;
598 }
599 if (IsEnvDefined(env, prefix, "ScaleMin", print))
600 {
601 SetScaleMin(GetEnvValue(env, prefix, "ScaleMin", fScaleMin));
602 rc = kTRUE;
603 }
604 if (IsEnvDefined(env, prefix, "ScaleMax", print))
605 {
606 SetScaleMax(GetEnvValue(env, prefix, "ScaleMax", fScaleMax));
607 rc = kTRUE;
608 }
609 if (IsEnvDefined(env, prefix, "PolynomOrder", print))
610 {
611 SetPolynomOrder(GetEnvValue(env, prefix, "PolynomOrder", fPolynomOrder));
612 rc = kTRUE;
613 }
614
615 if (IsEnvDefined(env, prefix, "MinimizationStrategy", print))
616 {
617 TString txt = GetEnvValue(env, prefix, "MinimizationStrategy", "");
618 txt = txt.Strip(TString::kBoth);
619 txt.ToLower();
620 if (txt==(TString)"significance")
621 fStrategy = kSignificance;
622 if (txt==(TString)"significancechi2")
623 fStrategy = kSignificanceChi2;
624 if (txt==(TString)"significanceexcess")
625 fStrategy = kSignificanceExcess;
626 if (txt==(TString)"excess")
627 fStrategy = kExcess;
628 if (txt==(TString)"gausssigma" || txt==(TString)"gaussigma")
629 fStrategy = kGaussSigma;
630 rc = kTRUE;
631 }
632 if (IsEnvDefined(env, prefix, "Scale", print))
633 {
634 fScaleUser = GetEnvValue(env, prefix, "Scale", fScaleUser);
635 rc = kTRUE;
636 }
637 if (IsEnvDefined(env, prefix, "ScaleMode", print))
638 {
639 TString txt = GetEnvValue(env, prefix, "ScaleMode", "");
640 txt = txt.Strip(TString::kBoth);
641 txt.ToLower();
642 if (txt==(TString)"none")
643 fScaleMode = kNone;
644 if (txt==(TString)"entries")
645 fScaleMode = kEntries;
646 if (txt==(TString)"integral")
647 fScaleMode = kIntegral;
648 if (txt==(TString)"offregion")
649 fScaleMode = kOffRegion;
650 if (txt==(TString)"background")
651 fScaleMode = kBackground;
652 if (txt==(TString)"leastsquare")
653 fScaleMode = kLeastSquare;
654 if (txt==(TString)"userscale")
655 fScaleMode = kUserScale;
656 if (txt==(TString)"fixed")
657 {
658 fScaleMode = kUserScale;
659 fScaleUser = fScaleFactor;
660 cout << "---------> " << fScaleFactor << " <----------" << endl;
661 }
662 rc = kTRUE;
663 }
664 if (IsEnvDefined(env, prefix, "SignalFunction", print))
665 {
666 TString txt = GetEnvValue(env, prefix, "SignalFunction", "");
667 txt = txt.Strip(TString::kBoth);
668 txt.ToLower();
669 if (txt==(TString)"gauss" || txt==(TString)"gaus")
670 SetSignalFunction(kGauss);
671 if (txt==(TString)"thetasq")
672 SetSignalFunction(kThetaSq);
673 rc = kTRUE;
674 }
675
676 return rc;
677}
Note: See TracBrowser for help on using the repository browser.