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

Detailed Description

Nedelec H1 FE operator definitions.

Authors
Chris Hansen
Date
August 2011

Data Types

type  oft_h1_cinterp
 Interpolate \( \nabla \times \) of a H1 vector field. More...
 
type  oft_h1_dinterp
 Interpolate \( \nabla \times \) of a H1 vector field. More...
 
type  oft_h1_divout
 Clean the divergence from a H1 vector field. More...
 
type  oft_h1_rinterp
 Interpolate a H1 vector field. More...
 
type  oft_h1_zerograd
 Orthogonalize a H1 vector field by zeroing the gradient subspace. More...
 

Functions/Subroutines

subroutine h1_base_pop (acors, afine)
 Transfer a base level Lagrange scalar field to the next MPI level.
 
subroutine h1_base_push (afine, acors)
 Transfer a MPI level Lagrange scalar field to the base level.
 
subroutine h1_bmc (a, hcpc, hcpv, new_tol)
 Compute the 0-th order gradient due to a jump plane.
 
subroutine h1_cinterp (self, cell, f, gop, val)
 Reconstruct \( \nabla \times \) of a Nedelec H1 vector field.
 
subroutine h1_curltp (a, b)
 Apply the curl transpose operator to a H1 field.
 
subroutine h1_dinterp (self, cell, f, gop, val)
 Reconstruct \( \nabla \cdot \) of a Nedelec H1 vector field.
 
subroutine h1_div (a, b)
 Apply the divergence operator to a H1 field.
 
subroutine h1_divout_apply (self, u)
 Remove divergence from a H1 vector field by adding a gradient correction.
 
subroutine h1_divout_delete (self)
 Clean-up internal storage for a oft_h1_divout object.
 
subroutine h1_divout_setup (self, bc)
 Setup matrix and solver with default.
 
subroutine h1_getmop (mat, bc)
 Construct mass matrix for a H1 representation.
 
subroutine h1_getmop_pre (pre, mats, level, nlevels)
 Compute eigenvalues and smoothing coefficients for the operator H1::MOP.
 
subroutine h1_grad (a, b, keep_boundary)
 Add the gradient of a H0 scalar field to a H1 vector field.
 
subroutine h1_gradtp (a, b)
 Apply the transposed gradient operator to a Nedelec H1 vector field.
 
subroutine h1_inject (afine, acors)
 Inject a fine level Lagrange scalar field to the next coarsest level.
 
subroutine h1_interp (acors, afine)
 Interpolate a coarse level Lagrange scalar field to the next finest level.
 
real(r8) function h1_jump_error (u, quad_order)
 Evaluate the jump error in a field over internal faces.
 
subroutine h1_mc (a, hcpc, hcpv, new_tol)
 Compute the 0-th order gradient due to a jump plane.
 
subroutine h1_mloptions
 Read-in options for the basic Nedelec H1(Curl) ML preconditioners.
 
subroutine h1_mop_eigs (minlev)
 Compute eigenvalues and smoothing coefficients for the operator H1::MOP.
 
subroutine h1_rinterp (self, cell, f, gop, val)
 Reconstruct a Nedelec H1 vector field.
 
subroutine h1_rinterp_delete (self)
 Destroy temporary internal storage.
 
subroutine h1_rinterp_setup (self)
 Setup interpolator for H1 vector fields.
 
subroutine h1_setup_interp (create_full)
 Construct interpolation matrices on each MG level.
 
subroutine h1_zerob (a)
 Zero a Nedelec H1 vector field at all boundary nodes.
 
subroutine h1_zerograd_apply (self, u)
 
subroutine h1_zerograd_delete (self)
 
subroutine h1_zeroi (a)
 Zero a Nedelec H1 vector field at all interior nodes.
 
subroutine h1curl_zerob (a)
 Zero the curl components of a Nedelec H1 vector field at all boundary nodes.
 
subroutine h1grad_zerop (a)
 Zero the gradient components of a Nedelec H1 vector field at all boundary nodes.
 
subroutine hgrad_ginterpmatrix (mat)
 Construct interpolation matrix for polynomial levels.
 
subroutine oft_h1_bproject (field, x)
 Boundary projection of a vector field onto a H1 basis.
 
subroutine oft_h1_project (field, x)
 Project a vector field onto a H1 basis.
 

Variables

real(r8), dimension(fem_max_levels), private df_mop =-1.d99
 
integer(i4), dimension(fem_max_levels), private nu_mop =0
 

Function/Subroutine Documentation

◆ h1_base_pop()

subroutine h1_base_pop ( 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]acorsVector to transfer
[in,out]afineFine vector from transfer

◆ h1_base_push()

subroutine h1_base_push ( 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]afineVector to transfer
[in,out]acorsFine vector from transfer

◆ h1_bmc()

subroutine h1_bmc ( class(oft_vector), intent(inout)  a,
real(r8), dimension(3), intent(in)  hcpc,
real(r8), dimension(3), intent(in)  hcpv,
real(r8), intent(in), optional  new_tol 
)

Compute the 0-th order gradient due to a jump plane.

The jump is represented as a circular surface defined by a center possition and surface normal. This method requires that the mesh contain a matching internal surface, such that no edge crosses the jump plane.

Note
The radius of the surface is represented by \( \left| hcpv \right| \)
Parameters
[in,out]aJump field
[in]hcpcJump plane center possition [3]
[in]hcpvJump plane normal vector [3]

◆ h1_cinterp()

subroutine h1_cinterp ( class(oft_h1_cinterp), 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 \times \) of a Nedelec H1 vector 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 [3]

◆ h1_curltp()

subroutine h1_curltp ( class(oft_vector), intent(inout)  a,
class(oft_vector), intent(inout)  b 
)

Apply the curl transpose operator to a H1 field.

◆ h1_dinterp()

subroutine h1_dinterp ( class(oft_h1_dinterp), 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 \cdot \) of a Nedelec H1 vector 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]

◆ h1_div()

subroutine h1_div ( class(oft_vector), intent(inout)  a,
class(oft_vector), intent(inout)  b 
)

Apply the divergence operator to a H1 field.

◆ h1_divout_apply()

subroutine h1_divout_apply ( class(oft_h1_divout), intent(inout)  self,
class(oft_vector), intent(inout)  u 
)

Remove divergence from a H1 vector field by adding a gradient correction.

Note
Should only be used via class oft_h1_divout
Parameters
[in,out]uField for divergence cleaning

◆ h1_divout_delete()

subroutine h1_divout_delete ( class(oft_h1_divout), intent(inout)  self)

Clean-up internal storage for a oft_h1_divout object.

Note
Should only be used via class oft_h1_divout

◆ h1_divout_setup()

subroutine h1_divout_setup ( class(oft_h1_divout), intent(inout)  self,
character(len=*), intent(in)  bc 
)

Setup matrix and solver with default.

Note
Should only be used via class oft_h1_divout
Parameters
[in]bcBoundary condition

◆ h1_getmop()

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

Construct mass matrix for a H1 representation.

Supported boundary conditions

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

◆ h1_getmop_pre()

subroutine h1_getmop_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 H1::MOP.

◆ h1_grad()

subroutine h1_grad ( class(oft_vector), intent(inout)  a,
class(oft_vector), intent(inout)  b,
logical, intent(in), optional  keep_boundary 
)

Add the gradient of a H0 scalar field to a H1 vector field.

Note
By default the 0-th order gradient subspace is represented on the H1(Curl) DOF, use the keep_boundary flag otherwise.
Parameters
[in,out]aScalar field
[in,out]bVector field for gradient
[in]keep_boundaryFlag to keep 0-th order boundary component (optional)

◆ h1_gradtp()

subroutine h1_gradtp ( class(oft_vector), intent(inout)  a,
class(oft_vector), intent(inout)  b 
)

Apply the transposed gradient operator to a Nedelec H1 vector field.

Parameters
[in,out]aInput field
[in,out]b\( G^{T} a \)

◆ h1_inject()

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

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

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

◆ h1_interp()

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

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

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

◆ h1_jump_error()

real(r8) function h1_jump_error ( class(oft_vector), intent(inout)  u,
integer(i4), intent(in)  quad_order 
)

Evaluate the jump error in a field over internal faces.

Note
Currently faces on domain boundaries are skipped, this is due to the fact that evaluting the error would require costly communication.
Parameters
[in,out]uH1 vector field to evaluate
[in]quad_orderDesired quadrature order for integration
Returns
Jump error metric

◆ h1_mc()

subroutine h1_mc ( class(oft_vector), intent(inout)  a,
real(r8), dimension(3), intent(in)  hcpc,
real(r8), dimension(3), intent(in)  hcpv,
real(r8), intent(in), optional  new_tol 
)

Compute the 0-th order gradient due to a jump plane.

The jump is represented as a circular surface defined by a center possition and surface normal. This method requires that the mesh contain a matching internal surface, such that no edge crosses the jump plane.

Note
The radius of the surface is represented by \( \left| hcpv \right| \)
Parameters
[in,out]aJump field
[in]hcpcJump plane center possition [3]
[in]hcpvJump plane normal vector [3]

◆ h1_mloptions()

subroutine h1_mloptions

Read-in options for the basic Nedelec H1(Curl) ML preconditioners.

◆ h1_mop_eigs()

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

Compute eigenvalues and smoothing coefficients for the operator H1::MOP.

◆ h1_rinterp()

subroutine h1_rinterp ( 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 Nedelec H1 vector 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 [3]

◆ h1_rinterp_delete()

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

Destroy temporary internal storage.

Note
Should only be used via class oft_h1_rinterp or children

◆ h1_rinterp_setup()

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

Setup interpolator for H1 vector fields.

Fetches local representation used for interpolation from vector object

Note
Should only be used via class oft_h1_rinterp or children

◆ h1_setup_interp()

subroutine h1_setup_interp ( logical, intent(in), optional  create_full)

Construct interpolation matrices on each MG level.

◆ h1_zerob()

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

Zero a Nedelec H1 vector field at all boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ h1_zerograd_apply()

subroutine h1_zerograd_apply ( class(oft_h1_zerograd), intent(inout)  self,
class(oft_vector), intent(inout)  u 
)

◆ h1_zerograd_delete()

subroutine h1_zerograd_delete ( class(oft_h1_zerograd), intent(inout)  self)

◆ h1_zeroi()

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

Zero a Nedelec H1 vector field at all interior nodes.

Parameters
[in,out]aField to be zeroed

◆ h1curl_zerob()

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

Zero the curl components of a Nedelec H1 vector field at all boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ h1grad_zerop()

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

Zero the gradient components of a Nedelec H1 vector field at all boundary nodes.

Parameters
[in,out]aField to be zeroed

◆ hgrad_ginterpmatrix()

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

Construct interpolation matrix for polynomial levels.

◆ oft_h1_bproject()

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

Boundary projection of a vector field onto a H1 basis.

Note
This subroutine only performs the integration of the field with boundary test functions for a H1 basis.
Parameters
[in,out]fieldVector field for projection
[in,out]xField projected onto H1 basis

◆ oft_h1_project()

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

Project a vector field onto a H1 basis.

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

Variable Documentation

◆ df_mop

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

◆ nu_mop

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