LSSTApplications  16.0+42,16.0-1-gce273f5+8,16.0-10-g230e10e+1,16.0-11-g9fe0e56+17,16.0-11-gce733cf+17,16.0-12-g5ad1ebf+9,16.0-12-gc85596e+2,16.0-13-gde155d7+2,16.0-14-g9428de4d,16.0-14-gc1cf4a94+2,16.0-15-g8e16a51+14,16.0-2-g0febb12+7,16.0-2-g839ba83+32,16.0-2-g9d5294e+22,16.0-2-gab3db49+7,16.0-2-gf41ba6b+6,16.0-2-gf4e7cdd+5,16.0-3-g6923fb6+15,16.0-3-g8e51203+2,16.0-3-g9645794+6,16.0-3-gcfd6c53+20,16.0-35-g34c7dfe62+1,16.0-4-g03cf288+11,16.0-4-g32d12de,16.0-4-g5f3a788+7,16.0-4-g7690030+30,16.0-4-g8a0f11a+16,16.0-4-ga5d8928+16,16.0-5-g0da18be+7,16.0-5-g4940a70,16.0-5-g563880a+2,16.0-5-g7742071+2,16.0-5-gb3f8a4b+26,16.0-6-g3610b4f+5,16.0-6-gf0acd13+14,16.0-8-g4dec96c+7,16.0-8-gc315727+16,16.0-9-g1de645c+7,16.0-9-gcc4efb7+6,w.2018.36
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 
457  std::shared_ptr<FrameSet> convert(Frame const &to, std::string const &domainlist = "");
458 
472  double distance(PointD const &point1, PointD const &point2) const {
473  assertPointLength(point1, "point1");
474  assertPointLength(point2, "point2");
475  return detail::safeDouble(astDistance(getRawPtr(), point1.data(), point2.data()));
476  }
477 
767  std::shared_ptr<FrameSet> findFrame(Frame const &tmplt, std::string const &domainlist = "");
768 
780  std::string format(int axis, double value) const {
781  char const *rawstr = astFormat(getRawPtr(), axis, value);
782  assertOK();
783  return std::string(rawstr);
784  }
785 
790  bool getActiveUnit() const {
791  bool ret = astGetActiveUnit(getRawPtr());
792  assertOK();
793  return ret;
794  }
795 
800  std::string getAlignSystem() const { return getC("AlignSystem"); }
801 
805  double getBottom(int axis) const { return getD(detail::formatAxisAttr("Bottom", axis)); }
806 
810  int getDigits() const { return getI("Digits"); }
811 
815  int getDigits(int axis) const { return getI(detail::formatAxisAttr("Digits", axis)); }
816 
820  bool getDirection(int axis) const { return getB(detail::formatAxisAttr("Direction", axis)); }
821 
825  std::string getDomain() const { return getC("Domain"); }
826 
830  double getDut1() const { return getD("Dut1"); }
831 
835  double getEpoch() const { return getD("Epoch"); }
836 
840  std::string getFormat(int axis) const { return getC(detail::formatAxisAttr("Format", axis)); }
841 
846  std::string getInternalUnit(int axis) const { return getC(detail::formatAxisAttr("InternalUnit", axis)); }
847 
851  std::string getLabel(int axis) const { return getC(detail::formatAxisAttr("Label", axis)); }
852 
856  bool getMatchEnd() const { return getB("MatchEnd"); }
857 
862  int getMaxAxes() const { return getI("MaxAxes"); }
863 
868  int getMinAxes() const { return getI("MinAxes"); }
869 
874  int getNAxes() const { return getI("NAxes"); }
875 
880  std::string getNormUnit(int axis) const { return getC(detail::formatAxisAttr("NormUnit", axis)); }
881 
885  double getObsAlt() const { return getD("ObsAlt"); }
886 
890  std::string getObsLat() const { return getC("ObsLat"); }
891 
895  std::string getObsLon() const { return getC("ObsLon"); }
896 
900  bool getPermute() const { return getB("Permute"); }
901 
905  bool getPreserveAxes() const { return getB("PreserveAxes"); }
906 
910  std::string getSymbol(int axis) const { return getC(detail::formatAxisAttr("Symbol", axis)); }
911 
916  std::string getSystem() const { return getC("System"); }
917 
921  std::string getTitle() const { return getC("Title"); }
922 
926  double getTop(int axis) const { return getD(detail::formatAxisAttr("Top", axis)); }
927 
931  std::string getUnit(int axis) const { return getC(detail::formatAxisAttr("Unit", axis)); }
932 
964  std::vector<double> intersect(std::vector<double> const &a1, std::vector<double> const &a2,
965  std::vector<double> const &b1, std::vector<double> const &b2) const;
966 
984  std::vector<int> ret(other.getNIn());
985  astMatchAxes(getRawPtr(), other.getRawPtr(), ret.data());
986  assertOK();
987  return ret;
988  }
989 
1021  CmpFrame under(Frame const &next) const;
1022 
1043  PointD norm(PointD value) const {
1044  astNorm(getRawPtr(), value.data());
1045  assertOK();
1046  detail::astBadToNan(value);
1047  return value;
1048  }
1049 
1073  PointD offset(PointD point1, PointD point2, double offset) const {
1074  assertPointLength(point1, "point1");
1075  assertPointLength(point2, "point2");
1076  PointD ret(getNIn());
1077  astOffset(getRawPtr(), point1.data(), point2.data(), offset, ret.data());
1078  assertOK();
1079  detail::astBadToNan(ret);
1080  return ret;
1081  }
1082 
1115  DirectionPoint offset2(PointD const &point1, double angle, double offset) const {
1116  detail::assertEqual(getNIn(), "naxes", 2, " cannot call offset2");
1117  assertPointLength(point1, "point1");
1118  PointD point2(getNIn());
1119  double offsetAngle = astOffset2(getRawPtr(), point1.data(), angle, offset, point2.data());
1120  assertOK();
1121  detail::astBadToNan(point2);
1122  return DirectionPoint(detail::safeDouble(offsetAngle), point2);
1123  }
1124 
1137  detail::assertEqual(perm.size(), "perm.size()", static_cast<std::size_t>(getNAxes()), "naxes");
1138  astPermAxes(getRawPtr(), perm.data());
1139  assertOK();
1140  }
1141 
1161  FrameMapping pickAxes(std::vector<int> const &axes) const;
1162 
1194  ResolvedPoint resolve(std::vector<double> const &point1, std::vector<double> const &point2,
1195  std::vector<double> const &point3) const;
1196 
1201  void setAlignSystem(std::string const &system) { setC("AlignSystem", system); }
1202 
1206  void setBottom(int axis, double bottom) { setD(detail::formatAxisAttr("Bottom", axis), bottom); }
1207 
1211  void setDigits(int digits) { setI("Digits", digits); }
1212 
1216  void setDigits(int axis, int digits) { setD(detail::formatAxisAttr("Digits", axis), digits); }
1217 
1221  void setDirection(bool direction, int axis) {
1222  setB(detail::formatAxisAttr("Direction", axis), direction);
1223  }
1224 
1228  virtual void setDomain(std::string const &domain) { setC("Domain", domain); }
1229 
1233  void setDut1(double dut1) { setD("Dut1", dut1); }
1234 
1238  void setEpoch(double epoch) { setD("Epoch", epoch); }
1239 
1243  void setEpoch(std::string const &epoch) { setC("Epoch", epoch); }
1244 
1248  void setFormat(int axis, std::string const &format) {
1249  setC(detail::formatAxisAttr("Format", axis), format);
1250  }
1251 
1255  void setLabel(int axis, std::string const &label) { setC(detail::formatAxisAttr("Label", axis), label); }
1256 
1260  void setMatchEnd(bool match) { setB("MatchEnd", match); }
1261 
1266  void setMaxAxes(int maxAxes) { setI("MaxAxes", maxAxes); }
1267 
1272  void setMinAxes(int minAxes) { setI("MinAxes", minAxes); }
1273 
1277  void setObsAlt(double alt) { setD("ObsAlt", alt); }
1278 
1282  void setObsLat(std::string const &lat) { setC("ObsLat", lat); }
1283 
1287  void setObsLon(std::string const &lon) { setC("ObsLon", lon); }
1288 
1293  void setActiveUnit(bool enable) {
1294  astSetActiveUnit(getRawPtr(), enable);
1295  assertOK();
1296  }
1297 
1301  void setPermute(bool permute) { setB("Permute", permute); }
1302 
1306  void setPreserveAxes(bool preserve) { setB("PreserveAxes", preserve); }
1307 
1311  void setSymbol(int axis, std::string const &symbol) {
1312  setC(detail::formatAxisAttr("Symbol", axis), symbol);
1313  }
1314 
1319  void setSystem(std::string const &system) { setC("System", system); }
1320 
1324  void setTitle(std::string const &title) { setC("Title", title); }
1325 
1329  void setTop(int axis, double top) { setD(detail::formatAxisAttr("Top", axis), top); }
1330 
1334  void setUnit(int axis, std::string const &unit) { setC(detail::formatAxisAttr("Unit", axis), unit); }
1335 
1500  NReadValue unformat(int axis, std::string const &str) const {
1501  double value;
1502  int nread = astUnformat(getRawPtr(), axis, str.c_str(), &value);
1503  assertOK();
1504  return NReadValue(nread, detail::safeDouble(value));
1505  }
1506 
1507 protected:
1517  explicit Frame(AstFrame *rawPtr) : Mapping(reinterpret_cast<AstMapping *>(rawPtr)) {
1518  if (!astIsAFrame(getRawPtr())) {
1520  os << "This is a " << getClassName() << ", which is not a Frame";
1521  throw std::invalid_argument(os.str());
1522  }
1523  }
1524 
1525  virtual std::shared_ptr<Object> copyPolymorphic() const override { return copyImpl<Frame, AstFrame>(); }
1526 
1527 private:
1534  template <typename T>
1535  void assertPointLength(T const &p, char const *name) const {
1536  if (static_cast<int>(p.size()) != getNIn()) {
1538  os << "point " << name << " has " << p.size() << " axes, but " << getNIn() << " required";
1539  throw std::invalid_argument(os.str());
1540  }
1541  }
1542 };
1543 
1544 } // namespace ast
1545 
1546 #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:825
void setDigits(int digits)
Set Digits for all axes: number of digits of precision.
Definition: Frame.h:1211
void setObsLat(std::string const &lat)
Set ObsLat: frame title.
Definition: Frame.h:1282
double getTop(int axis) const
Get Top: the highest axis value to display.
Definition: Frame.h:926
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
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:1243
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:1306
double distance(PointD const &point1, PointD const &point2) const
Find the distance between two points whose Frame coordinates are given.
Definition: Frame.h:472
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:1266
table::Key< int > b
void setDut1(double dut1)
Set Dut1: difference between the UT1 and UTC timescale (sec)
Definition: Frame.h:1233
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
table::Key< int > a
void setSymbol(int axis, std::string const &symbol)
Set Symbol(axis) for one axis: axis symbol.
Definition: Frame.h:1311
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:805
table::Key< double > angle
tuple options
Definition: lsstimport.py:47
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:1248
void permAxes(std::vector< int > perm)
Permute the order in which a Frame&#39;s axes occur.
Definition: Frame.h:1136
void astBadToNan(std::vector< double > &p)
Replace AST__BAD with a quiet NaN in a vector.
Definition: utils.h:58
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:1329
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:856
std::string getFormat(int axis) const
Get Format for one axis: format specification for axis values.
Definition: Frame.h:840
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:1238
std::string getSystem() const
Get System: coordinate system used to describe positions within the domain.
Definition: Frame.h:916
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:874
std::shared_ptr< Mapping > mapping
Mapping.
Definition: Frame.h:102
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:129
std::string getUnit(int axis) const
Get Unit(axis) for one axis: physical units for formatted axis values.
Definition: Frame.h:931
void setBottom(int axis, double bottom)
Set Bottom: the lowest axis value to display.
Definition: Frame.h:1206
T str(T... args)
Frame(AstFrame *rawPtr)
Construct a Frame from a pointer to a raw AstFrame.
Definition: Frame.h:1517
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:1255
std::shared_ptr< Frame > frame
Frame.
Definition: Frame.h:101
void setObsAlt(double alt)
Set ObsAlt: Geodetic altitude of observer (m).
Definition: Frame.h:1277
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:1073
int getDigits() const
Get Digits: the default used if no specific value specified for an axis.
Definition: Frame.h:810
virtual void setDomain(std::string const &domain)
Set Domain: coordinate system domain.
Definition: Frame.h:1228
bool getDirection(int axis) const
Get Direction for one axis: display axis in conventional direction?
Definition: Frame.h:820
std::string getObsLon() const
Get ObsLon: Geodetic longitude of observer.
Definition: Frame.h:895
int getDigits(int axis) const
Get Digits for one axis.
Definition: Frame.h:815
double getObsAlt() const
Get ObsAlt: Geodetic altitude of observer (m).
Definition: Frame.h:885
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:862
std::string getInternalUnit(int axis) const
Get InternalUnit(axis) read-only attribute for one axis: physical units for unformated axis values...
Definition: Frame.h:846
double getEpoch() const
Get Epoch: Epoch of observation.
Definition: Frame.h:835
std::string getTitle() const
Get Title: frame title.
Definition: Frame.h:921
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:905
std::string getLabel(int axis) const
Get Label(axis) for one axis: axis label.
Definition: Frame.h:851
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:830
void setObsLon(std::string const &lon)
Set ObsLon: Geodetic longitude of observer.
Definition: Frame.h:1287
std::string formatAxisAttr(std::string const &name, int axis)
Format an axis-specific attribute by appending the axis index.
Definition: utils.h:78
void setMatchEnd(bool match)
Set MatchEnd: match trailing axes?
Definition: Frame.h:1260
bool getPermute() const
Get Permute: allow axis permutation when used as a template?
Definition: Frame.h:900
double safeDouble(double val)
Return a double value after checking status and replacing AST__BAD with nan
Definition: utils.h:99
ItemVariant const * other
Definition: Schema.cc:55
std::string getAlignSystem() const
Get AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition: Frame.h:800
void assertEqual(T1 val1, std::string const &descr1, T2 val2, std::string const &descr2)
Definition: utils.h:47
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:1334
void setSystem(std::string const &system)
Set System: coordinate system used to describe positions within the domain.
Definition: Frame.h:1319
std::string getObsLat() const
Get ObsLat: Geodetic latitude of observer.
Definition: Frame.h:890
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:1221
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:1115
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:1500
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:780
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:910
void setMinAxes(int minAxes)
Get MinAxes: the minimum number of axes a frame found by findFrame may have.
Definition: Frame.h:1272
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:1293
void setAlignSystem(std::string const &system)
Set AlignSystem: the coordinate system used by convert and findFrame to align Frames.
Definition: Frame.h:1201
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:790
std::ostream * os
Definition: Schema.cc:737
virtual std::shared_ptr< Object > copyPolymorphic() const override
Return a deep copy of this object.
Definition: Frame.h:1525
void setPermute(bool permute)
Set Permute: allow axis permutation when used as a template?
Definition: Frame.h:1301
void setDigits(int axis, int digits)
Set Digits for one axis: number of digits of precision.
Definition: Frame.h:1216
void setTitle(std::string const &title)
Set Title: frame title.
Definition: Frame.h:1324
std::vector< int > matchAxes(Frame const &other) const
Look for corresponding axes between this frame and another.
Definition: Frame.h:983
PointD norm(PointD value) const
Normalise a set of Frame coordinate values which might be unsuitable for display (e.g.
Definition: Frame.h:1043
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:880
int getMinAxes() const
Get MinAxes: the maximum axes a frame found by findFrame may have.
Definition: Frame.h:868