34 procedure :: final => hexhedralelement_final
35 procedure :: genintgausslegendreintrpmat => hexhedralelement_gen_intgausslegendreintrpmat
42 private :: construct_element
53 elem, elemOrder_h, elemOrder_v, &
58 integer,
intent(in) :: elemOrder_h
59 integer,
intent(in) :: elemOrder_v
60 logical,
intent(in) :: LumpedMassMatFlag
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
72 elem%Nfaces = elem%Nfaces_h + elem%Nfaces_v
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
78 elem%Np = elem%Nfp_v * elem%Nnode_v
81 call construct_element(elem)
90 subroutine hexhedralelement_final(elem)
99 end subroutine hexhedralelement_final
102 subroutine construct_element(elem)
118 integer :: nodes_ijk(elem%Nnode_h1D, elem%Nnode_h1D, elem%Nnode_v)
120 real(RP) :: lglPts1D_h(elem%Nnode_h1D)
121 real(RP) :: lglPts1D_v(elem%Nnode_v)
123 real(DP) :: intWeight_lgl1DPts_h(elem%Nnode_h1D)
124 real(DP) :: intWeight_lgl1DPts_v(elem%Nnode_v)
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)
139 integer :: p1, p2, p3
145 integer :: fp, fp_h1, fp_h2, fp_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
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 /))
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 /))
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,:)
192 elem%Hslice(:,k) = reshape(nodes_ijk(:,:,k), (/ elem%Nfp_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)
208 call elem2d%Init( elem%PolyOrder_h, .false. )
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)
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
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
244 elem%Dx1(:,:) = 0.0_rp
245 elem%Dx2(:,:) = 0.0_rp
246 elem%Dx3(:,:) = 0.0_rp
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
254 elem%x1(n) = lglpts1d_h(i)
255 elem%x2(n) = lglpts1d_h(j)
256 elem%x3(n) = lglpts1d_v(k)
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))
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)
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)
294 if (elem%IsLumpedMatrix())
then
295 elem%invM(:,:) = 0.0_rp
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)
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) )
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) )
350 do f=1, elem%Nfaces_h
351 if (elem%IsLumpedMatrix())
then
352 massedge_h(:,:) = 0.0_rp
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)
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
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)
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
392 end subroutine construct_element
395 function hexhedralelement_gen_intgausslegendreintrpmat( this, IntrpPolyOrder, &
396 intw_intrp, x_intrp, y_intrp, z_intrp )
result(IntrpMat)
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)
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)
419 integer :: p1, p2, p3, p1_, p2_, p3_
420 integer :: n_, l_, m_
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_)
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)
450 intrpmat(:,:) = matmul(vint, this%invV)
453 end function hexhedralelement_gen_intgausslegendreintrpmat
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.