FE-Project
Loading...
Searching...
No Matches
scale_element_hexahedral.F90
Go to the documentation of this file.
1
9#include "scaleFElib.h"
11
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17
18 use scale_element_base, only: &
21 !-----------------------------------------------------------------------------
22 implicit none
23 private
24
25 !-----------------------------------------------------------------------------
26 !
27 !++ Public type & procedure
28 !
29 type, public, extends(elementbase3d) :: hexahedralelement
30 contains
31 procedure :: init => hexhedralelement_init
32 procedure :: final => hexhedralelement_final
33 procedure :: genintgausslegendreintrpmat => hexhedralelement_gen_intgausslegendreintrpmat
34 end type hexahedralelement
35
36contains
37!OCL SERIAL
39 elem, elemOrder_h, elemOrder_v, &
40 LumpedMassMatFlag )
41 implicit none
42
43 class(hexahedralelement), intent(inout) :: elem
44 integer, intent(in) :: elemOrder_h
45 integer, intent(in) :: elemOrder_v
46 logical, intent(in) :: LumpedMassMatFlag
47
48 !-----------------------------------------------------------------------------
49
50 elem%PolyOrder_h = elemorder_h
51 elem%PolyOrder_v = elemorder_v
52 elem%Nnode_h1D = elemorder_h + 1
53 elem%Nnode_v = elemorder_v + 1
54
55 elem%Nv = 8
56 elem%Nfaces_h = 4
57 elem%Nfaces_v = 2
58 elem%Nfaces = elem%Nfaces_h + elem%Nfaces_v
59
60 elem%Nfp_h = elem%Nnode_h1D*elem%Nnode_v
61 elem%Nfp_v = elem%Nnode_h1D**2
62 elem%NfpTot = elem%Nfp_h*elem%Nfaces_h + elem%Nfp_v*elem%Nfaces_v
63
64 elem%Np = elem%Nfp_v * elem%Nnode_v
65
66 call elementbase3d_init(elem, lumpedmassmatflag)
67 call construct_element(elem)
68
69 return
70 end subroutine hexhedralelement_init
71
72!OCL SERIAL
73 subroutine hexhedralelement_final(elem)
74 implicit none
75
76 class(hexahedralelement), intent(inout) :: elem
77 !-----------------------------------------------------------------------------
78
79 call elementbase3d_final(elem)
80
81 return
82 end subroutine hexhedralelement_final
83
84!OCL SERIAL
85 subroutine construct_element(elem)
86
88 use scale_polynominal, only: &
92 use scale_element_base, only: &
96
97 implicit none
98
99 type(hexahedralelement), intent(inout) :: elem
100
101 integer :: nodes_ijk(elem%Nnode_h1D, elem%Nnode_h1D, elem%Nnode_v)
102
103 real(RP) :: lglPts1D_h(elem%Nnode_h1D)
104 real(RP) :: lglPts1D_v(elem%Nnode_v)
105
106 real(DP) :: intWeight_lgl1DPts_h(elem%Nnode_h1D)
107 real(DP) :: intWeight_lgl1DPts_v(elem%Nnode_v)
108
109 real(RP) :: P1D_ori_h(elem%Nnode_h1D, elem%Nnode_h1D)
110 real(RP) :: P1D_ori_v(elem%Nnode_v, elem%Nnode_v)
111 real(RP) :: DP1D_ori_h(elem%Nnode_h1D, elem%Nnode_h1D)
112 real(RP) :: DP1D_ori_v(elem%Nnode_v, elem%Nnode_v)
113 real(RP) :: DLagr1D_h(elem%Nnode_h1D, elem%Nnode_h1D)
114 real(RP) :: DLagr1D_v(elem%Nnode_v, elem%Nnode_v)
115 real(RP) :: V2D_h(elem%Nfp_h, elem%Nfp_h)
116 real(RP) :: V2D_v(elem%Nfp_v, elem%Nfp_v)
117 real(RP) :: Emat(elem%Np, elem%NfpTot)
118 real(RP) :: MassEdge_h(elem%Nfp_h, elem%Nfp_h)
119 real(RP) :: MassEdge_v(elem%Nfp_v, elem%Nfp_v)
120
121 integer :: i, j, k
122 integer :: p1, p2, p3
123 integer :: n, l, f
124 integer :: Nord
125 integer :: is, ie
126
127 integer :: f_h, f_v
128 integer :: fp, fp_h1, fp_h2, fp_v
129
130 type(quadrilateralelement) :: elem2D
131 !-----------------------------------------------------------------------------
132
133 lglpts1d_h(:) = polynominal_gengausslobattopt( elem%PolyOrder_h )
134 p1d_ori_h(:,:) = polynominal_genlegendrepoly( elem%PolyOrder_h, lglpts1d_h )
135 dp1d_ori_h(:,:) = polynominal_gendlegendrepoly( elem%PolyOrder_h, lglpts1d_h, p1d_ori_h )
136 dlagr1d_h(:,:) = polynominal_gendlagrangepoly_lglpt( elem%PolyOrder_h, lglpts1d_h )
137
138 lglpts1d_v(:) = polynominal_gengausslobattopt( elem%PolyOrder_v )
139 p1d_ori_v(:,:) = polynominal_genlegendrepoly( elem%PolyOrder_v, lglpts1d_v )
140 dp1d_ori_v(:,:) = polynominal_gendlegendrepoly( elem%PolyOrder_v, lglpts1d_v, p1d_ori_v )
141 dlagr1d_v(:,:) = polynominal_gendlagrangepoly_lglpt( elem%PolyOrder_v, lglpts1d_v )
142
143 !* Preparation
144
145 do k=1, elem%Nnode_v
146 do j=1, elem%Nnode_h1D
147 do i=1, elem%Nnode_h1D
148 nodes_ijk(i,j,k) = i + (j-1)*elem%Nnode_h1D + (k-1)*elem%Nnode_h1D**2
149 end do
150 end do
151 end do
152
153 ! Set the mask to extract the values at faces
154
155 elem%Fmask_h(:,1) = reshape(nodes_ijk(:,1,:), (/ elem%Nfp_h /))
156 elem%Fmask_h(:,2) = reshape(nodes_ijk(elem%Nnode_h1D,:,:), (/ elem%Nfp_h /))
157 elem%Fmask_h(:,3) = reshape(nodes_ijk(:,elem%Nnode_h1D,:), (/ elem%Nfp_h /))
158 elem%Fmask_h(:,4) = reshape(nodes_ijk(1,:,:), (/ elem%Nfp_h /))
159
160 elem%Fmask_v(:,1) = reshape(nodes_ijk(:,:,1), (/ elem%Nfp_v /))
161 elem%Fmask_v(:,2) = reshape(nodes_ijk(:,:,elem%Nnode_v), (/ elem%Nfp_v /))
162
163 !- ColMask
164
165 do j=1, elem%Nnode_h1D
166 do i=1, elem%Nnode_h1D
167 n = i + (j-1)*elem%Nnode_h1D
168 elem%Colmask(:,n) = nodes_ijk(i,j,:)
169 end do
170 end do
171
172 != Hslice
173
174 do k=1, elem%Nnode_v
175 elem%Hslice(:,k) = reshape(nodes_ijk(:,:,k), (/ elem%Nfp_v /))
176 end do
177
178 !- IndexH2Dto3D
179
180 do k=1, elem%Nnode_v
181 do j=1, elem%Nnode_h1D
182 do i=1, elem%Nnode_h1D
183 n = i + (j-1)*elem%Nnode_h1D + (k-1)*elem%Nnode_h1D**2
184 elem%IndexH2Dto3D(n) = nodes_ijk(i,j,1)
185 end do
186 end do
187 end do
188
189 !- IndexH2Dto3D_bnd
190
191 call elem2d%Init( elem%PolyOrder_h, .false. )
192
193 do f_h=1, 4
194 do fp_v=1, elem%Nnode_v
195 do fp_h1=1, elem%Nnode_h1D
196 fp = fp_h1 + (fp_v-1)*elem%Nnode_h1D + (f_h-1)*elem%Nfp_h
197 elem%IndexH2Dto3D_bnd(fp) = elem2d%Fmask(fp_h1,f_h)
198 end do
199 end do
200 end do
201 do f_v=1, 2
202 do fp_h2=1, elem%Nnode_h1D
203 do fp_h1=1, elem%Nnode_h1D
204 fp = fp_h1 + (fp_h2-1)*elem%Nnode_h1D &
205 + (f_v-1) * elem%Nfp_v &
206 + 4 * elem%Nnode_h1D * elem%Nnode_v
207 elem%IndexH2Dto3D_bnd(fp) = fp_h1 + (fp_h2-1)*elem%Nnode_h1D
208 end do
209 end do
210 end do
211
212 call elem2d%Final()
213
214 !- IndexZ1Dto3D
215
216 do k=1, elem%Nnode_v
217 do j=1, elem%Nnode_h1D
218 do i=1, elem%Nnode_h1D
219 n = i + (j-1)*elem%Nnode_h1D + (k-1)*elem%Nnode_h1D**2
220 elem%IndexZ1Dto3D(n) = k
221 end do
222 end do
223 end do
224
225 !* Set the coordinates of LGL points, and the Vandermonde and differential matricies
226
227 elem%Dx1(:,:) = 0.0_rp
228 elem%Dx2(:,:) = 0.0_rp
229 elem%Dx3(:,:) = 0.0_rp
230
231 do k=1, elem%Nnode_v
232 do j=1, elem%Nnode_h1D
233 do i=1, elem%Nnode_h1D
234 n = i + (j-1)*elem%Nnode_h1D + (k-1)*elem%Nnode_h1D**2
235
236 !* Set the coordinates of LGL points
237 elem%x1(n) = lglpts1d_h(i)
238 elem%x2(n) = lglpts1d_h(j)
239 elem%x3(n) = lglpts1d_v(k)
240
241 !* Set the Vandermonde and differential matricies
242 do p3=1, elem%Nnode_v
243 do p2=1, elem%Nnode_h1D
244 do p1=1, elem%Nnode_h1D
245 l = p1 + (p2-1)*elem%Nnode_h1D + (p3-1)*elem%Nnode_h1D**2
246 elem%V(n,l) = (p1d_ori_h(i,p1)*p1d_ori_h(j,p2)*p1d_ori_v(k,p3)) &
247 * sqrt((dble(p1-1) + 0.5_dp)*(dble(p2-1) + 0.5_dp)*(dble(p3-1) + 0.5_dp))
248
249 if(p2==j .and. p3==k) elem%Dx1(n,l) = dlagr1d_h(p1,i)
250 if(p1==i .and. p3==k) elem%Dx2(n,l) = dlagr1d_h(p2,j)
251 if(p1==i .and. p2==j) elem%Dx3(n,l) = dlagr1d_v(p3,k)
252 end do
253 end do
254 end do
255 end do
256 end do
257 end do
258 elem%invV(:,:) = linalgebra_inv(elem%V)
259
260 !* Set the weights at LGL points to integrate over element
261
262 intweight_lgl1dpts_h(:) = polynominal_gengausslobattoptintweight(elem%PolyOrder_h)
263 intweight_lgl1dpts_v(:) = polynominal_gengausslobattoptintweight(elem%PolyOrder_v)
264
265 do k=1, elem%Nnode_v
266 do j=1, elem%Nnode_h1D
267 do i=1, elem%Nnode_h1D
268 l = i + (j-1)*elem%Nnode_h1D + (k-1)*elem%Nnode_h1D**2
269 elem%IntWeight_lgl(l) = &
270 intweight_lgl1dpts_h(i) * intweight_lgl1dpts_h(j) * intweight_lgl1dpts_v(k)
271 end do
272 end do
273 end do
274
275 !* Set the mass matrix
276
277 if (elem%IsLumpedMatrix()) then
278 elem%invM(:,:) = 0.0_rp
279 elem%M(:,:) = 0.0_rp
280 do k=1, elem%Nnode_v
281 do j=1, elem%Nnode_h1D
282 do i=1, elem%Nnode_h1D
283 l = i + (j-1)*elem%Nnode_h1D + (k-1)*elem%Nnode_h1D**2
284 elem%M(l,l) = elem%IntWeight_lgl(l)
285 elem%invM(l,l) = 1.0_dp/elem%IntWeight_lgl(l)
286 end do
287 end do
288 end do
289 else
290 call elementbase_construct_massmat( elem%V, elem%Np, & ! (in)
291 elem%M, elem%invM ) ! (out)
292 end if
293
294 !* Set the stiffness matrix
295
296 call elementbase_construct_stiffmat( elem%M, elem%invM, elem%Dx1, elem%Np, & ! (in)
297 elem%Sx1 ) ! (out)
298 call elementbase_construct_stiffmat( elem%M, elem%invM, elem%Dx2, elem%Np, & ! (in)
299 elem%Sx2 ) ! (out)
300 call elementbase_construct_stiffmat( elem%M, elem%invM, elem%Dx3, elem%Np, & ! (in)
301 elem%Sx3 ) ! (out)
302
303 !* Set the lift matrix
304
305 do k=1, elem%Nnode_v
306 do i=1, elem%Nnode_h1D
307 n = i + (k-1)*elem%Nnode_h1D
308 do p3=1, elem%Nnode_v
309 do p1=1, elem%Nnode_h1D
310 l = p1 + (p3-1)*elem%Nnode_h1D
311 v2d_h(n,l) = p1d_ori_h(i,p1)*p1d_ori_v(k,p3) &
312 * sqrt( (dble(p1-1) + 0.5_dp)*(dble(p3-1) + 0.5_dp) )
313 end do
314 end do
315 end do
316 end do
317 do j=1, elem%Nnode_h1D
318 do i=1, elem%Nnode_h1D
319 n = i + (j-1)*elem%Nnode_h1D
320 do p2=1, elem%Nnode_h1D
321 do p1=1, elem%Nnode_h1D
322 l = p1 + (p2-1)*elem%Nnode_h1D
323 v2d_v(n,l) = p1d_ori_h(i,p1)*p1d_ori_h(j,p2) &
324 * sqrt( (dble(p1-1) + 0.5_dp)*(dble(p2-1) + 0.5_dp) )
325 end do
326 end do
327 end do
328 end do
329
330 !--
331
332 emat(:,:) = 0.0_rp
333 do f=1, elem%Nfaces_h
334 if (elem%IsLumpedMatrix()) then
335 massedge_h(:,:) = 0.0_rp
336 do k=1, elem%Nnode_v
337 do i=1, elem%Nnode_h1D
338 l = i + (k-1)*elem%Nnode_h1D
339 massedge_h(l,l) = intweight_lgl1dpts_h(i) * intweight_lgl1dpts_v(k)
340 end do
341 end do
342 else
343 call elementbase_construct_massmat( v2d_h, elem%Nfp_h, & ! (in)
344 massedge_h ) ! (out)
345 end if
346
347 is = (f-1)*elem%Nfp_h + 1
348 ie = is + elem%Nfp_h - 1
349 emat(elem%Fmask_h(:,f), is:ie) = massedge_h
350 end do
351
352 do f=1, elem%Nfaces_v
353 if (elem%IsLumpedMatrix()) then
354 massedge_v(:,:) = 0.0_rp
355 do j=1, elem%Nnode_h1D
356 do i=1, elem%Nnode_h1D
357 l = i + (j-1)*elem%Nnode_h1D
358 massedge_v(l,l) = intweight_lgl1dpts_h(i) * intweight_lgl1dpts_h(j)
359 end do
360 end do
361 else
362 call elementbase_construct_massmat( v2d_v, elem%Nfp_v, & ! (in)
363 massedge_v ) ! (out)
364 end if
365
366 is = elem%Nfaces_h*elem%Nfp_h + (f-1)*elem%Nfp_v + 1
367 ie = is + elem%Nfp_v - 1
368 emat(elem%Fmask_v(:,f), is:ie) = massedge_v
369 end do
370
371 call elementbase_construct_liftmat( elem%invM, emat, elem%Np, elem%NfpTot, & ! (in)
372 elem%Lift ) ! (out)
373
374 return
375 end subroutine construct_element
376
377!OCL SERIAL
378 function hexhedralelement_gen_intgausslegendreintrpmat( this, IntrpPolyOrder, &
379 intw_intrp, x_intrp, y_intrp, z_intrp ) result(IntrpMat)
380
381 use scale_polynominal, only: &
385
386 implicit none
387
388 class(hexahedralelement), intent(in) :: this
389 integer, intent(in) :: IntrpPolyOrder
390 real(RP), intent(out), optional :: intw_intrp(IntrpPolyOrder**3)
391 real(RP), intent(out), optional :: x_intrp(IntrpPolyOrder**3)
392 real(RP), intent(out), optional :: y_intrp(IntrpPolyOrder**3)
393 real(RP), intent(out), optional :: z_intrp(IntrpPolyOrder**3)
394 real(RP) :: IntrpMat(IntrpPolyOrder**3,this%Np)
395
396 real(RP) :: r_int1D_i(IntrpPolyOrder)
397 real(RP) :: r_int1Dw_i(IntrpPolyOrder)
398 real(RP) :: P_int1D_ori_h(IntrpPolyOrder,this%Nnode_h1D)
399 real(RP) :: P_int1D_ori_v(IntrpPolyOrder,this%Nnode_v)
400 real(RP) :: Vint(IntrpPolyOrder**3,this%Np)
401
402 integer :: p1, p2, p3, p1_, p2_, p3_
403 integer :: n_, l_, m_
404 !-----------------------------------------------------
405
406 r_int1d_i(:) = polynominal_gengausslegendrept( intrppolyorder )
407 r_int1dw_i(:) = polynominal_gengausslegendreptintweight( intrppolyorder )
408 p_int1d_ori_h(:,:) = polynominal_genlegendrepoly( this%PolyOrder_h, r_int1d_i)
409 p_int1d_ori_v(:,:) = polynominal_genlegendrepoly( this%PolyOrder_v, r_int1d_i)
410
411 do p3_=1, intrppolyorder
412 do p2_=1, intrppolyorder
413 do p1_=1, intrppolyorder
414 n_= p1_ + (p2_-1)*intrppolyorder + (p3_-1)*intrppolyorder**2
415 if (present(intw_intrp)) intw_intrp(n_) = r_int1dw_i(p1_) * r_int1dw_i(p2_) * r_int1dw_i(p3_)
416 if (present(x_intrp)) x_intrp(n_) = r_int1d_i(p1_)
417 if (present(y_intrp)) y_intrp(n_) = r_int1d_i(p2_)
418 if (present(z_intrp)) z_intrp(n_) = r_int1d_i(p3_)
419
420 do p3=1, this%Nnode_v
421 do p2=1, this%Nnode_h1D
422 do p1=1, this%Nnode_h1D
423 l_ = p1 + (p2-1)*this%Nnode_h1D + (p3-1)*this%Nnode_h1D**2
424 vint(n_,l_) = p_int1d_ori_h(p1_,p1) * sqrt(dble(p1-1) + 0.5_dp) &
425 * p_int1d_ori_h(p2_,p2) * sqrt(dble(p2-1) + 0.5_dp) &
426 * p_int1d_ori_v(p3_,p3) * sqrt(dble(p3-1) + 0.5_dp)
427 end do
428 end do
429 end do
430 end do
431 end do
432 end do
433 intrpmat(:,:) = matmul(vint, this%invV)
434
435 return
436 end function hexhedralelement_gen_intgausslegendreintrpmat
437
438end module scale_element_hexahedral
module FElib / Element / Base
subroutine, public elementbase3d_init(elem, lumpedmat_flag)
subroutine, public elementbase_construct_massmat(v, np, massmat, invmassmat)
Construct mass matrix M^-1 = V V^T M = ( M^-1 )^-1.
subroutine, public elementbase_construct_stiffmat(massmat, invmassmat, dmat, np, stiffmat)
Construct stiffness matrix StiffMat_i = M^-1 ( M D_xi )^T.
subroutine, public elementbase_construct_liftmat(invm, emat, np, nfptot, liftmat)
Construct stiffness matrix StiffMat_i = M^-1 ( M D_xi )^T.
subroutine, public elementbase3d_final(elem)
module FElib / Element / hexahedron
subroutine hexhedralelement_init(elem, elemorder_h, elemorder_v, lumpedmassmatflag)
module FElib / Element / Quadrilateral
module common / Linear algebra
real(rp) function, dimension(size(a, 1), size(a, 2)), public linalgebra_inv(a)
module common / Polynominal
real(rp) function, dimension(nord), public polynominal_gengausslegendreptintweight(nord)
A function to calcuate the Gauss-Legendre weights.
real(rp) function, dimension(nord), public polynominal_gengausslegendrept(nord)
A function to calcuate the Gauss-Legendre points.
real(rp) function, dimension(size(x), nord+1), public polynominal_genlagrangepoly(nord, x_lgl, x)
A function to obtain the values of Lagrange basis functions which are evaluated over aribitary points...
real(rp) function, dimension(nord+1), public polynominal_gengausslobattopt(nord)
A function to calcuate the Legendre-Gauss-Lobtatto (LGL) points.
real(rp) function, dimension(size(x), nord+1), public polynominal_genlegendrepoly(nord, x)
A function to obtain the values of Legendre polynominals which are evaluated at aribitary points.
real(rp) function, dimension(size(x), nord+1), public polynominal_gendlegendrepoly(nord, x, p)
A function to obtain differential values of Legendre polynominals which are evaluated at aribitary po...
real(rp) function, dimension(nord+1), public polynominal_gengausslobattoptintweight(nord)
A function to calcuate the Gauss-Lobbato weights.
real(rp) function, dimension(nord+1, nord+1), public polynominal_gendlagrangepoly_lglpt(nord, x_lgl)
A function to obtain the differential values of Lagrange basis functions which are evaluated over ari...