isofit.core.sunposition

Classes

_sp

Sunposition

Compute sun position parameters given the time and location.

Functions

julian_day(dt)

Convert UTC datetimes or UTC timestamps to Julian days.

arcdist(p0, p1[, radians])

Angular distance between azimuth, zenith pairs.

observed_sunpos(dt, latitude, longitude, elevation[, ...])

Compute the observed coordinates of the sun as viewed at the given time and location.

topocentric_sunpos(dt, latitude, longitude[, ...])

Compute the topocentric coordinates of the sun as viewed at the given time and location.

sunpos(dt, latitude, longitude, elevation[, ...])

Compute the observed and topocentric coordinates of the sun as viewed at the given time and location.

Module Contents

class _sp[source]

.

static calendar_time(dt)[source]

.

static julian_day(dt)[source]

Calculate the Julian Day from a datetime.datetime object in UTC.

static julian_ephemeris_day(jd, deltat)[source]

Calculate the Julian Ephemeris Day from the Julian Day and delta-time = (terrestrial time - universal time) in seconds.

static julian_century(jd)[source]

Caluclate the Julian Century from Julian Day or Julian Ephemeris Day.

static julian_millennium(jc)[source]

Calculate the Julian Millennium from Julian Ephemeris Century.

_EHL_ = [[(175347046, 0.0, 0.0), (3341656, 4.6692568, 6283.07585), (34894, 4.6261, 12566.1517), (3497,...[source]
_EHB_ = [[(280, 3.199, 84334.662), (102, 5.422, 5507.553), (80, 3.88, 5223.69), (44, 3.7, 2352.87), (32,...[source]
_EHR_ = [[(100013989, 0.0, 0.0), (1670700, 3.0984635, 6283.07585), (13956, 3.05525, 12566.1517), (3084,...[source]
static heliocentric_longitude(jme)[source]

Compute the Earth Heliocentric Longitude (L) in degrees given the Julian Ephemeris Millennium.

static heliocentric_latitude(jme)[source]

Compute the Earth Heliocentric Latitude (B) in degrees given the Julian Ephemeris Millennium.

static heliocentric_radius(jme)[source]

Compute the Earth Heliocentric Radius (R) in astronimical units given the Julian Ephemeris Millennium.

static heliocentric_position(jme)[source]

Compute the Earth Heliocentric Longitude, Latitude, and Radius given the Julian Ephemeris Millennium.

Returns (L, B, R) where L = longitude in degrees, B = latitude in degrees, and R = radius in astronimical units.

static geocentric_position(helio_pos)[source]

Compute the geocentric latitude (Theta) and longitude (beta) (in degrees) of the sun given Earth’s heliocentric position (L, B, R).

_NLOY_[source]
_NLOab_[source]
_NLOcd_[source]
static ecliptic_obliquity(jme, delta_epsilon)[source]

Calculate the true obliquity of the ecliptic (epsilon, in degrees) given the Julian Ephemeris Millennium and the obliquity.

static nutation_obliquity(jce)[source]

Compute the nutation in longitude (delta_psi) and the true obliquity (epsilon) given the Julian Ephemeris Century.

static abberation_correction(R)[source]

Calculate the abberation correction (delta_tau, in degrees) given the Earth Heliocentric Radius (in AU).

static sun_longitude(helio_pos, delta_psi)[source]

Calculate the apparent sun longitude (lambda, in degrees) and geocentric longitude (beta, in degrees) given the earth heliocentric position and delta_psi.

static greenwich_sidereal_time(jd, delta_psi, epsilon)[source]

Calculate the apparent Greenwich sidereal time (v, in degrees) given the Julian Day.

static sun_ra_decl(llambda, epsilon, beta)[source]

Calculate the sun’s geocentric right ascension (alpha, in degrees) and declination (delta, in degrees).

static sun_topo_ra_decl_hour(latitude, longitude, elevation, jd, delta_t=0)[source]

Calculate the sun’s topocentric right ascension (alpha’), declination (delta’), and hour angle (H’).

static sun_topo_azimuth_zenith(latitude, delta_prime, H_prime, temperature=14.6, pressure=1013)[source]

Compute the sun’s topocentric azimuth and zenith angles.

Azimuth is measured eastward from north, zenith from vertical. Temperature = average temperature in C (default is 14.6 = global average in 2013). Pressure = average pressure in mBar (default 1013 = global average).

static norm_lat_lon(lat, lon)[source]

.

static topo_pos(t, lat, lon, elev, temp, press, dt)[source]

Compute RA,dec,H, all in degrees.

static pos(t, lat, lon, elev, temp, press, dt)[source]

Compute azimute,zenith,RA,dec,H all in degree.

julian_day(dt)[source]

Convert UTC datetimes or UTC timestamps to Julian days.

Parameters:

dt (array_like) – UTC datetime objects or UTC timestamps (as per datetime.utcfromtimestamp)

Returns:

jd – datetimes converted to fractional Julian days

Return type:

ndarray

arcdist(p0, p1, radians=False)[source]

Angular distance between azimuth, zenith pairs.

Parameters:
  • p0 (array_like, shape (..., 2))

  • p1 (array_like, shape (..., 2)) – p[…,0] = azimuth angles, p[…,1] = zenith angles

  • radians (boolean (default False)) – If False, angles are in degrees, otherwise in radians

Returns:

ad – Arcdistances between corresponding pairs in p0,p1 In degrees by default, in radians if radians=True

Return type:

array_like, shape is broadcast(p0,p1).shape

observed_sunpos(dt, latitude, longitude, elevation, temperature=None, pressure=None, delta_t=0, radians=False)[source]

Compute the observed coordinates of the sun as viewed at the given time and location.

Parameters:
  • dt (array_like) – UTC datetime objects or UTC timestamps (as per datetime.utcfromtimestamp) representing the times of observations

  • latitude (array_like) – decimal degrees, positive for north of the equator and east of Greenwich

  • longitude (array_like) – decimal degrees, positive for north of the equator and east of Greenwich

  • elevation (array_like) – meters, relative to the WGS-84 ellipsoid

  • temperature (array_like or None, optional) – celcius, default is 14.6 (global average in 2013)

  • pressure (array_like or None, optional) – millibar, default is 1013 (global average in ??)

  • delta_t (array_like, optional) – seconds, default is 0, difference between the earth’s rotation time (TT) and universal time (UT)

  • radians ({True, False}, optional) – return results in radians if True, degrees if False (default)

Returns:

coords – The shape of the array is parameters broadcast together, plus a final dimension for the coordinates. coords[…,0] = observed azimuth angle, measured eastward from north coords[…,1] = observed zenith angle, measured down from vertical

Return type:

ndarray, (…,2)

topocentric_sunpos(dt, latitude, longitude, temperature=None, pressure=None, delta_t=0, radians=False)[source]

Compute the topocentric coordinates of the sun as viewed at the given time and location.

Parameters:
  • dt (array_like) – UTC datetime objects or UTC timestamps (as per datetime.utcfromtimestamp) representing the times of observations

  • latitude (array_like) – decimal degrees, positive for north of the equator and east of Greenwich

  • longitude (array_like) – decimal degrees, positive for north of the equator and east of Greenwich

  • elevation (array_like) – meters, relative to the WGS-84 ellipsoid

  • temperature (array_like or None, optional) – celcius, default is 14.6 (global average in 2013)

  • pressure (array_like or None, optional) – millibar, default is 1013 (global average in ??)

  • delta_t (array_like, optional) – seconds, default is 0, difference between the earth’s rotation time (TT) and universal time (UT)

  • radians ({True, False}, optional) – return results in radians if True, degrees if False (default)

Returns:

coords – The shape of the array is parameters broadcast together, plus a final dimension for the coordinates. coords[…,0] = topocentric right ascension coords[…,1] = topocentric declination coords[…,2] = topocentric hour angle

Return type:

ndarray, (…,3)

sunpos(dt, latitude, longitude, elevation, temperature=None, pressure=None, delta_t=0, radians=False)[source]

Compute the observed and topocentric coordinates of the sun as viewed at the given time and location.

Parameters:
  • dt (array_like) – UTC datetime objects or UTC timestamps (as per datetime.utcfromtimestamp) representing the times of observations

  • latitude (array_like) – decimal degrees, positive for north of the equator and east of Greenwich

  • longitude (array_like) – decimal degrees, positive for north of the equator and east of Greenwich

  • elevation (array_like) – meters, relative to the WGS-84 ellipsoid

  • temperature (array_like or None, optional) – celcius, default is 14.6 (global average in 2013)

  • pressure (array_like or None, optional) – millibar, default is 1013 (global average in ??)

  • delta_t (array_like, optional) – seconds, default is 0, difference between the earth’s rotation time (TT) and universal time (UT)

  • radians ({True, False}, optional) – return results in radians if True, degrees if False (default)

Returns:

coords – The shape of the array is parameters broadcast together, plus a final dimension for the coordinates. coords[…,0] = observed azimuth angle, measured eastward from north coords[…,1] = observed zenith angle, measured down from vertical coords[…,2] = topocentric right ascension coords[…,3] = topocentric declination coords[…,4] = topocentric hour angle

Return type:

ndarray, (…,5)

class Sunposition(t, lat, lon, elev, temp, p, dt, rad, csv=False)[source]

Compute sun position parameters given the time and location.

t = None[source]
elev = None[source]
temp = None[source]
p = None[source]
dt = None[source]
rad = None[source]
az = None[source]
zen = None[source]
ra = None[source]
dec = None[source]
h = None[source]
lat[source]
lon[source]
property citation[source]

Print the citation.