<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd"> <html xmlns="http://www.w3.org/1999/xhtml"> <head> <meta http-equiv="Content-Type" content="text/html; charset=utf-8" /> <title>astropy.coordinates.builtin_systems — Astropy v0.2.4</title> <link rel="stylesheet" href="../../../_static/bootstrap-astropy.css" type="text/css" /> <link rel="stylesheet" href="../../../_static/pygments.css" type="text/css" /> <script type="text/javascript"> var DOCUMENTATION_OPTIONS = { URL_ROOT: '../../../', VERSION: '0.2.4', COLLAPSE_INDEX: false, FILE_SUFFIX: '.html', HAS_SOURCE: true }; </script> <script type="text/javascript" src="../../../_static/jquery.js"></script> <script type="text/javascript" src="../../../_static/underscore.js"></script> <script type="text/javascript" src="../../../_static/doctools.js"></script> <script type="text/javascript" src="../../../_static/sidebar.js"></script> <link rel="shortcut icon" href="../../../_static/astropy_logo.ico"/> <link rel="top" title="Astropy v0.2.4" href="../../../index.html" /> <link rel="up" title="Module code" href="../../index.html" /> </head> <body> <div class="topbar"> <a class="brand" title="Documentation Home" href="../../../index.html"></a> <ul> <li><a class="homelink" title="AstroPy Homepage" href="http://www.astropy.org"></a></li> <li><a title="General Index" href="../../../genindex.html">Index</a></li> <li><a title="Python Module Index" href="../../../py-modindex.html">Modules</a></li> <li> <form action="../../../search.html" method="get"> <input type="text" name="q" placeholder="Search" /> <input type="hidden" name="check_keywords" value="yes" /> <input type="hidden" name="area" value="default" /> </form> </li> </ul> </div> <div class="related"> <h3>Navigation</h3> <ul> <li> <a href="../../../index.html">Astropy v0.2.4</a> » </li> <li><a href="../../index.html" accesskey="U">Module code</a> »</li> </ul> </div> <div class="document"> <div class="documentwrapper"> <div class="bodywrapper"> <div class="body"> <h1>Source code for astropy.coordinates.builtin_systems</h1><div class="highlight"><pre> <span class="c"># Licensed under a 3-clause BSD style license - see LICENSE.rst</span> <span class="sd">"""</span> <span class="sd">This module contains the implementations of specific coordinate systems</span> <span class="sd">and the conversions between them.</span> <span class="sd">"""</span> <span class="kn">import</span> <span class="nn">numpy</span> <span class="kn">as</span> <span class="nn">np</span> <span class="kn">from</span> <span class="nn">.angles</span> <span class="kn">import</span> <span class="n">Angle</span> <span class="kn">from</span> <span class="nn">.coordsystems</span> <span class="kn">import</span> <span class="n">SphericalCoordinatesBase</span> <span class="kn">from</span> <span class="nn">..time</span> <span class="kn">import</span> <span class="n">Time</span> <span class="kn">from</span> <span class="nn">.</span> <span class="kn">import</span> <span class="n">transformations</span> <span class="kn">from</span> <span class="nn">..</span> <span class="kn">import</span> <span class="n">units</span> <span class="k">as</span> <span class="n">u</span> <span class="n">__all__</span> <span class="o">=</span> <span class="p">[</span><span class="s">'ICRSCoordinates'</span><span class="p">,</span> <span class="s">'FK5Coordinates'</span><span class="p">,</span> <span class="s">'FK4Coordinates'</span><span class="p">,</span> <span class="s">'FK4NoETermCoordinates'</span><span class="p">,</span> <span class="s">'GalacticCoordinates'</span><span class="p">,</span> <span class="s">'HorizontalCoordinates'</span> <span class="p">]</span> <span class="c"># The UTC time scale is not properly defined prior to 1960, so Time('B1950',</span> <span class="c"># scale='utc') will emit a warning. Instead, we use Time('B1950', scale='tai')</span> <span class="c"># which is equivalent, but does not emit a warning.</span> <span class="n">_EQUINOX_J2000</span> <span class="o">=</span> <span class="n">Time</span><span class="p">(</span><span class="s">'J2000'</span><span class="p">,</span> <span class="n">scale</span><span class="o">=</span><span class="s">'utc'</span><span class="p">)</span> <span class="n">_EQUINOX_B1950</span> <span class="o">=</span> <span class="n">Time</span><span class="p">(</span><span class="s">'B1950'</span><span class="p">,</span> <span class="n">scale</span><span class="o">=</span><span class="s">'tai'</span><span class="p">)</span> <span class="c">#<--------------Coordinate definitions; transformations are below--------------></span> <span class="nd">@transformations.coordinate_alias</span><span class="p">(</span><span class="s">'icrs'</span><span class="p">)</span> <div class="viewcode-block" id="ICRSCoordinates"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.ICRSCoordinates.html#astropy.coordinates.builtin_systems.ICRSCoordinates">[docs]</a><span class="k">class</span> <span class="nc">ICRSCoordinates</span><span class="p">(</span><span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> A coordinate in the ICRS.</span> <span class="sd"> If you're looking for "J2000" coordinates, and aren't sure if you</span> <span class="sd"> want to use this or `FK5Coordinates`, you probably want to use ICRS.</span> <span class="sd"> It's more well-defined as a catalog coordinate and is an inertial</span> <span class="sd"> system.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> {params}</span> <span class="sd"> obstime : `~astropy.time.Time` or None</span> <span class="sd"> The time of observation for this coordinate. If None, it will be taken</span> <span class="sd"> to be the same as the `equinox`.</span> <span class="sd"> Alternatively, a single argument that is any kind of spherical coordinate</span> <span class="sd"> can be provided, and will be converted to ICRSCoordinates and used as this</span> <span class="sd"> coordinate.</span> <span class="sd"> """</span> <span class="n">__doc__</span> <span class="o">=</span> <span class="n">__doc__</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">params</span><span class="o">=</span><span class="n">SphericalCoordinatesBase</span><span class="o">.</span><span class="n">_init_docstring_param_templ</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">lonnm</span><span class="o">=</span><span class="s">'ra'</span><span class="p">,</span> <span class="n">latnm</span><span class="o">=</span><span class="s">'dec'</span><span class="p">))</span> <span class="k">def</span> <span class="nf">__init__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="o">*</span><span class="n">args</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span> <span class="nb">super</span><span class="p">(</span><span class="n">ICRSCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">__init__</span><span class="p">()</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'obstime'</span><span class="p">,</span> <span class="bp">None</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified obstime is not None or a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="o">==</span> <span class="mi">1</span> <span class="ow">and</span> <span class="nb">len</span><span class="p">(</span><span class="n">kwargs</span><span class="p">)</span> <span class="o">==</span> <span class="mi">0</span> <span class="ow">and</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="n">newcoord</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">.</span><span class="n">transform_to</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">ra</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">dec</span> <span class="bp">self</span><span class="o">.</span><span class="n">_distance</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">_distance</span> <span class="k">else</span><span class="p">:</span> <span class="nb">super</span><span class="p">(</span><span class="n">ICRSCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">_initialize_latlon</span><span class="p">(</span><span class="s">'ra'</span><span class="p">,</span> <span class="s">'dec'</span><span class="p">,</span> <span class="bp">True</span><span class="p">,</span> <span class="n">args</span><span class="p">,</span> <span class="n">kwargs</span><span class="p">)</span> <span class="k">def</span> <span class="nf">__repr__</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">', Distance={0:.2g} {1!s}'</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_value</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">''</span> <span class="n">msg</span> <span class="o">=</span> <span class="s">"<{0} RA={1:.5f} deg, Dec={2:.5f} deg{3}>"</span> <span class="k">return</span> <span class="n">msg</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="o">.</span><span class="n">__name__</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="n">diststr</span><span class="p">)</span> <span class="nd">@property</span> <div class="viewcode-block" id="ICRSCoordinates.lonangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.ICRSCoordinates.html#astropy.coordinates.builtin_systems.ICRSCoordinates.lonangle">[docs]</a> <span class="k">def</span> <span class="nf">lonangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="ICRSCoordinates.latangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.ICRSCoordinates.html#astropy.coordinates.builtin_systems.ICRSCoordinates.latangle">[docs]</a> <span class="k">def</span> <span class="nf">latangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="ICRSCoordinates.equinox"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.ICRSCoordinates.html#astropy.coordinates.builtin_systems.ICRSCoordinates.equinox">[docs]</a> <span class="k">def</span> <span class="nf">equinox</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="n">_EQUINOX_J2000</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="ICRSCoordinates.obstime"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.ICRSCoordinates.html#astropy.coordinates.builtin_systems.ICRSCoordinates.obstime">[docs]</a> <span class="k">def</span> <span class="nf">obstime</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">equinox</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> </div></div> <span class="nd">@transformations.coordinate_alias</span><span class="p">(</span><span class="s">'fk5'</span><span class="p">)</span> <div class="viewcode-block" id="FK5Coordinates"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK5Coordinates.html#astropy.coordinates.builtin_systems.FK5Coordinates">[docs]</a><span class="k">class</span> <span class="nc">FK5Coordinates</span><span class="p">(</span><span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> A coordinate in the FK5 system.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> {params}</span> <span class="sd"> equinox : `~astropy.time.Time`, optional</span> <span class="sd"> The equinox for these coordinates. Defaults to J2000.</span> <span class="sd"> obstime : `~astropy.time.Time` or None</span> <span class="sd"> The time of observation for this coordinate. If None, it will be taken</span> <span class="sd"> to be the same as the `equinox`.</span> <span class="sd"> Alternatively, a single argument that is any kind of spherical coordinate</span> <span class="sd"> can be provided, and will be converted to `FK5Coordinates` and used as this</span> <span class="sd"> coordinate.</span> <span class="sd"> """</span> <span class="n">__doc__</span> <span class="o">=</span> <span class="n">__doc__</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">params</span><span class="o">=</span><span class="n">SphericalCoordinatesBase</span><span class="o">.</span><span class="n">_init_docstring_param_templ</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">lonnm</span><span class="o">=</span><span class="s">'ra'</span><span class="p">,</span> <span class="n">latnm</span><span class="o">=</span><span class="s">'dec'</span><span class="p">))</span> <span class="k">def</span> <span class="nf">__init__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="o">*</span><span class="n">args</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span> <span class="nb">super</span><span class="p">(</span><span class="n">FK5Coordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">__init__</span><span class="p">()</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'equinox'</span><span class="p">,</span> <span class="n">_EQUINOX_J2000</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'obstime'</span><span class="p">,</span> <span class="bp">None</span><span class="p">)</span> <span class="k">if</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified equinox is not a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified obstime is not None or a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="o">==</span> <span class="mi">1</span> <span class="ow">and</span> <span class="nb">len</span><span class="p">(</span><span class="n">kwargs</span><span class="p">)</span> <span class="o">==</span> <span class="mi">0</span> <span class="ow">and</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="n">newcoord</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">.</span><span class="n">transform_to</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">ra</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">dec</span> <span class="bp">self</span><span class="o">.</span><span class="n">_distance</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">_distance</span> <span class="k">else</span><span class="p">:</span> <span class="nb">super</span><span class="p">(</span><span class="n">FK5Coordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">_initialize_latlon</span><span class="p">(</span><span class="s">'ra'</span><span class="p">,</span> <span class="s">'dec'</span><span class="p">,</span> <span class="bp">True</span><span class="p">,</span> <span class="n">args</span><span class="p">,</span> <span class="n">kwargs</span><span class="p">)</span> <span class="k">def</span> <span class="nf">__repr__</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">', Distance={0:.2g} {1!s}'</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_value</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">''</span> <span class="n">msg</span> <span class="o">=</span> <span class="s">"<{0} RA={1:.5f} deg, Dec={2:.5f} deg{3}>"</span> <span class="k">return</span> <span class="n">msg</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="o">.</span><span class="n">__name__</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="n">diststr</span><span class="p">)</span> <span class="nd">@property</span> <div class="viewcode-block" id="FK5Coordinates.lonangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK5Coordinates.html#astropy.coordinates.builtin_systems.FK5Coordinates.lonangle">[docs]</a> <span class="k">def</span> <span class="nf">lonangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK5Coordinates.latangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK5Coordinates.html#astropy.coordinates.builtin_systems.FK5Coordinates.latangle">[docs]</a> <span class="k">def</span> <span class="nf">latangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK5Coordinates.equinox"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK5Coordinates.html#astropy.coordinates.builtin_systems.FK5Coordinates.equinox">[docs]</a> <span class="k">def</span> <span class="nf">equinox</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK5Coordinates.obstime"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK5Coordinates.html#astropy.coordinates.builtin_systems.FK5Coordinates.obstime">[docs]</a> <span class="k">def</span> <span class="nf">obstime</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">equinox</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> </div> <div class="viewcode-block" id="FK5Coordinates.precess_to"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK5Coordinates.html#astropy.coordinates.builtin_systems.FK5Coordinates.precess_to">[docs]</a> <span class="k">def</span> <span class="nf">precess_to</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">newequinox</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> Precesses the coordinates from their current `equinox` to a new equinox and</span> <span class="sd"> returns the resulting coordinate.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> newequinox : `~astropy.time.Time`</span> <span class="sd"> The equinox to precess these coordinates to.</span> <span class="sd"> Returns</span> <span class="sd"> -------</span> <span class="sd"> newcoord : FK5Coordinates</span> <span class="sd"> The new coordinate</span> <span class="sd"> """</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">precession_matrix_Capitaine</span> <span class="n">pmat</span> <span class="o">=</span> <span class="n">precession_matrix_Capitaine</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span><span class="p">,</span> <span class="n">newequinox</span><span class="p">)</span> <span class="n">v</span> <span class="o">=</span> <span class="p">[</span><span class="bp">self</span><span class="o">.</span><span class="n">x</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">y</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">z</span><span class="p">]</span> <span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="p">,</span> <span class="n">z</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">dot</span><span class="p">(</span><span class="n">pmat</span><span class="o">.</span><span class="n">A</span><span class="p">,</span> <span class="n">v</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="o">=</span><span class="n">y</span><span class="p">,</span> <span class="n">z</span><span class="o">=</span><span class="n">z</span><span class="p">,</span> <span class="n">unit</span><span class="o">=</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">newequinox</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="o">=</span><span class="n">y</span><span class="p">,</span> <span class="n">z</span><span class="o">=</span><span class="n">z</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">newequinox</span><span class="p">)</span> </div></div> <span class="nd">@transformations.coordinate_alias</span><span class="p">(</span><span class="s">'fk4'</span><span class="p">)</span> <div class="viewcode-block" id="FK4Coordinates"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4Coordinates.html#astropy.coordinates.builtin_systems.FK4Coordinates">[docs]</a><span class="k">class</span> <span class="nc">FK4Coordinates</span><span class="p">(</span><span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> A coordinate in the FK4 system.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> {params}</span> <span class="sd"> equinox : `~astropy.time.Time`, optional</span> <span class="sd"> The equinox for these coordinates. Defaults to B1950.</span> <span class="sd"> obstime : `~astropy.time.Time` or None</span> <span class="sd"> The time of observation for this coordinate. If None, it will be taken</span> <span class="sd"> to be the same as the `equinox`.</span> <span class="sd"> Alternatively, a single argument that is any kind of spherical coordinate</span> <span class="sd"> can be provided, and will be converted to `FK4Coordinates` and used as this</span> <span class="sd"> coordinate.</span> <span class="sd"> """</span> <span class="n">__doc__</span> <span class="o">=</span> <span class="n">__doc__</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">params</span><span class="o">=</span><span class="n">SphericalCoordinatesBase</span><span class="o">.</span><span class="n">_init_docstring_param_templ</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">lonnm</span><span class="o">=</span><span class="s">'ra'</span><span class="p">,</span> <span class="n">latnm</span><span class="o">=</span><span class="s">'dec'</span><span class="p">))</span> <span class="k">def</span> <span class="nf">__init__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="o">*</span><span class="n">args</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span> <span class="nb">super</span><span class="p">(</span><span class="n">FK4Coordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">__init__</span><span class="p">()</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'equinox'</span><span class="p">,</span> <span class="n">_EQUINOX_B1950</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'obstime'</span><span class="p">,</span> <span class="bp">None</span><span class="p">)</span> <span class="k">if</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified equinox is not a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified obstime is not None or a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="o">==</span> <span class="mi">1</span> <span class="ow">and</span> <span class="nb">len</span><span class="p">(</span><span class="n">kwargs</span><span class="p">)</span> <span class="o">==</span> <span class="mi">0</span> <span class="ow">and</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="n">newcoord</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">.</span><span class="n">transform_to</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">ra</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">dec</span> <span class="bp">self</span><span class="o">.</span><span class="n">_distance</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">_distance</span> <span class="k">else</span><span class="p">:</span> <span class="nb">super</span><span class="p">(</span><span class="n">FK4Coordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">_initialize_latlon</span><span class="p">(</span><span class="s">'ra'</span><span class="p">,</span> <span class="s">'dec'</span><span class="p">,</span> <span class="bp">True</span><span class="p">,</span> <span class="n">args</span><span class="p">,</span> <span class="n">kwargs</span><span class="p">)</span> <span class="k">def</span> <span class="nf">__repr__</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">', Distance={0:.2g} {1!s}'</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_value</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">''</span> <span class="n">msg</span> <span class="o">=</span> <span class="s">"<{0} RA={1:.5f} deg, Dec={2:.5f} deg{3}>"</span> <span class="k">return</span> <span class="n">msg</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="o">.</span><span class="n">__name__</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="n">diststr</span><span class="p">)</span> <span class="nd">@property</span> <div class="viewcode-block" id="FK4Coordinates.lonangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4Coordinates.html#astropy.coordinates.builtin_systems.FK4Coordinates.lonangle">[docs]</a> <span class="k">def</span> <span class="nf">lonangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK4Coordinates.latangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4Coordinates.html#astropy.coordinates.builtin_systems.FK4Coordinates.latangle">[docs]</a> <span class="k">def</span> <span class="nf">latangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK4Coordinates.equinox"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4Coordinates.html#astropy.coordinates.builtin_systems.FK4Coordinates.equinox">[docs]</a> <span class="k">def</span> <span class="nf">equinox</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK4Coordinates.obstime"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4Coordinates.html#astropy.coordinates.builtin_systems.FK4Coordinates.obstime">[docs]</a> <span class="k">def</span> <span class="nf">obstime</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">equinox</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> </div> <div class="viewcode-block" id="FK4Coordinates.precess_to"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4Coordinates.html#astropy.coordinates.builtin_systems.FK4Coordinates.precess_to">[docs]</a> <span class="k">def</span> <span class="nf">precess_to</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">newequinox</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> Precesses the coordinates from their current `equinox` to a new equinox.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> newequinox : `~astropy.time.Time`</span> <span class="sd"> The equinox to precess these coordinates to.</span> <span class="sd"> Returns</span> <span class="sd"> -------</span> <span class="sd"> newcoord : FK4Coordinates</span> <span class="sd"> The new coordinate</span> <span class="sd"> """</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">_precession_matrix_besselian</span> <span class="n">pmat</span> <span class="o">=</span> <span class="n">_precession_matrix_besselian</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span><span class="o">.</span><span class="n">byear</span><span class="p">,</span> <span class="n">newequinox</span><span class="o">.</span><span class="n">byear</span><span class="p">)</span> <span class="n">v</span> <span class="o">=</span> <span class="p">[</span><span class="bp">self</span><span class="o">.</span><span class="n">x</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">y</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">z</span><span class="p">]</span> <span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="p">,</span> <span class="n">z</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">dot</span><span class="p">(</span><span class="n">pmat</span><span class="o">.</span><span class="n">A</span><span class="p">,</span> <span class="n">v</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="o">=</span><span class="n">y</span><span class="p">,</span> <span class="n">z</span><span class="o">=</span><span class="n">z</span><span class="p">,</span> <span class="n">unit</span><span class="o">=</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">newequinox</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="o">=</span><span class="n">y</span><span class="p">,</span> <span class="n">z</span><span class="o">=</span><span class="n">z</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">newequinox</span><span class="p">)</span> </div></div> <span class="nd">@transformations.coordinate_alias</span><span class="p">(</span><span class="s">'fk4_no_e'</span><span class="p">)</span> <div class="viewcode-block" id="FK4NoETermCoordinates"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4NoETermCoordinates.html#astropy.coordinates.builtin_systems.FK4NoETermCoordinates">[docs]</a><span class="k">class</span> <span class="nc">FK4NoETermCoordinates</span><span class="p">(</span><span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> A coordinate in the FK4 system.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> {params}</span> <span class="sd"> equinox : `~astropy.time.Time`, optional</span> <span class="sd"> The equinox for these coordinates. Defaults to B1950.</span> <span class="sd"> obstime : `~astropy.time.Time` or None</span> <span class="sd"> The time of observation for this coordinate. If None, it will be taken</span> <span class="sd"> to be the same as the `equinox`.</span> <span class="sd"> Alternatively, a single argument that is any kind of spherical coordinate</span> <span class="sd"> can be provided, and will be converted to `FK4NoETermCoordinates` and used as this</span> <span class="sd"> coordinate.</span> <span class="sd"> """</span> <span class="n">__doc__</span> <span class="o">=</span> <span class="n">__doc__</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">params</span><span class="o">=</span><span class="n">SphericalCoordinatesBase</span><span class="o">.</span><span class="n">_init_docstring_param_templ</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">lonnm</span><span class="o">=</span><span class="s">'ra'</span><span class="p">,</span> <span class="n">latnm</span><span class="o">=</span><span class="s">'dec'</span><span class="p">))</span> <span class="k">def</span> <span class="nf">__init__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="o">*</span><span class="n">args</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span> <span class="nb">super</span><span class="p">(</span><span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">__init__</span><span class="p">()</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'equinox'</span><span class="p">,</span> <span class="n">_EQUINOX_B1950</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'obstime'</span><span class="p">,</span> <span class="bp">None</span><span class="p">)</span> <span class="k">if</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified equinox is not a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified obstime is not None or a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="o">==</span> <span class="mi">1</span> <span class="ow">and</span> <span class="nb">len</span><span class="p">(</span><span class="n">kwargs</span><span class="p">)</span> <span class="o">==</span> <span class="mi">0</span> <span class="ow">and</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="n">newcoord</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">.</span><span class="n">transform_to</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">ra</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">dec</span> <span class="bp">self</span><span class="o">.</span><span class="n">_distance</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">_distance</span> <span class="k">else</span><span class="p">:</span> <span class="nb">super</span><span class="p">(</span><span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">_initialize_latlon</span><span class="p">(</span><span class="s">'ra'</span><span class="p">,</span> <span class="s">'dec'</span><span class="p">,</span> <span class="bp">True</span><span class="p">,</span> <span class="n">args</span><span class="p">,</span> <span class="n">kwargs</span><span class="p">)</span> <span class="k">def</span> <span class="nf">__repr__</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">', Distance={0:.2g} {1!s}'</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_value</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">''</span> <span class="n">msg</span> <span class="o">=</span> <span class="s">"<{0} RA={1:.5f} deg, Dec={2:.5f} deg{3}>"</span> <span class="k">return</span> <span class="n">msg</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="o">.</span><span class="n">__name__</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="n">diststr</span><span class="p">)</span> <span class="nd">@property</span> <div class="viewcode-block" id="FK4NoETermCoordinates.lonangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4NoETermCoordinates.html#astropy.coordinates.builtin_systems.FK4NoETermCoordinates.lonangle">[docs]</a> <span class="k">def</span> <span class="nf">lonangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">ra</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK4NoETermCoordinates.latangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4NoETermCoordinates.html#astropy.coordinates.builtin_systems.FK4NoETermCoordinates.latangle">[docs]</a> <span class="k">def</span> <span class="nf">latangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">dec</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK4NoETermCoordinates.equinox"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4NoETermCoordinates.html#astropy.coordinates.builtin_systems.FK4NoETermCoordinates.equinox">[docs]</a> <span class="k">def</span> <span class="nf">equinox</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="FK4NoETermCoordinates.obstime"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4NoETermCoordinates.html#astropy.coordinates.builtin_systems.FK4NoETermCoordinates.obstime">[docs]</a> <span class="k">def</span> <span class="nf">obstime</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">equinox</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> </div> <div class="viewcode-block" id="FK4NoETermCoordinates.precess_to"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.FK4NoETermCoordinates.html#astropy.coordinates.builtin_systems.FK4NoETermCoordinates.precess_to">[docs]</a> <span class="k">def</span> <span class="nf">precess_to</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="n">newequinox</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> Precesses the coordinates from their current `equinox` to a new equinox.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> newequinox : `~astropy.time.Time`</span> <span class="sd"> The equinox to precess these coordinates to.</span> <span class="sd"> Returns</span> <span class="sd"> -------</span> <span class="sd"> newcoord : FK4NoETermCoordinates</span> <span class="sd"> The new coordinate</span> <span class="sd"> """</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">_precession_matrix_besselian</span> <span class="n">pmat</span> <span class="o">=</span> <span class="n">_precession_matrix_besselian</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span><span class="o">.</span><span class="n">byear</span><span class="p">,</span> <span class="n">newequinox</span><span class="o">.</span><span class="n">byear</span><span class="p">)</span> <span class="n">v</span> <span class="o">=</span> <span class="p">[</span><span class="bp">self</span><span class="o">.</span><span class="n">x</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">y</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">z</span><span class="p">]</span> <span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="p">,</span> <span class="n">z</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">dot</span><span class="p">(</span><span class="n">pmat</span><span class="o">.</span><span class="n">A</span><span class="p">,</span> <span class="n">v</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="o">=</span><span class="n">y</span><span class="p">,</span> <span class="n">z</span><span class="o">=</span><span class="n">z</span><span class="p">,</span> <span class="n">unit</span><span class="o">=</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">newequinox</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">x</span><span class="p">,</span> <span class="n">y</span><span class="o">=</span><span class="n">y</span><span class="p">,</span> <span class="n">z</span><span class="o">=</span><span class="n">z</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">newequinox</span><span class="p">)</span> </div></div> <span class="nd">@transformations.coordinate_alias</span><span class="p">(</span><span class="s">'galactic'</span><span class="p">)</span> <div class="viewcode-block" id="GalacticCoordinates"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.GalacticCoordinates.html#astropy.coordinates.builtin_systems.GalacticCoordinates">[docs]</a><span class="k">class</span> <span class="nc">GalacticCoordinates</span><span class="p">(</span><span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> A coordinate in Galactic Coordinates.</span> <span class="sd"> .. note::</span> <span class="sd"> Transformations from Galactic Coordinates to other systems are</span> <span class="sd"> not well-defined because of ambiguities in the definition of</span> <span class="sd"> Galactic Coordinates. See</span> <span class="sd"> `Lie et al. 2011 <http://dx.doi.org/10.1051/0004-6361/201014961>`</span> <span class="sd"> for more details on this. Here, we use the</span> <span class="sd"> `Reid & Brunthaler 2004 <http://dx.doi.org/10.1086/424960>`</span> <span class="sd"> definition for converting to/from FK5, and assume the IAU</span> <span class="sd"> definition applies for converting to FK4 *without* e-terms.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> {params}</span> <span class="sd"> obstime : `~astropy.time.Time` or None</span> <span class="sd"> The time of observation for this coordinate. If None, it will be taken</span> <span class="sd"> to be the same as the `equinox`.</span> <span class="sd"> Alternatively, a single argument that is any kind of spherical coordinate</span> <span class="sd"> can be provided, and will be converted to `GalacticCoordinates` and</span> <span class="sd"> used as this coordinate.</span> <span class="sd"> """</span> <span class="n">__doc__</span> <span class="o">=</span> <span class="n">__doc__</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">params</span><span class="o">=</span><span class="n">SphericalCoordinatesBase</span><span class="o">.</span><span class="n">_init_docstring_param_templ</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">lonnm</span><span class="o">=</span><span class="s">'l'</span><span class="p">,</span> <span class="n">latnm</span><span class="o">=</span><span class="s">'b'</span><span class="p">))</span> <span class="c"># North galactic pole and zeropoint of l in FK4/FK5 coordinates. Needed for</span> <span class="c"># transformations to/from FK4/5</span> <span class="c"># These are from Reid & Brunthaler 2004</span> <span class="n">_ngp_J2000</span> <span class="o">=</span> <span class="n">FK5Coordinates</span><span class="p">(</span><span class="mf">192.859508</span><span class="p">,</span> <span class="mf">27.128336</span><span class="p">,</span> <span class="n">unit</span><span class="o">=</span><span class="p">(</span><span class="n">u</span><span class="o">.</span><span class="n">degree</span><span class="p">,</span> <span class="n">u</span><span class="o">.</span><span class="n">degree</span><span class="p">))</span> <span class="n">_lon0_J2000</span> <span class="o">=</span> <span class="n">Angle</span><span class="p">(</span><span class="mf">122.932</span><span class="p">,</span> <span class="n">unit</span><span class="o">=</span><span class="n">u</span><span class="o">.</span><span class="n">degree</span><span class="p">)</span> <span class="c"># These are from the IAU's definition of galactic coordinates</span> <span class="n">_ngp_B1950</span> <span class="o">=</span> <span class="n">FK4Coordinates</span><span class="p">(</span><span class="mf">192.25</span><span class="p">,</span> <span class="mf">27.4</span><span class="p">,</span> <span class="n">unit</span><span class="o">=</span><span class="p">(</span><span class="n">u</span><span class="o">.</span><span class="n">degree</span><span class="p">,</span> <span class="n">u</span><span class="o">.</span><span class="n">degree</span><span class="p">))</span> <span class="n">_lon0_B1950</span> <span class="o">=</span> <span class="n">Angle</span><span class="p">(</span><span class="mi">123</span><span class="p">,</span> <span class="n">unit</span><span class="o">=</span><span class="n">u</span><span class="o">.</span><span class="n">degree</span><span class="p">)</span> <span class="k">def</span> <span class="nf">__init__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="o">*</span><span class="n">args</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span> <span class="nb">super</span><span class="p">(</span><span class="n">GalacticCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">__init__</span><span class="p">()</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'obstime'</span><span class="p">,</span> <span class="bp">None</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified obstime is not None or a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="o">==</span> <span class="mi">1</span> <span class="ow">and</span> <span class="nb">len</span><span class="p">(</span><span class="n">kwargs</span><span class="p">)</span> <span class="o">==</span> <span class="mi">0</span> <span class="ow">and</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="n">newcoord</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">.</span><span class="n">transform_to</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">l</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">l</span> <span class="bp">self</span><span class="o">.</span><span class="n">b</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">b</span> <span class="bp">self</span><span class="o">.</span><span class="n">_distance</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">_distance</span> <span class="k">else</span><span class="p">:</span> <span class="nb">super</span><span class="p">(</span><span class="n">GalacticCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">_initialize_latlon</span><span class="p">(</span><span class="s">'l'</span><span class="p">,</span> <span class="s">'b'</span><span class="p">,</span> <span class="bp">False</span><span class="p">,</span> <span class="n">args</span><span class="p">,</span> <span class="n">kwargs</span><span class="p">)</span> <span class="k">def</span> <span class="nf">__repr__</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">', Distance={0:.2g} {1!s}'</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_value</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">''</span> <span class="n">msg</span> <span class="o">=</span> <span class="s">"<{0} l={1:.5f} deg, b={2:.5f} deg{3}>"</span> <span class="k">return</span> <span class="n">msg</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="o">.</span><span class="n">__name__</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">l</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">b</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="n">diststr</span><span class="p">)</span> <span class="nd">@property</span> <div class="viewcode-block" id="GalacticCoordinates.lonangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.GalacticCoordinates.html#astropy.coordinates.builtin_systems.GalacticCoordinates.lonangle">[docs]</a> <span class="k">def</span> <span class="nf">lonangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">l</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="GalacticCoordinates.latangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.GalacticCoordinates.html#astropy.coordinates.builtin_systems.GalacticCoordinates.latangle">[docs]</a> <span class="k">def</span> <span class="nf">latangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">b</span> </div></div> <span class="nd">@transformations.coordinate_alias</span><span class="p">(</span><span class="s">'horizontal'</span><span class="p">)</span> <div class="viewcode-block" id="HorizontalCoordinates"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.HorizontalCoordinates.html#astropy.coordinates.builtin_systems.HorizontalCoordinates">[docs]</a><span class="k">class</span> <span class="nc">HorizontalCoordinates</span><span class="p">(</span><span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> A coordinate in the Horizontal or "az/el" system.</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> {params}</span> <span class="sd"> equinox : `~astropy.time.Time`, optional</span> <span class="sd"> The equinox for these coordinates. Defaults to J2000.</span> <span class="sd"> obstime : `~astropy.time.Time` or None</span> <span class="sd"> The time of observation for this coordinate. If None, it will be taken</span> <span class="sd"> to be the same as the `equinox`.</span> <span class="sd"> Alternatively, a single argument that is any kind of spherical coordinate</span> <span class="sd"> can be provided, and will be converted to `HorizontalCoordinates` and used</span> <span class="sd"> as this coordinate.</span> <span class="sd"> """</span> <span class="n">__doc__</span> <span class="o">=</span> <span class="n">__doc__</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">params</span><span class="o">=</span><span class="n">SphericalCoordinatesBase</span><span class="o">.</span><span class="n">_init_docstring_param_templ</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="n">lonnm</span><span class="o">=</span><span class="s">'az'</span><span class="p">,</span> <span class="n">latnm</span><span class="o">=</span><span class="s">'el'</span><span class="p">))</span> <span class="k">def</span> <span class="nf">__init__</span><span class="p">(</span><span class="bp">self</span><span class="p">,</span> <span class="o">*</span><span class="n">args</span><span class="p">,</span> <span class="o">**</span><span class="n">kwargs</span><span class="p">):</span> <span class="nb">super</span><span class="p">(</span><span class="n">HorizontalCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">__init__</span><span class="p">()</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'equinox'</span><span class="p">,</span> <span class="n">_EQUINOX_J2000</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="o">=</span> <span class="n">kwargs</span><span class="o">.</span><span class="n">pop</span><span class="p">(</span><span class="s">'obstime'</span><span class="p">,</span> <span class="bp">None</span><span class="p">)</span> <span class="k">if</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified equinox is not a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span> <span class="ow">and</span> <span class="ow">not</span> <span class="nb">isinstance</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span><span class="p">,</span> <span class="n">Time</span><span class="p">):</span> <span class="k">raise</span> <span class="ne">TypeError</span><span class="p">(</span><span class="s">'specified obstime is not None or a Time object'</span><span class="p">)</span> <span class="k">if</span> <span class="nb">len</span><span class="p">(</span><span class="n">args</span><span class="p">)</span> <span class="o">==</span> <span class="mi">1</span> <span class="ow">and</span> <span class="nb">len</span><span class="p">(</span><span class="n">kwargs</span><span class="p">)</span> <span class="o">==</span> <span class="mi">0</span> <span class="ow">and</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">SphericalCoordinatesBase</span><span class="p">):</span> <span class="n">newcoord</span> <span class="o">=</span> <span class="n">args</span><span class="p">[</span><span class="mi">0</span><span class="p">]</span><span class="o">.</span><span class="n">transform_to</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="p">)</span> <span class="bp">self</span><span class="o">.</span><span class="n">az</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">az</span> <span class="bp">self</span><span class="o">.</span><span class="n">el</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">el</span> <span class="bp">self</span><span class="o">.</span><span class="n">_distance</span> <span class="o">=</span> <span class="n">newcoord</span><span class="o">.</span><span class="n">_distance</span> <span class="k">else</span><span class="p">:</span> <span class="nb">super</span><span class="p">(</span><span class="n">HorizontalCoordinates</span><span class="p">,</span> <span class="bp">self</span><span class="p">)</span><span class="o">.</span><span class="n">_initialize_latlon</span><span class="p">(</span><span class="s">'az'</span><span class="p">,</span> <span class="s">'el'</span><span class="p">,</span> <span class="bp">False</span><span class="p">,</span> <span class="n">args</span><span class="p">,</span> <span class="n">kwargs</span><span class="p">)</span> <span class="k">def</span> <span class="nf">__repr__</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="ow">not</span> <span class="bp">None</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">', Distance={0:.2g} {1!s}'</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_value</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span><span class="p">)</span> <span class="k">else</span><span class="p">:</span> <span class="n">diststr</span> <span class="o">=</span> <span class="s">''</span> <span class="n">msg</span> <span class="o">=</span> <span class="s">"<{0} az={1:.5f} deg, el={2:.5f} deg{3}>"</span> <span class="k">return</span> <span class="n">msg</span><span class="o">.</span><span class="n">format</span><span class="p">(</span><span class="bp">self</span><span class="o">.</span><span class="n">__class__</span><span class="o">.</span><span class="n">__name__</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">az</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="bp">self</span><span class="o">.</span><span class="n">el</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="n">diststr</span><span class="p">)</span> <span class="nd">@property</span> <div class="viewcode-block" id="HorizontalCoordinates.lonangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.HorizontalCoordinates.html#astropy.coordinates.builtin_systems.HorizontalCoordinates.lonangle">[docs]</a> <span class="k">def</span> <span class="nf">lonangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">az</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="HorizontalCoordinates.latangle"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.HorizontalCoordinates.html#astropy.coordinates.builtin_systems.HorizontalCoordinates.latangle">[docs]</a> <span class="k">def</span> <span class="nf">latangle</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">el</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="HorizontalCoordinates.equinox"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.HorizontalCoordinates.html#astropy.coordinates.builtin_systems.HorizontalCoordinates.equinox">[docs]</a> <span class="k">def</span> <span class="nf">equinox</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_equinox</span> </div> <span class="nd">@property</span> <div class="viewcode-block" id="HorizontalCoordinates.obstime"><a class="viewcode-back" href="../../../_generated/astropy.coordinates.builtin_systems.HorizontalCoordinates.html#astropy.coordinates.builtin_systems.HorizontalCoordinates.obstime">[docs]</a> <span class="k">def</span> <span class="nf">obstime</span><span class="p">(</span><span class="bp">self</span><span class="p">):</span> <span class="k">if</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="ow">is</span> <span class="bp">None</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">equinox</span> <span class="k">else</span><span class="p">:</span> <span class="k">return</span> <span class="bp">self</span><span class="o">.</span><span class="n">_obstime</span> <span class="c">#<--------------------------------transformations------------------------------></span> <span class="c"># ICRS to/from FK5</span></div></div> <span class="nd">@transformations.static_transform_matrix</span><span class="p">(</span><span class="n">ICRSCoordinates</span><span class="p">,</span> <span class="n">FK5Coordinates</span><span class="p">)</span> <span class="k">def</span> <span class="nf">icrs_to_fk5</span><span class="p">():</span> <span class="sd">"""</span> <span class="sd"> B-matrix from USNO circular 179</span> <span class="sd"> """</span> <span class="kn">from</span> <span class="nn">.angles</span> <span class="kn">import</span> <span class="n">rotation_matrix</span> <span class="n">eta0</span> <span class="o">=</span> <span class="o">-</span><span class="mf">19.9</span> <span class="o">/</span> <span class="mf">3600000.</span> <span class="n">xi0</span> <span class="o">=</span> <span class="mf">9.1</span> <span class="o">/</span> <span class="mf">3600000.</span> <span class="n">da0</span> <span class="o">=</span> <span class="o">-</span><span class="mf">22.9</span> <span class="o">/</span> <span class="mf">3600000.</span> <span class="n">m1</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="o">-</span><span class="n">eta0</span><span class="p">,</span> <span class="s">'x'</span><span class="p">)</span> <span class="n">m2</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="n">xi0</span><span class="p">,</span> <span class="s">'y'</span><span class="p">)</span> <span class="n">m3</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="n">da0</span><span class="p">,</span> <span class="s">'z'</span><span class="p">)</span> <span class="k">return</span> <span class="n">m1</span> <span class="o">*</span> <span class="n">m2</span> <span class="o">*</span> <span class="n">m3</span> <span class="c"># can't be static because the equinox is needed</span> <span class="nd">@transformations.dynamic_transform_matrix</span><span class="p">(</span><span class="n">FK5Coordinates</span><span class="p">,</span> <span class="n">ICRSCoordinates</span><span class="p">)</span> <span class="k">def</span> <span class="nf">fk5_to_icrs</span><span class="p">(</span><span class="n">fk5c</span><span class="p">):</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">_precess_from_J2000_Capitaine</span> <span class="n">pmat</span> <span class="o">=</span> <span class="n">_precess_from_J2000_Capitaine</span><span class="p">(</span><span class="n">fk5c</span><span class="o">.</span><span class="n">equinox</span><span class="o">.</span><span class="n">jyear</span><span class="p">)</span><span class="o">.</span><span class="n">T</span> <span class="c"># transpose gets equinox -> J2000</span> <span class="n">fk5toicrsmat</span> <span class="o">=</span> <span class="n">icrs_to_fk5</span><span class="p">()</span><span class="o">.</span><span class="n">T</span> <span class="k">return</span> <span class="n">fk5toicrsmat</span> <span class="o">*</span> <span class="n">pmat</span> <span class="c"># FK4-NO-E to/from FK4</span> <span class="c"># In the present framework, we include two coordinate classes for FK4</span> <span class="c"># coordinates - one including the E-terms of aberration (FK4Coordinates), and</span> <span class="c"># one not including them (FK4NoETermCoordinates). In the following functions,</span> <span class="c"># we describe the transformation between these two.</span> <span class="k">def</span> <span class="nf">fk4_e_terms</span><span class="p">(</span><span class="n">equinox</span><span class="p">):</span> <span class="sd">"""</span> <span class="sd"> Return the e-terms of aberation vector</span> <span class="sd"> Parameters</span> <span class="sd"> ----------</span> <span class="sd"> equinox : Time object</span> <span class="sd"> The equinox for which to compute the e-terms</span> <span class="sd"> """</span> <span class="kn">from</span> <span class="nn">.</span> <span class="kn">import</span> <span class="n">earth_orientation</span> <span class="k">as</span> <span class="n">earth</span> <span class="c"># Constant of aberration at J2000</span> <span class="n">k</span> <span class="o">=</span> <span class="mf">0.0056932</span> <span class="c"># Eccentricity of the Earth's orbit</span> <span class="n">e</span> <span class="o">=</span> <span class="n">earth</span><span class="o">.</span><span class="n">eccentricity</span><span class="p">(</span><span class="n">equinox</span><span class="o">.</span><span class="n">jd</span><span class="p">)</span> <span class="n">e</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">radians</span><span class="p">(</span><span class="n">e</span><span class="p">)</span> <span class="c"># Mean longitude of perigee of the solar orbit</span> <span class="n">g</span> <span class="o">=</span> <span class="n">earth</span><span class="o">.</span><span class="n">mean_lon_of_perigee</span><span class="p">(</span><span class="n">equinox</span><span class="o">.</span><span class="n">jd</span><span class="p">)</span> <span class="n">g</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">radians</span><span class="p">(</span><span class="n">g</span><span class="p">)</span> <span class="c"># Obliquity of the ecliptic</span> <span class="n">o</span> <span class="o">=</span> <span class="n">earth</span><span class="o">.</span><span class="n">obliquity</span><span class="p">(</span><span class="n">equinox</span><span class="o">.</span><span class="n">jd</span><span class="p">,</span> <span class="n">algorithm</span><span class="o">=</span><span class="mi">1980</span><span class="p">)</span> <span class="n">o</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">radians</span><span class="p">(</span><span class="n">o</span><span class="p">)</span> <span class="k">return</span> <span class="n">e</span> <span class="o">*</span> <span class="n">k</span> <span class="o">*</span> <span class="n">np</span><span class="o">.</span><span class="n">sin</span><span class="p">(</span><span class="n">g</span><span class="p">),</span> \ <span class="o">-</span><span class="n">e</span> <span class="o">*</span> <span class="n">k</span> <span class="o">*</span> <span class="n">np</span><span class="o">.</span><span class="n">cos</span><span class="p">(</span><span class="n">g</span><span class="p">)</span> <span class="o">*</span> <span class="n">np</span><span class="o">.</span><span class="n">cos</span><span class="p">(</span><span class="n">o</span><span class="p">),</span> \ <span class="o">-</span><span class="n">e</span> <span class="o">*</span> <span class="n">k</span> <span class="o">*</span> <span class="n">np</span><span class="o">.</span><span class="n">cos</span><span class="p">(</span><span class="n">g</span><span class="p">)</span> <span class="o">*</span> <span class="n">np</span><span class="o">.</span><span class="n">sin</span><span class="p">(</span><span class="n">o</span><span class="p">)</span> <span class="nd">@transformations.transform_function</span><span class="p">(</span><span class="n">FK4Coordinates</span><span class="p">,</span> <span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="n">priority</span><span class="o">=</span><span class="mi">1</span><span class="p">)</span> <span class="k">def</span> <span class="nf">fk4_to_fk4_no_e</span><span class="p">(</span><span class="n">fk4c</span><span class="p">):</span> <span class="c"># Extract cartesian vector</span> <span class="n">r</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">([</span><span class="n">fk4c</span><span class="o">.</span><span class="n">x</span><span class="p">,</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">y</span><span class="p">,</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">z</span><span class="p">])</span> <span class="c"># Find distance (for re-normalization)</span> <span class="n">d_orig</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">sum</span><span class="p">(</span><span class="n">r</span> <span class="o">**</span> <span class="mi">2</span><span class="p">))</span> <span class="c"># Apply E-terms of aberration</span> <span class="n">eterms_a</span> <span class="o">=</span> <span class="n">fk4_e_terms</span><span class="p">(</span><span class="n">fk4c</span><span class="o">.</span><span class="n">equinox</span><span class="p">)</span> <span class="n">r</span> <span class="o">=</span> <span class="n">r</span> <span class="o">-</span> <span class="n">eterms_a</span> <span class="o">+</span> <span class="n">np</span><span class="o">.</span><span class="n">dot</span><span class="p">(</span><span class="n">r</span><span class="p">,</span> <span class="n">eterms_a</span><span class="p">)</span> <span class="o">*</span> <span class="n">r</span> <span class="c"># Find new distance (for re-normalization)</span> <span class="n">d_new</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">sum</span><span class="p">(</span><span class="n">r</span> <span class="o">**</span> <span class="mi">2</span><span class="p">))</span> <span class="c"># Renormalize</span> <span class="n">r</span> <span class="o">=</span> <span class="n">r</span> <span class="o">*</span> <span class="n">d_orig</span> <span class="o">/</span> <span class="n">d_new</span> <span class="n">unit</span> <span class="o">=</span> <span class="bp">None</span> <span class="k">if</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="bp">None</span> <span class="k">else</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span> <span class="n">result</span> <span class="o">=</span> <span class="n">FK4NoETermCoordinates</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">r</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">y</span><span class="o">=</span><span class="n">r</span><span class="p">[</span><span class="mi">1</span><span class="p">],</span> <span class="n">z</span><span class="o">=</span><span class="n">r</span><span class="p">[</span><span class="mi">2</span><span class="p">],</span> <span class="n">unit</span><span class="o">=</span><span class="n">unit</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">fk4c</span><span class="o">.</span><span class="n">equinox</span><span class="p">)</span> <span class="k">return</span> <span class="n">result</span> <span class="nd">@transformations.transform_function</span><span class="p">(</span><span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="n">FK4Coordinates</span><span class="p">,</span> <span class="n">priority</span><span class="o">=</span><span class="mi">1</span><span class="p">)</span> <span class="k">def</span> <span class="nf">fk4_no_e_to_fk4</span><span class="p">(</span><span class="n">fk4c</span><span class="p">):</span> <span class="c"># Extract cartesian vector</span> <span class="n">r</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">array</span><span class="p">([</span><span class="n">fk4c</span><span class="o">.</span><span class="n">x</span><span class="p">,</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">y</span><span class="p">,</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">z</span><span class="p">])</span> <span class="c"># Find distance (for re-normalization)</span> <span class="n">d_orig</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">sum</span><span class="p">(</span><span class="n">r</span> <span class="o">**</span> <span class="mi">2</span><span class="p">))</span> <span class="c"># Apply E-terms of aberration</span> <span class="n">eterms_a</span> <span class="o">=</span> <span class="n">fk4_e_terms</span><span class="p">(</span><span class="n">fk4c</span><span class="o">.</span><span class="n">equinox</span><span class="p">)</span> <span class="n">r0</span> <span class="o">=</span> <span class="n">r</span><span class="o">.</span><span class="n">copy</span><span class="p">()</span> <span class="k">for</span> <span class="n">j</span> <span class="ow">in</span> <span class="nb">range</span><span class="p">(</span><span class="mi">10</span><span class="p">):</span> <span class="n">r</span> <span class="o">=</span> <span class="p">(</span><span class="n">r0</span> <span class="o">+</span> <span class="n">eterms_a</span><span class="p">)</span> <span class="o">/</span> <span class="p">(</span><span class="mf">1.</span> <span class="o">+</span> <span class="n">np</span><span class="o">.</span><span class="n">dot</span><span class="p">(</span><span class="n">r</span><span class="p">,</span> <span class="n">eterms_a</span><span class="p">))</span> <span class="c"># Find new distance (for re-normalization)</span> <span class="n">d_new</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">sqrt</span><span class="p">(</span><span class="n">np</span><span class="o">.</span><span class="n">sum</span><span class="p">(</span><span class="n">r</span> <span class="o">**</span> <span class="mi">2</span><span class="p">))</span> <span class="c"># Renormalize</span> <span class="n">r</span> <span class="o">=</span> <span class="n">r</span> <span class="o">*</span> <span class="n">d_orig</span> <span class="o">/</span> <span class="n">d_new</span> <span class="n">unit</span> <span class="o">=</span> <span class="bp">None</span> <span class="k">if</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">distance</span> <span class="ow">is</span> <span class="bp">None</span> <span class="k">else</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">distance</span><span class="o">.</span><span class="n">_unit</span> <span class="n">result</span> <span class="o">=</span> <span class="n">FK4Coordinates</span><span class="p">(</span><span class="n">x</span><span class="o">=</span><span class="n">r</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="n">y</span><span class="o">=</span><span class="n">r</span><span class="p">[</span><span class="mi">1</span><span class="p">],</span> <span class="n">z</span><span class="o">=</span><span class="n">r</span><span class="p">[</span><span class="mi">2</span><span class="p">],</span> <span class="n">unit</span><span class="o">=</span><span class="n">unit</span><span class="p">,</span> <span class="n">equinox</span><span class="o">=</span><span class="n">fk4c</span><span class="o">.</span><span class="n">equinox</span><span class="p">)</span> <span class="k">return</span> <span class="n">result</span> <span class="c"># FK5 to/from FK4</span> <span class="c"># These transformations are very slightly prioritized >1 (lower priority number means</span> <span class="c"># better path) to prefer the FK5 path over FK4 when possible</span> <span class="c">#can't be static because the equinox is needed</span> <span class="c"># B1950->J2000 matrix from Murray 1989 A&A 218,325 eqn 28</span> <span class="n">B1950_TO_J2000_M</span> <span class="o">=</span> \ <span class="n">np</span><span class="o">.</span><span class="n">mat</span><span class="p">([[</span><span class="mf">0.9999256794956877</span><span class="p">,</span> <span class="o">-</span><span class="mf">0.0111814832204662</span><span class="p">,</span> <span class="o">-</span><span class="mf">0.0048590038153592</span><span class="p">],</span> <span class="p">[</span><span class="mf">0.0111814832391717</span><span class="p">,</span> <span class="mf">0.9999374848933135</span><span class="p">,</span> <span class="o">-</span><span class="mf">0.0000271625947142</span><span class="p">],</span> <span class="p">[</span><span class="mf">0.0048590037723143</span><span class="p">,</span> <span class="o">-</span><span class="mf">0.0000271702937440</span><span class="p">,</span> <span class="mf">0.9999881946023742</span><span class="p">]])</span> <span class="n">FK4_CORR</span> <span class="o">=</span> \ <span class="n">np</span><span class="o">.</span><span class="n">mat</span><span class="p">([[</span><span class="o">-</span><span class="mf">0.0026455262</span><span class="p">,</span> <span class="o">-</span><span class="mf">1.1539918689</span><span class="p">,</span> <span class="o">+</span><span class="mf">2.1111346190</span><span class="p">],</span> <span class="p">[</span><span class="o">+</span><span class="mf">1.1540628161</span><span class="p">,</span> <span class="o">-</span><span class="mf">0.0129042997</span><span class="p">,</span> <span class="o">+</span><span class="mf">0.0236021478</span><span class="p">],</span> <span class="p">[</span><span class="o">-</span><span class="mf">2.1112979048</span><span class="p">,</span> <span class="o">-</span><span class="mf">0.0056024448</span><span class="p">,</span> <span class="o">+</span><span class="mf">0.0102587734</span><span class="p">]])</span> <span class="o">*</span> <span class="mf">1.e-6</span> <span class="c"># This transformation can't be static because the observation date is needed.</span> <span class="nd">@transformations.dynamic_transform_matrix</span><span class="p">(</span><span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="n">FK5Coordinates</span><span class="p">,</span> <span class="n">priority</span><span class="o">=</span><span class="mi">1</span><span class="p">)</span> <span class="k">def</span> <span class="nf">fk4_no_e_to_fk5</span><span class="p">(</span><span class="n">fk4c</span><span class="p">,</span> <span class="n">skip_precession</span><span class="o">=</span><span class="bp">False</span><span class="p">):</span> <span class="c"># Add in correction terms for FK4 rotating system - Murray 89 eqn 29</span> <span class="c"># Note this is *julian century*, not besselian</span> <span class="n">T</span> <span class="o">=</span> <span class="p">(</span><span class="n">fk4c</span><span class="o">.</span><span class="n">obstime</span><span class="o">.</span><span class="n">jyear</span> <span class="o">-</span> <span class="mf">1950.</span><span class="p">)</span> <span class="o">/</span> <span class="mf">100.</span> <span class="n">B</span> <span class="o">=</span> <span class="n">B1950_TO_J2000_M</span> <span class="o">+</span> <span class="n">FK4_CORR</span> <span class="o">*</span> <span class="n">T</span> <span class="c"># If equinox is not B1950, need to precess first</span> <span class="k">if</span> <span class="n">skip_precession</span> <span class="ow">or</span> <span class="n">fk4c</span><span class="o">.</span><span class="n">equinox</span> <span class="o">==</span> <span class="n">_EQUINOX_B1950</span><span class="p">:</span> <span class="k">return</span> <span class="n">B</span> <span class="k">else</span><span class="p">:</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">_precession_matrix_besselian</span> <span class="k">return</span> <span class="n">B</span> <span class="o">*</span> <span class="n">_precession_matrix_besselian</span><span class="p">(</span><span class="n">fk4c</span><span class="o">.</span><span class="n">equinox</span><span class="o">.</span><span class="n">byear</span><span class="p">,</span> <span class="mi">1950</span><span class="p">)</span> <span class="c"># This transformation can't be static because the observation date is needed.</span> <span class="nd">@transformations.dynamic_transform_matrix</span><span class="p">(</span><span class="n">FK5Coordinates</span><span class="p">,</span> <span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="n">priority</span><span class="o">=</span><span class="mi">1</span><span class="p">)</span> <span class="k">def</span> <span class="nf">fk5_to_fk4_no_e</span><span class="p">(</span><span class="n">fk5c</span><span class="p">):</span> <span class="c"># Get transposed matrix from FK4 -> FK5 assuming equinox B1950 -> J2000</span> <span class="n">B</span> <span class="o">=</span> <span class="n">fk4_no_e_to_fk5</span><span class="p">(</span><span class="n">fk5c</span><span class="p">,</span> <span class="n">skip_precession</span><span class="o">=</span><span class="bp">True</span><span class="p">)</span><span class="o">.</span><span class="n">T</span> <span class="c"># If equinox is not B1950, need to precess first</span> <span class="k">if</span> <span class="n">fk5c</span><span class="o">.</span><span class="n">equinox</span> <span class="o">==</span> <span class="n">_EQUINOX_J2000</span><span class="p">:</span> <span class="k">return</span> <span class="n">B</span> <span class="k">else</span><span class="p">:</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">precession_matrix_Capitaine</span> <span class="k">return</span> <span class="n">B</span> <span class="o">*</span> <span class="n">precession_matrix_Capitaine</span><span class="p">(</span><span class="n">fk5c</span><span class="o">.</span><span class="n">equinox</span><span class="p">,</span> <span class="n">_EQUINOX_J2000</span><span class="p">)</span> <span class="c"># GalacticCoordinates to/from FK4/FK5</span> <span class="c"># can't be static because the equinox is needed</span> <span class="nd">@transformations.dynamic_transform_matrix</span><span class="p">(</span><span class="n">FK5Coordinates</span><span class="p">,</span> <span class="n">GalacticCoordinates</span><span class="p">)</span> <span class="k">def</span> <span class="nf">_fk5_to_gal</span><span class="p">(</span><span class="n">fk5coords</span><span class="p">):</span> <span class="kn">from</span> <span class="nn">.angles</span> <span class="kn">import</span> <span class="n">rotation_matrix</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">_precess_from_J2000_Capitaine</span> <span class="c"># needed mainly to support inverse from galactic</span> <span class="n">jequinox</span> <span class="o">=</span> <span class="mi">2000</span> <span class="k">if</span> <span class="n">fk5coords</span><span class="o">.</span><span class="n">equinox</span> <span class="ow">is</span> <span class="bp">None</span> <span class="k">else</span> <span class="n">fk5coords</span><span class="o">.</span><span class="n">equinox</span><span class="o">.</span><span class="n">jyear</span> <span class="n">mat1</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="mi">180</span> <span class="o">-</span> <span class="n">GalacticCoordinates</span><span class="o">.</span><span class="n">_lon0_J2000</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="s">'z'</span><span class="p">)</span> <span class="n">mat2</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="mi">90</span> <span class="o">-</span> <span class="n">GalacticCoordinates</span><span class="o">.</span><span class="n">_ngp_J2000</span><span class="o">.</span><span class="n">dec</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="s">'y'</span><span class="p">)</span> <span class="n">mat3</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="n">GalacticCoordinates</span><span class="o">.</span><span class="n">_ngp_J2000</span><span class="o">.</span><span class="n">ra</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="s">'z'</span><span class="p">)</span> <span class="c"># transpose gets equinox -> J2000</span> <span class="n">matprec</span> <span class="o">=</span> <span class="n">_precess_from_J2000_Capitaine</span><span class="p">(</span><span class="n">jequinox</span><span class="p">)</span><span class="o">.</span><span class="n">T</span> <span class="k">return</span> <span class="n">mat1</span> <span class="o">*</span> <span class="n">mat2</span> <span class="o">*</span> <span class="n">mat3</span> <span class="o">*</span> <span class="n">matprec</span> <span class="nd">@transformations.dynamic_transform_matrix</span><span class="p">(</span><span class="n">GalacticCoordinates</span><span class="p">,</span> <span class="n">FK5Coordinates</span><span class="p">)</span> <span class="k">def</span> <span class="nf">_gal_to_fk5</span><span class="p">(</span><span class="n">galcoords</span><span class="p">):</span> <span class="k">return</span> <span class="n">_fk5_to_gal</span><span class="p">(</span><span class="n">galcoords</span><span class="p">)</span><span class="o">.</span><span class="n">T</span> <span class="nd">@transformations.dynamic_transform_matrix</span><span class="p">(</span><span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="n">GalacticCoordinates</span><span class="p">,</span> <span class="n">priority</span><span class="o">=</span><span class="mi">1</span><span class="p">)</span> <span class="k">def</span> <span class="nf">_fk4_to_gal</span><span class="p">(</span><span class="n">fk4coords</span><span class="p">):</span> <span class="kn">from</span> <span class="nn">.angles</span> <span class="kn">import</span> <span class="n">rotation_matrix</span> <span class="kn">from</span> <span class="nn">.earth_orientation</span> <span class="kn">import</span> <span class="n">_precession_matrix_besselian</span> <span class="c"># needed mainly to support inverse from galactic</span> <span class="n">bequinox</span> <span class="o">=</span> <span class="mi">1950</span> <span class="k">if</span> <span class="n">fk4coords</span><span class="o">.</span><span class="n">equinox</span> <span class="ow">is</span> <span class="bp">None</span> <span class="k">else</span> <span class="n">fk4coords</span><span class="o">.</span><span class="n">equinox</span><span class="o">.</span><span class="n">byear</span> <span class="n">mat1</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="mi">180</span> <span class="o">-</span> <span class="n">GalacticCoordinates</span><span class="o">.</span><span class="n">_lon0_B1950</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="s">'z'</span><span class="p">)</span> <span class="n">mat2</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="mi">90</span> <span class="o">-</span> <span class="n">GalacticCoordinates</span><span class="o">.</span><span class="n">_ngp_B1950</span><span class="o">.</span><span class="n">dec</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="s">'y'</span><span class="p">)</span> <span class="n">mat3</span> <span class="o">=</span> <span class="n">rotation_matrix</span><span class="p">(</span><span class="n">GalacticCoordinates</span><span class="o">.</span><span class="n">_ngp_B1950</span><span class="o">.</span><span class="n">ra</span><span class="o">.</span><span class="n">degrees</span><span class="p">,</span> <span class="s">'z'</span><span class="p">)</span> <span class="n">matprec</span> <span class="o">=</span> <span class="n">_precession_matrix_besselian</span><span class="p">(</span><span class="n">bequinox</span><span class="p">,</span> <span class="mi">1950</span><span class="p">)</span> <span class="k">return</span> <span class="n">mat1</span> <span class="o">*</span> <span class="n">mat2</span> <span class="o">*</span> <span class="n">mat3</span> <span class="o">*</span> <span class="n">matprec</span> <span class="nd">@transformations.dynamic_transform_matrix</span><span class="p">(</span><span class="n">GalacticCoordinates</span><span class="p">,</span> <span class="n">FK4NoETermCoordinates</span><span class="p">,</span> <span class="n">priority</span><span class="o">=</span><span class="mi">1</span><span class="p">)</span> <span class="k">def</span> <span class="nf">_gal_to_fk4</span><span class="p">(</span><span class="n">galcoords</span><span class="p">):</span> <span class="k">return</span> <span class="n">_fk4_to_gal</span><span class="p">(</span><span class="n">galcoords</span><span class="p">)</span><span class="o">.</span><span class="n">T</span> <span class="k">def</span> <span class="nf">_make_transform_graph_docs</span><span class="p">():</span> <span class="sd">"""</span> <span class="sd"> Generates a string for use with the coordinate package's docstring</span> <span class="sd"> to show the available transforms and coordinate systems</span> <span class="sd"> """</span> <span class="kn">from</span> <span class="nn">inspect</span> <span class="kn">import</span> <span class="n">isclass</span> <span class="kn">from</span> <span class="nn">textwrap</span> <span class="kn">import</span> <span class="n">dedent</span> <span class="kn">from</span> <span class="nn">.transformations</span> <span class="kn">import</span> <span class="n">master_transform_graph</span> <span class="n">coosys</span> <span class="o">=</span> <span class="p">[</span><span class="n">item</span> <span class="k">for</span> <span class="n">item</span> <span class="ow">in</span> <span class="nb">globals</span><span class="p">()</span><span class="o">.</span><span class="n">values</span><span class="p">()</span> <span class="k">if</span> <span class="n">isclass</span><span class="p">(</span><span class="n">item</span><span class="p">)</span> <span class="ow">and</span> <span class="nb">issubclass</span><span class="p">(</span><span class="n">item</span><span class="p">,</span> <span class="n">SphericalCoordinatesBase</span><span class="p">)]</span> <span class="n">coosys</span><span class="o">.</span><span class="n">remove</span><span class="p">(</span><span class="n">SphericalCoordinatesBase</span><span class="p">)</span> <span class="n">graphstr</span> <span class="o">=</span> <span class="n">master_transform_graph</span><span class="o">.</span><span class="n">to_dot_graph</span><span class="p">(</span><span class="n">addnodes</span><span class="o">=</span><span class="n">coosys</span><span class="p">)</span> <span class="n">docstr</span> <span class="o">=</span> <span class="s">"""</span> <span class="s"> The diagram below shows all of the coordinate systems built into the</span> <span class="s"> `~astropy.coordinates` package, their aliases (usable for converting</span> <span class="s"> other coordinates to them using attribute-style access) and the</span> <span class="s"> pre-defined transformations between them. The user is free to</span> <span class="s"> override any of these transformations by defining new trasnformation</span> <span class="s"> between these systems, but the pre-defined transformations should be</span> <span class="s"> sufficient for typical usage.</span> <span class="s"> The graph also indicates the priority for each transformation as a</span> <span class="s"> number next to the arrow. These priorities are used to decide the</span> <span class="s"> preferred order when two trasnformation paths have the same number</span> <span class="s"> of steps. These priorities are defined such that path with a</span> <span class="s"> *smaller* total priority are favored over larger.</span> <span class="s"> E.g., the path from `ICRSCoordinates` to `GalacticCoordinates` goes</span> <span class="s"> through `FK5Coordinates` because the total path length is 2 instead</span> <span class="s"> of 2.03.</span> <span class="s"> .. graphviz::</span> <span class="s"> """</span> <span class="k">return</span> <span class="n">dedent</span><span class="p">(</span><span class="n">docstr</span><span class="p">)</span> <span class="o">+</span> <span class="s">' '</span> <span class="o">+</span> <span class="n">graphstr</span><span class="o">.</span><span class="n">replace</span><span class="p">(</span><span class="s">'</span><span class="se">\n</span><span class="s">'</span><span class="p">,</span> <span class="s">'</span><span class="se">\n</span><span class="s"> '</span><span class="p">)</span> <span class="n">_transform_graph_docs</span> <span class="o">=</span> <span class="n">_make_transform_graph_docs</span><span class="p">()</span> </pre></div> </div> </div> </div> <div class="sphinxsidebar"> <div class="sphinxsidebarwrapper"><h3>Page Contents</h3> </div> </div> <div class="clearer"></div> </div> <footer class="footer"> <p class="pull-right"> <a href="#">Back to Top</a></p> <p> © Copyright 2011-2013, The Astropy Developers.<br/> Created using <a href="http://sphinx.pocoo.org/">Sphinx</a> 1.1.3. Last built 22 Oct 2013. <br/> </p> </footer> </body> </html>