FE-Project
Loading...
Searching...
No Matches
scale_element_line.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(elementbase1d) :: lineelement
31 contains
32 procedure :: init => lineelement_init
33 procedure :: final => lineelement_final
34 procedure :: genintgausslegendreintrpmat => lineelement_gen_intgausslegendreintrpmat
35 end type lineelement
36
37contains
38!OCL SERIAL
39 subroutine lineelement_init( &
40 elem, elemOrder, &
41 LumpedMassMatFlag )
42
43 implicit none
44
45 class(lineelement), intent(inout) :: elem
46 integer, intent(in) :: elemOrder
47 logical, intent(in) :: LumpedMassMatFlag
48
49 !-----------------------------------------------------------------------------
50
51 elem%PolyOrder = elemorder
52 elem%Nv = 2
53 elem%Np = elemorder + 1
54 elem%Nfp = 1
55 elem%Nfaces = 2
56 elem%NfpTot = elem%Nfp*elem%Nfaces
57
58 call elementbase1d_init(elem, lumpedmassmatflag)
59 call construct_element(elem)
60
61 return
62 end subroutine lineelement_init
63
64!OCL SERIAL
65 subroutine lineelement_final(elem)
66 implicit none
67
68 class(lineelement), intent(inout) :: elem
69 !-----------------------------------------------------------------------------
70
71 call elementbase1d_final(elem)
72
73 return
74 end subroutine lineelement_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(lineelement), intent(inout) :: elem
91
92 integer :: nodes(elem%Np)
93
94 real(RP) :: lglPts1D(elem%Np)
95 real(DP) :: intWeight_lgl1DPts(elem%Np)
96
97 real(RP) :: P1D_ori(elem%Np, elem%Np)
98 real(RP) :: DP1D_ori(elem%Np, elem%Np)
99 real(RP) :: DLagr1D(elem%Np, elem%Np)
100 real(RP) :: Emat(elem%Np, elem%Nfp*elem%Nfaces)
101 real(RP) :: MassEdge(elem%Nfp, elem%Nfp)
102
103 integer :: i
104 integer :: p1
105 integer :: n, l, f
106 integer :: Nord
107
108 !-----------------------------------------------------------------------------
109
110 lglpts1d(:) = polynominal_gengausslobattopt( elem%PolyOrder )
111
112 p1d_ori(:,:) = polynominal_genlegendrepoly( elem%PolyOrder, lglpts1d )
113 dp1d_ori(:,:) = polynominal_gendlegendrepoly( elem%PolyOrder, lglpts1d, p1d_ori )
114 dlagr1d(:,:) = polynominal_gendlagrangepoly_lglpt(elem%PolyOrder, lglpts1d)
115
116 !* Preparation
117
118 do i=1, elem%Np
119 nodes(i) = i
120 end do
121
122 ! Set the mask to extract the values at faces
123
124 elem%Fmask(:,1) = 1
125 elem%Fmask(:,2) = elem%Np
126
127 !* Set the coordinates of LGL points, and the Vandermonde and differential matricies
128
129 elem%Dx1(:,:) = 0.0_rp
130
131 do n=1, elem%Np
132 !* Set the coordinates of LGL points
133 elem%x1(n) = lglpts1d(n)
134
135 !* Set the Vandermonde and differential matricies
136 do l=1, elem%Np
137 elem%V(n,l) = p1d_ori(n,l) * sqrt(dble(l-1) + 0.5_rp)
138 elem%Dx1(n,l) = dlagr1d(l,n)
139 end do
140 end do
141 elem%invV(:,:) = linalgebra_inv(elem%V)
142
143 !* Set the weights at LGL points to integrate over element
144 elem%IntWeight_lgl(:) = polynominal_gengausslobattoptintweight(elem%PolyOrder)
145
146 !* Set the mass matrix
147
148 if (elem%IsLumpedMatrix()) then
149 elem%invM(:,:) = 0.0_rp
150 elem%M(:,:) = 0.0_rp
151 do i=1, elem%Np
152 elem%M(i,i) = elem%IntWeight_lgl(i)
153 elem%invM(i,i) = 1.0_rp/elem%IntWeight_lgl(i)
154 end do
155 else
156 call elementbase_construct_massmat( elem%V, elem%Np, & ! (in)
157 elem%M, elem%invM ) ! (out)
158 end if
159
160 !* Set the stiffness matrix
161 call elementbase_construct_stiffmat( elem%M, elem%invM, elem%Dx1, elem%Np, & ! (in)
162 elem%Sx1 ) ! (out)
163
164 !* Set the lift matrix
165
166 emat(:,:) = 0.0_rp
167 do f=1, elem%Nfaces
168 massedge(:,:) = 0.0_rp
169 do l=1, elem%Nfp
170 massedge(l,l) = 1.0_rp
171 end do
172 emat(elem%Fmask(:,f), (f-1)*elem%Nfp+1:f*elem%Nfp) = massedge
173 end do
174 call elementbase_construct_liftmat( elem%invM, emat, elem%Np, elem%NfpTot, & ! (in)
175 elem%Lift ) ! (out)
176
177 return
178 end subroutine construct_element
179
180!OCL SERIAL
181 function lineelement_gen_intgausslegendreintrpmat( this, IntrpPolyOrder, &
182 intw_intrp, x_intrp ) result(IntrpMat)
183
184 use scale_polynominal, only: &
188
189 implicit none
190
191 class(lineelement), intent(in) :: this
192 integer, intent(in) :: IntrpPolyOrder
193 real(RP), intent(out), optional :: intw_intrp(IntrpPolyOrder)
194 real(RP), intent(out), optional :: x_intrp(IntrpPolyOrder)
195 real(RP) :: IntrpMat(IntrpPolyOrder,this%Np)
196
197 real(RP) :: r_int1D_i(IntrpPolyOrder)
198 real(RP) :: r_int1Dw_i(IntrpPolyOrder)
199 real(RP) :: P_int1D_ori(IntrpPolyOrder,this%PolyOrder+1)
200 real(RP) :: Vint(IntrpPolyOrder,this%PolyOrder+1)
201
202 integer :: p1, p1_
203 !-----------------------------------------------------
204
205 r_int1d_i(:) = polynominal_gengausslegendrept( intrppolyorder )
206 r_int1dw_i(:) = polynominal_gengausslegendreptintweight( intrppolyorder )
207 p_int1d_ori(:,:) = polynominal_genlegendrepoly( this%PolyOrder, r_int1d_i)
208
209 do p1_=1, intrppolyorder
210 if (present(intw_intrp)) intw_intrp(p1_) = r_int1dw_i(p1_)
211 if (present(x_intrp)) x_intrp(p1_) = r_int1d_i(p1_)
212 do p1=1, this%Np
213 vint(p1_,p1) = p_int1d_ori(p1_,p1) * sqrt(real(p1-1,kind=rp) + 0.5_rp)
214 end do
215 end do
216 intrpmat(:,:) = matmul(vint, this%invV)
217
218 return
219 end function lineelement_gen_intgausslegendreintrpmat
220
221end 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)
subroutine, public elementbase1d_final(elem)
module FElib / Element / line
subroutine lineelement_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...