SPEC 3.20
Stepped Pressure Equilibrium Code
|
The namelist locallist
controls the construction of the Beltrami fields in each volume.
More...
Variables | |
integer | inputlist::lbeltrami = 4 |
Control flag for solution of Beltrami equation. More... | |
integer | inputlist::linitgues = 1 |
controls how initial guess for Beltrami field is constructed More... | |
integer | inputlist::lposdef = 0 |
redundant; | |
real | inputlist::maxrndgues = 1.0 |
the maximum random number of the Beltrami field if Linitgues = 3 | |
integer | inputlist::lmatsolver = 3 |
1 for LU factorization, 2 for GMRES, 3 for GMRES matrix-free | |
integer | inputlist::nitergmres = 200 |
number of max iteration for GMRES | |
real | inputlist::epsgmres = 1e-14 |
the precision of GMRES | |
integer | inputlist::lgmresprec = 1 |
type of preconditioner for GMRES, 1 for ILU sparse matrix | |
real | inputlist::epsilu = 1e-12 |
the precision of incomplete LU factorization for preconditioning | |
The namelist locallist
controls the construction of the Beltrami fields in each volume.
The transformation to straight-fieldline coordinates is singular when the rotational-transform of the interfaces is rational; however, the rotational-transform is still well defined.
integer inputlist::lbeltrami = 4 |
Control flag for solution of Beltrami equation.
LBeltrami
= 1,3,5 or 7, (SQP) then the Beltrami field in each volume is constructed by minimizing the magnetic energy with the constraint of fixed helicity; this is achieved by using sequential quadratic programming as provided by E04UFF
. This approach has the benefit (in theory) of robustly constructing minimum energy solutions when multiple, i.e. bifurcated, solutions exist. LBeltrami
= 2,3,6 or 7, (Newton) then the Beltrami fields are constructed by employing a standard Newton method for locating an extremum of \(F\equiv \int B^2 dv - \mu (\int {\bf A}\cdot{\bf B}dv-{\cal K})\), where \(\mu\) is treated as an independent degree of freedom similar to the parameters describing the vector potential and \({\cal K}\) is the required value of the helicity; this is the standard Lagrange multipler approach for locating the constrained minimum; this method cannot distinguish saddle-type extrema from minima, and which solution that will be obtained depends on the initial guess; LBeltrami
= 4,5,6 or 7, (linear) it is assumed that the Beltrami fields are parameterized by \(\mu\); in this case, it is only required to solve \(\nabla \times {\bf B} = \mu {\bf B}\) which reduces to a system of linear equations; \(\mu\) may or may not be adjusted iteratively, depending on Lconstraint
, to satisfy either rotational-transform or helicity constraints; LBeltrami
= 1, only the SQP method will be employed; LBeltrami
= 2, only the Newton method will be employed; LBeltrami
= 4, only the linear method will be employed; LBeltrami
= 3, the SQP and the Newton method are used; LBeltrami
= 5, the SQP and the linear method are used; LBeltrami
= 6, the Newton and the linear method are used; LBeltrami
= 7, all three methods will be employed; Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().
integer inputlist::linitgues = 1 |
controls how initial guess for Beltrami field is constructed
Linitgues
= 0, the initial guess for the Beltrami field is trivial Linitgues
= 1, the initial guess for the Beltrami field is an integrable approximation Linitgues
= 2, the initial guess for the Beltrami field is read from file Linitgues
= 3, the initial guess for the Beltrami field will be randomized with the maximum maxrndgues
Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), sphdf5::mirror_input_to_outfile(), preset(), and allglobal::wrtend().