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

Detailed Description

Nedelec H1(Curl) FE operator definitions.

Authors
Chris Hansen
Date
August 2011

Data Types

type  oft_hcurl_cinterp
 Interpolate \( \nabla \times \) of a H1(Curl) field. More...
 
type  oft_hcurl_divout
 Clean the divergence from a H1(Curl) vector field. More...
 
type  oft_hcurl_orthog
 Orthogonalize a H1(Curl) vector against a library of given modes. More...
 
type  oft_hcurl_rinterp
 Interpolate a H1(Curl) field. More...
 

Functions/Subroutines

subroutine hcurl_base_pop (acors, afine)
 Transfer a base level H1(Curl) vector field to the next MPI level.
 
subroutine hcurl_base_push (afine, acors)
 Transfer a MPI level H1(Curl) vector field to the base level.
 
subroutine hcurl_cinterp (self, cell, f, gop, val)
 Reconstruct the curl of a Nedelec H1(Curl) vector field.
 
subroutine hcurl_div (a, b)
 Apply the divergence operator to a H1(Curl) field.
 
subroutine hcurl_divout_apply (self, u)
 Remove divergence from a H1(Curl) vector field by adding a gradient correction.
 
subroutine hcurl_divout_delete (self)
 Clean-up internal storage for a oft_hcurl_divout object.
 
subroutine hcurl_divout_setup (self, bc)
 Setup matrix and solver with default.
 
subroutine hcurl_getjmlb_pre (pre, mats, alam, level, nlevels)
 Construct default MG preconditioner for H1(Curl)::JMLB.
 
subroutine hcurl_getwop_pre (pre, mats, level, nlevels)
 Construct default MG preconditioner for H1(Curl)::WOP.
 
subroutine hcurl_ginterpmatrix (mat)
 Construct interpolation matrix for polynomial levels.
 
subroutine hcurl_grad (a, b)
 Add the gradient of a linear Lagrange scalar field to a H1(Curl) field.
 
subroutine hcurl_gradtp (a, b)
 Apply the transposed gradient operator to a Nedelec H1(Curl) vector field.
 
subroutine hcurl_inject (afine, acors)
 Inject a fine level H1(Curl) vector field to the next coarsest level.
 
subroutine hcurl_interp (acors, afine)
 Interpolate a coarse level H1(Curl) vector field to the next finest level.
 
subroutine hcurl_mloptions
 Read-in options for the basic Nedelec H1(Curl) ML preconditioners.
 
subroutine hcurl_orthog_apply (self, u)
 Orthogonalize a H1(Curl) vector against a library of given modes.
 
subroutine hcurl_orthog_delete (self)
 Clean-up internal storage for a oft_hcurl_orthog object.
 
subroutine hcurl_pinterpmatrix (mat)
 Construct interpolation matrix for polynomial levels.
 
subroutine hcurl_rinterp (self, cell, f, gop, val)
 Reconstruct a Nedelec H1(Curl) vector field.
 
subroutine hcurl_rinterp_delete (self)
 Destroy temporary internal storage.
 
subroutine hcurl_rinterp_setup (self)
 Setup interpolator for H1(Curl) fields.
 
subroutine hcurl_setup_interp
 Construct interpolation matrices on each MG level.
 
subroutine hcurl_wop_eigs (minlev)
 Compute eigenvalues and smoothing coefficients for the operator H1(Curl)::WOP.
 
subroutine hcurl_zerob (a)
 Zero the tangential component of a Nedelec H1(Curl) vector field on the boundary.
 
subroutine oft_hcurl_bcurl (field, x)
 Compute the boundary term for integration by parts of the curl operator using a HCurl basis.
 
subroutine oft_hcurl_getjmlb (mat, alam, bc)
 Construct force-free response matrix for a H1(Curl) representation.
 
subroutine oft_hcurl_getkop (mat, bc)
 Construct helicity matrix for a H1(Curl) representation.
 
subroutine oft_hcurl_getmop (mat, bc)
 Construct mass matrix for a H1(Curl) representation.
 
subroutine oft_hcurl_getwop (mat, bc)
 Construct energy matrix for a H1(Curl) representation.
 
subroutine oft_hcurl_project (field, x)
 Project a vector field onto a H1(Curl) basis.
 

Variables

real(r8), dimension(fem_max_levels) df_wop =-1.d99
 
integer(i4), dimension(fem_max_levels) nu_jmlb =20
 
integer(i4), dimension(fem_max_levels) nu_wop =0
 
real(r8), dimension(:,:), pointer oft_hcurl_cop => NULL()
 
real(r8), dimension(:,:), pointer oft_hcurl_rop => NULL()
 

Function/Subroutine Documentation

◆ hcurl_base_pop()

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

Transfer a base level H1(Curl) vector field to the next MPI level.

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

◆ hcurl_base_push()

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

Transfer a MPI level H1(Curl) vector field to the base level.

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

◆ hcurl_cinterp()

subroutine hcurl_cinterp ( class(oft_hcurl_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 the curl of a Nedelec H1(Curl) 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 gradient at f [3]

◆ hcurl_div()

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

Apply the divergence operator to a H1(Curl) field.

Note
This subroutine computes the divergence of a H1(Curl) field as projected on to a linear Lagrange scalar basis.

◆ hcurl_divout_apply()

subroutine hcurl_divout_apply ( class(oft_hcurl_divout), intent(inout)  self,
class(oft_vector), intent(inout)  u 
)

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

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

◆ hcurl_divout_delete()

subroutine hcurl_divout_delete ( class(oft_hcurl_divout), intent(inout)  self)

Clean-up internal storage for a oft_hcurl_divout object.

Note
Should only be used via class oft_hcurl_divout

◆ hcurl_divout_setup()

subroutine hcurl_divout_setup ( class(oft_hcurl_divout), intent(inout)  self,
character(len=*), intent(in)  bc 
)

Setup matrix and solver with default.

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

◆ hcurl_getjmlb_pre()

subroutine hcurl_getjmlb_pre ( class(oft_solver), intent(out), pointer  pre,
type(oft_matrix_ptr), dimension(:), intent(inout), pointer  mats,
real(r8), intent(in)  alam,
integer(i4), intent(in), optional  level,
integer(i4), intent(in), optional  nlevels 
)

Construct default MG preconditioner for H1(Curl)::JMLB.

◆ hcurl_getwop_pre()

subroutine hcurl_getwop_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 
)

Construct default MG preconditioner for H1(Curl)::WOP.

◆ hcurl_ginterpmatrix()

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

Construct interpolation matrix for polynomial levels.

◆ hcurl_grad()

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

Add the gradient of a linear Lagrange scalar field to a H1(Curl) field.

◆ hcurl_gradtp()

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

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

Note
Only the 0th order component is computed from the discrete gradient operator
Parameters
[in,out]aInput field
[in,out]b\( G^{T} a \)

◆ hcurl_inject()

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

Inject a fine level H1(Curl) vector field to the next coarsest level.

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

◆ hcurl_interp()

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

Interpolate a coarse level H1(Curl) vector field to the next finest level.

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

◆ hcurl_mloptions()

subroutine hcurl_mloptions

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

◆ hcurl_orthog_apply()

subroutine hcurl_orthog_apply ( class(oft_hcurl_orthog), intent(inout)  self,
class(oft_vector), intent(inout)  u 
)

Orthogonalize a H1(Curl) vector against a library of given modes.

Note
Used as a member function of oft_hcurl_orthog only
Parameters
[in,out]uField to orthogonalize

◆ hcurl_orthog_delete()

subroutine hcurl_orthog_delete ( class(oft_hcurl_orthog), intent(inout)  self)

Clean-up internal storage for a oft_hcurl_orthog object.

Note
Used as a member function of oft_hcurl_orthog only

◆ hcurl_pinterpmatrix()

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

Construct interpolation matrix for polynomial levels.

◆ hcurl_rinterp()

subroutine hcurl_rinterp ( class(oft_hcurl_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(Curl) 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]

◆ hcurl_rinterp_delete()

subroutine hcurl_rinterp_delete ( class(oft_hcurl_rinterp), intent(inout)  self)

Destroy temporary internal storage.

Note
Should only be used via class oft_hcurl_rinterp or children

◆ hcurl_rinterp_setup()

subroutine hcurl_rinterp_setup ( class(oft_hcurl_rinterp), intent(inout)  self)

Setup interpolator for H1(Curl) fields.

Fetches local representation used for interpolation from vector object

Note
Should only be used via class oft_hcurl_rinterp or children

◆ hcurl_setup_interp()

subroutine hcurl_setup_interp

Construct interpolation matrices on each MG level.

◆ hcurl_wop_eigs()

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

Compute eigenvalues and smoothing coefficients for the operator H1(Curl)::WOP.

◆ hcurl_zerob()

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

Zero the tangential component of a Nedelec H1(Curl) vector field on the boundary.

Parameters
[in,out]aField to be zeroed

◆ oft_hcurl_bcurl()

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

Compute the boundary term for integration by parts of the curl operator using a HCurl basis.

Parameters
[in,out]fieldVector field for projection
[in,out]x\( \int \left( \textbf{u}^T \times \textbf{F} \right) \cdot \textbf{dS} \)

◆ oft_hcurl_getjmlb()

subroutine oft_hcurl_getjmlb ( class(oft_matrix), intent(inout), pointer  mat,
real(r8), intent(in)  alam,
character(len=*), intent(in)  bc 
)

Construct force-free response matrix for a H1(Curl) representation.

Supported boundary conditions

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

◆ oft_hcurl_getkop()

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

Construct helicity matrix for a H1(Curl) representation.

Supported boundary conditions

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

◆ oft_hcurl_getmop()

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

Construct mass matrix for a H1(Curl) representation.

Supported boundary conditions

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

◆ oft_hcurl_getwop()

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

Construct energy matrix for a H1(Curl) representation.

Supported boundary conditions

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

◆ oft_hcurl_project()

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

Project a vector field onto a H1(Curl) basis.

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

Variable Documentation

◆ df_wop

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

◆ nu_jmlb

integer(i4), dimension(fem_max_levels) nu_jmlb =20

◆ nu_wop

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

◆ oft_hcurl_cop

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

◆ oft_hcurl_rop

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