LSSTApplications  12.1-5-gbdcc3ab+2,15.0+13,15.0+26,15.0-1-g19261fa+17,15.0-1-g60afb23+26,15.0-1-g615e0bb+18,15.0-1-g788a293+26,15.0-1-ga91101e+26,15.0-1-gae1598d+12,15.0-1-gd076f1f+24,15.0-1-gdf18595+5,15.0-1-gf4f1c34+12,15.0-11-g7db6e543+4,15.0-12-g3681e7a+4,15.0-15-gc15de322,15.0-16-g83b84f4,15.0-2-g100d730+19,15.0-2-g1f9c9cf+4,15.0-2-g8aea5f4+1,15.0-2-gf38729e+21,15.0-29-ga12a2b06e,15.0-3-g11fe1a0+14,15.0-3-g707930d+3,15.0-3-g9103c06+12,15.0-3-gd3cbb57+3,15.0-4-g2d82b59,15.0-4-g535e784+10,15.0-4-g92ca6c3+4,15.0-4-gf906033+2,15.0-5-g23e394c+14,15.0-5-g4be42a9,15.0-6-g69628aa,15.0-6-g86e3f3d+1,15.0-6-gfa9b38f+4,15.0-7-g949993c+3,15.0-8-g67a62d3+1,15.0-8-gcf05001+1,15.0-9-g1e7c341+1,w.2018.21
LSSTDataManagementBasePackage
Public Member Functions | Related Functions | List of all members
lsst::afw::geom::SpherePoint Class Referencefinal

Point in an unspecified spherical coordinate system. More...

#include <SpherePoint.h>

Public Member Functions

 SpherePoint (Angle const &longitude, Angle const &latitude)
 Construct a SpherePoint from a longitude and latitude. More...
 
 SpherePoint (double longitude, double latitude, AngleUnit units)
 Construct a SpherePoint from double longitude and latitude. More...
 
 SpherePoint (sphgeom::Vector3d const &vector)
 Construct a SpherePoint from a vector representing a direction. More...
 
 SpherePoint (sphgeom::LonLat const &lonLat)
 Construct a SpherePoint from a sphgeom::LonLat. More...
 
 SpherePoint ()
 Construct a SpherePoint with "nan" for longitude and latitude. More...
 
 SpherePoint (SpherePoint const &other) noexcept
 Create a copy of a SpherePoint. More...
 
 SpherePoint (SpherePoint &&other) noexcept
 Create a copy of a SpherePoint. More...
 
 ~SpherePoint ()
 
SpherePointoperator= (SpherePoint const &other) noexcept
 Overwrite this object with the value of another SpherePoint. More...
 
SpherePointoperator= (SpherePoint &&other) noexcept
 Overwrite this object with the value of another SpherePoint. More...
 
 operator sphgeom::LonLat () const
 Make SpherePoint implicitly convertible to LonLat. More...
 
Angle getLongitude () const noexcept
 The longitude of this point. More...
 
Angle getRa () const noexcept
 Synonym for getLongitude. More...
 
Angle getLatitude () const noexcept
 The latitude of this point. More...
 
Angle getDec () const noexcept
 Synonym for getLatitude. More...
 
Point2D getPosition (AngleUnit unit) const
 Return longitude, latitude as a Point2D object. More...
 
sphgeom::UnitVector3d getVector () const noexcept
 A unit vector representation of this point. More...
 
Angle operator[] (size_t index) const
 Longitude and latitude by index. More...
 
bool atPole () const noexcept
 true if this point is either coordinate pole. More...
 
bool isFinite () const noexcept
 true if this point is a well-defined position. More...
 
bool operator== (SpherePoint const &other) const noexcept
 true if two points have the same longitude and latitude. More...
 
bool operator!= (SpherePoint const &other) const noexcept
 false if two points represent the same position. More...
 
Angle bearingTo (SpherePoint const &other) const
 Orientation at this point of the great circle arc to another point. More...
 
Angle separation (SpherePoint const &other) const noexcept
 Angular distance between two points. More...
 
SpherePoint rotated (SpherePoint const &axis, Angle const &amount) const noexcept
 Return a point rotated from this one around an axis. More...
 
SpherePoint offset (Angle const &bearing, Angle const &amount) const
 Return a point offset from this one along a great circle. More...
 
std::pair< Angle, AnglegetTangentPlaneOffset (SpherePoint const &other) const
 Get the offset from a tangent plane centered at this point to another point. More...
 

Related Functions

(Note that these are not member functions.)

std::ostreamoperator<< (std::ostream &os, SpherePoint const &point)
 Print the value of a point to a stream. More...
 

Detailed Description

Point in an unspecified spherical coordinate system.

This class represents a point on a sphere in the mathematical rather than the astronomical sense. It does not refer to any astronomical reference frame, nor does it have any concept of (radial) distance.

Points can be represented either as spherical coordinates or as a unit vector. The adopted convention for converting between these two systems is that (0, 0) corresponds to <1, 0, 0>, that (π/2, 0) corresponds to <0, 1, 0>, and that (0, π/2) corresponds to <0, 0, 1>.

To keep usage simple, SpherePoint does not support modification of existing values; transformations of SpherePoints should be expressed as a new object. To support cases where a SpherePoint must be updated in-place, SpherePoint supports assignment to an existing object. One example is in containers of SpherePoints; no STL container has an atomic element-replacement method, so complicated constructions would need to be used if you couldn't overwrite an existing element.

Definition at line 61 of file SpherePoint.h.

Constructor & Destructor Documentation

◆ SpherePoint() [1/7]

lsst::afw::geom::SpherePoint::SpherePoint ( Angle const &  longitude,
Angle const &  latitude 
)

Construct a SpherePoint from a longitude and latitude.

Parameters
longitudeThe longitude of the point. +/-inf is treated as nan (because longitude is wrapped to a standard range and infinity cannot be wrapped)
latitudeThe latitude of the point. Must be in the interval [-π/2, π/2] radians.
Exceptions
pex::exceptions::InvalidParameterErrorThrown if latitude isout of range.
Exception Safety
Provides strong exception guarantee.

Definition at line 70 of file SpherePoint.cc.

71  : _longitude(longitude.wrap().asRadians()), _latitude(latitude.asRadians()) {
72  if (fabs(_latitude) > HALFPI) {
73  throw pexExcept::InvalidParameterError("Angle " + to_string(latitude.asDegrees()) +
74  " is not a valid latitude.");
75  }
76 }
table::Key< geom::Angle > longitude
Definition: VisitInfo.cc:166
double constexpr HALFPI
Definition: Angle.h:23
T to_string(T... args)
table::Key< geom::Angle > latitude
Definition: VisitInfo.cc:165
T fabs(T... args)
Reports invalid arguments.
Definition: Runtime.h:66

◆ SpherePoint() [2/7]

lsst::afw::geom::SpherePoint::SpherePoint ( double  longitude,
double  latitude,
AngleUnit  units 
)

Construct a SpherePoint from double longitude and latitude.

Parameters
longitudeThe longitude of the point. +/-inf is treated as nan (because longitude is wrapped to a standard range and infinity cannot be wrapped)
latitudeThe latitude of the point. Must be in the interval [-π/2, π/2] radians.
unitsThe units of longitude and latitude
Exceptions
pex::exceptions::InvalidParameterErrorThrown if latitude isout of range.
Exception Safety
Provides strong exception guarantee.

Definition at line 67 of file SpherePoint.cc.

68  : SpherePoint(longitude * units, latitude * units) {}
table::Key< geom::Angle > longitude
Definition: VisitInfo.cc:166
table::Key< geom::Angle > latitude
Definition: VisitInfo.cc:165
SpherePoint()
Construct a SpherePoint with "nan" for longitude and latitude.
Definition: SpherePoint.cc:110

◆ SpherePoint() [3/7]

lsst::afw::geom::SpherePoint::SpherePoint ( sphgeom::Vector3d const &  vector)
explicit

Construct a SpherePoint from a vector representing a direction.

Parameters
vectorA position whose direction will be stored as a SpherePoint. Must not be the zero vector. Need not be normalized, and the norm will not affect the value of the point.
Exceptions
pex::exceptions::InvalidParameterErrorThrown if vector is the zero vector.
Exception Safety
Provides strong exception guarantee.
Note
If the SpherePoint is at a pole then longitude is set to 0. That provides predictable behavior for bearingTo and offset.

Definition at line 78 of file SpherePoint.cc.

78  {
79  // sphgeom Vector3d has its own normalization,
80  // but its behavior is not documented for non-finite values
81  double norm = vector.getNorm();
82  if (norm <= 0.0) {
83  stringstream buffer;
84  buffer << "Vector " << vector << " has zero norm and cannot be normalized.";
85  throw pexExcept::InvalidParameterError(buffer.str());
86  }
87  // To avoid unexpected behavior from mixing finite elements and infinite norm
88  if (!isfinite(norm)) {
89  norm = NAN;
90  }
91  auto unitVector =
92  sphgeom::UnitVector3d::fromNormalized(vector.x() / norm, vector.y() / norm, vector.z() / norm);
93  _set(unitVector);
94 }
T norm(const T &x)
Definition: Integrate.h:194
T str(T... args)
T isfinite(T... args)
STL class.
Reports invalid arguments.
Definition: Runtime.h:66
static UnitVector3d fromNormalized(Vector3d const &v)
fromNormalized returns the unit vector equal to v, which is assumed to be normalized.
Definition: UnitVector3d.h:82

◆ SpherePoint() [4/7]

lsst::afw::geom::SpherePoint::SpherePoint ( sphgeom::LonLat const &  lonLat)

Construct a SpherePoint from a sphgeom::LonLat.

Parameters
lonLatThe lonLat
Exception Safety
Provides strong exception guarantee.

Definition at line 96 of file SpherePoint.cc.

97  : SpherePoint(lonLat.getLon().asRadians(), lonLat.getLat().asRadians(), radians) {}
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:87
SpherePoint()
Construct a SpherePoint with "nan" for longitude and latitude.
Definition: SpherePoint.cc:110

◆ SpherePoint() [5/7]

lsst::afw::geom::SpherePoint::SpherePoint ( )

Construct a SpherePoint with "nan" for longitude and latitude.

Definition at line 110 of file SpherePoint.cc.

110 : _longitude(nan("")), _latitude(nan("")) {}
T nan(T... args)

◆ SpherePoint() [6/7]

lsst::afw::geom::SpherePoint::SpherePoint ( SpherePoint const &  other)
defaultnoexcept

Create a copy of a SpherePoint.

Parameters
otherThe point to copy.
Exception Safety
Shall not throw exceptions.

◆ SpherePoint() [7/7]

lsst::afw::geom::SpherePoint::SpherePoint ( SpherePoint &&  other)
defaultnoexcept

Create a copy of a SpherePoint.

As SpherePoint(SpherePoint const&), except that other may be left in an unspecified but valid state.

◆ ~SpherePoint()

lsst::afw::geom::SpherePoint::~SpherePoint ( )
default

Member Function Documentation

◆ atPole()

bool lsst::afw::geom::SpherePoint::atPole ( ) const
inlinenoexcept

true if this point is either coordinate pole.

Returns
true if this point is at the north or south pole, false otherwise
Exception Safety
Shall not throw exceptions.

Definition at line 237 of file SpherePoint.h.

237  {
238  // Unit tests indicate we don't need to worry about epsilon-errors, in that
239  // Objects constructed from lat=90*degrees, <0,0,1>, etc. all have atPole() = true.
240  return fabs(_latitude) >= HALFPI;
241  }
double constexpr HALFPI
Definition: Angle.h:23
T fabs(T... args)

◆ bearingTo()

Angle lsst::afw::geom::SpherePoint::bearingTo ( SpherePoint const &  other) const

Orientation at this point of the great circle arc to another point.

The orientation is measured in a plane tangent to the sphere at this point, with 0 degrees along the direction of increasing longitude and 90 degrees along the direction of increasing latitude. Thus for any point except the pole, the arc is due east at 0 degrees and due north at 90 degrees. To understand the behavior at the poles it may help to imagine the point has been shifted slightly by decreasing the absolute value of its latitude.

If the points are on opposite sides of the sphere then the returned bearing may be any value.

Parameters
otherthe point to which to measure the bearing
Returns
the direction, as defined above, in the interval [0, 2π).
Exception Safety
Provides strong exception guarantee.
Note
For two points A and B, A.bearingTo(B) will in general not be 180 degrees away from B.bearingTo(A).

Definition at line 156 of file SpherePoint.cc.

156  {
157  Angle const deltaLon = other.getLongitude() - this->getLongitude();
158 
159  double const sinDelta1 = sin(getLatitude().asRadians());
160  double const sinDelta2 = sin(other.getLatitude().asRadians());
161  double const cosDelta1 = cos(getLatitude().asRadians());
162  double const cosDelta2 = cos(other.getLatitude().asRadians());
163 
164  // Adapted from http://www.movable-type.co.uk/scripts/latlong.html
165  double const y = sin(deltaLon) * cosDelta2;
166  double const x = cosDelta1 * sinDelta2 - sinDelta1 * cosDelta2 * cos(deltaLon);
167  return (90.0 * degrees - atan2(y, x) * radians).wrap();
168 }
lsst::afw::geom::Angle Angle
Definition: misc.h:33
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:87
Angle getLongitude() const noexcept
The longitude of this point.
Definition: SpherePoint.h:180
AngleUnit constexpr degrees
constant with units of degrees
Definition: Angle.h:88
int y
Definition: SpanSet.cc:44
double sin(Angle const &a)
Definition: Angle.h:102
double cos(Angle const &a)
Definition: Angle.h:103
T atan2(T... args)
double x
Angle getLatitude() const noexcept
The latitude of this point.
Definition: SpherePoint.h:192
ItemVariant const * other
Definition: Schema.cc:55

◆ getDec()

Angle lsst::afw::geom::SpherePoint::getDec ( ) const
inlinenoexcept

Synonym for getLatitude.

Definition at line 195 of file SpherePoint.h.

195 { return _latitude * radians; };
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:87

◆ getLatitude()

Angle lsst::afw::geom::SpherePoint::getLatitude ( ) const
inlinenoexcept

The latitude of this point.

Returns
the latitude, in the interval [-π/2, π/2] radians.
Exception Safety
Shall not throw exceptions.

Definition at line 192 of file SpherePoint.h.

192 { return _latitude * radians; };
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:87

◆ getLongitude()

Angle lsst::afw::geom::SpherePoint::getLongitude ( ) const
inlinenoexcept

The longitude of this point.

Returns
the longitude, in the interval [0, 2π) radians.
Exception Safety
Shall not throw exceptions.

Definition at line 180 of file SpherePoint.h.

180 { return _longitude * radians; };
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:87

◆ getPosition()

Point2D lsst::afw::geom::SpherePoint::getPosition ( AngleUnit  unit) const

Return longitude, latitude as a Point2D object.

Parameters
[in]unitUnits of returned data.

Definition at line 131 of file SpherePoint.cc.

131  {
132  return Point2D(getLongitude().asAngularUnits(unit), getLatitude().asAngularUnits(unit));
133 }
Angle getLongitude() const noexcept
The longitude of this point.
Definition: SpherePoint.h:180
Point< double, 2 > Point2D
Definition: Point.h:304
Angle getLatitude() const noexcept
The latitude of this point.
Definition: SpherePoint.h:192

◆ getRa()

Angle lsst::afw::geom::SpherePoint::getRa ( ) const
inlinenoexcept

Synonym for getLongitude.

Definition at line 183 of file SpherePoint.h.

183 { return _longitude * radians; };
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:87

◆ getTangentPlaneOffset()

std::pair< geom::Angle, geom::Angle > lsst::afw::geom::SpherePoint::getTangentPlaneOffset ( SpherePoint const &  other) const

Get the offset from a tangent plane centered at this point to another point.

This is suitable only for small angles.

Parameters
[in]otherCoordinate to which to compute offset
Returns
pair of Angles: Longitude and Latitude offsets

Definition at line 211 of file SpherePoint.cc.

211  {
212  geom::Angle const alpha1 = this->getLongitude();
213  geom::Angle const delta1 = this->getLatitude();
214  geom::Angle const alpha2 = other.getLongitude();
215  geom::Angle const delta2 = other.getLatitude();
216 
217  // Compute the projection of "other" on a tangent plane centered at this point
218  double const sinDelta1 = std::sin(delta1);
219  double const cosDelta1 = std::cos(delta1);
220  double const sinDelta2 = std::sin(delta2);
221  double const cosDelta2 = std::cos(delta2);
222  double const cosAlphaDiff = std::cos(alpha2 - alpha1);
223  double const sinAlphaDiff = std::sin(alpha2 - alpha1);
224 
225  double const div = cosDelta1 * cosAlphaDiff * cosDelta2 + sinDelta1 * sinDelta2;
226  double const xi = cosDelta1 * sinAlphaDiff / div;
227  double const eta = (cosDelta1 * cosAlphaDiff * sinDelta2 - sinDelta1 * cosDelta2) / div;
228 
229  return std::make_pair(xi * geom::radians, eta * geom::radians);
230 }
AngleUnit constexpr radians
constant with units of radians
Definition: Angle.h:87
Angle getLongitude() const noexcept
The longitude of this point.
Definition: SpherePoint.h:180
T sin(T... args)
A class representing an angle.
Definition: Angle.h:102
T make_pair(T... args)
T cos(T... args)
T div(T... args)
Angle getLatitude() const noexcept
The latitude of this point.
Definition: SpherePoint.h:192
ItemVariant const * other
Definition: Schema.cc:55

◆ getVector()

sphgeom::UnitVector3d lsst::afw::geom::SpherePoint::getVector ( ) const
noexcept

A unit vector representation of this point.

Returns
a unit vector whose direction corresponds to this point
Exception Safety
Shall not throw exceptions.

Definition at line 126 of file SpherePoint.cc.

126  {
127  return sphgeom::UnitVector3d::fromNormalized(cos(_longitude) * cos(_latitude),
128  sin(_longitude) * cos(_latitude), sin(_latitude));
129 }
double sin(Angle const &a)
Definition: Angle.h:102
double cos(Angle const &a)
Definition: Angle.h:103
static UnitVector3d fromNormalized(Vector3d const &v)
fromNormalized returns the unit vector equal to v, which is assumed to be normalized.
Definition: UnitVector3d.h:82

◆ isFinite()

bool lsst::afw::geom::SpherePoint::isFinite ( ) const
noexcept

true if this point is a well-defined position.

Returns
true if getLongitude(), getLatitude(), and getVector() return finite floating-point values; false if any return NaN or infinity.
Exception Safety
Shall not throw exceptions.

Definition at line 146 of file SpherePoint.cc.

146 { return isfinite(_longitude) && isfinite(_latitude); }
T isfinite(T... args)

◆ offset()

SpherePoint lsst::afw::geom::SpherePoint::offset ( Angle const &  bearing,
Angle const &  amount 
) const

Return a point offset from this one along a great circle.

For any point A, any angle bearing and any non-negative angle amount, if B = A.offset(bearing, amount)thenA.bearingTo(B) = amountandA.separationTo(B) = amount. Negative values ofamountare supported in the obvious manner: A.offset(bearing, delta) = A.offset(bearing + 180*degrees, -delta)`

Parameters
bearingthe direction in which to move this point, following the conventions described in bearingTo.
amountthe distance by which to move along the great circle defined by bearing
Returns
a new point created by rotating this point. If the new point is at the pole then its longitude will be 0.
Exception Safety
Provides strong exception guarantee.

Definition at line 182 of file SpherePoint.cc.

182  {
183  double const phi = bearing.asRadians();
184 
185  // let v = vector in the direction bearing points (tangent to surface of sphere)
186  // To do the rotation, use rotated() method.
187  // - must provide an axis of rotation: take the cross product r x v to get that axis (pole)
188 
189  auto r = getVector();
190 
191  // Get the vector v:
192  // let u = unit vector lying on a parallel of declination
193  // let w = unit vector along line of longitude = r x u
194  // the vector v must satisfy the following:
195  // r . v = 0
196  // u . v = cos(phi)
197  // w . v = sin(phi)
198 
199  // v is a linear combination of u and w
200  // v = cos(phi)*u + sin(phi)*w
201  sphgeom::Vector3d const u(-sin(_longitude), cos(_longitude), 0.0);
202  auto w = r.cross(u);
203  auto v = cos(phi) * u + sin(phi) * w;
204 
205  // take r x v to get the axis
206  SpherePoint axis = SpherePoint(r.cross(v));
207 
208  return rotated(axis, amount);
209 }
double sin(Angle const &a)
Definition: Angle.h:102
double cos(Angle const &a)
Definition: Angle.h:103
SpherePoint()
Construct a SpherePoint with "nan" for longitude and latitude.
Definition: SpherePoint.cc:110
SpherePoint rotated(SpherePoint const &axis, Angle const &amount) const noexcept
Return a point rotated from this one around an axis.
Definition: SpherePoint.cc:175
double w
Definition: CoaddPsf.cc:57
sphgeom::UnitVector3d getVector() const noexcept
A unit vector representation of this point.
Definition: SpherePoint.cc:126

◆ operator sphgeom::LonLat()

lsst::afw::geom::SpherePoint::operator sphgeom::LonLat ( ) const

Make SpherePoint implicitly convertible to LonLat.

Definition at line 120 of file SpherePoint.cc.

120  {
121  return sphgeom::LonLat::fromRadians(getLongitude().asRadians(), getLatitude().asRadians());
122 }
Angle getLongitude() const noexcept
The longitude of this point.
Definition: SpherePoint.h:180
Angle getLatitude() const noexcept
The latitude of this point.
Definition: SpherePoint.h:192
static LonLat fromRadians(double lon, double lat)
Definition: LonLat.h:55

◆ operator!=()

bool lsst::afw::geom::SpherePoint::operator!= ( SpherePoint const &  other) const
noexcept

false if two points represent the same position.

This operator shall always return the logical negation of operator==(); see its documentation for a detailed specification.

Definition at line 154 of file SpherePoint.cc.

154 { return !(*this == other); }
ItemVariant const * other
Definition: Schema.cc:55

◆ operator=() [1/2]

SpherePoint & lsst::afw::geom::SpherePoint::operator= ( SpherePoint const &  other)
defaultnoexcept

Overwrite this object with the value of another SpherePoint.

This is the only method that can alter the state of a SpherePoint after its construction, and is provided to allow in-place replacement of SpherePoints where swapping is not possible.

Parameters
otherThe object with which to overwrite this one.
Returns
a reference to this object.
Exception Safety
Shall not throw exceptions.

◆ operator=() [2/2]

SpherePoint & lsst::afw::geom::SpherePoint::operator= ( SpherePoint &&  other)
defaultnoexcept

Overwrite this object with the value of another SpherePoint.

As operator=(SpherePoint const&), except that other may be left in an unspecified but valid state.

◆ operator==()

bool lsst::afw::geom::SpherePoint::operator== ( SpherePoint const &  other) const
noexcept

true if two points have the same longitude and latitude.

Parameters
otherthe point to test for equality
Returns
true if this point has exactly the same values of getLongitude() and getLatitude() as other, false otherwise
Exception Safety
Shall not throw exceptions.
Note
Two points at the same pole will compare unequal if they have different longitudes, despite representing the same point on the unit sphere. This is important because the behavior of offset and bearingTo depend on longitude, even at the pole.
Warning
Points whose getLongitude(), getLatitude(), or getVector() methods return NaN shall not compare equal to any point, including themselves. This may break algorithms that assume object equality is reflexive; use isFinite() to filter objects if necessary.

Definition at line 148 of file SpherePoint.cc.

148  {
149  // Deliberate override of Style Guide 5-12
150  // Approximate FP comparison would make object equality intransitive
151  return _longitude == other._longitude && _latitude == other._latitude;
152 }
ItemVariant const * other
Definition: Schema.cc:55

◆ operator[]()

Angle lsst::afw::geom::SpherePoint::operator[] ( size_t  index) const

Longitude and latitude by index.

Parameters
indexthe index of the spherical coordinate to return. Must be either 0 or 1.
Returns
getLongitude() if index = 0, or getLatitude() if index = 1
Exceptions
pex::exceptions::OutOfRangeErrorThrown if index is neither 0 nor 1.
Exception Safety
Provides strong exception guarantee.

Definition at line 135 of file SpherePoint.cc.

135  {
136  switch (index) {
137  case 0:
138  return getLongitude();
139  case 1:
140  return getLatitude();
141  default:
142  throw pexExcept::OutOfRangeError("Index " + to_string(index) + " must be 0 or 1.");
143  }
144 }
Angle getLongitude() const noexcept
The longitude of this point.
Definition: SpherePoint.h:180
T to_string(T... args)
Angle getLatitude() const noexcept
The latitude of this point.
Definition: SpherePoint.h:192
Reports attempts to access elements outside a valid range of indices.
Definition: Runtime.h:89

◆ rotated()

SpherePoint lsst::afw::geom::SpherePoint::rotated ( SpherePoint const &  axis,
Angle const &  amount 
) const
noexcept

Return a point rotated from this one around an axis.

Parameters
axisa point defining the north pole of the rotation axis.
amountthe amount of rotation, where positive values represent right-handed rotations around axis.
Returns
a new point created by rotating this point. If the new point is at the pole then its longitude will be 0.
Exception Safety
Shall not throw exceptions.

Definition at line 175 of file SpherePoint.cc.

175  {
176  auto const rotation = Eigen::AngleAxisd(amount.asRadians(), asEigen(axis.getVector())).matrix();
177  auto const x = asEigen(getVector());
178  auto const xprime = rotation * x;
179  return SpherePoint(sphgeom::Vector3d(xprime[0], xprime[1], xprime[2]));
180 }
Eigen::Vector3d asEigen(sphgeom::Vector3d const &vector)
Definition: sphgeomUtils.h:38
double x
SpherePoint()
Construct a SpherePoint with "nan" for longitude and latitude.
Definition: SpherePoint.cc:110
sphgeom::UnitVector3d getVector() const noexcept
A unit vector representation of this point.
Definition: SpherePoint.cc:126

◆ separation()

Angle lsst::afw::geom::SpherePoint::separation ( SpherePoint const &  other) const
noexcept

Angular distance between two points.

Parameters
otherthe point to which to measure the separation
Returns
the length of the shortest (great circle) arc between the two points. Shall not be negative.
Exception Safety
Shall not throw exceptions.

Definition at line 170 of file SpherePoint.cc.

170  {
171  return haversine(getLongitude() - other.getLongitude(), getLatitude() - other.getLatitude(),
172  cos(getLatitude().asRadians()), cos(other.getLatitude().asRadians()));
173 }
Angle getLongitude() const noexcept
The longitude of this point.
Definition: SpherePoint.h:180
double cos(Angle const &a)
Definition: Angle.h:103
Angle getLatitude() const noexcept
The latitude of this point.
Definition: SpherePoint.h:192
ItemVariant const * other
Definition: Schema.cc:55

Friends And Related Function Documentation

◆ operator<<()

std::ostream & operator<< ( std::ostream os,
SpherePoint const &  point 
)
related

Print the value of a point to a stream.

The exact details of the string representation are unspecified and subject to change, but the following may be regarded as typical: "(10.543250, +32.830583)".

Parameters
osthe stream to which to print point
pointthe point to print to the stream
Returns
a reference to os
Exceptions
pex::exceptions::std::ostream::failureThrown if an I/O state flag was set that was registered with os.exceptions(). See the documentation of std::ostream for more details.
Exception Safety
Provides basic exception guarantee.

Definition at line 250 of file SpherePoint.cc.

250  {
251  // Can't provide atomic guarantee anyway for I/O, so ok to be sloppy.
252  auto oldFlags = os.setf(ostream::fixed);
253  auto oldPrecision = os.precision(6);
254 
255  os << "(" << point.getLongitude().asDegrees() << ", ";
256  os.setf(ostream::showpos);
257  os << point.getLatitude().asDegrees() << ")";
258 
259  os.flags(oldFlags);
260  os.precision(oldPrecision);
261  return os;
262 }
T setf(T... args)
T precision(T... args)
T flags(T... args)
std::ostream * os
Definition: Schema.cc:736

The documentation for this class was generated from the following files: