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