| 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 | // MHFalseSource | 
|---|
| 28 | // | 
|---|
| 29 | // Create a false source plot. For the Binning in x,y the object MBinning | 
|---|
| 30 | // "BinningFalseSource" is taken from the paremeter list. | 
|---|
| 31 | // | 
|---|
| 32 | // The binning in alpha is currently fixed as 18bins from 0 to 90deg. | 
|---|
| 33 | // | 
|---|
| 34 | // If MTime, MObservatory and MPointingPos is available in the paremeter | 
|---|
| 35 | // list each image is derotated. | 
|---|
| 36 | // | 
|---|
| 37 | // MHFalseSource fills a 3D histogram with alpha distribution for | 
|---|
| 38 | // each (derotated) x and y position. | 
|---|
| 39 | // | 
|---|
| 40 | // Only a radius of 1.2deg around the camera center is taken into account. | 
|---|
| 41 | // | 
|---|
| 42 | // The displayed histogram is the projection of alpha<fAlphaCut into | 
|---|
| 43 | // the x,y plain and the projection of alpha>90-fAlphaCut | 
|---|
| 44 | // | 
|---|
| 45 | // By using the context menu (find it in a small region between the single | 
|---|
| 46 | // pads) you can change the AlphaCut 'online' | 
|---|
| 47 | // | 
|---|
| 48 | // Here is a slightly simplified version of the algorithm: | 
|---|
| 49 | // ------------------------------------------------------ | 
|---|
| 50 | //    MHillas hil; // Taken as second argument in MFillH | 
|---|
| 51 | // | 
|---|
| 52 | //    MSrcPosCam src; | 
|---|
| 53 | //    MHillasSrc hsrc; | 
|---|
| 54 | //    hsrc.SetSrcPos(&src); | 
|---|
| 55 | // | 
|---|
| 56 | //    for (int ix=0; ix<nx; ix++) | 
|---|
| 57 | //        for (int iy=0; iy<ny; iy++) | 
|---|
| 58 | //        { | 
|---|
| 59 | //            TVector2 v(cx[ix], cy[iy]); //cx center of x-bin ix | 
|---|
| 60 | //            if (rho!=0)                 //cy center of y-bin iy | 
|---|
| 61 | //                v.Rotate(-rho);         //image rotation angle | 
|---|
| 62 | // | 
|---|
| 63 | //            src.SetXY(v);               //source position in the camera | 
|---|
| 64 | // | 
|---|
| 65 | //            if (!hsrc.Calc(hil))        //calc source dependant hillas | 
|---|
| 66 | //                return; | 
|---|
| 67 | // | 
|---|
| 68 | //            //fill absolute alpha into histogram | 
|---|
| 69 | //            const Double_t alpha = hsrc.GetAlpha(); | 
|---|
| 70 | //            fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w); | 
|---|
| 71 | //        } | 
|---|
| 72 | //    } | 
|---|
| 73 | // | 
|---|
| 74 | // Use MHFalse Source like this: | 
|---|
| 75 | // ----------------------------- | 
|---|
| 76 | //    MFillH fill("MHFalseSource", "MHillas"); | 
|---|
| 77 | // or | 
|---|
| 78 | //    MHFalseSource hist; | 
|---|
| 79 | //    hist.SetAlphaCut(12.5);  // The default value | 
|---|
| 80 | //    hist.SetBgmean(55);      // The default value | 
|---|
| 81 | //    MFillH fill(&hist, "MHillas"); | 
|---|
| 82 | // | 
|---|
| 83 | // To be implemented: | 
|---|
| 84 | // ------------------ | 
|---|
| 85 | //  - a switch to use alpha or |alpha| | 
|---|
| 86 | //  - taking the binning for alpha from the parlist (binning is | 
|---|
| 87 | //    currently fixed) | 
|---|
| 88 | //  - a possibility to change the fit-function (eg using a different TF1) | 
|---|
| 89 | //  - a possibility to change the fit algorithm (eg which paremeters | 
|---|
| 90 | //    are fixed at which time) | 
|---|
| 91 | //  - use a different varaible than alpha | 
|---|
| 92 | //  - a possiblity to change the algorithm which is used to calculate | 
|---|
| 93 | //    alpha (or another parameter) is missing (this could be done using | 
|---|
| 94 | //    a tasklist somewhere...) | 
|---|
| 95 | //  - a more clever (and faster) algorithm to fill the histogram, eg by | 
|---|
| 96 | //    calculating alpha once and fill the two coils around the mean | 
|---|
| 97 | //  - more drawing options... | 
|---|
| 98 | //  - Move Significance() to a more 'general' place and implement | 
|---|
| 99 | //    also different algorithms like (Li/Ma) | 
|---|
| 100 | //  - implement fit for best alpha distribution -- online (Paint) | 
|---|
| 101 | // | 
|---|
| 102 | ////////////////////////////////////////////////////////////////////////////// | 
|---|
| 103 | #include "MHFalseSource.h" | 
|---|
| 104 |  | 
|---|
| 105 | #include <TF1.h> | 
|---|
| 106 | #include <TH2.h> | 
|---|
| 107 | #include <TGraph.h> | 
|---|
| 108 | #include <TStyle.h> | 
|---|
| 109 | #include <TCanvas.h> | 
|---|
| 110 | #include <TPaveText.h> | 
|---|
| 111 | #include <TStopwatch.h> | 
|---|
| 112 |  | 
|---|
| 113 | #include "MGeomCam.h" | 
|---|
| 114 | #include "MSrcPosCam.h" | 
|---|
| 115 | #include "MHillasSrc.h" | 
|---|
| 116 | #include "MTime.h" | 
|---|
| 117 | #include "MObservatory.h" | 
|---|
| 118 | #include "MPointingPos.h" | 
|---|
| 119 | #include "MAstroSky2Local.h" | 
|---|
| 120 | #include "MStatusDisplay.h" | 
|---|
| 121 | #include "MMath.h" | 
|---|
| 122 |  | 
|---|
| 123 | #include "MBinning.h" | 
|---|
| 124 | #include "MParList.h" | 
|---|
| 125 |  | 
|---|
| 126 | #include "MLog.h" | 
|---|
| 127 | #include "MLogManip.h" | 
|---|
| 128 |  | 
|---|
| 129 | ClassImp(MHFalseSource); | 
|---|
| 130 |  | 
|---|
| 131 | using namespace std; | 
|---|
| 132 |  | 
|---|
| 133 | // -------------------------------------------------------------------------- | 
|---|
| 134 | // | 
|---|
| 135 | // Default Constructor | 
|---|
| 136 | // | 
|---|
| 137 | MHFalseSource::MHFalseSource(const char *name, const char *title) | 
|---|
| 138 | : fMm2Deg(-1), fUseMmScale(kTRUE), fAlphaCut(12.5), fBgMean(55) | 
|---|
| 139 | { | 
|---|
| 140 | // | 
|---|
| 141 | //   set the name and title of this object | 
|---|
| 142 | // | 
|---|
| 143 | fName  = name  ? name  : "MHFalseSource"; | 
|---|
| 144 | fTitle = title ? title : "3D-plot of Alpha vs x, y"; | 
|---|
| 145 |  | 
|---|
| 146 | fHist.SetDirectory(NULL); | 
|---|
| 147 |  | 
|---|
| 148 | fHist.SetName("Alpha"); | 
|---|
| 149 | fHist.SetTitle("3D-plot of Alpha vs x, y"); | 
|---|
| 150 | fHist.SetXTitle("x [\\circ]"); | 
|---|
| 151 | fHist.SetYTitle("y [\\circ]"); | 
|---|
| 152 | fHist.SetZTitle("\\alpha [\\circ]"); | 
|---|
| 153 | } | 
|---|
| 154 |  | 
|---|
| 155 | // -------------------------------------------------------------------------- | 
|---|
| 156 | // | 
|---|
| 157 | // Set the alpha cut (|alpha|<fAlphaCut) which is use to estimate the | 
|---|
| 158 | // number of excess events | 
|---|
| 159 | // | 
|---|
| 160 | void MHFalseSource::SetAlphaCut(Float_t alpha) | 
|---|
| 161 | { | 
|---|
| 162 | if (alpha<0) | 
|---|
| 163 | *fLog << warn << "Alpha<0... taking absolute value." << endl; | 
|---|
| 164 |  | 
|---|
| 165 | fAlphaCut = TMath::Abs(alpha); | 
|---|
| 166 |  | 
|---|
| 167 | Modified(); | 
|---|
| 168 | } | 
|---|
| 169 |  | 
|---|
| 170 | // -------------------------------------------------------------------------- | 
|---|
| 171 | // | 
|---|
| 172 | // Set mean alpha around which the off sample is determined | 
|---|
| 173 | // (fBgMean-fAlphaCut/2<|fAlpha|<fBgMean+fAlphaCut/2) which is use | 
|---|
| 174 | // to estimate the number of off events | 
|---|
| 175 | // | 
|---|
| 176 | void MHFalseSource::SetBgMean(Float_t alpha) | 
|---|
| 177 | { | 
|---|
| 178 | if (alpha<0) | 
|---|
| 179 | *fLog << warn << "Alpha<0... taking absolute value." << endl; | 
|---|
| 180 |  | 
|---|
| 181 | fBgMean = TMath::Abs(alpha); | 
|---|
| 182 |  | 
|---|
| 183 | Modified(); | 
|---|
| 184 | } | 
|---|
| 185 |  | 
|---|
| 186 | // -------------------------------------------------------------------------- | 
|---|
| 187 | // | 
|---|
| 188 | // Use this function to setup your own conversion factor between degrees | 
|---|
| 189 | // and millimeters. The conversion factor should be the one calculated in | 
|---|
| 190 | // MGeomCam. Use this function with Caution: You could create wrong values | 
|---|
| 191 | // by setting up your own scale factor. | 
|---|
| 192 | // | 
|---|
| 193 | void MHFalseSource::SetMm2Deg(Float_t mmdeg) | 
|---|
| 194 | { | 
|---|
| 195 | if (mmdeg<0) | 
|---|
| 196 | { | 
|---|
| 197 | *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl; | 
|---|
| 198 | return; | 
|---|
| 199 | } | 
|---|
| 200 |  | 
|---|
| 201 | if (fMm2Deg>=0) | 
|---|
| 202 | *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl; | 
|---|
| 203 |  | 
|---|
| 204 | fMm2Deg = mmdeg; | 
|---|
| 205 | } | 
|---|
| 206 |  | 
|---|
| 207 | // -------------------------------------------------------------------------- | 
|---|
| 208 | // | 
|---|
| 209 | // With this function you can convert the histogram ('on the fly') between | 
|---|
| 210 | // degrees and millimeters. | 
|---|
| 211 | // | 
|---|
| 212 | void MHFalseSource::SetMmScale(Bool_t mmscale) | 
|---|
| 213 | { | 
|---|
| 214 | if (fUseMmScale == mmscale) | 
|---|
| 215 | return; | 
|---|
| 216 |  | 
|---|
| 217 | if (fMm2Deg<0) | 
|---|
| 218 | { | 
|---|
| 219 | *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl; | 
|---|
| 220 | return; | 
|---|
| 221 | } | 
|---|
| 222 |  | 
|---|
| 223 | if (fUseMmScale) | 
|---|
| 224 | { | 
|---|
| 225 | fHist.SetXTitle("x [mm]"); | 
|---|
| 226 | fHist.SetYTitle("y [mm]"); | 
|---|
| 227 |  | 
|---|
| 228 | fHist.Scale(1./fMm2Deg); | 
|---|
| 229 | } | 
|---|
| 230 | else | 
|---|
| 231 | { | 
|---|
| 232 | fHist.SetXTitle("x [\\circ]"); | 
|---|
| 233 | fHist.SetYTitle("y [\\circ]"); | 
|---|
| 234 |  | 
|---|
| 235 | fHist.Scale(1./fMm2Deg); | 
|---|
| 236 | } | 
|---|
| 237 |  | 
|---|
| 238 | fUseMmScale = mmscale; | 
|---|
| 239 | } | 
|---|
| 240 |  | 
|---|
| 241 | // -------------------------------------------------------------------------- | 
|---|
| 242 | // | 
|---|
| 243 | // Calculate Significance as | 
|---|
| 244 | // significance = (s-b)/sqrt(s+k*k*b) mit k=s/b | 
|---|
| 245 | // | 
|---|
| 246 | // s: total number of events in signal region | 
|---|
| 247 | // b: number of background events in signal region | 
|---|
| 248 | // | 
|---|
| 249 | Double_t MHFalseSource::Significance(Double_t s, Double_t b) | 
|---|
| 250 | { | 
|---|
| 251 | return MMath::Significance(s, b); | 
|---|
| 252 | } | 
|---|
| 253 |  | 
|---|
| 254 | // -------------------------------------------------------------------------- | 
|---|
| 255 | // | 
|---|
| 256 | //  calculates the significance according to Li & Ma | 
|---|
| 257 | //  ApJ 272 (1983) 317, Formula 17 | 
|---|
| 258 | // | 
|---|
| 259 | //  s                    // s: number of on events | 
|---|
| 260 | //  b                    // b: number of off events | 
|---|
| 261 | //  alpha = t_on/t_off;  // t: observation time | 
|---|
| 262 | // | 
|---|
| 263 | Double_t MHFalseSource::SignificanceLiMa(Double_t s, Double_t b, Double_t alpha) | 
|---|
| 264 | { | 
|---|
| 265 | return MMath::SignificanceLiMa(s, b); | 
|---|
| 266 | } | 
|---|
| 267 |  | 
|---|
| 268 | // -------------------------------------------------------------------------- | 
|---|
| 269 | // | 
|---|
| 270 | // Set binnings (takes BinningFalseSource) and prepare filling of the | 
|---|
| 271 | // histogram. | 
|---|
| 272 | // | 
|---|
| 273 | // Also search for MTime, MObservatory and MPointingPos | 
|---|
| 274 | // | 
|---|
| 275 | Bool_t MHFalseSource::SetupFill(const MParList *plist) | 
|---|
| 276 | { | 
|---|
| 277 | const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam"); | 
|---|
| 278 | if (geom) | 
|---|
| 279 | { | 
|---|
| 280 | fMm2Deg = geom->GetConvMm2Deg(); | 
|---|
| 281 | fUseMmScale = kFALSE; | 
|---|
| 282 |  | 
|---|
| 283 | fHist.SetXTitle("x [\\circ]"); | 
|---|
| 284 | fHist.SetYTitle("y [\\circ]"); | 
|---|
| 285 | } | 
|---|
| 286 |  | 
|---|
| 287 | MBinning binsa; | 
|---|
| 288 | binsa.SetEdges(18, 0, 90); | 
|---|
| 289 |  | 
|---|
| 290 | const MBinning *bins = (MBinning*)plist->FindObject("BinningFalseSource"); | 
|---|
| 291 | if (!bins) | 
|---|
| 292 | { | 
|---|
| 293 | Float_t r = geom ? geom->GetMaxRadius() : 600; | 
|---|
| 294 | r /= 3; | 
|---|
| 295 | if (!fUseMmScale) | 
|---|
| 296 | r *= fMm2Deg; | 
|---|
| 297 |  | 
|---|
| 298 | MBinning b; | 
|---|
| 299 | b.SetEdges(20, -r, r); | 
|---|
| 300 | SetBinning(&fHist, &b, &b, &binsa); | 
|---|
| 301 | } | 
|---|
| 302 | else | 
|---|
| 303 | SetBinning(&fHist, bins, bins, &binsa); | 
|---|
| 304 |  | 
|---|
| 305 | fPointPos = (MPointingPos*)plist->FindObject(AddSerialNumber("MPointingPos")); | 
|---|
| 306 | if (!fPointPos) | 
|---|
| 307 | *fLog << warn << "MPointingPos not found... no derotation." << endl; | 
|---|
| 308 |  | 
|---|
| 309 | fTime = (MTime*)plist->FindObject(AddSerialNumber("MTime")); | 
|---|
| 310 | if (!fTime) | 
|---|
| 311 | *fLog << warn << "MTime not found... no derotation." << endl; | 
|---|
| 312 |  | 
|---|
| 313 | fObservatory = (MObservatory*)plist->FindObject(AddSerialNumber("MObservatory")); | 
|---|
| 314 | if (!fObservatory) | 
|---|
| 315 | *fLog << err << "MObservatory not found...  no derotation." << endl; | 
|---|
| 316 |  | 
|---|
| 317 |  | 
|---|
| 318 | return kTRUE; | 
|---|
| 319 | } | 
|---|
| 320 |  | 
|---|
| 321 | // -------------------------------------------------------------------------- | 
|---|
| 322 | // | 
|---|
| 323 | // Fill the histogram. For details see the code or the class description | 
|---|
| 324 | // | 
|---|
| 325 | Bool_t MHFalseSource::Fill(const MParContainer *par, const Stat_t w) | 
|---|
| 326 | { | 
|---|
| 327 | MHillas *hil = (MHillas*)par; | 
|---|
| 328 |  | 
|---|
| 329 | const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1)); | 
|---|
| 330 |  | 
|---|
| 331 | // Get camera rotation angle | 
|---|
| 332 | Double_t rho = 0; | 
|---|
| 333 | if (fTime && fObservatory && fPointPos) | 
|---|
| 334 | rho = fPointPos->RotationAngle(*fObservatory, *fTime); | 
|---|
| 335 | //if (fPointPos) | 
|---|
| 336 | //    rho = fPointPos->RotationAngle(*fObservatory); | 
|---|
| 337 |  | 
|---|
| 338 | MSrcPosCam src; | 
|---|
| 339 | MHillasSrc hsrc; | 
|---|
| 340 | hsrc.SetSrcPos(&src); | 
|---|
| 341 |  | 
|---|
| 342 | const Int_t nx = fHist.GetNbinsX(); | 
|---|
| 343 | const Int_t ny = fHist.GetNbinsY(); | 
|---|
| 344 |  | 
|---|
| 345 | Axis_t cx[nx]; | 
|---|
| 346 | Axis_t cy[ny]; | 
|---|
| 347 | fHist.GetXaxis()->GetCenter(cx); | 
|---|
| 348 | fHist.GetYaxis()->GetCenter(cy); | 
|---|
| 349 |  | 
|---|
| 350 | for (int ix=0; ix<nx; ix++) | 
|---|
| 351 | { | 
|---|
| 352 | for (int iy=0; iy<ny; iy++) | 
|---|
| 353 | { | 
|---|
| 354 | if (TMath::Hypot(cx[ix], cy[iy])>maxr) | 
|---|
| 355 | continue; | 
|---|
| 356 |  | 
|---|
| 357 | TVector2 v(cx[ix], cy[iy]); | 
|---|
| 358 | if (rho!=0) | 
|---|
| 359 | v.Rotate(-rho); | 
|---|
| 360 |  | 
|---|
| 361 | if (!fUseMmScale) | 
|---|
| 362 | v *= 1./fMm2Deg; | 
|---|
| 363 |  | 
|---|
| 364 | src.SetXY(v); | 
|---|
| 365 |  | 
|---|
| 366 | if (!hsrc.Calc(hil)) | 
|---|
| 367 | { | 
|---|
| 368 | *fLog << warn << "Calculation of MHillasSrc failed for x=" << cx[ix] << " y=" << cy[iy] << endl; | 
|---|
| 369 | return kFALSE; | 
|---|
| 370 | } | 
|---|
| 371 |  | 
|---|
| 372 | const Double_t alpha = hsrc.GetAlpha(); | 
|---|
| 373 | fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w); | 
|---|
| 374 | } | 
|---|
| 375 | } | 
|---|
| 376 |  | 
|---|
| 377 | return kTRUE; | 
|---|
| 378 | } | 
|---|
| 379 |  | 
|---|
| 380 | // -------------------------------------------------------------------------- | 
|---|
| 381 | // | 
|---|
| 382 | // Create projection for off data, taking sum of bin contents of | 
|---|
| 383 | // range (fBgMean-fAlphaCut/2, fBgMean+fAlphaCut) Making sure to take | 
|---|
| 384 | // the same number of bins than for on-data | 
|---|
| 385 | // | 
|---|
| 386 | void MHFalseSource::ProjectOff(TH2D *h2) | 
|---|
| 387 | { | 
|---|
| 388 | TAxis &axe = *fHist.GetZaxis(); | 
|---|
| 389 |  | 
|---|
| 390 | // Find range to cut (left edge and width) | 
|---|
| 391 | const Int_t f = axe.FindBin(fBgMean-fAlphaCut/2); | 
|---|
| 392 | const Int_t l = axe.FindBin(fAlphaCut)+f-1; | 
|---|
| 393 |  | 
|---|
| 394 | axe.SetRange(f, l); | 
|---|
| 395 | const Float_t cut1 = axe.GetBinLowEdge(f); | 
|---|
| 396 | const Float_t cut2 = axe.GetBinUpEdge(l); | 
|---|
| 397 | h2->SetTitle(Form("Distribution of %.1f\\circ<|\\alpha|<%.1f\\circ in x,y", cut1, cut2)); | 
|---|
| 398 |  | 
|---|
| 399 | // Get projection for range | 
|---|
| 400 | TH2D *p = (TH2D*)fHist.Project3D("xy_off"); | 
|---|
| 401 |  | 
|---|
| 402 | // Reset range | 
|---|
| 403 | axe.SetRange(0,9999); | 
|---|
| 404 |  | 
|---|
| 405 | // Move contents from projection to h2 | 
|---|
| 406 | h2->Reset(); | 
|---|
| 407 | h2->Add(p); | 
|---|
| 408 |  | 
|---|
| 409 | // Delete p | 
|---|
| 410 | delete p; | 
|---|
| 411 |  | 
|---|
| 412 | // Set Minimum as minimum value Greater Than 0 | 
|---|
| 413 | h2->SetMinimum(GetMinimumGT(*h2)); | 
|---|
| 414 | } | 
|---|
| 415 |  | 
|---|
| 416 | // -------------------------------------------------------------------------- | 
|---|
| 417 | // | 
|---|
| 418 | // Create projection for on data, taking sum of bin contents of | 
|---|
| 419 | // range (0, fAlphaCut) | 
|---|
| 420 | // | 
|---|
| 421 | void MHFalseSource::ProjectOn(TH2D *h3) | 
|---|
| 422 | { | 
|---|
| 423 | TAxis &axe = *fHist.GetZaxis(); | 
|---|
| 424 |  | 
|---|
| 425 | // Find range to cut | 
|---|
| 426 | axe.SetRangeUser(0, fAlphaCut); | 
|---|
| 427 | const Float_t cut = axe.GetBinUpEdge(axe.GetLast()); | 
|---|
| 428 | h3->SetTitle(Form("Distribution of |\\alpha|<%.1f\\circ in x,y", cut)); | 
|---|
| 429 |  | 
|---|
| 430 | // Get projection for range | 
|---|
| 431 | TH2D *p = (TH2D*)fHist.Project3D("xy_on"); | 
|---|
| 432 |  | 
|---|
| 433 | // Reset range | 
|---|
| 434 | axe.SetRange(0,9999); | 
|---|
| 435 |  | 
|---|
| 436 | // Move contents from projection to h3 | 
|---|
| 437 | h3->Reset(); | 
|---|
| 438 | h3->Add(p); | 
|---|
| 439 | delete p; | 
|---|
| 440 |  | 
|---|
| 441 | // Set Minimum as minimum value Greater Than 0 | 
|---|
| 442 | h3->SetMinimum(GetMinimumGT(*h3)); | 
|---|
| 443 | } | 
|---|
| 444 |  | 
|---|
| 445 | // -------------------------------------------------------------------------- | 
|---|
| 446 | // | 
|---|
| 447 | // Update the projections | 
|---|
| 448 | // | 
|---|
| 449 | void MHFalseSource::Paint(Option_t *opt) | 
|---|
| 450 | { | 
|---|
| 451 | // sigma = (s-b)/sqrt(s+k*k*b) mit k=s/b | 
|---|
| 452 |  | 
|---|
| 453 | gStyle->SetPalette(1, 0); | 
|---|
| 454 |  | 
|---|
| 455 | TVirtualPad *padsave = gPad; | 
|---|
| 456 |  | 
|---|
| 457 | TH1D* h1; | 
|---|
| 458 | TH2D* h2; | 
|---|
| 459 | TH2D* h3; | 
|---|
| 460 | TH2D* h4; | 
|---|
| 461 |  | 
|---|
| 462 | padsave->cd(3); | 
|---|
| 463 | if ((h3 = (TH2D*)gPad->FindObject("Alpha_xy_on"))) | 
|---|
| 464 | ProjectOn(h3); | 
|---|
| 465 |  | 
|---|
| 466 | padsave->cd(4); | 
|---|
| 467 | if ((h2 = (TH2D*)gPad->FindObject("Alpha_xy_off"))) | 
|---|
| 468 | ProjectOff(h2); | 
|---|
| 469 |  | 
|---|
| 470 | padsave->cd(2); | 
|---|
| 471 | if (h2 && h3 && (h4 = (TH2D*)gPad->FindObject("Alpha_xy_sig"))) | 
|---|
| 472 | { | 
|---|
| 473 | const Int_t nx = h4->GetXaxis()->GetNbins(); | 
|---|
| 474 | const Int_t ny = h4->GetYaxis()->GetNbins(); | 
|---|
| 475 |  | 
|---|
| 476 | Int_t maxx=nx/2; | 
|---|
| 477 | Int_t maxy=ny/2; | 
|---|
| 478 |  | 
|---|
| 479 | Int_t max = h4->GetBin(maxx, maxy); | 
|---|
| 480 |  | 
|---|
| 481 | for (int ix=0; ix<nx; ix++) | 
|---|
| 482 | for (int iy=0; iy<ny; iy++) | 
|---|
| 483 | { | 
|---|
| 484 | const Int_t n = h4->GetBin(ix+1, iy+1); | 
|---|
| 485 |  | 
|---|
| 486 | const Double_t s = h3->GetBinContent(n); | 
|---|
| 487 | const Double_t b = h2->GetBinContent(n); | 
|---|
| 488 |  | 
|---|
| 489 | const Double_t sig = Significance(s, b); | 
|---|
| 490 |  | 
|---|
| 491 | h4->SetBinContent(n, sig); | 
|---|
| 492 |  | 
|---|
| 493 | if (sig>h4->GetBinContent(max) && sig!=0) | 
|---|
| 494 | { | 
|---|
| 495 | max = n; | 
|---|
| 496 | maxx=ix; | 
|---|
| 497 | maxy=iy; | 
|---|
| 498 | } | 
|---|
| 499 | } | 
|---|
| 500 |  | 
|---|
| 501 | padsave->cd(1); | 
|---|
| 502 | if ((h1 = (TH1D*)gPad->FindObject("Alpha"))) | 
|---|
| 503 | { | 
|---|
| 504 | //h1->Reset(); | 
|---|
| 505 |  | 
|---|
| 506 | const Double_t x = fHist.GetXaxis()->GetBinCenter(maxx); | 
|---|
| 507 | const Double_t y = fHist.GetYaxis()->GetBinCenter(maxy); | 
|---|
| 508 | const Double_t s = h4->GetBinContent(max); | 
|---|
| 509 |  | 
|---|
| 510 | TH1 *h = fHist.ProjectionZ("Alpha", maxx, maxx, maxy, maxy); | 
|---|
| 511 | h->SetTitle(Form("Distribution of \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f)", x, y, s)); | 
|---|
| 512 | } | 
|---|
| 513 | } | 
|---|
| 514 |  | 
|---|
| 515 | gPad = padsave; | 
|---|
| 516 | } | 
|---|
| 517 |  | 
|---|
| 518 | // -------------------------------------------------------------------------- | 
|---|
| 519 | // | 
|---|
| 520 | // Draw the histogram | 
|---|
| 521 | // | 
|---|
| 522 | void MHFalseSource::Draw(Option_t *opt) | 
|---|
| 523 | { | 
|---|
| 524 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this); | 
|---|
| 525 | pad->SetBorderMode(0); | 
|---|
| 526 |  | 
|---|
| 527 | AppendPad(""); | 
|---|
| 528 |  | 
|---|
| 529 | pad->Divide(2, 2); | 
|---|
| 530 |  | 
|---|
| 531 | // draw the 2D histogram Sigmabar versus Theta | 
|---|
| 532 | pad->cd(1); | 
|---|
| 533 | gPad->SetBorderMode(0); | 
|---|
| 534 | TH1 *h1 = fHist.ProjectionZ("Alpha"); | 
|---|
| 535 | h1->SetDirectory(NULL); | 
|---|
| 536 | h1->SetTitle("Distribution of \\alpha"); | 
|---|
| 537 | h1->SetXTitle(fHist.GetZaxis()->GetTitle()); | 
|---|
| 538 | h1->SetYTitle("Counts"); | 
|---|
| 539 | h1->Draw(opt); | 
|---|
| 540 | h1->SetBit(kCanDelete); | 
|---|
| 541 |  | 
|---|
| 542 | pad->cd(4); | 
|---|
| 543 | gPad->SetBorderMode(0); | 
|---|
| 544 | fHist.GetZaxis()->SetRangeUser(fBgMean-fAlphaCut/2, fBgMean+fAlphaCut/2); | 
|---|
| 545 | TH1 *h2 = fHist.Project3D("xy_off"); | 
|---|
| 546 | h2->SetDirectory(NULL); | 
|---|
| 547 | h2->SetXTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 548 | h2->SetYTitle(fHist.GetYaxis()->GetTitle()); | 
|---|
| 549 | h2->Draw("colz"); | 
|---|
| 550 | h2->SetBit(kCanDelete); | 
|---|
| 551 |  | 
|---|
| 552 | pad->cd(3); | 
|---|
| 553 | gPad->SetBorderMode(0); | 
|---|
| 554 | fHist.GetZaxis()->SetRangeUser(0,fAlphaCut); | 
|---|
| 555 | TH1 *h3 = fHist.Project3D("xy_on"); | 
|---|
| 556 | fHist.GetZaxis()->SetRange(0,9999); | 
|---|
| 557 | h3->SetDirectory(NULL); | 
|---|
| 558 | h3->SetXTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 559 | h3->SetYTitle(fHist.GetYaxis()->GetTitle()); | 
|---|
| 560 | h3->Draw("colz"); | 
|---|
| 561 | h3->SetBit(kCanDelete); | 
|---|
| 562 |  | 
|---|
| 563 | pad->cd(2); | 
|---|
| 564 | gPad->SetBorderMode(0); | 
|---|
| 565 | fHist.GetZaxis()->SetRange(0,0); | 
|---|
| 566 | TH1 *h4 = fHist.Project3D("xy_sig"); // Do this to get the correct binning.... | 
|---|
| 567 | fHist.GetZaxis()->SetRange(0,9999); | 
|---|
| 568 | h4->SetTitle("Significance"); | 
|---|
| 569 | h4->SetDirectory(NULL); | 
|---|
| 570 | h4->SetXTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 571 | h4->SetYTitle(fHist.GetYaxis()->GetTitle()); | 
|---|
| 572 | h4->Draw("colz"); | 
|---|
| 573 | h4->SetBit(kCanDelete); | 
|---|
| 574 | } | 
|---|
| 575 |  | 
|---|
| 576 | // -------------------------------------------------------------------------- | 
|---|
| 577 | // | 
|---|
| 578 | // Everything which is in the main pad belongs to this class! | 
|---|
| 579 | // | 
|---|
| 580 | Int_t MHFalseSource::DistancetoPrimitive(Int_t px, Int_t py) | 
|---|
| 581 | { | 
|---|
| 582 | return 0; | 
|---|
| 583 | } | 
|---|
| 584 |  | 
|---|
| 585 | // -------------------------------------------------------------------------- | 
|---|
| 586 | // | 
|---|
| 587 | // Set all sub-pads to Modified() | 
|---|
| 588 | // | 
|---|
| 589 | void MHFalseSource::Modified() | 
|---|
| 590 | { | 
|---|
| 591 | if (!gPad) | 
|---|
| 592 | return; | 
|---|
| 593 |  | 
|---|
| 594 | TVirtualPad *padsave = gPad; | 
|---|
| 595 | padsave->Modified(); | 
|---|
| 596 | padsave->cd(1); | 
|---|
| 597 | gPad->Modified(); | 
|---|
| 598 | padsave->cd(2); | 
|---|
| 599 | gPad->Modified(); | 
|---|
| 600 | padsave->cd(3); | 
|---|
| 601 | gPad->Modified(); | 
|---|
| 602 | padsave->cd(4); | 
|---|
| 603 | gPad->Modified(); | 
|---|
| 604 | gPad->cd(); | 
|---|
| 605 | } | 
|---|
| 606 |  | 
|---|
| 607 | // -------------------------------------------------------------------------- | 
|---|
| 608 | // | 
|---|
| 609 | // This is a preliminary implementation of a alpha-fit procedure for | 
|---|
| 610 | // all possible source positions. It will be moved into its own | 
|---|
| 611 | // more powerfull class soon. | 
|---|
| 612 | // | 
|---|
| 613 | // The fit function is "gaus(0)+pol2(3)" which is equivalent to: | 
|---|
| 614 | //   [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2 | 
|---|
| 615 | // or | 
|---|
| 616 | //   A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2 | 
|---|
| 617 | // | 
|---|
| 618 | // Parameter [1] is fixed to 0 while the alpha peak should be | 
|---|
| 619 | // symmetric around alpha=0. | 
|---|
| 620 | // | 
|---|
| 621 | // Parameter [4] is fixed to 0 because the first derivative at | 
|---|
| 622 | // alpha=0 should be 0, too. | 
|---|
| 623 | // | 
|---|
| 624 | // In a first step the background is fitted between bgmin and bgmax, | 
|---|
| 625 | // while the parameters [0]=0 and [2]=1 are fixed. | 
|---|
| 626 | // | 
|---|
| 627 | // In a second step the signal region (alpha<sigmax) is fittet using | 
|---|
| 628 | // the whole function with parameters [1], [3], [4] and [5] fixed. | 
|---|
| 629 | // | 
|---|
| 630 | // The number of excess and background events are calculated as | 
|---|
| 631 | //   s = int(0, sigint, gaus(0)+pol2(3)) | 
|---|
| 632 | //   b = int(0, sigint,         pol2(3)) | 
|---|
| 633 | // | 
|---|
| 634 | // The Significance is calculated using the Significance() member | 
|---|
| 635 | // function. | 
|---|
| 636 | // | 
|---|
| 637 | void MHFalseSource::FitSignificance(Float_t sigint, Float_t sigmax, Float_t bgmin, Float_t bgmax, Byte_t polynom) | 
|---|
| 638 | { | 
|---|
| 639 | TH1D h0a("A",          "", 50,   0, 4000); | 
|---|
| 640 | TH1D h4a("chisq1",     "", 50,   0,   35); | 
|---|
| 641 | //TH1D h5a("prob1",      "", 50,   0,  1.1); | 
|---|
| 642 | TH1D h6("signifcance", "", 50, -20,   20); | 
|---|
| 643 |  | 
|---|
| 644 | TH1D h1("mu",    "Parameter \\mu",    50,   -1,    1); | 
|---|
| 645 | TH1D h2("sigma", "Parameter \\sigma", 50,    0,   90); | 
|---|
| 646 | TH1D h3("b",     "Parameter b",       50, -0.1,  0.1); | 
|---|
| 647 |  | 
|---|
| 648 | TH1D h0b("a",         "Parameter a (red), A (blue)", 50, 0, 4000); | 
|---|
| 649 | TH1D h4b("\\chi^{2}", "\\chi^{2} (red, green) / significance (black)", 50, 0, 35); | 
|---|
| 650 | //TH1D h5b("prob",      "Fit probability: Bg(red), F(blue)", 50, 0, 1.1); | 
|---|
| 651 |  | 
|---|
| 652 | h0a.SetLineColor(kBlue); | 
|---|
| 653 | h4a.SetLineColor(kBlue); | 
|---|
| 654 | //h5a.SetLineColor(kBlue); | 
|---|
| 655 | h0b.SetLineColor(kRed); | 
|---|
| 656 | h4b.SetLineColor(kRed); | 
|---|
| 657 | //h5b.SetLineColor(kRed); | 
|---|
| 658 |  | 
|---|
| 659 | TH1 *hist  = fHist.Project3D("xy_fit"); | 
|---|
| 660 | hist->SetDirectory(0); | 
|---|
| 661 | TH1 *hists = fHist.Project3D("xy_fit"); | 
|---|
| 662 | hists->SetDirectory(0); | 
|---|
| 663 | TH1 *histb = fHist.Project3D("xy_fit"); | 
|---|
| 664 | histb->SetDirectory(0); | 
|---|
| 665 | hist->Reset(); | 
|---|
| 666 | hists->Reset(); | 
|---|
| 667 | histb->Reset(); | 
|---|
| 668 | hist->SetNameTitle("Significance", | 
|---|
| 669 | Form("Fit Region: Signal<%.1f\\circ, %.1f\\circ<Bg<%.1f\\circ", | 
|---|
| 670 | sigmax, bgmin, bgmax)); | 
|---|
| 671 | hists->SetNameTitle("Excess",     Form("Number of excess events for \\alpha<%.0f\\circ", sigint)); | 
|---|
| 672 | histb->SetNameTitle("Background", Form("Number of background events for \\alpha<%.0f\\circ", sigint)); | 
|---|
| 673 | hist->SetXTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 674 | hists->SetXTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 675 | histb->SetXTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 676 | hist->SetYTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 677 | hists->SetYTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 678 | histb->SetYTitle(fHist.GetXaxis()->GetTitle()); | 
|---|
| 679 |  | 
|---|
| 680 | const Double_t w = fHist.GetZaxis()->GetBinWidth(1); | 
|---|
| 681 |  | 
|---|
| 682 | //                      xmin, xmax, npar | 
|---|
| 683 | //TF1 func("MyFunc", fcn, 0, 90, 6); | 
|---|
| 684 | // Implementing the function yourself is only about 5% faster | 
|---|
| 685 | TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90); | 
|---|
| 686 | TArrayD maxpar(func.GetNpar()); | 
|---|
| 687 |  | 
|---|
| 688 | /* | 
|---|
| 689 | func.SetParName(0, "A"); | 
|---|
| 690 | func.SetParName(1, "mu"); | 
|---|
| 691 | func.SetParName(2, "sigma"); | 
|---|
| 692 | func.SetParName(3, "a"); | 
|---|
| 693 | func.SetParName(4, "b"); | 
|---|
| 694 | func.SetParName(5, "c"); | 
|---|
| 695 | */ | 
|---|
| 696 |  | 
|---|
| 697 | func.FixParameter(1, 0); | 
|---|
| 698 | func.FixParameter(4, 0); | 
|---|
| 699 | func.SetParLimits(2, 0, 90); | 
|---|
| 700 | func.SetParLimits(3, -1, 1); | 
|---|
| 701 |  | 
|---|
| 702 | const Int_t nx = fHist.GetXaxis()->GetNbins(); | 
|---|
| 703 | const Int_t ny = fHist.GetYaxis()->GetNbins(); | 
|---|
| 704 | const Int_t nr = nx*nx+ny*ny; | 
|---|
| 705 |  | 
|---|
| 706 | Double_t maxalpha0=0; | 
|---|
| 707 | Double_t maxs=3; | 
|---|
| 708 |  | 
|---|
| 709 | Int_t maxx=0; | 
|---|
| 710 | Int_t maxy=0; | 
|---|
| 711 |  | 
|---|
| 712 | TStopwatch clk; | 
|---|
| 713 | clk.Start(); | 
|---|
| 714 |  | 
|---|
| 715 | *fLog << inf; | 
|---|
| 716 | *fLog << "Signal fit:     alpha < " << sigmax << endl; | 
|---|
| 717 | *fLog << "Integration:    alpha < " << sigint << endl; | 
|---|
| 718 | *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl; | 
|---|
| 719 | *fLog << "Polynom order:  " << (int)polynom << endl; | 
|---|
| 720 | *fLog << "Fitting False Source Plot..." << flush; | 
|---|
| 721 |  | 
|---|
| 722 | TH1 *h=0; | 
|---|
| 723 | for (int ix=1; ix<=nx; ix++) | 
|---|
| 724 | for (int iy=1; iy<=ny; iy++) | 
|---|
| 725 | { | 
|---|
| 726 | h = fHist.ProjectionZ("AlphaFit", ix, ix, iy, iy); | 
|---|
| 727 |  | 
|---|
| 728 | const Double_t alpha0 = h->GetBinContent(1); | 
|---|
| 729 |  | 
|---|
| 730 | // Check for the regios which is not filled... | 
|---|
| 731 | if (alpha0==0) | 
|---|
| 732 | continue; | 
|---|
| 733 |  | 
|---|
| 734 | if (alpha0>maxalpha0) | 
|---|
| 735 | maxalpha0=alpha0; | 
|---|
| 736 |  | 
|---|
| 737 | // First fit a polynom in the off region | 
|---|
| 738 | func.FixParameter(0, 0); | 
|---|
| 739 | func.FixParameter(2, 1); | 
|---|
| 740 | func.ReleaseParameter(3); | 
|---|
| 741 |  | 
|---|
| 742 | for (int i=5; i<func.GetNpar(); i++) | 
|---|
| 743 | func.ReleaseParameter(i); | 
|---|
| 744 |  | 
|---|
| 745 | h->Fit(&func, "N0Q", "", bgmin, bgmax); | 
|---|
| 746 | //*fLog << dbg << ix << "/" << iy << ":  " << func.GetParameter(3) << "    " << func.GetParameter(4) << endl; | 
|---|
| 747 |  | 
|---|
| 748 | h4a.Fill(func.GetChisquare()); | 
|---|
| 749 | //h5a.Fill(func.GetProb()); | 
|---|
| 750 |  | 
|---|
| 751 | //func.SetParLimits(0, 0.5*h->GetBinContent(1), 1.5*h->GetBinContent(1)); | 
|---|
| 752 | //func.SetParLimits(2, 5, 90); | 
|---|
| 753 |  | 
|---|
| 754 | func.ReleaseParameter(0); | 
|---|
| 755 | //func.ReleaseParameter(1); | 
|---|
| 756 | func.ReleaseParameter(2); | 
|---|
| 757 | func.FixParameter(3, func.GetParameter(3)); | 
|---|
| 758 | for (int i=5; i<func.GetNpar(); i++) | 
|---|
| 759 | func.FixParameter(i, func.GetParameter(i)); | 
|---|
| 760 |  | 
|---|
| 761 | // Do not allow signals smaller than the background | 
|---|
| 762 | const Double_t A  = alpha0-func.GetParameter(3); | 
|---|
| 763 | const Double_t dA = TMath::Abs(A); | 
|---|
| 764 | func.SetParLimits(0, -dA*4, dA*4); | 
|---|
| 765 | func.SetParLimits(2, 0, 90); | 
|---|
| 766 |  | 
|---|
| 767 | // Now fit a gaus in the on region on top of the polynom | 
|---|
| 768 | func.SetParameter(0, A); | 
|---|
| 769 | func.SetParameter(2, sigmax*0.75); | 
|---|
| 770 |  | 
|---|
| 771 | h->Fit(&func, "N0Q", "", 0, sigmax); | 
|---|
| 772 | //*fLog << dbg << "     " << func.GetParameter(0) << "    " << func.GetParameter(1) << "    " << func.GetParameter(2) << endl; | 
|---|
| 773 |  | 
|---|
| 774 | TArrayD p(func.GetNpar(), func.GetParameters()); | 
|---|
| 775 |  | 
|---|
| 776 | // Fill results into some histograms | 
|---|
| 777 | h0a.Fill(p[0]); | 
|---|
| 778 | h0b.Fill(p[3]); | 
|---|
| 779 | h1.Fill(p[1]); | 
|---|
| 780 | h2.Fill(p[2]); | 
|---|
| 781 | if (polynom>1) | 
|---|
| 782 | h3.Fill(p[5]); | 
|---|
| 783 | h4b.Fill(func.GetChisquare()); | 
|---|
| 784 | //h5b.Fill(func.GetProb()); | 
|---|
| 785 |  | 
|---|
| 786 | // Implementing the integral as analytical function | 
|---|
| 787 | // gives the same result in the order of 10e-5 | 
|---|
| 788 | // and it is not faster at all... | 
|---|
| 789 | //const Bool_t ok = /*p[0]>=0 && /*p[0]<alpha0*2 &&*/ p[2]>1.75;// && p[2]<88.5; | 
|---|
| 790 | /* | 
|---|
| 791 | if (p[0]<-fabs(alpha0-p[3])*7 && p[2]<alphaw/3) | 
|---|
| 792 | { | 
|---|
| 793 | func.SetParameter(0, alpha0-p[3]); | 
|---|
| 794 | cout << p[2] << " " << p[0] << " " << alpha0-p[3] << endl; | 
|---|
| 795 | } | 
|---|
| 796 | */ | 
|---|
| 797 |  | 
|---|
| 798 | // The fitted function returned units of | 
|---|
| 799 | // counts bin binwidth. To get the correct number | 
|---|
| 800 | // of events we must adapt the functions by dividing | 
|---|
| 801 | // the result of the integration by the bin-width | 
|---|
| 802 | const Double_t s = func.Integral(0, sigint)/w; | 
|---|
| 803 |  | 
|---|
| 804 | func.SetParameter(0, 0); | 
|---|
| 805 | func.SetParameter(2, 1); | 
|---|
| 806 |  | 
|---|
| 807 | const Double_t b   = func.Integral(0, sigint)/w; | 
|---|
| 808 | const Double_t sig = Significance(s, b); | 
|---|
| 809 |  | 
|---|
| 810 | const Int_t n = hist->GetBin(ix, iy); | 
|---|
| 811 | hists->SetBinContent(n, s-b); | 
|---|
| 812 | histb->SetBinContent(n, b); | 
|---|
| 813 |  | 
|---|
| 814 | hist->SetBinContent(n, sig); | 
|---|
| 815 | if (sig!=0) | 
|---|
| 816 | h6.Fill(sig); | 
|---|
| 817 |  | 
|---|
| 818 | if (sig>maxs && ix*ix+iy*iy<nr*nr/9) | 
|---|
| 819 | { | 
|---|
| 820 | maxs = sig; | 
|---|
| 821 | maxx = ix; | 
|---|
| 822 | maxy = iy; | 
|---|
| 823 | maxpar = p; | 
|---|
| 824 | } | 
|---|
| 825 | } | 
|---|
| 826 |  | 
|---|
| 827 | *fLog << inf << "done." << endl; | 
|---|
| 828 |  | 
|---|
| 829 | h0a.GetXaxis()->SetRangeUser(0, maxalpha0*1.5); | 
|---|
| 830 | h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5); | 
|---|
| 831 |  | 
|---|
| 832 | //hists->SetMinimum(GetMinimumGT(*hists)); | 
|---|
| 833 | histb->SetMinimum(GetMinimumGT(*histb)); | 
|---|
| 834 |  | 
|---|
| 835 | clk.Stop(); | 
|---|
| 836 | clk.Print(); | 
|---|
| 837 |  | 
|---|
| 838 | TCanvas *c=new TCanvas; | 
|---|
| 839 |  | 
|---|
| 840 | c->Divide(3,2, 0, 0); | 
|---|
| 841 | c->cd(1); | 
|---|
| 842 | gPad->SetBorderMode(0); | 
|---|
| 843 | hists->Draw("colz"); | 
|---|
| 844 | hists->SetBit(kCanDelete); | 
|---|
| 845 | c->cd(2); | 
|---|
| 846 | gPad->SetBorderMode(0); | 
|---|
| 847 | hist->Draw("colz"); | 
|---|
| 848 | hist->SetBit(kCanDelete); | 
|---|
| 849 | c->cd(3); | 
|---|
| 850 | gPad->SetBorderMode(0); | 
|---|
| 851 | histb->Draw("colz"); | 
|---|
| 852 | histb->SetBit(kCanDelete); | 
|---|
| 853 | c->cd(4); | 
|---|
| 854 | gPad->Divide(1,3, 0, 0); | 
|---|
| 855 | TVirtualPad *p=gPad; | 
|---|
| 856 | p->SetBorderMode(0); | 
|---|
| 857 | p->cd(1); | 
|---|
| 858 | gPad->SetBorderMode(0); | 
|---|
| 859 | h0b.DrawCopy(); | 
|---|
| 860 | h0a.DrawCopy("same"); | 
|---|
| 861 | p->cd(2); | 
|---|
| 862 | gPad->SetBorderMode(0); | 
|---|
| 863 | h3.DrawCopy(); | 
|---|
| 864 | p->cd(3); | 
|---|
| 865 | gPad->SetBorderMode(0); | 
|---|
| 866 | h2.DrawCopy(); | 
|---|
| 867 | c->cd(6); | 
|---|
| 868 | gPad->Divide(1,2, 0, 0); | 
|---|
| 869 | TVirtualPad *q=gPad; | 
|---|
| 870 | q->SetBorderMode(0); | 
|---|
| 871 | q->cd(1); | 
|---|
| 872 | gPad->SetBorderMode(0); | 
|---|
| 873 | gPad->SetBorderMode(0); | 
|---|
| 874 | h4b.DrawCopy(); | 
|---|
| 875 | h4a.DrawCopy("same"); | 
|---|
| 876 | h6.DrawCopy("same"); | 
|---|
| 877 | q->cd(2); | 
|---|
| 878 | gPad->SetBorderMode(0); | 
|---|
| 879 | //h5b.DrawCopy(); | 
|---|
| 880 | //h5a.DrawCopy("same"); | 
|---|
| 881 | c->cd(5); | 
|---|
| 882 | gPad->SetBorderMode(0); | 
|---|
| 883 | if (maxx>0 && maxy>0) | 
|---|
| 884 | { | 
|---|
| 885 | const char *title = Form(" \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f) ", | 
|---|
| 886 | fHist.GetXaxis()->GetBinCenter(maxx), | 
|---|
| 887 | fHist.GetYaxis()->GetBinCenter(maxy), maxs); | 
|---|
| 888 |  | 
|---|
| 889 | TH1 *result = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy); | 
|---|
| 890 | result->SetDirectory(NULL); | 
|---|
| 891 | result->SetNameTitle("Result \\alpha", title); | 
|---|
| 892 | result->SetBit(kCanDelete); | 
|---|
| 893 | result->SetXTitle("\\alpha [\\circ]"); | 
|---|
| 894 | result->SetYTitle("Counts"); | 
|---|
| 895 | result->Draw(); | 
|---|
| 896 |  | 
|---|
| 897 | TF1 f1("", func.GetExpFormula(), 0, 90); | 
|---|
| 898 | TF1 f2("", func.GetExpFormula(), 0, 90); | 
|---|
| 899 | f1.SetParameters(maxpar.GetArray()); | 
|---|
| 900 | f2.SetParameters(maxpar.GetArray()); | 
|---|
| 901 | f2.FixParameter(0, 0); | 
|---|
| 902 | f2.FixParameter(1, 0); | 
|---|
| 903 | f2.FixParameter(2, 1); | 
|---|
| 904 | f1.SetLineColor(kGreen); | 
|---|
| 905 | f2.SetLineColor(kRed); | 
|---|
| 906 |  | 
|---|
| 907 | f2.DrawCopy("same"); | 
|---|
| 908 | f1.DrawCopy("same"); | 
|---|
| 909 |  | 
|---|
| 910 | TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC"); | 
|---|
| 911 | leg->SetBorderSize(1); | 
|---|
| 912 | leg->SetTextSize(0.04); | 
|---|
| 913 | leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22); | 
|---|
| 914 | //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22); | 
|---|
| 915 | leg->AddLine(0, 0.65, 0, 0.65); | 
|---|
| 916 | leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11); | 
|---|
| 917 | leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11); | 
|---|
| 918 | leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11); | 
|---|
| 919 | leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11); | 
|---|
| 920 | leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11); | 
|---|
| 921 | if (polynom>1) | 
|---|
| 922 | leg->AddText(0.60, 0.14, Form("c=%.2f", maxpar[5]))->SetTextAlign(11); | 
|---|
| 923 | leg->SetBit(kCanDelete); | 
|---|
| 924 | leg->Draw(); | 
|---|
| 925 |  | 
|---|
| 926 | q->cd(2); | 
|---|
| 927 |  | 
|---|
| 928 | TGraph *g = new TGraph; | 
|---|
| 929 | g->SetPoint(0, 0, 0); | 
|---|
| 930 |  | 
|---|
| 931 | Int_t max=0; | 
|---|
| 932 | Float_t maxsig=0; | 
|---|
| 933 | for (int i=1; i<89; i++) | 
|---|
| 934 | { | 
|---|
| 935 | const Double_t s = f1.Integral(0, (float)i); | 
|---|
| 936 | const Double_t b = f2.Integral(0, (float)i); | 
|---|
| 937 |  | 
|---|
| 938 | const Double_t sig = Significance(s, b); | 
|---|
| 939 |  | 
|---|
| 940 | g->SetPoint(g->GetN(), i, sig); | 
|---|
| 941 |  | 
|---|
| 942 | if (sig>maxsig) | 
|---|
| 943 | { | 
|---|
| 944 | max = i; | 
|---|
| 945 | maxsig = sig; | 
|---|
| 946 | } | 
|---|
| 947 | } | 
|---|
| 948 |  | 
|---|
| 949 | g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha"); | 
|---|
| 950 | g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]"); | 
|---|
| 951 | g->GetHistogram()->SetYTitle("Significance"); | 
|---|
| 952 | g->SetBit(kCanDelete); | 
|---|
| 953 | g->Draw("AP"); | 
|---|
| 954 |  | 
|---|
| 955 | leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC"); | 
|---|
| 956 | leg->SetBorderSize(1); | 
|---|
| 957 | leg->SetTextSize(0.1); | 
|---|
| 958 | leg->AddText(Form("\\sigma_{max}=%.1f at \\alpha_{max}=%d\\circ", maxsig, max)); | 
|---|
| 959 | leg->SetBit(kCanDelete); | 
|---|
| 960 | leg->Draw(); | 
|---|
| 961 | } | 
|---|
| 962 | } | 
|---|