SPEC 3.20
Stepped Pressure Equilibrium Code

The namelist physicslist controls the geometry, profiles, and numerical resolution. More...

Collaboration diagram for physicslist:

Variables

integer inputlist::igeometry = 3
 selects Cartesian, cylindrical or toroidal geometry; More...
 
integer inputlist::istellsym = 1
 stellarator symmetry is enforced if Istellsym==1
 
integer inputlist::lfreebound = 0
 compute vacuum field surrounding plasma
 
real inputlist::phiedge = 1.0
 total enclosed toroidal magnetic flux;
 
real inputlist::curtor = 0.0
 total enclosed (toroidal) plasma current;
 
real inputlist::curpol = 0.0
 total enclosed (poloidal) linking current;
 
real inputlist::gamma = 0.0
 adiabatic index; cannot set \(|\gamma| = 1\)
 
integer inputlist::nfp = 1
 field periodicity More...
 
integer inputlist::nvol = 1
 number of volumes More...
 
integer inputlist::mpol = 0
 number of poloidal Fourier harmonics More...
 
integer inputlist::ntor = 0
 number of toroidal Fourier harmonics More...
 
integer, dimension(1:mnvol+1) inputlist::lrad = 4
 Chebyshev resolution in each volume. More...
 
integer inputlist::lconstraint = -1
 selects constraints; primarily used in ma02aa() and mp00ac(). More...
 
real, dimension(1:mnvol+1) inputlist::tflux = 0.0
 toroidal flux, \(\psi_t\), enclosed by each interface More...
 
real, dimension(1:mnvol+1) inputlist::pflux = 0.0
 poloidal flux, \(\psi_p\), enclosed by each interface
 
real, dimension(1:mnvol) inputlist::helicity = 0.0
 helicity, \({\cal K}\), in each volume, \({\cal V}_i\) More...
 
real inputlist::pscale = 0.0
 pressure scale factor More...
 
real, dimension(1:mnvol+1) inputlist::pressure = 0.0
 pressure in each volume More...
 
integer inputlist::ladiabatic = 0
 logical flag More...
 
real, dimension(1:mnvol+1) inputlist::adiabatic = 0.0
 adiabatic constants in each volume More...
 
real, dimension(1:mnvol+1) inputlist::mu = 0.0
 helicity-multiplier, \(\mu\), in each volume
 
real, dimension(1:mnvol+1) inputlist::ivolume = 0.0
 Toroidal current constraint normalized by \(\mu_0\) ( \(I_{volume} = \mu_0\cdot [A]\)), in each volume. This is a cumulative quantity: \(I_{\mathcal{V},i} = \int_0^{\psi_{t,i}} \mathbf{J}\cdot\mathbf{dS}\). Physically, it represents the sum of all non-pressure driven currents.
 
real, dimension(1:mnvol) inputlist::isurf = 0.0
 Toroidal current normalized by \(\mu_0\) at each interface (cumulative). This is the sum of all pressure driven currents.
 
integer, dimension(0:mnvol) inputlist::pl = 0
 "inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \). More...
 
integer, dimension(0:mnvol) inputlist::ql = 0
 "inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \). More...
 
integer, dimension(0:mnvol) inputlist::pr = 0
 "inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \). More...
 
integer, dimension(0:mnvol) inputlist::qr = 0
 "inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \). More...
 
real, dimension(0:mnvol) inputlist::iota = 0.0
 rotational-transform, \(\mbox{$\,\iota\!\!$-}\), on inner side of each interface More...
 
integer, dimension(0:mnvol) inputlist::lp = 0
 "outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\). More...
 
integer, dimension(0:mnvol) inputlist::lq = 0
 "outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\). More...
 
integer, dimension(0:mnvol) inputlist::rp = 0
 "outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\). More...
 
integer, dimension(0:mnvol) inputlist::rq = 0
 "outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\). More...
 
real, dimension(0:mnvol) inputlist::oita = 0.0
 rotational-transform, \(\mbox{$\,\iota\!\!$-}\), on outer side of each interface More...
 
real inputlist::mupftol = 1.0e-14
 accuracy to which \(\mu\) and \(\Delta\psi_p\) are required More...
 
integer inputlist::mupfits = 8
 an upper limit on the transform/helicity constraint iterations; More...
 
real inputlist::rpol = 1.0
 poloidal extent of slab (effective radius) More...
 
real inputlist::rtor = 1.0
 toroidal extent of slab (effective radius) More...
 
integer inputlist::lreflect = 0
 =1 reflect the upper and lower bound in slab, =0 do not reflect
 
real, dimension( 0:mntor) inputlist::rac = 0.0
 stellarator symmetric coordinate axis;
 
real, dimension( 0:mntor) inputlist::zas = 0.0
 stellarator symmetric coordinate axis;
 
real, dimension( 0:mntor) inputlist::ras = 0.0
 non-stellarator symmetric coordinate axis;
 
real, dimension( 0:mntor) inputlist::zac = 0.0
 non-stellarator symmetric coordinate axis;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::rbc = 0.0
 stellarator symmetric boundary components;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::zbs = 0.0
 stellarator symmetric boundary components;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::rbs = 0.0
 non-stellarator symmetric boundary components;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::zbc = 0.0
 non-stellarator symmetric boundary components;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::rwc = 0.0
 stellarator symmetric boundary components of wall;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::zws = 0.0
 stellarator symmetric boundary components of wall;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::rws = 0.0
 non-stellarator symmetric boundary components of wall;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::zwc = 0.0
 non-stellarator symmetric boundary components of wall;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::vns = 0.0
 stellarator symmetric normal field at boundary; vacuum component;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::bns = 0.0
 stellarator symmetric normal field at boundary; plasma component;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::vnc = 0.0
 non-stellarator symmetric normal field at boundary; vacuum component;
 
real, dimension(-mntor:mntor,-mmpol:mmpol) inputlist::bnc = 0.0
 non-stellarator symmetric normal field at boundary; plasma component;
 

Detailed Description

The namelist physicslist controls the geometry, profiles, and numerical resolution.

Variable Documentation

◆ igeometry

integer inputlist::igeometry = 3

selects Cartesian, cylindrical or toroidal geometry;

  • Igeometry=1 : Cartesian; geometry determined by \(R\);
  • Igeometry=2 : cylindrical; geometry determined by \(R\);
  • Igeometry=3 : toroidal; geometry determined by \(R\) and \(Z\);

Referenced by bnorml(), allglobal::broadcast_inputs(), allglobal::check_inputs(), coords(), dforce(), dfp100(), dfp200(), dvcfield(), evaluate_dbb(), evaluate_dmupfdx(), fcn1(), fcn2(), final_diagnostics(), hesian(), jo00aa(), lforce(), sphdf5::mirror_input_to_outfile(), newton(), packxi(), pc00ab(), pp00aa(), preset(), rzaxis(), spec(), stzxyz(), volume(), sphdf5::write_grid(), writereadgf(), and allglobal::wrtend().

◆ nfp

integer inputlist::nfp = 1

field periodicity

  • all Fourier representations are of the form \(\cos(m\theta-n N \zeta)\), \(\sin(m\theta-n N \zeta)\), where \(N\equiv\)Nfp
  • constraint: Nfp >= 1

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

◆ nvol

integer inputlist::nvol = 1

number of volumes

  • each volume \({\cal V}_l\) is bounded by the \({\cal I}_{l-1}\) and \({\cal I}_{l}\) interfaces
  • note that in cylindrical or toroidal geometry, \({\cal I}_{0}\) is the degenerate coordinate axis
  • constraint: Nvol<=MNvol

Referenced by brcast(), allglobal::broadcast_inputs(), allglobal::check_inputs(), df00ab(), dforce(), dfp100(), dfp200(), dvcfield(), evaluate_dbb(), evaluate_dmupfdx(), fcn1(), fcn2(), final_diagnostics(), sphdf5::hdfint(), hesian(), jo00aa(), lforce(), sphdf5::mirror_input_to_outfile(), newton(), packxi(), pc00ab(), pp00aa(), pp00ab(), preset(), spec(), stzxyz(), tr00ab(), volume(), wa00aa(), sphdf5::write_grid(), writereadgf(), and allglobal::wrtend().

◆ mpol

integer inputlist::mpol = 0

number of poloidal Fourier harmonics

  • all Fourier representations of doubly-periodic functions are of the form

    \begin{eqnarray} f(\theta,\zeta) & = & \sum_{n=0}^{\texttt{Ntor}} f_{0,n}\cos(-n \, \texttt{Nfp} \, \zeta) + \sum_{m=1}^{\texttt{Mpol}}\sum_{n=\texttt{-Ntor}}^{\texttt{Ntor}} f_{m,n}\cos(m\theta-n \, \texttt{Nfp} \, \zeta), \end{eqnarray}

    Internally these "double" summations are written as a "single" summation, e.g. \(f(\theta,\zeta) = \sum_j f_j \cos(m_j\theta-n_j\zeta)\).

Referenced by allocate_geometry_matrices(), bfield(), bfield_tangent(), allglobal::broadcast_inputs(), allglobal::check_inputs(), dfp200(), intghs(), intghs_workspace_init(), jo00aa(), ma00aa(), matrix(), sphdf5::mirror_input_to_outfile(), mtrxhs(), preset(), ra00aa(), spsint(), spsmat(), tr00ab(), wa00aa(), writereadgf(), and allglobal::wrtend().

◆ ntor

integer inputlist::ntor = 0

number of toroidal Fourier harmonics

  • all Fourier representations of doubly-periodic functions are of the form

    \begin{eqnarray} f(\theta,\zeta) & = & \sum_{n=0}^{\texttt{Ntor}} f_{0,n}\cos(-n \, \texttt{Nfp} \, \zeta) + \sum_{m=1}^{\texttt{Mpol}}\sum_{n=\texttt{-Ntor}}^{\texttt{Ntor}} f_{m,n}\cos(m\theta-n \, \texttt{Nfp} \, \zeta), \end{eqnarray}

    Internally these "double" summations are written as a "single" summation, e.g. \(f(\theta,\zeta) = \sum_j f_j \cos(m_j\theta-n_j\zeta)\).

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), coords(), dforce(), dfp200(), evaluate_dbb(), sphdf5::mirror_input_to_outfile(), mp00ac(), packxi(), preset(), ra00aa(), rzaxis(), stzxyz(), tr00ab(), wa00aa(), writereadgf(), and allglobal::wrtend().

◆ lrad

◆ lconstraint

integer inputlist::lconstraint = -1

selects constraints; primarily used in ma02aa() and mp00ac().

  • if Lconstraint==-1, then in the plasma regions \(\Delta\psi_t\), \(\mu\) and \(\Delta \psi_p\) are not varied and in the vacuum region (only for free-boundary) \(\Delta\psi_t\) and \(\Delta \psi_p\) are not varied, and \(\mu = 0\).
  • if Lconstraint==0, then in the plasma regions \(\Delta\psi_t\), \(\mu\) and \(\Delta \psi_p\) are not varied and in the vacuum region (only for free-boundary) \(\Delta\psi_t\) and \(\Delta \psi_p\) are varied to match the prescribed plasma current, curtor, and the "linking" current, curpol, and \(\mu = 0\)
  • if Lconstraint==1, then in the plasma regions \(\mu\) and \(\Delta\psi_p\) are adjusted in order to satisfy the inner and outer interface transform constraints (except in the simple torus, where the enclosed poloidal flux is irrelevant, and only \(\mu\) is varied to satisfy the outer interface transform constraint); and in the vacuum region \(\Delta\psi_t\) and \(\Delta \psi_p\) are varied to match the transform constraint on the boundary and to obtain the prescribed linking current, curpol, and \(\mu = 0\).
  • Todo:
    if Lconstraint==2, under reconstruction.

  • if Lconstraint.eq.3 , then the \(\mu\) and \(\psi_p\) variables are adjusted in order to satisfy the volume and surface toroidal current computed with lbpol() (excepted in the inner most volume, where the volume current is irrelevant). Not implemented yet in free boundary.

Referenced by brcast(), allglobal::broadcast_inputs(), allglobal::check_inputs(), dforce(), dfp100(), dfp200(), evaluate_dbb(), evaluate_dmupfdx(), get_lu_beltrami_matrices(), get_perturbed_solution(), ma02aa(), sphdf5::mirror_input_to_outfile(), mp00ac(), pp00aa(), preset(), spec(), and allglobal::wrtend().

◆ tflux

real, dimension(1:mnvol+1) inputlist::tflux = 0.0

toroidal flux, \(\psi_t\), enclosed by each interface

  • For each of the plasma volumes, this is a constraint: tflux is not varied
  • For the vacuum region (only if Lfreebound==1), tflux may be allowed to vary to match constraints
  • Note that tflux will be normalized so that tflux(Nvol) = 1.0, so that tflux is arbitrary up to a scale factor
  • See also
    phiedge

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

◆ helicity

real, dimension(1:mnvol) inputlist::helicity = 0.0

helicity, \({\cal K}\), in each volume, \({\cal V}_i\)

  • on exit, helicity is set to the computed values of \({\cal K} \equiv \int {\bf A}\cdot{\bf B}\;dv\)

Referenced by brcast(), allglobal::broadcast_inputs(), allglobal::check_inputs(), df00ab(), sphdf5::hdfint(), hesian(), ma02aa(), sphdf5::mirror_input_to_outfile(), mp00ac(), preset(), spec(), and allglobal::wrtend().

◆ pscale

real inputlist::pscale = 0.0

pressure scale factor

  • the initial pressure profile is given by pscale \(*\) pressure

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), dfp200(), evaluate_dbb(), lforce(), sphdf5::mirror_input_to_outfile(), spec(), volume(), and allglobal::wrtend().

◆ pressure

real, dimension(1:mnvol+1) inputlist::pressure = 0.0

pressure in each volume

  • The pressure is not held constant, but \(p_l V_l^\gamma = P_l\) is held constant, where \(P_l\) is determined by the initial pressures and the initial volumes, \(V_l\).
  • Note that if gamma==0.0, then \(p_l \equiv P_l\).
  • On output, the pressure is given by \(p_l = P_l/V_l^\gamma\), where \(V_l\) is the final volume.
  • pressure is only used in calculation of interface force-balance.

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

◆ ladiabatic

integer inputlist::ladiabatic = 0

logical flag

  • If Ladiabatic==0, the adiabatic constants are determined by the initial pressure and volume.
  • If Ladiabatic==1, the adiabatic constants are determined by the given input adiabatic.

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

◆ adiabatic

real, dimension(1:mnvol+1) inputlist::adiabatic = 0.0

adiabatic constants in each volume

  • The pressure is not held constant, but \(p_l V_l^\gamma = P_l \equiv\)adiabatic is constant.
  • Note that if gamma==0.0, then pressure==adiabatic.
  • pressure is only used in calculation of interface force-balance.

Referenced by allglobal::broadcast_inputs(), dfp200(), evaluate_dbb(), sphdf5::hdfint(), lforce(), sphdf5::mirror_input_to_outfile(), spec(), and allglobal::wrtend().

◆ pl

integer, dimension(0:mnvol) inputlist::pl = 0

"inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \).

If both \(q_l = 0\) and \(q_r = 0\), then the (inside) interface rotational-transform is defined by iota .

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

◆ ql

integer, dimension(0:mnvol) inputlist::ql = 0

"inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \).

If both \(q_l = 0\) and \(q_r = 0\), then the (inside) interface rotational-transform is defined by iota .

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

◆ pr

integer, dimension(0:mnvol) inputlist::pr = 0

"inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \).

If both \(q_l = 0\) and \(q_r = 0\), then the (inside) interface rotational-transform is defined by iota .

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

◆ qr

integer, dimension(0:mnvol) inputlist::qr = 0

"inside" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2 \).

If both \(q_l = 0\) and \(q_r = 0\), then the (inside) interface rotational-transform is defined by iota .

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

◆ iota

real, dimension(0:mnvol) inputlist::iota = 0.0

rotational-transform, \(\mbox{$\,\iota\!\!$-}\), on inner side of each interface

  • only relevant if illogical input for ql and qr are provided

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

◆ lp

integer, dimension(0:mnvol) inputlist::lp = 0

"outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\).

If both \(q_l = 0\) and \(q_r = 0\), then the (outer) interface rotational-transform is defined by oita .

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

◆ lq

integer, dimension(0:mnvol) inputlist::lq = 0

"outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\).

If both \(q_l = 0\) and \(q_r = 0\), then the (outer) interface rotational-transform is defined by oita .

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

◆ rp

integer, dimension(0:mnvol) inputlist::rp = 0

"outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\).

If both \(q_l = 0\) and \(q_r = 0\), then the (outer) interface rotational-transform is defined by oita .

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

◆ rq

integer, dimension(0:mnvol) inputlist::rq = 0

"outer" interface rotational-transform is \(\mbox{$\,\iota\!\!$-} = (p_l+\gamma p_r)/(q_l+\gamma q_r)\), where \(\gamma\) is the golden mean, \(\gamma = (1 + \sqrt 5 ) / 2\).

If both \(q_l = 0\) and \(q_r = 0\), then the (outer) interface rotational-transform is defined by oita .

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

◆ oita

real, dimension(0:mnvol) inputlist::oita = 0.0

rotational-transform, \(\mbox{$\,\iota\!\!$-}\), on outer side of each interface

  • only relevant if illogical input for ql and qr are provided

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

◆ mupftol

real inputlist::mupftol = 1.0e-14

accuracy to which \(\mu\) and \(\Delta\psi_p\) are required

  • only relevant if constraints on transform, enclosed currents etc. are to be satisfied iteratively, see Lconstraint

Referenced by allglobal::broadcast_inputs(), allglobal::check_inputs(), dforce(), evaluate_dmupfdx(), ma02aa(), sphdf5::mirror_input_to_outfile(), mp00ac(), and allglobal::wrtend().

◆ mupfits

integer inputlist::mupfits = 8

an upper limit on the transform/helicity constraint iterations;

  • only relevant if constraints on transform, enclosed currents etc. are to be satisfied iteratively, see Lconstraint
  • constraint: mupfits > 0

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

◆ rpol

real inputlist::rpol = 1.0

poloidal extent of slab (effective radius)

  • only relevant if Igeometry==1
  • poloidal size is \(L = 2\pi*\)rpol

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

◆ rtor

real inputlist::rtor = 1.0

toroidal extent of slab (effective radius)

  • only relevant if Igeometry==1
  • toroidal size is \(L = 2\pi*\)rtor

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