33 procedure :: final => quadrilateralelement_final
34 procedure :: genintgausslegendreintrpmat => quadrilateralelement_gen_intgausslegendreintrpmat
46 integer,
intent(in) :: elemOrder
47 logical,
intent(in) :: LumpedMassMatFlag
51 elem%PolyOrder = elemorder
53 elem%Np = (elemorder + 1)**2
54 elem%Nfp = elemorder + 1
56 elem%NfpTot = elem%Nfp*elem%Nfaces
59 call construct_element(elem)
65 subroutine quadrilateralelement_final(elem)
74 end subroutine quadrilateralelement_final
77 subroutine construct_element(elem)
92 integer :: nodes_ij(elem%Nfp, elem%Nfp)
94 real(RP) :: lglPts1D(elem%Nfp)
95 real(DP) :: intWeight_lgl1DPts(elem%Nfp)
97 real(RP) :: P1D_ori(elem%Nfp, elem%Nfp)
98 real(RP) :: DP1D_ori(elem%Nfp, elem%Nfp)
99 real(RP) :: DLagr1D(elem%Nfp, elem%Nfp)
100 real(RP) :: V1D(elem%Nfp, elem%Nfp)
101 real(RP) :: Emat(elem%Np, elem%NfpTot)
102 real(RP) :: MassEdge(elem%Nfp, elem%Nfp)
104 real(RP) :: eta, etac
105 real(RP) :: filter1D(elem%Nfp), filter2D(elem%Np)
122 nodes_ij(i,j) = i + (j-1)*elem%Nfp
128 elem%Fmask(:,1) = nodes_ij(:,1)
129 elem%Fmask(:,2) = nodes_ij(elem%Nfp,:)
130 elem%Fmask(:,3) = nodes_ij(:,elem%Nfp)
131 elem%Fmask(:,4) = nodes_ij(1,:)
135 elem%Dx1(:,:) = 0.0_rp
136 elem%Dx2(:,:) = 0.0_rp
140 n = i + (j-1)*elem%Nfp
143 elem%x1(n) = lglpts1d(i)
144 elem%x2(n) = lglpts1d(j)
149 l = p1 + (p2-1)*elem%Nfp
150 elem%V(n,l) = (p1d_ori(i,p1)*p1d_ori(j,p2)) &
151 * sqrt((dble(p1-1) + 0.5_dp)*(dble(p2-1) + 0.5_dp))
153 if(p2==j) elem%Dx1(n,l) = dlagr1d(p1,i)
154 if(p1==i) elem%Dx2(n,l) = dlagr1d(p2,j)
167 l = i + (j - 1)*elem%Nfp
168 elem%IntWeight_lgl(l) = &
169 intweight_lgl1dpts(i) * intweight_lgl1dpts(j)
175 if (elem%IsLumpedMatrix())
then
176 elem%invM(:,:) = 0.0_rp
180 l = i + (j - 1)*elem%Nfp
181 elem%M(l,l) = elem%IntWeight_lgl(l)
182 elem%invM(l,l) = 1.0_rp/elem%IntWeight_lgl(l)
200 v1d(:,p1) = p1d_ori(:,p1)*sqrt(dble(p1-1) + 0.5_dp)
206 if (elem%IsLumpedMatrix())
then
209 massedge(l,l) = intweight_lgl1dpts(l)
216 emat(elem%Fmask(:,f), (f-1)*elem%Nfp+1:f*elem%Nfp) = massedge
224 end subroutine construct_element
227 function quadrilateralelement_gen_intgausslegendreintrpmat( this, IntrpPolyOrder, &
228 intw_intrp, x_intrp, y_intrp )
result(IntrpMat)
238 integer,
intent(in) :: IntrpPolyOrder
239 real(RP),
intent(out),
optional :: intw_intrp(IntrpPolyOrder**2)
240 real(RP),
intent(out),
optional :: x_intrp(IntrpPolyOrder**2)
241 real(RP),
intent(out),
optional :: y_intrp(IntrpPolyOrder**2)
242 real(RP) :: IntrpMat(IntrpPolyOrder**2,this%Np)
244 real(RP) :: r_int1D_i(IntrpPolyOrder)
245 real(RP) :: r_int1Dw_i(IntrpPolyOrder)
246 real(RP) :: P_int1D_ori(IntrpPolyOrder,this%PolyOrder+1)
247 real(RP) :: Vint(IntrpPolyOrder**2,(this%PolyOrder+1)**2)
249 integer :: p1, p2, p1_, p2_
257 do p2_=1, intrppolyorder
258 do p1_=1, intrppolyorder
259 n_= p1_ + (p2_-1)*intrppolyorder
260 if (
present(intw_intrp)) intw_intrp(n_) = r_int1dw_i(p1_) * r_int1dw_i(p2_)
261 if (
present(x_intrp)) x_intrp(n_) = r_int1d_i(p1_)
262 if (
present(y_intrp)) y_intrp(n_) = r_int1d_i(p2_)
266 l_ = p1 + (p2-1)*this%Nfp
267 vint(n_,l_) = p_int1d_ori(p1_,p1) * sqrt(dble(p1-1) + 0.5_dp) &
268 * p_int1d_ori(p2_,p2) * sqrt(dble(p2-1) + 0.5_dp)
273 intrpmat(:,:) = matmul(vint, this%invV)
276 end function quadrilateralelement_gen_intgausslegendreintrpmat
280end module scale_element_quadrilateral
module FElib / Element / Base
subroutine, public elementbase2d_final(elem)
subroutine, public elementbase2d_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.
module FElib / Element / Quadrilateral
subroutine quadrilateralelement_init(elem, elemorder, lumpedmassmatflag)
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...