FE-Project
Loading...
Searching...
No Matches
scale_element_base.F90
Go to the documentation of this file.
1
9!-------------------------------------------------------------------------------
10#include "scaleFElib.h"
12
13 !-----------------------------------------------------------------------------
14 !
15 !++ used modules
16 !
17 use scale_precision
18
19 !-----------------------------------------------------------------------------
20 implicit none
21 private
22 !-----------------------------------------------------------------------------
23 !
24 !++ Public type & procedure
25 !
26
27 !- Base
28
29 type, public :: elementbase
30 integer :: np
31 integer :: nfaces
32 integer :: nfptot
33 integer :: nv
34 logical, private :: lumpedmatflag
35
36 real(rp), allocatable :: v(:,:)
37 real(rp), allocatable :: invv(:,:)
38 real(rp), allocatable :: m(:,:)
39 real(rp), allocatable :: invm(:,:)
40 real(rp), allocatable :: lift(:,:)
41 real(rp), allocatable :: intweight_lgl(:)
42 contains
43 procedure :: islumpedmatrix => elementbase_islumpedmatrix
44 end type elementbase
48
49 !- 1D
50
51 type, public, extends(elementbase) :: elementbase1d
52 integer :: polyorder
53 integer :: nfp
54 integer, allocatable :: fmask(:,:)
55
56 real(rp), allocatable :: x1(:)
57
58 real(rp), allocatable :: dx1(:,:)
59
60 real(rp), allocatable :: sx1(:,:)
61 end type elementbase1d
62
63 public :: elementbase1d_init
64 public :: elementbase1d_final
65
66 !- 2D
67
68 type, public, extends(elementbase) :: elementbase2d
69 integer :: polyorder
70 integer :: nfp
71 integer, allocatable :: fmask(:,:)
72
73 real(rp), allocatable :: x1(:)
74 real(rp), allocatable :: x2(:)
75
76 real(rp), allocatable :: dx1(:,:)
77 real(rp), allocatable :: dx2(:,:)
78
79 real(rp), allocatable :: sx1(:,:)
80 real(rp), allocatable :: sx2(:,:)
81 end type elementbase2d
82
83 public :: elementbase2d_init
84 public :: elementbase2d_final
85
86 interface
87 function elementbase2d_genintgausslegendreintrpmat(this, IntrpPolyOrder, &
88 intw_intrp, x_intrp, y_intrp ) result(IntrpMat)
89
90 import elementbase2d
91 import rp
92 class(elementbase2d), intent(in) :: this
93 integer, intent(in) :: intrppolyorder
94 real(rp), intent(out), optional :: intw_intrp(intrppolyorder**2)
95 real(rp), intent(out), optional :: x_intrp(intrppolyorder**2)
96 real(rp), intent(out), optional :: y_intrp(intrppolyorder**2)
97 real(rp) :: intrpmat(intrppolyorder**2,this%np)
99 end interface
100
101 !- 3D
102
103 type, public, extends(elementbase) :: elementbase3d
104 integer :: polyorder_h
105 integer :: nnode_h1d
106 integer :: nfaces_h
107 integer :: nfp_h
108 integer, allocatable :: fmask_h(:,:)
109
110 integer :: polyorder_v
111 integer :: nnode_v
112 integer :: nfaces_v
113 integer :: nfp_v
114 integer, allocatable :: fmask_v(:,:)
115
116 integer, allocatable :: colmask(:,:)
117 integer, allocatable :: hslice(:,:)
118 integer, allocatable :: indexh2dto3d(:)
119 integer, allocatable :: indexh2dto3d_bnd(:)
120 integer, allocatable :: indexz1dto3d(:)
121
122 real(rp), allocatable :: x1(:)
123 real(rp), allocatable :: x2(:)
124 real(rp), allocatable :: x3(:)
125
126 real(rp), allocatable :: dx1(:,:)
127 real(rp), allocatable :: dx2(:,:)
128 real(rp), allocatable :: dx3(:,:)
129
130 real(rp), allocatable :: sx1(:,:)
131 real(rp), allocatable :: sx2(:,:)
132 real(rp), allocatable :: sx3(:,:)
133 end type elementbase3d
134
135 public :: elementbase3d_init
136 public :: elementbase3d_final
137
138 !-----------------------------------------------------------------------------
139 !
140 !++ Private type & procedure
141 !
142
143 private :: elementbase_init
144 private :: elementbase_final
145
146contains
147 !-- Base Element ------------------------------------------------------------------------------
148
149!OCL SERIAL
150 subroutine elementbase_init( elem, lumpedmat_flag )
151 implicit none
152
153 class(elementbase), intent(inout) :: elem
154 logical, intent(in) :: lumpedmat_flag
155 !-----------------------------------------------------------------------------
156
157 allocate( elem%M(elem%Np, elem%Np) )
158 allocate( elem%invM(elem%Np, elem%Np) )
159 allocate( elem%V(elem%Np, elem%Np) )
160 allocate( elem%invV(elem%Np, elem%Np) )
161 allocate( elem%Lift(elem%Np, elem%NfpTot) )
162
163 allocate( elem%IntWeight_lgl(elem%Np) )
164
165 elem%LumpedMatFlag = lumpedmat_flag
166
167 return
168 end subroutine elementbase_init
169
170!OCL SERIAL
171 subroutine elementbase_final( elem )
172 implicit none
173
174 class(elementbase), intent(inout) :: elem
175 !-----------------------------------------------------------------------------
176
177 if ( allocated(elem%M) ) then
178 deallocate( elem%M )
179 deallocate( elem%invM )
180 deallocate( elem%V )
181 deallocate( elem%invV )
182 deallocate( elem%Lift )
183 deallocate( elem%IntWeight_lgl )
184 end if
185
186 return
187 end subroutine elementbase_final
188
189!OCL SERIAL
190 function elementbase_islumpedmatrix( elem ) result(lumpedmat_flag)
191 implicit none
192 class(elementbase), intent(in) :: elem
193 logical :: lumpedmat_flag
194 !---------------------------------------------
195
196 lumpedmat_flag = elem%LumpedMatFlag
197 return
198 end function elementbase_islumpedmatrix
199
200
204!OCL SERIAL
205 subroutine elementbase_construct_massmat( V, Np, &
206 MassMat, invMassMat )
208 implicit none
209 integer, intent(in) :: np
210 real(rp), intent(in) :: v(np,np)
211 real(rp), intent(out) :: massmat(np,np)
212 real(rp), intent(out), optional :: invmassmat(np,np)
213
214 real(rp) :: tmpmat(np,np)
215 real(rp) :: invm(np,np)
216 !------------------------------------
217
218 tmpmat(:,:) = transpose(v)
219 invm(:,:) = matmul( v, tmpmat )
220 massmat(:,:) = linalgebra_inv( invm )
221
222 if ( present(invmassmat) ) invmassmat(:,:) = invm(:,:)
223 return
224 end subroutine elementbase_construct_massmat
225
228!OCL SERIAL
229 subroutine elementbase_construct_stiffmat( MassMat, invMassMat, DMat, Np, &
230 StiffMat )
231 implicit none
232 integer, intent(in) :: np
233 real(rp), intent(in) :: massmat(np,np)
234 real(rp), intent(in) :: invmassmat(np,np)
235 real(rp), intent(in) :: dmat(np,np)
236 real(rp), intent(out) :: stiffmat(np,np)
237
238 real(rp) :: tmpmat1(np,np)
239 real(rp) :: tmpmat2(np,np)
240 !------------------------------------
241
242 tmpmat1(:,:) = matmul( massmat, dmat )
243 tmpmat2(:,:) = transpose( tmpmat1 )
244 stiffmat(:,:) = matmul( invmassmat, tmpmat2 )
245
246 return
247 end subroutine elementbase_construct_stiffmat
248
251!OCL SERIAL
252 subroutine elementbase_construct_liftmat( invM, EMat, Np, NfpTot, &
253 LiftMat )
254 implicit none
255 integer, intent(in) :: np
256 integer, intent(in) :: nfptot
257 real(rp), intent(in) :: invm(np,np)
258 real(rp), intent(in) :: emat(np,nfptot)
259 real(rp), intent(out) :: liftmat(np,nfptot)
260 !------------------------------------
261
262 liftmat(:,:) = matmul( invm, emat )
263 return
264 end subroutine elementbase_construct_liftmat
265
266 !-- 1D Element ------------------------------------------------------------------------------
267
268!OCL SERIAL
269 subroutine elementbase1d_init( elem, lumpedmat_flag )
270 implicit none
271
272 class(elementbase1d), intent(inout) :: elem
273 logical, intent(in) :: lumpedmat_flag
274 !-----------------------------------------------------------------------------
275
276 call elementbase_init( elem, lumpedmat_flag )
277
278 allocate( elem%x1(elem%Np) )
279 allocate( elem%Fmask(elem%Nfp, elem%Nfaces) )
280
281 allocate( elem%Dx1(elem%Np, elem%Np) )
282 allocate( elem%Sx1(elem%Np, elem%Np) )
283
284 return
285 end subroutine elementbase1d_init
286
287!OCL SERIAL
288 subroutine elementbase1d_final( elem )
289 implicit none
290
291 class(elementbase1d), intent(inout) :: elem
292 !-----------------------------------------------------------------------------
293
294 if ( allocated( elem%x1 ) ) then
295 deallocate( elem%x1 )
296 deallocate( elem%Dx1 )
297 deallocate( elem%Sx1 )
298 deallocate( elem%Fmask )
299 end if
300
301 call elementbase_final( elem )
302
303 return
304 end subroutine elementbase1d_final
305
306 !-- 2D Element ------------------------------------------------------------------------------
307
308!OCL SERIAL
309 subroutine elementbase2d_init( elem, lumpedmat_flag )
310 implicit none
311
312 class(elementbase2d), intent(inout) :: elem
313 logical, intent(in) :: lumpedmat_flag
314 !-----------------------------------------------------------------------------
315
316 call elementbase_init( elem, lumpedmat_flag )
317
318 allocate( elem%x1(elem%Np), elem%x2(elem%Np) )
319 allocate( elem%Fmask(elem%Nfp, elem%Nfaces) )
320
321 allocate( elem%Dx1(elem%Np, elem%Np), elem%Dx2(elem%Np, elem%Np) )
322 allocate( elem%Sx1(elem%Np, elem%Np), elem%Sx2(elem%Np, elem%Np) )
323
324 return
325 end subroutine elementbase2d_init
326
327!OCL SERIAL
328 subroutine elementbase2d_final( elem )
329 implicit none
330
331 class(elementbase2d), intent(inout) :: elem
332 !-----------------------------------------------------------------------------
333
334 if ( allocated( elem%x1 ) ) then
335 deallocate( elem%x1, elem%x2 )
336 deallocate( elem%Fmask )
337
338 deallocate( elem%Dx1, elem%Dx2 )
339 deallocate( elem%Sx1, elem%Sx2 )
340 end if
341
342 call elementbase_final( elem )
343
344 return
345 end subroutine elementbase2d_final
346
347 !-- 3D Element ------------------------------------------------------------------------------
348
349!OCL SERIAL
350 subroutine elementbase3d_init( elem, lumpedmat_flag )
351 implicit none
352
353 class(elementbase3d), intent(inout) :: elem
354 logical, intent(in) :: lumpedmat_flag
355 !-----------------------------------------------------------------------------
356
357 call elementbase_init( elem, lumpedmat_flag )
358
359 allocate( elem%x1(elem%Np), elem%x2(elem%Np), elem%x3(elem%Np) )
360 allocate( elem%Dx1(elem%Np, elem%Np), elem%Dx2(elem%Np, elem%Np), elem%Dx3(elem%Np, elem%Np) )
361 allocate( elem%Sx1(elem%Np, elem%Np), elem%Sx2(elem%Np, elem%Np), elem%Sx3(elem%Np, elem%Np) )
362 allocate( elem%Fmask_h(elem%Nfp_h, elem%Nfaces_h), elem%Fmask_v(elem%Nfp_v, elem%Nfaces_v) )
363 allocate( elem%Colmask(elem%Nnode_v,elem%Nfp_v))
364 allocate( elem%Hslice(elem%Nfp_v,elem%Nnode_v) )
365 allocate( elem%IndexH2Dto3D(elem%Np) )
366 allocate( elem%IndexH2Dto3D_bnd(elem%NfpTot) )
367 allocate( elem%IndexZ1Dto3D(elem%Np) )
368
369 return
370 end subroutine elementbase3d_init
371
372!OCL SERIAL
373 subroutine elementbase3d_final( elem )
374 implicit none
375
376 class(elementbase3d), intent(inout) :: elem
377 !-----------------------------------------------------------------------------
378
379 if ( allocated( elem%x1 ) ) then
380 deallocate( elem%x1, elem%x2, elem%x3 )
381 deallocate( elem%Dx1, elem%Dx2, elem%Dx3 )
382 deallocate( elem%Sx1, elem%Sx2, elem%Sx3 )
383 deallocate( elem%Fmask_h, elem%Fmask_v )
384 deallocate( elem%Colmask, elem%Hslice )
385 deallocate( elem%IndexH2Dto3D, elem%IndexH2Dto3D_bnd )
386 deallocate( elem%IndexZ1Dto3D )
387 end if
388
389 call elementbase_final( elem )
390
391 return
392 end subroutine elementbase3d_final
393
394end module scale_element_base
module FElib / Element / Base
subroutine, public elementbase2d_final(elem)
subroutine, public elementbase3d_init(elem, lumpedmat_flag)
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.
subroutine, public elementbase3d_final(elem)
subroutine, public elementbase1d_init(elem, lumpedmat_flag)
subroutine, public elementbase1d_final(elem)
module common / Linear algebra
real(rp) function, dimension(size(a, 1), size(a, 2)), public linalgebra_inv(a)