libTheSky
Routines to compute sky positions of Sun, Moon, planets and more
All Namespaces Files Functions Variables Pages
data.f90
Go to the documentation of this file.
1!> \file data.f90 Procedures to define constants and read data files for libTheSky
2
3
4! Copyright (c) 2002-2025 Marc van der Sluys - Nikhef/Utrecht University - marc.vandersluys.nl
5!
6! This file is part of the libTheSky package,
7! see: https://libthesky.sf.net/
8!
9! This is free software: you can redistribute it and/or modify it under the terms of the
10! European Union Public Licence 1.2 (EUPL 1.2).
11!
12! This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13! without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14! See the EU Public Licence for more details.
15!
16! You should have received a copy of the European Union Public Licence along with this code.
17! If not, see <https://www.eupl.eu/1.2/en/>.
18
19
20!***********************************************************************************************************************************
21!> \brief Procedures to set constants and read data files
22
24 implicit none
25 save
26
27
28contains
29
30
31 !*********************************************************************************************************************************
32 !> \brief Set the initial values of the variables and constants used in this package
33
35 use sufr_constants, only: set_sufr_constants, homedir
36 use sufr_system, only: error,quit_program_error
38
39 implicit none
40 integer, parameter :: ntrials = 11
41 integer :: try
42 character :: tryPaths(ntrials)*(99)
43 logical :: ex
44
45 call set_sufr_constants() ! Set libSUFR constants
46
47 library_name = 'libTheSky (libthesky.sourceforge.net)'
48
49 thesky_verbosity = 99 ! 99: print all messages - 0: be quiet
50
51 ! Find the libTheSky path to the data files:
52 trypaths(1) = trim(theskydatadir) ! Predefined?
53 if(len_trim(theskydatadir).eq.0 .or. len_trim(theskydatadir).eq.len(theskydatadir)) trypaths(1) = './' ! Current working dir
54
55 ! Try the user's home directory first:
56 trypaths(2) = trim(homedir)//'/share/libTheSky/'
57 trypaths(3) = trim(homedir)//'/local/share/libTheSky/'
58 trypaths(4) = trim(homedir)//'/usr/share/libTheSky/'
59 trypaths(5) = trim(homedir)//'/usr/local/share/libTheSky/'
60 trypaths(6) = trim(homedir)//'/opt/share/libTheSky/'
61 trypaths(7) = trim(homedir)//'/opt/local/share/libTheSky/'
62
63 ! Try default system paths:
64 trypaths(8) = '/usr/share/libTheSky/'
65 trypaths(9) = '/usr/local/share/libTheSky/'
66 trypaths(10) = '/opt/share/libTheSky/'
67 trypaths(11) = '/opt/local/share/libTheSky/'
68
69
70 do try = 1,ntrials
71 theskydatadir = trypaths(try)
72 inquire(file=trim(theskydatadir)//'planets.dat', exist=ex)
73 if(ex) exit
74 end do
75
76 if(.not.ex) then
77 call error('I could not find the libTheSky data directory. I tried the following possibilities:', 0)
78 do try= 1,ntrials
79 write(0,'(I4,2x,A)') try, trim(trypaths(try))
80 end do
81 write(0,'(A)') 'Please make sure you have the libTheSky data files installed. They can be found in the'
82 write(0,'(A)') 'libthesky-data-<YYYYMMDD>.tar.bz2 tarball, which can be downloaded from libthesky.sf.net'
83 write(0,'(A)') 'or github.com/astronomy/libTheSky. See the libTheSky INSTALL file for details.'
84 call quit_program_error('libTheSky: Data directory not found.', 1)
85 end if
86
87 deltat_forced = -huge(deltat_forced) ! Use this to force the DeltaT to a given value (if > -9999)
88 jd1820 = 2385800.5d0 ! JD for 1820-01-01 0h UT, used for DeltaT approximation
89
90 end subroutine set_thesky_constants
91 !*********************************************************************************************************************************
92
93
94
95 !*********************************************************************************************************************************
96 !> \brief Read input data files
97
98 subroutine thesky_readdata()
99 implicit none
100
101 call read_deltat()
102
103 call readplanetdata()
104 call readpluto()
105
106 call readmoondata_la()
107
108 call readplanetelements()
109
111 call readcometelements()
112 call readnutation()
113 call readconstellations()
114
115 end subroutine thesky_readdata
116 !*********************************************************************************************************************************
117
118
119
120
121
122 !*********************************************************************************************************************************
123 !> \brief Read deltat.dat with historical values for DeltaT
124 !!
125 !! \note
126 !! \see http://hemel.waarnemen.com/Computing/deltat.html
127
128 subroutine read_deltat()
129 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_error_quit
130 use sufr_dummy, only: dumstr
133
134 implicit none
135 integer :: i,status, dn, ip
136 character :: infile*(199)
137
138 call find_free_io_unit(ip)
139 infile = trim(theskydatadir)//'deltat.dat'
140 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
141 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
142
143 do i=1,7
144 read(ip,*) dumstr
145 end do
146
147 do i=1,deltat_nmax
148 read(ip,*, iostat=status) deltat_years(i),deltat_values(i)
149 if(status.lt.0) exit
150 if(status.gt.0) call file_read_error_quit(trim(infile), i, 1) ! Exit status 1
151 end do
152
153 close(ip)
154
155 deltat_n = i-1
156 deltat_minyr = nint(deltat_years(1))
158
159
160 ! Set the average acceleration of DeltaT, from -1500 to 2010 CE:
161 ! see http://hemel.waarnemen.com/Computing/deltat.html
162 deltat_accel = 3.975539d-3 ! Seconds / yr^2: 39.755 +- 0.048 s/cty^2 (fit.f90)
163
164 ! Compute the average change in DeltaT in the last dn years:
165 dn = 100
167
168 ! Take the most recent tabulated value of DeltaT:
170
171 end subroutine read_deltat
172 !*********************************************************************************************************************************
173
174
175
176
177
178
179
180
181
182
183
184 !*********************************************************************************************************************************
185 !> \brief Read constellation names: abbreviation, Latin, genitive, Dutch and English
186
188 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_end_error
189 use sufr_dummy, only: dumstr9
190
194
195 implicit none
196 integer :: i,j, ip,status
197 character :: infile*(199)
198
199 ! Read constellation names (abbr., Lat, Lat.Gen., Nl, En):
200 call find_free_io_unit(ip)
201 infile = trim(theskydatadir)//'const_names.dat'
202 open(ip, status='old', action='read', file=trim(infile), iostat=status)
203 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
204
205 read(ip,*) dumstr9
206 read(ip,*) dumstr9
207
208 do i=1,nconstel
209 read(ip,'(A3,2x,A19,2x,A19,2x,A17,2x,A18)', iostat=status) conabr(i),latconnames(i),genconnames(i), &
211
212 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
213 end do
214
215 close(ip)
216
217
218 ! Read constellation ID data:
219 infile = trim(theskydatadir)//'const_id.dat'
220 open(ip, status='old', action='read', file=trim(infile), iostat=status)
221 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
222
223 do i=1,nconid
224 read(ip,'(2F8.4,F9.4,1x,A3)', iostat=status) conidral(i),conidrau(i),coniddecl(i),conidabr(i)
225 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
226
227 do j=1,nconstel
228 if(conabr(j).eq.conidabr(i)) conid(i) = j ! Assign a constellation number ID rather than abbrev.
229 end do
230 end do
231
232 close(ip)
233 end subroutine readconstellations
234 !*********************************************************************************************************************************
235
236
237
238
239 !*********************************************************************************************************************************
240 !> \brief Read nutation input files
241
242 subroutine readnutation()
243 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_end_error
244 use sufr_kinds, only: double
246
247 implicit none
248 integer :: i,nu1(6),nu3, ip,status
249 real(double) :: nu2,nu4
250 character :: infile*(199)
251
252 call find_free_io_unit(ip)
253 infile = trim(theskydatadir)//'nutation.dat'
254 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
255 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
256
257 do i=1,63
258 read(ip,'(5(I2),2x,I7,2x,F6.1,2x,I6,2x,F4.1)', iostat=status) nu1,nu2,nu3,nu4
259 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
260
261 nutationdat(1:6,i) = nu1
262 nutationdat(7,i) = nu2
263 nutationdat(8,i) = nu3
264 nutationdat(9,i) = nu4
265 end do
266 close(ip)
267
268 end subroutine readnutation
269 !*********************************************************************************************************************************
270
271
272
273 !*********************************************************************************************************************************
274 !> \brief Reads VSOP planet data from planets.dat
275
276 subroutine readplanetdata()
277 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_end_error
278 use sufr_constants, only: plana,au
281
282 implicit none
283 integer :: pl, li, ip, status, powr,opowr, var, TOTnls(8)
284 character :: infile*(199)
285
286 vsopnls = reshape( (/ 2808,1620,2399, 671,426,585, 1080,349,997, & ! Number of lines in VSOP input files (l,b,r x 8 pl)
287 2393,915,2175, 1484,530,1469, 2358,966,2435, 1578,516,1897, 681,290,959/), (/3,8/))
288
289 totnls(1:8) = sum(vsopnls(:,1:8), 1) ! Total number of lines per planet; sum along dimension 1
290
291 ! VSOP87 terms are truncated at these accuracies:
292 vsoptruncs(1,:) = (/0.05d0, 0.5d0, 0.5d0, 0.5d0, 1.d0, 1.d0, 1.6d0, 1.d0/) * 1.d-9 ! L (rad)
293 vsoptruncs(2,:) = vsoptruncs(1,:) ! B (rad)
294 vsoptruncs(3,:) = vsoptruncs(1,:) * plana(1:8)/au ! R (AU)
295
296 ! Read planets.dat file:
297 call find_free_io_unit(ip)
298 infile = trim(theskydatadir)//'planets.dat'
299 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
300 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
301
302 vsopdat = 0.d0
303 opowr = 9 ! should be > powr = 0-5
304 do pl=1,8 ! Planet
305 var = 0 ! L,B,R
306 do li=1,totnls(pl)
307 read(ip,'(I1,F18.11,F14.11,F20.11)', iostat=status) powr, vsopdat(2:4,li,pl)
308 if(status.ne.0) call file_read_end_error(trim(infile), li, status, 1, 1) ! stopcode=1, exitstatus=1
309 vsopdat(1,li,pl) = dble(powr)
310
311 if(powr.ne.opowr) then ! Save the line number where the next block of (Planet, Variable (LBR), Power) starts
312 if(powr.lt.opowr) var = var+1 ! L,B,R
313 vsopnblk(powr,var,pl) = li
314 end if
315 opowr = powr
316 end do ! li
317 end do ! pl
318
319 close(ip)
320
321 end subroutine readplanetdata
322 !*********************************************************************************************************************************
323
324
325
326
327 !*********************************************************************************************************************************
328 !> \brief Read periodic terms for the position of Pluto
329
330 subroutine readpluto()
331 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_end_error
334
335 implicit none
336 integer :: i, ip, status
337 character :: infile*(199)
338
339 call find_free_io_unit(ip)
340 infile = trim(theskydatadir)//'pluto.dat'
341 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
342 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
343
344 do i=1,43
345 read(ip,'(3I3,6I10)', iostat=status) pluc(i,1), pluc(i,2), pluc(i,3), plul(i,1), plul(i,2), plub(i,1), plub(i,2), &
346 plur(i,1), plur(i,2)
347
348 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
349 end do
350 close(ip)
351
352 end subroutine readpluto
353 !*********************************************************************************************************************************
354
355
356 !*********************************************************************************************************************************
357 !> \brief Read the data needed to compute the orbital elements of the planets
358 !!
359 !! \note
360 !! - use planetelements() to compute the orbital elements
361 !!
362 !! - L, a, e, i, Omega, pi for Mer-Nep, mean equinox of date and J2000.0 (=2x6x8)
363 !! - Terms in one row are a_0, a_1, a_2, a_3 and should be used as: Sum(i=0,3) a_i t^i, with t the time in centuries since 2000.0
364 !! - Angles are still expressed in degrees
365
367 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_end_error
370
371 implicit none
372 integer :: eq,pl,el, ip,status
373 character :: infile*(199)
374
375 call find_free_io_unit(ip)
376 infile = trim(theskydatadir)//'planetelements.dat'
377 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
378 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
379
380 outerloop: do eq=1,2 ! Equinoctes
381 do pl=1,8 ! Planets
382 do el=1,6 ! Elements
383 if(eq.eq.2.and.(el.eq.2.or.el.eq.3)) cycle ! a, e are not in the second part (J2000.0?)
384
385 read(ip,'(F13.9, F18.10, 2F15.11)', iostat=status) plelemdata(eq, pl, el, 0:3)
386
387 if(status.ne.0) then
388 call file_read_end_error(trim(infile), 0, status, 0, 0) ! line=0, stopcode=0, exitstatus=0
389 exit outerloop
390 end if
391
392 end do ! el
393 end do ! pl
394 end do outerloop ! eq
395 close(ip)
396
397 plelemdata(2,1:8,2,0:3) = plelemdata(1,1:8,2,0:3) ! a_J2000 = a
398 plelemdata(2,1:8,3,0:3) = plelemdata(1,1:8,3,0:3) ! e_J2000 = e
399
400 end subroutine readplanetelements
401 !*********************************************************************************************************************************
402
403
404 !*********************************************************************************************************************************
405 !> \brief Read ELP-82B periodic terms from moondata.dat
406 !!
407 !! \note
408 !! Constants are defined in the module moondata
409 !!
410 !! \see
411 !! - Chapront-Touzé & Chapront, A&A, 124, 50 (1983)
412 !! - Chapront-Touzé & Chapront, A&A, 190, 342 (1988)
413 !! - ftp://cdsarc.u-strasbg.fr/pub/cats/VI/79/
414
415 subroutine readmoondata()
416 use sufr_kinds, only: double, dbl
417 use sufr_constants, only: pi, pio2, d2r, as2r, r2as
418 use sufr_system, only: find_free_io_unit,file_open_error_quit, file_read_end_error
419 use sufr_numerics, only: dne
420
424 use thesky_moondata, only: p1,p2,p3,p4,p5,q1,q2,q3,q4,q5, prec0
425
426 implicit none
427 integer :: i,ific,ir,itab,iv,j,k,iz, ip,status, filen(36),fl,il
428 real(double) :: prec,am,alpha,dtasm,precess, delnu,dele,delg,delnp,delep, xx,tgv,y,pha
429 character :: infile*(199)
430
431 filen = (/1023,918,704,347,316,237,14,11,8,14328,5233,6631,4384, &
432 833,1715,170,150,114,226,188,169,3,2,2,6,4,5,20,12,14,11,4,10,28,13,19/)
433 prec = 0.0_dbl
434
435 ! Make sure these are always defined:
436 delep = 0.0_dbl
437 delnp = 0.0_dbl
438 delg = 0.0_dbl
439 dele = 0.0_dbl
440 delnu = 0.0_dbl
441 dtasm = 0.0_dbl
442 am = 0.0_dbl
443
444 !*** Parameters:
445 if(ideb.eq.0) then ! then it's the first time you run this routine
446
447 ideb = 1
448 am = 0.074801329518d0 ! Ratio of the mean motions (EMB / Moon)
449 alpha = 0.002571881335d0 ! Ratio of the semi-major axis (Moon / EMB)
450 dtasm = 2.d0*alpha/(3.d0*am)
451
452 ! Lunar arguments:
453 ! W1: mean longitude of the Moon:
454 w(1,0) = (218.d0+18.d0*c1 + 59.95571d0*c2) * d2r
455 w(1,1) = 1732559343.73604d0 * as2r
456 ! w(1,2) = -5.8883d0 * as2r ! Original, but erroneous!
457 w(1,2) = -6.8084d0 * as2r ! CHECK - was -5.8883
458 ! Also -5.8883 in elp82b.pdf (originally PS?), Table C.
459 ! Should this be -6.8883, -6.8083 or -6.8084? All constants W(1,*) are identical to MPP02, except this one!
460 ! Effect w.r.t. MPP02 in 3000 BCE (one test calculation per century):
461 ! - -5.8883: Delta longitude gradually drifts to +0.65°, dlat ~ +-0.06°, ddist ~ +-300km
462 ! - -6.8883: dlon drifts to -0.065°, dlat ~ +-0.005°, ddist ~ +-13km
463 ! - -6.8083: dlon no longer drifts, ~ +-0.012°, dlat ~ +-0.006°, ddist ~ +-20km
464 ! - -6.8084: same as above
465 ! - This value comes from DE405 fits(?)
466 w(1,3) = 0.6604d-2 * as2r
467 w(1,4) = -0.3169d-4 * as2r
468
469 ! W2: mean longitude of the lunar perigee:
470 w(2,0) = (83.d0+21.d0*c1 + 11.67475d0*c2) * d2r
471 w(2,1) = 14643420.2632d0 * as2r
472 w(2,2) = -38.2776d0 * as2r
473 w(2,3) = -0.45047d-1 * as2r
474 w(2,4) = 0.21301d-3 * as2r
475
476 ! W3: mean longitude of the lunar ascending node:
477 w(3,0) = (125.d0+2.d0*c1 + 40.39816d0*c2) * d2r
478 w(3,1) = -6967919.3622d0 * as2r
479 w(3,2) = 6.3622d0 * as2r
480 w(3,3) = 0.7625d-2 * as2r
481 w(3,4) = -0.3586d-4 * as2r
482
483 ! Earth-Moon (EMB) elements:
484 ! Te: mean longitude of EMB:
485 eart(0) = (100.d0+27.d0*c1 + 59.22059d0*c2) * d2r
486 eart(1) = 129597742.2758d0 * as2r
487 eart(2) = -0.0202d0 * as2r
488 eart(3) = 0.9d-5 * as2r
489 eart(4) = 0.15d-6 * as2r
490
491 ! Pip: mean longitude of the perihelion of EMB:
492 peri(0) = (102.d0+56.d0*c1 + 14.42753d0*c2) * d2r
493 peri(1) = 1161.2283d0 * as2r
494 peri(2) = 0.5327d0 * as2r
495 peri(3) = -0.138d-3 * as2r
496 peri(4) = 0.d0
497
498
499 ! Precession constant:
500 precess = 5029.0966d0 * as2r
501
502 ! Planetary arguments: mean longitudes:
503 p(1,0) = (252.d0+15.d0*c1+3.25986d0*c2) * d2r
504 p(2,0) = (181.d0+58.d0*c1+47.28305d0*c2) * d2r
505 p(3,0) = eart(0)
506 p(4,0) = (355.d0+25.d0*c1+59.78866d0*c2) * d2r
507 p(5,0) = (34.d0+21.d0*c1+5.34212d0*c2) * d2r
508 p(6,0) = (50.d0+4.d0*c1+38.89694d0*c2) * d2r
509 p(7,0) = (314.d0+3.d0*c1+18.01841d0*c2) * d2r
510 p(8,0) = (304.d0+20.d0*c1+55.19575d0*c2) * d2r
511
512 ! Planetary arguments: mean motions:
513 p(1,1) = 538101628.68898d0 * as2r
514 p(2,1) = 210664136.43355d0 * as2r
515 p(3,1) = eart(1)
516 p(4,1) = 68905077.59284d0 * as2r
517 p(5,1) = 10925660.42861d0 * as2r
518 p(6,1) = 4399609.65932d0 * as2r
519 p(7,1) = 1542481.19393d0 * as2r
520 p(8,1) = 786550.32074d0 * as2r
521
522 ! Corrections of the constants (fit to DE200/LE200):
523 delnu = +0.55604d0/w(1,1) * as2r
524 dele = +0.01789d0 * as2r
525 delg = -0.08066d0 * as2r
526 delnp = -0.06424d0/w(1,1) * as2r
527 delep = -0.12879d0 * as2r
528
529 ! Delaunay's arguments (https://en.wikipedia.org/wiki/Orbital_elements#Delaunay_variables):
530 do i=0,4
531 del(1,i) = w(1,i) - eart(i)
532 del(4,i) = w(1,i) - w(3,i)
533 del(3,i) = w(1,i) - w(2,i)
534 del(2,i) = eart(i) - peri(i)
535 end do
536
537 del(1,0) = del(1,0) + pi
538 zeta(0) = w(1,0)
539 zeta(1) = w(1,1) + precess
540
541 ! Precession of the longitude of the ascending node of the mean ecliptic of date on fixed ecliptic J2000 (Laskar, 1986):
542 ! P: sine coefficients:
543 p1 = 0.10180391d-4
544 p2 = 0.47020439d-6
545 p3 = -0.5417367d-9
546 p4 = -0.2507948d-11
547 p5 = 0.463486d-14
548
549 ! Q: cosine coefficients:
550 q1 = -0.113469002d-3
551 q2 = 0.12372674d-6
552 q3 = 0.1265417d-8
553 q4 = -0.1371808d-11
554 q5 = -0.320334d-14
555
556 end if ! if(ideb.eq.0)
557
558
559
560 ! Read moondata file:
561 call find_free_io_unit(ip)
562 infile = trim(theskydatadir)//'moondata.dat'
563 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
564 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
565
566 if(dne(prec,prec0)) then
567 prec0 = prec
568 pre(1) = prec * r2as - 1.d-12
569 pre(2) = prec * r2as - 1.d-12
570 pre(3) = prec * ath
571
572 do ific=1,36 ! Loop over all original ELP files
573 fl = filen(ific)
574 ir = 0
575 itab = (ific+2)/3
576 iv = mod(ific-1,3) + 1
577
578 select case(ific)
579 case(1:3) ! Files: Main problem
580 do il=1,fl
581 read(ip,'(4I3,2x,F13.5,6(2x,F10.2))', iostat=status) ilu,coef
582 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
583
584 ir = ir+1
585 tgv = coef(2)+dtasm*coef(6)
586 if(ific.eq.3) coef(1) = coef(1)-2.d0*coef(1)*delnu/3.d0
587 xx = coef(1) + tgv*(delnp-am*delnu) + coef(3)*delg + coef(4)*dele + coef(5)*delep
588 zone(1) = xx
589 do k=0,4
590 y = 0.d0
591 do i=1,4
592 y = y + ilu(i)*del(i,k)
593 end do
594 zone(k+2) = y
595 end do
596 if(iv.eq.3) zone(2) = zone(2) + pio2
597 do i=1,6
598 if(iv.eq.1) pc1(i,ir) = zone(i)
599 if(iv.eq.2) pc2(i,ir) = zone(i)
600 if(iv.eq.3) pc3(i,ir) = zone(i)
601 end do
602 end do ! do il=1,fl
603
604
605 case(4:9,22:36) ! Files: Tides - Relativity - Solar eccentricity
606 do il=1,fl
607 read(ip,'(5I3,1x,F9.5,1x,F9.5)', iostat=status) iz,ilu,pha,xx
608 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
609
610 ir = ir + 1
611 zone(1) = xx
612 do k=0,1
613 if(k.eq.0) y = pha*d2r
614 if(k.ne.0) y = 0.d0
615 y = y+iz*zeta(k)
616 do i=1,4
617 y = y+ilu(i)*del(i,k)
618 end do
619 zone(k+2) = y
620 end do
621 j = nrang(iv,itab-1)+ir
622 do i=1,3
623 if(iv.eq.1) per1(i,j) = zone(i)
624 if(iv.eq.2) per2(i,j) = zone(i)
625 if(iv.eq.3) per3(i,j) = zone(i)
626 end do
627 end do ! do il=1,fl
628
629
630 case(10:21) ! Files: Planetary perturbations
631 do il=1,fl
632 read(ip,'(11I3,1x,F9.5,1x,F9.5)', iostat=status) ipla,pha,xx
633 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
634
635 ir = ir+1
636 zone(1) = xx
637 if(ific.lt.16) then
638 do k=0,1
639 if(k.eq.0) y = pha*d2r
640 if(k.ne.0) y = 0.d0
641 y = y + ipla(9)*del(1,k) + ipla(10)*del(3,k) + ipla(11)*del(4,k)
642 do i=1,8
643 y = y+ipla(i)*p(i,k)
644 end do
645 zone(k+2) = y
646 end do
647 else
648 do k=0,1
649 if(k.eq.0) y = pha*d2r
650 if(k.ne.0) y = 0.d0
651 do i=1,4
652 y = y+ipla(i+7)*del(i,k)
653 end do
654 do i=1,7
655 y = y+ipla(i)*p(i,k)
656 end do
657 zone(k+2) = y
658 end do
659 end if
660 j = nrang(iv,itab-1) + ir
661 do i=1,3
662 if(iv.eq.1) per1(i,j) = zone(i)
663 if(iv.eq.2) per2(i,j) = zone(i)
664 if(iv.eq.3) per3(i,j) = zone(i)
665 end do
666 end do ! do il=1,fl
667
668 end select
669
670 nterm(iv,itab) = ir
671 if(itab.eq.1) then
672 nrang(iv,itab) = 0
673 else
674 nrang(iv,itab) = nrang(iv,itab-1) + nterm(iv,itab)
675 end if
676
677 end do ! do ific=1,36 - original ELP files
678
679 end if ! if(prec.ne.prec0)
680 close(ip)
681
682 end subroutine readmoondata
683 !*********************************************************************************************************************************
684
685
686
687
688
689
690 !*********************************************************************************************************************************
691 !> \brief Read low-accuracy moon-position data from Meeus (moon_la.dat)
692
693 subroutine readmoondata_la()
694 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_end_error
695 use sufr_dummy, only: dumstr
698
699 implicit none
700 integer :: i, ip, status
701 character :: infile*(199)
702
703 call find_free_io_unit(ip)
704 infile = trim(theskydatadir)//'moon_la.dat'
705 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
706 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
707
708 read(ip,*) dumstr
709 read(ip,*) dumstr
710 read(ip,*) dumstr
711
712 do i=1,60
713 read(ip,'(I1,3I3,2I10)', iostat=status) moonla_arg(1:4,i),moonla_lrb(1:2,i) ! Read terms for L and R
714 if(status.ne.0) call file_read_end_error(trim(infile), 3+i, status, 1, 1) ! stopcode=1, exitstatus=1
715 end do
716
717 read(ip,*) dumstr
718 do i=1,60
719 read(ip,'(I1,3I3,2I10)', iostat=status) moonla_arg(5:8,i),moonla_lrb(3,i) ! Read terms for B
720 if(status.ne.0) call file_read_end_error(trim(infile), 64+i, status, 1, 1) ! stopcode=1, exitstatus=1
721 end do
722 close(ip)
723
724 end subroutine readmoondata_la
725 !*********************************************************************************************************************************
726
727
728
729 !*********************************************************************************************************************************
730 !> \brief Initialise ELP-MPP02 data and read the data files.
731 !!
732 !! \param mode Index of the corrections to the constants: 0: LLR observations for 1970-2000, 1: DE405 ephemeris for 1950-2060 (input)
733 !! \param ierr File error index: ierr=0: no error, ierr=1: file error (output)
734
736 implicit none
737 integer, intent(in) :: mode
738 integer, intent(out) :: ierr
739
740 integer, save :: modeInit = 999 ! Test whether data have been initialised
741
742 ! Initializing of constants and reading the files:
743 ierr = 0
744 ! print*,mode,modeInit
745 if(mode.ne.modeinit) then
746 call elp_mpp02_initialise(mode)
747 call elp_mpp02_read_files(ierr)
748 if(ierr.ne.0) return
749
750 modeinit = mode
751 end if
752
754 !*********************************************************************************************************************************
755
756
757
758 !*********************************************************************************************************************************
759 !> \brief Initialization of the constants and parameters used for the evaluation of the ELP/MPP02 series
760 !!
761 !! \param mode Index of the corrections to the constants: 0: LLR observations for 1970-2000, 1: DE405 ephemeris for 1950-2060
762 !!
763 !! \note
764 !!
765 !! - Returns:
766 !! Set of the constants of ELPMPP02 solution in module TheSky_elp_mpp02_constants.
767 !!
768 !! - Remarks:
769 !! The nominal values of some constants have to be corrected. There are two sets of corrections, one of which can be chosen
770 !! using the parameter 'mode' (used in elp_mpp02_initialise()):
771 !! - mode=0, the constants are fitted to LLR observations provided from 1970 to 2001; it is the default value;
772 !! - mode=1, the constants are fitted to DE405 ephemeris over one century (1950-2060); the lunar angles W1, W2, W3 receive also additive corrections to the secular coefficients.
773 !!
774 !! - Moon constants:
775 !! nu : mean motion of the Moon (W1(1,1)) (Nu)
776 !! g : half coefficient of sin(F) in latitude (Gamma)
777 !! e : half coefficient of sin(l) in longitude (E)
778 !! np : mean motion of EMB (eart(1)) (n')
779 !! ep : eccentricity of EMB (e')
780 !!
781 !! p is the precession rate and t is the time
782 !!
783 !! \see
784 !! - ELPdoc: Lunar solution ELP, version ELP/MPP02, Jean Chapront and Gerard Francou, October 2002
785
786 subroutine elp_mpp02_initialise(mode)
787 use sufr_kinds, only: double
788 use sufr_constants, only: r2as, pi
789 use sufr_system, only: quit_program_error
790 use thesky_elp_mpp02_constants, only: w,eart,peri, zeta,del, p,delnu,dele,delg,delnp,delep,dtasm,am, p1,p2,p3,p4,p5, q1,q2,q3,q4,q5
791
792 implicit none
793 integer, intent(in) :: mode
794 integer :: iD
795 real(double) :: bp(5,2), x2,x3,y2,y3, d21,d22,d23,d24,d25, d31,d32,d33,d34,d35, Cw2_1,Cw3_1
796 real(double) :: alpha,xa, Deart_0,Dperi, Dw1_1,Dgam,De, Deart_1,Dep, Dw1_0,Dw1_2,Dw2_0,Dw2_1,Dw3_0,Dw3_1
797
798 ! Value of the correction to the constant of precession:
799 real(double), parameter :: Dprec = -0.29965d0 ! Source: IAU 2000A
800
801
802 data bp/ 0.311079095d0, -0.4482398d-2, -0.110248500d-2, 0.1056062d-2, 0.50928d-4, -0.103837907d0, 0.6682870d-3, -0.129807200d-2, -0.1780280d-3, -0.37342d-4 /
803
804 if(mode.lt.0 .or. mode.gt.1) call quit_program_error('elp_mpp02_initialise(): mode must have value 0 or 1', 1)
805
806 ! Constants for the evaluation of the partial derivatives:
807 am = 0.074801329d0 ! Ratio of the mean motions (EMB / Moon)
808 alpha = 0.002571881d0 ! Ratio of the semi-major axis (Moon / EMB)
809 dtasm = (2*alpha)/(3*am) ! (2*alpha) / (3*am)
810 xa = (2*alpha)/3.d0
811
812
813 ! Corrections to constants:
814 if(mode.eq.0) then ! Default - LLR
815 ! Values of the corrections to the constants fitted to LLR. Fit 13-05-02 (2 iterations) except Phi and eps w2_1 and w3_1
816 ! See ELPdoc, Table 3 and paper, Table 1
817 dw1_0 = -0.10525d0
818 dw2_0 = 0.16826d0
819 dw3_0 = -0.10760d0
820 deart_0 = -0.04012d0
821 dperi = -0.04854d0
822 dw1_1 = -0.32311d0
823 dgam = 0.00069d0
824 de = 0.00005d0
825 deart_1 = 0.01442d0
826 dep = 0.00226d0
827 dw2_1 = 0.08017d0
828 dw3_1 = -0.04317d0
829 dw1_2 = -0.03794d0
830 else ! DE 405
831 ! Values of the corrections to the constants fitted to DE405 over the time interval (1950-2060)
832 dw1_0 = -0.07008d0
833 dw2_0 = 0.20794d0
834 dw3_0 = -0.07215d0
835 deart_0 = -0.00033d0
836 dperi = -0.00749d0
837 dw1_1 = -0.35106d0
838 dgam = 0.00085d0
839 de = -0.00006d0
840 deart_1 = 0.00732d0
841 dep = 0.00224d0
842 dw2_1 = 0.08017d0
843 dw3_1 = -0.04317d0
844 dw1_2 = -0.03743d0
845 end if
846
847 ! Fundamental arguments (Moon and EMB - ELPdoc, Table 1):
848 ! W1: mean longitude of the Moon:
849 w(1,0) = elp_dms2rad(218,18,59.95571d0+dw1_0) ! Source: ELP
850 w(1,1) = (1732559343.73604d0+dw1_1)/r2as ! Source: ELP
851 w(1,2) = ( -6.8084d0 +dw1_2)/r2as ! Source: DE405
852 w(1,3) = 0.66040d-2/r2as ! Source: ELP
853 w(1,4) = -0.31690d-4/r2as ! Source: ELP
854
855 ! W2: mean longitude of the lunar perigee:
856 w(2,0) = elp_dms2rad( 83,21,11.67475d0+dw2_0) ! Source: ELP
857 w(2,1) = ( 14643420.3171d0 +dw2_1)/r2as ! Source: DE405
858 w(2,2) = ( -38.2631d0)/r2as ! Source: DE405
859 w(2,3) = -0.45047d-1/r2as ! Source: ELP
860 w(2,4) = 0.21301d-3/r2as ! Source: ELP
861
862 ! W3: mean longitude of the lunar ascending node:
863 w(3,0) = elp_dms2rad(125, 2,40.39816d0+dw3_0) ! Source: ELP
864 w(3,1) = ( -6967919.5383d0 +dw3_1)/r2as ! Source: DE405
865 w(3,2) = 6.3590d0/r2as ! Source: DE405
866 w(3,3) = 0.76250d-2/r2as ! Source: ELP
867 w(3,4) = -0.35860d-4/r2as ! Source: ELP
868
869 ! Earth-Moon (EMB) elements:
870 ! Te: mean longitude of EMB:
871 eart(0) = elp_dms2rad(100,27,59.13885d0+deart_0) ! Source: VSOP2000
872 eart(1) = (129597742.29300d0 +deart_1)/r2as ! Source: VSOP2000
873 eart(2) = -0.020200d0/r2as ! Source: ELP
874 eart(3) = 0.90000d-5/r2as ! Source: ELP
875 eart(4) = 0.15000d-6/r2as ! Source: ELP
876
877 ! Pip: mean longitude of the perihelion of EMB:
878 peri(0) = elp_dms2rad(102,56,14.45766d0+dperi) ! Source: VSOP2000
879 peri(1) = 1161.24342d0/r2as ! Source: VSOP2000
880 peri(2) = 0.529265d0/r2as ! Source: VSOP2000
881 peri(3) = -0.11814d-3/r2as ! Source: VSOP2000
882 peri(4) = 0.11379d-4/r2as ! Source: VSOP2000
883
884 ! Corrections to the secular terms of Moon angles. This gives a better (long-term?) fit
885 ! to DE 406. See ELPdoc, Table 6/paper, Table 4, line 2:
886 if(mode.eq.1) then ! DE 405 / DE 406
887 w(1,3) = w(1,3) - 0.00018865d0/r2as
888 w(1,4) = w(1,4) - 0.00001024d0/r2as
889
890 w(2,2) = w(2,2) + 0.00470602d0/r2as
891 w(2,3) = w(2,3) - 0.00025213d0/r2as
892
893 w(3,2) = w(3,2) - 0.00261070d0/r2as
894 w(3,3) = w(3,3) - 0.00010712d0/r2as
895 end if
896
897 ! Corrections to the mean motions of the Moon angles W2 and W3, infered from the modifications of the constants:
898 x2 = w(2,1)/w(1,1)
899 x3 = w(3,1)/w(1,1)
900 y2 = am*bp(1,1)+xa*bp(5,1)
901 y3 = am*bp(1,2)+xa*bp(5,2)
902
903 d21 = x2-y2
904 d22 = w(1,1)*bp(2,1)
905 d23 = w(1,1)*bp(3,1)
906 d24 = w(1,1)*bp(4,1)
907 d25 = y2/am
908
909 d31 = x3-y3
910 d32 = w(1,1)*bp(2,2)
911 d33 = w(1,1)*bp(3,2)
912 d34 = w(1,1)*bp(4,2)
913 d35 = y3/am
914
915 cw2_1 = d21*dw1_1+d25*deart_1+d22*dgam+d23*de+d24*dep
916 cw3_1 = d31*dw1_1+d35*deart_1+d32*dgam+d33*de+d34*dep
917
918 w(2,1) = w(2,1)+cw2_1/r2as
919 w(3,1) = w(3,1)+cw3_1/r2as
920
921 ! Arguments of Delaunay:
922 do id=0,4
923 del(1,id) = w(1,id) - eart(id) ! D = W1 - Te + 180 degrees
924 del(2,id) = w(1,id) - w(3,id) ! F = W1 - W3
925 del(3,id) = w(1,id) - w(2,id) ! l = W1 - W2 mean anomaly of the Moon
926 del(4,id) = eart(id) - peri(id) ! l' = Te - Pip mean anomaly of EMB
927 end do
928 del(1,0) = del(1,0) + pi
929
930 ! Planetary arguments: mean longitudes for J2000 (from VSOP2000):
931 p(1,0) = elp_dms2rad(252, 15, 3.216919d0) ! Mercury
932 p(2,0) = elp_dms2rad(181, 58, 44.758419d0) ! Venus
933 p(3,0) = elp_dms2rad(100, 27, 59.138850d0) ! EMB (eart(0))
934 p(4,0) = elp_dms2rad(355, 26, 3.642778d0) ! Mars
935 p(5,0) = elp_dms2rad( 34, 21, 5.379392d0) ! Jupiter
936 p(6,0) = elp_dms2rad( 50, 4, 38.902495d0) ! Saturn
937 p(7,0) = elp_dms2rad(314, 3, 4.354234d0) ! Uranus
938 p(8,0) = elp_dms2rad(304, 20, 56.808371d0) ! Neptune
939
940 ! Planetary arguments: mean motions (from VSOP2000):
941 p(1,1) = 538101628.66888d0/r2as ! Mercury
942 p(2,1) = 210664136.45777d0/r2as ! Venus
943 p(3,1) = 129597742.29300d0/r2as ! EMB (eart(1))
944 p(4,1) = 68905077.65936d0/r2as ! Mars
945 p(5,1) = 10925660.57335d0/r2as ! Jupiter
946 p(6,1) = 4399609.33632d0/r2as ! Saturn
947 p(7,1) = 1542482.57845d0/r2as ! Uranus
948 p(8,1) = 786547.89700d0/r2as ! Neptune
949
950 p(1:8,2:4) = 0.d0
951
952
953 ! Zeta: Mean longitude of the Moon W1 + Rate of precession (pt):
954 zeta(0) = w(1,0)
955 zeta(1) = w(1,1) + (5029.0966d0+dprec)/r2as ! Include precession
956 zeta(2) = w(1,2)
957 zeta(3) = w(1,3)
958 zeta(4) = w(1,4)
959
960 ! Corrections to the parameters: Nu, E, Gamma, n' et e' (Source: ELP):
961 delnu = (+0.55604d0+dw1_1)/r2as/w(1,1) ! Correction to the mean motion of the Moon
962 dele = (+0.01789d0+de)/r2as ! Correction to the half coefficient of sin(l) in longitude
963 delg = (-0.08066d0+dgam)/r2as ! Correction to the half coefficient of sin(F) in latitude
964 delnp = (-0.06424d0+deart_1)/r2as/w(1,1) ! Correction to the mean motion of EMB
965 delep = (-0.12879d0+dep)/r2as ! Correction to the eccentricity of EMB
966
967 ! Precession of the longitude of the ascending node of the mean ecliptic of date on fixed ecliptic J2000 (Laskar, 1986):
968 ! P: sine coefficients:
969 p1 = 0.10180391d-04
970 p2 = 0.47020439d-06
971 p3 = -0.5417367d-09
972 p4 = -0.2507948d-11
973 p5 = 0.463486d-14
974
975 ! Q: cosine coefficients:
976 q1 = -0.113469002d-03
977 q2 = 0.12372674d-06
978 q3 = 0.1265417d-08
979 q4 = -0.1371808d-11
980 q5 = -0.320334d-14
981
982 end subroutine elp_mpp02_initialise
983 !*********************************************************************************************************************************
984
985
986
987 !*********************************************************************************************************************************
988 !> \brief Read the six data files containing the ELP/MPP02 series
989 !!
990 !! \param ierr File error index: ierr=0: no error, ierr=1: file error (output)
991 !!
992 !! \note
993 !! - module elp_mpp02_constants: Set of the constants of ELP/MPP02 solution (input)
994 !! - module elp_mpp02_series: Series of the ELP/MPP02 solution (output)
995 !!
996 !! \see
997 !! - Data files: ftp://cyrano-se.obspm.fr/pub/2_lunar_solutions/2_elpmpp02/
998
999 subroutine elp_mpp02_read_files(ierr)
1000 use sufr_kinds, only: double
1001 use sufr_constants, only: pi2,pio2
1002 use sufr_system, only: find_free_io_unit, file_open_error_quit
1004
1005 use thesky_elp_mpp02_series, only: cmpb,fmpb,nmpb, cper,fper,nper
1006 use thesky_elp_mpp02_constants, only: zeta,del, p,delnu,dele,delg,delnp,delep,dtasm,am
1007
1008 implicit none
1009 integer, intent(out) :: ierr
1010
1011 integer :: ip, i, ilu(4),ifi(16), icount,ipt,ir,it,iFile, k,iLine,nerr
1012 real(double) :: a,b(5),c, pha,s,tgv
1013 logical :: fexist
1014 character :: filename*(199)
1015
1016 ! Read the Main Problem series:
1017 ir=0
1018 nmpb=0
1019 ierr=1
1020
1021 ! Name of the (here single) ELPMPP02 file:
1022 filename = trim(theskydatadir)//'elp_mpp02.dat'
1023 inquire(file=trim(filename), exist=fexist)
1024 if(.not.fexist) call file_open_error_quit(trim(filename), 1, 1) ! 1: input file
1025
1026
1027 call find_free_io_unit(ip)
1028 open(ip,file=trim(filename),status='old',iostat=nerr)
1029
1030 do ifile=1,3 ! These used to be three files
1031 ierr=3
1032 read(ip,'(25x,I10)', iostat=nerr,end=100) nmpb(iFile,1)
1033 if(nerr.ne.0) return
1034
1035 nmpb(ifile,2) = ir+1
1036 nmpb(ifile,3) = nmpb(ifile,1) + nmpb(ifile,2) - 1
1037
1038 do iline=1,nmpb(ifile,1)
1039 ierr=4
1040 read(ip,'(4I3,2x,F13.5,5F12.2)', iostat=nerr,end=100) ilu,a,b
1041 if(nerr.ne.0) return
1042
1043 ir=ir+1
1044 tgv = b(1) + dtasm*b(5)
1045 if(ifile.eq.3) a = a - 2*a*delnu/3.d0
1046 cmpb(ir) = a + tgv*(delnp-am*delnu) + b(2)*delg + b(3)*dele + b(4)*delep
1047
1048 do k=0,4
1049 fmpb(k,ir)=0.d0
1050 do i=1,4
1051 fmpb(k,ir) = fmpb(k,ir) + ilu(i)*del(i,k)
1052 end do
1053 end do ! k
1054
1055 if(ifile.eq.3) fmpb(0,ir)=fmpb(0,ir)+pio2
1056 end do ! iLine
1057
1058 end do ! iFile
1059
1060
1061 ! Read the Perturbations series:
1062 ir=0
1063 nper = 0
1064
1065 do ifile=1,3 ! These used to be three files
1066 do it=0,3
1067 read(ip,'(25x,2I10)', iostat=nerr,end=100) nper(iFile,it,1),ipt
1068 ierr = 6
1069 if(nerr.ne.0) return
1070
1071 nper(ifile,it,2) = ir+1
1072 nper(ifile,it,3) = nper(ifile,it,1) + nper(ifile,it,2) - 1
1073 if(nper(ifile,it,1).eq.0) cycle
1074
1075 do iline=1,nper(ifile,it,1)
1076 read(ip,'(I5,2D20.13,16I3)', iostat=nerr,end=100) icount,s,c,ifi
1077 ierr = 7
1078 if(nerr.ne.0) return
1079
1080 ir = ir+1
1081 cper(ir) = sqrt(c**2+s**2)
1082 pha = atan2(c,s)
1083 if(pha.lt.0.d0) pha = pha+pi2
1084
1085 do k=0,4
1086 fper(k,ir)=0.d0
1087 if(k.eq.0) fper(k,ir) = pha
1088 do i=1,4
1089 fper(k,ir) = fper(k,ir) + ifi(i)*del(i,k)
1090 end do
1091 do i=5,12
1092 fper(k,ir) = fper(k,ir) + ifi(i)*p(i-4,k)
1093 end do
1094 fper(k,ir) = fper(k,ir) + ifi(13)*zeta(k)
1095 end do ! k
1096
1097 end do ! iLine
1098
1099 end do ! it
1100 end do ! iFile
1101
1102 close(ip)
1103
1104 ! Exit:
1105 ierr=0
1106 return
1107
1108 ! End of file error:
1109100 continue
1110 ierr=9
1111
1112 icount = icount ! Suppress 'variable set but not used' compiler warnings
1113 ipt = ipt ! Suppress 'variable set but not used' compiler warnings
1114
1115 end subroutine elp_mpp02_read_files
1116 !*********************************************************************************************************************************
1117
1118
1119 !*********************************************************************************************************************************
1120 !> \brief Function for the conversion: sexagesimal degrees -> radians
1121 !!
1122 !! \param deg Degrees
1123 !! \param min Minutes or arc
1124 !! \param sec Seconds of arc
1125 !!
1126 !!
1127 !! \retval elp_dms2rad The angle in radians
1128
1129 function elp_dms2rad(deg,min,sec)
1130 use sufr_kinds, only: double
1131 use sufr_constants, only: d2r
1132
1133 implicit none
1134 integer, intent(in) :: deg, min
1135 real(double), intent(in) :: sec
1136 real(double) :: elp_dms2rad
1137
1138 elp_dms2rad = (deg+min/60.d0+sec/3600.d0) * d2r
1139 end function elp_dms2rad
1140 !*********************************************************************************************************************************
1141
1142
1143
1144 !*********************************************************************************************************************************
1145 !> \brief Read orbital-element data for the asteroids
1146 !!
1147 !! \note
1148 !! - data are passed via the module planetdata, in asterElems(*,1:9):
1149 !! - Epoch (JD), a, e, i, omega, Omega, M, H, G, for J2000.0
1150 !!
1151 !! \see
1152 !! - asteroids.dat: http://sf.net/projects/libthesky/files/asteroids.dat.bz2 (selection: H<15, a<100 - 25% of all bodies)
1153 !! - original: http://ssd.jpl.nasa.gov/dat/ELEMENTS.NUMBR.gz
1154
1156 use sufr_constants, only: d2r
1157 use sufr_system, only: find_free_io_unit, file_read_end_error
1158 use sufr_dummy, only: dumstr9
1161
1162 implicit none
1163 integer :: i,epoch, ip,status
1164 character :: infile*(199)
1165
1166 call find_free_io_unit(ip)
1167 infile = trim(theskydatadir)//'asteroids.dat' ! (Reduced) copy of http://ssd.jpl.nasa.gov/dat/ELEMENTS.NUMBR.gz
1168 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
1169 if(status.ne.0) then
1170 write(0,'(/,A)') ' I could not open the file '//trim(infile)//' - you will not be able to compute asteroid data.'
1171 write(0,'(A,/)') ' Download http://sf.net/projects/libthesky/files/asteroids.dat.bz2 - unzip it and move it to '// &
1172 trim(infile)//' to fix this.'
1173 return
1174 end if
1175
1176 read(ip,*) dumstr9
1177 read(ip,*) dumstr9
1178
1179 do i=1,nasteroids
1180 read(ip,'(7x,A18,I5,F11.7,F11.8,3F10.5,F12.7,F6.2,F5.2)', iostat=status) asternames(i),epoch,asterelems(i,2:9)
1181 if(status.ne.0) call file_read_end_error(trim(infile), 2+i, status, 1, 1) ! stopcode=1, exitstatus=1
1182
1183 if(status.ne.0) then
1184 call file_read_end_error(trim(infile), 0, status, 0, 0) ! line=0, stopcode=0, exitstatus=0
1185 exit
1186 end if
1187
1188 asterelems(i,1) = dble(epoch) + 2400000.5d0
1189 end do
1190
1191 close(ip)
1192
1193 asterelems(1:nasteroids,4:7) = asterelems(1:nasteroids,4:7)*d2r
1194
1195 end subroutine readasteroidelements
1196 !*********************************************************************************************************************************
1197
1198
1199 !*********************************************************************************************************************************
1200 !> \brief Read orbital-element data for the comets
1201 !!
1202 !! \note
1203 !! - data are passed via the module comet_data, in cometElems(*,1:9)
1204 !! - colums: 1: Epoch (JD), 2: q, 3: e, 4: i, 5: omega, 6: Omega, 7: Tp, for J2000.0, 8: H, 9: G
1205 !!
1206 !! \see
1207 !! - comets.dat: http://ssd.jpl.nasa.gov/dat/ELEMENTS.COMET
1208 !! - comets_mpc.dat: http://www.minorplanetcenter.net/iau/Ephemerides/Comets/Soft00Cmt.txt
1209
1211 use sufr_kinds, only: double
1212 use sufr_constants, only: d2r, currentjd
1213 use sufr_system, only: find_free_io_unit, quit_program_error, file_open_error_quit, file_read_error_quit
1214 use sufr_date_and_time, only: cal2jd
1215 use sufr_dummy, only: dumstr9
1216
1219
1220 implicit none
1221 integer :: ci,epoch, yr,mnt,idy, status, ip
1222 real(double) :: dy
1223 character :: tpstr*(14),epstr*(8), infile*(199)
1224
1225 cometdatfile = 2 ! 1: comets.dat (MANY comets, no magnitude info), 2: comets_mpc.dat (currently visible comets + magn. info)
1226
1227 cometelems = 0.d0
1228 cometdiedatp = .false. ! By default, comets do not die at perihelion - set to .true. for a given comet to change this
1229 ncomets = 0
1230
1231 call find_free_io_unit(ip)
1232 select case(cometdatfile)
1233 case(1) ! comets.dat
1234 infile = trim(theskydatadir)//'comets.dat' ! Copy of http://ssd.jpl.nasa.gov/dat/ELEMENTS.COMET
1235 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
1236 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
1237
1238 read(ip,*) dumstr9
1239 read(ip,*) dumstr9
1240
1241 do ci=11,ncometsmax
1242 read(ip,'(A43,I8,F12.8,F11.8,3F10.5,1x,A14)', iostat=status) cometnames(ci), epoch, cometelems(ci,2), cometelems(ci,3), &
1243 cometelems(ci,4), cometelems(ci,5), cometelems(ci,6), tpstr
1244
1245 if(status.lt.0) exit
1246 if(status.gt.0) call file_read_error_quit(trim(infile), ci-10+2, 1)
1247
1248 ncomets = ncomets + 1
1249 cometelems(ci,1) = dble(epoch) + 2400000.5d0 ! JD of epoch
1250 read(tpstr,'(I4,I2,F8.5)', iostat=status) yr,mnt,dy
1251 if(status.ne.0) call quit_program_error('readCometElements(): Error reading perihelion date from string', 1)
1252 cometelems(ci,7) = cal2jd(yr,mnt,dy) ! JD of perihelion
1253 end do
1254
1255 close(ip)
1256
1257 case(2) ! comets_mpc.dat
1258 infile = trim(theskydatadir)//'comets_mpc.dat' ! Copy of http://www.minorplanetcenter.net/iau/Ephemerides/Comets/Soft00Cmt.txt
1259 open(ip, form='formatted', status='old', action='read', file=trim(infile), iostat=status)
1260 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
1261
1262 do ci=11,ncometsmax
1263 read(ip,'(A12, I6,I3,F8.4, 2F10.6,3F10.4, 2x,A8, F6.1,F5.1,2x,A56)', iostat=status) dumstr9, yr,mnt,dy, &
1264 cometelems(ci,2:3), cometelems(ci,5:6),cometelems(ci,4), epstr, cometelems(ci,8:9), cometnames(ci)
1265
1266 if(status.lt.0) exit
1267 if(status.gt.0) call file_read_error_quit(trim(infile), ci-10, 1)
1268
1269 ncomets = ncomets + 1
1270 cometelems(ci,7) = cal2jd(yr,mnt,dy) ! JD of perihelion
1271
1272 read(epstr,'(I4,I2,I2)', iostat=status) yr,mnt,idy
1273 if(status.ne.0) call quit_program_error('readCometElements(): Error reading epoch from string', 1)
1274 cometelems(ci,1) = cal2jd(yr,mnt,dble(idy)) ! JD of epoch
1275 if(len_trim(epstr).ne.8) cometelems(ci,1) = currentjd ! epstr is sometimes empty
1276 end do
1277
1278 close(ip)
1279
1280 case default
1281 call quit_program_error('readCometElements(): non-exisiting value of cometDatFile',0)
1282
1283 end select
1284
1285
1286 cometelems(1:ncomets,4:6) = cometelems(1:ncomets,4:6)*d2r ! omega, Omega, i
1287
1288 end subroutine readcometelements
1289 !*********************************************************************************************************************************
1290
1291
1292
1293
1294
1295
1296 !*********************************************************************************************************************************
1297 !> \brief Read the Bright Star Catalogue
1298
1299 subroutine read_bsc
1300 use sufr_system, only: find_free_io_unit, file_open_error_quit, file_read_error_quit, file_read_end_error
1301 use sufr_sorting, only: sorted_index_list
1305
1306
1307 implicit none
1308 integer :: i,rv,nsn,snr, ip,status
1309 character :: snam*(10), infile*(199)
1310
1311 call find_free_io_unit(ip)
1312 infile = trim(theskydatadir)//'bsc_data.dat'
1313 open(ip, action='read', form='formatted', status='old', file=trim(infile), iostat=status)
1314 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
1315
1316 do i=1,n_bsc
1317 read(ip,'(A10,1x,2F10.6,1x,2F7.3,I5,F6.2,F6.3,A2,A10,2x,A20,I8,1x,3F6.2)', iostat=status) &
1318 bsc_abbr(i),bsc_ra(i),bsc_dec(i),bsc_pma(i),bsc_pmd(i),rv,bsc_vm(i),bsc_par(i),bsc_mult(i),bsc_var(i), &
1319 bsc_sptype(i),bsc_sao(i),bsc_bv(i),bsc_ub(i),bsc_ri(i)
1320
1321 if(status.ne.0) call file_read_end_error(trim(infile), i, status, 1, 1) ! stopcode=1, exitstatus=1
1322
1323 if(abs(bsc_vm(i)).lt.1.d-5) bsc_vm(i) = 99.d0 ! There are some non-stellar objects in the catalogue, with vm=0.0
1324 bsc_rv(i) = dble(rv)
1325 end do
1326 close(ip)
1327
1328 ! Create an index sorted to visual magnitude:
1329 call sorted_index_list(bsc_vm(1:n_bsc),bsc_vm_indx(1:n_bsc)) ! size n_bsc
1330
1331
1332 ! Star names:
1333 bsc_name = ' '
1334 nsn = 79 ! Number of names in bsc_names.dat: actually 80, but Polaris is listed twice
1335 infile = trim(theskydatadir)//'bsc_names.dat'
1336 open(ip, action='read', form='formatted', status='old', file=trim(infile), iostat=status)
1337 if(status.ne.0) call file_open_error_quit(trim(infile), 1, 1) ! 1-input file, 1-exit status
1338
1339 do i=1,nsn
1340 read(ip,'(I4,2x,A10)', iostat=status) snr, snam
1341 if(status.lt.0) exit
1342 if(status.gt.0) call file_read_error_quit(trim(infile), i, 1)
1343
1344 bsc_name(snr) = snam
1345 end do
1346 close(ip)
1347
1348 end subroutine read_bsc
1349 !*********************************************************************************************************************************
1350
1351
1352
1353end module thesky_data
1354!***********************************************************************************************************************************
Data from the Bright Star Catalogue (BSC)
Definition modules.f90:390
integer, parameter n_bsc
Size of the Bright Star Catalogue (BSC)
Definition modules.f90:396
real(double), dimension(n_bsc) bsc_bv
B-V colours of the BSC stars.
Definition modules.f90:408
real(double), dimension(n_bsc) bsc_ub
U-B colours of the BSC stars.
Definition modules.f90:409
real(double), dimension(n_bsc) bsc_dec
Declinations of the BSC stars.
Definition modules.f90:401
real(double), dimension(n_bsc) bsc_par
Parallaxes of the BSC stars.
Definition modules.f90:407
character, dimension(10) bsc_abbr
Abbreviated names/codes for the BSC stars.
Definition modules.f90:413
real(double), dimension(n_bsc) bsc_vm
Visual magnitudes of the BSC stars.
Definition modules.f90:405
character, dimension(n_bsc) bsc_mult
Multiplicity codes for the BSC stars.
Definition modules.f90:414
real(double), dimension(n_bsc) bsc_pmd
Proper motions in declination of the BSC stars.
Definition modules.f90:403
character, dimension(20) bsc_sptype
Spectral types of the BSC stars.
Definition modules.f90:416
character, dimension(10) bsc_name
Proper names of the BSC stars.
Definition modules.f90:412
real(double), dimension(n_bsc) bsc_rv
Radial velocities of the BSC stars.
Definition modules.f90:404
real(double), dimension(n_bsc) bsc_pma
Proper motions in RA of the BSC stars.
Definition modules.f90:402
real(double), dimension(n_bsc) bsc_ra
Right ascensions of the BSC stars.
Definition modules.f90:400
integer, dimension(n_bsc) bsc_sao
SAO numbers of the BSC stars.
Definition modules.f90:397
real(double), dimension(n_bsc) bsc_ri
R-I colours of the BSC stars.
Definition modules.f90:410
integer, dimension(n_bsc) bsc_vm_indx
Index array for sorting to visual magnitude.
Definition modules.f90:398
character, dimension(n_bsc) bsc_var
Variability codes for the BSC stars.
Definition modules.f90:415
Data to compute comet positions.
Definition modules.f90:322
integer, parameter ncometsmax
Size of comet database.
Definition modules.f90:328
logical, dimension(ncometsmax) cometdiedatp
This comet died at perihelion (true/false)
Definition modules.f90:331
real(double), dimension(ncometsmax, 9) cometelems
Orbital elements of the comets: 1: JD of epoch (often J2000), 2: Perihelion distance (AU?...
Definition modules.f90:332
character, dimension(60) cometnames
Names of the comets.
Definition modules.f90:334
integer ncomets
Actual number of comets in database.
Definition modules.f90:329
integer cometdatfile
Data file to use 1: comets.dat (MANY comets, no magnitude info), 2: comets_mpc.dat (currently visible...
Definition modules.f90:330
Constants used in libTheSky.
Definition modules.f90:23
character, dimension(99) library_name
Name of this library.
Definition modules.f90:47
integer thesky_verbosity
Verbosity of libTheSky output.
Definition modules.f90:42
real(double) deltat_forced
Forced value for DeltaT, overriding computation.
Definition modules.f90:35
integer, parameter deltat_nmax
Maximum number of Delta-T measurements. Need ~430 until 2000.
Definition modules.f90:29
real(double) deltat_accel
Acceleration for DeltaT parabola.
Definition modules.f90:32
integer deltat_minyr
Start year of DeltaT measurements.
Definition modules.f90:39
real(double) deltat_change
Change for DeltaT parabola.
Definition modules.f90:33
real(double) jd1820
JD of 1820.0, for DeltaT.
Definition modules.f90:36
character, dimension(99) theskydatadir
Directory containing data files for libTheSky.
Definition modules.f90:46
real(double), dimension(deltat_nmax) deltat_values
Values of DeltaT.
Definition modules.f90:30
integer deltat_n
Actual number of DeltaT measurements.
Definition modules.f90:38
real(double) deltat_0
Zero point for DeltaT parabola.
Definition modules.f90:34
real(double), dimension(9, 63) nutationdat
Data for simple nutation function.
Definition modules.f90:44
real(double), dimension(deltat_nmax) deltat_years
Years for DeltaT values.
Definition modules.f90:31
integer deltat_maxyr
End year of DeltaT measurements.
Definition modules.f90:40
Procedures to set constants and read data files.
Definition data.f90:23
subroutine readcometelements()
Read orbital-element data for the comets.
Definition data.f90:1211
subroutine readmoondata()
Read ELP-82B periodic terms from moondata.dat.
Definition data.f90:416
subroutine readnutation()
Read nutation input files.
Definition data.f90:243
subroutine elp_mpp02_initialise_and_read_files(mode, ierr)
Initialise ELP-MPP02 data and read the data files.
Definition data.f90:736
subroutine readplanetelements()
Read the data needed to compute the orbital elements of the planets.
Definition data.f90:367
subroutine readconstellations()
Read constellation names: abbreviation, Latin, genitive, Dutch and English.
Definition data.f90:188
subroutine elp_mpp02_read_files(ierr)
Read the six data files containing the ELP/MPP02 series.
Definition data.f90:1000
subroutine read_bsc
Read the Bright Star Catalogue.
Definition data.f90:1300
subroutine elp_mpp02_initialise(mode)
Initialization of the constants and parameters used for the evaluation of the ELP/MPP02 series.
Definition data.f90:787
subroutine read_deltat()
Read deltat.dat with historical values for DeltaT.
Definition data.f90:129
subroutine readmoondata_la()
Read low-accuracy moon-position data from Meeus (moon_la.dat)
Definition data.f90:694
real(double) function elp_dms2rad(deg, min, sec)
Function for the conversion: sexagesimal degrees -> radians.
Definition data.f90:1130
subroutine readpluto()
Read periodic terms for the position of Pluto.
Definition data.f90:331
subroutine set_thesky_constants()
Set the initial values of the variables and constants used in this package.
Definition data.f90:35
subroutine readasteroidelements()
Read orbital-element data for the asteroids.
Definition data.f90:1156
subroutine thesky_readdata()
Read input data files.
Definition data.f90:99
subroutine readplanetdata()
Reads VSOP planet data from planets.dat.
Definition data.f90:277
ELP 2000-82B Moon data, needed to compute Moon positions.
Definition modules.f90:202
real(double), dimension(4, 0:4) del
Delaunay's variables (https://en.wikipedia.org/wiki/Orbital_elements#Delaunay_variables)
Definition modules.f90:230
real(double) q5
Precession cosine coefficient for ELP 2000-82B theory.
Definition modules.f90:223
real(double) q3
Precession cosine coefficient for ELP 2000-82B theory.
Definition modules.f90:221
real(double), dimension(3, 0:4) w
Constants for mean longitude.
Definition modules.f90:225
integer, dimension(3, 0:12) nrang
CHECK: Number of terms? in ELP 2000-82B data file.
Definition modules.f90:247
real(double), parameter c1
Constant for ELP 2000-82B theory (arcminutes to degrees)
Definition modules.f90:208
real(double) p3
Precession sine coefficient for ELP 2000-82B theory.
Definition modules.f90:215
real(double) q4
Precession cosine coefficient for ELP 2000-82B theory.
Definition modules.f90:222
real(double), dimension(0:4) peri
Mean longitude of the perihelion of the Earth-Moon barycentre (EMB)
Definition modules.f90:227
real(double), dimension(6, 704) pc3
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:239
real(double), dimension(0:1) zeta
Mean longitude (w) + rate precession (?)
Definition modules.f90:231
real(double) q2
Precession cosine coefficient for ELP 2000-82B theory.
Definition modules.f90:220
real(double), dimension(6, 918) pc2
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:238
real(double) prec0
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:249
real(double), dimension(6, 1023) pc1
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:237
real(double), dimension(3, 19537) per1
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:240
integer ideb
Memorise whether this routine has been run before.
Definition modules.f90:250
real(double), parameter ath
Constant for ELP 2000-82B theory (orbital separation?)
Definition modules.f90:210
real(double), dimension(8, 0:1) p
Planetary arguments: mean longitudes and mean motions.
Definition modules.f90:228
integer, dimension(3, 12) nterm
CHECK: Number of terms? in ELP 2000-82B data file.
Definition modules.f90:246
real(double), dimension(3, 6766) per2
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:241
integer, dimension(4) ilu
CHECK: Coefficients in ELP 2000-82B data file (integer)
Definition modules.f90:244
real(double) p5
Precession sine coefficient for ELP 2000-82B theory.
Definition modules.f90:217
real(double) p2
Precession sine coefficient for ELP 2000-82B theory.
Definition modules.f90:214
real(double) p1
Precession sine coefficient for ELP 2000-82B theory.
Definition modules.f90:213
real(double), dimension(0:4) eart
Earth-Moon barycentre (EMB) elements.
Definition modules.f90:226
integer, dimension(11) ipla
CHECK: Coefficients in ELP 2000-82B data file (integer)
Definition modules.f90:245
real(double), parameter c2
Constant for ELP 2000-82B theory (arcseconds to degrees)
Definition modules.f90:209
real(double), dimension(6) zone
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:236
real(double) q1
Precession cosine coefficient for ELP 2000-82B theory.
Definition modules.f90:219
real(double), dimension(7) coef
CHECK: Coefficients in ELP 2000-82B data file (float)
Definition modules.f90:235
real(double) p4
Precession sine coefficient for ELP 2000-82B theory.
Definition modules.f90:216
real(double), dimension(3, 8924) per3
CHECK: Something in ELP 2000-82B theory.
Definition modules.f90:242
real(double), dimension(3) pre
CHECK: does this actually do anything?
Definition modules.f90:234
Planet data, needed to compute planet positions.
Definition modules.f90:85
real(double), dimension(4, 6827, 10) vsopdat
Periodic terms for VSOP87.
Definition modules.f90:106
integer(long), dimension(3, 60) moonla_lrb
L,B,R data for the low-accuracy (la) Moon position.
Definition modules.f90:94
character, dimension(18) asternames
Names of the asteroids.
Definition modules.f90:114
integer, dimension(43, 2) plub
Constants for the latitude of Pluto.
Definition modules.f90:98
integer, dimension(3, 8) vsopnls
Numbers of lines in the VSOP input files (l,b,r x 8 pl)
Definition modules.f90:103
integer, dimension(43, 2) plur
Constants for the distance of Pluto.
Definition modules.f90:99
real(double), dimension(3, 8) vsoptruncs
Truncuate VSOP87 terms at these accuracies.
Definition modules.f90:107
real(double), dimension(2, 8, 6, 0:3) plelemdata
Data to compute planet orbital elements.
Definition modules.f90:111
integer, dimension(43, 3) pluc
Constants for the periodic terms for the position of Pluto.
Definition modules.f90:96
integer, dimension(43, 2) plul
Constants for the longitude of Pluto.
Definition modules.f90:97
integer, dimension(8, 60) moonla_arg
Arguments for the low-accuracy (la) Moon position.
Definition modules.f90:93
integer, parameter nasteroids
Number of entries in the asteroids array. Nasteroids is actually much larger; look at the first Naste...
Definition modules.f90:91
integer, dimension(0:5, 3, 8) vsopnblk
Line number in the VSOP data where the next block of (Planet, Variable (LBR), Power) starts.
Definition modules.f90:104
real(double), dimension(nasteroids, 9) asterelems
Asteroid orbital elements.
Definition modules.f90:113
Star and basic constellation data.
Definition modules.f90:343
character, dimension(19) latconnames
Latin constellation names.
Definition modules.f90:378
character, dimension(3) conidabr
Abbreviations of the constellation IDs.
Definition modules.f90:376
character, dimension(17) nlconnames
Dutch constellation names.
Definition modules.f90:380
integer, parameter nconid
Number of data points for constellation ID.
Definition modules.f90:368
character, dimension(19) genconnames
Genitives of the Latin constellation names.
Definition modules.f90:379
integer, parameter nconstel
Number of constellations.
Definition modules.f90:367
real(double), dimension(nconid) coniddecl
Constellation lower declination boundary for ID.
Definition modules.f90:373
character, dimension(18) enconnames
English constellation names.
Definition modules.f90:381
real(double), dimension(nconid) conidrau
Constellation uppwer RA boundary for ID.
Definition modules.f90:372
integer, dimension(nconid) conid
Constellation ID.
Definition modules.f90:369
character, dimension(3) conabr
Abbreviations of the constellations.
Definition modules.f90:375
real(double), dimension(nconid) conidral
Constellation lower RA boundary for ID.
Definition modules.f90:371