FE-Project
Loading...
Searching...
No Matches
scale_element_base Module Reference

module FElib / Element / Base More...

Data Types

type  elementbase
 Derived type representing an arbitrary finite element. More...
type  elementbase1d
 Derived type representing a 1D reference element. More...
type  elementbase2d
 Derived type representing a 2D reference element. More...
interface  elementbase2d_genintgausslegendreintrpmat
type  elementbase3d
 Derived type representing a 3D reference element. More...

Functions/Subroutines

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.
subroutine, public elementbase2d_init (elem, lumpedmat_flag)
 Initialize an object to manage a 2D reference element.
subroutine, public elementbase2d_final (elem)
 Finalize an object to manage a 2D reference element.
subroutine, public elementbase3d_init (elem, lumpedmat_flag)
 Initialize an object to manage a 3D reference element.
subroutine, public elementbase3d_final (elem)
 Finalize an object to manage a 3D reference element.

Detailed Description

module FElib / Element / Base

Description
A base module for finite element
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ elementbase_construct_massmat()

subroutine, public scale_element_base::elementbase_construct_massmat ( real(rp), dimension(np,np), intent(in) v,
integer, intent(in) np,
real(rp), dimension(np,np), intent(out) massmat,
real(rp), dimension(np,np), intent(out), optional invmassmat )

Construct mass matrix M^-1 = V V^T M = ( M^-1 )^-1.

Definition at line 212 of file scale_element_base.F90.

215 implicit none
216 integer, intent(in) :: Np
217 real(RP), intent(in) :: V(Np,Np)
218 real(RP), intent(out) :: MassMat(Np,Np)
219 real(RP), intent(out), optional :: invMassMat(Np,Np)
220
221 real(RP) :: tmpMat(Np,Np)
222 real(RP) :: invM(Np,Np)
223 !------------------------------------
224
225 tmpmat(:,:) = transpose(v)
226 invm(:,:) = matmul( v, tmpmat )
227 massmat(:,:) = linalgebra_inv( invm )
228
229 if ( present(invmassmat) ) invmassmat(:,:) = invm(:,:)
230 return
Module common / Linear algebra.
real(rp) function, dimension(size(a, 1), size(a, 2)), public linalgebra_inv(a)
Calculate a inversion of matrix A.

References scale_linalgebra::linalgebra_inv().

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_base::elementbase::islumpedmatrix(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ elementbase_construct_stiffmat()

subroutine, public scale_element_base::elementbase_construct_stiffmat ( real(rp), dimension(np,np), intent(in) massmat,
real(rp), dimension(np,np), intent(in) invmassmat,
real(rp), dimension(np,np), intent(in) dmat,
integer, intent(in) np,
real(rp), dimension(np,np), intent(out) stiffmat )

Construct stiffness matrix StiffMat_i = M^-1 ( M D_xi )^T.

Definition at line 236 of file scale_element_base.F90.

238 implicit none
239 integer, intent(in) :: Np
240 real(RP), intent(in) :: MassMat(Np,Np)
241 real(RP), intent(in) :: invMassMat(Np,Np)
242 real(RP), intent(in) :: DMat(Np,Np)
243 real(RP), intent(out) :: StiffMat(Np,Np)
244
245 real(RP) :: tmpMat1(Np,Np)
246 real(RP) :: tmpMat2(Np,Np)
247 !------------------------------------
248
249 tmpmat1(:,:) = matmul( massmat, dmat )
250 tmpmat2(:,:) = transpose( tmpmat1 )
251 stiffmat(:,:) = matmul( invmassmat, tmpmat2 )
252
253 return

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_base::elementbase::islumpedmatrix(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ elementbase_construct_liftmat()

subroutine, public scale_element_base::elementbase_construct_liftmat ( real(rp), dimension(np,np), intent(in) invm,
real(rp), dimension(np,nfptot), intent(in) emat,
integer, intent(in) np,
integer, intent(in) nfptot,
real(rp), dimension(np,nfptot), intent(out) liftmat )

Construct stiffness matrix StiffMat_i = M^-1 ( M D_xi )^T.

Definition at line 259 of file scale_element_base.F90.

261 implicit none
262 integer, intent(in) :: Np
263 integer, intent(in) :: NfpTot
264 real(RP), intent(in) :: invM(Np,Np)
265 real(RP), intent(in) :: EMat(Np,NfpTot)
266 real(RP), intent(out) :: LiftMat(Np,NfpTot)
267 !------------------------------------
268
269 liftmat(:,:) = matmul( invm, emat )
270 return

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_base::elementbase::islumpedmatrix(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ elementbase1d_init()

subroutine, public scale_element_base::elementbase1d_init ( class(elementbase1d), intent(inout) elem,
logical, intent(in) lumpedmat_flag )

Initialize an object to manage a 1D reference element.

Parameters
elemObject of finite element
elemFlag whether mass lumping is considered

Definition at line 280 of file scale_element_base.F90.

281 implicit none
282
283 class(ElementBase1D), intent(inout) :: elem
284 logical, intent(in) :: lumpedmat_flag
285 !-----------------------------------------------------------------------------
286
287 call elementbase_init( elem, lumpedmat_flag )
288
289 allocate( elem%x1(elem%Np) )
290 allocate( elem%Fmask(elem%Nfp, elem%Nfaces) )
291
292 allocate( elem%Dx1(elem%Np, elem%Np) )
293 allocate( elem%Sx1(elem%Np, elem%Np) )
294
295 return

Referenced by scale_element_line::lineelement_init().

◆ elementbase1d_final()

subroutine, public scale_element_base::elementbase1d_final ( class(elementbase1d), intent(inout) elem)

Finalize an object to manage a 1D reference element.

Definition at line 300 of file scale_element_base.F90.

301 implicit none
302
303 class(ElementBase1D), intent(inout) :: elem
304 !-----------------------------------------------------------------------------
305
306 if ( allocated( elem%x1 ) ) then
307 deallocate( elem%x1 )
308 deallocate( elem%Dx1 )
309 deallocate( elem%Sx1 )
310 deallocate( elem%Fmask )
311 end if
312
313 call elementbase_final( elem )
314
315 return

Referenced by scale_element_line::lineelement_init().

◆ elementbase2d_init()

subroutine, public scale_element_base::elementbase2d_init ( class(elementbase2d), intent(inout) elem,
logical, intent(in) lumpedmat_flag )

Initialize an object to manage a 2D reference element.

Parameters
elemObject of finite element
elemFlag whether mass lumping is considered

Definition at line 325 of file scale_element_base.F90.

326 implicit none
327
328 class(ElementBase2D), intent(inout) :: elem
329 logical, intent(in) :: lumpedmat_flag
330 !-----------------------------------------------------------------------------
331
332 call elementbase_init( elem, lumpedmat_flag )
333
334 allocate( elem%x1(elem%Np), elem%x2(elem%Np) )
335 allocate( elem%Fmask(elem%Nfp, elem%Nfaces) )
336
337 allocate( elem%Dx1(elem%Np, elem%Np), elem%Dx2(elem%Np, elem%Np) )
338 allocate( elem%Sx1(elem%Np, elem%Np), elem%Sx2(elem%Np, elem%Np) )
339
340 return

Referenced by scale_element_quadrilateral::quadrilateralelement_init().

◆ elementbase2d_final()

subroutine, public scale_element_base::elementbase2d_final ( class(elementbase2d), intent(inout) elem)

Finalize an object to manage a 2D reference element.

Definition at line 345 of file scale_element_base.F90.

346 implicit none
347
348 class(ElementBase2D), intent(inout) :: elem
349 !-----------------------------------------------------------------------------
350
351 if ( allocated( elem%x1 ) ) then
352 deallocate( elem%x1, elem%x2 )
353 deallocate( elem%Fmask )
354
355 deallocate( elem%Dx1, elem%Dx2 )
356 deallocate( elem%Sx1, elem%Sx2 )
357 end if
358
359 call elementbase_final( elem )
360
361 return

Referenced by scale_element_quadrilateral::quadrilateralelement_init().

◆ elementbase3d_init()

subroutine, public scale_element_base::elementbase3d_init ( class(elementbase3d), intent(inout) elem,
logical, intent(in) lumpedmat_flag )

Initialize an object to manage a 3D reference element.

Parameters
elemObject of finite element
elemFlag whether mass lumping is considered

Definition at line 371 of file scale_element_base.F90.

372 implicit none
373
374 class(ElementBase3D), intent(inout) :: elem
375 logical, intent(in) :: lumpedmat_flag
376 !-----------------------------------------------------------------------------
377
378 call elementbase_init( elem, lumpedmat_flag )
379
380 allocate( elem%x1(elem%Np), elem%x2(elem%Np), elem%x3(elem%Np) )
381 allocate( elem%Dx1(elem%Np, elem%Np), elem%Dx2(elem%Np, elem%Np), elem%Dx3(elem%Np, elem%Np) )
382 allocate( elem%Sx1(elem%Np, elem%Np), elem%Sx2(elem%Np, elem%Np), elem%Sx3(elem%Np, elem%Np) )
383 allocate( elem%Fmask_h(elem%Nfp_h, elem%Nfaces_h), elem%Fmask_v(elem%Nfp_v, elem%Nfaces_v) )
384 allocate( elem%Colmask(elem%Nnode_v,elem%Nfp_v))
385 allocate( elem%Hslice(elem%Nfp_v,elem%Nnode_v) )
386 allocate( elem%IndexH2Dto3D(elem%Np) )
387 allocate( elem%IndexH2Dto3D_bnd(elem%NfpTot) )
388 allocate( elem%IndexZ1Dto3D(elem%Np) )
389
390 return

Referenced by scale_element_hexahedral::hexhedralelement_init().

◆ elementbase3d_final()

subroutine, public scale_element_base::elementbase3d_final ( class(elementbase3d), intent(inout) elem)

Finalize an object to manage a 3D reference element.

Definition at line 395 of file scale_element_base.F90.

396 implicit none
397
398 class(ElementBase3D), intent(inout) :: elem
399 !-----------------------------------------------------------------------------
400
401 if ( allocated( elem%x1 ) ) then
402 deallocate( elem%x1, elem%x2, elem%x3 )
403 deallocate( elem%Dx1, elem%Dx2, elem%Dx3 )
404 deallocate( elem%Sx1, elem%Sx2, elem%Sx3 )
405 deallocate( elem%Fmask_h, elem%Fmask_v )
406 deallocate( elem%Colmask, elem%Hslice )
407 deallocate( elem%IndexH2Dto3D, elem%IndexH2Dto3D_bnd )
408 deallocate( elem%IndexZ1Dto3D )
409 end if
410
411 call elementbase_final( elem )
412
413 return

Referenced by scale_element_hexahedral::hexhedralelement_init().