SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
Loading...
Searching...
No Matches
Interferometric applications

Various tools supporting interferometric applications, such as u, v, w projections of observing stations and delay calculations along a line of sight. More...

Classes

class  supernovas::Interferometric
 u, v, w projections of an interferometric station along a line of sight. More...

Functions

Interferometric supernovas::AstrometricPosition::to_interferometric (const Equatorial &phase_center, const Coordinate &distance=Coordinate::at_Gpc(), const Velocity &relative_motion=Velocity::stationary()) const
 Returns u,v,w coordinates for a space-based interferometer station represented by this astrometric position.
Interferometric supernovas::Observer::to_interferometric (const Apparent &phase_center, enum novas_reference_system system=NOVAS_ICRS) const
 Returns projected u,v,w coordinates for an interferometer station.

Detailed Description

Various tools supporting interferometric applications, such as u, v, w projections of observing stations and delay calculations along a line of sight.

Function Documentation

◆ to_interferometric() [1/2]

Interferometric supernovas::AstrometricPosition::to_interferometric ( const Equatorial & phase_center,
const Coordinate & distance = Coordinate::at_Gpc(),
const Velocity & relative_motion = Velocity::stationary() ) const

Returns u,v,w coordinates for a space-based interferometer station represented by this astrometric position.

That is, it returns the u,v,w coordinates of this astrometric place (of a station), measured relative to the reference point (array reference location), for a given apparent reference line-of-sight on the sky (source) at the time as when the station location is defined. The u and v coordinates are the orthogonal projections of the site, relative to the array reference, in the direction of the local East and North respectively, as seen from the source; while w is the distance from the array reference along the line of sight.

You could also use Observer::to_interferometric() instead. However, using relative astrometric positions can overcome numerical precision issues for interferometers located far from the geocenter or the SSB, such as at L2. Specifically, Observer::to_interferometric() uses absolute geocentric and SSB station positions, and interferometric baselines are obtained from differencing these with the reference position used in defining the phase center direction. In contrast, astrometric positions are always defined as relative positions from the array reference, and will be more accurate, in general.

// The momentary position of a station relative to the array reference
AstrometricPosition station = ...;
// the location of the phase center
Equatorial eq = ...;
// u,v,w coordinates for a source observed at the time the station position is defined
// the geometric delay of the station, relative to the array reference
Interval& dt = uvw.geometric_delay();
Equatorial coordinates (RA, Dec / α, δ), representing the direction on the sky, for a particular type...
Definition supernovas.h:995
u, v, w projections of an interferometric station along a line of sight.
Definition supernovas.h:839
A signed time interval between two instants of time, in the astronomical timescale of choice.
Definition supernovas.h:502
Interferometric to_interferometric(const Equatorial &phase_center, const Coordinate &distance=Coordinate::at_Gpc(), const Velocity &relative_motion=Velocity::stationary()) const
Returns u,v,w coordinates for a space-based interferometer station represented by this astrometric po...
Definition AstrometricPosition.cpp:235

For ground-based interferometric application, see GeodeticObserver::to_interferometric() instead.

Parameters
phase_centerApparent equatorial coordinates of the interferometric phase center, as seen from the array reference position, at the time at which this position is defined.
distance(optional) Apparent distance to phase center, from the array reference, at the time of observation (default: 1 Gpc).
relative_motion(optional station's velocity vector relative to the reference position (default: stationary).
Returns
interferometric uvw projection this astrometric place viewed from the direction of the source at the time for which this position was defined. The u and v directions are aligned with the local East and Noth in the coordinate system in which the phase center was specified.
Since
1.6
See also
Observer::to_interferometric()

References supernovas::Vector::_array(), supernovas::Unit::AU, supernovas::Unit::day, supernovas::Position::distance(), supernovas::Equinox::from_system_type(), supernovas::Validating::is_valid(), novas_uvw(), obs_time(), supernovas::Vector::scaled(), supernovas::Equatorial::to_system(), and supernovas::Spherical::xyz().

◆ to_interferometric() [2/2]

Interferometric supernovas::Observer::to_interferometric ( const Apparent & phase_center,
enum novas_reference_system system = NOVAS_ICRS ) const

Returns projected u,v,w coordinates for an interferometer station.

That is, it returns u,v,w coordinates, for this observer (station), for a given apparent line-of-sight on the sky, at the time of observation and relative to the location it was calculated for. The u and v coordinates are the orthogonal projections of the site, relative to the array reference, in the directions of local East and North respectively (w.r.t. the coordinate system of choice), as seen from the source; while w is the distance from the array reference along the line of sight.

The supplied phase center defines the time of observation and the array reference position in its observing frame. For example:

// Define the array reference location (e.g. on Earth) as an observer place
auto reference = Observer::on_earth(...);
// Astrometric time of observation
Time time = ...;
// The observing frame of the array reference
Frame frame = reference.frame_at(time);
// E.g. an astronomical source at the interferometric phase center
Source source = ...;
// the apparent location of the phase center for an observer at the geocenter
Apparent app = source.apparent_in(frame);
// Define an interferometer station (e.g. on Earth)
auto station = Observer::on_earth(...);
// u,v,w coordinates of the station relative to the array reference (w.r.t. ICRS)
Interferometric uwv = station.to_interferometric(app);
// the geometric delay of the stanovas_site_uvwtion, relative to the array reference
Interval delay = uvw.geometric_delay()
Apparent position on sky as seen by an observer at a specific time of observation.
Definition supernovas.h:2485
An observing frame, defined by an observer location and precise time of observation.
Definition supernovas.h:1901
static GeodeticObserver on_earth(const Site &site, const EOP &eop)
Returns a new observer located at a fixed observing site.
Definition Observer.cpp:307
An abstract superclass for an astronomical source or target of observation.
Definition supernovas.h:1964
Precise astronomical time specification, supporting all relevant astronomical timescales (UT1,...
Definition supernovas.h:1750

NOTES:

  • This method does not take atmospheric refraction into account. Refraction is usually included as a separate atmospheric delay term on top of the geometric delays returned here, and the various other delay terms that account for optics, cabling, and digital electronics etc.
  • This method supports down to nanometer precision for sites on Earth or in Low Earth Orbit (LEO), and sub micron (<μm) precision even at the distance of the Moon. However, because its calculations are based on geocentric positions, the precision degrades with increasing distance from Earth, and may not be suitable for interferometers far from Earth. When precision is a concern, you might use AstrometricPosition::to_interferometric() instead, with positions (and velocities) of the stations defined relative to the array center – enabling higher precision projections than this method for interferometers far from Earth.
Parameters
phase_centerApparent place on sky from the array reference place, at the time of observation.
system(optional) Coordinate reference system type in which to return result. (default: ICRS). Specifically, the u and v directions of the returned projection will be aligned to the local East and North directions of the specified coordinate reference system type at the time of observation.
Returns
interferometric uvw _projection of this site in the ICRS, relative to the array reference, viewed from the direction of the source at the specified time. The _u and v directions are aligned with the local East and North w.r.t. the coordinate system of choice (default: ICRS).
Since
1.6
See also
AstrometricPosition::to_interferometric()

References supernovas::Vector::_array(), supernovas::Frame::accuracy(), supernovas::Unit::AU, supernovas::Unit::day, supernovas::Apparent::distance(), supernovas::Apparent::equatorial(), supernovas::Apparent::frame(), supernovas::Equinox::from_system_type(), gcrs_position_at(), gcrs_velocity_at(), supernovas::Validating::is_valid(), NOVAS_ICRS, NOVAS_REFERENCE_SYSTEMS, novas_uvw(), supernovas::Frame::observer(), supernovas::Geometric::position(), supernovas::Vector::scaled(), supernovas::Frame::time(), to_interferometric(), supernovas::Equatorial::to_system(), supernovas::Interferometric::undefined(), supernovas::Geometric::velocity(), and supernovas::Spherical::xyz().

Referenced by to_interferometric().