source: trunk/MagicSoft/Mars/mhflux/MAlphaFitter.cc@ 7462

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