SuperNOVAS C API v1.6
High-precision C/C++ astrometry library
Loading...
Searching...
No Matches
SuperNOVAS C99 vs. astropy

Nowadays astropy is widely used in the astronomy community and it is known for its simplicity and elegance, but it is rather slow (putting it mildly). In contrast, C is known to be fast, but it also has a bad reputation for being tedious and 'ugly'. However, below is a side-by-side comparison of equivalent program snippets for calculating CIRS apparent positions in astropy vs. SuperNOVAS for Antares for a given date and observer location:

astropysupernovas (C99)

from astropy import units as u
from astropy.coordinates import SkyCoord,
EarthLocation, Longitude, Latitude,
CIRS
# Define ICRS coordinates
source = SkyCoord(
'16h 29m 24.45970s', '−26d 25m 55.2094s',
d = u.AU / 5.89 * u.mas,
pmra = -12.11 * u.mas / u.yr,
pmdec = -23.30 * u.mas / u.yr,
rv = -3.4 * u.km / u.s)
# Observer location
loc = EarthLocation.from_geodetic(
Longitude(50.7374),
Latitude(7.0982),
height=60.0)
# Set time of observation
time = astropy.time.Time(
"2025-02-27T19:57:00.728+0200"
scale='tai')
# Observer frame & system
frame = CIRS(obstime=time, location=loc)
# apparent coordinates
app = source.transform_to(frame);

#include <novas.h>
cat_entry star;
object source;
sky_pos app;
// IERS Earth Orientation Parameters...
int leap_seconds = 37;
double dt1 = 0.06256; // [s]
double dx = 103.4; // [mas]
double dy = 396.2; // [mas]
// Define ICRS coordinates
make_cat_entry("Antares", "HIP", 80763,
novas_hms_hours("16h 29m 24.45970s"),
novas_dms_degrees("−26d 25m 55.2094s"),
-12.11, -23.30, 5.89, -3.4, &star);
make_cat_object(&star, &source);
// Observer location
make_gps_observer(50.7374, 7.0982, 60.0,
&loc);
// Set time of observation
"2025-02-27T19:57:00.728+0200",
leap_seconds, dut1, &time);
// Observer frame
&loc, &time, dx, dy, &frame);
// apparent coordinates in system
novas_sky_pos(&source, &frame, NOVAS_CIRS, &app);
struct novas_sky_pos sky_pos
Celestial object's place on the sky; contains the output from place().
int novas_sky_pos(const object *restrict object, const novas_frame *restrict frame, enum novas_reference_system sys, sky_pos *restrict out)
Calculates an apparent location on sky for the source.
Definition frames.c:824
@ NOVAS_CIRS
Celestial Intermediate Reference System: dynamical system of the true equator, with its origin at the...
Definition novas.h:872
int novas_make_frame(enum novas_accuracy accuracy, const observer *obs, const novas_timespec *time, double xp, double yp, novas_frame *frame)
Sets up a observing frame for a specific observer location, time of observation, and accuracy require...
Definition frames.c:509
@ NOVAS_FULL_ACCURACY
Use full precision calculations to micro-arcsecond accuracy.
Definition novas.h:969
int make_gps_observer(double latitude, double longitude, double height, observer *obs)
Initializes an observer data structure for a ground-based observer with the specified GPS / WGS84 loc...
Definition observer.c:278
struct novas_observer observer
Observer location.
struct novas_cat_entry cat_entry
Basic astrometric data for any sidereal object located outside the solar system.
short make_cat_entry(const char *restrict name, const char *restrict catalog, long cat_num, double ra, double dec, double pm_ra, double pm_dec, double parallax, double rad_vel, cat_entry *source)
Fully populates the data structure for a catalog source, such as a star, with all parameters provided...
Definition target.c:331
int make_cat_object(const cat_entry *star, object *source)
Populates and object data structure with the data for a catalog source.
Definition target.c:487
int novas_set_str_time(enum novas_timescale timescale, const char *restrict str, int leap, double dut1, novas_timespec *restrict time)
Sets astronomical time in a specific timescale using a string specification.
Definition timescale.c:501
@ NOVAS_TAI
Innternational Atomic Time (TAI).
Definition novas.h:1780
double novas_hms_hours(const char *restrict hms)
Returns the decimal hours for a HMS string specification.
Definition parse.c:247
double novas_dms_degrees(const char *restrict dms)
Returns the decimal degrees for a DMS string specification.
Definition parse.c:457
SuperNOVAS types, definitions, and function prototypes.
A set of parameters that uniquely define the place and time of observation.
Definition novas.h:1900
A structure, which defines a precise instant of time that can be extpressed in any of the astronomica...
Definition novas.h:1806

Yes, SuperNOVAS is a bit more verbose, but not painfully so. While astropy will automatically fetch IERS data (leap seconds, and Earth-orientation parameters), in SuperNOVAS you will have to set these explicitly. And, of course, the above comparison ignores error checking also. In Python, you can simply surround the above code in a try block, and then catch errors using except. In C, there is no catch-all solution like that. Instead, you will have to check the return values for each line, and decide if you need to bail early, e.g.:

if(make_cat_object(&star, &source) != 0) {
// oops, the above failed, bail if we must...
return -1;
}

Regardless of the slight omissions, the difference is not night and day. If you can handle astropy, chances are you can handle SuperNOVAS too. And the reward for dealing with slightly 'uglier' code, is orders of magnitude gain in the speed of execution. It is a trade-off that worth considering, at the least.


Copyright (C) 2025 Attila Kovács