libTheSky
Routines to compute sky positions of Sun, Moon, planets and more
All Namespaces Files Functions Variables Pages
riset.f90
Go to the documentation of this file.
1!> \file riset.f90 Compute rise, transit and set times, or beginning and end of twilight 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!***********************************************************************************************************************************
22!> \brief Procedures to determine rise, transit and set times
23
25 implicit none
26 save
27
28contains
29
30 !*********************************************************************************************************************************
31 !> \brief Rise, transit and set times routine for Sun, Moon, planets and asteroids.
32 !! This version recomputes positions (low accuracy first, then full accuracy) during the convergence process.
33 !! See riset_ipol() for a version using interpolation, to be used for planets only(!)
34 !!
35 !! \param jd Julian day number. The rise/transit and set data are computed for the (calendar) DAY of this JD (midnight-midnight),
36 !! unless for_night=.true. in which case this happens for the NIGHT starting at JD.
37 !! \param pl Planet/object number (-1 - 9: ES, Moon, Mer-Plu; >10000: asteroids)
38 !!
39 !! \param rt Rise time (hours) (output)
40 !! \param tt Transit time (hours) (output)
41 !! \param st Set time (hours) (output)
42 !!
43 !! \param rh Rising wind direction (rad) (output)
44 !! \param ta Transit altitude (rad) (output)
45 !! \param sh Setting wind direction (rad) (output)
46 !!
47 !! \param rsAlt Altitude to return rise/set data for (degrees; 0. is actual rise/set). rsAlt>90: compute transit only
48 !! \param for_night Compute rise/tranit/set data for the NIGHT starting at the indicated JD, rather than the calendar day.
49 !! Times between 0h and 12h are for the NEXT calendar day! (optional; default: false)
50 !!
51 !! \param ltime Passed to planet_position(). If .true., include light time, doubling the CPU time while gaining a bit of
52 !! accuracy (optional; default: false)
53 !! \param cWarn Warn upon convergence failure (optional; default: true)
54 !!
55 !! \param converge Number of iterations needed to converge (optional) (output)
56 !!
57 !! \note
58 !! - for rsAlt = 0.d0, rise and set times are computed
59 !! - for (rsAlt.ne.0), the routine calculates when alt=rsAlt is reached
60 !! - use rsAlt=-6,-12,-18 for the sun for twilight calculations (rsAlt is expressed in degrees).
61 !! - Returns times, rise/set azimuth and transit altitude
62 !! - Moon transit error ~0.15s (?)
63 !!
64 !! \todo
65 !! - This version sometimes finds the answer tmRad(i)>2pi, which is the answer for the next day...
66 !! - (only?) solution: use a solver on JD for az,alt
67 !!
68 !! \note Speed wrt riset_ipol(), transit only, when using low-accuracy approximations:
69 !! - Moon: 270x faster (used to be slowest by factor 3-5 wrt planets, factor 6.5 wrt Sun)
70 !! - Sun: 340x faster - full accuracy (VSOP): ~2-3x SLOWER!
71 !! - Mer: 20% slower
72 !! - Ven: 4% slower
73 !! - Mar: 10% slower
74 !! - Jup: 25% slower
75 !! - Nep: 30% slower
76 !!
77 !! \see
78 !! - riset_ipol()
79 !! - Meeus, Astronomical algorithms, Ch.15, but with geographic longitude east of Greenwich defined as > 0
80
81
82 subroutine riset(jd,pl, rt,tt,st, rh,ta,sh, rsAlt, for_night, ltime, cWarn, converge)
83 use sufr_kinds, only: double
84 use sufr_constants, only: pi2, d2r,am2r, r2h, enpname, earthr,au
85 use sufr_angles, only: rev, rev2
86 use sufr_date_and_time, only: cal2jd,jd2cal
87 use sufr_time2string, only: hm
88 use sufr_numerics, only: deq0, deq
89
90 use thesky_local, only: tz, lat0,lon0,deltat
92 use thesky_planetdata, only: planpos
94
95 implicit none
96 integer, intent(in) :: pl
97 real(double), intent(in) :: jd, rsAlt
98 real(double), intent(out) :: rt,tt,st,rh,ta,sh
99 logical, intent(in), optional :: for_night, ltime, cWarn
100 integer, intent(out), optional :: converge(3)
101
102 integer :: evi,iter,yr,mnt,tc,evMax, lconverge(3)
103 real(double) :: dy,day0, jd0,jd1,tmRad(3),tmhr(3), ra,dec
104 real(double) :: rsa,cosH0,h0,agst0,th0,dTmRad,accur, ha,alt,azAlt(3)
105 character :: event(3)*(13)
106 logical :: use_vsop, lfor_night, lltime, lcWarn
107
108
109 ! Take care of optional input variables:
110 lfor_night = .false. ! Compute rise/transit/set times for a given NIGHT, rather than a calendar day
111 if(present(for_night)) lfor_night = for_night
112
113 lltime = .false. ! Call planet_position() IGNORING light time by default (faster, lower accuracy)
114 if(present(ltime)) lltime = ltime
115
116 lcwarn = .true. ! Warn upon convergence failure
117 if(present(cwarn)) lcwarn = cwarn
118
119 ! Use the old interpolation routine for all but Moon and Sun:
120 if(pl.ne.0.and.pl.ne.3) then ! Not the Moon or the Sun, but a planet, comet or asteroid
121 call riset_ipol(jd,pl, rt,tt,st, rh,ta,sh, rsalt, for_night=lfor_night, ltime=lltime, cwarn=lcwarn)
122 return
123 end if
124
125
126 alt=0.d0; ha=0.d0; h0=0.d0; azalt=0.d0; tmrad=0.d0
127 tc = 1 ! 0: geocentric (unused), 1: topocentric; see also different rsa for the Moon below
128 event = ['Transit time ','Rise time ','Set time ']
129
130 if(abs(rsalt).gt.1.d-9) then
131 rsa = rsalt*d2r ! Use a user-specified altitude
132 else
133 rsa = -0.5667d0*d2r ! Standard rise/set altitude for planets
134 if(pl.eq.3) rsa = rsa - 16.d0*am2r ! Correct for radius of Sun, -16 arcminutes ! rsa = -0.8333d0*d2r - default for Sun
135 if(pl.eq.0) rsa = 0.125d0*d2r ! Approximate standard rise/set altitude for Moon, including parallax
136 end if
137
138 evmax = 3 ! 'Maximum' event to compute: transit only: evMax=1, +rise/set: evMax=3
139 if(rsalt.gt.90.d0) evmax = 1 ! Compute transit time and altitude only
140
141 call jd2cal(jd, yr,mnt,dy)
142 day0 = dble(floor(dy))-tz/24.d0 ! Midnight local time at the start of the desired DAY, needed for agst0
143 if(lfor_night) day0 = day0 + 0.5d0 ! Noon local time at the start of the desired NIGHT, needed for agst0
144 jd0 = cal2jd(yr,mnt,day0) ! Midnight local time at the start of the desired day, needed for agst0
145
146 call planet_position_la(jd0, pl, 3, 60) ! Compute low-accuracy positions - calc=2 computes ra,dec, calc=3 computes AGST
147
148 ra = planpos(5+tc*20)
149 dec = planpos(6+tc*20)
150 agst0 = planpos(45) ! AGST for midnight
151
152 if(pl.eq.0) then ! Moon
153 if(tc.eq.0) then ! Geocentric (unused)
154 rsa = asin(earthr/(planpos(4)*au))*0.7275d0-0.5667d0*d2r ! Exact altitude for Moon, ~0.1-0.2deg, for geocentric coord.
155 else ! tc.eq.1 - topocentric:
156 rsa = -0.8333d0*d2r ! For Moon, in combination with topocentric coordinates
157 end if
158 end if
159
160
161 if(evmax.eq.3) then ! Computing transit, rise and set
162 cosh0 = (sin(rsa)-sin(lat0)*sin(dec)) / (cos(lat0)*cos(dec)) ! Cosine of the hour angle of rise/set; Meeus, Eq.15.1
163 if(abs(cosh0).gt.1.d0) then ! Body never rises/sets
164 evmax = 1 ! Compute transit time and altitude only
165 else
166 h0 = rev(2*acos(cosh0))/2.d0 ! Hour angle of rise/set
167 end if
168 end if
169
170
171 tmrad(1) = rev(ra - lon0 - agst0) ! Transit time in radians; tmRad(1)=m0, tmRad(2)=m1, tmRad(3)=m2 in Meeus, but lon0 > 0 for E
172 if(evmax.eq.3) then
173 tmrad(2) = rev(tmrad(1) - h0) ! Rise time in radians
174 tmrad(3) = rev(tmrad(1) + h0) ! Set time in radians
175 end if
176
177
178 do evi=1,evmax ! Transit, rise, set
179 iter = 0
180 accur = 1.d-3 ! Accuracy (time in rad). Initially 1d-3 (~14s), later 1d-5 (~0.14s). Don't make this smaller than 1d-16
181 use_vsop = .false. ! Initially
182
183 dtmrad = huge(dtmrad)
184 do while(abs(dtmrad).ge.accur .or. .not.use_vsop)
185 th0 = agst0 + 1.002737909350795d0*tmrad(evi) ! Solar day in sidereal days in 2000; Expl.Suppl.tt.Astr.Almanac 3rdEd
186 jd1 = jd0 + tmrad(evi)/pi2 + deltat/86400.d0 ! ra, s -> days Eq.3.17 (removed '...37...' typo)
187
188 if(abs(dtmrad).le.accur) then
189 use_vsop = .true.
190 accur = 1.d-5 ! 1d-5 rad ~ 0.14s. Changing this to 1.d-4 (~1s) speeds the code yearly_moon_table code up by ~30%. Don't make this smaller than 1d-16.
191 end if
192
193 if(use_vsop) then
194 ! Accuracy of 1 min = 0.25 deg = 4e-3 rad. 1.d-6,1.d-2 saves ~43% CPU time for the Sun, <1 round-off errors/year
195 call planet_position(jd1,pl, lbaccur=1.d-6,raccur=1.d-2, ltime=lltime) ! Uses truncated VSOP87
196 else
197 call planet_position_la(jd1,pl, 2,60) ! Computes low-accuracy positions - calc=2 computes ra,dec - use 60 terms
198 end if
199
200 ra = planpos(5+tc*20) ! Right ascension
201 dec = planpos(6+tc*20) ! Declination
202
203 ha = rev2(th0 + lon0 - ra) ! Hour angle; Meeus p.103
204 ! ha = planpos(8+tc*20) ! Hour angle -/- DeltaT
205 alt = asin(sin(lat0)*sin(dec) + cos(lat0)*cos(dec)*cos(ha)) ! Altitude; Meeus, Eq.13.6
206
207 ! Correction to transit/rise/set times:
208 if(evi.eq.1) then ! Transit
209 dtmrad = -ha
210 else ! Rise/set
211 dtmrad = (alt-rsa)/(cos(dec)*cos(lat0)*sin(ha))
212 end if
213 tmrad(evi) = tmrad(evi) + dtmrad
214
215 iter = iter+1
216 if(iter.gt.30) exit ! do-while loop
217 end do ! do while(abs(dTmRad).ge.accur .or. .not.use_vsop)
218
219
220 if(iter.gt.30) then ! Convergence failed
221 if(lcwarn .and. (pl.ne.3.or.nint(rsalt).ne.-18)) &
222 write(0,'(A,F10.3,A)') ' * WARNING: riset(): Riset failed to '// &
223 'converge: '//trim(enpname(pl))//' '//trim(event(evi)),rsalt,'d *'
224
225 tmrad(evi) = 0.d0
226 azalt(evi) = 0.d0
227 else ! Result converged, store it
228 if(evi.eq.1) then
229 azalt(evi) = alt ! Transit altitude
230 else
231 azalt(evi) = atan2(sin(ha),(cos(ha)*sin(lat0)-tan(dec)*cos(lat0))) ! Rise,set hour angle -> azimuth
232 end if
233 end if
234
235 ! The Moon may rise/transit/set AFTER MIDNIGHT (t>24h; tmRad > 2pi):
236 if(pl.eq.0 .and. tmrad(evi).gt.pi2) then ! Moon; in case after m_i = m_i+1, m_f > 1
237 tmrad(evi) = 0.d0
238 azalt(evi) = 0.d0
239 end if
240
241 ! - If the time is < 0, the Moon's SET times (only?) may repeat the previous value on days where
242 ! the Moon does not set.
243 ! - Don't do this for planets, which may rise/transit/set twice a day!
244 ! - Don't do this for twilight (Sun with rsAlt<0), which can cross midnight in summer!
245 ! - This may be needed for planets/the Sun near the poles...
246 !
247 ! if(tmRad(evi).lt.0.d0 .and. deq0(rsAlt)) then ! Pre-2025
248 ! if((pl.eq.0.or.pl.eq.3) .and. tmRad(evi).lt.0.d0 .and. deq0(rsAlt)) then ! Sun and Moon only
249 ! if(pl.eq.0 .and. tmRad(evi).lt.0.d0 .and. deq0(rsAlt)) then ! Moon only
250 if(pl.eq.0 .and. tmrad(evi).lt.0.d0) then ! Moon only - all altitudes
251 tmrad(evi) = 0.d0
252 azalt(evi) = 0.d0
253 end if
254
255 lconverge(evi) = iter ! Number of iterations needed to converge
256
257 end do ! evi: 1-evMax to cycle through transit, rise, set
258
259
260 ! Store results:
261 tmhr = tmrad * r2h ! Times radians -> hours
262 if(lfor_night) then
263 tmhr = tmhr + 12 ! Times relative to noon -> relative to midnight
264 where(deq(tmhr,12.d0)) tmhr = 0.d0
265 end if
266
267 tt = tmhr(1) ! Transit time
268 rt = tmhr(2) ! Rise time
269 st = tmhr(3) ! Set time
270
271 ta = azalt(1) + refract(azalt(1)) ! Transit altitude + refraction
272 rh = azalt(2) ! Rise azimuth
273 sh = azalt(3) ! Set azimuth
274
275
276 ! Take care of optional output variables:
277 if(present(converge)) converge = lconverge ! Number of iterations needed to converge
278
279 end subroutine riset
280 !*********************************************************************************************************************************
281
282
283
284
285 !*********************************************************************************************************************************
286 !> \brief Old routine for rise, transit and set times for planets and asteroids - don't use this for Sun and Moon
287 !! Computes three sets of planet positions, and interpolates between them.
288 !!
289 !! \param jd Julian day number. The rise/transit and set data are computed for the (calendar) DAY of this JD (midnight-midnight),
290 !! unless for_night=.true. in which case this happens for the NIGHT starting at JD.
291 !! \param pl Planet/object number (-1 - 9: ES, Moon, Mer-Plu; >10000: asteroids)
292 !!
293 !! \param rt Rise time (hours) (output)
294 !! \param tt Transit time (hours) (output)
295 !! \param st Set time (hours) (output)
296 !!
297 !! \param rh Rising wind direction (rad) (output)
298 !! \param ta Transit altitude (rad) (output)
299 !! \param sh Setting wind direction (rad) (output)
300 !!
301 !! \param rsAlt Altitude to return rise/set data for (degrees; 0. is actual rise/set). rsAlt>90: compute transit only
302 !! \param for_night Compute rise/tranit/set data for the NIGHT starting at the indicated JD, rather than the calendar day.
303 !! Times between 0h and 12h are for the NEXT calendar day! (optional; default: false)
304 !!
305 !! \param ltime Passed to planet_position(). If .true., include light time, doubling the CPU time while gaining a bit of accur.
306 !! \param cWarn Warn upon convergence failure (optional; default: true)
307 !!
308 !! \note
309 !! - for rsAlt = 0.d0, rise and set times are computed
310 !! - for (rsAlt.ne.0), the routine calculates when alt=rsAlt is reached
311 !! - use rsAlt=-6,-12,-18 for the sun for twilight calculations (rsAlt is expressed in degrees).
312 !! - transit error for Moon <= ~10s (due to interpolation); typically ~4s?
313 !! - code gets confused around midnight - it sometimes returns a small negative tmdy(i), which is on the day before(!)
314 !!
315 !! \see
316 !! - Meeus, Astronomical algorithms, Ch.15, but with geographic longitude east of Greenwich defined as > 0
317
318 subroutine riset_ipol(jd,pl, rt,tt,st, rh,ta,sh, rsAlt, for_night, ltime, cWarn)
319 use sufr_kinds, only: double
320 use sufr_constants, only: pi,pi2, d2r,am2r, enpname, earthr,au
321 use sufr_system, only: warn
322 use sufr_angles, only: rev, rev2
323 use sufr_date_and_time, only: cal2jd,jd2cal
324 use sufr_numerics, only: deq0, deq
325
326 use thesky_local, only: tz, lat0,lon0,deltat
328 use thesky_planetdata, only: planpos
329 use thesky_coordinates, only: refract
330
331 implicit none
332 integer, intent(in) :: pl
333 real(double), intent(in) :: jd, rsAlt
334 real(double), intent(out) :: rt,tt,st, rh,ta,sh
335 logical, intent(in), optional :: for_night, ltime, cWarn
336
337 integer :: evi,iter,yr,mnt,tc,evMax, indic
338 real(double) :: dy,day0, jd0,jd1,jd2, tmdy(3),dTmdy(3), tmhr(3)
339 real(double) :: ra0,dec0,ra1,dec1,ra2,dec2,ra,dec, rsa,cosH0,h0,agst0,th0,n,dtm,accur, ha,alt,azAlt(3)
340 character :: event(3)*(13)
341 logical :: lfor_night, lltime, lcWarn
342 save :: indic
343
344
345 ! Take care of optional input variables:
346 lfor_night = .true.
347 if(present(for_night)) lfor_night = for_night
348
349 lltime = .true.
350 if(present(ltime)) lltime = ltime
351
352 lcwarn = .true.
353 if(present(cwarn)) lcwarn = cwarn
354
355 if((pl.eq.0.or.pl.eq.3) .and. indic.ne.12345) then
356 call warn("riset_ipol(): Don't use this routine for Sun or Moon - use riset() instead.", 0)
357 indic = 12345
358 end if
359
360 rt=0.d0; tt=0.d0; st=0.d0; rh=0.d0; ta=0.d0; sh=0.d0
361 alt = 0.d0; ha = 0.d0; dtmdy = 0.d0; h0 = 0.d0; dec = 0.d0
362 tc = 1 ! 0: geocentric (unused), 1: topocentricl; see also different rsa for the Moon below
363 event = ['Transit time ','Rise time ','Set time ']
364 accur = 1.d-6 ! Accuracy. 1d-6 ~ 0.1s. Don't make this smaller than 1d-16
365
366
367 rsa = -0.5667d0*d2r ! Standard rise/set altitude for planets
368 ! NOTE: Sun and Moon should not use this routine(!)
369 if(pl.eq.3) rsa = rsa - 16.0d0*am2r ! Compensate for radius of Sun, -16 arcminutes ! rsa = -0.8333d0*d2r - default for Sun
370 if(pl.eq.0) rsa = 0.125d0*d2r ! Approximate standard altitude for Moon
371
372 if(abs(rsalt).gt.1.d-9) rsa = rsalt*d2r
373
374 call jd2cal(jd,yr,mnt,dy)
375 day0 = dble(floor(dy))-tz/24.d0 ! Midnight local time of the desired day, needed for agst0
376 if(lfor_night) day0 = day0 + 0.5d0 ! Noon local time at the start of the desired NIGHT, needed for agst0
377 jd0 = cal2jd(yr,mnt,day0-1.d0) ! Midnight local time, day before
378 jd1 = cal2jd(yr,mnt,day0) ! Midnight local time, needed for agst0
379 jd2 = cal2jd(yr,mnt,day0+1.d0) ! Midnight local time, day after
380
381 call planet_position(jd0,pl, lbaccur=1.d-6,raccur=1.d-2, ltime=lltime)
382 ra0 = planpos(5+tc*20)
383 dec0 = planpos(6+tc*20)
384
385 call planet_position(jd1,pl, lbaccur=1.d-6,raccur=1.d-2, ltime=lltime)
386 ra1 = planpos(5+tc*20)
387 dec1 = planpos(6+tc*20)
388 agst0 = planpos(45) ! AGST for midnight
389
390 ! NOTE: Sun and Moon should not use this routine(!)
391 if(pl.eq.0) then
392 if(tc.eq.0) then ! Geocentric (unused)
393 rsa = asin(earthr/(planpos(4)*au))*0.7275d0-0.5667d0*d2r ! Exact altitude for Moon, ~0.1-0.2deg, for geocentric coord.
394 else ! tc.eq.1 - topocentric:
395 rsa = -0.8333d0*d2r ! For Moon, in combination with topocentric coordinates
396 end if
397 end if
398
399 call planet_position(jd2,pl, lbaccur=1.d-6,raccur=1.d-2, ltime=lltime)
400 ra2 = planpos(5+tc*20)
401 dec2 = planpos(6+tc*20)
402
403 if(abs(ra1-ra0).gt.2..or.abs(ra2-ra1).gt.2.) then ! ra flips 2pi <-> 0
404 if(ra0.ge.4.5.and.ra1.ge.4.5.and.ra2.le.2.) ra2 = ra2+2*pi
405 if(ra0.ge.4.5.and.ra1.le.2.and.ra2.le.2.) ra0 = ra0-2*pi
406 if(ra0.le.2.and.ra1.ge.4.5.and.ra2.ge.4.5) ra0 = ra0+2*pi
407 if(ra0.le.4.5.and.ra1.le.2.and.ra2.ge.4.5) ra2 = ra2-2*pi
408 ! write(6,'(A)') 'RA flipped'
409 end if
410
411 cosh0 = (sin(rsa)-sin(lat0)*sin(dec2)) / (cos(lat0)*cos(dec2)) ! Cosine of the hour angle of rise/set; Meeus, Eq.15.1
412
413 tmdy = 0.d0
414 azalt = 0.d0
415
416 evmax = 3 ! 'Maximum' event to compute: transit only: evMax=1, +rise/set: evMax=3
417 if(abs(cosh0).gt.1.d0) then ! Body never rises/sets
418 evmax = 1 ! Compute transit time and altitude only
419 else
420 h0 = rev(2*acos(cosh0))/2.d0 ! Hour angle of rise/set
421 end if
422
423
424 tmdy(1) = rev(ra2 - lon0 - agst0)/pi2 + dtmdy(1) ! Transit time in days; Meeus Eq. 15.2a, but lon0 > 0 for E
425 if(evmax.eq.3) then
426 tmdy(2) = rev(tmdy(1)*pi2 - h0)/pi2 + dtmdy(2) ! Rise time in days; Meeus Eq.15.2b
427 tmdy(3) = rev(tmdy(1)*pi2 + h0)/pi2 + dtmdy(3) ! Set time in days; Meeus Eq.15.2c
428 end if
429
430 evi = 1 ! i=1,evMax
431 do while(evi.le.evmax)
432 iter = 0
433
434 dtm = huge(dtm)
435 do while(abs(dtm).ge.accur)
436 th0 = agst0 + 6.300388092591991d0*tmdy(evi) ! Siderial day in radians; Meeus, p.103
437 n = tmdy(evi) + deltat/86400.d0
438 ra = rsipol( ra0, ra1, ra2, n) ! Interpolate right ascension
439 dec = rsipol(dec0,dec1,dec2, n) ! Interpolate declination
440
441 ha = rev2(th0 + lon0 - ra) ! Hour angle; Meeus p.103
442 alt = asin(sin(lat0)*sin(dec) + cos(lat0)*cos(dec)*cos(ha)) ! Altitude; Meeus, Eq.13.6
443
444 ! Correction to transit/rise/set times:
445 if(evi.eq.1) then ! Transit
446 dtm = -ha/pi2
447 else ! Rise/set
448 dtm = (alt-rsa)/(pi2*cos(dec)*cos(lat0)*sin(ha))
449 end if
450 tmdy(evi) = tmdy(evi) + dtm
451
452 iter = iter+1
453 if(iter.gt.30) exit ! do-while loop
454 end do ! do while(abs(dtm).ge.accur)
455
456
457 if(iter.gt.30) then ! Convergence failed
458 if(lcwarn .and. (pl.ne.3.or.nint(rsalt).ne.-18)) &
459 write(0,'(A,F10.3,A)') ' * WARNING: riset_ipol(): Riset failed '// &
460 'to converge: '//trim(enpname(min(pl,19)))//' '//trim(event(evi)),rsalt,'d'
461
462 tmdy(evi) = 0.d0
463 azalt(evi) = 0.d0
464 else ! Result converged, store it
465 if(evi.eq.1) then
466 azalt(evi) = alt ! Transit altitude
467 else
468 azalt(evi) = atan2(sin(ha),(cos(ha)*sin(lat0)-tan(dec)*cos(lat0))) ! Rise,set hour angle -> azimuth
469 end if
470 end if
471
472
473 ! NOTE: Sun and Moon should not use this routine(!)
474 if(pl.eq.0) then
475 if(tmdy(evi).lt.0.) then ! Sometimes, (at least) the Moon doesn't r,t,s twice, because m_i slightly > 1, and rev(m_i) << 1,
476 dtmdy(evi) = dtmdy(evi) + 1.d0 ! then initially m_f < 0, though if m_i = m_i + 1, then 0 < m_f < 1
477
478 ! Restart the do loop:
479 evi = 1
480 tmdy(1) = rev(ra2 - lon0 - agst0)/pi2 + dtmdy(1) ! Transit; tmdy(1)=m0, tmdy(2)=m1, tmdy(3)=m2 in Meeus, but lon0 > 0 for E
481 if(evmax.eq.3) then
482 tmdy(2) = rev(tmdy(1)*pi2 - h0)/pi2 + dtmdy(2) ! Rise
483 tmdy(3) = rev(tmdy(1)*pi2 + h0)/pi2 + dtmdy(3) ! Set
484 end if
485 cycle ! i
486 end if ! if(tmdy(evi).lt.0.)
487
488 if(tmdy(evi).gt.1.d0) then ! in case after m_i = m_i+1, m_f > 1
489 tmdy(evi) = 0.d0
490 azalt(evi) = 0.d0
491 end if ! if(tmdy(evi).gt.1.d0)
492 end if ! if(pl.eq.0)
493
494
495 ! NOTE: Sun and Moon should not use this routine(!)
496 ! - If the time is < 0, the Moon's SET times (only?) may repeat the previous value on days where
497 ! the Moon does not set.
498 ! - Don't do this for planets, which may rise/transit/set twice a day!
499 ! - Don't do this for twilight (Sun with rsAlt<0), which can cross midnight in summer!
500 ! - This may be needed for planets/the Sun near the poles...
501 !
502 ! if(tmdy(evi).lt.0.d0 .and. deq0(rsAlt)) then ! Pre-2025
503 ! if((pl.eq.0.or.pl.eq.3) .and. tmdy(evi).lt.0.d0 .and. deq0(rsAlt)) then
504 ! if(pl.eq.0 .and. tmdy(evi).lt.0.d0 .and. deq0(rsAlt)) then ! Moon only
505 if(pl.eq.0 .and. tmdy(evi).lt.0.d0) then ! Moon only - all altitudes
506 tmdy(evi) = 0.d0
507 azalt(evi) = 0.d0
508 end if
509
510 evi = evi+1 ! evi=1,evMax: cycle through transit, rise, set
511 end do ! do while(evi.le.evMax)
512
513
514 ! Store results:
515 tmhr = tmdy * 24 ! Times days -> hours
516 if(lfor_night) then
517 tmhr = tmhr + 12 ! Times relative to noon -> relative to midnight
518 where(deq(tmhr,12.d0)) tmhr = 0.d0
519 end if
520
521 rt = tmhr(2) ! Rise time
522 tt = tmhr(1) ! Transit time
523 st = tmhr(3) ! Set time
524
525 rh = azalt(2) ! Rise azimuth
526 ta = azalt(1) + refract(azalt(1)) ! Transit altitude + refraction
527 sh = azalt(3) ! Set azimuth
528
529 end subroutine riset_ipol
530 !*********************************************************************************************************************************
531
532
533
534
535 !*********************************************************************************************************************************
536 !> \brief Rise, transit and set times routine for an object with fixed ra & dec
537 !!
538 !! \param jd Julian day number
539 !! \param ra Right ascension (rad)
540 !! \param dec Declination (rad)
541 !!
542 !! \param rt Rise time (hours) (output)
543 !! \param tt Transit time (hours) (output)
544 !! \param st Set time (hours) (output)
545 !!
546 !! \param rh Rising wind direction (rad) (output)
547 !! \param ta Transit altitude (rad) (output)
548 !! \param sh Setting wind direction (rad) (output)
549 !!
550 !! \param rsAlt Altitude to return rise/set data for (degrees; 0. is actual rise/set). rsAlt>90: compute transit only
551 !! \param for_night Compute rise/tranit/set data for the NIGHT starting at the indicated JD, rather than the calendar day.
552 !! Times between 0h and 12h are for the NEXT calendar day! (optional; default: false)
553 !!
554 !! \param cWarn Warn upon convergence failure (optional; default: true)
555 !!
556 !! for rsAlt = 0., rise and set times are computed
557 !! for rsAlt.ne.0, the routine calculates when alt=rsAlt is reached
558 !! use rsAlt=-6,-12,-18 for the sun for twilight calculations
559 !! (rsAlt is expressed in degrees).
560 !!
561 !! \see
562 !! - Meeus, Astronomical algorithms, Ch.15, but with geographic longitude east of Greenwich defined as > 0
563
564 subroutine riset_ad(jd, ra,dec, rt,tt,st,rh,ta,sh, rsAlt, for_night, cWarn)
565 use sufr_kinds, only: double
566 use sufr_constants, only: pi2, d2r
567 use sufr_angles, only: rev, rev2
568 use sufr_date_and_time, only: cal2jd, jd2cal
569 use sufr_numerics, only: deq0, deq
570
572 use thesky_planetdata, only: planpos
573 use thesky_local, only: tz, lat0,lon0
574
575 implicit none
576 real(double), intent(in) :: jd, ra,dec, rsAlt
577 real(double), intent(out) :: rt,tt,st,rh,ta,sh
578 logical, intent(in), optional :: for_night, cWarn
579
580 integer :: evi,iter,yr,mnt,evMax
581 real(double) :: dy,day0,jd0, tmdy(3),tmhr(3), rsa,cosH0,h0,agst0,th0,dtm,accur, ha,alt,azAlt(3) !,n
582 character :: event(3)*(13)
583 logical :: lfor_night, lcWarn
584
585 ! Take care of optional input variables:
586 lfor_night = .true.
587 if(present(for_night)) lfor_night = for_night
588
589 lcwarn = .true.
590 if(present(cwarn)) lcwarn = cwarn
591
592 alt = 0.d0; ha = 0.d0; h0 = 0.d0
593 accur = 1.d-6 ! Accuracy. 1d-6 ~ 0.1s. Don't make this smaller than 1d-16
594 event = ['Transit time ','Rise time ','Set time ']
595
596 rsa = rsalt*d2r
597 evmax = 3 ! 'Maximum' event to compute: transit only: evMax=1, +rise/set: evMax=3
598 if(rsalt.gt.90.d0) evmax = 1 ! Compute transit time and altitude only
599
600 call jd2cal(jd, yr,mnt,dy)
601 day0 = dble(floor(dy))-tz/24.d0 ! Midnight local time of the desired day, needed for agst0
602 if(lfor_night) day0 = day0 + 0.5d0 ! Noon local time at the start of the desired NIGHT, needed for agst0
603 jd0 = cal2jd(yr,mnt,day0)
604
605 call planet_position(jd0,3)
606 agst0 = planpos(45) ! AGST for midnight (or noon)
607
608 tmdy = 0.d0
609 azalt = 0.d0
610 if(evmax.eq.3) then
611 cosh0 = (sin(rsa)-sin(lat0)*sin(dec))/(cos(lat0)*cos(dec)) ! Cosine of the hour angle of rise/set; Meeus, Eq.15.1
612 if(abs(cosh0).gt.1.d0) then ! Body never rises/sets
613 evmax = 1 ! Compute transit time and altitude only
614 else
615 h0 = rev(2*acos(cosh0))/2.d0 ! Hour angle of rise/set
616 end if
617 end if
618
619 tmdy(1) = rev(ra - lon0 - agst0)/pi2 ! Transit; tmdy(1)=m0, tmdy(2)=m1, tmdy(3)=m2 in Meeus, but lon0 > 0 for E
620 if(evmax.eq.3) then
621 tmdy(2) = rev(tmdy(1)*pi2 - h0)/pi2 ! Rise
622 tmdy(3) = rev(tmdy(1)*pi2 + h0)/pi2 ! Set
623 end if
624
625
626 do evi=1,evmax
627 iter = 0
628 dtm = huge(dtm)
629 do while(abs(dtm).ge.accur)
630 th0 = agst0 + 6.300388092591991d0*tmdy(evi) ! Siderial day in radians; Meeus, p.103
631
632 ha = rev2(th0 + lon0 - ra) ! Hour angle; Meeus p.103
633 alt = asin(sin(lat0)*sin(dec) + cos(lat0)*cos(dec)*cos(ha)) ! Altitude; Meeus Eq. 13.6
634
635 ! Correction on transit/rise/set times:
636 if(evi.eq.1) then ! Transit
637 dtm = -ha/pi2
638 else ! Rise/set
639 dtm = (alt-rsa)/(pi2*cos(dec)*cos(lat0)*sin(ha))
640 end if
641 tmdy(evi) = tmdy(evi) + dtm
642
643 iter = iter+1
644 if(iter.gt.30) exit ! do-while loop
645 end do ! do while(abs(dtm).ge.accur)
646
647
648 if(iter.gt.30) then ! Convergence failed
649 if(lcwarn) write(0,'(A,F10.3,A)') ' * WARNING: riset_ad(): Riset failed to converge: '//trim(event(evi)),rsalt,'d *'
650
651 tmdy(evi) = 0.d0
652 azalt(evi) = 0.d0
653 else ! Result converged, store it
654 if(evi.eq.1) then
655 azalt(evi) = alt ! Transit altitude
656 else
657 azalt(evi) = atan2(sin(ha),(cos(ha)*sin(lat0)-tan(dec)*cos(lat0))) ! Rise,set hour angle -> azimuth
658 end if
659 end if
660
661 if(tmdy(evi).lt.0. .and. deq0(rsalt)) then
662 tmdy(evi) = 0.d0
663 azalt(evi) = 0.d0
664 end if
665
666 end do ! evi
667
668
669 ! Store results:
670 tmhr = tmdy * 24 ! Times days -> hours
671 if(lfor_night) then
672 tmhr = tmhr + 12 ! Times relative to noon -> relative to midnight
673 where(deq(tmhr,12.d0)) tmhr = 0.d0
674 end if
675
676 rt = tmhr(2) ! Rise time
677 tt = tmhr(1) ! Transit time
678 st = tmhr(3) ! Set time
679
680 rh = azalt(2) ! Rise azimuth
681 ta = azalt(1) ! Transit altitude
682 sh = azalt(3) ! Set azimuth
683
684 end subroutine riset_ad
685 !*********************************************************************************************************************************
686
687
688
689 !*********************************************************************************************************************************
690 !> \brief Compute the date (m,d) at which an object with ra can be observed best, i.e., it transits at midnight
691 !!
692 !! \param jd00 Julian day?
693 !! \param ra0 Right ascension
694 !!
695 !! \param m Month of best observation (output)
696 !! \param d Day of month of best observation (output)
697 !!
698 !! \note
699 !! - VERY SLOW!!! -> use best_obs_date_ra instead!!! (but not (yet) in libTheSky...)
700 !! - riset() should be faster now (04/2010)
701
702 subroutine best_obs_date(jd00,ra0, m,d)
703 use sufr_kinds, only: double
704 use sufr_constants, only: r2h
705 use sufr_angles, only: rev,rv12
706 use sufr_date_and_time, only: jd2cal
707 use sufr_numerics, only: deq
708
711 use thesky_local, only: tz
712 use thesky_datetime, only: gettz
713
714 implicit none
715 real(double), intent(in) :: jd00,ra0
716 integer, intent(out) :: m,d
717 real(double) :: jd0,oldjd,z(nplanpos)
718 real(double) :: rt,stt,ott,st,rh,ta,sh,midnight,hh,dt,dd
719 integer :: y,iter
720
721 oldjd = -99.d0
722 jd0 = jd00
723
724 dd = 9999.d0
725 iter = 0
726 do while(abs(dd).gt.0.5d0)
727 call riset(jd0,3, rt,stt,st,rh,ta,sh, 0.d0) ! Sun transit time
728 midnight = mod(stt+36.d0,24.d0)
729
730 call planet_position(jd0,3) ! Position of the Sun
731 z = planpos
732
733 hh = rev(z(44) - ra0) ! Hour angle is LST - RA
734 ott = -hh*r2h ! Time of object transit = -ha
735 dt = rv12(ott-midnight) ! Difference between transit time and midnight on this day
736 dd = dt/24.d0*365.25d0 ! Estimate for difference in days between today and day we're looking for
737 jd0 = dble(nint(jd0 + dd + 0.5d0)) - 0.5d0 ! Force midnight UT
738 tz = gettz(jd0)
739 jd0 = jd0 - tz/24.d0 ! Local midnight
740
741 ! Consider result 'converged':
742 if(deq(jd0,oldjd)) exit
743 if(iter.gt.5 .and. abs(dd).lt.0.7d0) exit
744 if(iter.gt.10 .and. abs(dd).lt.1.d0) exit
745 if(iter.gt.100) exit
746
747 oldjd = jd0
748 iter = iter + 1
749 end do
750
751 call jd2cal(jd0, y,m,dd)
752 d = floor(dd)
753
754 end subroutine best_obs_date
755 !*********************************************************************************************************************************
756
757
758
759
760
761
762
763 !*********************************************************************************************************************************
764 !> \brief Quadratic interpolation of rise/set times for riset_ipol().
765 !!
766 !! \param y1 Coordinate for day 1
767 !! \param y2 Coordinate for day 2
768 !! \param y3 Coordinate for day 3
769 !! \param n Fraction of day
770 !!
771 !! \retval rsIpol Interpolated coordinate
772
773 function rsipol(y1,y2,y3, n)
774 use sufr_kinds, only: double
775
776 implicit none
777 real(double), intent(in) :: y1,y2,y3,n
778 real(double) :: rsipol, a,b,c
779
780 a = y2-y1
781 b = y3-y2
782 c = b-a
783
784 rsipol = y2 + n/2.d0 * (a + b + c*n) ! 1 + (1+y1/2)*n + (y1/2 + y2 + y3/2)*n^2
785
786 end function rsipol
787 !*********************************************************************************************************************************
788
789
790end module thesky_riset
791!***********************************************************************************************************************************
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...
Date and time procedures.
Definition date_time.f90:49
real(double) function gettz(jd, tz0, dsttp)
Returns time zone: tz0 (standard time) or tz0+1 (daylight-savings time)
Local parameters for libTheSky: location, date, time.
Definition modules.f90:56
real(double) tz
Current value of time zone, taking into account DST (hours; >0 is east of Greenwich)
Definition modules.f90:67
real(double) lon0
Longitude of the observer (rad)
Definition modules.f90:64
real(double) deltat
Current value of DeltaT (s)
Definition modules.f90:66
real(double) lat0
Latitude of the observer (rad)
Definition modules.f90:63
Planet data, needed to compute planet positions.
Definition modules.f90:85
integer, parameter nplanpos
Number of entries in the planpos array.
Definition modules.f90:116
real(double), dimension(nplanpos) planpos
Planpos[] is an array with many different types of coordinates and related variables:
Definition modules.f90:192
Procedures for planets.
Definition planets.f90:23
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:63
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:1227
Procedures to determine rise, transit and set times.
Definition riset.f90:24
subroutine riset_ipol(jd, pl, rt, tt, st, rh, ta, sh, rsalt, for_night, ltime, cwarn)
Old routine for rise, transit and set times for planets and asteroids - don't use this for Sun and Mo...
Definition riset.f90:319
subroutine riset_ad(jd, ra, dec, rt, tt, st, rh, ta, sh, rsalt, for_night, cwarn)
Rise, transit and set times routine for an object with fixed ra & dec.
Definition riset.f90:565
subroutine best_obs_date(jd00, ra0, m, d)
Compute the date (m,d) at which an object with ra can be observed best, i.e., it transits at midnight
Definition riset.f90:703
subroutine riset(jd, pl, rt, tt, st, rh, ta, sh, rsalt, for_night, ltime, cwarn, converge)
Rise, transit and set times routine for Sun, Moon, planets and asteroids. This version recomputes pos...
Definition riset.f90:83
real(double) function rsipol(y1, y2, y3, n)
Quadratic interpolation of rise/set times for riset_ipol().
Definition riset.f90:774