source: trunk/MagicSoft/Mars/mimage/MHillasSrc.cc@ 5032

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