SPEED
GET_NODE_DEPTH_FROM_SIMPLE.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine get_node_depth_from_simple (loc_n_num, nn_elev, nn_elem,
 Computes elevation from topography (XYZ.out).
 

Function/Subroutine Documentation

◆ get_node_depth_from_simple()

subroutine get_node_depth_from_simple ( integer*4, dimension(nn_s)  loc_n_num,
integer*4  nn_elev,
integer*4  nn_elem 
)

Computes elevation from topography (XYZ.out).

Note
See MAKE_NOT_HONORING.f90
Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]nn_snumber of local nodes
[in]loc_n_numlocal numeration vector
[in]nn_elevnumber nodes in the triangular grid
[in]x_elevelevation values of local nodes
[in]y_elevelevation values of local nodes
[in]z_elevelevation values of local nodes
[in]nn_elemnumber of triangular elements
[in]node1_elemindex triangle vertex
[in]node2_elemindex triangle vertex
[in]node3_elemindex triangle vertex
[in]cs_nnx_loclength cs_loc
[in]cs_loclocal connectivity vector
[in]xx_svertex x- coordinate of local nodes
[in]yy_svertex y- coordinate of local nodes
[in]zz_svertex z- coordinate of local nodes
[in]nmnumber of materials
[in]tmlabels for material vector
[in]sdpolynomial degree vector
[in]tagmatspecific material tag given in CASE option
[in]max_esmax topography spacing
[in]toltolerance given in CASE option
[out]zz_elevationelevation of the nodes from simple topography

Definition at line 47 of file GET_NODE_DEPTH_FROM_SIMPLE.f90.

48 xx_elev, yy_elev, zz_elev, &
49 node1_elem, node2_elem, node3_elem, &
50 cs_nnz_loc, cs_loc, nm, tm, sd, &
51 nn_s, xx_s, yy_s, zz_s, zz_elevation, &
52 tagmat, max_es,tol)
53
54
55
56 implicit none
57
58 integer*4 :: nn_elev,nn_elem,cs_nnz_loc,nm,ne,nn_s
59 integer*4 :: im,ie,i,j,k,nn,ip,isn,ic
60 integer*4 :: h
61 integer*4 :: tagmat
62
63 integer*4, dimension(nn_s) :: loc_n_num
64 integer*4, dimension(nn_elem) :: node1_elem,node2_elem,node3_elem
65 integer*4, dimension(0:cs_nnz_loc) :: cs_loc
66 integer*4, dimension(nm) :: tm
67 integer*4, dimension(nm) :: sd
68
69 real*8 :: dx,dy,dz,tol
70 real*8 :: x1,y1,z1
71 real*8 :: x2,y2,z2
72 real*8 :: x3,y3,z3
73 real*8 :: ux,uy,uz,vx,vy,vz
74 real*8 :: a,b,c
75 real*8 :: max_es
76 real*8 :: zz_interp
77 real*8 :: v0x,v0y,v1x,v1y,v2x,v2y
78 real*8 :: dot00,dot01,dot02,dot11,dot12
79 real*8 :: invdenom,u,v
80 real*8 :: d2min
81 real*8 :: zz_elev_min
82
83 real*8, dimension(:), allocatable :: ct,ww
84 real*8, dimension(nn_elev) :: xx_elev,yy_elev,zz_elev
85 real*8, dimension(nn_s) :: xx_s, yy_s, zz_s
86 real*8, dimension(nn_s) :: zz_elevation
87
88 real*8, dimension(:,:), allocatable :: dd
89
90 d2min = (5 * max_es)**2
91
92 zz_elev_min = zz_elev(1)
93 do i = 1,nn_elev
94 if (zz_elev(i).lt.zz_elev_min) then
95 zz_elev_min = zz_elev(i)
96 endif
97 enddo
98
99
100 nn = 2
101 allocate(ct(nn),ww(nn),dd(nn,nn))
102 call make_lgl_nw(nn,ct,ww,dd)
103
104 ne = cs_loc(0) - 1
105
106 do im = 1,nm
107 if ((sd(im) +1).ne.nn) then
108 deallocate(ct,ww,dd)
109 nn = sd(im) +1
110 allocate(ct(nn),ww(nn),dd(nn,nn))
111 call make_lgl_nw(nn,ct,ww,dd)
112 endif
113
114 do ie = 1,ne
115 if (cs_loc(cs_loc(ie -1) +0).eq.tagmat) then
116 do k = 1,nn
117 do j = 1,nn
118 do i = 1,nn
119
120 ip = nn*nn*(k -1) +nn*(j -1) +i
121 isn = cs_loc(cs_loc(ie -1) + ip)
122 ic = isn
123
124
125
126 if (zz_elevation(ic) .eq. -1.0e+30) then
127
128 do h = 1, nn_elem
129
130 x1 = xx_elev(node1_elem(h))
131 y1 = yy_elev(node1_elem(h))
132 z1 = zz_elev(node1_elem(h))
133
134 if (((x1 - xx_s(ic))**2 + (y1 - yy_s(ic))**2).le.d2min) then
135 x2 = xx_elev(node2_elem(h))
136 y2 = yy_elev(node2_elem(h))
137 z2 = zz_elev(node2_elem(h))
138
139 x3 = xx_elev(node3_elem(h))
140 y3 = yy_elev(node3_elem(h))
141 z3 = zz_elev(node3_elem(h))
142
143 !Point in triangle test
144 ! P = (xx_s(isn) yy_s(isn))
145 ! A = (X1 Y1)
146 ! B = (X2 Y2)
147 ! C = (X3 Y3)
148 ! Compute vectors
149 ! v0 = C - A
150 v0x=(x3 - x1)
151 v0y=(y3 - y1)
152 ! v1 = B - A
153 v1x=(x2 - x1)
154 v1y=(y2 - y1)
155 ! v2 = P - A
156 v2x=(xx_s(ic) - x1)
157 v2y=(yy_s(ic) - y1)
158
159 ! Compute dot products
160 ! [u].[v] = ux * vx + uy * vy
161 ! dot([u],[v])
162 !dot00 = dot(v0, v0)
163 dot00 = v0x * v0x + v0y * v0y
164 !dot01 = dot(v0, v1)
165 dot01 = v0x * v1x + v0y * v1y
166 !dot02 = dot(v0, v2)
167 dot02 = v0x * v2x + v0y * v2y
168 !dot11 = dot(v1, v1)
169 dot11 = v1x * v1x + v1y * v1y
170 !dot12 = dot(v1, v2)
171 dot12 = v1x * v2x + v1y * v2y
172
173 ! Compute barycentric coordinates
174 invdenom = 1 / (dot00 * dot11 - dot01 * dot01)
175 u = (dot11 * dot02 - dot01 * dot12) * invdenom
176 v = (dot00 * dot12 - dot01 * dot02) * invdenom
177 !Point in triangle test
178 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
179
180 if ( (u.ge.0.0d0).and.(v.ge.0.0d0).and.((u + v).le.1.0d0) ) then
181
182 ! Build up the plane passing through the points P1, P2 and P3
183
184 ux=(x1-x2)/sqrt((x1-x2)**2+(y1-y2)**2+(z1-z2)**2)
185 uy=(y1-y2)/sqrt((x1-x2)**2+(y1-y2)**2+(z1-z2)**2)
186 uz=(z1-z2)/sqrt((x1-x2)**2+(y1-y2)**2+(z1-z2)**2)
187 vx=(x3-x2)/sqrt((x3-x2)**2+(y3-y2)**2+(z3-z2)**2)
188 vy=(y3-y2)/sqrt((x3-x2)**2+(y3-y2)**2+(z3-z2)**2)
189 vz=(z3-z2)/sqrt((x3-x2)**2+(y3-y2)**2+(z3-z2)**2)
190
191 a = uy * vz - uz * vy
192 b = uz * vx - ux * vz
193 c = ux * vy - uy * vx
194
195 zz_interp = -a/c * (xx_s(ic)-x1) -b/c * (yy_s(ic)-y1) + z1
196 zz_elevation(ic) = ( zz_interp - zz_s(ic) )
197
198 if (abs(zz_elevation(ic)).lt.tol) then
199 zz_elevation(ic) = 0.0d0
200 endif
201
202 endif !if ( (u.ge.0.0d0).and.(v.ge.0.0d0).and.((u + v).le.1.0d0) ) then
203
204 if ( (u.ge.0.0d0).and.(v.ge.0.0d0).and.((u + v).le.1.0d0) ) exit
205
206 endif !if (((X1 - xx_s(isn))**2 + (Y1 - yy_s(isn))**2).le.d2min) then
207
208 enddo !do h = 1,nn_elem
209
210 endif !if (zz_elevation(isn).le.0.0d0) then
211
212
213
214 enddo !do k = 1,nn
215 enddo !do j = 1,nn
216 enddo !do i = 1,nn
217
218
219 endif
220 enddo
221 enddo
222
223 return
224
subroutine make_lgl_nw(nb_pnt, xq, wq, dd)
Makes Gauss-Legendre-Lobatto nodes, weigths and spectral derivatives.

References make_lgl_nw().

Referenced by make_nothonoring().

Here is the call graph for this function:
Here is the caller graph for this function: