SPEED
MAKE_NOTHONORING.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
19
43
44 subroutine make_nothonoring(loc_n_num, nn_loc,&
45 n_case, tag_case, val_case, tol_case, &
46 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
47 xs_loc, ys_loc, zs_loc, zs_elev, zs_all, vs, thick, sub_tag_all,mpi_id)
48
49! SUBROUTINE FOR CASEs
50! COMPUTATION OF THE ELEVATION
51!
52! (FILE -> 'XYZ.out')
53! z_elev
54! +--------------------x----------------------+ +
55! | - vvvvvvvvvvvvvv|vvvvvvvvvvvvvv- | | zz_spx_elevation (depth) == zs_elev
56! | --vvvvvvvvvvvvv| z_spx vvvvv-- | |
57! | ---vvvvvvvv x vvvvvvvv--- | +
58! | -----vvvv|vvvv----- | | zz_spx_alluvial (all_depth) == zs_all
59! | ----x---- | +
60! | z_alluvial |
61! | (FILE -> 'ALL.out') |
62! | |
63! | tag_mat = val_case |
64! +-------------------------------------------+
65!
66! v = ALLUVIAL BASIN MATERIAL
67! depth = z_elev - z_spx -> zz_spx_elevation == zs_elev
68! all_depth = z_spx - z_alluvial -> zz_spx_alluvial == zs_all
69!
70
71
72 character*70 :: file_case_xyz
73 character*70 :: file_case_all
74 character*70 :: file_case_vs
75
76 !character*70 :: file_nhe_proc, file_nhe_new
77 !integer*4 :: inode, unit_mpi
78 !real*8 :: rho,lambda,mu,gamma,qs,qp
79
80 integer*4 :: n_case, nn_loc, cs_nnz_loc, nm, mpi_id
81 integer*4 :: ncase,vcase,tcase
82 integer*4 :: n_elev,n_tria_elev
83 integer*4 :: start,finish
84 integer*4 :: n_all,n_tria_all, ival, icase
85 !integer*4 :: tag_case, val_case
86
87 integer*4, dimension (:), allocatable :: node1_all,node2_all,node3_all
88 integer*4, dimension (:), allocatable :: node1_elev,node2_elev,node3_elev
89 integer*4, dimension(3) :: clock
90 integer*4, dimension(0:cs_nnz_loc) :: cs_loc
91 integer*4, dimension(nm) :: tag_mat, sdeg_mat
92 integer*4, dimension(nn_loc) :: loc_n_num
93 integer*4, dimension(n_case) :: tag_case
94 integer*4, dimension(n_case) :: val_case
95 integer*4, dimension(nn_loc), intent(inout) :: sub_tag_all
96
97 real*8 :: tolerance
98 real*8 :: max_elev_spacing,max_all_spacing
99
100 real*8, dimension(n_case) :: tol_case
101 real*8, dimension (:), allocatable :: x_elev,y_elev,z_elev, vs_elev, sedim
102 real*8, dimension (:), allocatable :: x_all,y_all,z_all
103 real*8, dimension(nn_loc) :: xs_loc, ys_loc, zs_loc
104 real*8, dimension(nn_loc), intent(inout) :: zs_elev, zs_all, vs, thick
105
106
107 !Initialization for vs for the new not-honoring strategy
108 vs = 0.d0
109 thick = 0.d0
110
111
112!*************************************************************************************************
113! GRENOBLE - honoring
114!*************************************************************************************************
115 ncase = n_case
116
117
118 tcase = tag_case(1) !tag_case(ncase)
119
120
121
122 if (tcase.eq.1 .or. tcase .eq. 91) then
123
124 if (mpi_id.eq.0 .and. tcase .eq. 1) then
125 write(*,'(A)')
126 write(*,'(A)')'CASE 1: GRENOBLE honoring'
127 write(*,'(A)')'Reading Topography...'
128 elseif (mpi_id.eq.0 .and. tcase .eq. 91) then
129 write(*,'(A)')
130 write(*,'(A)')'CASE 91: IRPINIA'
131 write(*,'(A)')'Reading Topography...'
132 endif
133
134 file_case_xyz ='XYZ.out'
135
136 zs_elev = -1.0e+30
137
138
139 call read_dime_filexyz(file_case_xyz,n_elev,n_tria_elev)
140
141 allocate(x_elev(n_elev),y_elev(n_elev),z_elev(n_elev))
142 allocate(node1_elev(n_tria_elev))
143 allocate(node2_elev(n_tria_elev))
144 allocate(node3_elev(n_tria_elev))
145
146 call read_filexyz(file_case_xyz,n_elev,n_tria_elev,&
147 x_elev,y_elev,z_elev,&
148 node1_elev,node2_elev,node3_elev,&
149 max_elev_spacing)
150
151 do icase = 1, ncase
152 call get_node_depth_from_simple(loc_n_num, n_elev, n_tria_elev,&
153 x_elev, y_elev, z_elev,&
154 node1_elev, node2_elev, node3_elev,&
155 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat,&
156 nn_loc, xs_loc, ys_loc, zs_loc,&
157 zs_elev, val_case(icase), max_elev_spacing, tol_case(icase))
158 enddo
159 deallocate(x_elev,y_elev,z_elev,node1_elev,node2_elev,node3_elev)
160
161 if (mpi_id.eq.0) then
162 write(*,'(A)')'Done'
163 write(*,'(A)')
164 endif
165
166
167!*************************************************************************************************
168! Calssical not honoring (1 ALL and 1 XYZ file)
169!*************************************************************************************************
170
171
172 elseif (tcase .eq. 2 .or. tcase .eq. 3 .or. tcase .eq. 4 .or. tcase .eq. 6 &
173 .or. tcase .eq. 7 .or. tcase .eq. 8 .or. tcase .eq. 11 .or. tcase .eq. 12 &
174 .or. tcase .eq. 13 .or. tcase .eq. 14 .or. tcase .eq. 15 .or. tcase .eq. 18 &
175 .or. tcase .eq. 22 .or. tcase .eq. 27 .or. tcase .eq. 28 .or. tcase .eq. 40 &
176 .or. tcase .eq. 33 .or. tcase .eq. 38 .or. tcase .eq. 46 .or. tcase .eq. 60) then
177
178 if (mpi_id.eq. 0 .and. tcase .eq. 2) then
179 write(*,'(A)')
180 write(*,'(A)')'CASE 2: GRENOBLE'
181
182 elseif(mpi_id .eq. 0 .and. tcase .eq. 3) then
183 write(*,'(A)')
184 write(*,'(A)')'CASE 3: GUBBIO'
185
186 elseif(mpi_id .eq. 0 .and. tcase .eq. 4) then
187 write(*,'(A)')
188 write(*,'(A)')'CASE 4: SULMONA'
189
190 elseif(mpi_id .eq. 0 .and. tcase .eq. 6) then
191 write(*,'(A)')
192 write(*,'(A)')'CASE 6: FRIULI'
193
194 elseif(mpi_id .eq. 0 .and. tcase .eq. 7) then
195 write(*,'(A)')
196 write(*,'(A)')'CASE 7: AQUILA'
197
198 elseif(mpi_id .eq. 0 .and. tcase .eq. 8) then
199 write(*,'(A)')
200 write(*,'(A)')'CASE 8: SANTIAGO'
201
202 elseif(mpi_id .eq. 0 .and. tcase .eq. 11) then
203 write(*,'(A)')
204 write(*,'(A)')'CASE 11: CHRISTCHURCH'
205
206 elseif(mpi_id .eq. 0 .and. tcase .eq. 12) then
207 write(*,'(A)')
208 write(*,'(A)')'CASE 12: PO PLAIN'
209
210 elseif(mpi_id .eq. 0 .and. tcase .eq. 13) then
211 write(*,'(A)')
212 write(*,'(A)')'CASE 13: PO PLAIN-BEDROCK'
213
214 elseif(mpi_id .eq. 0 .and. tcase .eq. 14) then
215 write(*,'(A)')
216 write(*,'(A)')'CASE 14: WELLINGTON'
217
218 elseif(mpi_id .eq. 0 .and. tcase .eq. 15) then
219 write(*,'(A)')
220 write(*,'(A)')'CASE 15: MARSICA'
221
222 elseif(mpi_id .eq. 0 .and. tcase .eq. 18) then
223 write(*,'(A)')
224 write(*,'(A)')'CASE 18: BEIJING-TUTORIAL'
225
226 elseif(mpi_id .eq. 0 .and. tcase .eq. 22) then
227 write(*,'(A)')
228 write(*,'(A)')'CASE 22: NORCIA'
229
230 elseif(mpi_id .eq. 0 .and. tcase .eq. 23) then
231 write(*,'(A)')
232 write(*,'(A)')'CASE 33: GRONINGEN-ZE'
233
234 elseif(mpi_id .eq. 0 .and. tcase .eq. 27) then
235 write(*,'(A)')
236 write(*,'(A)')'CASE 27: AQUILA-OB'
237
238 elseif(mpi_id .eq. 0 .and. tcase .eq. 28) then
239 write(*,'(A)')
240 write(*,'(A)')'CASE 28: NORCIA-OB'
241
242 elseif(mpi_id .eq. 0 .and. tcase .eq. 38) then
243 write(*,'(A)')
244 write(*,'(A)')'CASE 38: MONTELIMAR'
245
246 elseif(mpi_id .eq. 0 .and. tcase .eq. 40) then
247 write(*,'(A)')
248 write(*,'(A)')'CASE 40: KUTCH'
249
250 elseif(mpi_id .eq. 0 .and. tcase .eq. 46) then
251 write(*,'(A)')
252 write(*,'(A)')'CASE 46: KUMAMOTO'
253
254 elseif(mpi_id .eq. 0 .and. tcase .eq. 50) then
255 write(*,'(A)')
256 write(*,'(A)')'CASE 60: JAKARTA'
257 endif
258
259
260 if(mpi_id .eq. 0) write(*,'(A)')'Reading Topography&Alluvial...'
261 file_case_xyz ='XYZ.out'
262 file_case_all ='ALL.out'
263
264 zs_elev = -1.0e+30
265 zs_all = -1.0e+30
266
267 call read_dime_filexyz(file_case_xyz,n_elev,n_tria_elev)
268 call read_dime_filexyz(file_case_all,n_all,n_tria_all)
269
270 allocate(x_elev(n_elev),y_elev(n_elev),z_elev(n_elev))
271 allocate(node1_elev(n_tria_elev), node2_elev(n_tria_elev), node3_elev(n_tria_elev))
272
273 allocate(x_all(n_all),y_all(n_all),z_all(n_all))
274 allocate(node1_all(n_tria_all),node2_all(n_tria_all),node3_all(n_tria_all))
275
276 call read_filexyz(file_case_xyz,n_elev,n_tria_elev,&
277 x_elev,y_elev,z_elev,&
278 node1_elev,node2_elev,node3_elev,&
279 max_elev_spacing)
280
281 call read_filexyz(file_case_all,n_all,n_tria_all,&
282 x_all,y_all,z_all,&
283 node1_all,node2_all,node3_all,&
284 max_all_spacing)
285
286 do icase = 1, ncase
287
288 call get_node_depth_from_alluvial(loc_n_num, n_all, n_tria_all, &
289 x_all, y_all, z_all, &
290 node1_all, node2_all, node3_all,&
291 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
292 nn_loc, xs_loc, ys_loc, zs_loc, &
293 zs_all, val_case(icase), max_all_spacing, tol_case(icase))
294
295 call get_node_depth_from_cmplx(loc_n_num, n_elev, n_tria_elev, &
296 x_elev, y_elev, z_elev, &
297 node1_elev, node2_elev, node3_elev,&
298 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
299 nn_loc, xs_loc, ys_loc, zs_loc, &
300 zs_elev, zs_all, &
301 val_case(icase), max_elev_spacing, tol_case(icase))
302
303 enddo
304
305 deallocate(x_elev, y_elev, z_elev, node1_elev, node2_elev, node3_elev)
306 deallocate(x_all, y_all, z_all, node1_all, node2_all, node3_all)
307
308 if (mpi_id.eq.0) then
309 write(*,'(A)')'Done'
310 write(*,'(A)')
311 endif
312
313!*************************************************************************************************
314! L'AQUILA, MEXICO-CITY MULTI BASIN
315!*************************************************************************************************
316
317 elseif (tcase.eq. 70 .or. tcase .eq. 45) then
318 if (mpi_id .eq. 0 .and. tcase .eq. 70 ) then
319 write(*,'(A)')
320 write(*,'(A)')'CASE 70: Aquila-multibasin'
321
322 elseif (mpi_id .eq. 0 .and. tcase .eq. 45) then
323 write(*,'(A)')
324 write(*,'(A)')'CASE 45: MEXICO-CITY'
325
326 endif
327
328
329 if(mpi_id .eq. 0) write(*,'(A)')'Reading Topography&Alluvial...'
330
331 file_case_xyz ='XYZ.out'
332
333 zs_elev = -1.0e+30
334 zs_all = 0.d0;
335
336 call read_dime_filexyz(file_case_xyz,n_elev,n_tria_elev)
337 allocate(x_elev(n_elev),y_elev(n_elev),z_elev(n_elev))
338 allocate(node1_elev(n_tria_elev), node2_elev(n_tria_elev), node3_elev(n_tria_elev))
339
340 call read_filexyz(file_case_xyz,n_elev,n_tria_elev,&
341 x_elev,y_elev,z_elev,&
342 node1_elev,node2_elev,node3_elev,&
343 max_elev_spacing)
344
345 call get_node_depth_from_cmplx(loc_n_num, n_elev, n_tria_elev, &
346 x_elev, y_elev, z_elev, &
347 node1_elev, node2_elev, node3_elev,&
348 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
349 nn_loc, xs_loc, ys_loc, zs_loc, &
350 zs_elev, zs_all, &
351 val_case(1), max_elev_spacing, tol_case(1))
352
353 deallocate(x_elev, y_elev, z_elev, node1_elev, node2_elev, node3_elev)
354
355 sub_tag_all = 3
356 ival = 3
357
358 do j = 1,2
359 if (j.eq.1) then
360 file_case_all ='ALL1.out'
361 else
362 file_case_all ='ALL2.out'
363 endif
364
365 zs_all = -1.0e+30
366
367 call read_dime_filexyz(file_case_all,n_all,n_tria_all)
368
369
370 allocate(x_all(n_all), y_all(n_all), z_all(n_all))
371 allocate(node1_all(n_tria_all), node2_all(n_tria_all), node3_all(n_tria_all))
372
373 call read_filexyz(file_case_all,n_all,n_tria_all,&
374 x_all,y_all,z_all,&
375 node1_all,node2_all,node3_all,&
376 max_all_spacing)
377
378
379 do icase = 1, ncase
380
381 call get_node_depth_from_alluvial(loc_n_num, n_all, n_tria_all, &
382 x_all, y_all, z_all, &
383 node1_all, node2_all, node3_all,&
384 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
385 nn_loc, xs_loc, ys_loc, zs_loc, &
386 zs_all, val_case(icase), max_all_spacing, tol_case(icase))
387 enddo
388
389 call make_subtag_alluvial(nn_loc, zs_all, j, sub_tag_all, xs_loc, ival)
390
391 deallocate(x_all, y_all, z_all, node1_all, node2_all, node3_all)
392
393 if (mpi_id.eq.0) then
394 write(*,'(A)')
395 write(*,'(A,I8)') 'ALLUVIAL Layer # ',j
396 endif
397
398 enddo
399
400 if (mpi_id.eq.0) then
401 write(*,'(A)') 'Done'
402 write(*,'(A)')
403 endif
404
405!*************************************************************************************************
406! MULTI - NOT honoring
407!*************************************************************************************************
408
409 elseif (tcase.eq. 5 .or. tcase .eq. 50) then
410 if (mpi_id.eq.0 .and. tcase .eq. 5) then
411 write(*,'(A)')
412 write(*,'(A)')'CASE 5: VOLVI for CASHIMA benchmark'
413
414 elseif (mpi_id.eq.0 .and. tcase .eq. 50) then
415 write(*,'(A)')
416 write(*,'(A)')'CASE 50: PLANE-WAVE benchmark'
417
418 endif
419
420 write(*,'(A)')'Reading Topography&Alluvial...'
421 sub_tag_all = 4
422 ival = 4
423
424 do j = 1,3
425 if (j.eq.1) then
426 file_case_all ='ALL1.out'
427 elseif (j.eq.2) then
428 file_case_all ='ALL2.out'
429 else
430 file_case_all ='ALL3.out'
431 endif
432
433 zs_all = -1.0e+30
434
435 call read_dime_filexyz(file_case_all,n_all,n_tria_all)
436
437 allocate(x_all(n_all), y_all(n_all), z_all(n_all))
438 allocate(node1_all(n_tria_all), node2_all(n_tria_all), node3_all(n_tria_all))
439
440 call read_filexyz(file_case_all,n_all,n_tria_all,&
441 x_all,y_all,z_all,&
442 node1_all,node2_all,node3_all,&
443 max_all_spacing)
444
445 do icase = 1, ncase
446
447 call get_node_depth_from_alluvial(loc_n_num, n_all, n_tria_all, &
448 x_all, y_all, z_all, &
449 node1_all, node2_all, node3_all,&
450 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
451 nn_loc, xs_loc, ys_loc, zs_loc, &
452 zs_all, val_case(icase), max_all_spacing, tol_case(icase))
453 enddo
454
455 call make_subtag_alluvial(nn_loc, zs_all, j, sub_tag_all, xs_loc, ival)
456
457 deallocate(x_all, y_all, z_all, node1_all, node2_all, node3_all)
458
459 if (mpi_id.eq.0) then
460 write(*,'(A)')
461 write(*,'(A,I8)') 'ALLUVIAL Layer # ',j
462 endif
463
464 enddo
465
466 if (mpi_id.eq.0) then
467 write(*,'(A)') 'Done'
468 write(*,'(A)')
469 endif
470
471
472!*************************************************************************************************
473! XYZ map - ALL map - VS30 map
474!*************************************************************************************************
475
476 elseif (tcase.eq. 16 .or. tcase.eq. 19 .or. tcase .eq. 20 .or. tcase .eq. 21 &
477 .or. tcase .eq. 29 .or. tcase .eq. 35) then
478
479 if (mpi_id.eq.0 .and. tcase .eq. 16) then
480 write(*,'(A)')
481 write(*,'(A)')'CASE 16: ISTANBUL'
482 endif
483 if (mpi_id.eq.0 .and. tcase .eq. 19) then
484 write(*,'(A)')
485 write(*,'(A)')'CASE 19: THESSALONIKI'
486 endif
487 if (mpi_id.eq.0 .and. tcase .eq. 20) then
488 write(*,'(A)')
489 write(*,'(A)')'CASE 20: ATHENS'
490 endif
491 if (mpi_id.eq.0 .and. tcase .eq. 21) then
492 write(*,'(A)')
493 write(*,'(A)')'CASE 21: BEIJING '
494 endif
495 if (mpi_id.eq.0 .and. tcase .eq. 29) then
496 write(*,'(A)')
497 write(*,'(A)')'CASE 29: THESS-BEDROCK'
498 endif
499 if (mpi_id.eq.0 .and. tcase .eq. 35) then
500 write(*,'(A)')
501 write(*,'(A)')'CASE 35: THESS+MYGD-FINAL'
502 endif
503
504 if (mpi_id.eq.0) write(*,'(A)')'Reading Topography&Alluvial...'
505
506
507 file_case_xyz ='XYZ.out'
508 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) file_case_all ='ALL.out'
509 file_case_vs = 'VS_RS.out'
510
511 zs_elev = 0.d0
512 zs_all = 1.d0
513 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) zs_all = -1.0e+30
514
515 call read_dime_filexyz(file_case_xyz,n_elev,n_tria_elev)
516
517 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) call read_dime_filexyz(file_case_all,n_all,n_tria_all)
518
519 allocate(x_elev(n_elev),y_elev(n_elev),z_elev(n_elev),&
520 vs_elev(n_tria_elev),sedim(n_tria_elev))
521 allocate(node1_elev(n_tria_elev), node2_elev(n_tria_elev), node3_elev(n_tria_elev))
522
523 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) allocate(x_all(n_all),y_all(n_all),z_all(n_all))
524 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) allocate(node1_all(n_tria_all),node2_all(n_tria_all),node3_all(n_tria_all))
525
526 call read_filexyz(file_case_xyz,n_elev,n_tria_elev,&
527 x_elev,y_elev,z_elev,&
528 node1_elev,node2_elev,node3_elev,&
529 max_elev_spacing)
530
531 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) &
532 call read_filexyz(file_case_all,n_all,n_tria_all,&
533 x_all,y_all,z_all,&
534 node1_all,node2_all,node3_all,&
535 max_all_spacing)
536
537
538 call read_filevs(file_case_vs, n_tria_elev, vs_elev, sedim)
539
540
541 do icase = 1, ncase
542 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) call get_node_depth_from_alluvial(loc_n_num, n_all, n_tria_all, &
543 x_all, y_all, z_all, &
544 node1_all, node2_all, node3_all,&
545 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
546 nn_loc, xs_loc, ys_loc, zs_loc, &
547 zs_all, val_case(icase), max_all_spacing, tol_case(icase))
548
549 call get_node_depth_and_vs(loc_n_num, n_elev, n_tria_elev, &
550 x_elev, y_elev, z_elev, vs_elev, sedim,&
551 node1_elev, node2_elev, node3_elev,&
552 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
553 nn_loc, xs_loc, ys_loc, zs_loc, &
554 zs_elev, zs_all, vs, thick, &
555 val_case(icase), max_elev_spacing, tol_case(icase))
556 enddo
557
558
559 deallocate(x_elev, y_elev, z_elev,vs_elev,sedim, node1_elev, node2_elev, node3_elev)
560 if(tcase .eq. 19 .or. tcase .eq. 21 .or. tcase .eq. 29) deallocate(x_all, y_all, z_all, node1_all, node2_all, node3_all)
561
562 if (mpi_id.eq.0) then
563 write(*,'(A)')'Done'
564 write(*,'(A)')
565 endif
566
567
568
569
570!*************************************************************************************************
571! Atene - Partenone
572!*************************************************************************************************
573
574 elseif (tcase.eq. 30) then
575 if (mpi_id.eq.0) then
576 write(*,'(A)')
577 write(*,'(A)')'CASE 30: ATHENS-PARTHENON'
578 write(*,'(A)')'Reading Topography&Alluvial...'
579 endif
580
581
582 sub_tag_all = 4
583 ival = 4
584
585 do j = 1,3
586 if (j.eq.1) then
587 file_case_all ='ALL1.out'
588 elseif (j.eq.2) then
589 file_case_all ='ALL2.out'
590 else
591 file_case_all ='ALL3.out'
592 endif
593
594 zs_all = -1.0e+30
595
596 call read_dime_filexyz(file_case_all,n_all,n_tria_all)
597
598 allocate(x_all(n_all), y_all(n_all), z_all(n_all))
599 allocate(node1_all(n_tria_all), node2_all(n_tria_all), node3_all(n_tria_all))
600
601 call read_filexyz(file_case_all,n_all,n_tria_all,&
602 x_all,y_all,z_all,&
603 node1_all,node2_all,node3_all,&
604 max_all_spacing)
605
606 do icase = 1, ncase
607 call get_node_depth_from_alluvial(loc_n_num, n_all, n_tria_all, &
608 x_all, y_all, z_all, &
609 node1_all, node2_all, node3_all,&
610 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
611 nn_loc, xs_loc, ys_loc, zs_loc, &
612 zs_all, val_case(icase), max_all_spacing, tol_case(icase))
613 enddo
614
615 call make_subtag_alluvial(nn_loc, zs_all, j, sub_tag_all, xs_loc, ival)
616
617 deallocate(x_all, y_all, z_all, node1_all, node2_all, node3_all)
618
619 if (mpi_id.eq.0) then
620 write(*,'(A)')
621 write(*,'(A,I8)') 'ALLUVIAL Layer # ',j
622 endif
623
624 enddo
625
626
627
628
629 if (mpi_id.eq.0) then
630 write(*,'(A)') 'Done'
631 write(*,'(A)')
632 endif
633
634
635!*************************************************************************************************
636! Groningen
637!*************************************************************************************************
638
639 elseif (tcase.eq. 31 .or. tcase .eq. 32) then
640 if (mpi_id.eq.0) then
641 write(*,'(A)')
642 if(tcase.eq. 31) write(*,'(A)')'CASE 31: GRONINGEN'
643 if(tcase.eq. 32) write(*,'(A)')'CASE 32: GRONINGEN'
644 write(*,'(A)')'Reading Topography&Alluvial...'
645 endif
646
647
648 sub_tag_all = 7
649 ival = 7
650
651 do j = 1,6
652 if (j.eq.1) then
653 file_case_all ='ALL1.out'
654 elseif(j.eq.2) then
655 file_case_all ='ALL2.out'
656 elseif(j.eq.3) then
657 file_case_all ='ALL3.out'
658 elseif(j.eq.4) then
659 file_case_all ='ALL4.out'
660 elseif(j.eq.5) then
661 file_case_all ='ALL5.out'
662 else
663 file_case_all ='ALL6.out'
664 endif
665
666 zs_all = -1.0e+30
667
668 call read_dime_filexyz(file_case_all,n_all,n_tria_all)
669
670 allocate(x_all(n_all), y_all(n_all), z_all(n_all))
671 allocate(node1_all(n_tria_all), node2_all(n_tria_all), node3_all(n_tria_all))
672
673 call read_filexyz(file_case_all,n_all,n_tria_all,&
674 x_all,y_all,z_all,&
675 node1_all,node2_all,node3_all,&
676 max_all_spacing)
677
678 do icase = 1, ncase
679 call get_node_depth_from_alluvial(loc_n_num, n_all, n_tria_all, &
680 x_all, y_all, z_all, &
681 node1_all, node2_all, node3_all,&
682 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
683 nn_loc, xs_loc, ys_loc, zs_loc, &
684 zs_all, val_case(icase), max_all_spacing, tol_case(icase))
685 enddo
686
687 call make_subtag_alluvial(nn_loc, zs_all, j, sub_tag_all, xs_loc, ival)
688
689 deallocate(x_all, y_all, z_all, node1_all, node2_all, node3_all)
690
691 if (mpi_id.eq.0) then
692 write(*,'(A)')
693 write(*,'(A,I8)') 'ALLUVIAL Layer # ',j
694 endif
695
696 enddo !do j = 1,3
697
698 !do i = 1, nn_loc
699 ! write(*,*) zs_loc(i), sub_tag_all(i)
700 !enddo
701 !read(*,*)
702
703
704 if (mpi_id.eq.0) then
705 write(*,'(A)') 'Done'
706 write(*,'(A)')
707 endif
708
709!*************************************************************************************************
710! Groningen layered model
711!*************************************************************************************************
712
713 !elseif (tcase.eq. 32) then
714
715 ! zs_elev = -1.0e+30
716 ! zs_all = -1.0e+30
717
718! if (mpi_id.eq.0) then
719! write(*,'(A)')
720! write(*,'(A)')'CASE 32: GRONINGEN'
721! write(*,'(A)')'Layered Model...'
722! endif
723! if (mpi_id.eq.0) then
724! write(*,'(A)') 'Done'
725! write(*,'(A)')
726! endif
727
728
729!*************************************************************************************************
730! TEST honoring -- only topgraphy surface
731!*************************************************************************************************
732
733 elseif (tcase.eq.98) then
734 if (mpi_id.eq.0) then
735 write(*,'(A)')
736 write(*,'(A)')'CASE 98: TEST honoring'
737 write(*,'(A)')'Reading Topography...'
738 endif
739
740 file_case_xyz ='XYZ.out'
741 zs_elev = -1.0e+30
742
743
744 call read_dime_filexyz(file_case_xyz,n_elev,n_tria_elev)
745
746 allocate(x_elev(n_elev),y_elev(n_elev),z_elev(n_elev))
747 allocate(node1_elev(n_tria_elev))
748 allocate(node2_elev(n_tria_elev))
749 allocate(node3_elev(n_tria_elev))
750
751 call read_filexyz(file_case_xyz,n_elev,n_tria_elev,&
752 x_elev,y_elev,z_elev,&
753 node1_elev,node2_elev,node3_elev,&
754 max_elev_spacing)
755
756 do icase = 1, ncase
757
758 call get_node_depth_from_simple(loc_n_num, n_elev, n_tria_elev,&
759 x_elev, y_elev, z_elev,&
760 node1_elev, node2_elev, node3_elev,&
761 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat,&
762 nn_loc, xs_loc, ys_loc, zs_loc,&
763 zs_elev, val_case(iacse), max_elev_spacing, tol_case(icase))
764
765 enddo
766
767 deallocate(x_elev,y_elev,z_elev,node1_elev,node2_elev,node3_elev)
768
769 if (mpi_id.eq.0) then
770 write(*,'(A)')'Done'
771 write(*,'(A)')
772 endif
773
774!*************************************************************************************************
775! TEST - NOT honoring - topo and alluvial surfaces
776!*************************************************************************************************
777
778 elseif (tcase.eq.99) then
779 if (mpi_id.eq.0) then
780 write(*,'(A)')
781 write(*,'(A)')'CASE 99: TEST not honoring'
782 write(*,'(A)')'Reading Topography&Alluvial...'
783 endif
784
785 file_case_xyz ='XYZ.out'
786 file_case_all ='ALL.out'
787
788 zs_elev = -1.0e+30
789 zs_all = -1.0e+30
790
791 call read_dime_filexyz(file_case_xyz,n_elev,n_tria_elev)
792 call read_dime_filexyz(file_case_all,n_all,n_tria_all)
793
794 allocate(x_elev(n_elev),y_elev(n_elev),z_elev(n_elev))
795 allocate(node1_elev(n_tria_elev), node2_elev(n_tria_elev), node3_elev(n_tria_elev))
796
797 allocate(x_all(n_all),y_all(n_all),z_all(n_all))
798 allocate(node1_all(n_tria_all),node2_all(n_tria_all),node3_all(n_tria_all))
799
800 call read_filexyz(file_case_xyz,n_elev,n_tria_elev,&
801 x_elev,y_elev,z_elev,&
802 node1_elev,node2_elev,node3_elev,&
803 max_elev_spacing)
804
805 call read_filexyz(file_case_all,n_all,n_tria_all,&
806 x_all,y_all,z_all,&
807 node1_all,node2_all,node3_all,&
808 max_all_spacing)
809
810
811 do icase = 1, ncase
812
813 call get_node_depth_from_alluvial(loc_n_num, n_all, n_tria_all, &
814 x_all, y_all, z_all, &
815 node1_all, node2_all, node3_all,&
816 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
817 nn_loc, xs_loc, ys_loc, zs_loc, &
818 zs_all, val_case(icase), max_all_spacing, tol_case(icase))
819
820 call get_node_depth_from_cmplx(loc_n_num, n_elev, n_tria_elev, &
821 x_elev, y_elev, z_elev, &
822 node1_elev, node2_elev, node3_elev,&
823 cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, &
824 nn_loc, xs_loc, ys_loc, zs_loc, &
825 zs_elev, zs_all, &
826 val_case(icase), max_elev_spacing, tol_case(icase))
827 enddo
828
829
830 deallocate(x_elev, y_elev, z_elev, node1_elev, node2_elev, node3_elev)
831 deallocate(x_all, y_all, z_all, node1_all, node2_all, node3_all)
832
833 if (mpi_id.eq.0) then
834 write(*,'(A)')'Done'
835 write(*,'(A)')
836 endif
837
838
839 endif ! TCASE
840
841
842
843 end subroutine make_nothonoring
844
subroutine get_node_depth_and_vs(loc_n_num, nn_elev, nn_elem, xx_elev, yy_elev, zz_elev, v
Computes elevation from topography (XYZ.out), vs30 velocity and depth of sediments according to the i...
subroutine get_node_depth_from_alluvial(loc_n_num, n_elev, nn_elem
Computes elevation from alluvial basin (ALL.out).
subroutine get_node_depth_from_cmplx(loc_n_num, nn_elev, nn_elem,
Computes elevation from topography (XYZ.out).
subroutine get_node_depth_from_simple(loc_n_num, nn_elev, nn_elem,
Computes elevation from topography (XYZ.out).
subroutine make_nothonoring(loc_n_num, nn_loc, n_case, tag_case, val_case, tol_case, cs_nnz_loc, cs_loc, nm, tag_mat, sdeg_mat, xs_loc, ys_loc, zs_loc, zs_elev, zs_all, vs, thick, sub_tag_all, mpi_id)
Makes not-honoring technique.
subroutine make_subtag_alluvial(nn_s, zz_elevation, j, sub_tag_all
Assignes labels for multi not-honoring technique.
subroutine read_dime_filexyz(filec, num_nodes, num_tria)
Reads dimension of triangular grids.
subroutine read_filevs(filec, num_tria, velocity, sediments)
Reads vs30 and thickness from VS_RS.out.
subroutine read_filexyz(filec, num_nodes, num_tria, x, y, z, node1, node2, node3, max_length)
Reads mesh of triangular elements (topography or alluvial basin)