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_lag_operators Module Reference

Detailed Description

Lagrange FE operator definitions.

Authors
Chris Hansen
Date
August 2011

Data Types

type  oft_lag_ginterp
 Interpolate \( \nabla \) of a Lagrange field. More...
 
type  oft_lag_rinterp
 Interpolate a Lagrange field. More...
 
type  oft_lag_vcinterp
 Interpolate \( \nabla \times \) of a Lagrange vector field. More...
 
type  oft_lag_vdinterp
 Interpolate \( \nabla \) of a Lagrange vector field. More...
 
type  oft_lag_vrinterp
 Interpolate a Lagrange vector field. More...
 
type  oft_lag_zerob
 Zero Lagrange FE vector at all global boundary nodes. More...
 
type  oft_lag_zerogrnd
 Zero Lagrange FE vector at all global grounding node(s) More...
 
type  oft_vlag_zerob
 Zero all components of vector Lagrange FE vector at all global boundary nodes. More...
 
type  oft_vlag_zeron
 Zero normal component of vector Lagrange FE vector at all global boundary nodes. More...
 
type  oft_vlag_zerot
 Zero tangential component of vector Lagrange FE vector at all global boundary nodes. More...
 

Functions/Subroutines

subroutine lag_base_pop (self, acors, afine)
 Transfer a base level Lagrange scalar field to the next MPI level.
 
subroutine lag_base_push (self, afine, acors)
 Transfer a MPI level Lagrange scalar field to the base level.
 
subroutine lag_div (fe_rep, a, reg)
 Compute the divergence of a Lagrange vector field.
 
subroutine lag_getlop_pre (ml_lag_rep, pre, mats, level, nlevels)
 Construct default MG preconditioner for LAG::LOP.
 
subroutine lag_ginterp_apply (self, cell, f, gop, val)
 Reconstruct the gradient of a Lagrange scalar field.
 
subroutine lag_lop_eigs (ml_lag_rep, minlev)
 Compute eigenvalues and smoothing coefficients for the operator LAG::LOP.
 
subroutine lag_mloptions ()
 Read-in options for the basic Lagrange ML preconditioners.
 
subroutine lag_rinterp (self, cell, f, gop, val)
 Reconstruct a Lagrange scalar field.
 
subroutine lag_rinterp_delete (self)
 Destroy temporary internal storage and nullify references.
 
subroutine lag_rinterp_setup (self, lag_rep)
 Setup interpolator for Lagrange scalar fields.
 
subroutine lag_setup_interp (ml_lag_rep, ml_vlag_rep)
 Construct interpolation matrices on each MG level.
 
subroutine lag_vbc_diag (lag_rep, j_lag, bc_type, nn)
 Get diagonal entries for a given boundary condition and desired node in a Lagrange vector field.
 
subroutine lag_vbc_tensor (lag_rep, j_lag, bc_type, nn)
 Get boundary condition tensor for desired node in a Lagrange vector field.
 
subroutine lag_vcinterp (self, cell, f, gop, val)
 Reconstruct the curl of a Lagrange vector field.
 
subroutine lag_vdinterp (self, cell, f, gop, val)
 Reconstruct \( \nabla v \) for a Lagrange vector field.
 
subroutine lag_vrinterp (self, cell, f, gop, val)
 Reconstruct a Lagrange vector field.
 
subroutine lag_vrinterp_delete (self)
 Destroy temporary internal storage and nullify references.
 
subroutine lag_vrinterp_setup (self, lag_rep)
 Setup interpolator for Lagrange vector fields.
 
subroutine oft_lag_getlop (fe_rep, mat, bc)
 Construct laplacian matrix for Lagrange scalar representation.
 
subroutine oft_lag_getmop (fe_rep, mat, bc)
 Construct mass matrix for Lagrange scalar representation.
 
subroutine oft_lag_getpdop (fe_rep, mat, field, bc, perp, be_flag)
 Construct parallel diffusion matrix for Lagrange scalar representation.
 
subroutine oft_lag_project (fe_rep, field, x)
 Project a scalar field onto a lagrange basis.
 
subroutine oft_lag_project_div (fe_rep, field, x)
 Project the divergence of a scalar field onto a lagrange basis.
 
subroutine oft_lag_vgetmop (vlag_rep, mat, bc)
 Construct mass matrix for Lagrange vector representation.
 
subroutine oft_lag_vproject (fe_rep, field, x)
 Project a vector field onto a lagrange basis.
 
subroutine vzerob_apply (self, a)
 Zero a surface Lagrange vector field at all edge nodes.
 
subroutine vzerob_delete (self)
 Destroy temporary internal storage and nullify references.
 
subroutine vzeron_apply (self, a)
 Zero normal component of a Lagrange vector field at boundary nodes.
 
subroutine vzerot_apply (self, a)
 Zero tangential component of a Lagrange vector field at boundary nodes.
 
subroutine zerob_apply (self, a)
 Zero a Lagrange scalar field at all boundary nodes.
 
subroutine zerob_delete (self)
 Destroy temporary internal storage and nullify references.
 
subroutine zerogrnd_apply (self, a)
 Zero a Lagrange scalar field at the global grounding node.
 

Variables

real(r8), dimension(fem_max_levels) df_lop =-1.d99
 
real(r8), dimension(fem_max_levels), private df_pdop =-1.d99
 
integer(i4), dimension(fem_max_levels) nu_lop =0
 
integer(i4), dimension(fem_max_levels), private nu_pdop =0
 
real(r8), dimension(:,:), pointer oft_lag_gop => NULL()
 
real(r8), dimension(:), pointer oft_lag_rop => NULL()
 

Function/Subroutine Documentation

◆ lag_base_pop()

subroutine lag_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 Lagrange scalar field to the next MPI level.

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

◆ lag_base_push()

subroutine lag_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 Lagrange scalar field to the base level.

Parameters
[in,out]afineNeeds docs
[in,out]acorsNeeds docs

◆ lag_div()

subroutine lag_div ( class(oft_afem_type), intent(inout), target  fe_rep,
class(oft_vector), intent(inout)  a,
real(r8), intent(out)  reg 
)

Compute the divergence of a Lagrange vector field.

Parameters
[in,out]aInput field
[out]reg\( \int_v \nabla \cdot a \; dV \)

◆ lag_getlop_pre()

subroutine lag_getlop_pre ( type(oft_ml_fem_type), intent(inout), target  ml_lag_rep,
class(oft_solver), intent(out), pointer  pre,
type(oft_matrix_ptr), dimension(:), intent(inout), pointer  mats,
integer(i4), intent(in), optional  level,
integer(i4), intent(in), optional  nlevels 
)

Construct default MG preconditioner for LAG::LOP.

Parameters
[out]preNeeds docs
[in,out]matsNeeds docs
[in]levelNeeds docs
[in]nlevelsNeeds docs

◆ lag_ginterp_apply()

subroutine lag_ginterp_apply ( class(oft_lag_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 Lagrange 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 gradient at f [3]

◆ lag_lop_eigs()

subroutine lag_lop_eigs ( type(oft_ml_fem_type), intent(inout), target  ml_lag_rep,
integer(i4), intent(in)  minlev 
)

Compute eigenvalues and smoothing coefficients for the operator LAG::LOP.

Parameters
[in]minlevNeeds docs

◆ lag_mloptions()

subroutine lag_mloptions

Read-in options for the basic Lagrange ML preconditioners.

◆ lag_rinterp()

subroutine lag_rinterp ( class(oft_lag_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 Lagrange 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]

◆ lag_rinterp_delete()

subroutine lag_rinterp_delete ( class(oft_lag_rinterp), intent(inout)  self)

Destroy temporary internal storage and nullify references.

◆ lag_rinterp_setup()

subroutine lag_rinterp_setup ( class(oft_lag_rinterp), intent(inout)  self,
class(oft_afem_type), intent(inout), target  lag_rep 
)

Setup interpolator for Lagrange scalar fields.

Fetches local representation used for interpolation from vector object

◆ lag_setup_interp()

subroutine lag_setup_interp ( class(oft_ml_fem_type), intent(inout)  ml_lag_rep,
class(oft_ml_fem_comp_type), intent(inout), optional  ml_vlag_rep 
)

Construct interpolation matrices on each MG level.

◆ lag_vbc_diag()

subroutine lag_vbc_diag ( class(oft_scalar_fem), intent(inout), target  lag_rep,
integer(i4), intent(in)  j_lag,
integer(i4), intent(in)  bc_type,
real(r8), dimension(3,3), intent(out)  nn 
)

Get diagonal entries for a given boundary condition and desired node in a Lagrange vector field.

Parameters
[in]j_lagLocal index of Lagrange node
[in]bc_typeDesired BC (1 -> norm, 2-> tang)
[out]nnDiagonal entries (3,3)

◆ lag_vbc_tensor()

subroutine lag_vbc_tensor ( class(oft_scalar_fem), intent(inout), target  lag_rep,
integer(i4), intent(in)  j_lag,
integer(i4), intent(in)  bc_type,
real(r8), dimension(3,3), intent(out)  nn 
)

Get boundary condition tensor for desired node in a Lagrange vector field.

Parameters
[in]j_lagLocal index of Lagrange node
[in]bc_typeDesired BC (1 -> norm, 2-> tang)
[out]nnBC tensor (3,3)

◆ lag_vcinterp()

subroutine lag_vcinterp ( class(oft_lag_vcinterp), 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 curl of a Lagrange vector 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]

◆ lag_vdinterp()

subroutine lag_vdinterp ( class(oft_lag_vdinterp), 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 \( \nabla v \) for a Lagrange vector field.

The tensor is packed using reshape, to retrieve use

dv = reshape(val,(/3,3/))
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]

◆ lag_vrinterp()

subroutine lag_vrinterp ( class(oft_lag_vrinterp), 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 Lagrange vector 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]

◆ lag_vrinterp_delete()

subroutine lag_vrinterp_delete ( class(oft_lag_vrinterp), intent(inout)  self)

Destroy temporary internal storage and nullify references.

◆ lag_vrinterp_setup()

subroutine lag_vrinterp_setup ( class(oft_lag_vrinterp), intent(inout)  self,
class(oft_afem_type), intent(inout), target  lag_rep 
)

Setup interpolator for Lagrange vector fields.

Fetches local representation used for interpolation from vector object

◆ oft_lag_getlop()

subroutine oft_lag_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 Lagrange 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_lag_getmop()

subroutine oft_lag_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 Lagrange scalar representation.

Supported boundary conditions

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

◆ oft_lag_getpdop()

subroutine oft_lag_getpdop ( class(oft_afem_type), intent(inout), target  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
class(fem_interp), intent(inout)  field,
character(len=*), intent(in)  bc,
real(r8), intent(in), optional  perp,
logical, dimension(:), intent(in), optional  be_flag 
)

Construct parallel diffusion matrix for Lagrange 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,out]fieldVector field defining \( \hat{b} \)
    [in]bcBoundary condition
    [in]perpValue of perpendicular conductivity (optional)
    [in]be_flagFlag for dirichlet nodes if different from boundary [ne] (optional)

◆ oft_lag_project()

subroutine oft_lag_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 lagrange basis.

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

◆ oft_lag_project_div()

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

Project the divergence of a scalar field onto a lagrange basis.

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

◆ oft_lag_vgetmop()

subroutine oft_lag_vgetmop ( class(oft_fem_comp_type), intent(inout), target  vlag_rep,
class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct mass matrix for Lagrange vector representation.

Supported boundary conditions

  • ‘'none’Full matrix -'all'Dirichlet for all components at boundary -'norm'Dirichlet for normal component at boundary -'tang'` Dirichlet for tangential component at boundary
    Parameters
    [in,out]matMatrix object
    [in]bcBoundary condition

◆ oft_lag_vproject()

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

Project a vector field onto a lagrange basis.

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

◆ vzerob_apply()

subroutine vzerob_apply ( class(oft_vlag_zerob), intent(inout)  self,
class(oft_vector), intent(inout)  a 
)

Zero a surface Lagrange vector field at all edge nodes.

Parameters
[in,out]aField to be zeroed

◆ vzerob_delete()

subroutine vzerob_delete ( class(oft_vlag_zerob), intent(inout)  self)

Destroy temporary internal storage and nullify references.

◆ vzeron_apply()

subroutine vzeron_apply ( class(oft_vlag_zeron), intent(inout)  self,
class(oft_vector), intent(inout)  a 
)

Zero normal component of a Lagrange vector field at boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ vzerot_apply()

subroutine vzerot_apply ( class(oft_vlag_zerot), intent(inout)  self,
class(oft_vector), intent(inout)  a 
)

Zero tangential component of a Lagrange vector field at boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ zerob_apply()

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

Zero a Lagrange scalar field at all boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ zerob_delete()

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

Destroy temporary internal storage and nullify references.

◆ zerogrnd_apply()

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

Zero a Lagrange 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

Variable Documentation

◆ df_lop

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

◆ df_pdop

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

◆ nu_lop

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

◆ nu_pdop

integer(i4), dimension(fem_max_levels), private nu_pdop =0
private

◆ oft_lag_gop

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

◆ oft_lag_rop

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