LSSTApplications  11.0-24-g0a022a1,14.0+77,15.0,15.0+1
LSSTDataManagementBasePackage
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"
33 #include "astshim/detail/utils.h"
34 
35 namespace ast {
36 
37 class CmpFrame;
38 class Frame;
39 
44 public:
51  DirectionPoint(double direction, PointD const &point) : direction(direction), point(point){};
52  double direction;
54 };
55 
59 class NReadValue {
60 public:
67  NReadValue(int nread, double value) : nread(nread), value(value){};
68  int nread;
69  double value;
70 };
71 
76 public:
82  explicit ResolvedPoint(int naxes) : point(naxes), d1(), d2() {}
84  double d1;
85  double d2;
86 };
87 
91 class FrameMapping {
92 public:
100  : frame(frame), mapping(mapping) {}
103 };
104 
105 class FrameSet;
106 
157 class Frame : public Mapping {
158  friend class Object;
159 
160 public:
168  explicit Frame(int naxes, std::string const &options = "")
169  : Mapping(reinterpret_cast<AstMapping *>(astFrame(naxes, "%s", options.c_str()))) {}
170 
171  virtual ~Frame() {}
172 
174  Frame(Frame const &) = default;
175  Frame(Frame &&) = default;
176  Frame &operator=(Frame const &) = delete;
177  Frame &operator=(Frame &&) = default;
178 
180  std::shared_ptr<Frame> copy() const { return std::static_pointer_cast<Frame>(copyPolymorphic()); }
181 
199  double angle(PointD const &a, PointD const &b, PointD const &c) const {
200  assertPointLength(a, "a");
201  assertPointLength(b, "b");
202  assertPointLength(c, "c");
203  return detail::safeDouble(astAngle(getRawPtr(), a.data(), b.data(), c.data()));
204  }
205 
227  double axAngle(PointD const &a, PointD const &b, int axis) const {
228  assertPointLength(a, "a");
229  assertPointLength(b, "b");
230  return detail::safeDouble(astAxAngle(getRawPtr(), a.data(), b.data(), axis));
231  }
232 
245  double axDistance(int axis, double v1, double v2) const {
246  return detail::safeDouble(astAxDistance(getRawPtr(), axis, v1, v2));
247  }
248 
261  double axOffset(int axis, double v1, double dist) const {
262  return detail::safeDouble(astAxOffset(getRawPtr(), axis, v1, dist));
263  }
264 
458  std::shared_ptr<FrameSet> convert(Frame const &to, std::string const &domainlist = "");
459 
473  double distance(PointD const &point1, PointD const &point2) const {
474  assertPointLength(point1, "point1");
475  assertPointLength(point2, "point2");
476  return detail::safeDouble(astDistance(getRawPtr(), point1.data(), point2.data()));
477  }
478 
768  std::shared_ptr<FrameSet> findFrame(Frame const &tmplt, std::string const &domainlist = "");
769 
781  std::string format(int axis, double value) const {
782  char const *rawstr = astFormat(getRawPtr(), axis, value);
783  assertOK();
784  return std::string(rawstr);
785  }
786 
791  bool getActiveUnit() const {
792  bool ret = astGetActiveUnit(getRawPtr());
793  assertOK();
794  return ret;
795  }
796 
801  std::string getAlignSystem() const { return getC("AlignSystem"); }
802 
806  double getBottom(int axis) const { return getD(detail::formatAxisAttr("Bottom", axis)); }
807 
811  int getDigits() const { return getI("Digits"); }
812 
816  int getDigits(int axis) const { return getI(detail::formatAxisAttr("Digits", axis)); }
817 
821  bool getDirection(int axis) const { return getB(detail::formatAxisAttr("Direction", axis)); }
822 
826  std::string getDomain() const { return getC("Domain"); }
827 
831  double getDut1() const { return getD("Dut1"); }
832 
836  double getEpoch() const { return getD("Epoch"); }
837 
841  std::string getFormat(int axis) const { return getC(detail::formatAxisAttr("Format", axis)); }
842 
847  std::string getInternalUnit(int axis) const { return getC(detail::formatAxisAttr("InternalUnit", axis)); }
848 
852  std::string getLabel(int axis) const { return getC(detail::formatAxisAttr("Label", axis)); }
853 
857  bool getMatchEnd() const { return getB("MatchEnd"); }
858 
863  int getMaxAxes() const { return getI("MaxAxes"); }
864 
869  int getMinAxes() const { return getI("MinAxes"); }
870 
875  int getNAxes() const { return getI("NAxes"); }
876 
881  std::string getNormUnit(int axis) const { return getC(detail::formatAxisAttr("NormUnit", axis)); }
882 
886  double getObsAlt() const { return getD("ObsAlt"); }
887 
891  std::string getObsLat() const { return getC("ObsLat"); }
892 
896  std::string getObsLon() const { return getC("ObsLon"); }
897 
901  bool getPermute() const { return getB("Permute"); }
902 
906  bool getPreserveAxes() const { return getB("PreserveAxes"); }
907 
911  std::string getSymbol(int axis) const { return getC(detail::formatAxisAttr("Symbol", axis)); }
912 
917  std::string getSystem() const { return getC("System"); }
918 
922  std::string getTitle() const { return getC("Title"); }
923 
927  double getTop(int axis) const { return getD(detail::formatAxisAttr("Top", axis)); }
928 
932  std::string getUnit(int axis) const { return getC(detail::formatAxisAttr("Unit", axis)); }
933 
965  std::vector<double> intersect(std::vector<double> const &a1, std::vector<double> const &a2,
966  std::vector<double> const &b1, std::vector<double> const &b2) const;
967 
984  std::vector<int> matchAxes(Frame const &other) const {
985  std::vector<int> ret(other.getNIn());
986  astMatchAxes(getRawPtr(), other.getRawPtr(), ret.data());
987  assertOK();
988  return ret;
989  }
990 
1022  CmpFrame under(Frame const &next) const;
1023 
1044  PointD norm(PointD value) const {
1045  astNorm(getRawPtr(), value.data());
1046  assertOK();
1047  detail::astBadToNan(value);
1048  return value;
1049  }
1050 
1074  PointD offset(PointD point1, PointD point2, double offset) const {
1075  assertPointLength(point1, "point1");
1076  assertPointLength(point2, "point2");
1077  PointD ret(getNIn());
1078  astOffset(getRawPtr(), point1.data(), point2.data(), offset, ret.data());
1079  assertOK();
1080  detail::astBadToNan(ret);
1081  return ret;
1082  }
1083 
1116  DirectionPoint offset2(PointD const &point1, double angle, double offset) const {
1117  detail::assertEqual(getNIn(), "naxes", 2, " cannot call offset2");
1118  assertPointLength(point1, "point1");
1119  PointD point2(getNIn());
1120  double offsetAngle = astOffset2(getRawPtr(), point1.data(), angle, offset, point2.data());
1121  assertOK();
1122  detail::astBadToNan(point2);
1123  return DirectionPoint(detail::safeDouble(offsetAngle), point2);
1124  }
1125 
1138  detail::assertEqual(perm.size(), "perm.size()", static_cast<std::size_t>(getNAxes()), "naxes");
1139  astPermAxes(getRawPtr(), perm.data());
1140  assertOK();
1141  }
1142 
1162  FrameMapping pickAxes(std::vector<int> const &axes) const;
1163 
1195  ResolvedPoint resolve(std::vector<double> const &point1, std::vector<double> const &point2,
1196  std::vector<double> const &point3) const;
1197 
1202  void setAlignSystem(std::string const &system) { setC("AlignSystem", system); }
1203 
1207  void setBottom(int axis, double bottom) { setD(detail::formatAxisAttr("Bottom", axis), bottom); }
1208 
1212  void setDigits(int digits) { setI("Digits", digits); }
1213 
1217  void setDigits(int axis, int digits) { setD(detail::formatAxisAttr("Digits", axis), digits); }
1218 
1222  void setDirection(bool direction, int axis) {
1223  setB(detail::formatAxisAttr("Direction", axis), direction);
1224  }
1225 
1229  virtual void setDomain(std::string const &domain) { setC("Domain", domain); }
1230 
1234  void setDut1(double dut1) { setD("Dut1", dut1); }
1235 
1239  void setEpoch(double epoch) { setD("Epoch", epoch); }
1240 
1244  void setEpoch(std::string const &epoch) { setC("Epoch", epoch); }
1245 
1249  void setFormat(int axis, std::string const &format) {
1250  setC(detail::formatAxisAttr("Format", axis), format);
1251  }
1252 
1256  void setLabel(int axis, std::string const &label) { setC(detail::formatAxisAttr("Label", axis), label); }
1257 
1261  void setMatchEnd(bool match) { setB("MatchEnd", match); }
1262 
1267  void setMaxAxes(int maxAxes) { setI("MaxAxes", maxAxes); }
1268 
1273  void setMinAxes(int minAxes) { setI("MinAxes", minAxes); }
1274 
1278  void setObsAlt(double alt) { setD("ObsAlt", alt); }
1279 
1283  void setObsLat(std::string const &lat) { setC("ObsLat", lat); }
1284 
1288  void setObsLon(std::string const &lon) { setC("ObsLon", lon); }
1289 
1294  void setActiveUnit(bool enable) {
1295  astSetActiveUnit(getRawPtr(), enable);
1296  assertOK();
1297  }
1298 
1302  void setPermute(bool permute) { setB("Permute", permute); }
1303 
1307  void setPreserveAxes(bool preserve) { setB("PreserveAxes", preserve); }
1308 
1312  void setSymbol(int axis, std::string const &symbol) {
1313  setC(detail::formatAxisAttr("Symbol", axis), symbol);
1314  }
1315 
1320  void setSystem(std::string const &system) { setC("System", system); }
1321 
1325  void setTitle(std::string const &title) { setC("Title", title); }
1326 
1330  void setTop(int axis, double top) { setD(detail::formatAxisAttr("Top", axis), top); }
1331 
1335  void setUnit(int axis, std::string const &unit) { setC(detail::formatAxisAttr("Unit", axis), unit); }
1336 
1501  NReadValue unformat(int axis, std::string const &str) const {
1502  double value;
1503  int nread = astUnformat(getRawPtr(), axis, str.c_str(), &value);
1504  assertOK();
1505  return NReadValue(nread, detail::safeDouble(value));
1506  }
1507 
1508 protected:
1518  explicit Frame(AstFrame *rawPtr) : Mapping(reinterpret_cast<AstMapping *>(rawPtr)) {
1519  if (!astIsAFrame(getRawPtr())) {
1520  std::ostringstream os;
1521  os << "This is a " << getClassName() << ", which is not a Frame";
1522  throw std::invalid_argument(os.str());
1523  }
1524  }
1525 
1526  virtual std::shared_ptr<Object> copyPolymorphic() const override { return copyImpl<Frame, AstFrame>(); }
1527 
1528 private:
1535  template <typename T>
1536  void assertPointLength(T const &p, char const *name) const {
1537  if (static_cast<int>(p.size()) != getNIn()) {
1538  std::ostringstream os;
1539  os << "point " << name << " has " << p.size() << " axes, but " << getNIn() << " required";
1540  throw std::invalid_argument(os.str());
1541  }
1542  }
1543 };
1544 
1545 } // namespace ast
1546 
1547 #endif
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:261
A CmpFrame is a compound Frame which allows two component Frames (of any class) to be merged together...
Definition: CmpFrame.h:60
Frame(int naxes, std::string const &options="")
Construct a Frame.
Definition: Frame.h:168
std::string getDomain() const
Get Domain: coordinate system domain.
Definition: Frame.h:826
void setDigits(int digits)
Set Digits for all axes: number of digits of precision.
Definition: Frame.h:1212
void setObsLat(std::string const &lat)
Set ObsLat: frame title.
Definition: Frame.h:1283
double getTop(int axis) const
Get Top: the highest axis value to display.
Definition: Frame.h:927
Struct returned by Frame::resolve containing a point and the resolved vector components.
Definition: Frame.h:75
AstObject const * getRawPtr() const
Get the raw AST pointer.
Definition: Object.h:291
table::Key< std::string > name
Definition: ApCorrMap.cc:71
double value
Value that was read.
Definition: Frame.h:69
void setEpoch(std::string const &epoch)
Set Epoch: Epoch of observation as a string.
Definition: Frame.h:1244
Struct returned by Frame::pickAxes containing a frame and a mapping.
Definition: Frame.h:91
void setPreserveAxes(bool preserve)
Set PreserveAxes: preserve axes?
Definition: Frame.h:1307
double distance(PointD const &point1, PointD const &point2) const
Find the distance between two points whose Frame coordinates are given.
Definition: Frame.h:473
AST wrapper classes and functions.
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:245
void setMaxAxes(int maxAxes)
Get MaxAxes: the maximum number of axes a frame found by findFrame may have.
Definition: Frame.h:1267
void setDut1(double dut1)
Set Dut1: difference between the UT1 and UTC timescale (sec)
Definition: Frame.h:1234
std::shared_ptr< Frame > copy() const
Return a deep copy of this object.
Definition: Frame.h:180
double d2
Resolved vector component 2.
Definition: Frame.h:85
void assertOK(AstObject *rawPtr1=nullptr, AstObject *rawPtr2=nullptr)
Throw std::runtime_error if AST&#39;s state is bad.
Definition: base.cc:49
void setSymbol(int axis, std::string const &symbol)
Set Symbol(axis) for one axis: axis symbol.
Definition: Frame.h:1312
DirectionPoint(double direction, PointD const &point)
Construct a DirectionPoint.
Definition: Frame.h:51
double getBottom(int axis) const
Get Bottom for one axis: the lowest axis value to display.
Definition: Frame.h:806
tuple options
Definition: lsstimport.py:53
std::vector< double > point
Point.
Definition: Frame.h:83
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:227
void setFormat(int axis, std::string const &format)
Set Format for one axis: format specification for axis values.
Definition: Frame.h:1249
void permAxes(std::vector< int > perm)
Permute the order in which a Frame&#39;s axes occur.
Definition: Frame.h:1137
void astBadToNan(std::vector< double > &p)
Replace AST__BAD with a quiet NaN in a vector.
Definition: utils.h:57
double direction
Direction, an angle in radians.
Definition: Frame.h:51
void setTop(int axis, double top)
Set Top for one axis: the highest axis value to display.
Definition: Frame.h:1330
An abstract base class for objects which transform one set of coordinates to another.
Definition: Mapping.h:59
bool getMatchEnd() const
Get MatchEnd: match trailing axes?
Definition: Frame.h:857
std::string getFormat(int axis) const
Get Format for one axis: format specification for axis values.
Definition: Frame.h:841
STL class.
PointD point
Point.
Definition: Frame.h:53
T data(T... args)
void setEpoch(double epoch)
Set Epoch: Epoch of observation as a double (years)
Definition: Frame.h:1239
std::string getSystem() const
Get System: coordinate system used to describe positions within the domain.
Definition: Frame.h:917
Struct returned by Frame::unformat containing the number of characters read and corresponding value...
Definition: Frame.h:59
Frame is used to represent a coordinate system.
Definition: Frame.h:157
int getNAxes() const
Get NAxes: the number of axes in the frame (i.e.
Definition: Frame.h:875
std::shared_ptr< Mapping > mapping
Mapping.
Definition: Frame.h:102
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:134
std::string getUnit(int axis) const
Get Unit(axis) for one axis: physical units for formatted axis values.
Definition: Frame.h:932
table::Key< int > b
void setBottom(int axis, double bottom)
Set Bottom: the lowest axis value to display.
Definition: Frame.h:1207
T str(T... args)
Frame(AstFrame *rawPtr)
Construct a Frame from a pointer to a raw AstFrame.
Definition: Frame.h:1518
T static_pointer_cast(T... args)
double d1
Resolved vector component 1.
Definition: Frame.h:84
void setLabel(int axis, std::string const &label)
Set Label(axis) for one axis: axis label.
Definition: Frame.h:1256
std::shared_ptr< Frame > frame
Frame.
Definition: Frame.h:101
void setObsAlt(double alt)
Set ObsAlt: Geodetic altitude of observer (m).
Definition: Frame.h:1278
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:1074
int getDigits() const
Get Digits: the default used if no specific value specified for an axis.
Definition: Frame.h:811
virtual void setDomain(std::string const &domain)
Set Domain: coordinate system domain.
Definition: Frame.h:1229
bool getDirection(int axis) const
Get Direction for one axis: display axis in conventional direction?
Definition: Frame.h:821
std::string getObsLon() const
Get ObsLon: Geodetic longitude of observer.
Definition: Frame.h:896
int getDigits(int axis) const
Get Digits for one axis.
Definition: Frame.h:816
double getObsAlt() const
Get ObsAlt: Geodetic altitude of observer (m).
Definition: Frame.h:886
std::string getClassName(AstObject const *rawObj)
Get the AST class name, changing CmpMap to SeriesMap or ParallelMap as appropriate.
Definition: utils.cc:37
T size(T... args)
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:199
int getMaxAxes() const
Get MaxAxes: the maximum axes a frame found by findFrame may have.
Definition: Frame.h:863
std::string getInternalUnit(int axis) const
Get InternalUnit(axis) read-only attribute for one axis: physical units for unformated axis values...
Definition: Frame.h:847
double getEpoch() const
Get Epoch: Epoch of observation.
Definition: Frame.h:836
std::string getTitle() const
Get Title: frame title.
Definition: Frame.h:922
int getNIn() const
Get NIn: the number of input axes.
Definition: Mapping.h:77
T c_str(T... args)
bool getPreserveAxes() const
Get PreserveAxes: preserve axes?
Definition: Frame.h:906
std::string getLabel(int axis) const
Get Label(axis) for one axis: axis label.
Definition: Frame.h:852
ResolvedPoint(int naxes)
Construct an empty ResolvedPoint.
Definition: Frame.h:82
double getDut1() const
Get Dut1: difference between the UT1 and UTC timescale (sec)
Definition: Frame.h:831
void setObsLon(std::string const &lon)
Set ObsLon: Geodetic longitude of observer.
Definition: Frame.h:1288
std::string formatAxisAttr(std::string const &name, int axis)
Format an axis-specific attribute by appending the axis index.
Definition: utils.h:77
void setMatchEnd(bool match)
Set MatchEnd: match trailing axes?
Definition: Frame.h:1261
bool getPermute() const
Get Permute: allow axis permutation when used as a template?
Definition: Frame.h:901
table::Key< int > a
double safeDouble(double val)
Return a double value after checking status and replacing AST__BAD with nan
Definition: utils.h:98
std::string getAlignSystem() const
Get AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition: Frame.h:801
void assertEqual(T1 val1, std::string const &descr1, T2 val2, std::string const &descr2)
Definition: utils.h:46
virtual ~Frame()
Definition: Frame.h:171
void setUnit(int axis, std::string const &unit)
Set Unit(axis) for one axis: physical units for formatted axis values.
Definition: Frame.h:1335
void setSystem(std::string const &system)
Set System: coordinate system used to describe positions within the domain.
Definition: Frame.h:1320
std::string getObsLat() const
Get ObsLat: Geodetic latitude of observer.
Definition: Frame.h:891
FrameMapping(std::shared_ptr< Frame > frame, std::shared_ptr< Mapping > mapping)
Construct a FrameMapping.
Definition: Frame.h:99
Struct returned by Frame::offset2 containing a direction and a point.
Definition: Frame.h:43
void setDirection(bool direction, int axis)
Set Direction for one axis: display axis in conventional direction?
Definition: Frame.h:1222
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:1116
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:1501
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:781
Abstract base class for all AST objects.
Definition: Object.h:51
std::string getSymbol(int axis) const
Get Symbol(axis) for one axis: axis symbol.
Definition: Frame.h:911
void setMinAxes(int minAxes)
Get MinAxes: the minimum number of axes a frame found by findFrame may have.
Definition: Frame.h:1273
A FrameSet consists of a set of one or more Frames (which describe coordinate systems), connected together by Mappings (which describe how the coordinate systems are inter-related).
Definition: FrameSet.h:99
void setActiveUnit(bool enable)
Set ActiveUnit: pay attention to units when one Frame is used to match another?
Definition: Frame.h:1294
void setAlignSystem(std::string const &system)
Set AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition: Frame.h:1202
NReadValue(int nread, double value)
Construct an NReadValue.
Definition: Frame.h:67
bool getActiveUnit() const
Get ActiveUnit: pay attention to units when one Frame is used to match another?
Definition: Frame.h:791
virtual std::shared_ptr< Object > copyPolymorphic() const override
Return a deep copy of this object.
Definition: Frame.h:1526
void setPermute(bool permute)
Set Permute: allow axis permutation when used as a template?
Definition: Frame.h:1302
void setDigits(int axis, int digits)
Set Digits for one axis: number of digits of precision.
Definition: Frame.h:1217
void assertPointLength(T const &p, char const *name) const
Assert that a point has the correct length.
Definition: Frame.h:1536
void setTitle(std::string const &title)
Set Title: frame title.
Definition: Frame.h:1325
std::vector< int > matchAxes(Frame const &other) const
Look for corresponding axes between this frame and another.
Definition: Frame.h:984
PointD norm(PointD value) const
Normalise a set of Frame coordinate values which might be unsuitable for display (e.g.
Definition: Frame.h:1044
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:881
int getMinAxes() const
Get MinAxes: the maximum axes a frame found by findFrame may have.
Definition: Frame.h:869