def longitude_true_ascending_node(epoch): """This method computes the longitude of the true ascending node of the Moon in degrees, for a given instant, measured from the mean equinox of the date. :param epoch: Instant to compute the Moon's true ascending node, as an py:class:`Epoch` object :type epoch: :py:class:`Epoch` :returns: The longitude of the true ascending node. :rtype: py:class:`Angle` :raises: TypeError if input value is of wrong type. >>> epoch = Epoch(1913, 5, 27.0) >>> Omega = Moon.longitude_true_ascending_node(epoch) >>> print(round(Omega, 4)) 0.8763 """ # First check that input values are of correct types if not (isinstance(epoch, Epoch)): raise TypeError("Invalid input type") # Let's start computing the longitude of the MEAN ascending node Omega = Moon.longitude_mean_ascending_node(epoch) # Get the time from J2000.0 in Julian centuries t = (epoch - JDE2000) / 36525.0 # Mean elongation of the Moon D = 297.8501921 + (445267.1114034 + (-0.0018819 + (1.0/545868.0 - t/113065000.0) * t) * t) * t # Sun's mean anomaly M = 357.5291092 + (35999.0502909 + (-0.0001536 + t/24490000.0) * t) * t # Moon's mean anomaly Mprime = 134.9633964 + (477198.8675055 + (0.0087414 + (1.0/69699.9 + t/14712000.0) * t) * t) * t # Moon's argument of latitude F = 93.2720950 + (483202.0175233 + (-0.0036539 + (-1.0/3526000.0 + t/863310000.0) * t) * t) * t # Reduce the angles to a [0 360] range D = Angle(Angle.reduce_deg(D)).to_positive() Dr = D.rad() M = Angle(Angle.reduce_deg(M)).to_positive() Mr = M.rad() Mprime = Angle(Angle.reduce_deg(Mprime)).to_positive() Mprimer = Mprime.rad() F = Angle(Angle.reduce_deg(F)).to_positive() Fr = F.rad() # Compute the periodic terms corr = (-1.4979 * sin(2.0 * (Dr - Fr)) - 0.15 * sin(Mr) - 0.1226 * sin(2.0 * Dr) + 0.1176 * sin(2.0 * Fr) - 0.0801 * sin(2.0 * (Mprimer - Fr))) Omega += Angle(corr) return Omega
def ephemeris_physical_observations(epoch): """This method uses Carrington's formulas to compute the following quantities: - P : position angle of the northern extremity of the axis of rotation - B0 : heliographic latitude of the center of the solar disk - L0 : heliographic longitude of the center of the solar disk :param epoch: Epoch to compute the parameters :type epoch: :py:class:`Epoch` :returns: Parameters P, B0 and L0, in a tuple :rtype: tuple :raises: TypeError if input value is of wrong type. >>> epoch = Epoch(1992, 10, 13) >>> p, b0, l0 = Sun.ephemeris_physical_observations(epoch) >>> print(round(p, 2)) 26.27 >>> print(round(b0, 2)) 5.99 >>> print(round(l0, 2)) 238.63 """ # First check that input values are of correct types if not isinstance(epoch, Epoch): raise TypeError("Invalid input type") # Compute the auxiliary parameters epoch += 0.00068 theta = (epoch() - 2398220.0) * 360.0 / 25.38 theta = Angle(theta) i = Angle(7.25) k = 73.6667 + 1.3958333 * (epoch() - 2396758.0) / 36525.0 k = Angle(k) lon, lat, r = Sun.apparent_geocentric_position(epoch, nutation=False) eps = true_obliquity(epoch) dpsi = nutation_longitude(epoch) lonp = lon + dpsi x = atan(-cos(lonp.rad()) * tan(eps.rad())) x = Angle(x, radians=True) delta = lon - k y = atan(-cos(delta.rad()) * tan(i.rad())) y = Angle(y, radians=True) p = x + y b0 = asin(sin(delta.rad()) * sin(i.rad())) b0 = Angle(b0, radians=True) eta = atan(tan(delta.rad()) * cos(i.rad())) eta = Angle(eta, radians=True) l0 = eta - theta return p, b0, l0.to_positive()
def apparent_longitude_coarse(epoch): """This method provides the Sun's apparent longitude with a relatively low accuracy of about 0.01 degree. :param epoch: Epoch to compute the position of the Sun :type epoch: :py:class:`Epoch` :returns: A tuple containing the sun_apparent (ecliptical) longitude (as an Angle object) and the radius vector in astronomical units. :rtype: tuple :raises: TypeError if input value is of wrong type. >>> epoch = Epoch(1992, 10, 13) >>> app_lon, r = Sun.apparent_longitude_coarse(epoch) >>> print(app_lon.dms_str(n_dec=0)) 199d 54' 32.0'' >>> print(round(r, 5)) 0.99766 """ # First check that input values are of correct types if not isinstance(epoch, Epoch): raise TypeError("Invalid input type") # First find the true longitude sun = Sun() true_lon, r = sun.true_longitude_coarse(epoch) # Compute the time in Julian centuries t = (epoch - JDE2000) / 36525.0 # Then correct for nutation and aberration omega = 125.04 - 1934.136 * t omega = Angle(omega) lambd = true_lon - 0.00569 - 0.00478 * sin(omega.rad()) return (lambd, r)
def beginning_synodic_rotation(number): """This method calculates the epoch when the Carrington's synodic rotation No. 'number' starts. :param number: Number of Carrington's synodic rotation :type number: int :returns: Epoch when the provided rotation starts :rtype: :py:class:`Epoch` :raises: TypeError if input value is of wrong type. >>> epoch = Sun.beginning_synodic_rotation(1699) >>> print(round(epoch(), 3)) 2444480.723 """ # First check that input values are of correct types if not isinstance(number, int): raise TypeError("Invalid input type") # Apply formula (29.1) jde = 2398140.227 + 27.2752316 * number # Now, find the correction using formula (29.2) m = 281.96 + 26.882476 * number m = Angle(m) m = m.rad() delta = 0.1454 * sin(m) - 0.0085 * sin(2.0 * m) - 0.0141 * cos(2.0 * m) # Apply the correction jde += delta return Epoch(jde)
def geometric_heliocentric_position(epoch): """This method computes the geometric heliocentric position of planet Pluto for a given epoch. :param epoch: Epoch to compute Pluto position, as an Epoch object :type epoch: :py:class:`Epoch` :returns: A tuple with the heliocentric longitude and latitude (as :py:class:`Angle` objects), and the radius vector (as a float, in astronomical units), in that order :rtype: tuple :raises: TypeError if input value is of wrong type. :raises: ValueError if input epoch outside the 1885-2099 range. >>> epoch = Epoch(1992, 10, 13.0) >>> l, b, r = Pluto.geometric_heliocentric_position(epoch) >>> print(round(l, 5)) 232.74071 >>> print(round(b, 5)) 14.58782 >>> print(round(r, 6)) 29.711111 """ # First check that input value is of correct types if not isinstance(epoch, Epoch): raise TypeError("Invalid input type") # Check that the input epoch is within valid range y = epoch.year() if y < 1885.0 or y > 2099.0: raise ValueError("Epoch outside the 1885-2099 range") t = (epoch - JDE2000) / 36525.0 jj = 34.35 + 3034.9057 * t ss = 50.08 + 1222.1138 * t pp = 238.96 + 144.96 * t # Compute the arguments corr_lon = 0.0 corr_lat = 0.0 corr_rad = 0.0 for n in range(len(PLUTO_ARGUMENT)): iii, jjj, kkk = PLUTO_ARGUMENT[n] alpha = Angle(iii * jj + jjj * ss + kkk * pp).to_positive() alpha = alpha.rad() sin_a = sin(alpha) cos_a = cos(alpha) a_lon, b_lon = PLUTO_LONGITUDE[n] corr_lon += a_lon * sin_a + b_lon * cos_a a_lat, b_lat = PLUTO_LATITUDE[n] corr_lat += a_lat * sin_a + b_lat * cos_a a_rad, b_rad = PLUTO_RADIUS_VECTOR[n] corr_rad += a_rad * sin_a + b_rad * cos_a # The coefficients in the tables were scaled up. Let's scale them down corr_lon /= 1000000.0 corr_lat /= 1000000.0 corr_rad /= 10000000.0 lon = Angle(238.958116 + 144.96 * t + corr_lon) lat = Angle(-3.908239 + corr_lat) radius = 40.7241346 + corr_rad return lon, lat, radius
def true_longitude_coarse(epoch): """This method provides the Sun's true longitude with a relatively low accuracy of about 0.01 degree. :param epoch: Epoch to compute the position of the Sun :type epoch: :py:class:`Epoch` :returns: A tuple containing the true (ecliptical) longitude (as an Angle object) and the radius vector in astronomical units. :rtype: tuple :raises: TypeError if input value is of wrong type. >>> epoch = Epoch(1992, 10, 13) >>> true_lon, r = Sun.true_longitude_coarse(epoch) >>> print(true_lon.dms_str(n_dec=0)) 199d 54' 36.0'' >>> print(round(r, 5)) 0.99766 """ # First check that input values are of correct types if not (isinstance(epoch, Epoch)): raise TypeError("Invalid input type") # Compute the time in Julian centuries t = (epoch - JDE2000) / 36525.0 # Compute the geometric mean longitude of the Sun l0 = 280.46646 + t * (36000.76983 + t * 0.0003032) l0 = Angle(l0) l0.to_positive() # Now, compute the mean anomaly of the Sun m = 357.52911 + t * (35999.05029 - t * 0.0001537) m = Angle(m) mrad = m.rad() # The eccentricity of the Earth's orbit e = 0.016708634 - t * (0.000042037 + t * 0.0000001267) # Equation of the center c = ((1.914602 - t * (0.004817 + t * 0.000014)) * sin(mrad) + (0.019993 - t * 0.000101) * sin(2.0 * mrad) + 0.000289 * sin(3.0 * mrad)) c = Angle(c) true_lon = l0 + c true_anom = m + c # Sun's radius vector r = (1.000001018 * (1.0 - e * e)) / (1.0 + e * cos(true_anom.rad())) return (true_lon, r)
def test_coordinates_equatorial2horizontal(): """Tests the equatorial2horizontal() method of Coordinates module""" lon = Angle(77, 3, 56) lat = Angle(38, 55, 17) ra = Angle(23, 9, 16.641, ra=True) dec = Angle(-6, 43, 11.61) theta0 = Angle(8, 34, 57.0896, ra=True) eps = Angle(23, 26, 36.87) delta = Angle(0, 0, ((-3.868 * cos(eps.rad())) / 15.0), ra=True) theta0 += delta h = theta0 - lon - ra azi, ele = equatorial2horizontal(h, dec, lat) assert abs(round(azi, 3) - 68.034) < TOL, \ "ERROR: 1st equatorial2horizontal() test, 'azimuth' doesn't match" assert abs(round(ele, 3) - 15.125) < TOL, \ "ERROR: 2nd equatorial2horizontal() test, 'elevation' doesn't match"
def apparent_rightascension_declination_coarse(epoch): """This method provides the Sun's apparent right ascension and declination with a relatively low accuracy of about 0.01 degree. :param epoch: Epoch to compute the position of the Sun :type epoch: :py:class:`Epoch` :returns: A tuple containing the right ascension and the declination (as Angle objects) and the radius vector in astronomical units. :rtype: tuple :raises: TypeError if input value is of wrong type. >>> epo = Epoch(1992, 10, 13) >>> ra, delta, r = Sun.apparent_rightascension_declination_coarse(epo) >>> print(ra.ra_str(n_dec=1)) 13h 13' 31.4'' >>> print(delta.dms_str(n_dec=0)) -7d 47' 6.0'' >>> print(round(r, 5)) 0.99766 """ # First check that input values are of correct types if not isinstance(epoch, Epoch): raise TypeError("Invalid input type") # Second, find the apparent longitude sun = Sun() app_lon, r = sun.apparent_longitude_coarse(epoch) # Compute the obliquity of the ecliptic e0 = mean_obliquity(epoch) # Compute the time in Julian centuries t = (epoch - JDE2000) / 36525.0 # Then correct for nutation and aberration omega = 125.04 - 1934.136 * t omega = Angle(omega) # Correct the obliquity e = e0 + 0.00256 * cos(omega.rad()) alpha = atan2(cos(e.rad()) * sin(app_lon.rad()), cos(app_lon.rad())) alpha = Angle(alpha, radians=True) alpha.to_positive() delta = asin(sin(e.rad()) * sin(app_lon.rad())) delta = Angle(delta, radians=True) return (alpha, delta, r)
def heliocentric_ecliptical_position(self, epoch): """This method computes the heliocentric position of a minor celestial body, providing the result in ecliptical coordinates. :param epoch: Epoch to compute geocentric position, as an Epoch object :type epoch: :py:class:`Epoch` :returns: A tuple containing longitude and latitude, as Angle objects :rtype: tuple :raises: TypeError if input value is of wrong type. >>> a = 2.2091404 >>> e = 0.8502196 >>> q = a * (1.0 - e) >>> i = Angle(11.94524) >>> omega = Angle(334.75006) >>> w = Angle(186.23352) >>> t = Epoch(1990, 10, 28.54502) >>> epoch = Epoch(1990, 10, 6.0) >>> minor = Minor(q, e, i, omega, w, t) >>> lon, lat = minor.heliocentric_ecliptical_position(epoch) >>> print(lon.dms_str(n_dec=1)) 66d 51' 57.8'' >>> print(lat.dms_str(n_dec=1)) 11d 56' 14.4'' """ # First check that input value is of correct types if not isinstance(epoch, Epoch): raise TypeError("Invalid input type") # Get the mean motion and other orbital parameters n = self._n a = self._a e = self._e i = self._i omega = self._omega w = self._w t = self._t # Time since perihelion t_peri = epoch - t # Now, compute the mean anomaly, in degrees m = t_peri * n m = Angle(m) # With the mean anomaly, use Kepler's equation to find E and v ee, v = kepler_equation(e, m) ee = Angle(ee).to_positive() # Get r er = ee.rad() r = a * (1.0 - e * cos(er)) # Compute the heliocentric rectangular ecliptical coordinates wr = w.rad() vr = Angle(v).rad() ur = wr + vr omer = omega.rad() ir = i.rad() x = r * (cos(omer) * cos(ur) - sin(omer) * sin(ur) * cos(ir)) y = r * (sin(omer) * cos(ur) + cos(omer) * sin(ur) * cos(ir)) z = r * sin(ir) * sin(ur) lon = atan2(y, x) lat = atan2(z, sqrt(x * x + y * y)) return Angle(lon, radians=True), Angle(lat, radians=True)
def geocentric_position(self, epoch): """This method computes the geocentric position of a minor celestial body (right ascension and declination) for the given epoch, and referred to the standard equinox J2000.0. Additionally, it also computes the elongation angle to the Sun. :param epoch: Epoch to compute geocentric position, as an Epoch object :type epoch: :py:class:`Epoch` :returns: A tuple containing the right ascension, the declination and the elongation angle to the Sun, as Angle objects :rtype: tuple :raises: TypeError if input value is of wrong type. >>> a = 2.2091404 >>> e = 0.8502196 >>> q = a * (1.0 - e) >>> i = Angle(11.94524) >>> omega = Angle(334.75006) >>> w = Angle(186.23352) >>> t = Epoch(1990, 10, 28.54502) >>> minor = Minor(q, e, i, omega, w, t) >>> epoch = Epoch(1990, 10, 6.0) >>> ra, dec, p = minor.geocentric_position(epoch) >>> print(ra.ra_str(n_dec=1)) 10h 34' 13.7'' >>> print(dec.dms_str(n_dec=0)) 19d 9' 32.0'' >>> print(round(p, 2)) 40.51 >>> t = Epoch(1998, 4, 14.4358) >>> q = 1.487469 >>> e = 1.0 >>> i = Angle(0.0) >>> omega = Angle(0.0) >>> w = Angle(0.0) >>> minor = Minor(q, e, i, omega, w, t) >>> epoch = Epoch(1998, 8, 5.0) >>> ra, dec, p = minor.geocentric_position(epoch) >>> print(ra.ra_str(n_dec=1)) 5h 45' 34.5'' >>> print(dec.dms_str(n_dec=0)) 23d 23' 53.0'' >>> print(round(p, 2)) 45.73 """ # First check that input value is of correct types if not isinstance(epoch, Epoch): raise TypeError("Invalid input type") # Get internal parameters aa, bb, cc = self._aa, self._bb, self._cc am, bm, cm = self._am, self._bm, self._cm # Get the mean motion and other orbital parameters n = self._n a = self._a e = self._e w = self._w t = self._t # Time since perihelion t_peri = epoch - t # Now, compute the mean anomaly, in degrees m = t_peri * n m = Angle(m) if e < 0.98: # Elliptic case # With the mean anomaly, use Kepler's equation to find E and v ee, v = kepler_equation(e, m) ee = Angle(ee).to_positive() # Get r er = ee.rad() rr = a * (1.0 - e * cos(er)) elif abs(e - 1.0) < self._tol: # Parabolic case q = self._q ww = (0.03649116245 * (epoch - self._t)) / (q * sqrt(q)) sp = ww / 3.0 iterate = True while iterate: s = (2.0 * sp * sp * sp + ww) / (3.0 * (sp * sp + 1.0)) iterate = abs(s - sp) > self._tol sp = s v = 2.0 * atan(s) v = Angle(v, radians=True) rr = q * (1.0 + s * s) else: # We are in the near-parabolic case v, rr = self._near_parabolic(t_peri) # Compute the heliocentric rectangular equatorial coordinates wr = w.rad() vr = Angle(v).rad() x = rr * am * sin(aa + wr + vr) y = rr * bm * sin(bb + wr + vr) z = rr * cm * sin(cc + wr + vr) # Now let's compute Sun's rectangular coordinates xs, ys, zs = Sun.rectangular_coordinates_j2000(epoch) xi = x + xs eta = y + ys zeta = z + zs delta = sqrt(xi * xi + eta * eta + zeta * zeta) # We need to correct for the effect of light-time. Compute delay tau tau = 0.0057755183 * delta # Recompute some critical parameters t_peri = epoch - t - tau # Now, compute the mean anomaly, in degrees m = t_peri * n m = Angle(m) if e < 0.98: # Elliptic case # With the mean anomaly, use Kepler's equation to find E and v ee, v = kepler_equation(e, m) ee = Angle(ee).to_positive() # Get r er = ee.rad() rr = a * (1.0 - e * cos(er)) elif abs(e - 1.0) < self._tol: # Parabolic case q = self._q ww = (0.03649116245 * (epoch - self._t)) / (q * sqrt(q)) sp = ww / 3.0 iterate = True while iterate: s = (2.0 * sp * sp * sp + ww) / (3.0 * (sp * sp + 1.0)) iterate = abs(s - sp) > self._tol sp = s v = 2.0 * atan(s) v = Angle(v, radians=True) rr = q * (1.0 + s * s) else: # We are in the near-parabolic case v, rr = self._near_parabolic(t_peri) # Compute the heliocentric rectangular equatorial coordinates wr = w.rad() vr = Angle(v).rad() x = rr * am * sin(aa + wr + vr) y = rr * bm * sin(bb + wr + vr) z = rr * cm * sin(cc + wr + vr) xi = x + xs eta = y + ys zeta = z + zs ra = Angle(atan2(eta, xi), radians=True) dec = Angle(atan2(zeta, sqrt(xi * xi + eta * eta)), radians=True) r_sun = sqrt(xs * xs + ys * ys + zs * zs) psi = acos((xi * xs + eta * ys + zeta * zs) / (r_sun * delta)) psi = Angle(psi, radians=True) return ra, dec, psi
def geocentric_ecliptical_pos(epoch): """This method computes the geocentric ecliptical position (longitude, latitude) of the Moon for a given instant, referred to the mean equinox of the date, as well as the Moon-Earth distance in kilometers and the equatorial horizontal parallax. :param epoch: Instant to compute the Moon's position, as an py:class:`Epoch` object :type epoch: :py:class:`Epoch` :returns: Tuple containing: * Geocentric longitude of the center of the Moon, as an py:class:`Epoch` object. * Geocentric latitude of the center of the Moon, as an py:class:`Epoch` object. * Distance in kilometers between the centers of Earth and Moon, in kilometers (float) * Equatorial horizontal parallax of the Moon, as an py:class:`Epoch` object. :rtype: tuple :raises: TypeError if input value is of wrong type. >>> epoch = Epoch(1992, 4, 12.0) >>> Lambda, Beta, Delta, ppi = Moon.geocentric_ecliptical_pos(epoch) >>> print(round(Lambda, 6)) 133.162655 >>> print(round(Beta, 6)) -3.229126 >>> print(round(Delta, 1)) 368409.7 >>> print(round(ppi, 5)) 0.99199 """ # First check that input values are of correct types if not (isinstance(epoch, Epoch)): raise TypeError("Invalid input type") # Get the time from J2000.0 in Julian centuries t = (epoch - JDE2000) / 36525.0 # Compute Moon's mean longitude, referred to mean equinox of date Lprime = 218.3164477 + (481267.88123421 + (-0.0015786 + (1.0/538841.0 - t/65194000.0) * t) * t) * t # Mean elongation of the Moon D = 297.8501921 + (445267.1114034 + (-0.0018819 + (1.0/545868.0 - t/113065000.0) * t) * t) * t # Sun's mean anomaly M = 357.5291092 + (35999.0502909 + (-0.0001536 + t/24490000.0) * t) * t # Moon's mean anomaly Mprime = 134.9633964 + (477198.8675055 + (0.0087414 + (1.0/69699.9 + t/14712000.0) * t) * t) * t # Moon's argument of latitude F = 93.2720950 + (483202.0175233 + (-0.0036539 + (-1.0/3526000.0 + t/863310000.0) * t) * t) * t # Let's compute some additional arguments A1 = 119.75 + 131.849 * t A2 = 53.09 + 479264.290 * t A3 = 313.45 + 481266.484 * t # Eccentricity of Earth's orbit around the Sun E = 1.0 + (-0.002516 - 0.0000074 * t) * t E2 = E * E # Reduce the angles to a [0 360] range Lprime = Angle(Angle.reduce_deg(Lprime)).to_positive() Lprimer = Lprime.rad() D = Angle(Angle.reduce_deg(D)).to_positive() Dr = D.rad() M = Angle(Angle.reduce_deg(M)).to_positive() Mr = M.rad() Mprime = Angle(Angle.reduce_deg(Mprime)).to_positive() Mprimer = Mprime.rad() F = Angle(Angle.reduce_deg(F)).to_positive() Fr = F.rad() A1 = Angle(Angle.reduce_deg(A1)).to_positive() A1r = A1.rad() A2 = Angle(Angle.reduce_deg(A2)).to_positive() A2r = A2.rad() A3 = Angle(Angle.reduce_deg(A3)).to_positive() A3r = A3.rad() # Let's store this results in a list, in preparation for using tables arguments = [Dr, Mr, Mprimer, Fr] # Now we use the tables of periodic terms. First for sigmal and sigmar sigmal = 0.0 sigmar = 0.0 for i, value in enumerate(PERIODIC_TERMS_LR_TABLE): argument = 0.0 for j in range(4): if PERIODIC_TERMS_LR_TABLE[i][j]: # Avoid multiply by zero argument += PERIODIC_TERMS_LR_TABLE[i][j] * arguments[j] coeffl = value[4] coeffr = value[5] if abs(value[1]) == 1: coeffl = coeffl * E coeffr = coeffr * E elif abs(value[1]) == 2: coeffl = coeffl * E2 coeffr = coeffr * E2 sigmal += coeffl * sin(argument) sigmar += coeffr * cos(argument) # Add the additive terms to sigmal sigmal += (3958.0 * sin(A1r) + 1962.0 * sin(Lprimer - Fr) + 318.0 * sin(A2r)) # Now use the tabla for sigmab sigmab = 0.0 for i, value in enumerate(PERIODIC_TERMS_B_TABLE): argument = 0.0 for j in range(4): if PERIODIC_TERMS_B_TABLE[i][j]: # Avoid multiply by zero argument += PERIODIC_TERMS_B_TABLE[i][j] * arguments[j] coeffb = value[4] if abs(value[1]) == 1: coeffb = coeffb * E elif abs(value[1]) == 2: coeffb = coeffb * E2 sigmab += coeffb * sin(argument) # Add the additive terms to sigmab sigmab += (-2235.0 * sin(Lprimer) + 382.0 * sin(A3r) + 175.0 * sin(A1r - Fr) + 175.0 * sin(A1r + Fr) + 127.0 * sin(Lprimer - Mprimer) - 115.0 * sin(Lprimer + Mprimer)) Lambda = Lprime + (sigmal / 1000000.0) Beta = Angle(sigmab / 1000000.0) Delta = 385000.56 + (sigmar / 1000.0) ppii = asin(6378.14 / Delta) ppi = Angle(ppii, radians=True) return Lambda, Beta, Delta, ppi
def get_equinox_solstice(year, target="spring"): """This method computes the times of the equinoxes or the solstices. :param year: Year we want to compute the equinox or solstice for :type year: int :param target: Corresponding equinox or solstice. It can be "spring", "summer", "autumn", "winter" :type target: str :returns: The instant of time when the equinox or solstice happens :rtype: :py:class:`Epoch` :raises: TypeError if input values are of wrong type. :raises: ValueError if 'target' value is invalid. >>> epoch = Sun.get_equinox_solstice(1962, target="summer") >>> y, m, d, h, mi, s = epoch.get_full_date() >>> print("{}/{}/{} {}:{}:{}".format(y, m, d, h, mi, round(s, 0))) 1962/6/21 21:24:42.0 """ # NOTE: The results from the previous example are computed using the # complete VSOP87 theory. The results provided by Meeus in exercises # 27.a and 27.b are computed with lower accuracy # First check that input values are of correct types if not (isinstance(year, int) and isinstance(target, str)): raise TypeError("Invalid input types") # Second, check that the target is correct if ((target != "spring") and (target != "summer") and (target != "autumn") and (target != "winter")): raise ValueError("'target' value is invalid") # Now we can start computing an approximate value (Tables 27.A, 27.B) if (year >= -1000) and (year < 1000): y = year / 1000.0 if target == "spring": jde0 = 1721139.29189 + y * (365242.1374 + y * (0.06134 + y * (0.00111 - y * 0.00071))) elif target == "summer": jde0 = 1721233.25401 + y * (365241.72562 + y * (-0.05323 + y * (0.00907 + y * 0.00025))) elif target == "autumn": jde0 = (1721325.70455 + y * (365242.49558 + y * (-0.11677 + y * (-0.00297 + y * 0.00074)))) elif target == "winter": jde0 = (1721414.39987 + y * (363242.88257 + y * (-0.00769 + y * (-0.00933 - y * 0.00006)))) elif (year >= 1000) and (year <= 3000): y = (year - 2000.0) / 1000.0 if target == "spring": jde0 = 2451623.80984 + y * (365242.37404 + y * (0.05169 + y * (-0.00411 - y * 0.00057))) elif target == "summer": jde0 = 2451716.56767 + y * (365241.62603 + y * (0.00325 + y * (0.00888 - y * 0.0003))) elif target == "autumn": jde0 = 2451810.21715 + y * (365242.01767 + y * (-0.11575 + y * (0.00337 + y * 0.00078))) elif target == "winter": jde0 = (2451900.05952 + y * (365242.74049 + y * (-0.06223 + y * (-0.00823 + y * 0.00032)))) else: raise ValueError("'year' value out of range") k = ["spring", "summer", "autumn", "winter"].index(target) epoch = Epoch(jde0) corr = 1.0 while abs(corr) > 0.0000025: lon, lat, r = Sun.apparent_geocentric_position(epoch) arg = k * 90.0 - lon.to_positive() arg = Angle(arg) corr = 58.0 * sin(arg.rad()) epoch += corr epoch -= corr return epoch
def rectangular_coordinates_equinox(epoch, equinox_epoch): """This method computes the rectangular geocentric equatorial coordinates (X, Y, Z) of the Sun, referred to an arbitrary mean equinox. The X axis is directed towards the vernal equinox (longitude 0), the Y axis lies in the plane of the equator and is directed towards longitude 90, and the Z axis is directed towards the north celestial pole. :param epoch: Epoch to compute Sun position, as an Epoch object :type epoch: :py:class:`Epoch` :param equinox_epoch: Epoch corresponding to the mean equinox :type equinox_epoch: :py:class:`Epoch` :returns: A tuple with the X, Y, Z values in astronomical units :rtype: tuple :raises: TypeError if input values are of wrong type. >>> epoch = Epoch(1992, 10, 13.0) >>> e_equinox = Epoch(2467616.0) >>> x, y, z = Sun.rectangular_coordinates_equinox(epoch, e_equinox) >>> print(round(x, 8)) -0.93368986 >>> print(round(y, 8)) -0.32235085 >>> print(round(z, 8)) -0.13977098 """ # First check that input values are of correct types if (not isinstance(epoch, Epoch) and not isinstance(equinox_epoch, Epoch)): raise TypeError("Invalid input types") # Second, compute Sun's rectangular coordinates w.r.t. J2000.0 x0, y0, z0 = Sun.rectangular_coordinates_j2000(epoch) # Third, computed auxiliary angles t = (equinox_epoch - JDE2000) / 36525.0 tt = (epoch - equinox_epoch) / 36525.0 # Compute the conversion parameters zeta = t * ((2306.2181 + tt * (1.39656 - 0.000139 * tt)) + t * ((0.30188 - 0.000344 * tt) + 0.017998 * t)) z = t * ((2306.2181 + tt * (1.39656 - 0.000139 * tt)) + t * ((1.09468 + 0.000066 * tt) + 0.018203 * t)) theta = t * (2004.3109 + tt * (-0.85330 - 0.000217 * tt) + t * (-(0.42665 + 0.000217 * tt) - 0.041833 * t)) # Redefine the former values as Angles, and compute them in radians zeta = Angle(0, 0, zeta) zetar = zeta.rad() z = Angle(0, 0, z) zr = z.rad() theta = Angle(0, 0, theta) thetar = theta.rad() xx = cos(zetar) * cos(zr) * cos(thetar) - sin(zetar) * sin(zr) xy = sin(zetar) * cos(zr) + cos(zetar) * sin(zr) * cos(thetar) xz = cos(zetar) * sin(thetar) yx = -cos(zetar) * sin(zr) - sin(zetar) * cos(zr) * cos(thetar) yy = cos(zetar) * cos(zr) - sin(zetar) * sin(zr) * cos(thetar) yz = -sin(zetar) * sin(thetar) zx = -cos(zr) * sin(thetar) zy = -sin(zr) * sin(thetar) zz = cos(thetar) xp = xx * x0 + yx * y0 + zx * z0 yp = xy * x0 + yy * y0 + zy * z0 zp = xz * x0 + yz * y0 + zz * z0 return xp, yp, zp
def test_angle_rad(): """Tests the rad() method of Angle class""" a = Angle(180.0) assert abs(a.rad() - pi) < TOL, \ "ERROR: In 1st rad() test, radians value doesn't match"