SPEC 3.20
Stepped Pressure Equilibrium Code

The namelist globallist controls the search for global force-balance. More...

Collaboration diagram for globallist:

Variables

integer inputlist::lfindzero = 0
 use Newton methods to find zero of force-balance, which is computed by dforce() More...
 
real inputlist::escale = 0.0
 controls the weight factor, BBweight, in the force-imbalance harmonics More...
 
real inputlist::opsilon = 1.0
 weighting of force-imbalance More...
 
real inputlist::pcondense = 2.0
 spectral condensation parameter More...
 
real inputlist::epsilon = 0.0
 weighting of spectral-width constraint More...
 
real inputlist::wpoloidal = 1.0
 "star-like" poloidal angle constraint radial exponential factor used in preset() to construct sweight
 
real inputlist::upsilon = 1.0
 weighting of "star-like" poloidal angle constraint used in preset() to construct sweight
 
real inputlist::forcetol = 1.0e-10
 required tolerance in force-balance error; only used as an initial check More...
 
real inputlist::c05xmax = 1.0e-06
 required tolerance in position, \({\bf x} \equiv \{ R_{i,v}, Z_{i,v}\}\)
 
real inputlist::c05xtol = 1.0e-12
 required tolerance in position, \({\bf x} \equiv \{ R_{i,v}, Z_{i,v}\}\) More...
 
real inputlist::c05factor = 1.0e-02
 used to control initial step size in C05NDF and C05PDF More...
 
logical inputlist::lreadgf = .true.
 read \(\nabla_{\bf x} {\bf F}\) from file ext.GF More...
 
integer inputlist::mfreeits = 0
 maximum allowed free-boundary iterations More...
 
real inputlist::bnstol = 1.0e-06
 redundant;
 
real inputlist::bnsblend = 0.666
 redundant;
 
real inputlist::gbntol = 1.0e-06
 required tolerance in free-boundary iterations More...
 
real inputlist::gbnbld = 0.666
 normal blend More...
 
real inputlist::vcasingeps = 1.e-12
 regularization of Biot-Savart; see bnorml(), casing()
 
real inputlist::vcasingtol = 1.e-08
 accuracy on virtual casing integral; see bnorml(), casing()
 
integer inputlist::vcasingits = 8
 minimum number of calls to adaptive virtual casing routine; see casing()
 
integer inputlist::vcasingper = 1
 periods of integragion in adaptive virtual casing routine; see casing()
 
integer inputlist::mcasingcal = 8
 minimum number of calls to adaptive virtual casing routine; see casing(); redundant;
 

Detailed Description

The namelist globallist controls the search for global force-balance.

Comments:

Variable Documentation

◆ lfindzero

integer inputlist::lfindzero = 0

use Newton methods to find zero of force-balance, which is computed by dforce()

  • if Lfindzero = 0, then dforce() is called once to compute the Beltrami fields consistent with the given geometry and constraints
  • if Lfindzero = 1, then call C05NDF (uses function values only), which iteratively calls dforce()
  • if Lfindzero = 2, then call C05PDF (uses derivative information), which iteratively calls dforce()

Referenced by brcast(), allglobal::broadcast_inputs(), allglobal::check_inputs(), dfp200(), fcn1(), fcn2(), hesian(), sphdf5::mirror_input_to_outfile(), newton(), packxi(), preset(), spec(), and allglobal::wrtend().

◆ escale

real inputlist::escale = 0.0

controls the weight factor, BBweight, in the force-imbalance harmonics

  • BBweight(i) \(\displaystyle \equiv \texttt{opsilon} \times \exp\left[-\texttt{escale} \times (m_i^2+n_i^2) \right]\)
  • defined in preset() ; used in dforce()
  • also see Eqn. \((\ref{eq:forcebalancemn_global})\)

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().

◆ opsilon

real inputlist::opsilon = 1.0

weighting of force-imbalance

  • used in dforce(); also see Eqn. \((\ref{eq:forcebalancemn_global})\)

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().

◆ pcondense

real inputlist::pcondense = 2.0

spectral condensation parameter

  • used in preset() to define mmpp(i) \(\equiv m_i^p\), where \(p\equiv \) pcondense
  • the angle freedom is exploited to minimize \(\displaystyle \texttt{epsilon} \sum_{i} m_i^p (R_{i}^2+Z_{i}^2)\) with respect to tangential variations in the interface geometry
  • also see Eqn. \((\ref{eq:spectralbalancemn_global})\)

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().

◆ epsilon

real inputlist::epsilon = 0.0

weighting of spectral-width constraint

  • used in dforce(); also see Eqn. \((\ref{eq:spectralbalancemn_global})\)

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), dforce(), dfp200(), evaluate_dbb(), sphdf5::mirror_input_to_outfile(), pc00ab(), and allglobal::wrtend().

◆ forcetol

real inputlist::forcetol = 1.0e-10

required tolerance in force-balance error; only used as an initial check

  • if the initially supplied interfaces are consistent with force-balance to within forcetol then the geometry of the interfaces is not altered
  • if not, then the geometry of the interfaces is changed in order to bring the configuration into force balance so that the geometry of interfaces is within c05xtol, defined below, of the true solution
  • to force execution of either C05NDF or C05PDF, regardless of the initial force imbalance, set forcetol < 0

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), fcn1(), fcn2(), sphdf5::mirror_input_to_outfile(), newton(), pc00aa(), pc00ab(), preset(), and allglobal::wrtend().

◆ c05xtol

real inputlist::c05xtol = 1.0e-12

required tolerance in position, \({\bf x} \equiv \{ R_{i,v}, Z_{i,v}\}\)

  • used by both C05NDF and C05PDF; see the NAG documents for further details on how the error is defined
  • constraint c05xtol > 0.0

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), fcn1(), fcn2(), sphdf5::mirror_input_to_outfile(), newton(), and allglobal::wrtend().

◆ c05factor

real inputlist::c05factor = 1.0e-02

used to control initial step size in C05NDF and C05PDF

  • constraint c05factor > 0.0
  • only relevant if Lfindzero > 0

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), fcn1(), fcn2(), sphdf5::mirror_input_to_outfile(), newton(), and allglobal::wrtend().

◆ lreadgf

logical inputlist::lreadgf = .true.

read \(\nabla_{\bf x} {\bf F}\) from file ext.GF

  • only used if Lfindzero = 2
  • only used in newton()

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), fcn1(), fcn2(), sphdf5::mirror_input_to_outfile(), newton(), and allglobal::wrtend().

◆ mfreeits

integer inputlist::mfreeits = 0

maximum allowed free-boundary iterations

  • only used if Lfreebound = 1
  • only used in xspech()

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), spec(), and allglobal::wrtend().

◆ gbntol

real inputlist::gbntol = 1.0e-06

required tolerance in free-boundary iterations

  • only used if Lfreebound = 1
  • only used in xspech()

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), spec(), and allglobal::wrtend().

◆ gbnbld

real inputlist::gbnbld = 0.666

normal blend

  • The "new" magnetic field at the computational boundary produced by the plasma currents is updated using a Picard scheme:

    \begin{eqnarray} ({\bf B}\cdot{\bf n})^{j+1} = \texttt{gBnbld} \times ({\bf B}\cdot{\bf n})^{j} + (1-\texttt{gBnbld}) \times ({\bf B}\cdot{\bf n})^{*}, \end{eqnarray}

    where \(j\) labels free-boundary iterations, and \(({\bf B}\cdot{\bf n})^{*}\) is computed by virtual casing.
  • only used if Lfreebound = 1
  • only used in xspech()

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), spec(), and allglobal::wrtend().