spacepy.pybats.bats.Bats2d

class spacepy.pybats.bats.Bats2d(filename, *args, blocksize=8, **kwargs)[source]

A child class of IdlFile tailored to 2D BATS-R-US output.

__init__(filename, *args, blocksize=8, **kwargs)[source]

Base class for “Data Model” representation data

Abstract method, reimplement

Attributes:
attrsdict

dictionary of the attributes of the SpaceData object

Methods

add_b_magsphere([target, loc, style, ...])

Create an array of field lines closed to the central body in the domain.

add_body([ax, facecolor, DoPlanet, ang])

Creates a circle of radius=self.attrs['rbody'] and returns the MatPlotLib Ellipse patch object for plotting.

add_contour(dim1, dim2, value[, nlev, ...])

Create a contour plot of variable value against dim1 on the x-axis and dim2 on the y-axis.

add_grid_plot([target, loc, do_label, ...])

Create a plot of the grid resolution by coloring regions of constant resolution.

add_pcolor(dim1, dim2, value[, zlim, ...])

Create a pcolor plot of variable value against dim1 on the x-axis and dim2 on the y-axis.

add_planet([ax, rad, ang])

Creates a circle of radius=self.attrs['rbody'] and returns the MatPlotLib Ellipse patch object for plotting.

add_stream_scatter(xcomp, ycomp[, nlines, ...])

Add a set of stream traces to a figure or axes that are distributed evenly but randomly throughout the plot domain.

calc_E()

Calculates the MHD electric field, -UxB.

calc_alfven()

Calculate the Alfven speed, B/(mu*rho)^(1/2) in km/s.

calc_all([exclude])

Perform all variable calculations (e.g. calculations that begin with "calc").

calc_b()

Calculates total B-field strength using all three B components.

calc_beta()

Calculates plasma beta (ratio of plasma to magnetic pressure, indicative of who - plasma or B-field - is "in charge" with regards to flow.

calc_gradP()

Calculate the pressure gradient force.

calc_j()

Calculates total current density strength using all three J components.

calc_jxb()

Calculates the JxB force assuming: -Units of J are uA/m2, units of B are nT.

calc_ndens()

Calculate number densities for each fluid.

calc_temp([units])

Calculate plasma temperature for each fluid.

calc_upar()

Calculate $vec{U} cdot vec{B}$ and store as 'upar' in self.

calc_uperp()

Calculate the magnitude of the velocity perpendicular to the magnetic field: $vec{U} times hat{b}$.

calc_utotal()

Calculate bulk velocity magnitude: $u^2 = u_X^2 + u_Y^2 + u_Z^2$.

calc_vort([conv])

Calculate the vorticity (curl of bulk velocity) for the direction orthogonal to the cut plane.

cfl(dt[, xcoord, ycoord])

Calculate the CFL number in each cell given time step dt.

extract(x, y, **kwargs)

For x, y of a 1D curve, extract values along that curve and return slice as a new Extraction object.

find_earth_lastclosed([tol, method, ...])

For Y=0 cuts, attempt to locate the last-closed magnetic field line for both day- and night-sides.

get_stream(x, y, xvar, yvar[, method, ...])

Trace a 2D streamline through the domain, returning a Stream object to the caller.

gradP_regular([cellsize, dim1range, dim2range])

Calculate pressure gradient on a regular grid.

gyroradius([velocities, m_avg])

Calculate the ion gyroradius in each cell.

inertial_length([m_avg])

Calculate the ion inertial length.

plasma_freq([m_avg])

Calculate the ion plasma frequency.

regrid([cellsize, dim1range, dim2range, debug])

Re-bin data to regular grid of spacing cellsize.

switch_frame(*args, **kwargs)

For files that have more than one data frame (i.e., *.outs files), load data from the iframe-th frame into the object replacing what is currently loaded.

vth([m_avg])

Calculate the thermal velocity.

Attributes

add_b_magsphere(target=None, loc=111, style='mag', DoLast=True, DoOpen=True, compX='bx', compY='bz', narrow=0, arrsize=12, method='rk4', tol=0.004363323129985824, DoClosed=True, colors=None, linestyles=None, maxPoints=1000000.0, nOpen=5, nClosed=15, arrstyle='->', add_body=False, **kwargs)[source]

Create an array of field lines closed to the central body in the domain. Add these lines to Matplotlib target object target. If no target is specified, a new figure and axis are created.

Note that this should currently only be used for GSM y=0 cuts of the magnetosphere.

A tuple containing the figure, axes, and LineCollection object is returned.

Basic styling (color and linestyle) can be handled with the style, colors, and linestyles kwargs. style can accept style names as defined in Stream, which colors and styles lines based on characteristics (e.g., open, closed). The default is ‘mag’, which colors open lines black and closed lines white. Alternatively, this kwarg works in a similar manner as it does in plot(), i.e., a string code such as “b-” (a solid blue line) or ‘r:’ (a dotted red line), etc. Both colors and linestyles work much as they do for LineCollection, but only a single value (not a list or tuple) should be provided. colors can be a CSS4 color name, an RGB tuple, or a string hex code. linestyles can be the name of the style (e.g., “dashed”) or a shortcut compatable with the style kwarg (e.g., “–“). See the documentation for the associated Matplotlib classes & functions to see all options. Note that linestyles and colors override style.

If the styling kwargs are used, they will set the colors for all lines except last-closed boundaries. Users may control groups individually using multiple calls and plotting one group at a time. Note that colors and linestyles kwargs will override style; colors allows for more flexibility concerning color choice.

Algorithm: This method, unlike its predecessor, starts by finding the last closed field lines via find_earth_lastclosed(). It then fills the regions between the open and closed regions. Currently, it does not treat purely IMF field lines.

Kwarg

Description

target

The figure or axes to place the resulting lines.

style

The color coding system for field lines. Defaults to ‘mag’. See spacepy.pybats.bats.Stream. Because lines are added as a LineCollection, only certain styles are allowed (i.e., line styles only, no marker styles).

loc

The location of the subplot on which to place the lines.

DoLast

Plot last-closed lines as red lines. Defaults to True.

DoOpen

Plot open field lines. Defaults to True.

DoClosed

Plot closed field lines. Defaults to True.

nOpen

Number of closed field lines to trace per hemisphere. Defaults to 5.

nClosed

Number of open field lines to trace per hemisphere. Defaults to 15.

narrow

Add “n” arrows to each line to indicate direction. Default is zero, or no arrows.

arrstyle

Set the arrow style in the same manner as Matplotlib’s annotate function. Default is ‘->’.

arrsize

Set the size, in points, of each directional arrow. Default is 12.

method

The tracing method; defaults to ‘rk4’. See spacepy.pybats.bats.Stream.

tol

Tolerance for finding open-closed boundary; see find_earth_lastclosed().

compX

Name of x-variable through which to trace, defaults to ‘bx’.

compY

Name of y-variable through which to trace, defaults to ‘bz’.

colors

Matplotlib-compatable color name (single) to apply to lines.

maxPoints

Set the maximum number of points in an field line integration. Defaults to one million.

linestyles

A single line style indicator, defaults to ‘-‘; see LineCollection for possible options.

add_body

Add an planet/body to the figure. Defaults to False.

Extra kwargs are passed to Matplotlib’s LineCollection class as described above.

Returns:
figmatplotlib Figure object
axmatplotlib Axes object
collectmatplotlib Collection object of trace results

Examples

>>> import matplotlib.pyplot as plt
>>> from spacepy.pybats import bats
>>> # Open a 2D slice, add a pressure contour.
>>> # Example file in spacepy/tests/data/pybats_test/:
>>> mhd = bats.Bats2d('./y0_binary.out')
>>> mhd.add_contour('x','z','p')
>>> # Add field lines using default styling:
>>> mhd.add_b_magsphere(target=plt.gca())
>>> # Add a subset of lines using custom styling:
>>> mhd.add_b_magsphere(target=plt.gca(), DoLast=False, DoOpen=False, style='g--')
add_body(ax=None, facecolor='lightgrey', DoPlanet=True, ang=0.0, **extra_kwargs)[source]

Creates a circle of radius=self.attrs[‘rbody’] and returns the MatPlotLib Ellipse patch object for plotting. If an axis is specified using the “ax” keyword, the patch is added to the plot. Default color is light grey; extra keywords are handed to the Ellipse generator function.

Because the body is rarely the size of the planet at the center of the modeling domain, add_planet is automatically called. This can be negated by using the DoPlanet kwarg.

add_contour(dim1, dim2, value, nlev=30, target=None, loc=111, title=None, xlabel=None, ylabel=None, ylim=None, xlim=None, add_cbar=False, clabel=None, filled=True, add_body=True, dolog=False, zlim=None, *args, **kwargs)[source]

Create a contour plot of variable value against dim1 on the x-axis and dim2 on the y-axis.

Simple example:

>>> self.add_contour('x', 'y', 'rho')

If kwarg target is None (default), a new figure is generated from scratch. If target is a matplotlib Figure object, a new axis is created to fill that figure at subplot location loc. If target is a matplotlib Axes object, the plot is placed into that axis.

Four values are returned: the matplotlib Figure and Axes objects, the matplotlib contour object, and the matplotlib colorbar object (defaults to False if not used.)

Kwarg

Description

target

Set plot destination. Defaults to new figure.

loc

Set subplot location. Defaults to 111.

nlev

Number of contour levels. Defaults to 30.

title

Sets title of axes. Default is no title.

xlabel

Sets x label of axes. Defaults to dim1 and units.

ylabel

Sets y label of axes. Defaults to dim2 and units.

xlim

Sets limits of x-axes. Defaults to whole domain.

ylim

Sets limits of y-axes. Defaults to whole domain.

zlim

Sets contour range. Defaults to variable max/min.

add_cbar

Adds colorbar to plot. Defaults to False.

clabel

Sets colorbar label. Defaults to var and units.

add_body

Places planetary body in plot. Defaults to True.

dolog

Sets use of logarithmic scale. Defaults to False.

add_grid_plot(target=None, loc=111, do_label=True, do_fill=True, show_nums=False, show_borders=True, cmap='jet_r', title='BATS-R-US Grid Layout')[source]

Create a plot of the grid resolution by coloring regions of constant resolution. Kwarg “target” specifies where to place the plot and can be a figure, an axis, or None. If target is a figure, a new subplot is created. The subplot location is set by the kwarg “loc”, which defaults to 111. If target is an axis, the plot is placed into that axis object. If target is None, a new figure and axis are created and used to display the plot.

Resolution labels can be disabled by setting kwarg do_label to False.

Plot title is set using the ‘title’ kwarg, defaults to ‘BATS-R-US Grid Layout’.

Note that if target is not an axis object, the axis will automatically flip the direction of positive X GSM and turn on equal aspect ratio. In other words, if you want a customized axis, it’s best to create one yourself.

Figure and axis, even if none given, are returned.

Parameters:
targetMatplotlib Figure or Axes object

Set plot destination. Defaults to new figure.

loc3-digit integer

Set subplot location. Defaults to 111.

do_labelbool

Adds resolution legend to righthand margin. Defaults to True.

do_fillbool

If true, fill blocks with color indicating grid resolution.

show_numsbool

Adds quadtree values to each region. Useful for debugging. Defaults to False.

show_bordersTrue

Show black borders around each quadtree leaf. Defaults to True.

cmapstring

Sets the color map used to color each region. Must be a Matplotlib named colormap. Defaults to ‘jet_r’.

titlestring

Sets the title at the top of the plot. Defaults to ‘BATS-R-US Grid Layout’.

add_pcolor(dim1, dim2, value, zlim=None, target=None, loc=111, title=None, xlabel=None, ylabel=None, ylim=None, xlim=None, add_cbar=False, clabel=None, add_body=True, dolog=False, *args, **kwargs)[source]

Create a pcolor plot of variable value against dim1 on the x-axis and dim2 on the y-axis. Pcolor plots shade each computational cell with the value at the cell center. Because no interpolation or smoothing is used in the visualization, pcolor plots are excellent for examining the raw output.

Simple example:

>>> self.add_pcolor('x', 'y', 'rho')

If kwarg target is None (default), a new figure is generated from scratch. If target is a matplotlib Figure object, a new axis is created to fill that figure at subplot location loc. If target is a matplotlib Axes object, the plot is placed into that axis.

Four values are returned: the matplotlib Figure and Axes objects, the matplotlib contour object, and the matplotlib colorbar object (defaults to False if not used.)

Kwarg

Description

target

Set plot destination. Defaults to new figure.

loc

Set subplot location. Defaults to 111.

title

Sets title of axes. Default is no title.

xlabel

Sets x label of axes. Defaults to dim1 and units.

ylabel

Sets y label of axes. Defaults to dim2 and units.

xlim

Sets limits of x-axes. Defaults to whole domain.

ylim

Sets limits of y-axes. Defaults to whole domain.

zlim

Sets color bar range. Defaults to variable max/min.

add_cbar

Adds colorbar to plot. Defaults to False.

clabel

Sets colorbar label. Defaults to var and units.

add_body

Places planetary body in plot. Defaults to True.

dolog

Sets use of logarithmic scale. Defaults to False.

add_planet(ax=None, rad=1.0, ang=0.0, **extra_kwargs)[source]

Creates a circle of radius=self.attrs[‘rbody’] and returns the MatPlotLib Ellipse patch object for plotting. If an axis is specified using the “ax” keyword, the patch is added to the plot.

Unlike the add_body method, the circle is colored half white (dayside) and half black (nightside) to coincide with the direction of the sun. Additionally, because the size of the planet is not intrinsically known to the MHD file, the kwarg “rad”, defaulting to 1.0, sets the size of the planet.

Extra keywords are handed to the Ellipse generator function.

add_stream_scatter(xcomp, ycomp, nlines=100, target=None, loc=111, method='rk4', maxPoints=1000000.0, xlim=None, ylim=None, narrow=0, arrsize=12, arrstyle='->', start_points=None, **kwargs)[source]

Add a set of stream traces to a figure or axes that are distributed evenly but randomly throughout the plot domain.

Lines will be seeded randomly over a given spatial range given by xlim and ylim OR the range of the axes (if target is set to a non-empty axes object) OR over the entire object domain (in that order).

Extra keyword args are handed to matplotlib’s LineCollection object: matplotlib.collections.LineCollection.

Parameters:
xcompstring

The first component of the vector field to trace (e.g., ‘bx’).

ycompstring

The second component of the vector field to trace (e.g., ‘bz’).

Returns:
figmatplotlib Figure object
axmatplotlib Axes object
collectmatplotlib Collection object of trace results
start_pointsnlines x 2 array of line starting points
Other Parameters:
targetMatplotlib Figure or Axes object

Set plot destination. Defaults to new figure.

loc3-digit integer

Set subplot location. Defaults to 111.

xlimTwo-element list/tuple

Set the range in 1st dimension over which lines will be seeded.

ylimTwo-element list/tuple

Set the range in 2nd dimension over which lines will be seeded.

nlinesint

Number of stream lines to create; default is 100.

start_pointsnlinesx2 array

Set start_points to define starting location of traces instead of using random points. This is useful for creating timeseries of plots.

maxPointsint

Set the maximum number of points in a single trace. Defaults to 1E6.

narrowint

Add “n” arrows to each line to indicate direction. Default is zero, or no lines. If narrow=1, arrows will be placed at start_points.

arrstylestring

Set the arrow style in the same manner as Matplotlib’s annotate function. Default is ‘->’.

arrsizeint

Set the size, in points, of each directional arrow. Default is 12.

calc_E()[source]

Calculates the MHD electric field, -UxB. Works for default MHD units of nT and km/s; if these units are not correct, an exception will be raised. Stores E in mV/m.

Values are saved as self[‘Ex’], self[‘Ey’], self[‘Ez’], and self[‘E’].

calc_alfven()[source]

Calculate the Alfven speed, B/(mu*rho)^(1/2) in km/s. This is performed for each species and the total fluid. The variable is saved under key “alfven” in self.data.

calc_all(exclude=[])[source]

Perform all variable calculations (e.g. calculations that begin with “calc”). Any exceptions raised by functions that could not be peformed (typicaly from missing variables) are discarded.

calc_b()[source]

Calculates total B-field strength using all three B components. Retains units of components. Additionally, the unit vector b-hat is calculated as well.

calc_beta()[source]

Calculates plasma beta (ratio of plasma to magnetic pressure, indicative of who - plasma or B-field - is “in charge” with regards to flow. Assumes: -pressure in units of nPa -B in units of nT. Resulting value is unitless. Values where b_total = zero are set to -1.0 in the final array.

Values are stored as self[‘beta’].

calc_gradP()[source]

Calculate the pressure gradient force.

calc_j()[source]

Calculates total current density strength using all three J components. Retains units of components, stores in self[‘j’]

calc_jxb()[source]

Calculates the JxB force assuming: -Units of J are uA/m2, units of B are nT. Under these assumptions, the value calculated is force density with units of nN/m^3.

Values are stored as ‘jb’ (magnitude) and ‘jbx’, ‘jby’, ‘jbz’.

calc_ndens()[source]

Calculate number densities for each fluid. Species mass is ascertained via recognition of fluid name (e.g. OpRho is clearly oxygen). A full list of recognized fluids/species can be found by exploring the dictionary mass found in bats. Composition is also calculated as percent of total number density.

New values are saved using the keys speciesN (e.g. opN) and speciesFrac (e.g. opFrac).

calc_temp(units='eV')[source]

Calculate plasma temperature for each fluid. Number density is calculated using calc_ndens if it hasn’t been done so already.

Temperature is obtained via density and pressure through the simple relationship P=nkT.

Use the units kwarg to set output units. Current choices are KeV, eV, and K. Default is eV.

calc_upar()[source]

Calculate $vec{U} cdot vec{B}$ and store as ‘upar’ in self. Result maintains units of velocity. Values are calculated for each fluid.

Values arestored as self[‘u_par’]; species prefixes are used for multi-species/fluid files.

Notes

Added in version 0.5.0.

calc_uperp()[source]

Calculate the magnitude of the velocity perpendicular to the magnetic field: $vec{U} times hat{b}$. Result maintains units of velocity.

Values are calculated for each fluid and stored as self[‘u_perp’] (using species prefixes as needed).

Notes

Added in version 0.5.0.

calc_utotal()[source]

Calculate bulk velocity magnitude: $u^2 = u_X^2 + u_Y^2 + u_Z^2$. This is done on a per-fluid basis.

calc_vort(conv=0.0001569612305760477)[source]

Calculate the vorticity (curl of bulk velocity) for the direction orthogonal to the cut plane. For example, if output file is a cut in the equatorial plane (GSM X-Y plane), only the z-component of the curl is calculated.

Output is saved as self[‘wD’], where D is the resulting dimension (e.g., ‘wz’ for the z-component of vorticity in the X-Y plane).

Vorticity is calculated for total fluid and by-species for multi-fluid output files.

Parameters:
None
Other Parameters:
convfloat

Required unit conversion such that output units are 1/s. Defaults to 1/RE (in km), which assumes grid is in RE and velocity is in km/s.

cfl(dt, xcoord='x', ycoord='z')[source]

Calculate the CFL number in each cell given time step dt.

Result is stored in self[‘cfl’].

extract(x, y, **kwargs)[source]

For x, y of a 1D curve, extract values along that curve and return slice as a new Extraction object. Valid keyword arguments are the same as for Extraction.

find_earth_lastclosed(tol=0.008726646259971648, method='rk4', max_iter=100, debug=False)[source]

For Y=0 cuts, attempt to locate the last-closed magnetic field line for both day- and night-sides. This is done using a bisection approach to precisely locate the transition between open and closed geometries. The method stops once this transition is found within a latitudinal tolerance of tol, which defaults to $pi/360.$, or one-half degree. The tracing method can be set via keyword and defaults to ‘rk4’ (4th order Runge Kutta, see Stream for more information). The maximum number of iterations the algorithm will take is set by max_iter, which defaults to 100. Latitudinal footprints of the last closed field lines at the inner boundary (not the ionosphere!) are also returned.

This method returns 5 objects:

  • The dipole tilt in radians

  • A tuple of the northern/southern hemisphere polar angle of footpoints for the dayside last-closed field line.

  • A tuple of the northern/southern hemisphere polar angle of footpoints for the nightside last-closed field line.

  • The dayside last-closed field line as a Stream object.

  • The nightside last-closed field line as a Stream object.

In each case, the angle is defined as elevation from the positive x-axis, in radians.

get_stream(x, y, xvar, yvar, method='rk4', style='mag', maxPoints=40000, extract=False)[source]

Trace a 2D streamline through the domain, returning a Stream object to the caller.

x and y set the starting point for the tracing.

xvar and yvar are string keys to self.data that define the vector field through which this function traces.

The method kwarg sets the numerical method to use for the tracing. Default is Runge-Kutta 4 (rk4).

gradP_regular(cellsize=None, dim1range=-1, dim2range=-1)[source]

Calculate pressure gradient on a regular grid. Note that if the Bats2d object is not on a regular grid, one of two things will happen. If kwarg cellsize is set, the value of cellsize will be used to call self.regrid and the object will be regridded using a cellsize of cellsize. Kwargs dim1range and dim2range can be used in the same way they are used in self.regrid to restrict the regridding to a smaller domain. If cellsize is not set and the object is on an irregular grid, an exception is raised.

The gradient is calculated using numpy.gradient. The output units are force density (N/m^3). Three variables are added to self.data: gradp(dim1), gradp(dim2), gradp. For example, if the object is an equatorial cut, the variables gradpx, gradpy, and gradp would be added representing the gradient in each direction and then the magnitude of the vector.

gyroradius(velocities=('u', 'vth'), m_avg=3.1)[source]

Calculate the ion gyroradius in each cell.

velocities to use in calculating the gyroradius are listed by name in the sequence argument velocities. If more than one variable is given, they are summed in quadrature.

m_avg denotes the average ion mass in AMU.

Result is stored in self[‘gyroradius’]

inertial_length(m_avg=3.1)[source]

Calculate the ion inertial length.

m_avg denotes the average ion mass in AMU.

Result is stored in self[‘inertial_length’].

plasma_freq(m_avg=3.1)[source]

Calculate the ion plasma frequency.

m_avg denotes the average ion mass in AMU.

Result is stored in self[‘plasm_freq’].

regrid(cellsize=1.0, dim1range=-1, dim2range=-1, debug=False)[source]

Re-bin data to regular grid of spacing cellsize. Action is performed on all data entries in the bats2d object.

switch_frame(*args, **kwargs)[source]

For files that have more than one data frame (i.e., *.outs files), load data from the iframe-th frame into the object replacing what is currently loaded.

vth(m_avg=3.1)[source]

Calculate the thermal velocity. m_avg denotes the average ion mass in AMU.

Result is stored in self[‘vth’].

find_block
qtree
attrs: collections.abc.Mapping