The Open FUSION Toolkit 1.0.0-beta6
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
hcurl_grad_operators.F90 File Reference
#include "local.h"

Data Types

type  oft_hcurl_grad_cinterp
 Interpolate \( \nabla \times \) of a H(Curl) + Grad(H^1) vector field. More...
type  oft_hcurl_grad_czerob
 Needs docs. More...
type  oft_hcurl_grad_dinterp
 Interpolate \( \nabla \times \) of a H(Curl) + Grad(H^1) vector field. More...
type  oft_hcurl_grad_divout
 Clean the divergence from a H(Curl) + Grad(H^1) vector field. More...
type  oft_hcurl_grad_gzero
 Orthogonalize a H(Curl) + Grad(H^1) vector field by zeroing the gradient subspace. More...
type  oft_hcurl_grad_gzerop
 Needs docs. More...
type  oft_hcurl_grad_rinterp
 Interpolate a H(Curl) + Grad(H^1) vector field. More...
type  oft_hcurl_grad_zerob
 Needs docs. More...
type  oft_hcurl_grad_zeroi
 Needs docs. More...
type  oft_ml_hcurl_grad_vecspace
 Needs docs. More...

Modules

module  oft_hcurl_grad_operators
 FE operators for full H(Curl) + Grad(H^1) space.

Functions/Subroutines

subroutine base_pop (self, acors, afine)
 Transfer a base level Lagrange scalar field to the next MPI level.
subroutine base_push (self, afine, acors)
 Transfer a MPI level Lagrange scalar field to the base level.
subroutine cinterp_apply (self, cell, f, gop, val)
 Reconstruct \( \nabla \times \) of a H(Curl) + Grad(H^1) vector field.
subroutine curl_zerob_apply (self, a)
 Zero the curl components of a H(Curl) + Grad(H^1) vector field at all boundary nodes.
subroutine dinterp_apply (self, cell, f, gop, val)
 Reconstruct \( \nabla \cdot \) of a H(Curl) + Grad(H^1) vector field.
subroutine divout_apply (self, a)
 Remove divergence from a H(Curl) + Grad(H^1) vector field by adding a gradient correction.
subroutine divout_delete (self)
 Clean-up internal storage for a oft_hcurl_grad_divout object.
subroutine divout_setup (self, ml_hcurl_grad_rep, bc, solver)
 Setup matrix and solver with default.
subroutine grad_zerop_apply (self, a)
 Zero the gradient components of a H(Curl) + Grad(H^1) vector field at all boundary nodes.
subroutine hcurl_grad_bmc (mesh, a, hcpc, hcpv, new_tol)
 Compute the 0-th order gradient due to a jump plane.
subroutine hcurl_grad_curltp (hcurl_grad_fe, a, b)
 Apply the curl transpose operator to a H(Curl) + Grad(H^1) field.
subroutine hcurl_grad_div (hcurl_grad_fe, a, b)
 Apply the divergence operator to a H(Curl) + Grad(H^1) field.
subroutine hcurl_grad_getmop (hcurl_grad_rep, mat, bc)
 Construct mass matrix for a H(Curl) + Grad(H^1) representation.
subroutine hcurl_grad_getmop_pre (ml_hcurl_aug_obj, pre, mats, level, nlevels)
 Compute eigenvalues and smoothing coefficients for the operator H(Curl) + Grad(H^1) mass matrix.
subroutine hcurl_grad_grad (hcurl_grad_fe, a, b, keep_boundary)
 Add the gradient of a H^1 scalar field to a H(Curl) + Grad(H^1) vector field.
subroutine hcurl_grad_gradtp (h1_fe, a, b)
 Apply the transposed gradient operator to a H(Curl) + Grad(H^1) vector field.
real(r8) function hcurl_grad_jump_error (hcurl_grad_fe, u, quad_order)
 Evaluate the jump error in a field over internal faces.
subroutine hcurl_grad_mc (mesh, a, hcpc, hcpv, new_tol)
 Compute the 0-th order gradient due to a jump plane.
subroutine hcurl_grad_mloptions
 Read-in options for the basic H(Curl) + Grad(H^1) space ML preconditioners.
subroutine hcurl_grad_mop_eigs (ml_hcurl_aug_obj, minlev)
 Compute eigenvalues and smoothing coefficients for the mass matrix.
subroutine hcurl_grad_setup_interp (ml_hcurl_aug_rep, ml_h1_rep, create_full)
 Construct interpolation matrices on each MG level.
subroutine hgrad_ginterpmatrix (mat)
 Construct interpolation matrix for polynomial levels.
subroutine ml_vecspace_inject (self, afine, acors)
 Interpolate a coarse level Lagrange scalar field to the next finest level.
subroutine ml_vecspace_interp (self, acors, afine)
 Interpolate a coarse level Lagrange scalar field to the next finest level.
subroutine oft_hcurl_grad_bproject (hcurl_grad_rep, field, x)
 Boundary projection of a vector field onto a H(Curl) + Grad(H^1) basis.
subroutine oft_hcurl_grad_project (hcurl_grad_rep, field, x)
 Project a vector field onto a H(Curl) + Grad(H^1) basis.
subroutine rinterp_apply (self, cell, f, gop, val)
 Reconstruct a H(Curl) + Grad(H^1) vector field.
subroutine rinterp_delete (self)
 Destroy temporary internal storage.
subroutine rinterp_setup1 (self, hcurl_grad_rep)
 Setup interpolator for H(Curl) + Grad(H^1) vector fields.
subroutine rinterp_setup2 (self, hcurl_rep, hgrad_rep)
 Setup interpolator for H(Curl) + Grad(H^1) vector fields.
subroutine zerob_apply (self, a)
 Zero a H(Curl) + Grad(H^1) vector field at all boundary nodes.
subroutine zerob_delete (self)
 Zero a H(Curl) + Grad(H^1) vector field at all boundary nodes.
subroutine zerograd_apply (self, a)
 Needs docs.
subroutine zerograd_delete (self)
 Needs docs.
subroutine zeroi_apply (self, a)
 Zero a H(Curl) + Grad(H^1) vector field at all interior nodes.

Variables

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

Function/Subroutine Documentation

◆ hgrad_ginterpmatrix()

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

Construct interpolation matrix for polynomial levels.

Parameters
[in,out]matInterpolation matrix