libTheSky
Routines to compute sky positions of Sun, Moon, planets and more
Loading...
Searching...
No Matches
planets.f90
Go to the documentation of this file.
1!> \file planets.f90 Procedures to compute planet positions and more for libTheSky
2
3
4! Copyright (c) 2002-2023 Marc van der Sluys - marc.vandersluys.nl
5!
6! This file is part of the libTheSky package,
7! see: http://libthesky.sf.net/
8!
9! This is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License
10! as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
11!
12! This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied
13! warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
14!
15! You should have received a copy of the GNU General Public License along with this code. If not, see
16! <http://www.gnu.org/licenses/>.
17
18
19!***********************************************************************************************************************************
20!> \brief Procedures for planets
21
23 implicit none
24 save
25
26
27contains
28
29
30 !*********************************************************************************************************************************
31 !> \brief Compute the position, distance, etc of a planet
32 !!
33 !! \param jd Julian date of computation
34 !! \param pl Planet number: Moon=0, Merc=1,Nep=8,Plu=9 3=Sun, 0=Moon, >10 for other objects
35 !!
36 !! \param lat Latitude of the observer (rad, optional)
37 !! \param lon Longitude of the observer (rad, optional)
38 !! \param hgt Altitude/elevation of the observer above sea level (metres, optional)
39 !!
40 !! \param LBaccur Desired accuracy of the heliocentric L,B in VSOP87 (rad, optional; defaults to full accuracy)
41 !! \param Raccur Desired accuracy of the heliocentric R in VSOP87 (AU, optional; defaults to full accuracy)
42 !!
43 !! \param ltime Correct for light time (optional; defaults to true). .false. saves ~50% in CPU time at the cost of some accuracy.
44 !! \param aber Correct for aberration (optional; defaults to true).
45 !! \param to_fk5 Convert coordinates to the FK5 system (optional; defaults to true).
46 !! \param assume_jde Assume JD provided is actually JDE (optional; defaults to false).
47 !!
48 !! \param lunar_theory Choose Lunar theory: 1: ELP82b, 2: ELP-MPP02/LLR, 3: ELP-MPP02/DE405 ('historical' - default)
49 !! \param nutat IAU nutation model to use: 0 (no nutation!), 1980 or 2000 (default).
50 !! \param magmdl Planet-magnitude model: 1: Müller (1893), 2: Meeus p.286, 3: AA 1992, 4: AA 2012 (optional, defaults to 4).
51 !! \param verbosity Verbosity for debug output (0-3). Defaults to 0: silent.
52 !!
53 !!
54 !! \note
55 !! - lat0 and lon0 can be provided through the module TheSky_local (rad, rad and m), or through the optional arguments.
56 !! Note that using the latter will update the former!
57 !! - results are returned in the array planpos() in the module TheSky_planetdata
58 !!
59
60 subroutine planet_position(jd,pl, lat,lon,hgt, LBaccur,Raccur, ltime,aber,to_fk5,assume_jde, &
61 lunar_theory,nutat,magmdl, verbosity)
62
63 use sufr_kinds, only: double
64 use sufr_constants, only: pi,pi2, r2d,r2h, au,earthr,pland, enpname, jd2000
65 use sufr_system, only: warn, quit_program_error
66 use sufr_angles, only: rev, rev2
67 use sufr_dummy, only: dumdbl1,dumdbl2
68 use sufr_text, only: d2s
69 use sufr_numerics, only: deq0
70
71 use thesky_vsop, only: vsop87d_lbr
76 use thesky_sun, only: sunmagn
78 use thesky_comets, only: cometgc
79 use thesky_planetdata, only: planpos, pl0
83
84 implicit none
85 real(double), intent(in) :: jd
86 integer, intent(in) :: pl
87 real(double), intent(in), optional :: lat,lon,hgt, LBaccur,Raccur
88 logical, intent(in), optional :: ltime,aber,to_fk5, assume_jde
89 integer, intent(in), optional :: lunar_theory, nutat,magmdl, verbosity
90
91 integer :: j, llunar_theory, lnutat,lmagmdl, lverb
92 real(double) :: jdl,jde,jde_lt, tjm,tjm0, llat,llon,lhgt, lLBaccur,lRaccur, dpsi,eps0,deps,eps,tau,tau1
93 real(double) :: hcl0,hcb0,hcr0, hcl,hcb,hcr, hcl00,hcb00,hcr00, sun_gcl,sun_gcb, gcx,gcy,gcz, gcx0,gcy0,gcz0, dhcr
94 real(double) :: gcl,gcb,delta,gcl0,gcb0,delta0
95 real(double) :: ra,dec,gmst,agst,lst,hh,az,alt,elon, topra,topdec,topl,topb,topdiam,topdelta,tophh,topaz,topalt
96 real(double) :: diam,illfr,pa,magn, parang,parang_ecl,hp, rES1,rES2
97 logical :: lltime,laber,lto_fk5
98
99 ! Make sure these are always defined:
100 gcl0 = 0.d0; gcb0 = 0.d0; hcr00 = 0.d0
101 gcx = 0.d0; gcy = 0.d0; gcz = 0.d0
102 gcx0 = 0.d0; gcy0 = 0.d0; gcz0 = 0.d0
103 topdelta = 0.d0; delta0 = 0.d0
104 diam = 0.d0
105 magn = 0.d0 ! Make sure magn is always defined
106
107 pl0 = pl ! Remember which planet was computed last
108
109
110 ! Handle optional variables:
111 ! Geographic location:
112 llat = lat0
113 llon = lon0
114 lhgt = height
115 if(present(lat)) llat = lat
116 if(present(lon)) llon = lon
117 if(present(hgt)) lhgt = hgt
118
119 ! Coordinate accuracy:
120 llbaccur = 0.d0 ! Seting L,B,R accuracy equal to VSOP87 truncation still skips some terms!
121 lraccur = 0.d0
122 if(present(lbaccur)) llbaccur = lbaccur
123 if(present(raccur)) lraccur = raccur
124
125 ! Light time, aberration and FK5 switches:
126 lltime = .true. ! Correct for light time by default
127 if(present(ltime)) lltime = ltime
128
129 laber = .true. ! Correct for aberration by default
130 if(present(aber)) laber = aber
131
132 lto_fk5 = .true. ! Convert to FK5 by default
133 if(present(to_fk5)) lto_fk5 = to_fk5
134
135 ! Select lunar theory:
136 ! llunar_theory = 1 ! Use (corrected) ELP82b as default
137 llunar_theory = 3 ! Use ELP-MPP02/DE405 ('historical') as default. This should be the "default default" :-)
138 if(present(lunar_theory)) llunar_theory = lunar_theory
139 if(llunar_theory.lt.1 .or. llunar_theory.gt.3) &
140 call quit_program_error('planet_position(): lunar_theory must be 1, 2 or 3', 1)
141
142 ! Nutation model:
143 lnutat = 2000 ! 2000 IAU nutation model
144 if(present(nutat)) then
145 lnutat = nutat
146 if(lnutat.ne.1980 .and. lnutat.ne.2000 .and. lnutat.ne.0) &
147 call quit_program_error('planet_position(): nutat must be 2000, 1980 or 0', 1)
148 end if
149
150 ! Planet-magnitude model:
151 lmagmdl = 4 ! AA 2012
152 if(present(magmdl)) then
153 lmagmdl = magmdl
154 if(lmagmdl.lt.1 .or. lmagmdl.gt.4) &
155 call quit_program_error('planet_position(): magmdl must be 1-4', 1)
156 end if
157
158 ! Verbosity:
159 lverb = 0 ! Silent
160 if(present(verbosity)) lverb = verbosity
161
162
163 ! Calc JDE and tjm:
164 jdl = jd ! Local version of JD
165 deltat = calc_deltat(jdl)
166 ! deltat = 0.d0
167 ! write(*,'(/,A,/)') '*** WARNING: libTheSky/planets.f90: setting DeltaT to 0 ***'
168
169 if(present(assume_jde) .and. assume_jde) then
170 jde = jdl ! JD provided is actually JDE
171 jdl = jde - deltat/86400.d0
172 else
173 jde = jdl + deltat/86400.d0
174 end if
175
176 tjm = (jde-jd2000)/365250.d0 ! Julian Millennia after 2000.0 in dyn. time, tau in Meeus
177 tjm0 = tjm
178 if(lverb.gt.0) then
179 print*,'JD: ', d2s(jdl,8) ! ms accuracy
180 print*,'DeltaT: ', d2s(deltat,3) ! ms accuracy
181 print*,'JDE: ', d2s(jde,8) ! ms accuracy
182 print*,'t_jm: ', d2s(tjm,10) ! ~s accuracy
183 end if
184
185 call vsop87d_lbr(tjm,3, hcl0,hcb0,hcr0, llbaccur,lraccur) ! Calculate the Earth's true heliocentric l,b,r
186
187 if(lverb.gt.3) then
188 print*
189 ! print*, 'Heliocentric longitude Earth, FK4: ', d2s(modulo(hcl0,pi2)*r2d, 9), ' deg'
190 ! print*, 'Heliocentric latitude Earth, FK4: ', d2s(hcb0*r2d, 9), ' deg'
191 ! print*, 'Heliocentric distance Earth: ', d2s(hcr0, 9), ' au'
192 print*
193 print*, 'Earth:'
194 print*, 't_jm: ', d2s(tjm, 10)
195 print*, 'Heliocentric longitude Earth: ', d2s(modulo(hcl0,pi2)*r2d, 9), ' deg'
196 print*, 'Heliocentric latitude Earth: ', d2s(hcb0*r2d, 9), ' deg'
197 print*, 'Heliocentric distance Earth: ', d2s(hcr0, 9), ' au'
198 end if
199
200
201 ! Iterate to get the light time tau, hence apparent positions:
202 tau = 0.d0; tau1 = 0.d0 ! Initial values for light time in days
203 j = 0 ! Count iterations; allows escape in case of infinite loop
204 do while(deq0(tau1) .or. abs((tau1-tau)/tau1).gt.1.d-10) ! Moon's tau ~10^-5; 1.d-10~10^-5 sec (1.d-7~10^-2 sec)
205
206 tau = tau1 ! On first iteration, tau=0 to get true positions
207 tjm = tjm0 - tau/365250.d0 ! Iterate to calculate light time
208 jde_lt = jd2000 + tjm*365250 ! JDE, corrected for light time
209 hcl = 0.d0; hcb = 0.d0; hcr = 0.d0 ! Heliocentric coordinates
210
211
212 ! Compute l,b,r for planet, Pluto, asteroid, Earth's shadow:
213 if(pl.gt.0.and.pl.lt.9) call vsop87d_lbr(tjm,pl, hcl,hcb,hcr, llbaccur,lraccur) ! Heliocentric l,b,r
214 if(pl.eq.9) call plutolbr(tjm*10.d0, hcl,hcb,hcr) ! This is for 2000.0, precess 10 lines below
215 if(pl.gt.10000) call asteroid_lbr(tjm,pl-10000, hcl,hcb,hcr) ! Heliocentric lbr for asteroids
216 if(pl.eq.-1) then ! Centre of Earth's shadow
217 call vsop87d_lbr(tjm,3, hcl,hcb,hcr, llbaccur,lraccur) ! = heliocentric coordinates of the Earth...
218
219 select case(llunar_theory)
220 case(1)
221 call elp82b_lbr(tjm, dumdbl1,dumdbl2,dhcr) ! only want dhcr
222 case(2)
223 call elp_mpp02_lbr(jde_lt, 0, dumdbl1,dumdbl2,dhcr) ! only want dhcr - 0: LLR mode
224 case(3)
225 call elp_mpp02_lbr(jde_lt, 1, dumdbl1,dumdbl2,dhcr) ! only want dhcr - 1: DE405 ('historical') mode
226 end select
227
228 hcr = hcr + dhcr ! + distance Earth-Moon
229 end if
230
231
232 ! Compute geocentric ecliptical position for the Moon:
233 if(pl.eq.0) then ! Moon
234 select case(llunar_theory)
235 case(1)
236 call elp82b_lbr(tjm, gcl,gcb,delta) ! Get apparent geocentric coordinates of the Moona
237 case(2)
238 call elp_mpp02_lbr(jde_lt, 0, gcl,gcb,delta) ! Get apparent geocentric coordinates of the Moon. Mode=0: LLR mode
239 case(3)
240 call elp_mpp02_lbr(jde_lt, 1, gcl,gcb,delta) ! Get apparent geocentric coordinates of the Moon. Mode=1: DE405 ('historical') mode
241 end select
242 end if
243
244 if(pl.ne.0.and.pl.lt.10 .or. pl.gt.10000) then ! Planet, asteroid, Earth's shadow
245 ! For Neptune's birthday:
246 ! print*,'TESTING!!!'
247 ! call precess_ecl(jde,jd2000,hcl,hcb) ! from JoD to J2000.0
248
249 call hc_spher_2_gc_rect(hcl,hcb,hcr, hcl0,hcb0,hcr0, gcx,gcy,gcz) ! Convert heliocentric l,b,r to geocentric x,y,z
250 call rect_2_spher(gcx,gcy,gcz, gcl,gcb,delta) ! Convert geocentric x,y,z to geocentric l,b,r
251 if(pl.eq.3) delta = hcr
252 if(pl.eq.9) call precess_ecl(jd2000,jde, gcl,gcb) ! Pluto: from J2000.0 to JoD
253 end if
254
255 if(pl.gt.10.and.pl.lt.10000) &
256 call cometgc(tjm,tjm0, pl, hcr,gcl,gcb,delta) ! Calc geocentric l,b,r coordinates for a comet
257
258
259 ! Store 'true' position, for tau=0:
260 if(j.eq.0) then
261 if(abs(tau).gt.1.d-10) call warn('planet_position(): tau != 0 on first light-time iteration', 0)
262
263 hcl00 = hcl ! Heliocentric l,b,r
264 hcb00 = hcb
265 hcr00 = hcr
266
267 gcx0 = gcx ! Geocentric x,y,z
268 gcy0 = gcy
269 gcz0 = gcz
270
271 delta0 = delta ! Geocentric l,b,r
272 gcl0 = gcl
273 gcb0 = gcb
274 end if
275
276 tau1 = 5.77551830441d-3 * delta ! Light time in days
277
278 if(lverb.gt.3) then
279 print*
280 print*
281 print*, 'Iteration: ', j
282 ! print*,'t_jm: ', d2s(tjm, 10)
283 ! print*
284 ! print*, 'Heliocentric longitude planet, FK4: ', d2s(modulo(hcl00,pi2)*r2d, 9), ' deg'
285 ! print*, 'Heliocentric latitude planet, FK4: ', d2s(hcb00*r2d, 9), ' deg'
286 ! print*, 'Heliocentric distance planet: ', d2s(hcr00, 9), ' au'
287 print*
288 print*, 'Planet:'
289 print*, 't_jm: ', d2s(tjm, 10)
290 print*, 'Heliocentric longitude: ', d2s(modulo(hcl,pi2)*r2d, 9), ' deg'
291 print*, 'Heliocentric latitude: ', d2s(hcb*r2d, 9), ' deg'
292 print*, 'Heliocentric distance: ', d2s(hcr, 9), ' au'
293 print*
294 print*, 'Geocentric x: ', d2s(gcx, 9), ' au'
295 print*, 'Geocentric y: ', d2s(gcy, 9), ' au'
296 print*, 'Geocentric z: ', d2s(gcz, 9), ' au'
297 print*
298 print*, 'Geocentric longitude: ', d2s(modulo(gcl,pi2)*r2d, 9), ' deg'
299 print*, 'Geocentric latitude: ', d2s(gcb*r2d, 9), ' deg'
300 print*, 'True geocentric distance: ', d2s(delta0, 9), ' au'
301 print*, 'Apparent geocentric distance: ', d2s(delta, 9), ' au'
302 print*
303 print*, 'JDE_i: ', d2s(jde_lt, 9)
304 print*, 't_jm: ', d2s(tjm, 10)
305 print*, 'Light time: ', d2s(tau,9), ' day'
306 print*, 'Light time1: ', d2s(tau1,9), ' day'
307 print*
308 end if
309
310 if(.not.lltime) exit ! Do not take into account light time
311 if(pl.eq.3) exit ! Want true position for Sun(!)
312
313 j = j+1
314 if(j.ge.30) then
315 call warn('planet_position(): Light time failed to converge for '//trim(enpname(pl)), 0)
316 exit
317 end if
318 end do
319
320 tau = tau1
321 tjm = tjm0 ! Still Julian Millennia since 2000.0 in dynamical time
322 delta = delta0 ! Want the TRUE distance, not APPARENT!
323
324 if(lverb.gt.3) then
325 print*
326 print*, 'Number of iterations: ', j
327 print*, 'Final light time: ', d2s(tau, 9)
328 print*, 'Final t_Jm: ', d2s(tjm, 10)
329 print*
330 end if
331
332 if(pl.eq.-1) then ! Earth's shadow
333 select case(llunar_theory)
334 case(1)
335 call elp82b_lbr(tjm, dumdbl1,dumdbl2,delta) ! Geocentric distance of the Moon for Earth shadow; only need delta
336 case(2)
337 call elp_mpp02_lbr(jde_lt, 0, dumdbl1,dumdbl2,delta) ! Geocentric distance of the Moon for Earth shadow; only need delta - LLR mode
338 case(3)
339 call elp_mpp02_lbr(jde_lt, 1, dumdbl1,dumdbl2,delta) ! Geocentric distance of the Moon for Earth shadow; only need delta - DE405 ('historical') mode
340 end select
341 end if
342
343 ! Earth -> Sun. Note that we want the TRUE position!
344
345 ! (As for other planets, we'd use the true position for the Earth (when the light arrives) and the
346 ! APPARENT position of the Sun (when the light left). However, since the positions are heliocentric, the
347 ! latter is always (0,0,0).
348 if(pl.eq.3) then
349 gcl = rev(hcl0+pi)
350 gcb = -hcb0
351 delta = hcr0
352 end if
353
354
355 if(lverb.gt.2) then
356 print*
357 print*
358 ! print*,'t_jm: ', d2s(tjm, 10)
359 ! print*
360 ! print*, 'Heliocentric longitude planet, FK4: ', d2s(modulo(hcl00,pi2)*r2d, 9), ' deg'
361 ! print*, 'Heliocentric latitude planet, FK4: ', d2s(hcb00*r2d, 9), ' deg'
362 ! print*, 'Heliocentric distance planet: ', d2s(hcr00, 9), ' au'
363 print*
364 print*, 'Planet after light-time iterations, before correction for nutation, aberration, FK5:'
365 print*, 'Final t_jm: ', d2s(tjm, 10)
366 print*, 'Final heliocentric longitude: ', d2s(modulo(hcl,pi2)*r2d, 9), ' deg'
367 print*, 'Final heliocentric latitude: ', d2s(hcb*r2d, 9), ' deg'
368 print*, 'Final heliocentric distance: ', d2s(hcr, 9), ' au'
369 print*
370 print*, 'Final geocentric x: ', d2s(gcx, 9), ' au'
371 print*, 'Final geocentric y: ', d2s(gcy, 9), ' au'
372 print*, 'Final geocentric z: ', d2s(gcz, 9), ' au'
373 print*
374 print*, 'Final geocentric longitude: ', d2s(modulo(gcl,pi2)*r2d, 9), ' deg'
375 print*, 'Final geocentric latitude: ', d2s(gcb*r2d, 9), ' deg'
376 print*, 'Final (true) geocentric distance: ', d2s(delta, 9), ' au'
377 print*
378 print*, 'Final jde_i: ', d2s(jde_lt, 9)
379 print*, 'Final t_jm: ', d2s(tjm, 10)
380 print*, 'Final light time: ', d2s(tau,9), ' day'
381 print*, 'Final light time1: ', d2s(tau1,9), ' day'
382 print*
383 end if
384
385 ! Compute and correct for nutation:
386 if(lnutat.eq.1980) then
387 call nutation(tjm, dpsi,eps0,deps) ! IAU 1980 nutation model: dpsi: nutation in longitude, deps: in obliquity
388 else ! if(lnutat.eq.2000) then
389 call nutation2000(jdl, dpsi,deps, eps0) ! IAU 2000 nutation model (originally doesn't provide eps0)
390
391 if(lnutat.eq.0) then ! No nutation(!)
392 dpsi = 0.d0
393 deps = 0.d0
394 end if
395 end if
396 eps = eps0 + deps ! Correct for nutation in obliquity: mean -> true obliquity of the ecliptic
397
398
399 ! Correct for aberration and nutation, and convert to FK5:
400 if(pl.lt.10) then ! I suppose this should also happen for comets, but this compares better...
401 if(laber .and. (pl.ne.0)) call aberration_ecl(tjm,hcl0, gcl,gcb) ! Aberration - not for the Moon
402 if(lto_fk5) call fk5(tjm, gcl,gcb)
403 gcl = gcl + dpsi ! Nutation in longitude
404 end if
405
406 ! sun_gcl,sun_gcb give gcl,gcb as if pl.eq.3
407 sun_gcl = rev(hcl0+pi)
408 sun_gcb = -hcb0
409
410 ! Aberration, FK5 and nutation for sun_gcl,sun_gcb:
411 if(laber) call aberration_ecl(tjm,hcl0, sun_gcl,sun_gcb)
412 if(lto_fk5) call fk5(tjm, sun_gcl,sun_gcb)
413 sun_gcl = sun_gcl + dpsi ! Nutation in longitude
414
415 ! FK5 conversion for l,b, hcl0,hcb0, hcl00,hcb00:
416 if(lto_fk5) then
417 call fk5(tjm, hcl,hcb)
418 call fk5(tjm, hcl0,hcb0)
419 call fk5(tjm, hcl00,hcb00)
420 end if
421
422 if(lverb.gt.3) then
423 print*
424 print*, 'Convert to the FK5 system: ', lto_fk5
425 print*, 'True heliocentric longitude, FK5: ', d2s(modulo(hcl00,pi2)*r2d, 9), ' deg'
426 print*, 'True heliocentric latitude, FK5: ', d2s(hcb00*r2d, 9), ' deg'
427 print*, 'True heliocentric distance: ', d2s(hcr00, 9), ' au'
428 print*
429 end if
430
431
432 ! Convert geocentric ecliptical to equatorial coordinates:
433 call ecl_2_eq(gcl,gcb,eps, ra,dec) ! RA, Dec
434
435
436 ! Sidereal time:
437 gmst = calc_gmst(jdl) ! Greenwich mean sidereal time
438 agst = rev(gmst + dpsi*cos(eps)) ! Correction for equation of the equinoxes -> Greenwich apparent sidereal time
439 lst = rev(agst + llon) ! Local apparent sidereal time, llon > 0 for E
440
441
442 ! Apparent diameter:
443 if(pl.ge.0.and.pl.lt.10) diam = atan(pland(pl)/(delta*au))
444 !if(pl.lt.10.and.pl.ge.0) diam = atan(planr(pl)/(2*delta*au))*2
445 if(pl.eq.-1) then
446 call earthshadow(delta,hcr0, res1,res2) ! Calc Earth shadow radii at Moon distance: (pen)umbra radius - geocentric
447 diam = 2*res1 ! Umbra diameter to planpos(12)
448 end if
449
450 ! Convert heliocentric to topocentric coordinates:
451 call geoc2topoc_ecl(gcl,gcb, delta,diam/2.d0, eps,lst, topl,topb,topdiam, lat=llat,hgt=lhgt) ! Geocentric to topoc: l, b, diam
452 if(pl.eq.-1) topdiam = res2 ! Earth penumbra radius at Moon distance
453 topdiam = 2*topdiam ! Was radius, now diameter
454 if(pl.ge.0.and.pl.lt.10) topdelta = pland(pl)/tan(topdiam)/au
455
456
457 ! Convert equatorial to horizontal coordinates:
458 call eq2horiz(ra,dec,agst, hh,az,alt, lat=llat,lon=llon)
459
460
461 ! Elongation:
462 elon = acos(cos(gcb)*cos(hcb0)*cos(gcl-rev(hcl0+pi)))
463 if(pl.gt.10) elon = acos((hcr0**2 + delta**2 - hcr**2)/(2*hcr0*delta)) ! For comets (only?) - CHECK: looks like phase angle...
464
465 ! Convert topocentric coordinates: ecliptical -> equatorial -> horizontal:
466 ! call geoc2topoc_eq(ra,dec,delta,hh,topra,topdec) ! Geocentric to topocentric
467 call ecl_2_eq(topl,topb,eps, topra,topdec) ! Topocentric l,b -> RA,dec - identical result (AA 1999-01-02/22); probably cheaper
468 call eq2horiz(topra,topdec,agst, tophh,topaz,topalt, lat=llat,lon=llon)
469
470 ! Phase angle:
471 pa = 0.d0 ! Make sure pa is always defined
472 if(pl.eq.0) then
473 pa = atan2( hcr0*sin(elon) , delta - hcr0*cos(elon) ) ! Moon
474 else if(pl.gt.0) then
475 pa = acos( (hcr**2 + delta**2 - hcr0**2) / (2*hcr*delta) ) ! Planets
476 end if
477 illfr = 0.5d0*(1.d0 + cos(pa)) ! Illuminated fraction of the disc
478
479
480 ! Compute magnitude:
481 if(pl.gt.0.and.pl.lt.10) then
482 magn = planet_magnitude(pl,hcr,delta,pa, lmagmdl)
483 if(pl.eq.6) magn = magn + dsatmagn(tjm*10., gcl,gcb) ! Correct Saturn's magnitude for rings
484 end if
485
486 if(pl.eq.0) magn = moonmagn(pa,delta) ! Moon
487 if(pl.gt.10000) magn = asteroid_magn(pl-10000,delta,hcr,pa) ! Asteroid magnitude (valid for |pa|<120°)
488
489 ! Comet magnitude:
490 if(pl.gt.10 .and. pl.lt.10000) then
491 diam = 0.d0
492 if(cometdiedatp(pl).ne.0 .and. jdl.gt.cometelems(pl,7)) then
493 magn = 99.9d0 ! Comet died at perihelion
494 else
495 magn = cometelems(pl,8) + 5*log10(delta) + 2.5d0*cometelems(pl,9)*log10(hcr) ! m = H + 5log(d) + 2.5*G*log(r)
496 magn = magn + comet_scatter_magnitude_correction(pa)
497 end if
498 topdelta = delta
499 end if
500
501
502 ! Parallactic angle:
503 !parang = atan2(sin(hh),tan(llat)*cos(dec) - sin(dec)*cos(hh)) ! Parallactic angle (between zenith and celestial north pole): geocentric - Meeus Eq.14.1
504 parang = atan2(sin(tophh),tan(llat)*cos(topdec) - sin(topdec)*cos(tophh)) ! Parallactic angle (between zenith and celestial north pole): topocentric
505 parang_ecl = atan( (cos(topl)*tan(eps)) / (sin(topl)*sin(topb)*tan(eps) - cos(topb)) ) ! "Parallactic angle" (between ECLIPTICAL and celestial north pole): topocentric - Meeus p.100
506
507 ! Horizontal parallax:
508 hp = asin(earthr/(delta*au)) ! Horizontal parallax
509
510
511 ! Overrule some values:
512 if(pl.eq.0) hcr = 0.d0 ! No heliocentric distance for the Moon
513 if(pl.eq.3) then ! Sun:
514 elon = 0.d0 ! Elongation
515 pa = 0.d0 ! Phase angle
516 magn = sunmagn(delta) ! Magnitude; changes ~0.05 over half a year...
517 illfr = 1.d0 ! Illuminated fraction
518 end if
519
520
521
522
523 ! Save variables:
524 ! Geocentric:
525 planpos = 0.d0 ! 1-12 are geocentric, repeated as topocentric in 21-32
526
527 planpos(1) = rev(gcl) ! Ecliptical longitude
528 planpos(2) = rev2(gcb) ! Ecliptical latitude
529 planpos(3) = hcr ! Distance to the Sun (heliocentric!)
530
531 planpos(4) = delta ! True(!) geocentric distance
532 planpos(5) = rev(ra) ! R.A.
533 planpos(6) = dec ! Declination
534
535 planpos(7) = tau ! Light time in days
536
537 planpos(8) = rev(hh) ! Hour angle
538 planpos(9) = rev(az) ! Azimuth
539 planpos(10) = alt ! Altitude
540
541 planpos(11) = rev(elon) ! Elongation
542 planpos(12) = diam ! Apparent diameter
543
544 if(lverb.gt.1) then
545 print*
546 print*, 'Final apparent geocentric longitude: ', d2s(rev(gcl)*r2d, 9), ' deg'
547 print*, 'Final apparent geocentric latitude: ', d2s(rev2(gcb)*r2d, 9), ' deg'
548 print*, 'Final apparent geocentric distance: ', d2s(delta, 9), ' au'
549 print*, 'Final apparent heliocentric distance: ', d2s(hcr, 9), ' au'
550 print*
551 print*, 'Final apparent geocentric right ascension: ', d2s(rev(ra)*r2h, 9), ' hr'
552 print*, 'Final apparent geocentric declination: ', d2s(dec*r2d, 9), ' deg'
553 print*
554 end if
555
556
557 planpos(13) = magn ! Apparent visual magnitude
558 planpos(14) = illfr ! Illuminated fraction
559 planpos(15) = rev(pa) ! Phase angle
560
561 planpos(16) = parang ! Topocentric parallactic angle (between celestial pole and zenith)
562 planpos(17) = hp ! Horizontal parallax
563 planpos(18) = parang_ecl ! Topocentric ecliptical "parallactic angle" (between celestial and ecliptical pole)
564
565 ! Topocentric:
566 planpos(21) = rev(topl)
567 planpos(22) = rev2(topb)
568 planpos(23) = delta0 ! True geocentric distance (== delta?)
569
570 planpos(24) = topdelta
571 planpos(25) = rev(topra)
572 planpos(26) = topdec
573
574 planpos(28) = rev(tophh) ! Hour angle
575 planpos(29) = rev(topaz)
576 planpos(30) = topalt
577 planpos(31) = topalt + refract(topalt) ! Topocentric altitude, corrected for (cheap) atmospheric refraction
578
579 planpos(32) = topdiam
580
581 ! Heliocentric:
582 planpos(33) = rev(hcl) ! Heliocentric apparent l,b,r
583 planpos(34) = rev2(hcb)
584 planpos(35) = hcr
585
586 planpos(36) = rev(hcl00) ! Heliocentric true l,b,r, converted to FK5
587 planpos(37) = rev2(hcb00)
588 planpos(38) = hcr00
589
590 ! Other variables:
591 planpos(39) = dble(pl) ! Remember for which planet data were computed
592 planpos(40) = jde ! JDE for these data
593
594 planpos(41) = rev(hcl0+pi+dpsi) ! Geocentric, true L,B,R for the Sun, in FK5, corrected for nutation
595 planpos(42) = rev2(-hcb0)
596 planpos(43) = hcr0
597
598 planpos(44) = lst ! Local APPARENT sidereal time
599 planpos(45) = rev(agst) ! Greenwich APPARENT sidereal time (in radians)
600 planpos(46) = tjm0 * 10 ! Apparent dynamical time in Julian Centuries since 2000.0
601
602 planpos(47) = dpsi ! Nutation in longitude
603 planpos(48) = eps ! True obliquity of the ecliptic; corrected for nutation
604 planpos(49) = rev(gmst) ! Greenwich MEAN sidereal time (in radians)
605 planpos(50) = eps0 ! Mean obliquity of the ecliptic; without nutation
606
607 ! Geocentric:
608 planpos(51) = rev(sun_gcl) ! Apparent geocentric longitude of the Sun (variable was treated as if pl.eq.3)
609 planpos(52) = rev2(sun_gcb) ! Apparent geocentric latitude of the Sun (variable was treated as if pl.eq.3)
610
611 planpos(58) = jdl ! JD used for calculation
612 planpos(59) = deltat ! DeltaT used for calculation
613
614 planpos(61) = gcx ! Apparent geocentric x,y,z
615 planpos(62) = gcy
616 planpos(63) = gcz
617
618 planpos(64) = gcx0 ! True geocentric x,y,z
619 planpos(65) = gcy0
620 planpos(66) = gcz0
621
622 planpos(67) = gcl0 ! True geocentric l,b,r
623 planpos(68) = gcb0
624 planpos(69) = delta0 ! (==delta?)
625
626 end subroutine planet_position
627 !*********************************************************************************************************************************
628
629
630
631
632
633
634 !*********************************************************************************************************************************
635 !> \brief Calculate Pluto's position l,b,r at time t
636 !!
637 !! \param t Dynamical time in Julian Centuries after 2000.0
638 !! \param l Heliocentric longitude (rad) (output)
639 !! \param b Heliocentric latitude (rad) (output)
640 !! \param r Heliocentric distance (AU?) (output)
641
642 subroutine plutolbr(t,l,b,r)
643 use sufr_kinds, only: double
644 use sufr_constants, only: d2r
646
647 implicit none
648 real(double), intent(in) :: t
649 real(double), intent(out) :: l,b,r
650 integer :: i
651 real(double) :: j,s,p,a, cosa,sina
652
653 j = 34.35d0 + 3034.9057d0 * t
654 s = 50.08d0 + 1222.1138 * t
655 p = 238.96d0 + 144.96d0 * t
656
657 l = 0.d0
658 b = 0.d0
659 r = 0.d0
660
661 do i=1,43
662 a = (pluc(i,1)*j + pluc(i,2)*s + pluc(i,3)*p) * d2r
663 cosa = cos(a)
664 sina = sin(a)
665 l = l + plul(i,1)*sina + plul(i,2)*cosa
666 b = b + plub(i,1)*sina + plub(i,2)*cosa
667 r = r + plur(i,1)*sina + plur(i,2)*cosa
668 end do
669
670 l = ( l*1.d-6 + 238.958116d0 + 144.96d0*t ) * d2r
671 b = ( b*1.d-6 - 3.908239d0 ) * d2r
672 r = r*1.d-7 + 40.7241346d0
673
674 end subroutine plutolbr
675 !*********************************************************************************************************************************
676
677
678
679
680 !*********************************************************************************************************************************
681 !> \brief Calculate orbital elements for the planets
682 !!
683 !! \param jd Julian day of calculation
684 !!
685 !! \note
686 !! Results are returned through the module planetdata, for planet pl and element el:
687 !! - plelems(pl,el): EoD
688 !! - plelems2000(pl,el): J2000.0
689 !!
690 !! \note
691 !! Elements (el):
692 !! - 1: L - mean longitude
693 !! - 2: a - semi-major axis
694 !! - 3: e - eccentricity
695 !! - 4: i - inclination
696 !! - 5: Omega - longitude of ascending node
697 !! - 6: pi - longitude of perihelion
698
699 subroutine planetelements(jd)
700 use sufr_kinds, only: double
701 use sufr_constants, only: d2r, jd2000
702 use sufr_system, only: quit_program_error
703 use sufr_angles, only: rev
704 use sufr_numerics, only: deq0
705
706 use thesky_local, only: deltat
709
710 implicit none
711 real(double), intent(in) :: jd
712 integer :: el,pl,te
713 real(double) :: jde,tm,tms(0:3)
714
715 if(deq0(plelemdata(1,1,1,1))) call quit_program_error('planetelements(): did you forget to call readplanetelements()?', 0)
716
717 deltat = calc_deltat(jd)
718 jde = jd + deltat/86400.d0
719 tm = (jde-jd2000)/36525.d0 ! Julian Centuries after 2000.0 in dynamical time
720
721
722 ! Powers of tm:
723 do te=0,3
724 tms(te) = tm**te
725 end do
726
727
728 do pl=1,8
729 do el=1,6
730
731 plelems(pl,el) = 0.d0
732 plelems2000(pl,el) = 0.d0
733
734 do te=0,3
735 plelems(pl,el) = plelems(pl,el) + plelemdata(1,pl,el,te) * tms(te)
736 plelems2000(pl,el) = plelems2000(pl,el) + plelemdata(2,pl,el,te) * tms(te)
737 end do
738
739 if(el.ne.2.and.el.ne.3) then
740 plelems(pl,el) = rev(plelems(pl,el) * d2r)
741 plelems2000(pl,el) = rev(plelems2000(pl,el) * d2r)
742 end if
743
744 end do
745 end do
746
747 end subroutine planetelements
748 !*********************************************************************************************************************************
749
750
751
752
753 !*********************************************************************************************************************************
754 !> \brief Calculate planet magnitude
755 !!
756 !! \param pl Planet ID
757 !! \param hc_dist Distance from the Sun (AU)
758 !! \param gc_dist Distance from the Earth (AU)
759 !! \param phang Phase angle (rad)
760 !! \param model Model to use: 1: Müller (1893), 2: Meeus p.286, 3: AA 1992, 4: AA 2012 (optional, defaults to 4).
761 !!
762 !! \retval planet_magnitude Apparent visual planet magnitiude
763 !!
764 !! \see
765 !! - Müller, POPot 8, 193 (1893) - https://ui.adsabs.harvard.edu/abs/1893POPot...8..193M/ - p.341/149 - Model 1.
766 !! - Meeus, Astronomical Algorithms, 1998, Ch.41 - Model 2 (origin unsure).
767 !! - Expl.Supl.tt.Astr.Almanac, 2nd Ed. (1992) - Model 3.
768 !! - Expl.Supl.tt.Astr.Almanac, 3rd Ed. (2012) - Model 4 (default).
769
770 function planet_magnitude(pl, hc_dist,gc_dist, phang, model)
771 use sufr_kinds, only: double
772 use sufr_constants, only: r2d
773 use sufr_system, only: quit_program_error
774
775 implicit none
776 integer, intent(in) :: pl
777 integer, intent(in), optional :: model
778 integer :: lpl, lmodel
779 real(double), intent(in) :: hc_dist,gc_dist,phang
780 real(double) :: planet_magnitude,pa,pa2,pa3,a0(9),a1(9),a2(9),a3(9)
781
782 ! Optional argument:
783 lmodel = 4 ! Model 4: Expl.Supl.tt.Astr.Almanac, 3rd Ed. (2012)
784 if(present(model)) lmodel = model
785
786 ! PA must be in degrees:
787 pa = phang*r2d
788 if(lmodel.eq.1 .and. pl.eq.1) pa = pa - 50.d0 ! For Mercury (deg) -- ONLY FOR MODEL 1 !!!
789 pa2 = pa*pa
790 pa3 = pa*pa2
791
792 lpl = pl ! Local pl
793
794 if(lmodel.eq.1) then ! Model 1: Müller (1893)
795 ! Mer Ven Ear Mars Jup Sat Ur Nep Pl
796 a0 = [-1.16d0, 4.d0, 0.d0, 1.3d0, 8.93d0, 8.68d0, 6.85d0, 7.05d0, 1.d0] * (-1)
797 a1 = [ 2.838d0, 1.322d0, 0.d0, 1.486d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-2
798 a2 = [ 1.023d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-4
799 a3 = [ 0.d0, 0.4247d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-6
800
801 else if(lmodel.eq.2) then ! Model 2: Meeus p.286
802 ! Mer Ven Ear Mars Jup Sat Ur Nep Pl
803 a0 = [ 0.42d0, 4.40d0, 0.d0, 1.52d0, 9.40d0, 8.88d0, 7.19d0, 6.87d0, 1.00d0] * (-1)
804 a1 = [ 3.80d0, 0.09d0, 0.d0, 1.6d0, 0.5d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-2
805 a2 = [-2.73d0, 2.39d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-4
806 a3 = [ 2.d0, -0.65d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-6
807
808 else if(lmodel.eq.3) then ! Model 3: Expl.Supl.tt.Astr.Almanac, 2nd Ed. (1992)
809 ! Mer Ven Ear Mars Jup Sat Ur Nep Pl
810 a0 = [ 0.36d0, 4.29d0, 0.d0, 1.52d0, 9.25d0, 8.88d0, 7.19d0, 6.87d0, 1.00d0] * (-1)
811 a1 = [ 3.80d0, 0.09d0, 0.d0, 1.6d0, 0.5d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-2
812 a2 = [-2.73d0, 2.39d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-4
813 a3 = [ 2.d0, -0.65d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-6
814
815 else if(lmodel.eq.4) then ! Model 4: Expl.Supl.tt.Astr.Almanac, 3rd Ed. (2012)
816 ! Mer Ven Ven2 Mars Jup Sat Ur Nep Pl
817 a0 = [ 0.60d0, 4.47d0, -0.98d0, 1.52d0, 9.40d0, 8.88d0, 7.19d0, 6.87d0, 1.01d0] * (-1)
818 a1 = [ 4.98d0, 1.03d0, -1.02d0, 1.6d0, 0.5d0, 4.4d0, 0.2d0, 0.d0, 0.d0] * 1.d-2
819 a2 = [-4.88d0, 0.57d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-4
820 a3 = [ 3.02d0, 0.13d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0, 0.d0] * 1.d-6
821
822 if(pl.eq.2 .and. pa.gt.163.6d0) lpl = 3 ! Venus 2
823 else
824 call quit_program_error('planet_magnitude(): model must be 1-4', 0)
825 end if
826
827 planet_magnitude = 5*log10(hc_dist*gc_dist) + a0(lpl) + a1(lpl)*pa + a2(lpl)*pa2 + a3(lpl)*pa3
828
829 end function planet_magnitude
830 !*********************************************************************************************************************************
831
832
833
834
835 !*********************************************************************************************************************************
836 !> \brief Calculate Saturn's magnitude
837 !!
838 !! \param t Dynamical time in Julian centuries after J2000.0
839 !!
840 !! \param gl Geocentric longitude (rad)
841 !! \param gb Geocentric latitude (rad)
842 !! \param d Geocentric distance (AU)
843 !!
844 !! \param l Heliocentric longitude (rad)
845 !! \param b Heliocentric latitude (rad)
846 !! \param r Heliocentric distance (AU)
847 !!
848 !! \param magmdl Model to use: 1: Müller (1893), 2: Meeus p.286; optional, default: 2.
849 !!
850 !! \see
851 !! - Müller, POPot 8, 193 (1893) - https://ui.adsabs.harvard.edu/abs/1893POPot...8..193M/ - p.341/149
852 !! - Meeus, Astronomical Algorithms, 1998, Ch.41
853
854 function satmagn(t, gl,gb,d, l,b,r, magmdl)
855 use sufr_kinds, only: double
856 use sufr_constants, only: d2r,r2d
857 use sufr_angles, only: rev2
858
859 implicit none
860 real(double), intent(in) :: t, gl,gb,d, l,b,r
861 integer, intent(in), optional :: magmdl
862 real(double) :: satmagn, in,om,bbb,n,ll,bb,u1,u2,du
863 integer :: lmagmdl
864
865 lmagmdl = 2
866 if(present(magmdl)) lmagmdl = magmdl
867
868 in = 0.490005d0 - 2.2686d-4*t + 7.d-8*t**2 ! Radians
869 om = 2.9584809d0 + 0.0243418d0*t + 7.19d-6*t**2 ! Radians
870
871 bbb = asin( sin(in)*cos(gb)*sin(gl-om) - cos(in)*sin(gb) )
872 n = (113.6655d0 + 0.8771d0*t)*d2r
873 ll = l - 0.01759d0*d2r/r
874 bb = b - 7.64d-4*d2r * cos(l-n)/r
875
876 u1 = atan2( sin(in)*sin(bb) + cos(in)*cos(bb)*sin(ll-om) , cos(bb)*cos(ll-om) )
877 u2 = atan2( sin(in)*sin(gb) + cos(in)*cos(gb)*sin(gl-om) , cos(gb)*cos(gl-om) )
878 du = abs(rev2(u1-u2)) ! rev2() is needed because u1 and u2 jump from pi to -pi at different moments
879
880 if(lmagmdl.eq.1) then
881 ! Model 1 - Müller (1983):
882 satmagn = -8.68d0 + 5.d0*log10(d*r) + 0.044d0*abs(du)*r2d - 2.60d0*sin(abs(bbb)) + 1.25d0*(sin(bbb))**2
883 else
884 ! Method 2 - Meeus p.286: - constant difference of 0.2m...(?)
885 satmagn = -8.88d0 + 5.d0*log10(d*r) + 0.044d0*abs(du)*r2d - 2.60d0*sin(abs(bbb)) + 1.25d0*(sin(bbb))**2
886 end if
887
888 end function satmagn
889 !*********************************************************************************************************************************
890
891
892 !*********************************************************************************************************************************
893 !> \brief Calculate a correction to Saturn's magnitude due to its rings.
894 !!
895 !! \param tjc Dynamical time in Julian centuries after J2000.0.
896 !!
897 !! \param glon Geocentric longitude (rad).
898 !! \param glat Geocentric latitude (rad).
899 !!
900 !! \see Müller, POPot 8, 193 (1893) - https://ui.adsabs.harvard.edu/abs/1893POPot...8..193M/ - p.341/149
901
902 function dsatmagn(tjc, glon,glat)
903 use sufr_kinds, only: double
904 use sufr_angles, only: rev2
905
906 implicit none
907 real(double), intent(in) :: tjc, glon,glat
908 real(double) :: dsatmagn, inc,omg, bbb ! , n,ll,bb, u1,u2,du, bbb
909
910 inc = 0.490005d0 - 2.2686d-4*tjc + 7.d-8*tjc**2 ! Radians
911 omg = 2.9584809d0 + 0.0243418d0*tjc + 7.19d-6*tjc**2 ! Radians
912
913 bbb = asin( sin(inc)*cos(glat)*sin(glon-omg) - cos(inc)*sin(glat) )
914
915 dsatmagn = -2.60d0*sin(abs(bbb)) + 1.25d0*(sin(bbb))**2
916
917 end function dsatmagn
918 !*********************************************************************************************************************************
919
920
921 !*********************************************************************************************************************************
922 !> \brief Calculate the umbra and penumbra geocentric radii of the Earth's shadow
923 !!
924 !! \param dm0 Distance of the Moon (AU)
925 !! \param ds0 Distance of the Sun (AU)
926 !!
927 !! \param r1 Umbra radius at distance of the Moon (rad) (output)
928 !! \param r2 Penumbra radius at distance of the Moon (rad) (output)
929 !!
930 !! \see Expl.sup. to the Astronomical Almanac, p.428
931
932 subroutine earthshadow(dm0,ds0, r1,r2)
933 use sufr_kinds, only: double
934 use sufr_constants, only: au, earthr,planr
935
936 implicit none
937 real(double), intent(in) :: dm0,ds0
938 real(double), intent(out) :: r1,r2
939 real(double) :: rs,re,ds,dm,p1,ps,ss
940
941 dm = dm0*au ! AU -> cm
942 ds = ds0*au ! AU -> cm
943
944 re = earthr
945 rs = planr(3)
946
947 p1 = asin(re/dm) * 0.998340d0 ! Moon parallax, compensated for Earth's atmosphere
948 ps = asin(re/ds) ! Solar parallax
949 ss = asin(rs/ds) ! Solar semi-diameter
950
951 r1 = 1.02d0 * (p1 + ps - ss)
952 r2 = 1.02d0 * (p1 + ps + ss)
953
954 end subroutine earthshadow
955 !*********************************************************************************************************************************
956
957
958 !*********************************************************************************************************************************
959 !> \brief Compute physical data for Jupiter
960 !!
961 !! \param jd Julian day for computation
962 !!
963 !! \param de Jovocentric latitude of the Earth = inclination of planet axis to plane of the sky as seen from Earth (output)
964 !! \param ds Jovocentric latitude of the Sun = inclination of planet axis to plane of the sky as seen from Earth (output)
965 !! \param omg1 Longitude of central meridian of System I (equator+-10deg) (output)
966 !! \param omg2 Longitude of central meridian of System II (output)
967 !!
968 !! \param dphase Phase correction (output)
969 !! \param pa Position angle of Jupiter's north pole (from N to E) (output)
970 !! \param in Inclination of Jupiter's rotation axis to the orbital plane (output)
971 !! \param om Longitude of node of Jupiter's equator on ecliptic (output)
972 !!
973 !! \see Meeus, Astronomical Algorithms, 1998, Ch.43, p.293-295
974
975 subroutine jupiterphys(jd, de,ds, omg1,omg2, dphase, pa, in,om)
976 use sufr_kinds, only: double
977 use sufr_constants, only: d2r, pi, jd1900
978 use sufr_angles, only: rev
979
980 use thesky_local, only: deltat
981 use thesky_planetdata, only: planpos
984
985 implicit none
986 real(double), intent(in) :: jd
987 real(double), intent(out) :: de,ds, omg1,omg2, dphase, pa, in,om
988 real(double) :: jde,d,t1,t2,t3
989 real(double) :: alp0,del0,w1,w2,x,y,z,eps0,eps,dpsi,l,b,alps,dels, r,r0,l0
990 real(double) :: u,v,alp,del,dze,delta,l1,b1
991
992 ! To get data for the exact (rounded off) moment of Meeus' example:
993 !d = 15690.00068d0
994 !jde = d + 2433282.5d0
995 !jd = jde - deltat/86400.d0
996
997 call planet_position(jd,5)
998 call planetelements(jd)
999
1000 deltat = calc_deltat(jd)
1001 jde = jd + deltat/86400.d0
1002
1003 ! Meeus, step 1, p.293:
1004 d = jde - 2433282.5d0 ! Since 1950 (!)
1005 t1 = d/36525.d0 ! Time since 1950 (!) in Julian centuries
1006
1007 t2 = (jde - jd1900)/36525.d0 ! Time since J1900.0 (!), in Julian centuries
1008 in = (3.120262d0 + 0.0006d0 * t2) * d2r ! Inclination of Jupiter's rotation axis to the orbital plane (Meeus, p.311)
1009
1010 t3 = jde - 2443000.5d0 - planpos(7) ! Days since 1976-08-10 (Meeus, p.305), minus light-travel time
1011 om = (316.518203d0 - 2.08362d-6*t3) * d2r ! Long. of node of Jupiter's equator on ecliptic, psi in Meeus, p.306
1012 om = om + (1.3966626d0*t1 + 3.088d-4*t1**2)*d2r ! Correct for precession (since B1950, but J1950 will do), Meeus p.311
1013
1014 alp0 = rev((268.d0 + 0.1061d0*t1)*d2r) ! Right ascension of Jupiter's north pole
1015 del0 = rev((64.5d0 - 0.0164d0*t1)*d2r) ! Declination of Jupiter's north pole
1016
1017 ! Meeus, step 2, p.294:
1018 w1 = rev((17.710d0 + 877.90003539d0*d)*d2r) ! Longitude system I
1019 w2 = rev((16.838d0 + 870.27003539d0*d)*d2r) ! Longitude system II
1020
1021 x = planpos(61) ! Apparent geocentric x
1022 y = planpos(62) ! Apparent geocentric y
1023 z = planpos(63) ! Apparent geocentric z
1024
1025 l = planpos(33) ! Heliocentric apparent l
1026 b = planpos(34) ! Heliocentric apparent b
1027 r = planpos(35) ! Heliocentric apparent r
1028 delta = planpos(4) ! Apparent geocentric distance
1029
1030 l0 = rev(planpos(41)+pi) ! Heliocentric, true coordinates of the Earth
1031 r0 = planpos(43)
1032
1033 ! Meeus, step 8, p.294:
1034 eps0 = planpos(50) ! Mean obliquity of the ecliptic; without nutation
1035 eps = planpos(48) ! True obliquity of the ecliptic; corrected for nutation
1036 dpsi = planpos(47) ! Nutation in longitude
1037
1038 ! Meeus, step 9, p.294:
1039 call ecl_2_eq(l,b,eps, alps,dels)
1040
1041 ! Meeus, step 10, p.294:
1042 ds = -sin(del0)*sin(dels) - cos(del0)*cos(dels)*cos(alp0-alps) ! Planetocentric declination of the Sun
1043
1044 ! Meeus, step 11, p.294:
1045 u = y*cos(eps0) - z*sin(eps0)
1046 v = y*sin(eps0) + z*cos(eps0)
1047
1048 alp = rev(atan2(u,x))
1049 del = atan2(v,sqrt(x**2+u**2))
1050 dze = rev( atan2( sin(del0)*cos(del)*cos(alp0-alp) - sin(del)*cos(del0) , cos(del)*sin(alp0-alp) ) )
1051
1052 ! Meeus, step 12, p.294:
1053 de = -sin(del0)*sin(del) - cos(del0)*cos(del)*cos(alp0-alp) ! Planetocentric declination of the Earth
1054
1055 ! Meeus, step 14, p.295 - phase correction - same sign as sin(l-l0):
1056 dphase = abs( (2*r*delta + r0**2 - r**2 - delta**2)/(4*r*delta) ) * sin(l-l0)/abs(sin(l-l0))
1057
1058 ! Meeus, step 13, p.295:
1059 omg1 = rev(w1 - dze - 5.07033d0*d2r * delta + dphase) ! Longitude of central meridian of System I (equator +- 10deg)
1060 omg2 = rev(w2 - dze - 5.02626d0*d2r * delta + dphase) ! Longitude of central meridian of System II (higher lat)
1061
1062 ! Meeus, steps 16, 17, p.295:
1063 call eq_2_ecl(alp0,del0,eps0, l1,b1) ! Has to be eps0
1064 call ecl_2_eq(l1+dpsi,b1,eps, alp0,del0) ! Has to be eps
1065
1066 ! Meeus, step 18, p.295:
1067 pa = atan2( cos(del0)*sin(alp0-alp) , sin(del0)*cos(del) - cos(del0)*sin(del)*cos(alp0-alp) ) ! PA of Jupiter's north pole
1068
1069 end subroutine jupiterphys
1070 !*********************************************************************************************************************************
1071
1072
1073 !*********************************************************************************************************************************
1074 !> \brief Compute physical data for Saturn
1075 !!
1076 !! \param jd Julian day of computation
1077 !!
1078 !! \param be Saturnicentric latitude of the Earth = inclination of planet axis to plane of the sky as seen from Earth (output)
1079 !! \param bs Saturnicentric latitude of the Sun = inclination of planet axis to plane of the sky as seen from the Sun (output)
1080 !!
1081 !! \param pa Position angle of Saturn's north pole (from N to E) (output)
1082 !! \param in Inclination of Saturn's rotation axis to the orbital plane (output)
1083 !! \param om Longitude of node of Saturn's equator on ecliptic (output)
1084 !!
1085 !! \param ar Projected major axis of ring (NOT semi-!!!) (output)
1086 !! \param br Projected minor axes of ring (NOT semi-!!!), = ar * sin(be) (output)
1087 !! \param du Difference between saturnicentric longitudes of Sun and Earth (Sun - Earth, NOT abs!!!) (output)
1088 !!
1089 !! \param pa_s Position angle of Saturn's north pole (from N to E), as seen from the SUN (output)
1090 !! \param ar_s Projected major axis of ring (NOT semi-!!!), as seen from the SUN (output)
1091 !! \param br_s Projected minor axes of ring (NOT semi-!!!), as seen from the SUN, = ar_s * sin(bs) (output)
1092 !!
1093 !!
1094 !! \see Meeus, Astronomical Algorithms, 1998, Ch.45 (The ring of Saturn)
1095 !!
1096 !! \todo
1097 !! - check true/apparent coordinates (see CHECK)
1098
1099 subroutine saturnphys(jd, be,bs, pa, in,om, ar,br, du, pa_s, ar_s,br_s)
1100 use sufr_kinds, only: double
1101 use sufr_constants, only: pi, as2r,d2r
1102 use sufr_angles, only: rev
1103
1106
1107 implicit none
1108 real(double), intent(in) :: jd
1109 real(double), intent(out) :: be,bs, pa, in,om, ar,br, du, pa_s, ar_s,br_s
1110 integer :: loop
1111 real(double) :: t, lam,bet, d,l,b,r, l0,b0,r0, x,y,z, dpsi,eps, ascn
1112 real(double) :: l2,b2, u1,u2, lam0,bet0, dl,db, ra,dec, ra0,dec0, sinbe
1113
1114 call planet_position(jd,6) ! Position of Saturn
1115
1116 ! Meeus, step 3:
1117 l = planpos(33) ! Heliocentric longitude of Saturn
1118 b = planpos(34) ! Heliocentric latitude of Saturn
1119 r = planpos(35) ! Heliocentric distance of Saturn
1120
1121 t = planpos(46) ! App. dyn. time in Julian Centuries since 2000.0
1122
1123 ! Meeus, step 10:
1124 dpsi = planpos(47) ! Nutation in longitude
1125 eps = planpos(48) ! True obliquity of the ecliptic; corrected for nutation
1126
1127 do loop=1,2 ! From Sun, Earth:
1128
1129 ! Meeus, step 2:
1130 if(loop.eq.1) then ! Seen from the Sun:
1131 lam = l
1132 bet = b
1133 d = r
1134 else ! Seen from the Earth:
1135 l0 = rev(planpos(41)+pi) ! Geocentric, true L,B,R for the Earth, in FK5 - CHECK - need apparent?
1136 b0 = -planpos(42)
1137 r0 = planpos(43)
1138
1139 ! Meeus, Ch. 45, step 5:
1140 call hc_spher_2_gc_rect(l,b,r, l0,b0,r0, x,y,z)
1141 call rect_2_spher(x,y,z, lam,bet,d)
1142 end if
1143
1144 ! Meeus, step 1:
1145 in = (28.075216d0 - 0.012998d0*t + 4.d-6 * t**2) * d2r ! Inclination of rotation axis, Eq. 45.1a
1146 om = (169.508470d0 + 1.394681d0*t + 4.12d-4 * t**2) * d2r ! Longitude of asc. node of equator - Omega, Eq. 45.1b
1147
1148 ! Meeus, step 6:
1149 sinbe = sin(in)*cos(bet)*sin(lam-om) - cos(in)*sin(bet)
1150 ar = 375.35d0*as2r/d ! Major axis of outer edge of outer ring
1151 br = ar*sinbe ! Minor axis of outer edge of outer ring
1152 be = asin(sinbe) ! B
1153
1154 if(loop.eq.2) then ! Seen from Earth:
1155 ! Meeus, step 7:
1156 call planetelements(jd)
1157 ascn = plelems(6,5) ! Longitude of Saturn's ascending node
1158 l2 = l - 0.01759d0*d2r / r ! Correct for the Sun's aberration on Saturn - lon
1159 b2 = b - 7.64d-4*d2r * cos(l-ascn) / r ! Correct for the Sun's aberration on Saturn - lat
1160
1161 ! Meeus, step 8:
1162 bs = asin(sin(in)*cos(b2)*sin(l2-om) - cos(in)*sin(b2)) ! B'
1163
1164 ! Meeus, step 9:
1165 u1 = atan2(sin(in)*sin(b2) + cos(in)*cos(b2)*sin(l2-om), cos(b2)*cos(l2-om)) ! Saturnicentric longitude of the Sun
1166 u2 = atan2(sin(in)*sin(bet) + cos(in)*cos(bet)*sin(lam-om), cos(bet)*cos(lam-om)) ! Saturnicentric longitude of the Earth
1167 du = u1-u2
1168 end if
1169
1170 ! Meeus, step 11:
1171 lam0 = om - pi/2.d0 ! Ecliptical longitude of Saturn's north pole
1172 bet0 = pi/2.d0 - in ! Ecliptical latitude of Saturn's north pole
1173
1174 if(loop.eq.2) then ! Seen from the Earth:
1175 ! Meeus, step 12 - correct for the aberration of Saturn:
1176 dl = 0.005693d0*d2r * cos(l0-lam)/cos(bet)
1177 db = 0.005693d0*d2r * sin(l0-lam)*sin(bet)
1178 lam = lam + dl
1179 bet = bet + db
1180
1181 ! Meeus, step 13 - correct for nutation in longitude:
1182 lam0 = lam0 + dpsi
1183 lam = lam + dpsi
1184 end if
1185
1186 ! Meeus, step 14:
1187 call ecl_2_eq(lam,bet,eps, ra,dec)
1188 call ecl_2_eq(lam0,bet0,eps, ra0,dec0)
1189
1190 ! Meeus, step 15:
1191 pa = atan2(cos(dec0)*sin(ra0-ra),sin(dec0)*cos(dec) - cos(dec0)*sin(dec)*cos(ra0-ra)) ! Position angle
1192
1193 if(loop.eq.1) then ! Save data for the Sun's pov:
1194 pa_s = pa
1195 ar_s = ar
1196 br_s = br
1197 end if
1198 end do ! loop=1,2 - from Sun, Earth
1199
1200 end subroutine saturnphys
1201 !*********************************************************************************************************************************
1202
1203
1204
1205
1206
1207 !*********************************************************************************************************************************
1208 !> \brief Compute low-accuracy planet positions. Sun and Moon have dedicated routines, for planets abridged VSOP87 is used.
1209 !!
1210 !! \param jd Julian Day of computation
1211 !! \param pl Planet number: 0: Moon, 1-2: Mer-Ven, 3: Sun, 4-9: Mar-Plu
1212 !!
1213 !! \param calc Calculate 1: l,b,r,diam, 2: + ra,dec, 3: + gmst,agst, 4: + az,alt, 5: + elon, mag, k, pa, parang, hp, GC LBR,
1214 !! 6: + topocentric positions (Sun and Moon only, optional - default = 5)
1215 !! \param nt Number of terms to use for the calculation (has an effect for Moon only; nt<=60);
1216 !! a smaller nt gives faster, but less accurate results (optional, default=60)
1217 !!
1218 !! \param lat Latitude of the observer (rad, optional)
1219 !! \param lon Longitude of the observer (rad, optional)
1220
1221 subroutine planet_position_la(jd, pl, calc,nt, lat,lon)
1222 use sufr_kinds, only: double
1223 use sufr_angles, only: rev,rev2
1224
1225 use thesky_planetdata, only: planpos
1226 use thesky_moon, only: moonpos_la
1227 use thesky_sun, only: sunpos_la
1229 use thesky_local, only: lat0,lon0
1230
1231 implicit none
1232 real(double), intent(in) :: jd
1233 integer, intent(in) :: pl
1234 integer, intent(in), optional :: calc,nt
1235 real(double), intent(in), optional :: lat,lon
1236
1237 integer :: lcalc,lnt
1238 real(double) :: dh_ref, parAng, dRA_ref
1239
1240 lcalc = 5
1241 if(present(calc)) lcalc = calc
1242
1243 if(present(lat)) lat0 = lat
1244 if(present(lon)) lon0 = lon
1245
1246
1247 planpos = 0.d0
1248 planpos(39) = dble(pl)
1249
1250 select case(pl)
1251 case(0) ! Moon
1252
1253 lnt = 60
1254 if(present(nt)) lnt = nt
1255 call moonpos_la(jd, lcalc,lnt)
1256
1257 case(3) ! Sun
1258
1259 call sunpos_la(jd, lcalc)
1260
1261 case default ! Planets
1262
1263 call planet_position(jd,pl, lbaccur=1.d-6,raccur=1.d-4, ltime=.false.) ! Truncated VSOP87, no light-time correction
1264
1265 end select
1266
1267
1268 ! Fill topocentric planpos elements (21-32) with geocentric values (1-12) for Moon and Sun:
1269 if(pl.eq.0 .or. pl.eq.3) then
1270 planpos(21:32) = planpos(1:12)
1271 planpos(31) = planpos(30) ! Altitude, "corrected" for refraction - except that it is not corrected in the cheap version
1272
1273
1274 ! Simple correction for horizontal parallax, Meeus p.281, assuming rho=1:
1275 if(lcalc.eq.5) planpos(30) = planpos(10) - asin(sin(planpos(17))*cos(planpos(10))) ! alt' = alt - Hp*cos(alt)
1276
1277 ! More precise conversion to topocentric coordinates:
1278 if(lcalc.ge.6) then ! This takes 15% more CPU time for the Moon (full accuracy; nt=60) and 110% more for the Sun!
1279 call geoc2topoc_ecl(planpos(1),planpos(2), planpos(4),planpos(12)/2.d0, planpos(48),planpos(44), & ! l,b, del,r, eps,lst
1280 planpos(21),planpos(22),planpos(32)) ! L,B, r
1281 planpos(32) = planpos(32)*2.d0 ! Apparent diameter
1282 planpos(24) = planpos(4)*planpos(12)/planpos(32) ! Distance
1283
1284 call geoc2topoc_eq(planpos(5),planpos(6), planpos(4),planpos(8), planpos(25),planpos(26)) ! ra,dec, del,hh, topo ra,dec
1285
1286 call eq2horiz(planpos(25),planpos(26),planpos(45), planpos(28),planpos(29),planpos(30)) ! topo ra,dec,agst, -> hh,az,alt
1287 end if
1288
1289 if(lcalc.ge.5) then
1290 dh_ref = refract(planpos(30)) ! Corrected for atmospheric refraction in altitude
1291 planpos(31) = planpos(30) + dh_ref ! Topocentric altitude, corrected for atmospheric refraction
1292
1293 parang = atan2(sin(planpos(28)),tan(lat0)*cos(planpos(26)) - sin(planpos(26))*cos(planpos(28))) ! Topoc. paral. angle
1294 dra_ref = dh_ref / cos(planpos(26)) * sin(parang)
1295 planpos(25) = rev(planpos(25) + dra_ref) ! Correct topocentric RA for refraction
1296 planpos(28) = rev(planpos(28) - dra_ref) ! Correct topocentric HA for refraction
1297 planpos(26) = planpos(26) + dh_ref*cos(parang)
1298 end if
1299 end if
1300
1301 end subroutine planet_position_la
1302 !*********************************************************************************************************************************
1303
1304 !*********************************************************************************************************************************
1305 !> \brief Correction for comet magnitude due to forward and backscattering
1306 !!
1307 !! \param phaseAng Phase angle of the comet (rad)
1308 !!
1309 !! \retval Magnitude correction, to be added to the magnitude (usually negative)
1310 !!
1311 !! \see https://ui.adsabs.harvard.edu/abs/2007ICQ....29..119M
1312 !!
1313
1315 use sufr_kinds, only: double
1316 use sufr_constants, only: pi ! , r2d
1317 use sufr_angles, only: rev2
1318
1319 implicit none
1320 real(double), intent(in) :: phaseang
1321 real(double) :: comet_scatter_magnitude_correction, theta,d90,gf,gb,k,phi
1322
1323 theta = rev2(pi - phaseang) ! Scattering angle = 180° - phase angle
1324
1325 d90 = 0.3d0 ! Dust-to-gas light ratio in the coma at theta=90°. 1 for a "normal" comet, ~10 for a "dusty" one. 3.3 ~ half way. 0.3 taken from http://astro.vanbuitenen.nl/comet/2020F8
1326 gf = 0.9d0 ! Forward-scattering asymmetry factor
1327 gb = -0.6d0 ! Backscattering asymmetry factor
1328 k = 0.95d0 ! Partition coefficient between forward and backscattering
1329
1330 phi = d90/(1+d90) * (k*((1+gf**2)/(1+gf**2 - 2*gf*cos(theta)))**1.5d0 + &
1331 k*((1+gb**2)/(1+gb**2 - 2*gb*cos(theta)))**1.5d0 + 1.d0/d90) ! Eq.1
1332
1333 comet_scatter_magnitude_correction = -2.5d0*log10(phi)
1334
1335 ! write(*,'(99F9.3)') phaseAng*r2d,theta*r2d,d90,gf,gb,k,phi,comet_scatter_magnitude_correction
1337 !*********************************************************************************************************************************
1338
1339
1340end module thesky_planets
1341!***********************************************************************************************************************************
Procedures for asteroids.
Definition asteroids.f90:22
subroutine asteroid_lbr(t, as, l, b, r)
Calculate heliocentric ecliptical coordinates l,b,r for asteroid as at time t (Julian Millennia after...
real(double) function asteroid_magn(as, delta, r, pa)
Calculate asteroid magnitude for asteroid as.
Data to compute comet positions.
Definition modules.f90:220
real(double), dimension(ncometsmax, 9) cometelems
Definition modules.f90:228
integer, dimension(ncometsmax) cometdiedatp
Definition modules.f90:227
Procedures for comets.
Definition comets.f90:22
subroutine cometgc(t, t0, comid, r, l, b, d)
Calc geocentric lbr coordinates for comet com at t1 (in Julian millennia Dynamical Time)
Definition comets.f90:265
Procedures for coordinates.
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,...
subroutine eq_2_ecl(ra, dec, eps, l, b)
Convert (geocentric) spherical equatorial coordinates RA, Dec (and eps) to spherical ecliptical coord...
subroutine rect_2_spher(x, y, z, l, b, r)
Convert rectangular coordinates x,y,z to spherical coordinates l,b,r.
subroutine ecl_2_eq(l, b, eps, ra, dec)
Convert (geocentric) spherical ecliptical coordinates l,b (and eps) to spherical equatorial coordinat...
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 aberration_ecl(t, l0, l, b)
Correct ecliptic longitude and latitiude for annual aberration.
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.
Date and time procedures.
Definition date_time.f90:48
real(double) function calc_gmst(jd, deltat)
Calculate Greenwich Mean Sidereal Time for any instant, in radians.
real(double) function calc_deltat(jd, force_recompute)
Compute DeltaT for a given JD.
Local parameters for libTheSky: location, date, time.
Definition modules.f90:44
real(double) height
Definition modules.f90:51
real(double) lon0
Definition modules.f90:51
real(double) deltat
Definition modules.f90:51
real(double) lat0
Definition modules.f90:51
Procedures for the Moon.
subroutine elp82b_lbr(tjj, ll, bb, rr)
Calculate the apparent geocentric ecliptical position of the Moon, using the Lunar Solution ELP 2000-...
real(double) function moonmagn(pa, delta)
Calculate the magnitude of the Moon.
subroutine elp_mpp02_lbr(jd, mode, lon, lat, rad)
Compute the spherical lunar coordinates using the ELP2000/MPP02 lunar theory in the dynamical mean ec...
subroutine moonpos_la(jd, calc, nt)
Quick, lower-accuracy lunar coordinates; ~600x faster than ELP.
Procedures for nutation.
Definition nutation.f90:24
subroutine nutation(t, dpsi, eps0, deps)
Calculate nutation - cheap routine from Meeus - as well as the mean obliquity of the ecliptic.
Definition nutation.f90:45
subroutine nutation2000(jd, dpsi_tot, deps_tot, eps0)
Compute nutation using the IAU 2000 model. Add the mean obliquity of the ecliptic by Laskar (1986).
Definition nutation.f90:101
Planet data, needed to compute planet positions.
Definition modules.f90:61
integer, dimension(43, 2) plub
Definition modules.f90:71
real(double), dimension(8, 6) plelems
Definition modules.f90:75
real(double), dimension(nplanpos) planpos
Planpos[] is an array with many different types of coordinates and related variables:
Definition modules.f90:153
real(double), dimension(8, 6) plelems2000
Definition modules.f90:75
integer, dimension(43, 2) plur
Definition modules.f90:71
real(double), dimension(2, 8, 6, 0:3) plelemdata
Definition modules.f90:75
integer, dimension(43, 3) pluc
Definition modules.f90:71
integer, dimension(43, 2) plul
Definition modules.f90:71
Procedures for planets.
Definition planets.f90:22
subroutine earthshadow(dm0, ds0, r1, r2)
Calculate the umbra and penumbra geocentric radii of the Earth's shadow.
Definition planets.f90:933
subroutine planet_position(jd, pl, lat, lon, hgt, lbaccur, raccur, ltime, aber, to_fk5, assume_jde, lunar_theory, nutat, magmdl, verbosity)
Compute the position, distance, etc of a planet.
Definition planets.f90:62
real(double) function satmagn(t, gl, gb, d, l, b, r, magmdl)
Calculate Saturn's magnitude.
Definition planets.f90:855
subroutine planet_position_la(jd, pl, calc, nt, lat, lon)
Compute low-accuracy planet positions. Sun and Moon have dedicated routines, for planets abridged VSO...
Definition planets.f90:1222
subroutine planetelements(jd)
Calculate orbital elements for the planets.
Definition planets.f90:700
subroutine jupiterphys(jd, de, ds, omg1, omg2, dphase, pa, in, om)
Compute physical data for Jupiter.
Definition planets.f90:976
real(double) function dsatmagn(tjc, glon, glat)
Calculate a correction to Saturn's magnitude due to its rings.
Definition planets.f90:903
subroutine saturnphys(jd, be, bs, pa, in, om, ar, br, du, pa_s, ar_s, br_s)
Compute physical data for Saturn.
Definition planets.f90:1100
real(double) function comet_scatter_magnitude_correction(phaseang)
Correction for comet magnitude due to forward and backscattering.
Definition planets.f90:1315
real(double) function planet_magnitude(pl, hc_dist, gc_dist, phang, model)
Calculate planet magnitude.
Definition planets.f90:771
subroutine plutolbr(t, l, b, r)
Calculate Pluto's position l,b,r at time t.
Definition planets.f90:643
Low-accuracy procedures for the Sun.
Definition sun.f90:22
real(double) function sunmagn(dist)
Calculate Sun magnitude.
Definition sun.f90:282
subroutine sunpos_la(jd, calc, lat, lon)
Low-accuracy solar coordinates.
Definition sun.f90:57
Procedures for VSOP.
Definition vsop.f90:33
subroutine vsop87d_lbr(tm, pl, lon, lat, rad, lbaccur, raccur)
Calculate true heliocentric ecliptic coordinates l,b,r for planet pl and the mean ecliptic and equino...
Definition vsop.f90:58