SPEC 3.20
Stepped Pressure Equilibrium Code
Initialization of the code

More...

Functions/Subroutines

subroutine preset
 Allocates and initializes internal arrays. More...
 

Detailed Description

Function/Subroutine Documentation

◆ preset()

subroutine preset

Allocates and initializes internal arrays.

LGdof and NGdof : number of geometrical degrees-of-freedom

  • LGdof \(\equiv\) the number of degrees-of-freedom in the geometry (i.e. Fourier harmonics) of each interface
  • NGdof \(\equiv\) total number of degrees-of-freedom in geometry, i.e. of all interfaces

iota and oita: rotational transform on interfaces

  • The input variables iota and oita are the rotational transform on "inner-side" and on the "outer-side" of each interface.
  • These quantities are formally inputs.
  • Note that if \(q_l+\gamma q_r \ne 0\), then iota is given by

    \begin{eqnarray} {{\,\iota\!\!\!}-} \equiv \frac{p_l + \gamma p_r}{q_l + \gamma q_r}, \end{eqnarray}

    where \(p_l \equiv\,\)pl, \(q_l \equiv\,\)ql , etc.; and similarly for oita .

dtflux(1:Mvol) and dpflux(1:Mvol): enclosed fluxes

  • dtflux \(\equiv \Delta \psi_{tor} / 2\pi\) and dpflux \(\equiv \Delta \psi_{pol} / 2\pi\) in each volume.
  • Note that the total toroidal flux enclosed by the plasma boundary is \(\Phi_{edge} \equiv\,\)phiedge .
  • \(\psi_{tor} \equiv\,\)tflux and \(\psi_{pol} \equiv\,\)pflux are immediately normalized (in readin() ) according to \(\psi_{tor,i} \rightarrow \psi_{tor,i} / \psi_{0}\) and \(\psi_{pol,i} \rightarrow \psi_{pol,i} / \psi_{0}\), where \(\psi_{0} \equiv \psi_{tor,N}\) on input.

sweight(1:Mvol): star-like angle constraint weight

  • the "star-like" poloidal angle constraint weights (only required for toroidal geometry, i.e. Igeometry=3) are given by

    \begin{eqnarray} \texttt{sweight}_v \equiv \texttt{upsilon} \times (l_v / N_{vol})^w, \end{eqnarray}

    where \(l_v\) is the volume number, and \(w \equiv\,\)wpoloidal.

TT(0:Mrad,0:1,0:1): Chebyshev polynomials at inner/outer interface

  • TT(0:Lrad,0:1,0:1) gives the Chebyshev polynomials, and their first derivative, evaluated at \(s=-1\) and \(s=+1\).
  • Precisely, TT(l,i,d) \(\equiv T_l^{(d)}(s_i)\) for \(s_0=-1\) and \(s_1=+1\).
  • Note that \(T_l^{(0)}(s)=s^l\) and \(T_l^{(1)}(s)=s^{l+1} l^2\) for \(s=\pm 1\).
  • Note that

    \begin{eqnarray} T_l(-1) = \left\{ \begin{array}{ccccccccccccccc}+1,& \textrm{ if $l$ is even,} \\ -1,& \textrm{ if $l$ is odd;} \end{array} \right. & \; \;& T_l(+1) = \left\{ \begin{array}{ccccccccccccccc}+1,& \textrm{ if $l$ is even,} \\ +1,& \textrm{ if $l$ is odd;} \end{array} \right. \\ T_l^\prime(-1) = \left\{ \begin{array}{ccccccccccccccc}-l^2,& \textrm{ if $l$ is even,} \\ +l^2,& \textrm{ if $l$ is odd;} \end{array} \right. &\; \;& T_l^\prime(+1) = \left\{ \begin{array}{ccccccccccccccc}+l^2,& \textrm{ if $l$ is even,} \\ +l^2,& \textrm{ if $l$ is odd.} \end{array} \right. \end{eqnarray}

  • TT(0:Mrad,0:1,0:1) is used in routines that explicity require interface information, such as
    • the interface force-balance routine, lforce()
    • the virtual casing routine, casing()
    • computing the rotational-transform on the interfaces, tr00ab()
    • computing the covariant components of the interface magnetic field, sc00aa()
    • enforcing the constraints on the Beltrami fields, matrix() and
    • computing the enclosed currents of the vacuum field, curent().

ImagneticOK(1:Mvol): Beltrami/vacuum error flag

  • error flags that indicate if the magnetic field in each volume has been successfully constructed
  • ImagneticOK is initialized to .false. in dforce() before the Beltrami solver routines are called. If the construction of the Beltrami field is successful (in either ma02aa() or mp00ac() ) then ImagneticOK is set to .true. .

Lhessianallocated

  • The internal logical variable, Lhessianallocated, indicates whether the `‘Hessian’' matrix of second-partial derivatives (really, the first derivatives of the force-vector) has been allocated, or not!

ki(1:mn,0:1): Fourier identification

  • Consider the "abbreviated" representation for a double Fourier series,

    \begin{eqnarray} \sum_i f_i \cos( m_i \theta - n_i \zeta) \equiv \sum_{n= 0 }^{ N_0} f_{0,n} \cos( -n\zeta) + \sum_{m=1}^{ M_0} \sum_{n=- N_0}^{ N_0} f_{m,n} \cos(m\theta-n\zeta), \end{eqnarray}

    and the same representation but with enhanced resolution,

    \begin{eqnarray} \sum_k \bar f_k \cos(\bar m_k \theta - \bar n_k \zeta) \equiv \sum_{n= 0 }^{ N_1} f_{0,n} \cos( -n\zeta) + \sum_{m=1}^{ M_1} \sum_{n=- N_1}^{ N_1} f_{m,n} \cos(m\theta-n\zeta), \label{eq:enhancedFourierrepresentation_preset} \end{eqnarray}

    with \(M_1 \ge M_0\) and \(N_1 \ge N_0\); then \(k_i\equiv\,\)ki(i,0) is defined such that \(\bar m_{k_i} = m_i\) and \(\bar n_{k_i} = n_i\).

kija(1:mn,1:mn,0:1), kijs(1:mn,1:mn,0:1): Fourier identification

  • Consider the following quantities, which are computed in ma00aa(), where \(\bar g^{\mu\nu} = \sum_k \bar g^{\mu\nu}_k \cos \alpha_k\) for \(\alpha_k \equiv m_k \theta - n_k \zeta\),

    \begin{eqnarray} \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} \cos\alpha_i \; \cos\alpha_j & = & \frac{1}{2} \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} ( + \cos\alpha_{k_{ij+}} + \cos\alpha_{k_{ij-}} ), \\ \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} \cos\alpha_i \; \sin\alpha_j & = & \frac{1}{2} \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} ( + \sin\alpha_{k_{ij+}} - \sin\alpha_{k_{ij-}} ), \\ \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} \sin\alpha_i \; \cos\alpha_j & = & \frac{1}{2} \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} ( + \sin\alpha_{k_{ij+}} + \sin\alpha_{k_{ij-}} ), \\ \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} \sin\alpha_i \; \sin\alpha_j & = & \frac{1}{2} \oint\!\!\!\oint \!d\theta d\zeta \,\,\, \bar g^{\mu\nu} ( - \cos\alpha_{k_{ij+}} + \cos\alpha_{k_{ij-}} ), \end{eqnarray}

    where \((m_{k_{ij+}}, n_{k_{ij+}}) = (m_i + m_j, n_i + n_j)\) and \((m_{k_{ij-}}, n_{k_{ij-}}) = (m_i - m_j, n_i - n_j)\); then kija(i,j,0) \(\equiv k_{ij+}\) and kijs(i,j,0) \(\equiv k_{ij-}\).
  • Note that Eqn. \((\ref{eq:enhancedFourierrepresentation_preset})\) does not include \(m<0\); so, if \(m_i - m_j < 0\) then \(k_{ij-}\) is re-defined such that \((m_{k_{ij-}}, n_{k_{ij-}}) = (m_j - m_i, n_j - n_i)\); and similarly for the case \(m=0\) and \(n<0\). Also, take care that the sign of the sine harmonics in the above expressions will change for these cases.

djkp

iotakki

cheby(0:Lrad,0:2): Chebyshev polynomial workspace

  • cheby(0:Lrad,0:2) is global workspace for computing the Chebyshev polynomials, and their derivatives, using the recurrence relations \(T_0(s) = 1\), \(T_1(s) = s\) and \(T_l(s) = 2 \, s \,T_{l-1}(s) - T_{l-2}(s)\).
  • These are computed as required, i.e. for arbitrary \(s\), in bfield(), jo00aa() and ma00aa().
  • Note that the quantities required for ma00aa() are for fixed \(s\), and so these quantities should be precomputed.

Iquad, gaussianweight, gaussianabscissae: Gauss-Legendre quadrature

  • The volume integrals are computed using a "Fourier" integration over the angles and by Gauss-Legendre quadrature over the radial, i.e. \(\displaystyle \int \!\! f(s) ds = \sum_k \omega_k f(s_k)\).
  • The quadrature resolution in each volume is give by Iquad(1:Mvol) which is determined as follows:
    • if Nquad.gt.0 , then Iquad(vvol)=Nquad
    • if Nquad.le.0 and .not.Lcoordinatesingularity , then Iquad(vvol)=2*Lrad(vvol)-Nquad
    • if Nquad.le.0 and Lcoordinatesingularity , then Iquad(vvol)=2*Lrad(vvol)-Nquad+Mpol
  • The Gaussian weights and abscissae are given by gaussianweight(1:maxIquad,1:Mvol) and gaussianabscissae(1:maxIquad,1:Mvol), which are computed using modified Numerical Recipes routine gauleg() .
  • Iquad \(_v\) is passed through to ma00aa() to compute the volume integrals of the metric elements; also see jo00aa(), where Iquad \(_v\) is used to compute the volume integrals of \(||\nabla\times{\bf B} - \mu {\bf B}||\).

LBsequad, LBnewton and LBlinear

  • LBsequad, LBnewton and LBlinear depend simply on LBeltrami , which is described in global.f90 .

BBweight(1:mn): weighting of force-imbalance harmonics

  • weight on force-imbalance harmonics;

    \begin{eqnarray} \texttt{BBweight}_i \equiv \texttt{opsilon} \times \exp\left[ - \texttt{escale} \times (m_i^2 + n_i^2) \right] \end{eqnarray}

  • this is only used in dforce() in constructing the force-imbalance vector

mmpp(1:mn): spectral condensation weight factors

  • spectral condensation weight factors;

    \begin{eqnarray} \texttt{mmpp(i)} \equiv m_i^p, \end{eqnarray}

    where \(p \equiv\,\)pcondense .

NAdof, Ate, Aze, Ato and Azo: degrees-of-freedom in magnetic vector potential

  • NAdof(1:Mvol) \(\equiv\) total number of degrees-of-freedom in magnetic vector potential, including Lagrange multipliers, in each volume. This can de deduced from matrix().
  • The components of the vector potential, \({\bf A}=A_\theta \nabla + A_\zeta \nabla \zeta\), are

    \begin{eqnarray} A_\theta(s,\theta,\zeta) &=& \sum_{i,l} {\color{red} A_{\theta,e,i,l}} \; {\overline T}_{l,i}(s) \cos\alpha_i + \sum_{i,l} {\color{Orange} A_{\theta,o,i,l}} \; {\overline T}_{l,i}(s) \sin\alpha_i, \label{eq:At_preset} \\ A_\zeta( s,\theta,\zeta) &=& \sum_{i,l} {\color{blue} A_{\zeta, e,i,l}} \; {\overline T}_{l,i}(s) \cos\alpha_i + \sum_{i,l} {\color{Cerulean}A_{\zeta ,o,i,l}} \; {\overline T}_{l,i}(s) \sin\alpha_i, \label{eq:Az_preset} \end{eqnarray}

    where \({\overline T}_{l,i}(s) \equiv \bar s^{m_i/2} \, T_l(s)\), \(T_l(s)\) is the Chebyshev polynomial, and \(\alpha_j \equiv m_j\theta-n_j\zeta\). The regularity factor, \(\bar s^{m_i/2}\), where \(\bar s \equiv (1+s)/2\), is only included if there is a coordinate singularity in the domain (i.e. only in the innermost volume, and only in cylindrical and toroidal geometry.)
  • The Chebyshev-Fourier harmonics of the covariant components of the magnetic vector potential are kept in

    \begin{eqnarray} {\color{red} A_{\theta,e,i,l}} &\equiv& \texttt{Ate(v,0,j)}\%\texttt{s(l)} , \\ {\color{blue} A_{\zeta, e,i,l}} &\equiv& \texttt{Aze(v,0,j)}\%\texttt{s(l)} , \\ {\color{Orange} A_{\theta,o,i,l}} &\equiv& \texttt{Ato(v,0,j)}\%\texttt{s(l)} , \mathrm{and} \\ {\color{Cerulean}A_{\zeta ,o,i,l}} &\equiv& \texttt{Azo(v,0,j)}\%\texttt{s(l)} ; \end{eqnarray}

    where \(v=1,\texttt{Mvol}\) labels volume, \(j=1,\texttt{mn}\) labels Fourier harmonic, and \(l=0,\,\)Lrad \((v)\) labels Chebyshev polynomial. (These arrays also contains derivative information.)
  • If Linitguess=1 , a guess for the initial state for the Beltrami fields is constructed. An initial state is required for iterative solvers of the Beltrami fields, see LBeltrami .
  • If Linitguess=2 , the initial state for the Beltrami fields is read from file (see ra00aa() ). An initial state is required for iterative solvers of the Beltrami fields, see LBeltrami .

workspace

goomne, goomno: metric information These are defined in metrix() , and used in ma00aa().

gssmne, gssmno: metric information These are defined in metrix() , and used in ma00aa().

gstmne, gstmno: metric information These are defined in metrix() , and used in ma00aa().

gszmne, gszmno: metric information These are defined in metrix() , and used in ma00aa().

gttmne, gttmno: metric information These are defined in metrix() , and used in ma00aa().

gtzmne, gtzmno: metric information These are defined in metrix() , and used in ma00aa().

gzzmne, gzzmno: metric information These are defined in metrix() , and used in ma00aa().

cosi(1:Ntz,1:mn) and sini(1:Ntz,1:mn)

  • Trigonometric factors used in various Fast Fourier transforms, where

    \begin{eqnarray} \texttt{cosi}_{j,i} & = & \cos( m_i \theta_j - n_i \zeta_j ), \\ \texttt{sini}_{j,i} & = & \sin( m_i \theta_j - n_i \zeta_j ). \end{eqnarray}

psifactor(1:mn,1:Mvol): coordinate "pre-conditioning" factor

  • In toroidal geometry, the coordinate "pre-conditioning" factor is

    \begin{eqnarray} f_{j,v} \equiv \left\{ \begin{array}{lcccccc}\psi_{t,v}^{0 }&,&\mbox{for $m_j=0$}, \\ \psi_{t,v}^{m_j/2}&,&\mbox{otherwise}. \end{array}\right. \end{eqnarray}

    where \(\psi_{t,v} \equiv\,\)tflux is the (normalized?) toroidal flux enclosed by the \(v\)-th interface.
  • psifactor is used in packxi(), dforce() and hesian().
  • inifactor is similarly constructed, with

    \begin{eqnarray} f_{j,v} \equiv \left\{ \begin{array}{lcccccc}\psi_{t,v}^{ 1 /2}&,&\mbox{for $m_j=0$}, \\ \psi_{t,v}^{m_j/2}&,&\mbox{otherwise}. \end{array}\right. \end{eqnarray}

    and used only for the initialization of the surfaces taking into account axis information if provided.

Bsupumn and Bsupvmn

diotadxup and glambda: transformation to straight fieldline angle

  • Given the Beltrami fields in any volume, the rotational-transform on the adjacent interfaces may be determined (in tr00ab()) by constructing the straight fieldline angle on the interfaces.
  • The rotational transform on the inner or outer interface of a given volume depends on the magnetic field in that volume, i.e. \({{\,\iota\!\!\!}-}_\pm = {{\,\iota\!\!\!}-}({\bf B}_\pm)\), so that

    \begin{eqnarray} \delta {{\,\iota\!\!\!}-}_\pm = \frac{\partial {{\,\iota\!\!\!}-}_\pm}{\partial {\bf B}_\pm} \cdot \delta {\bf B_\pm}. \end{eqnarray}

  • The magnetic field depends on the Fourier harmonics of both the inner and outer interface geometry (represented here as \(x_j\)), the helicity multiplier, and the enclosed poloidal flux, i.e. \({\bf B_\pm} = {\bf B_\pm}(x_j, \mu, \Delta \psi_p)\), so that

    \begin{eqnarray} \delta {\bf B_\pm} = \frac{\partial {\bf B}_\pm}{\partial x_j } \delta x_j + \frac{\partial {\bf B}_\pm}{\partial \mu } \delta \mu + \frac{\partial {\bf B}_\pm}{\partial \Delta \psi_p} \delta \Delta \psi_p. \end{eqnarray}

  • The rotational-transforms, thus, can be considered to be functions of the geometry, the helicity-multiplier and the enclosed poloidal flux, \({{\,\iota\!\!\!}-}_{\pm} = {{\,\iota\!\!\!}-}_{\pm}(x_j,\mu,\Delta\psi_p)\).
  • The rotational-transform, and its derivatives, on the inner and outer interfaces of each volume is stored in diotadxup(0:1,-1:2,1:Mvol) . The indices label:
    • the first index labels the inner or outer interface,
    • the the second one labels derivative, with
      • -1 : indicating the derivative with respect to the interface geometry, i.e. \(\displaystyle \frac{\partial {{\,\iota\!\!\!}-}_{\pm}}{\partial x_j}\),
      • 0 : the rotational-transform itself,
      • 1,2 : the derivatives with respec to \(\mu\) and \(\Delta \psi_p\), i.e. \(\displaystyle \frac{\partial {{\,\iota\!\!\!}-}_{\pm}}{\partial \mu}\) and \(\displaystyle \frac{\partial {{\,\iota\!\!\!}-}_{\pm}}{\partial \Delta \psi_p}\);
    • The third index labels volume.
  • The values of diotadxup are assigned in mp00aa() after calling tr00ab().

vvolume, lBBintegral and lABintegral

  • volume integrals

    \begin{eqnarray} \texttt{vvolume(i)} &=& \int_{{\cal V}_i} \, dv\\ \texttt{lBBintegral(i)} &=& \int_{{\cal V}_i} {\bf B} \cdot {\bf B} \, dv\\ \texttt{lABintegral(i)} &=& \int_{{\cal V}_i} {\bf A} \cdot {\bf B} \, dv \end{eqnarray}

References allglobal::ajk, allglobal::ate, allglobal::ato, allglobal::aze, allglobal::azo, allglobal::bbe, allglobal::bbo, allglobal::bbweight, allglobal::beltramierror, allglobal::bemn, allglobal::bloweremn, allglobal::bloweromn, inputlist::bnc, bnorml(), inputlist::bns, allglobal::bomn, allglobal::bsupumn, allglobal::bsupvmn, allglobal::btemn, allglobal::btomn, allglobal::bzemn, allglobal::bzomn, allglobal::cfmn, allglobal::cheby, allglobal::comn, coords(), allglobal::cosi, fftw_interface::cplxin, fftw_interface::cplxout, allglobal::cpus, allglobal::diotadxup, allglobal::ditgpdxtp, allglobal::djkm, allglobal::djkp, allglobal::dpflux, allglobal::dradr, allglobal::dradz, allglobal::drbc, allglobal::drbs, allglobal::drij, allglobal::drodr, allglobal::drodz, allglobal::dtflux, allglobal::dxyz, allglobal::dzadr, allglobal::dzadz, allglobal::dzbc, allglobal::dzbs, allglobal::dzij, allglobal::dzodr, allglobal::dzodz, allglobal::efmn, inputlist::escale, allglobal::evmn, inputlist::forcetol, allglobal::fse, allglobal::fso, allglobal::gaussianabscissae, allglobal::gaussianweight, get_cheby(), get_zernike(), get_zernike_rm(), allglobal::glambda, allglobal::gmreslastsolution, constants::goldenmean, allglobal::goomne, allglobal::goomno, allglobal::gssmne, allglobal::gssmno, allglobal::gstmne, allglobal::gstmno, allglobal::gszmne, allglobal::gszmno, allglobal::gteta, allglobal::gttmne, allglobal::gttmno, allglobal::gtzmne, allglobal::gtzmno, allglobal::guvij, allglobal::gvuij, allglobal::gzeta, allglobal::gzzmne, allglobal::gzzmno, constants::half, allglobal::halfmm, inputlist::helicity, allglobal::hnt, allglobal::hnz, allglobal::ibnc, allglobal::ibns, allglobal::iemn, inputlist::igeometry, allglobal::iie, allglobal::iio, allglobal::ijimag, allglobal::ijreal, allglobal::im, allglobal::imagneticok, allglobal::ime, inputlist::impol, allglobal::ims, allglobal::in, allglobal::ine, allglobal::inifactor, allglobal::ins, inputlist::intor, allglobal::iomn, inputlist::iota, allglobal::iotakadd, allglobal::iotakkii, allglobal::iotaksgn, allglobal::iotaksub, allglobal::ipdtdpf, allglobal::iquad, allglobal::irbc, allglobal::irbs, allglobal::irij, inputlist::istellsym, allglobal::ivnc, allglobal::ivns, inputlist::ivolume, allglobal::izbc, allglobal::izbs, allglobal::izij, allglobal::jiimag, allglobal::jireal, allglobal::jkimag, allglobal::jkreal, allglobal::jxyz, allglobal::ki, allglobal::kija, allglobal::kijs, allglobal::kjimag, allglobal::kjreal, allglobal::labintegral, allglobal::lbbintegral, inputlist::lbeltrami, allglobal::lblinear, allglobal::lbnewton, allglobal::lbsequad, inputlist::lcheck, inputlist::lconstraint, allglobal::lcoordinatesingularity, inputlist::lfindzero, inputlist::lfreebound, allglobal::lgdof, inputlist::lgmresprec, allglobal::lhessianallocated, inputlist::lhevalues, inputlist::lhevectors, inputlist::lhmatrix, allglobal::liluprecond, inputlist::linitgues, inputlist::linitialize, allglobal::lma, inputlist::lmatsolver, allglobal::lmavalue, allglobal::lmb, allglobal::lmbvalue, allglobal::lmc, allglobal::lmcvalue, allglobal::lmd, allglobal::lmdvalue, allglobal::lme, allglobal::lmevalue, allglobal::lmf, allglobal::lmfvalue, allglobal::lmg, allglobal::lmgvalue, allglobal::lmh, allglobal::lmhvalue, allglobal::lmns, allglobal::lmpol, allglobal::lntor, allglobal::localconstraint, inputlist::lp, inputlist::lperturbed, inputlist::lq, inputlist::lrad, inputlist::lreflect, matrix(), inputlist::maxrndgues, allglobal::mmpp, allglobal::mn, allglobal::mne, allglobal::mns, inputlist::mpol, inputlist::mregular, inputlist::mu, constants::mu0, allglobal::myid, allglobal::nadof, inputlist::ndiscrete, allglobal::ndmas, allglobal::ndmasmax, allglobal::nfielddof, inputlist::nfp, allglobal::ngdof, allglobal::notmatrixfree, allglobal::notstellsym, inputlist::nppts, inputlist::nquad, allglobal::nt, inputlist::ntor, allglobal::ntz, inputlist::nvol, allglobal::nxyz, allglobal::nz, allglobal::odmn, allglobal::ofmn, inputlist::oita, constants::one, inputlist::opsilon, fileunits::ounit, inputlist::pcondense, allglobal::pemn, inputlist::pflux, inputlist::phiedge, constants::pi2, inputlist::pl, fftw_interface::planb, fftw_interface::planf, allglobal::pomn, inputlist::pr, preset(), allglobal::psifactor, inputlist::ql, inputlist::qr, constants::quart, ra00aa(), inputlist::rac, inputlist::ras, inputlist::rbc, inputlist::rbs, allglobal::regumm, allglobal::rij, inputlist::rp, inputlist::rq, allglobal::rscale, allglobal::rtm, allglobal::rtt, inputlist::rwc, inputlist::rws, rzaxis(), allglobal::semn, allglobal::sfmn, allglobal::sg, allglobal::simn, allglobal::sini, numerical::small, allglobal::smpol, allglobal::sntor, allglobal::somn, allglobal::sontz, numerical::sqrtmachprec, allglobal::sweight, inputlist::tflux, allglobal::trij, allglobal::tt, allglobal::tzij, inputlist::upsilon, inputlist::vnc, inputlist::vns, numerical::vsmall, allglobal::vvolume, inputlist::wpoloidal, allglobal::yesstellsym, inputlist::zac, inputlist::zas, inputlist::zbc, inputlist::zbs, allglobal::zernike, constants::zero, allglobal::zij, inputlist::zwc, and inputlist::zws.

Referenced by hesian(), preset(), and xspech().

Here is the call graph for this function:
Here is the caller graph for this function: