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 "MSrcPosCam.h"
|
---|
81 |
|
---|
82 | ClassImp(MHillasSrc);
|
---|
83 |
|
---|
84 | using namespace std;
|
---|
85 |
|
---|
86 | // --------------------------------------------------------------------------
|
---|
87 | //
|
---|
88 | // Default constructor.
|
---|
89 | //
|
---|
90 | MHillasSrc::MHillasSrc(const char *name, const char *title)
|
---|
91 | {
|
---|
92 | fName = name ? name : "MHillasSrc";
|
---|
93 | fTitle = title ? title : "Parameters depending in source position";
|
---|
94 | }
|
---|
95 |
|
---|
96 | void MHillasSrc::Reset()
|
---|
97 | {
|
---|
98 | fDist = -1;
|
---|
99 | fAlpha = 0;
|
---|
100 | fCosDeltaAlpha = 0;
|
---|
101 |
|
---|
102 | fDCA = -1;
|
---|
103 | fDCADelta = 0;
|
---|
104 | }
|
---|
105 |
|
---|
106 | // --------------------------------------------------------------------------
|
---|
107 | //
|
---|
108 | // Calculation of source-dependent parameters
|
---|
109 | // In case you don't call Calc from within an eventloop make sure, that
|
---|
110 | // you call the Reset member function before.
|
---|
111 | //
|
---|
112 | Int_t MHillasSrc::Calc(const MHillas &hillas)
|
---|
113 | {
|
---|
114 | const Double_t mx = hillas.GetMeanX(); // [mm]
|
---|
115 | const Double_t my = hillas.GetMeanY(); // [mm]
|
---|
116 |
|
---|
117 | const Double_t sx = mx - fSrcPos->GetX(); // [mm]
|
---|
118 | const Double_t sy = my - fSrcPos->GetY(); // [mm]
|
---|
119 |
|
---|
120 | const Double_t sd = hillas.GetSinDelta(); // [1]
|
---|
121 | const Double_t cd = hillas.GetCosDelta(); // [1]
|
---|
122 |
|
---|
123 | //
|
---|
124 | // Distance from source position to center of ellipse.
|
---|
125 | // If the distance is 0 distance, Alpha is not specified.
|
---|
126 | // The calculation has failed and returnes kFALSE.
|
---|
127 | //
|
---|
128 | const Double_t dist = TMath::Sqrt(sx*sx + sy*sy); // [mm]
|
---|
129 | if (dist==0)
|
---|
130 | return 1;
|
---|
131 |
|
---|
132 | //
|
---|
133 | // Calculate Alpha and Cosda = cos(d,a)
|
---|
134 | // The sign of Cosda will be used for quantities containing
|
---|
135 | // a head-tail information
|
---|
136 | //
|
---|
137 | // *OLD* const Double_t arg = (sy-tand*sx) / (dist*sqrt(tand*tand+1));
|
---|
138 | // *OLD* fAlpha = asin(arg)*kRad2Deg;
|
---|
139 | //
|
---|
140 | const Double_t arg2 = cd*sx + sd*sy; // [mm]
|
---|
141 | if (arg2==0)
|
---|
142 | return 2;
|
---|
143 |
|
---|
144 | const Double_t arg1 = cd*sy - sd*sx; // [mm]
|
---|
145 |
|
---|
146 | //
|
---|
147 | // Due to numerical uncertanties in the calculation of the
|
---|
148 | // square root (dist) and arg1 it can happen (in less than 1e-5 cases)
|
---|
149 | // that the absolute value of arg exceeds 1. Because this uncertainty
|
---|
150 | // results in an Delta Alpha which is still less than 1e-3 we don't care
|
---|
151 | // about this uncertainty in general and simply set values which exceed
|
---|
152 | // to 1 saving its sign.
|
---|
153 | //
|
---|
154 | const Double_t arg = arg1/dist;
|
---|
155 | fAlpha = TMath::Abs(arg)>1 ? TMath::Sign(90., arg) : TMath::ASin(arg)*TMath::RadToDeg(); // [deg]
|
---|
156 |
|
---|
157 | fCosDeltaAlpha = arg2/dist; // [1]
|
---|
158 | fDist = dist; // [mm]
|
---|
159 |
|
---|
160 | // ----- Calculation of Distance to closest approach 'DCA' -----
|
---|
161 |
|
---|
162 | // Components of DCA vector
|
---|
163 | const Double_t fpd1 = sx - arg2*cd; // [mm]
|
---|
164 | const Double_t fpd2 = sy - arg2*sd; // [mm]
|
---|
165 |
|
---|
166 | // Determine the correct sign of the DCA (cross product of DCA vector and the
|
---|
167 | // vector going from the intersection point of the DCA vector with the shower axis
|
---|
168 | // to the COG)
|
---|
169 | const Double_t sign = arg2*cd*fpd2 - arg2*sd*fpd1;
|
---|
170 | fDCA = TMath::Sign(TMath::Sqrt(fpd1*fpd1 + fpd2*fpd2), sign); // [mm]
|
---|
171 |
|
---|
172 | // Calculate angle of the shower axis with respect to the x-axis
|
---|
173 | fDCADelta = TMath::ACos((sx-fpd1)/TMath::Abs(arg2))*TMath::RadToDeg(); // [deg]
|
---|
174 |
|
---|
175 | // Enlarge the interval of fDdca to [0, 2pi]
|
---|
176 | if (sy < fpd2)
|
---|
177 | fDCADelta = 360 - fDCADelta;
|
---|
178 |
|
---|
179 | SetReadyToSave();
|
---|
180 |
|
---|
181 | return 0;
|
---|
182 | }
|
---|
183 |
|
---|
184 | void MHillasSrc::Paint(Option_t *opt)
|
---|
185 | {
|
---|
186 | const Float_t x = fSrcPos ? fSrcPos->GetX() : 0;
|
---|
187 | const Float_t y = fSrcPos ? fSrcPos->GetY() : 0;
|
---|
188 |
|
---|
189 | TMarker m;
|
---|
190 | m.SetMarkerColor(kYellow);
|
---|
191 | m.SetMarkerStyle(kStar);
|
---|
192 | m.PaintMarker(x, y);
|
---|
193 | /*
|
---|
194 | m.SetMarkerColor(kMagenta);
|
---|
195 | m.PaintMarker(fDist*cos((fDCADelta-fAlpha)*TMath::DegToRad()),
|
---|
196 | fDist*sin((fDCADelta-fAlpha)*TMath::DegToRad()));
|
---|
197 |
|
---|
198 | TLine line;
|
---|
199 | line.SetLineWidth(1);
|
---|
200 | line.SetLineColor(108);
|
---|
201 |
|
---|
202 | line.PaintLine(-radius, y, radius, y);
|
---|
203 | line.PaintLine(x, -radius, x, radius);
|
---|
204 |
|
---|
205 | // COG line
|
---|
206 | TLine line(x, y, meanX, meanY);
|
---|
207 | line.SetLineWidth(1);
|
---|
208 | line.SetLineColor(2);
|
---|
209 | line.Paint();*/
|
---|
210 | }
|
---|
211 |
|
---|
212 | // --------------------------------------------------------------------------
|
---|
213 | //
|
---|
214 | // Print contents of MHillasSrc to *fLog
|
---|
215 | //
|
---|
216 | void MHillasSrc::Print(Option_t *) const
|
---|
217 | {
|
---|
218 | *fLog << all;
|
---|
219 | *fLog << "Source dependant Image Parameters (" << GetName() << ")" << endl;
|
---|
220 | *fLog << " - Dist [mm] = " << fDist << endl;
|
---|
221 | *fLog << " - Alpha [deg] = " << fAlpha << endl;
|
---|
222 | *fLog << " - CosDeltaAlpha = " << fCosDeltaAlpha << endl;
|
---|
223 | *fLog << " - DCA [mm] = " << fDCA << endl;
|
---|
224 | *fLog << " - DCA delta [deg] = " << fDCADelta << endl;
|
---|
225 | }
|
---|
226 |
|
---|
227 | // --------------------------------------------------------------------------
|
---|
228 | //
|
---|
229 | // Print contents of MHillasSrc to *fLog depending on the geometry in
|
---|
230 | // units of deg.
|
---|
231 | //
|
---|
232 | void MHillasSrc::Print(const MGeomCam &geom) const
|
---|
233 | {
|
---|
234 | *fLog << all;
|
---|
235 | *fLog << "Source dependant Image Parameters (" << GetName() << ")" << endl;
|
---|
236 | *fLog << " - Dist [deg] = " << fDist*geom.GetConvMm2Deg() << endl;
|
---|
237 | *fLog << " - Alpha [deg] = " << fAlpha << endl;
|
---|
238 | *fLog << " - CosDeltaAlpha = " << fCosDeltaAlpha << endl;
|
---|
239 | *fLog << " - DCA [deg] = " << fDCA*geom.GetConvMm2Deg() << endl;
|
---|
240 | *fLog << " - DCA delta [deg] = " << fDCADelta << endl;
|
---|
241 | }
|
---|
242 |
|
---|
243 | // --------------------------------------------------------------------------
|
---|
244 | //
|
---|
245 | // This function is ment for special usage, please never try to set
|
---|
246 | // values via this function
|
---|
247 | //
|
---|
248 | void MHillasSrc::Set(const TArrayF &arr)
|
---|
249 | {
|
---|
250 | if (arr.GetSize() != 5)
|
---|
251 | return;
|
---|
252 |
|
---|
253 | fAlpha = arr.At(0); // [deg] angle of major axis with vector to src
|
---|
254 | fDist = arr.At(1); // [mm] distance between src and center of ellipse
|
---|
255 | fCosDeltaAlpha = arr.At(2); // [1] cosine of angle between d and a
|
---|
256 | fDCA = arr.At(3); // [mm]
|
---|
257 | fDCADelta = arr.At(4); // [mm]
|
---|
258 | }
|
---|