![]() |
SPEC 3.20
Stepped Pressure Equilibrium Code
|
The namelist globallist controls the search for global force-balance.
More...
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::vcnt = 256 |
| theta resolution of the real space grid on which the surface current is computed during virtual casing [0, 2pi] | |
| integer | inputlist::vcnz = 256 |
| zeta resolution of the real space grid on which the surface current is computed during virtual casing [0, 2pi], over one field period | |
| integer | inputlist::mcasingcal = 8 |
| minimum number of calls to adaptive virtual casing routine; see casing(); redundant; | |
The namelist globallist controls the search for global force-balance.
Comments:
\begin{eqnarray} F_{i,v} \equiv [[ p+B^2/2 ]]_{i,v} \times \exp\left[-\texttt{escale}(m_i^2+n_i^2) \right] \times \texttt{opsilon}, \label{eq:forcebalancemn_global} \end{eqnarray}
and spectral-condensation constraints, \(I_{i,v}\), and the "star-like" angle constraints, \(S_{i,v,}\), (see lforce() for details)\begin{eqnarray} F_{i,v} \equiv \texttt{epsilon} \times I_{i,v} + \texttt{upsilon} \times \left( \psi_v^\omega S_{i,v,1} - \psi_{v+1}^\omega S_{i,v+1,0} \right), \label{eq:spectralbalancemn_global} \end{eqnarray}
where \(\psi_v\equiv\) normalized toroidal flux,tflux, and \(\omega\equiv\) wpoloidal. | integer inputlist::lfindzero = 0 |
use Newton methods to find zero of force-balance, which is computed by dforce()
Lfindzero = 0, then dforce() is called once to compute the Beltrami fields consistent with the given geometry and constraints Lfindzero = 1, then call C05NDF (uses function values only), which iteratively calls dforce() 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().
| 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]\) Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().
| real inputlist::opsilon = 1.0 |
weighting of force-imbalance
Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().
| real inputlist::pcondense = 2.0 |
spectral condensation parameter
mmpp(i) \(\equiv m_i^p\), where \(p\equiv \) pcondense Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().
| real inputlist::epsilon = 0.0 |
weighting of spectral-width constraint
Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), dforce(), dfp200(), evaluate_dbb(), sphdf5::mirror_input_to_outfile(), pc00ab(), and allglobal::wrtend().
| real inputlist::forcetol = 1.0e-10 |
required tolerance in force-balance error; only used as an initial check
forcetol then the geometry of the interfaces is not altered c05xtol, defined below, of the true solution 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().
| real inputlist::c05xtol = 1.0e-12 |
required tolerance in position, \({\bf x} \equiv \{ R_{i,v}, Z_{i,v}\}\)
C05NDF and C05PDF; see the NAG documents for further details on how the error is defined c05xtol > 0.0 Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), fcn1(), fcn2(), sphdf5::mirror_input_to_outfile(), newton(), and allglobal::wrtend().
| real inputlist::c05factor = 1.0e-02 |
used to control initial step size in C05NDF and C05PDF
c05factor > 0.0 Lfindzero > 0 Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), fcn1(), fcn2(), sphdf5::mirror_input_to_outfile(), newton(), and allglobal::wrtend().
| logical inputlist::lreadgf = .true. |
read \(\nabla_{\bf x} {\bf F}\) from file ext.GF
Lfindzero = 2 Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), fcn1(), fcn2(), sphdf5::mirror_input_to_outfile(), newton(), and allglobal::wrtend().
| integer inputlist::mfreeits = 0 |
maximum allowed free-boundary iterations
Lfreebound = 1 Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), spec(), and allglobal::wrtend().
| real inputlist::gbntol = 1.0e-06 |
required tolerance in free-boundary iterations
Lfreebound = 1 Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), spec(), and allglobal::wrtend().
| real inputlist::gbnbld = 0.666 |
normal blend
\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.Lfreebound = 1 Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), spec(), and allglobal::wrtend().