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

Detailed Description

Nedelec H0 FE operator definitions.

Authors
Chris Hansen
Date
August 2011

Data Types

type  oft_h0_ginterp
 Interpolate \( \nabla \) of a H0 field. More...
 
type  oft_h0_rinterp
 Interpolate a H0 field. More...
 

Functions/Subroutines

subroutine h0_base_pop (acors, afine)
 Transfer a base level H0 scalar field to the next MPI level.
 
subroutine h0_base_push (afine, acors)
 Transfer a MPI level H0 scalar field to the base level.
 
subroutine h0_getlop_pre (pre, mats, level, nlevels)
 Compute eigenvalues and smoothing coefficients for the operator H0::LOP.
 
subroutine h0_ginterp_apply (self, cell, f, gop, val)
 Reconstruct the gradient of a Nedelec H0 scalar field.
 
subroutine h0_ginterpmatrix (mat)
 Construct interpolation matrix for transfer between geometric levels of H0 finite element space.
 
subroutine h0_inject (afine, acors)
 Inject a fine level H0 scalar field to the next coarsest level.
 
subroutine h0_interp (acors, afine)
 Interpolate a coarse level H0 scalar field to the next finest level.
 
subroutine h0_lop_eigs (minlev)
 Compute eigenvalues and smoothing coefficients for the operator H0::LOP.
 
subroutine h0_mloptions
 Read-in options for the basic Nedelec H0 ML preconditioners.
 
subroutine h0_pinterpmatrix (mat)
 Construct interpolation matrix for transfer between polynomial levels of H0 finite element space.
 
subroutine h0_rinterp (self, cell, f, gop, val)
 Reconstruct a Nedelec H0 scalar field.
 
subroutine h0_rinterp_delete (self)
 Destroy temporary internal storage.
 
subroutine h0_rinterp_setup (self)
 Setup interpolator for H0 scalar fields.
 
subroutine h0_setup_interp
 Construct interpolation matrices for transfer between H0 finite element spaces.
 
subroutine h0_zerob (a)
 Zero a Nedelec H0 scalar field at all boundary nodes.
 
subroutine h0_zerogrnd (a)
 Zero a Nedelec H0 scalar field at the global grounding node.
 
subroutine h0_zeroi (a)
 Zero a Nedelec H0 scalar field at all interior nodes.
 
subroutine oft_h0_getlop (mat, bc)
 Construct laplacian matrix for H0 scalar representation.
 
subroutine oft_h0_getmop (mat, bc)
 Construct mass matrix for H0 scalar representation.
 
subroutine oft_h0_project (field, x)
 Project a scalar field onto a H0 basis.
 

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_h0_gop => NULL()
 
real(r8), dimension(:), pointer oft_h0_rop => NULL()
 

Function/Subroutine Documentation

◆ h0_base_pop()

subroutine h0_base_pop ( class(oft_vector), intent(inout)  acors,
class(oft_vector), intent(inout)  afine 
)

Transfer a base level H0 scalar field to the next MPI level.

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

◆ h0_base_push()

subroutine h0_base_push ( class(oft_vector), intent(inout)  afine,
class(oft_vector), intent(inout)  acors 
)

Transfer a MPI level H0 scalar field to the base level.

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

◆ h0_getlop_pre()

subroutine h0_getlop_pre ( 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 
)

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

◆ h0_ginterp_apply()

subroutine h0_ginterp_apply ( class(oft_h0_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 Nedelec H0 scalar field.

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

◆ h0_ginterpmatrix()

subroutine h0_ginterpmatrix ( class(oft_matrix), intent(inout), pointer  mat)

Construct interpolation matrix for transfer between geometric levels of H0 finite element space.

◆ h0_inject()

subroutine h0_inject ( class(oft_vector), intent(inout)  afine,
class(oft_vector), intent(inout)  acors 
)

Inject a fine level H0 scalar field to the next coarsest level.

Note
The global H0 level in decremented by one in this subroutine
Parameters
[in]afineVector to inject
[in,out]acorsCoarse vector from injection

◆ h0_interp()

subroutine h0_interp ( class(oft_vector), intent(inout)  acors,
class(oft_vector), intent(inout)  afine 
)

Interpolate a coarse level H0 scalar field to the next finest level.

Note
The global H0 level in incremented by one in this subroutine
Parameters
[in]acorsVector to interpolate
[in,out]afineFine vector from interpolation

◆ h0_lop_eigs()

subroutine h0_lop_eigs ( integer(i4), intent(in)  minlev)

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

◆ h0_mloptions()

subroutine h0_mloptions

Read-in options for the basic Nedelec H0 ML preconditioners.

◆ h0_pinterpmatrix()

subroutine h0_pinterpmatrix ( class(oft_matrix), intent(inout), pointer  mat)

Construct interpolation matrix for transfer between polynomial levels of H0 finite element space.

◆ h0_rinterp()

subroutine h0_rinterp ( class(oft_h0_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 Nedelec H0 scalar field.

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

◆ h0_rinterp_delete()

subroutine h0_rinterp_delete ( class(oft_h0_rinterp), intent(inout)  self)

Destroy temporary internal storage.

Note
Should only be used via class oft_h0_rinterp or children

◆ h0_rinterp_setup()

subroutine h0_rinterp_setup ( class(oft_h0_rinterp), intent(inout)  self)

Setup interpolator for H0 scalar fields.

Fetches local representation used for interpolation from vector object

Note
Should only be used via class oft_h0_rinterp or children

◆ h0_setup_interp()

subroutine h0_setup_interp

Construct interpolation matrices for transfer between H0 finite element spaces.

◆ h0_zerob()

subroutine h0_zerob ( class(oft_vector), intent(inout)  a)

Zero a Nedelec H0 scalar field at all boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ h0_zerogrnd()

subroutine h0_zerogrnd ( class(oft_vector), intent(inout)  a)

Zero a Nedelec H0 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

◆ h0_zeroi()

subroutine h0_zeroi ( class(oft_vector), intent(inout)  a)

Zero a Nedelec H0 scalar field at all interior nodes.

Parameters
[in,out]aField to be zeroed

◆ oft_h0_getlop()

subroutine oft_h0_getlop ( class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct laplacian matrix for H0 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_h0_getmop()

subroutine oft_h0_getmop ( class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct mass matrix for H0 scalar representation.

Supported boundary conditions

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

◆ oft_h0_project()

subroutine oft_h0_project ( class(fem_interp), intent(inout)  field,
class(oft_vector), intent(inout)  x 
)

Project a scalar field onto a H0 basis.

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

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_h0_gop

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

◆ oft_h0_rop

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