siderust-cpp
Header-only C++ wrapper for siderust
Loading...
Searching...
No Matches
coordinates_examples.cpp
Go to the documentation of this file.
1
11#include <siderust/siderust.hpp>
12
13#include <cstdio>
14
15using namespace siderust;
16using namespace qtty::literals;
17
19 std::printf("1) Geodetic -> ECEF cartesian\n");
20
21 Geodetic obs(-17.8890, 28.7610, 2396.0);
22 auto ecef = obs.to_cartesian();
23 auto ecef_km = obs.to_cartesian<qtty::Kilometer>();
24
25 std::printf(" Geodetic lon=%.4f deg lat=%.4f deg h=%.1f m\n",
26 obs.lon.value(), obs.lat.value(), obs.height.value());
27 std::printf(" ECEF x=%.2f m y=%.2f m z=%.2f m\n\n",
28 ecef.x().value(), ecef.y().value(), ecef.z().value());
29 std::printf(" ECEF x=%.2f km y=%.2f km z=%.2f km\n\n",
30 ecef_km.x().value(), ecef_km.y().value(), ecef_km.z().value());
31}
32
34 std::printf("2) Spherical direction frame conversions\n");
35
36 spherical::direction::ICRS vega_icrs(279.23473, 38.78369);
37 auto jd = JulianDate::from_utc({2026, 7, 15, 22, 0, 0});
38
39 auto ecl = vega_icrs.to<frames::EclipticMeanJ2000>(jd);
40 auto eq_mod = vega_icrs.to<frames::EquatorialMeanOfDate>(jd);
41 auto hor = vega_icrs.to_horizontal(jd, ROQUE_DE_LOS_MUCHACHOS);
42
43 std::printf(" ICRS RA=%.5f Dec=%.5f\n", vega_icrs.ra().value(), vega_icrs.dec().value());
44 std::printf(" Ecliptic lon=%.5f lat=%.5f\n", ecl.lon().value(), ecl.lat().value());
45 std::printf(" Equatorial(MOD) RA=%.5f Dec=%.5f\n", eq_mod.ra().value(), eq_mod.dec().value());
46 std::printf(" Horizontal az=%.5f alt=%.5f\n\n", hor.az().value(), hor.alt().value());
47}
48
50 std::printf("3) Spherical position + extracting direction\n");
51
53 120.0_deg, -25.0_deg, 2.0e17_m);
54 auto dir = target.direction();
55
56 std::printf(" Position RA=%.2f Dec=%.2f dist=%.3e m\n",
57 target.ra().value(), target.dec().value(), target.distance().value());
58 std::printf(" Direction-only RA=%.2f Dec=%.2f\n\n",
59 dir.ra().value(), dir.dec().value());
60}
61
63 std::printf("4) Cartesian coordinate creation + unit conversion\n");
64
65 cartesian::Direction<frames::ICRS> axis_x(1.0, 0.0, 0.0);
67
68 auto x_km = sample_helio_au.x().to<qtty::Kilometer>();
69 auto y_km = sample_helio_au.y().to<qtty::Kilometer>();
70
71 std::printf(" Direction<ICRS> x=%.1f y=%.1f z=%.1f\n", axis_x.x, axis_x.y, axis_x.z);
72 std::printf(" Position<Heliocentric,Ecl,AU> x=%.3f AU y=%.3f AU z=%.3f AU\n",
73 sample_helio_au.x().value(), sample_helio_au.y().value(), sample_helio_au.z().value());
74 std::printf(" Same position in km x=%.2f y=%.2f\n\n", x_km.value(), y_km.value());
75}
76
78 std::printf("5) Typed ephemeris coordinates\n");
79
80 auto jd = JulianDate::J2000();
81 auto earth = ephemeris::earth_heliocentric(jd); // cartesian::position::EclipticMeanJ2000<qtty::AstronomicalUnit>
82 auto moon = ephemeris::moon_geocentric(jd); // cartesian::position::MoonGeocentric<qtty::Kilometer>
83
84 std::printf(" Earth heliocentric (AU) x=%.8f y=%.8f z=%.8f\n",
85 earth.x().value(), earth.y().value(), earth.z().value());
86 std::printf(" Moon geocentric (km) x=%.3f y=%.3f z=%.3f\n\n",
87 moon.x().value(), moon.y().value(), moon.z().value());
88}
89
90int main() {
91 std::printf("=== Coordinate Creation & Conversion Examples ===\n\n");
92
98
99 std::printf("Done.\n");
100 return 0;
101}
static void geodetic_and_ecef_example()
static void spherical_direction_example()
static void cartesian_and_units_example()
static void ephemeris_typed_example()
static void spherical_position_example()
int main()
cartesian::position::EclipticMeanJ2000< qtty::AstronomicalUnit > earth_heliocentric(const JulianDate &jd)
Earth's heliocentric position (EclipticMeanJ2000, AU) via VSOP87.
Definition ephemeris.hpp:46
cartesian::position::MoonGeocentric< qtty::Kilometer > moon_geocentric(const JulianDate &jd)
Moon's geocentric position (EclipticMeanJ2000, km) via ELP2000.
Definition ephemeris.hpp:56
const Geodetic ROQUE_DE_LOS_MUCHACHOS
Roque de los Muchachos Observatory (La Palma, Spain).
Planet earth()
Definition bodies.hpp:142
Umbrella header for the siderust C++ wrapper library.
Geodetic position (WGS84 ellipsoid).
Definition geodetic.hpp:28
qtty::Degree lon
Longitude (east positive).
Definition geodetic.hpp:29
cartesian::Position< centers::Geocentric, frames::ECEF, U > to_cartesian() const
Convert geodetic (WGS84/ECEF) to cartesian position.
qtty::Degree lat
Latitude (north positive).
Definition geodetic.hpp:30
qtty::Meter height
Height above ellipsoid.
Definition geodetic.hpp:31
A unit-vector direction in Cartesian form, compile-time frame-tagged.
Definition cartesian.hpp:27
double x
X component (unitless).
Definition cartesian.hpp:30
double y
Y component (unitless).
Definition cartesian.hpp:31
double z
Z component (unitless).
Definition cartesian.hpp:32
A 3D Cartesian position, compile-time tagged by center, frame, unit.
Definition cartesian.hpp:53
Mean ecliptic & equinox of J2000.0.
Definition frames.hpp:51
Mean equatorial of date (precessed, no nutation).
Definition frames.hpp:59
A direction on the celestial sphere, compile-time tagged by frame.
Definition spherical.hpp:36
std::enable_if_t< frames::has_horizontal_transform_v< F_ >, Direction< frames::Horizontal > > to_horizontal(const JulianDate &jd, const Geodetic &observer) const
Transform to the horizontal (alt-az) frame.
auto to(const JulianDate &jd) const -> decltype(this->template to_frame< Target >(jd))
Shorthand: .to<Target>(jd) (calls to_frame).
qtty::Degree dec() const
Definition spherical.hpp:69
qtty::Degree ra() const
Definition spherical.hpp:66
A spherical position (direction + distance), compile-time tagged.
Direction< F > direction() const
Extract the direction component.
qtty::Degree ra() const
qtty::Degree dec() const