Reference

Data types

Observatory

AstroLib.jl defines a new Observatory type. This can be used to define a new object holding information about an observing site. It is a [composite type] whose fields are

  • name (String type): the name of the site
  • latitude (Float64 type): North-ward latitude of the site in degrees
  • longitude (Float64 type): East-ward longitude of the site in degrees
  • altitude (Float64 type): altitude of the site in meters
  • tz (Float64 type): the number of hours of offset from UTC

The type constructor Observatory can be used to create a new Observatory object. Its syntax is

Observatory(name, lat, long, alt, tz)

name should be a string; lat, long, and tz should be anything that can be converted to a floating number with ten function; alt should be a real number.

A predefined list of some observing sites is provided with AstroLib.observatories constant. It is a dictionary whose keys are the abbreviated names of the observatories. For example, you can access information of the European Southern Observatory with

julia> obs = AstroLib.observatories["eso"]
Observatory: European Southern Observatory
latitude:    -29.256666666666668°N
longitude:   -70.73°E
altitude:    2347.0 m
time zone:   UTC-4

julia> obs.longitude
-70.73

You can list all keys of the dictionary with

keys(AstroLib.observatories)

Feel free to contribute new sites or adjust information of already present ones.

Planet

The package provides Planet type to hold information about Solar System planets. Its fields are

  • Designation:
    • name: the name
  • Physical characteristics:
    • radius: mean radius in meters
    • eqradius: equatorial radius in meters
    • polradius: polar radius in meters
    • mass: mass in kilogram
  • Orbital characteristics (epoch J2000):
    • ecc: eccentricity of the orbit
    • axis: semi-major axis of the orbit in meters
    • period: sidereal orbital period in seconds

The constructor has this syntax:

Planet(name, radius, eqradius, polradius, mass, ecc, axis, period)

The list of Solar System planets, from Mercury to Pluto, is available with AstroLib.planets dictionary. The keys of this dictionary are the lowercase names of the planets. For example:

julia> AstroLib.planets["mercury"]
Planet:            Mercury
mean radius:       2.4397e6 m
equatorial radius: 2.4397e6 m
polar radius:      2.4397e6 m
mass:              3.3011e23 kg
eccentricity:      0.20563069
semi-major axis:   5.790905e10 m
period:            5.790905e10 s

julia> AstroLib.planets["mars"].eqradius
3.3962e6

julia> AstroLib.planets["saturn"].mass
5.6834e26

Functions organized by category

Coordinates and positions

adstring(), aitoff(), altaz2hadec(), baryvel(), bprecess(), co_aberration(), co_nutate(), co_refract(), eci2geo(), eq2hor(), eqpole(), euler(), gcirc(), geo2eci(), geo2geodetic(), geo2mag(), geodetic2geo(), hadec2altaz(), helio_rv(), helio(), hor2eq(), jprecess(), mag2geo(), mean_obliquity(), planet_coords(), polrec(), posang(), precess(), precess_cd(), precess_xyz(), premat(), radec(), recpol(), true_obliquity(), zenpos()

Time and date

ct2lst(), daycnv(), get_date(), get_juldate(), helio_jd(), jdcnv(), juldate(), month_cnv(), nutate(), ydn2md(), ymd2dn()

Moon and sun

moonpos(), mphase(), sunpos(), xyz()

Utilities

airtovac(), calz_unred(), deredd(), flux2mag(), gal_uvw(), imf(), ismeuv(), kepler_solver(), lsf_rotate(), mag2flux(), paczynski(), planck_freq(), planck_wave(), rad2sec(), rhotheta(), sec2rad(), sixty(), sphdist(), ten(), tic_one(), ticpos(), tics(), trueanom(), uvbybeta(), vactoair()

Miscellaneous (non-astronomy) functions

ordinal()

Types and functions organized alphabetically

AstroLib.POLELATLONGConstant

List of locations of North Magnetic Pole since 1590.

This is provided by World Magnetic Model (https://www.ngdc.noaa.gov/geomag/data/poles/NP.xy).

source
AstroLib.planetsConstant

List of planets of the Solar System, from Mercury to Pluto. The elements of the list have Planet type.

Reference for most quantities is the Planetary Fact Sheet: http://nssdc.gsfc.nasa.gov/planetary/factsheet/index.html and the Keplerian Elements for Approximate Positions of the Major Planets: https://ssd.jpl.nasa.gov/txt/pelemt1.txt

source
AstroLib.ObservatoryType

Type holding information about an observing site. Its fields are:

  • name: the name of the site
  • latitude: North-ward latitude of the site in degrees
  • longitude: East-ward longitude of the site in degrees
  • altitude: altitude of the site in meters
  • tz: the number of hours of offset from UTC
source
AstroLib.PlanetType

Type holding information about a planet. Its fields are:

Designation:

  • name: the name

Physical characteristics:

  • radius: mean radius in meters
  • eqradius: equatorial radius in meters
  • polradius: polar radius in meters
  • mass: mass in kilogram

Orbital characteristics (epoch J2000):

  • ecc: eccentricity of the orbit
  • axis: semi-major axis of the orbit in meters
  • period: sidereal orbital period in seconds

Position characteristics (epoch J2000):

  • inc: inclination in degrees
  • asc_long: longitude of the ascending node in degrees
  • per_long: longitude of perihelion in degrees
  • mean_long: mean longitude in degrees
source
AstroLib.adstringMethod
adstring(ra::Real, dec::Real[, precision::Int=2, truncate::Bool=true]) -> string
adstring([ra, dec]) -> string
adstring(dec) -> string
adstring([ra], [dec]) -> ["string1", "string2", ...]

Purpose

Returns right ascension and declination as string(s) in sexagesimal format.

Explanation

Takes right ascension and declination expressed in decimal format, converts them to sexagesimal and return a formatted string. The precision of right ascension and declination can be specified.

Arguments

Arguments of this function are:

  • ra: right ascension in decimal degrees. It is converted to hours before printing.
  • dec: declination in decimal degrees.

The function can be called in different ways:

  • Two numeric arguments: first is ra, the second is dec.
  • An iterable (array, tuple) of two elements: (ra, dec).
  • One numeric argument: it is assumed only dec is provided.

Optional keywords affecting the output format are always available:

  • precision (optional integer keyword): specifies the number of digits of declination seconds. The number of digits for right ascension seconds is always assumed to be one more precision. If the function is called with only dec as input, precision default to 1, in any other case defaults to 0.
  • truncate (optional boolean keyword): if true, then the last displayed digit in the output is truncated in precision rather than rounded. This option is useful if adstring is used to form an official IAU name (see http://vizier.u-strasbg.fr/Dic/iau-spec.htx) with coordinate specification.

Output

The function returns one string. The format of strings can be specified with precision and truncate keywords, see above.

Example

julia> using AstroLib

julia> adstring(30.4, -1.23, truncate=true)
" 02 01 35.9  -01 13 48"

julia> adstring.([30.4, -15.63], [-1.23, 48.41], precision=1)
2-element Array{String,1}:
 " 02 01 36.00  -01 13 48.0"
 " 22 57 28.80  +48 24 36.0"
source
AstroLib.airtovacMethod
airtovac(wave_air) -> wave_vacuum

Purpose

Converts air wavelengths to vacuum wavelengths.

Explanation

Wavelengths are corrected for the index of refraction of air under standard conditions. Wavelength values below $2000 Å$ will not be altered, take care within $[1 Å, 2000 Å]$. Uses relation of Ciddor (1996).

Arguments

  • wave_air: the wavelength in air.

Output

Vacuum wavelength in angstroms.

Method

Uses relation of Ciddor (1996), Applied Optics 62, 958.

Example

If the air wavelength is w = 6056.125 (a Krypton line), then airtovac(w) yields a vacuum wavelength of 6057.8019.

julia> using AstroLib

julia> airtovac(6056.125)
6057.801930991426

Notes

vactoair converts vacuum wavelengths to air wavelengths.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.aitoffMethod
aitoff(l, b) -> x, y

Purpose

Convert longitude l and latitude b to (x, y) using an Aitoff projection.

Explanation

This function can be used to create an all-sky map in Galactic coordinates with an equal-area Aitoff projection. Output map coordinates are zero longitude centered.

Arguments

  • l: longitude, scalar or vector, in degrees.
  • b: latitude, number of elements as l, in degrees.

Coordinates can be given also as a 2-tuple (l, b).

Output

2-tuple (x, y).

  • x: x coordinate, same number of elements as l. x is normalized to be in $[-180, 180]$.
  • y: y coordinate, same number of elements as l. y is normalized to be in $[-90, 90]$.

Example

Get $(x ,y)$ Aitoff coordinates of Sirius, whose Galactic coordinates are $(227.23, -8.890)$.

julia> using AstroLib

julia> x, y = aitoff(227.23, -8.890)
(-137.92196683723276, -11.772527357473054)

Notes

See AIPS memo No. 46 (ftp://ftp.aoc.nrao.edu/pub/software/aips/TEXT/PUBL/AIPSMEMO46.PS), page 4, for details of the algorithm. This version of aitoff assumes the projection is centered at b=0 degrees.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.altaz2hadecMethod
altaz2hadec(alt, az, lat) -> ha, dec

Purpose

Convert Horizon (Alt-Az) coordinates to Hour Angle and Declination.

Explanation

Can deal with the NCP singularity. Intended mainly to be used by program hor2eq.

Arguments

Input coordinates may be either a scalar or an array, of the same dimension.

  • alt: local apparent altitude, in degrees, scalar or array.
  • az: the local apparent azimuth, in degrees, scalar or vector, measured east of north!!! If you have measured azimuth west-of-south (like the book Meeus does), convert it to east of north via: az = (az + 180) % 360.
  • lat: the local geodetic latitude, in degrees, scalar or array.

alt and az can be given as a 2-tuple (alt, az).

Output

2-tuple (ha, dec)

  • ha: the local apparent hour angle, in degrees. The hour angle is the time that right ascension of 0 hours crosses the local meridian. It is unambiguously defined.
  • dec: the local apparent declination, in degrees.

The output coordinates are always floating points and have the same type (scalar or array) as the input coordinates.

Example

Arcturus is observed at an apparent altitude of 59d,05m,10s and an azimuth (measured east of north) of 133d,18m,29s while at the latitude of +43.07833 degrees. What are the local hour angle and declination of this object?

julia> using AstroLib

julia> ha, dec = altaz2hadec(ten(59,05,10), ten(133,18,29), 43.07833)
(336.6828582472844, 19.182450965120402)

The widely available XEPHEM code gets:

Hour Angle = 336.683
Declination = 19.1824

Notes

hadec2altaz converts Hour Angle and Declination to Horizon (Alt-Az) coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.baryvelMethod
baryvel(dje, deq) -> dvelh, dvelb

Purpose

Calculates heliocentric and barycentric velocity components of Earth.

Explanation

Baryvel takes into account the Earth-Moon motion, and is useful for radial velocity work to an accuracy of ~1 m/s.

Arguments

  • dje: julian ephemeris date
  • deq (optional): epoch of mean equinox of dvelh and dvelb. If deq is not provided, then it is assumed to be equal to dje.

Output

  • dvelh: heliocentric velocity component. in km/s
  • dvelb: barycentric velocity component. in km/s

Example

Compute the radial velocity of the Earth toward Altair on 15-Feb-1994 using both the original Stumpf algorithm.

julia> using AstroLib

julia> jd = jdcnv(1994, 2, 15, 0)
2.4493985e6

julia> baryvel(jd, 2000)
([-17.0724258266945, -22.81120895274765, -9.889315408506354], [-17.080834081384847, -22.80470807516409, -9.886258269159352])

Notes

The 3-vectors outputs dvelh and dvelb are given in a right-handed coordinate system with the +X axis toward the Vernal Equinox, and +Z axis toward the celestial pole.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.bprecessFunction
bprecess(ra, dec[, epoch]) -> ra1950, dec1950
bprecess(ra, dec, muradec[, parallax=parallax, radvel=radvel]) -> ra1950, dec1950

Purpose

Precess positions from J2000.0 (FK5) to B1950.0 (FK4).

Explanation

Calculates the mean place of a star at B1950.0 on the FK4 system from the mean place at J2000.0 on the FK5 system.

bprecess function has two methods, one for each of the following cases:

  • the proper motion is known and non-zero
  • the proper motion is unknown or known to be exactly zero (i.e. extragalactic radio sources). Better precision can be achieved in this case by inputting the epoch of the original observations.

Arguments

The function has 2 methods. The common mandatory arguments are:

  • ra: input J2000 right ascension, in degrees.
  • dec: input J2000 declination, in degrees.

The two methods have a different third argument (see "Explanation" section for more details). It can be one of the following:

  • muradec: 2-element vector containing the proper motion in seconds of arc per tropical century in right ascension and declination.
  • epoch: scalar giving epoch of original observations.

If none of these two arguments is provided (so bprecess is fed only with right ascension and declination), it is assumed that proper motion is exactly zero and epoch = 2000.

If it is used the method involving muradec argument, the following keywords are available:

  • parallax (optional numerical keyword): stellar parallax, in seconds of arc.
  • radvel (optional numerical keyword): radial velocity in km/s.

Right ascension and declination can be passed as the 2-tuple (ra, dec). You can also pass ra, dec, parallax, and radvel as arrays, all of the same length N. In that case, muradec should be a matrix 2×N.

Output

The 2-tuple of right ascension and declination in 1950, in degrees, of input coordinates is returned. If ra and dec (and other possible optional arguments) are arrays, the 2-tuple of arrays (ra1950, dec1950) of the same length as the input coordinates is returned.

Method

The algorithm is taken from the Explanatory Supplement to the Astronomical Almanac 1992, page 186. See also Aoki et al (1983), A&A, 128, 263. URL: http://adsabs.harvard.edu/abs/1983A%26A...128..263A.

Example

The SAO2000 catalogue gives the J2000 position and proper motion for the star HD

  1. Find the B1950 position.
  • RA(2000) = 13h 42m 12.740s
  • Dec(2000) = 8d 23' 17.69''
  • Mu(RA) = -.0257 s/yr
  • Mu(Dec) = -.090 ''/yr
julia> using AstroLib

julia> muradec = 100*[-15*0.0257, -0.090]; # convert to century proper motion

julia> ra = ten(13, 42, 12.74)*15;

julia> decl = ten(8, 23, 17.69);

julia> adstring(bprecess(ra, decl, muradec), precision=2)
" 13 39 44.526  +08 38 28.63"

Notes

"When transferring individual observations, as opposed to catalog mean place, the safest method is to transform the observations back to the epoch of the observation, on the FK4 system (or in the system that was used to to produce the observed mean place), convert to the FK5 system, and transform to the the epoch and equinox of J2000.0" – from the Explanatory Supplement (1992), p. 180

jprecess performs the precession to J2000 coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.calz_unredFunction
calz_unred(wave, flux, ebv[, r_v]) -> deredden_wave

Purpose

Deredden a galaxy spectrum using the Calzetti et al. (2000) recipe.

Explanation

Calzetti et al. (2000, ApJ 533, 682; http://adsabs.harvard.edu/abs/2000ApJ...533..682C) developed a recipe for dereddening the spectra of galaxies where massive stars dominate the radiation output, valid between $0.12$ to $2.2$ microns. (calz_unred extrapolates between $0.12$ and $0.0912$ microns.)

Arguments

  • wave: wavelength (Angstroms)
  • flux: calibrated flux.
  • ebv: color excess E(B-V). If a negative ebv is supplied, then fluxes will be reddened rather than deredenned. Note that the supplied color excess should be that derived for the stellar continuum, EBV(stars), which is related to the reddening derived from the gas, EBV(gas), via the Balmer decrement by EBV(stars) = 0.44*EBV(gas).
  • r_v (optional): ratio of total to selective extinction, default is 4.05. Calzetti et al. (2000) estimate $r_v = 4.05 \pm 0.80$ from optical-IR observations of 4 starbursts.

Output

Unreddened flux, same units as flux. Flux values will be left unchanged outside valid domain ($0.0912$ - $2.2$ microns).

Example

Estimate how a flat galaxy spectrum (in wavelength) between $1200 Å$ and $3200 Å$ is altered by a reddening of E(B-V) = 0.1.

wave = collect(1200:50:3150);
flux = ones(wave);
flux_new = calz_unred.(wave, flux, -0.1);

Using a plotting tool you can visualize the unreddend flux. For example, with PyPlot.jl

using PyPlot
plot(wave, flux_new)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.co_aberrationFunction
co_aberration(jd, ra, dec[, eps=NaN]) -> d_ra, d_dec

Purpose

Calculate changes to right ascension and declination due to the effect of annual aberration

Explanation

With reference to Meeus, Chapter 23

Arguments

  • jd: julian date, scalar or vector
  • ra: right ascension in degrees, scalar or vector
  • dec: declination in degrees, scalar or vector
  • eps (optional): true obliquity of the ecliptic (in radians). It will be calculated if no argument is specified.

Output

The 2-tuple (d_ra, d_dec):

  • d_ra: correction to right ascension due to aberration, in arc seconds
  • d_dec: correction to declination due to aberration, in arc seconds

Example

Compute the change in RA and Dec of Theta Persei (RA = 2h46m,11.331s, Dec = 49d20',54.5'') due to aberration on 2028 Nov 13.19 TD

julia> using AstroLib

julia> jd = jdcnv(2028,11,13,4, 56)
2.4620887055555554e6

julia> co_aberration(jd,ten(2,46,11.331)*15,ten(49,20,54.54))
(30.04404628365077, 6.699400463119431)

dra = 30.04404628365103'' (≈ 2.003s) ddec = 6.699400463118504''

Notes

Code of this function is based on IDL Astronomy User's Library.

The output dra is not multiplied by cos(dec), so that apparentra = ra + d_ra/3600.

These formula are from Meeus, Chapters 23. Accuracy is much better than 1 arcsecond. The maximum deviation due to annual aberration is 20.49'' and occurs when the Earth's velocity is perpendicular to the direction of the star.

This function calls true_obliquity and sunpos.

source
AstroLib.co_nutateMethod
co_nutate(jd, ra, dec) -> d_ra, d_dec, eps, d_psi, d_eps

Purpose

Calculate changes in RA and Dec due to nutation of the Earth's rotation

Explanation

Calculates necessary changes to ra and dec due to the nutation of the Earth's rotation axis, as described in Meeus, Chap 23. Uses formulae from Astronomical Almanac, 1984, and does the calculations in equatorial rectangular coordinates to avoid singularities at the celestial poles.

Arguments

  • jd: julian date, scalar or vector
  • ra: right ascension in degrees, scalar or vector
  • dec: declination in degrees, scalar or vector

Output

The 5-tuple (d_ra, d_dec, eps, d_psi, d_eps):

  • d_ra: correction to right ascension due to nutation, in degrees
  • d_dec: correction to declination due to nutation, in degrees
  • eps: the true obliquity of the ecliptic
  • d_psi: nutation in the longitude of the ecliptic
  • d_eps: nutation in the obliquity of the ecliptic

Example

Example 23a in Meeus: On 2028 Nov 13.19 TD the mean position of Theta Persei is 2h 46m 11.331s 49d 20' 54.54''. Determine the shift in position due to the Earth's nutation.

julia> using AstroLib

julia> jd = jdcnv(2028,11,13,4,56)
2.4620887055555554e6

julia> co_nutate(jd,ten(2,46,11.331) * 15,ten(49,20,54.54))
(0.004400660977140092, 0.00172668646508356, 0.40904016038217555, 14.859389427896472, 2.703809037235057)

Notes

Code of this function is based on IDL Astronomy User's Library.

The output of d_ra and d_dec in IDL AstroLib is in arcseconds, however it is in degrees here.

This function calls mean_obliquity and nutate.

source
AstroLib.co_refractFunction
co_refract(old_alt[, altitude=0, pressure=NaN, temperature=NaN,
           epsilon=0.25, to_observe=false]) -> aout

Purpose

Calculate correction to altitude due to atmospheric refraction.

Explanation

Because the index of refraction of air is not precisely 1.0, the atmosphere bends all incoming light, making a star or other celestial object appear at a slightly different altitude (or elevation) than it really is. It is important to understand the following definitions:

  • Observed Altitude: The altitude that a star is seen to be, with a telescope. This is where it appears in the sky. This is should be always greater than the apparent altitude.

  • Apparent Altitude: The altitude that a star would be at, if ~there were no atmosphere~ (sometimes called the "true" altitude). This is usually calculated from an object's celestial coordinates. Apparent altitude should always be smaller than the observed altitude.

Thus, for example, the Sun's apparent altitude when you see it right on the horizon is actually -34 arcminutes.

This program uses a couple of simple formulae to estimate the effect for most optical and radio wavelengths. Typically, you know your observed altitude (from an observation), and want the apparent altitude. To go the other way, this program uses an iterative approach.

Arguments

  • old_alt: observed altitude in degrees. If to_observe is set to true, this should be apparent altitude
  • altitude (optional): the height of the observing location, in meters. This is only used to determine an approximate temperature and pressure, if these are not specified separately. Default is 0 i.e. sea level
  • pressure (optional): the pressure at the observing location, in millibars. Default is NaN
  • temperature (optional): the temperature at the observing location, in Kelvins. Default is NaN
  • epsilon (optional): the accuracy to obtain, in arcseconds. If to_observe is true, then it will be calculated. Default is 0.25 arcseconds
  • to_observe (optional boolean keyword): if set to true, it is assumed that old_alt has apparent altitude as its input and the observed altitude will be found

Output

  • aout: apparent altitude, in degrees. Observed altitude is returned if to_observe is set to true

Example

The lower limb of the Sun is observed to have altitude of 0d 30'. Calculate the the true (i.e. apparent) altitude of the Sun's lower limb using mean conditions of air pressure and temperature.

julia> using AstroLib

julia> co_refract(0.5)
0.02584736873098442

Notes

If altitude is set but the temperature or pressure is not, the program will make an intelligent guess for the temperature and pressure.

Wavelength Dependence

This correction is 0 at zenith, about 1 arcminute at 45 degrees, and 34 arcminutes at the horizon for optical wavelengths. The correction is non-negligible at all wavelengths, but is not very easily calculable. These formulae assume a wavelength of 550 nm, and will be accurate to about 4 arcseconds for all visible wavelengths, for elevations of 10 degrees and higher. Amazingly, they are also accurate for radio frequencies less than ~ 100 GHz.

References

  • Meeus, Astronomical Algorithms, Chapter 15.
  • Explanatory Supplement to the Astronomical Almanac, 1992.
  • Methods of Experimental Physics, Vol 12 Part B, Astrophysics, Radio Telescopes, Chapter 2.5, "Refraction Effects in the Neutral Atmosphere", by R.K. Crane.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.co_refract_forwardMethod
co_refract_forward(alt, pre, temp) -> ref

Purpose

A function used by co_refract to find apparent (or observed) altitude

Arguments

  • alt: the observed (or apparent) altitude, in degrees
  • pre: pressure, in millibars
  • temp: temperature, in Kelvins

Output

  • ref: the atmospheric refraction, in minutes of arc

Notes

The atmospheric refraction is calculated by Saemundsson's formula

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.ct2lstMethod
ct2lst(longitude, jd) -> local_sidereal_time
ct2lst(longitude, tz, date) -> local_sidereal_time

Purpose

Convert from Local Civil Time to Local Mean Sidereal Time.

Arguments

The function can be called in two different ways. The only argument common to both methods is longitude:

  • longitude: the longitude in degrees (east of Greenwich) of the place for which the local sidereal time is desired. The Greenwich mean sidereal time (GMST) can be found by setting longitude = 0.

The civil date to be converted to mean sidereal time can be specified either by providing the Julian days:

  • jd: this is number of Julian days for the date to be converted.

or the time zone and the date:

  • tz: the time zone of the site in hours, positive East of the Greenwich meridian (ahead of GMT). Use this parameter to easily account for Daylight Savings time (e.g. -4=EDT, -5 = EST/CDT).
  • date: this is the local civil time with type DateTime.

Output

The local sidereal time for the date/time specified in hours.

Method

The Julian days of the day and time is question is used to determine the number of days to have passed since 2000-01-01. This is used in conjunction with the GST of that date to extrapolate to the current GST; this is then used to get the LST. See Astronomical Algorithms by Jean Meeus, p. 84 (Eq. 11-4) for the constants used.

Example

Find the Greenwich mean sidereal time (GMST) on 2008-07-30 at 15:53 in Baltimore, Maryland (longitude=-76.72 degrees). The timezone is EDT or tz=-4

julia> using AstroLib, Dates

julia> lst = ct2lst(-76.72, -4, DateTime(2008, 7, 30, 15, 53))
11.356505172312609

julia> sixty(lst)
3-element StaticArrays.SArray{Tuple{3},Float64,1,3} with indices SOneTo(3):
 11.0
 21.0
 23.418620325392112

Find the Greenwich mean sidereal time (GMST) on 2015-11-24 at 13:21 in Heidelberg, Germany (longitude=08° 43' E). The timezone is CET or tz=1. Provide ct2lst only with the longitude of the place and the number of Julian days.

julia> using AstroLib, Dates

julia> longitude=ten(8, 43); # Convert longitude to decimals.

julia> jd = jdcnv(DateTime(2015, 11, 24, 13, 21) - Dates.Hour(1));
# Get number of Julian days. Remember to subtract the time zone in
# order to convert local time to UTC.

julia> lst = ct2lst(longitude, jd) # Calculate Greenwich Mean Sidereal Time.
17.140685171005316

julia> sixty(lst)
3-element StaticArrays.SArray{Tuple{3},Float64,1,3} with indices SOneTo(3):
 17.0
  8.0
 26.466615619137883

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.daycnvFunction
daycnv(julian_days) -> DateTime

Purpose

Converts Julian days number to Gregorian calendar dates.

Explanation

Takes the number of Julian calendar days since epoch -4713-11-24T12:00:00 and returns the corresponding proleptic Gregorian Calendar date.

Argument

  • julian_days: Julian days number.

Output

Proleptic Gregorian Calendar date, of type DateTime, corresponding to the given Julian days number.

Example

julia> using AstroLib

julia> daycnv(2440000)
1968-05-23T12:00:00

Notes

jdcnv is the inverse of this function.

source
AstroLib.dereddMethod
deredd(Eby, by, m1, c1, ub) -> by0, m0, c0, ub0

Purpose

Deredden stellar Stromgren parameters given for a value of E(b-y)

Arguments

  • Eby: color index E(b-y), scalar (E(b-y) = 0.73*E(B-V))
  • by: b-y color (observed)
  • m1: Stromgren line blanketing parameter (observed)
  • c1: Stromgren Balmer discontinuity parameter (observed)
  • ub: u-b color (observed)

All arguments can be either scalars or arrays all of the same length.

Output

The 4-tuple (by0, m0, c0, ub0).

  • by0: b-y color (dereddened)
  • m0: line blanketing index (dereddened)
  • c0: Balmer discontinuity parameter (dereddened)
  • ub0: u-b color (dereddened)

These are scalars or arrays of the same length as the input arguments.

Example

julia> using AstroLib

julia> deredd(0.5, 0.2, 1.0, 1.0, 0.1)
(-0.3, 1.165, 0.905, -0.665)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.eci2geoMethod
eci2geo(x, y, z, jd) -> latitude, longitude, altitude

Purpose

Convert Earth-centered inertial coordinates to geographic spherical coordinates.

Explanation

Converts from ECI (Earth-Centered Inertial) (x, y, z) rectangular coordinates to geographic spherical coordinates (latitude, longitude, altitude). Julian day is also needed as input.

ECI coordinates are in km from Earth center at the supplied time (True of Date). Geographic coordinates assume the Earth is a perfect sphere, with radius equal to its equatorial radius.

Arguments

  • x: ECI x coordinate at jd, in kilometers.
  • y: ECI y coordinate at jd, in kilometers.
  • z: ECI z coordinate at jd, in kilometers.
  • jd: Julian days.

The three coordinates can be passed as a 3-tuple (x, y, z). In addition, x, y, z, and jd can be given as arrays of the same length.

Output

The 3-tuple of geographical coordinate (latitude, longitude, altitude).

  • latitude: latitude, in degrees.
  • longitude: longitude, in degrees.
  • altitude: altitude, in kilometers.

If ECI coordinates are given as arrays, a 3-tuple of arrays of the same length is returned.

Example

Obtain the geographic direction of the vernal point on 2015-06-30T14:03:12.857, in geographic coordinates, at altitude 600 km. Note: equatorial radii of Solar System planets in meters are stored into AstroLib.planets dictionary.

julia> using AstroLib

julia> x = AstroLib.planets["earth"].eqradius*1e-3 + 600;

julia> lat, long, alt = eci2geo(x, 0, 0, jdcnv("2015-06-30T14:03:12.857"))
(0.0, 230.87301833205856, 600.0)

These coordinates can be further transformed into geodetic coordinates using geo2geodetic or into geomagnetic coordinates using geo2mag.

Notes

geo2eci converts geographic spherical coordinates to Earth-centered inertial coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.eq2horFunction
eq2hor(ra, dec, jd[, obsname; ws=false, B1950=false, precession=true, nutate=true,
       aberration=true, refract=true, pressure=NaN, temperature=NaN]) -> alt, az, ha

eq2hor(ra, dec, jd, lat, lon[, altitude=0; ws=false, B1950=false,
       precession=true, nutate=true, aberration=true, refract=true,
       pressure=NaN, temperature=NaN]) -> alt, az, ha

Purpose

Convert celestial (ra-dec) coords to local horizon coords (alt-az).

Explanation

This code calculates horizon (alt,az) coordinates from equatorial (ra,dec) coords. It performs precession, nutation, aberration, and refraction corrections.

Arguments

This function has two base methods. With one you can specify the name of the observatory, if present in AstroLib.observatories, with the other one you can provide the coordinates of the observing site and, optionally, the altitude.

Common mandatory arguments:

  • ra: right ascension of object, in degrees
  • dec: declination of object, in degrees
  • jd: julian date

Other positional arguments:

  • obsname: set this to a valid observatory name in AstroLib.observatories.

or

  • lat: north geodetic latitude of location, in degrees.
  • lon: AST longitude of location, in degrees. You can specify west longitude with a negative sign.
  • altitude: the altitude of the observing location, in meters. It is 0 by default

Optional keyword arguments:

  • ws (optional boolean keyword): set this to true to get the azimuth measured westward from south (not East of North)
  • B1950 (optional boolean keyword): set this to true if the ra and dec are specified in B1950 (FK4 coordinates) instead of J2000 (FK5). This is false by default
  • precession (optional boolean keyword): set this to false for no precession correction, true by default
  • nutate (optional boolean keyword): set this to false for no nutation, true by default
  • aberration (optional boolean keyword): set this to false for no aberration correction, true by default
  • refract (optional boolean keyword): set this to false for no refraction correction, true by default
  • pressure (optional keyword): the pressure at the observing location, in millibars. Default value is NaN
  • temperature (optional keyword): the temperature at the observing location, in Kelvins. Default value is NaN

Output

  • alt: altitude of horizon coords, in degrees
  • az: azimuth angle measured East from North (unless ws is true), in degrees
  • ha: hour angle, in degrees

Example

julia> using AstroLib

julia> alt_o, az_o = eq2hor(ten(6,40,58.2)*15, ten(9,53,44), 2460107.25, ten(50,31,36),
                            ten(6,51,18), 369, pressure = 980, temperature=283)
(16.423991509721567, 265.60656932130564, 76.11502253130612)

julia> adstring(az_o, alt_o)
" 17 42 25.6  +16 25 26"

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.eqpoleMethod
eqpole(l, b[; southpole = false]) -> x, y

Purpose

Convert right ascension $l$ and declination $b$ to coordinate $(x, y)$ using an equal-area polar projection.

Explanation

The output $x$ and $y$ coordinates are scaled to be in the range $[-90, 90]$ and to go from equator to pole to equator. Output map points can be centered on the north pole or south pole.

Arguments

  • l: longitude, scalar or vector, in degrees
  • b: latitude, same number of elements as right ascension, in degrees
  • southpole (optional boolean keyword): keyword to indicate that the plot is to be centered on the south pole instead of the north pole. Default is false.

Output

The 2-tuple $(x, y)$:

  • $x$ coordinate, same number of elements as right ascension, normalized to be in the range $[-90, 90]$.
  • $y$ coordinate, same number of elements as declination, normalized to be in the range $[-90, 90]$.

Example

julia> using AstroLib

julia> eqpole(100, 35, southpole=true)
(-111.18287262822456, -19.604540237028665)

julia> eqpole(80, 19)
(72.78853915267848, 12.83458333897169)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.eulerMethod
euler(ai, bi, select[, FK4=true, radians=true])

Purpose

Transform between Galactic, celestial, and ecliptic coordinates.

Explanation

The function is used by the astro procedure.

Arguments

  • ai: input longitude, scalar or vector.

  • bi: input latitude, scalar or vector.

  • select : integer input specifying type of coordinate transformation. SELECT From To | SELECT From To 1 RA-Dec (2000) Galactic | 4 Ecliptic RA-Dec 2 Galactic RA-DEC | 5 Ecliptic Galactic 3 RA-Dec Ecliptic | 6 Galactic Ecliptic

  • FK4 (optional boolean keyword) : if this keyword is set to true, then input and output celestial and ecliptic coordinates should be given in equinox B1950. When false, by default, they should be given in equinox J2000.

  • radians (optional boolean keyword) : if this keyword is set to true, all input and output angles are in radians rather than degrees.

Output

a 2-tuple (ao, bo):

  • ao: output longitude in degrees.
  • bo: output latitude in degrees.

Example

Find the Galactic coordinates of Cyg X-1 (ra=299.590315, dec=35.201604)

julia> using AstroLib

julia> euler(299.590315, 35.201604, 1)
(71.33498957116959, 3.0668335310640984)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.flux2magFunction
flux2mag(flux[, zero_point, ABwave=number]) -> magnitude

Purpose

Convert from flux expressed in erg/(s cm² Å) to magnitudes.

Explanation

This is the reverse of mag2flux.

Arguments

  • flux: the flux to be converted in magnitude, expressed in erg/(s cm² Å).
  • zero_point: the zero point level of the magnitude. If not

supplied then defaults to 21.1 (Code et al 1976). Ignored if the ABwave keyword is supplied

  • ABwave (optional numeric keyword): wavelength in Angstroms.

If supplied, then returns Oke AB magnitudes (Oke & Gunn 1983, ApJ, 266, 713; http://adsabs.harvard.edu/abs/1983ApJ...266..713O).

Output

The magnitude.

If the ABwave keyword is set then magnitude is given by the expression

\[\text{ABmag} = -2.5\log_{10}(f) - 5\log_{10}(\text{ABwave}) - 2.406\]

Otherwise, magnitude is given by the expression

\[\text{mag} = -2.5\log_{10}(\text{flux}) - \text{zero point}\]

Example

julia> using AstroLib

julia> flux2mag(5.2e-15)
14.609991640913002

julia> flux2mag(5.2e-15, 15)
20.709991640913003

julia> flux2mag(5.2e-15, ABwave=15)
27.423535345634598

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.gal_uvwMethod
gal_uvw(ra, dec, pmra, pmdec, vrad, plx[, lsr=true]) -> u, v, w

Purpose

Calculate the Galactic space velocity $(u, v, w)$ of a star.

Explanation

Calculates the Galactic space velocity $(u, v, w)$ of a star given its (1) coordinates, (2) proper motion, (3) parallax, and (4) radial velocity.

Arguments

User must supply a position, proper motion, radial velocity and parallax. Either scalars or arrays all of the same length can be supplied.

(1) Position:

  • ra: right ascension, in degrees
  • dec: declination, in degrees

(2) Proper Motion

  • pmra: proper motion in right ascension in arc units (typically milli-arcseconds/yr). If given $\mu_\alpha$ – proper motion in seconds of time/year – then this is equal to $15 \mu_\alpha \cos(\text{dec})$.
  • pmdec: proper motion in declination (typically mas/yr).

(3) Radial Velocity

  • vrad: velocity in km/s

(4) Parallax

  • plx: parallax with same distance units as proper motion measurements typically milliarcseconds (mas)

If you know the distance in parsecs, then set plx to $1000/\text{distance}$, if proper motion measurements are given in milli-arcseconds/yr.

There is an additional optional keyword:

  • lsr (optional boolean keyword): if this keyword is set to true, then the output velocities will be corrected for the solar motion $(u, v, w)_\odot = (-8.5, 13.38, 6.49)$ (Coşkunoǧlu et al. 2011 MNRAS, 412, 1237; DOI:10.1111/j.1365-2966.2010.17983.x) to the local standard of rest (LSR). Note that the value of the solar motion through the LSR remains poorly determined.

Output

The 3-tuple $(u, v, w)$

  • $u$: velocity (km/s) positive toward the Galactic anticenter
  • $v$: velocity (km/s) positive in the direction of Galactic rotation
  • $w$: velocity (km/s) positive toward the North Galactic Pole

Method

Follows the general outline of Johnson & Soderblom (1987, AJ, 93, 864; DOI:10.1086/114370) except that $u$ is positive outward toward the Galactic anticenter, and the J2000 transformation matrix to Galactic coordinates is taken from the introduction to the Hipparcos catalog.

Example

Compute the U,V,W coordinates for the halo star HD 6755. Use values from Hipparcos catalog, and correct to the LSR.

julia> using AstroLib

julia> ra=ten(1,9,42.3)*15.; dec = ten(61,32,49.5);

julia> pmra = 627.89;  pmdec = 77.84; # mas/yr

julia> vrad = -321.4; dis = 129; # distance in parsecs

julia> u, v, w = gal_uvw(ra, dec, pmra, pmdec, vrad, 1e3/dis, lsr=true)
(118.2110474553902, -466.4828898385057, 88.16573278565097)

Notes

This function does not take distance as input. See "Arguments" section above for how to provide it using parallax argument.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.gcircMethod
gcirc(units, ra1, dec1, ra2, dec2) -> angular_distance

Purpose

Computes rigorous great circle arc distances.

Explanation

Input position can be either radians, sexagesimal right ascension and declination, or degrees.

Arguments

  • units: integer, can be either 0, or 1, or 2. Describes units of inputs and output:
    • 0: everything (input right ascensions and declinations, and output distance) is radians
    • 1: right ascensions are in decimal hours, declinations in decimal degrees, output distance in arc seconds
    • 2: right ascensions and declinations are in degrees, output distance in arc seconds
  • ra1: right ascension or longitude of point 1
  • dec1: declination or latitude of point 1
  • ra2: right ascension or longitude of point 2
  • dec2: declination or latitude of point 2

Both ra1 and dec1, and ra2 and dec2 can be given as 2-tuples (ra1, dec1) and (ra2, dec2).

Output

Angular distance on the sky between points 1 and 2, as a AbstractFloat. See units argument above for the units.

Method

"Haversine formula" see http://en.wikipedia.org/wiki/Great-circle_distance.

Example

julia> using AstroLib

julia> gcirc(0, 120, -43, 175, +22)
1.590442261600714

Notes

  • The function sphdist provides an alternate method of computing a spherical

distance.

  • The Haversine formula can give rounding errors for antipodal points.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.geo2eciMethod
geo2eci(latitude, longitude, altitude, jd) -> x, y, z

Purpose

Convert geographic spherical coordinates to Earth-centered inertial coordinates.

Explanation

Converts from geographic spherical coordinates (latitude, longitude, altitude) to ECI (Earth-Centered Inertial) (x, y, z) rectangular coordinates. Julian days is also needed.

Geographic coordinates assume the Earth is a perfect sphere, with radius equal to its equatorial radius. ECI coordinates are in km from Earth center at epoch TOD (True of Date).

Arguments

  • latitude: geographic latitude, in degrees.
  • longitude: geographic longitude, in degrees.
  • altitude: geographic altitude, in kilometers.
  • jd: Julian days.

The three coordinates can be passed as a 3-tuple (latitude, longitude, altitude). In addition, latitude, longitude, altitude, and jd can be given as arrays of the same length.

Output

The 3-tuple of ECI (x, y, z) coordinates, in kilometers. The TOD epoch is the supplied jd time.

If geographical coordinates are given as arrays, a 3-tuple of arrays of the same length is returned.

Example

Obtain the ECI coordinates of the intersection of the equator and Greenwich's meridian on 2015-06-30T14:03:12.857

julia> using AstroLib

julia> geo2eci(0, 0, 0, jdcnv("2015-06-30T14:03:12.857"))
(-4024.8671780315185, 4947.835465127513, 0.0)

Notes

eci2geo converts Earth-centered inertial coordinates to geographic spherical coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.geo2geodeticMethod
geo2geodetic(latitude, longitude, altitude) -> latitude, longitude, altitude
geo2geodetic(latitude, longitude, altitude, planet) -> latitude, longitude, altitude
geo2geodetic(latitude, longitude, altitude, equatorial_radius, polar_radius) -> latitude, longitude, altitude

Purpose

Convert from geographic (or planetographic) to geodetic coordinates.

Explanation

Converts from geographic (latitude, longitude, altitude) to geodetic (latitude, longitude, altitude). In geographic coordinates, the Earth is assumed a perfect sphere with a radius equal to its equatorial radius. The geodetic (or ellipsoidal) coordinate system takes into account the Earth's oblateness.

Geographic and geodetic longitudes are identical. Geodetic latitude is the angle between local zenith and the equatorial plane. Geographic and geodetic altitudes are both the closest distance between the satellite and the ground.

Arguments

The function has two base methods. The arguments common to all methods and always mandatory are latitude, longitude, and altitude:

  • latitude: geographic latitude, in degrees.
  • longitude: geographic longitude, in degrees.
  • altitude: geographic altitude, in kilometers.

In order to convert to geodetic coordinates, you can either provide custom equatorial and polar radii of the planet or use the values of one of the planets of Solar System (Pluto included).

If you want to use the method with explicit equatorial and polar radii the additional mandatory arguments are:

  • equatorial_radius: value of the equatorial radius of the body, in kilometers.
  • polar_radius: value of the polar radius of the body, in kilometers.

Instead, if you want to use the method with the selection of a planet, the only additional argument is the planet name:

  • planet (optional string argument): string with the name of the Solar System planet, from "Mercury" to "Pluto". If omitted (so, when only latitude, longitude, and altitude are provided), the default is "Earth".

In all cases, the three coordinates can be passed as a 3-tuple (latitude, longitude, altitude). In addition, geographical latitude, longitude, and altitude can be given as arrays of the same length.

Output

The 3-tuple (latitude, longitude, altitude) in geodetic coordinates, for the body with specified equatorial and polar radii (Earth by default).

If geographical coordinates are given as arrays, a 3-tuple of arrays of the same length is returned.

Method

Stephen P. Keeler and Yves Nievergelt, "Computing geodetic coordinates", SIAM Rev. Vol. 40, No. 2, pp. 300-309, June 1998 (DOI:10.1137/S0036144597323921).

Planetary constants are from Planetary Fact Sheet (http://nssdc.gsfc.nasa.gov/planetary/factsheet/index.html).

Example

Locate the Earth geographic North pole (latitude: 90°, longitude: 0°, altitude 0 km), in geodetic coordinates:

julia> using AstroLib

julia> geo2geodetic(90, 0, 0)
(90.0, 0.0, 21.38499999999931)

The same for Jupiter:

julia> using AstroLib

julia> geo2geodetic(90, 0, 0, "Jupiter")
(90.0, 0.0, 4638.0)

Find geodetic coordinates for point of geographic coordinates (latitude, longitude, altitude) = (43.16°, -24.32°, 3.87 km) on a planet with equatorial radius 8724.32 km and polar radius 8619.19 km:

julia> using AstroLib

julia> geo2geodetic(43.16, -24.32, 3.87, 8724.32, 8619.19)
(43.849399515234516, -24.32, 53.53354478670965)

Notes

Whereas the conversion from geodetic to geographic coordinates is given by an exact, analytical formula, the conversion from geographic to geodetic isn't. Approximative iterations (as used here) exist, but tend to become less good with increasing eccentricity and altitude. The formula used in this routine should give correct results within six digits for all spatial locations, for an ellipsoid (planet) with an eccentricity similar to or less than Earth's. More accurate results can be obtained via calculus, needing a non-determined amount of iterations.

In any case, the function geodetic2geo, which converts from geodetic (or planetodetic) to geographic coordinates, can be used to estimate the accuracy of geo2geodetic.

julia> using AstroLib

julia> collect(geodetic2geo(geo2geodetic(67.2, 13.4, 1.2))) - [67.2, 13.4, 1.2]
3-element Array{Float64,1}:
 -3.5672513831741526e-9
  0.0
  9.484211194177306e-10

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.geo2magFunction
geo2mag(latitude, longitude[, year]) -> geomagnetic_latitude, geomagnetic_longitude

Purpose

Convert from geographic to geomagnetic coordinates.

Explanation

Converts from geographic (latitude, longitude) to geomagnetic (latitude, longitude). Altitude is not involved in this function.

Arguments

  • latitude: geographic latitude (North), in degrees.
  • longitude: geographic longitude (East), in degrees.
  • year (optional numerical argument): the year in which to perform conversion. If omitted, defaults to current year.

The coordinates can be passed as arrays of the same length.

Output

The 2-tuple of magnetic (latitude, longitude) coordinates, in degrees.

If geographical coordinates are given as arrays, a 2-tuple of arrays of the same length is returned.

Example

Kyoto has geographic coordinates 35° 00' 42'' N, 135° 46' 06'' E, find its geomagnetic coordinates in 2016:

julia> using AstroLib

julia> geo2mag(ten(35,0,42), ten(135,46,6), 2016)
(36.86579228937769, -60.184060536651614)

Notes

This function uses list of North Magnetic Pole positions provided by World Magnetic Model (https://www.ngdc.noaa.gov/geomag/data/poles/NP.xy).

mag2geo converts geomagnetical coordinates to geographic coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.geodetic2geoMethod
geodetic2geo(latitude, longitude, altitude) -> latitude, longitude, altitude
geodetic2geo(latitude, longitude, altitude, planet) -> latitude, longitude, altitude
geodetic2geo(latitude, longitude, altitude, equatorial_radius, polar_radius) -> latitude, longitude, altitude

Purpose

Convert from geodetic (or planetodetic) to geographic coordinates.

Explanation

Converts from geodetic (latitude, longitude, altitude) to geographic (latitude, longitude, altitude). In geographic coordinates, the Earth is assumed a perfect sphere with a radius equal to its equatorial radius. The geodetic (or ellipsoidal) coordinate system takes into account the Earth's oblateness.

Geographic and geodetic longitudes are identical. Geodetic latitude is the angle between local zenith and the equatorial plane. Geographic and geodetic altitudes are both the closest distance between the satellite and the ground.

Arguments

The function has two base methods. The arguments common to all methods and always mandatory are latitude, longitude, and altitude:

  • latitude: geodetic latitude, in degrees.
  • longitude: geodetic longitude, in degrees.
  • altitude: geodetic altitude, in kilometers.

In order to convert to geographic coordinates, you can either provide custom equatorial and polar radii of the planet or use the values of one of the planets of Solar System (Pluto included).

If you want to use the method with explicit equatorial and polar radii the additional mandatory arguments are:

  • equatorial_radius: value of the equatorial radius of the body, in kilometers.
  • polar_radius: value of the polar radius of the body, in kilometers.

Instead, if you want to use the method with the selection of a planet, the only additional argument is the planet name:

  • planet (optional string argument): string with the name of the Solar System planet, from "Mercury" to "Pluto". If omitted (so, when only latitude, longitude, and altitude are provided), the default is "Earth".

In all cases, the three coordinates can be passed as a 3-tuple (latitude, longitude, altitude). In addition, geodetic latitude, longitude, and altitude can be given as arrays of the same length.

Output

The 3-tuple (latitude, longitude, altitude) in geographic coordinates, for the body with specified equatorial and polar radii (Earth by default).

If geodetic coordinates are given as arrays, a 3-tuple of arrays of the same length is returned.

Method

Stephen P. Keeler and Yves Nievergelt, "Computing geodetic coordinates", SIAM Rev. Vol. 40, No. 2, pp. 300-309, June 1998 (DOI:10.1137/S0036144597323921).

Planetary constants from "Allen's Astrophysical Quantities", Fourth Ed., (2000).

Example

Find geographic coordinates of geodetic North pole (latitude: 90°, longitude: 0°, altitude 0 km) of the Earth:

julia> using AstroLib

julia> geodetic2geo(90, 0, 0)
(90.0, 0.0, -21.38499999999931)

The same for Jupiter:

julia> using AstroLib

julia> geodetic2geo(90, 0, 0, "Jupiter")
(90.0, 0.0, -4638.0)

Find geographic coordinates for point of geodetic coordinates (latitude, longitude, altitude) = (43.16°, -24.32°, 3.87 km) on a planet with equatorial radius 8724.32 km and polar radius 8619.19 km:

julia> using AstroLib

julia> geodetic2geo(43.16, -24.32, 3.87, 8724.32, 8619.19)
(42.46772711708433, -24.32, -44.52902080669082)

Notes

geo2geodetic converts from geographic (or planetographic) to geodetic coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.get_dateMethod
get_date([date, old=true, timetag=true]) -> string

Purpose

Returns the UTC date in "CCYY-MM-DD" format for FITS headers.

Explanation

This is the format required by the DATE and DATE-OBS keywords in a FITS header.

Argument

  • date (optional): the date in UTC standard. If omitted, defaults to the current UTC time. Each element can be either a DateTime type or anything that can be converted to that type. In the case of vectorial input, each element is considered as a date, so you cannot provide a date by parts.
  • old (optional boolean keyword): see below.
  • timetag (optional boolean keyword): see below.

Output

A string with the date formatted according to the given optional keywords.

  • When no optional keywords (timetag and old) are supplied, the format of the output string is "CCYY-MM-DD" (year-month-day part of the date), where CCYY represents a 4-digit calendar year, MM the 2-digit ordinal number of a calendar month within the calendar year, and DD the 2-digit ordinal number of a day within the calendar month.
  • If the boolean keyword old is true (default: false), the year-month-day part of date has "DD/MM/YY" format. This is the formerly (pre-1997) recommended for FITS. Note that this format is now deprecated because it uses only a 2-digit representation of the year.
  • If the boolean keyword timetag is true (default: false), "Thh:mm:ss" is appended to the year-month-day part of the date, where <hh> represents the hour in the day, <mm> the minutes, <ss> the seconds, and the literal 'T' the ISO 8601 time designator.

Note that old and timetag keywords can be used together, so that the output string will have "DD/MM/YYThh:mm:ss" format.

Example

julia> using AstroLib, Dates

julia> get_date(DateTime(21937, 05, 30, 09, 59, 00), timetag=true)
"21937-05-30T09:59:00"

Notes

  1. A discussion of the DATExxx syntax in FITS headers can be found in

http://www.cv.nrao.edu/fits/documents/standards/year2000.txt

  1. Those who wish to use need further flexibility in their date formats (e.g. to

use TAI time) should look at Bill Thompson's time routines in http://sohowww.nascom.nasa.gov/solarsoft/gen/idl/time

source
AstroLib.get_juldateMethod
get_juldate() -> julian_days

Purpose

Return the number of Julian days for current time.

Explanation

Return for current time the number of Julian calendar days since epoch -4713-11-24T12:00:00 as a floating point.

Example

get_juldate()
daycnv(get_juldate())

Notes

Use jdcnv to get the number of Julian days for a different date.

source
AstroLib.hadec2altazMethod
hadec2altaz(ha, dec, lat[, ws=true]) -> alt, az

Purpose

Convert Hour Angle and Declination to Horizon (Alt-Az) coordinates.

Explanation

Can deal with the NCP singularity. Intended mainly to be used by program eq2hor.

Arguments

Input coordinates may be either a scalar or an array, of the same dimension.

  • ha: the local apparent hour angle, in degrees. The hour angle is the time that right ascension of 0 hours crosses the local meridian. It is unambiguously defined.
  • dec: the local apparent declination, in degrees.
  • lat: the local geodetic latitude, in degrees, scalar or array.
  • ws (optional boolean keyword): if true, the output azimuth is measured West from South. The default is to measure azimuth East from North.

ha and dec can be given as a 2-tuple (ha, dec).

Output

2-tuple (alt, az)

  • alt: local apparent altitude, in degrees.
  • az: the local apparent azimuth, in degrees.

The output coordinates are always floating points and have the same type (scalar or array) as the input coordinates.

Example

Arcturus is observed at an apparent hour angle of 336.6829 and a declination of 19.1825 while at the latitude of +43° 4' 42''. What are the local altitude and azimuth of this object?

julia> using AstroLib

julia> alt, az = hadec2altaz(336.6829, 19.1825, ten(43, 4, 42))
(59.08617155005685, 133.3080693440254)

Notes

altaz2hadec converts Horizon (Alt-Az) coordinates to Hour Angle and Declination.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.helioFunction
helio(jd, list[, radians=true]) -> hrad, hlong, hlat

Purpose

Compute heliocentric coordinates for the planets.

Explanation

The mean orbital elements for epoch J2000 are used. These are derived from a 250 yr least squares fit of the DE 200 planetary ephemeris to a Keplerian orbit where each element is allowed to vary linearly with time. Useful mainly for dates between 1800 and 2050, this solution fits the terrestrial planet orbits to ~25'' or better, but achieves only ~600'' for Saturn.

Arguments

  • jd: julian date, scalar or vector
  • num: integer denoting planet number, scalar or vector 1 = Mercury, 2 = Venus, ... 9 = Pluto
  • radians(optional): if this keyword is set to true, than the longitude and latitude output are in radians rather than degrees.

Output

  • hrad: the heliocentric radii, in astronomical units.
  • hlong: the heliocentric (ecliptic) longitudes, in degrees.
  • hlat: the heliocentric latitudes in degrees.

Example

(1) Find heliocentric position of Venus on August 23, 2000

julia> using AstroLib

julia> helio(jdcnv(2000,08,23,0), 2)
(0.7213758288364316, 198.39093251916148, 2.887355631705488)

(2) Find the current heliocentric positions of all the planets

julia> using AstroLib

julia> helio.([jdcnv(1900)], 1:9)
9-element Array{Tuple{Float64,Float64,Float64},1}:
 (0.4207394142180803, 202.60972662618906, 3.0503005607270532)
 (0.7274605731764012, 344.5381482401048, -3.3924346961624785)
 (0.9832446886519147, 101.54969268801035, 0.012669354526696368)
 (1.4212659241051142, 287.8531100442217, -1.5754626002228043)
 (5.386813769590955, 235.91306092135062, 0.9131692817310215)
 (10.054339927304339, 268.04069870870387, 1.0851704598594278)
 (18.984683376211326, 250.0555468087738, 0.05297087029604253)
 (29.87722677219009, 87.07244903504716, -1.245060583142733)
 (46.9647515992327, 75.94692594417324, -9.576681044165511)

Notes

This program is based on the two-body model and thus neglects interactions between the planets.

The coordinates are given for equinox 2000 and not the equinox of the supplied date.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.helio_jdMethod
helio_jd(date, ra, dec[, B1950=true, diff=false]) -> jd_helio
helio_jd(date, ra, dec[, B1950=true, diff=true]) -> time_diff

Purpose

Convert geocentric (reduced) Julian date to heliocentric Julian date.

Explanation

This procedure corrects for the extra light travel time between the Earth and the Sun.

An online calculator for this quantity is available at http://www.physics.sfasu.edu/astro/javascript/hjd.html

Users requiring more precise calculations and documentation should look at the IDL code available at http://astroutils.astronomy.ohio-state.edu/time/

Arguments

  • date: reduced Julian date (= JD - 2400000). You can use juldate() to calculate the reduced Julian date.
  • ra and dec: right ascension and declination in degrees. Default equinox is J2000.
  • B1950 (optional boolean keyword): if set to true, then input coordinates are assumed to be in equinox B1950 coordinates. Default is false.
  • diff (optional boolean keyword): if set to true, the function returns the time difference (heliocentric JD - geocentric JD) in seconds. Default is false.

Output

The return value depends on the value of diff optional keywords:

  • if diff is false (default), then the heliocentric reduced Julian date is returned.
  • if diff is true, then the time difference in seconds between the geocentric and heliocentric Julian date is returned.

Example

What is the heliocentric Julian date of an observation of V402 Cygni (J2000: RA = 20 9 7.8, Dec = 37 09 07) taken on June 15, 2016 at 11:40 UT?

julia> using AstroLib

julia> jd = juldate(2016, 6, 15, 11, 40);

julia> helio_jd(jd, ten(20, 9, 7.8) * 15, ten(37, 9, 7))
57554.98808289718

Notes

Wayne Warren (Raytheon ITSS) has compared the results of this algorithm with the FORTRAN subroutines in the STARLINK SLALIB library (see http://star-www.rl.ac.uk/).

                                                 Time Diff (sec)
     Date               RA(2000)   Dec(2000)  STARLINK      IDL

1999-10-29T00:00:00.0  21 08 25.  -67 22 00.  -59.0        -59.0
1999-10-29T00:00:00.0  02 56 33.4 +00 26 55.  474.1        474.1
1940-12-11T06:55:00.0  07 34 41.9 -00 30 42.  366.3        370.2
1992-02-29T03:15:56.2  12 56 27.4 +42 10 17.  350.8        350.9
2000-03-01T10:26:31.8  14 28 36.7 -20 42 11.  243.7        243.7
2100-02-26T09:18:24.2  08 26 51.7 +85 47 28.  104.0        108.8

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.helio_rvFunction
helio_rv(jd, T, P, V_0, K[, e, ω]) -> rv

Purpose

Return the heliocentric radial velocity of a spectroscopic binary.

Explanation

This function will return the heliocentric radial velocity of a spectroscopic binary star at a given heliocentric date given its orbit.

Arguments

  • jd: time of observation, as number of Julian days.
  • T: time of periastron passage (max. +ve velocity for circular orbits), same time system as jd
  • P: the orbital period in same units as jd
  • V_0: systemic velocity
  • K: velocity semi-amplitude in the same units as V_0
  • e: eccentricity of the orbit. It defaults to 0 if omitted
  • ω: longitude of periastron in degrees. It defaults to 0 if omitted

Output

The predicted heliocentric radial velocity in the same units as Gamma for the date(s) specified by jd.

Example

(1) What was the heliocentric radial velocity of the primary component of HU Tau at 1730 UT 25 Oct 1994?

julia> using AstroLib

julia> jd = juldate(94, 10, 25, 17, 30); # Obtain Geocentric Julian days

julia> hjd = helio_jd(jd, ten(04, 38, 16) * 15, ten(20, 41, 05)); # Convert to HJD

julia> helio_rv(hjd, 46487.5303, 2.0563056, -6, 59.3)
-62.965570107789475

NB: the functions juldate and helio_jd return a reduced HJD (HJD - 2400000) and so T and P must be specified in the same fashion.

(2) Plot two cycles of an eccentric orbit, $e=0.6$, $\omega=45\degree$ for both components of a binary star. Use PyPlot.jl for plotting.

using PyPlot
φ = range(0, stop=2, length=1000); # Generate 1000 phase points
plot(φ ,helio_rv.(φ, 0, 1, 0, 100, 0.6, 45)) # Plot 1st component
plot(φ ,helio_rv.(φ, 0, 1, 0, 100, 0.6, 45+180)) # Plot 2nd component

Notes

The user should ensure consistency with all time systems being used (i.e. jd and t should be in the same units and time system). Generally, users should reduce large time values by subtracting a large constant offset, which may improve numerical accuracy.

If using the the function juldate and helio_jd, the reduced HJD time system must be used throughtout.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.hor2eqFunction
hor2eq(alt, az, jd[, obsname; ws=false, B1950=false, precession=true, nutate=true,
       aberration=true, refract=true, lat=NaN, lon=NaN, altitude=0, pressure=NaN,
       temperature=NaN]) -> ra, dec, ha

hor2eq(alt, az, jd, lat, lon[, altitude=0; ws=false, B1950=false,
       precession=true, nutate=true, aberration=true, refract=true, pressure=NaN,
       temperature=NaN]) -> ra, dec, ha

Purpose

Converts local horizon coordinates (alt-az) to equatorial (ra-dec) coordinates.

Explanation

This is a function to calculate equatorial (ra,dec) coordinates from horizon (alt,az) coords. It is accurate to about 1 arcsecond or better. It performs precession, nutation, aberration, and refraction corrections.

Arguments

This function has two base methods. With one you can specify the name of the observatory, if present in AstroLib.observatories, with the other one you can provide the coordinates of the observing site and, optionally, the altitude.

Common mandatory arguments:

  • alt: altitude of horizon coords, in degrees
  • az: azimuth angle measured East from North (unless ws is true), in degrees
  • jd: julian date

Other positional arguments:

  • obsname: set this to a valid observatory name in AstroLib.observatories.

or

  • lat: north geodetic latitude of location, in degrees.
  • lon: AST longitude of location, in degrees. You can specify west longitude with a negative sign.
  • altitude: the altitude of the observing location, in meters. It is 0 by default

Optional keyword arguments:

  • ws (optional boolean keyword): set this to true to get the azimuth measured westward from south. This is false by default
  • B1950 (optional boolean keyword): set this to true if the ra and dec are specified in B1950 (FK4 coordinates) instead of J2000 (FK5). This is false by default
  • precession (optional boolean keyword): set this to false for no precession, true by default
  • nutate (optional boolean keyword): set this to false for no nutation, true by default
  • aberration (optional boolean keyword): set this to false for no aberration correction, true by default
  • refract (optional boolean keyword): set this to false for no refraction correction, true by default
  • pressure (optional keyword): the pressure at the observing location, in millibars. Default value is NaN
  • temperature (optional keyword): the temperature at the observing location, in Kelvins. Default value is NaN

Output

  • ra: right ascension of object, in degrees (FK5)
  • dec: declination of the object, in degrees (FK5)
  • ha: hour angle, in degrees

Example

You are at Kitt Peak National Observatory, looking at a star at azimuth angle 264d 55m 06s and elevation 37d 54m 41s (in the visible). Today is Dec 25, 2041 and the local time is 10 PM precisely. What is the right ascension and declination (J2000) of the star you're looking at? The temperature here is about 0 Celsius, and the pressure is 781 millibars. The Julian date for this time is 2466879.7083333

julia> using AstroLib

julia> ra_o, dec_o = hor2eq(ten(37,54,41), ten(264,55,06), 2466879.7083333,
                            "kpno", pressure = 781, temperature = 273)
(3.3224480269254717, 15.19061543702944, 54.61174536229464)

julia> adstring(ra_o, dec_o)
" 00 13 17.4  +15 11 26"

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.imfMethod
imf(mass, expon, mass_range) -> psi

Purpose

Compute an N-component power-law logarithmic initial mass function (IMF).

Explanation

The function is normalized so that the total mass distribution equals one solar mass.

Arguments

  • mass: mass in units of solar mass, vector.
  • expon: power law exponent, vector. The number of values in expon equals the number of different power-law components in the IMF.
  • mass_range: vector containing the mass upper and lower limits of the IMF and masses where the IMF exponent changes. The number of values in massrange should be one more than in expon. The values in massrange should be monotonically increasing and positive.

Output

  • psi: mass function, number of stars per unit logarithmic mass interval evaluated for supplied masses.

Example

Show the number of stars per unit mass interval at 3 Msun for a Salpeter (expon = -1.35) IMF, with a mass range from 0.1 MSun to 110 Msun.

julia> using AstroLib

julia> imf([3], [-1.35], [0.1, 110]) / 3
1-element Array{Float64,1}:
 0.01294143518151214

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.ismeuvFunction
ismeuv(wave, hcol[, he1col=hcol*0.1, he2col=0, fano=false]) -> tau

Purpose

Compute the continuum interstellar EUV optical depth

Explanation

The EUV optical depth is computed from the photoionization of hydrogen and helium.

Arguments

  • wave: wavelength value (in Angstroms). Useful range is 40 - 912 A; at shorter wavelength metal opacity should be considered, at longer wavelengths there is no photoionization.
  • hcol: interstellar hydrogen column density in cm-2.
  • he1col (optional): neutral helium column density in cm-2. Default is 0.1*hcol (10% of hydrogen column)
  • he2col (optional): ionized helium column density in cm-2 Default is 0.
  • fano (optional boolean keyword): If this keyword is true, then the 4 strongest auto-ionizing resonances of He I are included. The shape of these resonances is given by a Fano profile - see Rumph, Bowyer, & Vennes 1994, AJ, 107, 2108. If these resonances are included then the input wavelength vector should have a fine (>~0.01 A) grid between 190 A and 210 A, since the resonances are very narrow.

Output

  • tau: Vector giving resulting optical depth, non-negative values.

Example

One has a model EUV spectrum with wavelength, w (in Angstroms). Find the EUV optical depth by 1e18 cm-2 of HI, with N(HeI)/N(HI) = N(HeII)/N(HI) = 0.05.

julia> using AstroLib

julia> ismeuv.([670, 910], 1e19, 5e17, 5e17)
2-element Array{Float64,1}:
 27.35393320556168
 62.683796028917286

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.jdcnvFunction
jdcnv(date) -> julian_days

Purpose

Convert proleptic Gregorian Calendar date in UTC standard to number of Julian days.

Explanation

Takes the given proleptic Gregorian date in UTC standard and returns the number of Julian calendar days since epoch -4713-11-24T12:00:00.

Argument

  • date: date in proleptic Gregorian Calendar. Each element can be either a DateTime or anything that can be converted directly to DateTime.

Output

Number of Julian days, as a floating point.

Example

Find the Julian days number at 2016 August 23, 03:39:06.

julia> using AstroLib, Dates

julia> jdcnv(DateTime(2016, 08, 23, 03, 39, 06))
2.4576236521527776e6

julia> jdcnv(2016, 08, 23, 03, 39, 06)
2.4576236521527776e6

julia> jdcnv("2016-08-23T03:39:06")
2.4576236521527776e6

Notes

This is the inverse of daycnv.

get_juldate returns the number of Julian days for current time. It is equivalent to jdcnv(now(Dates.UTC)).

For the conversion of Julian date to number of Julian days, use juldate.

source
AstroLib.jprecessFunction
jprecess(ra, dec[, epoch]) -> ra2000, dec2000
jprecess(ra, dec, muradec[, parallax=parallax, radvel=radvel]) -> ra2000, dec2000

Purpose

Precess positions from B1950.0 (FK4) to J2000.0 (FK5).

Explanation

Calculate the mean place of a star at J2000.0 on the FK5 system from the mean place at B1950.0 on the FK4 system.

jprecess function has two methods, one for each of the following cases:

  • the proper motion is known and non-zero
  • the proper motion is unknown or known to be exactly zero (i.e. extragalactic radio sources). Better precision can be achieved in this case by inputting the epoch of the original observations.

Arguments

The function has 2 methods. The common mandatory arguments are:

  • ra: input B1950 right ascension, in degrees.
  • dec: input B1950 declination, in degrees.

The two methods have a different third argument (see "Explanation" section for more details). It can be one of the following:

  • muradec: 2-element vector containing the proper motion in seconds of arc per tropical century in right ascension and declination.
  • epoch: scalar giving epoch of original observations.

If none of these two arguments is provided (so jprecess is fed only with right ascension and declination), it is assumed that proper motion is exactly zero and epoch = 1950.

If it is used the method involving muradec argument, the following keywords are available:

  • parallax (optional numerical keyword): stellar parallax, in seconds of arc.
  • radvel (optional numerical keyword): radial velocity in km/s.

Right ascension and declination can be passed as the 2-tuple (ra, dec). You can also pass ra, dec, parallax, and radvel as arrays, all of the same length N. In that case, muradec should be a matrix 2×N.

Output

The 2-tuple of right ascension and declination in 2000, in degrees, of input coordinates is returned. If ra and dec (and other possible optional arguments) are arrays, the 2-tuple of arrays (ra2000, dec2000) of the same length as the input coordinates is returned.

Method

The algorithm is taken from the Explanatory Supplement to the Astronomical Almanac 1992, page 184. See also Aoki et al (1983), A&A, 128, 263. URL: http://adsabs.harvard.edu/abs/1983A%26A...128..263A.

Example

The SAO catalogue gives the B1950 position and proper motion for the star HD 119288. Find the J2000 position.

  • RA(1950) = 13h 39m 44.526s
  • Dec(1950) = 8d 38' 28.63''
  • Mu(RA) = -.0259 s/yr
  • Mu(Dec) = -.093 ''/yr
julia> using AstroLib

julia> muradec = 100 * [-15*0.0259, -0.093]; # convert to century proper motion

julia> ra = ten(13, 39, 44.526)*15;

julia> decl = ten(8, 38, 28.63);

julia> adstring(jprecess(ra, decl, muradec), precision=2)
" 13 42 12.740  +08 23 17.69"

Notes

"When transferring individual observations, as opposed to catalog mean place, the safest method is to tranform the observations back to the epoch of the observation, on the FK4 system (or in the system that was used to to produce the observed mean place), convert to the FK5 system, and transform to the the epoch and equinox of J2000.0" – from the Explanatory Supplement (1992), p. 180

bprecess performs the precession to B1950 coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.juldateMethod
juldate(date::DateTime) -> reduced_julia_days

Purpose

Convert from calendar to Reduced Julian Days.

Explanation

Julian Day Number is a count of days elapsed since Greenwich mean noon on 1 January 4713 B.C. Julian Days are the number of Julian days followed by the fraction of the day elapsed since the preceding noon.

This function takes the given date and returns the number of Julian calendar days since epoch 1858-11-16T12:00:00 (Reduced Julian Days = Julian Days - 2400000).

Argument

  • date: date in Julian Calendar, UTC standard. Each element can be given in DateTime type or anything that can be converted to that type.

Output

The number of Reduced Julian Days is returned.

Example

Get number of Reduced Julian Days at 2016-03-20T15:24:00.

julia> using AstroLib, Dates

julia> juldate(DateTime(2016, 03, 20, 15, 24))
57468.14166666667

julia> juldate(2016, 03, 20, 15, 24)
57468.14166666667

julia> juldate("2016-03-20T15:24")
57468.14166666667

Notes

Julian Calendar is assumed, thus before 1582-10-15T00:00:00 this function is not the inverse of daycnv. For the conversion proleptic Gregorian date to number of Julian days, use jdcnv, which is the inverse of daycnv.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.kepler_solverFunction
kepler_solver(M, e) -> E

Purpose

Solve Kepler's equation in the elliptic motion regime ($0 \leq e \leq 1$) and return eccentric anomaly $E$.

Explanation

In order to find the position of a body in elliptic motion (e.g., in the two-body problem) at a given time $t$, one has to solve the Kepler's equation

$M(t) = E(t) - e\sin E(t)$

where $M(t) = (t - t_{0})/P$ is the mean anomaly, $E(t)$ the eccentric anomaly, $e$ the eccentricity of the orbit, $t_0$ is the time of periapsis passage, and $P$ is the period of the orbit. Usually the eccentricity is given and one wants to find the eccentric anomaly $E(t)$ at a specific time $t$, so that also the mean anomaly $M(t)$ is known.

Arguments

  • M: mean anomaly.
  • e: eccentricity, in the elliptic motion regime ($0 \leq e \leq 1$)

Output

The eccentric anomaly $E$, restricted to the range $[-\pi, \pi]$.

Method

Many different numerical methods exist to solve Kepler's equation. This function implements the algorithm proposed in Markley (1995) Celestial Mechanics and Dynamical Astronomy, 63, 101 (DOI:10.1007/BF00691917). This method is not iterative, requires only four transcendental function evaluations, and has been proved to be fast and efficient over the entire range of elliptic motion $0 \leq e \leq 1$.

Example

(1) Find the eccentric anomaly for an orbit with eccentricity $e = 0.7$ and for $M(t) = 8\pi/3$.

julia> using AstroLib

julia> ecc = 0.7;

julia> E = kepler_solver(8pi/3, ecc)
2.5085279492864223

(2) Plot the eccentric anomaly as a function of mean anomaly for eccentricity $e = 0$, $0.5$, $0.9$. Recall that kepler_solver gives $E \in [-\pi, \pi]$, use mod2pi to have it in $[0, 2\pi]$. Use PyPlot.jl for plotting.

using AstroLib, PyPlot
M = range(0, stop=2pi, length=1001)[1:end-1];
for ecc in (0, 0.5, 0.9); plot(M, mod2pi.(kepler_solver.(M, ecc))); end

Notes

The true anomaly can be calculated with trueanom function.

source
AstroLib.lsf_rotateFunction
lsf_rotate(delta_v, v_sin_i[, epsilon = 0.3]) -> velocity_grid, lsf

Purpose

Create a 1-d convolution kernel to broaden a spectrum from a rotating star.

Explanation

Can be used to derive the broadening effect (LSF, line spread function) due to rotation on a synthetic stellar spectrum. Assumes constant limb darkening across the disk.

Arguments

  • delta_v: the step increment (in km/s) in the output rotation kernel

  • v_sin_i: the rotational velocity projected along the line of sight (km/s)

  • epsilon (optional numeric argument): the limb-darkening coefficient, default = 0.6 which is typical for photospheric lines. The specific intensity $I$ at any angle $\theta$ from the specific intensity $I_{\text{cen}}$ at the center of the disk is given by:

    $I = I_{\text{cen}}\cdot(1 - \varepsilon\cdot(1 - \cos(\theta)))$

Output

The 2-tuple (velocity_grid, lsf):

  • velocity_grid: vector of velocity grid with the same number of elements as lsf (see below)
  • lsf: the convolution kernel vector for the specified rotational velocity. The number of points in lsf will be always be odd (the kernel is symmetric) and equal to either ceil(2*v_sin_i/delta_v) or ceil(2*v_sin_i/delta_v) + 1, whichever number is odd. Elements of lsf will always be of type AbstractFloat. To actually compute the broadening, the spectrum should be convolved with the rotational lsf

Example

Plot the line spread function for a star rotating at 90 km/s in velocity space every 3 km/s. Use PyPlot.jl for plotting.

using PyPlot
plot(lsf_rotate(3, 90)...)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.mag2fluxFunction
mag2flux(mag[, zero_point, ABwave=number]) -> flux

Purpose

Convert from magnitudes to flux expressed in erg/(s cm² Å).

Explanation

This is the reverse of flux2mag.

Arguments

  • mag: the magnitude to be converted in flux.
  • zero_point: the zero point level of the magnitude. If not supplied then defaults to

21.1 (Code et al 1976). Ignored if the ABwave keyword is supplied

  • ABwave (optional numeric keyword): wavelength, in Angstroms. If supplied, then the

input mag is assumed to contain Oke AB magnitudes (Oke & Gunn 1983, ApJ, 266, 713; http://adsabs.harvard.edu/abs/1983ApJ...266..713O).

Output

The flux.

If the ABwave keyword is set, then the flux is given by the expression

\[\text{flux} = 10^{-0.4(\text{mag} +2.406 + 4\log_{10}(\text{ABwave}))}\]

Otherwise the flux is given by

\[\text{flux} = 10^{-0.4(\text{mag} + \text{zero point})}\]

Example

julia> using AstroLib

julia> mag2flux(8.3)
1.7378008287493692e-12

julia> mag2flux(8.3, 12)
7.58577575029182e-9

julia> mag2flux(8.3, ABwave=12)
3.6244115683017193e-7

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.mag2geoFunction
mag2geo(latitude, longitude[, year]) -> geographic_latitude, geographic_longitude

Purpose

Convert from geomagnetic to geographic coordinates.

Explanation

Converts from geomagnetic (latitude, longitude) to geographic (latitude, longitude). Altitude is not involved in this function.

Arguments

  • latitude: geomagnetic latitude (North), in degrees.
  • longitude: geomagnetic longitude (East), in degrees.
  • year (optional numerical argument): the year in which to perform conversion. If omitted, defaults to current year.

The coordinates can be passed as arrays of the same length.

Output

The 2-tuple of geographic (latitude, longitude) coordinates, in degrees.

If geomagnetic coordinates are given as arrays, a 2-tuple of arrays of the same length is returned.

Example

Find position of North Magnetic Pole in 2016

julia> using AstroLib

julia> mag2geo(90, 0, 2016)
(86.395, -166.29000000000002)

Notes

This function uses list of North Magnetic Pole positions provided by World Magnetic Model (https://www.ngdc.noaa.gov/geomag/data/poles/NP.xy).

geo2mag converts geographic coordinates to geomagnetic coordinates.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.mean_obliquityMethod
mean_obliquity(jd) -> m_eps

Purpose

Return the mean obliquity of the ecliptic for a given Julian date

Explanation

The function is used by the co_nutate procedure.

Arguments

  • jd: julian date

Output

  • m_eps: mean obliquity of the ecliptic, in radians

Example

julia> using AstroLib

julia> mean_obliquity(jdcnv(1978,01,7,11, 01))
0.4091425159336512

Notes

The algorithm used to find the mean obliquity(eps0) is mentioned in USNO Circular 179, but the canonical reference for the IAU adoption is apparently Hilton et al., 2006, Celest.Mech.Dyn.Astron. 94, 351. 2000

source
AstroLib.month_cnvMethod
month_cnv(number[, shor=true, up=true, low=true]) -> month_name
month_cnv(name) -> number

Purpose

Convert between a month English name and the equivalent number.

Explanation

For example, converts from "January" to 1 or vice-versa.

Arguments

The functions has two methods, one with numeric input (and three possible boolean keywords) and the other one with string input.

Numeric input arguments:

  • number: the number of the month to be converted to month name.
  • short (optional boolean keyword): if true, the abbreviated (3-character) name of the month will be returned, e.g. "Apr" or "Oct". Default is false.
  • up (optional boolean keyword): if true, the name of the month will be all in upper case, e.g. "APRIL" or "OCTOBER". Default is false.
  • low (optional boolean keyword): if true, the name of the month will be all in lower case, e.g. "april" or "october". Default is false.

String input argument:

  • name: month name to be converted to month number.

Output

The month name or month number, depending on the input. For numeric input, the format of the month name is influenced by the optional keywords.

Example

julia> using AstroLib

julia> month_cnv.(["janua", "SEP", "aUgUsT"])
3-element Array{Int64,1}:
 1
 9
 8

julia> month_cnv.([2, 12, 6], short=true, low=true)
3-element Array{String,1}:
 "feb"
 "dec"
 "jun"
source
AstroLib.moonposMethod
moonpos(jd[, radians=true]) -> ra, dec, dis, geolong, geolat

Purpose

Compute the right ascension and declination of the Moon at specified Julian date.

Arguments

  • jd: the Julian ephemeris date. It can be either a scalar or an array
  • radians (optional boolean keyword): if set to true, then all output angular quantities are given in radians rather than degrees. The default is false

Output

The 5-tuple (ra, dec, dis, geolong, geolat):

  • ra: apparent right ascension of the Moon in degrees, referred to the true equator of the specified date(s)
  • dec: the declination of the Moon in degrees
  • dis: the distance between the centre of the Earth and the centre of the Moon in kilometers
  • geolong: apparent longitude of the moon in degrees, referred to the ecliptic of the specified date(s)
  • geolat: apparent longitude of the moon in degrees, referred to the ecliptic of the specified date(s)

If jd is an array, then all output quantities are arrays of the same length as jd.

Method

Derived from the Chapront ELP2000/82 Lunar Theory (Chapront-Touze' and Chapront, 1983, 124, 50), as described by Jean Meeus in Chapter 47 of ``Astronomical Algorithms'' (Willmann-Bell, Richmond), 2nd edition, 1998. Meeus quotes an approximate accuracy of 10" in longitude and 4" in latitude, but he does not give the time range for this accuracy.

Comparison of the IDL procedure with the example in ``Astronomical Algorithms'' reveals a very small discrepancy (~1 km) in the distance computation, but no difference in the position calculation.

Example

(1) Find the position of the moon on April 12, 1992

julia> using AstroLib

julia> jd = jdcnv(1992, 4, 12);

julia> adstring(moonpos(jd)[1:2],precision=1)
" 08 58 45.23  +13 46 06.1"

This is within 1" from the position given in the Astronomical Almanac.

(2) Plot the Earth-moon distance during 2016 with sampling of 6 hours. Use PyPlot.jl for plotting

using PyPlot
points = DateTime(2016):Dates.Hour(6):DateTime(2017);
plot(points, moonpos(jdcnv.(points))[3])

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.mphaseMethod
mphase(jd) -> k

Purpose

Return the illuminated fraction of the Moon at given Julian date(s).

Arguments

  • jd: the Julian ephemeris date.

Output

The illuminated fraction $k$ of Moon's disk, with $0 \leq k \leq 1$. $k = 0$ indicates a new moon, while $k = 1$ stands for a full moon.

Method

Algorithm from Chapter 46 of "Astronomical Algorithms" by Jean Meeus (Willmann-Bell, Richmond) 1991. sunpos and moonpos are used to get positions of the Sun and the Moon, and the Moon distance. The selenocentric elongation of the Earth from the Sun (phase angle) is then computed, and used to determine the illuminated fraction.

Example

Plot the illuminated fraction of the Moon for every day in January 2018 with a hourly sampling. Use PyPlot.jl for plotting

using PyPlot
points = DateTime(2018,01,01):Dates.Hour(1):DateTime(2018,01,31,23,59,59);
plot(points, mphase.(jdcnv.(points)))

Note that in this calendar month there are two full moons, this event is called blue moon.

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.nutateMethod
nutate(jd) -> long, obliq

Purpose

Return the nutation in longitude and obliquity for a given Julian date.

Arguments

  • jd: Julian ephemeris date, it can be either a scalar or a vector

Output

The 2-tuple (long, obliq), where

  • long: the nutation in longitude
  • obl: the nutation in latitude

If jd is an array, long and obl are arrays of the same length.

Method

Uses the formula in Chapter 22 of ``Astronomical Algorithms'' by Jean Meeus (1998, 2nd ed.) which is based on the 1980 IAU Theory of Nutation and includes all terms larger than 0.0003".

Example

(1) Find the nutation in longitude and obliquity 1987 on Apr 10 at Oh. This is example 22.a from Meeus

julia> using AstroLib

julia> jd = jdcnv(1987, 4, 10);

julia> nutate(jd)
(-3.787931077110494, 9.44252069864449)

(2) Plot the daily nutation in longitude and obliquity during the 21st century. Use PyPlot.jl for plotting.

using PyPlot
years = DateTime(2000):DateTime(2100);
long, obl = nutate(jdcnv.(years));
plot(years, long); plot(years, obl)

You can see both the dominant large scale period of nutation, of 18.6 years, and smaller oscillations with shorter periods.

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.ordinalMethod
ordinal(num) -> result

Purpose

Convert an integer to a correct English ordinal string.

Explanation

The first four ordinal strings are "1st", "2nd", "3rd", "4th" ....

Arguments

  • num: number to be made ordinal. It should be of type int.

Output

  • result: ordinal string, such as '1st' '3rd '164th' '87th' etc

Example

julia> using AstroLib

julia> ordinal.(1:5)
5-element Array{String,1}:
 "1st"
 "2nd"
 "3rd"
 "4th"
 "5th"

Notes

This function does not support float arguments, unlike the IDL implementation. Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.paczynskiMethod
paczynski(u) -> amplification

Purpose

Calculate gravitational microlensing amplification of a point-like source by a single point-like lens.

Explanation

Return the gravitational microlensing amplification of a point-like source by a single point-like lens, using Paczyński formula

\[A(u) = \frac{u^2 + 2}{u\sqrt{u^2 + 4}}\]

where $u$ is the projected distance between the lens and the source in units of Einstein radii.

In order to speed up calculations for extreme values of $u$, the following asyntotic expressions for $A(u)$ are used:

\[A(u) = \begin{cases} 1/u & |u| \ll 1 \\ \text{sgn}(u) & |u| \gg 1 \end{cases}\]

Arguments

  • u: projected distance between the lens and the source, in units of Einstein radii

Output

The microlensing amplification for the given distance.

Example

Calculate the microlensing amplification for $u = 10^{-10}$, $10^{-1}$, $1$, $10$, $10^{10}$:

julia> using AstroLib

julia> paczynski.([1e-10, 1e-1, 1, 10, 1e10])
5-element Array{Float64,1}:
  1.0e10
 10.037461005722337
  1.3416407864998738
  1.0001922892047386
  1.0

Notes

The expression of $A(u)$ of microlensing amplification has been given by Bohdan Paczyński in

The same expression was actually found by Albert Einstein half a century earlier:

source
AstroLib.planck_freqMethod
planck_freq(frequency, temperature) -> black_body_flux

Purpose

Calculate the flux of a black body per unit frequency.

Explanation

Return the spectral radiance of a black body per unit frequency using Planck's law

$B_\nu(\nu, T) = \frac{2h\nu ^3}{c^2} \frac{1}{e^\frac{h\nu}{k_\mathrm{B}T} - 1}$

Arguments

  • frequency: frequency at which the flux is to be calculated, in Hertz.
  • temperature: the equilibrium temperature of the black body, in Kelvin.

Output

The spectral radiance of the black body, in units of W/(sr·m²·Hz).

Example

Plot the spectrum of a black body in $[10^{12}, 10^{15.4}]$ Hz at $8000$ K. Use PyPlot.jl for plotting.

using PyPlot
frequency = exp10.(range(12, stop=15.4, length=1000));
temperature = ones(frequency)*8000;
flux = planck_freq.(frequency, temperature);
plot(frequency, flux)

Notes

planck_wave calculates the flux of a black body per unit wavelength.

source
AstroLib.planck_waveMethod
planck_wave(wavelength, temperature) -> black_body_flux

Purpose

Calculate the flux of a black body per unit wavelength.

Explanation

Return the spectral radiance of a black body per unit wavelength using Planck's law

$B_\lambda(\lambda, T) =\frac{2hc^2}{\lambda^5}\frac{1}{e^{\frac{hc}{\lambda k_\mathrm{B}T}} - 1}$

Arguments

  • wavelength: wavelength at which the flux is to be calculated, in meters.
  • temperature: the equilibrium temperature of the black body, in Kelvin.

Output

The spectral radiance of the black body, in units of W/(sr·m³).

Example

Plot the spectrum of a black body in $[0, 3]$ µm at $5000$ K. Use PyPlot.jl for plotting.

using PyPlot
wavelength = range(0, stop=3e-6, length=1000);
temperature = ones(wavelength)*5000;
flux = planck_wave.(wavelength, temperature);
plot(wavelength, flux)

Notes

planck_freq calculates the flux of a black body per unit frequency.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.planet_coordsMethod
planet_coords(date, num)

Purpose

Find right ascension and declination for the planets when provided a date as input.

Explanation

This function uses the helio to get the heliocentric ecliptic coordinates of the planets at the given date which it then converts these to geocentric ecliptic coordinates ala "Astronomical Algorithms" by Jean Meeus (1991, p 209). These are then converted to right ascension and declination using euler.

The accuracy between the years 1800 and 2050 is better than 1 arcminute for the terrestial planets, but reaches 10 arcminutes for Saturn. Before 1850 or after 2050 the accuracy can get much worse.

Arguments

  • date: Can be either a single date or an array of dates. Each element can be either a DateTime type or Julian Date. It can be a scalar or vector.
  • num: integer denoting planet number, scalar or vector 1 = Mercury, 2 = Venus, ... 9 = Pluto. If not in that change, then the program will throw an error.

Output

  • ra: right ascension of planet(J2000), in degrees
  • dec: declination of the planet(J2000), in degrees

Example

Find the RA, Dec of Venus on 1992 Dec 20

julia> using AstroLib, Dates

julia> adstring(planet_coords(DateTime(1992,12,20),2))
" 21 05 02.8  -18 51 41"

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.polrecMethod
polrec(radius, angle[, degrees=true]) -> x, y

Purpose

Convert 2D polar coordinates to rectangular coordinates.

Explanation

This is the partial inverse function of recpol.

Arguments

  • radius: radial coordinate of the point. It may be a scalar or an array.
  • angle: the angular coordinate of the point. It may be a scalar or an array of the same lenth as radius.
  • degrees (optional boolean keyword): if true, the angle is assumed to be in degrees, otherwise in radians. It defaults to false.

Mandatory arguments can also be passed as the 2-tuple (radius, angle), so that it is possible to execute recpol(polrec(radius, angle)).

Output

A 2-tuple (x, y) with the rectangular coordinate of the input. If radius and angle are arrays, x and y are arrays of the same length as radius and angle.

Example

Get rectangular coordinates $(x, y)$ of the point with polar coordinates $(r, \varphi) = (1.7, 227)$, with angle $\varphi$ expressed in degrees.

julia> using AstroLib

julia> x, y = polrec(1.7, 227, degrees=true)
(-1.1593972121062475, -1.2433012927525897)
source
AstroLib.posangMethod
posang(units, ra1, dec1, ra2, dec2) -> angular_distance

Purpose

Compute rigorous position angle of point 2 relative to point 1.

Explanation

Computes the rigorous position angle of point 2 (with given right ascension and declination) using point 1 (with given right ascension and declination) as the center.

Arguments

  • units: integer, can be either 0, or 1, or 2. Describes units of inputs and

output: * 0: everything (input right ascensions and declinations, and output distance) is radians * 1: right ascensions are in decimal hours, declinations in decimal degrees, output distance in degrees * 2: right ascensions and declinations are in degrees, output distance in degrees

  • ra1: right ascension or longitude of point 1
  • dec1: declination or latitude of point 1
  • ra2: right ascension or longitude of point 2
  • dec2: declination or latitude of point 2

Both ra1 and dec1, and ra2 and dec2 can be given as 2-tuples (ra1, dec1) and (ra2, dec2).

Output

Angle of the great circle containing [ra2, dec2] from the meridian containing [ra1, dec1], in the sense north through east rotating about [ra1, dec1]. See units argument above for units.

Method

The "four-parts formula" from spherical trigonometry (p. 12 of Smart's Spherical Astronomy or p. 12 of Green' Spherical Astronomy).

Example

Mizar has coordinates (ra, dec) = (13h 23m 55.5s, +54° 55' 31''). Its companion, Alcor, has coordinates (ra, dec) = (13h 25m 13.5s, +54° 59' 17''). Find the position angle of Alcor with respect to Mizar.

julia> using AstroLib

julia> posang(1, ten(13, 25, 13.5), ten(54, 59, 17), ten(13, 23, 55.5), ten(54, 55, 31))
-108.46011246802047

Notes

  • The function sphdist provides an alternate method of computing a spherical

distance.

  • Note that posang is not commutative: the position angle between A and B is $\theta$, then the position angle between B and A is $180 + \theta$.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.precessMethod
precess(ra, dec, equinox1, equinox2[, FK4=true, radians=true]) -> prec_ra, prec_dec

Purpose

Precess coordinates from equinox1 to equinox2.

Explanation

The default (ra, dec) system is FK5 based on epoch J2000.0 but FK4 based on B1950.0 is available via the FK4 boolean keyword.

Arguments

  • ra: input right ascension, scalar or vector, in degrees, unless the radians keyword is set to true
  • dec: input declination, scalar or vector, in degrees, unless the radians keyword is set to true
  • equinox1: original equinox of coordinates, numeric scalar.
  • equinox2: equinox of precessed coordinates.
  • FK4 (optional boolean keyword): if this keyword is set to true, the FK4 (B1950.0) system precession angles are used to compute the precession matrix. When it is false, the default, use FK5 (J2000.0) precession angles.
  • radians (optional boolean keyword): if this keyword is set to true, then the input and output right ascension and declination vectors are in radians rather than degrees.

Output

The 2-tuple (ra, dec) of coordinates modified by precession.

Example

The Pole Star has J2000.0 coordinates (2h, 31m, 46.3s, 89d 15' 50.6"); compute its coordinates at J1985.0

julia> using AstroLib

julia> ra, dec = ten(2,31,46.3)*15, ten(89,15,50.6)
(37.94291666666666, 89.26405555555556)

julia> adstring(precess(ra, dec, 2000, 1985), precision=1)
" 02 16 22.73  +89 11 47.3"

Precess the B1950 coordinates of Eps Ind (RA = 21h 59m,33.053s, DEC = (-56d, 59', 33.053") to equinox B1975.

julia> using AstroLib

julia> ra, dec = ten(21, 59, 33.053) * 15, ten(-56, 59, 33.053)
(329.88772083333333, -56.992514722222225)

julia> adstring(precess(ra, dec, 1950, 1975, FK4=true), precision=1)
" 22 01 15.46  -56 52 18.7"

Method

Algorithm from "Computational Spherical Astronomy" by Taff (1983), p. 24. (FK4). FK5 constants from "Explanatory Supplement To The Astronomical Almanac" 1992, page 104 Table 3.211.1 (https://archive.org/details/131123ExplanatorySupplementAstronomicalAlmanac).

Notes

Accuracy of precession decreases for declination values near 90 degrees. precess should not be used more than 2.5 centuries from 2000 on the FK5 system (1950.0 on the FK4 system). If you need better accuracy, use bprecess or jprecess as needed.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.precess_cdFunction
precess_cd(cd, epoch1, epoch2, crval_old, crval_new[, FK4=true]) -> cd

Purpose

Precess the coordinate description matrix.

Explanation

The coordinate matrix is precessed from epoch1 to epoch2.

Arguments

  • cd: 2 x 2 coordinate description matrix in degrees
  • epoch1: original equinox of coordinates, scalar
  • epoch2: equinox of precessed coordinates, scalar
  • crval_old: 2 element vector containing right ascension and declination in degrees of the reference pixel in the original equinox
  • crval_new: 2 element vector giving crval in the new equinox
  • FK4 (optional boolean keyword): if this keyword is set to true, then the precession constants are taken in the FK4 reference frame. When it is false, the default is the FK5 frame

Output

  • cd: coordinate description containing precessed values

Example

julia> using AstroLib

julia> precess_cd([20 60; 45 45], 1950, 2000, [34, 58], [12, 83])
2×2 Array{Float64,2}:
  48.8944  147.075
 110.188   110.365

Notes

Code of this function is based on IDL Astronomy User's Library. This function should not be used for values more than 2.5 centuries from the year 1900. This function calls sec2rad, precess and bprecess.

source
AstroLib.precess_xyzMethod
precess_xyz(x, y, z, equinox1, equinox2) -> prec_x, prec_y, prec_z

Purpose

Precess equatorial geocentric rectangular coordinates.

Arguments

  • x, y, z: scalars or vectors giving heliocentric rectangular coordinates.
  • equinox1: original equinox of coordinates, numeric scalar.
  • equinox2: equinox of precessed coordinates, numeric scalar.

Input coordinates can be given also a 3-tuple (x, y, z).

Output

The 3-tuple (x, y, z) of coordinates modified by precession.

Example

Precess 2000 equinox coordinates (1, 1, 1) to 2050.

julia> using AstroLib

julia> precess_xyz(1, 1, 1, 2000, 2050)
(0.9838854500981734, 1.0110925876508692, 1.0048189888146941)

Method

The equatorial geocentric rectangular coordinates are converted to right ascension and declination, precessed in the normal way, then changed back to x, y and z using unit vectors.

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.prematMethod
premat(equinox1, equinox2[, FK4=true]) -> precession_matrix

Purpose

Return the precession matrix needed to go from equinox1 to equinox2.

Explanation

This matrix is used by precess and baryvel to precess astronomical coordinates.

Arguments

  • equinox1: original equinox of coordinates.
  • equinox2: equinox of precessed coordinates.
  • FK4 (optional boolean keyword): if this keyword is set to true, the FK4 (B1950.0) system precession angles are used to compute the precession matrix. When it is false, the default, use FK5 (J2000.0) precession angles.

Output

A 3×3 matrix, used to precess equatorial rectangular coordinates.

Example

Return the precession matrix from 1950.0 to 1975.0 in the FK4 system

julia> using AstroLib

julia> premat(1950,1975,FK4=true)
3×3 StaticArrays.SArray{Tuple{3,3},Float64,2,9} with indices SOneTo(3)×SOneTo(3):
 0.999981    -0.00558775  -0.00242909
 0.00558775   0.999984    -6.78691e-6
 0.00242909  -6.78633e-6   0.999997

Method

FK4 constants from "Computational Spherical Astronomy" by Taff (1983), p. 24. (FK4). FK5 constants from "Explanatory Supplement To The Astronomical Almanac" 1992, page 104 Table 3.211.1 (https://archive.org/details/131123ExplanatorySupplementAstronomicalAlmanac).

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.rad2secMethod
rad2sec(rad) -> seconds

Purpose

Convert from radians to seconds.

Argument

  • rad: number of radians.

Output

The number of seconds corresponding to rad.

Example

julia> using AstroLib

julia> rad2sec(1)
206264.80624709636

Notes

Use sec2rad to convert seconds to radians.

source
AstroLib.radecMethod
radec(ra::Real, dec::Real[, hours=true]) -> ra_hours, ra_minutes, ra_seconds, dec_degrees, dec_minutes, dec_seconds

Purpose

Convert right ascension and declination from decimal to sexagesimal units.

Explanation

The conversion is to sexagesimal hours for right ascension, and sexagesimal degrees for declination.

Arguments

  • ra: decimal right ascension, scalar or array. It is expressed in degrees, unless the optional keyword hours is set to true.
  • dec: declination in decimal degrees, scalar or array, same number of elements as ra.
  • hours (optional boolean keyword): if false (the default), ra is assumed to be given in degrees, otherwise ra is assumed to be expressed in hours.

Output

A 6-tuple of AbstractFloat:

(ra_hours, ra_minutes, ra_seconds, dec_degrees, dec_minutes, dec_seconds)

If ra and dec are arrays, also each element of the output 6-tuple are arrays of the same dimension.

Example

Position of Sirius in the sky is (ra, dec) = (6.7525, -16.7161), with right ascension expressed in hours. Its sexagesimal representation is given by

julia> using AstroLib

julia> radec(6.7525, -16.7161, hours=true)
(6.0, 45.0, 9.0, -16.0, 42.0, 57.9600000000064)
source
AstroLib.recpolMethod
recpol(x, y[, degrees=true]) -> radius, angle

Purpose

Convert 2D rectangular coordinates to polar coordinates.

Explanation

This is the partial inverse function of polrec.

Arguments

  • x: the abscissa coordinate of the point. It may be a scalar or an array.
  • y: the ordinate coordinate of the point. It may be a scalar or an array of the same lenth as x.
  • degrees (optional boolean keyword): if true, the output angle is given

in degrees, otherwise in radians. It defaults to false.

Mandatory arguments may also be passed as the 2-tuple (x, y), so that it is possible to execute polrec(recpol(x, y)).

Output

A 2-tuple (radius, angle) with the polar coordinates of the input. The coordinate angle coordinate lies in the range $[-\pi, \pi]$ if degrees=false, or $[-180, 180]$ when degrees=true.

If x and y are arrays, radius and angle are arrays of the same length as radius and angle.

Example

Calculate polar coordinates $(r, \varphi)$ of point with rectangular coordinates $(x, y) = (2.24, -1.87)$.

julia> using AstroLib

julia> r, phi = recpol(2.24, -1.87)
(2.9179616172938263, -0.6956158538564537)

Angle $\varphi$ is given in radians.

source
AstroLib.rhothetaMethod
rhotheta(period, periastron, eccentricity, semimajor_axis, inclination, omega, omega2, epoch) -> rho, theta

Purpose

Calculate the separation and position angle of a binary star.

Explanation

This function will return the separation $\rho$ and position angle $\theta$ of a visual binary star derived from its orbital elements. The algorithms described in the following book will be used: Meeus J., 1992, Astronomische Algorithmen, Barth. Compared to the examples given at page 400 and no discrepancy found.

Arguments

  • period: period [year]
  • periastro: time of periastron passage [year]
  • eccentricity: eccentricity of the orbit
  • semimajor_axis: semi-major axis [arc second]
  • inclination: inclination angle [degree]
  • omega: node [degree]
  • omega2: longitude of periastron [degree]
  • epoch: epoch of observation [year]

All input parameters have to be scalars.

Output

The 2-tuple $(\rho, \theta)$, where

  • $\rho$ is separation [arc second], and
  • $\theta$ is position angle (degree).

Example

Find the position of Eta Coronae Borealis at the epoch 2016

julia> using AstroLib

julia> ρ, θ = rhotheta(41.623, 1934.008, 0.2763, 0.907, 59.025, 23.717, 219.907, 2016)
(0.6351167848659552, 214.42513387396497)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.sec2radMethod
sec2rad(sec) -> radians

Purpose

Convert from seconds to radians.

Argument

  • sec: number of seconds.

Output

The number of radians corresponding to sec.

Example

julia> using AstroLib

julia> sec2rad(3600 * 30)
0.5235987755982988

Notes

Use rad2sec to convert radians to seconds.

source
AstroLib.sixtyMethod
sixty(number) -> [deg, min, sec]

Purpose

Converts a decimal number to sexagesimal.

Explanation

The reverse of ten function.

Argument

  • number: decimal number to be converted to sexagesimal.

Output

An array of three AbstractFloat, that are the sexagesimal counterpart (degrees, minutes, seconds) of number.

Example

julia> using AstroLib

julia> sixty(-0.615)
3-element StaticArrays.SArray{Tuple{3},Float64,1,3} with indices SOneTo(3):
 -0.0
 36.0
 54.0

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.sphdistMethod
sphdist(long1, lat1, long2, lat2[, degrees=true]) -> angular_distance

Purpose

Angular distance between points on a sphere.

Arguments

  • long1: longitude of point 1
  • lat1: latitude of point 1
  • long2: longitude of point 2
  • lat2: latitude of point 2
  • degrees (optional boolean keyword): if true, all angles, including the output distance, are assumed to be in degrees, otherwise they are all in radians. It defaults to false.

Output

Angular distance on a sphere between points 1 and 2, as an AbstractFloat. It is expressed in radians unless degrees keyword is set to true.

Example

julia> using AstroLib

julia> sphdist(120, -43, 175, +22)
1.5904422616007134

Notes

  • gcirc function is similar to sphdist, but may be more suitable for astronomical applications.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.sunposMethod
sunpos(jd[, radians=true]) -> ra, dec, elong, obliquity

Purpose

Compute the right ascension and declination of the Sun at a given date.

Arguments

  • jd: the Julian date of when you want to calculate Sun position. It can be either a scalar or a vector. Use jdcnv to get the Julian date for a given date and time.
  • radians (optional boolean keyword): if set to true, all output quantities are given in radians. The default is false, so all quantities are given in degrees.

Output

The 4-tuple (ra, dec, elong, obliquity):

  • ra: the right ascension of the Sun at that date
  • dec: the declination of the Sun at that date
  • elong: ecliptic longitude of the Sun at that date
  • obliquity: the obliquity of the ecliptic

All quantities are given in degrees, unless radians keyword is set to true (see "Arguments" section). If jd is an array, arrays of the same given as jd are returned.

Method

Uses a truncated version of Newcomb's Sun. Adapted from the IDL routine SUN_POS by CD Pike, which was adapted from a FORTRAN routine by B. Emerson (RGO).

Example

(1) Find the apparent right ascension and declination of the Sun on May 1, 1982

julia> using AstroLib

julia> adstring(sunpos(jdcnv(1982, 5, 1))[1:2], precision=2)
" 02 31 32.614  +14 54 34.92"

The Astronomical Almanac gives 02 31 32.58 +14 54 34.9 so the error for this case is < 0.5".

(2) Plot the apparent right ascension, in hours, and declination of the Sun, in degrees, for every day in 2016. Use PyPlot.jl for plotting.

using PyPlot
using Dates

days = DateTime(2016):Day(1):DateTime(2016, 12, 31);
ra, declin = sunpos(jdcnv.(days));
plot(days, ra/15); plot(days, declin)

Notes

Patrick Wallace (Rutherford Appleton Laboratory, UK) has tested the accuracy of a C adaptation of the present algorithm and found the following results. From 1900-2100 sunpos gave 7.3 arcsec maximum error, 2.6 arcsec RMS. Over the shorter interval 1950-2050 the figures were 6.4 arcsec max, 2.2 arcsec RMS.

The returned ra and dec are in the given date's equinox.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.tenFunction
ten(deg[, min, sec]) -> decimal
ten("deg:min:sec") -> decimal

Purpose

Converts a sexagesimal number or string to decimal.

Explanation

ten is the inverse of the sixty function.

Arguments

ten takes as argument either three scalars (deg, min, sec) or a string. The string should have the form "deg:min:sec" or "deg min sec". Also any iterable like (deg, min, sec) or [deg, min, sec] is accepted as argument.

If minutes and seconds are not specified they default to zero.

Output

The decimal conversion of the sexagesimal numbers provided is returned.

Method

The formula used for the conversion is

\[\mathrm{sign}(\mathrm{deg})·\left(|\mathrm{deg}| + \frac{\mathrm{min}}{60} + \frac{\mathrm{sec}}{3600}\right)\]

Example

julia> using AstroLib

julia> ten(-0.0, 19, 47)
-0.3297222222222222

julia> ten("+5:14:58")
5.249444444444444

julia> ten("-10 26")
-10.433333333333334

julia> ten((-10, 26))
-10.433333333333334

Notes

These functions cannot deal with -0 (negative integer zero) in numeric input. If it is important to give sense to negative zero, you can either make sure to pass a floating point negative zero -0.0 (this is the best option), or use negative minutes and seconds, or non-integer negative degrees and minutes.

source
AstroLib.tic_oneFunction
tic_one(zmin, pixx, incr[, ra=true]) -> min2, tic1

Purpose

Determine the position of the first tic mark for astronomical images.

Explanation

For use in labelling images with right ascension and declination axes. This routine determines the position in pixels of the first tic.

Arguments

  • zmin: astronomical coordinate value at axis zero point (degrees or hours).
  • pixx: distance in pixels between tic marks (usually obtained from tics).
  • incr - increment in minutes for labels (usually an even number obtained from the procedure tics).
  • ra (optional boolean keyword): if true, incremental value being entered is in minutes of time, else it is assumed that value is in else it's in minutes of arc. Default is false.

Output

The 2 tuple (min2, tic1):

  • min2: astronomical coordinate value at first tic mark
  • tic1: position in pixels of first tic mark

Example

Suppose a declination axis has a value of 30.2345 degrees at its zero point. A tic mark is desired every 10 arc minutes, which corresponds to 12.74 pixels, with increment for labels being 10 minutes. Then

julia> using AstroLib

julia> tic_one(30.2345, 12.74, 10)
(30.333333333333332, 7.554820000000081)

yields values of min2 ≈ 30.333 and tic1 ≈ 7.55482, i.e. the first tic mark should be labeled 30 deg 20 minutes and be placed at pixel value 7.55482.

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.ticposMethod
ticpos(deglen, pixlen, ticsize) -> ticsize, incr, units

Purpose

Specify distance between tic marks for astronomical coordinate overlays.

Explanation

User inputs number an approximate distance between tic marks, and the axis length in degrees. ticpos will return a distance between tic marks such that the separation is a round multiple in arc seconds, arc minutes, or degrees.

Arguments

  • deglen: length of axis in degrees, positive scalar
  • pixlen: length of axis in plotting units (pixels), postive scalar
  • ticsize: distance between tic marks (pixels). This value will be adjusted by ticpos such that the distance corresponds to a round multiple in the astronomical coordinate.

Output

The 3-tuple (ticsize, incr, units):

  • ticsize: distance between tic marks (pixels), positive scalar
  • incr: incremental value for tic marks in round units given by the units parameter
  • units: string giving units of ticsize, either 'Arc Seconds', 'Arc Minutes', or 'Degrees'

Example

Suppose a 512 x 512 image array corresponds to 0.2 x 0.2 degrees on the sky. A tic mark is desired in round angular units, approximately every 75 pixels. Then

julia> using AstroLib

julia> ticpos(0.2, 512, 75)
(85.33333333333333, 2.0, "Arc Minutes")

i.e. a good tic mark spacing is every 2 arc minutes, corresponding to 85.333 pixels.

Notes

All the arguments taken as input are assumed to be positive in nature.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.ticsFunction
tics(radec_min, radec_max, numx, ticsize[, ra=true]) -> ticsize, incr

Purpose

Compute a nice increment between tic marks for astronomical images.

Explanation

For use in labelling a displayed image with right ascension or declination axes. An approximate distance between tic marks is input, and a new value is computed such that the distance between tic marks is in simple increments of the tic label values.

Arguements

  • radec_min : minimum axis value (degrees).
  • radec_min : maximum axis value (degrees).
  • numx : number of pixels in x direction.
  • ticsize : distance between tic marks (pixels).
  • ra (optional boolean keyword): if true, incremental value would be in minutes of time. Default is false.

Output

A 2-tuple (ticsize, incr):

  • ticsize : distance between tic marks (pixels).
  • incr : incremental value for tic labels. The format is dependent on the optional keyword. If true (i.e for right ascension), it's in minutes of time, else it's in minutes of arc (for declination).

Example

julia> using AstroLib

julia> tics(55, 60, 100.0, 1/2)
(0.66, 2.0)

julia> tics(30, 60, 12, 2, true)
(2.75, 30.0)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.true_obliquityMethod
true_obliquity(jd) -> t_eps

Purpose

Return the true obliquity of the ecliptic for a given Julian date

Explanation

The function is used by the co_aberration procedure.

Arguments

  • jd: Julian date.

Output

  • t_eps: true obliquity of the ecliptic, in radians

Example

julia> using AstroLib

julia> true_obliquity(jdcnv(1978,01,7,11, 01))
0.4090953896211926

Notes

The function calls mean_obliquity.

source
AstroLib.trueanomMethod
trueanom(E, e) -> true anomaly

Purpose

Calculate true anomaly for a particle in elliptic orbit with eccentric anomaly $E$ and eccentricity $e$.

Explanation

In the two-body problem, once that the Kepler's equation is solved and $E(t)$ is determined, the polar coordinates $(r(t), \theta(t))$ of the body at time $t$ in the elliptic orbit are given by

$\theta(t) = 2\arctan \left(\sqrt{\frac{1 + e}{1 - e}} \tan\frac{E(t)}{2} \right)$

$r(t) = \frac{a(1 - e^{2})}{1 + e\cos(\theta(t) - \theta_{0})}$

in which $a$ is the semi-major axis of the orbit, and $\theta_0$ the value of angular coordinate at time $t = t_{0}$.

Arguments

  • E: eccentric anomaly.
  • e: eccentricity, in the elliptic motion regime ($0 \leq e \leq 1$)

Output

The true anomaly.

Example

Plot the true anomaly as a function of mean anomaly for eccentricity $e = 0$, $0.5$, $0.9$. Use PyPlot.jl for plotting.

using PyPlot
M = range(0, stop=2pi, length=1001)[1:end-1];
for ecc in (0, 0.5, 0.9)
    plot(M, mod2pi.(trueanom.(kepler_solver.(M, ecc), ecc)))
end

Notes

The eccentric anomaly can be calculated with kepler_solver function.

source
AstroLib.uvbybetaFunction
uvbybeta(by, m1, c1, n[, hbeta=NaN, eby_in=NaN]) -> te, mv, eby, delm0, radius

Purpose

Derive dereddened colors, metallicity, and Teff from Stromgren colors.

Arguments

  • by: Stromgren b-y color
  • m1: Stromgren line-blanketing parameter
  • c1: Stromgren Balmer discontinuity parameter
  • n: Integer which can be any value between 1 to 8, giving approximate stellar classification. (1) B0 - A0, classes III - V, 2.59 < Hbeta < 2.88,-0.20 < c0 < 1.00 (2) B0 - A0, class Ia , 2.52 < Hbeta < 2.59,-0.15 < c0 < 0.40 (3) B0 - A0, class Ib , 2.56 < Hbeta < 2.61,-0.10 < c0 < 0.50 (4) B0 - A0, class II , 2.58 < Hbeta < 2.63,-0.10 < c0 < 0.10 (5) A0 - A3, classes III - V, 2.87 < Hbeta < 2.93,-0.01 < (b-y)o < 0.06 (6) A3 - F0, classes III - V, 2.72 < Hbeta < 2.88, 0.05 < (b-y)o < 0.22 (7) F1 - G2, classes III - V, 2.60 < Hbeta < 2.72, 0.22 < (b-y)o < 0.39 (8) G2 - M2, classes IV - V, 0.20 < m0 < 0.76, 0.39 < (b-y)o < 1.00
  • hbeta (optional): H-beta line strength index. If it is not supplied, then by default its value will be NaN and the code will estimate a value based on by, m1,and c1. It is not used for stars in group 8.
  • eby_in (optional): specifies the E(b-y) color to use. If not supplied, then by default its value will be NaN and E(b-y) will be estimated from the Stromgren colors.

Output

  • te: approximate effective temperature
  • mv: absolute visible magnitude
  • eby: color excess E(b-y)
  • delm0: metallicity index, delta m0, may not be calculable for early B stars and so returns NaN.
  • radius: stellar radius (R/R(solar))

Example

Suppose 5 stars have the following Stromgren parameters

by = [-0.001 ,0.403, 0.244, 0.216, 0.394] m1 = [0.105, -0.074, -0.053, 0.167, 0.186] c1 = [0.647, 0.215, 0.051, 0.785, 0.362] hbeta = [2.75, 2.552, 2.568, 2.743, 0] nn = [1,2,3,7,8]

Determine the stellar parameters

julia> using AstroLib

julia> by = [-0.001 ,0.403, 0.244, 0.216, 0.394];

julia> m1 = [0.105, -0.074, -0.053, 0.167, 0.186];

julia> c1 = [0.647, 0.215, 0.051, 0.785, 0.362];

julia> hbeta = [2.75, 2.552, 2.568, 2.743, 0];

julia> nn = [1,2,3,7,8];

julia> uvbybeta.(by, m1, c1, nn, hbeta)
5-element Array{NTuple{5,Float64},1}:
 (13057.535222326893, -0.27375469585031265, 0.04954396423248884, -0.008292894218734928, 2.7136529525371897)
 (14025.053834219656, -6.907050783073221, 0.4140562248995983, NaN, 73.50771722263974)
 (18423.76405400214, -5.935816553877892, 0.2828247876690783, NaN, 39.84106215808709)
 (7210.507090112837, 2.2180408083364167, 0.018404079180028038, 0.018750927360588615, 2.0459018065648165)
 (5755.671513413262, 3.9449408311022, -0.025062997393370458, 0.03241423718769865, 1.5339239690774464)

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.vactoairMethod
vactoair(wave_vacuum) -> wave_air

Purpose

Converts vacuum wavelengths to air wavelengths.

Explanation

Corrects for the index of refraction of air under standard conditions. Wavelength values below $2000 Å$ will not be altered. Uses relation of Ciddor (1996).

Arguments

  • wave_vacuum: vacuum wavelength in angstroms. Wavelengths are corrected for the index of refraction of air under standard conditions. Wavelength values below $2000 Å$ will not be altered, take care within $[1 Å, 2000 Å]$.

Output

Air wavelength in angstroms.

Method

Uses relation of Ciddor (1996), Applied Optics 35, 1566 (http://adsabs.harvard.edu/abs/1996ApOpt..35.1566C).

Example

If the vacuum wavelength is w = 2000, then vactoair(w) yields an air wavelength of 1999.353.

julia> using AstroLib

julia> vactoair(2000)
1999.3526230448367

Notes

airtovac converts air wavelengths to vacuum wavelengths.

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.xyzFunction
xyz(jd[, equinox]) -> x, y, z, v_x, v_y, v_z

Purpose

Calculate geocentric $x$, $y$, and $z$ and velocity coordinates of the Sun.

Explanation

Calculates geocentric $x$, $y$, and $z$ vectors and velocity coordinates ($dx$, $dy$ and $dz$) of the Sun. (The positive $x$ axis is directed towards the equinox, the $y$-axis, towards the point on the equator at right ascension 6h, and the $z$ axis toward the north pole of the equator). Typical position accuracy is $<10^{-4}$ AU (15000 km).

Arguments

  • jd: number of Reduced Julian Days for the wanted date. It can be either a scalar or a vector.
  • equinox (optional numeric argument): equinox of output. Default is 1950.

You can use juldate to get the number of Reduced Julian Days for the selected dates.

Output

The 6-tuple $(x, y, z, v_x, v_y, v_z)$, where

  • $x, y, z$: scalars or vectors giving heliocentric rectangular coordinates (in AU) for each date supplied. Note that $\sqrt{x^2 + y^2 + z^2}$ gives the Earth-Sun distance for the given date.
  • $v_x, v_y, v_z$: velocity vectors corresponding to $x, y$, and $z$.

Example

What were the rectangular coordinates and velocities of the Sun on 1999-01-22T00:00:00 (= JD 2451200.5) in J2000 coords? Note: Astronomical Almanac (AA) is in TDT, so add 64 seconds to UT to convert.

julia> using AstroLib, Dates

julia> jd = juldate(DateTime(1999, 1, 22))
51200.5

julia> xyz(jd + 64/86400, 2000)
(0.514568709240398, -0.7696326261820209, -0.33376880143023935, 0.014947267514079971, 0.008314838205477328, 0.003606857607575486)

Compare to Astronomical Almanac (1999 page C20)

            x  (AU)        y  (AU)     z (AU)
xyz:      0.51456871   -0.76963263  -0.33376880
AA:       0.51453130   -0.7697110   -0.3337152
abs(err): 0.00003739    0.00007839   0.00005360
abs(err)
    (km):   5609          11759         8040

NOTE: Velocities in AA are for Earth/Moon barycenter (a very minor offset) see AA 1999 page E3

           x vel (AU/day) y vel (AU/day)   z vel (AU/day)
xyz:      -0.014947268   -0.0083148382    -0.0036068576
AA:       -0.01494574    -0.00831185      -0.00360365
abs(err):  0.000001583    0.0000029886     0.0000032076
abs(err)
 (km/sec): 0.00265        0.00519          0.00557

Notes

Code of this function is based on IDL Astronomy User's Library.

source
AstroLib.ydn2mdMethod
ydn2md(year, day) -> date

Purpose

Convert from year and day number of year to a date.

Explanation

Returns the date corresponding to the day of year.

Arguments

  • year: the year, as an integer.
  • day: the day of year, as an integer.

Output

The date, of Date type, of $\text{day} - 1$ days after January 1st of year.

Example

Find the date of the 60th and 234th days of the year 2016.

julia> using AstroLib

julia> ydn2md.(2016, [60, 234])
2-element Array{Dates.Date,1}:
 2016-02-29
 2016-08-21

Note

ymd2dn converts from a date to day of the year.

source
AstroLib.ymd2dnFunction
ymd2dn(date) -> number_of_days

Purpose

Convert from a date to day of the year.

Explanation

Returns the day of the year for date with January 1st being day 1.

Arguments

  • date: the date with Date type. Can be a single date or an array of dates.

Output

The day of the year for the given date. If date is an array, returns an array of days.

Example

Find the days of the year for March 5 in the years 2015 and 2016 (this is a leap year).

julia> using AstroLib, Dates

julia> ymd2dn.([Date(2015, 3, 5), Date(2016, 3, 5)])
2-element Array{Int64,1}:
 64
 65

Note

ydn2md converts from year and day number of year to a date.

source
AstroLib.zenposFunction
zenpos(jd, latitude, longitude) -> zenith_right_ascension, declination
zenpos(date, latitude, longitude, tz) -> zenith_right_ascension, declination

Purpose

Return the zenith right ascension and declination in radians for a given Julian date or a local civil time and timezone.

Explanation

The local sidereal time is computed with the help of ct2lst, which is the right ascension of the zenith. This and the observatories latitude (corresponding to the declination) are converted to radians and returned as the zenith direction.

Arguments

The function can be called in two different ways. The arguments common to both methods are latitude and longitude:

  • latitude : latitude of the desired location.
  • longitude : longitude of the desired location.

The zenith direction can be computed either by providing the Julian date:

  • jd : the Julian date of the date and time for which the zenith position is desired.

or the time zone and the date:

  • tz: the time zone (in hours) of the desired location (e.g. 4 = EDT, 5 = EST)
  • date: the local civil time with type DateTime.

Output

A 2-tuple (ra, dec):

  • ra : the right ascension (in radians) of the zenith.
  • dec : the declination (in radians) of the zenith.

Example

julia> using AstroLib, Dates

julia> zenpos(DateTime(2017, 04, 25, 18, 59), 43.16, -24.32, 4)
(0.946790432684706, 0.7532841051607526)

julia> zenpos(jdcnv(2016, 05, 05, 13, 41), ten(35,0,42), ten(135,46,6))
(3.5757821152779536, 0.6110688599440813)

Notes

Code of this function is based on IDL Astronomy User's Library.

source