libTheSky
Routines to compute sky positions of Sun, Moon, planets and more
All Namespaces Files Functions Variables Pages
thesky_riset Module Reference

Procedures to determine rise, transit and set times. More...

Functions/Subroutines

subroutine riset (jd, pl, rt, tt, st, rh, ta, sh, rsalt, for_night, ltime, cwarn, converge)
 Rise, transit and set times routine for Sun, Moon, planets and asteroids. This version recomputes positions (low accuracy first, then full accuracy) during the convergence process. See riset_ipol() for a version using interpolation, to be used for planets only(!)
 
subroutine riset_ipol (jd, pl, rt, tt, st, rh, ta, sh, rsalt, for_night, ltime, cwarn)
 Old routine for rise, transit and set times for planets and asteroids - don't use this for Sun and Moon Computes three sets of planet positions, and interpolates between them.
 
subroutine riset_ad (jd, ra, dec, rt, tt, st, rh, ta, sh, rsalt, for_night, cwarn)
 Rise, transit and set times routine for an object with fixed ra & dec.
 
subroutine best_obs_date (jd00, ra0, m, d)
 Compute the date (m,d) at which an object with ra can be observed best, i.e., it transits at midnight

 
real(double) function rsipol (y1, y2, y3, n)
 Quadratic interpolation of rise/set times for riset_ipol().
 

Detailed Description

Procedures to determine rise, transit and set times.

Function/Subroutine Documentation

◆ best_obs_date()

subroutine thesky_riset::best_obs_date ( real(double), intent(in) jd00,
real(double), intent(in) ra0,
integer, intent(out) m,
integer, intent(out) d )

Compute the date (m,d) at which an object with ra can be observed best, i.e., it transits at midnight

Parameters
jd00Julian day?
ra0Right ascension
mMonth of best observation (output)
dDay of month of best observation (output)
Note
  • VERY SLOW!!! -> use best_obs_date_ra instead!!! (but not (yet) in libTheSky...)
  • riset() should be faster now (04/2010)

Definition at line 702 of file riset.f90.

References thesky_datetime::gettz(), thesky_planetdata::nplanpos, thesky_planets::planet_position(), thesky_planetdata::planpos, riset(), and thesky_local::tz.

Here is the call graph for this function:

◆ riset()

subroutine thesky_riset::riset ( real(double), intent(in) jd,
integer, intent(in) pl,
real(double), intent(out) rt,
real(double), intent(out) tt,
real(double), intent(out) st,
real(double), intent(out) rh,
real(double), intent(out) ta,
real(double), intent(out) sh,
real(double), intent(in) rsalt,
logical, intent(in), optional for_night,
logical, intent(in), optional ltime,
logical, intent(in), optional cwarn,
integer, dimension(3), intent(out), optional converge )

Rise, transit and set times routine for Sun, Moon, planets and asteroids. This version recomputes positions (low accuracy first, then full accuracy) during the convergence process. See riset_ipol() for a version using interpolation, to be used for planets only(!)

Parameters
jdJulian day number. The rise/transit and set data are computed for the (calendar) DAY of this JD (midnight-midnight), unless for_night=.true. in which case this happens for the NIGHT starting at JD.
plPlanet/object number (-1 - 9: ES, Moon, Mer-Plu; >10000: asteroids)
rtRise time (hours) (output)
ttTransit time (hours) (output)
stSet time (hours) (output)
rhRising wind direction (rad) (output)
taTransit altitude (rad) (output)
shSetting wind direction (rad) (output)
rsAltAltitude to return rise/set data for (degrees; 0. is actual rise/set). rsAlt>90: compute transit only
for_nightCompute rise/tranit/set data for the NIGHT starting at the indicated JD, rather than the calendar day. Times between 0h and 12h are for the NEXT calendar day! (optional; default: false)
ltimePassed to planet_position(). If .true., include light time, doubling the CPU time while gaining a bit of accuracy (optional; default: false)
cWarnWarn upon convergence failure (optional; default: true)
convergeNumber of iterations needed to converge (optional) (output)
Note
  • for rsAlt = 0.d0, rise and set times are computed
  • for (rsAlt.ne.0), the routine calculates when alt=rsAlt is reached
  • use rsAlt=-6,-12,-18 for the sun for twilight calculations (rsAlt is expressed in degrees).
  • Returns times, rise/set azimuth and transit altitude
  • Moon transit error ~0.15s (?)
Todo
  • This version sometimes finds the answer tmRad(i)>2pi, which is the answer for the next day...
    • (only?) solution: use a solver on JD for az,alt
Note
Speed wrt riset_ipol(), transit only, when using low-accuracy approximations:
  • Moon: 270x faster (used to be slowest by factor 3-5 wrt planets, factor 6.5 wrt Sun)
  • Sun: 340x faster - full accuracy (VSOP): ~2-3x SLOWER!
  • Mer: 20% slower
  • Ven: 4% slower
  • Mar: 10% slower
  • Jup: 25% slower
  • Nep: 30% slower
See also
  • riset_ipol()
  • Meeus, Astronomical algorithms, Ch.15, but with geographic longitude east of Greenwich defined as > 0

Definition at line 82 of file riset.f90.

References thesky_local::deltat, thesky_local::lat0, thesky_local::lon0, thesky_planets::planet_position(), thesky_planets::planet_position_la(), thesky_planetdata::planpos, thesky_coordinates::refract(), riset_ipol(), and thesky_local::tz.

Referenced by best_obs_date(), thesky_visibility::best_planet_visibility(), and thesky_visibility::planet_visibility_tonight().

Here is the call graph for this function:

◆ riset_ad()

subroutine thesky_riset::riset_ad ( real(double), intent(in) jd,
real(double), intent(in) ra,
real(double), intent(in) dec,
real(double), intent(out) rt,
real(double), intent(out) tt,
real(double), intent(out) st,
real(double), intent(out) rh,
real(double), intent(out) ta,
real(double), intent(out) sh,
real(double), intent(in) rsalt,
logical, intent(in), optional for_night,
logical, intent(in), optional cwarn )

Rise, transit and set times routine for an object with fixed ra & dec.

Parameters
jdJulian day number
raRight ascension (rad)
decDeclination (rad)
rtRise time (hours) (output)
ttTransit time (hours) (output)
stSet time (hours) (output)
rhRising wind direction (rad) (output)
taTransit altitude (rad) (output)
shSetting wind direction (rad) (output)
rsAltAltitude to return rise/set data for (degrees; 0. is actual rise/set). rsAlt>90: compute transit only
for_nightCompute rise/tranit/set data for the NIGHT starting at the indicated JD, rather than the calendar day. Times between 0h and 12h are for the NEXT calendar day! (optional; default: false)
cWarnWarn upon convergence failure (optional; default: true)

for rsAlt = 0., rise and set times are computed for rsAlt.ne.0, the routine calculates when alt=rsAlt is reached use rsAlt=-6,-12,-18 for the sun for twilight calculations (rsAlt is expressed in degrees).

See also
  • Meeus, Astronomical algorithms, Ch.15, but with geographic longitude east of Greenwich defined as > 0

Definition at line 564 of file riset.f90.

References thesky_local::lat0, thesky_local::lon0, thesky_planets::planet_position(), thesky_planetdata::planpos, and thesky_local::tz.

Here is the call graph for this function:

◆ riset_ipol()

subroutine thesky_riset::riset_ipol ( real(double), intent(in) jd,
integer, intent(in) pl,
real(double), intent(out) rt,
real(double), intent(out) tt,
real(double), intent(out) st,
real(double), intent(out) rh,
real(double), intent(out) ta,
real(double), intent(out) sh,
real(double), intent(in) rsalt,
logical, intent(in), optional for_night,
logical, intent(in), optional ltime,
logical, intent(in), optional cwarn )

Old routine for rise, transit and set times for planets and asteroids - don't use this for Sun and Moon Computes three sets of planet positions, and interpolates between them.

Parameters
jdJulian day number. The rise/transit and set data are computed for the (calendar) DAY of this JD (midnight-midnight), unless for_night=.true. in which case this happens for the NIGHT starting at JD.
plPlanet/object number (-1 - 9: ES, Moon, Mer-Plu; >10000: asteroids)
rtRise time (hours) (output)
ttTransit time (hours) (output)
stSet time (hours) (output)
rhRising wind direction (rad) (output)
taTransit altitude (rad) (output)
shSetting wind direction (rad) (output)
rsAltAltitude to return rise/set data for (degrees; 0. is actual rise/set). rsAlt>90: compute transit only
for_nightCompute rise/tranit/set data for the NIGHT starting at the indicated JD, rather than the calendar day. Times between 0h and 12h are for the NEXT calendar day! (optional; default: false)
ltimePassed to planet_position(). If .true., include light time, doubling the CPU time while gaining a bit of accur.
cWarnWarn upon convergence failure (optional; default: true)
Note
  • for rsAlt = 0.d0, rise and set times are computed
  • for (rsAlt.ne.0), the routine calculates when alt=rsAlt is reached
  • use rsAlt=-6,-12,-18 for the sun for twilight calculations (rsAlt is expressed in degrees).
  • transit error for Moon <= ~10s (due to interpolation); typically ~4s?
  • code gets confused around midnight - it sometimes returns a small negative tmdy(i), which is on the day before(!)
See also
  • Meeus, Astronomical algorithms, Ch.15, but with geographic longitude east of Greenwich defined as > 0

Definition at line 318 of file riset.f90.

References thesky_local::deltat, thesky_local::lat0, thesky_local::lon0, thesky_planets::planet_position(), thesky_planetdata::planpos, thesky_coordinates::refract(), rsipol(), and thesky_local::tz.

Referenced by riset().

Here is the call graph for this function:

◆ rsipol()

real(double) function thesky_riset::rsipol ( real(double), intent(in) y1,
real(double), intent(in) y2,
real(double), intent(in) y3,
real(double), intent(in) n )

Quadratic interpolation of rise/set times for riset_ipol().

Parameters
y1Coordinate for day 1
y2Coordinate for day 2
y3Coordinate for day 3
nFraction of day
Return values
rsIpolInterpolated coordinate

Definition at line 773 of file riset.f90.

References rsipol().

Referenced by riset_ipol(), and rsipol().

Here is the call graph for this function: