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