xbout.boutdataset.BoutDatasetAccessor#
- class xbout.boutdataset.BoutDatasetAccessor(ds)[source]#
Bases:
object
Contains BOUT-specific methods to use on BOUT++ datasets opened using
open_boutdataset()
.These BOUT-specific methods and attributes are accessed via the bout accessor, e.g.
ds.bout.options
returns aBoutOptionsFile
instance.Methods
__init__
(ds)Add Cartesian (X,Y,Z) coordinates.
animate_list
(variables[, animate_over, ...])- type variables:
Create a new Dataset with all 3d variables transformed to non-field-aligned coordinates
from_region
(name[, with_guards])Get a logically-rectangular section of data from a certain region.
get_bounding_surfaces
([coords])Get bounding surfaces.
get_field_aligned
(name[, caching])Get a field-aligned version of a variable, calculating (and caching in the Dataset) if necessary
integrate_midpoints
(variable, *[, dims, ...])Integrate using the midpoint rule for spatial dimensions, and trapezium rule for time.
interpolate_from_unstructured
(variables, *)Interpolate Dataset onto new grids of some existing coordinates
interpolate_parallel
(variables, **kwargs)Interpolate in the parallel direction to get a higher resolution version of a subset of variables.
interpolate_to_cartesian
([nX, nY, nZ, ...])Interpolate the Dataset to a regular Cartesian grid.
remove_yboundaries
(**kwargs)Remove y-boundary points, if present, from the Dataset
save
([savepath, filetype, variables, ...])Save data variables to a netCDF file.
Create a new Dataset with all 3d variables transformed to field-aligned coordinates, which are shifted with respect to the base coordinates by an angle zShift
to_restart
([variables, savepath, nxpe, ...])Write out a timestep as a set of netCDF BOUT.restart files.
Attributes
The default factor to increase resolution when doing parallel interpolation
- add_cartesian_coordinates()[source]#
Add Cartesian (X,Y,Z) coordinates.
- Returns:
Dataset with new coordinates added
, which are named'X_cartesian'
,
'Y_cartesian'
, and'Z_cartesian'
- animate_list(variables, animate_over=None, save_as=None, show=False, fps=10, nrows=None, ncols=None, poloidal_plot=False, axis_coords=None, subplots_adjust=None, vmin=None, vmax=None, logscale=None, titles=None, aspect=None, extend=None, controls='both', tight_layout=True, **kwargs)[source]#
- Parameters:
variables (
list
ofstr
orBoutDataArray
) – The variables to plot. For any string passed, the corresponding variable in this DataSet is used - then the calling DataSet must have only 3 dimensions. It is possible to pass BoutDataArrays to allow more flexible plots, e.g. with different variables being plotted against different axes.animate_over (
str
, optional) – Dimension over which to animate, defaults to the time dimensionsave_as (
str
, optional) – If passed, a gif is created with this filenameshow (
bool
, optional) – Call pyplot.show() to display the animationfps (
float
, optional) – Indicates the number of frames per second to playnrows (
int
, optional) – Specify the number of rows of plotsncols (
int
, optional) – Specify the number of columns of plotspoloidal_plot (
bool
orsequence
ofbool
, optional) – If set to True, make all 2D animations in the poloidal plane instead of using grid coordinates, per variable if sequence is givenaxis_coords (
None
,str
,dict
orlist
ofNone
,str
ordict
) –Coordinates to use for axis labelling.
None: Use the dimension coordinate for each axis, if it exists.
”index”: Use the integer index values.
dict: keys are dimension names, values set axis_coords for each axis separately. Values can be: None, “index”, the name of a 1d variable or coordinate (which must have the dimension given by ‘key’), or a 1d numpy array, dask array or DataArray whose length matches the length of the dimension given by ‘key’.
Only affects time coordinate for plots with poloidal_plot=True. If a list is passed, it must have the same length as ‘variables’ and gives the axis_coords setting for each plot individually. The setting to use for the ‘animate_over’ coordinate can be passed in one or more dict values, but must be the same in all dicts if given more than once.
subplots_adjust (
dict
, optional) – Arguments passed to fig.subplots_adjust()()vmin (
float
orsequence
offloats
) – Minimum value for color scale, per variable if a sequence is givenvmax (
float
orsequence
offloats
) – Maximum value for color scale, per variable if a sequence is givenlogscale (
bool
orfloat
,sequence
ofbool
orfloat
, optional) – If True, default to a logarithmic color scale instead of a linear one. If a non-bool type is passed it is treated as a float used to set the linear threshold of a symmetric logarithmic scale as linthresh=min(abs(vmin),abs(vmax))*logscale, defaults to 1e-5 if True is passed. Per variable if sequence is given.titles (
sequence
ofstr
orNone
, optional) – Custom titles for each plot. Pass None in the sequence to use the default for a certain variableaspect (
str
orNone
, orsequence
ofstr
orNone
, optional) – Argument to set_aspect() for each plot. Defaults to “equal” for poloidal plots and “auto” for others.controls (
string
orNone
, default"both"
) – By default, add both the timeline and play/pause toggle to the animation. If “timeline” is passed add only the timeline, if “toggle” is passed add only the play/pause toggle. If None or an empty string is passed, add neither.tight_layout (
bool
ordict
, optional) – If set to False, don’t call tight_layout() on the figure. If a dict is passed, the dict entries are passed as arguments to tight_layout()**kwargs (
dict
, optional) – Additional keyword arguments are passed on to each animation function, per variable if a sequence is given.
- Returns:
An animatplot.Animation object.
- Return type:
animation
- property fine_interpolation_factor#
The default factor to increase resolution when doing parallel interpolation
- from_field_aligned()[source]#
Create a new Dataset with all 3d variables transformed to non-field-aligned coordinates
- from_region(name, with_guards=None)[source]#
Get a logically-rectangular section of data from a certain region. Includes guard cells from neighbouring regions.
- get_bounding_surfaces(coords=('R', 'Z'))[source]#
Get bounding surfaces. Surfaces are returned as arrays of points describing a polygon, assuming the third spatial dimension is a symmetry direction.
- Parameters:
coords (
(str
,str)
, default(``
”R”, ``"Z"
)
) – Pair of names of coordinates whose values are used to give the positions of the points in the result- Returns:
result – Each DataArray in the list contains points on a boundary, with size (<number of points in the bounding polygon>, 2). Points wind clockwise around the outside domain, and anti-clockwise around the inside (if there is an inner boundary).
- Return type:
list
ofDataArrays
- get_field_aligned(name, caching=True)[source]#
Get a field-aligned version of a variable, calculating (and caching in the Dataset) if necessary
- integrate_midpoints(variable, *, dims=None, cumulative_t=False)[source]#
Integrate using the midpoint rule for spatial dimensions, and trapezium rule for time.
The quantity being integrated is assumed to be a scalar variable.
When doing a 1d integral in the ‘y’ dimension, the integral is calculated as a poloidal integral if the variable is on the standard grid (
direction_y
attribute is “Standard”), or as a parallel-to-B integral if the variable is on the field-aligned grid (direction_y
attribute is “Aligned”).When doing a 2d integral over ‘x’ and ‘y’ dimensions, the integral will be over poloidal cross-sections if the variable is not field-aligned (
direction_y == "Standard"
) and over field-aligned surfaces if the variable is field-aligned (direction_ == "Aligned"
). The latter seems unlikely to be useful as the surfaces depend on the arbitrary origin used for zShift.Is a method of BoutDataset accessor rather than of BoutDataArray so we can use other variables like \(J\), \(g^{11}\), \(g_{22}\) for the integration.
Note the xarray.DataArray.integrate() method uses the trapezium rule, which is not consistent with the way BOUT++ defines grid spacings as cell widths. Also, this way for example:
inner = da.isel(x=slice(i)).bout.integrate_midpoints() outer = da.isel(x=slice(i, None).bout.integrate_midpoints() total = da.bout.integrate_midpoints() inner + outer == total
while with the trapezium rule you would have to select
radial=slice(i+1)
for inner to get a similar relation to be true.- Parameters:
variable (
str
orDataArray
) – Name of the variable to integrate, or the variable itself as a DataArray.dims (
str
,list
ofstr
or...
) – Dimensions to integrate over. Can be any combination of of the dimensions of the Dataset. Defaults to integration over all spatial dimensions. If...
is passed, integrate over all dimensions including time.cumulative_t (
bool
, defaultFalse
) – If integrating in time, return the cumulative integral (integral from the beginning up to each point in the time dimension) instead of the definite integral.
- interpolate_from_unstructured(variables, *, fill_value=nan, structured_output=True, unstructured_dim_name='unstructured_dim', **kwargs)[source]#
Interpolate Dataset onto new grids of some existing coordinates
- Parameters:
variables (
str
orsequence
ofstr
or...
) – The names of the variables to interpolate. If ‘variables=…’ is passed explicitly, then interpolate all variables in the Dataset.**kwargs (
(str
,array)
) – Each keyword is the name of a coordinate in the DataArray, the argument is a 1d array giving the values of that coordinate on the output gridfill_value (
float
) – fill_value passed through to scipy.interpolation.griddatastructured_output (
bool
, defaultTrue
) – If True, treat output coordinates values as a structured grid. If False, output coordinate values must all have the same length and are not broadcast together.unstructured_dim_name (
str
, default"unstructured_dim"
) – Name used for the dimension in the output that replaces the dimensions of the interpolated coordinates. Only used if structured_output=False.
- Returns:
Dataset interpolated onto a new, structured grid
- Return type:
Dataset
- interpolate_parallel(variables, **kwargs)[source]#
Interpolate in the parallel direction to get a higher resolution version of a subset of variables.
Note that the high-resolution variables are all loaded into memory, so most likely it is necessary to select only a small number. The toroidal_points argument can also be used to reduce the memory demand.
- Parameters:
variables (
str
orsequence
ofstr
or...
) – The names of the variables to interpolate. If ‘variables=…’ is passed explicitly, then interpolate all variables in the Dataset.n (
int
, optional) – The factor to increase the resolution by. Defaults to the value set by BoutDataset.setupParallelInterp(), or 10 if that has not been called.toroidal_points (
int
orsequence
ofint
, optional) – If int, number of toroidal points to output, applies a stride to toroidal direction to save memory usage. If sequence of int, the indexes of toroidal points for the output.method (
str
, optional) – The interpolation method to use. Options from xarray.DataArray.interp(), currently: linear, nearest, zero, slinear, quadratic, cubic. Default is ‘cubic’.
- Returns:
A new Dataset containing a high-resolution versions
ofthe variables. The new
Dataset is a valid BoutDataset
,although containing only the specified variables.
- interpolate_to_cartesian(nX=300, nY=300, nZ=100, *, use_float32=True, fill_value=nan)[source]#
Interpolate the Dataset to a regular Cartesian grid.
This method is intended to be used to produce data for visualisation, which normally does not require double-precision values, so by default the data is converted to
numpy.float32
. Passuse_float32=False
to retain the original precision.- Parameters:
nX (
int (default 300)
) – Number of grid points in the X directionnY (
int (default 300)
) – Number of grid points in the Y directionnZ (
int (default 100)
) – Number of grid points in the Z directionuse_float32 (
bool (default True)
) – Downgrade precision tonumpy.float32
?fill_value (
float (default np.nan)
) – Value to use for points outside the interpolation domain (passed toscipy.interpolate.RegularGridInterpolator
)
See also
BoutDataArray.interpolate_to_cartesian
- save(savepath='./boutdata.nc', filetype='NETCDF4', variables=None, save_dtype=None, separate_vars=False, pre_load=False)[source]#
Save data variables to a netCDF file.
- Parameters:
savepath (
str
, optional) –filetype (
str
, optional) –variables (
list
ofstr
, optional) – Variables from the dataset to save. Default is to save all of them.separate_vars (
bool
, optional) – If this is true then every variable which depends on time (but not solely on time) will be saved into a different output file. The files are labelled by the name of the variable. Variables which don’t meet this criterion will be present in every output file.pre_load (
bool
, optional) – When saving separate variables, will load each variable into memory before saving to file, which can be considerably faster.
Examples
If
separate_vars=True
, then multiple files will be created. These can all be opened and merged in one go using a call of the form:ds = xr.open_mfdataset('boutdata_*.nc', combine='nested', concat_dim=None)
- to_field_aligned()[source]#
Create a new Dataset with all 3d variables transformed to field-aligned coordinates, which are shifted with respect to the base coordinates by an angle zShift
- to_restart(variables=None, *, savepath='.', nxpe=None, nype=None, tind=-1, prefix='BOUT.restart', overwrite=False)[source]#
Write out a timestep as a set of netCDF BOUT.restart files.
If processor decomposition is not specified then data will be saved using the decomposition it had when loaded.
- Parameters:
variables (
str
orsequence
ofstr
, optional) – The evolving variables needed in the restart files. If not given explicitly, all time-evolving variables in the Dataset will be used, which may result in larger restart files than necessary. If there is no time-dimension in the Dataset (e.g. if it was loaded from restart files), then all variables will be added if this argument is not given explicitly.savepath (
str
, default'.'
) – Directory to save the created restart files undernxpe (
int
, optional) – Number of processors in the x-direction. If not given, keep the number used for the original simulationnype (
int
, optional) – Number of processors in the y-direction. If not given, keep the number used for the original simulationtind (
int
, default-1
) – Time-index of the slice to write to the restart files. Note, when creating restart files from ‘dump’ files it is recommended to open the Dataset using the full time range and use thetind
argument here, rather than selecting a time point manually, so that the calculation ofhist_hi
in the output can be correct (which requires knowing the existing value ofhist_hi
(output step count at the end of the simulation),tind
and the total number of time points in the current output data).prefix (
str
, default"BOUT.restart"
) – Prefix to use for names of restart filesoverwrite (
bool
, defaultFalse
) – By default, raises if restart file already exists. Set to True to overwrite existing files