-
Notifications
You must be signed in to change notification settings - Fork 93
/
setup_sources_receivers.f90
2104 lines (1763 loc) · 77.3 KB
/
setup_sources_receivers.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 setup_sources_receivers()
use specfem_par, only: myrank,IMAIN,NSOURCES,NSTEP, &
theta_source,phi_source, &
TOPOGRAPHY,ibathy_topo, &
USE_DISTANCE_CRITERION,xyz_midpoints,xadj,adjncy
use kdtree_search, only: kdtree_delete,kdtree_nodes_location,kdtree_nodes_index
implicit none
! setup for point search
call setup_point_search_arrays()
! locates sources and determines simulation start time t0
call setup_sources()
! reads in stations file and locates receivers
call setup_receivers()
! write source and receiver VTK files for Paraview
call setup_sources_receivers_VTKfile()
! pre-compute source arrays
call setup_sources_precompute_arrays()
! pre-compute receiver interpolation factors
call setup_receivers_precompute_intp()
! user output
if (myrank == 0) then
write(IMAIN,*)
write(IMAIN,*) 'Total number of samples for seismograms = ',NSTEP
write(IMAIN,*)
if (NSOURCES > 1) write(IMAIN,*) 'Using ',NSOURCES,' point sources'
call flush_IMAIN()
endif
call synchronize_all()
! frees arrays
deallocate(theta_source,phi_source)
! topography array no more needed
if (TOPOGRAPHY) then
if (allocated(ibathy_topo) ) deallocate(ibathy_topo)
endif
! frees memory
if (USE_DISTANCE_CRITERION) deallocate(xyz_midpoints)
deallocate(xadj,adjncy)
! deletes tree arrays
deallocate(kdtree_nodes_location)
deallocate(kdtree_nodes_index)
! deletes search tree nodes
call kdtree_delete()
end subroutine setup_sources_receivers
!
!-------------------------------------------------------------------------------------------------
!
subroutine setup_point_search_arrays()
use constants, only: &
NDIM,NGLLX,NGLLY,NGLLZ,MIDX,MIDY,MIDZ,IMAIN,TWO_PI,R_UNIT_SPHERE,DEGREES_TO_RADIANS, &
USE_DISTANCE_CRITERION
use specfem_par, only: &
NCHUNKS_VAL,NEX_XI_VAL,NEX_ETA_VAL,ANGULAR_WIDTH_XI_IN_DEGREES_VAL,ANGULAR_WIDTH_ETA_IN_DEGREES_VAL, &
LAT_LON_MARGIN,myrank
use specfem_par, only: &
nspec => NSPEC_CRUST_MANTLE,nglob => NGLOB_CRUST_MANTLE
use specfem_par_crustmantle, only: &
ibool => ibool_crust_mantle, &
xstore => xstore_crust_mantle,ystore => ystore_crust_mantle,zstore => zstore_crust_mantle
! for point search
use specfem_par, only: &
element_size,typical_size_squared, &
anchor_iax,anchor_iay,anchor_iaz, &
lat_min,lat_max,lon_min,lon_max,xyz_midpoints
use kdtree_search, only: kdtree_setup,kdtree_set_verbose, &
kdtree_num_nodes,kdtree_nodes_location,kdtree_nodes_index
implicit none
! local parameters
integer :: ispec,iglob,ier
integer :: i,j,k,inodes
double precision ANGULAR_WIDTH_XI_RAD,ANGULAR_WIDTH_ETA_RAD
! determines tree points
logical :: use_midpoints_only
! compute typical size of elements at the surface
! (normalized)
if (NCHUNKS_VAL == 6) then
! estimation for global meshes (assuming 90-degree chunks)
element_size = TWO_PI * R_UNIT_SPHERE / (4.d0 * NEX_XI_VAL)
else
! estimation for 1-chunk meshes
ANGULAR_WIDTH_XI_RAD = ANGULAR_WIDTH_XI_IN_DEGREES_VAL * DEGREES_TO_RADIANS
ANGULAR_WIDTH_ETA_RAD = ANGULAR_WIDTH_ETA_IN_DEGREES_VAL * DEGREES_TO_RADIANS
element_size = max( ANGULAR_WIDTH_XI_RAD/NEX_XI_VAL,ANGULAR_WIDTH_ETA_RAD/NEX_ETA_VAL ) * R_UNIT_SPHERE
endif
! use 10 times the distance as a criterion for source detection
typical_size_squared = (10.d0 * element_size)**2
! limits receiver search
if (USE_DISTANCE_CRITERION) then
! retrieves latitude/longitude range of this slice
call xyz_2_latlon_minmax(nspec,nglob,ibool,xstore,ystore,zstore,lat_min,lat_max,lon_min,lon_max)
! adds search margin
lat_min = lat_min - LAT_LON_MARGIN
lat_max = lat_max + LAT_LON_MARGIN
lon_min = lon_min - LAT_LON_MARGIN
lon_max = lon_max + LAT_LON_MARGIN
! limits latitude to [-90.0,90.0]
if (lat_min < -90.d0 ) lat_min = -90.d0
if (lat_max > 90.d0 ) lat_max = 90.d0
! limits longitude to [0.0,360.0]
if (lon_min < 0.d0 ) lon_min = 0.d0
if (lon_min > 360.d0 ) lon_min = 360.d0
if (lon_max < 0.d0 ) lon_max = 0.d0
if (lon_max > 360.d0 ) lon_max = 360.d0
! prepares midpoints coordinates
allocate(xyz_midpoints(NDIM,nspec),stat=ier)
if (ier /= 0 ) call exit_MPI(myrank,'Error allocating array xyz_midpoints')
! store x/y/z coordinates of center point
do ispec = 1,nspec
iglob = ibool(MIDX,MIDY,MIDZ,ispec)
xyz_midpoints(1,ispec) = dble(xstore(iglob))
xyz_midpoints(2,ispec) = dble(ystore(iglob))
xyz_midpoints(3,ispec) = dble(zstore(iglob))
enddo
endif
! define (i,j,k) indices of the control/anchor points
call hex_nodes_anchor_ijk(anchor_iax,anchor_iay,anchor_iaz)
! setups adjacency array to search neighbors
call setup_adjacency_neighbors()
! kd-tree setup for point localization
!
! determines tree size
if (NEX_ETA_VAL > 100 .and. NEX_XI_VAL > 100) then
! high-resolution mesh
! only midpoints for search, should be sufficient to get accurate location
use_midpoints_only = .true.
else
! low-resolution mesh
! uses element's inner points
use_midpoints_only = .false.
endif
! sets total number of tree points
if (use_midpoints_only) then
! small tree size
kdtree_num_nodes = nspec
else
! uses all internal GLL points for search tree
! internal GLL points ( 2 to NGLLX-1 )
kdtree_num_nodes = nspec * (NGLLX-2)*(NGLLY-2)*(NGLLZ-2)
endif
! allocates tree arrays
allocate(kdtree_nodes_location(NDIM,kdtree_num_nodes),stat=ier)
if (ier /= 0) stop 'Error allocating kdtree_nodes_location arrays'
allocate(kdtree_nodes_index(kdtree_num_nodes),stat=ier)
if (ier /= 0) stop 'Error allocating kdtree_nodes_index arrays'
! prepares search arrays, each element takes its internal GLL points for tree search
kdtree_nodes_index(:) = 0
kdtree_nodes_location(:,:) = 0.0
! adds tree nodes
inodes = 0
if (use_midpoints_only) then
! sets up tree nodes
do ispec = 1,nspec
iglob = ibool(MIDX,MIDY,MIDZ,ispec)
! counts nodes
inodes = inodes + 1
if (inodes > kdtree_num_nodes) stop 'Error index inodes bigger than kdtree_num_nodes'
! adds node index (index points to same ispec for all internal GLL points)
kdtree_nodes_index(inodes) = ispec
! adds node location
kdtree_nodes_location(1,inodes) = xstore(iglob)
kdtree_nodes_location(2,inodes) = ystore(iglob)
kdtree_nodes_location(3,inodes) = zstore(iglob)
enddo
else
! all internal GLL points
do ispec = 1,nspec
do k = 2,NGLLZ-1
do j = 2,NGLLY-1
do i = 2,NGLLX-1
iglob = ibool(i,j,k,ispec)
! counts nodes
inodes = inodes + 1
if (inodes > kdtree_num_nodes) stop 'Error index inodes bigger than kdtree_num_nodes'
! adds node index (index points to same ispec for all internal GLL points)
kdtree_nodes_index(inodes) = ispec
! adds node location
kdtree_nodes_location(1,inodes) = xstore(iglob)
kdtree_nodes_location(2,inodes) = ystore(iglob)
kdtree_nodes_location(3,inodes) = zstore(iglob)
enddo
enddo
enddo
enddo
endif
if (inodes /= kdtree_num_nodes) stop 'Error index inodes does not match kdtree_num_nodes'
! tree verbosity
if (myrank == 0) call kdtree_set_verbose(IMAIN)
! creates kd-tree for searching point locations in locate_point() routine
call kdtree_setup()
end subroutine setup_point_search_arrays
!
!-------------------------------------------------------------------------------------------------
!
subroutine setup_adjacency_neighbors()
use constants, only: &
NDIM,NGLLX,NGLLY,NGLLZ,MIDX,MIDY,MIDZ,IMAIN, &
USE_DISTANCE_CRITERION
use shared_parameters, only: R_PLANET_KM
use specfem_par, only: &
myrank, &
nspec => NSPEC_CRUST_MANTLE
use specfem_par_crustmantle, only: &
ibool => ibool_crust_mantle, &
xstore => xstore_crust_mantle,ystore => ystore_crust_mantle,zstore => zstore_crust_mantle
! for point search
use specfem_par, only: element_size,typical_size_squared,xyz_midpoints, &
xadj,adjncy,num_neighbors_all
use kdtree_search, only: kdtree_setup,kdtree_delete, &
kdtree_nodes_location,kdtree_nodes_index,kdtree_num_nodes, &
kdtree_count_nearest_n_neighbors,kdtree_get_nearest_n_neighbors, &
kdtree_search_index,kdtree_search_num_nodes
implicit none
integer :: num_neighbors,num_neighbors_max
integer,dimension(8) :: iglob_corner,iglob_corner2
integer :: ispec_ref,ispec,iglob,icorner,ier !,jj
! temporary
integer,parameter :: MAX_NEIGHBORS = 50 ! maximum number of neighbors (around 37 should be sufficient for crust/mantle)
integer,dimension(:),allocatable :: tmp_adjncy ! temporary adjacency
integer :: inum_neighbor
! timer MPI
double precision :: time1,tCPU
double precision, external :: wtime
! kd-tree search
integer :: ielem,inodes
integer :: nsearch_points
integer :: num_elements,num_elements_max
integer :: ielem_counter,num_elements_actual_max
!integer, parameter :: max_search_points = 2000
double precision :: r_search
double precision :: xyz_target(NDIM)
double precision :: dist_squared,dist_squared_max
logical :: is_neighbor
logical,parameter :: DO_BRUTE_FORCE_SEARCH = .false.
if (myrank == 0) then
write(IMAIN,*) 'adjacency:'
write(IMAIN,*) ' total number of elements in this slice = ',nspec
write(IMAIN,*)
call flush_IMAIN()
endif
! get MPI starting time
time1 = wtime()
! adjacency arrays
!
! how to use:
! num_neighbors = xadj(ispec+1)-xadj(ispec)
! do i = 1,num_neighbors
! ! get neighbor
! ispec_neighbor = adjncy(xadj(ispec) + i)
! ..
! enddo
allocate(xadj(nspec + 1),stat=ier)
if (ier /= 0) stop 'Error allocating xadj'
allocate(tmp_adjncy(MAX_NEIGHBORS*nspec),stat=ier)
if (ier /= 0) stop 'Error allocating tmp_adjncy'
xadj(:) = 0
tmp_adjncy(:) = 0
! kd-tree search
if (.not. DO_BRUTE_FORCE_SEARCH) then
! kd-tree search
! search radius around element midpoints
!
! note: typical search size is using 10 times the size of a surface element;
! we take here 6 times the surface element size
! - since at low resolutions NEX < 64 and large element sizes, this search radius needs to be large enough, and,
! - due to doubling layers (elements at depth will become bigger, however radius shrinks)
!
! the search radius r_search given as routine argument must be non-squared
r_search = 6.0 * element_size
! user output
if (myrank == 0) then
write(IMAIN,*) ' using kd-tree search radius = ',r_search * R_PLANET_KM,'(km)'
write(IMAIN,*)
call flush_IMAIN()
endif
! kd-tree setup for adjacency search
!
! uses only element midpoint location
kdtree_num_nodes = nspec
! allocates tree arrays
allocate(kdtree_nodes_location(NDIM,kdtree_num_nodes),stat=ier)
if (ier /= 0) stop 'Error allocating kdtree_nodes_location arrays'
allocate(kdtree_nodes_index(kdtree_num_nodes),stat=ier)
if (ier /= 0) stop 'Error allocating kdtree_nodes_index arrays'
! prepares search arrays, each element takes its internal GLL points for tree search
kdtree_nodes_index(:) = 0
kdtree_nodes_location(:,:) = 0.0
! adds tree nodes
inodes = 0
do ispec = 1,nspec
! sets up tree nodes
iglob = ibool(MIDX,MIDY,MIDZ,ispec)
! counts nodes
inodes = inodes + 1
if (inodes > kdtree_num_nodes ) stop 'Error index inodes bigger than kdtree_num_nodes'
! adds node index (index points to same ispec for all internal GLL points)
kdtree_nodes_index(inodes) = ispec
! adds node location
kdtree_nodes_location(1,inodes) = xstore(iglob)
kdtree_nodes_location(2,inodes) = ystore(iglob)
kdtree_nodes_location(3,inodes) = zstore(iglob)
enddo
if (inodes /= kdtree_num_nodes ) stop 'Error index inodes does not match nnodes_local'
! alternative: to avoid allocating/deallocating search index arrays, though there is hardly a speedup
!allocate(kdtree_search_index(max_search_points),stat=ier)
!if (ier /= 0) stop 'Error allocating array kdtree_search_index'
! creates kd-tree for searching
call kdtree_setup()
endif
! gets maximum number of neighbors
inum_neighbor = 0
num_neighbors_max = 0
num_neighbors_all = 0
num_elements_max = 0
num_elements_actual_max = 0
dist_squared_max = 0.d0
do ispec_ref = 1,nspec
! the eight corners of the current element
iglob_corner(1) = ibool(1,1,1,ispec_ref)
iglob_corner(2) = ibool(NGLLX,1,1,ispec_ref)
iglob_corner(3) = ibool(NGLLX,NGLLY,1,ispec_ref)
iglob_corner(4) = ibool(1,NGLLY,1,ispec_ref)
iglob_corner(5) = ibool(1,1,NGLLZ,ispec_ref)
iglob_corner(6) = ibool(NGLLX,1,NGLLZ,ispec_ref)
iglob_corner(7) = ibool(NGLLX,NGLLY,NGLLZ,ispec_ref)
iglob_corner(8) = ibool(1,NGLLY,NGLLZ,ispec_ref)
! midpoint for search radius
iglob = ibool(MIDX,MIDY,MIDZ,ispec_ref)
xyz_target(1) = xstore(iglob)
xyz_target(2) = ystore(iglob)
xyz_target(3) = zstore(iglob)
if (DO_BRUTE_FORCE_SEARCH) then
! loops over all other elements to find closest neighbors
num_elements = nspec
else
! looks only at elements in kd-tree search radius
! gets number of tree points within search radius
! (within search sphere)
call kdtree_count_nearest_n_neighbors(xyz_target,r_search,nsearch_points)
! debug
!print *,' total number of search elements: ',nsearch_points,'ispec',ispec_ref
! alternative: limits search results
!if (nsearch_points > max_search_points) nsearch_points = max_search_points
! sets number of search nodes to get
kdtree_search_num_nodes = nsearch_points
! allocates search index
allocate(kdtree_search_index(kdtree_search_num_nodes),stat=ier)
if (ier /= 0) stop 'Error allocating array kdtree_search_index'
! gets closest n points around target (within search sphere)
call kdtree_get_nearest_n_neighbors(xyz_target,r_search,nsearch_points)
! loops over search radius
num_elements = nsearch_points
endif
! statistics
if (num_elements > num_elements_max) num_elements_max = num_elements
ielem_counter = 0
! counts number of neighbors
num_neighbors = 0
do ielem = 1,num_elements
! gets element index
if (DO_BRUTE_FORCE_SEARCH) then
ispec = ielem
else
! kd-tree search radius
! gets search point/element index
ispec = kdtree_search_index(ielem)
! checks index
if (ispec < 1 .or. ispec > nspec) stop 'Error element index is invalid'
endif
! skip reference element
if (ispec == ispec_ref) cycle
! distance to reference element
dist_squared = (xyz_target(1) - xyz_midpoints(1,ispec))*(xyz_target(1) - xyz_midpoints(1,ispec)) &
+ (xyz_target(2) - xyz_midpoints(2,ispec))*(xyz_target(2) - xyz_midpoints(2,ispec)) &
+ (xyz_target(3) - xyz_midpoints(3,ispec))*(xyz_target(3) - xyz_midpoints(3,ispec))
! exclude elements that are too far from target
if (USE_DISTANCE_CRITERION) then
! we compare squared distances instead of distances themselves to significantly speed up calculations
if (dist_squared > typical_size_squared) cycle
endif
ielem_counter = ielem_counter + 1
! checks if element has a corner iglob from reference element
is_neighbor = .false.
iglob_corner2(1) = ibool(1,1,1,ispec)
iglob_corner2(2) = ibool(NGLLX,1,1,ispec)
iglob_corner2(3) = ibool(NGLLX,NGLLY,1,ispec)
iglob_corner2(4) = ibool(1,NGLLY,1,ispec)
iglob_corner2(5) = ibool(1,1,NGLLZ,ispec)
iglob_corner2(6) = ibool(NGLLX,1,NGLLZ,ispec)
iglob_corner2(7) = ibool(NGLLX,NGLLY,NGLLZ,ispec)
iglob_corner2(8) = ibool(1,NGLLY,NGLLZ,ispec)
do icorner = 1,8
! checks if corner also has reference element
if (any(iglob_corner(:) == iglob_corner2(icorner))) then
is_neighbor = .true.
exit
endif
! alternative: (slightly slower with 12.4s compared to 11.4s with any() intrinsic function)
!do jj = 1,8
! if (iglob == iglob_corner(jj)) then
! is_neighbor = .true.
! exit
! endif
!enddo
!if (is_neighbor) exit
enddo
! counts neighbors to reference element
if (is_neighbor) then
! adds to adjacency
inum_neighbor = inum_neighbor + 1
! checks
if (inum_neighbor > MAX_NEIGHBORS*nspec) stop 'Error maximum neighbors exceeded'
! adds element
tmp_adjncy(inum_neighbor) = ispec
! for statistics
num_neighbors = num_neighbors + 1
! maximum distance to reference element
if (dist_squared > dist_squared_max) dist_squared_max = dist_squared
endif
enddo ! ielem
! checks if neighbors were found
if (num_neighbors == 0) then
print *,'Error: rank ',myrank,' - element ',ispec_ref,'has no neighbors!'
print *,' element midpoint location: ',xyz_target(:)*R_PLANET_KM
print *,' typical element size : ',element_size*R_PLANET_KM,'(km)'
print *,' brute force search : ',DO_BRUTE_FORCE_SEARCH
print *,' distance criteria : ',USE_DISTANCE_CRITERION
print *,' typical search distance : ',typical_size_squared*R_PLANET_KM,'(km)'
print *,' kd-tree r_search : ',r_search*R_PLANET_KM,'(km)'
print *,' search elements : ',num_elements
call exit_MPI(myrank,'Error adjacency invalid')
endif
! statistics
if (num_neighbors > num_neighbors_max) num_neighbors_max = num_neighbors
if (ielem_counter > num_elements_actual_max) num_elements_actual_max = ielem_counter
! adjacency indexing
xadj(ispec_ref + 1) = inum_neighbor
! how to use:
!num_neighbors = xadj(ispec+1)-xadj(ispec)
!do i = 1,num_neighbors
! ! get neighbor
! ispec_neighbor = adjncy(xadj(ispec) + i)
!enddo
! frees kdtree search array
if (.not. DO_BRUTE_FORCE_SEARCH) then
deallocate(kdtree_search_index)
endif
enddo ! ispec_ref
! total number of neighbors
num_neighbors_all = inum_neighbor
! allocates compacted array
allocate(adjncy(num_neighbors_all),stat=ier)
if (ier /= 0) stop 'Error allocating tmp_adjncy'
adjncy(1:num_neighbors_all) = tmp_adjncy(1:num_neighbors_all)
! checks
if (minval(adjncy(:)) < 1 .or. maxval(adjncy(:)) > nspec) stop 'Invalid adjncy array'
! frees temporary array
deallocate(tmp_adjncy)
if (.not. DO_BRUTE_FORCE_SEARCH) then
! frees current tree memory
! deletes tree arrays
deallocate(kdtree_nodes_location)
deallocate(kdtree_nodes_index)
! deletes search tree nodes
call kdtree_delete()
endif
if (myrank == 0) then
! elapsed time since beginning of neighbor detection
tCPU = wtime() - time1
write(IMAIN,*) ' maximum search elements = ',num_elements_max
write(IMAIN,*) ' maximum of actual search elements (after distance criterion) = ',num_elements_actual_max
write(IMAIN,*)
write(IMAIN,*) ' estimated typical element size at surface = ',element_size*R_PLANET_KM,'(km)'
write(IMAIN,*) ' maximum distance between neighbor centers = ',sqrt(dist_squared_max)*R_PLANET_KM,'(km)'
if (.not. DO_BRUTE_FORCE_SEARCH) then
if (sqrt(dist_squared_max) > r_search - 0.5*element_size) then
write(IMAIN,*) '***'
write(IMAIN,*) '*** Warning: consider increasing the kd-tree search radius to improve this neighbor setup ***'
write(IMAIN,*) '***'
endif
endif
write(IMAIN,*)
write(IMAIN,*) ' maximum neighbors found per element = ',num_neighbors_max,'(should be 37 for globe meshes)'
write(IMAIN,*) ' total number of neighbors = ',num_neighbors_all
write(IMAIN,*)
write(IMAIN,*) ' Elapsed time for detection of neighbors in seconds = ',tCPU
write(IMAIN,*)
call flush_IMAIN()
endif
end subroutine setup_adjacency_neighbors
!
!-------------------------------------------------------------------------------------------------
!
subroutine setup_sources()
use specfem_par
use specfem_par_crustmantle
use specfem_par_movie
implicit none
! local parameters
integer :: isource,ier
character(len=MAX_STRING_LEN) :: filename
! user output
if (myrank == 0) then
write(IMAIN,*)
write(IMAIN,*) 'sources:',NSOURCES
call flush_IMAIN()
endif
! allocate arrays for source
allocate(islice_selected_source(NSOURCES), &
ispec_selected_source(NSOURCES), &
Mxx(NSOURCES), &
Myy(NSOURCES), &
Mzz(NSOURCES), &
Mxy(NSOURCES), &
Mxz(NSOURCES), &
Myz(NSOURCES),stat=ier)
if (ier /= 0 ) call exit_MPI(myrank,'Error allocating source arrays')
! initializes arrays
islice_selected_source(:) = -1
ispec_selected_source(:) = 0
Mxx(:) = 0.d0; Myy(:) = 0.d0; Mzz(:) = 0.d0
Mxy(:) = 0.d0; Mxz(:) = 0.d0; Myz(:) = 0.d0
allocate(xi_source(NSOURCES), &
eta_source(NSOURCES), &
gamma_source(NSOURCES),stat=ier)
if (ier /= 0 ) call exit_MPI(myrank,'Error allocating source arrays')
xi_source(:) = 0.d0; eta_source(:) = 0.d0; gamma_source(:) = 0.d0
allocate(tshift_src(NSOURCES), &
hdur(NSOURCES), &
hdur_Gaussian(NSOURCES),stat=ier)
if (ier /= 0 ) call exit_MPI(myrank,'Error allocating source arrays')
tshift_src(:) = 0.d0; hdur(:) = 0.d0; hdur_Gaussian(:) = 0.d0
allocate(theta_source(NSOURCES), &
phi_source(NSOURCES),stat=ier)
if (ier /= 0 ) call exit_MPI(myrank,'Error allocating source arrays')
theta_source(:) = 0.d0; phi_source(:) = 0.d0
allocate(nu_source(NDIM,NDIM,NSOURCES),stat=ier)
if (ier /= 0 ) call exit_MPI(myrank,'Error allocating source arrays')
nu_source(:,:,:) = 0.d0
if (USE_FORCE_POINT_SOURCE) then
allocate(force_stf(NSOURCES), &
factor_force_source(NSOURCES), &
comp_dir_vect_source_E(NSOURCES), &
comp_dir_vect_source_N(NSOURCES), &
comp_dir_vect_source_Z_UP(NSOURCES),stat=ier)
if (ier /= 0) stop 'error allocating arrays for force point sources'
force_stf(:) = 0
factor_force_source(:) = 0.d0
comp_dir_vect_source_E(:) = 0.d0
comp_dir_vect_source_N(:) = 0.d0
comp_dir_vect_source_Z_UP(:) = 0.d0
endif
! sources
! moved open statement and writing of first lines into sr.vtk before the
! call to locate_sources, where further write statements to that file follow
if (myrank == 0) then
! write source and receiver VTK files for Paraview
filename = trim(OUTPUT_FILES)//'/sr_tmp.vtk'
open(IOUT_VTK,file=trim(filename),status='unknown',iostat=ier)
if (ier /= 0 ) call exit_MPI(myrank,'Error opening temporary file sr_temp.vtk')
write(IOUT_VTK,'(a)') '# vtk DataFile Version 2.0'
write(IOUT_VTK,'(a)') 'Source and Receiver VTK file'
write(IOUT_VTK,'(a)') 'ASCII'
write(IOUT_VTK,'(a)') 'DATASET UNSTRUCTURED_GRID'
! LQY -- won't be able to know NSOURCES+nrec at this point...
write(IOUT_VTK, '(a,i6,a)') 'POINTS ', NSOURCES, ' float'
! closing file, rest of information will be appended later on
close(IOUT_VTK)
endif
! locate sources in the mesh
call locate_sources()
! determines onset time
call setup_stf_constants()
! count number of sources located in this slice
nsources_local = 0
if (SIMULATION_TYPE == 1 .or. SIMULATION_TYPE == 3) then
do isource = 1,NSOURCES
if (myrank == islice_selected_source(isource)) nsources_local = nsources_local + 1
enddo
endif
! determines number of times steps for simulation
call setup_timesteps()
! prints source time functions and spectrum to output files
if (PRINT_SOURCE_TIME_FUNCTION) call print_stf_file()
! get information about event name and location
! (e.g. needed for SAC seismograms)
! The following line is added for get_event_info subroutine.
! Because the way NSOURCES_SAC was declared has been changed.
! The rest of the changes in this program is just the updates of the subroutines that
! I did changes, e.g., adding/removing parameters. by Ebru Bozdag
call get_event_info_parallel(yr_SAC,jda_SAC,mo_SAC, da_SAC, ho_SAC,mi_SAC,sec_SAC, &
event_name_SAC,t_cmt_SAC,t_shift_SAC, &
elat_SAC,elon_SAC,depth_SAC,mb_SAC,ms_SAC,cmt_lat_SAC, &
cmt_lon_SAC,cmt_depth_SAC,cmt_hdur_SAC,NSOURCES, &
Mrr,Mtt,Mpp,Mrt,Mrp,Mtp)
! noise simulations ignore the CMTSOLUTIONS sources but employ a noise-spectrum source S_squared instead
! checks if anything to do for noise simulations
if (NOISE_TOMOGRAPHY /= 0) then
if (myrank == 0) then
write(IMAIN,*) 'noise simulation will ignore CMT sources'
call flush_IMAIN()
endif
endif
! syncs after source setup
call synchronize_all()
end subroutine setup_sources
!
!-------------------------------------------------------------------------------------------------
!
subroutine setup_stf_constants()
use specfem_par
use specfem_par_movie
implicit none
! local parameters
integer :: isource
! makes smaller hdur for movies
logical,parameter :: USE_SMALLER_HDUR_MOVIE = .false. ! by default off, to use same HDUR_MOVIE as specified in Par_file
if (abs(minval(tshift_src)) > TINYVAL) &
call exit_MPI(myrank,'one tshift_src must be zero, others must be positive')
! filter source time function by Gaussian with hdur = HDUR_MOVIE when writing movies or shakemaps
if (MOVIE_SURFACE .or. MOVIE_VOLUME) then
! smaller hdur_movie will do
if (USE_SMALLER_HDUR_MOVIE) then
! hdur_movie gets assigned an automatic value based on the simulation resolution
! this will make that a bit smaller to have a higher-frequency movie output
HDUR_MOVIE = 0.5 * HDUR_MOVIE
endif
! new hdur for simulation
hdur = sqrt(hdur**2 + HDUR_MOVIE**2)
if (myrank == 0) then
write(IMAIN,*)
write(IMAIN,*) 'Each source is being convolved with HDUR_MOVIE = ',sngl(HDUR_MOVIE)
write(IMAIN,*) 'Total source hdur = ',sngl(hdur),'(s)'
write(IMAIN,*)
call flush_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
if (USE_FORCE_POINT_SOURCE) then
! point force sources
! (might start depending on the frequency given by hdur)
! note: point force sources will give the dominant frequency in hdur, thus the main period is 1/hdur.
! also, these sources might use a Ricker source time function instead of a Gaussian.
! For a Ricker source time function, a start time ~1.2 * main_period is a good choice.
t0 = 0.d0
do isource = 1,NSOURCES
select case(force_stf(isource))
case (0)
! Gaussian source time function
t0 = min(t0,1.5d0 * (tshift_src(isource) - hdur(isource)))
case (1)
! Ricker source time function
t0 = min(t0,1.2d0 * (tshift_src(isource) - 1.0d0/hdur(isource)))
case (2)
! Heaviside
t0 = min(t0,1.5d0 * (tshift_src(isource) - hdur(isource)))
case (3)
! Monochromatic
t0 = 0.d0
case (4)
! Gaussian source time function by Meschede et al. (2011)
t0 = min(t0,1.5d0 * (tshift_src(isource) - hdur(isource)))
case default
stop 'unsupported force_stf value!'
end select
enddo
! start time defined as positive value, will be subtracted
t0 = - t0
else
! moment tensors
if (USE_MONOCHROMATIC_CMT_SOURCE) then
! (based on monochromatic functions)
t0 = 0.d0
else
! (based on Heaviside functions)
t0 = - 1.5d0 * minval( tshift_src(:) - hdur(:) )
endif
endif
! uses an external file for source time function, which starts at time 0.0
if (EXTERNAL_SOURCE_TIME_FUNCTION) then
hdur(:) = 0.d0
t0 = 0.d0
endif
! checks if user set USER_T0 to fix simulation start time
! note: USER_T0 has to be positive
if (USER_T0 > 0.d0) then
! user cares about origin time and time shifts of the CMTSOLUTION
! and wants to fix simulation start time to a constant start time
! time 0 on time axis will correspond to given origin time
! notifies user
if (myrank == 0) then
write(IMAIN,*) 'USER_T0: ',USER_T0
write(IMAIN,*) 't0: ',t0,'min_tshift_src_original: ',min_tshift_src_original
write(IMAIN,*)
endif
! checks if automatically set t0 is too small
! note: min_tshift_src_original can be a positive or negative time shift (minimum from all tshift)
if (t0 <= USER_T0 + min_tshift_src_original) then
! by default, tshift_src(:) holds relative time shifts with a minimum time shift set to zero
! re-adds (minimum) original time shift such that sources will kick in
! according to their absolute time shift
tshift_src(:) = tshift_src(:) + min_tshift_src_original
! sets new simulation start time such that
! simulation starts at t = - t0 = - USER_T0
t0 = USER_T0
! notifies user
if (myrank == 0) then
write(IMAIN,*) ' set new simulation start time: ', - t0
write(IMAIN,*)
endif
else
! start time needs to be at least t0 for numerical stability
! notifies user
if (myrank == 0) then
write(IMAIN,*) 'Error: USER_T0 is too small'
write(IMAIN,*) ' must make one of three adjustments:'
write(IMAIN,*) ' - increase USER_T0 to be at least: ',t0-min_tshift_src_original
write(IMAIN,*) ' - decrease time shift in CMTSOLUTION file'
write(IMAIN,*) ' - decrease hdur in CMTSOLUTION file'
call flush_IMAIN()
endif
call exit_mpi(myrank,'Error USER_T0 is set but too small')
endif
else if (USER_T0 < 0.d0) then
if (myrank == 0) then
write(IMAIN,*) 'Error: USER_T0 is negative, must be set zero or positive!'
endif
call exit_mpi(myrank,'Error negative USER_T0 parameter in constants.h')
endif
end subroutine setup_stf_constants
!
!-------------------------------------------------------------------------------------------------
!
subroutine setup_timesteps()
use specfem_par
implicit none
! local parameters
logical :: is_initial_guess
! checks if set by initial guess from read_compute_parameters() routine
if (NTSTEP_BETWEEN_OUTPUT_SEISMOS == NSTEP) then
is_initial_guess = .true.
else
is_initial_guess = .false.
endif
! from initial guess in read_compute_parameters:
! compute total number of time steps, rounded to next multiple of 100
! NSTEP = 100 * (int(RECORD_LENGTH_IN_MINUTES * 60.d0 / (100.d0*DT)) + 1)
!
! adds initial t0 time to update number of time steps and reach full record length
if (abs(t0) > 0.d0) then
! note: for zero length, nstep has minimal of 5 timesteps for testing
! we won't extend this
!
! careful: do not use RECORD_LENGTH_IN_MINUTES here, as it is only read by the main process
! when reading the parameter file, but it is not broadcasted to all other processes
! NSTEP gets broadcasted, so we work with this values
if (NSTEP /= 5) then
! extend by bulk of 100 steps to account for half-duration rise time
NSTEP = NSTEP + 100 * (int( abs(t0) / (100.d0*DT)) + 1)
endif
endif
! if doing benchmark runs to measure scaling of the code for a limited number of time steps only
if (DO_BENCHMARK_RUN_ONLY) NSTEP = NSTEP_FOR_BENCHMARK
! overrides NSTEP in case specified in Par_file
if (USER_NSTEP > 0) then
! overrides NSTEP
NSTEP = USER_NSTEP
endif
! checks length for symmetry in case of noise simulations
if (NOISE_TOMOGRAPHY /= 0) then
if (mod(NSTEP+1,2) /= 0) then
print *,'Error noise simulation: invalid time steps = ',NSTEP,', NSTEP + 1 must be a multiple of 2 due to branch symmetry'
call exit_MPI(myrank,'Error noise simulation: number of timesteps must be symmetric, due to +/- branches')
endif
endif
! make sure NSTEP is a multiple of NTSTEP_BETWEEN_OUTPUT_SAMPLE
! if not, increase it a little bit, to the next multiple
if (mod(NSTEP,NTSTEP_BETWEEN_OUTPUT_SAMPLE) /= 0) then
if (NOISE_TOMOGRAPHY /= 0) then
if (myrank == 0) then
print *,'Noise simulation: Invalid number of NSTEP = ',NSTEP
print *,'Must be a multiple of NTSTEP_BETWEEN_OUTPUT_SAMPLE = ',NTSTEP_BETWEEN_OUTPUT_SAMPLE
endif
stop 'Error: NSTEP must be a multiple of NTSTEP_BETWEEN_OUTPUT_SAMPLE'
else
NSTEP = (NSTEP/NTSTEP_BETWEEN_OUTPUT_SAMPLE + 1)*NTSTEP_BETWEEN_OUTPUT_SAMPLE
! user output
if (myrank == 0) then
print *
print *,'NSTEP is not a multiple of NTSTEP_BETWEEN_OUTPUT_SAMPLE'
print *,'thus increasing it automatically to the next multiple, which is ',NSTEP
print *
endif
endif
endif
! output seismograms at least once at the end of the simulation
NTSTEP_BETWEEN_OUTPUT_SEISMOS = min(NSTEP,NTSTEP_BETWEEN_OUTPUT_SEISMOS)
! make sure NSTEP_BETWEEN_OUTPUT_SEISMOS is a multiple of NTSTEP_BETWEEN_OUTPUT_SAMPLE
if (mod(NTSTEP_BETWEEN_OUTPUT_SEISMOS,NTSTEP_BETWEEN_OUTPUT_SAMPLE) /= 0) then
if (myrank == 0) then
print *,'Invalid number of NTSTEP_BETWEEN_OUTPUT_SEISMOS = ',NTSTEP_BETWEEN_OUTPUT_SEISMOS
print *,'Must be a multiple of NTSTEP_BETWEEN_OUTPUT_SAMPLE = ',NTSTEP_BETWEEN_OUTPUT_SAMPLE
endif