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

Detailed Description

H(Curl) FE operator definitions.

Authors
Chris Hansen
Date
August 2011

Data Types

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

Functions/Subroutines

subroutine hcurl_base_pop (self, acors, afine)
 Transfer a base level H(Curl) vector field to the next MPI level.
 
subroutine hcurl_base_push (self, afine, acors)
 Transfer a MPI level H(Curl) vector field to the base level.
 
subroutine hcurl_cinterp (self, cell, f, gop, val)
 Reconstruct the curl of a H(Curl) vector field.
 
subroutine hcurl_div (hcurl_fe, lag_fe, a, b)
 Apply the divergence operator to a H(Curl) field.
 
subroutine hcurl_divout_apply (self, a)
 Remove divergence from a H(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, ml_hcurl_rep, ml_lag_rep, bc, solver)
 Setup matrix and solver with default.
 
subroutine hcurl_getjmlb_pre (ml_hcurl_rep, pre, mats, alam, level, nlevels)
 Construct default MG preconditioner for H(Curl)::JMLB.
 
subroutine hcurl_getwop_pre (ml_hcurl_rep, pre, mats, level, nlevels)
 Construct default MG preconditioner for H(Curl)::WOP.
 
subroutine hcurl_grad (hcurl_fe, a, b)
 Add the gradient of a linear Lagrange scalar field to a H(Curl) field.
 
subroutine hcurl_gradtp (hcurl_fe, a, b)
 Apply the transposed gradient operator to a H(Curl) vector field.
 
subroutine hcurl_mloptions (ml_hcurl_obj)
 Read-in options for the basic H(Curl) ML preconditioners.
 
subroutine hcurl_orthog_apply (self, a)
 Orthogonalize a H(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_rinterp (self, cell, f, gop, val)
 Reconstruct a H(Curl) vector field.
 
subroutine hcurl_rinterp_delete (self)
 Destroy temporary internal storage.
 
subroutine hcurl_rinterp_setup (self, hcurl_rep)
 Setup interpolator for H(Curl) fields.
 
subroutine hcurl_setup_interp (ml_hcurl_rep)
 Construct interpolation matrices on each MG level.
 
subroutine hcurl_wop_eigs (ml_hcurl_rep, minlev)
 Compute eigenvalues and smoothing coefficients for the operator H(Curl)::WOP.
 
subroutine oft_hcurl_bcurl (fe_rep, bfe_rep, field, x)
 Compute the boundary term for integration by parts of the curl operator using a HCurl basis.
 
subroutine oft_hcurl_getjmlb (fe_rep, mat, alam, bc)
 Construct force-free response matrix for a H(Curl) representation.
 
subroutine oft_hcurl_getkop (fe_rep, mat, bc)
 Construct helicity matrix for a H(Curl) representation.
 
subroutine oft_hcurl_getmop (fe_rep, mat, bc)
 Construct mass matrix for a H(Curl) representation.
 
subroutine oft_hcurl_getwop (fe_rep, mat, bc)
 Construct energy matrix for a H(Curl) representation.
 
subroutine oft_hcurl_project (fe_rep, field, x)
 Project a vector field onto a H(Curl) basis.
 
subroutine zerob_apply (self, a)
 Zero the tangential component of a H(Curl) vector field on the boundary.
 
subroutine zerob_delete (self)
 Zero the tangential component of a H(Curl) vector field on the boundary.
 

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_ml_fe_vecspace), intent(inout)  self,
class(oft_vector), intent(inout)  acors,
class(oft_vector), intent(inout)  afine 
)

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

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

◆ hcurl_base_push()

subroutine hcurl_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(Curl) vector field to the base level.

Parameters
[in,out]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 H(Curl) 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]

◆ hcurl_div()

subroutine hcurl_div ( class(oft_afem_type), intent(inout)  hcurl_fe,
class(oft_afem_type), intent(inout)  lag_fe,
class(oft_vector), intent(inout)  a,
class(oft_vector), intent(inout)  b 
)

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

Note
This subroutine computes the divergence of a H(Curl) field as projected on to a linear Lagrange scalar basis.
Parameters
[in,out]aNeeds docs
[in,out]bNeeds docs

◆ hcurl_divout_apply()

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

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

Parameters
[in,out]aField 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.

◆ hcurl_divout_setup()

subroutine hcurl_divout_setup ( class(oft_hcurl_divout), intent(inout)  self,
class(oft_ml_fem_type), intent(inout), target  ml_hcurl_rep,
class(oft_ml_fem_type), intent(inout), target  ml_lag_rep,
character(len=*), intent(in)  bc,
class(oft_solver), intent(in), optional, target  solver 
)

Setup matrix and solver with default.

Parameters
[in]bcBoundary condition

◆ hcurl_getjmlb_pre()

subroutine hcurl_getjmlb_pre ( type(oft_ml_fem_type), intent(inout), target  ml_hcurl_rep,
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 H(Curl)::JMLB.

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

◆ hcurl_getwop_pre()

subroutine hcurl_getwop_pre ( type(oft_ml_fem_type), intent(inout), target  ml_hcurl_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 H(Curl)::WOP.

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

◆ hcurl_grad()

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

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

Parameters
[in,out]aNeeds docs
[in,out]bNeeds docs

◆ hcurl_gradtp()

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

Apply the transposed gradient operator to a H(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_mloptions()

subroutine hcurl_mloptions ( class(oft_ml_fem_type), intent(inout)  ml_hcurl_obj)

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

◆ hcurl_orthog_apply()

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

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

Note
Used as a member function of oft_hcurl_orthog only
Parameters
[in,out]aField 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_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 H(Curl) 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]

◆ hcurl_rinterp_delete()

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

Destroy temporary internal storage.

◆ hcurl_rinterp_setup()

subroutine hcurl_rinterp_setup ( class(oft_hcurl_rinterp), intent(inout)  self,
class(oft_afem_type), intent(inout), target  hcurl_rep 
)

Setup interpolator for H(Curl) fields.

Fetches local representation used for interpolation from vector object

◆ hcurl_setup_interp()

subroutine hcurl_setup_interp ( class(oft_ml_fem_type), intent(inout)  ml_hcurl_rep)

Construct interpolation matrices on each MG level.

◆ hcurl_wop_eigs()

subroutine hcurl_wop_eigs ( type(oft_ml_fem_type), intent(inout), target  ml_hcurl_rep,
integer(i4), intent(in)  minlev 
)

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

Parameters
[in]minlevNeeds docs

◆ oft_hcurl_bcurl()

subroutine oft_hcurl_bcurl ( class(oft_afem_type), intent(inout), target  fe_rep,
class(oft_afem_type), intent(inout), target  bfe_rep,
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_afem_type), intent(inout), target  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
real(r8), intent(in)  alam,
character(len=*), intent(in)  bc 
)

Construct force-free response matrix for a H(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_afem_type), intent(inout), target  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct helicity matrix for a H(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_afem_type), intent(inout), target  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct mass matrix for a H(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_afem_type), intent(inout), target  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct energy matrix for a H(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(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 H(Curl) basis.

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

◆ zerob_apply()

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

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

Parameters
[in,out]aField to be zeroed

◆ zerob_delete()

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

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

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