| 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): Wolfgang Wittek 2003 <mailto:wittek@mppmu.mpg.de>
|
|---|
| 19 | !
|
|---|
| 20 | ! Copyright: MAGIC Software Development, 2000-2003
|
|---|
| 21 | !
|
|---|
| 22 | !
|
|---|
| 23 | \* ======================================================================== */
|
|---|
| 24 |
|
|---|
| 25 | ///////////////////////////////////////////////////////////////////////
|
|---|
| 26 | //
|
|---|
| 27 | // MHCT1Supercuts
|
|---|
| 28 | //
|
|---|
| 29 | // This class contains histograms for the supercuts
|
|---|
| 30 | //
|
|---|
| 31 | // the histograms are filled during the optimization of the supercuts
|
|---|
| 32 | //
|
|---|
| 33 | ///////////////////////////////////////////////////////////////////////
|
|---|
| 34 | #include "MHCT1Supercuts.h"
|
|---|
| 35 |
|
|---|
| 36 | #include <math.h>
|
|---|
| 37 |
|
|---|
| 38 | #include <TH1.h>
|
|---|
| 39 | #include <TH2.h>
|
|---|
| 40 | #include <TPad.h>
|
|---|
| 41 | #include <TCanvas.h>
|
|---|
| 42 |
|
|---|
| 43 | #include "MParList.h"
|
|---|
| 44 | #include "MHFindSignificance.h"
|
|---|
| 45 |
|
|---|
| 46 | #include "MLog.h"
|
|---|
| 47 | #include "MLogManip.h"
|
|---|
| 48 |
|
|---|
| 49 |
|
|---|
| 50 |
|
|---|
| 51 | ClassImp(MHCT1Supercuts);
|
|---|
| 52 |
|
|---|
| 53 | using namespace std;
|
|---|
| 54 |
|
|---|
| 55 | // --------------------------------------------------------------------------
|
|---|
| 56 | //
|
|---|
| 57 | // Setup four histograms for Alpha, and Dist
|
|---|
| 58 | //
|
|---|
| 59 | MHCT1Supercuts::MHCT1Supercuts(const char *name, const char *title)
|
|---|
| 60 | {
|
|---|
| 61 | //
|
|---|
| 62 | // set the name and title of this object
|
|---|
| 63 | //
|
|---|
| 64 | fName = name ? name : "MHCT1Supercuts";
|
|---|
| 65 | fTitle = title ? title : "Container for supercuts";
|
|---|
| 66 |
|
|---|
| 67 |
|
|---|
| 68 | fDegree = new TH1F("Degree", "Degree of polynomial", 5, -0.5, 4.5);
|
|---|
| 69 | fProb = new TH1F("Prob", "chi2 probability", 40, 0, 1.0);
|
|---|
| 70 | fNdf = new TH1F("NDF", "NDF of polynomial fit", 60, -0.5, 59.5);
|
|---|
| 71 | fGamma = new TH1F("gamma", "gamma", 40, 0.0, 8.0);
|
|---|
| 72 | fNexNon = new TH1F("NexNon", "Nex / Non", 50, 0.0, 1.0);
|
|---|
| 73 | fSigLiMa= new TH1F("Significance", "significance of gamma signal",
|
|---|
| 74 | 60, 0.0, 120.0);
|
|---|
| 75 | fSigtoBackg= new TH2F("SigtoBackg", "Significance vs signal/backg ratio",
|
|---|
| 76 | 50, 0.0, 10.0, 60, 0.0, 120.0);
|
|---|
| 77 | fSigDegree = new TH2F("SigDegree", "Significance vs Degree of polynomial",
|
|---|
| 78 | 5, -0.5, 4.5, 60, 0.0, 120.0);
|
|---|
| 79 | fSigNbins = new TH2F("SigNbins", "Significance vs number of bins",
|
|---|
| 80 | 40, -0.5, 79.5, 60, 0.0, 120.0);
|
|---|
| 81 |
|
|---|
| 82 | fDegree->SetDirectory(NULL);
|
|---|
| 83 | fProb->SetDirectory(NULL);
|
|---|
| 84 | fNdf->SetDirectory(NULL);
|
|---|
| 85 | fGamma->SetDirectory(NULL);
|
|---|
| 86 | fNexNon->SetDirectory(NULL);
|
|---|
| 87 | fSigLiMa->SetDirectory(NULL);
|
|---|
| 88 | fSigtoBackg->SetDirectory(NULL);
|
|---|
| 89 | fSigDegree->SetDirectory(NULL);
|
|---|
| 90 | fSigNbins->SetDirectory(NULL);
|
|---|
| 91 |
|
|---|
| 92 | fDegree->SetXTitle("order of polynomial");
|
|---|
| 93 | fProb->SetXTitle("chi2 probability of polynomial fit");
|
|---|
| 94 | fNdf->SetXTitle("NDF of polynomial fit");
|
|---|
| 95 | fGamma->SetXTitle("gamma");
|
|---|
| 96 | fNexNon->SetXTitle("Nex / Non");
|
|---|
| 97 | fSigLiMa->SetXTitle("significance");
|
|---|
| 98 |
|
|---|
| 99 | fSigtoBackg->SetXTitle("signa./background ratio");
|
|---|
| 100 | fSigtoBackg->SetYTitle("significance");
|
|---|
| 101 |
|
|---|
| 102 | fSigDegree->SetXTitle("order of polynomial");
|
|---|
| 103 | fSigDegree->SetYTitle("significance");
|
|---|
| 104 |
|
|---|
| 105 | fSigNbins->SetXTitle("number of bins");
|
|---|
| 106 | fSigNbins->SetYTitle("significance");
|
|---|
| 107 |
|
|---|
| 108 | fDegree->SetYTitle("Counts");
|
|---|
| 109 | fProb->SetYTitle("Counts");
|
|---|
| 110 | fNdf->SetYTitle("Counts");
|
|---|
| 111 | fGamma->SetYTitle("Counts");
|
|---|
| 112 | fNexNon->SetYTitle("Counts");
|
|---|
| 113 | fSigLiMa->SetYTitle("Counts");
|
|---|
| 114 |
|
|---|
| 115 | fSigtoBackg->SetZTitle("Counts");
|
|---|
| 116 | fSigDegree->SetZTitle("Counts");
|
|---|
| 117 | fSigNbins->SetZTitle("Counts");
|
|---|
| 118 | }
|
|---|
| 119 |
|
|---|
| 120 | // --------------------------------------------------------------------------
|
|---|
| 121 | //
|
|---|
| 122 | // Delete the histograms
|
|---|
| 123 | //
|
|---|
| 124 | MHCT1Supercuts::~MHCT1Supercuts()
|
|---|
| 125 | {
|
|---|
| 126 | delete fDegree;
|
|---|
| 127 | delete fProb;
|
|---|
| 128 | delete fNdf;
|
|---|
| 129 | delete fGamma;
|
|---|
| 130 | delete fNexNon;
|
|---|
| 131 | delete fSigLiMa;
|
|---|
| 132 |
|
|---|
| 133 | delete fSigtoBackg;
|
|---|
| 134 | delete fSigDegree;
|
|---|
| 135 | delete fSigNbins;
|
|---|
| 136 | }
|
|---|
| 137 |
|
|---|
| 138 | // --------------------------------------------------------------------------
|
|---|
| 139 | //
|
|---|
| 140 | // Set the pointers (if there are any)
|
|---|
| 141 | //
|
|---|
| 142 |
|
|---|
| 143 | Bool_t MHCT1Supercuts::SetupFill(const MParList *plist)
|
|---|
| 144 | {
|
|---|
| 145 | return kTRUE;
|
|---|
| 146 | }
|
|---|
| 147 |
|
|---|
| 148 | // --------------------------------------------------------------------------
|
|---|
| 149 | //
|
|---|
| 150 | // Fill the histograms from the 'MHFindSignificance' container
|
|---|
| 151 | //
|
|---|
| 152 |
|
|---|
| 153 | Bool_t MHCT1Supercuts::Fill(const MParContainer *par, const Stat_t w)
|
|---|
| 154 | {
|
|---|
| 155 | if (!par)
|
|---|
| 156 | {
|
|---|
| 157 | *fLog << err << "MHCT1Supercuts::Fill: Pointer (!=NULL) expected." << endl;
|
|---|
| 158 | return kFALSE;
|
|---|
| 159 | }
|
|---|
| 160 |
|
|---|
| 161 | MHFindSignificance &h = *(MHFindSignificance*)par;
|
|---|
| 162 |
|
|---|
| 163 | fDegree ->Fill(h.GetDegree( ), w);
|
|---|
| 164 | fProb ->Fill(h.GetProb(), w);
|
|---|
| 165 | fNdf ->Fill(h.GetNdf(), w);
|
|---|
| 166 | fGamma ->Fill(h.GetGamma(), w);
|
|---|
| 167 |
|
|---|
| 168 | Double_t ratio = h.GetNon()>0.0 ? h.GetNex()/h.GetNon() : 0.0;
|
|---|
| 169 | fNexNon ->Fill(ratio, w);
|
|---|
| 170 |
|
|---|
| 171 | fSigLiMa ->Fill(h.GetSigLiMa(), w);
|
|---|
| 172 |
|
|---|
| 173 | Double_t sigtobackg = h.GetNbg()!=0.0 ? h.GetNex() / h.GetNbg() : 0.0;
|
|---|
| 174 | fSigtoBackg->Fill(sigtobackg, h.GetSigLiMa(), w);
|
|---|
| 175 |
|
|---|
| 176 | fSigDegree ->Fill(h.GetDegree(), h.GetSigLiMa(), w);
|
|---|
| 177 | fSigNbins ->Fill(h.GetMbins(), h.GetSigLiMa(), w);
|
|---|
| 178 |
|
|---|
| 179 | return kTRUE;
|
|---|
| 180 | }
|
|---|
| 181 |
|
|---|
| 182 |
|
|---|
| 183 | // --------------------------------------------------------------------------
|
|---|
| 184 | //
|
|---|
| 185 | // Creates a new canvas and draws the two histograms into it.
|
|---|
| 186 | // Be careful: The histograms belongs to this object and won't get deleted
|
|---|
| 187 | // together with the canvas.
|
|---|
| 188 | //
|
|---|
| 189 | void MHCT1Supercuts::Draw(Option_t *)
|
|---|
| 190 | {
|
|---|
| 191 | //TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
|---|
| 192 | //pad->SetBorderMode(0);
|
|---|
| 193 | //AppendPad("");
|
|---|
| 194 |
|
|---|
| 195 | TCanvas *pad = new TCanvas("Supercuts", "Supercut plots", 900, 900);
|
|---|
| 196 | gROOT->SetSelectedPad(NULL);
|
|---|
| 197 |
|
|---|
| 198 | pad->Divide(3, 3);
|
|---|
| 199 |
|
|---|
| 200 | pad->cd(1);
|
|---|
| 201 | gPad->SetBorderMode(0);
|
|---|
| 202 | fDegree->Draw();
|
|---|
| 203 |
|
|---|
| 204 | pad->cd(2);
|
|---|
| 205 | gPad->SetBorderMode(0);
|
|---|
| 206 | fProb->Draw();
|
|---|
| 207 |
|
|---|
| 208 | pad->cd(3);
|
|---|
| 209 | gPad->SetBorderMode(0);
|
|---|
| 210 | fNdf->Draw();
|
|---|
| 211 |
|
|---|
| 212 | pad->cd(4);
|
|---|
| 213 | gPad->SetBorderMode(0);
|
|---|
| 214 | fGamma->Draw();
|
|---|
| 215 |
|
|---|
| 216 | pad->cd(5);
|
|---|
| 217 | gPad->SetBorderMode(0);
|
|---|
| 218 | fNexNon->Draw();
|
|---|
| 219 |
|
|---|
| 220 | pad->cd(6);
|
|---|
| 221 | gPad->SetBorderMode(0);
|
|---|
| 222 | fSigLiMa->Draw();
|
|---|
| 223 |
|
|---|
| 224 | pad->cd(7);
|
|---|
| 225 | gPad->SetBorderMode(0);
|
|---|
| 226 | fSigtoBackg->Draw();
|
|---|
| 227 |
|
|---|
| 228 | pad->cd(8);
|
|---|
| 229 | gPad->SetBorderMode(0);
|
|---|
| 230 | fSigDegree->Draw();
|
|---|
| 231 |
|
|---|
| 232 | pad->cd(9);
|
|---|
| 233 | gPad->SetBorderMode(0);
|
|---|
| 234 | fSigNbins->Draw();
|
|---|
| 235 |
|
|---|
| 236 | pad->Modified();
|
|---|
| 237 | pad->Update();
|
|---|
| 238 | }
|
|---|
| 239 | //=========================================================================
|
|---|
| 240 |
|
|---|
| 241 |
|
|---|
| 242 |
|
|---|
| 243 |
|
|---|
| 244 |
|
|---|
| 245 |
|
|---|
| 246 |
|
|---|
| 247 |
|
|---|
| 248 |
|
|---|