| 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, 1/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
|
|---|
| 19 | ! Markus Gaug, 02/2004 <mailto:markus@ifae.es>
|
|---|
| 20 | !
|
|---|
| 21 | ! Copyright: MAGIC Software Development, 2000-2004
|
|---|
| 22 | !
|
|---|
| 23 | !
|
|---|
| 24 | \* ======================================================================== */
|
|---|
| 25 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 26 | //
|
|---|
| 27 | // MGCamDisplays
|
|---|
| 28 | //
|
|---|
| 29 | // Graphical interfaces to display the camera with fits and projections
|
|---|
| 30 | //
|
|---|
| 31 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 32 | #include "MGCamDisplays.h"
|
|---|
| 33 |
|
|---|
| 34 | #include <TStyle.h>
|
|---|
| 35 | #include <TCanvas.h>
|
|---|
| 36 |
|
|---|
| 37 | #include "MHCamera.h"
|
|---|
| 38 | #include "MGeomCam.h"
|
|---|
| 39 | #include "TPad.h"
|
|---|
| 40 | #include "TProfile.h"
|
|---|
| 41 | #include "TF1.h"
|
|---|
| 42 |
|
|---|
| 43 | #include "MStatusDisplay.h"
|
|---|
| 44 |
|
|---|
| 45 | ClassImp(MGCamDisplays);
|
|---|
| 46 |
|
|---|
| 47 | using namespace std;
|
|---|
| 48 |
|
|---|
| 49 | // --------------------------------------------------------------------------
|
|---|
| 50 | //
|
|---|
| 51 | // Default constructor.
|
|---|
| 52 | //
|
|---|
| 53 | MGCamDisplays::MGCamDisplays()
|
|---|
| 54 | {
|
|---|
| 55 | }
|
|---|
| 56 |
|
|---|
| 57 | // --------------------------------------------------------------------------
|
|---|
| 58 | //
|
|---|
| 59 | // Draw the MHCamera into the MStatusDisplay:
|
|---|
| 60 | //
|
|---|
| 61 | // 1) Draw it as histogram (MHCamera::DrawCopy("hist")
|
|---|
| 62 | // 2) Draw it as a camera, with MHCamera::SetPrettyPalette() set.
|
|---|
| 63 | // 3) If "rad" is not zero, draw its values vs. the radius from the camera center.
|
|---|
| 64 | // (DrawRadialProfile())
|
|---|
| 65 | // 4) Depending on the variable "fit", draw the values projection on the y-axis
|
|---|
| 66 | // (DrawProjection()):
|
|---|
| 67 | // 0: don't draw
|
|---|
| 68 | // 1: Draw fit to Single Gauss (for distributions flat-fielded over the whole camera)
|
|---|
| 69 | // 2: Draw and fit to Double Gauss (for distributions different for inner and outer pixels)
|
|---|
| 70 | // 3: Draw and fit to Triple Gauss (for distributions with inner, outer pixels and outliers)
|
|---|
| 71 | // 4: Draw and fit to Polynomial grade 0: (for the probability distributions)
|
|---|
| 72 | // >4: Draw and don;t fit.
|
|---|
| 73 | //
|
|---|
| 74 | void MGCamDisplays::CamDraw(TCanvas &c, const Int_t x, const Int_t y, const MHCamera &cam1,
|
|---|
| 75 | const Int_t fit, const Int_t rad)
|
|---|
| 76 | {
|
|---|
| 77 |
|
|---|
| 78 | c.cd(x);
|
|---|
| 79 | gPad->SetBorderMode(0);
|
|---|
| 80 | gPad->SetTicks();
|
|---|
| 81 | MHCamera *obj1=(MHCamera*)cam1.DrawCopy("hist");
|
|---|
| 82 | obj1->SetDirectory(NULL);
|
|---|
| 83 |
|
|---|
| 84 | c.cd(x+y);
|
|---|
| 85 | gPad->SetBorderMode(0);
|
|---|
| 86 | obj1->SetPrettyPalette();
|
|---|
| 87 | obj1->Draw();
|
|---|
| 88 |
|
|---|
| 89 | if (rad)
|
|---|
| 90 | {
|
|---|
| 91 | c.cd(x+2*y);
|
|---|
| 92 | gPad->SetBorderMode(0);
|
|---|
| 93 | gPad->SetTicks();
|
|---|
| 94 | DrawRadialProfile(obj1);
|
|---|
| 95 | }
|
|---|
| 96 |
|
|---|
| 97 |
|
|---|
| 98 | if (!fit)
|
|---|
| 99 | return;
|
|---|
| 100 |
|
|---|
| 101 | c.cd(rad ? x+3*y : x+2*y);
|
|---|
| 102 | gPad->SetBorderMode(0);
|
|---|
| 103 | gPad->SetTicks();
|
|---|
| 104 | DrawProjection(obj1, fit);
|
|---|
| 105 | }
|
|---|
| 106 |
|
|---|
| 107 | // --------------------------------------------------------------------------
|
|---|
| 108 | //
|
|---|
| 109 | // Draw a projection of MHCamera onto the y-axis values. Depending on the
|
|---|
| 110 | // variable fit, the following fits are performed:
|
|---|
| 111 | //
|
|---|
| 112 | // 1: Single Gauss (for distributions flat-fielded over the whole camera)
|
|---|
| 113 | // 2: Double Gauss (for distributions different for inner and outer pixels)
|
|---|
| 114 | // 3: Triple Gauss (for distributions with inner, outer pixels and outliers)
|
|---|
| 115 | // 4: flat (for the probability distributions)
|
|---|
| 116 | //
|
|---|
| 117 | // Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are
|
|---|
| 118 | // drawn separately, for inner and outer pixels.
|
|---|
| 119 | //
|
|---|
| 120 | void MGCamDisplays::DrawProjection(MHCamera *obj, Int_t fit) const
|
|---|
| 121 | {
|
|---|
| 122 |
|
|---|
| 123 | TH1D *obj2 = (TH1D*)obj->Projection(obj->GetName());
|
|---|
| 124 | obj2->SetDirectory(0);
|
|---|
| 125 | obj2->Draw();
|
|---|
| 126 | obj2->SetBit(kCanDelete);
|
|---|
| 127 |
|
|---|
| 128 | if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
|
|---|
| 129 | {
|
|---|
| 130 | TArrayI s0(3);
|
|---|
| 131 | s0[0] = 6;
|
|---|
| 132 | s0[1] = 1;
|
|---|
| 133 | s0[2] = 2;
|
|---|
| 134 |
|
|---|
| 135 | TArrayI s1(3);
|
|---|
| 136 | s1[0] = 3;
|
|---|
| 137 | s1[1] = 4;
|
|---|
| 138 | s1[2] = 5;
|
|---|
| 139 |
|
|---|
| 140 | TArrayI inner(1);
|
|---|
| 141 | inner[0] = 0;
|
|---|
| 142 |
|
|---|
| 143 | TArrayI outer(1);
|
|---|
| 144 | outer[0] = 1;
|
|---|
| 145 |
|
|---|
| 146 | // Just to get the right (maximum) binning
|
|---|
| 147 | TH1D *half[4];
|
|---|
| 148 | half[0] = obj->ProjectionS(s0, inner, "Sector 6-1-2 Inner");
|
|---|
| 149 | half[1] = obj->ProjectionS(s1, inner, "Sector 3-4-5 Inner");
|
|---|
| 150 | half[2] = obj->ProjectionS(s0, outer, "Sector 6-1-2 Outer");
|
|---|
| 151 | half[3] = obj->ProjectionS(s1, outer, "Sector 3-4-5 Outer");
|
|---|
| 152 |
|
|---|
| 153 | for (int i=0; i<4; i++)
|
|---|
| 154 | {
|
|---|
| 155 | half[i]->SetLineColor(kRed+i);
|
|---|
| 156 | half[i]->SetDirectory(0);
|
|---|
| 157 | half[i]->SetBit(kCanDelete);
|
|---|
| 158 | half[i]->Draw("same");
|
|---|
| 159 | }
|
|---|
| 160 | }
|
|---|
| 161 |
|
|---|
| 162 | const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
|
|---|
| 163 | const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
|
|---|
| 164 | const Double_t integ = obj2->Integral("width")/2.5;
|
|---|
| 165 | const Double_t mean = obj2->GetMean();
|
|---|
| 166 | const Double_t rms = obj2->GetRMS();
|
|---|
| 167 | const Double_t width = max-min;
|
|---|
| 168 |
|
|---|
| 169 | const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
|
|---|
| 170 | "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
|
|---|
| 171 |
|
|---|
| 172 | const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
|
|---|
| 173 | "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"
|
|---|
| 174 | "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
|
|---|
| 175 | TF1 *f=0;
|
|---|
| 176 | switch (fit)
|
|---|
| 177 | {
|
|---|
| 178 | case 1:
|
|---|
| 179 | f = new TF1("sgaus", "gaus(0)", min, max);
|
|---|
| 180 | f->SetLineColor(kYellow);
|
|---|
| 181 | f->SetBit(kCanDelete);
|
|---|
| 182 | f->SetParNames("Area", "#mu", "#sigma");
|
|---|
| 183 | f->SetParameters(integ/rms, mean, rms);
|
|---|
| 184 | f->SetParLimits(0, 0, integ);
|
|---|
| 185 | f->SetParLimits(1, min, max);
|
|---|
| 186 | f->SetParLimits(2, 0, width/1.5);
|
|---|
| 187 |
|
|---|
| 188 | obj2->Fit(f, "QLR");
|
|---|
| 189 | break;
|
|---|
| 190 |
|
|---|
| 191 | case 2:
|
|---|
| 192 | f = new TF1("dgaus",dgausformula.Data(),min,max);
|
|---|
| 193 | f->SetLineColor(kYellow);
|
|---|
| 194 | f->SetBit(kCanDelete);
|
|---|
| 195 | f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2");
|
|---|
| 196 | f->SetParameters(integ,(min+mean)/2.,width/4.,
|
|---|
| 197 | integ/width/2.,(max+mean)/2.,width/4.);
|
|---|
| 198 | // The left-sided Gauss
|
|---|
| 199 | f->SetParLimits(0,integ-1.5 , integ+1.5);
|
|---|
| 200 | f->SetParLimits(1,min+(width/10.), mean);
|
|---|
| 201 | f->SetParLimits(2,0 , width/2.);
|
|---|
| 202 | // The right-sided Gauss
|
|---|
| 203 | f->SetParLimits(3,0 , integ);
|
|---|
| 204 | f->SetParLimits(4,mean, max-(width/10.));
|
|---|
| 205 | f->SetParLimits(5,0 , width/2.);
|
|---|
| 206 | obj2->Fit(f,"QLRM");
|
|---|
| 207 | break;
|
|---|
| 208 |
|
|---|
| 209 | case 3:
|
|---|
| 210 | f = new TF1("tgaus",tgausformula.Data(),min,max);
|
|---|
| 211 | f->SetLineColor(kYellow);
|
|---|
| 212 | f->SetBit(kCanDelete);
|
|---|
| 213 | f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
|
|---|
| 214 | "A_{2}","#mu_{2}","#sigma_{2}",
|
|---|
| 215 | "A_{3}","#mu_{3}","#sigma_{3}");
|
|---|
| 216 | f->SetParameters(integ,(min+mean)/2,width/4.,
|
|---|
| 217 | integ/width/3.,(max+mean)/2.,width/4.,
|
|---|
| 218 | integ/width/3.,mean,width/2.);
|
|---|
| 219 | // The left-sided Gauss
|
|---|
| 220 | f->SetParLimits(0,integ-1.5,integ+1.5);
|
|---|
| 221 | f->SetParLimits(1,min+(width/10.),mean);
|
|---|
| 222 | f->SetParLimits(2,width/15.,width/2.);
|
|---|
| 223 | // The right-sided Gauss
|
|---|
| 224 | f->SetParLimits(3,0.,integ);
|
|---|
| 225 | f->SetParLimits(4,mean,max-(width/10.));
|
|---|
| 226 | f->SetParLimits(5,width/15.,width/2.);
|
|---|
| 227 | // The Gauss describing the outliers
|
|---|
| 228 | f->SetParLimits(6,0.,integ);
|
|---|
| 229 | f->SetParLimits(7,min,max);
|
|---|
| 230 | f->SetParLimits(8,width/4.,width/1.5);
|
|---|
| 231 | obj2->Fit(f,"QLRM");
|
|---|
| 232 | break;
|
|---|
| 233 |
|
|---|
| 234 | case 4:
|
|---|
| 235 | obj2->Fit("pol0", "Q");
|
|---|
| 236 | obj2->GetFunction("pol0")->SetLineColor(kYellow);
|
|---|
| 237 | break;
|
|---|
| 238 |
|
|---|
| 239 | case 9:
|
|---|
| 240 | break;
|
|---|
| 241 |
|
|---|
| 242 | default:
|
|---|
| 243 | obj2->Fit("gaus", "Q");
|
|---|
| 244 | obj2->GetFunction("gaus")->SetLineColor(kYellow);
|
|---|
| 245 | break;
|
|---|
| 246 | }
|
|---|
| 247 | }
|
|---|
| 248 |
|
|---|
| 249 | // --------------------------------------------------------------------------
|
|---|
| 250 | //
|
|---|
| 251 | // Draw a projection of MHCamera vs. the radius from the central pixel.
|
|---|
| 252 | //
|
|---|
| 253 | // The inner and outer pixels are drawn separately, both fitted by a polynomial
|
|---|
| 254 | // of grade 1.
|
|---|
| 255 | //
|
|---|
| 256 | void MGCamDisplays::DrawRadialProfile(MHCamera *obj) const
|
|---|
| 257 | {
|
|---|
| 258 |
|
|---|
| 259 | TProfile *obj2 = (TProfile*)obj->RadialProfile(obj->GetName());
|
|---|
| 260 | obj2->SetDirectory(0);
|
|---|
| 261 | obj2->Draw();
|
|---|
| 262 | obj2->SetBit(kCanDelete);
|
|---|
| 263 |
|
|---|
| 264 | if (obj->GetGeomCam().InheritsFrom("MGeomCamMagic"))
|
|---|
| 265 | {
|
|---|
| 266 |
|
|---|
| 267 | TArrayI s0(6);
|
|---|
| 268 | s0[0] = 1;
|
|---|
| 269 | s0[1] = 2;
|
|---|
| 270 | s0[2] = 3;
|
|---|
| 271 | s0[3] = 4;
|
|---|
| 272 | s0[4] = 5;
|
|---|
| 273 | s0[5] = 6;
|
|---|
| 274 |
|
|---|
| 275 | TArrayI inner(1);
|
|---|
| 276 | inner[0] = 0;
|
|---|
| 277 |
|
|---|
| 278 | TArrayI outer(1);
|
|---|
| 279 | outer[0] = 1;
|
|---|
| 280 |
|
|---|
| 281 | // Just to get the right (maximum) binning
|
|---|
| 282 | TProfile *half[2];
|
|---|
| 283 | half[0] = obj->RadialProfileS(s0, inner,Form("%s%s",obj->GetName(),"Inner"));
|
|---|
| 284 | half[1] = obj->RadialProfileS(s0, outer,Form("%s%s",obj->GetName(),"Outer"));
|
|---|
| 285 |
|
|---|
| 286 | for (Int_t i=0; i<2; i++)
|
|---|
| 287 | {
|
|---|
| 288 | Double_t min = obj->GetGeomCam().GetMinRadius(i);
|
|---|
| 289 | Double_t max = obj->GetGeomCam().GetMaxRadius(i);
|
|---|
| 290 |
|
|---|
| 291 | half[i]->SetLineColor(kRed+i);
|
|---|
| 292 | half[i]->SetDirectory(0);
|
|---|
| 293 | half[i]->SetBit(kCanDelete);
|
|---|
| 294 | half[i]->Draw("same");
|
|---|
| 295 | half[i]->Fit("pol1","Q","",min,max);
|
|---|
| 296 | half[i]->GetFunction("pol1")->SetLineColor(kRed+i);
|
|---|
| 297 | half[i]->GetFunction("pol1")->SetLineWidth(1);
|
|---|
| 298 | }
|
|---|
| 299 | }
|
|---|
| 300 | }
|
|---|
| 301 |
|
|---|