SPEED
MAKE_SEISMIC_MOMENT_OR_EXPLOSIVE_SOURCE.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
24
25 subroutine make_seismic_moment_or_explosive_source()
26
27 use max_var
28 use speed_par
29
30
31 implicit none
32
33 include 'SPEED.MPI'
34
35!*****************************************************************************************
36! SEISMIC MOMENT
37!*****************************************************************************************
38
39 if (mpi_id.eq.0 .and. nload_sism_el.gt.0) &
40 write(*,'(A)') '-----------Building the Seismic Moment vector----------'
41
42
43! Dimensioning vector 'num_node_sism'(nodes number along each faults) - begin
44
45 if (nload_sism_el.gt.0) then
46
49
50 do i = 1,nload_sism_el
51 if ((((val_sism_el(i,1).eq.val_sism_el(i,4)).and.(val_sism_el(i,4).eq.val_sism_el(i,7))) .and. &
52 (val_sism_el(i,7).eq.val_sism_el(i,10))).and.(((val_sism_el(i,2).eq.val_sism_el(i,5)).and. &
53 (val_sism_el(i,5).eq.val_sism_el(i,8))) .and. (val_sism_el(i,8).eq.val_sism_el(i,11))) &
54 .and.(((val_sism_el(i,3).eq.val_sism_el(i,6)).and.(val_sism_el(i,6).eq.val_sism_el(i,9))) &
55 .and. (val_sism_el(i,9).eq.val_sism_el(i,12)))) then
56
57 num_node_sism(i) = 1
58
59 else
60
67
68
69 endif
70 enddo !i = 1,nload_sism_el
71
73 call mpi_barrier(mpi_comm, mpi_ierr)
74 call mpi_allgather(num_node_sism, nload_sism_el, speed_integer, sism_el_glo, &
75 nload_sism_el, speed_integer, mpi_comm, mpi_ierr)
76
77
78 if(nload_sism_el .eq. 1) then
79 num_node_sism(1) = 1
80 else
81 do i = 1, nload_sism_el
82 do j = 1, mpi_np
83 vec(j) = sism_el_glo(i + (j-1)*nload_sism_el)
84 enddo
85 num_node_sism(i) = sum(vec)
86 enddo
87 endif
88
89 deallocate(vec, sism_el_glo)
90
91 !Checking the maximum number of nodes for single seismic faults
92
94
95 do i = 1,nload_sism_el
97 enddo
98
99
105
108
109 !Searching the node 'id' in the global numeration for each seismic faults.
110 !sour_node_sism = node id (global numeration) generating the 'i'th seismic fault
111 do i = 1,nload_sism_el
112 if (num_node_sism(i) .eq. 1 .and. nload_sism_el .eq. 1) then
113
114
118
119
120 !global numeration
121 !write(*,*) mpi_id, 'coord', xx_spx_loc(sour_node_sism(1,i)), yy_spx_loc(sour_node_sism(1,i)), zz_spx_loc(sour_node_sism(1,i))
122 !write(*,*) mpi_id, 'local index', sour_node_sism(1,i)
123 !write(*,*) mpi_id, 'local dist', dist_sour_node_sism(1,i)
124
128
129
130
131
133 j = 1;
134
135
138
139 call mpi_barrier(mpi_comm, mpi_ierr)
140 call mpi_allgather(sour_node_sism(1,i), 1, speed_integer, sism_el_glo, 1, speed_integer, mpi_comm, mpi_ierr)
141 call mpi_allgather(dist_sour_node_sism(1,i), 1, speed_double, dist_el_glo, 1, speed_double, mpi_comm, mpi_ierr)
142 call mpi_allgather(pos_sour_node_x(1,i), 1, speed_double, posx_el_glo, 1, speed_double, mpi_comm, mpi_ierr)
143 call mpi_allgather(pos_sour_node_y(1,i), 1, speed_double, posy_el_glo, 1, speed_double, mpi_comm, mpi_ierr)
144 call mpi_allgather(pos_sour_node_z(1,i), 1, speed_double, posz_el_glo, 1, speed_double, mpi_comm, mpi_ierr)
145
146
147
148 !if(mpi_id .eq. 0) then
149 ! write(*,*) 'global index' , sism_el_glo
150 ! write(*,*) 'global distances', dist_el_glo
151 !endif
152
153
155
156
162
163
164
165 !write(*,*) mpi_id, 'min ind', sour_node_sism(1,i), 'min dist', dist_sour_node_sism(1,i)
166
167
169
170
171
172
173 else !if(num_node_sism(i) .gt. 1) then
174
178 val_sism_el(i,10),val_sism_el(i,11),val_sism_el(i,12),&
184
188
189 sism_el_glo = 0
190 dist_el_glo = 0.d0
191 posx_el_glo = 0.d0; posy_el_glo = 0.d0; posz_el_glo = 0.d0;
192
193 call mpi_barrier(mpi_comm, mpi_ierr)
194 call mpi_allgather(sour_node_sism(:,i), max_num_node_sism, speed_integer, &
196
197 call mpi_allgather(dist_sour_node_sism(:,i), max_num_node_sism, speed_double, &
199
200 call mpi_allgather(pos_sour_node_x(:,i), max_num_node_sism, speed_double, &
202
203 call mpi_allgather(pos_sour_node_y(:,i), max_num_node_sism, speed_double, &
205
206 call mpi_allgather(pos_sour_node_z(:,i), max_num_node_sism, speed_double, &
208
209
210 j = 0
212 if(sism_el_glo(k).ne. 0) then
213
214 j = j + 1
220
221 endif
222 enddo
223 deallocate(sism_el_glo, dist_el_glo)
225
227 num_node_sism(i) = j;
229
230 endif
231
232 ! if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')
233 ! if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')'Seismic faults'
234 ! if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')
235 ! if (mpi_id.eq.0) write(*,'(A,I6,A,I6,A)')'Seismic Faults - ',i, ' is generated by ',num_node_sism(i),' nodes'
236 ! if (mpi_id.eq.0) write(*,'(I10,I12)')(j,sour_node_sism(j,i), j=1,num_node_sism(i))
237 ! if (mpi_id.eq.0) write(*,'(A)')
238 if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')
239 if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')'Seismic faults'
240 if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')
241 if (mpi_id.eq.0) write(*,'(A,I6,A,I6,A)')'Seismic Faults - ',i, ' is generated by ',j,' nodes'
242 do k = 1, j
243 if (mpi_id.eq.0) write(*,*) 'Node :', sour_node_sism(k,i)
244 if (mpi_id.eq.0) write(*,*) 'Position :', pos_sour_node_x(k,i), pos_sour_node_y(k,i), pos_sour_node_z(k,i)
245 if (mpi_id.eq.0) write(*,*) 'distance node-Hypo in [m] = ', dist_sour_node_sism(k,i)
246 enddo
247 if (mpi_id.eq.0) write(*,'(A)')
248
249
250
251 enddo !i = 1,nload_sism_el
252
253
254 endif !if (nload_sism_el.gt.0) then
255
256 ! Dimensioning vector 'num_node_sism'(nodes number along each faults) - end
257
258
259 if (nfunc.le.0) nfunc = 1
260
261 if (nload_sism_el.gt.0) then
263 allocate (tau_seismic_moment(nload_sism_el,1))
264 endif
265
266 if ((mpi_id.eq.0).and.(nload_sism_el.gt.0)) write(*,'(A)')'Seismic Moment vector built'
267 call mpi_barrier(mpi_comm, mpi_ierr)
268 !stop
269
270!*****************************************************************************************
271! EXPLOSIVE SOURCE
272!*****************************************************************************************
273
274 if (mpi_id.eq.0 .and. nload_expl_el.gt.0 ) write(*,'(A)')'----------Building the Explosive Source vector---------'
275
276 ! Dimensioning vector 'num_node_expl'(nodes number along each explosion) - begin
277 !
278 if (nload_expl_el.gt.0) then
279
280 allocate (num_node_expl(nload_expl_el))
281
282 do i = 1,nload_expl_el
283 if ((((val_expl_el(i,1).eq.val_expl_el(i,4)).and.(val_expl_el(i,4).eq.val_expl_el(i,7))) &
284 .and. (val_expl_el(i,7).eq.val_expl_el(i,10))).and.(((val_expl_el(i,2).eq.val_expl_el(i,5)) &
285 .and.(val_expl_el(i,5).eq.val_expl_el(i,8))) .and. (val_expl_el(i,8).eq.val_expl_el(i,11))) &
286 .and.(((val_expl_el(i,3).eq.val_expl_el(i,6)).and.(val_expl_el(i,6).eq.val_expl_el(i,9))) &
287 .and. (val_expl_el(i,9).eq.val_expl_el(i,12)))) then
288
290
291 else
295 val_expl_el(i,10),val_expl_el(i,11),val_expl_el(i,12),&
298 endif
299 enddo !i = 1,nload_expl_el
300
301
303 call mpi_barrier(mpi_comm, mpi_ierr)
304 call mpi_allgather(num_node_expl, nload_expl_el, speed_integer, expl_el_glo, &
305 nload_expl_el, speed_integer, mpi_comm, mpi_ierr)
306
307 if(nload_expl_el .eq. 1) then
308 num_node_expl(1) = 1
309 else
310 do i = 1, nload_expl_el
311 do j = 1, mpi_np
312 vec(j) = expl_el_glo(i + (j-1)*nload_expl_el)
313 enddo
314
315 num_node_expl(i) = sum(vec)
316 enddo
317 endif
318
319
320 deallocate(vec, expl_el_glo)
321
322 !Checking the maximum number of nodes for single seismic faults
323
325
326 do i = 1,nload_expl_el
327 if (num_node_expl(i).gt.max_num_node_expl) then
329 endif
330 enddo
331
336
337 !Searching the node 'id' in the global numeration for each Explosive Source.
338 !sour_node_expl = node id (global numeration) generating the 'i'th Explosion
339 do i = 1,nload_expl_el
340 if (num_node_expl(i).eq.1 .and. nload_expl_el .eq. 1 ) then
341
345
346
348
350 call mpi_barrier(mpi_comm, mpi_ierr)
351 call mpi_allgather(sour_node_expl(1,i), 1, speed_integer, expl_el_glo, 1, speed_integer, mpi_comm, mpi_ierr)
352 call mpi_allgather(dist_sour_node_expl(1,i), 1, speed_double, dist_el_glo, 1, speed_double, mpi_comm, mpi_ierr)
353
356
357 deallocate(expl_el_glo, dist_el_glo)
358
359
360 else
364 val_expl_el(i,10),val_expl_el(i,11),val_expl_el(i,12),&
369
370
372 call mpi_barrier(mpi_comm, mpi_ierr)
373 call mpi_allgather(sour_node_expl(:,i), max_num_node_expl, speed_integer, &
375 call mpi_allgather(dist_sour_node_expl(:,i), max_num_node_expl, speed_double, &
377 j = 0
379 if(expl_el_glo(k).ne. 0) then
380
381 j = j + 1
384
385 endif
386 enddo
387 deallocate(expl_el_glo, dist_el_glo)
388
389 endif
390 if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')
391 if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')'Explosive Source'
392 if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')
393 if (mpi_id.eq.0) write(*,'(A,I6,A,I6,A)')'Explosive Source - ',i, ' is generated by ',num_node_expl(i),' nodes'
394 if (mpi_id.eq.0) write(*,'(I6,I6)')(j,sour_node_expl(j,i), j=1,num_node_expl(i))
395 if (mpi_id.eq.0 .and. i .eq. 1) write(*,'(A)')
396 enddo !i = 1,nload_expl_el
397
398
399 endif !if (nload_expl_el.gt.0) then
400
401 ! Dimensioning vector 'num_node_expl'(nodes number along each explosion area) - end
402 if (nfunc.le.0) nfunc = 1
403
405 if ((mpi_id.eq.0).and.(nload_expl_el.gt.0)) write(*,'(A)')'Explosive Source vector built'
406
407
408
409 end subroutine make_seismic_moment_or_explosive_source
410
subroutine get_dime_expl(xipo, yipo, zipo, x1, y1, z1, x2, y2, z2, x3, y3, z3, nnod, xs, ys, zs, node_expl, nnloc)
Computes local number of nodes where explosive source is imposed.
subroutine get_dime_sism(xipo, yipo, zipo, x1, y1, z1, x2, y2, z2, x3, y3, z3, nnod, xs, ys, zs, node_sism, nn_loc, mpi_id, ind, loc_n_n
Computes local number of nodes where seismic source is imposed.
subroutine get_minvalues(i_glo, v_glo, n_glo, i_loc, n_loc, np)
Computes positions of minimum values.
subroutine get_nearest_node(n, xs, ys, zs, xt, yt, zt, nt, dist_min)
Computes the nearest node with respet to (xt,yt,zt).
subroutine read_expl(xipo, yipo, zipo, x1, y1, z1, x2, y2, z2, x3, y3, z3, nnod, xs, ys, zs, num_ne, sour_ne, i, dist_sour_ne, nl_expl, max_num_ne, loc_n_num, nn_loc)
Generates explosive triangular faults.
Definition READ_EXPL.f90:56
subroutine read_sism(xipo, yipo, zipo, x1, y1, z1, x2, y2, z2, x3, y3, z3, nnod, xs, ys, zs, num_ns, sour_ns, i_sism, dist_sour_ns, nl_sism, max_num_ns, loc_n_num, nn_loc, pos_sour_nx, pos_sour_ny, pos_sour_nz)
Generates seismic triangular faults.
Definition READ_SISM.f90:57
Set maximal bounds.
Definition MODULES.f90:54
Contains SPEED PARAMETERS used in (SPEED, READ_INPUT_FILES, MAKE_PARTION_AND_MPI_FILES,...
Definition MODULES.f90:207
integer *4 nload_sism_el
Definition MODULES.f90:269
real *8, dimension(:,:), allocatable pos_sour_node_y
Definition MODULES.f90:474
real *8, dimension(:,:), allocatable pos_sour_node_x
Definition MODULES.f90:474
real *8, dimension(:,:), allocatable factor_seismic_moment
Definition MODULES.f90:474
integer *4, dimension(:,:), allocatable sour_node_expl
Definition MODULES.f90:383
integer *4 nfunc
Definition MODULES.f90:269
integer *4 j
Definition MODULES.f90:263
integer *4, dimension(:), allocatable sism_el_glo
Definition MODULES.f90:363
integer *4, dimension(:), allocatable num_node_expl
Definition MODULES.f90:363
real *8, dimension(:,:), allocatable factor_explosive_source
Definition MODULES.f90:474
real *8, dimension(:,:), allocatable val_sism_el
Definition MODULES.f90:463
real *8, dimension(:), allocatable zz_spx_loc
Definition MODULES.f90:408
integer *4 i
Definition MODULES.f90:263
real *8, dimension(:,:), allocatable dist_sour_node_expl
Definition MODULES.f90:474
integer *4 nload_expl_el
Definition MODULES.f90:269
real *8, dimension(:), allocatable yy_spx_loc
Definition MODULES.f90:408
real *8, dimension(:,:), allocatable dist_sour_node_sism
Definition MODULES.f90:474
integer *4, dimension(:), allocatable vec
Definition MODULES.f90:375
integer *4 k
Definition MODULES.f90:263
integer *4 mpi_comm
Definition MODULES.f90:308
integer *4, dimension(:), allocatable local_node_num
Definition MODULES.f90:322
integer *4 mpi_ierr
Definition MODULES.f90:298
integer *4, dimension(:), allocatable expl_el_glo
Definition MODULES.f90:363
real *8, dimension(:), allocatable posy_el_glo
Definition MODULES.f90:417
integer *4, dimension(:), allocatable num_node_sism
Definition MODULES.f90:363
integer *4 max_num_node_expl
Definition MODULES.f90:269
integer *4 nnode_dom
Definition MODULES.f90:269
integer *4 max_num_node_sism
Definition MODULES.f90:269
real *8, dimension(:), allocatable posx_el_glo
Definition MODULES.f90:417
integer *4 mpi_np
Definition MODULES.f90:308
real *8, dimension(:), allocatable posz_el_glo
Definition MODULES.f90:417
real *8, dimension(:,:), allocatable val_expl_el
Definition MODULES.f90:463
integer *4 mpi_id
Definition MODULES.f90:308
integer *4 nnod_loc
Definition MODULES.f90:269
real *8, dimension(:), allocatable xx_spx_loc
Definition MODULES.f90:408
real *8, dimension(:,:), allocatable pos_sour_node_z
Definition MODULES.f90:474
integer *4, dimension(:,:), allocatable sour_node_sism
Definition MODULES.f90:383
real *8, dimension(:,:), allocatable tau_seismic_moment
Definition MODULES.f90:474
real *8, dimension(:), allocatable dist_el_glo
Definition MODULES.f90:417