-
Notifications
You must be signed in to change notification settings - Fork 0
/
meshfem3D_models.F90
executable file
·2163 lines (1783 loc) · 70.2 KB
/
meshfem3D_models.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 G l o b e
! ----------------------------
!
! Main historical authors: Dimitri Komatitsch and Jeroen Tromp
! Princeton University, USA
! and CNRS / University of Marseille, France
! (there are currently many more authors!)
! (c) Princeton University and CNRS / University of Marseille, April 2014
!
! 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 3 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.
!
!=====================================================================
subroutine meshfem3D_models_broadcast()
! preparing model parameter coefficients on all processes
use shared_parameters, only: LOCAL_PATH,SAVE_MESH_FILES,ADD_SCATTERING_PERTURBATIONS
use meshfem_models_par
implicit none
! local parameters
integer :: ier
! sets up spline coefficients for ellipticity
if (ELLIPTICITY) call make_ellipticity(nspl_ellip,rspl_ellip,ellipicity_spline,ellipicity_spline2)
! read topography and bathymetry file
if (TOPOGRAPHY) then
! arrays for elevations
allocate(ibathy_topo(NX_BATHY,NY_BATHY),stat=ier)
if (ier /= 0) stop 'Error allocating ibathy_topo array'
! initializes
ibathy_topo(:,:) = 0
! sets up topo/bathy
call model_topo_bathy_broadcast(ibathy_topo,LOCAL_PATH)
! saves VTK output
if (SAVE_MESH_FILES) call meshfem3D_plot_VTK_topo_bathy()
endif
!---
!
! ADD YOUR MODEL HERE
!
!---
! reads 1D reference models
call meshfem3D_reference_model_broadcast()
! reads in 3D mantle models
call meshfem3D_mantle_broadcast()
! reads in crustal model
call meshfem3D_crust_broadcast()
! attenuation
if (ATTENUATION) then
call model_attenuation_broadcast()
! 3D attenuation
if (ATTENUATION_GLL) then
call model_attenuation_gll_broadcast()
else if (ATTENUATION_3D) then
! Colleen's model defined originally between 24.4km and 650km
call model_atten3D_QRFSI12_broadcast()
else
! sets up attenuation coefficients according to the chosen, "pure" 1D model
! (including their 1D-crustal profiles)
call model_attenuation_setup(REFERENCE_1D_MODEL,CRUSTAL)
endif
endif
! sets up scattering perturbations
if (ADD_SCATTERING_PERTURBATIONS) call model_scattering_broadcast()
end subroutine meshfem3D_models_broadcast
!
!-------------------------------------------------------------------------------------------------
!
subroutine meshfem3D_reference_model_broadcast()
! preparing model parameter coefficients on all processes for reference models
use meshfem_models_par
implicit none
!---
!
! ADD YOUR MODEL HERE
!
!---
! re-defines/initializes models 1066a and ak135 and ref
! ( with possible external crustal model: if CRUSTAL is set to true
! it strips the 1-D crustal profile and replaces it with mantle properties)
select case (REFERENCE_1D_MODEL)
case (REFERENCE_MODEL_1066A)
call model_1066a_broadcast(CRUSTAL)
case (REFERENCE_MODEL_AK135F_NO_MUD)
call model_ak135_broadcast(CRUSTAL)
case (REFERENCE_MODEL_1DREF)
call model_1dref_broadcast(CRUSTAL)
case (REFERENCE_MODEL_SEA1D)
call model_sea1d_broadcast(CRUSTAL)
case (REFERENCE_MODEL_CASE65TAY)
call model_case65TAY_broadcast(CRUSTAL)
case (REFERENCE_MODEL_MARS_1D)
call model_mars_1D_broadcast(CRUSTAL)
case (REFERENCE_MODEL_CCREM)
call model_ccrem_broadcast(CRUSTAL)
case(REFERENCE_MODEL_SEMUCB) ! Berkeley model always has external crust
call model_1dberkeley_broadcast()
case(REFERENCE_MODEL_PREM_A3D) ! Berkeley model always has external crust
call model_1dberkeley_broadcast()
case (REFERENCE_MODEL_VPREMOON)
call model_vpremoon_broadcast(CRUSTAL)
end select
end subroutine meshfem3D_reference_model_broadcast
!
!-------------------------------------------------------------------------------------------------
!
subroutine meshfem3D_mantle_broadcast()
! preparing model parameter coefficients on all processes for mantle models
use constants, only: myrank,IMAIN
use meshfem_models_par
implicit none
!---
!
! ADD YOUR MODEL HERE
!
!---
! GLL model uses S29EA as reference 3D model
if (THREE_D_MODEL == THREE_D_MODEL_GLL) then
! sets to initial reference model from which iterations started
THREE_D_MODEL = GLL_REFERENCE_MODEL
endif
! 3D mantle models
if (MODEL_3D_MANTLE_PERTUBATIONS) then
select case (THREE_D_MODEL)
case (THREE_D_MODEL_S20RTS)
call model_s20rts_broadcast()
case (THREE_D_MODEL_S40RTS)
call model_s40rts_broadcast()
case(THREE_D_MODEL_MANTLE_SH)
call model_mantle_sh_broadcast()
case (THREE_D_MODEL_SEA99_JP3D)
! the variables read are declared and stored in structure model_sea99_s_par and model_jp3d_par
call model_sea99_s_broadcast()
call model_jp3d_broadcast()
case (THREE_D_MODEL_SEA99)
! the variables read are declared and stored in structure model_sea99_s_par
call model_sea99_s_broadcast()
case (THREE_D_MODEL_JP3D)
! the variables read are declared and stored in structure model_jp3d_par
call model_jp3d_broadcast()
case (THREE_D_MODEL_S362ANI,THREE_D_MODEL_S362WMANI, &
THREE_D_MODEL_S362ANI_PREM,THREE_D_MODEL_S29EA)
! the variables read are declared and stored in structure model_s362ani_par
call model_s362ani_broadcast(THREE_D_MODEL)
case (THREE_D_MODEL_PPM)
! Point Profile Models
! the variables read are declared and stored in structure model_ppm_par
call model_ppm_broadcast()
case (THREE_D_MODEL_GAPP2)
! GAP model
call model_gapp2_broadcast()
case (THREE_D_MODEL_SGLOBE,THREE_D_MODEL_SGLOBE_ISO)
! SGLOBE-rani model
call model_sglobe_broadcast()
case (THREE_D_MODEL_ANISO_MANTLE)
! Montagner anisotropic model
call model_aniso_mantle_broadcast()
case (THREE_D_MODEL_BKMNS_GLAD)
! GLAD model expansion on block-mantle-spherical-harmonics
if (myrank == 0) write(IMAIN,*) 'setting up both s362ani and bkmns models...'
! needs s362ani topo for 410/660
call model_s362ani_broadcast(THREE_D_MODEL_S362ANI)
! block-mantle-spherical-harmonics model
call model_bkmns_mantle_broadcast()
case (THREE_D_MODEL_SPIRAL)
! SPiRal model
call model_mantle_spiral_broadcast()
case(THREE_D_MODEL_BERKELEY)
! Berkeley SEMUCB model (should be used along with the berkeley crust)
call model_berkeley_broadcast()
case (THREE_D_MODEL_HETEROGEN_PREM)
!chris modif checker 02/20/21
call model_heterogen_mntl_broadcast()
case (THREE_D_MODEL_SH_MARS)
! Mars spherical harmonics model
call model_SH_mars_broadcast()
case default
call exit_MPI(myrank,'3D model not defined')
end select
endif
! adds additional perturbations on top of a reference 3D model
! using "arbitrary" mantle models (for example density models derived from geodynamics modeling)
if (HETEROGEN_3D_MANTLE) &
call model_heterogen_mntl_broadcast()
! Enclose this in an ifdef so we don't link to netcdf
! if we don't need it.
! Collaborative Earth Model (CEM)
#ifdef USE_CEM
if (CEM_REQUEST .or. CEM_ACCEPT) &
call model_cem_broadcast()
#endif
! IRIS EMC models
#ifdef USE_EMC
if (EMC_MODEL) &
call model_EMC_broadcast()
#endif
! GLL model
if (MODEL_GLL) &
call model_gll_broadcast()
end subroutine meshfem3D_mantle_broadcast
!
!-------------------------------------------------------------------------------------------------
!
subroutine meshfem3D_crust_broadcast()
! preparing model parameter coefficients on all processes for crustal models
use constants, only: myrank,IMAIN
use shared_parameters, only: SAVE_MESH_FILES
use meshfem_models_par
implicit none
! checks if anything to do
if (.not. CRUSTAL) return
!---
!
! ADD YOUR MODEL HERE
!
!---
select case (REFERENCE_CRUSTAL_MODEL)
case (ICRUST_CRUST1)
! crust 1.0
call model_crust_1_0_broadcast()
case (ICRUST_CRUST2)
! default
! crust 2.0
call model_crust_2_0_broadcast()
case (ICRUST_CRUSTMAPS)
! general crustmaps
call model_crustmaps_broadcast()
case (ICRUST_EPCRUST)
! EPcrust (regional crustal model for Europe)
call model_epcrust_broadcast()
! by default crust 1.0 (global coverage)
call model_crust_1_0_broadcast()
case (ICRUST_CRUST_SH)
! SH crustmaps
call model_crust_sh_broadcast()
case (ICRUST_EUCRUST)
! EUcrust07 Vp crustal structure (regional crustal model)
call model_eucrust_broadcast()
! by default (vs,rho,eta,moho) from crust 1.0 (global coverage)
call model_crust_1_0_broadcast()
case (ICRUST_SGLOBECRUST)
! modified crust2.0 for SGLOBE-rani
call model_sglobecrust_broadcast()
case (ICRUST_BKMNS_GLAD)
! GLAD model expansion on block-mantle-spherical-harmonics
if (myrank == 0) write(IMAIN,*) 'setting up both crust 2.0 and bkmns crust models ...'
! needs crust 2.0 for moho depths
call model_crust_2_0_broadcast()
! crustal structure in blocks between topography and 80km depth
call model_bkmns_crust_broadcast()
case (ICRUST_SPIRAL)
! anisotropic crust from SPiRaL
call model_crust_spiral_broadcast()
case (ICRUST_SH_MARS)
! Mars SH model (defines both crust & mantle)
call model_SH_mars_broadcast()
case (ICRUST_BERKELEY)
! EPcrust
call model_berkeley_crust_broadcast(myrank)
case default
stop 'crustal model type not defined'
end select
! plot moho as VTK file output
if (SAVE_MESH_FILES) call meshfem3D_plot_VTK_crust_moho()
end subroutine meshfem3D_crust_broadcast
!
!-------------------------------------------------------------------------------------------------
!
subroutine meshfem3D_models_get1D_val(iregion_code,idoubling, &
r_prem,rho,vpv,vph,vsv,vsh,eta_aniso, &
Qkappa,Qmu,RICB,RCMB, &
RTOPDDOUBLEPRIME,R80,R120,R220,R400,R670,R771, &
RMOHO,RMIDDLE_CRUST)
! reference model values
!
! for a given location radius (r_prem, which is the point's radius with tolerance factor),
! this calculates density and velocities
!
! note: if CRUSTAL is set, it strips the 1-D crustal profile and mantle gets expanded
! up to the surface.
! only exception is JP1D...
!
! routine returns: rho,vpv,vph,vsv,vsh,eta_aniso,Qkappa,Qmu
use shared_parameters, only: REGIONAL_MESH_CUTOFF
use meshfem_models_par
implicit none
integer,intent(in) :: iregion_code,idoubling
double precision,intent(in) :: r_prem
double precision,intent(inout) :: vpv,vph,vsv,vsh,eta_aniso,rho
double precision,intent(inout) :: Qkappa,Qmu
double precision,intent(in) :: RICB,RCMB,RTOPDDOUBLEPRIME,R80,R120,R220,R400, &
R670,R771,RMOHO,RMIDDLE_CRUST
! local parameters
double precision :: drhodr,vp,vs
logical :: check_doubling_flag
! initializes
if (REGIONAL_MESH_CUTOFF) then
! no need to check doubling flags
check_doubling_flag = .false.
else
! check if mesh doubling flags and layering are consistent
check_doubling_flag = .true.
endif
vpv = ZERO
vph = ZERO
vsv = ZERO
vsh = ZERO
eta_aniso = ONE
rho = ZERO
Qmu = ZERO
Qkappa = ZERO
!---
!
! ADD YOUR MODEL HERE
!
!---
! gets 1-D reference model parameters
select case (REFERENCE_1D_MODEL)
case (REFERENCE_MODEL_PREM,REFERENCE_MODEL_PREM2)
! PREM (by Dziewonski & Anderson) - used also as background for 3D models
if (TRANSVERSE_ISOTROPY) then
! default PREM:
! gets anisotropic PREM parameters, with radial anisotropic extension (from moho to surface for crustal model)
call model_prem_aniso(r_prem,rho,vpv,vph,vsv,vsh,eta_aniso,Qkappa,Qmu,idoubling,CRUSTAL,check_doubling_flag)
! calculates isotropic values
vp = sqrt(((8.d0+4.d0*eta_aniso)*vph*vph + 3.d0*vpv*vpv &
+ (8.d0 - 8.d0*eta_aniso)*vsv*vsv)/15.d0)
if (iregion_code == IREGION_OUTER_CORE) then
! fluid with zero shear speed
vs = 0.d0
else
vs = sqrt(((1.d0-2.d0*eta_aniso)*vph*vph + vpv*vpv &
+ 5.d0*vsh*vsh + (6.d0+4.d0*eta_aniso)*vsv*vsv)/15.d0)
endif
!daniel todo:
! specific 3D models with PREM references which would become too fast at shorter periods ( < 40s Love waves)
!
!select case (THREE_D_MODEL)
!
! eventually sgloberani, check...
!case (THREE_D_MODEL_SGLOBE,THREE_D_MODEL_SGLOBE_ISO)
! ! gets anisotropic PREM parameters, with isotropic extension (from moho to surface for crustal model)
! call model_prem_aniso_extended_isotropic(r_prem,rho,vpv,vph,vsv,vsh,eta_aniso,Qkappa,Qmu, &
! idoubling,CRUSTAL,check_doubling_flag)
!
! eventually also Ritsema models, check...
!case (THREE_D_MODEL_S20RTS,THREE_D_MODEL_S40RTS)
! ! gets anisotropic PREM parameters, with isotropic extension (from moho to surface for crustal model)
! call model_prem_aniso_extended_isotropic(r_prem,rho,vpv,vph,vsv,vsh,eta_aniso,Qkappa,Qmu, &
! idoubling,CRUSTAL,check_doubling_flag)
!
!case default
! continue
!end select
else
! isotropic PREM model
call model_prem_iso(r_prem,rho,drhodr,vp,vs,Qkappa,Qmu,idoubling,CRUSTAL,check_doubling_flag)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
endif
case (REFERENCE_MODEL_1DREF)
! 1D-REF also known as STW105 (by Kustowski et al.) - used also as background for 3D models
call model_1dref(r_prem,rho,vpv,vph,vsv,vsh,eta_aniso,Qkappa,Qmu,iregion_code,CRUSTAL)
! calculates isotropic values
vp = sqrt(((8.d0+4.d0*eta_aniso)*vph*vph + 3.d0*vpv*vpv &
+ (8.d0 - 8.d0*eta_aniso)*vsv*vsv)/15.d0)
if (iregion_code == IREGION_OUTER_CORE) then
! fluid with zero shear speed
vs = 0.d0
else
vs = sqrt(((1.d0-2.d0*eta_aniso)*vph*vph + vpv*vpv &
+ 5.d0*vsh*vsh + (6.d0+4.d0*eta_aniso)*vsv*vsv)/15.d0)
endif
if (.not. TRANSVERSE_ISOTROPY) then
if (.not. (MODEL_3D_MANTLE_PERTUBATIONS .and. iregion_code == IREGION_CRUST_MANTLE)) then
! this case here is only executed for 1D_ref_iso
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
endif
endif
case (REFERENCE_MODEL_1066A)
! 1066A (by Gilbert & Dziewonski) - pure isotropic model, used in 1D model mode only
call model_1066a(r_prem,rho,vp,vs,Qkappa,Qmu,iregion_code)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case (REFERENCE_MODEL_AK135F_NO_MUD)
! AK135 (by Kennett et al. ) - pure isotropic model, used in 1D model mode only
call model_ak135(r_prem,rho,vp,vs,Qkappa,Qmu,iregion_code)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case (REFERENCE_MODEL_IASP91)
! IASP91 (by Kennett & Engdahl) - pure isotropic model, used in 1D model mode only
call model_iasp91(r_prem,rho,vp,vs,Qkappa,Qmu,idoubling, &
ONE_CRUST,check_doubling_flag,RICB,RCMB,RTOPDDOUBLEPRIME, &
R771,R670,R400,R220,R120,RMOHO,RMIDDLE_CRUST)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case(REFERENCE_MODEL_SEMUCB)
! BERKELEY MODEL SEMUCB
call model_1dberkeley(r_prem,rho,vpv,vph,vsv,vsh,eta_aniso,Qkappa,Qmu,iregion_code,CRUSTAL)
case(REFERENCE_MODEL_PREM_A3D)
! BERKELEY MODEL SEMUCB
call model_1dberkeley(r_prem,rho,vpv,vph,vsv,vsh,eta_aniso,Qkappa,Qmu,iregion_code,CRUSTAL)
case (REFERENCE_MODEL_JP1D)
!JP1D (by Zhao et al.) - pure isotropic model, used also as background for 3D models
call model_jp1d(r_prem,rho,vp,vs,Qkappa,Qmu,idoubling, &
check_doubling_flag,RICB,RCMB,RTOPDDOUBLEPRIME, &
R670,R220,R771,R400,R80,RMOHO,RMIDDLE_CRUST)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case (REFERENCE_MODEL_SEA1D)
! SEA1D (by Lebedev & Nolet) - pure isotropic model, used also as background for 3D models
call model_sea1d(r_prem,rho,vp,vs,Qkappa,Qmu,iregion_code)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case (REFERENCE_MODEL_CCREM)
! CCREM (by Ma & Tkalcic) - pure isotropic model
call model_ccrem(r_prem,rho,vp,vs,Qkappa,Qmu,idoubling,iregion_code)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
! Mars 1D models
case (REFERENCE_MODEL_SOHL)
! Mars
call model_Sohl(r_prem,rho,drhodr,vp,vs,Qkappa,Qmu,idoubling,CRUSTAL,check_doubling_flag)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case (REFERENCE_MODEL_CASE65TAY)
! Mars
call model_case65TAY(r_prem,rho,vp,vs,Qkappa,Qmu,iregion_code)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case (REFERENCE_MODEL_MARS_1D)
! Mars
call model_mars_1D(r_prem,rho,drhodr,vp,vs,Qkappa,Qmu,idoubling,CRUSTAL,.true.,iregion_code)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
! Moon 1D models
case (REFERENCE_MODEL_VPREMOON)
! Moon
call model_vpremoon(r_prem,rho,drhodr,vp,vs,Qkappa,Qmu,idoubling,CRUSTAL,.true.,iregion_code)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
case default
stop 'unknown 1D reference Earth model in meshfem3D_models_get1D_val()'
end select
! needs to set vpv,vph,vsv,vsh and eta_aniso for isotropic models
if (.not. TRANSVERSE_ISOTROPY) then
! in the case of s362iso we want to save the anisotropic constants for the Voigt average
if (.not. (REFERENCE_1D_MODEL == REFERENCE_MODEL_1DREF &
.and. MODEL_3D_MANTLE_PERTUBATIONS &
.and. iregion_code == IREGION_CRUST_MANTLE)) then
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.d0
endif
endif ! TRANSVERSE_ISOTROPY
end subroutine meshfem3D_models_get1D_val
!
!-------------------------------------------------------------------------------------------------
!
subroutine meshfem3D_models_get3Dmntl_val(iregion_code,r_prem,rho, &
vpv,vph,vsv,vsh,eta_aniso, &
RCMB,RMOHO, &
r,theta,phi, &
c11,c12,c13,c14,c15,c16,c22,c23,c24,c25,c26, &
c33,c34,c35,c36,c44,c45,c46,c55,c56,c66, &
ispec,i,j,k)
use meshfem_models_par
use constants,only: EARTH_R_KM
implicit none
integer, intent(in) :: iregion_code
double precision, intent(in) :: r_prem
double precision, intent(inout) :: rho
double precision, intent(inout) :: vpv,vph,vsv,vsh,eta_aniso
double precision,intent(in) :: RCMB,RMOHO
double precision,intent(in) :: r,theta,phi
! the 21 coefficients for an anisotropic medium in reduced notation
double precision,intent(inout) :: c11,c12,c13,c14,c15,c16,c22,c23,c24,c25,c26, &
c33,c34,c35,c36,c44,c45,c46,c55,c56,c66
! heterogen model and CEM needs these (CEM to determine iglob)
integer, intent(in) :: ispec, i, j, k
! local parameters
double precision :: r_used
double precision :: dvp,dvs,drho,vp,vs,moho,sediment
double precision :: dvpv,dvph,dvsv,dvsh,deta
double precision :: lat,lon
double precision :: A,C,L,N,F
real(kind=4) :: xcolat,xlon,xrad
real(kind=4) :: xdvpv,xdvph,xdvsv,xdvsh
logical :: found_crust,suppress_mantle_extension,is_inside_region,convert_tiso_to_cij
! initializes perturbation values
dvs = ZERO
dvp = ZERO
drho = ZERO
dvpv = ZERO
dvph = ZERO
dvsv = ZERO
dvsh = ZERO
deta = ZERO
xdvpv = 0.d0
xdvph = 0.d0
xdvsv = 0.d0
xdvsh = 0.d0
r_used = ZERO
suppress_mantle_extension = .false.
! note: global mantle models are in general defined with respect to a spherical Earth.
! furthermore, geocentric and geographic colatitude (theta) is the same for a spherical Earth.
!
! the only model using point locations at "true" Earth positions is the EMC model, where the mesh
! can be either elliptical or not, depending on the Par_file 'ELLIPTICITY' flag.
!
! TODO: in case, we will have more 3D (tomographic) models in future that need "true" positions, we can re-visit this.
!---
!
! ADD YOUR MODEL HERE
!
!---
if (iregion_code == IREGION_CRUST_MANTLE) then
! crust/mantle
! sets flag when mantle should not be extended to surface
if (r_prem >= RMOHO/R_PLANET .and. .not. CRUSTAL) then
suppress_mantle_extension = .true.
endif
! gets parameters for isotropic 3D mantle model
!
! note: there can be transverse isotropy in the mantle, but only Lame parameters
! like kappav,kappah,muv,muh and eta_aniso are used for these simulations
!
! note: in general, models here make use of perturbation values with respect to their
! corresponding 1-D reference models
if (MODEL_3D_MANTLE_PERTUBATIONS .and. r_prem > RCMB/R_PLANET .and. .not. suppress_mantle_extension) then
! extend 3-D mantle model above the Moho to the surface before adding the crust
if (r_prem > RCMB/R_PLANET .and. r_prem < RMOHO/R_PLANET) then
! GLL point is in mantle region, takes exact location
r_used = r
else ! else if (r_prem >= RMOHO/R_PLANET) then
if (CRUSTAL) then
! GLL point is above moho
! takes radius slightly below moho radius, this will then "extend the mantle up to the surface";
! crustal values will be superimposed later on
!
! note: this assumes that all the following mantle models are defined below RMOHO.
! this is in general true, almost all mantle models fix the moho depth at 24.4 km (set as RMOHO reference)
! and invert their velocity models for structures below.
!
! however, in case a mantle model extends to shallower depths, those velocity regions will be ignored.
!
! to do in future: we might want to let the mantle routines decide where to this upper cut-off radius value
r_used = 0.999999d0*RMOHO/R_PLANET
endif
endif
! gets model parameters
select case (THREE_D_MODEL)
case (THREE_D_MODEL_S20RTS)
! s20rts
call mantle_s20rts(r_used,theta,phi,dvs,dvp,drho)
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
rho = rho*(1.0d0+drho)
case (THREE_D_MODEL_S40RTS)
! s40rts
call mantle_s40rts(r_used,theta,phi,dvs,dvp,drho)
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
rho = rho*(1.0d0+drho)
case(THREE_D_MODEL_MANTLE_SH)
! full_sh model
lat = (PI/2.0d0-theta)*180.0d0/PI
lon = phi*180.0d0/PI
if (lon > 180.0d0) lon = lon - 360.0d0
call mantle_sh(lat,lon,r_used,dvpv,dvph,dvsv,dvsh,deta,drho)
vpv = vpv*(1.0d0+dvpv)
vph = vph*(1.0d0+dvph)
vsv = vsv*(1.0d0+dvsv)
vsh = vsh*(1.0d0+dvsh)
eta_aniso = eta_aniso*(1.0d0+deta)
rho = rho*(1.0d0+drho)
case (THREE_D_MODEL_SEA99_JP3D)
! sea99 + jp3d1994
call model_sea99_s(r_used,theta,phi,dvs)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
! use Lebedev model sea99 as background and add vp & vs perturbation from Zhao 1994 model jp3d
call model_jp3d_iso_zhao(r_used,theta,phi,vp,vs,dvp,dvs,rho,moho,sediment,found_crust,is_inside_region)
if (is_inside_region) then
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
endif
case (THREE_D_MODEL_SEA99)
! sea99 Vs-only
call model_sea99_s(r_used,theta,phi,dvs)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
case (THREE_D_MODEL_JP3D)
! jp3d1994
call model_jp3d_iso_zhao(r_used,theta,phi,vp,vs,dvp,dvs,rho,moho,sediment,found_crust,is_inside_region)
if (is_inside_region) then
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
endif
case (THREE_D_MODEL_S362ANI,THREE_D_MODEL_S362WMANI, &
THREE_D_MODEL_S362ANI_PREM,THREE_D_MODEL_S29EA)
! 3D Harvard models s362ani, s362wmani, s362ani_prem and s2.9ea
xcolat = sngl(theta*180.0d0/PI)
xlon = sngl(phi*180.0d0/PI)
xrad = sngl(r_used*R_PLANET_KM)
call model_s362ani_subshsv(xcolat,xlon,xrad,xdvsh,xdvsv,xdvph,xdvpv)
! to use speed values from the 1D reference model but with 3D mesh variations
if (USE_1D_REFERENCE) then
! sets all 3D variations in the mantle to zero
xdvpv = 0.d0
xdvph = 0.d0
xdvsv = 0.d0
xdvsh = 0.d0
endif
if (TRANSVERSE_ISOTROPY) then
! tiso perturbation
vpv = vpv*(1.0d0+dble(xdvpv))
vph = vph*(1.0d0+dble(xdvph))
vsv = vsv*(1.0d0+dble(xdvsv))
vsh = vsh*(1.0d0+dble(xdvsh))
else
! isotropic model
vpv = vpv+xdvpv
vph = vph+xdvph
vsv = vsv+xdvsv
vsh = vsh+xdvsh
! isotropic average (considers anisotropic parameterization eta,vsv,vsh,vpv,vph)
vp = sqrt(((8.d0+4.d0*eta_aniso)*vph*vph + 3.d0*vpv*vpv &
+ (8.d0 - 8.d0*eta_aniso)*vsv*vsv)/15.d0)
vs = sqrt(((1.d0-2.d0*eta_aniso)*vph*vph + vpv*vpv &
+ 5.d0*vsh*vsh + (6.d0+4.d0*eta_aniso)*vsv*vsv)/15.d0)
vpv = vp
vph = vp
vsv = vs
vsh = vs
eta_aniso = 1.0d0
endif
case (THREE_D_MODEL_PPM)
! point profile model
call model_PPM(r_used,theta,phi,dvs,dvp,drho)
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
rho = rho*(1.0d0+drho)
case (THREE_D_MODEL_GAPP2)
! 3D GAP model (Obayashi)
call mantle_gapmodel(r_used,theta,phi,dvs,dvp,drho)
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvs)
vsh = vsh*(1.0d0+dvs)
rho = rho*(1.0d0+drho)
case (THREE_D_MODEL_SGLOBE,THREE_D_MODEL_SGLOBE_ISO)
! 3D SGLOBE-rani model (Chang)
! normally mantle perturbations are taken from 24.4km (R_MOHO) up.
! we need to add the if statement for sgloberani_iso or sgloberani_aniso to take from 50km up:
if (r_prem > RCMB/R_PLANET .and. r_prem < 6321000.d0/R_PLANET) then
r_used = r
else ! if (r_prem >= 6321000.d0/R_PLANET) then
! this will then "extend the mantle up to the surface" from 50km depth
r_used = 6321000.d0/R_PLANET
endif
call mantle_sglobe(r_used,theta,phi,vsv,vsh,dvsv,dvsh,dvp,drho)
! updates velocities from reference model
if (TRANSVERSE_ISOTROPY) then
! tiso perturbation
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvsv)
vsh = vsh*(1.0d0+dvsh)
else
! isotropic model
vpv = vpv*(1.0d0+dvp)
vph = vph*(1.0d0+dvp)
vsv = vsv*(1.0d0+dvsv)
vsh = vsh*(1.0d0+dvsh)
! Voigt average
vp = sqrt( (2.d0*vpv**2 + vph**2)/3.d0 )
vs = sqrt( (2.d0*vsv**2 + vsh**2)/3.d0 )
vph = vp
vpv = vp
vsh = vs
vsv = vs
eta_aniso = 1.d0
endif
rho = rho*(1.0d0+drho)
case (THREE_D_MODEL_ANISO_MANTLE)
! Montagner anisotropic model
call model_aniso_mantle(r_used,theta,phi,vpv,vph,vsv,vsh,eta_aniso,rho, &
c11,c12,c13,c14,c15,c16,c22,c23,c24,c25,c26, &
c33,c34,c35,c36,c44,c45,c46,c55,c56,c66)
case (THREE_D_MODEL_BKMNS_GLAD)
! GLAD model expansion on block-mantle-spherical-harmonics
! model takes actual position (includes topography between 80km depth and surface topo)
r_used = r
call model_bkmns_mantle(r_used,theta,phi,vpv,vph,vsv,vsh,eta_aniso,rho)
case (THREE_D_MODEL_SPIRAL)
! SPiRaL v1.4 anisotropic model
lat = (PI/2.0d0-theta)*180.0d0/PI
lon = phi*180.0d0/PI
! input: lat in range [-90,90]; lon in range [-180,180]
if (lon > 180.0d0) lon = lon - 360.0d0
call model_mantle_spiral(r_used,lat,lon,vpv,vph,vsv,vsh,eta_aniso,rho, &
c11,c12,c13,c14,c15,c16,c22,c23,c24,c25,c26, &
c33,c34,c35,c36,c44,c45,c46,c55,c56,c66)
case (THREE_D_MODEL_BERKELEY)
! 3D Berkeley Model SEMUCB
call model_berkeley_shsv(dble(r_used*EARTH_R_KM),theta,phi,dvsh,dvsv, &
dvph,dvpv,drho,eta_aniso,iregion_code,CRUSTAL)
!
! call model_prem_iso(myrank,r_prem,rho,drhodr,vp,vs,Qkappa,&
! Qmu,idoubling,CRUSTAL, ONE_CRUST,.true.,RICB,RCMB, &
! RTOPDDOUBLEPRIME, R600,R670,R220,R771,R400,R80,RMOHO,RMIDDLE_CRUST,ROCEAN)
! to use speed values from the 1D reference model but with 3D mesh variations
if( USE_1D_REFERENCE ) then
! sets all 3D variations in the mantle to zero
dvpv = 0.d0
dvph = 0.d0
dvsv = 0.d0
dvsh = 0.d0
endif
if(TRANSVERSE_ISOTROPY) then
vpv=vpv*(1.0d0+dble(dvpv))
vph=vph*(1.0d0+dble(dvph))
vsv=vsv*(1.0d0+dble(dvsv))
vsh=vsh*(1.0d0+dble(dvsh))
rho=rho*(1.0d0+drho)
else
vpv=vpv+dvpv
vph=vph+dvph
vsv=vsv+dvsv
vsh=vsh+dvsh
vp = sqrt(((8.d0+4.d0*eta_aniso)*vph*vph + 3.d0*vpv*vpv &
+ (8.d0 - 8.d0*eta_aniso)*vsv*vsv)/15.d0)
vs = sqrt(((1.d0-2.d0*eta_aniso)*vph*vph + vpv*vpv &
+ 5.d0*vsh*vsh + (6.d0+4.d0*eta_aniso)*vsv*vsv)/15.d0)
vpv=vp
vph=vp
vsv=vs
vsh=vs
eta_aniso=1.0d0
rho=rho*(1.0d0+drho)
endif
case (THREE_D_MODEL_HETEROGEN_PREM)
! chris modif checkers 02/20/21
call model_heterogen_mantle(ispec,i,j,k,r_prem,theta,phi,dvs,dvp,drho)
! vpv = vpv*(1.0d0+dvp) ! correct format with initial heterogenous model
vpv = dvp
vph = dvp
vsv = dvs
vsh = dvs
rho = drho
case (THREE_D_MODEL_SH_MARS)
! Mars model expansion on spherical harmonics
r_used = r ! takes actual position (between CMB and surface)
call model_SH_mars_crustmantle(r_used,theta,phi,vpv,vph,vsv,vsh,eta_aniso,rho)
case default
print *,'Error: do not recognize value for THREE_D_MODEL ',THREE_D_MODEL
stop 'unknown 3D Earth model in meshfem3D_models_get3Dmntl_val(), please check... '
end select ! THREE_D_MODEL
endif ! MODEL_3D_MANTLE_PERTUBATIONS