-
Notifications
You must be signed in to change notification settings - Fork 146
/
interpolate.f90
670 lines (547 loc) · 22.7 KB
/
interpolate.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
! this program creates first a specfem mesh based on the SEP-format mesh, then interpolates the wave speeds
! onto all GLL points and runs a new simulation with the interpolated mesh.
program main_interpolate
implicit none
!-----------------------------------------------------------------------
! SEP format model files
character(len=512),parameter :: sep_directory = './SEG_2D_SALT/'
character(len=512),parameter :: sep_header_file_vp = 'vp.H'
character(len=512),parameter :: sep_header_file_vs = 'vs.H' ! not used if scaling from vp
character(len=512),parameter :: sep_header_file_rho = 'rho.H' ! not used if scaling from vp
! Option to scale VS & RHO from VP
logical,parameter :: SCALE_FROM_VP = .true.
! Interpolation: closest_point ('closet_point=.true.') or bilinear interpolation ('closet_point=.false.')
logical,parameter :: INTERPOLATE_FROM_CLOSEST_POINT = .true.
! Option for interpolation points: 1 == corners only / 2 == corner & midpoints used for interpolation
integer,parameter :: INTERPOLATION_POINTS = 1
! lets user choose number of MPI processes
logical,parameter :: USER_INPUT = .false.
!-----------------------------------------------------------------------
integer :: NX,NY,NZ,esize,nproc,ier
real :: OX,OY,OZ,DX,DY,DZ
real,dimension(:,:),allocatable :: vp_SEP,vs_SEP,rho_SEP
character(len=512) :: system_command,sep_file,data_format,line
character(len=128) :: tmpstring
character(len=128) :: numstring
character(len=3) ::num
! user output
print *
print *, '***********************************************************'
print *, 'Model Interpolation'
print *, '***********************************************************'
print *
! user input
if (USER_INPUT) then
print *, 'please input the number of processes:'
read(*,*) nproc
print *
else
! reads from default Par_file
nproc = 0
open(unit=15,file='./DATA/Par_file',status='old',iostat=ier)
if (ier /= 0) stop 'Error opening DATA/Par_file'
do while (ier == 0)
read(15,'(a512)',iostat=ier) line
if (ier == 0) then
! left trim
line = adjustl(line)
! skip comment lines
if (line(1:1) == '#') cycle
! suppress trailing comment
if (index(line,'#') > 0) line = line(1:index(line,'#')-1)
! checks parameter name
!print *,'line: ',line(1:5)
if (line(1:5) == 'NPROC') then
read(line(index(line,'=')+1:len_trim(line)),'(i3)') nproc
exit
endif
endif
enddo
close(15)
if (nproc < 1) stop 'Error could not read nproc from Par_file'
endif
! user output
print *, 'setting:'
print *, ' SEP model from directory : ',trim(sep_directory)
print *, ' using number of processes: ',nproc
if (INTERPOLATE_FROM_CLOSEST_POINT) then
print *, ' interpolation type : from closest point'
else
print *, ' interpolation type : bilinear'
endif
print *, ' interpolation points set : ',INTERPOLATION_POINTS
print *, ' using scaling from Vp : ',SCALE_FROM_VP
print *
! reading SEP models
print *, 'reading SEP models...'
! VP model
! gets header info for VP model in SEP format
call READ_SEP_HEADER(sep_directory,sep_header_file_vp, &
sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ, &
esize,data_format)
! assuming that VP,VS,RHO models have the same dimension
allocate(vp_SEP(NX,NZ),vs_SEP(NX,NZ),rho_SEP(NX,NZ))
! reads in VP model
call READ_SEP_MODEL_2D(sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ,esize,data_format,vp_SEP)
! VS model
if (.not. SCALE_FROM_VP) then
call READ_SEP_HEADER(sep_directory,sep_header_file_vs, &
sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ, &
esize,data_format)
call READ_SEP_MODEL_2D(sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ,esize,data_format,vs_SEP)
endif
! RHO model
if (.not. SCALE_FROM_VP) then
call READ_SEP_HEADER(sep_directory,sep_header_file_rho, &
sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ, &
esize,data_format)
call READ_SEP_MODEL_2D(sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ,esize,data_format,rho_SEP)
endif
! creates interface file
print *
print *,'Writing interface data to file DATA/interface_industry.dat'
print *
open(unit=15,file='./DATA/interface_industry.dat',status='unknown')
write(15,'(a)')'#'
write(15,'(a)')'# number of interfaces'
write(15,'(a)')'#'
write(15,'(a)')' 2'
write(15,'(a)')'#'
write(15,'(a)')'# for each interface below, we give the number of points and then x,z for each point'
write(15,'(a)')'#'
write(15,'(a)')'#'
write(15,'(a)')'# interface number 1 (bottom of the mesh)'
write(15,'(a)')'#'
write(15,'(a)')' 2'
write(15,*) OX,OZ
write(15,*) OX+(NX-1)*DX,OZ
write(15,'(a)')'#'
write(15,'(a)')'# interface number 2 (topography, top of the mesh)'
write(15,'(a)')'#'
write(15,'(a)')' 2'
write(15,*) OX,OZ+(NZ-1)*DZ
write(15,*) OX+(NX-1)*DX,OZ+(NZ-1)*DZ
write(15,'(a)')'#'
write(15,'(a)')'# for each layer, we give the number of spectral elements in the vertical direction'
write(15,'(a)')'#'
write(15,'(a)')'#'
write(15,'(a)')'# layer number 1 (bottom layer)'
write(15,'(a)')'#'
write(15,*) (NZ-1)/INTERPOLATION_POINTS
close(15)
! model geometry
print *
print *, '************ Preparing Par_file ************ '
print *
!print *, 'OX = ',OX
write(numstring,*) OX
write(tmpstring,*) 'xmin =',trim(numstring)
tmpstring = adjustl(tmpstring)
!write(*,*) trim(tmpstring)
write(system_command,*) 'sed -i "s/^xmin .*/',trim(tmpstring),'/" ./DATA/Par_file'
write (*,*) trim(system_command)
call execute_command_line(trim(system_command),exitstat=ier)
if (ier /= 0) stop
write(numstring,*) OX+(NX-1)*DX
write(tmpstring,*) 'xmax =',trim(numstring)
tmpstring = adjustl(tmpstring)
!write(*,*) trim(tmpstring)
write(system_command,*) 'sed -i "s/^xmax .*/',trim(tmpstring),'/" ./DATA/Par_file'
write (*,*) trim(system_command)
call execute_command_line(system_command,exitstat=ier)
if (ier /= 0) stop
write(numstring,*) (NX-1)/INTERPOLATION_POINTS
write(tmpstring,*) 'nx =',trim(numstring)
tmpstring = adjustl(tmpstring)
!write(*,*) trim(tmpstring)
write(system_command,*) 'sed -i "s/^nx .*/',trim(tmpstring),'/" ./DATA/Par_file'
write (*,*) trim(system_command)
call execute_command_line(trim(system_command),exitstat=ier)
if (ier /= 0) stop
! nbregions
!call system('sed -i "$ d" ./DATA/Par_file')
write(tmpstring,*) '1',(NX-1)/INTERPOLATION_POINTS,'1',(NZ-1)/INTERPOLATION_POINTS,'1'
tmpstring = adjustl(tmpstring)
!write (*,*) trim(tmpstring)
!write(system_command,*) 'echo ${tmpstring}'
!write(system_command,*) 'sed -i "/nbregions /a',trim(tmpstring),'" ./DATA/Par_file'
! original Par_file line for regions to replace is: 1 20 1 20 1
write(system_command,*) 'sed -i "s/^1 20 .*/',trim(tmpstring),'/" ./DATA/Par_file'
write (*,*) trim(system_command)
call execute_command_line(trim(system_command),exitstat=ier)
if (ier /= 0) stop
! number of processes
if (USER_INPUT) then
write(tmpstring,*) 'NPROC =',nproc
tmpstring = adjustl(tmpstring)
!write(*,*) trim(tmpstring)
write(system_command,*) 'sed -i "s/^NPROC .*/',trim(tmpstring),'/" ./DATA/Par_file'
write (*,*) trim(system_command)
call execute_command_line(trim(system_command),exitstat=ier)
if (ier /= 0) stop
endif
print *
! backup Par_file
call execute_command_line('cp ./DATA/Par_file ./DATA/Par_file.org',exitstat=ier)
if (ier /= 0) stop
call sleep(1)
! note: first simulation is used to create a mesh and initial model files,
! which will be read in a second run for interpolation
! user output
print *
print *, '************ Setting up first simulation ************ '
print *
! Par_file flags to create new mesh files
call execute_command_line('cp ./DATA/Par_file.org ./DATA/Par_file',exitstat=ier)
if (ier /= 0) stop
call execute_command_line('sed -i "s/^MODEL .*/MODEL = default/" ./DATA/Par_file',exitstat=ier)
if (ier /= 0) stop
call execute_command_line('sed -i "s/^SAVE_MODEL .*/SAVE_MODEL = legacy/" ./DATA/Par_file',exitstat=ier)
if (ier /= 0) stop
call execute_command_line('sed -i "s/^NSTEP .*/NSTEP = 5 /" ./DATA/Par_file',exitstat=ier)
if (ier /= 0) stop
call execute_command_line('cp ./DATA/Par_file ./OUTPUT_FILES/Par_file.step_1',exitstat=ier)
if (ier /= 0) stop
write(num,'(i2.2)') nproc
! runs mesher
print *,'call xmeshfem2d'
call execute_command_line('mpirun -np ' //num//' ./bin/xmeshfem2D > OUTPUT_FILES/output_mesher.step_1.txt',exitstat=ier)
if (ier /= 0) stop
! runs first forward simulation to generate new specfem files
print *,'call xspecfem2d'
call execute_command_line('mpirun -np ' //num//' ./bin/xspecfem2D > OUTPUT_FILES/output_solver.step_1.txt',exitstat=ier)
if (ier /= 0) stop
print *
print *, '************ Setting up new simulation ************ '
print *
! now uses default specfem file format
! backup Par_file
call execute_command_line('cp ./DATA/Par_file.org ./DATA/Par_file',exitstat=ier)
if (ier /= 0) stop
call execute_command_line('sed -i "s/^MODEL .*/MODEL = legacy/" ./DATA/Par_file',exitstat=ier)
if (ier /= 0) stop
call execute_command_line('sed -i "s/^SAVE_MODEL .*/SAVE_MODEL = default/" ./DATA/Par_file',exitstat=ier)
if (ier /= 0) stop
call execute_command_line('cp ./DATA/Par_file ./OUTPUT_FILES/Par_file.step_2',exitstat=ier)
if (ier /= 0) stop
! runs mesher
print *,'call xmeshfem2d'
call execute_command_line('mpirun -np ' //num//' ./bin/xmeshfem2D',exitstat=ier)
if (ier /= 0) stop
! scaling model values
if (SCALE_FROM_VP) then
! VS model, scaled or in SEP format
vs_SEP(:,:) = vp_SEP(:,:) / 1.732
! RHO model, scaled or in SEP format
! note: by default, just assumes constant density
! constant density
rho_SEP(:,:) = 1000.0
! scaling rule
!rho_SEP = (1.6612 * (vp_SEP/14500*4500 / 1000.0) &
! -0.4720 * (vp_SEP/14500*4500 / 1000.0)**2 &
! +0.0671 * (vp_SEP/14500*4500 / 1000.0)**3 &
! -0.0043 * (vp_SEP/14500*4500 / 1000.0)**4 &
! +0.000106*(vp_SEP/14500*4500 / 1000.0)**5 )*1000.0
endif
! user output
print *
print *, '************ Interpolating model ************'
print *
print *, 'input model:'
print *, ' Vp min/max = ',minval(vp_SEP(:,:)),maxval(vp_SEP(:,:))
print *, ' Vs min/max = ',minval(vs_SEP(:,:)),maxval(vs_SEP(:,:))
print *, ' rho min/max = ',minval(rho_SEP(:,:)),maxval(rho_SEP(:,:))
call interpolate_slowness_gll(vp_SEP,vs_SEP,rho_SEP,NX,NY,NZ,DX,DY,DZ,OX,OY,OZ,INTERPOLATE_FROM_CLOSEST_POINT,nproc)
! user output
print *
print *, 'Interpolation done'
print *
! runs forward simulation
print *,'call xspecfem2d'
call execute_command_line('mpirun -np ' //num//' ./bin/xspecfem2D ',exitstat=ier)
if (ier /= 0) stop
end program main_interpolate
!
!-----------------------------------------------------------------------------------
!
subroutine READ_SEP_HEADER(sep_directory,sep_header_file, &
sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ, &
esize,data_format)
implicit none
!!! input
character(len=512),intent(in) :: sep_directory ! where the sep header file and sep file are saved
character(len=512),intent(in) :: sep_header_file ! file name of the sep header file
!!! output
character(len=512),intent(out) :: sep_file ! sep file specified in sep header file
integer,intent(out) :: NX,NY,NZ ! n1,n2,n3 in sep header file
real,intent(out) :: OX,OY,OZ ! o1,o2,o3 in sep header file
real,intent(out) :: DX,DY,DZ ! d1,d2,d3 in sep header file
integer,intent(out) :: esize ! esize in sep header file
character(len=512),intent(out) :: data_format ! data_format in sep header file
!!! local
integer :: ier
character(len=512) :: junk,sep_header_file_complete
sep_header_file_complete = trim(adjustl(sep_directory))//trim(adjustl(sep_header_file))
open(unit=13,file=trim(adjustl(sep_header_file_complete)),status='old',iostat=ier)
print *
print *, '*******************************************************************************'
print *, 'reading sep header file: '
print *, trim(adjustl(sep_header_file_complete))
if (ier /= 0) stop 'ERROR: cannot open sep header file'
read(13,'(a3a)') junk, sep_file
read(13,'(a3i10)') junk, NX
read(13,'(a3i10)') junk, NY
read(13,'(a3i10)') junk, NZ
read(13,'(a3f20.0)') junk, OX
read(13,'(a3f20.0)') junk, OY
read(13,'(a3f20.0)') junk, OZ
read(13,'(a3f20.0)') junk, DX
read(13,'(a3f20.0)') junk, DY
read(13,'(a3f20.0)') junk, DZ
read(13,'(a6i10)') junk, esize
read(13,'(a13a)') junk, data_format
close(13)
sep_file = trim(adjustl(sep_directory))//trim(adjustl(sep_file))
data_format = data_format(1:len_trim(adjustl(data_format))-1)
print *
print *, 'sep file specified in the header file is: ', trim(adjustl(sep_file))
print *, 'NX,NY,NZ = ', NX,NY,NZ
print *, 'OX,OY,OZ = ', OX,OY,OZ
print *, 'DX,DY,DZ = ', DX,DY,DZ
print *, 'esize = ', esize
print *, 'data_format = ', trim(adjustl(data_format))
print *, '*******************************************************************************'
print *
end subroutine READ_SEP_HEADER
!
!-----------------------------------------------------------------------------------
!
subroutine READ_SEP_MODEL_2D(sep_file,NX,NY,NZ,OX,OY,OZ,DX,DY,DZ,esize,data_format, &
model)
implicit none
!!! input
character(len=512),intent(in) :: sep_file ! sep file specified in sep header file
integer,intent(in) :: NX,NY,NZ ! n1,n2,n3 in sep header file
real,intent(in) :: OX,OY,OZ ! o1,o2,o3 in sep header file
real,intent(in) :: DX,DY,DZ ! d1,d2,d3 in sep header file
integer,intent(in) :: esize ! esize in sep header file
character(len=512),intent(in) :: data_format ! data_format in sep header file
!!! output
real,dimension(NX,NZ),intent(out) :: model
!!! local
integer :: ier,idummy
character(len=1) :: strdummy
! to avoid compiler warning
idummy = int(OX * OY * OZ)
idummy = int(DX * DY * DZ)
idummy = esize
strdummy(1:1) = data_format(1:1)
! check whether the model is 2D (NY==1)
if (NY /= 1) stop 'ERROR: this only works for 2D problems (NY/n2 must be 1)'
! note that we keep NY as general in the following (for 3D problems in the future)
open(unit=14,file=trim(adjustl(sep_file)),access='direct',status='old',recl=4*NX*NY*NZ,iostat=ier)
print *, '*******************************************************************************'
print *, 'reading sep file: '
print *, trim(adjustl(sep_file))
if (ier /= 0) stop 'ERROR: cannot open sep file'
read(14,rec=1,iostat=ier) model(:,:)
close(14)
if (ier /= 0) stop 'ERROR: reading sep file'
print *, 'done reading sucessfully'
print *, '*******************************************************************************'
print *
end subroutine READ_SEP_MODEL_2D
!
!-----------------------------------------------------------------------------------
!
subroutine interpolate_slowness_gll(vp,vs,rho,NX,NY,NZ,DX,DY,DZ,OX,OY,OZ,closest_point,nproc)
implicit none
!!! constants
integer, parameter :: NGLLX = 5, NGLLZ = 5
!double precision, parameter :: vp_water = 5500 ! upper bound
!!! input
integer,intent(in) :: NX,NY,NZ
real,intent(in) :: OX,OY,OZ
real,intent(in) :: DX,DY,DZ
real,dimension(NX,NZ),intent(in) :: vp,vs,rho
logical,intent(in) :: closest_point
integer,intent(in) :: nproc
!!! output
! currently no output, all information is saved in output file
! local parameters
integer :: iproc
integer :: ier,ix,iz,ix1,iz1,i,npoints
integer, dimension(NGLLX*NGLLZ) :: iglob
double precision, dimension(NGLLX*NGLLZ) :: rho_new,vp_new,vs_new,x,z
double precision :: rho_temp,vp_temp,vs_temp,a1,b1,tmp1,tmp2,tmp3
character(len=512) :: input_file,output_file
character(len=512) :: system_command
double precision :: vp_min,vp_max,vs_min,vs_max,rho_min,rho_max
real :: rdummy
! to avoid compiler warning
rdummy = OY * DY * NY
print *
print *, 'Interpolating GLL model values'
print *
! interpolating model values onto all GLL points
vp_min = 1.d30
vp_max = -1.d30
vs_min = 1.d30
vs_max = - 1.d30
rho_min = 1.d30
rho_max = - 1.d30
! for each process
do iproc = 1,nproc
! user output
print *, 'iproc = ',iproc
! input file
write(input_file,'(a,i6.6,a)') 'DATA/proc',iproc-1,'_model_velocity.dat_input'
open(unit=15,file=trim(adjustl(input_file)),status='old',iostat=ier)
if (ier /= 0) stop 'Error opening DATA/proc*****_model_velocity.dat_input file.'
! counts number of points
npoints = 0
do while (ier == 0)
! loops over GLL points of every element
do i = 1,NGLLX*NGLLZ
! reads model point
read(15,'(I10,5e15.5e4)',iostat=ier) iglob(i),x(i),z(i),rho_temp,vp_temp,vs_temp
if (ier /= 0) exit
npoints = npoints + 1
enddo
enddo
rewind(15)
print *, 'number of points = ',npoints
! number of lines must be a multiple of NGLLX * NGLLZ
if (mod(npoints,(NGLLX*NGLLZ)) /= 0) then
print *,'Error: invalid number of lines in input file '//trim(input_file)
print *,' number of points = ',npoints,' should be a multiple of NGLLX * NGLLZ = ',NGLLX,' * ',NGLLZ
stop 'Error invalid number of lines in input file'
endif
! output file
write(output_file,'(a,i6.6,a)') 'DATA/proc',iproc-1,'_model_velocity.dat_tmp'
open(unit=16,file=trim(adjustl(output_file)),status='unknown',iostat=ier)
if (ier /= 0) stop 'Error opening output file DATA/proc****_model_velocity.dat_tmp for interpolation.'
! interpolation
if (closest_point) then
! takes value from closest point
do while (ier == 0)
! loops over GLL points of every element
do i = 1,NGLLX*NGLLZ
! reads model point
read(15,'(I10,5e15.5e4)',iostat=ier) iglob(i),x(i),z(i),rho_temp,vp_temp,vs_temp
if (ier /= 0) exit
ix = NINT((x(i)-OX)/DX) + 1
iz = NINT((z(i)-OZ)/DZ) + 1
if (ix > NX) ix = NX
if (ix < 1) ix = 1
if (iz > NZ) iz = NZ
if (iz < 1) iz = 1
! closest point values
rho_new(i) = rho(ix,iz)
vp_new(i) = vp(ix,iz)
vs_new(i) = vs(ix,iz)
! statistics
if (rho_new(i) < rho_min) rho_min = rho_new(i)
if (rho_new(i) > rho_max) rho_max = rho_new(i)
if (vp_new(i) < vp_min) vp_min = vp_new(i)
if (vp_new(i) > vp_max) vp_max = vp_new(i)
if (vs_new(i) < vs_min) vs_min = vs_new(i)
if (vs_new(i) > vs_max) vs_max = vs_new(i)
enddo
! writes out interpolated values
! converting from slowness to wave speed values again & buoyancy to density
if (ier == 0) then
do i = 1,NGLLX*NGLLZ
write(16,'(I10,5e15.5e4)') iglob(i),x(i),z(i), rho_new(i), vp_new(i), vs_new(i)
enddo
endif
enddo
else
! bilinear interpolation
!
! (ix,iz) (ix+1,iz)
! x---------- x
! | * |
! | |
! | |
! x-----------x
! (ix,iz+1) (ix+1,iz+1)
!
do while (ier == 0)
! loops over GLL points of every element
do i = 1,NGLLX*NGLLZ
! reads model point
read(15,'(I10,5e15.5e4)',iostat=ier) iglob(i),x(i),z(i),rho_temp,vp_temp,vs_temp
if (ier /= 0) exit
ix = INT((x(i)-OX)/DX) + 1
iz = INT((z(i)-OZ)/DZ) + 1
ix1 = ix+1
iz1 = iz+1
if (ix1 > NX) ix1 = NX
if (iz1 > NZ) iz1 = NZ
a1 = x(i)-OX-(ix-1)*DX
a1 = a1/DX
b1 = z(i)-OZ-(iz-1)*DZ
b1 = b1/DZ
! note:
! - Vp: we use slowness values for interpolation. this garantees better interpolation results for traveltimes.
! - Vs: shear speeds might be zero in a coupled fluid-elastic model, so it is inverted directly.
! - density: it will be inverted to have "lightness" or "buoyancy" for interpolation.
! this is also done in a staggered-grid approach, see e.g.,
! Virieux 1987, "P-SV wave propagation in heterogeneous media: Velocity-stress finite-difference method",
! Geophysics, 51, p. 889 - 901.
! slowness
tmp1 = (1.d0-a1)*1.d0/vp(ix,iz) + a1*1.d0/vp(ix1,iz)
tmp2 = (1.d0-a1)*1.d0/vp(ix,iz1) + a1*1.d0/vp(ix1,iz1)
tmp3 = (1.d0-b1)*tmp1 + b1*tmp2
! converts back to Vp
vp_new(i) = 1.d0/tmp3
! Vs
tmp1 = (1.d0-a1)*vs(ix,iz) + a1*vs(ix1,iz)
tmp2 = (1.d0-a1)*vs(ix,iz1) + a1*vs(ix1,iz1)
tmp3 = (1.d0-b1)*tmp1 + b1*tmp2
vs_new(i) = tmp3
! buoyancy
tmp1 = (1.d0-a1)*1.d0/rho(ix,iz) + a1*1.d0/rho(ix1,iz)
tmp2 = (1.d0-a1)*1.d0/rho(ix,iz1) + a1*1.d0/rho(ix1,iz1)
tmp3 = (1.d0-b1)*tmp1 + b1*tmp2
! converts back to density
rho_new(i) = 1.d0/tmp3
! statistics
if (rho_new(i) < rho_min) rho_min = rho_new(i)
if (rho_new(i) > rho_max) rho_max = rho_new(i)
if (vp_new(i) < vp_min) vp_min = vp_new(i)
if (vp_new(i) > vp_max) vp_max = vp_new(i)
if (vs_new(i) < vs_min) vs_min = vs_new(i)
if (vs_new(i) > vs_max) vs_max = vs_new(i)
enddo
! writes out interpolated values
if (ier == 0) then
do i = 1,NGLLX*NGLLZ
write(16,'(I10,5e15.5e4)') iglob(i),x(i),z(i), rho_new(i), vp_new(i), vs_new(i)
enddo
endif
enddo
endif
! closes files
close(15)
close(16)
! replaces model file
print *, 'replacing model file'
write(system_command,*) 'mv -f ',trim(adjustl(input_file)),' ',trim(adjustl(input_file))//'.org'
write(*,*) trim(system_command)
call execute_command_line(trim(system_command),exitstat=ier)
if (ier /= 0) stop
write(system_command,*) 'mv -f ',trim(adjustl(output_file)),' ',trim(adjustl(input_file))
write(*,*) trim(system_command)
call execute_command_line(trim(system_command),exitstat=ier)
if (ier /= 0) stop
enddo ! iproc
! user output
print *
print *, 'interpolated output model:'
print *, ' Vp min/max = ',sngl(vp_min),sngl(vp_max)
print *, ' Vs min/max = ',sngl(vs_min),sngl(vs_max)
print *, ' rho min/max = ',sngl(rho_min),sngl(rho_max)
print *
end subroutine interpolate_slowness_gll