-
Notifications
You must be signed in to change notification settings - Fork 15
/
synthetic.f90
2164 lines (2129 loc) · 61.8 KB
/
synthetic.f90
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
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: MODULE
! CODE: FORTRAN 90
! This module declares variable for global use, that is, for
! USE in any subroutine or function or other module.
! Variables whose values are SAVEd can have their most
! recent values reused in any routine.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
MODULE globalp
IMPLICIT NONE
INTEGER, PARAMETER :: i10=SELECTED_REAL_KIND(6)
INTEGER :: checkstat
INTEGER, SAVE :: nvx,nvz,nnx,nnz,fom,gdx,gdz
INTEGER, SAVE :: vnl,vnr,vnt,vnb,nrnx,nrnz,sgdl,rbint
INTEGER, SAVE :: nnxr,nnzr,asgr
INTEGER, DIMENSION (:,:), ALLOCATABLE :: nsts,nstsr,srs
REAL(KIND=i10), SAVE :: gox,goz,dnx,dnz,dvx,dvz,snb,earth
REAL(KIND=i10), SAVE :: goxd,gozd,dvxd,dvzd,dnxd,dnzd
REAL(KIND=i10), SAVE :: drnx,drnz,gorx,gorz
REAL(KIND=i10), SAVE :: dnxr,dnzr,goxr,gozr
REAL(KIND=i10), DIMENSION (:,:), ALLOCATABLE, SAVE :: velv,veln,velnb
REAL(KIND=i10), DIMENSION (:,:), ALLOCATABLE, SAVE :: ttn,ttnr
!REAL(KIND=i10), DIMENSION (:), ALLOCATABLE, SAVE :: rcx,rcz
REAL(KIND=i10), PARAMETER :: pi=3.1415926535898
!!!--------------------------------------------------------------
!! modified by Hongjian Fang @ USTC
! real,dimension(:),allocatable,save::rw
! integer,dimension(:),allocatable,save::iw,col
! real,dimension(:,:,:),allocatable::vpf,vsf
! real,dimension(:),allocatable,save::obst,cbst,wt,dtres
!! integer,dimension(:),allocatable,save::cbst_stat
! real,dimension(:,:,:),allocatable,save::sen_vs,sen_vp,sen_rho
!!! real,dimension(:,:,:),allocatable,save::sen_vsRc,sen_vpRc,sen_rhoRc
!!! real,dimension(:,:,:),allocatable,save::sen_vsRg,sen_vpRg,sen_rhoRg
!!! real,dimension(:,:,:),allocatable,save::sen_vsLc,sen_vpLc,sen_rhoLc
!!! real,dimension(:,:,:),allocatable,save::sen_vsLg,sen_vpLg,sen_rhoLg
!!! integer,save:: count1,count2
! integer*8,save:: nar
! integer,save:: iter,maxiter
!!!--------------------------------------------------------------
!
! nvx,nvz = B-spline vertex values
! dvx,dvz = B-spline vertex separation
! velv(i,j) = velocity values at control points
! nnx,nnz = Number of nodes of grid in x and z
! nnxr,nnzr = Number of nodes of refined grid in x and z
! gox,goz = Origin of grid (theta,phi)
! goxr, gozr = Origin of refined grid (theta,phi)
! dnx,dnz = Node separation of grid in x and z
! dnxr,dnzr = Node separation of refined grid in x and z
! veln(i,j) = velocity values on a refined grid of nodes
! velnb(i,j) = Backup of veln required for source grid refinement
! ttn(i,j) = traveltime field on the refined grid of nodes
! ttnr(i,j) = ttn for refined grid
! nsts(i,j) = node status (-1=far,0=alive,>0=close)
! nstsr(i,j) = nsts for refined grid
! checkstat = check status of memory allocation
! fom = use first-order(0) or mixed-order(1) scheme
! snb = Maximum size of narrow band as fraction of nnx*nnz
! nrc = number of receivers
! rcx(i),rcz(i) = (x,z) coordinates of receivers
! earth = radius of Earth (in km)
! goxd,gozd = gox,goz in degrees
! dvxd,dvzd = dvx,dvz in degrees
! dnzd,dnzd = dnx,dnz in degrees
! gdx,gdz = grid dicing in x and z
! vnl,vnr,vnb,vnt = Bounds of refined grid
! nrnx,nrnz = Number of nodes in x and z for refined grid
! gorx,gorz = Grid origin of refined grid
! sgdl = Source grid dicing level
! rbint = Ray-boundary intersection (0=no, 1=yes).
! asgr = Apply source grid refinement (0=no,1=yes)
! srs = Source-receiver status (0=no path, 1=path exists)
!
END MODULE globalp
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: MODULE
! CODE: FORTRAN 90
! This module contains all the subroutines used to calculate
! the first-arrival traveltime field through the grid.
! Subroutines are:
! (1) travel
! (2) fouds1
! (3) fouds2
! (4) addtree
! (5) downtree
! (6) updtree
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
MODULE traveltime
USE globalp
IMPLICIT NONE
INTEGER ntr
TYPE backpointer
INTEGER(KIND=2) :: px,pz
END TYPE backpointer
TYPE(backpointer), DIMENSION (:), ALLOCATABLE :: btg
!
! btg = backpointer to relate grid nodes to binary tree entries
! px = grid-point in x
! pz = grid-point in z
! ntr = number of entries in binary tree
!
CONTAINS
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: SUBROUTINE
! CODE: FORTRAN 90
! This subroutine is passed the location of a source, and from
! this point the first-arrival traveltime field through the
! velocity grid is determined.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE travel(scx,scz,urg)
IMPLICIT NONE
INTEGER :: isx,isz,sw,i,j,ix,iz,urg,swrg
REAL(KIND=i10) :: scx,scz,vsrc,dsx,dsz,ds
REAL(KIND=i10), DIMENSION (2,2) :: vss
! isx,isz = grid cell indices (i,j,k) which contains source
! scx,scz = (r,x,y) location of source
! sw = a switch (0=off,1=on)
! ix,iz = j,k position of "close" point with minimum traveltime
! maxbt = maximum size of narrow band binary tree
! rd2,rd3 = substitution variables
! vsrc = velocity at source
! vss = velocity at nodes surrounding source
! dsx, dsz = distance from source to cell boundary in x and z
! ds = distance from source to nearby node
! urg = use refined grid (0=no,1=yes,2=previously used)
! swrg = switch to end refined source grid computation
!
! The first step is to find out where the source resides
! in the grid of nodes. The cell in which it resides is
! identified by the "north-west" node of the cell. If the
! source lies on the edge or corner (a node) of the cell, then
! this scheme still applies.
!
isx=INT((scx-gox)/dnx)+1
isz=INT((scz-goz)/dnz)+1
sw=0
IF(isx.lt.1.or.isx.gt.nnx)sw=1
IF(isz.lt.1.or.isz.gt.nnz)sw=1
IF(sw.eq.1)then
isx=90.0-isx*180.0/pi
isz=isz*180.0/pi
WRITE(6,*)"Source lies outside bounds of model (lat,long)= ",isx,isz
WRITE(6,*)"TERMINATING PROGRAM!!!"
STOP
ENDIF
IF(isx.eq.nnx)isx=isx-1
IF(isz.eq.nnz)isz=isz-1
!
! Set all values of nsts to -1 if beginning from a source
! point.
!
IF(urg.NE.2)nsts=-1
!
! set initial size of binary tree to zero
!
ntr=0
IF(urg.EQ.2)THEN
!
! In this case, source grid refinement has been applied, so
! the initial narrow band will come from resampling the
! refined grid.
!
DO i=1,nnx
DO j=1,nnz
IF(nsts(j,i).GT.0)THEN
CALL addtree(j,i)
ENDIF
ENDDO
ENDDO
ELSE
!
! In general, the source point need not lie on a grid point.
! Bi-linear interpolation is used to find velocity at the
! source point.
!
nsts=-1
DO i=1,2
DO j=1,2
vss(i,j)=veln(isz-1+j,isx-1+i)
ENDDO
ENDDO
dsx=(scx-gox)-(isx-1)*dnx
dsz=(scz-goz)-(isz-1)*dnz
CALL bilinear(vss,dsx,dsz,vsrc)
!
! Now find the traveltime at the four surrounding grid points. This
! is calculated approximately by assuming the traveltime from the
! source point to each node is equal to the the distance between
! the two points divided by the average velocity of the points
!
DO i=1,2
DO j=1,2
ds=SQRT((dsx-(i-1)*dnx)**2+(dsz-(j-1)*dnz)**2)
ttn(isz-1+j,isx-1+i)=2.0*ds/(vss(i,j)+vsrc)
CALL addtree(isz-1+j,isx-1+i)
ENDDO
ENDDO
ENDIF
!
! Now calculate the first-arrival traveltimes at the
! remaining grid points. This is done via a loop which
! repeats the procedure of finding the first-arrival
! of all "close" points, adding it to the set of "alive"
! points and updating the points surrounding the new "alive"
! point. The process ceases when the binary tree is empty,
! in which case all grid points are "alive".
!
DO WHILE(ntr.gt.0)
!
! First, check whether source grid refinement is
! being applied; if so, then there is a special
! exit condition.
!
IF(urg.EQ.1)THEN
ix=btg(1)%px
iz=btg(1)%pz
swrg=0
IF(ix.EQ.1)THEN
IF(vnl.NE.1)swrg=1
ENDIF
IF(ix.EQ.nnx)THEN
IF(vnr.NE.nnx)swrg=1
ENDIF
IF(iz.EQ.1)THEN
IF(vnt.NE.1)swrg=1
ENDIF
IF(iz.EQ.nnz)THEN
IF(vnb.NE.nnz)swrg=1
ENDIF
IF(swrg.EQ.1)THEN
nsts(iz,ix)=0
EXIT
ENDIF
ENDIF
!
! Set the "close" point with minimum traveltime
! to "alive"
!
ix=btg(1)%px
iz=btg(1)%pz
nsts(iz,ix)=0
!
! Update the binary tree by removing the root and
! sweeping down the tree.
!
CALL downtree
!
! Now update or find values of up to four grid points
! that surround the new "alive" point.
!
! Test points that vary in x
!
DO i=ix-1,ix+1,2
IF(i.ge.1.and.i.le.nnx)THEN
IF(nsts(iz,i).eq.-1)THEN
!
! This option occurs when a far point is added to the list
! of "close" points
!
IF(fom.eq.0)THEN
CALL fouds1(iz,i)
ELSE
CALL fouds2(iz,i)
ENDIF
CALL addtree(iz,i)
ELSE IF(nsts(iz,i).gt.0)THEN
!
! This happens when a "close" point is updated
!
IF(fom.eq.0)THEN
CALL fouds1(iz,i)
ELSE
CALL fouds2(iz,i)
ENDIF
CALL updtree(iz,i)
ENDIF
ENDIF
ENDDO
!
! Test points that vary in z
!
DO i=iz-1,iz+1,2
IF(i.ge.1.and.i.le.nnz)THEN
IF(nsts(i,ix).eq.-1)THEN
!
! This option occurs when a far point is added to the list
! of "close" points
!
IF(fom.eq.0)THEN
CALL fouds1(i,ix)
ELSE
CALL fouds2(i,ix)
ENDIF
CALL addtree(i,ix)
ELSE IF(nsts(i,ix).gt.0)THEN
!
! This happens when a "close" point is updated
!
IF(fom.eq.0)THEN
CALL fouds1(i,ix)
ELSE
CALL fouds2(i,ix)
ENDIF
CALL updtree(i,ix)
ENDIF
ENDIF
ENDDO
ENDDO
END SUBROUTINE travel
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: SUBROUTINE
! CODE: FORTRAN 90
! This subroutine calculates a trial first-arrival traveltime
! at a given node from surrounding nodes using the
! First-Order Upwind Difference Scheme (FOUDS) of
! Sethian and Popovici (1999).
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE fouds1(iz,ix)
IMPLICIT NONE
INTEGER :: j,k,ix,iz,tsw1,swsol
REAL(KIND=i10) :: trav,travm,slown,tdsh,tref
REAL(KIND=i10) :: a,b,c,u,v,em,ri,risti
REAL(KIND=i10) :: rd1
!
! ix = NS position of node coordinate for determination
! iz = EW vertical position of node coordinate for determination
! trav = traveltime calculated for trial node
! travm = minimum traveltime calculated for trial node
! slown = slowness at (iz,ix)
! tsw1 = traveltime switch (0=first time,1=previously)
! a,b,c,u,v,em = Convenience variables for solving quadratic
! tdsh = local traveltime from neighbouring node
! tref = reference traveltime at neighbouring node
! ri = Radial distance
! risti = ri*sin(theta) at point (iz,ix)
! rd1 = dummy variable
! swsol = switch for solution (0=no solution, 1=solution)
!
! Inspect each of the four quadrants for the minimum time
! solution.
!
tsw1=0
slown=1.0/veln(iz,ix)
ri=earth
risti=ri*sin(gox+(ix-1)*dnx)
DO j=ix-1,ix+1,2
DO k=iz-1,iz+1,2
IF(j.GE.1.AND.j.LE.nnx)THEN
IF(k.GE.1.AND.k.LE.nnz)THEN
!
! There are seven solution options in
! each quadrant.
!
swsol=0
IF(nsts(iz,j).EQ.0)THEN
swsol=1
IF(nsts(k,ix).EQ.0)THEN
u=ri*dnx
v=risti*dnz
em=ttn(k,ix)-ttn(iz,j)
a=u**2+v**2
b=-2.0*u**2*em
c=u**2*(em**2-v**2*slown**2)
tref=ttn(iz,j)
ELSE
a=1.0
b=0.0
c=-slown**2*ri**2*dnx**2
tref=ttn(iz,j)
ENDIF
ELSE IF(nsts(k,ix).EQ.0)THEN
swsol=1
a=1.0
b=0.0
c=-(slown*risti*dnz)**2
tref=ttn(k,ix)
ENDIF
!
! Now find the solution of the quadratic equation
!
IF(swsol.EQ.1)THEN
rd1=b**2-4.0*a*c
IF(rd1.LT.0.0)rd1=0.0
tdsh=(-b+sqrt(rd1))/(2.0*a)
trav=tref+tdsh
IF(tsw1.EQ.1)THEN
travm=MIN(trav,travm)
ELSE
travm=trav
tsw1=1
ENDIF
ENDIF
ENDIF
ENDIF
ENDDO
ENDDO
ttn(iz,ix)=travm
END SUBROUTINE fouds1
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: SUBROUTINE
! CODE: FORTRAN 90
! This subroutine calculates a trial first-arrival traveltime
! at a given node from surrounding nodes using the
! Mixed-Order (2nd) Upwind Difference Scheme (FOUDS) of
! Popovici and Sethian (2002).
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE fouds2(iz,ix)
IMPLICIT NONE
INTEGER :: j,k,j2,k2,ix,iz,tsw1
INTEGER :: swj,swk,swsol
REAL(KIND=i10) :: trav,travm,slown,tdsh,tref,tdiv
REAL(KIND=i10) :: a,b,c,u,v,em,ri,risti,rd1
!
! ix = NS position of node coordinate for determination
! iz = EW vertical position of node coordinate for determination
! trav = traveltime calculated for trial node
! travm = minimum traveltime calculated for trial node
! slown = slowness at (iz,ix)
! tsw1 = traveltime switch (0=first time,1=previously)
! a,b,c,u,v,em = Convenience variables for solving quadratic
! tdsh = local traveltime from neighbouring node
! tref = reference traveltime at neighbouring node
! ri = Radial distance
! risti = ri*sin(theta) at point (iz,ix)
! swj,swk = switches for second order operators
! tdiv = term to divide tref by depending on operator order
! swsol = switch for solution (0=no solution, 1=solution)
!
! Inspect each of the four quadrants for the minimum time
! solution.
!
tsw1=0
slown=1.0/veln(iz,ix)
ri=earth
risti=ri*sin(gox+(ix-1)*dnx)
DO j=ix-1,ix+1,2
IF(j.GE.1.AND.j.LE.nnx)THEN
swj=-1
IF(j.eq.ix-1)THEN
j2=j-1
IF(j2.GE.1)THEN
IF(nsts(iz,j2).EQ.0)swj=0
ENDIF
ELSE
j2=j+1
IF(j2.LE.nnx)THEN
IF(nsts(iz,j2).EQ.0)swj=0
ENDIF
ENDIF
IF(nsts(iz,j).EQ.0.AND.swj.EQ.0)THEN
swj=-1
IF(ttn(iz,j).GT.ttn(iz,j2))THEN
swj=0
ENDIF
ELSE
swj=-1
ENDIF
DO k=iz-1,iz+1,2
IF(k.GE.1.AND.k.LE.nnz)THEN
swk=-1
IF(k.eq.iz-1)THEN
k2=k-1
IF(k2.GE.1)THEN
IF(nsts(k2,ix).EQ.0)swk=0
ENDIF
ELSE
k2=k+1
IF(k2.LE.nnz)THEN
IF(nsts(k2,ix).EQ.0)swk=0
ENDIF
ENDIF
IF(nsts(k,ix).EQ.0.AND.swk.EQ.0)THEN
swk=-1
IF(ttn(k,ix).GT.ttn(k2,ix))THEN
swk=0
ENDIF
ELSE
swk=-1
ENDIF
!
! There are 8 solution options in
! each quadrant.
!
swsol=0
IF(swj.EQ.0)THEN
swsol=1
IF(swk.EQ.0)THEN
u=2.0*ri*dnx
v=2.0*risti*dnz
em=4.0*ttn(iz,j)-ttn(iz,j2)-4.0*ttn(k,ix)
em=em+ttn(k2,ix)
a=v**2+u**2
b=2.0*em*u**2
c=u**2*(em**2-slown**2*v**2)
tref=4.0*ttn(iz,j)-ttn(iz,j2)
tdiv=3.0
ELSE IF(nsts(k,ix).EQ.0)THEN
u=risti*dnz
v=2.0*ri*dnx
em=3.0*ttn(k,ix)-4.0*ttn(iz,j)+ttn(iz,j2)
a=v**2+9.0*u**2
b=6.0*em*u**2
c=u**2*(em**2-slown**2*v**2)
tref=ttn(k,ix)
tdiv=1.0
ELSE
u=2.0*ri*dnx
a=1.0
b=0.0
c=-u**2*slown**2
tref=4.0*ttn(iz,j)-ttn(iz,j2)
tdiv=3.0
ENDIF
ELSE IF(nsts(iz,j).EQ.0)THEN
swsol=1
IF(swk.EQ.0)THEN
u=ri*dnx
v=2.0*risti*dnz
em=3.0*ttn(iz,j)-4.0*ttn(k,ix)+ttn(k2,ix)
a=v**2+9.0*u**2
b=6.0*em*u**2
c=u**2*(em**2-v**2*slown**2)
tref=ttn(iz,j)
tdiv=1.0
ELSE IF(nsts(k,ix).EQ.0)THEN
u=ri*dnx
v=risti*dnz
em=ttn(k,ix)-ttn(iz,j)
a=u**2+v**2
b=-2.0*u**2*em
c=u**2*(em**2-v**2*slown**2)
tref=ttn(iz,j)
tdiv=1.0
ELSE
a=1.0
b=0.0
c=-slown**2*ri**2*dnx**2
tref=ttn(iz,j)
tdiv=1.0
ENDIF
ELSE
IF(swk.EQ.0)THEN
swsol=1
u=2.0*risti*dnz
a=1.0
b=0.0
c=-u**2*slown**2
tref=4.0*ttn(k,ix)-ttn(k2,ix)
tdiv=3.0
ELSE IF(nsts(k,ix).EQ.0)THEN
swsol=1
a=1.0
b=0.0
c=-slown**2*risti**2*dnz**2
tref=ttn(k,ix)
tdiv=1.0
ENDIF
ENDIF
!
! Now find the solution of the quadratic equation
!
IF(swsol.EQ.1)THEN
rd1=b**2-4.0*a*c
IF(rd1.LT.0.0)rd1=0.0
tdsh=(-b+sqrt(rd1))/(2.0*a)
trav=(tref+tdsh)/tdiv
IF(tsw1.EQ.1)THEN
travm=MIN(trav,travm)
ELSE
travm=trav
tsw1=1
ENDIF
ENDIF
ENDIF
ENDDO
ENDIF
ENDDO
ttn(iz,ix)=travm
END SUBROUTINE fouds2
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: SUBROUTINE
! CODE: FORTRAN 90
! This subroutine adds a value to the binary tree by
! placing a value at the bottom and pushing it up
! to its correct position.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE addtree(iz,ix)
IMPLICIT NONE
INTEGER :: ix,iz,tpp,tpc
TYPE(backpointer) :: exch
!
! ix,iz = grid position of new addition to tree
! tpp = tree position of parent
! tpc = tree position of child
! exch = dummy to exchange btg values
!
! First, increase the size of the tree by one.
!
ntr=ntr+1
!
! Put new value at base of tree
!
nsts(iz,ix)=ntr
btg(ntr)%px=ix
btg(ntr)%pz=iz
!
! Now filter the new value up to its correct position
!
tpc=ntr
tpp=tpc/2
DO WHILE(tpp.gt.0)
IF(ttn(iz,ix).lt.ttn(btg(tpp)%pz,btg(tpp)%px))THEN
nsts(iz,ix)=tpp
nsts(btg(tpp)%pz,btg(tpp)%px)=tpc
exch=btg(tpc)
btg(tpc)=btg(tpp)
btg(tpp)=exch
tpc=tpp
tpp=tpc/2
ELSE
tpp=0
ENDIF
ENDDO
END SUBROUTINE addtree
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: SUBROUTINE
! CODE: FORTRAN 90
! This subroutine updates the binary tree after the root
! value has been used. The root is replaced by the value
! at the bottom of the tree, which is then filtered down
! to its correct position. This ensures that the tree remains
! balanced.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE downtree
IMPLICIT NONE
INTEGER :: tpp,tpc
REAL(KIND=i10) :: rd1,rd2
TYPE(backpointer) :: exch
!
! tpp = tree position of parent
! tpc = tree position of child
! exch = dummy to exchange btg values
! rd1,rd2 = substitution variables
!
! Replace root of tree with its last value
!
IF(ntr.EQ.1)THEN
ntr=ntr-1
RETURN
ENDIF
nsts(btg(ntr)%pz,btg(ntr)%px)=1
btg(1)=btg(ntr)
!
! Reduce size of tree by one
!
ntr=ntr-1
!
! Now filter new root down to its correct position
!
tpp=1
tpc=2*tpp
DO WHILE(tpc.lt.ntr)
!
! Check which of the two children is smallest - use the smallest
!
rd1=ttn(btg(tpc)%pz,btg(tpc)%px)
rd2=ttn(btg(tpc+1)%pz,btg(tpc+1)%px)
IF(rd1.gt.rd2)THEN
tpc=tpc+1
ENDIF
!
! Check whether the child is smaller than the parent; if so, then swap,
! if not, then we are done
!
rd1=ttn(btg(tpc)%pz,btg(tpc)%px)
rd2=ttn(btg(tpp)%pz,btg(tpp)%px)
IF(rd1.lt.rd2)THEN
nsts(btg(tpp)%pz,btg(tpp)%px)=tpc
nsts(btg(tpc)%pz,btg(tpc)%px)=tpp
exch=btg(tpc)
btg(tpc)=btg(tpp)
btg(tpp)=exch
tpp=tpc
tpc=2*tpp
ELSE
tpc=ntr+1
ENDIF
ENDDO
!
! If ntr is an even number, then we still have one more test to do
!
IF(tpc.eq.ntr)THEN
rd1=ttn(btg(tpc)%pz,btg(tpc)%px)
rd2=ttn(btg(tpp)%pz,btg(tpp)%px)
IF(rd1.lt.rd2)THEN
nsts(btg(tpp)%pz,btg(tpp)%px)=tpc
nsts(btg(tpc)%pz,btg(tpc)%px)=tpp
exch=btg(tpc)
btg(tpc)=btg(tpp)
btg(tpp)=exch
ENDIF
ENDIF
END SUBROUTINE downtree
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: SUBROUTINE
! CODE: FORTRAN 90
! This subroutine updates a value on the binary tree. The FMM
! should only produce updated values that are less than their
! prior values.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
SUBROUTINE updtree(iz,ix)
IMPLICIT NONE
INTEGER :: ix,iz,tpp,tpc
TYPE(backpointer) :: exch
!
! ix,iz = grid position of new addition to tree
! tpp = tree position of parent
! tpc = tree position of child
! exch = dummy to exchange btg values
!
! Filter the updated value to its correct position
!
tpc=nsts(iz,ix)
tpp=tpc/2
DO WHILE(tpp.gt.0)
IF(ttn(iz,ix).lt.ttn(btg(tpp)%pz,btg(tpp)%px))THEN
nsts(iz,ix)=tpp
nsts(btg(tpp)%pz,btg(tpp)%px)=tpc
exch=btg(tpc)
btg(tpc)=btg(tpp)
btg(tpp)=exch
tpc=tpp
tpp=tpc/2
ELSE
tpp=0
ENDIF
ENDDO
END SUBROUTINE updtree
END MODULE traveltime
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! MAIN PROGRAM
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! TYPE: PROGRAM
! CODE: FORTRAN 90
! This program is designed to implement the Fast Marching
! Method (FMM) for calculating first-arrival traveltimes
! through a 2-D continuous velocity medium in spherical shell
! coordinates (x=theta or latitude, z=phi or longitude).
! It is written in Fortran 90, although it is probably more
! accurately described as Fortran 77 with some of the Fortran 90
! extensions.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!PROGRAM tomo_surf
subroutine synthetic(nx,ny,nz,nparpi,vels,obst, &
goxdf,gozdf,dvxdf,dvzdf,kmaxRc,kmaxRg,kmaxLc,kmaxLg, &
tRc,tRg,tLc,tLg,wavetype,igrt,periods,depz,minthk, &
scxf,sczf,rcxf,rczf,nrc1,nsrcsurf1,knum1,kmax,nsrcsurf,nrcf,noiselevel)
USE globalp
USE traveltime
IMPLICIT NONE
!CHARACTER (LEN=30) ::grid,frechet
!CHARACTER (LEN=40) :: sources,receivers,otimes
!CHARACTER (LEN=30) :: travelt,rtravel,wrays,cdum
INTEGER :: i,j,k,l,nsrc,tnr,urg
INTEGER :: sgs,isx,isz,sw,idm1,idm2,nnxb,nnzb
INTEGER :: ogx,ogz,grdfx,grdfz,maxbt
REAL(KIND=i10) :: x,z,goxb,gozb,dnxb,dnzb
!REAL(KIND=i10), DIMENSION (:,:), ALLOCATABLE :: scxf,sczf
!REAL(KIND=i10), DIMENSION (:,:,:), ALLOCATABLE :: rcxf,rczf
!
! sources = File containing source locations
! receivers = File containing receiver locations
! grid = File containing grid of velocity vertices for
! resampling on a finer grid with cubic B-splines
! frechet = output file containing matrix of frechet derivatives
! travelt = File name for storage of traveltime field
! wttf = Write traveltimes to file? (0=no,>0=source id)
! fom = Use first-order(0) or mixed-order(1) scheme
! nsrc = number of sources
! scx,scz = source location in r,x,z
! scx,scz = source location in r,x,z
! x,z = temporary variables for source location
! fsrt = find source-receiver traveltimes? (0=no,1=yes)
! rtravel = output file for source-receiver traveltimes
! cdum = dummy character variable ! wrgf = write ray geometries to file? (<0=all,0=no,>0=source id.)
! wrays = file containing raypath geometries
! cfd = calculate Frechet derivatives? (0=no, 1=yes)
! tnr = total number of receivers
! sgs = Extent of refined source grid
! isx,isz = cell containing source
! nnxb,nnzb = Backup for nnz,nnx
! goxb,gozb = Backup for gox,goz
! dnxb,dnzb = Backup for dnx,dnz
! ogx,ogz = Location of refined grid origin
! gridfx,grdfz = Number of refined nodes per cell
! urg = use refined grid (0=no,1=yes,2=previously used)
! maxbt = maximum size of narrow band binary tree
! otimes = file containing source-receiver association information
!c-----------------------------------------------------------------
! variables defined by Hongjian Fang
integer nx,ny,nz
integer kmax,nsrcsurf,nrcf
real vels(nx,ny,nz)
real obst(*)
real goxdf,gozdf,dvxdf,dvzdf
integer kmaxRc,kmaxRg,kmaxLc,kmaxLg
real*8 tRc(*),tRg(*),tLc(*),tLg(*)
integer wavetype(nsrcsurf,kmax)
integer periods(nsrcsurf,kmax),nrc1(nsrcsurf,kmax),nsrcsurf1(kmax)
integer knum1(kmax),igrt(nsrcsurf,kmax)
real scxf(nsrcsurf,kmax),sczf(nsrcsurf,kmax),rcxf(nrcf,nsrcsurf,kmax),rczf(nrcf,nsrcsurf,kmax)
real minthk
integer nparpi
real vpz(nz),vsz(nz),rhoz(nz),depz(nz)
real*8 pvRc(nx*ny,kmaxRc),pvRg(nx*ny,kmaxRg),pvLc(nx*ny,kmaxLc),pvLg(nx*ny,kmaxLg)
real*8 velf(ny*nx)
integer kmax1,kmax2,kmax3,count1
integer igr
integer iwave
integer knumi,srcnum
real cbst1
real noiselevel
real gaussian
external gaussian
integer ii,jj,kk,nn,istep
gdx=3
gdz=3
asgr=1
sgdl=5
sgs=5
earth=6371.0
fom=1
snb=0.5
goxd=goxdf
gozd=gozdf
dvxd=dvxdf
dvzd=dvzdf
nvx=nx-2
nvz=ny-2
ALLOCATE(velv(0:nvz+1,0:nvx+1), STAT=checkstat)
IF(checkstat > 0)THEN
WRITE(6,*)'Error with ALLOCATE: SUBROUTINE gridder: REAL velv'
ENDIF
!
! Convert from degrees to radians
!
dvx=dvxd*pi/180.0
dvz=dvzd*pi/180.0
gox=(90.0-goxd)*pi/180.0
goz=gozd*pi/180.0
!
! Compute corresponding values for propagation grid.
!
nnx=(nvx-1)*gdx+1
nnz=(nvz-1)*gdz+1
dnx=dvx/gdx
dnz=dvz/gdz
dnxd=dvxd/gdx
dnzd=dvzd/gdz
ALLOCATE(veln(nnz,nnx), STAT=checkstat)
IF(checkstat > 0)THEN
WRITE(6,*)'Error with ALLOCATE: SUBROUTINE gridder: REAL veln'
ENDIF
!
! Call a subroutine which reads in the velocity grid
!
!CALL gridder(grid)
!
! Read in all source coordinates.
!
!
! Now work out, source by source, the first-arrival traveltime
! field plus source-receiver traveltimes
! and ray paths if required. First, allocate memory to the
! traveltime field array
!
ALLOCATE(ttn(nnz,nnx), STAT=checkstat)
IF(checkstat > 0)THEN
WRITE(6,*)'Error with ALLOCATE: PROGRAM fmmin2d: REAL ttn'
ENDIF
rbint=0
!
! Allocate memory for node status and binary trees
!
ALLOCATE(nsts(nnz,nnx))
maxbt=NINT(snb*nnx*nnz)
ALLOCATE(btg(maxbt))
!allocate(fdm(0:nvz+1,0:nvx+1))
! if(kmaxRc.gt.0) then
! iwave=2
! igr=0
! call depthkernel(nx,ny,nz,vels,pvRc,sen_vsRc,sen_vpRc, &
! sen_rhoRc,iwave,igr,kmaxRc,tRc,depz,minthk)
! endif
!
! if(kmaxRg.gt.0) then
! iwave=2
! igr=1
! call depthkernel(nx,ny,nz,vels,pvRg,sen_vsRg,sen_vpRg, &
! sen_rhoRg,iwave,igr,kmaxRg,tRg,depz,minthk)
! endif
!
! if(kmaxLc.gt.0) then
! iwave=1
! igr=0
! call depthkernel(nx,ny,nz,vels,pvLc,sen_vsLc,sen_vpLc, &
! sen_rhoLc,iwave,igr,kmaxLc,tLc,depz,minthk)
! endif
!
! if(kmaxLg.gt.0) then
! iwave=1
! igr=1
! call depthkernel(nx,ny,nz,vels,pvLg,sen_vsLg,sen_vpLg, &
! sen_rhoLg,iwave,igr,kmaxLg,tLg,depz,minthk)
! endif
!nar=0
count1=0
!sen_vs=0
!sen_vp=0
!sen_rho=0
kmax1=kmaxRc
kmax2=kmaxRc+kmaxRg
kmax3=kmaxRc+kmaxRg+kmaxLc
do knumi=1,kmax
do srcnum=1,nsrcsurf1(knum1(knumi))
if(wavetype(srcnum,knum1(knumi))==2.and.igrt(srcnum,knum1(knumi))==0) then
velf(1:nx*ny)=pvRc(1:nx*ny,periods(srcnum,knum1(knumi)))
! sen_vs(:,1:kmax1,:)=sen_vsRc(:,1:kmaxRc,:)!(:,nt(istep),:)
! sen_vp(:,1:kmax1,:)=sen_vpRc(:,1:kmaxRc,:)!(:,nt(istep),:)
! sen_rho(:,1:kmax1,:)=sen_rhoRc(:,1:kmaxRc,:)!(:,nt(istep),:)
endif
if(wavetype(srcnum,knum1(knumi))==2.and.igrt(srcnum,knum1(knumi))==1) then
velf(1:nx*ny)=pvRg(1:nx*ny,periods(srcnum,knum1(knumi)))
! sen_vs(:,kmax1+1:kmax2,:)=sen_vsRg(:,1:kmaxRg,:)!(:,nt,:)
! sen_vp(:,kmax1+1:kmax2,:)=sen_vpRg(:,1:kmaxRg,:)!(:,nt,:)
! sen_rho(:,kmax1+1:kmax2,:)=sen_rhoRg(:,1:kmaxRg,:)!(:,nt,:)
endif
if(wavetype(srcnum,knum1(knumi))==1.and.igrt(srcnum,knum1(knumi))==0) then
velf(1:nx*ny)=pvLc(1:nx*ny,periods(srcnum,knum1(knumi)))
! sen_vs(:,kmax2+1:kmax3,:)=sen_vsLc(:,1:kmaxLc,:)!(:,nt,:)
! sen_vp(:,kmax2+1:kmax3,:)=sen_vpLc(:,1:kmaxLc,:)!(:,nt,:)
! sen_rho(:,kmax2+1:kmax3,:)=sen_rhoLc(:,1:kmaxLc,:)!(:,nt,:)
endif
if(wavetype(srcnum,knum1(knumi))==1.and.igrt(srcnum,knum1(knumi))==1) then
velf(1:nx*ny)=pvLg(1:nx*ny,periods(srcnum,knum1(knumi)))
! sen_vs(:,kmax3+1:kmax,:)=sen_vsLg(:,1:kmaxLg,:)!(:,nt,:)
! sen_vp(:,kmax3+1:kmax,:)=sen_vpLg(:,1:kmaxLg,:)!(:,nt,:)
! sen_rho(:,kmax3+1:kmax,:)=sen_rhoLg(:,1:kmaxLg,:)!(:,nt,:)
endif
call gridder(velf)
x=scxf(srcnum,knum1(knumi))
z=sczf(srcnum,knum1(knumi))
!
! Begin by computing refined source grid if required
!
urg=0
IF(asgr.EQ.1)THEN
!
! Back up coarse velocity grid to a holding matrix
!
ALLOCATE(velnb(nnz,nnx))
velnb=veln
nnxb=nnx
nnzb=nnz
dnxb=dnx
dnzb=dnz
goxb=gox
gozb=goz
!
! Identify nearest neighbouring node to source
!
isx=INT((x-gox)/dnx)+1
isz=INT((z-goz)/dnz)+1
sw=0
IF(isx.lt.1.or.isx.gt.nnx)sw=1
IF(isz.lt.1.or.isz.gt.nnz)sw=1