-
Notifications
You must be signed in to change notification settings - Fork 220
/
specfem3D.f90
3211 lines (2674 loc) · 121 KB
/
specfem3D.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
!=====================================================================
!
! S p e c f e m 3 D V e r s i o n 1 . 4
! ---------------------------------------
!
! Dimitri Komatitsch and Jeroen Tromp
! Seismological Laboratory - California Institute of Technology
! (c) California Institute of Technology September 2006
!
! This program is free software; you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation; either version 2 of the License, or
! (at your option) any later version.
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License along
! with this program; if not, write to the Free Software Foundation, Inc.,
! 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
!
!=====================================================================
!
! United States Government Sponsorship Acknowledged.
subroutine specfem3D
implicit none
include "constants.h"
! include values created by the mesher
include "OUTPUT_FILES/values_from_mesher.h"
!=============================================================================!
! !
! specfem3D is a 3-D spectral-element solver for a local or regional model. !
! It uses a mesh generated by program meshfem3D !
! !
!=============================================================================!
!
! If you use this code for your own research, please cite some of these articles:
!
! @ARTICLE{KoLiTrSuStSh04,
! author={Dimitri Komatitsch and Qinya Liu and Jeroen Tromp and Peter S\"{u}ss
! and Christiane Stidham and John H. Shaw},
! year=2004,
! title={Simulations of Ground Motion in the {L}os {A}ngeles {B}asin
! based upon the Spectral-Element Method},
! journal={Bull. Seism. Soc. Am.},
! volume=94,
! number=1,
! pages={187-206}}
!
! @ARTICLE{KoTr99,
! author={D. Komatitsch and J. Tromp},
! year=1999,
! title={Introduction to the spectral-element method for 3-{D} seismic wave propagation},
! journal={Geophys. J. Int.},
! volume=139,
! number=3,
! pages={806-822},
! doi={10.1046/j.1365-246x.1999.00967.x}}
!
! @ARTICLE{KoVi98,
! author={D. Komatitsch and J. P. Vilotte},
! title={The spectral-element method: an efficient tool to simulate the seismic response of 2{D} and 3{D} geological structures},
! journal={Bull. Seismol. Soc. Am.},
! year=1998,
! volume=88,
! number=2,
! pages={368-392}}
!
! If you use the kernel capabilities of the code, please cite
!
! @ARTICLE{LiTr06,
! author={Qinya Liu and Jeroen Tromp},
! title={Finite-frequency kernels based on adjoint methods},
! journal={Bull. Seismol. Soc. Am.},
! year=2006,
! volume=96,
! number=6,
! pages={2383-2397},
! doi={10.1785/0120060041}}
!
! Reference frame - convention:
! ----------------------------
!
! The code uses the following convention for the reference frame:
!
! - X axis is East
! - Y axis is North
! - Z axis is up
!
! Note that this convention is different from both the Aki-Richards convention
! and the Harvard CMT convention.
!
! Let us recall that the Aki-Richards convention is:
!
! - X axis is North
! - Y axis is East
! - Z axis is down
!
! and that the Harvard CMT convention is:
!
! - X axis is South
! - Y axis is East
! - Z axis is up
!
! To report bugs or suggest improvements to the code, please send an email
! to Jeroen Tromp <jtromp AT caltech.edu> and/or use our online
! bug tracking system at http://www.geodynamics.org/roundup .
!
! Evolution of the code:
! ---------------------
!
! MPI v. 1.4 Dimitri Komatitsch, University of Pau, Qinya Liu and others, Caltech, September 2006:
! better adjoint and kernel calculations, faster and better I/Os
! on very large systems, new Pyre version, many small improvements and bug fixes
! MPI v. 1.3 Dimitri Komatitsch, University of Pau, and Qinya Liu, Caltech, July 2005:
! serial version, regular mesh, adjoint and kernel calculations, ParaView support
! MPI v. 1.2 Min Chen and Dimitri Komatitsch, Caltech, July 2004:
! full anisotropy, volume movie
! MPI v. 1.1 Dimitri Komatitsch, Caltech, October 2002: Zhu's Moho map, scaling
! of Vs with depth, Hauksson's regional model, attenuation, oceans, movies
! MPI v. 1.0 Dimitri Komatitsch, Caltech, May 2002: first MPI version
! based on global code
! memory variables and standard linear solids for attenuation
double precision, dimension(N_SLS) :: tau_mu_dble,tau_sigma_dble,beta_dble
double precision factor_scale_dble,one_minus_sum_beta_dble
real(kind=CUSTOM_REAL), dimension(NUM_REGIONS_ATTENUATION,N_SLS) :: tau_mu,tau_sigma,beta
real(kind=CUSTOM_REAL), dimension(NUM_REGIONS_ATTENUATION) :: factor_scale,one_minus_sum_beta
real(kind=CUSTOM_REAL), dimension(NUM_REGIONS_ATTENUATION,N_SLS) :: tauinv,factor_common, alphaval,betaval,gammaval
integer iattenuation
double precision scale_factor
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_ATTENUATION,N_SLS) :: &
R_xx,R_yy,R_xy,R_xz,R_yz
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_ATTENUATION) :: &
epsilondev_xx,epsilondev_yy,epsilondev_xy,epsilondev_xz,epsilondev_yz
! ADJOINT
real(kind=CUSTOM_REAL), dimension(NUM_REGIONS_ATTENUATION,N_SLS) :: b_alphaval, b_betaval, b_gammaval
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_ATT_AND_KERNEL,N_SLS) :: &
b_R_xx,b_R_yy,b_R_xy,b_R_xz,b_R_yz
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_ATT_AND_KERNEL) :: b_epsilondev_xx, &
b_epsilondev_yy,b_epsilondev_xy,b_epsilondev_xz,b_epsilondev_yz
! ADJOINT
integer NPOIN2DMAX_XY
! use integer array to store topography values
integer NX_TOPO,NY_TOPO
double precision ORIG_LAT_TOPO,ORIG_LONG_TOPO,DEGREES_PER_CELL_TOPO
character(len=100) topo_file
integer, dimension(:,:), allocatable :: itopo_bathy
integer, dimension(NSPEC2DMAX_XMIN_XMAX_VAL) :: ibelm_xmin,ibelm_xmax
integer, dimension(NSPEC2DMAX_YMIN_YMAX_VAL) :: ibelm_ymin,ibelm_ymax
integer, dimension(NSPEC2D_BOTTOM_VAL) :: ibelm_bottom
integer, dimension(NSPEC2D_TOP_VAL) :: ibelm_top
real(kind=CUSTOM_REAL), dimension(NGLLY,NGLLZ,NSPEC2DMAX_XMIN_XMAX_VAL) :: jacobian2D_xmin,jacobian2D_xmax
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLZ,NSPEC2DMAX_YMIN_YMAX_VAL) :: jacobian2D_ymin,jacobian2D_ymax
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NSPEC2D_BOTTOM_VAL) :: jacobian2D_bottom
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NSPEC2D_TOP_VAL) :: jacobian2D_top
real(kind=CUSTOM_REAL), dimension(NDIM,NGLLY,NGLLZ,NSPEC2DMAX_XMIN_XMAX_VAL) :: normal_xmin,normal_xmax
real(kind=CUSTOM_REAL), dimension(NDIM,NGLLX,NGLLZ,NSPEC2DMAX_YMIN_YMAX_VAL) :: normal_ymin,normal_ymax
real(kind=CUSTOM_REAL), dimension(NDIM,NGLLX,NGLLY,NSPEC2D_BOTTOM_VAL) :: normal_bottom
real(kind=CUSTOM_REAL), dimension(NDIM,NGLLX,NGLLY,NSPEC2D_TOP_VAL) :: normal_top
! Moho mesh
integer,dimension(NSPEC2D_MOHO_BOUN) :: ibelm_moho_top, ibelm_moho_bot
real(CUSTOM_REAL), dimension(NDIM,NGLLX,NGLLY,NSPEC2D_MOHO_BOUN) :: normal_moho
integer :: nspec2D_moho, njunk
logical, dimension(NSPEC_BOUN) :: is_moho_top, is_moho_bot
! buffers for send and receive between faces of the slices and the chunks
real(kind=CUSTOM_REAL), dimension(NPOIN2DMAX_XY_VAL) :: buffer_send_faces_scalar,buffer_received_faces_scalar
real(kind=CUSTOM_REAL), dimension(NDIM,NPOIN2DMAX_XY_VAL) :: buffer_send_faces_vector,buffer_received_faces_vector
! -----------------
! mesh parameters
integer, dimension(NGLLX,NGLLY,NGLLZ,NSPEC_AB_VAL) :: ibool
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_AB_VAL) :: &
xix,xiy,xiz,etax,etay,etaz,gammax,gammay,gammaz,jacobian
real(kind=CUSTOM_REAL), dimension(NGLOB_AB_VAL) :: xstore,ystore,zstore
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_AB_VAL) :: &
kappastore,mustore
! material properties in case of a fully anisotropic material
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_ANISO) :: &
c11store,c12store,c13store,c14store,c15store,c16store,c22store, &
c23store,c24store,c25store,c26store,c33store,c34store,c35store, &
c36store,c44store,c45store,c46store,c55store,c56store,c66store
! flag for sediments
logical not_fully_in_bedrock(NSPEC_AB_VAL)
logical, dimension(NGLLX,NGLLY,NGLLZ,NSPEC_AB_VAL) :: flag_sediments
! Stacey
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_AB_VAL) :: rho_vp,rho_vs
! local to global mapping
integer, dimension(NSPEC_AB_VAL) :: idoubling
! mass matrix
real(kind=CUSTOM_REAL), dimension(NGLOB_AB_VAL) :: rmass
! additional mass matrix for ocean load
! ocean load mass matrix is always allocated statically even if no oceans
real(kind=CUSTOM_REAL), dimension(NGLOB_AB_VAL) :: rmass_ocean_load
logical, dimension(NGLOB_AB_VAL) :: updated_dof_ocean_load
real(kind=CUSTOM_REAL) additional_term,force_normal_comp
! displacement, velocity, acceleration
real(kind=CUSTOM_REAL), dimension(NDIM,NGLOB_AB_VAL) :: displ,veloc,accel
real(kind=CUSTOM_REAL) xixl,xiyl,xizl,etaxl,etayl,etazl,gammaxl,gammayl,gammazl,jacobianl
real(kind=CUSTOM_REAL) duxdxl,duxdyl,duxdzl,duydxl,duydyl,duydzl,duzdxl,duzdyl,duzdzl
real(kind=CUSTOM_REAL) duxdxl_plus_duydyl,duxdxl_plus_duzdzl,duydyl_plus_duzdzl
real(kind=CUSTOM_REAL) duxdyl_plus_duydxl,duzdxl_plus_duxdzl,duzdyl_plus_duydzl
real(kind=CUSTOM_REAL) sigma_xx,sigma_yy,sigma_zz,sigma_xy,sigma_xz,sigma_yz
real(kind=CUSTOM_REAL) hp1,hp2,hp3
real(kind=CUSTOM_REAL) fac1,fac2,fac3
real(kind=CUSTOM_REAL) lambdal,kappal,mul,lambdalplus2mul
real(kind=CUSTOM_REAL) c11,c12,c13,c14,c15,c16,c22,c23,c24,c25,c26,c33,c34,c35,c36,c44,c45,c46,c55,c56,c66
real(kind=CUSTOM_REAL) tempx1l,tempx2l,tempx3l
real(kind=CUSTOM_REAL) tempy1l,tempy2l,tempy3l
real(kind=CUSTOM_REAL) tempz1l,tempz2l,tempz3l
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ) :: &
tempx1,tempx2,tempx3,tempy1,tempy2,tempy3,tempz1,tempz2,tempz3
! time scheme
real(kind=CUSTOM_REAL) deltat,deltatover2,deltatsqover2
! ADJOINT
real(kind=CUSTOM_REAL) b_additional_term,b_force_normal_comp, kappa_k, mu_k
real(kind=CUSTOM_REAL), dimension(NDIM,NGLOB_ADJOINT) :: b_displ, b_veloc, b_accel
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ,NSPEC_ADJOINT) :: rho_kl, mu_kl, kappa_kl, &
rhop_kl, beta_kl, alpha_kl
real(kind=CUSTOM_REAL) dsxx,dsxy,dsxz,dsyy,dsyz,dszz
real(kind=CUSTOM_REAL) b_duxdxl,b_duxdyl,b_duxdzl,b_duydxl,b_duydyl,b_duydzl,b_duzdxl,b_duzdyl,b_duzdzl
real(kind=CUSTOM_REAL) b_duxdxl_plus_duydyl,b_duxdxl_plus_duzdzl,b_duydyl_plus_duzdzl
real(kind=CUSTOM_REAL) b_duxdyl_plus_duydxl,b_duzdxl_plus_duxdzl,b_duzdyl_plus_duydzl
real(kind=CUSTOM_REAL) b_dsxx,b_dsxy,b_dsxz,b_dsyy,b_dsyz,b_dszz
real(kind=CUSTOM_REAL) b_sigma_xx,b_sigma_yy,b_sigma_zz,b_sigma_xy,b_sigma_xz,b_sigma_yz
real(kind=CUSTOM_REAL) b_tempx1l,b_tempx2l,b_tempx3l
real(kind=CUSTOM_REAL) b_tempy1l,b_tempy2l,b_tempy3l
real(kind=CUSTOM_REAL) b_tempz1l,b_tempz2l,b_tempz3l
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ) :: &
b_tempx1,b_tempx2,b_tempx3,b_tempy1,b_tempy2,b_tempy3,b_tempz1,b_tempz2,b_tempz3
real(kind=CUSTOM_REAL), dimension(:,:,:,:), allocatable :: absorb_xmin, absorb_xmax, &
absorb_ymin, absorb_ymax, absorb_zmin ! for absorbing b.c.
integer reclen_xmin, reclen_xmax, reclen_ymin, reclen_ymax, reclen_zmin, reclen1, reclen2
real(kind=CUSTOM_REAL) b_deltat, b_deltatover2, b_deltatsqover2
! ADJOINT
! for attenuation
real(kind=CUSTOM_REAL) R_xx_val,R_yy_val
real(kind=CUSTOM_REAL) factor_loc,alphaval_loc,betaval_loc,gammaval_loc,Sn,Snp1
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ) :: epsilondev_xx_loc, &
epsilondev_yy_loc, epsilondev_xy_loc, epsilondev_xz_loc, epsilondev_yz_loc
real(kind=CUSTOM_REAL) epsilon_trace_over_3
! ADJOINT
real(kind=CUSTOM_REAL) b_R_xx_val,b_R_yy_val
real(kind=CUSTOM_REAL) b_alphaval_loc,b_betaval_loc,b_gammaval_loc,b_Sn,b_Snp1
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ) :: b_epsilondev_xx_loc, &
b_epsilondev_yy_loc, b_epsilondev_xy_loc, b_epsilondev_xz_loc, b_epsilondev_yz_loc
real(kind=CUSTOM_REAL) b_epsilon_trace_over_3
! ADJOINT
integer l
integer i_SLS
real(kind=CUSTOM_REAL) one_minus_sum_beta_use,minus_sum_beta,vs_val,Q_mu
integer iselected,iattenuation_sediments,int_Q_mu
! Moho kernel
integer ispec2D_moho_top, ispec2D_moho_bot, k_top, k_bot, ispec_top, ispec_bot, iglob_top, iglob_bot
real(kind=CUSTOM_REAL), dimension(NDIM,NDIM,NGLLX,NGLLY,NGLLZ,NSPEC2D_MOHO_BOUN) :: dsdx_top, dsdx_bot, b_dsdx_top, b_dsdx_bot
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NSPEC2D_MOHO_BOUN) :: moho_kl
real(kind=CUSTOM_REAL) :: kernel_moho_top, kernel_moho_bot
! --------
! parameters for the source
integer it,isource
integer, dimension(:), allocatable :: islice_selected_source,ispec_selected_source
integer yr,jda,ho,mi
real(kind=CUSTOM_REAL) stf_used
real(kind=CUSTOM_REAL), dimension(:,:,:,:), allocatable :: sourcearray
real(kind=CUSTOM_REAL), dimension(:,:,:,:,:), allocatable :: sourcearrays
!ADJOINT
character(len=150) adj_source_file
real(kind=CUSTOM_REAL), dimension(:,:,:,:,:), allocatable :: adj_sourcearray
real(kind=CUSTOM_REAL), dimension(:,:,:,:,:,:), allocatable :: adj_sourcearrays
!ADJOINT
double precision sec,stf
double precision, dimension(:), allocatable :: Mxx,Myy,Mzz,Mxy,Mxz,Myz
double precision, dimension(:), allocatable :: xi_source,eta_source,gamma_source
double precision, dimension(:), allocatable :: t_cmt,hdur,hdur_gaussian
double precision, dimension(:), allocatable :: utm_x_source,utm_y_source
double precision, external :: comp_source_time_function
double precision :: t0
! receiver information
character(len=150) rec_filename,filtered_rec_filename,dummystring
integer nrec,nrec_local,nrec_tot_found,irec_local,ios
integer, allocatable, dimension(:) :: islice_selected_rec,ispec_selected_rec,number_receiver_global
double precision, allocatable, dimension(:) :: xi_receiver,eta_receiver,gamma_receiver
double precision hlagrange
! ADJOINT
integer nrec_simulation, nadj_rec_local
! source frechet derivatives
real(kind=CUSTOM_REAL) :: displ_s(NDIM,NGLLX,NGLLY,NGLLZ), eps_s(NDIM,NDIM), eps_m_s(NDIM), stf_deltat
real(kind=CUSTOM_REAL), dimension(:), allocatable :: Mxx_der,Myy_der,Mzz_der,Mxy_der,Mxz_der,Myz_der
real(kind=CUSTOM_REAL), dimension(:,:), allocatable :: sloc_der
double precision, dimension(:,:), allocatable :: hpxir_store,hpetar_store,hpgammar_store
! ADJOINT
! timing information for the stations
double precision, allocatable, dimension(:,:,:) :: nu
character(len=MAX_LENGTH_STATION_NAME), allocatable, dimension(:) :: station_name
character(len=MAX_LENGTH_NETWORK_NAME), allocatable, dimension(:) :: network_name
! seismograms
double precision dxd,dyd,dzd,vxd,vyd,vzd,axd,ayd,azd
real(kind=CUSTOM_REAL), dimension(:,:,:), allocatable :: seismograms_d,seismograms_v,seismograms_a
real(kind=CUSTOM_REAL), dimension(:,:,:,:), allocatable :: seismograms_eps
integer i,j,k,ispec,irec,iglob
! Gauss-Lobatto-Legendre points of integration and weights
double precision, dimension(NGLLX) :: xigll,wxgll
double precision, dimension(NGLLY) :: yigll,wygll
double precision, dimension(NGLLZ) :: zigll,wzgll
! array with derivatives of Lagrange polynomials and precalculated products
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLX) :: hprime_xx,hprimewgll_xx
real(kind=CUSTOM_REAL), dimension(NGLLY,NGLLY) :: hprime_yy,hprimewgll_yy
real(kind=CUSTOM_REAL), dimension(NGLLZ,NGLLZ) :: hprime_zz,hprimewgll_zz
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY) :: wgllwgll_xy
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLZ) :: wgllwgll_xz
real(kind=CUSTOM_REAL), dimension(NGLLY,NGLLZ) :: wgllwgll_yz
! Lagrange interpolators at receivers
double precision, dimension(:), allocatable :: hxir,hetar,hpxir,hpetar,hgammar,hpgammar
double precision, dimension(:,:), allocatable :: hxir_store,hetar_store,hgammar_store
! 2-D addressing and buffers for summation between slices
integer, dimension(NPOIN2DMAX_XMIN_XMAX_VAL) :: iboolleft_xi,iboolright_xi
integer, dimension(NPOIN2DMAX_YMIN_YMAX_VAL) :: iboolleft_eta,iboolright_eta
! for addressing of the slices
integer, dimension(0:NPROC_XI_VAL-1,0:NPROC_ETA_VAL) :: addressing
integer, dimension(0:NPROC_VAL-1) :: iproc_xi_slice,iproc_eta_slice
! proc numbers for MPI
integer myrank,sizeprocs
integer npoin2D_xi,npoin2D_eta
integer iproc_xi,iproc_eta,iproc,iproc_read
! maximum of the norm of the displacement
real(kind=CUSTOM_REAL) Usolidnorm,Usolidnorm_all
! ADJOINT
real(kind=CUSTOM_REAL) b_Usolidnorm, b_Usolidnorm_all
! ADJOINT
! timer MPI
double precision, external :: wtime
integer ihours,iminutes,iseconds,int_tCPU
double precision time_start,tCPU
! parameters read from parameter file
integer NER_SEDIM,NER_BASEMENT_SEDIM,NER_16_BASEMENT, &
NER_MOHO_16,NER_BOTTOM_MOHO,NEX_XI,NEX_ETA, &
NPROC_XI,NPROC_ETA,NTSTEP_BETWEEN_OUTPUT_SEISMOS,NSTEP,UTM_PROJECTION_ZONE,SIMULATION_TYPE
integer NSOURCES
double precision UTM_X_MIN,UTM_X_MAX,UTM_Y_MIN,UTM_Y_MAX,Z_DEPTH_BLOCK
double precision DT,LATITUDE_MIN,LATITUDE_MAX,LONGITUDE_MIN,LONGITUDE_MAX,HDUR_MOVIE
double precision THICKNESS_TAPER_BLOCK_HR,THICKNESS_TAPER_BLOCK_MR,VP_MIN_GOCAD,VP_VS_RATIO_GOCAD_TOP,VP_VS_RATIO_GOCAD_BOTTOM
logical HARVARD_3D_GOCAD_MODEL,TOPOGRAPHY,ATTENUATION,USE_OLSEN_ATTENUATION, &
OCEANS,IMPOSE_MINIMUM_VP_GOCAD,HAUKSSON_REGIONAL_MODEL, &
BASEMENT_MAP,MOHO_MAP_LUPEI,ABSORBING_CONDITIONS,SAVE_FORWARD
logical ANISOTROPY,SAVE_MESH_FILES,PRINT_SOURCE_TIME_FUNCTION
logical MOVIE_SURFACE,MOVIE_VOLUME,CREATE_SHAKEMAP,SAVE_DISPLACEMENT, &
USE_HIGHRES_FOR_MOVIES,SUPPRESS_UTM_PROJECTION,USE_REGULAR_MESH
integer NTSTEP_BETWEEN_FRAMES,NTSTEP_BETWEEN_OUTPUT_INFO
character(len=150) OUTPUT_FILES,LOCAL_PATH,prname,prname_Q,MODEL
! parameters deduced from parameters read from file
integer NPROC,NEX_PER_PROC_XI,NEX_PER_PROC_ETA
integer NER
integer NSPEC2D_A_XI,NSPEC2D_B_XI, &
NSPEC2D_A_ETA,NSPEC2D_B_ETA, &
NSPEC2DMAX_XMIN_XMAX,NSPEC2DMAX_YMIN_YMAX, &
NSPEC2D_BOTTOM,NSPEC2D_TOP, &
NPOIN2DMAX_XMIN_XMAX,NPOIN2DMAX_YMIN_YMAX, &
NSPEC_AB, NGLOB_AB
! names of the data files for all the processors in MPI
character(len=150) outputname
! Stacey conditions put back
integer nspec2D_xmin,nspec2D_xmax,nspec2D_ymin,nspec2D_ymax,ispec2D
integer, dimension(2,NSPEC2DMAX_YMIN_YMAX_VAL) :: nimin,nimax,nkmin_eta
integer, dimension(2,NSPEC2DMAX_XMIN_XMAX_VAL) :: njmin,njmax,nkmin_xi
real(kind=CUSTOM_REAL) vx,vy,vz,nx,ny,nz,tx,ty,tz,vn,weight
! to save movie frames
integer ipoin, nmovie_points, iloc, iorderi(NGNOD2D), iorderj(NGNOD2D)
real(kind=CUSTOM_REAL), dimension(:), allocatable :: &
store_val_x,store_val_y,store_val_z, &
store_val_ux,store_val_uy,store_val_uz, &
store_val_norm_displ,store_val_norm_veloc,store_val_norm_accel
real(kind=CUSTOM_REAL), dimension(:,:), allocatable :: &
store_val_x_all,store_val_y_all,store_val_z_all, &
store_val_ux_all,store_val_uy_all,store_val_uz_all
! to save full 3D snapshot of velocity
real(kind=CUSTOM_REAL), dimension(NGLLX,NGLLY,NGLLZ) :: dvxdxl,dvxdyl,dvxdzl,dvydxl,dvydyl,dvydzl,dvzdxl,dvzdyl,dvzdzl
real(kind=CUSTOM_REAL), dimension(:,:,:,:),allocatable:: div, curl_x, curl_y, curl_z
! ************** PROGRAM STARTS HERE **************
! sizeprocs returns number of processes started
! (should be equal to NPROC)
! myrank is the rank of each process, between 0 and sizeprocs-1.
! as usual in MPI, process 0 is in charge of coordinating everything
! and also takes care of the main output
call world_size(sizeprocs)
call world_rank(myrank)
! read the parameter file
call read_parameter_file(LATITUDE_MIN,LATITUDE_MAX,LONGITUDE_MIN,LONGITUDE_MAX, &
UTM_X_MIN,UTM_X_MAX,UTM_Y_MIN,UTM_Y_MAX,Z_DEPTH_BLOCK, &
NER_SEDIM,NER_BASEMENT_SEDIM,NER_16_BASEMENT,NER_MOHO_16,NER_BOTTOM_MOHO, &
NEX_XI,NEX_ETA,NPROC_XI,NPROC_ETA,NTSTEP_BETWEEN_OUTPUT_SEISMOS,NSTEP,UTM_PROJECTION_ZONE,DT, &
ATTENUATION,USE_OLSEN_ATTENUATION,HARVARD_3D_GOCAD_MODEL,TOPOGRAPHY,LOCAL_PATH,NSOURCES, &
THICKNESS_TAPER_BLOCK_HR,THICKNESS_TAPER_BLOCK_MR,VP_MIN_GOCAD,VP_VS_RATIO_GOCAD_TOP,VP_VS_RATIO_GOCAD_BOTTOM, &
OCEANS,IMPOSE_MINIMUM_VP_GOCAD,HAUKSSON_REGIONAL_MODEL,ANISOTROPY, &
BASEMENT_MAP,MOHO_MAP_LUPEI,ABSORBING_CONDITIONS, &
MOVIE_SURFACE,MOVIE_VOLUME,CREATE_SHAKEMAP,SAVE_DISPLACEMENT, &
NTSTEP_BETWEEN_FRAMES,USE_HIGHRES_FOR_MOVIES,HDUR_MOVIE, &
SAVE_MESH_FILES,PRINT_SOURCE_TIME_FUNCTION, &
NTSTEP_BETWEEN_OUTPUT_INFO,SUPPRESS_UTM_PROJECTION,MODEL,USE_REGULAR_MESH,SIMULATION_TYPE,SAVE_FORWARD)
if (sizeprocs == 1 .and. (NPROC_XI /= 1 .or. NPROC_ETA /= 1)) then
stop 'must have NPROC_XI = NPROC_ETA = 1 for a serial run'
endif
! check simulation type
if (SIMULATION_TYPE /= 1 .and. SIMULATION_TYPE /= 2 .and. SIMULATION_TYPE /= 3) &
call exit_mpi(myrank,'SIMULATION_TYPE can only be 1, 2, or 3')
! check simulation parameters
if (SIMULATION_TYPE /= 1 .and. NSOURCES > 1000) call exit_mpi(myrank, 'for adjoint simulations, NSOURCES <= 1000')
! LQY -- note: kernel simulations with attenuation turned on has been implemented
! compute other parameters based upon values read
call compute_parameters(NER,NEX_XI,NEX_ETA,NPROC_XI,NPROC_ETA, &
NPROC,NEX_PER_PROC_XI,NEX_PER_PROC_ETA, &
NER_BOTTOM_MOHO,NER_MOHO_16,NER_16_BASEMENT,NER_BASEMENT_SEDIM,NER_SEDIM, &
NSPEC_AB,NSPEC2D_A_XI,NSPEC2D_B_XI, &
NSPEC2D_A_ETA,NSPEC2D_B_ETA, &
NSPEC2DMAX_XMIN_XMAX,NSPEC2DMAX_YMIN_YMAX,NSPEC2D_BOTTOM,NSPEC2D_TOP, &
NPOIN2DMAX_XMIN_XMAX,NPOIN2DMAX_YMIN_YMAX,NGLOB_AB,USE_REGULAR_MESH)
! get the base pathname for output files
call get_value_string(OUTPUT_FILES, 'OUTPUT_FILES', 'OUTPUT_FILES')
! open main output file, only written to by process 0
if(myrank == 0 .and. IMAIN /= ISTANDARD_OUTPUT) &
open(unit=IMAIN,file=trim(OUTPUT_FILES)//'/output_solver.txt',status='unknown')
if(myrank == 0) then
write(IMAIN,*)
write(IMAIN,*) '**********************************************'
write(IMAIN,*) '**** Specfem 3-D Solver - MPI version f90 ****'
write(IMAIN,*) '**********************************************'
write(IMAIN,*)
write(IMAIN,*)
if(FIX_UNDERFLOW_PROBLEM) write(IMAIN,*) 'Fixing slow underflow trapping problem using small initial field'
write(IMAIN,*)
write(IMAIN,*) 'There are ',sizeprocs,' MPI processes'
write(IMAIN,*) 'Processes are numbered from 0 to ',sizeprocs-1
write(IMAIN,*)
write(IMAIN,*) 'There are ',NEX_XI,' elements along xi'
write(IMAIN,*) 'There are ',NEX_ETA,' elements along eta'
write(IMAIN,*)
write(IMAIN,*) 'There are ',NPROC_XI,' slices along xi'
write(IMAIN,*) 'There are ',NPROC_ETA,' slices along eta'
write(IMAIN,*) 'There is a total of ',NPROC,' slices'
write(IMAIN,*)
write(IMAIN,*) ' NDIM = ',NDIM
write(IMAIN,*)
write(IMAIN,*) ' NGLLX = ',NGLLX
write(IMAIN,*) ' NGLLY = ',NGLLY
write(IMAIN,*) ' NGLLZ = ',NGLLZ
write(IMAIN,*)
! write information about precision used for floating-point operations
if(CUSTOM_REAL == SIZE_REAL) then
write(IMAIN,*) 'using single precision for the calculations'
else
write(IMAIN,*) 'using double precision for the calculations'
endif
write(IMAIN,*)
write(IMAIN,*) 'smallest and largest possible floating-point numbers are: ',tiny(1._CUSTOM_REAL),huge(1._CUSTOM_REAL)
write(IMAIN,*)
endif
! check that the code is running with the requested nb of processes
if(sizeprocs /= NPROC) call exit_MPI(myrank,'wrong number of MPI processes')
! check that we have at least one source
if(NSOURCES < 1) call exit_MPI(myrank,'need at least one source')
! open file with global slice number addressing
if(myrank == 0) then
open(unit=IIN,file=trim(OUTPUT_FILES)//'/addressing.txt',status='old',action='read')
do iproc = 0,NPROC-1
read(IIN,*) iproc_read,iproc_xi,iproc_eta
if(iproc_read /= iproc) call exit_MPI(myrank,'incorrect slice number read')
addressing(iproc_xi,iproc_eta) = iproc
iproc_xi_slice(iproc) = iproc_xi
iproc_eta_slice(iproc) = iproc_eta
enddo
close(IIN)
endif
! broadcast the information read on the master to the nodes
call bcast_all_i(addressing,NPROC_XI*NPROC_ETA)
call bcast_all_i(iproc_xi_slice,NPROC)
call bcast_all_i(iproc_eta_slice,NPROC)
! determine local slice coordinates using addressing
iproc_xi = iproc_xi_slice(myrank)
iproc_eta = iproc_eta_slice(myrank)
! define maximum size for message buffers
NPOIN2DMAX_XY = max(NPOIN2DMAX_XMIN_XMAX,NPOIN2DMAX_YMIN_YMAX)
! start reading the databases
! read arrays created by the mesher
call read_arrays_solver(myrank,xstore,ystore,zstore, &
xix,xiy,xiz,etax,etay,etaz,gammax,gammay,gammaz,jacobian, &
flag_sediments,not_fully_in_bedrock,rho_vp,rho_vs,ANISOTROPY, &
c11store,c12store,c13store,c14store,c15store,c16store,c22store, &
c23store,c24store,c25store,c26store,c33store,c34store,c35store, &
c36store,c44store,c45store,c46store,c55store,c56store,c66store, &
kappastore,mustore,ibool,idoubling,rmass,rmass_ocean_load,LOCAL_PATH,OCEANS)
! check that the number of points in this slice is correct
if(minval(ibool(:,:,:,:)) /= 1 .or. maxval(ibool(:,:,:,:)) /= NGLOB_AB_VAL) &
call exit_MPI(myrank,'incorrect global numbering: iboolmax does not equal NGLOB')
! read 2-D addressing for summation between slices with MPI
call read_arrays_buffers_solver(myrank,iboolleft_xi, &
iboolright_xi,iboolleft_eta,iboolright_eta, &
npoin2D_xi,npoin2D_eta, &
NPOIN2DMAX_XMIN_XMAX,NPOIN2DMAX_YMIN_YMAX,LOCAL_PATH)
! $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
if(myrank == 0) then
write(IMAIN,*) '******************************************'
write(IMAIN,*) 'There are ',NEX_XI,' elements along xi'
write(IMAIN,*) 'There are ',NEX_ETA,' elements along eta'
write(IMAIN,*)
write(IMAIN,*) 'There are ',NPROC_XI,' slices along xi'
write(IMAIN,*) 'There are ',NPROC_ETA,' slices along eta'
write(IMAIN,*) 'There is a total of ',NPROC,' slices'
write(IMAIN,*) '******************************************'
write(IMAIN,*)
endif
! set up GLL points, weights and derivation matrices
call define_derivation_matrices(xigll,yigll,zigll,wxgll,wygll,wzgll, &
hprime_xx,hprime_yy,hprime_zz, &
hprimewgll_xx,hprimewgll_yy,hprimewgll_zz, &
wgllwgll_xy,wgllwgll_xz,wgllwgll_yz)
! allocate 1-D Lagrange interpolators and derivatives
allocate(hxir(NGLLX))
allocate(hpxir(NGLLX))
allocate(hetar(NGLLY))
allocate(hpetar(NGLLY))
allocate(hgammar(NGLLZ))
allocate(hpgammar(NGLLZ))
! create name of database
call create_name_database(prname,myrank,LOCAL_PATH)
if (ATTENUATION .and. ((SIMULATION_TYPE == 1 .and. SAVE_FORWARD) .or. SIMULATION_TYPE == 3)) &
call create_name_database(prname_Q,myrank,LOCAL_PATH_Q)
! boundary parameters
open(unit=27,file=prname(1:len_trim(prname))//'ibelm.bin',status='old',action='read',form='unformatted')
read(27) ibelm_xmin
read(27) ibelm_xmax
read(27) ibelm_ymin
read(27) ibelm_ymax
read(27) ibelm_bottom
read(27) ibelm_top
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'normal.bin',status='old',action='read',form='unformatted')
read(27) normal_xmin
read(27) normal_xmax
read(27) normal_ymin
read(27) normal_ymax
read(27) normal_bottom
read(27) normal_top
close(27)
! moho boundary
if (SAVE_MOHO_MESH .and. SIMULATION_TYPE == 3) then
moho_kl = ZERO
open(unit=27,file=prname(1:len_trim(prname))//'ibelm_moho.bin',status='unknown',form='unformatted')
read(27) nspec2D_moho
read(27) njunk
read(27) njunk
read(27) ibelm_moho_top
read(27) ibelm_moho_bot
close(27)
if (nspec2D_moho /= NSPEC2D_BOTTOM) call exit_mpi(myrank, "nspec2D_moho /= NSPEC2D_BOTTOM for Moho mesh")
open(unit=27,file=prname(1:len_trim(prname))//'normal_moho.bin',status='old',form='unformatted')
read(27) normal_moho
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'is_moho.bin',status='old',form='unformatted')
read(27) is_moho_top
read(27) is_moho_bot
close(27)
endif
! Stacey put back
open(unit=27,file=prname(1:len_trim(prname))//'nspec2D.bin',status='unknown',form='unformatted')
read(27) nspec2D_xmin
read(27) nspec2D_xmax
read(27) nspec2D_ymin
read(27) nspec2D_ymax
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'jacobian2D.bin',status='old',action='read',form='unformatted')
read(27) jacobian2D_xmin
read(27) jacobian2D_xmax
read(27) jacobian2D_ymin
read(27) jacobian2D_ymax
read(27) jacobian2D_bottom
read(27) jacobian2D_top
close(27)
! Stacey put back
! read arrays for Stacey conditions
if(ABSORBING_CONDITIONS) then
open(unit=27,file=prname(1:len_trim(prname))//'nimin.bin',status='unknown',form='unformatted')
read(27) nimin
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'nimax.bin',status='unknown',form='unformatted')
read(27) nimax
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'njmin.bin',status='unknown',form='unformatted')
read(27) njmin
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'njmax.bin',status='unknown',form='unformatted')
read(27) njmax
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'nkmin_xi.bin',status='unknown',form='unformatted')
read(27) nkmin_xi
close(27)
open(unit=27,file=prname(1:len_trim(prname))//'nkmin_eta.bin',status='unknown',form='unformatted')
read(27) nkmin_eta
close(27)
! read in absorbing wavefield saved by forward simulations
if (nspec2D_xmin > 0 .and. (SIMULATION_TYPE == 3 .or. (SIMULATION_TYPE == 1 .and. SAVE_FORWARD))) then
allocate(absorb_xmin(NDIM,NGLLY,NGLLZ,nspec2D_xmin))
reclen_xmin = CUSTOM_REAL * (NDIM * NGLLY * NGLLZ * nspec2D_xmin)
if (SIMULATION_TYPE == 3) then
open(unit=31,file=trim(prname)//'absorb_xmin.bin',status='old',action='read',form='unformatted',access='direct', &
recl=reclen_xmin+2*4)
else
open(unit=31,file=trim(prname)//'absorb_xmin.bin',status='unknown',form='unformatted',access='direct',&
recl=reclen_xmin+2*4)
endif
endif
if (nspec2D_xmax > 0 .and. (SIMULATION_TYPE == 3 .or. (SIMULATION_TYPE == 1 .and. SAVE_FORWARD))) then
allocate(absorb_xmax(NDIM,NGLLY,NGLLZ,nspec2D_xmax))
reclen_xmax = CUSTOM_REAL * (NDIM * NGLLY * NGLLZ * nspec2D_xmax)
if (SIMULATION_TYPE == 3) then
open(unit=32,file=trim(prname)//'absorb_xmax.bin',status='old',action='read',form='unformatted',access='direct', &
recl=reclen_xmax+2*4)
else
open(unit=32,file=trim(prname)//'absorb_xmax.bin',status='unknown',form='unformatted',access='direct', &
recl=reclen_xmax+2*4)
endif
endif
if (nspec2D_ymin > 0 .and. (SIMULATION_TYPE == 3 .or. (SIMULATION_TYPE == 1 .and. SAVE_FORWARD))) then
allocate(absorb_ymin(NDIM,NGLLX,NGLLZ,nspec2D_ymin))
reclen_ymin = CUSTOM_REAL * (NDIM * NGLLX * NGLLZ * nspec2D_ymin)
if (SIMULATION_TYPE == 3) then
open(unit=33,file=trim(prname)//'absorb_ymin.bin',status='old',action='read',form='unformatted',access='direct',&
recl=reclen_ymin+2*4)
else
open(unit=33,file=trim(prname)//'absorb_ymin.bin',status='unknown',form='unformatted',access='direct',&
recl=reclen_ymin+2*4)
endif
endif
if (nspec2D_ymax > 0 .and. (SIMULATION_TYPE == 3 .or. (SIMULATION_TYPE == 1 .and. SAVE_FORWARD))) then
allocate(absorb_ymax(NDIM,NGLLX,NGLLZ,nspec2D_ymax))
reclen_ymax = CUSTOM_REAL * (NDIM * NGLLX * NGLLZ * nspec2D_ymax)
if (SIMULATION_TYPE == 3) then
open(unit=34,file=trim(prname)//'absorb_ymax.bin',status='old',action='read',form='unformatted',access='direct',&
recl=reclen_ymax+2*4)
else
open(unit=34,file=trim(prname)//'absorb_ymax.bin',status='unknown',form='unformatted',access='direct',&
recl=reclen_ymax+2*4)
endif
endif
if (NSPEC2D_BOTTOM > 0 .and. (SIMULATION_TYPE == 3 .or. (SIMULATION_TYPE == 1 .and. SAVE_FORWARD))) then
allocate(absorb_zmin(NDIM,NGLLX,NGLLY,NSPEC2D_BOTTOM))
reclen_zmin = CUSTOM_REAL * (NDIM * NGLLX * NGLLY * NSPEC2D_BOTTOM)
if (SIMULATION_TYPE == 3) then
open(unit=35,file=trim(prname)//'absorb_zmin.bin',status='old',action='read',form='unformatted',access='direct',&
recl=reclen_zmin+2*4)
else
open(unit=35,file=trim(prname)//'absorb_zmin.bin',status='unknown',form='unformatted',access='direct',&
recl=reclen_zmin+2*4)
endif
endif
endif
! $$$$$$$$$$$$$$$$$$$$$$$$ SOURCES $$$$$$$$$$$$$$$$$
! read topography and bathymetry file
if(TOPOGRAPHY .or. OCEANS) then
NX_TOPO = NX_TOPO_SOCAL
NY_TOPO = NY_TOPO_SOCAL
ORIG_LAT_TOPO = ORIG_LAT_TOPO_SOCAL
ORIG_LONG_TOPO = ORIG_LONG_TOPO_SOCAL
DEGREES_PER_CELL_TOPO = DEGREES_PER_CELL_TOPO_SOCAL
topo_file = TOPO_FILE_SOCAL
allocate(itopo_bathy(NX_TOPO,NY_TOPO))
call read_topo_bathy_file(itopo_bathy,NX_TOPO,NY_TOPO,topo_file)
if(myrank == 0) then
write(IMAIN,*)
write(IMAIN,*) 'regional topography file read ranges in m from ', &
minval(itopo_bathy),' to ',maxval(itopo_bathy)
write(IMAIN,*)
endif
endif
! write source and receiver VTK files for Paraview
if (myrank == 0) then
open(IOVTK,file=trim(OUTPUT_FILES)//'/sr.vtk',status='unknown')
write(IOVTK,'(a)') '# vtk DataFile Version 2.0'
write(IOVTK,'(a)') 'Source and Receiver VTK file'
write(IOVTK,'(a)') 'ASCII'
write(IOVTK,'(a)') 'DATASET POLYDATA'
! LQY -- cannot figure out NSOURCES+nrec at this point
write(IOVTK, '(a,i6,a)') 'POINTS ', 2, ' float'
endif
! allocate arrays for source
allocate(islice_selected_source(NSOURCES))
allocate(ispec_selected_source(NSOURCES))
allocate(Mxx(NSOURCES))
allocate(Myy(NSOURCES))
allocate(Mzz(NSOURCES))
allocate(Mxy(NSOURCES))
allocate(Mxz(NSOURCES))
allocate(Myz(NSOURCES))
allocate(xi_source(NSOURCES))
allocate(eta_source(NSOURCES))
allocate(gamma_source(NSOURCES))
allocate(t_cmt(NSOURCES))
allocate(hdur(NSOURCES))
allocate(hdur_gaussian(NSOURCES))
allocate(utm_x_source(NSOURCES))
allocate(utm_y_source(NSOURCES))
! locate sources in the mesh
call locate_source(ibool,NSOURCES,myrank,NSPEC_AB,NGLOB_AB, &
xstore,ystore,zstore,xigll,yigll,zigll,NPROC, &
sec,t_cmt,yr,jda,ho,mi,utm_x_source,utm_y_source, &
NSTEP,DT,hdur,Mxx,Myy,Mzz,Mxy,Mxz,Myz, &
islice_selected_source,ispec_selected_source, &
xi_source,eta_source,gamma_source, &
LATITUDE_MIN,LATITUDE_MAX,LONGITUDE_MIN,LONGITUDE_MAX,Z_DEPTH_BLOCK, &
TOPOGRAPHY,itopo_bathy,UTM_PROJECTION_ZONE, &
PRINT_SOURCE_TIME_FUNCTION,SUPPRESS_UTM_PROJECTION, &
NX_TOPO,NY_TOPO,ORIG_LAT_TOPO,ORIG_LONG_TOPO,DEGREES_PER_CELL_TOPO)
if(minval(t_cmt) /= 0.) call exit_MPI(myrank,'one t_cmt must be zero, others must be positive')
! filter source time function by Gaussian with hdur = HDUR_MOVIE when outputing movies or shakemaps
if (MOVIE_SURFACE .or. MOVIE_VOLUME .or. CREATE_SHAKEMAP) then
hdur = sqrt(hdur**2 + HDUR_MOVIE**2)
if(myrank == 0) then
write(IMAIN,*)
write(IMAIN,*) 'Each source is being convolved with HDUR_MOVIE = ',HDUR_MOVIE
write(IMAIN,*)
endif
endif
! convert the half duration for triangle STF to the one for gaussian STF
hdur_gaussian = hdur/SOURCE_DECAY_MIMIC_TRIANGLE
! define t0 as the earliest start time
t0 = - 1.5d0 * minval(t_cmt-hdur)
!$$$$$$$$$$$$$$$$$$ RECEIVERS $$$$$$$$$$$$$$$$$$$$$
if (SIMULATION_TYPE == 1) then
call get_value_string(filtered_rec_filename, 'solver.STATIONS_FILTERED', 'DATA/STATIONS_FILTERED')
! get total number of stations
open(unit=IIN,file=filtered_rec_filename,iostat=ios,status='old',action='read')
nrec = 0
do while(ios == 0)
read(IIN,"(a)",iostat=ios) dummystring
if(ios == 0) nrec = nrec + 1
enddo
close(IIN)
if(nrec < 1) call exit_MPI(myrank,'need at least one receiver')
else
call get_value_string(rec_filename, 'solver.STATIONS', 'DATA/STATIONS_ADJOINT')
call get_value_string(filtered_rec_filename, 'solver.STATIONS_FILTERED', 'DATA/STATIONS_ADJOINT_FILTERED')
call station_filter(myrank,rec_filename,filtered_rec_filename,nrec, &
LATITUDE_MIN, LATITUDE_MAX, LONGITUDE_MIN, LONGITUDE_MAX)
if (nrec < 1) call exit_MPI(myrank, 'adjoint simulation needs at least one source')
call sync_all()
endif
if(myrank == 0) then
write(IMAIN,*)
if (SIMULATION_TYPE == 1 .or. SIMULATION_TYPE == 3) then
write(IMAIN,*) 'Total number of receivers = ', nrec
else
write(IMAIN,*) 'Total number of adjoint sources = ', nrec
endif
write(IMAIN,*)
endif
if(nrec < 1) call exit_MPI(myrank,'need at least one receiver')
! allocate memory for receiver arrays
allocate(islice_selected_rec(nrec))
allocate(ispec_selected_rec(nrec))
allocate(xi_receiver(nrec))
allocate(eta_receiver(nrec))
allocate(gamma_receiver(nrec))
allocate(station_name(nrec))
allocate(network_name(nrec))
allocate(nu(NDIM,NDIM,nrec))
! locate receivers in the mesh
call locate_receivers(ibool,myrank,NSPEC_AB,NGLOB_AB, &
xstore,ystore,zstore,xigll,yigll,zigll,filtered_rec_filename, &
nrec,islice_selected_rec,ispec_selected_rec, &
xi_receiver,eta_receiver,gamma_receiver,station_name,network_name,nu, &
NPROC,utm_x_source(1),utm_y_source(1), &
TOPOGRAPHY,itopo_bathy,UTM_PROJECTION_ZONE,SUPPRESS_UTM_PROJECTION, &
NX_TOPO,NY_TOPO,ORIG_LAT_TOPO,ORIG_LONG_TOPO,DEGREES_PER_CELL_TOPO)
!###################### SOURCE ARRAYS ################
if (SIMULATION_TYPE == 1 .or. SIMULATION_TYPE == 3) then
allocate(sourcearray(NDIM,NGLLX,NGLLY,NGLLZ))
allocate(sourcearrays(NSOURCES,NDIM,NGLLX,NGLLY,NGLLZ))
! compute source arrays
do isource = 1,NSOURCES
! check that the source slice number is okay
if(islice_selected_source(isource) < 0 .or. islice_selected_source(isource) > NPROC-1) &
call exit_MPI(myrank,'something is wrong with the source slice number')
! compute source arrays in source slice
if(myrank == islice_selected_source(isource)) then
call compute_arrays_source(ispec_selected_source(isource), &
xi_source(isource),eta_source(isource),gamma_source(isource),sourcearray, &
Mxx(isource),Myy(isource),Mzz(isource),Mxy(isource),Mxz(isource),Myz(isource), &
xix,xiy,xiz,etax,etay,etaz,gammax,gammay,gammaz, &
xigll,yigll,zigll,NSPEC_AB)
sourcearrays(isource,:,:,:,:) = sourcearray(:,:,:,:)
endif
enddo
endif
if (SIMULATION_TYPE == 2 .or. SIMULATION_TYPE == 3) then
nadj_rec_local = 0
do irec = 1,nrec
if(myrank == islice_selected_rec(irec))then
! check that the source slice number is okay
if(islice_selected_rec(irec) < 0 .or. islice_selected_rec(irec) > NPROC-1) &
call exit_MPI(myrank,'something is wrong with the source slice number in adjoint simulation')
nadj_rec_local = nadj_rec_local + 1
endif
enddo
allocate(adj_sourcearray(NSTEP,NDIM,NGLLX,NGLLY,NGLLZ))
if (nadj_rec_local > 0) allocate(adj_sourcearrays(nadj_rec_local,NSTEP,NDIM,NGLLX,NGLLY,NGLLZ))
irec_local = 0
do irec = 1, nrec
! compute only adjoint source arrays in the local slice
if(myrank == islice_selected_rec(irec)) then
irec_local = irec_local + 1
adj_source_file = trim(station_name(irec))//'.'//trim(network_name(irec))
call compute_arrays_adjoint_source(myrank, adj_source_file, &
xi_receiver(irec), eta_receiver(irec), gamma_receiver(irec), &
adj_sourcearray, xigll,yigll,zigll,NSTEP)
adj_sourcearrays(irec_local,:,:,:,:,:) = adj_sourcearray(:,:,:,:,:)
endif
enddo
endif
!--- select local receivers
! count number of receivers located in this slice
nrec_local = 0
if (SIMULATION_TYPE == 1 .or. SIMULATION_TYPE == 3) then
nrec_simulation = nrec
do irec = 1,nrec
if(myrank == islice_selected_rec(irec)) nrec_local = nrec_local + 1
enddo
else
nrec_simulation = NSOURCES
do isource = 1, NSOURCES
if(myrank == islice_selected_source(isource)) nrec_local = nrec_local + 1
enddo
endif
if (nrec_local > 0) then
! allocate Lagrange interpolators for receivers
allocate(hxir_store(nrec_local,NGLLX))
allocate(hetar_store(nrec_local,NGLLY))
allocate(hgammar_store(nrec_local,NGLLZ))
! define local to global receiver numbering mapping
allocate(number_receiver_global(nrec_local))
irec_local = 0
if (SIMULATION_TYPE == 1 .or. SIMULATION_TYPE == 3) then
do irec = 1,nrec
if(myrank == islice_selected_rec(irec)) then
irec_local = irec_local + 1
number_receiver_global(irec_local) = irec
endif
enddo
else