pyrokinetics.equilibrium.equilibrium.Equilibrium#
- class pyrokinetics.equilibrium.equilibrium.Equilibrium(R, Z, psi_RZ, psi, F, FF_prime, p, p_prime, q, R_major, r_minor, Z_mid, psi_lcfs, a_minor, B_0=None, I_p=None, clockwise_phi=False, cocos=None, eq_type=None, surface_interps=None)[source]#
Bases:
DatasetWrapper,ReadableFromFileContains a solution of the Grad-Shafranov equation, which defines a tokamak plasma equilibrium. Users are not expected to initialise
Equilibriumobjects directly, and in most cases should instead make use of the functionread_equilibrium, which can read popular plasma equilibrium file types such as GEQDSK, generated by tools such as FreeGS and EFIT, or netCDF4 files generated by TRANSP or VMEC.The user may specify the COCOS convention followed by their inputs. Alternatively, If supplied with
b_axis,current, and (optionally)clockwise_phi, inputs will be converted to the convention COCOS 11. Otherwise, it is assumed that the inputs are already COCOS 11.- Parameters:
R (
ArrayLike,units [meter]) – Linearly spaced and monotonically increasing 1D grid of major radius coordinates. This is the radius from the central column of a tokamak, and not the radial distance from the magnetic axis.Z (
ArrayLike,units [meter]) – Linearly spaced and monotonically increasing 1D grid of tokamak z-coordinates. This is usually the height above the plasma midplane, but z=0 may be set at any reference point.psi_RZ (
ArrayLike,units [weber]) – 2D grid defining the poloidal magnetic flux function :math::psiwith respect toRandZ. Should have the shape(len(r), len(z)). If supplied in units of Webers per radian (i.e. following COCOS 1 to 8), this will be converted.psi (
ArrayLike,units [weber]) – 1D grid defining the poloidal magnetic flux function. This grid defines magnetic flux surface coordinates, on which most other parameters are defined.psi[0]should be the value ofpsion the magnetic axis.psishould be monotonically increasing or decreasing. If supplied in units of Weber per radian (i.e. following COCOS 11 to 18), this will be converted.psi_tor (
ArrayLike,units [weber]) – 1D grid defining the toroidal magnetic flux functionpsi_tor.F (
ArrayLike,units [meter * tesla]) – 1D grid defining the poloidal current function with respect topsi. Should have the same length aspsi.FF_prime (
ArrayLike,units [meter**2 * tesla**2 / weber]) – 1D grid defining the poloidal current functionfmultiplied by its derivative with respect topsi. Should have the same length aspsi.p (
ArrayLike,units [pascal]) – 1D grid defining the plasma pressure with respect topsi. Should have the same length aspsi.p_prime (
ArrayLike,units [pascal / weber]) – 1D grid defining the derivative of the plasma pressure with respect topsi. Should have the same length aspsi.q (
ArrayLike,units [dimensionless]) – 1D grid defining the ‘safety factor’ with respect topsi. Should have same length aspsi.R_major (
ArrayLike,units [meter]) – 1D grid of the major radius positions of the center of each flux surface with respect topsi. Should have the same length aspsi. This should be given by the mean of the maximum and minimum major radii of the flux surface.R_major[0]should be the major radius of the magnetic axis.r_minor (
ArrayLike,units[meter]) – 1D grid of the minor radius of each flux surface with respect topsi. Should have the same length aspsi. This should be half of the difference between the maximum and minimum major radii of a flux surface.r_minor[0]should be 0.0.Z_mid (
ArrayLike,units [meter]) – 1D grid of the z-midpoint of each flux surface with respace topsi. Should have the same length aspsi. This should be the mean of the maximum and minimum z-positions of the flux surface.Z_mid[0]should be the z-position of the magnetic axis.psi_lcfs (
float,units [weber]) – The value ofpsion the last closed flux surface (LCFS).a_minor (
float,units [meter]) – The minor radius of the last closed flux surface (LCFS). The minor radius of a flux surface is half of the difference between its maximum and minimum major radii.B_axis (
Optional[float],units [tesla]) – Vacuum toroidal magnetic field. Used to determine COCOS convention of the inputs.current (
Optional[float],units [ampere]) – Plasma current. Used to determine COCOS convention of the inputs.clockwise_phi (
bool, defaultFalse) – Determines whether the \(\phi\) grid increases clockwise or anti-clockwise when viewed from above. Used to determine COCOS convention of the inputs.cocos (
Optional[int]) – Asserts that the inputs follow a particular COCOS convention. If set, prevents the automatic detection of COCOS conventions using the values ofb_axis,current,q,psi,r_minor, andclockwise_phi. Results will be converted from this convention to COCOS 11. Possible values are between 1 to 8, or between 11 to 18.eq_type (
Optional[str]) – A label denoting the type of Equilibrium, such as “GEQDSK”, “TRANSP”, “VMEC”, etc.surface_interps (
Optional[dict]) – A dictionary of interpolations of R, Z and Bpol as a function of psin and theta Useful for codes that directly spit out these surfacesB_0 (Optional[float])
I_p (Optional[float])
- data#
The internal representation of the
Equilibriumobject. The functions__getattr__and__getitem__redirect most attribute/indexing lookups here, but the Dataset itself may be accessed directly by the user if they wish to perform more complex manipulations.- Type:
- psi_lcfs#
Poloidal magnetic flux function \(\psi\) on the last closed flux surface.
- Type:
float,units [weber]
See also
read_equilibriumCreate
Equilibriumfrom a file.FluxSurfaceAn individual flux surface created from an
Equilibrium.
Notes
Here we make use of COCOS 11 [Sauter et al. 2013], which is the ITER standard for defining a tokamak coordinate system:
The cylindrical coordinate system representing the whole Tokamak is right-handed and follows the ordering \((R,\theta,Z)\).
When viewed from above the Tokamak, the toroidal angle \(\phi\) increases in the counter-clockwise direction.
The coordinate system in the poloidal plane is right-handed and follows the ordering \((r,\theta,\phi)\).
When viewing the poloidal plane on the right hand side of the major vertical axis, the poloidal angle
thetaincreases in the clockwise direction.
The Grad-Shafranov equation [1] may be written:
\[\frac{\partial^2 \psi}{\partial R^2} - \frac{1}{R}\frac{\partial\psi}{\partial R} + \frac{\partial^2 \psi}{\partial Z^2} = -\mu_0 r^2 p'(\psi) - F(\psi)F'(\psi)\]\(\psi\) is the poloidal magnetic flux function, defined as:
\[\vec{B_p} = \frac{1}{2\pi r} \nabla\psi \times \nabla\phi\]where \(\vec{B_p}\) is the poloidal magnetic field vector. \(\psi\) is arbitrary to within an additive constant. The choice to include the factor of \(2\pi\) is part of the COCOS 11 convention. \(F\) is the poloidal current function, and is defined similarly:
\[\mu_0 \vec{J_p} = \nabla F \times \nabla\phi\]where \(\vec{J_p}\) is the poloidal current density. Comparing this to Ampere’s law gives the relation:
\[F = R B_\phi\]It may be shown that \(F\) is a function of \(\psi\) only, as is the plasma pressure \(p\), and the ‘safety factor’ \(q(\psi)\):
\[q = \frac{1}{2\pi} \int\frac{\vec{B}\cdot\nabla\phi}{\vec{B}\cdot\nabla\theta}d\theta\]\(q\) describes the number of times a magnetic field line wraps around the toroidal direction before arriving at the same point in the poloidal plane.
\(\psi\) is often used to index flux surfaces, so one can effectively work in the poloidal coordinate system \((\psi,\theta,\phi)\). It is common for researchers to instead make use of the normalised \(\psi_n\), which is scaled so that \(\psi_n=0\) at the magnetic axis, and \(\psi_n=1\) at the Last Closed Flux Surface (LCFS).
- __init__(R, Z, psi_RZ, psi, F, FF_prime, p, p_prime, q, R_major, r_minor, Z_mid, psi_lcfs, a_minor, B_0=None, I_p=None, clockwise_phi=False, cocos=None, eq_type=None, surface_interps=None)[source]#
- Parameters:
R (numpy.typing.ArrayLike)
Z (numpy.typing.ArrayLike)
psi_RZ (numpy.typing.ArrayLike)
psi (numpy.typing.ArrayLike)
F (numpy.typing.ArrayLike)
FF_prime (numpy.typing.ArrayLike)
p (numpy.typing.ArrayLike)
p_prime (numpy.typing.ArrayLike)
q (numpy.typing.ArrayLike)
R_major (numpy.typing.ArrayLike)
r_minor (numpy.typing.ArrayLike)
Z_mid (numpy.typing.ArrayLike)
psi_lcfs (float)
a_minor (float)
B_0 (float | None)
I_p (float | None)
clockwise_phi (bool)
cocos (int | None)
eq_type (str | None)
surface_interps (dict | None)
- Return type:
None
Methods
B_poloidal(R, Z)Return the magnitude of the polooidal magnetic flux density at the position(s)
(R, Z).B_radial(R, Z)Return the radial magnetic flux density at the position(s)
(R, Z).B_toroidal(R, Z)Return the toroidal magnetic flux density at the position(s)
(R, Z).B_vertical(R, Z)Return the vertical magnetic flux density at the position(s)
(R, Z).F(psi_n)Return poloidal current function \(F\) at the normalised poloidal magnetic flux function \(\psi_n\).
FF_prime(psi_n)Return poloidal current function \(F\) multiplied by its derivative with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
F_prime(psi_n)Return derivative of the poloidal current function \(f\) with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
R_major(psi_n)Return the major radius position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
R_major_prime(psi_n)Return derivative with respect to \(\psi\) of the major radius position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
Z_mid(psi_n)Return the vertical position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
Z_mid_prime(psi_n)Return the derivative with respect to \(\psi\) of the vertical position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
__init__(R, Z, psi_RZ, psi, F, FF_prime, p, ...)contour([ax, cbar, aspect, show, x_label, ...])Plot \(\psi\) over the \((R, Z)\) grid using a 2D coloured contour plot.
flux_surface(psi_n[, surface_interps, ntheta])Generate a FluxSurface object representing the flux surface with normalised poloidal magnetic flux function \(\psi_n\).
from_file(path[, file_type])Read a file from disk, returning an instance of this class.
from_netcdf(path, *args[, ...])Initialise an Equilibrium from a Pyrokinetics netCDF file.
p(psi_n)Return plasma pressure for a given normalised poloidal magnetic flux function \(\psi_n\).
p_prime(psi_n)Return derivative of the plasma pressure with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
plot(quantity[, ax, psi_n, show, x_label, ...])Plot a quantity defined on the \(\psi\) grid.
psi(psi_n)Return actual poloidal magnetic flux function \(\psi\) for a given normalised \(\psi_n\).
psi_n(psi)Return normalised poloidal magnetic flux function \(\psi_n\) for a given actual \(\psi\).
psi_tor(psi_n)Return the toroidal of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
q(psi_n)Return the safety factor for a given normalised poloidal magnetic flux function \(\psi_n\).
q_prime(psi_n)Return the derivative of the safety factor with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
r_minor(psi_n)Return half of the width of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
r_minor_prime(psi_n)Return derivative with respect to \(\psi\) of the width of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
rho(psi_n)Return the normalised minor radius of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
rho_tor(psi_n)Return the square root of the normalised toroidal flux of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
supported_file_types()Returns a list of all registered file types.
to_netcdf(*args, **kwargs)Writes self.data to disk.
upscale(R_upscale, Z_upscale)Attributes
attrsRedirects to underlying Xarray Dataset attrs.
coordsRedirects to underlying Xarray Dataset coords.
Property for managing the underlying Xarray Dataset.
data_varsRedirects to underlying Xarray Dataset data_vars.
dimsRedirects to underlying Xarray Dataset dims.
sizesRedirects to underlying Xarray Dataset sizes.
- B_poloidal(R, Z)[source]#
Return the magnitude of the polooidal magnetic flux density at the position(s)
(R, Z).- Parameters:
R (
ArrayLike,units [meter]) – Major radius positions.Z (
ArrayLike,units [meter]) – Vertical positions. Should have the same shape asR, or be broadcastable toR.
- Return type:
np.ndarray,units [tesla]
- B_radial(R, Z)[source]#
Return the radial magnetic flux density at the position(s)
(R, Z).- Parameters:
R (
ArrayLike,units [meter]) – Major radius positions.Z (
ArrayLike,units [meter]) – Vertical positions. Should have the same shape asR, or be broadcastable toR.
- Return type:
np.ndarray,units [tesla]
- B_toroidal(R, Z)[source]#
Return the toroidal magnetic flux density at the position(s)
(R, Z).- Parameters:
R (
ArrayLike,units [meter]) – Major radius positions.Z (
ArrayLike,units [meter]) – Vertical positions. Should have the same shape asR, or be broadcastable toR.
- Return type:
np.ndarray,units [tesla]
- B_vertical(R, Z)[source]#
Return the vertical magnetic flux density at the position(s)
(R, Z).- Parameters:
R (
ArrayLike,units [meter]) – Major radius positions.Z (
ArrayLike,units [meter]) – Vertical positions. Should have the same shape asR, or be broadcastable toR.
- Return type:
np.ndarray,units [tesla]
- F(psi_n)[source]#
Return poloidal current function \(F\) at the normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter * tesla]
- FF_prime(psi_n)[source]#
Return poloidal current function \(F\) multiplied by its derivative with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter**2 * tesla**2 / weber]
- F_prime(psi_n)[source]#
Return derivative of the poloidal current function \(f\) with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter * tesla / weber]
- R_major(psi_n)[source]#
Return the major radius position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter]
- R_major_prime(psi_n)[source]#
Return derivative with respect to \(\psi\) of the major radius position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter * radian / weber]
- Z_mid(psi_n)[source]#
Return the vertical position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter]
- Z_mid_prime(psi_n)[source]#
Return the derivative with respect to \(\psi\) of the vertical position of the midpoint of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter * radian / weber]
- contour(ax=None, cbar=True, aspect=False, show=False, x_label=None, y_label=None, z_label=None, **kwargs)[source]#
Plot \(\psi\) over the \((R, Z)\) grid using a 2D coloured contour plot. Uses Matplotlib’s
contourf.- Parameters:
ax (
Optional[plt.Axes]) – Axes object on which to plot. If not provided, a new figure is created.cbar (
bool, defaultTrue) – If True, builds a colourbar next to the generated plot.aspect (
bool, defaultFalse) – If True, ensures the axes have the correct aspect ratio. If the user supplies their ownax, has no effect.show (
bool, defaultFalse) – Immediately show Figure after creation.x_label (
Optional[str], defaultNone) – Overwrite the default x label. Set to an empty string""to disable.y_label (
Optional[str], defaultNone) – Overwrite the default y label. Set to an empty string""to disable.z_label (
Optional[str], defaultNone) – Overwrite the default colorbar label. Set to an empty string""to disable. Ignored ifcbaris False.**kwargs – Additional arguments to pass to Matplotlib’s
contourfcall.
- Returns:
The Axes object created after plotting.
- Return type:
plt.Axes
- flux_surface(psi_n, surface_interps=False, ntheta=256)[source]#
Generate a FluxSurface object representing the flux surface with normalised poloidal magnetic flux function \(\psi_n\). This Dataset-like object contains information such as the path swept out by the flux surface in
(R, z)coordinates, the magnetic flux density along the path, and quantities such as pressure, safety factor, the poloidal current function \(f\), and their derivatives with respect to \(\psi\) on the flux surface. It also contains derivatives with respect to the minor radius of the flux surface, and normalised versions of \(\psi\)- Return type:
- Parameters:
psi_n (float)
- psi_n: float
Normalised poloidal flux of FluxSurface to load
- surface_interps: bool
Use direct surface interpolation instead of contour fitting
- classmethod from_netcdf(path, *args, overwrite_metadata=False, overwrite_title=None, **kwargs)[source]#
Initialise an Equilibrium from a Pyrokinetics netCDF file.
This behaves like DatasetWrapper.from_netcdf, but also recreates the internal spline objects so that methods depending on them (e.g. B_radial, q, flux_surface, etc.) continue to work after reload.
- Return type:
- Parameters:
- p(psi_n)[source]#
Return plasma pressure for a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [pascal]
- p_prime(psi_n)[source]#
Return derivative of the plasma pressure with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [pascal / weber]
- plot(quantity, ax=None, psi_n=False, show=False, x_label=None, y_label=None, **kwargs)[source]#
Plot a quantity defined on the \(\psi\) grid.
- Parameters:
quantity (
str) – Name of the quantity to plot. Must be defined over the gridpsi.ax (
Optional[plt.Axes]) – Axes object on which to plot. If not provided, a new figure is created.psi_n (
bool, defaultFalse) – If True, plot against normalised \(\psi_n\). Otherwise, plot against \(\psi\).show (
bool, defaultFalse) – Immediately show Figure after creation.x_label (
Optional[str], defaultNone) – Overwrite the default x label. Set to an empty string""to disable.y_label (
Optional[str], defaultNone) – Overwrite the default y label. Set to an empty string""to disable.**kwargs – Additional arguments to pass to Matplotlib’s
plotcall.
- Returns:
The Axes object created after plotting.
- Return type:
plt.Axes- Raises:
ValueError – If
quantityis not a quantity defined over the \(\psi\) grid, or is not the name of an Equilibrium quantity.
- psi(psi_n)[source]#
Return actual poloidal magnetic flux function \(\psi\) for a given normalised \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless]) – Normalised poloidal magnetic flux function, scaled such that \(\psi_n=0\) on the magnetic axis and \(\psi_n=1\) on the last closed flux surface. Results are undefined outside this range.- Returns:
Actual poloidal magnetic flux. Note that the result is arbitrary within an additive constant.
- Return type:
np.ndarray,units [weber]
- psi_n(psi)[source]#
Return normalised poloidal magnetic flux function \(\psi_n\) for a given actual \(\psi\).
- Parameters:
psi (
ArrayLike,units [weber]) – Actual poloidal magnetic flux. Note that the result is arbitrary within an additive constant.- Returns:
Normalised poloidal magnetic flux function, scaled such that \(\psi_n=0\) on the magnetic axis and \(\psi_n=1\) on the last closed flux surface. Results are undefined outside this range.
- Return type:
np.ndarray,units [dimensionless]
- psi_tor(psi_n)[source]#
Return the toroidal of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\). This is the same as \(\psi_{tor} = \int q d\psi\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [weber]
- q(psi_n)[source]#
Return the safety factor for a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [dimensionless]
- q_prime(psi_n)[source]#
Return the derivative of the safety factor with respect to \(\psi\) for a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [weber**-1]
- r_minor(psi_n)[source]#
Return half of the width of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter]
- r_minor_prime(psi_n)[source]#
Return derivative with respect to \(\psi\) of the width of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [meter * radian / weber]
- rho(psi_n)[source]#
Return the normalised minor radius of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\). This is the same as
r_minor/a_minor.- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [dimensionless]
- rho_tor(psi_n)[source]#
Return the square root of the normalised toroidal flux of the flux surface represented by a given normalised poloidal magnetic flux function \(\psi_n\). This is the same as \(\sqrt{psi_{tor}/\psi_{tor}^{LCFS}}\).
- Parameters:
psi_n (
ArrayLike,units [dimensionless])- Return type:
np.ndarray,units [dimensionless]