Skip to content
Permalink
master
Switch branches/tags

Name already in use

A tag already exists with the provided branch name. Many Git commands accept both tag and branch names, so creating this branch may cause unexpected behavior. Are you sure you want to create this branch?
Go to file
 
 
Cannot retrieve contributors at this time
# This program is distributed under the terms of the GNU General Purpose License (GPL).
# Refer to http://www.gnu.org/licenses/gpl.txt
#
# This file is part of eqtools.
#
# eqtools is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# eqtools is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with eqtools. If not, see <http://www.gnu.org/licenses/>.
"""This module provides the core classes for :py:mod:`eqtools`, including the
base :py:class:`Equilibrium` class.
"""
import scipy
import scipy.interpolate
import scipy.integrate
import scipy.constants
import re
import warnings
# Python 2/3 cross compatibility:
import sys
if sys.version_info.major >= 3:
long = int
# Constants to determine how plot labels are formatted:
B_LABEL = '$B$ [T]'
J_LABEL = '$j$ [MA/m$^2$]'
class ModuleWarning(Warning):
"""Warning class to notify the user of unavailable modules.
"""
pass
try:
from . import trispline
_has_trispline = True
except ImportError:
warnings.warn("trispline module could not be loaded -- tricubic spline "
"interpolation will not be available.",
ModuleWarning)
_has_trispline = False
try:
import matplotlib.pyplot as plt
import matplotlib.widgets as mplw
import matplotlib.gridspec as mplgs
import matplotlib.patches as mpatches
import matplotlib.path as mpath
from mpl_toolkits.mplot3d import Axes3D
from mpl_toolkits.axes_grid1 import make_axes_locatable
from matplotlib.colorbar import ColorbarBase
from matplotlib.colors import Normalize
from .filewriter import gfile
except Exception:
warnings.warn(
"matplotlib modules could not be loaded -- plotting and gfile"
" writing will not be available.",
ModuleWarning
)
class PropertyAccessMixin(object):
"""Mixin to implement access of getter methods through a property-type
interface without the need to apply a decorator to every property.
For any getter `obj.getSomething()`, the call `obj.Something` will do the
same thing.
This is accomplished by overriding :py:meth:`__getattribute__` such that if
an attribute `ATTR` does not exist it then attempts to call `self.getATTR()`.
If `self.getATTR()` does not exist, an :py:class:`AttributeError` will be
raised as usual.
Also overrides :py:meth:`__setattr__` such that it will raise an
:py:class:`AttributeError` when attempting to write an attribute `ATTR` for
which there is already a method `getATTR`.
"""
def __getattribute__(self, name):
"""Get an attribute.
Tries to get attribute as-written. If this fails, tries to call the
method `get<name>` with no arguments. If this fails, raises
:py:class:`AttributeError`. This effectively generates a Python
'property' for each getter method.
Args:
name (String): Name of the attribute to retrieve. If the instance
has an attribute with this name, the attribute is returned. If
the instance does not have an attribute with this name but does
have a method called 'get'+name, this method is called and the
result is returned.
Returns:
The value of the attribute requested.
Raises:
AttributeError: If neither attribute name or method 'get'+name exist.
"""
try:
return super(Equilibrium, self).__getattribute__(name)
except AttributeError:
try:
return super(Equilibrium, self).__getattribute__('get'+name)()
except AttributeError:
raise AttributeError(
"%(class)s object has no attribute '%(n)s' or method 'get%(n)s'"
% {'class': self.__class__.__name__, 'n': name}
)
def __setattr__(self, name, value):
"""Set an attribute.
Raises :py:class:`AttributeError` if the object already has a method
'get'+name, as creation of such an attribute would interfere with the
automatic property generation in :py:meth:`__getattribute__`.
Args:
name (String): Name of the attribute to set.
value (Object): Value to set the attribute to.
Raises:
AttributeError: If a method called 'get'+name already exists.
"""
if hasattr(self, 'get'+name):
raise AttributeError(
"%(class)s object already has getter method 'get%(n)s', creating "
"attribute '%(n)s' will conflict with automatic property "
"generation." % {'class': self.__class__.__name__, 'n': name}
)
else:
super(Equilibrium, self).__setattr__(name, value)
"""The following is a dictionary to implement length unit conversions. The first
key is the unit are converting FROM, the second the unit you are converting TO.
Supports: m, cm, mm, in, ft, yd, smoot, cubit, hand
"""
_length_conversion = {'m': {'m': 1.0,
'cm': 100.0,
'mm': 1000.0,
'in': 39.37,
'ft': 39.37 / 12.0,
'yd': 39.37 / (12.0 * 3.0),
'smoot': 39.37 / 67.0,
'cubit': 39.37 / 18.0,
'hand': 39.37 / 4.0},
'cm': {'m': 0.01,
'cm': 1.0,
'mm': 10.0,
'in': 39.37 / 100.0,
'ft': 39.37 / (100.0 * 12.0),
'yd': 39.37 / (100.0 * 12.0 * 3.0),
'smoot': 39.37 / (100.0 * 67.0),
'cubit': 39.37 / (100.0 * 18.0),
'hand': 39.37 / (100.0 * 4.0)},
'mm': {'m': 0.001,
'cm': 0.1,
'mm': 1.0,
'in': 39.37 / 1000.0,
'ft': 39.37 / (1000.0 * 12.0),
'yd': 39.37 / (1000.0 * 12.0 * 3.0),
'smoot': 39.37 / (1000.0 * 67.0),
'cubit': 39.37 / (1000.0 * 18.0),
'hand': 39.37 / (1000.0 * 4.0)},
'in': {'m': 1.0 / 39.37,
'cm': 100.0 / 39.37,
'mm': 1000.0 / 39.37,
'in': 1.0,
'ft': 1.0 / 12.0,
'yd': 1.0 / (12.0 * 3.0),
'smoot': 1.0 / 67.0,
'cubit': 1.0 / 18.0,
'hand': 1.0 / 4.0},
'ft': {'m': 12.0 / 39.37,
'cm': 12.0 * 100.0 / 39.37,
'mm': 12.0 * 1000.0 / 39.37,
'in': 12.0,
'ft': 1.0,
'yd': 1.0 / 3.0,
'smoot': 12.0 / 67.0,
'cubit': 12.0 / 18.0,
'hand': 12.0 / 4.0},
'yd': {'m': 3.0 * 12.0 / 39.37,
'cm': 3.0 * 12.0 * 100.0 / 39.37,
'mm': 3.0 * 12.0 * 1000.0 / 39.37,
'in': 3.0 * 12.0,
'ft': 3.0,
'yd': 1.0,
'smoot': 3.0 * 12.0 / 67.0,
'cubit': 3.0 * 12.0 / 18.0,
'hand': 3.0 * 12.0 / 4.0},
'smoot': {'m': 67.0 / 39.37,
'cm': 67.0 * 100.0 / 39.37,
'mm': 67.0 * 1000.0 / 39.37,
'in': 67.0,
'ft': 67.0 / 12.0,
'yd': 67.0 / (12.0 * 3.0),
'smoot': 1.0,
'cubit': 67.0 / 18.0,
'hand': 67.0 / 4.0},
'cubit': {'m': 18.0 / 39.37,
'cm': 18.0 * 100.0 / 39.37,
'mm': 18.0 * 1000.0 / 39.37,
'in': 18.0,
'ft': 18.0 / 12.0,
'yd': 18.0 / (12.0 * 3.0),
'smoot': 18.0 / 67.0,
'cubit': 1.0,
'hand': 18.0 / 4.0},
'hand': {'m': 4.0 / 39.37,
'cm': 4.0 * 100.0 / 39.37,
'mm': 4.0 * 1000.0 / 39.37,
'in': 4.0,
'ft': 4.0 / 12.0,
'yd': 4.0 / (12.0 * 3.0),
'smoot': 4.0 / 67.0,
'cubit': 4.0 / 18.0,
'hand': 1.0}}
def inPolygon(polyx, polyy, pointx, pointy):
"""Function calculating whether a given point is within a 2D polygon.
Given an array of X,Y coordinates describing a 2D polygon, checks whether a
point given by x,y coordinates lies within the polygon. Operates via a
ray-casting approach - the function projects a semi-infinite ray parallel to
the positive horizontal axis, and counts how many edges of the polygon this
ray intersects. For a simply-connected polygon, this determines whether the
point is inside (even number of crossings) or outside (odd number of
crossings) the polygon, by the Jordan Curve Theorem.
Args:
polyx (Array-like): Array of x-coordinates of the vertices of the polygon.
polyy (Array-like): Array of y-coordinates of the vertices of the polygon.
pointx (Int or float): x-coordinate of test point.
pointy (Int or float): y-coordinate of test point.
Returns:
result (Boolean): True/False result for whether the point is contained within the polygon.
"""
# generator function for "lines" - pairs of (x,y) coords describing each edge of the polygon.
def lines():
p0x = polyx[-1]
p0y = polyy[-1]
p0 = (p0x, p0y)
for i, x in enumerate(polyx):
y = polyy[i]
p1 = (x, y)
yield p0, p1
p0 = p1
result = False
for p0, p1 in lines():
if (
(p0[1] > pointy) != (p1[1] > pointy)
) and (
pointx < ((p1[0]-p0[0])*(pointy-p0[1])/(p1[1]-p0[1]) + p0[0])
):
result = not result
return result
class Equilibrium(object):
"""Abstract class of data handling object for magnetic reconstruction outputs.
Defines the mapping routines and method fingerprints necessary. Each
variable or set of variables is recovered with a corresponding getter method.
Essential data for mapping are pulled on initialization (psirz grid, for
example) to frontload overhead. Additional data are pulled at the first
request and stored for subsequent usage.
.. note:: This abstract class should not be used directly. Device- and code-
specific subclasses are set up to account for inter-device/-code
differences in data storage.
Keyword Args:
length_unit (String): Sets the base unit used for any quantity whose
dimensions are length to any power. Valid options are:
=========== ===========================================================================================
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' whatever the default in the tree is (no conversion is performed, units may be inconsistent)
=========== ===========================================================================================
Default is 'm' (all units taken and returned in meters).
tspline (Boolean): Sets whether or not interpolation in time is
performed using a tricubic spline or nearest-neighbor interpolation.
Tricubic spline interpolation requires at least four complete
equilibria at different times. It is also assumed that they are
functionally correlated, and that parameters do not vary out of
their boundaries (derivative = 0 boundary condition). Default is
False (use nearest-neighbor interpolation).
monotonic (Boolean): Sets whether or not the "monotonic" form of time
window finding is used. If True, the timebase must be monotonically
increasing. Default is False (use slower, safer method).
verbose (Boolean): Allows or blocks console readout during operation.
Defaults to True, displaying useful information for the user. Set to
False for quiet usage or to avoid console clutter for multiple
instances.
Raises:
ValueError: If `length_unit` is not a valid unit specifier.
ValueError: If `tspline` is True but module trispline did not load
successfully.
"""
def __init__(
self, length_unit='m', tspline=False, monotonic=True, verbose=True
):
if (
(length_unit != 'default') and
(not (length_unit in _length_conversion))
):
raise ValueError(
"Unit '%s' not a valid unit specifier!" % length_unit
)
else:
self._length_unit = length_unit
self._tricubic = bool(tspline)
self._monotonic = bool(monotonic)
self._verbose = bool(verbose)
# These are indexes of splines, and become higher dimensional splines
# with the setting of the tspline keyword.
self._psiOfRZSpline = {}
self._phiNormSpline = {}
self._volNormSpline = {}
self._RmidSpline = {}
self._magRSpline = {}
self._magZSpline = {}
self._RmidOutSpline = {}
self._psiOfPsi0Spline = {}
self._psiOfLCFSSpline = {}
self._RmidToPsiNormSpline = {}
self._phiNormToPsiNormSpline = {}
self._volNormToPsiNormSpline = {}
self._AOutSpline = {}
self._qSpline = {}
self._FSpline = {}
self._FToPsinormSpline = {}
self._FFPrimeSpline = {}
self._pSpline = {}
self._pPrimeSpline = {}
self._vSpline = {}
self._BtVacSpline = {}
def __str__(self):
"""String representation of this instance.
Returns:
string (String): String describing this object.
"""
return 'This is an abstract class. Please use machine-specific subclass.'
def __getstate__(self):
"""Deletes all of the stored splines, since they aren't pickleable.
"""
self._psiOfRZSpline = {}
self._phiNormSpline = {}
self._volNormSpline = {}
self._RmidSpline = {}
self._magRSpline = {}
self._magZSpline = {}
self._RmidOutSpline = {}
self._psiOfPsi0Spline = {}
self._psiOfLCFSSpline = {}
self._RmidToPsiNormSpline = {}
self._phiNormToPsiNormSpline = {}
self._volNormToPsiNormSpline = {}
self._AOutSpline = {}
self._qSpline = {}
self._FSpline = {}
self._FToPsinormSpline = {}
self._FFPrimeSpline = {}
self._pSpline = {}
self._pPrimeSpline = {}
self._vSpline = {}
self._BtVacSpline = {}
return self.__dict__
####################
# Mapping routines #
####################
def rho2rho(self, origin, destination, *args, **kwargs):
r"""Convert from one coordinate to another.
Args:
origin (String): Indicates which coordinates the data are given in.
Valid options are:
======= ========================
RZ R,Z coordinates
psinorm Normalized poloidal flux
phinorm Normalized toroidal flux
volnorm Normalized volume
Rmid Midplane major radius
r/a Normalized minor radius
======= ========================
Additionally, each valid option may be prepended with 'sqrt'
to specify the square root of the desired unit.
destination (String): Indicates which coordinates to convert to.
Valid options are:
======= =================================
psinorm Normalized poloidal flux
phinorm Normalized toroidal flux
volnorm Normalized volume
Rmid Midplane major radius
r/a Normalized minor radius
q Safety factor
F Flux function :math:`F=RB_{\phi}`
FFPrime Flux function :math:`FF'`
p Pressure
pprime Pressure gradient
v Flux surface volume
======= =================================
Additionally, each valid option may be prepended with 'sqrt'
to specify the square root of the desired unit.
rho (Array-like or scalar float): Values of the starting coordinate
to map to the new coordinate. Will be two arguments `R`, `Z` if
`origin` is 'RZ'.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`rho`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `rho` (or the meshgrid of `R`
and `Z` if `make_grid` is True).
Keyword Args:
sqrt (Boolean): Set to True to return the square root of `rho`. Only
the square root of positive values is taken. Negative values are
replaced with zeros, consistent with Steve Wolfe's IDL
implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `rho` are evaluated at
each value in `t`. If True, `t` must have only one dimension (or
be a scalar). If False, `t` must match the shape of `rho` or be
a scalar. Default is True (evaluate ALL `rho` at EACH element in
`t`).
make_grid (Boolean): Only applicable if `origin` is 'RZ'. Set to
True to pass `R` and `Z` through :py:func:`scipy.meshgrid`
before evaluating. If this is set to True, `R` and `Z` must each
only have a single dimension, but can have different lengths.
Default is False (do not form meshgrid).
rho (Boolean): Set to True to return r/a (normalized minor radius)
instead of Rmid when `destination` is Rmid. Default is False
(return major radius, Rmid).
length_unit (String or 1): Length unit that quantities are
given/returned in, as applicable. If a string is given, it must
be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`rho` or (`rho`, `time_idxs`)
* **rho** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `rho`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Raises:
ValueError: If `origin` is not one of the supported values.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value at r/a=0.6, t=0.26s::
psi_val = Eq_instance.rho2rho('r/a', 'psinorm', 0.6, 0.26)
Find psinorm values at r/a points 0.6 and 0.8 at the
single time t=0.26s::
psi_arr = Eq_instance.rho2rho('r/a', 'psinorm', [0.6, 0.8], 0.26)
Find psinorm values at r/a of 0.6 at times t=[0.2s, 0.3s]::
psi_arr = Eq_instance.rho2rho('r/a', 'psinorm', 0.6, [0.2, 0.3])
Find psinorm values at (r/a, t) points (0.6, 0.2s) and (0.5, 0.3s)::
psi_arr = Eq_instance.rho2rho('r/a', 'psinorm', [0.6, 0.5], [0.2, 0.3], each_t=False)
"""
if origin.startswith('sqrt'):
args = list(args)
args[0] = scipy.asarray(args[0])**2
origin = origin[4:]
if destination.startswith('sqrt'):
kwargs['sqrt'] = True
destination = destination[4:]
if origin == 'RZ':
return self.rz2rho(destination, *args, **kwargs)
elif origin == 'Rmid':
return self.rmid2rho(destination, *args, **kwargs)
elif origin == 'r/a':
return self.roa2rho(destination, *args, **kwargs)
elif origin == 'psinorm':
return self.psinorm2rho(destination, *args, **kwargs)
elif origin == 'phinorm':
return self.phinorm2rho(destination, *args, **kwargs)
elif origin == 'volnorm':
return self.volnorm2rho(destination, *args, **kwargs)
else:
raise ValueError(
"rho2rho: Unsupported origin coordinate method '%s'!" % origin
)
def rz2psi(
self, R, Z, t, return_t=False, make_grid=False, each_t=True,
length_unit=1
):
r"""Converts the passed R, Z, t arrays to psi (unnormalized poloidal flux) values.
What is usually returned by EFIT is the stream function,
:math:`\psi=\psi_p/(2\pi)` which has units of Wb/rad.
Args:
R (Array-like or scalar float): Values of the radial coordinate to
map to poloidal flux. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `Z` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `R` must
have exactly one dimension.
Z (Array-like or scalar float): Values of the vertical coordinate to
map to poloidal flux. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `R` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `Z` must
have exactly one dimension.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R`, `Z`. If the `each_t` keyword is True, then `t` must be
scalar or have exactly one dimension. If the `each_t` keyword is
False, `t` must have the same shape as `R` and `Z` (or their
meshgrid if `make_grid` is True).
Keyword Args:
each_t (Boolean): When True, the elements in `R`, `Z` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R` and
`Z` or be a scalar. Default is True (evaluate ALL `R`, `Z` at
EACH element in `t`).
make_grid (Boolean): Set to True to pass `R` and `Z` through
:py:func:`scipy.meshgrid` before evaluating. If this is set to
True, `R` and `Z` must each only have a single dimension, but
can have different lengths. Default is False (do not form
meshgrid).
length_unit (String or 1): Length unit that `R`, `Z` are given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`psi` or (`psi`, `time_idxs`)
* **psi** (`Array or scalar float`) - The unnormalized poloidal
flux. If all of the input arguments are scalar, then a scalar is
returned. Otherwise, a scipy Array is returned. If `R` and `Z`
both have the same shape then `psi` has this shape as well,
unless the `make_grid` keyword was True, in which case `psi` has
shape (len(`Z`), len(`R`)).
* **time_idxs** (Array with same shape as `psi`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psi value at R=0.6m, Z=0.0m, t=0.26s::
psi_val = Eq_instance.rz2psi(0.6, 0, 0.26)
Find psi values at (R, Z) points (0.6m, 0m) and (0.8m, 0m) at the
single time t=0.26s. Note that the `Z` vector must be fully
specified, even if the values are all the same::
psi_arr = Eq_instance.rz2psi([0.6, 0.8], [0, 0], 0.26)
Find psi values at (R, Z) points (0.6m, 0m) at times t=[0.2s, 0.3s]::
psi_arr = Eq_instance.rz2psi(0.6, 0, [0.2, 0.3])
Find psi values at (R, Z, t) points (0.6m, 0m, 0.2s) and
(0.5m, 0.2m, 0.3s)::
psi_arr = Eq_instance.rz2psi([0.6, 0.5], [0, 0.2], [0.2, 0.3], each_t=False)
Find psi values on grid defined by 1D vector of radial positions `R`
and 1D vector of vertical positions `Z` at time t=0.2s::
psi_mat = Eq_instance.rz2psi(R, Z, 0.2, make_grid=True)
"""
# Check inputs and process into flat arrays with units of meters:
(
R, Z, t, time_idxs, unique_idxs, single_time, single_val,
original_shape
) = self._processRZt(
R,
Z,
t,
make_grid=make_grid,
each_t=each_t,
length_unit=length_unit,
compute_unique=True
)
if self._tricubic:
out_vals = scipy.reshape(
self._getFluxTriSpline().ev(t, Z, R),
original_shape
)
else:
if single_time:
out_vals = self._getFluxBiSpline(time_idxs[0]).ev(Z, R)
if single_val:
out_vals = out_vals[0]
else:
out_vals = scipy.reshape(out_vals, original_shape)
elif each_t:
out_vals = scipy.zeros(
scipy.concatenate(
([len(time_idxs), ], original_shape)
).astype(int)
)
for idx, t_idx in enumerate(time_idxs):
out_vals[idx] = self._getFluxBiSpline(
t_idx
).ev(Z, R).reshape(original_shape)
else:
out_vals = scipy.zeros_like(t, dtype=float)
for t_idx in unique_idxs:
t_mask = (time_idxs == t_idx)
out_vals[t_mask] = self._getFluxBiSpline(
t_idx
).ev(Z[t_mask], R[t_mask])
out_vals = scipy.reshape(out_vals, original_shape)
# Correct for current sign:
out_vals = -1.0 * out_vals * self.getCurrentSign()
if return_t:
if self._tricubic:
return out_vals, (t, single_time, single_val, original_shape)
else:
return out_vals, (
time_idxs, unique_idxs, single_time, single_val,
original_shape
)
else:
return out_vals
def rz2psinorm(self, R, Z, t, return_t=False, sqrt=False, make_grid=False,
each_t=True, length_unit=1):
r"""Calculates the normalized poloidal flux at the given (R, Z, t).
Uses the definition:
.. math::
\texttt{psi\_norm} = \frac{\psi - \psi(0)}{\psi(a) - \psi(0)}
Args:
R (Array-like or scalar float): Values of the radial coordinate to
map to psinorm. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `Z` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `R` must
have exactly one dimension.
Z (Array-like or scalar float): Values of the vertical coordinate to
map to psinorm. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `R` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `Z` must
have exactly one dimension.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R`, `Z`. If the `each_t` keyword is True, then `t` must be
scalar or have exactly one dimension. If the `each_t` keyword is
False, `t` must have the same shape as `R` and `Z` (or their
meshgrid if `make_grid` is True).
Keyword Args:
sqrt (Boolean): Set to True to return the square root of psinorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R`, `Z` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R` and
`Z` or be a scalar. Default is True (evaluate ALL `R`, `Z` at
EACH element in `t`).
make_grid (Boolean): Set to True to pass `R` and `Z` through
:py:func:`scipy.meshgrid` before evaluating. If this is set to
True, `R` and `Z` must each only have a single dimension, but
can have different lengths. Default is False (do not form
meshgrid).
length_unit (String or 1): Length unit that `R`, `Z` are given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`psinorm` or (`psinorm`, `time_idxs`)
* **psinorm** (`Array or scalar float`) - The normalized poloidal
flux. If all of the input arguments are scalar, then a scalar is
returned. Otherwise, a scipy Array is returned. If `R` and `Z`
both have the same shape then `psinorm` has this shape as well,
unless the `make_grid` keyword was True, in which case `psinorm`
has shape (len(`Z`), len(`R`)).
* **time_idxs** (Array with same shape as `psinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the
appropriate extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value at R=0.6m, Z=0.0m, t=0.26s::
psi_val = Eq_instance.rz2psinorm(0.6, 0, 0.26)
Find psinorm values at (R, Z) points (0.6m, 0m) and (0.8m, 0m) at the
single time t=0.26s. Note that the `Z` vector must be fully specified,
even if the values are all the same::
psi_arr = Eq_instance.rz2psinorm([0.6, 0.8], [0, 0], 0.26)
Find psinorm values at (R, Z) points (0.6m, 0m) at times t=[0.2s, 0.3s]::
psi_arr = Eq_instance.rz2psinorm(0.6, 0, [0.2, 0.3])
Find psinorm values at (R, Z, t) points (0.6m, 0m, 0.2s) and (0.5m, 0.2m, 0.3s)::
psi_arr = Eq_instance.rz2psinorm([0.6, 0.5], [0, 0.2], [0.2, 0.3], each_t=False)
Find psinorm values on grid defined by 1D vector of radial positions `R`
and 1D vector of vertical positions `Z` at time t=0.2s::
psi_mat = Eq_instance.rz2psinorm(R, Z, 0.2, make_grid=True)
"""
psi, blob = self.rz2psi(
R,
Z,
t,
return_t=True,
make_grid=make_grid,
each_t=each_t,
length_unit=length_unit
)
if self._tricubic:
psi_boundary = self._getLCFSPsiSpline()(blob[0]).reshape(blob[-1])
psi_0 = self._getPsi0Spline()(blob[0]).reshape(blob[-1])
else:
psi_boundary = self.getFluxLCFS()[blob[0]]
psi_0 = self.getFluxAxis()[blob[0]]
# If there is more than one time point, we need to expand these
# arrays to be broadcastable:
if not blob[-3]:
if each_t:
for k in range(0, len(blob[-1])):
psi_boundary = scipy.expand_dims(psi_boundary, -1)
psi_0 = scipy.expand_dims(psi_0, -1)
else:
psi_boundary = psi_boundary.reshape(blob[-1])
psi_0 = psi_0.reshape(blob[-1])
psi_norm = (psi - psi_0) / (psi_boundary - psi_0)
if sqrt:
if psi_norm.ndim == 0:
if psi_norm < 0.0:
psi_norm = 0.0
else:
scipy.place(psi_norm, psi_norm < 0, 0)
out = scipy.sqrt(psi_norm)
else:
out = psi_norm
# Unwrap single values to ensure least surprise:
if blob[-2] and blob[-3] and not self._tricubic:
out = out[0]
if return_t:
return out, blob
else:
return out
def rz2phinorm(self, *args, **kwargs):
r"""Calculates the normalized toroidal flux.
Uses the definitions:
.. math::
\texttt{phi} &= \int q(\psi)\,d\psi\\
\texttt{phi\_norm} &= \frac{\phi}{\phi(a)}
This is based on the IDL version efit_rz2rho.pro by Steve Wolfe.
Args:
R (Array-like or scalar float): Values of the radial coordinate to
map to phinorm. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `Z` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `R` must
have exactly one dimension.
Z (Array-like or scalar float): Values of the vertical coordinate to
map to phinorm. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `R` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `Z` must
have exactly one dimension.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R`, `Z`. If the `each_t` keyword is True, then `t` must be
scalar or have exactly one dimension. If the `each_t` keyword is
False, `t` must have the same shape as `R` and `Z` (or their
meshgrid if `make_grid` is True).
Keyword Args:
sqrt (Boolean): Set to True to return the square root of phinorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R`, `Z` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R` and
`Z` or be a scalar. Default is True (evaluate ALL `R`, `Z` at
EACH element in `t`).
make_grid (Boolean): Set to True to pass `R` and `Z` through
:py:func:`scipy.meshgrid` before evaluating. If this is set to
True, `R` and `Z` must each only have a single dimension, but
can have different lengths. Default is False (do not form
meshgrid).
length_unit (String or 1): Length unit that `R`, `Z` are given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting psinorm to phinorm.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`phinorm` or (`phinorm`, `time_idxs`)
* **phinorm** (`Array or scalar float`) - The normalized toroidal
flux. If all of the input arguments are scalar, then a scalar is
returned. Otherwise, a scipy Array is returned. If `R` and `Z`
both have the same shape then `phinorm` has this shape as well,
unless the `make_grid` keyword was True, in which case `phinorm`
has shape (len(`Z`), len(`R`)).
* **time_idxs** (Array with same shape as `phinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single phinorm value at R=0.6m, Z=0.0m, t=0.26s::
phi_val = Eq_instance.rz2phinorm(0.6, 0, 0.26)
Find phinorm values at (R, Z) points (0.6m, 0m) and (0.8m, 0m) at the
single time t=0.26s. Note that the `Z` vector must be fully specified,
even if the values are all the same::
phi_arr = Eq_instance.rz2phinorm([0.6, 0.8], [0, 0], 0.26)
Find phinorm values at (R, Z) points (0.6m, 0m) at times t=[0.2s, 0.3s]::
phi_arr = Eq_instance.rz2phinorm(0.6, 0, [0.2, 0.3])
Find phinorm values at (R, Z, t) points (0.6m, 0m, 0.2s) and (0.5m, 0.2m, 0.3s)::
phi_arr = Eq_instance.rz2phinorm([0.6, 0.5], [0, 0.2], [0.2, 0.3], each_t=False)
Find phinorm values on grid defined by 1D vector of radial positions `R`
and 1D vector of vertical positions `Z` at time t=0.2s::
phi_mat = Eq_instance.rz2phinorm(R, Z, 0.2, make_grid=True)
"""
return self._RZ2Quan(self._getPhiNormSpline, *args, **kwargs)
def rz2volnorm(self, *args, **kwargs):
"""Calculates the normalized flux surface volume.
Based on the IDL version efit_rz2rho.pro by Steve Wolfe.
Args:
R (Array-like or scalar float): Values of the radial coordinate to
map to volnorm. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `Z` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `R` must
have exactly one dimension.
Z (Array-like or scalar float): Values of the vertical coordinate to
map to volnorm. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `R` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `Z` must
have exactly one dimension.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R`, `Z`. If the `each_t` keyword is True, then `t` must be
scalar or have exactly one dimension. If the `each_t` keyword is
False, `t` must have the same shape as `R` and `Z` (or their
meshgrid if `make_grid` is True).
Keyword Args:
sqrt (Boolean): Set to True to return the square root of volnorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R`, `Z` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R` and
`Z` or be a scalar. Default is True (evaluate ALL `R`, `Z` at
EACH element in `t`).
make_grid (Boolean): Set to True to pass `R` and `Z` through
:py:func:`scipy.meshgrid` before evaluating. If this is set to
True, `R` and `Z` must each only have a single dimension, but
can have different lengths. Default is False (do not form
meshgrid).
length_unit (String or 1): Length unit that `R`, `Z` are given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting psinorm to volnorm.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`volnorm` or (`volnorm`, `time_idxs`)
* **volnorm** (`Array or scalar float`) - The normalized volume.
If all of the input arguments are scalar, then a scalar is
returned. Otherwise, a scipy Array is returned. If `R` and `Z`
both have the same shape then `volnorm` has this shape as well,
unless the `make_grid` keyword was True, in which case `volnorm`
has shape (len(`Z`), len(`R`)).
* **time_idxs** (Array with same shape as `volnorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single volnorm value at R=0.6m, Z=0.0m, t=0.26s::
psi_val = Eq_instance.rz2volnorm(0.6, 0, 0.26)
Find volnorm values at (R, Z) points (0.6m, 0m) and (0.8m, 0m) at the
single time t=0.26s. Note that the `Z` vector must be fully specified,
even if the values are all the same::
vol_arr = Eq_instance.rz2volnorm([0.6, 0.8], [0, 0], 0.26)
Find volnorm values at (R, Z) points (0.6m, 0m) at times t=[0.2s, 0.3s]::
vol_arr = Eq_instance.rz2volnorm(0.6, 0, [0.2, 0.3])
Find volnorm values at (R, Z, t) points (0.6m, 0m, 0.2s) and (0.5m, 0.2m, 0.3s)::
vol_arr = Eq_instance.rz2volnorm([0.6, 0.5], [0, 0.2], [0.2, 0.3], each_t=False)
Find volnorm values on grid defined by 1D vector of radial positions `R`
and 1D vector of vertical positions `Z` at time t=0.2s::
vol_mat = Eq_instance.rz2volnorm(R, Z, 0.2, make_grid=True)
"""
return self._RZ2Quan(self._getVolNormSpline, *args, **kwargs)
def rz2rmid(self, *args, **kwargs):
"""Maps the given points to the outboard midplane major radius, Rmid.
Based on the IDL version efit_rz2rmid.pro by Steve Wolfe.
Args:
R (Array-like or scalar float): Values of the radial coordinate to
map to Rmid. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `Z` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `R` must
have exactly one dimension.
Z (Array-like or scalar float): Values of the vertical coordinate to
map to Rmid. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `R` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `Z` must
have exactly one dimension.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R`, `Z`. If the `each_t` keyword is True, then `t` must be
scalar or have exactly one dimension. If the `each_t` keyword is
False, `t` must have the same shape as `R` and `Z` (or their
meshgrid if `make_grid` is True).
Keyword Args:
sqrt (Boolean): Set to True to return the square root of Rmid.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R`, `Z` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R` and
`Z` or be a scalar. Default is True (evaluate ALL `R`, `Z` at
EACH element in `t`).
make_grid (Boolean): Set to True to pass `R` and `Z` through
:py:func:`scipy.meshgrid` before evaluating. If this is set to
True, `R` and `Z` must each only have a single dimension, but
can have different lengths. Default is False (do not form
meshgrid).
rho (Boolean): Set to True to return r/a (normalized minor radius)
instead of Rmid. Default is False (return major radius, Rmid).
length_unit (String or 1): Length unit that `R`, `Z` are given in,
AND that `Rmid` is returned in. If a string is given, it must
be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting psinorm to Rmid.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`Rmid` or (`Rmid`, `time_idxs`)
* **Rmid** (`Array or scalar float`) - The outboard midplan major
radius. If all of the input arguments are scalar, then a scalar
is returned. Otherwise, a scipy Array is returned. If `R` and `Z`
both have the same shape then `Rmid` has this shape as well,
unless the `make_grid` keyword was True, in which case `Rmid`
has shape (len(`Z`), len(`R`)).
* **time_idxs** (Array with same shape as `Rmid`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single Rmid value at R=0.6m, Z=0.0m, t=0.26s::
R_mid_val = Eq_instance.rz2rmid(0.6, 0, 0.26)
Find R_mid values at (R, Z) points (0.6m, 0m) and (0.8m, 0m) at the
single time t=0.26s. Note that the `Z` vector must be fully specified,
even if the values are all the same::
R_mid_arr = Eq_instance.rz2rmid([0.6, 0.8], [0, 0], 0.26)
Find Rmid values at (R, Z) points (0.6m, 0m) at times t=[0.2s, 0.3s]::
R_mid_arr = Eq_instance.rz2rmid(0.6, 0, [0.2, 0.3])
Find Rmid values at (R, Z, t) points (0.6m, 0m, 0.2s) and (0.5m, 0.2m, 0.3s)::
R_mid_arr = Eq_instance.rz2rmid([0.6, 0.5], [0, 0.2], [0.2, 0.3], each_t=False)
Find Rmid values on grid defined by 1D vector of radial positions `R`
and 1D vector of vertical positions `Z` at time t=0.2s::
R_mid_mat = Eq_instance.rz2rmid(R, Z, 0.2, make_grid=True)
"""
# Steve Wolfe's version has an extra (linear) interpolation step for
# small psi_norm. Should check to see if we need this still with the
# scipy spline. So far looks fine...
# Convert units from meters to desired target but keep units consistent
# with rho keyword:
if kwargs.get('rho', False):
unit_factor = 1
else:
unit_factor = self._getLengthConversionFactor(
'm', kwargs.get('length_unit', 1)
)
return unit_factor * self._RZ2Quan(
self._getRmidSpline, *args, **kwargs
)
def rz2roa(self, *args, **kwargs):
"""Maps the given points to the normalized minor radius, r/a.
Based on the IDL version efit_rz2rmid.pro by Steve Wolfe.
Args:
R (Array-like or scalar float): Values of the radial coordinate to
map to r/a. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `Z` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `R` must
have exactly one dimension.
Z (Array-like or scalar float): Values of the vertical coordinate to
map to r/a. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `R` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `Z` must
have exactly one dimension.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R`, `Z`. If the `each_t` keyword is True, then `t` must be
scalar or have exactly one dimension. If the `each_t` keyword is
False, `t` must have the same shape as `R` and `Z` (or their
meshgrid if `make_grid` is True).
Keyword Args:
sqrt (Boolean): Set to True to return the square root of r/a.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R`, `Z` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R` and
`Z` or be a scalar. Default is True (evaluate ALL `R`, `Z` at
EACH element in `t`).
make_grid (Boolean): Set to True to pass `R` and `Z` through
:py:func:`scipy.meshgrid` before evaluating. If this is set to
True, `R` and `Z` must each only have a single dimension, but
can have different lengths. Default is False (do not form
meshgrid).
length_unit (String or 1): Length unit that `R`, `Z` are given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting psinorm to Rmid.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`roa` or (`roa`, `time_idxs`)
* **roa** (`Array or scalar float`) - The normalized minor radius.
If all of the input arguments are scalar, then a scalar
is returned. Otherwise, a scipy Array is returned. If `R` and `Z`
both have the same shape then `roa` has this shape as well,
unless the `make_grid` keyword was True, in which case `roa`
has shape (len(`Z`), len(`R`)).
* **time_idxs** (Array with same shape as `roa`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single r/a value at R=0.6m, Z=0.0m, t=0.26s::
roa_val = Eq_instance.rz2roa(0.6, 0, 0.26)
Find r/a values at (R, Z) points (0.6m, 0m) and (0.8m, 0m) at the
single time t=0.26s. Note that the Z vector must be fully specified,
even if the values are all the same::
roa_arr = Eq_instance.rz2roa([0.6, 0.8], [0, 0], 0.26)
Find r/a values at (R, Z) points (0.6m, 0m) at times t=[0.2s, 0.3s]::
roa_arr = Eq_instance.rz2roa(0.6, 0, [0.2, 0.3])
Find r/a values at (R, Z, t) points (0.6m, 0m, 0.2s) and (0.5m, 0.2m, 0.3s)::
roa_arr = Eq_instance.rz2roa([0.6, 0.5], [0, 0.2], [0.2, 0.3], each_t=False)
Find r/a values on grid defined by 1D vector of radial positions `R`
and 1D vector of vertical positions `Z` at time t=0.2s::
roa_mat = Eq_instance.rz2roa(R, Z, 0.2, make_grid=True)
"""
# Steve Wolfe's version has an extra (linear) interpolation step for
# small psi_norm. Should check to see if we need this still with the
# scipy spline. So far looks fine...
kwargs['rho'] = True
return self._RZ2Quan(self._getRmidSpline, *args, **kwargs)
def rz2rho(self, method, *args, **kwargs):
r"""Convert the passed (R, Z, t) coordinates into one of several coordinates.
Args:
method (String): Indicates which coordinates to convert to. Valid
options are:
======= =================================
psinorm Normalized poloidal flux
phinorm Normalized toroidal flux
volnorm Normalized volume
Rmid Midplane major radius
r/a Normalized minor radius
q Safety factor
F Flux function :math:`F=RB_{\phi}`
FFPrime Flux function :math:`FF'`
p Pressure
pprime Pressure gradient
v Flux surface volume
======= =================================
Additionally, each valid option may be prepended with 'sqrt'
to specify the square root of the desired unit.
R (Array-like or scalar float): Values of the radial coordinate to
map to `rho`. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `Z` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `R` must
have exactly one dimension.
Z (Array-like or scalar float): Values of the vertical coordinate to
map to `rho`. If `R` and `Z` are both scalar values,
they are used as the coordinate pair for all of the values in
`t`. Must have the same shape as `R` unless the `make_grid`
keyword is set. If the `make_grid` keyword is True, `Z` must
have exactly one dimension.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R`, `Z`. If the `each_t` keyword is True, then `t` must be
scalar or have exactly one dimension. If the `each_t` keyword is
False, `t` must have the same shape as `R` and `Z` (or their
meshgrid if `make_grid` is True).
Keyword Args:
sqrt (Boolean): Set to True to return the square root of `rho`.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R`, `Z` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R` and
`Z` or be a scalar. Default is True (evaluate ALL `R`, `Z` at
EACH element in `t`).
make_grid (Boolean): Set to True to pass `R` and `Z` through
:py:func:`scipy.meshgrid` before evaluating. If this is set to
True, `R` and `Z` must each only have a single dimension, but
can have different lengths. Default is False (do not form
meshgrid).
rho (Boolean): Set to True to return r/a (normalized minor radius)
instead of Rmid when `destination` is Rmid. Default is False
(return major radius, Rmid).
length_unit (String or 1): Length unit that `R`, `Z` are given in,
AND that `Rmid` is returned in. If a string is given, it must
be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`rho` or (`rho`, `time_idxs`)
* **rho** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `rho`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Raises:
ValueError: If `method` is not one of the supported values.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value at R=0.6m, Z=0.0m, t=0.26s::
psi_val = Eq_instance.rz2rho('psinorm', 0.6, 0, 0.26)
Find psinorm values at (R, Z) points (0.6m, 0m) and (0.8m, 0m) at the
single time t=0.26s. Note that the `Z` vector must be fully specified,
even if the values are all the same::
psi_arr = Eq_instance.rz2rho('psinorm', [0.6, 0.8], [0, 0], 0.26)
Find psinorm values at (R, Z) points (0.6m, 0m) at times t=[0.2s, 0.3s]::
psi_arr = Eq_instance.rz2rho('psinorm', 0.6, 0, [0.2, 0.3])
Find psinorm values at (R, Z, t) points (0.6m, 0m, 0.2s) and (0.5m, 0.2m, 0.3s)::
psi_arr = Eq_instance.rz2rho('psinorm', [0.6, 0.5], [0, 0.2], [0.2, 0.3], each_t=False)
Find psinorm values on grid defined by 1D vector of radial positions `R`
and 1D vector of vertical positions `Z` at time t=0.2s::
psi_mat = Eq_instance.rz2rho('psinorm', R, Z, 0.2, make_grid=True)
"""
if method.startswith('sqrt'):
kwargs['sqrt'] = True
method = method[4:]
if method == 'psinorm':
return self.rz2psinorm(*args, **kwargs)
elif method == 'phinorm':
return self.rz2phinorm(*args, **kwargs)
elif method == 'volnorm':
return self.rz2volnorm(*args, **kwargs)
elif method == 'Rmid':
return self.rz2rmid(*args, **kwargs)
elif method == 'r/a':
kwargs['rho'] = True
return self.rz2rmid(*args, **kwargs)
elif method == 'q':
return self.rz2q(*args, **kwargs)
elif method == 'F':
return self.rz2F(*args, **kwargs)
elif method == 'FFPrime':
return self.rz2FFPrime(*args, **kwargs)
elif method == 'p':
return self.rz2p(*args, **kwargs)
elif method == 'pprime':
return self.rz2pprime(*args, **kwargs)
elif method == 'v':
return self.rz2v(*args, **kwargs)
else:
raise ValueError(
"rz2rho: Unsupported normalized coordinate method '%s'!" % method
)
def rmid2roa(
self, R_mid, t, each_t=True, return_t=False, sqrt=False, blob=None,
length_unit=1
):
"""Convert the passed (R_mid, t) coordinates into r/a.
Args:
R_mid (Array-like or scalar float): Values of the outboard midplane
major radius to map to r/a.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R_mid`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `R_mid`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of r/a.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R_mid` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R_mid`
or be a scalar. Default is True (evaluate ALL `R_mid` at EACH
element in `t`).
length_unit (String or 1): Length unit that `R_mid` is given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`roa` or (`roa`, `time_idxs`)
* **roa** (`Array or scalar float`) - Normalized midplane minor
radius. If all of the input arguments are scalar, then a scalar
is returned. Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `roa`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single r/a value at R_mid=0.6m, t=0.26s::
roa_val = Eq_instance.rmid2roa(0.6, 0.26)
Find roa values at R_mid points 0.6m and 0.8m at the
single time t=0.26s.::
roa_arr = Eq_instance.rmid2roa([0.6, 0.8], 0.26)
Find roa values at R_mid of 0.6m at times t=[0.2s, 0.3s]::
roa_arr = Eq_instance.rmid2roa(0.6, [0.2, 0.3])
Find r/a values at (R_mid, t) points (0.6m, 0.2s) and (0.5m, 0.3s)::
roa_arr = Eq_instance.rmid2roa([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
# TODO: Make this map inboard to outboard!
# It looks like this is never actually called with pre-computed time
# indices internally, so I am going to not support that functionality
# for now.
if blob is not None:
raise NotImplementedError("Passing of time indices not supported!")
(
R_mid,
dum,
t,
time_idxs,
unique_idxs,
single_time,
single_val,
original_shape
) = self._processRZt(
R_mid,
R_mid,
t,
make_grid=False,
each_t=each_t,
length_unit=length_unit,
compute_unique=False,
convert_only=True
)
if self._tricubic:
roa = self._rmid2roa(R_mid, t).reshape(original_shape)
else:
if single_time:
roa = self._rmid2roa(R_mid, time_idxs[0])
if single_val:
roa = roa[0]
else:
roa = roa.reshape(original_shape)
elif each_t:
roa = scipy.zeros(
scipy.concatenate(
([len(time_idxs), ], original_shape)
).astype(int)
)
for idx, t_idx in enumerate(time_idxs):
roa[idx] = self._rmid2roa(
R_mid, t_idx
).reshape(original_shape)
else:
roa = self._rmid2roa(R_mid, time_idxs).reshape(original_shape)
if sqrt:
if roa.ndim == 0:
if roa < 0:
roa = 0.0
else:
scipy.place(roa, roa < 0, 0.0)
roa = scipy.sqrt(roa)
if return_t:
if self._tricubic:
return roa, (t, single_time, single_val, original_shape)
else:
return roa, (
time_idxs, unique_idxs, single_time, single_val,
original_shape
)
else:
return roa
def rmid2psinorm(self, R_mid, t, **kwargs):
"""Calculates the normalized poloidal flux corresponding to the passed R_mid (mapped outboard midplane major radius) values.
Args:
R_mid (Array-like or scalar float): Values of the outboard midplane
major radius to map to psinorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R_mid`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `R_mid`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of psinorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R_mid` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R_mid`
or be a scalar. Default is True (evaluate ALL `R_mid` at EACH
element in `t`).
length_unit (String or 1): Length unit that `R_mid` is given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`psinorm` or (`psinorm`, `time_idxs`)
* **psinorm** (`Array or scalar float`) - Normalized poloidal flux.
If all of the input arguments are scalar, then a scalar is
returned. Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `psinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value for Rmid=0.7m, t=0.26s::
psinorm_val = Eq_instance.rmid2psinorm(0.7, 0.26)
Find psinorm values at R_mid values of 0.5m and 0.7m at the single time
t=0.26s::
psinorm_arr = Eq_instance.rmid2psinorm([0.5, 0.7], 0.26)
Find psinorm values at R_mid=0.5m at times t=[0.2s, 0.3s]::
psinorm_arr = Eq_instance.rmid2psinorm(0.5, [0.2, 0.3])
Find psinorm values at (R_mid, t) points (0.6m, 0.2s) and (0.5m, 0.3s)::
psinorm_arr = Eq_instance.rmid2psinorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self._psinorm2Quan(
self._getRmidToPsiNormSpline, R_mid, t, check_space=True, **kwargs
)
def rmid2phinorm(self, *args, **kwargs):
r"""Calculates the normalized toroidal flux.
Uses the definitions:
.. math::
\texttt{phi} &= \int q(\psi)\,d\psi
\texttt{phi\_norm} &= \frac{\phi}{\phi(a)}
This is based on the IDL version efit_rz2rho.pro by Steve Wolfe.
Args:
R_mid (Array-like or scalar float): Values of the outboard midplane
major radius to map to phinorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R_mid`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `R_mid`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of phinorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R_mid` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R_mid`
or be a scalar. Default is True (evaluate ALL `R_mid` at EACH
element in `t`).
length_unit (String or 1): Length unit that `R_mid` is given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`phinorm` or (`phinorm`, `time_idxs`)
* **phinorm** (`Array or scalar float`) - Normalized toroidal flux.
If all of the input arguments are scalar, then a scalar is
returned. Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `phinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single phinorm value at R_mid=0.6m, t=0.26s::
phi_val = Eq_instance.rmid2phinorm(0.6, 0.26)
Find phinorm values at R_mid points 0.6m and 0.8m at the single time
t=0.26s::
phi_arr = Eq_instance.rmid2phinorm([0.6, 0.8], 0.26)
Find phinorm values at R_mid point 0.6m at times t=[0.2s, 0.3s]::
phi_arr = Eq_instance.rmid2phinorm(0.6, [0.2, 0.3])
Find phinorm values at (R, t) points (0.6m, 0.2s) and (0.5m, 0.3s)::
phi_arr = Eq_instance.rmid2phinorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self._Rmid2Quan(self._getPhiNormSpline, *args, **kwargs)
def rmid2volnorm(self, *args, **kwargs):
"""Calculates the normalized flux surface volume.
Based on the IDL version efit_rz2rho.pro by Steve Wolfe.
Args:
R_mid (Array-like or scalar float): Values of the outboard midplane
major radius to map to volnorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R_mid`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `R_mid`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of volnorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R_mid` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R_mid`
or be a scalar. Default is True (evaluate ALL `R_mid` at EACH
element in `t`).
length_unit (String or 1): Length unit that `R_mid` is given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`volnorm` or (`volnorm`, `time_idxs`)
* **volnorm** (`Array or scalar float`) - Normalized volume.
If all of the input arguments are scalar, then a scalar is
returned. Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `volnorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single volnorm value at R_mid=0.6m, t=0.26s::
vol_val = Eq_instance.rmid2volnorm(0.6, 0.26)
Find volnorm values at R_mid points 0.6m and 0.8m at the single time
t=0.26s::
vol_arr = Eq_instance.rmid2volnorm([0.6, 0.8], 0.26)
Find volnorm values at R_mid points 0.6m at times t=[0.2s, 0.3s]::
vol_arr = Eq_instance.rmid2volnorm(0.6, [0.2, 0.3])
Find volnorm values at (R_mid, t) points (0.6m, 0.2s) and (0.5m, 0.3s)::
vol_arr = Eq_instance.rmid2volnorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self._Rmid2Quan(self._getVolNormSpline, *args, **kwargs)
def rmid2rho(self, method, R_mid, t, **kwargs):
r"""Convert the passed (R_mid, t) coordinates into one of several coordinates.
Args:
method (String): Indicates which coordinates to convert to. Valid
options are:
======= =================================
psinorm Normalized poloidal flux
phinorm Normalized toroidal flux
volnorm Normalized volume
r/a Normalized minor radius
F Flux function :math:`F=RB_{\phi}`
FFPrime Flux function :math:`FF'`
p Pressure
pprime Pressure gradient
v Flux surface volume
======= =================================
Additionally, each valid option may be prepended with 'sqrt'
to specify the square root of the desired unit.
R_mid (Array-like or scalar float): Values of the outboard midplane
major radius to map to rho.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`R_mid`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `R_mid`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of rho.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `R_mid` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `R_mid`
or be a scalar. Default is True (evaluate ALL `R_mid` at EACH
element in `t`).
length_unit (String or 1): Length unit that `R_mid` is given in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`rho` or (`rho`, `time_idxs`)
* **rho** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `rho`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value at R_mid=0.6m, t=0.26s::
psi_val = Eq_instance.rmid2rho('psinorm', 0.6, 0.26)
Find psinorm values at R_mid points 0.6m and 0.8m at the
single time t=0.26s.::
psi_arr = Eq_instance.rmid2rho('psinorm', [0.6, 0.8], 0.26)
Find psinorm values at R_mid of 0.6m at times t=[0.2s, 0.3s]::
psi_arr = Eq_instance.rmid2rho('psinorm', 0.6, [0.2, 0.3])
Find psinorm values at (R_mid, t) points (0.6m, 0.2s) and (0.5m, 0.3s)::
psi_arr = Eq_instance.rmid2rho('psinorm', [0.6, 0.5], [0.2, 0.3], each_t=False)
"""
if method.startswith('sqrt'):
kwargs['sqrt'] = True
method = method[4:]
if method == 'psinorm':
return self.rmid2psinorm(R_mid, t, **kwargs)
elif method == 'r/a':
return self.rmid2roa(R_mid, t, **kwargs)
elif method == 'phinorm':
return self.rmid2phinorm(R_mid, t, **kwargs)
elif method == 'volnorm':
return self.rmid2volnorm(R_mid, t, **kwargs)
elif method == 'q':
return self.rmid2q(R_mid, t, **kwargs)
elif method == 'F':
return self.rmid2F(R_mid, t, **kwargs)
elif method == 'FFPrime':
return self.rmid2FFPrime(R_mid, t, **kwargs)
elif method == 'p':
return self.rmid2p(R_mid, t, **kwargs)
elif method == 'pprime':
return self.rmid2pprime(R_mid, t, **kwargs)
elif method == 'v':
return self.rmid2v(R_mid, t, **kwargs)
else:
# Default back to the old kludge that wastes time in rz2psi:
# TODO: This doesn't handle length units properly!
Z_mid = self.getMagZSpline()(t)
if kwargs.get('each_t', True):
# Need to override the default in _processRZt, since we are doing
# the shaping here:
kwargs['each_t'] = False
try:
iter(t)
except TypeError:
# For a single t, there will only be a single value of Z_mid and
# we only need to make it have the same shape as R_mid. Note
# that ones_like appears to be clever enough to handle the case
# of a scalar R_mid.
Z_mid = Z_mid * scipy.ones_like(R_mid, dtype=float)
else:
# For multiple t, we need to repeat R_mid for every t, then
# repeat the corresponding Z_mid that many times for each such
# entry.
t = scipy.asarray(t)
if t.ndim != 1:
raise ValueError(
"rmid2rho: When using the each_t keyword, "
"t must have only one dimension."
)
R_mid = scipy.tile(
R_mid,
scipy.concatenate(
(
[len(t), ],
scipy.ones_like(
scipy.shape(R_mid), dtype=float
)
)
) # may need to be declared as ints
)
# TODO: Is there a clever way to do this without a loop?
Z_mid_temp = scipy.ones_like(R_mid, dtype=float)
t_temp = scipy.ones_like(R_mid, dtype=float)
for k in range(0, len(Z_mid)):
Z_mid_temp[k] *= Z_mid[k]
t_temp[k] *= t[k]
Z_mid = Z_mid_temp
t = t_temp
return self.rz2rho(method, R_mid, Z_mid, t, **kwargs)
def roa2rmid(
self, roa, t, each_t=True, return_t=False, blob=None, length_unit=1
):
"""Convert the passed (r/a, t) coordinates into Rmid.
Args:
roa (Array-like or scalar float): Values of the normalized minor
radius to map to Rmid.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`roa`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `roa`.
Keyword Args:
each_t (Boolean): When True, the elements in `roa` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `roa`
or be a scalar. Default is True (evaluate ALL `roa` at EACH
element in `t`).
length_unit (String or 1): Length unit that `Rmid` is returned in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`Rmid` or (`Rmid`, `time_idxs`)
* **Rmid** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `Rmid`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single R_mid value at r/a=0.6, t=0.26s::
R_mid_val = Eq_instance.roa2rmid(0.6, 0.26)
Find R_mid values at r/a points 0.6 and 0.8 at the
single time t=0.26s.::
R_mid_arr = Eq_instance.roa2rmid([0.6, 0.8], 0.26)
Find R_mid values at r/a of 0.6 at times t=[0.2s, 0.3s]::
R_mid_arr = Eq_instance.roa2rmid(0.6, [0.2, 0.3])
Find R_mid values at (roa, t) points (0.6, 0.2s) and (0.5, 0.3s)::
R_mid_arr = Eq_instance.roa2rmid([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
# It looks like this is never actually called with pre-computed time
# indices internally, so I am going to not support that functionality
# for now.
if blob is not None:
raise NotImplementedError("Passing of time indices not supported!")
(
roa,
dum,
t,
time_idxs,
unique_idxs,
single_time,
single_val,
original_shape
) = self._processRZt(
roa,
roa,
t,
make_grid=False,
each_t=each_t,
length_unit=length_unit,
compute_unique=False,
check_space=False
)
if self._tricubic:
R_mid = self._roa2rmid(roa, t).reshape(original_shape)
else:
if single_time:
R_mid = self._roa2rmid(roa, time_idxs[0])
if single_val:
R_mid = R_mid[0]
else:
R_mid = R_mid.reshape(original_shape)
elif each_t:
R_mid = scipy.zeros(
scipy.concatenate(
([len(time_idxs), ], original_shape)
).astype(int)
)
for idx, t_idx in enumerate(time_idxs):
R_mid[idx] = self._roa2rmid(
roa, t_idx
).reshape(original_shape)
else:
R_mid = self._roa2rmid(roa, time_idxs).reshape(original_shape)
R_mid *= self._getLengthConversionFactor('m', length_unit)
if return_t:
if self._tricubic:
return R_mid, (t, single_time, single_val, original_shape)
else:
return R_mid, (
time_idxs, unique_idxs, single_time, single_val,
original_shape
)
else:
return R_mid
def roa2psinorm(self, *args, **kwargs):
"""Convert the passed (r/a, t) coordinates into psinorm.
Args:
roa (Array-like or scalar float): Values of the normalized minor
radius to map to psinorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`roa`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `roa`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of psinorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `roa` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `roa`
or be a scalar. Default is True (evaluate ALL `roa` at EACH
element in `t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`psinorm` or (`psinorm`, `time_idxs`)
* **psinorm** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `psinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value at r/a=0.6, t=0.26s::
psinorm_val = Eq_instance.roa2psinorm(0.6, 0.26)
Find psinorm values at r/a points 0.6 and 0.8 at the
single time t=0.26s.::
psinorm_arr = Eq_instance.roa2psinorm([0.6, 0.8], 0.26)
Find psinorm values at r/a of 0.6 at times t=[0.2s, 0.3s]::
psinorm_arr = Eq_instance.roa2psinorm(0.6, [0.2, 0.3])
Find psinorm values at (roa, t) points (0.6, 0.2s) and (0.5, 0.3s)::
psinorm_arr = Eq_instance.roa2psinorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self.roa2rho('psinorm', *args, **kwargs)
def roa2phinorm(self, *args, **kwargs):
"""Convert the passed (r/a, t) coordinates into phinorm.
Args:
roa (Array-like or scalar float): Values of the normalized minor
radius to map to phinorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`roa`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `roa`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of phinorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `roa` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `roa`
or be a scalar. Default is True (evaluate ALL `roa` at EACH
element in `t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`phinorm` or (`phinorm`, `time_idxs`)
* **phinorm** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `phinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single phinorm value at r/a=0.6, t=0.26s::
phinorm_val = Eq_instance.roa2phinorm(0.6, 0.26)
Find phinorm values at r/a points 0.6 and 0.8 at the
single time t=0.26s.::
phinorm_arr = Eq_instance.roa2phinorm([0.6, 0.8], 0.26)
Find phinorm values at r/a of 0.6 at times t=[0.2s, 0.3s]::
phinorm_arr = Eq_instance.roa2phinorm(0.6, [0.2, 0.3])
Find phinorm values at (roa, t) points (0.6, 0.2s) and (0.5, 0.3s)::
phinorm_arr = Eq_instance.roa2phinorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self.roa2rho('phinorm', *args, **kwargs)
def roa2volnorm(self, *args, **kwargs):
"""Convert the passed (r/a, t) coordinates into volnorm.
Args:
roa (Array-like or scalar float): Values of the normalized minor
radius to map to volnorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`roa`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `roa`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of volnorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `roa` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `roa`
or be a scalar. Default is True (evaluate ALL `roa` at EACH
element in `t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`volnorm` or (`volnorm`, `time_idxs`)
* **volnorm** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `volnorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single volnorm value at r/a=0.6, t=0.26s::
volnorm_val = Eq_instance.roa2volnorm(0.6, 0.26)
Find volnorm values at r/a points 0.6 and 0.8 at the
single time t=0.26s.::
volnorm_arr = Eq_instance.roa2volnorm([0.6, 0.8], 0.26)
Find volnorm values at r/a of 0.6 at times t=[0.2s, 0.3s]::
volnorm_arr = Eq_instance.roa2volnorm(0.6, [0.2, 0.3])
Find volnorm values at (roa, t) points (0.6, 0.2s) and (0.5, 0.3s)::
volnorm_arr = Eq_instance.roa2volnorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self.roa2rho('volnorm', *args, **kwargs)
def roa2rho(self, method, *args, **kwargs):
r"""Convert the passed (r/a, t) coordinates into one of several coordinates.
Args:
method (String): Indicates which coordinates to convert to.
Valid options are:
======= =================================
psinorm Normalized poloidal flux
phinorm Normalized toroidal flux
volnorm Normalized volume
Rmid Midplane major radius
q Safety factor
F Flux function :math:`F=RB_{\phi}`
FFPrime Flux function :math:`FF'`
p Pressure
pprime Pressure gradient
v Flux surface volume
======= =================================
Additionally, each valid option may be prepended with 'sqrt'
to specify the square root of the desired unit.
roa (Array-like or scalar float): Values of the normalized minor
radius to map to rho.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`roa`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `roa`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of rho.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `roa` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `roa`
or be a scalar. Default is True (evaluate ALL `roa` at EACH
element in `t`).
length_unit (String or 1): Length unit that `Rmid` is returned in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`rho` or (`rho`, `time_idxs`)
* **rho** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `rho`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value at r/a=0.6, t=0.26s::
psi_val = Eq_instance.roa2rho('psinorm', 0.6, 0.26)
Find psinorm values at r/a points 0.6 and 0.8 at the
single time t=0.26s::
psi_arr = Eq_instance.roa2rho('psinorm', [0.6, 0.8], 0.26)
Find psinorm values at r/a of 0.6 at times t=[0.2s, 0.3s]::
psi_arr = Eq_instance.roa2rho('psinorm', 0.6, [0.2, 0.3])
Find psinorm values at (r/a, t) points (0.6, 0.2s) and (0.5, 0.3s)::
psi_arr = Eq_instance.roa2rho('psinorm', [0.6, 0.5], [0.2, 0.3], each_t=False)
"""
if method == 'Rmid':
return self.roa2rmid(*args, **kwargs)
else:
kwargs['convert_roa'] = True
return self.rmid2rho(method, *args, **kwargs)
def psinorm2rmid(self, psi_norm, t, **kwargs):
"""Calculates the outboard R_mid location corresponding to the passed psinorm (normalized poloidal flux) values.
Args:
psi_norm (Array-like or scalar float): Values of the normalized
poloidal flux to map to Rmid.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`psi_norm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `psi_norm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of Rmid. Only
the square root of positive values is taken. Negative values are
replaced with zeros, consistent with Steve Wolfe's IDL
implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `psi_norm` are evaluated at
each value in `t`. If True, `t` must have only one dimension (or
be a scalar). If False, `t` must match the shape of `psi_norm` or be
a scalar. Default is True (evaluate ALL `psi_norm` at EACH element in
`t`).
rho (Boolean): Set to True to return r/a (normalized minor radius)
instead of Rmid. Default is False (return major radius, Rmid).
length_unit (String or 1): Length unit that `Rmid` is returned in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`Rmid` or (`Rmid`, `time_idxs`)
* **Rmid** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `Rmid`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single R_mid value for psinorm=0.7, t=0.26s::
R_mid_val = Eq_instance.psinorm2rmid(0.7, 0.26)
Find R_mid values at psi_norm values of 0.5 and 0.7 at the single time
t=0.26s::
R_mid_arr = Eq_instance.psinorm2rmid([0.5, 0.7], 0.26)
Find R_mid values at psi_norm=0.5 at times t=[0.2s, 0.3s]::
R_mid_arr = Eq_instance.psinorm2rmid(0.5, [0.2, 0.3])
Find R_mid values at (psinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
R_mid_arr = Eq_instance.psinorm2rmid([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
# Convert units from meters to desired target but keep units consistent
# with rho keyword:
if kwargs.get('rho', False):
unit_factor = 1
else:
unit_factor = self._getLengthConversionFactor(
'm', kwargs.get('length_unit', 1)
)
return unit_factor * self._psinorm2Quan(
self._getRmidSpline,
psi_norm,
t,
**kwargs
)
def psinorm2roa(self, psi_norm, t, **kwargs):
"""Calculates the normalized minor radius location corresponding to the passed psi_norm (normalized poloidal flux) values.
Args:
psi_norm (Array-like or scalar float): Values of the normalized
poloidal flux to map to r/a.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`psi_norm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `psi_norm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of r/a. Only
the square root of positive values is taken. Negative values are
replaced with zeros, consistent with Steve Wolfe's IDL
implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `psi_norm` are evaluated at
each value in `t`. If True, `t` must have only one dimension (or
be a scalar). If False, `t` must match the shape of `psi_norm` or be
a scalar. Default is True (evaluate ALL `psi_norm` at EACH element in
`t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`roa` or (`roa`, `time_idxs`)
* **roa** (`Array or scalar float`) - Normalized midplane minor
radius. If all of the input arguments are scalar, then a scalar
is returned. Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `roa`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single r/a value for psinorm=0.7, t=0.26s::
roa_val = Eq_instance.psinorm2roa(0.7, 0.26)
Find r/a values at psi_norm values of 0.5 and 0.7 at the single time
t=0.26s::
roa_arr = Eq_instance.psinorm2roa([0.5, 0.7], 0.26)
Find r/a values at psi_norm=0.5 at times t=[0.2s, 0.3s]::
roa_arr = Eq_instance.psinorm2roa(0.5, [0.2, 0.3])
Find r/a values at (psinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
roa_arr = Eq_instance.psinorm2roa([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
kwargs['rho'] = True
return self._psinorm2Quan(self._getRmidSpline, psi_norm, t, **kwargs)
def psinorm2volnorm(self, psi_norm, t, **kwargs):
"""Calculates the normalized volume corresponding to the passed psi_norm (normalized poloidal flux) values.
Args:
psi_norm (Array-like or scalar float): Values of the normalized
poloidal flux to map to volnorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`psi_norm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `psi_norm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of volnorm. Only
the square root of positive values is taken. Negative values are
replaced with zeros, consistent with Steve Wolfe's IDL
implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `psi_norm` are evaluated at
each value in `t`. If True, `t` must have only one dimension (or
be a scalar). If False, `t` must match the shape of `psi_norm` or be
a scalar. Default is True (evaluate ALL `psi_norm` at EACH element in
`t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`volnorm` or (`volnorm`, `time_idxs`)
* **volnorm** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `volnorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single volnorm value for psinorm=0.7, t=0.26s::
volnorm_val = Eq_instance.psinorm2volnorm(0.7, 0.26)
Find volnorm values at psi_norm values of 0.5 and 0.7 at the single time
t=0.26s::
volnorm_arr = Eq_instance.psinorm2volnorm([0.5, 0.7], 0.26)
Find volnorm values at psi_norm=0.5 at times t=[0.2s, 0.3s]::
volnorm_arr = Eq_instance.psinorm2volnorm(0.5, [0.2, 0.3])
Find volnorm values at (psinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
volnorm_arr = Eq_instance.psinorm2volnorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self._psinorm2Quan(
self._getVolNormSpline, psi_norm, t, **kwargs
)
def psinorm2phinorm(self, psi_norm, t, **kwargs):
"""Calculates the normalized toroidal flux corresponding to the passed psi_norm (normalized poloidal flux) values.
Args:
psi_norm (Array-like or scalar float): Values of the normalized
poloidal flux to map to phinorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`psi_norm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `psi_norm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of phinorm. Only
the square root of positive values is taken. Negative values are
replaced with zeros, consistent with Steve Wolfe's IDL
implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `psi_norm` are evaluated at
each value in `t`. If True, `t` must have only one dimension (or
be a scalar). If False, `t` must match the shape of `psi_norm` or be
a scalar. Default is True (evaluate ALL `psi_norm` at EACH element in
`t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`phinorm` or (`phinorm`, `time_idxs`)
* **phinorm** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `phinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single phinorm value for psinorm=0.7, t=0.26s::
phinorm_val = Eq_instance.psinorm2phinorm(0.7, 0.26)
Find phinorm values at psi_norm values of 0.5 and 0.7 at the single time
t=0.26s::
phinorm_arr = Eq_instance.psinorm2phinorm([0.5, 0.7], 0.26)
Find phinorm values at psi_norm=0.5 at times t=[0.2s, 0.3s]::
phinorm_arr = Eq_instance.psinorm2phinorm(0.5, [0.2, 0.3])
Find phinorm values at (psinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
phinorm_arr = Eq_instance.psinorm2phinorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self._psinorm2Quan(self._getPhiNormSpline, psi_norm, t, **kwargs)
def psinorm2rho(self, method, *args, **kwargs):
r"""Convert the passed (psinorm, t) coordinates into one of several coordinates.
Args:
method (String): Indicates which coordinates to convert to.
Valid options are:
======= =================================
phinorm Normalized toroidal flux
volnorm Normalized volume
Rmid Midplane major radius
r/a Normalized minor radius
q Safety factor
F Flux function :math:`F=RB_{\phi}`
FFPrime Flux function :math:`FF'`
p Pressure
pprime Pressure gradient
v Flux surface volume
======= =================================
Additionally, each valid option may be prepended with 'sqrt'
to specify the square root of the desired unit.
psi_norm (Array-like or scalar float): Values of the normalized
poloidal flux to map to rho.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`psi_norm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `psi_norm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of rho. Only
the square root of positive values is taken. Negative values are
replaced with zeros, consistent with Steve Wolfe's IDL
implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `psi_norm` are evaluated at
each value in `t`. If True, `t` must have only one dimension (or
be a scalar). If False, `t` must match the shape of `psi_norm` or be
a scalar. Default is True (evaluate ALL `psi_norm` at EACH element in
`t`).
rho (Boolean): Set to True to return r/a (normalized minor radius)
instead of Rmid. Default is False (return major radius, Rmid).
length_unit (String or 1): Length unit that `Rmid` is returned in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`rho` or (`rho`, `time_idxs`)
* **rho** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `rho`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Raises:
ValueError: If `method` is not one of the supported values.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single phinorm value at psinorm=0.6, t=0.26s::
phi_val = Eq_instance.psinorm2rho('phinorm', 0.6, 0.26)
Find phinorm values at phinorm of 0.6 and 0.8 at the
single time t=0.26s::
phi_arr = Eq_instance.psinorm2rho('phinorm', [0.6, 0.8], 0.26)
Find phinorm values at psinorm of 0.6 at times t=[0.2s, 0.3s]::
phi_arr = Eq_instance.psinorm2rho('phinorm', 0.6, [0.2, 0.3])
Find phinorm values at (psinorm, t) points (0.6, 0.2s) and (0.5m, 0.3s)::
phi_arr = Eq_instance.psinorm2rho('phinorm', [0.6, 0.5], [0.2, 0.3], each_t=False)
"""
if method.startswith('sqrt'):
kwargs['sqrt'] = True
method = method[4:]
if method == 'phinorm':
return self.psinorm2phinorm(*args, **kwargs)
elif method == 'volnorm':
return self.psinorm2volnorm(*args, **kwargs)
elif method == 'Rmid':
return self.psinorm2rmid(*args, **kwargs)
elif method == 'r/a':
kwargs['rho'] = True
return self.psinorm2rmid(*args, **kwargs)
elif method == 'q':
return self.psinorm2q(*args, **kwargs)
elif method == 'F':
return self.psinorm2F(*args, **kwargs)
elif method == 'FFPrime':
return self.psinorm2FFPrime(*args, **kwargs)
elif method == 'p':
return self.psinorm2p(*args, **kwargs)
elif method == 'pprime':
return self.psinorm2pprime(*args, **kwargs)
elif method == 'v':
return self.psinorm2v(*args, **kwargs)
else:
raise ValueError(
"psinorm2rho: Unsupported normalized coordinate method '%s'!" % method
)
def phinorm2psinorm(self, phinorm, t, **kwargs):
"""Calculates the normalized poloidal flux corresponding to the passed phinorm (normalized toroidal flux) values.
Args:
phinorm (Array-like or scalar float): Values of the normalized
toroidal flux to map to psinorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`phinorm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `phinorm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of psinorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `phinorm` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `phinorm`
or be a scalar. Default is True (evaluate ALL `phinorm` at EACH
element in `t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`psinorm` or (`psinorm`, `time_idxs`)
* **psinorm** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `psinorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single psinorm value for phinorm=0.7, t=0.26s::
psinorm_val = Eq_instance.phinorm2psinorm(0.7, 0.26)
Find psinorm values at phinorm values of 0.5 and 0.7 at the single time
t=0.26s::
psinorm_arr = Eq_instance.phinorm2psinorm([0.5, 0.7], 0.26)
Find psinorm values at phinorm=0.5 at times t=[0.2s, 0.3s]::
psinorm_arr = Eq_instance.phinorm2psinorm(0.5, [0.2, 0.3])
Find psinorm values at (phinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
psinorm_arr = Eq_instance.phinorm2psinorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self._psinorm2Quan(
self._getPhiNormToPsiNormSpline, phinorm, t, **kwargs
)
def phinorm2volnorm(self, *args, **kwargs):
"""Calculates the normalized flux surface volume corresponding to the passed phinorm (normalized toroidal flux) values.
Args:
phinorm (Array-like or scalar float): Values of the normalized
toroidal flux to map to volnorm.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`phinorm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `phinorm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of volnorm.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `phinorm` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `phinorm`
or be a scalar. Default is True (evaluate ALL `phinorm` at EACH
element in `t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`volnorm` or (`volnorm`, `time_idxs`)
* **volnorm** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `volnorm`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single volnorm value for phinorm=0.7, t=0.26s::
volnorm_val = Eq_instance.phinorm2volnorm(0.7, 0.26)
Find volnorm values at phinorm values of 0.5 and 0.7 at the single time
t=0.26s::
volnorm_arr = Eq_instance.phinorm2volnorm([0.5, 0.7], 0.26)
Find volnorm values at phinorm=0.5 at times t=[0.2s, 0.3s]::
volnorm_arr = Eq_instance.phinorm2volnorm(0.5, [0.2, 0.3])
Find volnorm values at (phinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
volnorm_arr = Eq_instance.phinorm2volnorm([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
return self._phinorm2Quan(self._getVolNormSpline, *args, **kwargs)
def phinorm2rmid(self, *args, **kwargs):
"""Calculates the mapped outboard midplane major radius corresponding to the passed phinorm (normalized toroidal flux) values.
Args:
phinorm (Array-like or scalar float): Values of the normalized
toroidal flux to map to Rmid.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`phinorm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `phinorm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of Rmid.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `phinorm` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `phinorm`
or be a scalar. Default is True (evaluate ALL `phinorm` at EACH
element in `t`).
rho (Boolean): Set to True to return r/a (normalized minor radius)
instead of Rmid. Default is False (return major radius, Rmid).
length_unit (String or 1): Length unit that `Rmid` is returned in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
value is 1 (use meters).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`Rmid` or (`Rmid`, `time_idxs`)
* **Rmid** (`Array or scalar float`) - The converted coordinates. If
all of the input arguments are scalar, then a scalar is returned.
Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `Rmid`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single Rmid value for phinorm=0.7, t=0.26s::
Rmid_val = Eq_instance.phinorm2rmid(0.7, 0.26)
Find Rmid values at phinorm values of 0.5 and 0.7 at the single time
t=0.26s::
Rmid_arr = Eq_instance.phinorm2rmid([0.5, 0.7], 0.26)
Find Rmid values at phinorm=0.5 at times t=[0.2s, 0.3s]::
Rmid_arr = Eq_instance.phinorm2rmid(0.5, [0.2, 0.3])
Find Rmid values at (phinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
Rmid_arr = Eq_instance.phinorm2rmid([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
# Convert units from meters to desired target but keep units consistent
# with rho keyword:
if kwargs.get('rho', False):
unit_factor = 1
else:
unit_factor = self._getLengthConversionFactor(
'm', kwargs.get('length_unit', 1)
)
return unit_factor * self._phinorm2Quan(
self._getRmidSpline,
*args,
**kwargs
)
def phinorm2roa(self, phi_norm, t, **kwargs):
"""Calculates the normalized minor radius corresponding to the passed phinorm (normalized toroidal flux) values.
Args:
phinorm (Array-like or scalar float): Values of the normalized
toroidal flux to map to r/a.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`phinorm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `phinorm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of r/a.
Only the square root of positive values is taken. Negative
values are replaced with zeros, consistent with Steve Wolfe's
IDL implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `phinorm` are evaluated
at each value in `t`. If True, `t` must have only one dimension
(or be a scalar). If False, `t` must match the shape of `phinorm`
or be a scalar. Default is True (evaluate ALL `phinorm` at EACH
element in `t`).
k (positive int): The degree of polynomial spline interpolation to
use in converting coordinates.
return_t (Boolean): Set to True to return a tuple of (`rho`,
`time_idxs`), where `time_idxs` is the array of time indices
actually used in evaluating `rho` with nearest-neighbor
interpolation. (This is mostly present as an internal helper.)
Default is False (only return `rho`).
Returns:
`roa` or (`roa`, `time_idxs`)
* **roa** (`Array or scalar float`) - Normalized midplane minor
radius. If all of the input arguments are scalar, then a scalar
is returned. Otherwise, a scipy Array is returned.
* **time_idxs** (Array with same shape as `roa`) - The indices
(in :py:meth:`self.getTimeBase`) that were used for
nearest-neighbor interpolation. Only returned if `return_t` is
True.
Examples:
All assume that `Eq_instance` is a valid instance of the appropriate
extension of the :py:class:`Equilibrium` abstract class.
Find single r/a value for phinorm=0.7, t=0.26s::
roa_val = Eq_instance.phinorm2roa(0.7, 0.26)
Find r/a values at phinorm values of 0.5 and 0.7 at the single time
t=0.26s::
roa_arr = Eq_instance.phinorm2roa([0.5, 0.7], 0.26)
Find r/a values at phinorm=0.5 at times t=[0.2s, 0.3s]::
roa_arr = Eq_instance.phinorm2roa(0.5, [0.2, 0.3])
Find r/a values at (phinorm, t) points (0.6, 0.2s) and (0.5, 0.3s)::
roa_arr = Eq_instance.phinorm2roa([0.6, 0.5], [0.2, 0.3], each_t=False)
"""
kwargs['rho'] = True
return self._phinorm2Quan(self._getRmidSpline, phi_norm, t, **kwargs)
def phinorm2rho(self, method, *args, **kwargs):
r"""Convert the passed (phinorm, t) coordinates into one of several coordinates.
Args:
method (String): Indicates which coordinates to convert to.
Valid options are:
======= =================================
psinorm Normalized poloidal flux
volnorm Normalized volume
Rmid Midplane major radius
r/a Normalized minor radius
q Safety factor
F Flux function :math:`F=RB_{\phi}`
FFPrime Flux function :math:`FF'`
p Pressure
pprime Pressure gradient
v Flux surface volume
======= =================================
Additionally, each valid option may be prepended with 'sqrt'
to specify the square root of the desired unit.
phinorm (Array-like or scalar float): Values of the normalized
toroidal flux to map to rho.
t (Array-like or scalar float): Times to perform the conversion at.
If `t` is a single value, it is used for all of the elements of
`phinorm`. If the `each_t` keyword is True, then `t` must be scalar
or have exactly one dimension. If the `each_t` keyword is False,
`t` must have the same shape as `phinorm`.
Keyword Args:
sqrt (Boolean): Set to True to return the square root of rho. Only
the square root of positive values is taken. Negative values are
replaced with zeros, consistent with Steve Wolfe's IDL
implementation efit_rz2rho.pro. Default is False.
each_t (Boolean): When True, the elements in `phinorm` are evaluated at
each value in `t`. If True, `t` must have only one dimension (or
be a scalar). If False, `t` must match the shape of `phinorm` or be
a scalar. Default is True (evaluate ALL `phinorm` at EACH element in
`t`).
rho (Boolean): Set to True to return r/a (normalized minor radius)
instead of Rmid. Default is False (return major radius, Rmid).
length_unit (String or 1): Length unit that `Rmid` is returned in.
If a string is given, it must be a valid unit specifier:
=========== ===========
'm' meters
'cm' centimeters
'mm' millimeters
'in' inches
'ft' feet
'yd' yards
'smoot' smoots
'cubit' cubits
'hand' hands
'default' meters
=========== ===========
If length_unit is 1 or None, meters are assumed. The default
va