FE-Project
All Classes Namespaces Files Functions Variables Pages
scale_element_quadrilateral.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 !-----------------------------------------------------------------------------
23 implicit none
24 private
25
26 !-----------------------------------------------------------------------------
27 !
28 !++ Public type & procedure
29 !
30 type, public, extends(elementbase2d) :: quadrilateralelement
31 contains
32 procedure :: init => quadrilateralelement_init
33 procedure :: final => quadrilateralelement_final
34 procedure :: genintgausslegendreintrpmat => quadrilateralelement_gen_intgausslegendreintrpmat
36
37contains
38!OCL SERIAL
40 elem, elemOrder, &
41 LumpedMassMatFlag )
42
43 implicit none
44
45 class(quadrilateralelement), intent(inout) :: elem
46 integer, intent(in) :: elemOrder
47 logical, intent(in) :: LumpedMassMatFlag
48
49 !-----------------------------------------------------------------------------
50
51 elem%PolyOrder = elemorder
52 elem%Nv = 4
53 elem%Np = (elemorder + 1)**2
54 elem%Nfp = elemorder + 1
55 elem%Nfaces = 4
56 elem%NfpTot = elem%Nfp*elem%Nfaces
57
58 call elementbase2d_init(elem, lumpedmassmatflag)
59 call construct_element(elem)
60
61 return
62 end subroutine quadrilateralelement_init
63
64!OCL SERIAL
65 subroutine quadrilateralelement_final(elem)
66 implicit none
67
68 class(quadrilateralelement), intent(inout) :: elem
69 !-----------------------------------------------------------------------------
70
71 call elementbase2d_final(elem)
72
73 return
74 end subroutine quadrilateralelement_final
75
76!OCL SERIAL
77 subroutine construct_element(elem)
78
80 use scale_polynominal, only: &
84 use scale_element_base, only: &
87
88 implicit none
89
90 type(quadrilateralelement), intent(inout) :: elem
91
92 integer :: nodes_ij(elem%Nfp, elem%Nfp)
93
94 real(RP) :: lglPts1D(elem%Nfp)
95 real(DP) :: intWeight_lgl1DPts(elem%Nfp)
96
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)
103
104 real(RP) :: eta, etac
105 real(RP) :: filter1D(elem%Nfp), filter2D(elem%Np)
106
107 integer :: i, j
108 integer :: p1, p2
109 integer :: n, l, f
110 integer :: Nord
111 !-----------------------------------------------------------------------------
112
113 lglpts1d(:) = polynominal_gengausslobattopt( elem%PolyOrder )
114 p1d_ori(:,:) = polynominal_genlegendrepoly( elem%PolyOrder, lglpts1d )
115 dp1d_ori(:,:) = polynominal_gendlegendrepoly( elem%PolyOrder, lglpts1d, p1d_ori )
116 dlagr1d(:,:) = polynominal_gendlagrangepoly_lglpt(elem%PolyOrder, lglpts1d)
117
118 !* Preparation
119
120 do j=1, elem%Nfp
121 do i=1, elem%Nfp
122 nodes_ij(i,j) = i + (j-1)*elem%Nfp
123 end do
124 end do
125
126 ! Set the mask to extract the values at faces
127
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,:)
132
133 !* Set the coordinates of LGL points, and the Vandermonde and differential matricies
134
135 elem%Dx1(:,:) = 0.0_rp
136 elem%Dx2(:,:) = 0.0_rp
137
138 do j=1, elem%Nfp
139 do i=1, elem%Nfp
140 n = i + (j-1)*elem%Nfp
141
142 !* Set the coordinates of LGL points
143 elem%x1(n) = lglpts1d(i)
144 elem%x2(n) = lglpts1d(j)
145
146 !* Set the Vandermonde and differential matricies
147 do p2=1, elem%Nfp
148 do p1=1, elem%Nfp
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))
152
153 if(p2==j) elem%Dx1(n,l) = dlagr1d(p1,i)
154 if(p1==i) elem%Dx2(n,l) = dlagr1d(p2,j)
155 end do
156 end do
157 end do
158 end do
159 elem%invV(:,:) = linalgebra_inv(elem%V)
160
161 !* Set the weights at LGL points to integrate over element
162
163 intweight_lgl1dpts(:) = polynominal_gengausslobattoptintweight(elem%PolyOrder)
164
165 do j=1, elem%Nfp
166 do i=1, elem%Nfp
167 l = i + (j - 1)*elem%Nfp
168 elem%IntWeight_lgl(l) = &
169 intweight_lgl1dpts(i) * intweight_lgl1dpts(j)
170 end do
171 end do
172
173 !* Set the mass matrix
174
175 if (elem%IsLumpedMatrix()) then
176 elem%invM(:,:) = 0.0_rp
177 elem%M(:,:) = 0.0_rp
178 do j=1, elem%Nfp
179 do i=1, elem%Nfp
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)
183 end do
184 end do
185 else
186 call elementbase_construct_massmat( elem%V, elem%Np, & ! (in)
187 elem%M, elem%invM ) ! (out)
188 end if
189
190 !* Set the stiffness matrix
191
192 call elementbase_construct_stiffmat( elem%M, elem%invM, elem%Dx1, elem%Np, & ! (in)
193 elem%Sx1 ) ! (out)
194 call elementbase_construct_stiffmat( elem%M, elem%invM, elem%Dx2, elem%Np, & ! (in)
195 elem%Sx2 ) ! (out)
196
197 !* Set the lift matrix
198
199 do p1=1, elem%Nfp
200 v1d(:,p1) = p1d_ori(:,p1)*sqrt(dble(p1-1) + 0.5_dp)
201 end do
202
203 emat(:,:) = 0.0_rp
204 do f=1, elem%Nfaces
205
206 if (elem%IsLumpedMatrix()) then
207 massedge = 0.0_rp
208 do l=1, elem%Nfp
209 massedge(l,l) = intweight_lgl1dpts(l)
210 end do
211 else
212 call elementbase_construct_massmat( v1d, elem%Nfp, & ! (in)
213 massedge ) ! (out)
214 end if
215
216 emat(elem%Fmask(:,f), (f-1)*elem%Nfp+1:f*elem%Nfp) = massedge
217 end do
218 call elementbase_construct_liftmat( elem%invM, emat, elem%Np, elem%NfpTot, & ! (in)
219 elem%Lift ) ! (out)
220
221 !* Construct filter matrix
222
223 return
224 end subroutine construct_element
225
226!OCL SERIAL
227 function quadrilateralelement_gen_intgausslegendreintrpmat( this, IntrpPolyOrder, &
228 intw_intrp, x_intrp, y_intrp ) result(IntrpMat)
229
230 use scale_polynominal, only: &
234
235 implicit none
236
237 class(quadrilateralelement), intent(in) :: this
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)
243
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)
248
249 integer :: p1, p2, p1_, p2_
250 integer :: n_, l_
251 !-----------------------------------------------------
252
253 r_int1d_i(:) = polynominal_gengausslegendrept( intrppolyorder )
254 r_int1dw_i(:) = polynominal_gengausslegendreptintweight( intrppolyorder )
255 p_int1d_ori(:,:) = polynominal_genlegendrepoly( this%PolyOrder, r_int1d_i)
256
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_)
263
264 do p2=1, this%Nfp
265 do p1=1, this%Nfp
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)
269 end do
270 end do
271 end do
272 end do
273 intrpmat(:,:) = matmul(vint, this%invV)
274
275 return
276 end function quadrilateralelement_gen_intgausslegendreintrpmat
277
278 !-------------------
279
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...