FE-Project
Loading...
Searching...
No Matches
scale_element_line.F90
Go to the documentation of this file.
1!> module FElib / Element / line
2!!
3!! @par Description
4!! A module for a line 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 !-----------------------------------------------------------------------------
23 implicit none
24 private
25
26 !-----------------------------------------------------------------------------
27 !
28 !++ Public type & procedure
29 !
30
31 !> Derived type representing a line element
32 type, public, extends(elementbase1d) :: lineelement
33 contains
34 procedure :: init => lineelement_init
35 procedure :: final => lineelement_final
36 procedure :: genintgausslegendreintrpmat => lineelement_gen_intgausslegendreintrpmat
37 end type lineelement
38
39 !-----------------------------------------------------------------------------
40 !
41 !++ Private procedure
42 !
43 private :: construct_element
44
45contains
46
47!> Initialize an object to manage a hexahedral element
48!!
49!! @param elem Object of finite element
50!! @param elemOrder Polynomial order
51!! @param LumpedMassMatFlag Flag whether mass lumping is considered
52!OCL SERIAL
53 subroutine lineelement_init( &
54 elem, elemOrder, &
55 LumpedMassMatFlag )
56
57 implicit none
58
59 class(lineelement), intent(inout) :: elem
60 integer, intent(in) :: elemOrder
61 logical, intent(in) :: LumpedMassMatFlag
62
63 !-----------------------------------------------------------------------------
64
65 elem%PolyOrder = elemorder
66 elem%Nv = 2
67 elem%Np = elemorder + 1
68 elem%Nfp = 1
69 elem%Nfaces = 2
70 elem%NfpTot = elem%Nfp*elem%Nfaces
71
72 call elementbase1d_init(elem, lumpedmassmatflag)
73 call construct_element(elem)
74
75 return
76 end subroutine lineelement_init
77
78!> Finalize an object to manage a line element
79!!
80!! @param elem Object of finite element
81!OCL SERIAL
82 subroutine lineelement_final(elem)
83 implicit none
84
85 class(lineelement), intent(inout) :: elem
86 !-----------------------------------------------------------------------------
87
88 call elementbase1d_final(elem)
89
90 return
91 end subroutine lineelement_final
92
93!OCL SERIAL
94 subroutine construct_element(elem)
95
97 use scale_polynominal, only: &
101 use scale_element_base, only: &
104
105 implicit none
106
107 type(lineelement), intent(inout) :: elem
108
109 integer :: nodes(elem%Np)
110
111 real(RP) :: lglPts1D(elem%Np)
112 real(DP) :: intWeight_lgl1DPts(elem%Np)
113
114 real(RP) :: P1D_ori(elem%Np, elem%Np)
115 real(RP) :: DP1D_ori(elem%Np, elem%Np)
116 real(RP) :: DLagr1D(elem%Np, elem%Np)
117 real(RP) :: Emat(elem%Np, elem%Nfp*elem%Nfaces)
118 real(RP) :: MassEdge(elem%Nfp, elem%Nfp)
119
120 integer :: i
121 integer :: p1
122 integer :: n, l, f
123 integer :: Nord
124
125 !-----------------------------------------------------------------------------
126
127 lglpts1d(:) = polynominal_gengausslobattopt( elem%PolyOrder )
128
129 p1d_ori(:,:) = polynominal_genlegendrepoly( elem%PolyOrder, lglpts1d )
130 dp1d_ori(:,:) = polynominal_gendlegendrepoly( elem%PolyOrder, lglpts1d, p1d_ori )
131 dlagr1d(:,:) = polynominal_gendlagrangepoly_lglpt(elem%PolyOrder, lglpts1d)
132
133 !* Preparation
134
135 do i=1, elem%Np
136 nodes(i) = i
137 end do
138
139 ! Set the mask to extract the values at faces
140
141 elem%Fmask(:,1) = 1
142 elem%Fmask(:,2) = elem%Np
143
144 !* Set the coordinates of LGL points, and the Vandermonde and differential matricies
145
146 elem%Dx1(:,:) = 0.0_rp
147
148 do n=1, elem%Np
149 !* Set the coordinates of LGL points
150 elem%x1(n) = lglpts1d(n)
151
152 !* Set the Vandermonde and differential matricies
153 do l=1, elem%Np
154 elem%V(n,l) = p1d_ori(n,l) * sqrt(dble(l-1) + 0.5_rp)
155 elem%Dx1(n,l) = dlagr1d(l,n)
156 end do
157 end do
158 elem%invV(:,:) = linalgebra_inv(elem%V)
159
160 !* Set the weights at LGL points to integrate over element
161 elem%IntWeight_lgl(:) = polynominal_gengausslobattoptintweight(elem%PolyOrder)
162
163 !* Set the mass matrix
164
165 if (elem%IsLumpedMatrix()) then
166 elem%invM(:,:) = 0.0_rp
167 elem%M(:,:) = 0.0_rp
168 do i=1, elem%Np
169 elem%M(i,i) = elem%IntWeight_lgl(i)
170 elem%invM(i,i) = 1.0_rp/elem%IntWeight_lgl(i)
171 end do
172 else
173 call elementbase_construct_massmat( elem%V, elem%Np, & ! (in)
174 elem%M, elem%invM ) ! (out)
175 end if
176
177 !* Set the stiffness matrix
178 call elementbase_construct_stiffmat( elem%M, elem%invM, elem%Dx1, elem%Np, & ! (in)
179 elem%Sx1 ) ! (out)
180
181 !* Set the lift matrix
182
183 emat(:,:) = 0.0_rp
184 do f=1, elem%Nfaces
185 massedge(:,:) = 0.0_rp
186 do l=1, elem%Nfp
187 massedge(l,l) = 1.0_rp
188 end do
189 emat(elem%Fmask(:,f), (f-1)*elem%Nfp+1:f*elem%Nfp) = massedge
190 end do
191 call elementbase_construct_liftmat( elem%invM, emat, elem%Np, elem%NfpTot, & ! (in)
192 elem%Lift ) ! (out)
193
194 return
195 end subroutine construct_element
196
197!OCL SERIAL
198 function lineelement_gen_intgausslegendreintrpmat( this, IntrpPolyOrder, &
199 intw_intrp, x_intrp ) result(IntrpMat)
200
201 use scale_polynominal, only: &
205
206 implicit none
207
208 class(lineelement), intent(in) :: this
209 integer, intent(in) :: IntrpPolyOrder
210 real(RP), intent(out), optional :: intw_intrp(IntrpPolyOrder)
211 real(RP), intent(out), optional :: x_intrp(IntrpPolyOrder)
212 real(RP) :: IntrpMat(IntrpPolyOrder,this%Np)
213
214 real(RP) :: r_int1D_i(IntrpPolyOrder)
215 real(RP) :: r_int1Dw_i(IntrpPolyOrder)
216 real(RP) :: P_int1D_ori(IntrpPolyOrder,this%PolyOrder+1)
217 real(RP) :: Vint(IntrpPolyOrder,this%PolyOrder+1)
218
219 integer :: p1, p1_
220 !-----------------------------------------------------
221
222 r_int1d_i(:) = polynominal_gengausslegendrept( intrppolyorder )
223 r_int1dw_i(:) = polynominal_gengausslegendreptintweight( intrppolyorder )
224 p_int1d_ori(:,:) = polynominal_genlegendrepoly( this%PolyOrder, r_int1d_i)
225
226 do p1_=1, intrppolyorder
227 if (present(intw_intrp)) intw_intrp(p1_) = r_int1dw_i(p1_)
228 if (present(x_intrp)) x_intrp(p1_) = r_int1d_i(p1_)
229 do p1=1, this%Np
230 vint(p1_,p1) = p_int1d_ori(p1_,p1) * sqrt(real(p1-1,kind=rp) + 0.5_rp)
231 end do
232 end do
233 intrpmat(:,:) = matmul(vint, this%invV)
234
235 return
236 end function lineelement_gen_intgausslegendreintrpmat
237
238end module scale_element_line
module FElib / Element / Base
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 elementbase1d_init(elem, lumpedmat_flag)
Initialize an object to manage a 1D reference element.
subroutine, public elementbase1d_final(elem)
Finalize an object to manage a 1D reference element.
module FElib / Element / line
subroutine lineelement_init(elem, elemorder, lumpedmassmatflag)
Initialize an object to manage a hexahedral element.
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 1D reference element.
Derived type representing a line element.