API Reference

Python

The booz_xform module provides a class Booz_xform, described in detail here. The same python module also includes several functions for plotting, described on the Plotting page.

The class Booz_xform provides a small number of functions to load input data, drive the calculation, and save output data. In addition, all the input and output data are available as properties of the class. The properties are scalars, 1D numpy arrays, and 2D numpy arrays.

Properties of the class are labeled below as being an input or output. The input quantities can be either read or written to in python, and should be set before calling run(). (Many of these properties can be set by calling read_wout()). The output quantities are read-only, and they are populated when run() is called.

class Booz_xform

A class holding the input and output information for a transformation to Boozer coordinates

init_from_vmec(self: booz_xform._booz_xform.Booz_xform, ns: int, iotas: numpy.ndarray[numpy.float64[m, 1]], rmnc0: numpy.ndarray[numpy.float64[m, n]], rmns0: numpy.ndarray[numpy.float64[m, n]], zmnc0: numpy.ndarray[numpy.float64[m, n]], zmns0: numpy.ndarray[numpy.float64[m, n]], lmnc0: numpy.ndarray[numpy.float64[m, n]], lmns0: numpy.ndarray[numpy.float64[m, n]], bmnc0: numpy.ndarray[numpy.float64[m, n]], bmns0: numpy.ndarray[numpy.float64[m, n]], bsubumnc0: numpy.ndarray[numpy.float64[m, n]], bsubumns0: numpy.ndarray[numpy.float64[m, n]], bsubvmnc0: numpy.ndarray[numpy.float64[m, n]], bsubvmns0: numpy.ndarray[numpy.float64[m, n]], phips: numpy.ndarray[numpy.float64[m, 1]] = array([], dtype=float64)) None

Handle conversion of the radial grids used in vmec to to the radial grid used by booz_xform. This involves truncation of the leading 0 for quantities on vmec’s half grid, and interpolation for quantities on vmec’s full grid. This function is called by the higher level function read_wout(), so usually users do not need to call this function directly.

Parameters:
  • ns – The number of radial vmec surfaces.

  • iotas – Iota on vmec’s half grid.

  • rmnc0 – Vmec’s original rmnc array, on the full grid.

  • rmns0 – Vmec’s original rmns array, on the full grid. For stellarator-symmetric configurations this array is ignored and need not be specified.

  • zmnc0 – Vmec’s original zmnc array, on the full grid. For stellarator-symmetric configurations this array is ignored and need not be specified.

  • zmns0 – Vmec’s original zmns array, on the full grid.

  • lmnc0 – Vmec’s original lmnc array, on the half grid. For stellarator-symmetric configurations this array is ignored and need not be specified.

  • lmns0 – Vmec’s original lmns array, on the half grid.

  • bmnc0 – Vmec’s original bmnc array, on the half grid.

  • bmns0 – Vmec’s original bmns array, on the half grid. For stellarator-symmetric configurations this array is ignored and need not be specified.

  • bsubumnc0 – Vmec’s original bsubumnc array, on the half grid.

  • bsubumns0 – Vmec’s original bsubumns array, on the half grid. For stellarator-symmetric configurations this array is ignored and need not be specified.

  • bsubvmnc0 – Vmec’s original bsubvmnc array, on the half grid.

  • bsubvmns0 – Vmec’s original bsubvmns array, on the half grid. For stellarator-symmetric configurations this array is ignored and need not be specified.

  • phip – Phip on vmec’s full grid. (Defaults to empty Vector)

read_boozmn(self: booz_xform._booz_xform.Booz_xform, filename: str) None

Read in the results of an earlier transformation from a classic boozmn_*.nc NetCDF file.

Parameters:

filename – The full name of the file to load.

read_wout(self: booz_xform._booz_xform.Booz_xform, filename: str, flux: bool = False) None

Read input information from a VMEC wout*.nc file.

Parameters:
  • filename – The full name of the file to load.

  • flux – If true, the poloidal flux is read. (Defaults to false)

run(self: booz_xform._booz_xform.Booz_xform) None

Run the transformation to Boozer coordinates on all the selected flux surfaces. The input arrays should be initialized before this step.

write_boozmn(self: booz_xform._booz_xform.Booz_xform, filename: str) None

Write the results of the transformation to a classic boozmn_*.nc NetCDF file. This function should only be called after booz_xform.Booz_xform.run().

Parameters:

filename – The full name of the file to save.

property Boozer_G

(1D float array of length ns_b, output) Coefficient of grad zeta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on the magnetic surfaces used for output quantities.

property Boozer_G_all

(1D float array of length ns_in, output) Coefficient of grad zeta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on all the magnetic surfaces for which input data was provided.

property Boozer_I

(1D float array of length ns_b, output) Coefficient of grad theta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on the magnetic surfaces used for output quantities.

property Boozer_I_all

(1D float array of length ns_in, output) Coefficient of grad theta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on all the magnetic surfaces for which input data was provided.

property aspect

(float, input) The aspect ratio of the configuration. This value is not used for anything by booz_xform, and does not need to be set. It is provided as a means to pass this value from the input equilibrium to booz_xform output files.

property asym

(bool, input) True if the configuration is not stellarator-symmetric, false if the configuration is stellarator-symmetric.

property bmnc

(2D float array of size mnmax_nyq x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of the magnetic field strength B.

property bmnc_b

(2D float array of size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes of the magnetic field strength in Boozer coordinates.

property bmns

(2D float array of size mnmax_nyq x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of the magnetic field strength B. For stellarator-symmetric configurations, this array is not used and need not be specified.

property bmns_b

(2D float array of size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes of the magnetic field strength in Boozer coordinates. If the configuration is stellarator-symmetric, this quantity is zero so the array will have size 0 x 0.

property bsubumnc

(2D float array of size mnmax_nyq x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d theta_0) where r is the position vector.

property bsubumns

(2D float array of size mnmax_nyq x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d theta_0) where r is the position vector. For stellarator-symmetric configurations, this array is not used and need not be specified.

property bsubvmnc

(2D float array of size mnmax_nyq x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d zeta_0) where r is the position vector.

property bsubvmns

(2D float array of size mnmax_nyq x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d zeta_0) where r is the position vector. For stellarator-symmetric configurations, this array is not used and need not be specified.

property compute_surfs

(1D integer array, input) Indices of ns_in-sized radial arrays, specifying the flux surfaces for which the transformation to Boozer coordinates will be performed. All values should be >= 0 and < ns_in. The array compute_surfs is similar to the array jlist in the earlier fortran booz_xform program, with compute_surfs = jlist - 2.

property gmnc_b

(2D float array of size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Jacobian of (psi, theta_B, zeta_B) coordinates.

property gmns_b

(2D float array of size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Jacobian of (psi, theta_B, zeta_B) coordinates. If the configuration is stellarator-symmetric, this quantity is zero so this array will have size 0 x 0.

property iota

(1D float array of length ns_in, input) Rotational transform on the input radial surfaces.

property lmnc

(2D float array of size mnmax_nyq x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of lambda = theta^* - theta_0, the difference between the original poloidal angle theta_0 and the straight field line poloidal angle theta^*. For stellarator-symmetric configurations, this array is not used and need not be specified.

property lmns

(2D float array of size mnmax_nyq x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of lambda = theta^* - theta_0, the difference between the original poloidal angle theta_0 and the straight field line poloidal angle theta^*.

property mboz

(int, input) Maximum poloidal mode number for representing output quantities in Boozer coordinates.

property mnboz

(int, output) Total number of Fourier modes for output data.

property mnmax

(int, input) Number of Fourier modes for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

property mnmax_nyq

(int, input) Total number of Fourier modes for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

property mpol

(int, input) Maximum poloidal mode number for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

property mpol_nyq

(int, input) Maximum poloidal mode number for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

property nboz

(int, input) Maximum toroidal mode number (divided by nfp) for representing output quantities in Boozer coordinates. For example, if nboz=2 and nfp=10, the toroidal modes used will be n=-20, -10, 0, 10, 20.

property nfp

(int, input) Number of field periods, i.e. the discrete toroidal rotation symmetry.

property ns_b

(int, output) Number of flux surfaces on which output data are available.

property ns_in

(int, input) Number of flux surfaces on which the input data is supplied. The transformation to Boozer coordinates is not necessarily run on all of these surfaces, only the ones indicated by compute_surfs.

property ntor

(int, input) Maximum toroidal mode number (divided by nfp) for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

property ntor_nyq

(int, input) Maximum toroidal mode number (divided by nfp) for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

property numnc_b

(2D float array of size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the toroidal angle difference nu = zeta_B - zeta_0. If the configuration is stellarator-symmetric, this quantity is zero so this array will have size 0 x 0.

property numns_b

(2D float array of size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the toroidal angle difference nu = zeta_B - zeta_0.

property phi

(1D float array of length ns_in, input) Toroidal flux normalized by 2*pi, evaluated on all the magnetic surfaces for which input data was provided.

property phip

(1D float array of length ns_in, input) Poloidal flux normalized by 2*pi, evaluated on all the magnetic surfaces for which input data was provided.

property rmnc

(2D float array of size mnmax x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of the major radius R.

property rmnc_b

(2D float array of size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the major radius R of the flux surfaces.

property rmns

(2D float array of size mnmax x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of the major radius R of the flux surfaces. For stellarator-symmetric configurations, this array is not used and need not be specified.

property rmns_b

(2D float array of size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the major radius R of the flux surfaces. If the configuration is stellarator-symmetric, this quantity is zero so this array will have size 0 x 0.

property s_b

(1D float array of length ns_b, output) Values of normalized toroidal flux s defining the magnetic surfaces for the output data.

property s_in

(1D float array of length ns_in, input) The values of normalized toroidal flux s for the input data. Here, s is the toroidal flux normalized to the value at the plasma boundary. This information is not needed for the coordinate transformation itself, but is useful for plotting output.

property toroidal_flux

(float, input) The boundary toroidal flux of the configuration. This value is not used for anything by booz_xform, and does not need to be set. It is provided as a means to pass this value from the input equilibrium to booz_xform output files.

property verbose

(int, input) Set this to 0 for no output to stdout, 1 for some output, 2 for lots of output

property xm

(1D integer array of length mnmax, input) The poloidal Fourier mode numbers for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

property xm_b

(1D int array of length mnboz, output) Poloidal mode numbers for the output data.

property xm_nyq

(1D integer array of length mnmax_nyq, input) The poloidal Fourier mode numbers for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

property xn

(1D integer array of length mnmax, input) The toroidal Fourier mode numbers for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns. The values should all be integer multiples of nfp.

property xn_b

(1D int array of length mnboz, output) Toroidal mode numbers for the output data. These values are all integer multiples of nfp.

property xn_nyq

(1D integer array of length mnmax_nyq, input) The toroidal Fourier mode numbers for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns. The values should all be integer multiples of nfp.

property zmnc

(2D float array of size mnmax x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of the Cartesian coordinate Z of the flux surfaces. For stellarator-symmetric configurations, this array is not used and need not be specified.

property zmnc_b

(2D float array of size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Cartesian coordinate Z of the flux surfaces. If the configuration is stellarator-symmetric, this quantity is zero so array will have size 0 x 0.

property zmns

(2D float array of size mnmax x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of the Cartesian coordinate Z of the flux surfaces.

property zmns_b

(2D float array of size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Cartesian coordinate Z of the flux surfaces.

Note

While 1D and 2D array parameters can appear on the left-hand side of assignment operations, individual array elements or array slices cannot. So for instance, if b is an instance of booz_xform.Booz_Xform, b.rmnc = [[1, 0.1],[1.1, 0.1]] will behave as expected, but b.rmnc[0,0] = 2.0 will not. This is due to the copying of data when interfacing python to C++.

C++

All definitions for the booz_xform code reside in the namespace booz_xform::, which will be omitted for brevity here.

1D and 2D arrays are handled using the Eigen library.

All floating point values have the type boozfloat, which is merely a typedef for double. All 1D arrays of integers (such as xm) have the type IntVector, which is merely a typedef for Eigen::ArrayXi. All 1D arrays of floating point values have the type Vector, which is merely a typedef for Eigen::ArrayXd. 2D arrays of floating point values have the type Matrix, which is a typedef for Eigen::ArrayXXd. If desired, the entire code could be switched to single precision by changing these typedef statements, which are defined in vector_matrix.hpp.

Public variables below are labeled below as being an input or output. The input quantities should be set before calling run(). (Many of these input quantities can be set by calling read_wout()). The output quantities are populated when run() is called.

The main functionality of the booz_xform code is contained in a single class, Booz_xform:

class Booz_xform

Public Functions

Booz_xform()

Constructor.

Create a Booz_xform object with no data.

void read_wout(std::string filename, bool flux = false)

Read input data from a VMEC wout_*.nc file.

This method also calls init_from_vmec().

Parameters:

filename[in] The name of the VMEC wout file to load

void init_from_vmec(int ns, Vector &iotas, Matrix &rmnc, Matrix &rmns, Matrix &zmnc, Matrix &zmns, Matrix &lmnc, Matrix &lmns, Matrix &bmnc, Matrix &bmns, Matrix &bsubumnc, Matrix &bsubumns, Matrix &bsubvmnc, Matrix &bsubvmns, Vector &phips = defaultInitPtr)

Handle radial dimension of vmec arrays.

This method also handles radial interpolation of the full-grid quantities rmnc, rmns, zmnc, and zmns onto the half-grid points. It also discards the first radial index for half-grid quantities, which is populated with zeros in vmec.

The input parameters to this function are the variables of the same name in a VMEC wout file.

void run()

Carry out the transformation calculation.

This method is the main computationally intensive step, and it should be called after the equilibrium data is loaded.

void write_boozmn(std::string filename)

Write results to a boozmn_*.nc output file.

This method writes a classic output file, and it should be called after run() has completed.

Parameters:

filename[in] The full name of the boozmn_*.nc file to write.

void read_boozmn(std::string filename)

Read previously calculated results from a boozmn_*.nc output file.

Parameters:

filename[in] The full name of the boozmn_*.nc file to read.

void init()

Initialize variables that are common to all surfaces that share the same (mpol, ntor) resolution.

In the original fortran code, this work was done in read_wout_booz.f, setup_booz.f, boozer_coords.f, and foranl.f.

void surface_solve(int js_b)

Execute the transformation to Boozer coordinates for a single flux surface.

This subroutine corresponds to boozer_coords.f in the fortran version.

Parameters:

js_b[in] The index into compute_surfs indicating the surface on which the Boozer transformation will be calculated.

Public Members

int verbose

(input) Set this to 0 for no output to stdout, 1 for some output, 2 for lots of output.

bool asym

(input) false if the configuration is stellarator-symmetric, true otherwise.

int nfp

(input) Number of field periods, i.e. the discrete toroidal rotation symmetry.

int mpol

(input) Maximum poloidal mode number for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

int ntor

(input) Maximum toroidal mode number (divided by nfp) for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

int mnmax

(input) mnmax Number of Fourier modes for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

int mpol_nyq

(input) Maximum poloidal mode number for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

int ntor_nyq

(input) Maximum toroidal mode number (divided by nfp) for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

int mnmax_nyq

(input) Total number of Fourier modes for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

IntVector xm

(size mnmax, input) The poloidal Fourier mode numbers for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns.

IntVector xn

(size mnmax, input) The toroidal Fourier mode numbers for the input arrays rmnc, rmns, zmnc, zmns, lmnc, and lmns. The values should all be integer multiples of nfp.

IntVector xm_nyq

(size mnmax_nyq, input) The poloidal Fourier mode numbers for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns.

IntVector xn_nyq

(size mnmax_nyq, input) The toroidal Fourier mode numbers for the input arrays bmnc, bmns, bsubumnc, bsubumns, bsubvmnc, and bsubvmns. The values should all be integer multiples of nfp.

int ns_in

(input) The number of flux surfaces on which input data is supplied. The transformation to Boozer coordinates is not necessarily run on all of these surfaces, only the ones indicated by compute_surfs.

Vector s_in

(size ns_in, input) Values of normalized toroidal flux for which the input data is stored. These numbers are used only for plotting.

Vector iota

(size ns_in, input) Rotational transform on the input radial surfaces.

Matrix rmnc

(size mnmax x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of the major radius R.

Matrix rmns

(size mnmax x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of the major radius R. For stellarator-symmetric configurations, this array is not used and need not be specified.

Matrix zmnc

(size mnmax x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of the Cartesian coordinate Z of the flux surfaces. For stellarator-symmetric configurations, this array is not used and need not be specified.

Matrix zmns

(size mnmax x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of the Cartesian coordinate Z of the flux surfaces.

Matrix lmnc

(size mnmax x ns_in, input) cos(m * theta_0 - n zeta_0) Fourier modes of lambda = theta^* - theta_0, the difference between the original poloidal angle theta_0 and the straight field line poloidal angle theta^*. For stellarator-symmetric configurations, this array is not used and need not be specified.

Matrix lmns

(size mnmax x ns_in, input) sin(m * theta_0 - n zeta_0) Fourier modes of lambda = theta^* - theta_0, the difference between the original poloidal angle theta_0 and the straight field line poloidal angle theta^*.

Matrix bmnc

(size mnmax_nyq x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of the magnetic field strength B.

Matrix bmns

(size mnmax_nyq x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of the magnetic field strength B. For stellarator-symmetric configurations, this array is not used and need not be specified.

Matrix bsubumnc

(size mnmax_nyq x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d theta_0) where r is the position vector.

Matrix bsubumns

(size mnmax_nyq x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d theta_0) where r is the position vector. For stellarator-symmetric configurations, this array is not used and need not be specified.

Matrix bsubvmnc

(size mnmax_nyq x ns_in, input) cos(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d zeta_0) where r is the position vector.

Matrix bsubvmns

(size mnmax_nyq x ns_in, input) sin(m * theta_0 - n * zeta_0) Fourier modes of B dot (d r / d zeta_0) where r is the position vector. For stellarator-symmetric configurations, this array is not used and need not be specified.

int mboz

(input) Maximum poloidal mode number for representing output quantities in Boozer coordinates.

int nboz

(input) Maximum toroidal mode number (divided by nfp) for representing output quantities in Boozer coordinates. For example, if nboz=2 and nfp=10, the toroidal modes used will be n=-20, -10, 0, 10, 20.

IntVector compute_surfs

(input) Indices of ns_in-sized radial arrays, specifying the flux surfaces for which the transformation to Boozer coordinates will be performed. All values should be >= 0 and < ns_in. The array compute_surfs is similar to the array jlist in the earlier fortran booz_xform program, with compute_surfs = jlist - 2.

boozfloat aspect

(input) The aspect ratio of the configuration. This value is not used for anything by booz_xform, and does not need to be set. It is provided as a means to pass this value from the input equilibrium to booz_xform output files.

boozfloat toroidal_flux

(input) The boundary toroidal flux of the configuration (not divided by (2*pi)). This value is not used for anything by booz_xform, and does not need to be set. It is provided as a means to pass this value from the input equilibrium to booz_xform output files.

int ns_b

(output) Number of surfaces on which the transformation is calculated.

Vector s_b

(size ns_b, output) Values of normalized toroidal flux for which the output data is stored. These numbers are used only for plotting.

int mnboz

(output) Total number of Fourier modes for output data.

IntVector xm_b

(size mnboz, output) Poloidal mode numbers used for the Fourier representation of output quantities, i.e. functions of the Boozer angles.

IntVector xn_b

(size mnboz, output) Toroidal mode numbers used for the Fourier representation of output quantities, i.e. functions of the Boozer angles.

Matrix bmnc_b

(size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes of the magnetic field strength in Boozer coordinates.

Matrix bmns_b

(size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes of the magnetic field strength in Boozer coordinates. If the configuration is stellarator-symmetric, this quantity is zero so the array will have size 0 x 0.

Matrix gmnc_b

(size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Jacobian of (psi, theta_B, zeta_B) coordinates.

Matrix gmns_b

(size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Jacobian of (psi, theta_B, zeta_B) coordinates. If the configuration is stellarator-symmetric, this quantity is zero so this array will have size 0 x 0.

Matrix rmnc_b

(size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the major radius R of the flux surfaces.

Matrix rmns_b

(size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the major radius R of the flux surfaces. If the configuration is stellarator-symmetric, this quantity is zero so this array will have size 0 x 0.

Matrix zmnc_b

(size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Cartesian coordinate Z of the flux surfaces. If the configuration is stellarator-symmetric, this quantity is zero so array will have size 0 x 0.

Matrix zmns_b

(size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the Cartesian coordinate Z of the flux surfaces.

Matrix numnc_b

(size mnboz x ns_b, output) cos(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the toroidal angle difference nu = zeta_B - zeta_0. If the configuration is stellarator-symmetric, this quantity is zero so this array will have size 0 x 0.

Matrix numns_b

(size mnboz x ns_b, output) sin(m * theta_B - n * zeta_B) Fourier modes (with respect to Boozer coordinates) of the toroidal angle difference nu = zeta_B - zeta_0.

Vector Boozer_G

(size ns_b, output) Coefficient of grad zeta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on the magnetic surfaces used for output quantities.

Vector Boozer_G_all

(size ns_in, output) Coefficient of grad zeta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on all the magnetic surfaces for which input data was provided.

Vector Boozer_I

(size ns_b, output) Coefficient of grad theta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on the magnetic surfaces used for output quantities.

Vector Boozer_I_all

(size ns_in, output) Coefficient of grad theta_B in the covariant representation of the magnetic field vector B in Boozer coordinates, evaluated on all the magnetic surfaces for which input data was provided.

Vector phi

(size ns_in + 1, output) Uniformly spaced grid going from 0 to the boundary toroidal flux (not divided by (2*pi)), evaluated on full vmec grid.

Vector phip

(size ns_in + 1, output) Uniformly spaced grid going from 0 to the boundary poloidal flux (not divided by (2*pi)), evaluated on full vmec grid.