-
Notifications
You must be signed in to change notification settings - Fork 19
/
minos_bran.f
3184 lines (3146 loc) · 93.8 KB
/
minos_bran.f
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
c
c MINEOS version 1.0 by Guy Masters, John Woodhouse, and Freeman Gilbert
c
c This program is free software; you can redistribute it and/or modify
c it under the terms of the GNU General Public License as published by
c the Free Software Foundation; either version 2 of the License, or
c (at your option) any later version.
c
c This program is distributed in the hope that it will be useful,
c but WITHOUT ANY WARRANTY; without even the implied warranty of
c MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
c GNU General Public License for more details.
c
c You should have received a copy of the GNU General Public License
c along with this program; if not, write to the Free Software
c Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
c
c****************************************************************************
c minos_bran program produces the solution of the seismic normal mode
c eigenvalue-eigenfunction problem.
c
c****************************************************************************
c this program uses one input and two output files
c the input file contains the model
c the output files are 1) a model listing + a summary of mode properties
c and 2) a file for the eigenfunctions
c the program then asks for some control info described below
c
c structure of model file
c card 1 : title (up to 256 chars long)
c card 2 : ifanis,tref,ifdeck (unformatted)
c ifanis=1 for anisotropic model, 0 for isotropic
c tref=ref period(secs) of model for dispersion correction.
c if tref is .le. 0. no correction is made.
c ifdeck=1 for card deck model, 0 for polynomial model.
c *** card deck model ***
c card 3 : n,nic,noc (unformatted)
c n=no of levels,nic=index of solid side of icb,noc=index of
c fluid side of mcb.note that n must be .le. 223.
c card 4-n: r,rho,vpv,vsv,qkappa,qshear,vph,vsh,eta (f8.0,3f9.2,2f9.1,2f9.2,
c f9.5)
c if isotropic vp=vpv,vs=vsv and vph,vsh,eta are unspecified.
c if the q model is not specified then no disp. correction is made.
c s.i. units are used,e.g. radius in metres and velocities in m/s.
c *** polynomial model ***
c card 3 : nreg,nic,noc,rx (unformatted)
c nreg is the number of regions in the model,nic and noc as before
c and rx is the normalising radius for the polynomials.
c rx is given in kms and is usually 6371.
c card 4 : nlay,r1,r2 (unformatted)
c nlay is the number of levels to be used in the region extending
c from radius r1 to r2 (in kms).
c card 5-9(iso) or card 5-12(ani) : coefs (5f9.5)
c 5 sets of coefficients are required for each region of an isotropic
c model and 8 sets for an anisotropic model.each polynomial can be
c up to a quartic in r/rx ( 5 coefficients) and are ordered thusly:
c rho,vpv,vsv,qkappa,qshear,vph,vsh,eta. the coeffs are given in the
c usual mixed seismological units (rho in g/cc,vel in km/s etc.)
c conversion to s.i. is done by the program.
c cards 4-9(iso) or 4-12(ani) are repeated for each region of the
c model
c
c control parameters (from screen)
c eps,wgrav
c eps controls the accuracy of the runge-kutta integration scheme.
c the relative accuracy of an eigen frequency will be 2-3 x eps.
c it also controls the precision with which a root is found and
c the minimum relative separation of two roots with the same angular
c order.it is safe to set eps=1.d-7.
c wgrav is the frequency in millihertz above which gravitational terms
c are neglected-this gives about a factor of 3 increase in speed.
c jcom
c jcom=1 radial modes, 2 toroidal modes, 3 spheroidal modes,
c 4 inner core toroidal modes.
c lmin,lmax,wmin,wmax,nmin,nmax
c lmin - lmax defines the range of angular orders to be computed.
c if jcom=1 this is ignored. wmin - wmax defines the frequency range
c to be computed (in millihertz)
c nmin-nmax specifies the branch numbers to be computed -- nmin=0
c is the fundamental mode
c
c model listing
c this is an ascii file which lists the model and mode properties
c i.e. phase velocity in km/s,frequency,period,group velocity in km/s,
c q and a parameter which is the ratio of kinetic to potential energy
c minus 1 which should be small ( of order eps )if the eigenfunction is
c accurate and if there are enough radial knots in the model to allow
c quodratures to be done accurately. (you will probably see some degradation
c in this parameter for strongly exponential modes such as stoneley modes ).
c
c eigenfunction file
c****** file name "none" will suppress calculation of eigenfunctions
c this is a fixed record length binary file with an entry for each mode
c written as
c write(ioeig) (abuf(i),i=1,nvec)
c where nvec is 5+6*npts for spheroidal modes and 5+2*npts for toroidal
c and radial modes. the first five words of abuf are n,l,frequecy,q and
c group velocity. the rest of abuf contains w(1..npts) and wp(1..npts)
c for toroidal modes and u(1..npts),up(1..npts),v(1..npts),vp(1..npts),
c p(1..npts) and pp(1..npts) for spheroidal modes. these are as in
c woodhouse and dahlen 1978 except that w,wp,v and vp must be divided
c by sqrt(l*(l+1)). the normalisation is such that
c frequency**2 times integral(rho*w*w*r*r) is 1 for toroidal modes
c frequency**2 times integral(rho*(u*u+l*(l+1)*v*v)*r*r) is 1 for
c spheroidal modes
c the model has been normalised such that a density of 5515 mg/m**3 = 1 ;
c pi*g=1 where g is the grav constant (6.6723e-11) and the radius of the
c earth (rn=6371000 m) is 1. these normalizations result in
c acceleration normalisation = pi*g*rhobar*rn
c velocity normalisation = rn*sqrt(pi*g*rhobar)= vn
c frequency normalization = vn/rn
c
c**************************************************************************
character*256 filnam
print *,'input model file:'
read(5,100) filnam
100 format(a256)
c MB added one line below
print *,filnam(1:lnblnk(filnam))
open(7,file=filnam,status='old',form='formatted',iostat=iret)
print *,'output file:'
read(5,100) filnam
c MB added one line below
print *,filnam(1:lnblnk(filnam))
open(8,file=filnam,form='formatted',iostat=iret)
call model(7,8)
close(7)
print *,'eigenfunction file (output):'
read(5,100) filnam
c MB added one line below
print *,filnam(1:lnblnk(filnam))
ifreq=1
if(filnam(1:4).eq.'none') ifreq=0
open(3,file=filnam,form='unformatted',iostat=iret)
call wtable(8,3,ifreq)
close(8)
close(3)
stop
end
subroutine wtable(iout,ioeig,ifreq)
c*** makes up table of frequencies ***
implicit real*8(a-h,o-z)
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
common/shanks/b(46),c(10),dx,step(8),stepf,maxo,in
common/mtab/we(2),de(2),ke(2),wtry,bm
dimension wt(2)
data inss/5/
cmhz=pi/500.d0
stepf=1.d0
print *,'enter eps and wgrav'
read(*,*) eps,wgrav
c MB added one line below
print *,eps,wgrav
eps=max(eps,1.d-12)
eps1=eps
eps2=eps
wgrav=wgrav*cmhz
write(iout,100) eps,eps1,wgrav
100 format(/,'integration precision =',g12.4,' root precision =',
+ g12.4,' gravity cut off =',g12.4,' rad/s',///,6x,'mode',
+ 8x,'phs vel',7x,'w(mhz)',10x,'t(secs)',6x,'grp vel(km/s)',
+ 8x,'q',13x,'raylquo',/)
call steps(eps)
print *,'enter jcom (1=rad;2=tor;3=sph;4=ictor)'
read(*,*) jcom
c MB added one line below
print *,jcom
if(jcom.lt.1.or.jcom.gt.4) jcom=3
print *,'enter lmin,lmax,wmin,wmax,nmin,nmax'
read(*,*) lmin,lmax,wmin,wmax,normin,normax
c MB added one line below
print *,lmin,lmax,wmin,wmax,normin,normax
wmin=max(wmin,0.1d0)
wmin=wmin*cmhz
wmax=wmax*cmhz
if(lmin.le.0) lmin=1
if(jcom.eq.1) then
lmin=0
lmax=0
end if
normin=max(normin,0)
normax=max(normax,normin)
ncall = 0
do 50 nor=normin,normax
wt(1)=wmin
wt(2)=wmax
wtry=0.d0
bm=0.d0
do 10 l=lmin,lmax
if(wt(1).ge.wmax) goto 50
nord=nor
nup=nor
ndn=nor-1
if(l.eq.1) then
nup=ndn
ndn=ndn-1
end if
knsw=1
maxo=inss
fl=l
fl1=fl+1.d0
fl2=fl+fl1
fl3=fl*fl1
sfl3=sqrt(fl3)
we(1)=wt(1)
we(2)=wt(2)
if(wtry.ne.0.d0) we(2)=wtry
call detqn(we(1),ke(1),de(1),0)
if(ke(1).gt.ndn) goto 10
call detqn(we(2),ke(2),de(2),0)
if(ke(2).lt.nup) then
we(2)=wt(2)
call detqn(we(2),ke(2),de(2),0)
if(ke(2).lt.nup) goto 50
end if
c*** bracket this mode ***
wx=0.5d0*(we(1)+we(2))
ktry=0
15 if(ke(1).eq.ndn.and.ke(2).eq.nup) goto 40
ktry=ktry+1
if(ktry.gt.50) goto 10
call detqn(wx,kx,dx,0)
if(kx.le.ndn) then
we(1)=wx
ke(1)=kx
de(1)=dx
else
we(2)=wx
ke(2)=kx
de(2)=dx
end if
wx=0.5d0*(we(1)+we(2))
goto 15
c*** find roots ***
40 knsw=0
maxo=8
c added argument for number of times called (mhr)
ncall = ncall + 1
call rotspl(eps1,wt,iout,ioeig,ifreq,ncall)
10 continue
50 continue
return
end
subroutine rotspl(eps1,wt,iout,ioeig,ifreq,ncall)
c*** find roots by spline interpolation ***
implicit real*8(a-h,o-z)
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
common/mtab/we(2),de(2),ke(2),wtry,bm
dimension x(20),det(20),qx(3,20),wrk(60),wt(*),kchar(4)
character*2 kchar
data tol/1.d-9/,itmax/15/,kchar/' s',' t',' s',' c'/
if(de(1)*de(2).gt.0.d0) return
nord=ke(2)
if(l.eq.1) nord=nord+1
det(1)=de(1)
det(2)=de(2)
x(1)=we(1)
x(2)=we(2)
5 c=x(1)
if(dabs(det(1)).gt.dabs(det(2))) c=x(2)
10 ntry=2
grad=(x(1)-x(2))/(det(1)-det(2))
b=x(1)-det(1)*grad
15 t=dabs(b*eps1)
if(dabs(b-c).lt.t) goto 65
call detqn(b,knt,fb,0)
ind=1
do 20 m=2,ntry
ind=ind+1
20 if(b.lt.x(m)) goto 25
25 ntry=ntry+1
j2=ntry
30 j1=j2-1
x(j2)=x(j1)
det(j2)=det(j1)
if(j1.eq.ind) goto 35
j2=j2-1
goto 30
35 x(ind)=b
det(ind)=fb
idn=0
do 40 m=2,ntry
idn=idn+1
iup=idn+1
40 if(det(idn)*det(iup).le.0.d0) goto 45
45 ind=iup
if(dabs(det(idn)).lt.dabs(det(iup))) ind=idn
c=x(ind)
if(ntry.ge.itmax) goto 60
call dsplin(ntry,x,det,qx,wrk)
del=-det(ind)/qx(1,ind)
50 delx=-det(ind)/(qx(1,ind)+del*qx(2,ind))
if(dabs(delx-del).lt.tol) goto 55
if(del*delx.lt.0.d0) goto 60
del=delx
goto 50
55 b=c+delx
if(b.ge.x(idn).and.b.le.x(iup)) goto 15
60 x(1)=x(idn)
x(2)=x(iup)
det(1)=det(idn)
det(2)=det(iup)
goto 10
c*** write out frequencies ***
65 call detqn(b,knt,fb,ifreq)
tcom=2.d0*pi/b
wmhz=1000.d0/tcom
cvel=b*rn/(l+0.5)/1000.
if(ifreq.eq.1) then
wdiff=(b-wray*wn)/b
gcom=vn*cg/1000.d0
qmod=0.d0
if(qinv.gt.0.d0) qmod=1.d0/qinv
c remove writing resulting to stnd output (mhr)
c print 200,nord,kchar(jcom),l,cvel,wmhz,tcom,gcom,qmod,wdiff
write(iout,200) nord,kchar(jcom),l,cvel,wmhz,tcom,gcom,qmod,
+ wdiff
200 format(i5,a2,i5,6g16.7)
call modout(b,qmod,gcom,ioeig,ncall)
else
cg=5000./vn
c*** we estimate group velocity using previous frequencies
if(bm.ne.0.d0) cg=(b-bm)/wn
gcom=vn*cg/1000.d0
c added print of gcom. qmod & wdiff not computed w/o egnfcns. (mhr)
c remove writing resulting to stnd output (mhr)
c print 200,nord,kchar(jcom),l,cvel,wmhz,tcom,gcom
write(iout,200) nord,kchar(jcom),l,cvel,wmhz,tcom,gcom
end if
wtry=b+2.d0*cg*wn
wt(1)=b
bm=b
return
end
subroutine modout(wcom,qmod,gcom,ioeig,ncall)
implicit real*8(a-h,o-z)
parameter (mk=350)
real*4 abuf,buf,ww,qq,gc
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
common/rindx/nic,noc,nsl,nicp1,nocp1,nslp1,n
common/eifx/a(14,mk),dum(mk)
c common block name changes from buf$ to c_buf (mhr)
common/c_buf/nn,ll,ww,qq,gc,buf(6*mk)
dimension abuf(6*mk+5)
equivalence (nn,abuf)
nn=nord
ll=l
ww=wcom
qq=qmod
gc=gcom
if(jcom.eq.3) then
nvec=6*n+5
do i=1,n
buf(i)=a(1,i)
buf(i+n)=a(2,i)
buf(i+n+n)=a(3,i)
buf(i+3*n)=a(4,i)
buf(i+4*n)=a(5,i)
buf(i+5*n)=a(6,i)
enddo
else
nvec=2*n+5
if(jcom.eq.2) then
do i=1,noc
a(1,i)=0.d0
a(2,i)=0.d0
enddo
end if
do i=1,n
buf(i)=a(1,i)
buf(i+n)=a(2,i)
enddo
end if
c include writing egfcns (mhr)
write(ioeig) (abuf(i),i=1,nvec)
return
end
subroutine model(iin,iout)
implicit real*8(a-h,o-z)
parameter (mk=350)
save
integer*4 ititle(20)
real*8 lcon,ncon,lspl,nspl
common r(mk),fmu(mk),flam(mk),qshear(mk),qkappa(mk),
+ xa2(mk),xlam(mk),rho(mk),qro(3,mk),g(mk),qg(3,mk),
+ fcon(mk),fspl(3,mk),lcon(mk),lspl(3,mk),ncon(mk),
+ nspl(3,mk),ccon(mk),cspl(3,mk),acon(mk),aspl(3,mk)
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
common/eifx/vpv(mk),vph(mk),vsv(mk),vsh(mk),eta(mk),wrk(mk*10)
common/rindx/nic,noc,nsl,nicp1,nocp1,nslp1,n
data bigg,tau/6.6723d-11,1.d3/,rhobar/5515.d0/
pi=3.14159265358979d0
read(iin,100) (ititle(i),i=1,20)
100 format(20a4)
read(iin,*) ifanis,tref,ifdeck
if(ifdeck.eq.0) go to 1000
c*** card deck model ***
read(iin,*) n,nic,noc
read(iin,*) (r(i),rho(i),vpv(i),vsv(i),
+ qkappa(i),qshear(i),vph(i),vsh(i),eta(i),i=1,n)
105 format(f8.0, 3f9.2, 2f9.1, 2f9.2, f9.5)
cc 105 format(f10.1,3f10.3,2f10.1,2f10.3,f10.5)
go to 2000
c*** polynomial model ***
1000 read(iin,*) nreg,nic,noc,rx
rx=rx*tau
n=0
knt=0
jj=5
if(ifanis.ne.0) jj=8
do nn=1,nreg
read(iin,*) nlay,r1,r2
r1=r1*tau
r2=r2*tau
dr=(r2-r1)/float(nlay-1)
do i=1,nlay
n=n+1
r(n)=r1+dr*float(i-1)
enddo
do j=1,jj
read(iin,110) (wrk(i),i=1,5)
110 format(5f9.5)
if(ifanis.eq.2) then
do i=1,nlay
ind=knt+i
rt=r(ind)/rx
val=wrk(1)+rt*(wrk(2)+rt*(wrk(3)+rt*(wrk(4)+rt*wrk(5))))
if(j.eq.1) rho(ind)=val*tau
if(j.eq.2) vph(ind)=val*tau
if(j.eq.3) vsv(ind)=val*tau
if(j.eq.4) qkappa(ind)=val
if(j.eq.5) qshear(ind)=val
if(j.eq.6) vpv(ind)=sqrt(val)*vph(ind)
if(j.eq.7) vsh(ind)=sqrt(val)*vsv(ind)
if(j.eq.8) eta(ind)=val
enddo
else
do i=1,nlay
ind=knt+i
rt=r(ind)/rx
val=wrk(1)+rt*(wrk(2)+rt*(wrk(3)+rt*(wrk(4)+rt*wrk(5))))
if(j.eq.1) rho(ind)=val*tau
if(j.eq.2) vpv(ind)=val*tau
if(j.eq.3) vsv(ind)=val*tau
if(j.eq.4) qkappa(ind)=val
if(j.eq.5) qshear(ind)=val
if(j.eq.6) vph(ind)=val*tau
if(j.eq.7) vsh(ind)=val*tau
if(j.eq.8) eta(ind)=val
enddo
end if
enddo
knt=knt+nlay
enddo
2000 if(ifanis.eq.0) then
do i=1,n
vph(i)=vpv(i)
vsh(i)=vsv(i)
eta(i)=1.d0
enddo
end if
if(iout.ge.0) then
c*** write out model ***
write(iout,900) (ititle(k),k=1,20),tref
900 format(1x,20a4,' ref per =',f6.1,' secs',///,2x,'level',
1 4x,'radius',8x,'rho',9x,'vpv',9x,'vph',9x,'vsv',
2 9x,'vsh',9x,'eta',9x,'qmu ',8x,'qkap',/)
write(iout,905) (i,r(i),rho(i),vpv(i),vph(i),vsv(i),vsh(i),
1 eta(i),qshear(i),qkappa(i),i=1,n)
905 format(3x,i3,f12.1,5f12.2,f12.5,2f12.2)
end if
c*** normalise and spline ***
rn=r(n)
gn=pi*bigg*rhobar*rn
vn2=gn*rn
vn=sqrt(vn2)
wn=vn/rn
do i=1,n
r(i)=r(i)/rn
if(i.gt.1.and.dabs(r(i)-r(i-1)).lt.1.d-7) r(i)=r(i-1)
if(qshear(i).gt.0.d0) qshear(i)=1.d0/qshear(i)
if(qkappa(i).gt.0.d0) qkappa(i)=1.d0/qkappa(i)
rho(i)=rho(i)/rhobar
acon(i)=rho(i)*vph(i)*vph(i)/vn2
ccon(i)=rho(i)*vpv(i)*vpv(i)/vn2
lcon(i)=rho(i)*vsv(i)*vsv(i)/vn2
ncon(i)=rho(i)*vsh(i)*vsh(i)/vn2
fcon(i)=eta(i)*(acon(i)-2.d0*lcon(i))
fmu(i)=(acon(i)+ccon(i)-2.d0*fcon(i)+5.d0*ncon(i)
+ +6.d0*lcon(i))/15.d0
flam(i)=(4.d0*(acon(i)+fcon(i)-ncon(i))+ccon(i))/9.d0
+ -2.d0*fmu(i)/3.d0
rat=4.d0*fmu(i)/(3.d0*(flam(i)+2.d0*fmu(i)))
xlam(i)=((1.d0-rat)*qkappa(i)-.5d0*rat*qshear(i))/(1.d0
+ -1.5d0*rat)
xa2(i)=(1.d0-rat)*qkappa(i)+rat*qshear(i)
enddo
call drspln(1,n,r,rho,qro,wrk)
call grav(g,rho,qro,r,n)
call drspln(1,n,r,g,qg,wrk)
call drspln(1,n,r,fcon,fspl,wrk)
call drspln(1,n,r,lcon,lspl,wrk)
if(ifanis.ne.0) then
call drspln(1,n,r,acon,aspl,wrk)
call drspln(1,n,r,ccon,cspl,wrk)
call drspln(1,n,r,ncon,nspl,wrk)
end if
nsl=n+1
65 nsl=nsl-1
if(vsv(nsl).le.0.d0) go to 65
nicp1=nic+1
nocp1=noc+1
nslp1=nsl+1
tref=0.5d0*tref/pi
return
end
subroutine detqn(wdim,knt,det,ifeif)
c**** supevises the integration of the equations,it returns the value
c**** of the secular determinant as det and the count of zero crossings.
implicit real*8(a-h,o-z)
parameter (mk=350)
save
real*8 lcon,ncon,lspl,nspl
common r(mk),fmu(mk),flam(mk),qshear(mk),qkappa(mk),
+ xa2(mk),xlam(mk),rho(mk),qro(3,mk),g(mk),qg(3,mk),
+ fcon(mk),fspl(3,mk),lcon(mk),lspl(3,mk),ncon(mk),
+ nspl(3,mk),ccon(mk),cspl(3,mk),acon(mk),aspl(3,mk)
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
common/eifx/a(14,mk),dum(mk)
common/rindx/nic,noc,nsl,nicp1,nocp1,nslp1,n
dimension ass(14),vf(mk),zi(4)
iback=0
w=wdim/wn
wsq=w*w
iexp=0
kount=0
kg=0
fct=0.d0
if(tref.gt.0.d0) fct=2.d0*dlog(tref*wdim)/pi
goto (2,3,1,3),jcom
1 if(wdim.le.wgrav) kg=1
nvefm=2+kg*3
nvesm=5+kg*9
call sdepth(wdim,ls)
if(ls.le.2) then
r10=4.5d-4*(fl+.5d0)/wdim
if(r10.lt.r(2)) then
r(1)=r10
g(1)=rho(1)*r(1)*1.333333333333333d0
ls=1
end if
end if
if(ls.le.nic) then
c*** propagate through inner core ***
call spsm(ls,nvesm,ass)
call sprpmn(ls,nic,ass,vf,nvesm,iexp)
r(1)=0.d0
g(1)=0.d0
call sfbm(ass,kg,iback)
end if
if(ls.le.noc) then
c*** propagate through outer core ***
is=max(ls,nicp1)
if(is.eq.ls) call fpsm(ls,nvefm,ass)
call fprpmn(is,noc,ass,vf,nvefm,iexp)
call fsbm(ass,kg,iback)
end if
is=max(ls,nocp1)
if(is.eq.ls) call spsm(ls,nvesm,ass)
c*** propagate through mantle ***
call sprpmn(is,nsl,ass,vf,nvesm,iexp)
if(nsl.eq.n) then
dnorm=a(1,nsl)*a(1,nsl)
do i=2,nvesm
dnorm=dnorm+a(i,nsl)*a(i,nsl)
enddo
det=a(5,nsl)/sqrt(dnorm)
else
call sfbm(ass,kg,iback)
c*** propagate through ocean ***
call fprpmn(nslp1,n,ass,vf,nvefm,iexp)
if(kg.eq.0) then
det=a(2,n)/sqrt(a(1,n)*a(1,n)+a(2,n)*a(2,n))
else
det=a(5,n)/sqrt(a(1,n)**2+a(2,n)**2+a(3,n)**2+
+ a(4,n)**2+a(5,n)**2)
end if
end if
if(ls.gt.noc) det=-det
if(knsw.eq.1) then
if(ls.gt.noc) kount=kount-2
irem=mod(kount,2)
if(irem.eq.0.and.det.lt.0.d0) kount=kount+1
if(irem.ne.0.and.det.gt.0.d0) kount=kount+1
knt=kount
end if
if(ifeif.eq.0) return
c*** this does eigenfunction calculation for spheroidal modes ***
iback=1
jexp=0
nbakf=1+kg*3
nbaks=4+kg*10
do i=1,nbaks
ass(i)=0.d0
enddo
if(n.ne.nsl) then
if(kg.eq.0) then
ass(1)=dsign(1.d0,a(1,n))
else
asi1=a(3,n)*a(3,n)
asi2=a(4,n)*a(4,n)
if(asi2.le.asi1) ass(1)=dsign(1.d0,a(3,n))
if(asi2.gt.asi1) ass(2)=dsign(1.d0,a(2,n))
end if
call fprpmn(n,nslp1,ass,vf,nbakf,jexp)
call fsbm(ass,kg,iback)
else
asi1=a(3,n)*a(3,n)
asi2=a(4,n)*a(4,n)
if(kg.ne.0) then
asi1=asi1+a(12,n)*a(12,n)
asi2=asi2+a(11,n)*a(11,n)
end if
if(asi2.le.asi1) ass(1)=dsign(1.d0,a(3,n))
if(asi2.gt.asi1) ass(2)=dsign(1.d0,a(2,n))
end if
nto=max(ls,nocp1)
call sprpmn(nsl,nto,ass,vf,nbaks,jexp)
if(nto.eq.ls) goto 90
call sfbm(ass,kg,iback)
nto=max(ls,nicp1)
call fprpmn(noc,nto,ass,vf,nbakf,jexp)
if(nto.eq.ls) goto 90
call fsbm(ass,kg,iback)
nto=max(ls,2)
call sprpmn(nic,nto,ass,vf,nbaks,jexp)
90 if(dabs(det).gt.1.d-4.and.ls.le.noc) call remedy(ls)
call eifout(ls)
return
c*** radial modes ***
2 ls=2
call rps(ls,ass)
call rprop(ls,n,ass)
det=a(2,n)/dsqrt(a(1,n)*a(1,n)+a(2,n)*a(2,n))
knt=kount-1
if(ifeif.eq.0) return
a(1,1)=0.d0
a(2,1)=0.d0
do i=ls,n
ff=fcon(i)*(1.d0+xlam(i)*fct)
cc=ccon(i)*(1.d0+xa2(i)*fct)
a(2,i)=(a(2,i)-2.d0*ff*a(1,i)/r(i))/cc
enddo
zi(1)=0.d0
zi(2)=0.d0
zi(3)=0.d0
do i=ls,n
im=i-1
if(r(i).ne.r(im)) call gauslv(r(im),r(i),im,zi,3)
enddo
rnrm=1.d0/(w*dsqrt(zi(1)))
cg=0.d0
qinv=zi(2)/(wsq*zi(1))
wray=dsqrt(zi(3)/zi(1))
do i=2,n
do j=1,2
a(j,i)=a(j,i)*rnrm
enddo
enddo
return
c*** toroidal modes ***
3 if(jcom.eq.2) then
nb=nocp1
n2=nsl
ass(1)=1.d0
ass(2)=0.d0
else
nb=2
a(1,1)=0.d0
a(2,1)=0.d0
n2=nic
end if
q=0.d0
ls=nb
call startl(ls,n2,fmu,ls,q)
if(ls.ne.nocp1) call tps(ls,ass)
call tprop(ls,n2,ass)
det=a(2,n2)/dsqrt(a(1,n2)*a(1,n2)+a(2,n2)*a(2,n2))
if(knsw.eq.1) then
knt=kount-2
if(jcom.ne.4) then
if(l.ne.1) then
knt=knt+1
irem=mod(knt,2)
if(irem.eq.0.and.det.gt.0.d0) knt=knt+1
if(irem.ne.0.and.det.lt.0.d0) knt=knt+1
end if
end if
end if
if(ifeif.eq.0) return
do i=ls,n2
a(2,i)=a(1,i)/r(i)+a(2,i)/(lcon(i)*(1.d0+qshear(i)*fct))
enddo
if(ls.ne.nb) then
do i=nb,ls-1
a(1,i)=0.d0
a(2,i)=0.d0
enddo
end if
do i=1,4
zi(i)=0.d0
enddo
do i=ls,n2
im=i-1
if(r(i).ne.r(im)) call gauslv(r(im),r(i),im,zi,4)
enddo
rnrm=1.d0/(w*dsqrt(zi(1)))
cg=(fl+0.5d0)*zi(2)/(w*zi(1))
qinv=zi(3)/(wsq*zi(1))
wray=dsqrt(zi(4)/zi(1))
do i=ls,n2
do j=1,2
a(j,i)=a(j,i)*rnrm
enddo
enddo
return
end
subroutine sprpmn(jf,jl,f,h,nvesm,iexp)
c*** propagate a minor vector in a solid region from level jf to jl ***
implicit real*8(a-h,o-z)
parameter (mk=350)
save
real*8 lcon,ncon,lspl,nspl
common r(mk),fmu(mk),flam(mk),qshear(mk),qkappa(mk),
+ xa2(mk),xlam(mk),rho(mk),qro(3,mk),g(mk),qg(3,mk),
+ fcon(mk),fspl(3,mk),lcon(mk),lspl(3,mk),ncon(mk),
+ nspl(3,mk),ccon(mk),cspl(3,mk),acon(mk),aspl(3,mk)
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
common/eifx/ar(14,mk),inorm(mk),jjj(mk)
common/shanks/b(46),c(10),dx,step(8),stepf,maxo,in
dimension f(*),h(nvesm,*),s(14),fp(14),rne(6)
data econst/1048576.d0/
maxo1=maxo-1
jud=1
if(jl.lt.jf) jud=-1
y=r(jf)
i=jf
go to 45
10 x=y
y=r(i)
if(y.eq.x) goto 45
iq=min(i,i-jud)
qff=1.d0+xlam(iq)*fct
qll=1.d0+qshear(iq)*fct
qaa=1.d0+xa2(iq)*fct
xi=g(i)/y
vpsq=(flam(i)+2.d0*fmu(i))/rho(i)
vssq=fmu(i)/rho(i)
alfsq=(wsq+4.d0*rho(i)+xi)/vpsq
betasq=wsq/vssq
delsq=sqrt((betasq-alfsq)**2+4.d0*fl3*xi*xi/(vssq*vpsq))
fksq=.5d0*(alfsq+betasq+delsq)-fl3/(x*x)
qt=sqrt(abs(fksq))+sqrt(abs(fksq-delsq))+2.d0/r(iq)
q=(qt+float(kg)*sfl3/x)/stepf
del=float(jud)*step(maxo)/q
dxs=0.d0
15 y=x+del
if(float(jud)*(y-r(i)).gt.0.d0) y=r(i)
dx=y-x
if(dx.ne.dxs) call baylis(q,maxo1)
dxs=dx
do j=1,nvesm
s(j)=f(j)
enddo
do ni=1,in
z=x+c(ni)
call derms(iq,z,f,h(1,ni),0,qff,qll,qaa)
call rkdot(f,s,h,nvesm,ni)
enddo
if(knsw.eq.1) then
call derms(iq,y,f,fp,1,qff,qll,qaa)
call zknt(s,h,f,fp,x,y,1)
end if
x=y
if(y.ne.r(i)) goto 15
45 size=abs(f(1))
do j=2,nvesm
size=max(size,dabs(f(j)))
enddo
55 if(size.lt.1024.d0) goto 65
do j=1,nvesm
f(j)=f(j)/econst
enddo
size=size/econst
iexp=iexp+20
goto 55
65 if(iback.ne.0) then
inorm(i)=inorm(i)+iexp
if(kg.ne.0) then
t1=f(4)+f(8)
t2=t1+f(4)
t1=t1+f(8)
t3=f(8)-f(4)
rne(1)=ar(6,i)*f(10)-ar(14,i)*f(9)+ar(13,i)*t3
1 -ar(1,i)*f(7)-ar(7,i)*f(6)+ar(8,i)*f(5)
2 +ar(12,i)*f(3)-ar(2,i)*f(2)+ar(3,i)*f(1)
rne(2)=ar(6,i)*f(13)+ar(14,i)*t2+ar(13,i)*f(12)
1 -ar(1,i)*f(11)-ar(9,i)*f(6)-ar(7,i)*f(5)
2 +ar(11,i)*f(3)-ar(4,i)*f(2)-ar(2,i)*f(1)
rne(3)=ar(6,i)*f(14)-ar(7,i)*t1-ar(8,i)*f(12)
1 +ar(13,i)*f(11)-ar(9,i)*f(9)+ar(14,i)*f(7)
2 +ar(10,i)*f(3)+ar(11,i)*f(2)+ar(12,i)*f(1)
rne(4)=ar(14,i)*f(14)+ar(7,i)*f(13)+ar(12,i)*f(12)
1 -ar(2,i)*f(11)-ar(9,i)*f(10)-ar(11,i)*t3
2 +ar(4,i)*f(7)+ar(10,i)*f(5)+ar(5,i)*f(1)
rne(5)=ar(13,i)*f(14)+ar(8,i)*f(13)-ar(12,i)*t2
1 -ar(3,i)*f(11)+ar(7,i)*f(10)-ar(11,i)*f(9)
2 -ar(2,i)*f(7)+ar(10,i)*f(6)+ar(5,i)*f(2)
rne(6)=ar(1,i)*f(14)+ar(13,i)*f(13)-ar(2,i)*t1
1 -ar(3,i)*f(12)+ar(14,i)*f(10)-ar(4,i)*f(9)
2 -ar(11,i)*f(6)-ar(12,i)*f(5)+ar(5,i)*f(3)
else
rne(1)=-ar(1,i)*f(3)+ar(2,i)*f(2)-ar(3,i)*f(1)
rne(2)=-ar(1,i)*f(4)+ar(4,i)*f(2)+ar(2,i)*f(1)
rne(3)=-ar(2,i)*f(4)+ar(4,i)*f(3)-ar(5,i)*f(1)
rne(4)=-ar(3,i)*f(4)-ar(2,i)*f(3)-ar(5,i)*f(2)
end if
do jj=1,6
ar(jj,i)=rne(jj)
enddo
else
inorm(i)=iexp
do j=1,nvesm
ar(j,i)=f(j)
enddo
end if
if(i.eq.jl) return
i=i+jud
go to 10
end
subroutine fprpmn(jf,jl,f,h,nvefm,iexp)
c*** propagate the minor vector in a fluid region from level jf to jl ***
implicit real*8(a-h,o-z)
parameter (mk=350)
save
real*8 lcon,ncon,lspl,nspl
common r(mk),fmu(mk),flam(mk),qshear(mk),qkappa(mk),
+ xa2(mk),xlam(mk),rho(mk),qro(3,mk),g(mk),qg(3,mk),
+ fcon(mk),fspl(3,mk),lcon(mk),lspl(3,mk),ncon(mk),
+ nspl(3,mk),ccon(mk),cspl(3,mk),acon(mk),aspl(3,mk)
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
common/eifx/ar(14,mk),inorm(mk),jjj(mk)
common/shanks/b(46),c(10),dx,step(8),stepf,maxo,in
dimension f(*),h(nvefm,1),s(5),fp(5)
data econst/1048576.d0/
if(nvefm.eq.1) then
do i=jl,jf
inorm(i)=inorm(i)+iexp
do j=1,2
ar(j,i)=ar(j,i)*f(1)
enddo
enddo
return
end if
maxo1=maxo-1
jud=1
if(jl.lt.jf) jud=-1
y=r(jf)
i=jf
go to 45
10 x=y
y=r(i)
if(y.eq.x) goto 45
iq=min(i,i-jud)
qff=1.d0+xlam(iq)*fct
xi=g(i)/y
alfsq=(wsq+4.d0*rho(i)+xi-fl3*xi*xi/wsq)*rho(i)/flam(i)
q=(sqrt(abs(alfsq-fl3/(x*x)))+1.d0/r(iq)+float(kg)*sfl3/x)/stepf
del=float(jud)*step(maxo)/q
dxs=0.d0
15 y=x+del
if(float(jud)*(y-r(i)).gt.0.d0) y=r(i)
dx=y-x
if(dx.ne.dxs) call baylis(q,maxo1)
dxs=dx
do j=1,nvefm
s(j)=f(j)
enddo
do ni=1,in
z=x+c(ni)
call dermf(iq,z,f,h(1,ni),0,qff)
call rkdot(f,s,h,nvefm,ni)
enddo
if(knsw.eq.1) then
call dermf(iq,y,f,fp,1,qff)
call zknt(s,h,f,fp,x,y,0)
end if
x=y
if(y.ne.r(i)) go to 15
45 size=dabs(f(1))
do j=2,nvefm
size=dmax1(size,dabs(f(j)))
enddo
55 if(size.lt.1024.d0) goto 65
do j=1,nvefm
f(j)=f(j)/econst
enddo
size=size/econst
iexp=iexp+20
goto 55
65 if(iback.ne.0) then
inorm(i)=inorm(i)+iexp
rne2 =-ar(1,i)*f(4)+ar(4,i)*f(2)+ar(2,i)*f(1)
ar(1,i)=-ar(1,i)*f(3)+ar(2,i)*f(2)-ar(3,i)*f(1)
rne3 =-ar(2,i)*f(4)+ar(4,i)*f(3)-ar(5,i)*f(1)
ar(4,i)=-ar(3,i)*f(4)-ar(2,i)*f(3)-ar(5,i)*f(2)
ar(2,i)=rne2
ar(3,i)=rne3
else
inorm(i)=iexp
do j=1,nvefm
ar(j,i)=f(j)
enddo
end if
if(i.eq.jl) return
i=i+jud
go to 10
end
subroutine derms(iq,z,f,fp,iknt,qff,qll,qaa)
c*** calculates minor vector derivative (fp) in a solid ***
implicit real*8(a-h,o-z)
parameter (mk=350)
save
real*8 nn,ll,lcon,ncon,lspl,nspl
common r(mk),fmu(mk),flam(mk),qshear(mk),qkappa(mk),
+ xa2(mk),xlam(mk),rho(mk),qro(3,mk),g(mk),qg(3,mk),
+ fcon(mk),fspl(3,mk),lcon(mk),lspl(3,mk),ncon(mk),
+ nspl(3,mk),ccon(mk),cspl(3,mk),acon(mk),aspl(3,mk)
common/bits/pi,rn,vn,wn,w,wsq,wray,qinv,cg,wgrav,tref,fct,eps,fl,
+ fl1,fl2,fl3,sfl3,jcom,nord,l,kg,kount,knsw,ifanis,iback
dimension f(*),fp(*)
if(iknt.ne.0) then
if(kg.eq.0) then
fp(3)=s22*f(1)-2.d0*t12*f(2)+b33*f(3)+c11*f(5)
fp(4)=-s11*f(1)+2.d0*t21*f(2)-b33*f(4)-c22*f(5)
fp(5)=-2.d0*s12*f(2)+s11*f(3)-s22*f(4)-b11*f(5)
else
fp(3)=s22*f(1)+b32*f(2)+b33*f(3)+c11*f(5)+b313*f(13)
fp(4)=-s11*f(1)+b42*f(2)+b44*f(4)-c22*f(5)+b414*f(14)
fp(5)=b52*f(2)+s11*f(3)-s22*f(4)+b55*f(5)-b313*f(11)
+ +b414*f(12)
end if
return
end if
t=z-r(iq)
if(t.eq.0.d0) then
ro=rho(iq)
gr=g(iq)
ff=fcon(iq)*qff
ll=lcon(iq)*qll
nn=ncon(iq)*qll
cc=ccon(iq)*qaa
aa=acon(iq)*qaa
else
ro=rho(iq)+t*(qro(1,iq)+t*(qro(2,iq)+t*qro(3,iq)))
gr=g(iq)+t*(qg(1,iq)+t*(qg(2,iq)+t*qg(3,iq)))
ff=(fcon(iq)+t*(fspl(1,iq)+t*(fspl(2,iq)+t*fspl(3,iq))))*qff
ll=(lcon(iq)+t*(lspl(1,iq)+t*(lspl(2,iq)+t*lspl(3,iq))))*qll
if(ifanis.eq.0) then
nn=ll
cc=ff+ll+ll
aa=cc
else
nn=(ncon(iq)+t*(nspl(1,iq)+t*(nspl(2,iq)+t*nspl(3,iq))))*qll
cc=(ccon(iq)+t*(cspl(1,iq)+t*(cspl(2,iq)+t*cspl(3,iq))))*qaa
aa=(acon(iq)+t*(aspl(1,iq)+t*(aspl(2,iq)+t*aspl(3,iq))))*qaa
end if
end if
zr=1.d0/z
sfl3z=sfl3*zr
rogr=ro*gr
c11=1.d0/cc
c22=1.d0/ll
dmg=aa-nn-ff*ff*c11
zdmg=zr*dmg
t11=-2.d0*ff*zr*c11+zr