SPEED
READ_EXPL.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

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.
 

Function/Subroutine Documentation

◆ read_expl()

subroutine read_expl ( 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  num_ne,
integer*4, dimension(max_num_ne,nl_expl)  sour_ne,
integer*4  i,
real*8, dimension(max_num_ne,nl_expl)  dist_sour_ne,
integer*4  nl_expl,
integer*4  max_num_ne,
integer*4, dimension(nn_loc)  loc_n_num,
integer*4  nn_loc 
)

Generates explosive triangular faults.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]XipoHipocenter x-coordinate
[in]YipoHipocenter y-coordinate
[in]ZipoHipocenter z-coordinate
[in]X1x-coordinate of the triangular fault 1-node
[in]Y1y-coordinate of the triangular fault 1-node
[in]Z1z-coordinate of the triangular fault 1-node
[in]X2x-coordinate of the triangular fault 2-node
[in]Y2y-coordinate of the triangular fault 2-node
[in]Z2z-coordinate of the triangular fault 2-node
[in]X3x-coordinate of the triangular fault 3-node
[in]Y3y-coordinate of the triangular fault 3-node
[in]Z3z-coordinate of the triangular fault 3-node
[in]xsx-coord local spectral nodes
[in]ysy-coord local spectral nodes
[in]zsz-coord local spectral nodes
[in]num_nenumber of explosive nodes
[in]icolumn index for filling sour_ne
[in]nl_explnumber of explosive loads
[in]nn_locnumber of local nodes
[in]nnodnumber of local nodes
[in]loc_n_numlocal node numeration
[in]max_num_nemax number of nodes for the specific triangular fault
[out]sour_neinfo about explosive nodes
[out]dist_sour_nedistance of explosive nodes in sour_ne

Definition at line 48 of file READ_EXPL.f90.

56
57
58 implicit none
59
60 integer*4 :: nnod,node_expl,num_ne,nl_expl,max_num_ne, nn_loc
61 integer*4 :: i,isn
62
63 integer*4, dimension(nn_loc) :: loc_n_num(nn_loc)
64
65 integer*4, dimension(max_num_ne,nl_expl) :: sour_ne
66
67 real*8 :: xipo,yipo,zipo,x1,y1,z1,x2,y2,z2,x3,y3,z3,ux,uy,uz,vx,vy,vz,tol
68 real*8 :: xmax,xmin,ymax,ymin,zmax,zmin,a,b,c
69 real*8 :: den
70
71 real*8, dimension(nn_loc) :: xs,ys,zs
72
73 real*8, dimension(max_num_ne,nl_expl) :: dist_sour_ne
74
75 node_expl = 0
76 tol = 5.0d0
77
78 ux=(x1-x2)/sqrt((x1-x2)**2+(y1-y2)**2+(z1-z2)**2)
79 uy=(y1-y2)/sqrt((x1-x2)**2+(y1-y2)**2+(z1-z2)**2)
80 uz=(z1-z2)/sqrt((x1-x2)**2+(y1-y2)**2+(z1-z2)**2)
81 vx=(x3-x2)/sqrt((x3-x2)**2+(y3-y2)**2+(z3-z2)**2)
82 vy=(y3-y2)/sqrt((x3-x2)**2+(y3-y2)**2+(z3-z2)**2)
83 vz=(z3-z2)/sqrt((x3-x2)**2+(y3-y2)**2+(z3-z2)**2)
84
85 a = uy * vz - uz * vy
86 b = uz * vx - ux * vz
87 c = ux * vy - uy * vx
88 den = sqrt(a**2 + b**2 + c**2)
89
90 !Check on the relative position of fault end points
91 !In this preliminary part we check if P1 point is the lower-left fault point or not
92 !If this check would fail the two points are swapped
93
94
95 xmax=max(x1,x2,x3)+tol
96 xmin=min(x1,x2,x3)-tol
97 ymax=max(y1,y2,y3)+tol
98 ymin=min(y1,y2,y3)-tol
99 zmax=max(z1,z2,z3)+tol
100 zmin=min(z1,z2,z3)-tol
101
102
103 do isn = 1,nnod
104
105 if ( dabs( a*(xs(isn)-x1) + b*(ys(isn)-y1) + c*(zs(isn)-z1)/den ).le.tol) then
106 if ((xs(isn).ge.xmin).and.(xs(isn).le.xmax)) then
107 if ((ys(isn).ge.ymin).and.(ys(isn).le.ymax)) then
108 if ((zs(isn).ge.zmin).and.(zs(isn).le.zmax)) then
109
110 node_expl = node_expl + 1
111 sour_ne(node_expl,i) = loc_n_num(isn)
112 dist_sour_ne(node_expl,i) = sqrt((xipo - xs(isn))**2 + (yipo - ys(isn))**2 + (zipo - zs(isn))**2)
113
114 endif
115 endif
116 endif
117 endif
118
119 enddo
120
121
122 return