SPEED
GET_DIME_SISM.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

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.
 

Function/Subroutine Documentation

◆ get_dime_sism()

subroutine get_dime_sism ( real*8  xipo,
real*8  yipo,
real*8  zipo,
real*8  x1,
real*8  y1,
real*8  z1,
real*8  x2,
real*8  y2,
real*8  z2,
real*8  x3,
real*8  y3,
real*8  z3,
integer*4  nnod,
real*8, dimension(nn_loc)  xs,
real*8, dimension(nn_loc)  ys,
real*8, dimension(nn_loc)  zs,
integer*4  node_sism,
integer*4  nn_loc,
integer*4  mpi_id,
integer*4  ind,
  loc_n_n 
)

Computes local number of nodes where seismic source is imposed.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]Xipox-coordinate of Hypocenter
[in]Yipoy-coordinate of Hypocenter
[in]Zipoz-coordinate of Hypocenter
[in]X1vertex coordinate (x) of the triangular fault
[in]Y1vertex coordinate (y) of the triangular fault
[in]Z1vertex coordinate (z) of the triangular fault
[in]X2vertex coordinate (x) of the triangular fault
[in]Y2vertex coordinate (y) of the triangular fault
[in]Z2vertex coordinate (z) of the triangular fault
[in]X3vertex coordinate (x) of the triangular fault
[in]Y3vertex coordinate (y) of the triangular fault
[in]Z3vertex coordinate (z) of the triangular fault
[in]nnodnumber of local nodes (according to nodedomain.mpi)
[in]xsx-coordinate of spectral nodes
[in]ysy-coordinate of spectral nodes
[in]zsz-coordinate of spectral nodes
[in]nn_locTotal number of local nodes
[in]loc_n_numnumeration of local nodes
[in]mpi_idid of MPI process
[in]indindex dummy
[out]node_sismlocal number of seismic source

Definition at line 47 of file GET_DIME_SISM.f90.

53
54 implicit none
55
56 integer*4 :: isn,i,nn_loc,mpi_id,ind
57 integer*4 :: nnod,node_sism
58
59 integer*4, dimension(nn_loc) :: loc_n_num
60
61 real*8 :: xipo,yipo,zipo,x1,y1,z1,x2,y2,z2,x3,y3,z3,x4,y4,z4,ux,uy,uz,vx,vy,vz,tol
62 real*8 :: p1x,p1y,p1z,p2x,p2y,p2z
63 real*8 :: epsilon, twopi, anglesum, costheta,m1,m2
64
65 real*8, dimension(4) :: x,y,z
66 real*8, dimension(nn_loc) :: xs,ys,zs
67
68 node_sism = 0
69
70 tol = 6.28
71 epsilon = 1.d-8
72 twopi = 6.283185307179586476925287
73
74
75 x(1) = x1
76 y(1) = y1
77 z(1) = z1
78
79 x(2) = x2
80 y(2) = y2
81 z(2) = z2
82
83 x(3) = x3
84 y(3) = y3
85 z(3) = z3
86
87 x(4) = x1
88 y(4) = y1
89 z(4) = z1
90
91 do isn = 1,nnod
92 anglesum = 0.0d0
93 do i = 1,3
94 p1x = x(i) - xs(isn)
95 p1y = y(i) - ys(isn)
96 p1z = z(i) - zs(isn)
97 p2x = x(i+1) - xs(isn)
98 p2y = y(i+1) - ys(isn)
99 p2z = z(i+1) - zs(isn)
100
101 m1 = dsqrt(p1x*p1x + p1y*p1y + p1z*p1z)
102 m2 = dsqrt(p2x*p2x + p2y*p2y + p2z*p2z)
103
104 if ((m1*m2).le.epsilon) then
105 anglesum = twopi
106 else
107 costheta = (p1x*p2x + p1y*p2y + p1z*p2z)/(m1*m2)
108 anglesum = anglesum + dacos(costheta)
109 endif
110
111
112 enddo
113
114
115 if (anglesum.ge.tol) then
116 node_sism = node_sism + 1
117 endif
118
119 enddo
120
121
122 return