Choral
Data Types | Functions/Subroutines
fespacexk_mod Module Reference

DERIVED TYPE feSpacexk: define \( Y = [ X_h ]^k \) for \( X_h \) a finite element space. More...

Data Types

interface  clear
 destructor More...
 
interface  fespacexk
 The type feSpacexk defines \( Y = [ X_h ]^k \) for \( X_h \) a finite element space. More...
 
interface  interp_vect_func
 
interface  l2_dist
 integral L2 distance More...
 
interface  l2_dist_grad
 
interface  print
 print a short description More...
 
interface  valid
 

Functions/Subroutines

subroutine fespacexk_clear (Y)
 Destructor for feSpacexk type. More...
 
type(fespacexk) function fespacexk_create (X_h, k)
  Constructor for the type feSpacexk More...
 
logical function fespacexk_valid (Y)
 Check if the structure content is correct. More...
 
subroutine fespacexk_interp_vect_func (u_h, Y, u_1, u_2, u_3)
 Interpolate a function u : R^3 –> R^d given by its components u_1, u_2 and u_3 on Y = [X_h]^d. More...
 
subroutine fespacexk_print (Y)
 Print a short description. More...
 
subroutine, public extract_component (u_c, u, Y, c)
  Extract the component \( c \) of a finite element function \( u \in Y = [X_h]^k \). More...
 
real(rp) function l2_dist_2 (uh, Y, qdm, u_1, u_2, u_3)
 Returns \( \left( \int_O \vert u - u_h \vert^2 \dx \right) ^{1/2} \). More...
 
real(rp) function l2_dist_grad_2 (uh, Y, qdm, grad_u1, grad_u2, grad_u3)
 Returns \( \left( \int_O \vert \nabla u - \nabla u_h \vert^2 dx \right)^{1/2} \). More...
 

Detailed Description

DERIVED TYPE feSpacexk: define \( Y = [ X_h ]^k \) for \( X_h \) a finite element space.

Given a finite element space \( X_h \) of type feSpace:

\(Y\) as \(k\) times more degrees of freedom than \(X_h\). They are ordered as follows:

Author
Charles Pierre

Function/Subroutine Documentation

◆ extract_component()

subroutine, public fespacexk_mod::extract_component ( real(rp), dimension(:), allocatable  u_c,
real(rp), dimension(:), intent(in)  u,
type(fespacexk), intent(in)  Y,
integer, intent(in)  c 
)

Extract the component \( c \) of a finite element function \( u \in Y = [X_h]^k \).

INPUT

  • \( u = (u_1, \dots, u_k) \in Y = [X_h]^k \) a finite element function given by the vector \( u \) corresponding to its decomposition on the basis of \( Y \)
  • \( Y = [X_h]^d \)
  • \( c \), the component to be extracted.

OUTPUT:

  • \( u_c \in X_h \) the component \( c \) of \( u \). It is also a finite element function given by the vector \( u_c \) corresponding to its decomposition on the basis of \( X_h \)

Definition at line 375 of file feSpacexk_mod.f90.

Here is the call graph for this function:

◆ fespacexk_clear()

subroutine fespacexk_mod::fespacexk_clear ( type(fespacexk), intent(inout)  Y)
private

Destructor for feSpacexk type.

Definition at line 163 of file feSpacexk_mod.f90.

◆ fespacexk_create()

type(fespacexk) function fespacexk_mod::fespacexk_create ( type(fespace), target  X_h,
integer, intent(in)  k 
)
private

Constructor for the type feSpacexk

  • OUTPUT:
    • Y = \( [X_h]^k \)
  • INPUT:
    • X_h = feSpace, finite element space
    • k = integer

Definition at line 184 of file feSpacexk_mod.f90.

Here is the call graph for this function:

◆ fespacexk_interp_vect_func()

subroutine fespacexk_mod::fespacexk_interp_vect_func ( real(rp), dimension(:), allocatable  u_h,
type(fespacexk), intent(in)  Y,
procedure(r3tor)  u_1,
procedure(r3tor)  u_2,
procedure(r3tor), optional  u_3 
)
private

Interpolate a function u : R^3 –> R^d given by its components u_1, u_2 and u_3 on Y = [X_h]^d.

If d = 2 : u(x) = [u_1(x), u_2(x)] If d = 3 : u(x) = [u_1(x), u_2(x), u_3(x)]

Definition at line 282 of file feSpacexk_mod.f90.

Here is the call graph for this function:

◆ fespacexk_print()

subroutine fespacexk_mod::fespacexk_print ( type(fespacexk), intent(in)  Y)
private

Print a short description.

Definition at line 342 of file feSpacexk_mod.f90.

◆ fespacexk_valid()

logical function fespacexk_mod::fespacexk_valid ( type (fespacexk), intent(in)  Y)
private

Check if the structure content is correct.

Definition at line 257 of file feSpacexk_mod.f90.

◆ l2_dist_2()

real(rp) function fespacexk_mod::l2_dist_2 ( real(rp), dimension(:), intent(in)  uh,
type(fespacexk), intent(in)  Y,
type(quadmesh), intent(in)  qdm,
procedure(r3tor)  u_1,
procedure(r3tor)  u_2,
procedure(r3tor), optional  u_3 
)
private

Returns \( \left( \int_O \vert u - u_h \vert^2 \dx \right) ^{1/2} \).

  • \( Y = [X_h]^d \) where \( X_h \) is a finite element space on the mesh \( \T \) which is of dimension \( d \).
  • \( u~: \R^3 \mapsto \R^d \) is a vector function, given by ite components:
    \( u = [u_1, u_2]\) if \( d=2\)
    \( u = [u_1, u_2, u_3]\) if \( d=3\)
  • \( u_h\in Y^k \), \( u_h~: \Omega \mapsto \R^d\) is a vector finite element function
  • 'qdm' = integration method on the mesh \( \T \),
  • \( O \) is the integration domain, see quadmesh_mod for a definition.

Definition at line 431 of file feSpacexk_mod.f90.

Here is the call graph for this function:

◆ l2_dist_grad_2()

real(rp) function fespacexk_mod::l2_dist_grad_2 ( real(rp), dimension(:), intent(in)  uh,
type(fespacexk), intent(in)  Y,
type(quadmesh), intent(in)  qdm,
procedure(r3tor3)  grad_u1,
procedure(r3tor3)  grad_u2,
procedure(r3tor3), optional  grad_u3 
)
private

Returns \( \left( \int_O \vert \nabla u - \nabla u_h \vert^2 dx \right)^{1/2} \).

  • \( Y = [X_h]^d \) where \( X_h \) is a finite element space on the mesh \( \T \) which is of dimension \( d \).
  • \( u~: \R^3 \mapsto \R^d \) is a vector function,
  • \( \nabla u~: \R^3 \mapsto \R^d \times \R^d\) is a matrix function, given by ite components:
    \( \nabla u = [\nabla u_1, \nabla u_2]\) if \( d=2\)
    \( \nabla u = (\nabla u_1, \nabla u_2, \nabla u_3]\) if \( d=3\)
    where
    \( u = [u_1, u_2]\) if \( d=2\)
    \( u = [u_1, u_2, u_3]\) if \( d=3\)
  • \( u_h\in Y^k \), \( u_h~: \Omega \mapsto \R^d\) is a vector finite element function
  • 'qdm' = integration method on the mesh \( \T \),
  • \( O \) is the integration domain, see quadmesh_mod for a definition.

Definition at line 514 of file feSpacexk_mod.f90.

Here is the call graph for this function: