The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines | Variables
oft_h1_operators Module Reference

Detailed Description

H^1 FE operator definitions.

Authors
Chris Hansen
Date
August 2011

Data Types

type  oft_h1_ginterp
 Interpolate \( \nabla \) of a H^1 field. More...
 
type  oft_h1_rinterp
 Interpolate a H^1 field. More...
 
type  oft_h1_zerob
 Needs docs. More...
 
type  oft_h1_zerogrnd
 Needs docs. More...
 
type  oft_h1_zeroi
 Needs docs. More...
 

Functions/Subroutines

subroutine ginterp_apply (self, cell, f, gop, val)
 Reconstruct the gradient of a H^1 scalar field.
 
subroutine h1_base_pop (self, acors, afine)
 Transfer a base level H^1 scalar field to the next MPI level.
 
subroutine h1_base_push (self, afine, acors)
 Transfer a MPI level H^1 scalar field to the base level.
 
subroutine h1_getlop_pre (ml_h1_rep, pre, mats, bc_type, level, nlevels)
 Compute eigenvalues and smoothing coefficients for the operator H^1::LOP.
 
subroutine h1_lop_eigs (ml_h1_rep, minlev)
 Compute eigenvalues and smoothing coefficients for the operator H^1::LOP.
 
subroutine h1_mloptions
 Read-in options for the basic H^1 ML preconditioners.
 
subroutine h1_setup_interp (ml_h1_rep)
 Construct interpolation matrices for transfer between H^1 finite element spaces.
 
subroutine oft_h1_getlop (fe_rep, mat, bc)
 Construct laplacian matrix for H^1 scalar representation.
 
subroutine oft_h1_getmop (fe_rep, mat, bc)
 Construct mass matrix for H^1 scalar representation.
 
subroutine oft_h1_project (fe_rep, field, x)
 Project a scalar field onto a H^1 basis.
 
subroutine rinterp_apply (self, cell, f, gop, val)
 Reconstruct a H^1 scalar field.
 
subroutine rinterp_delete (self)
 Destroy temporary internal storage.
 
subroutine rinterp_setup (self, h1_rep)
 Setup interpolator for H^1 scalar fields.
 
subroutine zerob_apply (self, a)
 Zero a H^1 scalar field at all boundary nodes.
 
subroutine zerob_delete (self)
 Zero a H^1 scalar field at all boundary nodes.
 
subroutine zerogrnd_apply (self, a)
 Zero a H^1 scalar field at the global grounding node.
 
subroutine zeroi_apply (self, a)
 Zero a H^1 scalar field at all interior nodes.
 

Variables

real(r8), dimension(fem_max_levels) df_lop =-1.d99
 
integer(i4), dimension(fem_max_levels) nu_lop =0
 
real(r8), dimension(:,:), pointer oft_h1_gop => NULL()
 
real(r8), dimension(:), pointer oft_h1_rop => NULL()
 

Function/Subroutine Documentation

◆ ginterp_apply()

subroutine ginterp_apply ( class(oft_h1_ginterp), intent(inout)  self,
integer(i4), intent(in)  cell,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(3,4), intent(in)  gop,
real(r8), dimension(:), intent(out)  val 
)

Reconstruct the gradient of a H^1 scalar field.

Parameters
[in]cellCell for interpolation
[in]fPosition in cell in logical coord [4]
[in]gopLogical gradient vectors at f [3,4]
[out]valReconstructed field at f [3]

◆ h1_base_pop()

subroutine h1_base_pop ( class(oft_ml_fe_vecspace), intent(inout)  self,
class(oft_vector), intent(inout)  acors,
class(oft_vector), intent(inout)  afine 
)

Transfer a base level H^1 scalar field to the next MPI level.

Parameters
[in,out]acorsVector to transfer
[in,out]afineFine vector from transfer

◆ h1_base_push()

subroutine h1_base_push ( class(oft_ml_fe_vecspace), intent(inout)  self,
class(oft_vector), intent(inout)  afine,
class(oft_vector), intent(inout)  acors 
)

Transfer a MPI level H^1 scalar field to the base level.

Parameters
[in,out]afineVector to transfer
[in,out]acorsFine vector from transfer

◆ h1_getlop_pre()

subroutine h1_getlop_pre ( type(oft_ml_fem_type), intent(inout), target  ml_h1_rep,
class(oft_solver), intent(out), pointer  pre,
type(oft_matrix_ptr), dimension(:), intent(inout), pointer  mats,
character(len=*), intent(in)  bc_type,
integer(i4), intent(in), optional  level,
integer(i4), intent(in), optional  nlevels 
)

Compute eigenvalues and smoothing coefficients for the operator H^1::LOP.

Parameters
[in]bc_typeBoundary condition

◆ h1_lop_eigs()

subroutine h1_lop_eigs ( type(oft_ml_fem_type), intent(inout), target  ml_h1_rep,
integer(i4), intent(in)  minlev 
)

Compute eigenvalues and smoothing coefficients for the operator H^1::LOP.

◆ h1_mloptions()

subroutine h1_mloptions

Read-in options for the basic H^1 ML preconditioners.

◆ h1_setup_interp()

subroutine h1_setup_interp ( class(oft_ml_fem_type), intent(inout)  ml_h1_rep)

Construct interpolation matrices for transfer between H^1 finite element spaces.

◆ oft_h1_getlop()

subroutine oft_h1_getlop ( class(oft_afem_type), intent(inout), target  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct laplacian matrix for H^1 scalar representation.

Supported boundary conditions

  • ‘'none’Full matrix -'zerob'Dirichlet for all boundary DOF -'grnd'` Dirichlet for only groundin point
    Parameters
    [in,out]matMatrix object
    [in]bcBoundary condition

◆ oft_h1_getmop()

subroutine oft_h1_getmop ( class(oft_afem_type), intent(inout), target  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct mass matrix for H^1 scalar representation.

Supported boundary conditions

  • ‘'none’Full matrix -'zerob'` Dirichlet for all boundary DOF
    Parameters
    [in,out]matMatrix object
    [in]bcBoundary condition

◆ oft_h1_project()

subroutine oft_h1_project ( class(oft_afem_type), intent(inout), target  fe_rep,
class(fem_interp), intent(inout)  field,
class(oft_vector), intent(inout)  x 
)

Project a scalar field onto a H^1 basis.

Note
This subroutine only performs the integration of the field with test functions for a H^1 basis. To retrieve the correct projection the result must be multiplied by the inverse of H^1::MOP.
Parameters
[in,out]fieldVector field for projection
[in,out]xField projected onto H^1 basis

◆ rinterp_apply()

subroutine rinterp_apply ( class(oft_h1_rinterp), intent(inout)  self,
integer(i4), intent(in)  cell,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(3,4), intent(in)  gop,
real(r8), dimension(:), intent(out)  val 
)

Reconstruct a H^1 scalar field.

Parameters
[in]cellCell for interpolation
[in]fPosition in cell in logical coord [4]
[in]gopLogical gradient vectors at f [3,4]
[out]valReconstructed field at f [1]

◆ rinterp_delete()

subroutine rinterp_delete ( class(oft_h1_rinterp), intent(inout)  self)

Destroy temporary internal storage.

◆ rinterp_setup()

subroutine rinterp_setup ( class(oft_h1_rinterp), intent(inout)  self,
class(oft_afem_type), intent(inout), target  h1_rep 
)

Setup interpolator for H^1 scalar fields.

Fetches local representation used for interpolation from vector object

◆ zerob_apply()

subroutine zerob_apply ( class(oft_h1_zerob), intent(inout)  self,
class(oft_vector), intent(inout)  a 
)

Zero a H^1 scalar field at all boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ zerob_delete()

subroutine zerob_delete ( class(oft_h1_zerob), intent(inout)  self)

Zero a H^1 scalar field at all boundary nodes.

◆ zerogrnd_apply()

subroutine zerogrnd_apply ( class(oft_h1_zerogrnd), intent(inout)  self,
class(oft_vector), intent(inout)  a 
)

Zero a H^1 scalar field at the global grounding node.

Note
The possition of this node is defined by the mesh pointer igrnd in mesh
Parameters
[in,out]aField to be zeroed

◆ zeroi_apply()

subroutine zeroi_apply ( class(oft_h1_zeroi), intent(inout)  self,
class(oft_vector), intent(inout)  a 
)

Zero a H^1 scalar field at all interior nodes.

Parameters
[in,out]aField to be zeroed

Variable Documentation

◆ df_lop

real(r8), dimension(fem_max_levels) df_lop =-1.d99

◆ nu_lop

integer(i4), dimension(fem_max_levels) nu_lop =0

◆ oft_h1_gop

real(r8), dimension(:,:), pointer oft_h1_gop => NULL()

◆ oft_h1_rop

real(r8), dimension(:), pointer oft_h1_rop => NULL()