| 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 | !
 | 
|---|
| 20 | !   Copyright: MAGIC Software Development, 2000-2005
 | 
|---|
| 21 | !
 | 
|---|
| 22 | !
 | 
|---|
| 23 | \* ======================================================================== */
 | 
|---|
| 24 | 
 | 
|---|
| 25 | //////////////////////////////////////////////////////////////////////////////
 | 
|---|
| 26 | //
 | 
|---|
| 27 | // MHDisp
 | 
|---|
| 28 | //
 | 
|---|
| 29 | // Create a false source plot using disp.
 | 
|---|
| 30 | //
 | 
|---|
| 31 | // Currently the use of this class requires to be after MFMagicCuts
 | 
|---|
| 32 | // in the tasklist. Switching of the M3L cut in MFMagicCuts is recommended.
 | 
|---|
| 33 | //
 | 
|---|
| 34 | // Class Version 3:
 | 
|---|
| 35 | // ----------------
 | 
|---|
| 36 | //  + Double_t fScaleMin;    // [deg] Minimum circle for integration of off-data for scaling
 | 
|---|
| 37 | //  + Double_t fScaleMax;    // [deg] Maximum circle for integration of off-data for scaling
 | 
|---|
| 38 | //
 | 
|---|
| 39 | //////////////////////////////////////////////////////////////////////////////
 | 
|---|
| 40 | #include "MHDisp.h"
 | 
|---|
| 41 | 
 | 
|---|
| 42 | #include <TStyle.h>
 | 
|---|
| 43 | #include <TCanvas.h>
 | 
|---|
| 44 | 
 | 
|---|
| 45 | #include <TF1.h>
 | 
|---|
| 46 | #include <TF2.h>
 | 
|---|
| 47 | #include <TProfile.h>
 | 
|---|
| 48 | 
 | 
|---|
| 49 | #include "MLog.h"
 | 
|---|
| 50 | #include "MLogManip.h"
 | 
|---|
| 51 | 
 | 
|---|
| 52 | #include "MMath.h"
 | 
|---|
| 53 | #include "MString.h"
 | 
|---|
| 54 | #include "MBinning.h"
 | 
|---|
| 55 | 
 | 
|---|
| 56 | #include "MParList.h"
 | 
|---|
| 57 | #include "MParameters.h"
 | 
|---|
| 58 | 
 | 
|---|
| 59 | #include "MHillas.h"
 | 
|---|
| 60 | #include "MSrcPosCam.h"
 | 
|---|
| 61 | #include "MPointingPos.h"
 | 
|---|
| 62 | #include "MPointingDev.h"
 | 
|---|
| 63 | 
 | 
|---|
| 64 | ClassImp(MHDisp);
 | 
|---|
| 65 | 
 | 
|---|
| 66 | using namespace std;
 | 
|---|
| 67 | 
 | 
|---|
| 68 | // --------------------------------------------------------------------------
 | 
|---|
| 69 | //
 | 
|---|
| 70 | // Default Constructor
 | 
|---|
| 71 | //
 | 
|---|
| 72 | MHDisp::MHDisp(const char *name, const char *title)
 | 
|---|
| 73 |     : fDisp(0), fDeviation(0), fSrcAnti(0), fHalf(kFALSE), fSmearing(-1),
 | 
|---|
| 74 |     fWobble(kFALSE), fScaleMin(0.325), fScaleMax(0.475)
 | 
|---|
| 75 | {
 | 
|---|
| 76 |     //
 | 
|---|
| 77 |     //   set the name and title of this object
 | 
|---|
| 78 |     //
 | 
|---|
| 79 |     fName  = name  ? name  : "MHDisp";
 | 
|---|
| 80 |     fTitle = title ? title : "3D-plot using Disp vs x, y";
 | 
|---|
| 81 | 
 | 
|---|
| 82 |     fHist.SetName("Alpha");
 | 
|---|
| 83 |     fHist.SetTitle("3D-plot of ThetaSq vs x, y");
 | 
|---|
| 84 |     fHist.SetXTitle("dx [\\circ]");
 | 
|---|
| 85 |     fHist.SetYTitle("dy [\\circ]");
 | 
|---|
| 86 |     fHist.SetZTitle("Eq.cts");
 | 
|---|
| 87 | 
 | 
|---|
| 88 |     fHist.SetDirectory(NULL);
 | 
|---|
| 89 |     fHistBg.SetDirectory(NULL);
 | 
|---|
| 90 |     fHistBg1.SetDirectory(NULL);
 | 
|---|
| 91 |     fHistBg2.SetDirectory(NULL);
 | 
|---|
| 92 | 
 | 
|---|
| 93 |     fHist.SetBit(TH1::kNoStats);
 | 
|---|
| 94 | }
 | 
|---|
| 95 | 
 | 
|---|
| 96 | // --------------------------------------------------------------------------
 | 
|---|
| 97 | //
 | 
|---|
| 98 | // Set binnings (takes BinningFalseSource) and prepare filling of the
 | 
|---|
| 99 | // histogram.
 | 
|---|
| 100 | //
 | 
|---|
| 101 | // Also search for MTime, MObservatory and MPointingPos
 | 
|---|
| 102 | //
 | 
|---|
| 103 | Bool_t MHDisp::SetupFill(const MParList *plist)
 | 
|---|
| 104 | {
 | 
|---|
| 105 |     if (!MHFalseSource::SetupFill(plist))
 | 
|---|
| 106 |         return kFALSE;
 | 
|---|
| 107 | 
 | 
|---|
| 108 |     fDisp = (MParameterD*)plist->FindObject("Disp", "MParameterD");
 | 
|---|
| 109 |     if (!fDisp)
 | 
|---|
| 110 |     {
 | 
|---|
| 111 |         *fLog << err << "Disp [MParameterD] not found... abort." << endl;
 | 
|---|
| 112 |         return kFALSE;
 | 
|---|
| 113 |     }
 | 
|---|
| 114 | 
 | 
|---|
| 115 |     if (fWobble)
 | 
|---|
| 116 |     {
 | 
|---|
| 117 |         fSrcAnti = (MSrcPosCam*)plist->FindObject("MSrcPosAnti", "MSrcPosCam");
 | 
|---|
| 118 |         if (!fSrcAnti)
 | 
|---|
| 119 |         {
 | 
|---|
| 120 |             *fLog << err << "MSrcPosAnti [MSrcPosCam] not found... abort." << endl;
 | 
|---|
| 121 |             return kFALSE;
 | 
|---|
| 122 |         }
 | 
|---|
| 123 | 
 | 
|---|
| 124 |         *fLog << inf << "Wobble mode initialized. " << endl;
 | 
|---|
| 125 |     }
 | 
|---|
| 126 | 
 | 
|---|
| 127 |     fDeviation = (MPointingDev*)plist->FindObject("MPointingDev");
 | 
|---|
| 128 |     if (!fDeviation)
 | 
|---|
| 129 |         *fLog << warn << "MPointingDev not found... ignored." << endl;
 | 
|---|
| 130 | 
 | 
|---|
| 131 |     MBinning binsx, binsy;
 | 
|---|
| 132 |     binsx.SetEdges(fHist, 'x');
 | 
|---|
| 133 |     binsy.SetEdges(fHist, 'y');
 | 
|---|
| 134 |     MH::SetBinning(fHistBg, binsx, binsy);
 | 
|---|
| 135 | 
 | 
|---|
| 136 |     if (!fHistOff)
 | 
|---|
| 137 |     {
 | 
|---|
| 138 |         MH::SetBinning(fHistBg1, binsx, binsy);
 | 
|---|
| 139 |         MH::SetBinning(fHistBg2, binsx, binsy);
 | 
|---|
| 140 |     }
 | 
|---|
| 141 | 
 | 
|---|
| 142 |     return kTRUE;
 | 
|---|
| 143 | }
 | 
|---|
| 144 | 
 | 
|---|
| 145 | // --------------------------------------------------------------------------
 | 
|---|
| 146 | //
 | 
|---|
| 147 | // Fill the histogram. For details see the code or the class description
 | 
|---|
| 148 | // 
 | 
|---|
| 149 | Int_t MHDisp::Fill(const MParContainer *par, const Stat_t w)
 | 
|---|
| 150 | {
 | 
|---|
| 151 |     const MHillas *hil = dynamic_cast<const MHillas*>(par);
 | 
|---|
| 152 |     if (!hil)
 | 
|---|
| 153 |     {
 | 
|---|
| 154 |         *fLog << err << "MHDisp::Fill: No container specified!" << endl;
 | 
|---|
| 155 |         return kERROR;
 | 
|---|
| 156 |     }
 | 
|---|
| 157 | 
 | 
|---|
| 158 |     // Get camera rotation angle
 | 
|---|
| 159 |     Double_t rho = 0;
 | 
|---|
| 160 |     if (fTime && fObservatory && fPointPos)
 | 
|---|
| 161 |         rho = fPointPos->RotationAngle(*fObservatory, *fTime, fDeviation);
 | 
|---|
| 162 | 
 | 
|---|
| 163 |     // Vector of length disp in direction of shower
 | 
|---|
| 164 |     // Move origin of vector to center-of-gravity of shower and derotate
 | 
|---|
| 165 |     TVector2 pos1 = hil->GetMean()*fMm2Deg + hil->GetNormAxis()*fDisp->GetVal();
 | 
|---|
| 166 | 
 | 
|---|
| 167 |     const TVector2 src = fSrcPos->GetXY()*fMm2Deg;
 | 
|---|
| 168 | 
 | 
|---|
| 169 |     Double_t w0 = 1;
 | 
|---|
| 170 |     if (fWobble)
 | 
|---|
| 171 |     {
 | 
|---|
| 172 |         const TVector2 anti = fSrcAnti->GetXY()*fMm2Deg;
 | 
|---|
| 173 | 
 | 
|---|
| 174 |         // Skip off-data not in the same half than the source (here: anti-source)
 | 
|---|
| 175 |         // Could be replaced by some kind of theta cut...
 | 
|---|
| 176 |         if (!fHistOff)
 | 
|---|
| 177 |         {
 | 
|---|
| 178 |             Double_t r = anti.Mod()>0.2*1.7 ? 0.2*1.7 : anti.Mod();
 | 
|---|
| 179 | 
 | 
|---|
| 180 |             // In wobble mode processing the off-data, the anti-source
 | 
|---|
| 181 |             // position is our source position. Check if this is a possible
 | 
|---|
| 182 |             // gamma. If it is, do not fill it into our off-data histogram
 | 
|---|
| 183 |             if ((pos1-anti).Mod()<r)
 | 
|---|
| 184 |                 return kTRUE;
 | 
|---|
| 185 | 
 | 
|---|
| 186 |             // Because afterwards the plots are normalized with the number
 | 
|---|
| 187 |             // of entries the  non-overlap  region corresponds to half the
 | 
|---|
| 188 |             // contents, the overlap region  to  full contents.  By taking
 | 
|---|
| 189 |             // the mean of both distributions  we get the same result than
 | 
|---|
| 190 |             // we would have gotten using full  off-events with a slightly
 | 
|---|
| 191 |             // increased uncertainty
 | 
|---|
| 192 |             // FIXME: The delta stuff could be replaced by a 2*antitheta cut...
 | 
|---|
| 193 |             //w0 = delta>25 ? 1 : 2;
 | 
|---|
| 194 | 
 | 
|---|
| 195 |             w0 = (pos1+anti).Mod()<r ? 2 : 1;
 | 
|---|
| 196 |         }
 | 
|---|
| 197 | 
 | 
|---|
| 198 |         // When processing off-data the anti-source is the real source
 | 
|---|
| 199 |         const TVector2 srcpos = fHistOff ? src : anti;
 | 
|---|
| 200 | 
 | 
|---|
| 201 |         // Define by the source position which histogram to fill
 | 
|---|
| 202 |         if (TMath::Abs(srcpos.DeltaPhi(fFormerSrc))*TMath::RadToDeg()>90)
 | 
|---|
| 203 |             fHalf = !fHalf;
 | 
|---|
| 204 |         fFormerSrc = srcpos;
 | 
|---|
| 205 |     }
 | 
|---|
| 206 | 
 | 
|---|
| 207 |     // If on-data: Derotate source and P. Translate P to center.
 | 
|---|
| 208 |     TVector2 m;
 | 
|---|
| 209 |     if (fHistOff)
 | 
|---|
| 210 |     {
 | 
|---|
| 211 |         // Derotate all position around camera center by -rho
 | 
|---|
| 212 |         // Move origin of vector to center-of-gravity of shower and derotate
 | 
|---|
| 213 |         pos1=pos1.Rotate(-rho);
 | 
|---|
| 214 | 
 | 
|---|
| 215 |         // Shift the source position to 0/0
 | 
|---|
| 216 |         if (fSrcPos)
 | 
|---|
| 217 |         {
 | 
|---|
| 218 |             // m: Position of the camera center in the FS plot
 | 
|---|
| 219 |             m = src.Rotate(-rho);
 | 
|---|
| 220 |             pos1 -= m;
 | 
|---|
| 221 |         }
 | 
|---|
| 222 |     }
 | 
|---|
| 223 | 
 | 
|---|
| 224 |     // -------------------------------------------------
 | 
|---|
| 225 |     //  The following algorithm may look complicated...
 | 
|---|
| 226 |     //            It is optimized for speed
 | 
|---|
| 227 |     // -------------------------------------------------
 | 
|---|
| 228 | 
 | 
|---|
| 229 |     // Define a vector used to calculate rotated x and y component
 | 
|---|
| 230 |     const TVector2 rot(TMath::Sin(rho), TMath::Cos(rho));
 | 
|---|
| 231 | 
 | 
|---|
| 232 |     // Fold event position with the PSF and fill it into histogram
 | 
|---|
| 233 |     const TAxis &axex = *fHist.GetXaxis();
 | 
|---|
| 234 |     const TAxis &axey = *fHist.GetYaxis();
 | 
|---|
| 235 | 
 | 
|---|
| 236 |     const Int_t nx = axex.GetNbins();
 | 
|---|
| 237 |     const Int_t ny = axey.GetNbins();
 | 
|---|
| 238 | 
 | 
|---|
| 239 |     // normg: Normalization for Gauss [sqrt(2*pi)*sigma]^2
 | 
|---|
| 240 |     const Double_t weight = axex.GetBinWidth(1)*axey.GetBinWidth(1)*w*w0;
 | 
|---|
| 241 |     const Double_t psf    = 2*fSmearing*fSmearing;
 | 
|---|
| 242 |     const Double_t pi2    = fSmearing * TMath::Pi()/2;
 | 
|---|
| 243 |     const Double_t normg  = pi2*pi2 * TMath::Sqrt(TMath::TwoPi()) / weight;
 | 
|---|
| 244 |     const Int_t    bz     = fHist.GetZaxis()->FindFixBin(0);
 | 
|---|
| 245 | 
 | 
|---|
| 246 |     TH2 &hbg = fHalf ? fHistBg1 : fHistBg2;
 | 
|---|
| 247 | 
 | 
|---|
| 248 |     const Bool_t smear = fSmearing>0;
 | 
|---|
| 249 |     if (!smear)
 | 
|---|
| 250 |     {
 | 
|---|
| 251 |         if (!fHistOff)
 | 
|---|
| 252 |             hbg.Fill(pos1.X(), pos1.Y(), w*w0);
 | 
|---|
| 253 |         else
 | 
|---|
| 254 |             fHist.Fill(pos1.X(), pos1.Y(), 0.0, w*w0);
 | 
|---|
| 255 |     }
 | 
|---|
| 256 | 
 | 
|---|
| 257 |     // To calculate significance map smear with 2*theta-cut and
 | 
|---|
| 258 |     // do not normalize gauss, for event map use theta-cut/2 instead
 | 
|---|
| 259 |     if (smear || fHistOff)
 | 
|---|
| 260 |     {
 | 
|---|
| 261 |         for (int x=1; x<=nx; x++)
 | 
|---|
| 262 |         {
 | 
|---|
| 263 |             const Double_t cx = axex.GetBinCenter(x);
 | 
|---|
| 264 |             const Double_t px = cx-pos1.X();
 | 
|---|
| 265 | 
 | 
|---|
| 266 |             for (int y=1; y<=ny; y++)
 | 
|---|
| 267 |             {
 | 
|---|
| 268 |                 const Double_t cy = axey.GetBinCenter(y);
 | 
|---|
| 269 |                 const Double_t sp = Sq(px, cy-pos1.Y());
 | 
|---|
| 270 |                 if (smear)
 | 
|---|
| 271 |                 {
 | 
|---|
| 272 |                     const Double_t dp = sp/psf;
 | 
|---|
| 273 | 
 | 
|---|
| 274 |                     // Values below 1e-3 (>3.5sigma) are not filled into the histogram
 | 
|---|
| 275 |                     if (dp<4)
 | 
|---|
| 276 |                     {
 | 
|---|
| 277 |                         const Double_t rc = TMath::Exp(-dp)/normg;
 | 
|---|
| 278 |                         if (!fHistOff)
 | 
|---|
| 279 |                             hbg.AddBinContent(hbg.GetBin(x, y), rc);
 | 
|---|
| 280 |                         else
 | 
|---|
| 281 |                             fHist.AddBinContent(fHist.GetBin(x, y, bz), rc);
 | 
|---|
| 282 |                     }
 | 
|---|
| 283 |                 }
 | 
|---|
| 284 | 
 | 
|---|
| 285 |                 if (!fHistOff)
 | 
|---|
| 286 |                     continue;
 | 
|---|
| 287 | 
 | 
|---|
| 288 |                 // If we are filling the signal already (fHistOff!=NULL)
 | 
|---|
| 289 |                 // we also fill the background by projecting the
 | 
|---|
| 290 |                 // background in the camera into the sky plot.
 | 
|---|
| 291 |                 const TVector2 v = TVector2(cx+m.X(), cy+m.Y());
 | 
|---|
| 292 | 
 | 
|---|
| 293 |                 // Speed up further: Xmax-fWobble
 | 
|---|
| 294 |                 if (v.Mod()>axex.GetXmax()) // Gains 10% speed
 | 
|---|
| 295 |                     continue;
 | 
|---|
| 296 | 
 | 
|---|
| 297 |                 const Int_t    bx = axex.FindFixBin(v^rot);
 | 
|---|
| 298 |                 const Int_t    by = axey.FindFixBin(v*rot);
 | 
|---|
| 299 |                 const Double_t bg = fHistOff->GetBinContent(bx, by, bz);
 | 
|---|
| 300 | 
 | 
|---|
| 301 |                 fHistBg.AddBinContent(fHistBg.GetBin(x, y), bg);
 | 
|---|
| 302 |             }
 | 
|---|
| 303 |         }
 | 
|---|
| 304 |     }
 | 
|---|
| 305 | 
 | 
|---|
| 306 |     if (fHistOff)
 | 
|---|
| 307 |         fHistBg.SetEntries(fHistBg.GetEntries()+1);
 | 
|---|
| 308 | 
 | 
|---|
| 309 |     if (!smear)
 | 
|---|
| 310 |         return kTRUE;
 | 
|---|
| 311 | 
 | 
|---|
| 312 |     if (!fHistOff)
 | 
|---|
| 313 |         hbg.SetEntries(hbg.GetEntries()+1);
 | 
|---|
| 314 |     else
 | 
|---|
| 315 |         fHist.SetEntries(fHist.GetEntries()+1);
 | 
|---|
| 316 | 
 | 
|---|
| 317 |     return kTRUE;
 | 
|---|
| 318 | }
 | 
|---|
| 319 | 
 | 
|---|
| 320 | // --------------------------------------------------------------------------
 | 
|---|
| 321 | //
 | 
|---|
| 322 | // Compile the background in camera coordinates from the two background
 | 
|---|
| 323 | // histograms
 | 
|---|
| 324 | //
 | 
|---|
| 325 | Bool_t MHDisp::Finalize()
 | 
|---|
| 326 | {
 | 
|---|
| 327 |     if (fHistOff)
 | 
|---|
| 328 |         return kTRUE;
 | 
|---|
| 329 | 
 | 
|---|
| 330 |     const Int_t bz = fHist.GetZaxis()->FindFixBin(0);
 | 
|---|
| 331 | 
 | 
|---|
| 332 |     const Double_t n1 = fHistBg1.GetEntries();
 | 
|---|
| 333 |     const Double_t n2 = fHistBg2.GetEntries();
 | 
|---|
| 334 | 
 | 
|---|
| 335 |     for (int x=1; x<=fHist.GetNbinsX(); x++)
 | 
|---|
| 336 |         for (int y=1; y<=fHist.GetNbinsY(); y++)
 | 
|---|
| 337 |         {
 | 
|---|
| 338 |             const Int_t bin1 = fHistBg1.GetBin(x, y);
 | 
|---|
| 339 | 
 | 
|---|
| 340 |             const Double_t rc =
 | 
|---|
| 341 |                 (n1==0?0:fHistBg1.GetBinContent(bin1)/n1)+
 | 
|---|
| 342 |                 (n2==0?0:fHistBg2.GetBinContent(bin1)/n2);
 | 
|---|
| 343 | 
 | 
|---|
| 344 |             fHist.SetBinContent(x, y, bz, rc/2);
 | 
|---|
| 345 |         }
 | 
|---|
| 346 | 
 | 
|---|
| 347 |     fHist.SetEntries(n1+n2);
 | 
|---|
| 348 | 
 | 
|---|
| 349 |     // Result corresponds to one smeared background event
 | 
|---|
| 350 | 
 | 
|---|
| 351 |     return kTRUE;
 | 
|---|
| 352 | }
 | 
|---|
| 353 | 
 | 
|---|
| 354 | // --------------------------------------------------------------------------
 | 
|---|
| 355 | //
 | 
|---|
| 356 | // Make sure that if the scale is changed by the context menu all subpads
 | 
|---|
| 357 | // are redrawn.
 | 
|---|
| 358 | //
 | 
|---|
| 359 | void MHDisp::Update()
 | 
|---|
| 360 | {
 | 
|---|
| 361 |     TVirtualPad *pad = gPad;
 | 
|---|
| 362 |     for (int i=1; i<=6; i++)
 | 
|---|
| 363 |     {
 | 
|---|
| 364 |         if (pad->GetPad(i))
 | 
|---|
| 365 |             pad->GetPad(i)->Modified();
 | 
|---|
| 366 |     }
 | 
|---|
| 367 | }
 | 
|---|
| 368 | 
 | 
|---|
| 369 | // --------------------------------------------------------------------------
 | 
|---|
| 370 | //
 | 
|---|
| 371 | // Return the mean signal in h around (0,0) in a distance between
 | 
|---|
| 372 | // 0.325 and 0.475deg
 | 
|---|
| 373 | //
 | 
|---|
| 374 | Double_t MHDisp::GetOffSignal(TH1 &h) const
 | 
|---|
| 375 | {
 | 
|---|
| 376 |     const TAxis &axex = *h.GetXaxis();
 | 
|---|
| 377 |     const TAxis &axey = *h.GetYaxis();
 | 
|---|
| 378 | 
 | 
|---|
| 379 |     Int_t    cnt = 0;
 | 
|---|
| 380 |     Double_t sum = 0;
 | 
|---|
| 381 |     for (int x=0; x<h.GetNbinsX(); x++)
 | 
|---|
| 382 |         for (int y=0; y<h.GetNbinsY(); y++)
 | 
|---|
| 383 |         {
 | 
|---|
| 384 |             const Double_t d = TMath::Hypot(axex.GetBinCenter(x+1), axey.GetBinCenter(y+1));
 | 
|---|
| 385 |             if (d>fScaleMin && d<fScaleMax)
 | 
|---|
| 386 |             {
 | 
|---|
| 387 |                 sum += h.GetBinContent(x+1,y+1);
 | 
|---|
| 388 |                 cnt++;
 | 
|---|
| 389 |             }
 | 
|---|
| 390 |         }
 | 
|---|
| 391 | 
 | 
|---|
| 392 |     return sum/cnt;
 | 
|---|
| 393 | }
 | 
|---|
| 394 | 
 | 
|---|
| 395 | // --------------------------------------------------------------------------
 | 
|---|
| 396 | //
 | 
|---|
| 397 | // Fill h2 with the radial profile of h1 around (x0, y0)
 | 
|---|
| 398 | //
 | 
|---|
| 399 | void MHDisp::Profile(TH1 &h2, const TH2 &h1, Axis_t x0, Axis_t y0) const
 | 
|---|
| 400 | {
 | 
|---|
| 401 |     const TAxis &axex = *h1.GetXaxis();
 | 
|---|
| 402 |     const TAxis &axey = *h1.GetYaxis();
 | 
|---|
| 403 | 
 | 
|---|
| 404 |     h2.Reset();
 | 
|---|
| 405 | 
 | 
|---|
| 406 |     for (int x=1; x<=axex.GetNbins(); x++)
 | 
|---|
| 407 |         for (int y=1; y<=axey.GetNbins(); y++)
 | 
|---|
| 408 |         {
 | 
|---|
| 409 |             const Double_t dx = axex.GetBinCenter(x)-x0;
 | 
|---|
| 410 |             const Double_t dy = axey.GetBinCenter(y)-y0;
 | 
|---|
| 411 | 
 | 
|---|
| 412 |             const Double_t r  = TMath::Hypot(dx, dy);
 | 
|---|
| 413 | 
 | 
|---|
| 414 |             h2.Fill(r, h1.GetBinContent(x, y));
 | 
|---|
| 415 |         }
 | 
|---|
| 416 | }
 | 
|---|
| 417 | 
 | 
|---|
| 418 | // --------------------------------------------------------------------------
 | 
|---|
| 419 | //
 | 
|---|
| 420 | // Remove contents of histogram around a circle.
 | 
|---|
| 421 | //
 | 
|---|
| 422 | void MHDisp::MakeDot(TH2 &h2) const
 | 
|---|
| 423 | {
 | 
|---|
| 424 |     const TAxis &axex = *h2.GetXaxis();
 | 
|---|
| 425 |     const TAxis &axey = *h2.GetYaxis();
 | 
|---|
| 426 | 
 | 
|---|
| 427 |     const Double_t rmax = (fWobble ? axex.GetXmax()/*-0.7*/ : axex.GetXmax()) - axex.GetBinWidth(1);
 | 
|---|
| 428 | 
 | 
|---|
| 429 |     for (int x=1; x<=axex.GetNbins(); x++)
 | 
|---|
| 430 |         for (int y=1; y<=axey.GetNbins(); y++)
 | 
|---|
| 431 |         {
 | 
|---|
| 432 |             const Int_t bin = h2.GetBin(x,y);
 | 
|---|
| 433 | 
 | 
|---|
| 434 |             const Double_t px = h2.GetBinCenter(x);
 | 
|---|
| 435 |             const Double_t py = h2.GetBinCenter(y);
 | 
|---|
| 436 | 
 | 
|---|
| 437 |             if (rmax<TMath::Hypot(px, py))
 | 
|---|
| 438 |                 h2.SetBinContent(bin, 0);
 | 
|---|
| 439 |             //else
 | 
|---|
| 440 |             //    if (h2.GetBinContent(bin)==0)
 | 
|---|
| 441 |             //        h2.SetBinContent(bin, 1e-10);
 | 
|---|
| 442 |         }
 | 
|---|
| 443 | }
 | 
|---|
| 444 | 
 | 
|---|
| 445 | // --------------------------------------------------------------------------
 | 
|---|
| 446 | //
 | 
|---|
| 447 | // Calculate from signal and background the significance map
 | 
|---|
| 448 | //
 | 
|---|
| 449 | void MHDisp::MakeSignificance(TH2 &s, const TH2 &h1, const TH2 &h2, const Double_t scale) const
 | 
|---|
| 450 | {
 | 
|---|
| 451 |     const TAxis &axex = *s.GetXaxis();
 | 
|---|
| 452 |     const TAxis &axey = *s.GetYaxis();
 | 
|---|
| 453 | 
 | 
|---|
| 454 |     const Int_t nx = axex.GetNbins();
 | 
|---|
| 455 |     const Int_t ny = axey.GetNbins();
 | 
|---|
| 456 | 
 | 
|---|
| 457 |     const Int_t n = TMath::Nint(0.2/axex.GetBinWidth(1));
 | 
|---|
| 458 | 
 | 
|---|
| 459 |     const Double_t sc = h2.GetEntries()/fHistOff->GetEntries();
 | 
|---|
| 460 | 
 | 
|---|
| 461 |     for (int x=1; x<=nx; x++)
 | 
|---|
| 462 |         for (int y=1; y<=ny; y++)
 | 
|---|
| 463 |         {
 | 
|---|
| 464 |             Double_t S=0;
 | 
|---|
| 465 | 
 | 
|---|
| 466 |             // Only calculate significances for pixels far enough
 | 
|---|
| 467 |             // from the border to get all expected pixels.
 | 
|---|
| 468 |             if (TMath::Hypot((float)x-0.5*nx, (float)y-0.5*ny)<0.5*axex.GetNbins()-n)
 | 
|---|
| 469 |             {
 | 
|---|
| 470 |                 Double_t sig=0;
 | 
|---|
| 471 |                 Double_t bg =0;
 | 
|---|
| 472 | 
 | 
|---|
| 473 |                 // Integral a region of n pixels around x/y
 | 
|---|
| 474 |                 for (int dx=-n; dx<=n; dx++)
 | 
|---|
| 475 |                     for (int dy=-n; dy<=n; dy++)
 | 
|---|
| 476 |                     {
 | 
|---|
| 477 |                         if (TMath::Hypot((float)dx, (float)dy)>n)
 | 
|---|
| 478 |                             continue;
 | 
|---|
| 479 | 
 | 
|---|
| 480 |                         const Int_t  bin = s.GetBin(x+dx,y+dy);
 | 
|---|
| 481 | 
 | 
|---|
| 482 |                         sig += h1.GetBinContent(bin);
 | 
|---|
| 483 |                         bg  += h2.GetBinContent(bin);
 | 
|---|
| 484 |                     }
 | 
|---|
| 485 | 
 | 
|---|
| 486 |                 // Scale such, that the statistical error corresponds to
 | 
|---|
| 487 |                 // the amount of off-data used to determin the background
 | 
|---|
| 488 |                 // model, not to the background itself. This gives
 | 
|---|
| 489 |                 // significances as calculated from the theta-sq plot.
 | 
|---|
| 490 |                 S = sig>0 ? MMath::SignificanceLiMaSigned(sig, bg*scale/sc, sc) : 0;
 | 
|---|
| 491 |             }
 | 
|---|
| 492 | 
 | 
|---|
| 493 |             s.SetBinContent(x, y, S);
 | 
|---|
| 494 |         }
 | 
|---|
| 495 | }
 | 
|---|
| 496 | 
 | 
|---|
| 497 | void MHDisp::Profile1D(const char *name, const TH2 &h) const
 | 
|---|
| 498 | {
 | 
|---|
| 499 |     if (!gPad)
 | 
|---|
| 500 |         return;
 | 
|---|
| 501 | 
 | 
|---|
| 502 |     TH1D *hrc = dynamic_cast<TH1D*>(gPad->FindObject(name));
 | 
|---|
| 503 |     if (!hrc)
 | 
|---|
| 504 |         return;
 | 
|---|
| 505 | 
 | 
|---|
| 506 |     hrc->Reset();
 | 
|---|
| 507 | 
 | 
|---|
| 508 |     //const_cast<TH2&>(h).SetMaximum();
 | 
|---|
| 509 |     const Double_t max = h.GetMaximum();
 | 
|---|
| 510 | 
 | 
|---|
| 511 |     MBinning(51, -max*1.1, max*1.1).Apply(*hrc);
 | 
|---|
| 512 | 
 | 
|---|
| 513 |     for (int x=1; x<=h.GetXaxis()->GetNbins(); x++)
 | 
|---|
| 514 |         for (int y=1; y<=h.GetXaxis()->GetNbins(); y++)
 | 
|---|
| 515 |         {
 | 
|---|
| 516 |             const Int_t   bin  = h.GetBin(x,y);
 | 
|---|
| 517 |             const Double_t sig = h.GetBinContent(bin);
 | 
|---|
| 518 |             if (sig!=0)
 | 
|---|
| 519 |                 hrc->Fill(sig);
 | 
|---|
| 520 |         }
 | 
|---|
| 521 | 
 | 
|---|
| 522 |     gPad->SetLogy(hrc->GetMaximum()>0);
 | 
|---|
| 523 | 
 | 
|---|
| 524 |     if (!fHistOff || hrc->GetRMS()<0.1)
 | 
|---|
| 525 |         return;
 | 
|---|
| 526 | 
 | 
|---|
| 527 |     // ---------- Fix mean ----------
 | 
|---|
| 528 |     // TF1 *g = (TF1*)gROOT->GetFunction("gaus");
 | 
|---|
| 529 |     // g->FixParameter(1, 0);
 | 
|---|
| 530 |     // hrc->Fit("gaus", "BQI");
 | 
|---|
| 531 | 
 | 
|---|
| 532 |     hrc->Fit("gaus", "QI");
 | 
|---|
| 533 | 
 | 
|---|
| 534 |     TF1 *f = hrc->GetFunction("gaus");
 | 
|---|
| 535 |     if (f)
 | 
|---|
| 536 |     {
 | 
|---|
| 537 |         f->SetLineWidth(1);
 | 
|---|
| 538 |         f->SetLineColor(kBlue);
 | 
|---|
| 539 |     }
 | 
|---|
| 540 | }
 | 
|---|
| 541 | 
 | 
|---|
| 542 | void MHDisp::Paint(Option_t *o)
 | 
|---|
| 543 | {
 | 
|---|
| 544 |     // Compile Background if necessary
 | 
|---|
| 545 |     Finalize();
 | 
|---|
| 546 | 
 | 
|---|
| 547 |     // Paint result
 | 
|---|
| 548 |     TVirtualPad *pad = gPad;
 | 
|---|
| 549 | 
 | 
|---|
| 550 |     pad->cd(1);
 | 
|---|
| 551 | 
 | 
|---|
| 552 |     // Project on data onto yx-plane
 | 
|---|
| 553 |     fHist.GetZaxis()->SetRange(0,9999);
 | 
|---|
| 554 |     TH2 *h1=(TH2*)fHist.Project3D("yx_on");
 | 
|---|
| 555 | 
 | 
|---|
| 556 |     // Set Glow-palette (PaintSurface doesn't allow more than 99 colors)
 | 
|---|
| 557 |     MH::SetPalette(fHistOff?"glowsym":"glow1", 99);
 | 
|---|
| 558 |     h1->SetContour(99);
 | 
|---|
| 559 | 
 | 
|---|
| 560 |     TH2 *hx=0;
 | 
|---|
| 561 | 
 | 
|---|
| 562 |     Double_t scale = 1;
 | 
|---|
| 563 |     if (fHistOff)
 | 
|---|
| 564 |     {
 | 
|---|
| 565 |         // Project off data onto yx-plane and subtract it from on-data
 | 
|---|
| 566 |         fHistOff->GetZaxis()->SetRange(0,9999);
 | 
|---|
| 567 |         TH2 *h=(TH2*)fHistOff->Project3D("yx_off");
 | 
|---|
| 568 | 
 | 
|---|
| 569 |         const Double_t h1off = GetOffSignal(*h1);
 | 
|---|
| 570 |         const Double_t hoff  = GetOffSignal(fHistBg);
 | 
|---|
| 571 | 
 | 
|---|
| 572 |         scale = hoff==0 ? -1 : -h1off/hoff;
 | 
|---|
| 573 | 
 | 
|---|
| 574 |         hx = (TH2*)pad->GetPad(4)->FindObject("Alpha_yx_sig");
 | 
|---|
| 575 |         if (hx)
 | 
|---|
| 576 |         {
 | 
|---|
| 577 |             hx->SetContour(99);
 | 
|---|
| 578 |             MakeSignificance(*hx, *h1, fHistBg, TMath::Abs(scale));
 | 
|---|
| 579 |             MakeDot(*hx);
 | 
|---|
| 580 |             MakeSymmetric(hx);
 | 
|---|
| 581 |         }
 | 
|---|
| 582 | 
 | 
|---|
| 583 |         h1->Add(&fHistBg, scale);
 | 
|---|
| 584 | 
 | 
|---|
| 585 |         MakeDot(*h1);
 | 
|---|
| 586 |         MakeSymmetric(h1);
 | 
|---|
| 587 | 
 | 
|---|
| 588 |         delete h;
 | 
|---|
| 589 |     }
 | 
|---|
| 590 | 
 | 
|---|
| 591 |     pad->cd(3);
 | 
|---|
| 592 |     TH1 *h2 = (TH1*)gPad->FindObject("RadProf");
 | 
|---|
| 593 | 
 | 
|---|
| 594 |     TString opt(o);
 | 
|---|
| 595 |     opt.ToLower();
 | 
|---|
| 596 | 
 | 
|---|
| 597 |     if (h1 && h2)
 | 
|---|
| 598 |     {
 | 
|---|
| 599 |         Int_t ix, iy, iz;
 | 
|---|
| 600 |         h1->GetMaximumBin(ix, iy, iz);
 | 
|---|
| 601 | 
 | 
|---|
| 602 |         const Double_t x0 = h1->GetXaxis()->GetBinCenter(ix);
 | 
|---|
| 603 |         const Double_t y0 = h1->GetYaxis()->GetBinCenter(iy);
 | 
|---|
| 604 |         //const Double_t w0 = h1->GetXaxis()->GetBinWidth(1);
 | 
|---|
| 605 | 
 | 
|---|
| 606 |         Profile(*h2, *h1, 0, 0);
 | 
|---|
| 607 |         //Profile(*h2, *h1, x0, y0);
 | 
|---|
| 608 | 
 | 
|---|
| 609 |         if (h2->GetEntries()>0)
 | 
|---|
| 610 |         {
 | 
|---|
| 611 |             // Replace with MAlphaFitter?
 | 
|---|
| 612 |             TF1 func("fcn", "gaus + [3]*x*x + [4]");
 | 
|---|
| 613 |             func.SetLineWidth(1);
 | 
|---|
| 614 |             func.SetLineColor(kBlue);
 | 
|---|
| 615 | 
 | 
|---|
| 616 |             func.SetParLimits(2, h2->GetBinWidth(1), 1.0);
 | 
|---|
| 617 | 
 | 
|---|
| 618 |             func.SetParameter(0, h2->GetBinContent(1));
 | 
|---|
| 619 |             func.FixParameter(1, 0);
 | 
|---|
| 620 |             func.SetParameter(2, 0.12);
 | 
|---|
| 621 |             if (fHistOff)
 | 
|---|
| 622 |                 func.FixParameter(3, 0);
 | 
|---|
| 623 |             func.SetParameter(4, 0);//h2->GetBinContent(20));
 | 
|---|
| 624 |             h2->Fit(&func, "IQ", "", 0, 1.0);
 | 
|---|
| 625 | 
 | 
|---|
| 626 |             h2->SetTitle(MString::Format("P=(%.2f\\circ/%.2f\\circ) \\omega=%.2f\\circ f=%.2f",
 | 
|---|
| 627 |                                          x0, y0, func.GetParameter(2), TMath::Abs(scale)));
 | 
|---|
| 628 |         }
 | 
|---|
| 629 |     }
 | 
|---|
| 630 | 
 | 
|---|
| 631 |     pad->cd(5);
 | 
|---|
| 632 |     if (h1)
 | 
|---|
| 633 |         Profile1D("Exc", *h1);
 | 
|---|
| 634 | 
 | 
|---|
| 635 |     pad->cd(6);
 | 
|---|
| 636 |     if (hx)
 | 
|---|
| 637 |         Profile1D("Sig", *hx);
 | 
|---|
| 638 | }
 | 
|---|
| 639 | 
 | 
|---|
| 640 | void MHDisp::Draw(Option_t *o)
 | 
|---|
| 641 | {
 | 
|---|
| 642 |     TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
 | 
|---|
| 643 |     const Int_t col = pad->GetFillColor();
 | 
|---|
| 644 |     pad->SetBorderMode(0);
 | 
|---|
| 645 | 
 | 
|---|
| 646 |     AppendPad(o);
 | 
|---|
| 647 | 
 | 
|---|
| 648 |     // ----- Pad number 4 -----
 | 
|---|
| 649 |     TString name = MString::Format("%s_4", pad->GetName());
 | 
|---|
| 650 |     TPad *p = new TPad(name,name, 0.525/*0.5025*/, 0.3355, 0.995, 0.995, col, 0, 0);
 | 
|---|
| 651 |     p->SetNumber(4);
 | 
|---|
| 652 |     p->Draw();
 | 
|---|
| 653 |     p->cd();
 | 
|---|
| 654 | 
 | 
|---|
| 655 |     TH1 *h3 = fHist.Project3D("yx_sig");
 | 
|---|
| 656 |     h3->SetTitle("Significance Map");
 | 
|---|
| 657 |     h3->SetDirectory(NULL);
 | 
|---|
| 658 |     h3->SetXTitle(fHist.GetXaxis()->GetTitle());
 | 
|---|
| 659 |     h3->SetYTitle(fHist.GetYaxis()->GetTitle());
 | 
|---|
| 660 |     h3->SetMinimum(0);
 | 
|---|
| 661 |     h3->Draw("colz");
 | 
|---|
| 662 |     h3->SetBit(kCanDelete);
 | 
|---|
| 663 | 
 | 
|---|
| 664 |     if (fHistOff)
 | 
|---|
| 665 |         GetCatalog()->Draw("mirror same *");
 | 
|---|
| 666 | 
 | 
|---|
| 667 |     // ----- Pad number 1 -----
 | 
|---|
| 668 |     pad->cd();
 | 
|---|
| 669 |     name = MString::Format("%s_1", pad->GetName());
 | 
|---|
| 670 |     p = new TPad(name,name, 0.005, 0.3355, 0.475/*0.4975*/, 0.995, col, 0, 0);
 | 
|---|
| 671 |     p->SetNumber(1);
 | 
|---|
| 672 |     p->Draw();
 | 
|---|
| 673 |     p->cd();
 | 
|---|
| 674 | 
 | 
|---|
| 675 |     h3 = fHist.Project3D("yx_on");
 | 
|---|
| 676 |     h3->SetTitle("Distribution of equivalent events");
 | 
|---|
| 677 |     h3->SetDirectory(NULL);
 | 
|---|
| 678 |     h3->SetXTitle(fHist.GetXaxis()->GetTitle());
 | 
|---|
| 679 |     h3->SetYTitle(fHist.GetYaxis()->GetTitle());
 | 
|---|
| 680 |     h3->SetMinimum(0);
 | 
|---|
| 681 |     h3->Draw("colz");
 | 
|---|
| 682 |     h3->SetBit(kCanDelete);
 | 
|---|
| 683 | 
 | 
|---|
| 684 |     if (fHistOff)
 | 
|---|
| 685 |         GetCatalog()->Draw("mirror same *");
 | 
|---|
| 686 | 
 | 
|---|
| 687 |     // ----- Pad number 2 -----
 | 
|---|
| 688 |     pad->cd();
 | 
|---|
| 689 |     name = MString::Format("%s_2", pad->GetName());
 | 
|---|
| 690 |     p = new TPad(name,name, 0.005, 0.005, 0.2485, 0.3315, col, 0, 0);
 | 
|---|
| 691 |     p->SetNumber(2);
 | 
|---|
| 692 |     p->Draw();
 | 
|---|
| 693 |     p->cd();
 | 
|---|
| 694 |     h3->Draw("surf3");
 | 
|---|
| 695 | 
 | 
|---|
| 696 |     // ----- Pad number 3 -----
 | 
|---|
| 697 |     pad->cd();
 | 
|---|
| 698 |     name = MString::Format("%s_3", pad->GetName());
 | 
|---|
| 699 |     p = new TPad(name,name, 0.2525, 0.005, 0.4985, 0.3315, col, 0, 0);
 | 
|---|
| 700 |     p->SetNumber(3);
 | 
|---|
| 701 |     p->SetGrid();
 | 
|---|
| 702 |     p->Draw();
 | 
|---|
| 703 |     p->cd();
 | 
|---|
| 704 | 
 | 
|---|
| 705 |     const Double_t maxr = TMath::Hypot(h3->GetXaxis()->GetXmax(), h3->GetYaxis()->GetXmax());
 | 
|---|
| 706 |     const Int_t    nbin = (h3->GetNbinsX()+h3->GetNbinsY())/2;
 | 
|---|
| 707 |     TProfile *h = new TProfile("RadProf", "Radial Profile", nbin, 0, maxr);
 | 
|---|
| 708 |     h->SetDirectory(0);
 | 
|---|
| 709 |     //TH1F *h = new TH1F("RadProf", "Radial Profile", nbin, 0, maxr);
 | 
|---|
| 710 |     //h->Sumw2();
 | 
|---|
| 711 |     h->SetXTitle("\\vartheta [\\circ]");
 | 
|---|
| 712 |     h->SetYTitle("<cts>/\\Delta R");
 | 
|---|
| 713 |     h->SetBit(kCanDelete);
 | 
|---|
| 714 |     h->Draw();
 | 
|---|
| 715 | 
 | 
|---|
| 716 |     // ----- Pad number 5 -----
 | 
|---|
| 717 |     pad->cd();
 | 
|---|
| 718 |     name = MString::Format("%s_5", pad->GetName());
 | 
|---|
| 719 |     p = new TPad(name,name, 0.5025, 0.005, 0.7485, 0.3315, col, 0, 0);
 | 
|---|
| 720 |     p->SetNumber(5);
 | 
|---|
| 721 |     p->SetGrid();
 | 
|---|
| 722 |     p->Draw();
 | 
|---|
| 723 |     p->cd();
 | 
|---|
| 724 | 
 | 
|---|
| 725 |     h3 = new TH1D("Exc", "Distribution of excess events", 10, -1, 1);
 | 
|---|
| 726 |     h3->SetDirectory(NULL);
 | 
|---|
| 727 |     h3->SetXTitle("N");
 | 
|---|
| 728 |     h3->SetYTitle("Counts");
 | 
|---|
| 729 |     h3->Draw();
 | 
|---|
| 730 |     h3->SetBit(kCanDelete);
 | 
|---|
| 731 | 
 | 
|---|
| 732 |     // ----- Pad number 6 -----
 | 
|---|
| 733 |     pad->cd();
 | 
|---|
| 734 |     name = MString::Format("%s_6", pad->GetName());
 | 
|---|
| 735 |     p = new TPad(name,name, 0.7525, 0.005, 0.9995, 0.3315, col, 0, 0);
 | 
|---|
| 736 |     p->SetNumber(6);
 | 
|---|
| 737 |     p->SetGrid();
 | 
|---|
| 738 |     p->Draw();
 | 
|---|
| 739 |     p->cd();
 | 
|---|
| 740 | 
 | 
|---|
| 741 |     h3 = new TH1D("Sig", "Distribution of significances", 10, -1, 1);
 | 
|---|
| 742 |     h3->SetDirectory(NULL);
 | 
|---|
| 743 |     h3->SetXTitle("N");
 | 
|---|
| 744 |     h3->SetYTitle("Counts");
 | 
|---|
| 745 |     h3->Draw();
 | 
|---|
| 746 |     h3->SetBit(kCanDelete);
 | 
|---|
| 747 | }
 | 
|---|
| 748 | 
 | 
|---|
| 749 | // --------------------------------------------------------------------------
 | 
|---|
| 750 | //
 | 
|---|
| 751 | // The following resources are available:
 | 
|---|
| 752 | //
 | 
|---|
| 753 | //    MHDisp.Smearing: 0.1
 | 
|---|
| 754 | //    MHDisp.Wobble:   on/off
 | 
|---|
| 755 | //
 | 
|---|
| 756 | Int_t MHDisp::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
 | 
|---|
| 757 | {
 | 
|---|
| 758 |     Int_t rc = MHFalseSource::ReadEnv(env, prefix, print);
 | 
|---|
| 759 |     if (rc==kERROR)
 | 
|---|
| 760 |         return kERROR;
 | 
|---|
| 761 | 
 | 
|---|
| 762 |     if (IsEnvDefined(env, prefix, "Smearing", print))
 | 
|---|
| 763 |     {
 | 
|---|
| 764 |         rc = kTRUE;
 | 
|---|
| 765 |         SetSmearing(GetEnvValue(env, prefix, "Smearing", fSmearing));
 | 
|---|
| 766 |     }
 | 
|---|
| 767 |     if (IsEnvDefined(env, prefix, "Wobble", print))
 | 
|---|
| 768 |     {
 | 
|---|
| 769 |         rc = kTRUE;
 | 
|---|
| 770 |         SetWobble(GetEnvValue(env, prefix, "Wobble", fWobble));
 | 
|---|
| 771 |     }
 | 
|---|
| 772 |     if (IsEnvDefined(env, prefix, "ScaleMin", print))
 | 
|---|
| 773 |     {
 | 
|---|
| 774 |         rc = kTRUE;
 | 
|---|
| 775 |         SetScaleMin(GetEnvValue(env, prefix, "ScaleMin", fScaleMin));
 | 
|---|
| 776 |     }
 | 
|---|
| 777 |     if (IsEnvDefined(env, prefix, "ScaleMax", print))
 | 
|---|
| 778 |     {
 | 
|---|
| 779 |         rc = kTRUE;
 | 
|---|
| 780 |         SetScaleMax(GetEnvValue(env, prefix, "ScaleMax", fScaleMax));
 | 
|---|
| 781 |     }
 | 
|---|
| 782 | 
 | 
|---|
| 783 |     return rc;
 | 
|---|
| 784 | }
 | 
|---|
| 785 | 
 | 
|---|