| 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 12/2000 <mailto:tbretz@uni-sw.gwdg.de>
|
|---|
| 19 | ! Author(s): Harald Kornmayer 1/2001
|
|---|
| 20 | ! Author(s): Rudolf Bock 10/2001 <mailto:Rudolf.Bock@cern.ch>
|
|---|
| 21 | ! Author(s): Wolfgang Wittek 06/2002 <mailto:wittek@mppmu.mpg.de>
|
|---|
| 22 | !
|
|---|
| 23 | ! Copyright: MAGIC Software Development, 2000-2002
|
|---|
| 24 | !
|
|---|
| 25 | !
|
|---|
| 26 | \* ======================================================================== */
|
|---|
| 27 |
|
|---|
| 28 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 29 | //
|
|---|
| 30 | // MHillasSrc
|
|---|
| 31 | //
|
|---|
| 32 | // Storage Container for image parameters
|
|---|
| 33 | //
|
|---|
| 34 | // source-dependent image parameters
|
|---|
| 35 | //
|
|---|
| 36 | //
|
|---|
| 37 | // Version 1:
|
|---|
| 38 | // ----------
|
|---|
| 39 | // fAlpha angle between major axis and line source-to-center
|
|---|
| 40 | // fDist distance from source to center of ellipse
|
|---|
| 41 | //
|
|---|
| 42 | //
|
|---|
| 43 | // Version 2:
|
|---|
| 44 | // ----------
|
|---|
| 45 | // fHeadTail added
|
|---|
| 46 | //
|
|---|
| 47 | //
|
|---|
| 48 | // Version 3:
|
|---|
| 49 | // ----------
|
|---|
| 50 | // fCosDeltaAlpha cosine of angle between d and a, where
|
|---|
| 51 | // - d is the vector from the source position to the
|
|---|
| 52 | // center of the ellipse
|
|---|
| 53 | // - a is a vector along the main axis of the ellipse,
|
|---|
| 54 | // defined with positive x-component
|
|---|
| 55 | //
|
|---|
| 56 | //
|
|---|
| 57 | // Version 4:
|
|---|
| 58 | // ----------
|
|---|
| 59 | // fHeadTail removed
|
|---|
| 60 | //
|
|---|
| 61 | //
|
|---|
| 62 | // Version 5:
|
|---|
| 63 | // ----------
|
|---|
| 64 | // - added Float_t fDCA; // [mm] Distance to closest approach 'DCA'
|
|---|
| 65 | // - added Float_t fDCADelta; // [deg] Angle of the shower axis with respect
|
|---|
| 66 | // to the x-axis [0,2pi]
|
|---|
| 67 | //
|
|---|
| 68 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 69 | #include "MHillasSrc.h"
|
|---|
| 70 |
|
|---|
| 71 | #include <TArrayF.h>
|
|---|
| 72 |
|
|---|
| 73 | #include <TLine.h>
|
|---|
| 74 | #include <TMarker.h>
|
|---|
| 75 |
|
|---|
| 76 | #include "MLog.h"
|
|---|
| 77 | #include "MLogManip.h"
|
|---|
| 78 |
|
|---|
| 79 | #include "MGeomCam.h"
|
|---|
| 80 | #include "MHillas.h"
|
|---|
| 81 | #include "MSrcPosCam.h"
|
|---|
| 82 |
|
|---|
| 83 | ClassImp(MHillasSrc);
|
|---|
| 84 |
|
|---|
| 85 | using namespace std;
|
|---|
| 86 |
|
|---|
| 87 | // --------------------------------------------------------------------------
|
|---|
| 88 | //
|
|---|
| 89 | // Default constructor.
|
|---|
| 90 | //
|
|---|
| 91 | MHillasSrc::MHillasSrc(const char *name, const char *title)
|
|---|
| 92 | {
|
|---|
| 93 | fName = name ? name : "MHillasSrc";
|
|---|
| 94 | fTitle = title ? title : "Parameters depending in source position";
|
|---|
| 95 | }
|
|---|
| 96 |
|
|---|
| 97 | void MHillasSrc::Reset()
|
|---|
| 98 | {
|
|---|
| 99 | fDist = -1;
|
|---|
| 100 | fAlpha = 0;
|
|---|
| 101 | fCosDeltaAlpha = 0;
|
|---|
| 102 |
|
|---|
| 103 | fDCA = -1;
|
|---|
| 104 | fDCADelta = 0;
|
|---|
| 105 | }
|
|---|
| 106 |
|
|---|
| 107 | // --------------------------------------------------------------------------
|
|---|
| 108 | //
|
|---|
| 109 | // Calculation of source-dependent parameters
|
|---|
| 110 | // In case you don't call Calc from within an eventloop make sure, that
|
|---|
| 111 | // you call the Reset member function before.
|
|---|
| 112 | //
|
|---|
| 113 | Int_t MHillasSrc::Calc(const MHillas &hillas)
|
|---|
| 114 | {
|
|---|
| 115 | const Double_t mx = hillas.GetMeanX(); // [mm]
|
|---|
| 116 | const Double_t my = hillas.GetMeanY(); // [mm]
|
|---|
| 117 |
|
|---|
| 118 | const Double_t sx = mx - fSrcPos->GetX(); // [mm]
|
|---|
| 119 | const Double_t sy = my - fSrcPos->GetY(); // [mm]
|
|---|
| 120 |
|
|---|
| 121 | const Double_t sd = hillas.GetSinDelta(); // [1]
|
|---|
| 122 | const Double_t cd = hillas.GetCosDelta(); // [1]
|
|---|
| 123 |
|
|---|
| 124 | //
|
|---|
| 125 | // Distance from source position to center of ellipse.
|
|---|
| 126 | // If the distance is 0 distance, Alpha is not specified.
|
|---|
| 127 | // The calculation has failed and returnes kFALSE.
|
|---|
| 128 | //
|
|---|
| 129 | const Double_t dist = TMath::Sqrt(sx*sx + sy*sy); // [mm]
|
|---|
| 130 | if (dist==0)
|
|---|
| 131 | return 1;
|
|---|
| 132 |
|
|---|
| 133 | //
|
|---|
| 134 | // Calculate Alpha and Cosda = cos(d,a)
|
|---|
| 135 | // The sign of Cosda will be used for quantities containing
|
|---|
| 136 | // a head-tail information
|
|---|
| 137 | //
|
|---|
| 138 | // *OLD* const Double_t arg = (sy-tand*sx) / (dist*sqrt(tand*tand+1));
|
|---|
| 139 | // *OLD* fAlpha = asin(arg)*kRad2Deg;
|
|---|
| 140 | //
|
|---|
| 141 | const Double_t arg2 = cd*sx + sd*sy; // [mm]
|
|---|
| 142 | if (arg2==0)
|
|---|
| 143 | return 2;
|
|---|
| 144 |
|
|---|
| 145 | const Double_t arg1 = cd*sy - sd*sx; // [mm]
|
|---|
| 146 |
|
|---|
| 147 | //
|
|---|
| 148 | // Due to numerical uncertanties in the calculation of the
|
|---|
| 149 | // square root (dist) and arg1 it can happen (in less than 1e-5 cases)
|
|---|
| 150 | // that the absolute value of arg exceeds 1. Because this uncertainty
|
|---|
| 151 | // results in an Delta Alpha which is still less than 1e-3 we don't care
|
|---|
| 152 | // about this uncertainty in general and simply set values which exceed
|
|---|
| 153 | // to 1 saving its sign.
|
|---|
| 154 | //
|
|---|
| 155 | const Double_t arg = arg1/dist;
|
|---|
| 156 | fAlpha = TMath::Abs(arg)>1 ? TMath::Sign(90., arg) : TMath::ASin(arg)*TMath::RadToDeg(); // [deg]
|
|---|
| 157 |
|
|---|
| 158 | fCosDeltaAlpha = arg2/dist; // [1]
|
|---|
| 159 | fDist = dist; // [mm]
|
|---|
| 160 |
|
|---|
| 161 | // ----- Calculation of Distance to closest approach 'DCA' -----
|
|---|
| 162 |
|
|---|
| 163 | // Components of DCA vector
|
|---|
| 164 | const Double_t fpd1 = sx - arg2*cd; // [mm]
|
|---|
| 165 | const Double_t fpd2 = sy - arg2*sd; // [mm]
|
|---|
| 166 |
|
|---|
| 167 | // Determine the correct sign of the DCA (cross product of DCA vector and the
|
|---|
| 168 | // vector going from the intersection point of the DCA vector with the shower axis
|
|---|
| 169 | // to the COG)
|
|---|
| 170 | const Double_t sign = arg2*cd*fpd2 - arg2*sd*fpd1;
|
|---|
| 171 | fDCA = TMath::Sign(TMath::Sqrt(fpd1*fpd1 + fpd2*fpd2), sign); // [mm]
|
|---|
| 172 |
|
|---|
| 173 | // Calculate angle of the shower axis with respect to the x-axis
|
|---|
| 174 | fDCADelta = TMath::ACos((sx-fpd1)/TMath::Abs(arg2))*TMath::RadToDeg(); // [deg]
|
|---|
| 175 |
|
|---|
| 176 | // Enlarge the interval of fDdca to [0, 2pi]
|
|---|
| 177 | if (sy < fpd2)
|
|---|
| 178 | fDCADelta = 360 - fDCADelta;
|
|---|
| 179 |
|
|---|
| 180 | SetReadyToSave();
|
|---|
| 181 |
|
|---|
| 182 | return 0;
|
|---|
| 183 | }
|
|---|
| 184 |
|
|---|
| 185 | void MHillasSrc::Paint(Option_t *opt)
|
|---|
| 186 | {
|
|---|
| 187 | const Float_t x = fSrcPos ? fSrcPos->GetX() : 0;
|
|---|
| 188 | const Float_t y = fSrcPos ? fSrcPos->GetY() : 0;
|
|---|
| 189 |
|
|---|
| 190 | TMarker m;
|
|---|
| 191 | m.SetMarkerColor(kYellow);
|
|---|
| 192 | m.SetMarkerStyle(kStar);
|
|---|
| 193 | m.PaintMarker(x, y);
|
|---|
| 194 | /*
|
|---|
| 195 | m.SetMarkerColor(kMagenta);
|
|---|
| 196 | m.PaintMarker(fDist*cos((fDCADelta-fAlpha)*TMath::DegToRad()),
|
|---|
| 197 | fDist*sin((fDCADelta-fAlpha)*TMath::DegToRad()));
|
|---|
| 198 |
|
|---|
| 199 | TLine line;
|
|---|
| 200 | line.SetLineWidth(1);
|
|---|
| 201 | line.SetLineColor(108);
|
|---|
| 202 |
|
|---|
| 203 | line.PaintLine(-radius, y, radius, y);
|
|---|
| 204 | line.PaintLine(x, -radius, x, radius);
|
|---|
| 205 |
|
|---|
| 206 | // COG line
|
|---|
| 207 | TLine line(x, y, meanX, meanY);
|
|---|
| 208 | line.SetLineWidth(1);
|
|---|
| 209 | line.SetLineColor(2);
|
|---|
| 210 | line.Paint();*/
|
|---|
| 211 | }
|
|---|
| 212 |
|
|---|
| 213 | // --------------------------------------------------------------------------
|
|---|
| 214 | //
|
|---|
| 215 | // Print contents of MHillasSrc to *fLog
|
|---|
| 216 | //
|
|---|
| 217 | void MHillasSrc::Print(Option_t *) const
|
|---|
| 218 | {
|
|---|
| 219 | *fLog << all;
|
|---|
| 220 | *fLog << GetDescriptor() << endl;
|
|---|
| 221 | *fLog << " - Dist [mm] = " << fDist << endl;
|
|---|
| 222 | *fLog << " - Alpha [deg] = " << fAlpha << endl;
|
|---|
| 223 | *fLog << " - CosDeltaAlpha = " << fCosDeltaAlpha << endl;
|
|---|
| 224 | *fLog << " - DCA [mm] = " << fDCA << endl;
|
|---|
| 225 | *fLog << " - DCA delta [deg] = " << fDCADelta << endl;
|
|---|
| 226 | }
|
|---|
| 227 |
|
|---|
| 228 | // --------------------------------------------------------------------------
|
|---|
| 229 | //
|
|---|
| 230 | // Print contents of MHillasSrc to *fLog depending on the geometry in
|
|---|
| 231 | // units of deg.
|
|---|
| 232 | //
|
|---|
| 233 | void MHillasSrc::Print(const MGeomCam &geom) const
|
|---|
| 234 | {
|
|---|
| 235 | *fLog << all;
|
|---|
| 236 | *fLog << GetDescriptor() << endl;
|
|---|
| 237 | *fLog << " - Dist [deg] = " << fDist*geom.GetConvMm2Deg() << endl;
|
|---|
| 238 | *fLog << " - Alpha [deg] = " << fAlpha << endl;
|
|---|
| 239 | *fLog << " - CosDeltaAlpha = " << fCosDeltaAlpha << endl;
|
|---|
| 240 | *fLog << " - DCA [deg] = " << fDCA*geom.GetConvMm2Deg() << endl;
|
|---|
| 241 | *fLog << " - DCA delta [deg] = " << fDCADelta << endl;
|
|---|
| 242 | }
|
|---|
| 243 |
|
|---|
| 244 | // --------------------------------------------------------------------------
|
|---|
| 245 | //
|
|---|
| 246 | // This function is ment for special usage, please never try to set
|
|---|
| 247 | // values via this function
|
|---|
| 248 | //
|
|---|
| 249 | void MHillasSrc::Set(const TArrayF &arr)
|
|---|
| 250 | {
|
|---|
| 251 | if (arr.GetSize() != 5)
|
|---|
| 252 | return;
|
|---|
| 253 |
|
|---|
| 254 | fAlpha = arr.At(0); // [deg] angle of major axis with vector to src
|
|---|
| 255 | fDist = arr.At(1); // [mm] distance between src and center of ellipse
|
|---|
| 256 | fCosDeltaAlpha = arr.At(2); // [1] cosine of angle between d and a
|
|---|
| 257 | fDCA = arr.At(3); // [mm]
|
|---|
| 258 | fDCADelta = arr.At(4); // [mm]
|
|---|
| 259 | }
|
|---|