bluemira.equilibria.equilibrium
Plasma MHD equilibrium and state objects
Classes
Stabilisation strategies for vertical position control |
|
Base class for magneto-hydrodynamic states |
|
Class for loading a fixed boundary plasma equilibrium. |
|
Base class for magneto-hydrodynamic states with a CoilSet |
|
Represents the breakdown state |
|
Modes for how to calculate qpsi |
|
Represents the equilibrium state, including plasma and coil currents |
Module Contents
- class bluemira.equilibria.equilibrium.VerticalPositionControlType(*args, **kwds)
Bases:
enum.EnumStabilisation strategies for vertical position control
- VIRTUAL
- FEEDBACK
- classmethod _missing_(value: str)
- Parameters:
value (str)
- class bluemira.equilibria.equilibrium.MHDState(grid: bluemira.equilibria.grid.Grid, *, o_point_fallback: bluemira.equilibria.profiles.OPointCalcOptions = OPointCalcOptions.GRID_CENTRE)
Base class for magneto-hydrodynamic states
- Parameters:
o_point_fallback (bluemira.equilibria.profiles.OPointCalcOptions)
- x: numpy.typing.NDArray[numpy.float64] | None = None
- z: numpy.typing.NDArray[numpy.float64] | None = None
- dx: float | None = None
- dz: float | None = None
- limiter: bluemira.equilibria.limiter.Limiter | None = None
- _o_point_fallback
- _label: str | None = None
- property label: str
A name used to identify the MHD state.
- Return type:
str
- set_grid(grid: bluemira.equilibria.grid.Grid)
Sets a Grid object for an Equilibrium, and sets the G-S operator and G-S solver on the grid.
- Parameters:
grid (bluemira.equilibria.grid.Grid) – The grid upon which to solve the Equilibrium
- classmethod _get_eqdsk(filename: pathlib.Path | str, from_cocos: int | None = 11, *, qpsi_positive: bool | None = None, full_coil: bool = False, regrid_nx_nz: tuple[int, int] | None = None, **kwargs) tuple[eqdsk.EQDSKInterface, bluemira.equilibria.grid.Grid, tuple[int, int] | None]
Get eqdsk data from file for read in
- Parameters:
filename (pathlib.Path | str) – Filename
from_cocos (int | None) – The COCOS index of the EQDSK file. Used when the determined COCOS is ambiguous. Will raise if given and not one of the determined COCOS indices.
qpsi_positive (bool | None) – Whether qpsi is positive or not, required for identification when qpsi is not present in the file.
full_coil (bool) – Whether the eqdsk dxc and dzc represents the full coil width or half coil width
regrid_nx_nz (tuple[int, int] | None) – Modify grid to use new nx and nz values
- Returns:
e – Instance if EQDSKInterface with the EQDSK file read in
psi – psi array
coilset – Coilset from eqdsk
grid – Grid from eqdsk
limiter – Limiter instance if any limiters are in file
- Return type:
tuple[eqdsk.EQDSKInterface, bluemira.equilibria.grid.Grid, tuple[int, int] | None]
- _prepare_eqdsk(data: dict[str, Any], filename: pathlib.Path | str, header: str = 'bluemira_equilibria', directory: str | None = None, filetype: str = 'json') eqdsk.EQDSKInterface
Prepares the Equilibrium Object eqdsk file
- Returns:
The initial EQDSKInterface data structure
- Return type:
EQDSKInterface
- Raises:
ValueError – Cant find data directory
- Parameters:
data (dict[str, Any])
filename (pathlib.Path | str)
header (str)
directory (str | None)
filetype (str)
- class bluemira.equilibria.equilibrium.FixedPlasmaEquilibrium(grid: bluemira.equilibria.grid.Grid, lcfs: bluemira.geometry.coordinates.Coordinates, profiles: bluemira.equilibria.profiles.Profile, psi: numpy.typing.NDArray[numpy.float64], psi_ax: float, psi_b: float, filename: pathlib.Path | str | None = None, *, label: str = 'Fixed Plasma Equilibrium', o_point_fallback: bluemira.equilibria.profiles.OPointCalcOptions = OPointCalcOptions.GRID_CENTRE)
Bases:
MHDStateClass for loading a fixed boundary plasma equilibrium.
- Parameters:
profiles (bluemira.equilibria.profiles.Profile)
psi (numpy.typing.NDArray[numpy.float64])
psi_ax (float)
psi_b (float)
filename (pathlib.Path | str | None)
label (str)
o_point_fallback (bluemira.equilibria.profiles.OPointCalcOptions)
- _psi
- _jtor
- profiles
- plasma
- _lcfs
- filename = None
- _label = 'Fixed Plasma Equilibrium'
- classmethod from_eqdsk(filename: pathlib.Path | str, from_cocos: int | None = 11, *, qpsi_positive: bool | None = None, full_coil: bool = False, regrid_nx_nz: tuple[int, int] | None = None, **kwargs)
Initialises a Breakdown Object from an eqdsk file. Note that this will involve recalculation of the magnetic flux.
- Parameters:
filename (pathlib.Path | str) – Filename
from_cocos (int | None) – The COCOS index of the EQDSK file. Used when the determined COCOS is ambiguous. Will raise if given and not one of the determined COCOS indices.
qpsi_positive (bool | None) – Whether qpsi is positive or not, required for identification when qpsi is not present in the file.
full_coil (bool) – Whether the eqdsk dxc and dzc represents the full coil width or half coil width
regrid_nx_nz (tuple[int, int] | None) – Modify grid to use new nx and nz values
- to_eqdsk(data, filename, header='bluemira_equilibria', directory=None, filetype='json', to_cocos=BLUEMIRA_DEFAULT_COCOS, **kwargs)
Writes the FixedPlasmaEquilibrium Object to an eqdsk file
- get_LCFS() bluemira.geometry.coordinates.Coordinates
Get the Last Closed FLux Surface (LCFS).
- Returns:
The Coordinates of the LCFS
- Return type:
- Bx(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total radial magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bx. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bx. If None, returns values at all grid points
- Returns:
Radial magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- Bz(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total vertical magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bz. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bz. If None, returns values at all grid points
- Returns:
Vertical magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- Bp(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total poloidal magnetic field at point(s) (x, z)
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bp. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bp. If None, returns values at all grid points
- Returns:
Poloidal magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- psi(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total poloidal magnetic flux, either for the whole grid, or for specified x, z coordinates.
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return psi. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return psi. If None, returns values at all grid points
- Returns:
Poloidal magnetic flux at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- plot(ax: matplotlib.axes.Axes | None = None, *, field: bluemira.equilibria.diagnostics.EqBPlotParam = EqBPlotParam.PSI)
Plots the FixedPlasmaEquilibrium object onto ax
- Returns:
the plot axis
- Parameters:
ax (matplotlib.axes.Axes | None)
- class bluemira.equilibria.equilibrium.CoilSetMHDState(grid: bluemira.equilibria.grid.Grid, coilset: bluemira.equilibria.coils.CoilSet, *, o_point_fallback: bluemira.equilibria.profiles.OPointCalcOptions = OPointCalcOptions.GRID_CENTRE)
Bases:
MHDStateBase class for magneto-hydrodynamic states with a CoilSet
- Parameters:
coilset (bluemira.equilibria.coils.CoilSet)
o_point_fallback (bluemira.equilibria.profiles.OPointCalcOptions)
- _psi_green = None
- _bx_green = None
- _bz_green = None
- _coilset
- property coilset
Equilibria coilset
- classmethod _get_eqdsk(filename: pathlib.Path | str, from_cocos: int | None = 11, *, qpsi_positive: bool | None = None, full_coil: bool = False, regrid_nx_nz: tuple[int, int] | None = None, user_coils: bluemira.equilibria.coils.CoilSet | None = None, force_symmetry: bool = False, force_symmetry_rtol: float | None = None, **kwargs) tuple[eqdsk.EQDSKInterface, bluemira.equilibria.grid.Grid, tuple[int, int] | None, bluemira.equilibria.coils.CoilSet, bluemira.equilibria.limiter.Limiter | None]
Get eqdsk data from file for read in
- Parameters:
filename (pathlib.Path | str) – Filename
from_cocos (int | None) – The COCOS index of the EQDSK file. Used when the determined COCOS is ambiguous. Will raise if given and not one of the determined COCOS indices.
qpsi_positive (bool | None) – Whether qpsi is positive or not, required for identification when qpsi is not present in the file.
full_coil (bool) – Whether the eqdsk dxc and dzc represents the full coil width or half coil width
regrid_nx_nz (tuple[int, int] | None) – Modify grid to use new nx and nz values
user_coils (bluemira.equilibria.coils.CoilSet | None) – Coilset provided by the user. Set current, j_max and b_max to zero in user_coils.
force_symmetry (bool) – Whether or not to force symmetrisation in the CoilSet
force_symmetry_rtol (float | None) – Relative tolerance used when comparing coil values. The values for the secondary coil in the pair will be set to be equal to the primary coil values if they are within force_symmetry_rtol.
- Returns:
e – Instance if EQDSKInterface with the EQDSK file read in
psi – psi array
regrid_nx_nz – modified regrid argument, possibly overwritten if previous grid states matched
coilset – Coilset from eqdsk
grid – Grid from eqdsk
limiter – Limiter instance if any limiters are in file
- Return type:
tuple[eqdsk.EQDSKInterface, bluemira.equilibria.grid.Grid, tuple[int, int] | None, bluemira.equilibria.coils.CoilSet, bluemira.equilibria.limiter.Limiter | None]
- _remap_greens()
Stores Green’s functions arrays in a dictionary of coils. Used upon initialisation and must be called after meshing of coils.
Notes
Modifies:
- ._pgreen:
Greens function coil mapping for psi
- ._bxgreen:
Greens function coil mapping for Bx
- .bzgreen:
Greens function coil mapping for Bz
- get_coil_forces() numpy.typing.NDArray[numpy.float64]
Returns the Fx and Fz force at the centre of the control coils
- Returns:
Fx, Fz array of forces on coils [N]
- Return type:
numpy.typing.NDArray[numpy.float64]
Notes
Will not work for symmetric circuits
- get_coil_fields() numpy.typing.NDArray[numpy.float64]
Returns the poloidal magnetic fields on the control coils (approximate peak at the middle inner radius of the coil)
- Returns:
The Bp array of fields on coils [T]
- Return type:
numpy.typing.NDArray[numpy.float64]
- reset_grid(grid: bluemira.equilibria.grid.Grid, psi: numpy.typing.NDArray[numpy.float64] | None = None)
Reset the grid for the MHDState.
- Parameters:
grid (bluemira.equilibria.grid.Grid) – The grid to set the MHDState on
psi (numpy.typing.NDArray[numpy.float64] | None) – Initial psi array to use
- _set_init_plasma(grid: bluemira.equilibria.grid.Grid, psi: numpy.typing.NDArray[numpy.float64] | None = None)
- Parameters:
psi (numpy.typing.NDArray[numpy.float64] | None)
- class bluemira.equilibria.equilibrium.Breakdown(coilset: bluemira.equilibria.coils.CoilSet, grid: bluemira.equilibria.grid.Grid, psi: numpy.typing.NDArray[numpy.float64] | None = None, filename: pathlib.Path | str | None = None, **kwargs)
Bases:
CoilSetMHDStateRepresents the breakdown state
- Parameters:
coilset (bluemira.equilibria.coils.CoilSet) – The set of coil objects which the equilibrium will be solved with
grid (bluemira.equilibria.grid.Grid) – The grid which to solve over
psi (numpy.typing.NDArray[numpy.float64] | None) – The initial psi array (optional)
filename (pathlib.Path | str | None) – The filename of the breakdown (optional)
- plasma
- limiter
- breakdown_point
- filename = None
- classmethod from_eqdsk(filename: pathlib.Path | str, from_cocos: int | None = 11, *, qpsi_positive: bool | None = None, full_coil: bool = False, regrid_nx_nz: tuple[int, int] | None = None, force_symmetry: bool, user_coils: bluemira.equilibria.coils.CoilSet | None = None, **kwargs)
Initialises a Breakdown Object from an eqdsk file. Note that this will involve recalculation of the magnetic flux.
- Parameters:
filename (pathlib.Path | str) – Filename
from_cocos (int | None) – The COCOS index of the EQDSK file. Used when the determined COCOS is ambiguous. Will raise if given and not one of the determined COCOS indices.
qpsi_positive (bool | None) – Whether qpsi is positive or not, required for identification when qpsi is not present in the file.
full_coil (bool) – Whether the eqdsk dxc and dzc represents the full coil width or half coil width
regrid_nx_nz (tuple[int, int] | None) – Modify grid to use new nx and nz values
force_symmetry (bool) – Whether or not to force symmetrisation in the CoilSet
user_coils (bluemira.equilibria.coils.CoilSet | None) – Coilset provided by the user. Set current, j_max and b_max to zero in user_coils.
- to_dict() dict[str, Any]
Creates a dictionary for a Breakdown object
- Returns:
A dictionary for the Breakdown object
- Return type:
dict[str, Any]
- to_eqdsk(filename: pathlib.Path | str, header: str = 'bluemira_breakdown', directory: str | None = None, filetype: str = 'json', to_cocos: int = BLUEMIRA_DEFAULT_COCOS, **kwargs)
Writes the Breakdown Object to an eqdsk file
- Parameters:
filename (pathlib.Path | str)
header (str)
directory (str | None)
filetype (str)
to_cocos (int)
- set_breakdown_point(x_bd: float, z_bd: float)
Set the point at which the centre of the breakdown region is defined.
- Parameters:
x_bd (float) – The x coordinate of the centre of the breakdown region
z_bd (float) – The z coordinate of the centre of the breakdown region
- property breakdown_psi: float
The poloidal magnetic flux at the centre of the breakdown region.
- Returns:
The minimum poloidal magnetic flux at the edge of the breakdown region [V.s/rad]
- Return type:
float
- Bx(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total radial magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bx. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bx. If None, returns values at all grid points
- Returns:
Radial magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- Bz(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total vertical magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bz. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bz. If None, returns values at all grid points
- Returns:
Vertical magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- Bp(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total poloidal magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bp. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bp. If None, returns values at all grid points
- Returns:
Poloidal magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- dBx_dz(x: numpy.typing.ArrayLike, z: numpy.typing.ArrayLike)
Total differential of the radial magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike) – Radial coordinates for which to return Bz. If None, returns values at all grid points
z (numpy.typing.ArrayLike) – Vertical coordinates for which to return Bz. If None, returns values at all grid points
- Returns:
Differential of the radial magnetic field at x, z
Notes
As there is not plasma this is equivalent to dBz_dx
- dBz_dx(x: numpy.typing.ArrayLike, z: numpy.typing.ArrayLike)
Total differential of the vertical magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike) – Radial coordinates for which to return Bz. If None, returns values at all grid points
z (numpy.typing.ArrayLike) – Vertical coordinates for which to return Bz. If None, returns values at all grid points
- Returns:
Differential of the vertical magnetic field at x, z
Notes
As there is not plasma this is equivalent to dBx_dz
- psi(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Returns the poloidal magnetic flux, either for the whole grid, or for specified x, z coordinates, including contributions from the coilset.
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return psi. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return psi. If None, returns values at all grid points
- Returns:
Poloidal magnetic flux at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- get_coil_Bp() numpy.typing.NDArray[numpy.float64]
- Returns:
The poloidal field within each coil
- Return type:
numpy.typing.NDArray[numpy.float64]
- plot(ax: matplotlib.axes.Axes | None = None, *, Bp: bool = False)
Plots the breakdown object onto ax
- Returns:
The plotting axis
- Parameters:
ax (matplotlib.axes.Axes | None)
Bp (bool)
- class bluemira.equilibria.equilibrium.QpsiCalcMode(*args, **kwds)
Bases:
enum.EnumModes for how to calculate qpsi
- Parameters:
0 – Don’t Calculate qpsi
1 – Calculate qpsi
2 – Fill qpsi grid with Zeros
- NO_CALC = 0
- CALC = 1
- ZEROS = 2
- class bluemira.equilibria.equilibrium.Equilibrium(coilset: bluemira.equilibria.coils.CoilSet, grid: bluemira.equilibria.grid.Grid, profiles: bluemira.equilibria.profiles.Profile, *, force_symmetry: bool = False, vcontrol: str | None = None, limiter: bluemira.equilibria.limiter.Limiter | None = None, psi: numpy.typing.NDArray[numpy.float64] | None = None, jtor: numpy.typing.NDArray[numpy.float64] | None = None, filename: pathlib.Path | str | None = None, label: str = 'Equilibrium', o_point_fallback: bluemira.equilibria.profiles.OPointCalcOptions = OPointCalcOptions.GRID_CENTRE)
Bases:
CoilSetMHDStateRepresents the equilibrium state, including plasma and coil currents
- Parameters:
coilset (bluemira.equilibria.coils.CoilSet) – The set of coil objects which the equilibrium will be solved with
grid (bluemira.equilibria.grid.Grid) – The grid on which to calculate the Equilibrium
profiles (bluemira.equilibria.profiles.Profile) – The plasma profiles to use in the Equilibrium
force_symmetry (bool) – Controls whether symmetry of the plasma contribution to psi across z=0 is strictly enforced in the linear system formed during solve step.
vcontrol (str | None) – Type of virtual plasma control to enact
limiter (bluemira.equilibria.limiter.Limiter | None) – Limiter conditions to apply to equilibrium
psi (numpy.typing.NDArray[numpy.float64] | None) – Magnetic flux [V.s] applied to X, Z grid
jtor (numpy.typing.NDArray[numpy.float64] | None) – The toroidal current density array of the plasma. Default = None will cause the jtor array to be constructed later as necessary.
filename (pathlib.Path | str | None) – The filename of the Equilibrium. Default = None (no file)
label (str) – The name used to identify this equilibrium
o_point_fallback (bluemira.equilibria.profiles.OPointCalcOptions)
- force_symmetry: bool = False
- _jtor = None
- profiles
- _o_points = None
- _x_points = None
- _eqdsk = None
- _li_flag: bool = False
- plasma = None
- controller = None
- limiter: bluemira.equilibria.limiter.Limiter | None = None
- filename: pathlib.Path | str | None = None
- _label: str = 'Equilibrium'
- _kwargs: dict[str, str | None]
- classmethod from_eqdsk(filename: pathlib.Path | str, from_cocos: int | None = 11, *, qpsi_positive: bool | None = None, full_coil: bool = False, regrid_nx_nz: tuple[int, int] | None = None, force_symmetry: bool = False, user_coils: bluemira.equilibria.coils.CoilSet | None = None, o_point_fallback: bluemira.equilibria.profiles.OPointCalcOptions = OPointCalcOptions.GRID_CENTRE, **kwargs)
Initialises an Equilibrium Object from an eqdsk file. Note that this will involve recalculation of the magnetic flux. Because of the nature of the (non-linear) Grad-Shafranov equation, values of psi may differ from those stored in eqdsk.
NOTE: Need to solve again with some profiles in order to re-find…
- Parameters:
filename (pathlib.Path | str) – Filename
from_cocos (int | None) – The COCOS index of the EQDSK file. Used when the determined COCOS is ambiguous. Will raise if given and not one of the determined COCOS indices.
qpsi_positive (bool | None) – Whether qpsi is positive or not, required for identification when qpsi is not present in the file.
full_coil (bool) – Whether the eqdsk dxc and dzc represents the full coil width or half coil width
regrid_nx_nz (tuple[int, int] | None) – Modify grid to use new nx and nz values
force_symmetry (bool) – Whether or not to force symmetrisation in the CoilSet
user_coils (bluemira.equilibria.coils.CoilSet | None) – Coilset provided by the user. Set current, j_max and b_max to zero in user_coils.
o_point_fallback (bluemira.equilibria.profiles.OPointCalcOptions)
- to_dict(qpsi_calcmode: int | QpsiCalcMode = QpsiCalcMode.NO_CALC, to_cocos: int = BLUEMIRA_DEFAULT_COCOS) dict[str, Any]
Creates dictionary for equilibrium object, in preparation for saving to a file format
- Parameters:
qpsi_calcmode (int | QpsiCalcMode) – don’t calculate: 0, calculate qpsi: 1, fill with zeros: 2
to_cocos (int)
- Returns:
A dictionary of the Equilibrium object values, sufficient for EQDSK
- Return type:
dict[str, Any]
- to_eqdsk(filename: pathlib.Path | str, header: str = 'bluemira_equilibria', directory: str | None = None, filetype: str = 'json', qpsi_calcmode: int | QpsiCalcMode = QpsiCalcMode.NO_CALC, to_cocos: int = BLUEMIRA_DEFAULT_COCOS, **kwargs)
Writes the Equilibrium Object to an eqdsk file
- Parameters:
filename (pathlib.Path | str)
header (str)
directory (str | None)
filetype (str)
qpsi_calcmode (int | QpsiCalcMode)
to_cocos (int)
- __getstate__()
Get the state of the Equilibrium object. Used in pickling.
- __setstate__(d)
Get the state of the Equilibrium object. Used in unpickling.
- set_grid(grid: bluemira.equilibria.grid.Grid)
Sets a Grid object for an Equilibrium, and sets the G-S operator and G-S solver on the grid.
- Parameters:
grid (bluemira.equilibria.grid.Grid) – The grid upon which to solve the Equilibrium
- _calculate_jtor(grid, psi_plasma)
Compute jtor from plasma psi.
- Returns:
Recalculated jtor.
- reset_grid(grid: bluemira.equilibria.grid.Grid, psi: numpy.typing.NDArray[numpy.float64] | None = None, **kwargs)
Reset the grid for the Equilibrium.
- Parameters:
grid (bluemira.equilibria.grid.Grid) – The grid to set the Equilibrium on
psi (numpy.typing.NDArray[numpy.float64] | None) – Psi array to use
- _set_init_plasma(grid: bluemira.equilibria.grid.Grid, psi: numpy.typing.NDArray[numpy.float64] | None = None, j_tor: numpy.typing.NDArray[numpy.float64] | None = None)
- Parameters:
psi (numpy.typing.NDArray[numpy.float64] | None)
j_tor (numpy.typing.NDArray[numpy.float64] | None)
- set_vcontrol(vcontrol_str: str | None = None)
Sets the vertical position controller
- Parameters:
vcontrol_str (str | None) – Vertical control strategy
- solve(jtor: numpy.ndarray | None = None, psi: numpy.ndarray | None = None)
Re-calculates the plasma equilibrium given new profiles
Linear Grad-Shafranov solve
- Parameters:
jtor (numpy.ndarray | None) – The toroidal current density on the finite difference grid [A/m^2]
psi (numpy.ndarray | None) – The poloidal magnetic flux on the finite difference grid [V.s/rad]
- Raises:
EquilibriaError – No O-point found
Notes
- Modifies the following in-place:
.plasma_psi .psi_func ._I_p ._Jtor
- solve_li(jtor: numpy.typing.NDArray[numpy.float64] | None = None, psi: numpy.typing.NDArray[numpy.float64] | None = None)
Optimises profiles to match input li Re-calculates the plasma equilibrium given new profiles
Linear Grad-Shafranov solve
- Parameters:
jtor (numpy.typing.NDArray[numpy.float64] | None) – The 2-D array toroidal current at each (x, z) point (optional)
psi (numpy.typing.NDArray[numpy.float64] | None) – The 2-D array of poloidal magnetic flux at each (x, z) point (optional)
- Raises:
EquilibriaError – No BetaLiIpProfile found
StopIteration – stop iterating
Notes
Modifies the following in-place:
.plasma_psi .psi_func ._I_p ._Jtor
jtor argument input is not used but kept for consistency with solve
- _update_plasma(plasma_psi: numpy.typing.NDArray[numpy.float64], j_tor: numpy.typing.NDArray[numpy.float64])
Update the plasma
- Parameters:
plasma_psi (numpy.typing.NDArray[numpy.float64])
j_tor (numpy.typing.NDArray[numpy.float64])
- _int_dxdz(func: numpy.typing.NDArray[numpy.float64]) float
Returns the double-integral of a function over the space
\(\int_Z\int_X f(x, z) dXdZ\)
- Parameters:
func (numpy.typing.NDArray[numpy.float64]) – a 2-D function map
- Returns:
The integral value of the field in 2-D
- Return type:
float
- effective_centre() tuple[float, float]
Jeon calculation for the effective current centre of the plasma
\(X_{cur}^{2}=\dfrac{1}{I_{p}}\int X^{2}J_{\phi,pl}(X, Z)d{\Omega}_{pl}\)
\(Z_{cur}=\dfrac{1}{I_{p}}\int ZJ_{\phi,pl}(X, Z)d{\Omega}_{pl}\)
- Returns:
xcur – The radial position of the effective current centre
zcur – The vertical position of the effective current centre
- Return type:
tuple[float, float]
- Bx(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total radial magnetic field at point (x, z) from coils and plasma
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bx. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bx. If None, returns values at all grid points
- Returns:
Radial magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- Bz(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total vertical magnetic field at point (x, z) from coils and plasma
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bz. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bz. If None, returns values at all grid points
- Returns:
Vertical magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- Bp(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Total poloidal magnetic field at point (x, z) from coils and plasma
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return Bp. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return Bp. If None, returns values at all grid points
- Returns:
Poloidal magnetic field at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- Bt(x: numpy.typing.ArrayLike) float | numpy.typing.NDArray[numpy.float64]
Toroidal magnetic field at point (x, z) from TF coils
- Parameters:
x (numpy.typing.ArrayLike) – Radial coordinates for which to return Bt.
- Returns:
Toroidal magnetic field at x
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- psi(x: numpy.typing.ArrayLike | None = None, z: numpy.typing.ArrayLike | None = None) float | numpy.typing.NDArray[numpy.float64]
Returns the poloidal magnetic flux, either for the whole grid, or for specified x, z coordinates, including contributions from: plasma, coilset, and vertical stabilisation controller (default None)
- Parameters:
x (numpy.typing.ArrayLike | None) – Radial coordinates for which to return psi. If None, returns values at all grid points
z (numpy.typing.ArrayLike | None) – Vertical coordinates for which to return psi. If None, returns values at all grid points
- Returns:
Poloidal magnetic flux at x, z
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- psi_norm() numpy.typing.NDArray[numpy.float64]
- Returns:
2-D x-z normalised poloidal flux map
- Return type:
numpy.typing.NDArray[numpy.float64]
- pressure_map(psi_n: float | None = None) numpy.typing.NDArray[numpy.float64]
- Parameters:
psi_n (float | None) – The normalised psi value for masking. Values outside the closed psi_n flux surface will be masked. Default is psi_n of the LCFS.
- Returns:
Plasma pressure map.
- Return type:
numpy.typing.NDArray[numpy.float64]
- dBx_dz(x: numpy.typing.ArrayLike, z: numpy.typing.ArrayLike)
Total differential of the radial magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike) – Radial coordinates for which to return Bz. If None, returns values at all grid points
z (numpy.typing.ArrayLike) – Vertical coordinates for which to return Bz. If None, returns values at all grid points
- Returns:
Differential of the radial magnetic field at x, z
- dBz_dx(x: numpy.typing.ArrayLike, z: numpy.typing.ArrayLike)
Total differential of the vertical magnetic field at point (x, z) from coils
- Parameters:
x (numpy.typing.ArrayLike) – Radial coordinates for which to return Bz. If None, returns values at all grid points
z (numpy.typing.ArrayLike) – Vertical coordinates for which to return Bz. If None, returns values at all grid points
- Returns:
Differential of the vertical magnetic field at x, z
- _get_core_mask(psi_n: float | None = None) numpy.typing.NDArray[numpy.float64]
- Parameters:
psi_n (float | None) – The normalised psi value for masking. Values outside the closed psi_n flux surface will be masked. Default is psi_n of the LCFS.
- Returns:
A 2-D masking array for the plasma core.
- Return type:
numpy.typing.NDArray[numpy.float64]
- q(psinorm: float | collections.abc.Iterable[float], o_points: collections.abc.Iterable | None = None, x_points: collections.abc.Iterable | None = None) float | numpy.typing.NDArray[numpy.float64]
- Returns:
The safety factor at given psinorm.
- Parameters:
psinorm (float | collections.abc.Iterable[float])
o_points (collections.abc.Iterable | None)
x_points (collections.abc.Iterable | None)
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- fRBpol(psinorm: numpy.typing.ArrayLike) float | numpy.typing.NDArray[numpy.float64]
- Returns:
f = R*Bt at specified values of normalised psi.
- Parameters:
psinorm (numpy.typing.ArrayLike)
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- fvac() numpy.typing.NDArray[numpy.float64]
- Returns:
The vacuum f = R*Bt.
- Return type:
numpy.typing.NDArray[numpy.float64]
- pprime(psinorm: numpy.typing.ArrayLike) float | numpy.typing.NDArray[numpy.float64]
- Returns:
p’ at given normalised psi
- Parameters:
psinorm (numpy.typing.ArrayLike)
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- ffprime(psinorm: numpy.typing.ArrayLike) float | numpy.typing.NDArray[numpy.float64]
- Returns:
ff’ at given normalised psi
- Parameters:
psinorm (numpy.typing.ArrayLike)
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- pressure(psinorm: numpy.typing.ArrayLike) float | numpy.typing.NDArray[numpy.float64]
- Returns:
Plasma pressure at specified values of normalised psi
- Parameters:
psinorm (numpy.typing.ArrayLike)
- Return type:
float | numpy.typing.NDArray[numpy.float64]
- get_flux_surface(psi_n: float, psi: numpy.typing.NDArray[numpy.float64] | None = None, o_points: collections.abc.Iterable | None = None, x_points: collections.abc.Iterable | None = None) bluemira.geometry.coordinates.Coordinates
Get a flux surface Coordinates. NOTE: Continuous surface (bridges grid)
- Parameters:
psi_n (float) – Normalised flux value of surface
psi (numpy.typing.NDArray[numpy.float64] | None) – Flux map
o_points (collections.abc.Iterable | None)
x_points (collections.abc.Iterable | None)
- Returns:
Flux surface Coordinates
- Return type:
- get_LCFS(psi: numpy.ndarray | None = None, psi_n_tol: float = 1e-06, delta_start=0.01) bluemira.geometry.coordinates.Coordinates
Get the Last Closed FLux Surface (LCFS).
- Parameters:
psi (numpy.ndarray | None) – The psi field on which to compute the LCFS. Will re-calculate if set to None
psi_n_tol (float) – The normalised psi tolerance to use when finding the LCFS
- Returns:
The Coordinates of the LCFS
- Return type:
- get_separatrix(psi: numpy.typing.NDArray[numpy.float64] | None = None, psi_n_tol: float = 1e-06) bluemira.geometry.coordinates.Coordinates | list[bluemira.geometry.coordinates.Coordinates]
Get the plasma separatrix(-ices).
- Parameters:
psi (numpy.typing.NDArray[numpy.float64] | None) – The flux array. Will re-calculate if set to None
psi_n_tol (float) – The normalised psi tolerance to use when finding the separatrix
- Returns:
The separatrix coordinates (Coordinates for SN, List[Coordinates]] for DN)
- Return type:
bluemira.geometry.coordinates.Coordinates | list[bluemira.geometry.coordinates.Coordinates]
- _clear_OX_points()
Speed optimisation for storing OX point searches in a single iteration of the solve. Large grids can cause OX finding to be expensive..
- get_OX_points(psi: numpy.typing.NDArray[numpy.float64] | None = None, *, force_update: bool = False, o_point_fallback: bluemira.equilibria.profiles.OPointCalcOptions | None = None) tuple[list[bluemira.equilibria.find.Opoint], list[bluemira.equilibria.find.Xpoint | bluemira.equilibria.find.Lpoint]]
- Returns:
The O-points
The X-points
- Parameters:
psi (numpy.typing.NDArray[numpy.float64] | None)
force_update (bool)
o_point_fallback (bluemira.equilibria.profiles.OPointCalcOptions | None)
- Return type:
tuple[list[bluemira.equilibria.find.Opoint], list[bluemira.equilibria.find.Xpoint | bluemira.equilibria.find.Lpoint]]
- get_OX_psis(psi: numpy.typing.NDArray[numpy.float64] | None = None) tuple[float, float]
The psi at the base
- Returns:
psi at O-point
psi at X-point
- Parameters:
psi (numpy.typing.NDArray[numpy.float64] | None)
- Return type:
tuple[float, float]
- get_midplane(x: float, z: float, x_psi: float) tuple[float, float]
Get the position at the midplane for a given psi value.
- Parameters:
x (float) – Starting x coordinate about which to search for a psi surface
z (float) – Starting z coordinate about which to search for a psi surface
x_psi (float) – Flux value
- Returns:
xMP – x coordinate of the midplane point with flux value Xpsi
zMP – z coordinate of the midplane point with flux value Xpsi
- Return type:
tuple[float, float]
- analyse_core(n_points: int = 50, *, plot: bool = True, ax=None) tuple[bluemira.equilibria.flux_surfaces.CoreResults, matplotlib.axes.Axes]
Analyse the shape and characteristics of the plasma core.
- Parameters:
n_points (int) – Number of points in normalised psi space to analyse
plot (bool)
- Returns:
Result dataclass
- Return type:
tuple[bluemira.equilibria.flux_surfaces.CoreResults, matplotlib.axes.Axes]
- analyse_plasma() bluemira.equilibria.physics.EqSummary
Analyse the energetic and magnetic characteristics of the plasma.
- Return type:
- analyse_coils(print_table=True) tuple[dict[str, Any], float, float]
Analyse and summarise the electro-magneto-mechanical characteristics of the equilibrium and coilset.
- Returns:
table – Tabulated currents fields and forces
fz_c_stot – The sum of the forces on the CS
fsep – CS separation force
- Return type:
tuple[dict[str, Any], float, float]
- property is_double_null: bool
Whether or not the Equilibrium is a double-null Equilibrium.
- Return type:
Whether or not the Equilibrium is a double-null Equilibrium.
- plot(ax: matplotlib.axes.Axes | None = None, *, plasma: bool = False, show_ox: bool = True)
Plot the equilibrium magnetic flux surfaces object onto ax.
- Returns:
The plot axis
- Parameters:
ax (matplotlib.axes.Axes | None)
plasma (bool)
show_ox (bool)
- plot_field(ax: matplotlib.axes.Axes | None = None, *, show_ox: bool = True)
Plot the equilibrium field structure onto ax.
- Returns:
The plot axis
- Parameters:
ax (matplotlib.axes.Axes | None)
show_ox (bool)
- plot_core(ax=None)
Plot a 1-D section through the magnetic axis.
- Returns:
The plot axis