Observed and apparent altitude

PyAstronomy.pyasl.co_refract_forward(alt, pressure=1010.0, temperature=10.0)

Converts the observed into the apparent (real) altitude.

The observed altitude is the altitude that a star is seen to be with a telescope. This is where it appears in the sky. The observed altitude is always greater than the the apparent altitude, which is the altitude that a star would be at, if there were no atmosphere (sometimes called “true” altitude).

Parameters
altfloat or array

Observed altitude of an object in DEGREES.

pressurefloat or array, optional

Atmospheric pressure in MILLIBAR. Default pressure is 1010 mbar. If a single value is given, it will be used for all given altitudes.

temperaturefloat or array, optional

Ground temperature in degrees Celsius. Default temperature is 10 Celsius. If a single value is given, it will be used for all given altitudes.

Returns
Altitude correctionarray

An array holding the altitude correction [deg]. To convert observed altitude into apparent (real) altitude, the correction needs to be subtracted from the observed altitude.

Notes

Note

This function was ported from the IDL Astronomy User’s Library.

IDL - Documentation

function co_refract_forward, a, P=P, T=T

INPUTS
a = The observed (apparent) altitude, in DEGREES.

May be scalar or vector.

INPUT KEYWORDS

P: Pressure [in millibars]. Default is 1010 millibars. [scalar or vector] T: Ground Temp [in Celsius]. Default is 0 Celsius. [scalar or vector]

PyAstronomy.pyasl.co_refract(alt, observer_alt=0.0, pressure=None, temperature=None, epsilon=0.25, convert_to_observed=False, full_output=True, maxiter=10000)

Convert between apparent (real) altitude and observed altitude.

This routine converts between the apparent (real) altitude of the object, which does not include the influence of the atmosphere, and the observed altitude, which is the altitude as seen through the atmosphere.

The convert_to_observed flag determines the direction of the conversion. By default, observed altitude is converted into apparent altitude.

Parameters
altfloat or array

Altitude of an object in DEGREES. Whether the value is interpreted as apparent or observed altitude depends on the convert_to_observed flag. By default, it refers to the apparent (real) altitude.

observer_altfloat or array, optional

Altitude of the observer in METER. Default is 0 (sea level).

pressurefloat or array, optional

Atmospheric pressure in MILLIBAR. Default pressure is 1010 mbar. If observer_alt is given, an estimate for the real atmospheric pressure is calculated and used.

temperaturefloat or array, optional

Atmospheric temperature at the observing location in Celsius. If not specified, the temperature will be calculated assuming a ground temperature of 10 degrees Celsius.

epsilonfloat, optional

If convert_to_observed is TRUE, it specifies the accuracy of the calculated altitude in ARCSECONDS that should be reached by the iteration process.

convert_to_observedboolean, optional

If set True, an iterative method is used to calculate the observed altitude of an object, which includes atmospheric refraction. If False (default), the given altitude will be interpreted as the observed altitude and the apparent (real) altitude will be calculated using co_refract_forward().

full_outputboolean, optional

If True (default), pressure and temperature used in the calculation will be returned as well.

maxiterint, optional

The maximum number of iterations in the calculation. Throws PyAAlgorithmFailure if exceeded. Prevents endless loop. Default is 10000.

Returns
Altitudearray

By default, this will be the observed altitude of the object in degrees. If convert_to_observed was set to False, the number refers to the apparent (real) altitude.

Pressurearray

The pressure [mbar] used in the calculations (only returned if full_output is True).

Temperaturearray

The temperature used in the calculations [K] (only returned if full_output is True).

Notes

Note

This function was ported from the IDL Astronomy User’s Library.

IDL - Documentation

NAME:

CO_REFRACT()

PURPOSE:

Calculate correction to altitude due to atmospheric refraction.

DESCRIPTION:

CO_REFRACT can calculate both apparent altitude from observed altitude and vice-versa.

CALLING SEQUENCE:
new_alt = CO_REFRACT(old_alt, [ ALTITUDE= , PRESSURE= , $

TEMPERATURE= , /TO_OBSERVED , EPSILON= ])

INPUT:
old_alt - Observed (apparent) altitude, in DEGREES. (apparent if keyword

/TO_OBSERVED set). May be scalar or vector.

OUTPUT:
Function returns apparent (observed) altitude, in DEGREES. (observed if

keyword /TO_OBSERVED set). Will be of same type as input altitude(s).

OPTIONAL KEYWORD INPUTS:
ALTITUDEThe 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=0, i.e. sea level]

PRESSURE : The pressure at the observing location, in millibars. TEMPERATURE: The temperature at the observing location, in Kelvin. EPSILON: When keyword /TO_OBSERVED has been set, this is the accuracy

to obtain via the iteration, in arcseconds [default = 0.25

arcseconds].

/TO_OBSERVED: Set this keyword to go from Apparent->Observed altitude,

using the iterative technique.

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

DESCRIPTION:

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 always GREATER than the apparent altitude.

Apparent Altitude: The altitude that a star would be at, if *there were no

atmosphere* (sometimes called “true” altitude). This is usually calculated from an object’s celestial coordinates. Apparent altitude is always LOWER 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 couple 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.

EXAMPLE:

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

IDL> print, co_refract(0.5) ===> 0.025degrees (1.55’)

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.

It is important to understand that these formulae really can’t do better than about 30 arcseconds of accuracy very close to the horizon, as variable atmospheric effects become very important.

REFERENCES:
  1. Meeus, Astronomical Algorithms, Chapter 15.

  2. Explanatory Supplement to the Astronomical Almanac, 1992.

  3. Methods of Experimental Physics, Vol 12 Part B, Astrophysics, Radio Telescopes, Chapter 2.5, “Refraction Effects in the Neutral Atmosphere”, by R.K. Crane.

DEPENDENCIES:

CO_REFRACT_FORWARD (contained in this file and automatically compiled).

AUTHOR:
Chris O’Dell

Univ. of Wisconsin-Madison

Observational Cosmology Laboratory Email: odell@cmb.physics.wisc.edu

REVISION HISTORY:

version 1 (May 31, 2002) Update iteration formula, W. Landsman June 2002 Corrected slight bug associated with scalar vs. vector temperature and

pressure inputs. 6/10/2002

Fixed problem with vector input when /TO_OBSERVED set W. Landsman Dec 2005 Allow arrays with more than 32767 elements W.Landsman/C.Dickinson Feb 2010

Example: co_refract_forward and co_refract

from __future__ import print_function, division
from PyAstronomy import pyasl
import datetime
import numpy as np

# Assume, a star is observed at an altitude of 50 degrees
alt = 50.
# Now one wants to know the real altitude of the star, i.e.,
# the altitude corrected for atmospheric refraction.
print()
print("Get apparent (real) altitude of a star with observed altitude of " +
      str(alt) + " degrees")
print("  ->  Apparent altitude = ", alt - pyasl.co_refract_forward(alt))

print()
print("You are not observing from sea level, but from an altitude of 5000 meter.")
print(("Apparent altitude = %9.5f, estimated pressure [mbar] = %9.5f, " +
       "estimated temperature [K] = %9.5f") %
      pyasl.co_refract(alt, observer_alt=5000, convert_to_observed=False))

print()
print("Convert apparent (real) altitude into observed altitude.")
print("Apparent altitude = " + str(alt) + " degrees", end=' ')
print(" -> Observed altitude = " + str(pyasl.co_refract(alt, full_output=False,
                                                        convert_to_observed=True)[0]))

print()
print("The same object observed from different observer altitudes")
apparentAltitudes = np.repeat(30.0, 10)
obsalts = np.linspace(0., 5000., len(apparentAltitudes))
r = pyasl.co_refract(apparentAltitudes, observer_alt=obsalts, convert_to_observed=True)
for i in range(len(r[0])):
    print("Observed altitude [deg] = %g, pressure [mbar] = %g, temperature [K] = %g"
          % (r[0][i], r[1][i], r[2][i]))