-
Notifications
You must be signed in to change notification settings - Fork 14
/
Copy pathsolid.for
1774 lines (1441 loc) · 55.5 KB
/
solid.for
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
*** Solid Earth tides (SET) prediction (section 7.1.2 in the 2003 IERS Conventions).
*** The used syntax is mainly f77 with MIL-STD-1753 extension.
*** Author: Dennis Milbert, 2018-06-01. The code is available at:
*** http://geodesyworld.github.io/SOFTS/solid.htm and can be downloaded as:
*** wget http://geodesyworld.github.io/SOFTS/solid.for.txt -O solid.for
*** The code is based on dehanttideinel.f provided by V. Dehant, S. Mathews,
*** J. Gipson and C. Bruyninx. The latest version of dehanttideinel.f and its
*** dependent subroutines can be download from IERS conventions website as:
*** wget -r -l1 --no-parent -R "index.html*" -nH --cut-dirs=3 https://iers-conventions.obspm.fr/content/chapter7/software/dehanttideinel
*** Sep 2020: modify solid() to solid_point/grid() as subroutines, Z. Yunjun and S. Sangha.
*** Apr 2023: return numpy arrays instead of writing txt file, S. Staniewicz.
subroutine solid_grid(iyr,imo,idy,ihh,imm,iss,
* glad0,steplat,nlat,glod0,steplon,nlon,tide_e,tide_n,tide_u)
*** calculate solid earth tides (SET) for one spatial grid given the date/time
*** Arguments: iyr/imo/idy/ihh/imm/iss - int, date/time for YYYY/MM/DD/HH/MM/SS
*** glad0/glad1/steplat - float, north(Y_FIRST)/south/step(negative) in deg
*** glod0/glod1/steplon - float, west (X_FIRST)/east /step(positive) in deg
*** Returns: tide_e/tide_n/tide_u - 2D np.ndarray, east/north/up component of SET in m
implicit double precision(a-h,o-z)
dimension rsun(3),rmoon(3),etide(3),xsta(3)
integer iyr,imo,idy,ihh,imm,iss
integer nlat,nlon
double precision glad0,steplat
double precision glod0,steplon
real(8), intent(out), dimension(nlat,nlon) :: tide_e
real(8), intent(out), dimension(nlat,nlon) :: tide_n
real(8), intent(out), dimension(nlat,nlon) :: tide_u
!***^ leap second table limit flag
logical lflag
common/stuff/rad,pi,pi2
common/comgrs/a,e2
!f2py intent(in) iyr,imo,idy,ihh,imm,iss,glad0,steplat,nlat,glod0,steplon,nlon
!f2py intent(out) tide_e,tide_n,tide_u
*** constants
pi=4.d0*datan(1.d0)
pi2=pi+pi
rad=180.d0/pi
*** grs80
a=6378137.d0
e2=6.69438002290341574957d-03
*** input section
if(iyr.lt.1901.or.iyr.gt.2099) then
print *, 'ERROR: year NOT in [1901-2099]:',iyr
return
endif
if(imo.lt.1.or.imo.gt.12) then
print *, 'ERROR: month NOT in [1-12]:',imo
return
endif
if(idy.lt.1.or.idy.gt.31) then
print *, 'ERROR: day NOT in [1-31]:',idy
return
endif
if(ihh.lt.0.or.ihh.gt.23) then
print *, 'ERROR: hour NOT in [0-23]:',ihh
return
endif
if(imm.lt.0.or.imm.gt.59) then
print *, 'ERROR: minute NOT in [0-59]:',imm
return
endif
if(iss.lt.0.or.iss.gt.59) then
print *, 'ERROR: second NOT in [0-59]:',iss
return
endif
if(glad0.lt.-90.d0.or.glad0.gt.90.d0) then
print *, 'ERROR: lat0 NOT in [-90,+90]:',glad0
return
endif
if(glod0.lt.-360.d0.or.glod0.gt.360.d0) then
print *, 'ERROR: lon0 NOT in [-360,+360]',glod0
return
endif
glad1=glad0+nlat*steplat
glod1=glod0+nlon*steplon
*** loop over the grid
do ilat=1,nlat
do ilon=1,nlon
glad = glad0 + (ilat-1)*steplat
glod = glod0 + (ilon-1)*steplon
*** position of observing point (positive East)
if(glod.lt. 0.d0) glod=glod+360.d0
if(glod.ge.360.d0) glod=glod-360.d0
gla0=glad/rad
glo0=glod/rad
eht0=0.d0
call geoxyz(gla0,glo0,eht0,x0,y0,z0)
xsta(1)=x0
xsta(2)=y0
xsta(3)=z0
*** here comes the sun (and the moon) (go, tide!)
!***^ UTC time system
ihr= ihh
imn= imm
sec= iss
call civmjd(iyr,imo,idy,ihr,imn,sec,mjd,fmjd)
!***^ normalize civil time
call mjdciv(mjd,fmjd,iyr,imo,idy,ihr,imn,sec)
call setjd0(iyr,imo,idy)
!***^ false means flag not raised
!***^ mjd/fmjd in UTC
!***^ mjd/fmjd in UTC
!***^ mjd/fmjd in UTC
lflag=.false.
call sunxyz (mjd,fmjd,rsun,lflag)
call moonxyz(mjd,fmjd,rmoon,lflag)
call detide (xsta,mjd,fmjd,rsun,rmoon,etide,lflag)
xt = etide(1)
yt = etide(2)
zt = etide(3)
*** determine local geodetic horizon components (topocentric)
!***^ tide vector
call rge(gla0,glo0,ut,vt,wt,xt, yt, zt)
call mjdciv(mjd,fmjd +0.001d0/86400.d0,
* iyr,imo,idy,ihr,imn,sec-0.001d0)
!*** write output respective arrays
tide_e(ilat, ilon) = vt
tide_n(ilat, ilon) = ut
tide_u(ilat, ilon) = wt
enddo
enddo
*** end of processing and flag for leap second
if(lflag) then
print *, 'Mild Warning -- time crossed leap second table'
print *, ' boundaries. Boundary edge value used instead'
endif
return
end
*-----------------------------------------------------------------------
subroutine solid_point(glad,glod,iyr,imo,idy,step_sec,
* secs,tide_e,tide_n,tide_u)
*** calculate SET at given location for one day with step_sec seconds resolution
*** Arguments: glad/glod - float, latitude/longitude in deg
*** iyr/imo/idy - int, start date/time in UTC
*** step_sec - int, time step in seconds
*** Returns: secs - 1D np.ndarray, seconds since start
*** tide_e/tide_n/tide_u - 1D np.ndarray, east/north/up component of SET in m
implicit double precision(a-h,o-z)
dimension rsun(3),rmoon(3),etide(3),xsta(3)
double precision glad,glod
integer iyr,imo,idy
integer nloop, step_sec
double precision tdel2
real(8), intent(out), dimension(60*60*24/step_sec) :: secs
real(8), intent(out), dimension(60*60*24/step_sec) :: tide_e
real(8), intent(out), dimension(60*60*24/step_sec) :: tide_n
real(8), intent(out), dimension(60*60*24/step_sec) :: tide_u
!*** leap second table limit flag
logical lflag
common/stuff/rad,pi,pi2
common/comgrs/a,e2
!f2py intent(in) glad,glod,iyr,imo,idy,step_sec
!f2py intent(out) secs,tide_e,tide_n,tide_u
*** constants
pi=4.d0*datan(1.d0)
pi2=pi+pi
rad=180.d0/pi
*** grs80
a=6378137.d0
e2=6.69438002290341574957d-03
*** check inputs section
if(glad.lt.-90.d0.or.glad.gt.90.d0) then
print *, 'ERROR: lat NOT in [-90,+90]:',glad
return
endif
if(glod.lt.-360.d0.or.glod.gt.360.d0) then
print *, 'ERROR: lon NOT in [-360,+360]:',glod
return
endif
if(iyr.lt.1901.or.iyr.gt.2099) then
print *, 'ERROR: year NOT in [1901-2099]:',iyr
return
endif
if(imo.lt.1.or.imo.gt.12) then
print *, 'ERROR: month NOT in [1-12]:',imo
return
endif
if(idy.lt.1.or.idy.gt.31) then
print *, 'ERROR: day NOT in [1-31]:',idy
return
endif
*** position of observing point (positive East)
if(glod.lt. 0.d0) glod=glod+360.d0
if(glod.ge.360.d0) glod=glod-360.d0
gla0=glad/rad
glo0=glod/rad
eht0=0.d0
call geoxyz(gla0,glo0,eht0,x0,y0,z0)
xsta(1)=x0
xsta(2)=y0
xsta(3)=z0
*** here comes the sun (and the moon) (go, tide!)
!*** UTC time system
ihr= 0
imn= 0
sec=0.d0
call civmjd(iyr,imo,idy,ihr,imn,sec,mjd,fmjd)
!*** normalize civil time
call mjdciv(mjd,fmjd,iyr,imo,idy,ihr,imn,sec)
call setjd0(iyr,imo,idy)
*** loop over time
nloop=60*60*24/step_sec
tdel2=1.d0/DFLOAT(nloop)
do iloop=1,nloop
!*** false means flag not raised
!*** mjd/fmjd in UTC
!*** mjd/fmjd in UTC
!*** mjd/fmjd in UTC
lflag=.false.
call sunxyz (mjd,fmjd,rsun,lflag)
call moonxyz(mjd,fmjd,rmoon,lflag)
call detide (xsta,mjd,fmjd,rsun,rmoon,etide,lflag)
xt = etide(1)
yt = etide(2)
zt = etide(3)
*** determine local geodetic horizon components (topocentric)
!*** tide vector
call rge(gla0,glo0,ut,vt,wt,xt, yt, zt)
call mjdciv(mjd,fmjd,iyr,imo,idy,ihr,imn,sec)
tsec=ihr*3600.d0+imn*60.d0+sec
secs(iloop) = tsec
tide_e(iloop) = vt
tide_n(iloop) = ut
tide_u(iloop) = wt
!*** update fmjd for the next round
fmjd=fmjd+tdel2
!*** force 1 sec. granularity
fmjd=(idnint(fmjd*86400.d0))/86400.d0
enddo
*** end of processing and flag of leap second
if(lflag) then
print *, 'Mild Warning -- time crossed leap second table'
print *, ' boundaries. Boundary edge value used instead'
endif
return
end
*@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
subroutine detide(xsta,mjd,fmjd,xsun,xmon,dxtide,lflag)
*** computation of tidal corrections of station displacements caused
*** by lunar and solar gravitational attraction
*** UTC version
*** step 1 (here general degree 2 and 3 corrections +
*** call st1idiu + call st1isem + call st1l1)
*** + step 2 (call step2diu + call step2lon + call step2idiu)
*** it has been decided that the step 3 un-correction for permanent tide
*** would *not* be applied in order to avoid jump in the reference frame
*** (this step 3 must added in order to get the mean tide station position
*** and to be conformed with the iag resolution.)
*** inputs
*** xsta(i),i=1,2,3 -- geocentric position of the station (ITRF/ECEF)
*** xsun(i),i=1,2,3 -- geoc. position of the sun (ECEF)
*** xmon(i),i=1,2,3 -- geoc. position of the moon (ECEF)
*** mjd,fmjd -- modified julian day (and fraction) (in UTC time)
****old calling sequence*****************************************************
*** dmjd -- time in mean julian date (including day fraction)
*** fhr=hr+zmin/60.+sec/3600. -- hr in the day
*** outputs
*** dxtide(i),i=1,2,3 -- displacement vector (ITRF)
*** lflag -- leap second table limit flag, false:flag not raised
*** author iers 1996 : v. dehant, s. mathews and j. gipson
*** (test between two subroutines)
*** author iers 2000 : v. dehant, c. bruyninx and s. mathews
*** (test in the bernese program by c. bruyninx)
*** created: 96/03/23 (see above)
*** modified from dehanttideinelMJD.f by Dennis Milbert 2006sep10
*** bug fix regarding fhr (changed calling sequence, too)
*** modified to reflect table 7.5a and b IERS Conventions 2003
*** modified to use TT time system to call step 2 functions
*** sign correction by V.Dehant to match eq.16b, p.81, Conventions
*** applied by Dennis Milbert 2007may05
*** UTC version by Dennis Milbert 2018june01
implicit double precision(a-h,o-z)
double precision xsta(3),xsun(3),xmon(3),dxtide(3),xcorsta(3)
double precision h20,l20,h3,l3,h2,l2,fmjd
integer mjd
double precision mass_ratio_sun,mass_ratio_moon
logical lflag,leapflag
!*** leap second table limit flag
!*** leap second table limit flag
save /limitflag/
common/limitflag/leapflag
*** nominal second degree and third degree love numbers and shida numbers
data h20/0.6078d0/,l20/0.0847d0/,h3/0.292d0/,l3/0.015d0/
*** internal support for new calling sequence
*** first, convert UTC time into TT time (and, bring leapflag into variable)
leapflag=lflag
!*** UTC time (sec of day)
!*** TT time (sec of day)
!*** TT time (fract. day)
tsecutc =fmjd*86400.d0
tsectt =utc2ttt(tsecutc)
fmjdtt =tsectt/86400.d0
lflag = leapflag
!*** float MJD in TT
dmjdtt=mjd+fmjdtt
*** commented line was live code in dehanttideinelMJD.f
*** changed on the suggestion of Dr. Don Kim, UNB -- 09mar21
*** Julian date for 2000 January 1 00:00:00.0 UT is JD 2451544.5
*** MJD for 2000 January 1 00:00:00.0 UT is MJD 51544.0
***** t=(dmjdtt-51545.d0)/36525.d0
!*** days to centuries, TT
!*** days to centuries, TT
t=(dmjdtt-51544.d0)/36525.d0
!*** hours in the day, TT
fhr=(dmjdtt-int(dmjdtt))*24.d0
*** scalar product of station vector with sun/moon vector
call sprod(xsta,xsun,scs,rsta,rsun)
call sprod(xsta,xmon,scm,rsta,rmon)
scsun=scs/rsta/rsun
scmon=scm/rsta/rmon
*** computation of new h2 and l2
cosphi=dsqrt(xsta(1)*xsta(1) + xsta(2)*xsta(2))/rsta
h2=h20-0.0006d0*(1.d0-3.d0/2.d0*cosphi*cosphi)
l2=l20+0.0002d0*(1.d0-3.d0/2.d0*cosphi*cosphi)
*** p2-term
p2sun=3.d0*(h2/2.d0-l2)*scsun*scsun-h2/2.d0
p2mon=3.d0*(h2/2.d0-l2)*scmon*scmon-h2/2.d0
*** p3-term
p3sun=5.d0/2.d0*(h3-3.d0*l3)*scsun**3+3.d0/2.d0*(l3-h3)*scsun
p3mon=5.d0/2.d0*(h3-3.d0*l3)*scmon**3+3.d0/2.d0*(l3-h3)*scmon
*** term in direction of sun/moon vector
x2sun=3.d0*l2*scsun
x2mon=3.d0*l2*scmon
x3sun=3.d0*l3/2.d0*(5.d0*scsun*scsun-1.d0)
x3mon=3.d0*l3/2.d0*(5.d0*scmon*scmon-1.d0)
*** factors for sun/moon
mass_ratio_sun=332945.943062d0
mass_ratio_moon=0.012300034d0
re =6378136.55d0
fac2sun=mass_ratio_sun*re*(re/rsun)**3
fac2mon=mass_ratio_moon*re*(re/rmon)**3
fac3sun=fac2sun*(re/rsun)
fac3mon=fac2mon*(re/rmon)
*** total displacement
do i=1,3
dxtide(i)=fac2sun*( x2sun*xsun(i)/rsun + p2sun*xsta(i)/rsta ) +
* fac2mon*( x2mon*xmon(i)/rmon + p2mon*xsta(i)/rsta ) +
* fac3sun*( x3sun*xsun(i)/rsun + p3sun*xsta(i)/rsta ) +
* fac3mon*( x3mon*xmon(i)/rmon + p3mon*xsta(i)/rsta )
enddo
call zero_vec8(xcorsta)
*** corrections for the out-of-phase part of love numbers
*** (part h_2^(0)i and l_2^(0)i )
*** first, for the diurnal band
call st1idiu(xsta,xsun,xmon,fac2sun,fac2mon,xcorsta)
dxtide(1)=dxtide(1)+xcorsta(1)
dxtide(2)=dxtide(2)+xcorsta(2)
dxtide(3)=dxtide(3)+xcorsta(3)
*** second, for the semi-diurnal band
call st1isem(xsta,xsun,xmon,fac2sun,fac2mon,xcorsta)
dxtide(1)=dxtide(1)+xcorsta(1)
dxtide(2)=dxtide(2)+xcorsta(2)
dxtide(3)=dxtide(3)+xcorsta(3)
*** corrections for the latitude dependence of love numbers (part l^(1) )
call st1l1(xsta,xsun,xmon,fac2sun,fac2mon,xcorsta)
dxtide(1)=dxtide(1)+xcorsta(1)
dxtide(2)=dxtide(2)+xcorsta(2)
dxtide(3)=dxtide(3)+xcorsta(3)
*** consider corrections for step 2
*** corrections for the diurnal band:
*** first, we need to know the date converted in julian centuries
*** this is now handled at top of code (also convert to TT time system)
***** t=(dmjd-51545.)/36525.
***** fhr=dmjd-int(dmjd) !*** this is/was a buggy line (day vs. hr)
*** second, the diurnal band corrections,
*** (in-phase and out-of-phase frequency dependence):
call step2diu(xsta,fhr,t,xcorsta)
dxtide(1)=dxtide(1)+xcorsta(1)
dxtide(2)=dxtide(2)+xcorsta(2)
dxtide(3)=dxtide(3)+xcorsta(3)
*** corrections for the long-period band,
*** (in-phase and out-of-phase frequency dependence):
call step2lon(xsta,fhr,t,xcorsta)
dxtide(1)=dxtide(1)+xcorsta(1)
dxtide(2)=dxtide(2)+xcorsta(2)
dxtide(3)=dxtide(3)+xcorsta(3)
*** consider corrections for step 3
*-----------------------------------------------------------------------
* The code below is commented to prevent restoring deformation
* due to permanent tide. All the code above removes
* total tidal deformation with conventional Love numbers.
* The code above realizes a conventional tide free crust (i.e. ITRF).
* This does NOT conform to Resolution 16 of the 18th General Assembly
* of the IAG (1983). This resolution has not been implemented by
* the space geodesy community in general (c.f. IERS Conventions 2003).
*-----------------------------------------------------------------------
*** uncorrect for the permanent tide (only if you want mean tide system)
*** pi=3.141592654
*** sinphi=xsta(3)/rsta
*** cosphi=dsqrt(xsta(1)**2+xsta(2)**2)/rsta
*** cosla=xsta(1)/cosphi/rsta
*** sinla=xsta(2)/cosphi/rsta
*** dr=-dsqrt(5./4./pi)*h2*0.31460*(3./2.*sinphi**2-0.5)
*** dn=-dsqrt(5./4./pi)*l2*0.31460*3.*cosphi*sinphi
*** dxtide(1)=dxtide(1)-dr*cosla*cosphi+dn*cosla*sinphi
*** dxtide(2)=dxtide(2)-dr*sinla*cosphi+dn*sinla*sinphi
*** dxtide(3)=dxtide(3)-dr*sinphi -dn*cosphi
return
end
*-----------------------------------------------------------------------
subroutine st1l1(xsta,xsun,xmon,fac2sun,fac2mon,xcorsta)
*** this subroutine gives the corrections induced by the latitude dependence
*** given by l^(1) in mahtews et al (1991)
*** input: xsta,xsun,xmon,fac3sun,fac3mon
*** output: xcorsta
implicit double precision (a-h,o-z)
double precision xsta(3),xsun(3),xmon(3),xcorsta(3)
double precision l1,l1d,l1sd,fac2sun,fac2mon
data l1d/0.0012d0/,l1sd/0.0024d0/
rsta=enorm8(xsta)
sinphi=xsta(3)/rsta
cosphi=dsqrt(xsta(1)**2+xsta(2)**2)/rsta
sinla=xsta(2)/cosphi/rsta
cosla=xsta(1)/cosphi/rsta
rmon=enorm8(xmon)
rsun=enorm8(xsun)
*** for the diurnal band
l1=l1d
dnsun=-l1*sinphi**2*fac2sun*xsun(3)*(xsun(1)*cosla+xsun(2)*sinla)
* /rsun**2
dnmon=-l1*sinphi**2*fac2mon*xmon(3)*(xmon(1)*cosla+xmon(2)*sinla)
* /rmon**2
desun=l1*sinphi*(cosphi**2-sinphi**2)*fac2sun*xsun(3)*
* (xsun(1)*sinla-xsun(2)*cosla)/rsun**2
demon=l1*sinphi*(cosphi**2-sinphi**2)*fac2mon*xmon(3)*
* (xmon(1)*sinla-xmon(2)*cosla)/rmon**2
de=3.d0*(desun+demon)
dn=3.d0*(dnsun+dnmon)
xcorsta(1)=-de*sinla-dn*sinphi*cosla
xcorsta(2)= de*cosla-dn*sinphi*sinla
xcorsta(3)= dn*cosphi
*** for the semi-diurnal band
l1=l1sd
costwola=cosla**2-sinla**2
sintwola=2.d0*cosla*sinla
dnsun=-l1/2.d0*sinphi*cosphi*fac2sun*((xsun(1)**2-xsun(2)**2)*
* costwola+2.d0*xsun(1)*xsun(2)*sintwola)/rsun**2
dnmon=-l1/2.d0*sinphi*cosphi*fac2mon*((xmon(1)**2-xmon(2)**2)*
* costwola+2.d0*xmon(1)*xmon(2)*sintwola)/rmon**2
desun=-l1/2.d0*sinphi**2*cosphi*fac2sun*((xsun(1)**2-xsun(2)**2)*
* sintwola-2.d0*xsun(1)*xsun(2)*costwola)/rsun**2
demon=-l1/2.d0*sinphi**2*cosphi*fac2mon*((xmon(1)**2-xmon(2)**2)*
* sintwola-2.d0*xmon(1)*xmon(2)*costwola)/rmon**2
de=3.d0*(desun+demon)
dn=3.d0*(dnsun+dnmon)
xcorsta(1)=xcorsta(1)-de*sinla-dn*sinphi*cosla
xcorsta(2)=xcorsta(2)+de*cosla-dn*sinphi*sinla
xcorsta(3)=xcorsta(3) +dn*cosphi
return
end
*-----------------------------------------------------------------------
subroutine step2diu(xsta,fhr,t,xcorsta)
*** last change: vd 17 may 00 1:20 pm
*** these are the subroutines for the step2 of the tidal corrections.
*** they are called to account for the frequency dependence
*** of the love numbers.
implicit double precision (a-h,o-z)
double precision xsta(3),xcorsta(3),datdi(9,31)
double precision deg2rad, fhr, t
data deg2rad/0.017453292519943295769d0/
*** note, following table is derived from dehanttideinelMJD.f (2000oct30 16:10)
*** has minor differences from that of dehanttideinel.f (2000apr17 14:10)
*** D.M. edited to strictly follow published table 7.5a (2006aug08 13:46)
*** cf. table 7.5a of IERS conventions 2003 (TN.32, pg.82)
*** columns are s,h,p,N',ps, dR(ip),dR(op),dT(ip),dT(op)
*** units of mm
data ((datdi(i,j),i=1,9),j=1,31)/
* -3., 0., 2., 0., 0.,-0.01,-0.01, 0.0 , 0.0,
* -3., 2., 0., 0., 0.,-0.01,-0.01, 0.0 , 0.0,
* -2., 0., 1.,-1., 0.,-0.02,-0.01, 0.0 , 0.0,
*****-----------------------------------------------------------------------
****** -2., 0., 1., 0., 0.,-0.08,-0.05, 0.01,-0.02,
!*** original entry
!*** table 7.5a
* -2., 0., 1., 0., 0.,-0.08, 0.00, 0.01, 0.01,
*****-----------------------------------------------------------------------
* -2., 2.,-1., 0., 0.,-0.02,-0.01, 0.0 , 0.0,
*****-----------------------------------------------------------------------
****** -1., 0., 0.,-1., 0.,-0.10,-0.05, 0.0 ,-0.02,
!*** original entry
!*** table 7.5a
* -1., 0., 0.,-1., 0.,-0.10, 0.00, 0.00, 0.00,
*****-----------------------------------------------------------------------
****** -1., 0., 0., 0., 0.,-0.51,-0.26,-0.02,-0.12,
!*** original entry
!*** table 7.5a
* -1., 0., 0., 0., 0.,-0.51, 0.00,-0.02, 0.03,
*****-----------------------------------------------------------------------
* -1., 2., 0., 0., 0., 0.01, 0.0 , 0.0 , 0.0,
* 0.,-2., 1., 0., 0., 0.01, 0.0 , 0.0 , 0.0,
* 0., 0.,-1., 0., 0., 0.02, 0.01, 0.0 , 0.0,
*****-----------------------------------------------------------------------
****** 0., 0., 1., 0., 0., 0.06, 0.02, 0.0 , 0.01,
!*** original entry
!*** table 7.5a
* 0., 0., 1., 0., 0., 0.06, 0.00, 0.00, 0.00,
*****-----------------------------------------------------------------------
* 0., 0., 1., 1., 0., 0.01, 0.0 , 0.0 , 0.0,
* 0., 2.,-1., 0., 0., 0.01, 0.0 , 0.0 , 0.0,
!*** table 7.5a
* 1.,-3., 0., 0., 1.,-0.06, 0.00, 0.00, 0.00,
* 1.,-2., 0., 1., 0., 0.01, 0.0 , 0.0 , 0.0,
*****-----------------------------------------------------------------------
****** 1.,-2., 0., 0., 0.,-1.23,-0.05, 0.06,-0.06,
!*** original entry
!*** table 7.5a
* 1.,-2., 0., 0., 0.,-1.23,-0.07, 0.06, 0.01,
*****-----------------------------------------------------------------------
* 1.,-1., 0., 0.,-1., 0.02, 0.0 , 0.0 , 0.0,
* 1.,-1., 0., 0., 1., 0.04, 0.0 , 0.0 , 0.0,
!*** table 7.5a
* 1., 0., 0.,-1., 0.,-0.22, 0.01, 0.01, 0.00,
*****-----------------------------------------------------------------------
****** 1., 0., 0., 0., 0.,12.02,-0.45,-0.66, 0.17,
!*** original entry
!*** table 7.5a
* 1., 0., 0., 0., 0.,12.00,-0.78,-0.67,-0.03,
*****-----------------------------------------------------------------------
****** 1., 0., 0., 1., 0., 1.73,-0.07,-0.10, 0.02,
!*** original entry
!*** table 7.5a
* 1., 0., 0., 1., 0., 1.73,-0.12,-0.10, 0.00,
*****-----------------------------------------------------------------------
* 1., 0., 0., 2., 0.,-0.04, 0.0 , 0.0 , 0.0,
*****-----------------------------------------------------------------------
****** 1., 1., 0., 0.,-1.,-0.50, 0.0 , 0.03, 0.0,
!*** original entry
!*** table 7.5a
* 1., 1., 0., 0.,-1.,-0.50,-0.01, 0.03, 0.00,
*****-----------------------------------------------------------------------
* 1., 1., 0., 0., 1., 0.01, 0.0 , 0.0 , 0.0,
*****-----------------------------------------------------------------------
****** 0., 1., 0., 1.,-1.,-0.01, 0.0 , 0.0 , 0.0,
!*** original entry
!*** v.dehant 2007
* 1., 1., 0., 1.,-1.,-0.01, 0.0 , 0.0 , 0.0,
*****-----------------------------------------------------------------------
* 1., 2.,-2., 0., 0.,-0.01, 0.0 , 0.0 , 0.0,
*****-----------------------------------------------------------------------
****** 1., 2., 0., 0., 0.,-0.12, 0.01, 0.01, 0.0,
!*** original entry
!*** table 7.5a
* 1., 2., 0., 0., 0.,-0.11, 0.01, 0.01, 0.00,
*****-----------------------------------------------------------------------
* 2.,-2., 1., 0., 0.,-0.01, 0.0 , 0.0 , 0.0,
* 2., 0.,-1., 0., 0.,-0.02, 0.02, 0.0 , 0.01,
* 3., 0., 0., 0., 0., 0.0 , 0.01, 0.0 , 0.01,
* 3., 0., 0., 1., 0., 0.0 , 0.01, 0.0 , 0.0/
s=218.31664563d0+481267.88194d0*t-0.0014663889d0*t*t
* +0.00000185139d0*t**3
tau=fhr*15.d0+280.4606184d0+36000.7700536d0*t+0.00038793d0*t*t
* -0.0000000258d0*t**3-s
pr=1.396971278*t+0.000308889*t*t+0.000000021*t**3
* +0.000000007*t**4
s=s+pr
h=280.46645d0+36000.7697489d0*t+0.00030322222d0*t*t
* +0.000000020*t**3-0.00000000654*t**4
p=83.35324312d0+4069.01363525d0*t-0.01032172222d0*t*t
* -0.0000124991d0*t**3+0.00000005263d0*t**4
zns=234.95544499d0 +1934.13626197d0*t-0.00207561111d0*t*t
* -0.00000213944d0*t**3+0.00000001650d0*t**4
ps=282.93734098d0+1.71945766667d0*t+0.00045688889d0*t*t
* -0.00000001778d0*t**3-0.00000000334d0*t**4
*** reduce angles to between 0 and 360
s= dmod( s,360.d0)
tau=dmod(tau,360.d0)
h= dmod( h,360.d0)
p= dmod( p,360.d0)
zns=dmod(zns,360.d0)
ps= dmod( ps,360.d0)
rsta=dsqrt(xsta(1)**2+xsta(2)**2+xsta(3)**2)
sinphi=xsta(3)/rsta
cosphi=dsqrt(xsta(1)**2+xsta(2)**2)/rsta
cosla=xsta(1)/cosphi/rsta
sinla=xsta(2)/cosphi/rsta
zla = datan2(xsta(2),xsta(1))
do i=1,3
xcorsta(i)=0.d0
enddo
do j=1,31
thetaf=(tau+datdi(1,j)*s+datdi(2,j)*h+datdi(3,j)*p+
* datdi(4,j)*zns+datdi(5,j)*ps)*deg2rad
dr=datdi(6,j)*2.d0*sinphi*cosphi*sin(thetaf+zla)+
* datdi(7,j)*2.d0*sinphi*cosphi*cos(thetaf+zla)
dn=datdi(8,j)*(cosphi**2-sinphi**2)*sin(thetaf+zla)+
* datdi(9,j)*(cosphi**2-sinphi**2)*cos(thetaf+zla)
***** following correction by V.Dehant to match eq.16b, p.81, 2003 Conventions
***** de=datdi(8,j)*sinphi*cos(thetaf+zla)+
de=datdi(8,j)*sinphi*cos(thetaf+zla)-
* datdi(9,j)*sinphi*sin(thetaf+zla)
xcorsta(1)=xcorsta(1)+dr*cosla*cosphi-de*sinla
* -dn*sinphi*cosla
xcorsta(2)=xcorsta(2)+dr*sinla*cosphi+de*cosla
* -dn*sinphi*sinla
xcorsta(3)=xcorsta(3)+dr*sinphi+dn*cosphi
enddo
do i=1,3
xcorsta(i)=xcorsta(i)/1000.d0
enddo
return
end
*-----------------------------------------------------------------------
subroutine step2lon(xsta,fhr,t,xcorsta)
implicit double precision (a-h,o-z)
double precision deg2rad,fhr,t
double precision xsta(3),xcorsta(3),datdi(9,5)
data deg2rad/0.017453292519943295769d0/
*** cf. table 7.5b of IERS conventions 2003 (TN.32, pg.82)
*** columns are s,h,p,N',ps, dR(ip),dT(ip),dR(op),dT(op)
*** IERS cols.= s,h,p,N',ps, dR(ip),dR(op),dT(ip),dT(op)
*** units of mm
data ((datdi(i,j),i=1,9),j=1,5)/
* 0, 0, 0, 1, 0, 0.47, 0.23, 0.16, 0.07,
* 0, 2, 0, 0, 0, -0.20,-0.12,-0.11,-0.05,
* 1, 0,-1, 0, 0, -0.11,-0.08,-0.09,-0.04,
* 2, 0, 0, 0, 0, -0.13,-0.11,-0.15,-0.07,
* 2, 0, 0, 1, 0, -0.05,-0.05,-0.06,-0.03/
s=218.31664563d0+481267.88194d0*t-0.0014663889d0*t*t
* +0.00000185139d0*t**3
pr=1.396971278*t+0.000308889*t*t+0.000000021*t**3
* +0.000000007*t**4
s=s+pr
h=280.46645d0+36000.7697489d0*t+0.00030322222d0*t*t
* +0.000000020*t**3-0.00000000654*t**4
p=83.35324312d0+4069.01363525d0*t-0.01032172222d0*t*t
* -0.0000124991d0*t**3+0.00000005263d0*t**4
zns=234.95544499d0 +1934.13626197d0*t-0.00207561111d0*t*t
* -0.00000213944d0*t**3+0.00000001650d0*t**4
ps=282.93734098d0+1.71945766667d0*t+0.00045688889d0*t*t
* -0.00000001778d0*t**3-0.00000000334d0*t**4
rsta=dsqrt(xsta(1)**2+xsta(2)**2+xsta(3)**2)
sinphi=xsta(3)/rsta
cosphi=dsqrt(xsta(1)**2+xsta(2)**2)/rsta
cosla=xsta(1)/cosphi/rsta
sinla=xsta(2)/cosphi/rsta
*** reduce angles to between 0 and 360
s= dmod( s,360.d0)
***** tau=dmod(tau,360.d0) !*** tau not used here--09jul28
h= dmod( h,360.d0)
p= dmod( p,360.d0)
zns=dmod(zns,360.d0)
ps= dmod( ps,360.d0)
dr_tot=0.d0
dn_tot=0.d0
do i=1,3
xcorsta(i)=0.d0
enddo
*** 1 2 3 4 5 6 7 8 9
*** columns are s,h,p,N',ps, dR(ip),dT(ip),dR(op),dT(op)
do j=1,5
thetaf=(datdi(1,j)*s+datdi(2,j)*h+datdi(3,j)*p+
* datdi(4,j)*zns+datdi(5,j)*ps)*deg2rad
dr=datdi(6,j)*(3.d0*sinphi**2-1.d0)/2.*cos(thetaf)+
* datdi(8,j)*(3.d0*sinphi**2-1.d0)/2.*sin(thetaf)
dn=datdi(7,j)*(cosphi*sinphi*2.d0)*cos(thetaf)+
* datdi(9,j)*(cosphi*sinphi*2.d0)*sin(thetaf)
de=0.d0
dr_tot=dr_tot+dr
dn_tot=dn_tot+dn
xcorsta(1)=xcorsta(1)+dr*cosla*cosphi-de*sinla
* -dn*sinphi*cosla
xcorsta(2)=xcorsta(2)+dr*sinla*cosphi+de*cosla
* -dn*sinphi*sinla
xcorsta(3)=xcorsta(3)+dr*sinphi+dn*cosphi
enddo
do i=1,3
xcorsta(i)=xcorsta(i)/1000.d0
enddo
return
end
*-----------------------------------------------------------------------
subroutine st1idiu(xsta,xsun,xmon,fac2sun,fac2mon,xcorsta)
*** this subroutine gives the out-of-phase corrections induced by
*** mantle inelasticity in the diurnal band
*** input: xsta,xsun,xmon,fac2sun,fac2mon
*** output: xcorsta
implicit double precision (a-h,o-z)
double precision xsta(3),xsun(3),xmon(3),xcorsta(3)
double precision fac2sun,fac2mon
data dhi/-0.0025d0/,dli/-0.0007d0/
rsta=enorm8(xsta)
sinphi=xsta(3)/rsta
cosphi=dsqrt(xsta(1)**2+xsta(2)**2)/rsta
cos2phi=cosphi**2-sinphi**2
sinla=xsta(2)/cosphi/rsta
cosla=xsta(1)/cosphi/rsta
rmon=enorm8(xmon)
rsun=enorm8(xsun)
drsun=-3.d0*dhi*sinphi*cosphi*fac2sun*xsun(3)*(xsun(1)*
* sinla-xsun(2)*cosla)/rsun**2
drmon=-3.d0*dhi*sinphi*cosphi*fac2mon*xmon(3)*(xmon(1)*
* sinla-xmon(2)*cosla)/rmon**2
dnsun=-3.d0*dli*cos2phi*fac2sun*xsun(3)*(xsun(1)*sinla-
* xsun(2)*cosla)/rsun**2
dnmon=-3.d0*dli*cos2phi*fac2mon*xmon(3)*(xmon(1)*sinla-
* xmon(2)*cosla)/rmon**2
desun=-3.d0*dli*sinphi*fac2sun*xsun(3)*
* (xsun(1)*cosla+xsun(2)*sinla)/rsun**2
demon=-3.d0*dli*sinphi*fac2mon*xmon(3)*
* (xmon(1)*cosla+xmon(2)*sinla)/rmon**2
dr=drsun+drmon
dn=dnsun+dnmon
de=desun+demon
xcorsta(1)=dr*cosla*cosphi-de*sinla-dn*sinphi*cosla
xcorsta(2)=dr*sinla*cosphi+de*cosla-dn*sinphi*sinla
xcorsta(3)=dr*sinphi +dn*cosphi
return
end
*-----------------------------------------------------------------------
subroutine st1isem(xsta,xsun,xmon,fac2sun,fac2mon,xcorsta)
*** this subroutine gives the out-of-phase corrections induced by
*** mantle inelasticity in the diurnal band
*** input: xsta,xsun,xmon,fac2sun,fac2mon
*** output: xcorsta
implicit double precision (a-h,o-z)
double precision xsta(3),xsun(3),xmon(3),xcorsta(3)
double precision fac2sun, fac2mon
data dhi/-0.0022d0/,dli/-0.0007d0/
rsta=enorm8(xsta)
sinphi=xsta(3)/rsta
cosphi=dsqrt(xsta(1)**2+xsta(2)**2)/rsta
sinla=xsta(2)/cosphi/rsta
cosla=xsta(1)/cosphi/rsta
costwola=cosla**2-sinla**2
sintwola=2.d0*cosla*sinla
rmon=enorm8(xmon)
rsun=enorm8(xsun)
drsun=-3.d0/4.d0*dhi*cosphi**2*fac2sun*((xsun(1)**2-xsun(2)**2)*
* sintwola-2.*xsun(1)*xsun(2)*costwola)/rsun**2
drmon=-3.d0/4.d0*dhi*cosphi**2*fac2mon*((xmon(1)**2-xmon(2)**2)*
* sintwola-2.*xmon(1)*xmon(2)*costwola)/rmon**2
dnsun=1.5d0*dli*sinphi*cosphi*fac2sun*((xsun(1)**2-xsun(2)**2)*
* sintwola-2.d0*xsun(1)*xsun(2)*costwola)/rsun**2
dnmon=1.5d0*dli*sinphi*cosphi*fac2mon*((xmon(1)**2-xmon(2)**2)*
* sintwola-2.d0*xmon(1)*xmon(2)*costwola)/rmon**2
desun=-3.d0/2.d0*dli*cosphi*fac2sun*((xsun(1)**2-xsun(2)**2)*
* costwola+2.*xsun(1)*xsun(2)*sintwola)/rsun**2
demon=-3.d0/2.d0*dli*cosphi*fac2mon*((xmon(1)**2-xmon(2)**2)*
* costwola+2.d0*xmon(1)*xmon(2)*sintwola)/rmon**2
dr=drsun+drmon
dn=dnsun+dnmon
de=desun+demon
xcorsta(1)=dr*cosla*cosphi-de*sinla-dn*sinphi*cosla
xcorsta(2)=dr*sinla*cosphi+de*cosla-dn*sinphi*sinla
xcorsta(3)=dr*sinphi+dn*cosphi
return
end
*-----------------------------------------------------------------------
subroutine sprod(x,y,scal,r1,r2)
*** computation of the scalar-product of two vectors and their norms
*** input: x(i),i=1,2,3 -- components of vector x
*** y(i),i=1,2,3 -- components of vector y
*** output: scal -- scalar product of x and y
*** r1,r2 -- lengths of the two vectors x and y
implicit double precision (a-h,o-z)
double precision x(3),y(3)
double precision r1,r2,scal
r1=dsqrt(x(1)*x(1) + x(2)*x(2) + x(3)*x(3))
r2=dsqrt(y(1)*y(1) + y(2)*y(2) + y(3)*y(3))
scal= x(1)*y(1) + x(2)*y(2) + x(3)*y(3)
return
end
*-----------------------------------------------------------------------
double precision function enorm8(a)
*** compute euclidian norm of a vector (of length 3)
double precision a(3)
enorm8=dsqrt(a(1)*a(1) + a(2)*a(2) + a(3)*a(3))
return
end
*-----------------------------------------------------------------------
subroutine zero_vec8(v)
*** initialize a vector (of length 3) to zero
double precision v(3)
v(1)=0.d0
v(2)=0.d0
v(3)=0.d0
return
end
*@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
subroutine moonxyz(mjd,fmjd,rm,lflag)
*** get low-precision, geocentric coordinates for moon (ECEF)
*** UTC version
*** input: mjd/fmjd, is Modified Julian Date (and fractional) in UTC time
*** output: rm, is geocentric lunar position vector [m] in ECEF
*** lflag -- leap second table limit flag, false:flag not raised
*** 1."satellite orbits: models, methods, applications" montenbruck & gill(2000)
*** section 3.3.2, pg. 72-73
*** 2."astronomy on the personal computer, 4th ed." montenbruck & pfleger (2005)
*** section 3.2, pg. 38-39 routine MiniMoon
implicit double precision(a-h,o-z)
double precision rm(3)
integer mjd
double precision fmjd
logical lflag,leapflag
!*** leap second table limit flag
!*** leap second table limit flag
save /limitflag/
common/limitflag/leapflag
common/stuff/rad,pi,pi2
*** use TT for lunar ephemerides
leapflag=lflag
!*** UTC time (sec of day)
!*** TT time (sec of day)
!*** TT time (fract. day)
tsecutc=fmjd*86400.d0
tsectt =utc2ttt(tsecutc)
fmjdtt =tsectt/86400.d0
lflag = leapflag
*** julian centuries since 1.5 january 2000 (J2000)
*** (note: also low precision use of mjd --> tjd)
!*** Julian Date, TT
!*** julian centuries, TT
tjdtt = mjd+fmjdtt+2400000.5d0
t = (tjdtt - 2451545.d0)/36525.d0
*** el0 -- mean longitude of Moon (deg)
*** el -- mean anomaly of Moon (deg)
*** elp -- mean anomaly of Sun (deg)
*** f -- mean angular distance of Moon from ascending node (deg)
*** d -- difference between mean longitudes of Sun and Moon (deg)
*** equations 3.47, p.72
el0=218.31617d0 + 481267.88088d0*t -1.3972*t
el =134.96292d0 + 477198.86753d0*t
elp=357.52543d0 + 35999.04944d0*t