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

Detailed Description

Grad-Shafranov implementation for TokaMaker.

Authors
Chris Hansen
Date
August 2011

Data Types

type  circular_curr
 Interpolation class for uniform source with simple tokamak representation. More...
 
type  coil_region
 Internal coil region structure. More...
 
type  cond_region
 Internal wall region structure. More...
 
interface  flux_cofs_get
 Needs Docs. More...
 
interface  flux_cofs_set
 Needs Docs. More...
 
type  flux_func
 Abstract flux function prototype. More...
 
interface  flux_func_eval
 Needs Docs. More...
 
interface  flux_func_update
 Needs Docs. More...
 
type  gs_b_interp
 Interpolate magnetic field for a G-S solution. More...
 
type  gs_curvature_interp
 Interpolate magnetic field for a G-S solution. More...
 
type  gs_eq
 Grad-Shafranov equilibrium object. More...
 
type  gs_j_interp
 Interpolate magnetic field for a G-S solution. More...
 
type  gs_prof_interp
 Interpolate G-S profiles at a specific point in space. More...
 
type  gs_region_info
 Information for non-continuous regions. More...
 
type  gsinv_interp
 Need docs. More...
 
type  oft_gs_zerob
 Zero scalar Lagrange FE field on all boundary nodes and all nodes outside the plasma. More...
 

Functions/Subroutines

subroutine build_dels (mat, self, bc, dt, scale)
 Construct \( \frac{1}{R} \Delta^* \) operator, with or without time-dependence.
 
subroutine build_mrop (fe_rep, mat, bc)
 Construct mass matrix for a boundary Lagrange scalar representation.
 
subroutine circle_interp (self, cell, f, gop, val)
 Evaluate uniform source over simple plasma cross-section.
 
subroutine compute_bcmat (self)
 Compute boundary condition matrix for free-boundary case.
 
real(r8) function dummy_fpp (self, psi)
 
subroutine get_olbp (smesh, olbp)
 Generate oriented loop from boundary points.
 
subroutine gs_analyze_saddles (self, o_point, o_psi, x_point, x_psi)
 Needs Docs.
 
subroutine gs_b_interp_apply (self, cell, f, gop, val)
 Reconstruct magnetic field from a Grad-Shafranov solution.
 
subroutine gs_bcrosskappa (self, bcross_kappa)
 Needs docs.
 
subroutine gs_coil_mutual (self, icoil, b, mutual)
 Compute inductance between coil and given poloidal flux.
 
subroutine gs_coil_mutual_distributed (self, icoil, b, curr_dist, mutual)
 Compute inductance between a coil with non-uniform current distribution and given poloidal flux.
 
subroutine gs_coil_source (self, icoil, b)
 Calculates coil current source for given coil.
 
subroutine gs_coil_source_distributed (self, icoil, b, curr_dist)
 Calculates coil current source for given coil with non-uniform current distribution.
 
subroutine gs_curvature (self, br, bt, bz, solver_in)
 Needs docs.
 
subroutine gs_curvature_apply (self, cell, f, gop, val)
 Reconstruct magnetic field from a Grad-Shafranov solution.
 
subroutine gs_destroy (self)
 Destroy internal storage and nullify references for Grad-Shafranov object.
 
real(8) function gs_dflux (self)
 Compute diamagentic flux.
 
character(len=40) function gs_err_reason (ierr)
 Compute Grad-Shafranov solution for current flux definitions.
 
subroutine gs_find_saddle (self, psi_scale_len, psi_x, pt, stype)
 Needs Docs.
 
subroutine gs_fit_isoflux (self, psi_full, ierr)
 Compute coil currents to best fit isoflux, flux, and saddle targets at current solution.
 
subroutine gs_fixed_vflux (self, pts, fluxes)
 Compute required vacuum flux for fixed boundary equilibrium.
 
subroutine gs_gen_source (self, source_fun, b)
 Compute RHS source from an arbitrary current distribution \( J_{\phi} \).
 
subroutine gs_get_chi (self)
 Compute toroidal flux potential from Grad-Shafranov solution.
 
subroutine gs_get_cond_scales (self, vals, skip_fixed)
 Needs Docs.
 
subroutine gs_get_cond_weights (self, vals, skip_fixed)
 Needs Docs.
 
subroutine gs_get_qprof (gseq, nr, psi_q, prof, dl, rbounds, zbounds, ravgs)
 Get q profile for equilibrium.
 
subroutine gs_helicity (self, ener, helic)
 Compute the magnetic energy and helicity of a fixed boundary equilibrium.
 
subroutine gs_init (self)
 Initialize Grad-Shafranov object by computing operators and allocating storage.
 
subroutine gs_init_psi (self, ierr, r0, a, kappa, delta, curr_source)
 Initialize \( \psi \).
 
real(8) function gs_itor (self, psi_vec)
 Compute toroidal current for Grad-Shafranov equilibrium.
 
subroutine gs_itor_nl (self, itor, centroid)
 Compute toroidal current for Grad-Shafranov equilibrium.
 
subroutine gs_j_interp_apply (self, cell, f, gop, val)
 Reconstruct magnetic field from a Grad-Shafranov solution.
 
subroutine gs_j_interp_delete (self)
 Destroy temporary internal storage and nullify references.
 
subroutine gs_j_interp_setup (self, gs)
 Needs Docs.
 
subroutine gs_lin_solve (self, adjust_r0, ierr)
 Compute solution to linearized Grad-Shafranov without updating \( \psi \) for RHS.
 
subroutine gs_load_limiters (self)
 Needs Docs.
 
subroutine gs_mat_create (fe_rep, new)
 Create G-S matrix with necessary augmentations for free-boundary BC.
 
subroutine gs_plasma_mutual (self, b, mutual, itor)
 Compute inductance between plasma current and given poloidal flux.
 
subroutine gs_prof_interp_apply (self, cell, f, gop, val)
 Reconstruct a component of a Grad-Shafranov solution.
 
subroutine gs_prof_interp_delete (self)
 Destroy temporary internal storage and nullify references.
 
subroutine gs_prof_interp_setup (self, gs)
 Needs Docs.
 
subroutine gs_project_b (self, br, bt, bz, solver_in, normalized)
 Needs Docs.
 
subroutine gs_psi2pt (self, psi_target, pt, pt_con, vec, psi_int)
 Find position of psi along a vector search direction.
 
subroutine gs_psi2r (self, psi_target, pt, psi_int)
 Find position of psi along a radial chord.
 
subroutine gs_save_fields (self, pts, npts, filename)
 Compute magnetic fields from Grad-Shafranov equilibrium.
 
subroutine gs_save_mug (self, filename, legacy)
 Needs Docs.
 
subroutine gs_save_prof (self, filename, mpsi_sample)
 Needs Docs.
 
subroutine gs_set_cond_weights (self, vals, skip_fixed)
 Needs Docs.
 
subroutine gs_setup (self, ml_lag_2d)
 Needs Docs.
 
subroutine gs_setup_walls (self, skip_load, make_plot)
 Needs Docs.
 
subroutine gs_solve (self, ierr)
 Compute Grad-Shafranov solution for current flux function definitions and targets.
 
subroutine gs_source (self, a, b, b2, b3, itor_alam, itor_press, estore)
 Compute plasma component of RHS source for Grad-Shafranov equation.
 
logical function gs_test_bounds (self, pt)
 Test whether a point is inside the LCFS.
 
subroutine gs_trace_surf (gseq, psi_in, points, npoints)
 Trace a single specified flux surface.
 
subroutine gs_update_bounds (self, track_opoint)
 Needs Docs.
 
subroutine gs_vac_solve (self, psi_sol, rhs_source, ierr)
 Compute Grad-Shafranov solution for vacuum (no plasma)
 
subroutine gs_vacuum_solve (self, pol_flux, source, ierr)
 Solve for \( \psi \) in vacuum for a given source \( \int \phi^T J_{\phi} dA \).
 
subroutine gs_wall_source (self, dpsi_dt, b)
 Calculate RHS source for quasi-static solve from previous \( \psi \).
 
subroutine gsinv_apply (self, cell, f, gop, val)
 Needs Docs.
 
subroutine gsinv_destroy (self)
 Needs Docs.
 
subroutine gsinv_setup (self, lag_rep)
 Needs Docs.
 
subroutine psi2pt_error (m, n, cofs, err, iflag)
 Needs docs.
 
subroutine psimax_error (m, n, cofs, err, iflag)
 Needs docs.
 
subroutine psimax_error_grad (m, n, cofs, err, jac_mat, ldjac_mat, iflag)
 Needs docs.
 
subroutine set_bcmat (self, mat)
 Add boundary condition terms for free-boundary case to matrix.
 
subroutine zerob_apply (self, a)
 Zero a 2D Lagrange scalar field at all nodes outside the plasma or on the global boundary.
 
subroutine zerob_delete (self)
 Destroy temporary internal storage and nullify references.
 

Variables

integer(i4cell_active = 0
 
real(r8), parameter gs_epsilon = 1.d-12
 Epsilon used for radial coordinate.
 
integer(4), parameter max_xpoints = 20
 
type(oft_lag_brinterp), pointer psi_eval_active => NULL()
 
type(oft_lag_bg2interp), pointer psi_g2eval_active => NULL()
 
type(oft_lag_bginterp), pointer psi_geval_active => NULL()
 
real(r8psi_target_active = 0.d0
 
real(r8), dimension(2) pt_con_active = [0.d0,0.d0]
 
real(r8qp_int_tol = 1.d-12
 
real(r8), dimension(2) vec_con_active = [0.d0,0.d0]
 

Function/Subroutine Documentation

◆ build_dels()

subroutine build_dels ( class(oft_matrix), intent(inout), pointer  mat,
class(gs_eq), intent(inout)  self,
character(len=*), intent(in)  bc,
real(8), intent(in), optional  dt,
real(8), intent(in), optional  scale 
)

Construct \( \frac{1}{R} \Delta^* \) operator, with or without time-dependence.

Supported boundary conditions

  • ‘'none’Full matrix -'zerob'Dirichlet for all boundary DOF -'grnd'` Dirichlet for only groundin point
    Parameters
    [in,out]matMatrix object
    [in,out]selfG-S object
    [in]bcBoundary condition
    [in]dtTimestep size for time-dependent version
    [in]scaleGlobal scale factor

◆ build_mrop()

subroutine build_mrop ( type(oft_scalar_bfem), intent(inout)  fe_rep,
class(oft_matrix), intent(inout), pointer  mat,
character(len=*), intent(in)  bc 
)

Construct mass matrix for a boundary Lagrange scalar representation.

Supported boundary conditions

  • ‘'none’Full matrix -'zerob'` Dirichlet for all boundary DOF
    Parameters
    [in,out]fe_rep2D Lagrange FE representation
    [in,out]matMatrix object
    [in]bcBoundary condition

◆ circle_interp()

subroutine circle_interp ( class(circular_curr), intent(inout)  self,
integer(4), intent(in)  cell,
real(8), dimension(:), intent(in)  f,
real(8), dimension(3,3), intent(in)  gop,
real(8), dimension(:), intent(out)  val 
)

Evaluate uniform source over simple plasma cross-section.

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

◆ compute_bcmat()

subroutine compute_bcmat ( class(gs_eq), intent(inout)  self)

Compute boundary condition matrix for free-boundary case.

◆ dummy_fpp()

real(r8) function dummy_fpp ( class(flux_func), intent(inout)  self,
real(r8), intent(in)  psi 
)

◆ get_olbp()

subroutine get_olbp ( class(oft_bmesh), intent(inout)  smesh,
integer(4), dimension(:), intent(out)  olbp 
)

Generate oriented loop from boundary points.

◆ gs_analyze_saddles()

subroutine gs_analyze_saddles ( class(gs_eq), intent(inout)  self,
real(8), dimension(2), intent(inout)  o_point,
real(8), intent(inout)  o_psi,
real(8), dimension(2,max_xpoints), intent(out)  x_point,
real(8), dimension(max_xpoints), intent(out)  x_psi 
)

Needs Docs.

◆ gs_b_interp_apply()

subroutine gs_b_interp_apply ( class(gs_b_interp), intent(inout)  self,
integer(4), intent(in)  cell,
real(8), dimension(:), intent(in)  f,
real(8), dimension(3,3), intent(in)  gop,
real(8), dimension(:), intent(out)  val 
)

Reconstruct magnetic field from a Grad-Shafranov solution.

Parameters
[in,out]selfInterpolation object
[in]cellCell for interpolation
[in]fPosition in cell in logical coord [3]
[in]gopLogical gradient vectors at f [3,3]
[out]valReconstructed field at f [3]

◆ gs_bcrosskappa()

subroutine gs_bcrosskappa ( class(gs_eq), intent(inout), target  self,
class(oft_vector), intent(inout), target  bcross_kappa 
)

Needs docs.

◆ gs_coil_mutual()

subroutine gs_coil_mutual ( class(gs_eq), intent(inout)  self,
integer(4), intent(in)  icoil,
class(oft_vector), intent(inout)  b,
real(8), intent(out)  mutual 
)

Compute inductance between coil and given poloidal flux.

Parameters
[in,out]selfG-S object
[in]icoilCoil index
[in,out]b\( \psi \) for mutual calculation
[out]mutualMutual inductance \( \int I_C \psi dV / I_C \)

◆ gs_coil_mutual_distributed()

subroutine gs_coil_mutual_distributed ( class(gs_eq), intent(inout)  self,
integer(4), intent(in)  icoil,
class(oft_vector), intent(inout)  b,
real(8), dimension(:), intent(in), pointer  curr_dist,
real(8), intent(out)  mutual 
)

Compute inductance between a coil with non-uniform current distribution and given poloidal flux.

Parameters
[in,out]selfG-S object
[in]icoilCoil index
[in,out]b\( \psi \) for mutual calculation
[in]curr_distRelative current density
[out]mutualMutual inductance \( \int I_C \psi dV / I_C \)

◆ gs_coil_source()

subroutine gs_coil_source ( class(gs_eq), intent(inout)  self,
integer(4), intent(in)  icoil,
class(oft_vector), intent(inout)  b 
)

Calculates coil current source for given coil.

Parameters
[in,out]selfG-S Object
[in]icoilCoil index
[in,out]bResulting source field

◆ gs_coil_source_distributed()

subroutine gs_coil_source_distributed ( class(gs_eq), intent(inout)  self,
integer(4), intent(in)  icoil,
class(oft_vector), intent(inout)  b,
real(8), dimension(:), intent(in), pointer  curr_dist 
)

Calculates coil current source for given coil with non-uniform current distribution.

Parameters
[in,out]selfG-S object
[in]icoilCoil index
[in,out]bResulting source field
[in]curr_distRelative current density

◆ gs_curvature()

subroutine gs_curvature ( class(gs_eq), intent(inout), target  self,
class(oft_vector), intent(inout), target  br,
class(oft_vector), intent(inout), target  bt,
class(oft_vector), intent(inout), target  bz,
class(oft_solver), intent(inout), optional, target  solver_in 
)

Needs docs.

◆ gs_curvature_apply()

subroutine gs_curvature_apply ( class(gs_curvature_interp), intent(inout)  self,
integer(4), intent(in)  cell,
real(8), dimension(:), intent(in)  f,
real(8), dimension(3,3), intent(in)  gop,
real(8), dimension(:), intent(out)  val 
)

Reconstruct magnetic field from a Grad-Shafranov solution.

Parameters
[in,out]selfInterpolation object
[in]cellCell for interpolation
[in]fPosition in cell in logical coord [3]
[in]gopLogical gradient vectors at f [3,3]
[out]valReconstructed field at f [3]

◆ gs_destroy()

subroutine gs_destroy ( class(gs_eq), intent(inout)  self)

Destroy internal storage and nullify references for Grad-Shafranov object.

Parameters
[in,out]selfG-S object

◆ gs_dflux()

real(8) function gs_dflux ( class(gs_eq), intent(inout)  self)

Compute diamagentic flux.

Parameters
[in,out]selfG-S object
Returns
Toroidal flux increment due to plasma

◆ gs_err_reason()

character(len=40) function gs_err_reason ( integer(4), intent(in)  ierr)

Compute Grad-Shafranov solution for current flux definitions.

Parameters
[in]ierrError flag
Returns
String representation of error

◆ gs_find_saddle()

subroutine gs_find_saddle ( class(gs_eq), intent(inout)  self,
real(8), intent(in)  psi_scale_len,
real(8), intent(inout)  psi_x,
real(8), dimension(2), intent(inout)  pt,
integer(4), intent(out)  stype 
)

Needs Docs.

◆ gs_fit_isoflux()

subroutine gs_fit_isoflux ( class(gs_eq), intent(inout)  self,
class(oft_vector), intent(inout), target  psi_full,
integer(4), intent(out)  ierr 
)

Compute coil currents to best fit isoflux, flux, and saddle targets at current solution.

Parameters
[in,out]selfG-S object
[in,out]psi_fullCurrent \( \psi \)
[out]ierrError flag

◆ gs_fixed_vflux()

subroutine gs_fixed_vflux ( class(gs_eq), intent(inout)  self,
real(8), dimension(:,:), intent(inout), pointer  pts,
real(8), dimension(:), intent(inout), pointer  fluxes 
)

Compute required vacuum flux for fixed boundary equilibrium.

Parameters
[in,out]selfG-S object
[in,out]ptsLocations of boundary points
[in,out]fluxesRequired flux at each point

◆ gs_gen_source()

subroutine gs_gen_source ( class(gs_eq), intent(inout)  self,
class(bfem_interp), intent(inout)  source_fun,
class(oft_vector), intent(inout)  b 
)

Compute RHS source from an arbitrary current distribution \( J_{\phi} \).

Parameters
[in,out]selfG-S object
[in,out]source_funInterpolation object for \( J_{\phi} \)
[in,out]bResulting source field

◆ gs_get_chi()

subroutine gs_get_chi ( class(gs_eq), intent(inout)  self)

Compute toroidal flux potential from Grad-Shafranov solution.

Parameters
[in,out]selfG-S object

◆ gs_get_cond_scales()

subroutine gs_get_cond_scales ( class(gs_eq), intent(inout)  self,
real(8), dimension(:), intent(inout)  vals,
logical, intent(in)  skip_fixed 
)

Needs Docs.

◆ gs_get_cond_weights()

subroutine gs_get_cond_weights ( class(gs_eq), intent(inout)  self,
real(8), dimension(:), intent(inout)  vals,
logical, intent(in)  skip_fixed 
)

Needs Docs.

◆ gs_get_qprof()

subroutine gs_get_qprof ( class(gs_eq), intent(inout)  gseq,
integer(4), intent(in)  nr,
real(8), dimension(nr), intent(in)  psi_q,
real(8), dimension(nr), intent(out)  prof,
real(8), intent(out), optional  dl,
real(8), dimension(2,2), intent(out), optional  rbounds,
real(8), dimension(2,2), intent(out), optional  zbounds,
real(8), dimension(nr,3), intent(out), optional  ravgs 
)

Get q profile for equilibrium.

Parameters
[in,out]gseqG-S object
[in]nrNumber of flux surfaces to sample
[in]psi_qLocations to sample in normalized flux
[out]profq value at each sampling location
[out]dlArc length of surface psi_q(1) (should be LCFS)
[out]rboundsRadial bounds of surface psi_q(1) (should be LCFS)
[out]zboundsVertical bounds of surface psi_q(1) (should be LCFS)
[out]ravgsFlux surface averages <R>, <1/R>, and dV/dPsi

◆ gs_helicity()

subroutine gs_helicity ( class(gs_eq), intent(inout)  self,
real(8), intent(out)  ener,
real(8), intent(out)  helic 
)

Compute the magnetic energy and helicity of a fixed boundary equilibrium.

Note
Helicity computed by this subroutine is only valid for equilibria with no normal field on the boundary
Parameters
[in,out]selfG-S object
[out]enerTotal magnetic energy
[out]helicTotal magnetic helicity

◆ gs_init()

subroutine gs_init ( class(gs_eq), intent(inout)  self)

Initialize Grad-Shafranov object by computing operators and allocating storage.

Parameters
[in,out]selfG-S object

◆ gs_init_psi()

subroutine gs_init_psi ( class(gs_eq), intent(inout)  self,
integer(4), intent(out)  ierr,
real(8), dimension(2), intent(in), optional  r0,
real(8), intent(in), optional  a,
real(8), intent(in), optional  kappa,
real(8), intent(in), optional  delta,
real(8), dimension(:), intent(in), optional  curr_source 
)

Initialize \( \psi \).

\( \psi \) can be initialized in one of four ways:

  1. If no optional arguments are passed a Taylor state is computed and used
  2. If r0 < 0.0 is passed a uniform current across the entire plasma region is used
  3. If r0 and a, and optionally kappa and delta, are passed then a uniform current across a cross-section defined by those parameters is used
  4. If curr_source is passed then the source defined by those values is used

In all cases the solution is scaled to match the target Ip value

Parameters
[in,out]selfG-S object
[out]ierrError flag
[in]r0Center for cross-section initialization
[in]aMinor radius for cross-section initialization
[in]kappaElongation for cross-section initialization
[in]deltaTriangularity for cross-section initialization
[in]curr_sourceExplicit current source

◆ gs_itor()

real(8) function gs_itor ( class(gs_eq), intent(inout)  self,
class(oft_vector), intent(inout), optional  psi_vec 
)

Compute toroidal current for Grad-Shafranov equilibrium.

Parameters
[in,out]selfG-S object
[in,out]psi_vecNeeds docs
Returns
Toroidal current

◆ gs_itor_nl()

subroutine gs_itor_nl ( class(gs_eq), intent(inout)  self,
real(8), intent(out)  itor,
real(8), dimension(2), intent(out), optional  centroid 
)

Compute toroidal current for Grad-Shafranov equilibrium.

Parameters
[in,out]selfG-S object
[out]itorToroidal current
[out]centroidCurrent centroid (optional) [2]

◆ gs_j_interp_apply()

subroutine gs_j_interp_apply ( class(gs_j_interp), intent(inout)  self,
integer(4), intent(in)  cell,
real(8), dimension(:), intent(in)  f,
real(8), dimension(3,3), intent(in)  gop,
real(8), dimension(:), intent(out)  val 
)

Reconstruct magnetic field from a Grad-Shafranov solution.

Parameters
[in,out]selfInterpolation object
[in]cellCell for interpolation
[in]fPosition in cell in logical coord [3]
[in]gopLogical gradient vectors at f [3,3]
[out]valReconstructed field at f [3]

◆ gs_j_interp_delete()

subroutine gs_j_interp_delete ( class(gs_j_interp), intent(inout)  self)

Destroy temporary internal storage and nullify references.

◆ gs_j_interp_setup()

subroutine gs_j_interp_setup ( class(gs_j_interp), intent(inout)  self,
class(gs_eq), intent(inout), target  gs 
)

Needs Docs.

◆ gs_lin_solve()

subroutine gs_lin_solve ( class(gs_eq), intent(inout)  self,
logical, intent(in)  adjust_r0,
integer(4), intent(out), optional  ierr 
)

Compute solution to linearized Grad-Shafranov without updating \( \psi \) for RHS.

Parameters
[in,out]selfG-S object
[in]adjust_r0Needs docs
[out]ierrError flag

◆ gs_load_limiters()

subroutine gs_load_limiters ( class(gs_eq), intent(inout)  self)

Needs Docs.

Parameters
[in,out]selfG-S object

◆ gs_mat_create()

subroutine gs_mat_create ( class(oft_scalar_bfem), intent(inout)  fe_rep,
class(oft_matrix), intent(out), pointer  new 
)

Create G-S matrix with necessary augmentations for free-boundary BC.

◆ gs_plasma_mutual()

subroutine gs_plasma_mutual ( class(gs_eq), intent(inout)  self,
class(oft_vector), intent(inout)  b,
real(8), intent(out)  mutual,
real(8), intent(out)  itor 
)

Compute inductance between plasma current and given poloidal flux.

Parameters
[in,out]selfG-S solver object
[in,out]b\( \psi \) for mutual calculation
[out]mutualMutual inductance \( \int J_p \psi dV / I_p \)
[out]itorPlasma toroidal current

◆ gs_prof_interp_apply()

subroutine gs_prof_interp_apply ( class(gs_prof_interp), intent(inout)  self,
integer(4), intent(in)  cell,
real(8), dimension(:), intent(in)  f,
real(8), dimension(3,3), intent(in)  gop,
real(8), dimension(:), intent(out)  val 
)

Reconstruct a component of a Grad-Shafranov solution.

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

◆ gs_prof_interp_delete()

subroutine gs_prof_interp_delete ( class(gs_prof_interp), intent(inout)  self)

Destroy temporary internal storage and nullify references.

◆ gs_prof_interp_setup()

subroutine gs_prof_interp_setup ( class(gs_prof_interp), intent(inout)  self,
class(gs_eq), intent(inout), target  gs 
)

Needs Docs.

◆ gs_project_b()

subroutine gs_project_b ( class(gs_eq), intent(inout), target  self,
class(oft_vector), intent(inout)  br,
class(oft_vector), intent(inout)  bt,
class(oft_vector), intent(inout)  bz,
class(oft_solver), intent(inout), optional, target  solver_in,
logical, intent(in), optional  normalized 
)

Needs Docs.

◆ gs_psi2pt()

subroutine gs_psi2pt ( class(gs_eq), intent(inout)  self,
real(8), intent(in)  psi_target,
real(8), dimension(2), intent(inout)  pt,
real(8), dimension(2), intent(in)  pt_con,
real(8), dimension(2), intent(in)  vec,
type(oft_lag_brinterp), intent(inout), optional, target  psi_int 
)

Find position of psi along a vector search direction.

Parameters
[in,out]selfG-S object
[in]psi_targetTarget \( \psi \) value to find
[in,out]ptGuess location (input); Closest point found (output) [2]
[in]pt_conLocation defining origin of search path
[in]vecVector defining direction of search path
[in,out]psi_intInterpolation object (created internally if not passed)

◆ gs_psi2r()

subroutine gs_psi2r ( class(gs_eq), intent(inout)  self,
real(8), intent(in)  psi_target,
real(8), dimension(2), intent(inout)  pt,
type(oft_lag_brinterp), intent(inout), optional, target  psi_int 
)

Find position of psi along a radial chord.

Parameters
[in,out]selfG-S object
[in]psi_targetTarget \( \psi \) value to find
[in,out]ptGuess location (input); Closest point found, by changing R only (output) [2]

◆ gs_save_fields()

subroutine gs_save_fields ( class(gs_eq), intent(inout)  self,
real(8), dimension(2,npts), intent(in)  pts,
integer(4), intent(in)  npts,
character(len=*), intent(in)  filename 
)

Compute magnetic fields from Grad-Shafranov equilibrium.

Parameters
[in,out]selfG-S object
[in]ptsSampling locations [2,npts]
[in]nptsNumber of points to sample
[in]filenameOutput file for field data

◆ gs_save_mug()

subroutine gs_save_mug ( class(gs_eq), intent(inout), target  self,
character(len=*), intent(in), optional  filename,
logical, intent(in), optional  legacy 
)

Needs Docs.

◆ gs_save_prof()

subroutine gs_save_prof ( class(gs_eq), intent(inout), target  self,
character(len=*), intent(in)  filename,
integer(4), intent(in), optional  mpsi_sample 
)

Needs Docs.

◆ gs_set_cond_weights()

subroutine gs_set_cond_weights ( class(gs_eq), intent(inout)  self,
real(8), dimension(:), intent(in)  vals,
logical, intent(in)  skip_fixed 
)

Needs Docs.

◆ gs_setup()

subroutine gs_setup ( class(gs_eq), intent(inout)  self,
class(oft_ml_fem_type), intent(inout), target  ml_lag_2d 
)

Needs Docs.

Parameters
[in,out]selfG-S object

◆ gs_setup_walls()

subroutine gs_setup_walls ( class(gs_eq), intent(inout)  self,
logical, intent(in), optional  skip_load,
logical, intent(in), optional  make_plot 
)

Needs Docs.

Parameters
[in,out]selfG-S object

◆ gs_solve()

subroutine gs_solve ( class(gs_eq), intent(inout)  self,
integer(4), intent(out), optional  ierr 
)

Compute Grad-Shafranov solution for current flux function definitions and targets.

Parameters
[in,out]selfG-S object
[out]ierrError flag

◆ gs_source()

subroutine gs_source ( class(gs_eq), intent(inout)  self,
class(oft_vector), intent(inout), target  a,
class(oft_vector), intent(inout)  b,
class(oft_vector), intent(inout)  b2,
class(oft_vector), intent(inout)  b3,
real(8), intent(out)  itor_alam,
real(8), intent(out)  itor_press,
real(8), intent(out)  estore 
)

Compute plasma component of RHS source for Grad-Shafranov equation.

Parameters
[in,out]selfG-S object
[in,out]a\( \psi \)
[in,out]bFull RHS source
[in,out]b2F*F' component of source (including alam)
[in,out]b3P' component of source (without pnorm)

◆ gs_test_bounds()

logical function gs_test_bounds ( class(gs_eq), intent(inout)  self,
real(8), dimension(2), intent(in)  pt 
)

Test whether a point is inside the LCFS.

Parameters
[in,out]selfG-S object
[in]ptLocation to test in/out of plasma

◆ gs_trace_surf()

subroutine gs_trace_surf ( class(gs_eq), intent(inout)  gseq,
real(8), intent(in)  psi_in,
real(8), dimension(:,:), intent(out), pointer  points,
integer(4), intent(out)  npoints 
)

Trace a single specified flux surface.

Parameters
[in,out]gseqG-S object
[in]psi_inLocations of surface to trace in normalized flux
[out]pointsTraced surface
[out]npointsNumber of points in traced surface

◆ gs_update_bounds()

subroutine gs_update_bounds ( class(gs_eq), intent(inout)  self,
logical, intent(in), optional  track_opoint 
)

Needs Docs.

◆ gs_vac_solve()

subroutine gs_vac_solve ( class(gs_eq), intent(inout)  self,
class(oft_vector), intent(inout)  psi_sol,
class(bfem_interp), intent(inout), optional  rhs_source,
integer(4), intent(out), optional  ierr 
)

Compute Grad-Shafranov solution for vacuum (no plasma)

Parameters
[in,out]selfG-S object
[in,out]psi_solInput: BCs for \( \psi \), Output: solution
[in,out]rhs_sourceSpecified current source (optional)
[out]ierrError flag

◆ gs_vacuum_solve()

subroutine gs_vacuum_solve ( class(gs_eq), intent(inout)  self,
class(oft_vector), intent(inout)  pol_flux,
class(oft_vector), intent(inout)  source,
integer(4), intent(out), optional  ierr 
)

Solve for \( \psi \) in vacuum for a given source \( \int \phi^T J_{\phi} dA \).

Source field must already be projected onto the appropriate Lagrange FE basis

Parameters
[in,out]selfG-S object
[in,out]pol_flux\( \psi \) (output)
[in,out]sourceCurrent source
[out]ierrError flag

◆ gs_wall_source()

subroutine gs_wall_source ( class(gs_eq), intent(inout)  self,
class(oft_vector), intent(inout)  dpsi_dt,
class(oft_vector), intent(inout)  b 
)

Calculate RHS source for quasi-static solve from previous \( \psi \).

Parameters
[in,out]selfG-S object
[in,out]dpsi_dt\( \psi \) at start of step
[in,out]bResulting source field

◆ gsinv_apply()

subroutine gsinv_apply ( class(gsinv_interp), intent(inout)  self,
integer(4), intent(in)  cell,
real(8), dimension(:), intent(in)  f,
real(8), dimension(3,3), intent(in)  gop,
real(8), dimension(:), intent(out)  val 
)

Needs Docs.

◆ gsinv_destroy()

subroutine gsinv_destroy ( class(gsinv_interp), intent(inout)  self)

Needs Docs.

◆ gsinv_setup()

subroutine gsinv_setup ( class(gsinv_interp), intent(inout)  self,
class(oft_afem_type), intent(inout), target  lag_rep 
)

Needs Docs.

◆ psi2pt_error()

subroutine psi2pt_error ( integer(4), intent(in)  m,
integer(4), intent(in)  n,
real(8), dimension(n), intent(in)  cofs,
real(8), dimension(m), intent(out)  err,
integer(4), intent(inout)  iflag 
)

Needs docs.

◆ psimax_error()

subroutine psimax_error ( integer(4), intent(in)  m,
integer(4), intent(in)  n,
real(8), dimension(n), intent(in)  cofs,
real(8), dimension(m), intent(out)  err,
integer(4), intent(inout)  iflag 
)

Needs docs.

◆ psimax_error_grad()

subroutine psimax_error_grad ( integer(4), intent(in)  m,
integer(4), intent(in)  n,
real(8), dimension(n), intent(in)  cofs,
real(8), dimension(m), intent(out)  err,
real(8), dimension(ldjac_mat,n), intent(out)  jac_mat,
integer(4), intent(in)  ldjac_mat,
integer(4), intent(in)  iflag 
)

Needs docs.

◆ set_bcmat()

subroutine set_bcmat ( class(gs_eq), intent(inout)  self,
class(oft_matrix), intent(inout)  mat 
)

Add boundary condition terms for free-boundary case to matrix.

Parameters
[in,out]selfG-S object
[in,out]matMatrix object

◆ zerob_apply()

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

Zero a 2D Lagrange scalar field at all nodes outside the plasma or on the global boundary.

Parameters
[in,out]selfBC object
[in,out]aField to be zeroed

◆ zerob_delete()

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

Destroy temporary internal storage and nullify references.

Parameters
[in,out]selfBC object

Variable Documentation

◆ cell_active

integer(i4) cell_active = 0

◆ gs_epsilon

real(r8), parameter gs_epsilon = 1.d-12

Epsilon used for radial coordinate.

◆ max_xpoints

integer(4), parameter max_xpoints = 20

◆ psi_eval_active

type(oft_lag_brinterp), pointer psi_eval_active => NULL()

◆ psi_g2eval_active

type(oft_lag_bg2interp), pointer psi_g2eval_active => NULL()

◆ psi_geval_active

type(oft_lag_bginterp), pointer psi_geval_active => NULL()

◆ psi_target_active

real(r8) psi_target_active = 0.d0

◆ pt_con_active

real(r8), dimension(2) pt_con_active = [0.d0,0.d0]

◆ qp_int_tol

real(r8) qp_int_tol = 1.d-12

◆ vec_con_active

real(r8), dimension(2) vec_con_active = [0.d0,0.d0]