libTheSky
Routines to compute sky positions of Sun, Moon, planets and more
All Namespaces Files Functions Variables Pages
coordinates.f90
Go to the documentation of this file.
1!> \file coordinates.f90 Procedures to perform coordinate transformations, apply precession, and more for libTheSky
2
3
4! Copyright (c) 2002-2025 Marc van der Sluys - Nikhef/Utrecht University - marc.vandersluys.nl
5!
6! This file is part of the libTheSky package,
7! see: https://libthesky.sf.net/
8!
9! This is free software: you can redistribute it and/or modify it under the terms of the
10! European Union Public Licence 1.2 (EUPL 1.2).
11!
12! This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13! without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14! See the EU Public Licence for more details.
15!
16! You should have received a copy of the European Union Public Licence along with this code.
17! If not, see <https://www.eupl.eu/1.2/en/>.
18
19
20!***********************************************************************************************************************************
21!> \brief Procedures for coordinates
22
24 implicit none
25 save
26
27
28contains
29
30
31 !*********************************************************************************************************************************
32 !> \brief Compute the geocentric rectangular coordinates of a planet, from its and the Earth's heliocentric spherical position
33 !!
34 !! \param l Heliocentric longitude of planet
35 !! \param b Heliocentric latitude of planet
36 !! \param r Heliocentric distance of planet
37 !!
38 !! \param l0 Heliocentric longitude of Earth
39 !! \param b0 Heliocentric latitude of Earth
40 !! \param r0 Heliocentric distance of Earth
41 !!
42 !! \param x Geocentric rectangular X of planet (output)
43 !! \param y Geocentric rectangular Y of planet (output)
44 !! \param z Geocentric rectangular Z of planet (output)
45
46 subroutine hc_spher_2_gc_rect(l,b,r, l0,b0,r0, x,y,z)
47 use sufr_kinds, only: double
48
49 implicit none
50 real(double), intent(in) :: l,b,r, l0,b0,r0
51 real(double), intent(out) :: x,y,z
52
53 x = r * cos(b) * cos(l) - r0 * cos(b0) * cos(l0)
54 y = r * cos(b) * sin(l) - r0 * cos(b0) * sin(l0)
55 z = r * sin(b) - r0 * sin(b0)
56
57 end subroutine hc_spher_2_gc_rect
58 !*********************************************************************************************************************************
59
60
61 !*********************************************************************************************************************************
62 !> \brief Convert spherical, ecliptical coordinates to rectangular, equatorial coordinates of an object - both geocentric
63 !!
64 !! \param l Heliocentric longitude of planet
65 !! \param b Heliocentric latitude of planet
66 !! \param r Heliocentric distance of planet
67 !! \param eps Obliquity of the ecliptic
68 !!
69 !! \param x Geocentric rectangular X of planet (output)
70 !! \param y Geocentric rectangular Y of planet (output)
71 !! \param z Geocentric rectangular Z of planet (output)
72
73 subroutine ecl_spher_2_eq_rect(l,b,r, eps, x,y,z)
74 use sufr_kinds, only: double
75
76 implicit none
77 real(double), intent(in) :: l,b,r,eps
78 real(double), intent(out) :: x,y,z
79
80 x = r * cos(b) * cos(l)
81 y = r * (cos(b) * sin(l) * cos(eps) - sin(b) * sin(eps))
82 z = r * (cos(b) * sin(l) * sin(eps) + sin(b) * cos(eps))
83
84 end subroutine ecl_spher_2_eq_rect
85 !*********************************************************************************************************************************
86
87
88
89 !*********************************************************************************************************************************
90 !> \brief Compute the geocentric equatorial rectangular coordinates of the Sun, from Earth's heliocentric, spherical position
91 !!
92 !! \param t1 Dynamical time in Julian Millennia after 2000.0
93 !! \param l0 Heliocentric longitude of Earth
94 !! \param b0 Heliocentric latitude of Earth
95 !! \param r0 Heliocentric distance of Earth
96 !!
97 !! \param x Geocentric rectangular X of planet (output)
98 !! \param y Geocentric rectangular Y of planet (output)
99 !! \param z Geocentric rectangular Z of planet (output)
100 !!
101 !! \todo To be disposed, has been replcaed by ecl_spher_2_eq_rect above
102
103 subroutine calcsunxyz(t1, l0,b0,r0, x,y,z)
104 use sufr_kinds, only: double
105 use sufr_constants, only: pi
106 use sufr_angles, only: rev
107
108 use thesky_nutation, only: nutation
109
110 implicit none
111 real(double), intent(in) :: t1, l0,b0,r0
112 real(double), intent(out) :: x,y,z
113 real(double) :: t, l,b,r, dpsi,eps,eps0,deps
114
115 t = t1
116 l = rev(l0+pi) ! Heliocentric longitude of the Earth -> geocentric longitude of the Sun
117 b = -b0 ! Heliocentric latitude of the Earth -> geocentric latitude of the Sun
118 r = r0 ! Heliocentric distance of the Earth = geocentric distance of the Sun
119
120 call fk5(t, l,b)
121 call nutation(t,dpsi,eps0,deps)
122 eps = eps0
123
124 x = r * cos(b) * cos(l)
125 y = r * (cos(b) * sin(l) * cos(eps) - sin(b) * sin(eps))
126 z = r * (cos(b) * sin(l) * sin(eps) + sin(b) * cos(eps))
127
128 !jd1 = t1*365250.d0 + jd2000 !From equinox of date...
129 !jd2 = comepoche ! ... to equinox of elements
130 !call precess_xyz(jd1,jd2,x,y,z)
131
132 end subroutine calcsunxyz
133 !*********************************************************************************************************************************
134
135
136 !*********************************************************************************************************************************
137 !> \brief Compute the precession of the equinoxes in rectangular coordinates, from jd1 to jd2
138 !!
139 !! \param jd1 Original Julian day
140 !! \param jd2 Target Julian day
141 !!
142 !! \param x Geocentric rectangular X (output)
143 !! \param y Geocentric rectangular Y (output)
144 !! \param z Geocentric rectangular Z (output)
145 !!
146 !! \see Meeus, Astronomical Algorithms, 1998, Ch. 21
147
148 subroutine precess_xyz(jd1,jd2, x,y,z)
149 use sufr_kinds, only: double
150 use sufr_constants, only: jd2000
151
152 implicit none
153 real(double), intent(in) :: jd1,jd2
154 real(double), intent(inout) :: x,y,z
155 real(double) :: t1,t2,t22,t23, x1,y1,z1
156 real(double) :: dz,ze,th,sd,cd,sz,cz,st,ct
157 real(double) :: xx,xy,xz,yx,yy,yz,zx,zy,zz
158
159 t1 = (jd1 - jd2000)/36525.d0 ! t since 2000.0 in Julian centuries
160 t2 = (jd2 - jd1)/36525.d0 ! dt in Julian centuries
161 t22 = t2*t2 ! t2^2
162 t23 = t22*t2 ! t2^3
163
164 ! Meeus, Eq. 21.2:
165 dz = (1.11808609d-2 + 6.770714d-6*t1 - 6.739d-10*t1*t1)*t2 + (1.463556d-6 - 1.668d-9*t1)*t22 + 8.725677d-8*t23
166 ze = (1.11808609d-2 + 6.770714d-6*t1 - 6.739d-10*t1*t1)*t2 + (5.307158d-6 + 3.20d-10*t1)*t22 + 8.825063d-8*t23
167 th = (9.71717346d-3 - 4.136915d-6*t1 - 1.0520d-9*t1*t1)*t2 - (2.068458d-6 + 1.052d-9*t1)*t22 - 2.028121d-7*t23
168
169 sd = sin(dz)
170 cd = cos(dz)
171 sz = sin(ze)
172 cz = cos(ze)
173 st = sin(th)
174 ct = cos(th)
175
176 ! Ref?
177 xx = cd*cz*ct - sd*sz
178 xy = sd*cz + cd*sz*ct
179 xz = cd*st
180 yx = -cd*sz - sd*cz*ct
181 yy = cd*cz - sd*sz*ct
182 yz = -sd*st
183 zx = -cz*st
184 zy = -sz*st
185 zz = ct
186
187 x1 = xx*x + yx*y + zx*z
188 y1 = xy*x + yy*y + zy*z
189 z1 = xz*x + yz*y + zz*z
190
191 x = x1
192 y = y1
193 z = z1
194
195 end subroutine precess_xyz
196 !*********************************************************************************************************************************
197
198
199
200 !*********************************************************************************************************************************
201 !> \brief Compute the precession of the equinoxes in equatorial coordinates, from jd1 to jd2
202 !!
203 !! \param jd1 Original Julian day
204 !! \param jd2 Target Julian day
205 !! \param a1 Right ascension (IO, rad) (output)
206 !! \param d1 Declination (IO, rad) (output)
207 !!
208 !! \see Meeus, Astronomical Algorithms, 1998, Ch.21, p.134
209
210 subroutine precess_eq(jd1,jd2, a1,d1)
211 use sufr_kinds, only: double
212 use sufr_constants, only: jd2000
213
214 implicit none
215 real(double), intent(in) :: jd1,jd2
216 real(double), intent(inout) :: a1,d1
217 real(double) :: t1,t12, t2,t22,t23, dz,ze,th,a,b,c
218
219 t1 = (jd1 - jd2000)/36525.d0 ! t since 2000.0 in Julian centuries
220 t12 = t1**2 ! t1^2
221 t2 = (jd2 - jd1)/36525.d0 ! dt in Julian centuries
222 t22 = t2**2 ! t2^2
223 t23 = t22*t2 ! t2^3
224
225 ! Meeus, Eq. 21.2:
226 dz = (1.11808609d-2 + 6.770714d-6*t1 - 6.739d-10*t12)*t2 + (1.463556d-6 - 1.668d-9*t1)*t22 + 8.725677d-8*t23
227 ze = (1.11808609d-2 + 6.770714d-6*t1 - 6.739d-10*t12)*t2 + (5.307158d-6 + 3.20d-10*t1)*t22 + 8.825063d-8*t23
228 th = (9.71717346d-3 - 4.136915d-6*t1 - 1.0520d-9*t12)*t2 - (2.068458d-6 + 1.052d-9*t1)*t22 - 2.028121d-7*t23
229
230 ! Meeus, Eq. 21.4:
231 a = cos(d1) * sin(a1+dz)
232 b = cos(th) * cos(d1) * cos(a1+dz) - sin(th) * sin(d1)
233 c = sin(th) * cos(d1) * cos(a1+dz) + cos(th) * sin(d1)
234
235 a1 = atan2(a,b) + ze
236 d1 = asin(c)
237
238 end subroutine precess_eq
239 !*********************************************************************************************************************************
240
241
242
243 !*********************************************************************************************************************************
244 !> \brief Compute the precession of the equinoxes in equatorial coordinates, from yr1 to yr2
245 !!
246 !! \param yr1 Original year
247 !! \param yr2 Target year
248 !! \param a1 Right ascension (IO, rad) (output)
249 !! \param d1 Declination (IO, rad) (output)
250 !!
251 !! \see Meeus, Astronomical Algorithms, 1998, Ch.21, p.134
252
253 subroutine precess_eq_yr(yr1,yr2,a1,d1)
254 use sufr_kinds, only: double
255
256 implicit none
257 real(double), intent(in) :: yr1,yr2
258 real(double), intent(inout) :: a1,d1
259 real(double) :: t1,t12, t2,t22,t23, dz,ze,th,a,b,c
260
261 t1 = (yr1 - 2000.d0)/100.d0 ! t since 2000.0 in Julian centuries
262 t12 = t1**2 ! t1^2
263 t2 = (yr2 - yr1)/100.d0 ! dt in Julian centuries
264 t22 = t2*t2 ! t2^2
265 t23 = t22*t2 ! t2^3
266
267 ! Meeus, Eq. 21.2:
268 dz = (1.11808609d-2 + 6.770714d-6*t1 - 6.739d-10*t12)*t2 + (1.463556d-6 - 1.668d-9*t1)*t22 + 8.725677d-8*t23
269 ze = (1.11808609d-2 + 6.770714d-6*t1 - 6.739d-10*t12)*t2 + (5.307158d-6 + 3.20d-10*t1)*t22 + 8.825063d-8*t23
270 th = (9.71717346d-3 - 4.136915d-6*t1 - 1.0520d-9*t12)*t2 - (2.068458d-6 + 1.052d-9*t1)*t22 - 2.028121d-7*t23
271
272 ! Meeus, Eq. 21.4:
273 a = cos(d1) * sin(a1+dz)
274 b = cos(th) * cos(d1) * cos(a1+dz) - sin(th) * sin(d1)
275 c = sin(th) * cos(d1) * cos(a1+dz) + cos(th) * sin(d1)
276
277 a1 = atan2(a,b) + ze
278 d1 = asin(c)
279
280 end subroutine precess_eq_yr
281 !*********************************************************************************************************************************
282
283
284
285 !*********************************************************************************************************************************
286 !> \brief Compute the precession of the equinoxes in geocentric ecliptical coordinates
287 !!
288 !! \param jd1 Original Julian day
289 !! \param jd2 Target Julian day
290 !! \param l Ecliptic longitude (I/O, rad) (output)
291 !! \param b Ecliptic latitude (I/O, rad) (output)
292 !!
293 !! \see Meeus, Astronomical Algorithms, 1998, Ch.21, p.136
294
295 subroutine precess_ecl(jd1,jd2, l,b)
296 use sufr_kinds, only: double
297 use sufr_constants, only: jd2000
298
299 implicit none
300 real(double), intent(in) :: jd1,jd2
301 real(double), intent(inout) :: l,b
302 real(double) :: t1,t12, t2,t22,t23, eta,pii,p,aa,bb,cc
303
304 t1 = (jd1 - jd2000)/36525.d0 ! t since 2000.0 in Julian centuries
305 t12 = t1**2 ! t1^2
306 t2 = (jd2 - jd1)/36525.d0 ! dt in Julian centuries
307 t22 = t2*t2 ! t2^2
308 t23 = t22*t2 ! t2^3
309
310 ! Meeus, Eq. 21.5 - the different powers of t1/t2 have been checked and are ok:
311 eta = (2.278765d-4 - 3.2012d-7 *t1 + 2.899d-9*t12) *t2 + (-1.60085d-7 + 2.899d-9*t1) *t22 + 2.909d-10*t23
312 pii = 3.0521686858d0 + 1.59478437d-2 *t1 + 2.939037d-6*t12 - (4.2169525d-3 + 2.44787d-6*t1) *t2 + 1.7143d-7*t22
313 p = (2.438174835d-2 + 1.077382d-5 *t1 - 2.036d-10*t12) *t2 + (5.38691d-6 - 2.036d-10*t1) *t22 - 2.91d-11*t23
314
315 ! Meeus, Eq. 21.7:
316 aa = cos(eta) * cos(b) * sin(pii-l) - sin(eta) * sin(b)
317 bb = cos(b) * cos(pii-l)
318 cc = cos(eta) * sin(b) + sin(eta) * cos(b) * sin(pii-l)
319
320 l = p + pii - atan2(aa,bb)
321 b = asin(cc)
322
323 end subroutine precess_ecl
324 !*********************************************************************************************************************************
325
326
327
328 !*********************************************************************************************************************************
329 !> \brief Compute the precession of the equinoxes in orbital elements
330 !!
331 !! \param jd1 Original Julian day
332 !! \param jd2 Target Julian day
333 !!
334 !! \param i Inclination (output)
335 !! \param o1 Argument of perihelion (output)
336 !! \param o2 Longitude of ascending node (output)
337 !!
338 !! \see Meeus, Astronomical Algorithms, 1998, Ch.24
339
340
341 subroutine precess_orb(jd1,jd2,i,o1,o2)
342 use sufr_kinds, only: double
343 use sufr_constants, only: jd2000
344
345 implicit none
346 real(double), intent(in) :: jd1,jd2
347 real(double), intent(inout) :: i,o1,o2
348 real(double) :: t1,t12,t2,t22,t23, eta,pii,p,psi,aa,bb,cc,dd
349
350 t1 = (jd1 - jd2000)/36525.d0 ! t since 2000.0 in Julian centuries
351 t12 = t1**2 ! t1^2
352 t2 = (jd2 - jd1)/36525.d0 ! dt in Julian centuries
353 t22 = t2*t2 ! t2^2
354 t23 = t22*t2 ! t2^3
355
356 ! Meeus, Eq. 21.5:
357 eta = (2.278765d-4 - 3.2012d-7 * t1 + 2.899d-9*t12) *t2 + (-1.60085d-7 + 2.899d-9*t1) * t22 + 2.909d-10 * t23
358 pii = 3.0521686858d0 + 1.59478437d-2 * t1 + 2.939037d-6*t12 - (4.2169525d-3 + 2.44787d-6*t1) * t2 + 1.7143d-7 * t22
359 p = (2.438174835d-2 + 1.077382d-5 * t1 - 2.036d-10*t12) *t2 + (5.38691d-6 - 2.036d-10*t1) * t22 - 2.91d-11 * t23
360 psi = pii + p
361
362 ! Meeus, Eq. 24.2:
363 aa = sin(i) * sin(o2-pii)
364 bb = -sin(eta) * cos(i) + cos(eta) * sin(i) * cos(o2-pii)
365 cc = -sin(eta) * sin(o2-pii)
366 dd = sin(i) * cos(eta) - cos(i) * sin(eta) * cos(o2-pii)
367
368 i = asin(sqrt(aa*aa+bb*bb)) ! Inclination
369 o2 = atan2(aa,bb) + psi ! Longitude of ascending node (Omega)
370 o1 = o1 + atan2(cc,dd) ! Argument of perihelion (omega)
371
372 end subroutine precess_orb
373 !*********************************************************************************************************************************
374
375
376
377 !*********************************************************************************************************************************
378 !> \brief Convert rectangular coordinates x,y,z to spherical coordinates l,b,r
379 !!
380 !! \param x Rectangular x coordinate (same unit as r)
381 !! \param y Rectangular y coordinate (same unit as r)
382 !! \param z Rectangular z coordinate (same unit as r)
383 !!
384 !! \param l Longitude (rad) (output)
385 !! \param b Latitude (rad) (output)
386 !! \param r Distance (same unit as x,y,z) (output)
387
388 subroutine rect_2_spher(x,y,z, l,b,r)
389 use sufr_kinds, only: double
390 use sufr_angles, only: rev
391 use sufr_numerics, only: deq0
392
393 implicit none
394 real(double), intent(in) :: x,y,z
395 real(double), intent(out) :: l,b,r
396 real(double) :: x2,y2
397
398 if(deq0(x) .and. deq0(y) .and. deq0(z)) then
399 l = 0.d0
400 b = 0.d0
401 r = 0.d0
402 else
403 x2 = x**2
404 y2 = y**2
405
406 l = rev( atan2(y,x) ) ! Longitude
407 b = atan2(z, sqrt(x2 + y2)) ! Latitude
408 r = sqrt(x2 + y2 + z**2) ! Distance
409 end if
410
411 end subroutine rect_2_spher
412 !*********************************************************************************************************************************
413
414
415 !*********************************************************************************************************************************
416 !> \brief Correct ecliptic longitude and latitiude for annual aberration
417 !!
418 !! \param t Dynamical time in Julian Millennia since 2000.0
419 !! \param l0 Earth longitude (rad)
420 !!
421 !! \param l Longitude of the object (rad, I/O) (output)
422 !! \param b Latitude of the object (rad, I/O) (output)
423 !!
424 !! \see Meeus, Astronomical Algorithms, 1998, Ch.23
425
426 subroutine aberration_ecl(t,l0, l,b)
427 use sufr_kinds, only: double
428 use sufr_constants, only: pi
429 use sufr_angles, only: rev
430
431 implicit none
432 real(double), intent(in) :: t,l0
433 real(double), intent(inout) :: l,b
434 real(double) :: tt,tt2, k,odot,e,pii,dl,db
435
436 tt = t*10 ! Time in Julian Centuries
437 tt2 = tt*tt
438
439 k = 9.93650849745d-5 ! Constant of aberration, radians
440 odot = rev(l0+pi) ! Longitude of the Sun
441
442 e = 0.016708634d0 - 4.2037d-5 * tt - 1.267d-7 * tt2 ! Eccentricity of the Earth's orbit
443 pii = 1.79659568d0 + 3.001024d-2 * tt + 8.0285d-6 * tt2 ! Longitude of perihelion of "
444
445 dl = k * ( -cos(odot-l) + e * cos(pii-l)) / cos(b) ! Meeus, Eq. 23.2
446 db = -k * sin(b) * (sin(odot-l) - e * sin(pii-l))
447
448 l = l + dl
449 b = b + db
450
451 end subroutine aberration_ecl
452 !*********************************************************************************************************************************
453
454 !*********************************************************************************************************************************
455 !> \brief Correct equatorial coordinates for annual aberration - moderate accuracy, use for stars
456 !!
457 !! \param jd Julian day of epoch
458 !!
459 !! \param ra Right ascension of the object (rad)
460 !! \param dec Declination of the object (rad)
461 !!
462 !! \param dra Correction in right ascension due to aberration (rad) (output)
463 !! \param ddec Correction in declination due to aberration(rad) (output)
464 !!
465 !! \param eps0 The mean obliquity of the ecliptic
466 !!
467 !!
468 !! \see Meeus, Astronomical Algorithms, 1998, Ch.23
469
470 subroutine aberration_eq(jd, ra,dec, dra, ddec, eps0)
471 use sufr_kinds, only: double
472 use sufr_constants, only: jd2000
473 use sufr_angles, only: rev
474 use sufr_dummy, only: dumdbl1,dumdbl2
475 use sufr_numerics, only: dne
476 use thesky_nutation, only: nutation
477
478 implicit none
479 real(double), intent(in) :: jd, ra,dec
480 real(double), intent(out) :: dra,ddec
481 real(double), intent(in), optional :: eps0
482 real(double) :: jdold, tjc,tjc2, l0,mas,sec,odot, k,ee,pii, eps
483 save jdold, tjc,odot,k,ee,pii, eps
484
485
486 ! Compute these constants (independent of ra,dec) only if jd has changed since the last call
487 if(dne(jd,jdold)) then
488
489 ! Meeus, Ch. 25:
490 tjc = (jd-jd2000)/36525.0d0 ! Julian centuries since 2000.0, Meeus Eq. 25.1
491 tjc2 = tjc**2
492 l0 = 4.8950631684d0 + 628.331966786d0 * tjc + 5.291838d-6 *tjc2 ! Mean longitude of the Sun, Meeus Eq. 25.2
493
494 ! Mean anomaly of the Sun, Meeus Eq. 25.3 -> p.144 (nutation):
495 mas = 6.240035881d0 + 628.301956024d0 * tjc - 2.79776d-6 * tjc2 - 5.817764d-8 * tjc**3
496
497 ! Sun's equation of the centre, Meeus p.164:
498 sec = (3.3416109d-2 - 8.40725d-5*tjc - 2.443d-7*tjc2)*sin(mas) + (3.489437d-4-1.763d-6*tjc)*sin(2*mas) + 5.044d-6*sin(3*mas)
499 odot = rev(l0 + sec) ! Sun's true longitude
500
501 ! Meeus, p.151:
502 k = 9.93650849745d-5 ! Constant of aberration, in radians (kappa in Meeus)
503 ee = 0.016708634d0 - 4.2037d-5 * tjc - 1.267d-7 * tjc2 ! Eccentricity of the Earth's orbit
504 pii = 1.79659568d0 + 3.001024d-2 * tjc + 8.0285d-6 * tjc2 ! Longitude of perihelion of "
505
506
507 ! If eps0 is provided, use it - otherwise compute it:
508 if(present(eps0)) then
509 eps = eps0
510 else
511 call nutation(tjc/10.d0, dumdbl1, eps, dumdbl2) ! Note: eps is actually eps0
512 end if
513
514 end if
515 jdold = jd
516
517
518 ! Note that eps in Meeus is actually eps0, the *mean* obliquity of the ecliptic
519 ! Meeus, Eq. 23.3:
520 dra = k * ( -(cos(ra) * cos(odot) * cos(eps) + sin(ra) * sin(odot)) / cos(dec) &
521 + ee * (cos(ra) * cos(pii) * cos(eps) + sin(ra) * sin(pii)) / cos(dec) )
522
523 ddec = k * ( -(cos(odot) * cos(eps) * (tan(eps) * cos(dec) - sin(ra) * sin(dec)) + cos(ra) * sin(dec) * sin(odot)) &
524 + ee * (cos(pii) * cos(eps) * (tan(eps) * cos(dec) - sin(ra) * sin(dec)) + cos(ra) * sin(dec) * sin(pii)) )
525
526 end subroutine aberration_eq
527 !*********************************************************************************************************************************
528
529
530
531 !*********************************************************************************************************************************
532 !> \brief Convert coordinates to the FK5 system
533 !!
534 !! \param t Dynamical time in Julian Millennia since 2000.0
535 !! \param l Longitude (rad, I/O) (output)
536 !! \param b Latitude (rad, I/O) (output)
537
538 subroutine fk5(t, l,b)
539 use sufr_kinds, only: double
540
541 implicit none
542 real(double), intent(in) :: t
543 real(double), intent(inout) :: l,b
544 real(double) :: tt,l2,dl,db
545
546 tt = t*10 ! Julian Centuries
547
548 l2 = l - 0.02438225d0*tt - 5.41052d-6*tt*tt
549 dl = -4.379322d-7 + 1.89853d-7*(cos(l2) + sin(l2))*tan(b)
550 db = 1.89853d-7*(cos(l2) - sin(l2))
551
552 l = l + dl
553 b = b + db
554
555 end subroutine fk5
556 !*********************************************************************************************************************************
557
558
559 !*********************************************************************************************************************************
560 !> \brief Convert (geocentric) spherical ecliptical coordinates l,b (and eps) to spherical equatorial coordinates RA, Dec
561 !!
562 !! \param l Longitude (rad)
563 !! \param b Latitude (rad)
564 !! \param eps Obliquity of the ecliptic
565 !!
566 !! \param ra Right ascension (rad) (output)
567 !! \param dec Declination (rad) (output)
568 !!
569 !! \see Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.43
570
571 subroutine ecl_2_eq(l,b,eps, ra,dec)
572 use sufr_kinds, only: double
573 use sufr_angles, only: rev
574
575 implicit none
576 real(double), intent(in) :: l,b,eps
577 real(double), intent(out) :: ra,dec
578
579 ra = rev(atan2( sin(l) * cos(eps) - tan(b) * sin(eps), cos(l) ))
580 dec = asin( sin(b) * cos(eps) + cos(b) * sin(eps) * sin(l) )
581
582 end subroutine ecl_2_eq
583 !*********************************************************************************************************************************
584
585
586 !*********************************************************************************************************************************
587 !> \brief Convert (geocentric) spherical equatorial coordinates RA, Dec (and eps) to spherical ecliptical coordinates l,b
588 !!
589 !! \param ra Right ascension (rad)
590 !! \param dec Declination (rad)
591 !! \param eps Obliquity of the ecliptic
592 !!
593 !! \param l Longitude (rad) (output)
594 !! \param b Latitude (rad) (output)
595 !!
596 !! \see Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.42
597
598 subroutine eq_2_ecl(ra,dec,eps, l,b)
599 use sufr_kinds, only: double
600 use sufr_angles, only: rev
601
602 implicit none
603 real(double), intent(in) :: ra,dec,eps
604 real(double), intent(out) :: l,b
605
606 l = rev(atan2( sin(ra) * cos(eps) + tan(dec) * sin(eps), cos(ra) ))
607 b = asin( sin(dec) * cos(eps) - cos(dec) * sin(eps) * sin(ra) )
608
609 end subroutine eq_2_ecl
610 !*********************************************************************************************************************************
611
612
613 !*********************************************************************************************************************************
614 !> \brief Convert spherical equatorial coordinates (RA, dec, agst) to spherical horizontal coordinates (hh, az, alt)
615 !!
616 !! \param ra Right ascension
617 !! \param dec Declination
618 !! \param agst Greenwich sidereal time
619 !!
620 !! \param hh Local hour angle (output)
621 !! \param az Azimuth (output)
622 !! \param alt Altitude (output)
623 !!
624 !! \param lat Geographical latitude (rad), optional
625 !! \param lon Geographical longitude (rad), optional
626 !!
627 !! \see Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.44, 14.45
628
629 subroutine eq2horiz(ra,dec,agst, hh,az,alt, lat,lon)
630 use sufr_kinds, only: double
631 use sufr_angles, only: rev,rev2
632 use thesky_local, only: lat0,lon0
633
634 implicit none
635 real(double), intent(in) :: ra,dec,agst
636 real(double), intent(in), optional :: lat,lon
637 real(double), intent(out) :: hh,az,alt
638 real(double) :: llat,llon, sinhh,coshh, sinlat,coslat, sindec,cosdec,tandec
639
640 ! Lat, lon added later as optional variables:
641 llat = lat0
642 llon = lon0
643 if(present(lat)) llat = lat
644 if(present(lon)) llon = lon
645
646 hh = rev( agst + llon - ra ) ! Local Hour Angle (agst since ra is also corrected for nutation?)
647
648 ! Some preparation, saves ~29%:
649 sinhh = sin(hh)
650 coshh = cos(hh)
651 sinlat = sin(llat)
652 coslat = sqrt(1.d0-sinlat**2) ! Cosine of a latitude is always positive
653 sindec = sin(dec)
654 cosdec = sqrt(1.d0-sindec**2) ! Cosine of a declination is always positive
655 tandec = sindec/cosdec
656
657 az = rev( atan2( sinhh, coshh * sinlat - tandec * coslat )) ! Azimuth
658 alt = rev2( asin( sinlat * sindec + coslat * cosdec * coshh )) ! Altitude
659
660 end subroutine eq2horiz
661 !*********************************************************************************************************************************
662
663
664 !*********************************************************************************************************************************
665 !> \brief Convert spherical horizontal coordinates (az, alt, agst) to spherical equatorial coordinates (hh, RA, dec)
666 !!
667 !! \param az Azimuth
668 !! \param alt Altitude
669 !! \param agst Greenwich sidereal time
670 !!
671 !! \param hh Local hour angle (output)
672 !! \param ra Right ascension (output)
673 !! \param dec Declination (output)
674 !!
675 !! \param lat Geographical latitude (rad), optional
676 !! \param lon Geographical longitude (rad), optional
677 !!
678 !! \see Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.46, 14.47
679
680 subroutine horiz2eq(az,alt, agst, hh,ra,dec, lat,lon)
681 use sufr_kinds, only: double
682 use sufr_angles, only: rev,rev2
683 use thesky_local, only: lat0,lon0
684
685 implicit none
686 real(double), intent(in) :: az,alt,agst
687 real(double), intent(out) :: ra,dec,hh
688 real(double), intent(in), optional :: lat,lon
689 real(double) :: sinlat,coslat, cosaz,sinalt,cosalt,tanalt, llat,llon
690
691 ! Lat, lon added later as optional variables:
692 llat = lat0
693 llon = lon0
694 if(present(lat)) llat = lat
695 if(present(lon)) llon = lon
696
697 ! Some preparation to save time:
698 sinlat = sin(llat)
699 coslat = sqrt(1.d0-sinlat**2) ! Cosine of a latitude is always positive
700 cosaz = cos(az)
701 sinalt = sin(alt)
702 cosalt = sqrt(1.d0-sinalt**2) ! Cosine of a altitude is always positive
703 tanalt = sinalt/cosalt
704
705 hh = rev( atan2( sin(az), cosaz * sinlat + tanalt * coslat )) ! Local Hour Angle
706 dec = rev2( asin( sinlat * sinalt - coslat * cosalt * cosaz )) ! Declination
707 ra = rev( agst + llon - hh ) ! Right ascension
708
709 end subroutine horiz2eq
710 !*********************************************************************************************************************************
711
712
713 !*********************************************************************************************************************************
714 !> \brief Convert spherical equatorial coordinates (RA, dec) to spherical galactic coordinates (l,b), for J2000.0!!!
715 !!
716 !! \param ra Right ascension
717 !! \param dec Declination
718 !!
719 !! \param l Longitude (output)
720 !! \param b Latitude (output)
721 !!
722 !! \see Meeus, Astronomical Algorithms, 1998, Ch.13
723
724 subroutine eq2gal(ra,dec, l,b)
725 use sufr_kinds, only: double
726 use sufr_angles, only: rev
727
728 implicit none
729 real(double), intent(in) :: ra,dec
730 real(double), intent(out) :: l,b
731 real(double) :: ra0,dec0,l0
732
733 ra0 = 3.36603472d0 ! RA of GNP, in rad, J2000.0 (12h51m26.3s ~ 192.8596deg)
734 dec0 = 0.473478737d0 ! Decl. of gal. NP in rad, J2000.0 (27d07'42"=27.1283deg)
735 l0 = 5.287194d0 ! J2000.0?, = 302.9339deg (the 303deg in Meeus, p.94) = l0+180deg
736
737 l = rev( l0 - atan2( sin(ra0-ra), cos(ra0-ra) * sin(dec0) - tan(dec) * cos(dec0) ))
738 b = asin( sin(dec) * sin(dec0) + cos(dec) * cos(dec0) * cos(ra0-ra) )
739
740 end subroutine eq2gal
741 !*********************************************************************************************************************************
742
743
744 !*********************************************************************************************************************************
745 !> \brief Convert spherical galactic coordinates (l,b) to spherical equatorial coordinates (RA, dec), for J2000.0!!!
746 !!
747 !! \param l Longitude
748 !! \param b Latitude
749 !!
750 !! \param ra Right ascension (output)
751 !! \param dec Declination (output)
752
753 subroutine gal2eq(l,b, ra,dec)
754 use sufr_kinds, only: double
755 use sufr_angles, only: rev
756
757 real(double), intent(in) :: l,b
758 real(double), intent(out) :: ra,dec
759 real(double) :: ra0,dec0,l0
760
761 ra0 = 0.22444207d0 ! RA of GNP - 12h, in rad, J2000.0 (12h51m26.3s - 12h ~ 12.8596deg)
762 dec0 = 0.473478737d0 ! Decl. of gal. NP in rad, J2000.0 (from Sterrengids?: 27d07'42"=27.1283deg)
763 l0 = 2.145601d0 ! J2000.0?, = 123.9339deg (the 123deg in Meeus, p.94)
764
765 ra = rev(atan2( sin(l-l0), cos(l-l0) * sin(dec0) - tan(b) * cos(dec0)) + ra0)
766 dec = asin( sin(b) * sin(dec0) + cos(b) * cos(dec0) * cos(l-l0))
767
768 end subroutine gal2eq
769 !*********************************************************************************************************************************
770
771
772
773 !*********************************************************************************************************************************
774 !> \brief Convert spherical ecliptical coordinates from the geocentric to the topocentric system
775 !!
776 !! \param gcl Geocentric ecliptic longitude of the object (rad)
777 !! \param gcb Geocentric ecliptic latitude of the object (rad)
778 !! \param gcr Geocentric distance of the object
779 !! \param gcs Geocentric semi-diameter of the object (rad)
780 !!
781 !! \param eps Obliquity of the ecliptic (rad)
782 !! \param lst Local sidereal time (rad)
783 !!
784 !! \param tcl Topocentric ecliptic longitude of the object (rad) (output)
785 !! \param tcb Topocentric ecliptic latitude of the object (rad) (output)
786 !! \param tcs Topocentric semi-diameter of the object (rad) (output)
787 !!
788 !!
789 !! \param lat Latitude of the observer (rad, optional)
790 !! \param hgt Altitude/elevation of the observer above sea level (metres, optional)
791 !!
792 !! \note lat/lat0 and hgt/height can be provided through the module TheSky_local (rad, and m), or through the optional arguments.
793 !! Note that using the latter will update the former!
794 !!
795 !! \see Meeus, Astronomical Algorithms, 1998, Ch. 11 and 40
796
797 subroutine geoc2topoc_ecl(gcl,gcb, gcr,gcs, eps,lst, tcl,tcb, tcs, lat,hgt)
798 use sufr_kinds, only: double
799 use sufr_constants, only: earthr,au
800 use sufr_angles, only: rev
801 use thesky_local, only: lat0, height
802
803 implicit none
804 real(double), intent(in) :: gcl,gcb,gcr,gcs,eps,lst
805 real(double), intent(in), optional :: lat,hgt
806 real(double), intent(out) :: tcl,tcb,tcs
807 real(double) :: llat,lhgt, ba,re,u,rs,rc,sinHp,n
808
809 ! Handle optional variables:
810 llat = lat0
811 lhgt = height
812 if(present(lat)) llat = lat
813 if(present(hgt)) lhgt = hgt
814
815 ! Meeus, Ch.11, p.82:
816 ba = 0.996647189335d0 ! b/a = 1-f: flattening of the Earth - WGS84 ellipsoid
817 re = earthr*1.d-2 ! Equatorial radius of the Earth in metres
818 ! (http://earth-info.nga.mil/GandG/publications/tr8350.2/wgs84fin.pdf)
819
820 u = atan(ba*tan(llat))
821 rs = ba*sin(u) + lhgt/re*sin(llat)
822 rc = cos(u) + lhgt/re*cos(llat)
823
824
825 sinhp = sin(earthr/au)/gcr ! Sine of the horizontal parallax, Meeus, Eq. 40.1
826
827 ! Meeus, Ch.40, p.282:
828 n = cos(gcl)*cos(gcb) - rc*sinhp*cos(lst)
829
830 tcl = rev( atan2( sin(gcl)*cos(gcb) - sinhp*(rs*sin(eps) + rc*cos(eps)*sin(lst)) , n ) ) ! Topocentric longitude
831 tcb = atan((cos(tcl)*(sin(gcb) - sinhp*(rs*cos(eps) - rc*sin(eps)*sin(lst))))/n) ! Topocentric latitude
832 tcs = asin(cos(tcl)*cos(tcb)*sin(gcs)/n) ! Topocentric semi-diameter
833
834 end subroutine geoc2topoc_ecl
835 !*********************************************************************************************************************************
836
837
838 !*********************************************************************************************************************************
839 !> \brief Convert geocentric equatorial coordinates to topocentric
840 !!
841 !! \param gcra Geocentric right ascension
842 !! \param gcd Geocentric declination
843 !! \param gcr Geocentric distance
844 !! \param gch Geocentric hour angle
845 !!
846 !! \param tcra Topocentric right ascension (output)
847 !! \param tcd Topocentric declination (output)
848 !!
849 !! \param lat Latitude of the observer (rad, optional)
850 !! \param hgt Altitude/elevation of the observer above sea level (metres, optional)
851 !!
852 !! \note lat/lat0 and hgt/height can be provided through the module TheSky_local (rad, and m), or through the optional arguments.
853 !! Note that using the latter will update the former!
854 !!
855 !!
856 !! \see Meeus, Astronomical Algorithms, 1998, Ch. 11 and 40
857
858 subroutine geoc2topoc_eq(gcra,gcd,gcr,gch, tcra,tcd, lat,hgt)
859 use sufr_kinds, only: double
860 use sufr_constants, only: earthr,au
861 use thesky_local, only: lat0, height
862
863 implicit none
864 real(double), intent(in) :: gcra,gcd,gcr,gch
865 real(double), intent(in), optional :: lat,hgt
866 real(double), intent(out) :: tcra,tcd
867 real(double) :: llat,lhgt, ba,re,u,rs,rc,sinHp,dra
868
869 ! Handle optional variables:
870 llat = lat0
871 lhgt = height
872 if(present(lat)) llat = lat
873 if(present(hgt)) lhgt = hgt
874
875 ! Meeus, Ch.11, p.82:
876 ba = 0.996647189335d0 ! b/a = 1-f: flattening of the Earth - WGS84 ellipsoid
877 re = earthr*1.d-2 ! Equatorial radius of the Earth in metres
878 ! (http://earth-info.nga.mil/GandG/publications/tr8350.2/wgs84fin.pdf)
879
880 u = atan(ba*tan(llat))
881 rs = ba*sin(u) + lhgt/re*sin(llat)
882 rc = cos(u) + lhgt/re*cos(llat)
883
884 ! Meeus Ch.40:
885 sinhp = sin(earthr/au)/gcr ! Sine of the horizontal parallax, Meeus, Eq. 40.1
886
887 dra = atan2( -rc*sinhp*sin(gch) , cos(gcd)-rc*sinhp*cos(gch) ) ! Meeus, Eq. 40.2
888 tcra = gcra + dra ! Topocentric right ascension
889 tcd = atan2( (sin(gcd)-rs*sinhp)*cos(dra) , cos(gcd)-rc*sinhp*cos(gch) ) ! Topocentric declination - Meeus, Eq. 40.3
890
891 end subroutine geoc2topoc_eq
892 !*********************************************************************************************************************************
893
894
895 !*********************************************************************************************************************************
896 !> \brief Compute the atmospheric refraction for a given true altitude. You should add the result to the uncorrected altitude
897 !! in order to obtain the observed altitude. Return 0 if alt + refract < -0.3°.
898 !!
899 !! \param alt True (computed) altitude (rad)
900 !! \param press Air pressure (hPa; optional)
901 !! \param temp Air temperature (degrees Celsius; optional)
902 !!
903 !! \retval refract Refraction in altitude (rad). You should *add* the result to the uncorrected altitude.
904 !!
905 !! \see
906 !! - T. Saemundsson, "Atmospheric refraction", Sky & Telescope vol.72, p.70 (1986), converted to radians;
907 !! - Meeus (1998), Eq. 16.4 ff.
908 !!
909
910 function refract(alt, press,temp)
911 use sufr_kinds, only: double
912 use sufr_constants, only: pio2, d2r
913
914 implicit none
915 real(double), intent(in) :: alt
916 real(double), intent(in), optional :: press,temp
917 real(double) :: refract
918
919 if(abs(alt).gt.pio2) then ! for |alt| > 90° refraction is meaningless
920 refract = 0.d0
921 else
922 if(alt.lt.-1.102d0*d2r) then ! Maximum possible refraction on Earth, for T=-90°C and P=1085 mbar
923 ! Without this, refract reaches a maximum around alt=-1.9° and becomes similar to the value for |alt| for lower alt (is this bad?).
924 ! See https://en.wikipedia.org/wiki/Lowest_temperature_recorded_on_Earth, https://en.wikipedia.org/wiki/Atmospheric_pressure
925 refract = 0.663455d0*d2r ! Refraction for alt=-1.102°, but P=1010 mbar and T=10°C.
926 else
927 refract = 2.9670597d-4/tan(alt + 3.137559d-3/(alt + 8.91863d-2)) ! Overstated accuracy in translation from degrees
928 end if
929
930 if(present(press)) refract = refract * press/1010.d0 ! Correct for pressure
931 if(present(temp)) refract = refract * 283.d0/(273.d0 + temp) ! Correct for temperature
932 end if
933
934 end function refract
935 !*********************************************************************************************************************************
936
937
938
939 !*********************************************************************************************************************************
940 !> \brief Compute the atmospheric refraction for a given true altitude, using single-precision values. This
941 !! is a wrapper for refract().
942 !!
943 !! \param alt True (computed) altitude (rad)
944 !! \param press Air pressure (hPa; optional)
945 !! \param temp Air temperature (degrees Celsius; optional)
946 !!
947 !! \retval refract Refraction in altitude (rad). You should add the result to the uncorrected altitude.
948 !!
949 !! \see refract() for details.
950
951 function refract_sp(alt, press,temp)
952 use sufr_kinds, only: double
953 use sufr_constants, only: rpio2
954
955 implicit none
956 real, intent(in) :: alt
957 real, intent(in), optional :: press,temp
958 real :: refract_sp
959 real(double) :: lalt, lpress,ltemp
960
961 if(abs(alt).ge.rpio2) then ! |alt| >= 90° refraction is meaningless
962 refract_sp = 0.0
963 else
964 lalt = dble(alt)
965 lpress = 1010.d0 ! Default value for refract()
966 if(present(press)) lpress = dble(press)
967 ltemp = 10.d0 ! Default value for refract()
968 if(present(temp)) ltemp = dble(temp)
969
970 refract_sp = real(refract(lalt, lpress,ltemp))
971 end if
972
973 end function refract_sp
974 !*********************************************************************************************************************************
975
976
977
978 !*********************************************************************************************************************************
979 !> \brief Compute the atmospheric refraction of light for a given *true* altitude. Return 0 for alt<-0.9°.
980 !! This is a wrapper for aref(), which does the opposite (compute refraction for an *apparent* zenith angle).
981 !! This is an expensive way to go about(!)
982 !!
983 !! \param alt0 The true (theoretical, computed) altitude of the object in radians
984 !!
985 !! \param h0 The height of the observer above sea level in metres
986 !! \param lat0 The latitude of the observer in radians
987 !!
988 !! \param t0 The temperature at the observer in degrees Celsius (e.g. 10)
989 !! \param p0 The pressure at the observer in millibars (e.g. 1010)
990 !! \param rh The relative humidity at the observer in percent (e.g. 50)
991 !!
992 !! \param lam The wavelength of the light at the observer in nanometres (e.g. 550)
993 !! \param dTdh The temperature lapse rate dT/dh in Kelvin/metre in the troposphere (only the absolute value is used; e.g. 0.0065)
994 !!
995 !! \param eps The desired precision in *arcseconds* (e.g. 1.d-3)
996 !!
997 !! \retval atmospheric_refraction The refraction at the observer in radians
998 !!
999 !!
1000 !! \todo Adapt aref() to compute the integral in the other direction for a direct method(?)
1001 !!
1002
1003 function atmospheric_refraction(alt0, h0,lat0, t0,p0,rh, lam,dTdh, eps)
1004 use sufr_kinds, only: double
1005 use sufr_constants, only: d2r,r2d, pio2
1006
1007 implicit none
1008 real(double), intent(in) :: alt0, h0,lat0, t0,p0, rh,lam, dtdh, eps
1009 real(double) :: atmospheric_refraction, z0,zi,zio,dz, ph,lamm,t0k,rhfr
1010
1012 if(abs(alt0).ge.pio2) return ! No refraction if |alt| >= 90°
1013 if(alt0.lt.-0.9d0*d2r) return ! No refraction if alt < -0.9° - allowing this can cause aref() to fail
1014
1015 ! Convert some variables/units:
1016 z0 = 90.d0 - alt0*r2d ! Altitude (rad) -> zenith angle (in *degrees!*)
1017 ph = lat0*r2d ! Latitude: rad -> *degrees!*
1018 lamm = lam / 1000.d0 ! Wavelength: nanometre -> micrometre
1019 t0k = t0 + 273.15d0 ! Temperature: degC -> K
1020 rhfr = rh/100.d0 ! Relative humidity: % -> fraction
1021
1022 ! Set initial values for the iteration:
1023 zi = z0
1024 zio = huge(zio)
1025 dz = huge(dz)
1026
1027 do while(dz.gt.eps/3600.d0) ! Iterate until the desired accuracy is achieved
1028 atmospheric_refraction = aref(zi, h0,ph, t0k,p0,rhfr, lamm,dtdh, eps) ! Do the *inverse* calculation: refraction for given *apparent* alt (in *degrees!*)
1029 zi = z0 - atmospheric_refraction ! Compute the new zenith angle
1030 dz = abs(zi-zio) ! Compute the absolute difference with the previous iteration
1031 zio = zi ! Store the current value as the 'old value' for the next iteration
1032 end do
1033
1034 atmospheric_refraction = atmospheric_refraction * d2r ! Convert result back to radians
1035
1036 end function atmospheric_refraction
1037 !*********************************************************************************************************************************
1038
1039
1040 !*********************************************************************************************************************************
1041 !> \brief Compute the atmospheric refraction of light for a given *apparent* (observed) altitude of a celestial object
1042 !! and a given observer. Note that you will usually want to use the *true* altitude as input instead.
1043 !! The method is based on N.A.O Technical Notes 59 and 63 and a paper by Standish and Auer 'Astronomical Refraction:
1044 !! Computational Method for all Zenith Angles'. Return 0 if alt<-1.102°.
1045 !!
1046 !! \param alt0 The observed (apparent) altitude of the object (rad).
1047 !!
1048 !! \param h0 The height of the observer above sea level (metres).
1049 !! \param ph The latitude of the observer (rad).
1050 !!
1051 !! \param t0 The temperature at the observer's site; optional, default 10 (°C).
1052 !! \param p0 The air pressure at the observer's site; optional, default: 1010 (mbar).
1053 !! \param rh The relative humidity at the observer's site; optional, default: 0.5 (fraction: 0-1).
1054 !!
1055 !! \param lam The wavelength of the light; optional, default: 550 (nm).
1056 !! \param dTdh The temperature gradient dT/dh in the troposphere; the absolute value is used; optional, default: 6.5e-3 (Kelvin/metre).
1057 !!
1058 !! \param eps The desired precision; optional, default: 1e-4; don't make this smaller than ~1e-8 (rad).
1059 !!
1060 !! \retval atmospheric_refraction_apparent The refraction at the observer (rad).
1061 !!
1062 !! \see Hohenkerk & Sinclair, HMNAO technical note 63 (1985).
1063
1064 function atmospheric_refraction_apparent(alt0, h0,ph, t0,p0,rh, lam,dTdh,eps)
1065 use sufr_kinds, only: double
1066 use sufr_constants, only: pio2, d2r,r2d, earthr
1067
1068 implicit none
1069 real(double), intent(in) :: alt0, h0,ph
1070 real(double), intent(in), optional :: t0,p0, rh,lam, dtdh, eps
1071
1072 real(double), parameter :: gcr=8314.36d0, md=28.966d0, mw=18.016d0, gamma=18.36d0
1073 real(double), parameter :: ht=11.d3, hs=80.d3, z2=11.2684d-06
1074
1075 integer :: i,in,is,istart, j,k
1076 real(double) :: atmospheric_refraction_apparent, z0, t0l,p0l,rhl, laml, dtdhl, epsl
1077 real(double) :: n,n0,nt,nts,ns,a(10), dndr,dndr0,dndrs,dndrt,dndrts, f,f0, fb,fe,ff,fo,fs,ft,fts,gb, h
1078 real(double) :: pw0,r,r0,refo,refp,reft,rg,rs,rt, sk0,step,t,t0o,tg,tt,z,z1,zs,zt,zts, earthrm
1079
1081 if(alt0.lt.-1.102d0*d2r) return ! Cannot *observe* object below horizon. However, for the *inverse* case, using atmospheric_refraction() to iterate on this function, we may start slightly below h=0.
1082 ! 1.102° is the maximum possible refraction on Earth, for T=-90°C and P=1085 mbar, according to refract().
1083
1084 ! Original algorithm was in degrees and kelvin and used the zenith angle as input iso altitude:
1085 z0 = (pio2 - alt0)*r2d ! alt -> za; rad -> deg
1086
1087 ! Deal with optional parameters:
1088 if(present(t0)) then ! Temperature
1089 t0l = t0 + 273.15d0 ! Use user-specified value; °C -> K
1090 else
1091 t0l = 283.15d0 ! 10°C -> K by default
1092 end if
1093 if(present(p0)) then ! Air pressure
1094 p0l = p0 ! Use user-specified value
1095 else
1096 p0l = 1010 ! 1010 mbar by default
1097 end if
1098 if(present(rh)) then ! Relative humidity
1099 rhl = rh ! Use user-specified value
1100 else
1101 rhl = 0.5 ! 50% by default
1102 end if
1103 if(present(lam)) then ! Wavelength
1104 laml = lam*1.d-3 ! Use user-specified value; nm -> μm
1105 else
1106 laml = 550*1.d-3 ! 550 nm -> μm by default
1107 end if
1108 if(present(dtdh)) then ! Temperature gradient
1109 dtdhl = dtdh ! Use user-specified value
1110 else
1111 dtdhl = 6.5d-3 ! Use default value: 6.5e-3 K/m
1112 end if
1113 if(present(eps)) then ! Desired accuracy
1114 epsl = eps ! Use user-specified value
1115 else
1116 epsl = 1d-4 ! Use default value: 1e-4 rad ~ 0.006°
1117 end if
1118
1119
1120 ! Always defined:
1121 z = 0.d0; reft=0.d0
1122
1123 ! Set up parameters defined at the observer for the atmosphere:
1124 gb = 9.784d0*(1.d0 - 0.0026d0*cos(2.d0*ph) - 0.00000028d0*h0)
1125 z1 = (287.604d0 + 1.6288d0/(laml**2) + 0.0136d0/(laml**4)) * (273.15d0/1013.25d0)*1.d-6
1126
1127 a(1) = abs(dtdhl)
1128 a(2) = (gb*md)/gcr
1129 a(3) = a(2)/a(1)
1130 a(4) = gamma
1131 pw0 = rhl*(t0l/247.1d0)**a(4)
1132 a(5) = pw0*(1.d0 - mw/md)*a(3)/(a(4)-a(3))
1133 a(6) = p0l + a(5)
1134 a(7) = z1*a(6)/t0l
1135 a(8) = ( z1*a(5) + z2*pw0)/t0l
1136 a(9) = (a(3) - 1.d0)*a(1)*a(7)/t0l
1137 a(10) = (a(4) - 1.d0)*a(1)*a(8)/t0l
1138
1139 ! At the Observer:
1140 earthrm = earthr*1.d-2 ! cm -> m
1141 r0 = earthrm + h0
1142 call troposphere_model(r0,t0l,a,r0,t0o,n0,dndr0)
1143 sk0 = n0 * r0 * sin(z0*d2r)
1144
1145
1146 f0 = refi(r0,n0,dndr0)
1147
1148 ! At the Tropopause in the Troposphere:
1149 rt = earthrm + ht
1150 call troposphere_model(r0,t0l,a,rt,tt,nt,dndrt)
1151 zt = asin(sk0/(rt*nt))*r2d
1152 ft = refi(rt,nt,dndrt)
1153
1154 ! At the Tropopause in the Stratosphere:
1155 call stratosphere_model(rt,tt,nt,a(2),rt,nts,dndrts)
1156 zts = asin(sk0/(rt*nts))*r2d
1157 fts = refi(rt,nts,dndrts)
1158
1159 ! At the stratosphere limit:
1160 rs = earthrm + hs
1161 call stratosphere_model(rt,tt,nt,a(2),rs,ns,dndrs)
1162 zs = asin(sk0/(rs*ns))*r2d
1163 fs = refi(rs,ns,dndrs)
1164
1165 ! Integrate the refraction integral in the troposhere and stratosphere,
1166 ! i.e. total refraction = refraction troposhere + refraction stratopshere
1167
1168 ! Initial step lengths etc:
1169 refo = -huge(refo)
1170 is = 16
1171 do k = 1,2
1172 istart = 0
1173 fe = 0.d0
1174 fo = 0.d0
1175
1176 if(k.eq.1) then
1177 h = (zt - z0)/dble(is)
1178 fb = f0
1179 ff = ft
1180 else if(k.eq.2) then
1181 h = (zs - zts )/dble(is)
1182 fb = fts
1183 ff = fs
1184 end if
1185
1186 in = is - 1
1187 is = is/2
1188 step = h
1189
1190200 continue
1191
1192 do i = 1,in
1193 if(i.eq.1.and.k.eq.1) then
1194 z = z0 + h
1195 r = r0
1196 else if(i.eq.1.and.k.eq.2) then
1197 z = zts + h
1198 r = rt
1199 else
1200 z = z + step
1201 end if
1202
1203 ! Given the Zenith distance (Z) find R:
1204 rg = r
1205 do j = 1,4
1206 if(k.eq.1) then
1207 call troposphere_model(r0,t0l,a,rg,tg,n,dndr)
1208 else if(k.eq.2) then
1209 call stratosphere_model(rt,tt,nt,a(2),rg,n,dndr)
1210 end if
1211 rg = rg - ( (rg*n - sk0/sin(z*d2r))/(n + rg*dndr) )
1212 end do
1213 r = rg
1214
1215 ! Find Refractive index and Integrand at R:
1216 if(k.eq.1) then
1217 call troposphere_model(r0,t0l,a,r,t,n,dndr)
1218 else if(k.eq.2) then
1219 call stratosphere_model(rt,tt,nt,a(2),r,n,dndr)
1220 end if
1221
1222 f = refi(r,n,dndr)
1223
1224 if(istart.eq.0.and.mod(i,2).eq.0) then
1225 fe = fe + f
1226 else
1227 fo = fo + f
1228 end if
1229 end do
1230
1231 ! Evaluate the integrand using Simpson's Rule:
1232 refp = h*(fb + 4.d0*fo + 2.d0*fe + ff)/3.d0
1233
1234 if(abs(refp-refo).gt.0.5d0*epsl*r2d) then
1235 is = 2*is
1236 in = is
1237 step = h
1238 h = h/2.d0
1239 fe = fe + fo
1240 fo = 0.d0
1241 refo = refp
1242 if(istart.eq.0) istart = 1
1243 goto 200
1244 end if
1245 if(k.eq.1) reft = refp
1246 end do
1247
1248 atmospheric_refraction_apparent = (reft + refp)*d2r ! Convert result from degrees to radians
1249
1251 !*********************************************************************************************************************************
1252
1253
1254 !*********************************************************************************************************************************
1255 !> \brief Compute the atmospheric refraction of light for a given *apparent* (observed) zenith angle.
1256 !! The method is based on N.A.O Technical Notes 59 and 63 and a paper by Standish and Auer 'Astronomical Refraction:
1257 !! Computational Method for all Zenith Angles'. Return 0 if z0>91.102°.
1258 !!
1259 !! \param z0 The observed zenith distance of the object in *degrees*(!)
1260 !!
1261 !! \param h0 The height of the observer above sea level in metres
1262 !! \param ph The latitude of the observer in *degrees*(!)
1263 !!
1264 !! \param t0 The temperature at the observer in Kelvin
1265 !! \param p0 The pressure at the observer in millibars
1266 !! \param rh The relative humidity at the observer as a fraction (0-1)
1267 !!
1268 !! \param lam The wavelength of the light at the observer in micrometres
1269 !! \param dTdh The temperature lapse rate dT/dh in Kelvin/metre in the troposphere; the absolute value is used
1270 !!
1271 !! \param eps The desired precision in *arcseconds*(!)
1272 !!
1273 !! \retval aref The refraction at the observer in *degrees*(!)
1274 !!
1275 !! \see Hohenkerk & Sinclair, HMNAO technical note 63 (1985)
1276
1277 function aref(z0, h0,ph, t0,p0,rh, lam,dTdh,eps)
1278 use sufr_kinds, only: double
1279 use sufr_constants, only: d2r,r2d
1280
1281 implicit none
1282 real(double), intent(in) :: z0, h0,ph, t0,p0, rh,lam, dtdh, eps
1283 real(double) :: aref, alt0
1284
1285 alt0 = (90.d0-z0)*d2r ! za -> alt; deg -> rad
1286 ! input: deg->rad, K->°C, μm->nm; output: rad -> deg:
1287 aref = atmospheric_refraction_apparent(alt0, h0,ph*d2r, t0-273.15d0,p0,rh, lam*1.d3, dtdh, eps*d2r)*r2d
1288
1289 end function aref
1290 !*********************************************************************************************************************************
1291
1292
1293 !*********************************************************************************************************************************
1294 !> \brief The refraction integrand for atmospheric_refraction_apparent()
1295 !!
1296 !! \param r The current distance from the centre of the Earth in metres
1297 !! \param n The refractive index at R
1298 !! \param dndr The rate the refractive index is changing at R
1299 !!
1300 !! \retval refi The integrand of the refraction function
1301
1302 function refi(r,n,dndr)
1303 use sufr_kinds, only: double
1304 implicit none
1305 real(double), intent(in) :: r,n,dndr
1306 real(double) :: refi
1307
1308 refi = r*dndr/(n + r*dndr)
1309 end function refi
1310 !*********************************************************************************************************************************
1311
1312
1313 !*********************************************************************************************************************************
1314 !> \brief Troposphere model for atmospheric_refraction_apparent()
1315 !!
1316 !! \param r0 The height of the observer from the centre of the Earth
1317 !! \param t0 The temperature at the observer in Kelvin
1318 !! \param a Constants defined at the observer
1319 !! \param r The current distance from the centre of the Earth in metres
1320 !!
1321 !! \param t The temperature at R in Kelvin (output)
1322 !! \param n The refractive index at R (output)
1323 !! \param dndr The rate the refractive index is changing at R (output)
1324
1325 subroutine troposphere_model(r0,t0, a,r, t,n,dndr)
1326 use sufr_kinds, only: double
1327 implicit none
1328 real(double), intent(in) :: r0,t0, a(10),r
1329 real(double), intent(out) :: t,n,dndr
1330 real(double) :: tt0, tt01,tt02
1331
1332 t = t0 - a(1)*(r-r0)
1333 tt0 = t/t0
1334 tt01 = tt0**(a(3)-2.d0)
1335 tt02 = tt0**(a(4)-2.d0)
1336
1337 n = 1.d0 + ( a(7)*tt01 - a(8)*tt02 )*tt0
1338 dndr = -a(9)*tt01 + a(10)*tt02
1339
1340 end subroutine troposphere_model
1341 !*********************************************************************************************************************************
1342
1343
1344 !*********************************************************************************************************************************
1345 !> \brief Stratosphere model for atmospheric_refraction_apparent()
1346 !!
1347 !! \param rt The height of the tropopause from the centre of the Earth in metres
1348 !! \param tt The temperature at the tropopause in Kelvin
1349 !! \param nt The refractive index at the tropopause
1350 !! \param a Constant of the atmospheric model = G*MD/R
1351 !! \param r The current distance from the centre of the Earth in metres
1352 !!
1353 !! \param n The refractive index at R (output)
1354 !! \param dndr The rate the refractive index is changing at R (output)
1355
1356 subroutine stratosphere_model(rt,tt,nt, a,r, n,dndr)
1357 use sufr_kinds, only: double
1358 implicit none
1359 real(double), intent(in) :: rt,tt,nt, a,r
1360 real(double), intent(out) :: n,dndr
1361 real(double) :: b
1362
1363 b = a/tt
1364 n = 1.d0 + (nt - 1.d0)*exp(-b*(r-rt))
1365 dndr = -b*(nt-1.d0)*exp(-b*(r-rt))
1366
1367 end subroutine stratosphere_model
1368 !*********************************************************************************************************************************
1369
1370end module thesky_coordinates
1371!***********************************************************************************************************************************
1372
Procedures for coordinates.
subroutine precess_xyz(jd1, jd2, x, y, z)
Compute the precession of the equinoxes in rectangular coordinates, from jd1 to jd2.
real(double) function refract(alt, press, temp)
Compute the atmospheric refraction for a given true altitude. You should add the result to the uncorr...
subroutine eq2horiz(ra, dec, agst, hh, az, alt, lat, lon)
Convert spherical equatorial coordinates (RA, dec, agst) to spherical horizontal coordinates (hh,...
real(double) function refi(r, n, dndr)
The refraction integrand for atmospheric_refraction_apparent()
subroutine horiz2eq(az, alt, agst, hh, ra, dec, lat, lon)
Convert spherical horizontal coordinates (az, alt, agst) to spherical equatorial coordinates (hh,...
subroutine eq_2_ecl(ra, dec, eps, l, b)
Convert (geocentric) spherical equatorial coordinates RA, Dec (and eps) to spherical ecliptical coord...
real(double) function atmospheric_refraction(alt0, h0, lat0, t0, p0, rh, lam, dtdh, eps)
Compute the atmospheric refraction of light for a given true altitude. Return 0 for alt<-0....
subroutine stratosphere_model(rt, tt, nt, a, r, n, dndr)
Stratosphere model for atmospheric_refraction_apparent()
subroutine gal2eq(l, b, ra, dec)
Convert spherical galactic coordinates (l,b) to spherical equatorial coordinates (RA,...
real(double) function atmospheric_refraction_apparent(alt0, h0, ph, t0, p0, rh, lam, dtdh, eps)
Compute the atmospheric refraction of light for a given apparent (observed) altitude of a celestial o...
subroutine aberration_eq(jd, ra, dec, dra, ddec, eps0)
Correct equatorial coordinates for annual aberration - moderate accuracy, use for stars.
subroutine precess_orb(jd1, jd2, i, o1, o2)
Compute the precession of the equinoxes in orbital elements.
real function refract_sp(alt, press, temp)
Compute the atmospheric refraction for a given true altitude, using single-precision values....
subroutine ecl_spher_2_eq_rect(l, b, r, eps, x, y, z)
Convert spherical, ecliptical coordinates to rectangular, equatorial coordinates of an object - both ...
subroutine eq2gal(ra, dec, l, b)
Convert spherical equatorial coordinates (RA, dec) to spherical galactic coordinates (l,...
subroutine rect_2_spher(x, y, z, l, b, r)
Convert rectangular coordinates x,y,z to spherical coordinates l,b,r.
subroutine precess_eq(jd1, jd2, a1, d1)
Compute the precession of the equinoxes in equatorial coordinates, from jd1 to jd2.
subroutine ecl_2_eq(l, b, eps, ra, dec)
Convert (geocentric) spherical ecliptical coordinates l,b (and eps) to spherical equatorial coordinat...
real(double) function aref(z0, h0, ph, t0, p0, rh, lam, dtdh, eps)
Compute the atmospheric refraction of light for a given apparent (observed) zenith angle....
subroutine fk5(t, l, b)
Convert coordinates to the FK5 system.
subroutine hc_spher_2_gc_rect(l, b, r, l0, b0, r0, x, y, z)
Compute the geocentric rectangular coordinates of a planet, from its and the Earth's heliocentric sph...
subroutine calcsunxyz(t1, l0, b0, r0, x, y, z)
Compute the geocentric equatorial rectangular coordinates of the Sun, from Earth's heliocentric,...
subroutine aberration_ecl(t, l0, l, b)
Correct ecliptic longitude and latitiude for annual aberration.
subroutine troposphere_model(r0, t0, a, r, t, n, dndr)
Troposphere model for atmospheric_refraction_apparent()
subroutine precess_eq_yr(yr1, yr2, a1, d1)
Compute the precession of the equinoxes in equatorial coordinates, from yr1 to yr2.
subroutine precess_ecl(jd1, jd2, l, b)
Compute the precession of the equinoxes in geocentric ecliptical coordinates.
subroutine geoc2topoc_eq(gcra, gcd, gcr, gch, tcra, tcd, lat, hgt)
Convert geocentric equatorial coordinates to topocentric.
subroutine geoc2topoc_ecl(gcl, gcb, gcr, gcs, eps, lst, tcl, tcb, tcs, lat, hgt)
Convert spherical ecliptical coordinates from the geocentric to the topocentric system.
Local parameters for libTheSky: location, date, time.
Definition modules.f90:56
real(double) height
Altitude of the observer above sea level (m)
Definition modules.f90:65
real(double) lon0
Longitude of the observer (rad)
Definition modules.f90:64
real(double) lat0
Latitude of the observer (rad)
Definition modules.f90:63
Procedures for nutation.
Definition nutation.f90:25
subroutine nutation(t, dpsi, eps0, deps)
Calculate nutation - cheap routine from Meeus - as well as the mean obliquity of the ecliptic.
Definition nutation.f90:46