Postprocessing#

In this notebook we take a deeper look at postprocessing, using the python interface of GVEC.

GVEC is parallelised with OpenMP and we can set the number of OpenMP threads before importing the gvec package.

import os

os.environ["OMP_NUM_THREADS"] = "2"

This tutorial requires matplotlib to be installed in addition to gvec. numpy and xarray are dependencies of gvec and should therefore be already installed.

import matplotlib.pyplot as plt
import numpy as np
import xarray as xr

import gvec

In this tutorial we will use the equilibrium of the elliptic stellarator, computed in the tutorial Elliptic Stellarator(๐ŸŒ).

To obtain a gvec.State(๐ŸŒ) object for this equilibrium, gvec.load_state(๐ŸŒ) can be used with specific paths for the parameter and statefile. Alternatively, when each equilibrium is saved in a separate directory, gvec.find_state(๐ŸŒ) can load the state with only the directory provided.

state = gvec.find_state("./run_ellipstell/")

Derived Quantities#

One of the central features of the GVEC python bindings is the ability to compute a variety of derived quantities. A list of computable quantities is returned by gvec.table_of_quantities(markdown=True) and also shown at the end of this notebook.

The requested quantities will be stored in an xarray.Dataset, which groups several varaibles, their coordinates and metadata, similar to a pandas Dataframe. Such a Dataset can also be stored directly as netCDF. GVEC will also automatically determine which quantities are necessary to compute the desired quantities and add them to the Dataset, thereby also caching them for future computations.

The grid on which the quantities are to be evaluated, can be specified explicitly (with an array or float) or automatically with either an integer for linear spacing (within one field period) or "int" for the integration points required by GVEC.

The state.evaluate(๐ŸŒ) method can be used to create a new Dataset and compute desired quantities in a single step:

ev = state.evaluate(
    "pos", "B", rho=[0.1, 0.5, 0.9], theta=20, zeta=0.0
)
ev
<xarray.Dataset> Size: 27kB
Dimensions:     (rad: 3, pol: 20, tor: 1, xyz: 3)
Coordinates:
  * rho         (rad) float64 24B 0.1 0.5 0.9
  * theta       (pol) float64 160B 0.0 0.3142 0.6283 ... 5.341 5.655 5.969
  * zeta        (tor) float64 8B 0.0
  * xyz         (xyz) <U1 12B 'x' 'y' 'z'
Dimensions without coordinates: rad, pol, tor
Data variables: (12/44)
    X1          (rad, pol, tor) float64 480B 3.092 3.086 3.068 ... 3.946 4.128
    dX1_dr      (rad, pol, tor) float64 480B 1.217 1.153 0.9674 ... 1.652 1.838
    dX1_dt      (rad, pol, tor) float64 480B 0.0 -0.03914 ... 0.754 0.3979
    dX1_dz      (rad, pol, tor) float64 480B 0.0 0.02456 ... -0.5363 -0.2756
    dX1_drr     (rad, pol, tor) float64 480B -0.02021 -0.06073 ... 3.7 3.349
    dX1_drt     (rad, pol, tor) float64 480B 0.0 -0.4048 ... 0.777 0.4011
    ...          ...
    Jac         (rad, pol, tor) float64 480B 0.2556 0.2549 0.2527 ... 3.559 4.34
    e_theta     (xyz, rad, pol, tor) float64 1kB 0.0 -0.03914 ... 0.4543 0.5512
    e_zeta      (xyz, rad, pol, tor) float64 1kB 0.0 0.02456 ... 1.495 1.643
    B_contra_t  (rad, pol, tor) float64 480B 0.08812 0.08847 ... 0.09179 0.08512
    B_contra_z  (rad, pol, tor) float64 480B 0.1215 0.122 ... 0.06608 0.05889
    B           (rad, pol, tor, xyz) float64 1kB 0.0 -0.3757 ... -0.2431 0.1437

An existing dataset can be extended using state.compute(๐ŸŒ):

state.compute(ev, "J", "V")

The individual quantities, within an xarray.Dataset can be accessed using ev.B or ev["B"] and are enriched with coordinate information which can be used to a variety of operations, e.g.:

mean_B = ev.B.mean(dim=("rad", "pol", "tor"))
B2 = xr.dot(ev.B, ev.B, dim="xyz")
JxB = xr.cross(ev.J, ev.B, dim="xyz")

Indexing is best done using ev.B.sel(rho=0.5) (by value) or ev.B.isel(rad=0) (by position) and the raw numpy array can be extracted with ev.B.values. To ensure a particular order, use ev.B.transpose("rad", "pol", "tor", "xyz").values. To convert a scalar (e.g. the volume V) into a python scalar use the ev.V.item() method. The .squeeze() method can be used to remove any dimensions of length 1.

For more information, see the xarray documentation.

Computing integrals#

Specifying "int" for rho, theta and zeta in the Evaluations function, chooses the grid to be the integration points used by GVEC internally. The functions radial_integral, fluxsurface_integral and volume_integral can then be used to perform integration with the apropriate weights. E.g. to compute the volume averaged plasma beta, one could use:

ev = state.evaluate(
    "mod_B",
    "mu0",
    "p",
    "Jac",
    "V",
    rho="int",
    theta="int",
    zeta="int",
)
beta = ev.p / (ev.mod_B**2 / (2 * ev.mu0))
beta_avg = gvec.volume_integral(beta * ev.Jac) / ev.V
beta_avg.item()
0.008802419491949296

For supported quantities which require integration (e.g. the volume V computed above), an auxiliary dataset is created to compute them if the grid does not conform to the integration points.

Straight-fieldline coordinates#

GVEC uses flux aligned coordinates (internally also called logical coordinates) \(\rho,\vartheta,\zeta\) to describe the equilibrium. Available straight fieldline coordinates are PEST- \((\rho,\vartheta_P,\zeta)\) and Boozer coordinates \((\rho,\vartheta_B,\zeta_B)\) for which their associated quantities are denoted with a suffix or subscript P and B respectively.

PEST transform#

The PEST coordinates can be computed directly from the GVEC solution with \(\vartheta_P = \vartheta + \lambda(\rho,\vartheta,\zeta)\). To obtain a grid in PEST coordinates, you can use state.evaluate_sfl(..., sfl="pest")(๐ŸŒ), which performs a Newton-Raphson search for the respective logical coordinates and then evaluates the requested quantities.

ev = state.evaluate_sfl(
    "mod_B",
    "B_contra_t_P",
    "B_contra_z",  # B_contra_z_P = B_contra_z
    rho=1.0,
    theta=40,
    zeta=50,
    sfl="pest",
).squeeze()

fig, axs = plt.subplots(1, 2, figsize=(15, 6), layout="constrained")
fig.suptitle(
    (
        r"Streamlines of $\mathbf{B}$ in PEST"
        + r" and logical coordinates at $\rho=1.0$."
    )
)

streamplot_kwargs = dict(
    color="black",
    broken_streamlines=False,
    density=1,
    integration_direction="both",
    start_points=np.vstack(
        [
            np.linspace(0, 2 * np.pi / state.nfp, 22)[1:-1],
            np.linspace(2 * np.pi, 0, 22)[1:-1],
        ]
    ).T,
)

ax = axs[0]
c = ax.contour(
    ev.zeta, ev.theta_P, ev.mod_B, levels=21, cmap="plasma"
)
fig.colorbar(
    c, ax=[axs[0], axs[1]], label=f"${ev.mod_B.attrs['symbol']}$"
)
ax.streamplot(
    ev.zeta.values,
    ev.theta_P.values,
    ev.B_contra_z,
    ev.B_contra_t_P,
    **streamplot_kwargs,
)
ax.set_xlabel(r"$\zeta$")
ax.set_ylabel(r"$\vartheta_P$")
ax.set_title("Grid in PEST coordinates")

ev = state.evaluate(
    "mod_B", "B_contra_t", "B_contra_z", rho=1.0, theta=40, zeta=50
).squeeze()

ax = axs[1]
ax.contour(
    ev.zeta, ev.theta, ev.mod_B, levels=21, cmap=c.cmap, norm=c.norm
)
ax.streamplot(
    ev.zeta.values,
    ev.theta.values,
    ev.B_contra_z.values,
    ev.B_contra_t.values,
    **streamplot_kwargs,
)
ax.set_xlabel(r"$\zeta$")
ax.set_ylabel(r"$\vartheta$")
ax.set_title("Grid in logical coordinates")
Text(0.5, 1.0, 'Grid in logical coordinates')
../../_images/8a5a0eded962b993474b96843390cc6284fe992a2b2adb1b5a5b3d4d3a85d6ad.png

Note

In GVEC the PEST coordinates are generalized also to non-cylindrical coordinates. If a configuration uses a G-Frame or other non-RZฯ† \(h\)-map, the PEST coordinates retain the toroidal angle \(\zeta\) as defined in the \(h\)-map, which doesnโ€™t necessarily coincide with the cylindrical angle \(\varphi\).

Boozer transform#

To evaluate the equilibrium in Boozer angles, you can use state.evaluate_sfl(..., sfl="boozer")(๐ŸŒ), which performs a Boozer transform to obtain a set of \(\vartheta,\zeta\) points which correspond to your desired grid in \(\vartheta_B,\zeta_B\). The evaluations with this new dataset work the same as above, note however that the suffixes t and z still refer to components/derivatives with respect to \(\vartheta,\zeta\). Some additional quantities, like B_contra_t_B or e_zeta_B are now also available.

The argument MNfactor specifies how many fourier modes should be used for the Boozer potential \(\nu_B\) (and recomputed straight fieldline potential \(\lambda\)) relative to maximum fourier modes used for computing the equilibrium solution. This parameter has a strong influence on the accuracy and required computational effort!

ev = state.evaluate_sfl(
    "mod_B",
    "B_contra_t_B",
    "B_contra_z_B",
    rho=1.0,
    theta=40,
    zeta=50,
    sfl="boozer",
).squeeze()

fig, axs = plt.subplots(1, 2, figsize=(15, 6), layout="constrained")
fig.suptitle(
    (
        r"Streamlines of $\mathbf{B}$ in boozer"
        + r" and logical coordinates at $\rho=1.0$."
    )
)

streamplot_kwargs = dict(
    color="black",
    broken_streamlines=False,
    density=1,
    integration_direction="both",
    start_points=np.vstack(
        [
            np.linspace(0, 2 * np.pi / state.nfp, 22)[1:-1],
            np.linspace(2 * np.pi, 0, 22)[1:-1],
        ]
    ).T,
)

ax = axs[0]
c = ax.contour(
    ev.zeta_B, ev.theta_B, ev.mod_B, levels=21, cmap="plasma"
)
fig.colorbar(
    c, ax=[axs[0], axs[1]], label=f"${ev.mod_B.attrs['symbol']}$"
)
ax.streamplot(
    ev.zeta_B.values,
    ev.theta_B.values,
    ev.B_contra_z_B,
    ev.B_contra_t_B,
    **streamplot_kwargs,
)
ax.set_xlabel(r"$\zeta_B$")
ax.set_ylabel(r"$\theta_B$")
ax.set_title("Grid in boozer coordinates")

ev = state.evaluate(
    "mod_B", "B_contra_t", "B_contra_z", rho=1.0, theta=40, zeta=50
).squeeze()

ax = axs[1]
ax.contour(
    ev.zeta, ev.theta, ev.mod_B, levels=21, cmap=c.cmap, norm=c.norm
)
ax.streamplot(
    ev.zeta.values,
    ev.theta.values,
    ev.B_contra_z.values,
    ev.B_contra_t.values,
    **streamplot_kwargs,
)
ax.set_xlabel(r"$\zeta$")
ax.set_ylabel(r"$\theta$")
ax.set_title("Grid in logical coordinates")
Text(0.5, 1.0, 'Grid in logical coordinates')
../../_images/e1ac7b1412c8ae39d7cfc316ae0592aba20444db2b34e878f95d1d7654678e49.png

Note

The Boozer transform recomputes \(\lambda\) with a higher resolution (to satisfy the integrability condition for \(\nu_B\))! Therefore some quantities will differ between the equilibrium evaluation and Boozer evaluation.

In particular \(\langle B_\vartheta \rangle, \langle B_\zeta \rangle\) will differ from \(B_{\vartheta_B},B_{\zeta_B}\) by an offset.

Field-aligned grid#

The state.evaluate_sfl method and EvaluationsBoozer(๐ŸŒ) ( or EvaluationsPEST(๐ŸŒ) ) factory function can also be used to generate a non-tensorproduct grid by providing (up to 3D) arrays for the values of \(\theta_B,\zeta_B\). This can be used for example to create a field-aligned grid:

rho = [0.5, 1.0]  # radial positions
alpha = np.linspace(
    0, 2 * np.pi, 100, endpoint=False
)  # fieldline label
phi = np.linspace(
    0, 2 * np.pi / state.nfp, 101
)  # angle along the fieldline

# evaluate the rotational transform (fieldline angle) on the desired surfaces
iota = state.evaluate("iota", rho=rho, theta=None, zeta=None).iota

# 3D toroidal and poloidal arrays that correspond to fieldline coordinates for each surface
theta_B = (
    alpha[None, :, None]
    + iota.data[:, None, None] * phi[None, None, :]
)

# create the grid
ev = gvec.EvaluationsBoozer(
    rho=rho, theta_B=theta_B, zeta_B=phi, state=state, MNfactor=5
)

# set the fiedline label as poloidal coordinate & index (not necessary, but good practice)
ev["alpha"] = ("pol", alpha)
ev["alpha"].attrs = dict(
    symbol=r"\alpha", long_name="fieldline label"
)
ev = ev.set_coords("alpha").set_xindex("alpha")

state.compute(ev, "B", "B_contra_t_B", "B_contra_z_B", "mod_B")
ev = ev.sel(rho=0.5)

fig, ax = plt.subplots(figsize=(8, 6), layout="constrained")

c = ax.contour(
    ev.zeta_B, ev.alpha, ev.mod_B, levels=21, cmap="plasma"
)
fig.colorbar(c, ax=ax, label=f"${ev.mod_B.attrs['symbol']}$")
ax.set_xlabel(f"${ev.zeta_B.attrs['symbol']}$")
ax.set_ylabel(f"${ev.alpha.attrs['symbol']}$")
ax.set_title(
    (
        f"${ev.mod_B.attrs['symbol']}$"
        + r" on a field-aligned grid at $\rho=0.5$"
    )
)
Text(0.5, 1.0, '$\\left|\\mathbf{B}\\right|$ on a field-aligned grid at $\\rho=0.5$')
../../_images/7250d18f31a1a42d692b0bc15ad60e839245268b8ca9c0ed47ebf907c50c5183.png

Available quantities for evaluation#

The following table contains the quantities that can be evaluated with the python bindings, it can be generated with

gvec.table_of_quantities(markdown=True)

label

long name

symbol

A_surface

surface area

\(A_\text{surface}\)

B

magnetic field

\(\mathbf{B}\)

B_contra_t

poloidal component of the magnetic field

\(B^\theta\)

B_contra_t_B

contravariant \(\theta\) component of the magnetic field in Boozer coordinates

\(B^{\theta_B}\)

B_contra_t_P

contravariant \(\theta\) component of the magnetic field in PEST coordinates

\(B^{\theta_P}\)

B_contra_z

toroidal component of the magnetic field

\(B^\zeta\)

B_contra_z_B

contravariant \(\zeta\) component of the magnetic field in Boozer coordinates

\(B^{\zeta_B}\)

B_rho_B

\(\rho\) component of the magnetic field in Boozer coordinates

\(B_{\rho_B}\)

B_rho_P

\(\rho\) component of the magnetic field in PEST coordinates

\(B_{\rho_P}\)

B_theta_B

\(\theta\) component of the magnetic field in Boozer coordinates

\(B_{\theta_B}\)

B_theta_P

\(\theta\) component of the magnetic field in PEST coordinates

\(B_{\theta_P}\)

B_theta_avg

flux-surface averaged poloidal magnetic field

\(\overline{B_\theta}\)

B_zeta_B

\(\zeta\) component of the magnetic field in Boozer coordinates

\(B_{\zeta_B}\)

B_zeta_P

\(\zeta\) component of the magnetic field in PEST coordinates

\(B_{\zeta_P}\)

B_zeta_avg

flux-surface averaged toroidal magnetic field

\(\overline{B_\zeta}\)

D_Merc

Mercier criterion

\(D_\text{Merc}\)

D_Merc_Curr

Current contribution to the Mercier criterion

\(D_\text{M,Curr}\)

D_Merc_Geod

Geodesic contribution to the Mercier criterion

\(D_\text{M,Geod}\)

D_Merc_Shear

Shear contribution to the Mercier criterion

\(D_\text{M,Shear}\)

D_Merc_Well

Magnetic well contribution to the Mercier criterion

\(D_\text{M,Well}\)

F

MHD force

\(F\)

F_r_avg

radial force balance

\(\overline{F_\rho}\)

II_tt

poloidal component of the second fundamental form

\(\mathrm{II}_{\theta\theta}\)

II_tt_B

t,t component of the second fundamental form in Boozer coordinates

\(\mathrm{II}_{\theta_B \theta_B}\)

II_tt_P

t,t component of the second fundamental form in PEST coordinates

\(\mathrm{II}_{\theta_P \theta_P}\)

II_tz

poloidal-toroidal component of the second fundamental form

\(\mathrm{II}_{\theta\zeta}\)

II_tz_B

t,z component of the second fundamental form in Boozer coordinates

\(\mathrm{II}_{\theta_B \zeta_B}\)

II_tz_P

t,z component of the second fundamental form in PEST coordinates

\(\mathrm{II}_{\theta_P \zeta_P}\)

II_zz

toroidal component of the second fundamental form

\(\mathrm{II}_{\zeta\zeta}\)

II_zz_B

z,z component of the second fundamental form in Boozer coordinates

\(\mathrm{II}_{\zeta_B \zeta_B}\)

II_zz_P

z,z component of the second fundamental form in PEST coordinates

\(\mathrm{II}_{\zeta_P \zeta_P}\)

I_pol

poloidal current, relative to the magnetic axis

\(I_\text{pol}\)

I_tor

toroidal current enclosed by flux surface

\(I_\text{tor}\)

J

current density

\(\mathbf{J}\)

J_contra_r

contravariant radial current density

\(J^{\rho}\)

J_contra_t

contravariant poloidal current density

\(J^{\theta}\)

J_contra_t_B

contravariant \(\theta\) component of the current density in Boozer coordinates

\(J^{\theta_B}\)

J_contra_t_P

contravariant \(\theta\) component of the current density in PEST coordinates

\(J^{\theta_P}\)

J_contra_z

contravariant toroidal current density

\(J^{\zeta}\)

J_contra_z_B

contravariant \(\zeta\) component of the current density in Boozer coordinates

\(J^{\zeta_B}\)

J_rho_B

\(\rho\) component of the current density in Boozer coordinates

\(J_{\rho_B}\)

J_rho_P

\(\rho\) component of the current density in PEST coordinates

\(J_{\rho_P}\)

J_theta_B

\(\theta\) component of the current density in Boozer coordinates

\(J_{\theta_B}\)

J_theta_P

\(\theta\) component of the current density in PEST coordinates

\(J_{\theta_P}\)

J_zeta_B

\(\zeta\) component of the current density in Boozer coordinates

\(J_{\zeta_B}\)

J_zeta_P

\(\zeta\) component of the current density in PEST coordinates

\(J_{\zeta_P}\)

Jac

Jacobian determinant

\(\mathcal{J}\)

Jac_B

Jacobian determinant in Boozer coordinates

\(\mathcal{J}_B\)

Jac_P

Jacobian determinant in PEST coordinates

\(\mathcal{J}_P\)

Jac_h

reference Jacobian determinant

\(\mathcal{J}_h\)

Jac_l

logical Jacobian determinant

\(\mathcal{J}_l\)

LA

straight field line potential

\(\lambda\)

L_axis

length of the magnetic axis

\(L_\text{axis}\)

L_gradB

magnetic gradient scale length

\(L_{\nabla\mathbf{B}}\)

N_FP

number of field periods

\(N_\text{FP}\)

Phi

toroidal magnetic flux

\(\Phi\)

Phi_edge

toroidal magnetic flux at the edge

\(\Phi_0\)

V

volume

\(V\)

W_MHD

total MHD energy

\(W_\text{MHD}\)

X1

first reference coordinate

\(X^1\)

X2

second reference coordinate

\(X^2\)

aspect_ratio

effective aspect ratio

\(a_\text{eff}\)

beta_avg

volume-averaged plasma beta

\(\overline{\beta}\)

chi

poloidal magnetic flux

\(\chi\)

dA

differential area element

\(dA\)

dB_contra_t_dr

radial derivative of the poloidal magnetic field

\(\frac{\partial B^\theta}{\partial \rho}\)

dB_contra_t_dt

poloidal derivative of the poloidal magnetic field

\(\frac{\partial B^\theta}{\partial \theta}\)

dB_contra_t_dz

toroidal derivative of the poloidal magnetic field

\(\frac{\partial B^\theta}{\partial \zeta}\)

dB_contra_z_dr

radial derivative of the toroidal magnetic field

\(\frac{\partial B^\zeta}{\partial \rho}\)

dB_contra_z_dt

poloidal derivative of the toroidal magnetic field

\(\frac{\partial B^\zeta}{\partial \theta}\)

dB_contra_z_dz

toroidal derivative of the toroidal magnetic field

\(\frac{\partial B^\zeta}{\partial \zeta}\)

dB_dr

radial derivative of the magnetic field

\(\frac{\partial \mathbf{B}}{\partial \rho}\)

dB_dt

poloidal derivative of the magnetic field

\(\frac{\partial \mathbf{B}}{\partial \theta}\)

dB_dz

toroidal derivative of the magnetic field

\(\frac{\partial \mathbf{B}}{\partial \zeta}\)

dB_theta_avg_dr

derivative of the flux-surface averaged poloidal magnetic field

\(\frac{d\overline{B_\theta}}{d\rho}\)

dI_tor_dr

derivative of the toroidal current enclosed by the flux surface

\(\frac{dI_\text{tor}}{d\rho}\)

dJac_dr

radial derivative of the Jacobian determinant

\(\frac{\partial \mathcal{J}}{\partial \rho}\)

dJac_dt

poloidal derivative of the Jacobian determinant

\(\frac{\partial \mathcal{J}}{\partial \theta}\)

dJac_dz

toroidal derivative of the Jacobian determinant

\(\frac{\partial \mathcal{J}}{\partial \zeta}\)

dJac_h_dr

radial derivative of the reference Jacobian determinant

\(\frac{\partial \mathcal{J}_h}{\partial \rho}\)

dJac_h_dt

poloidal derivative of the reference Jacobian determinant

\(\frac{\partial \mathcal{J}_h}{\partial \theta}\)

dJac_h_dz

toroidal derivative of the reference Jacobian determinant

\(\frac{\partial \mathcal{J}_h}{\partial \zeta}\)

dJac_l_dr

radial derivative of the logical Jacobian determinant

\(\frac{\partial \mathcal{J}_l}{\partial \rho}\)

dJac_l_dt

poloidal derivative of the logical Jacobian determinant

\(\frac{\partial \mathcal{J}_l}{\partial \theta}\)

dJac_l_dz

toroidal derivative of the logical Jacobian determinant

\(\frac{\partial \mathcal{J}_l}{\partial \zeta}\)

dLA_dr

radial derivative of the straight field line potential

\(\frac{\partial \lambda}{\partial \rho}\)

dLA_drr

second radial derivative of the straight field line potential

\(\frac{\partial^2 \lambda}{\partial \rho^2}\)

dLA_drt

radial-poloidal derivative of the straight field line potential

\(\frac{\partial^2 \lambda}{\partial \rho\partial \theta}\)

dLA_drz

radial-toroidal derivative of the straight field line potential

\(\frac{\partial^2 \lambda}{\partial \rho\partial \zeta}\)

dLA_dt

poloidal derivative of the straight field line potential

\(\frac{\partial \lambda}{\partial \theta}\)

dLA_dtt

second poloidal derivative of the straight field line potential

\(\frac{\partial^2 \lambda}{\partial \theta^2}\)

dLA_dtz

poloidal-toroidal derivative of the straight field line potential

\(\frac{\partial^2 \lambda}{\partial \theta\partial \zeta}\)

dLA_dz

toroidal derivative of the straight field line potential

\(\frac{\partial \lambda}{\partial \zeta}\)

dLA_dzz

second toroidal derivative of the straight field line potential

\(\frac{\partial^2 \lambda}{\partial \zeta^2}\)

dNU_B_dt

poloidal derivative of the Boozer potential computed from the magnetic field

\(\left.\frac{\partial \nu_B}{\partial \theta}\right|_\text{def.}\)

dNU_B_dz

toroidal derivative of the Boozer potential computed from the magnetic field

\(\left.\frac{\partial \nu_B}{\partial \zeta}\right|_\text{def.}\)

dPhi_dr

toroidal magnetic flux gradient

\(\frac{d\Phi}{d\rho}\)

dPhi_drr

toroidal magnetic flux curvature

\(\frac{d^2\Phi}{d\rho^2}\)

dV_dPhi_n

derivative of the plasma volume w.r.t. normalized toroidal magnetic flux

\(\frac{dV}{d\Phi_n}\)

dV_dPhi_n2

second derivative of the plasma volume w.r.t. normalized toroidal magnetic flux

\(\frac{d^2V}{d\Phi_n^2}\)

dX1_dr

radial derivative of the first reference coordinate

\(\frac{\partial X^1}{\partial \rho}\)

dX1_drr

second radial derivative of the first reference coordinate

\(\frac{\partial^2 X^1}{\partial \rho^2}\)

dX1_drt

radial-poloidal derivative of the first reference coordinate

\(\frac{\partial^2 X^1}{\partial \rho\partial \theta}\)

dX1_drz

radial-toroidal derivative of the first reference coordinate

\(\frac{\partial^2 X^1}{\partial \rho\partial \zeta}\)

dX1_dt

poloidal derivative of the first reference coordinate

\(\frac{\partial X^1}{\partial \theta}\)

dX1_dtt

second poloidal derivative of the first reference coordinate

\(\frac{\partial^2 X^1}{\partial \theta^2}\)

dX1_dtz

poloidal-toroidal derivative of the first reference coordinate

\(\frac{\partial^2 X^1}{\partial \theta\partial \zeta}\)

dX1_dz

toroidal derivative of the first reference coordinate

\(\frac{\partial X^1}{\partial \zeta}\)

dX1_dzz

second toroidal derivative of the first reference coordinate

\(\frac{\partial^2 X^1}{\partial \zeta^2}\)

dX2_dr

radial derivative of the second reference coordinate

\(\frac{\partial X^2}{\partial \rho}\)

dX2_drr

second radial derivative of the second reference coordinate

\(\frac{\partial^2 X^2}{\partial \rho^2}\)

dX2_drt

radial-poloidal derivative of the second reference coordinate

\(\frac{\partial^2 X^2}{\partial \rho\partial \theta}\)

dX2_drz

radial-toroidal derivative of the second reference coordinate

\(\frac{\partial^2 X^2}{\partial \rho\partial \zeta}\)

dX2_dt

poloidal derivative of the second reference coordinate

\(\frac{\partial X^2}{\partial \theta}\)

dX2_dtt

second poloidal derivative of the second reference coordinate

\(\frac{\partial^2 X^2}{\partial \theta^2}\)

dX2_dtz

poloidal-toroidal derivative of the second reference coordinate

\(\frac{\partial^2 X^2}{\partial \theta\partial \zeta}\)

dX2_dz

toroidal derivative of the second reference coordinate

\(\frac{\partial X^2}{\partial \zeta}\)

dX2_dzz

second toroidal derivative of the second reference coordinate

\(\frac{\partial^2 X^2}{\partial \zeta^2}\)

db_dr

radial derivative of the normalized magnetic field

\(\frac{\partial \mathbf{b}}{\partial \rho}\)

db_dt

poloidal derivative of the normalized magnetic field

\(\frac{\partial \mathbf{b}}{\partial \theta}\)

db_dz

toroidal derivative of the normalized magnetic field

\(\frac{\partial \mathbf{b}}{\partial \zeta}\)

dchi_dr

poloidal magnetic flux gradient

\(\frac{d\chi}{d\rho}\)

dchi_drr

poloidal magnetic flux curvature

\(\frac{d^2\chi}{d\rho^2}\)

dg_rr_dr

radial derivative of the rr component of the metric tensor

\(\frac{\partial g_{\rho\rho}}{\partial \rho}\)

dg_rr_dt

poloidal derivative of the rr component of the metric tensor

\(\frac{\partial g_{\rho\rho}}{\partial \theta}\)

dg_rr_dz

toroidal derivative of the rr component of the metric tensor

\(\frac{\partial g_{\rho\rho}}{\partial \zeta}\)

dg_rt_dr

radial derivative of the rt component of the metric tensor

\(\frac{\partial g_{\rho\theta}}{\partial \rho}\)

dg_rt_dt

poloidal derivative of the rt component of the metric tensor

\(\frac{\partial g_{\rho\theta}}{\partial \theta}\)

dg_rt_dz

toroidal derivative of the rt component of the metric tensor

\(\frac{\partial g_{\rho\theta}}{\partial \zeta}\)

dg_rz_dr

radial derivative of the rz component of the metric tensor

\(\frac{\partial g_{\rho\zeta}}{\partial \rho}\)

dg_rz_dt

poloidal derivative of the rz component of the metric tensor

\(\frac{\partial g_{\rho\zeta}}{\partial \theta}\)

dg_rz_dz

toroidal derivative of the rz component of the metric tensor

\(\frac{\partial g_{\rho\zeta}}{\partial \zeta}\)

dg_tt_dr

radial derivative of the tt component of the metric tensor

\(\frac{\partial g_{\theta\theta}}{\partial \rho}\)

dg_tt_dt

poloidal derivative of the tt component of the metric tensor

\(\frac{\partial g_{\theta\theta}}{\partial \theta}\)

dg_tt_dz

toroidal derivative of the tt component of the metric tensor

\(\frac{\partial g_{\theta\theta}}{\partial \zeta}\)

dg_tz_dr

radial derivative of the tz component of the metric tensor

\(\frac{\partial g_{\theta\zeta}}{\partial \rho}\)

dg_tz_dt

poloidal derivative of the tz component of the metric tensor

\(\frac{\partial g_{\theta\zeta}}{\partial \theta}\)

dg_tz_dz

toroidal derivative of the tz component of the metric tensor

\(\frac{\partial g_{\theta\zeta}}{\partial \zeta}\)

dg_zz_dr

radial derivative of the zz component of the metric tensor

\(\frac{\partial g_{\zeta\zeta}}{\partial \rho}\)

dg_zz_dt

poloidal derivative of the zz component of the metric tensor

\(\frac{\partial g_{\zeta\zeta}}{\partial \theta}\)

dg_zz_dz

toroidal derivative of the zz component of the metric tensor

\(\frac{\partial g_{\zeta\zeta}}{\partial \zeta}\)

diota_dr

rotational transform gradient

\(\frac{d\iota}{d\rho}\)

diota_drr

rotational transform curvature

\(\frac{d^2\iota}{d\rho^2}\)

dmod_B_dr

radial derivative of the modulus of the magnetic field

\(\frac{\partial \left|\mathbf{B}\right|}{\partial \rho}\)

dmod_B_dr_B

radial Boozer derivative of the modulus of the magnetic field

\(\frac{\partial\left|\mathbf{B}\right|}{\partial \rho_B}\)

dmod_B_dr_P

radial PEST-like derivative of the modulus of the magnetic field

\(\frac{\partial\left|\mathbf{B}\right|}{\partial \rho_P}\)

dmod_B_dt

poloidal derivative of the modulus of the magnetic field

\(\frac{\partial \left|\mathbf{B}\right|}{\partial \theta}\)

dmod_B_dt_B

poloidal Boozer derivative of the modulus of the magnetic field

\(\frac{\partial\left|\mathbf{B}\right|}{\partial \theta_B}\)

dmod_B_dt_P

poloidal PEST-like derivative of the modulus of the magnetic field

\(\frac{\partial\left|\mathbf{B}\right|}{\partial \vartheta_P}\)

dmod_B_dz

toroidal derivative of the modulus of the magnetic field

\(\frac{\partial \left|\mathbf{B}\right|}{\partial \zeta}\)

dmod_B_dz_B

toroidal Boozer derivative of the modulus of the magnetic field

\(\frac{\partial\left|\mathbf{B}\right|}{\partial \zeta_B}\)

dmod_B_dz_P

toroidal PEST-like derivative of the modulus of the magnetic field

\(\frac{\partial\left|\mathbf{B}\right|}{\partial \zeta_P}\)

dp_dr

pressure gradient

\(\frac{dp}{d\rho}\)

dp_drr

pressure curvature

\(\frac{d^2p}{d\rho^2}\)

e_q1

first reference tangent basis vector

\(\mathbf{e}_{q^1}\)

e_q2

second reference tangent basis vector

\(\mathbf{e}_{q^2}\)

e_q3

toroidal reference tangent basis vector

\(\mathbf{e}_{q^3}\)

e_rho

radial tangent basis vector

\(\mathbf{e}_\rho\)

e_rho_B

radial tangent basis vector in Boozer coordinates

\(\mathbf{e}_{\rho_B}\)

e_rho_P

poloidal tangent basis vector in PEST coordinates

\(\mathbf{e}_{\theta_P}\)

e_theta

poloidal tangent basis vector

\(\mathbf{e}_\theta\)

e_theta_B

poloidal tangent basis vector in Boozer coordinates

\(\mathbf{e}_{\theta_B}\)

e_theta_P

poloidal tangent basis vector in PEST coordinates

\(\mathbf{e}_{\theta_P}\)

e_zeta

toroidal tangent basis vector

\(\mathbf{e}_\zeta\)

e_zeta_B

toroidal tangent basis vector in Boozer coordinates

\(\mathbf{e}_{\zeta_B}\)

e_zeta_P

toroidal tangent basis vector in PEST coordinates

\(\mathbf{e}_{\zeta_P}\)

elongation

effective elongation

\(E_\text{eff}\)

g_rr

rr component of the metric tensor

\(g_{\rho\rho}\)

g_rr_B

rr component of the metric tensor in Boozer coordinates

\(g_{\rho_B \rho_B}\)

g_rr_P

rr component of the metric tensor in PEST coordinates

\(g_{\rho_P \rho_P}\)

g_rt

rt component of the metric tensor

\(g_{\rho\theta}\)

g_rt_B

rt component of the metric tensor in Boozer coordinates

\(g_{\rho_B \theta_B}\)

g_rt_P

rt component of the metric tensor in PEST coordinates

\(g_{\rho_P \theta_P}\)

g_rz

rz component of the metric tensor

\(g_{\rho\zeta}\)

g_rz_B

rz component of the metric tensor in Boozer coordinates

\(g_{\rho_B \zeta_B}\)

g_rz_P

rz component of the metric tensor in PEST coordinates

\(g_{\rho_P \zeta_P}\)

g_tt

tt component of the metric tensor

\(g_{\theta\theta}\)

g_tt_B

tt component of the metric tensor in Boozer coordinates

\(g_{\theta_B \theta_B}\)

g_tt_P

tt component of the metric tensor in PEST coordinates

\(g_{\theta_P \theta_P}\)

g_tz

tz component of the metric tensor

\(g_{\theta\zeta}\)

g_tz_B

tz component of the metric tensor in Boozer coordinates

\(g_{\theta_B \zeta_B}\)

g_tz_P

tz component of the metric tensor in PEST coordinates

\(g_{\theta_P \zeta_P}\)

g_zz

zz component of the metric tensor

\(g_{\zeta\zeta}\)

g_zz_B

zz component of the metric tensor in Boozer coordinates

\(g_{\zeta_B \zeta_B}\)

g_zz_P

zz component of the metric tensor in PEST coordinates

\(g_{\zeta_P \zeta_P}\)

gamma

adiabatic index

\(\gamma\)

grad_mod_B

gradient of the modulus of the magnetic field

\(\nabla\left|\mathbf{B}\right|\)

grad_rho

radial reciprocal basis vector

\(\nabla\rho\)

grad_theta

poloidal reciprocal basis vector

\(\nabla\theta\)

grad_theta_B

poloidal reciprocal basis vector in Boozer coordinates

\(\nabla\theta_B\)

grad_theta_P

poloidal reciprocal basis vector in PEST coordinates

\(\nabla \theta_P\)

grad_zeta

toroidal reciprocal basis vector

\(\nabla\zeta\)

grad_zeta_B

toroidal reciprocal basis vector in Boozer coordinates

\(\nabla\zeta_B\)

iota

rotational transform

\(\iota\)

iota_0

geometric contribution to the rotational transform

\(\iota_0\)

iota_avg

average rotational transform

\(\overline{\iota}\)

iota_avg2

rotational transform averaged over rho^2

\(\overline{\iota}_2\)

iota_curr

toroidal current contribution to the rotational transform

\(\iota_\text{curr}\)

iota_curr_0

factor to the toroidal current contribution to the rotational transform

\(\iota_{\text{curr},0}\)

k_q1q1

q1-q1 reference curvature vector

\(k_{q^1q^1}\)

k_q1q2

q1-q2 reference curvature vector

\(k_{q^1q^2}\)

k_q1q3

q1-q3 reference curvature vector

\(k_{q^1q^3}\)

k_q2q2

q2-q2 reference curvature vector

\(k_{q^2q^2}\)

k_q2q3

q2-q3 reference curvature vector

\(k_{q^2q^3}\)

k_q3q3

q3-q3 reference curvature vector

\(k_{q^3q^3}\)

k_rr

rr logical curvature vector

\(\mathbf{k}_{\rho\rho}\)

k_rt

rt logical curvature vector

\(\mathbf{k}_{\rho\theta}\)

k_rz

rz logical curvature vector

\(\mathbf{k}_{\rho\zeta}\)

k_tt

tt logical curvature vector

\(\mathbf{k}_{\theta\theta}\)

k_tt_B

tt boozer curvature vector

\(\mathbf{k}_{\theta_B \theta_B}\)

k_tt_P

tt PEST curvature vector

\(\mathbf{k}_{\theta_P \theta_P}\)

k_tz

tz logical curvature vector

\(\mathbf{k}_{\theta\zeta}\)

k_tz_B

tz boozer curvature vector

\(\mathbf{k}_{\theta_B \zeta_B}\)

k_tz_P

tz PEST curvature vector

\(\mathbf{k}_{\theta_P \zeta_P}\)

k_zz

zz logical curvature vector

\(\mathbf{k}_{\zeta\zeta}\)

k_zz_B

zz boozer curvature vector

\(\mathbf{k}_{\zeta_B \zeta_B}\)

k_zz_P

zz PEST curvature vector

\(\mathbf{k}_{\zeta_P \zeta_P}\)

kappa_B

field line curvature

\(\mathbf{\kappa}_B\)

kappa_G

geodesic curvature

\(\kappa_G\)

mirror_ratio

mirror ratio

\(\Delta_\text{mirror}\)

mod_B

modulus of the magnetic field

\(\left|\mathbf{B}\right|\)

mod_F

modulus of the MHD force

\(\left|F\right|\)

mod_J

modulus of the current density

\(\left|\mathbf{J}\right|\)

mod_e_rho

modulus of the radial tangent basis vector

\(\left|\mathbf{e}_\rho\right|\)

mod_e_theta

modulus of the poloidal tangent basis vector

\(\left|\mathbf{e}_\theta\right|\)

mod_e_zeta

modulus of the toroidal tangent basis vector

\(\left|\mathbf{e}_\zeta\right|\)

mod_grad_rho

modulus of the radial reciprocal basis vector

\(\left|\nabla\rho\right|\)

mod_grad_theta

modulus of the poloidal reciprocal basis vector

\(\left|\nabla\theta\right|\)

mod_grad_zeta

modulus of the toroidal reciprocal basis vector

\(\left|\nabla\zeta\right|\)

mu0

magnetic constant

\(\mu_0\)

normal

surface normal

\(\mathbf{n}\)

p

pressure

\(p\)

pos

position vector

\(\mathbf{x}\)

r_major

effective major radius

\(r_\text{major,eff}\)

r_minor

effective minor radius

\(r_\text{minor,eff}\)

shear

global magnetic shear

\(s_g\)

shear_avg

average global magnetic shear

\(\overline{s_g}\)

shear_avg2

global magnetic shear averaged over rho^2

\(\overline{s_g}_2\)

theta_P

poloidal angle in PEST coordinates

\(\theta_P\)

vacuum_magnetic_well_depth

vacuum magnetic well depth

\(d_\text{well}\)

xyz

cartesian vector components

\((x,y,z)\)