LSST Applications 24.1.6,g063fba187b+56b85ce14a,g0f08755f38+df8a265115,g12f32b3c4e+891a09f10d,g1524ad2192+7a5d7b3fbd,g1653933729+a8ce1bb630,g168dd56ebc+a8ce1bb630,g28da252d5a+07cb1400be,g2bbee38e9b+ae03bbfc84,g2bc492864f+ae03bbfc84,g3156d2b45e+6e55a43351,g347aa1857d+ae03bbfc84,g35bb328faa+a8ce1bb630,g3a166c0a6a+ae03bbfc84,g3e281a1b8c+c5dd892a6c,g414038480c+6b9177ef31,g41af890bb2+8f257c4c0b,g781aacb6e4+a8ce1bb630,g7af13505b9+7137b3b17d,g80478fca09+6df6903293,g82479be7b0+091ce1d07f,g858d7b2824+df8a265115,g89c8672015+f4add4ffd5,g9125e01d80+a8ce1bb630,g9726552aa6+414189b318,ga5288a1d22+4a2bca08d7,gacef1a1666+c9a8ff65f4,gb58c049af0+d64f4d3760,gbcfae0f0a0+de1d42d831,gc28159a63d+ae03bbfc84,gcf0d15dbbd+72117bf34e,gda6a2b7d83+72117bf34e,gdaeeff99f8+1711a396fd,ge500cccec5+c8c9c9af63,ge79ae78c31+ae03bbfc84,gf0baf85859+c1f95f4921,gfa517265be+df8a265115,gfa999e8aa5+17cd334064,gfb92a5be7c+df8a265115
LSST Data Management Base Package
Loading...
Searching...
No Matches
Frame.h
Go to the documentation of this file.
1/*
2 * LSST Data Management System
3 * Copyright 2017 AURA/LSST.
4 *
5 * This product includes software developed by the
6 * LSST Project (http://www.lsst.org/).
7 *
8 * This program is free software: you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License as published by
10 * the Free Software Foundation, either version 3 of the License, or
11 * (at your option) any later version.
12 *
13 * This program is distributed in the hope that it will be useful,
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 * GNU General Public License for more details.
17 *
18 * You should have received a copy of the LSST License Statement and
19 * the GNU General Public License along with this program. If not,
20 * see <https://www.lsstcorp.org/LegalNotices/>.
21 */
22#ifndef ASTSHIM_FRAME_H
23#define ASTSHIM_FRAME_H
24
25#include <memory>
26#include <sstream>
27#include <stdexcept>
28#include <utility>
29
30#include "astshim/Mapping.h"
31#include "astshim/Object.h"
32#include "astshim/base.h"
34
35namespace ast {
36
37class CmpFrame;
38class Frame;
39
55
60public:
68 int nread;
69 double value;
70};
71
76public:
82 explicit ResolvedPoint(int naxes) : point(naxes), d1(), d2() {}
84 double d1;
85 double d2;
86};
87
104
105class FrameSet;
106
157class Frame : public Mapping {
158 friend class Object;
159
160public:
168 explicit Frame(int naxes, std::string const &options = "")
169 : Mapping(reinterpret_cast<AstMapping *>(astFrame(naxes, "%s", options.c_str()))) {
170 assertOK();
171 }
172
173 virtual ~Frame() {}
174
176 Frame(Frame const &) = default;
177 Frame(Frame &&) = default;
178 Frame &operator=(Frame const &) = delete;
179 Frame &operator=(Frame &&) = default;
180
182 std::shared_ptr<Frame> copy() const { return std::static_pointer_cast<Frame>(copyPolymorphic()); }
183
201 double angle(PointD const &a, PointD const &b, PointD const &c) const {
202 assertPointLength(a, "a");
203 assertPointLength(b, "b");
204 assertPointLength(c, "c");
205 return detail::safeDouble(astAngle(getRawPtr(), a.data(), b.data(), c.data()));
206 }
207
229 double axAngle(PointD const &a, PointD const &b, int axis) const {
230 assertPointLength(a, "a");
231 assertPointLength(b, "b");
232 return detail::safeDouble(astAxAngle(getRawPtr(), a.data(), b.data(), axis));
233 }
234
247 double axDistance(int axis, double v1, double v2) const {
248 return detail::safeDouble(astAxDistance(getRawPtr(), axis, v1, v2));
249 }
250
263 double axOffset(int axis, double v1, double dist) const {
264 return detail::safeDouble(astAxOffset(getRawPtr(), axis, v1, dist));
265 }
266
459 std::shared_ptr<FrameSet> convert(Frame const &to, std::string const &domainlist = "");
460
474 double distance(PointD const &point1, PointD const &point2) const {
475 assertPointLength(point1, "point1");
476 assertPointLength(point2, "point2");
477 return detail::safeDouble(astDistance(getRawPtr(), point1.data(), point2.data()));
478 }
479
769 std::shared_ptr<FrameSet> findFrame(Frame const &tmplt, std::string const &domainlist = "");
770
782 std::string format(int axis, double value) const {
783 char const *rawstr = astFormat(getRawPtr(), axis, value);
784 assertOK();
785 return std::string(rawstr);
786 }
787
792 bool getActiveUnit() const {
793 bool ret = astGetActiveUnit(getRawPtr());
794 assertOK();
795 return ret;
796 }
797
802 std::string getAlignSystem() const { return getC("AlignSystem"); }
803
807 double getBottom(int axis) const { return getD(detail::formatAxisAttr("Bottom", axis)); }
808
812 int getDigits() const { return getI("Digits"); }
813
817 int getDigits(int axis) const { return getI(detail::formatAxisAttr("Digits", axis)); }
818
822 bool getDirection(int axis) const { return getB(detail::formatAxisAttr("Direction", axis)); }
823
827 std::string getDomain() const { return getC("Domain"); }
828
832 double getDut1() const { return getD("Dut1"); }
833
837 double getEpoch() const { return getD("Epoch"); }
838
842 std::string getFormat(int axis) const { return getC(detail::formatAxisAttr("Format", axis)); }
843
848 std::string getInternalUnit(int axis) const { return getC(detail::formatAxisAttr("InternalUnit", axis)); }
849
853 std::string getLabel(int axis) const { return getC(detail::formatAxisAttr("Label", axis)); }
854
858 bool getMatchEnd() const { return getB("MatchEnd"); }
859
864 int getMaxAxes() const { return getI("MaxAxes"); }
865
870 int getMinAxes() const { return getI("MinAxes"); }
871
876 int getNAxes() const { return getI("NAxes"); }
877
882 std::string getNormUnit(int axis) const { return getC(detail::formatAxisAttr("NormUnit", axis)); }
883
887 double getObsAlt() const { return getD("ObsAlt"); }
888
892 std::string getObsLat() const { return getC("ObsLat"); }
893
897 std::string getObsLon() const { return getC("ObsLon"); }
898
902 bool getPermute() const { return getB("Permute"); }
903
907 bool getPreserveAxes() const { return getB("PreserveAxes"); }
908
912 std::string getSymbol(int axis) const { return getC(detail::formatAxisAttr("Symbol", axis)); }
913
918 std::string getSystem() const { return getC("System"); }
919
923 std::string getTitle() const { return getC("Title"); }
924
928 double getTop(int axis) const { return getD(detail::formatAxisAttr("Top", axis)); }
929
933 std::string getUnit(int axis) const { return getC(detail::formatAxisAttr("Unit", axis)); }
934
967 std::vector<double> const &b1, std::vector<double> const &b2) const;
968
985 std::vector<int> matchAxes(Frame const &other) const {
986 std::vector<int> ret(other.getNIn());
987 astMatchAxes(getRawPtr(), other.getRawPtr(), ret.data());
988 assertOK();
989 return ret;
990 }
991
1023 CmpFrame under(Frame const &next) const;
1024
1045 PointD norm(PointD value) const {
1046 astNorm(getRawPtr(), value.data());
1047 assertOK();
1048 detail::astBadToNan(value);
1049 return value;
1050 }
1051
1075 PointD offset(PointD point1, PointD point2, double offset) const {
1076 assertPointLength(point1, "point1");
1077 assertPointLength(point2, "point2");
1078 PointD ret(getNIn());
1079 astOffset(getRawPtr(), point1.data(), point2.data(), offset, ret.data());
1080 assertOK();
1082 return ret;
1083 }
1084
1117 DirectionPoint offset2(PointD const &point1, double angle, double offset) const {
1118 detail::assertEqual(getNIn(), "naxes", 2, " cannot call offset2");
1119 assertPointLength(point1, "point1");
1120 PointD point2(getNIn());
1121 double offsetAngle = astOffset2(getRawPtr(), point1.data(), angle, offset, point2.data());
1122 assertOK();
1123 detail::astBadToNan(point2);
1124 return DirectionPoint(detail::safeDouble(offsetAngle), point2);
1125 }
1126
1139 detail::assertEqual(perm.size(), "perm.size()", static_cast<std::size_t>(getNAxes()), "naxes");
1140 astPermAxes(getRawPtr(), perm.data());
1141 assertOK();
1142 }
1143
1163 FrameMapping pickAxes(std::vector<int> const &axes) const;
1164
1196 ResolvedPoint resolve(std::vector<double> const &point1, std::vector<double> const &point2,
1197 std::vector<double> const &point3) const;
1198
1203 void setAlignSystem(std::string const &system) { setC("AlignSystem", system); }
1204
1208 void setBottom(int axis, double bottom) { setD(detail::formatAxisAttr("Bottom", axis), bottom); }
1209
1213 void setDigits(int digits) { setI("Digits", digits); }
1214
1218 void setDigits(int axis, int digits) { setD(detail::formatAxisAttr("Digits", axis), digits); }
1219
1223 void setDirection(bool direction, int axis) {
1224 setB(detail::formatAxisAttr("Direction", axis), direction);
1225 }
1226
1230 virtual void setDomain(std::string const &domain) { setC("Domain", domain); }
1231
1235 void setDut1(double dut1) { setD("Dut1", dut1); }
1236
1240 void setEpoch(double epoch) { setD("Epoch", epoch); }
1241
1245 void setEpoch(std::string const &epoch) { setC("Epoch", epoch); }
1246
1250 void setFormat(int axis, std::string const &format) {
1251 setC(detail::formatAxisAttr("Format", axis), format);
1252 }
1253
1257 void setLabel(int axis, std::string const &label) { setC(detail::formatAxisAttr("Label", axis), label); }
1258
1262 void setMatchEnd(bool match) { setB("MatchEnd", match); }
1263
1268 void setMaxAxes(int maxAxes) { setI("MaxAxes", maxAxes); }
1269
1274 void setMinAxes(int minAxes) { setI("MinAxes", minAxes); }
1275
1279 void setObsAlt(double alt) { setD("ObsAlt", alt); }
1280
1284 void setObsLat(std::string const &lat) { setC("ObsLat", lat); }
1285
1289 void setObsLon(std::string const &lon) { setC("ObsLon", lon); }
1290
1295 void setActiveUnit(bool enable) {
1296 astSetActiveUnit(getRawPtr(), enable);
1297 assertOK();
1298 }
1299
1303 void setPermute(bool permute) { setB("Permute", permute); }
1304
1308 void setPreserveAxes(bool preserve) { setB("PreserveAxes", preserve); }
1309
1313 void setSymbol(int axis, std::string const &symbol) {
1314 setC(detail::formatAxisAttr("Symbol", axis), symbol);
1315 }
1316
1321 void setSystem(std::string const &system) { setC("System", system); }
1322
1326 void setTitle(std::string const &title) { setC("Title", title); }
1327
1331 void setTop(int axis, double top) { setD(detail::formatAxisAttr("Top", axis), top); }
1332
1336 void setUnit(int axis, std::string const &unit) { setC(detail::formatAxisAttr("Unit", axis), unit); }
1337
1502 NReadValue unformat(int axis, std::string const &str) const {
1503 double value;
1504 int nread = astUnformat(getRawPtr(), axis, str.c_str(), &value);
1505 assertOK();
1506 return NReadValue(nread, detail::safeDouble(value));
1507 }
1508
1509protected:
1519 explicit Frame(AstFrame *rawPtr) : Mapping(reinterpret_cast<AstMapping *>(rawPtr)) {
1520 if (!astIsAFrame(getRawPtr())) {
1522 os << "This is a " << getClassName() << ", which is not a Frame";
1523 throw std::invalid_argument(os.str());
1524 }
1525 }
1526
1527 virtual std::shared_ptr<Object> copyPolymorphic() const override { return copyImpl<Frame, AstFrame>(); }
1528
1529private:
1536 template <typename T>
1537 void assertPointLength(T const &p, char const *name) const {
1538 if (static_cast<int>(p.size()) != getNIn()) {
1540 os << "point " << name << " has " << p.size() << " axes, but " << getNIn() << " required";
1541 throw std::invalid_argument(os.str());
1542 }
1543 }
1544};
1545
1546} // namespace ast
1547
1548#endif
table::Key< double > angle
std::ostream * os
Definition Schema.cc:557
table::Key< int > b
A CmpFrame is a compound Frame which allows two component Frames (of any class) to be merged together...
Definition CmpFrame.h:60
Struct returned by Frame::offset2 containing a direction and a point.
Definition Frame.h:43
DirectionPoint(double direction, PointD const &point)
Construct a DirectionPoint.
Definition Frame.h:51
double direction
Direction, an angle in radians.
Definition Frame.h:52
PointD point
Point.
Definition Frame.h:53
Frame is used to represent a coordinate system.
Definition Frame.h:157
void permAxes(std::vector< int > perm)
Permute the order in which a Frame's axes occur.
Definition Frame.h:1138
NReadValue unformat(int axis, std::string const &str) const
Read a formatted coordinate value (given as a character string) for a Frame axis and return the numbe...
Definition Frame.h:1502
double axDistance(int axis, double v1, double v2) const
Return a signed value representing the axis increment from axis value v1 to axis value v2.
Definition Frame.h:247
std::string getFormat(int axis) const
Get Format for one axis: format specification for axis values.
Definition Frame.h:842
int getMaxAxes() const
Get MaxAxes: the maximum axes a frame found by findFrame may have.
Definition Frame.h:864
double getDut1() const
Get Dut1: difference between the UT1 and UTC timescale (sec)
Definition Frame.h:832
void setSymbol(int axis, std::string const &symbol)
Set Symbol(axis) for one axis: axis symbol.
Definition Frame.h:1313
void setActiveUnit(bool enable)
Set ActiveUnit: pay attention to units when one Frame is used to match another?
Definition Frame.h:1295
std::string getObsLat() const
Get ObsLat: Geodetic latitude of observer.
Definition Frame.h:892
int getMinAxes() const
Get MinAxes: the maximum axes a frame found by findFrame may have.
Definition Frame.h:870
void setMaxAxes(int maxAxes)
Get MaxAxes: the maximum number of axes a frame found by findFrame may have.
Definition Frame.h:1268
virtual void setDomain(std::string const &domain)
Set Domain: coordinate system domain.
Definition Frame.h:1230
double axOffset(int axis, double v1, double dist) const
Return an axis value formed by adding a signed axis increment onto a supplied axis value.
Definition Frame.h:263
std::string getSymbol(int axis) const
Get Symbol(axis) for one axis: axis symbol.
Definition Frame.h:912
double axAngle(PointD const &a, PointD const &b, int axis) const
Find the angle, as seen from point A, between the positive direction of a specified axis,...
Definition Frame.h:229
double getBottom(int axis) const
Get Bottom for one axis: the lowest axis value to display.
Definition Frame.h:807
void setPermute(bool permute)
Set Permute: allow axis permutation when used as a template?
Definition Frame.h:1303
virtual std::shared_ptr< Object > copyPolymorphic() const override
Return a deep copy of this object.
Definition Frame.h:1527
bool getPermute() const
Get Permute: allow axis permutation when used as a template?
Definition Frame.h:902
virtual ~Frame()
Definition Frame.h:173
int getDigits() const
Get Digits: the default used if no specific value specified for an axis.
Definition Frame.h:812
Frame(Frame const &)=default
Copy constructor: make a deep copy.
int getNAxes() const
Get NAxes: the number of axes in the frame (i.e.
Definition Frame.h:876
std::string getAlignSystem() const
Get AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition Frame.h:802
double distance(PointD const &point1, PointD const &point2) const
Find the distance between two points whose Frame coordinates are given.
Definition Frame.h:474
void setDigits(int digits)
Set Digits for all axes: number of digits of precision.
Definition Frame.h:1213
std::string getDomain() const
Get Domain: coordinate system domain.
Definition Frame.h:827
void setFormat(int axis, std::string const &format)
Set Format for one axis: format specification for axis values.
Definition Frame.h:1250
std::string getSystem() const
Get System: coordinate system used to describe positions within the domain.
Definition Frame.h:918
void setDirection(bool direction, int axis)
Set Direction for one axis: display axis in conventional direction?
Definition Frame.h:1223
std::shared_ptr< FrameSet > findFrame(Frame const &tmplt, std::string const &domainlist="")
Find a coordinate system with specified characteristics.
Definition Frame.cc:41
PointD offset(PointD point1, PointD point2, double offset) const
Find the point which is offset a specified distance along the geodesic curve between two other points...
Definition Frame.h:1075
double angle(PointD const &a, PointD const &b, PointD const &c) const
Find the angle at point B between the line joining points A and B, and the line joining points C and ...
Definition Frame.h:201
std::string getNormUnit(int axis) const
Get NormUnit(axis) read-only attribute for one frame: normalised physical units for formatted axis va...
Definition Frame.h:882
std::shared_ptr< Frame > copy() const
Return a deep copy of this object.
Definition Frame.h:182
Frame & operator=(Frame const &)=delete
double getObsAlt() const
Get ObsAlt: Geodetic altitude of observer (m).
Definition Frame.h:887
std::string getInternalUnit(int axis) const
Get InternalUnit(axis) read-only attribute for one axis: physical units for unformated axis values.
Definition Frame.h:848
void setPreserveAxes(bool preserve)
Set PreserveAxes: preserve axes?
Definition Frame.h:1308
void setMatchEnd(bool match)
Set MatchEnd: match trailing axes?
Definition Frame.h:1262
void setObsAlt(double alt)
Set ObsAlt: Geodetic altitude of observer (m).
Definition Frame.h:1279
void setAlignSystem(std::string const &system)
Set AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition Frame.h:1203
void setBottom(int axis, double bottom)
Set Bottom: the lowest axis value to display.
Definition Frame.h:1208
std::vector< double > intersect(std::vector< double > const &a1, std::vector< double > const &a2, std::vector< double > const &b1, std::vector< double > const &b2) const
Find the point of intersection between two geodesic curves.
Definition Frame.cc:51
std::vector< int > matchAxes(Frame const &other) const
Look for corresponding axes between this frame and another.
Definition Frame.h:985
bool getPreserveAxes() const
Get PreserveAxes: preserve axes?
Definition Frame.h:907
CmpFrame under(Frame const &next) const
Combine this frame with another to form a compound frame (CmpFrame), with the axes of this frame foll...
Definition Frame.cc:66
void setTitle(std::string const &title)
Set Title: frame title.
Definition Frame.h:1326
bool getMatchEnd() const
Get MatchEnd: match trailing axes?
Definition Frame.h:858
void setObsLon(std::string const &lon)
Set ObsLon: Geodetic longitude of observer.
Definition Frame.h:1289
void setTop(int axis, double top)
Set Top for one axis: the highest axis value to display.
Definition Frame.h:1331
Frame & operator=(Frame &&)=default
PointD norm(PointD value) const
Normalise a set of Frame coordinate values which might be unsuitable for display (e....
Definition Frame.h:1045
std::string getObsLon() const
Get ObsLon: Geodetic longitude of observer.
Definition Frame.h:897
void setDut1(double dut1)
Set Dut1: difference between the UT1 and UTC timescale (sec)
Definition Frame.h:1235
std::string getUnit(int axis) const
Get Unit(axis) for one axis: physical units for formatted axis values.
Definition Frame.h:933
double getTop(int axis) const
Get Top: the highest axis value to display.
Definition Frame.h:928
bool getActiveUnit() const
Get ActiveUnit: pay attention to units when one Frame is used to match another?
Definition Frame.h:792
bool getDirection(int axis) const
Get Direction for one axis: display axis in conventional direction?
Definition Frame.h:822
Frame(Frame &&)=default
void setUnit(int axis, std::string const &unit)
Set Unit(axis) for one axis: physical units for formatted axis values.
Definition Frame.h:1336
void setObsLat(std::string const &lat)
Set ObsLat: frame title.
Definition Frame.h:1284
DirectionPoint offset2(PointD const &point1, double angle, double offset) const
Find the point which is offset a specified distance along the geodesic curve at a given angle from a ...
Definition Frame.h:1117
std::string getLabel(int axis) const
Get Label(axis) for one axis: axis label.
Definition Frame.h:853
FrameMapping pickAxes(std::vector< int > const &axes) const
Create a new Frame whose axes are copied from an existing Frame along with other Frame attributes,...
Definition Frame.cc:68
void setLabel(int axis, std::string const &label)
Set Label(axis) for one axis: axis label.
Definition Frame.h:1257
void setSystem(std::string const &system)
Set System: coordinate system used to describe positions within the domain.
Definition Frame.h:1321
Frame(int naxes, std::string const &options="")
Construct a Frame.
Definition Frame.h:168
std::string format(int axis, double value) const
Return a string containing the formatted (character) version of a coordinate value for a Frame axis.
Definition Frame.h:782
void setDigits(int axis, int digits)
Set Digits for one axis: number of digits of precision.
Definition Frame.h:1218
void setEpoch(double epoch)
Set Epoch: Epoch of observation as a double (years)
Definition Frame.h:1240
std::string getTitle() const
Get Title: frame title.
Definition Frame.h:923
void setMinAxes(int minAxes)
Get MinAxes: the minimum number of axes a frame found by findFrame may have.
Definition Frame.h:1274
Frame(AstFrame *rawPtr)
Construct a Frame from a pointer to a raw AstFrame.
Definition Frame.h:1519
ResolvedPoint resolve(std::vector< double > const &point1, std::vector< double > const &point2, std::vector< double > const &point3) const
Resolve a vector into two orthogonal components.
Definition Frame.cc:84
void setEpoch(std::string const &epoch)
Set Epoch: Epoch of observation as a string.
Definition Frame.h:1245
int getDigits(int axis) const
Get Digits for one axis.
Definition Frame.h:817
double getEpoch() const
Get Epoch: Epoch of observation.
Definition Frame.h:837
std::shared_ptr< FrameSet > convert(Frame const &to, std::string const &domainlist="")
Compute a frameset that describes the conversion between this frame and another frame.
Definition Frame.cc:31
Struct returned by Frame::pickAxes containing a frame and a mapping.
Definition Frame.h:91
FrameMapping(std::shared_ptr< Frame > frame, std::shared_ptr< Mapping > mapping)
Construct a FrameMapping.
Definition Frame.h:99
std::shared_ptr< Mapping > mapping
Mapping.
Definition Frame.h:102
std::shared_ptr< Frame > frame
Frame.
Definition Frame.h:101
A FrameSet consists of a set of one or more Frames (which describe coordinate systems),...
Definition FrameSet.h:99
An abstract base class for objects which transform one set of coordinates to another.
Definition Mapping.h:59
int getNIn() const
Get NIn: the number of input axes.
Definition Mapping.h:77
Struct returned by Frame::unformat containing the number of characters read and corresponding value.
Definition Frame.h:59
NReadValue(int nread, double value)
Construct an NReadValue.
Definition Frame.h:67
double value
Value that was read.
Definition Frame.h:69
int nread
Number of characters that was read.
Definition Frame.h:68
Abstract base class for all AST objects.
Definition Object.h:51
void setD(std::string const &attrib, double value)
Set the value of an attribute as a double.
Definition Object.h:472
double getD(std::string const &attrib) const
Get the value of an attribute as a double.
Definition Object.h:374
std::string const getC(std::string const &attrib) const
Get the value of an attribute as a string.
Definition Object.h:361
std::string getClassName() const
Get Class: the name of the class (e.g.
Definition Object.h:139
void setI(std::string const &attrib, int value)
Set the value of an attribute as an int.
Definition Object.h:496
bool getB(std::string const &attrib) const
Get the value of an attribute as a bool.
Definition Object.h:348
int getI(std::string const &attrib) const
Get the value of an attribute as an int.
Definition Object.cc:178
void setC(std::string const &attrib, std::string const &value)
Set the value of an attribute as a string.
Definition Object.h:460
void setB(std::string const &attrib, bool value)
Set the value of an attribute as a bool.
Definition Object.h:448
AstObject const * getRawPtr() const
Get the raw AST pointer.
Definition Object.h:292
Struct returned by Frame::resolve containing a point and the resolved vector components.
Definition Frame.h:75
std::vector< double > point
Point.
Definition Frame.h:83
ResolvedPoint(int naxes)
Construct an empty ResolvedPoint.
Definition Frame.h:82
double d2
Resolved vector component 2.
Definition Frame.h:85
double d1
Resolved vector component 1.
Definition Frame.h:84
T data(T... args)
std::string formatAxisAttr(std::string const &name, int axis)
Format an axis-specific attribute by appending the axis index.
Definition utils.h:79
void assertEqual(T1 val1, std::string const &descr1, T2 val2, std::string const &descr2)
Definition utils.h:48
double safeDouble(double val)
Return a double value after checking status and replacing AST__BAD with nan
Definition utils.h:100
void astBadToNan(std::vector< double > &p)
Replace AST__BAD with a quiet NaN in a vector.
Definition utils.h:59
AST wrapper classes and functions.
void assertOK(AstObject *rawPtr1=nullptr, AstObject *rawPtr2=nullptr)
Throw std::runtime_error if AST's state is bad.
Definition base.cc:49
T size(T... args)