LSSTApplications  19.0.0+10,19.0.0+80,19.0.0-1-g20d9b18+31,19.0.0-1-g49a97f9+4,19.0.0-1-g8c57eb9+31,19.0.0-1-g9a028c0+13,19.0.0-1-ga72da6b+3,19.0.0-1-gb77924a+15,19.0.0-1-gbfe0924+66,19.0.0-1-gc3516c3,19.0.0-1-gefe1d0d+49,19.0.0-1-gf8cb8b4+3,19.0.0-14-g7511ce4+6,19.0.0-16-g3dc1a33c+6,19.0.0-17-g59f0e8a+4,19.0.0-17-g9c22e3c+9,19.0.0-18-g35bb99870+2,19.0.0-19-g2772d4a+9,19.0.0-2-g260436e+53,19.0.0-2-g31cdcee,19.0.0-2-g9675b69+3,19.0.0-2-gaa2795f,19.0.0-2-gbcc4de1,19.0.0-2-gd6f004e+6,19.0.0-2-gde8e5e3+5,19.0.0-2-gff6972b+15,19.0.0-21-ge306cd8,19.0.0-21-gf122e69+2,19.0.0-3-g2d43a51+2,19.0.0-3-gf3b1435+6,19.0.0-4-g56feb96,19.0.0-4-gac56cce+17,19.0.0-49-gce872c1+1,19.0.0-5-g66946eb+6,19.0.0-5-gd8897ba+6,19.0.0-51-gfc4a647e,19.0.0-7-g686a884+5,19.0.0-7-g886f805+5,19.0.0-8-g305ff64,w.2020.17
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  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 
966  std::vector<double> intersect(std::vector<double> const &a1, std::vector<double> const &a2,
967  std::vector<double> const &b1, std::vector<double> const &b2) const;
968 
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();
1081  detail::astBadToNan(ret);
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 
1509 protected:
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 
1529 private:
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
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
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
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
std::string getDomain() const
Get Domain: coordinate system domain.
Definition: Frame.h:827
void setDigits(int digits)
Set Digits for all axes: number of digits of precision.
Definition: Frame.h:1213
void setObsLat(std::string const &lat)
Set ObsLat: frame title.
Definition: Frame.h:1284
double getTop(int axis) const
Get Top: the highest axis value to display.
Definition: Frame.h:928
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:292
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:1245
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:1308
double distance(PointD const &point1, PointD const &point2) const
Find the distance between two points whose Frame coordinates are given.
Definition: Frame.h:474
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:247
void setMaxAxes(int maxAxes)
Get MaxAxes: the maximum number of axes a frame found by findFrame may have.
Definition: Frame.h:1268
table::Key< int > b
void setDut1(double dut1)
Set Dut1: difference between the UT1 and UTC timescale (sec)
Definition: Frame.h:1235
std::shared_ptr< Frame > copy() const
Return a deep copy of this object.
Definition: Frame.h:182
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
table::Key< int > a
void setSymbol(int axis, std::string const &symbol)
Set Symbol(axis) for one axis: axis symbol.
Definition: Frame.h:1313
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:807
table::Key< double > angle
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:229
ItemVariant const * other
Definition: Schema.cc:56
void setFormat(int axis, std::string const &format)
Set Format for one axis: format specification for axis values.
Definition: Frame.h:1250
void permAxes(std::vector< int > perm)
Permute the order in which a Frame&#39;s axes occur.
Definition: Frame.h:1138
void astBadToNan(std::vector< double > &p)
Replace AST__BAD with a quiet NaN in a vector.
Definition: utils.h:59
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:1331
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:858
std::string getFormat(int axis) const
Get Format for one axis: format specification for axis values.
Definition: Frame.h:842
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:1240
table::Key< int > to
std::string getSystem() const
Get System: coordinate system used to describe positions within the domain.
Definition: Frame.h:918
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:876
std::shared_ptr< Mapping > mapping
Mapping.
Definition: Frame.h:102
std::string getUnit(int axis) const
Get Unit(axis) for one axis: physical units for formatted axis values.
Definition: Frame.h:933
void setBottom(int axis, double bottom)
Set Bottom: the lowest axis value to display.
Definition: Frame.h:1208
T str(T... args)
Frame(AstFrame *rawPtr)
Construct a Frame from a pointer to a raw AstFrame.
Definition: Frame.h:1519
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:1257
std::shared_ptr< Frame > frame
Frame.
Definition: Frame.h:101
void setObsAlt(double alt)
Set ObsAlt: Geodetic altitude of observer (m).
Definition: Frame.h:1279
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
int getDigits() const
Get Digits: the default used if no specific value specified for an axis.
Definition: Frame.h:812
virtual void setDomain(std::string const &domain)
Set Domain: coordinate system domain.
Definition: Frame.h:1230
bool getDirection(int axis) const
Get Direction for one axis: display axis in conventional direction?
Definition: Frame.h:822
std::string getObsLon() const
Get ObsLon: Geodetic longitude of observer.
Definition: Frame.h:897
int getDigits(int axis) const
Get Digits for one axis.
Definition: Frame.h:817
double getObsAlt() const
Get ObsAlt: Geodetic altitude of observer (m).
Definition: Frame.h:887
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:201
int getMaxAxes() const
Get MaxAxes: the maximum axes a frame found by findFrame may have.
Definition: Frame.h:864
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
double getEpoch() const
Get Epoch: Epoch of observation.
Definition: Frame.h:837
std::string getTitle() const
Get Title: frame title.
Definition: Frame.h:923
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:907
std::string getLabel(int axis) const
Get Label(axis) for one axis: axis label.
Definition: Frame.h:853
def convert(gen2root, gen3root, instrumentClass, calibFilterType, skymapName=None, skymapConfig=None, calibs=None, reruns=[], config=None, transferMode="auto")
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:832
void setObsLon(std::string const &lon)
Set ObsLon: Geodetic longitude of observer.
Definition: Frame.h:1289
std::string formatAxisAttr(std::string const &name, int axis)
Format an axis-specific attribute by appending the axis index.
Definition: utils.h:79
void setMatchEnd(bool match)
Set MatchEnd: match trailing axes?
Definition: Frame.h:1262
bool getPermute() const
Get Permute: allow axis permutation when used as a template?
Definition: Frame.h:902
double safeDouble(double val)
Return a double value after checking status and replacing AST__BAD with nan
Definition: utils.h:100
std::string getAlignSystem() const
Get AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition: Frame.h:802
void assertEqual(T1 val1, std::string const &descr1, T2 val2, std::string const &descr2)
Definition: utils.h:48
virtual ~Frame()
Definition: Frame.h:173
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 setSystem(std::string const &system)
Set System: coordinate system used to describe positions within the domain.
Definition: Frame.h:1321
std::string getObsLat() const
Get ObsLat: Geodetic latitude of observer.
Definition: Frame.h:892
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:1223
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::ostream * os
Definition: Schema.cc:746
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
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
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:912
void setMinAxes(int minAxes)
Get MinAxes: the minimum number of axes a frame found by findFrame may have.
Definition: Frame.h:1274
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:1295
void setAlignSystem(std::string const &system)
Set AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition: Frame.h:1203
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:792
virtual std::shared_ptr< Object > copyPolymorphic() const override
Return a deep copy of this object.
Definition: Frame.h:1527
void setPermute(bool permute)
Set Permute: allow axis permutation when used as a template?
Definition: Frame.h:1303
void setDigits(int axis, int digits)
Set Digits for one axis: number of digits of precision.
Definition: Frame.h:1218
void setTitle(std::string const &title)
Set Title: frame title.
Definition: Frame.h:1326
std::vector< int > matchAxes(Frame const &other) const
Look for corresponding axes between this frame and another.
Definition: Frame.h:985
PointD norm(PointD value) const
Normalise a set of Frame coordinate values which might be unsuitable for display (e.g.
Definition: Frame.h:1045
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
int getMinAxes() const
Get MinAxes: the maximum axes a frame found by findFrame may have.
Definition: Frame.h:870