summaryrefslogtreecommitdiffstats
path: root/ast/pal/palDe2h.c
blob: a250e9ed13f36568cefba3f3546cc1c4b05b4a1d (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
/*
*+
*  Name:
*     palDe2h

*  Purpose:
*     Equatorial to horizon coordinates: HA,Dec to Az,E

*  Language:
*     Starlink ANSI C

*  Type of Module:
*     Library routine

*  Invocation:
*     palDe2h( double ha, double dec, double phi, double * az, double * el );

*  Arguments:
*     ha = double * (Given)
*        Hour angle (radians)
*     dec = double * (Given)
*        Declination (radians)
*     phi = double (Given)
*        Observatory latitude (radians)
*     az = double * (Returned)
*        Azimuth (radians)
*     el = double * (Returned)
*        Elevation (radians)

*  Description:
*     Convert equatorial to horizon coordinates.

*  Authors:
*     PTW: Pat Wallace (STFC)
*     TIMJ: Tim Jenness (JAC, Hawaii)
*     {enter_new_authors_here}

*  Notes:
*     - All the arguments are angles in radians.
*     - Azimuth is returned in the range 0-2pi;  north is zero,
*       and east is +pi/2.  Elevation is returned in the range
*       +/-pi/2.
*     - The latitude must be geodetic.  In critical applications,
*       corrections for polar motion should be applied.
*     - In some applications it will be important to specify the
*       correct type of hour angle and declination in order to
*       produce the required type of azimuth and elevation.  In
*       particular, it may be important to distinguish between
*       elevation as affected by refraction, which would
*       require the "observed" HA,Dec, and the elevation
*       in vacuo, which would require the "topocentric" HA,Dec.
*       If the effects of diurnal aberration can be neglected, the
*       "apparent" HA,Dec may be used instead of the topocentric
*       HA,Dec.
*     - No range checking of arguments is carried out.
*     - In applications which involve many such calculations, rather
*       than calling the present routine it will be more efficient to
*       use inline code, having previously computed fixed terms such
*       as sine and cosine of latitude, and (for tracking a star)
*       sine and cosine of declination.

*  History:
*     2012-02-08 (TIMJ):
*        Initial version with documentation taken from Fortran SLA
*        Adapted with permission from the Fortran SLALIB library.
*     {enter_further_changes_here}

*  Copyright:
*     Copyright (C) 1995 Rutherford Appleton Laboratory
*     Copyright (C) 2012 Science and Technology Facilities Council.
*     All Rights Reserved.

*  Licence:
*     This program is free software: you can redistribute it and/or
*     modify it under the terms of the GNU Lesser General Public
*     License as published by the Free Software Foundation, either
*     version 3 of the License, or (at your option) any later
*     version.
*
*     This program is distributed in the hope that it will be useful,
*     but WITHOUT ANY WARRANTY; without even the implied warranty of
*     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
*     GNU Lesser General Public License for more details.
*
*     You should have received a copy of the GNU Lesser General
*     License along with this program.  If not, see
*     <http://www.gnu.org/licenses/>.

*  Bugs:
*     {note_any_bugs_here}
*-
*/

#include "pal.h"
#include "palmac.h"
#include <math.h>

void
palDe2h ( double ha, double dec, double phi, double *az, double *el) {

  double sh;
  double ch;
  double sd;
  double cd;
  double sp;
  double cp;

  double a;

  double x;
  double y;
  double z;
  double r;

  /*  Useful trig functions */
  sh = sin(ha);
  ch = cos(ha);
  sd = sin(dec);
  cd = cos(dec);
  sp = sin(phi);
  cp = cos(phi);

  /*  Az,El as x,y,z */
  x = -ch * cd * sp + sd * cp;
  y = -sh * cd;
  z = ch * cd * cp + sd * sp;

  /*  To spherical */
  r = sqrt(x * x + y * y);
  if (r == 0.) {
    a = 0.;
  } else {
    a = atan2(y, x);
  }
  if (a < 0.) {
    a += PAL__D2PI;
  }
  *az = a;
  *el = atan2(z, r);

  return;
}