xarray.DataArray
¶xarray
provides a convenient method to handle labeled multi-dimensional arrays. It integrates well with numpy
and pandas
for fast array manipulation, as well as matplotlib
for visualization. This makes it a good candidate to carry discretisedfield.Field
data and additional metadata such as field name and unit (which appear in the plots rendered by matplotlib
).
discretisedfield.Field
to xarray.DataArray
¶to_xarray
method of the discretisedfield.Field
object is utilised to export the field data to a DataArray.
import discretisedfield as df
p1 = (0, 0, 0)
p2 = (10e-9, 10e-9, 10e-9)
cell = (5e-9, 5e-9, 5e-9)
mesh = df.Mesh(p1=p1, p2=p2, cell=cell)
field = df.Field(mesh=mesh, dim=3, value=(0, 0, 1))
field
xarray = field.to_xarray()
xarray
<xarray.DataArray 'field' (x: 2, y: 2, z: 2, comp: 3)> array([[[[0., 0., 1.], [0., 0., 1.]], [[0., 0., 1.], [0., 0., 1.]]], [[[0., 0., 1.], [0., 0., 1.]], [[0., 0., 1.], [0., 0., 1.]]]]) Coordinates: * x (x) float64 2.5e-09 7.5e-09 * y (y) float64 2.5e-09 7.5e-09 * z (z) float64 2.5e-09 7.5e-09 * comp (comp) <U1 'x' 'y' 'z' Attributes: units: None cell: (5e-09, 5e-09, 5e-09) p1: (0, 0, 0) p2: (1e-08, 1e-08, 1e-08)
The DataArray
has four dimensions. Dimensions x
, y
, and z
represent the geometric axes and their respective coordinates store the midpoints of cell edges along respective axes, while comp
represents the components of the field. There will be no comp
dimension if the stored field is a scalar field.
By default, the DataArray
has four attributes, namely units
, cell
, p1
, and p2
. The last three attributes store the information about the mesh, and they prove useful for mesh reconstruction from a DataArray
object. The units
attribute stores the unit of the field.
By default, the DataArray name is 'field'
and units
attribute is None
. This can be set by passing the to_xarray
method with name
and units
arguments.
xarray_2 = field.to_xarray(name="m", units="A/m")
print(f"Name of the DataArray is {xarray_2.name}.")
print(f"Units of the DataArray is {xarray_2.attrs['units']}.")
Name of the DataArray is m. Units of the DataArray is A/m.
xarray.DataArray
methods¶The underlying array can be accessed through:
xarray.values
array([[[[0., 0., 1.], [0., 0., 1.]], [[0., 0., 1.], [0., 0., 1.]]], [[[0., 0., 1.], [0., 0., 1.]], [[0., 0., 1.], [0., 0., 1.]]]])
It is possible to select the cross-section of the magnetisation along a particular axis as:
xarray.sel(x=2e-9, method="nearest")
<xarray.DataArray 'field' (y: 2, z: 2, comp: 3)> array([[[0., 0., 1.], [0., 0., 1.]], [[0., 0., 1.], [0., 0., 1.]]]) Coordinates: x float64 2.5e-09 * y (y) float64 2.5e-09 7.5e-09 * z (z) float64 2.5e-09 7.5e-09 * comp (comp) <U1 'x' 'y' 'z' Attributes: units: None cell: (5e-09, 5e-09, 5e-09) p1: (0, 0, 0) p2: (1e-08, 1e-08, 1e-08)
The DataArray.sel
method returns a new DataArray
with only one coordinate along x axis and hence a cross-section. The method="nearest"
argument matches the value along x axis with the nearest midpoint.
Similarly, one can also obtain a desired component of the field as:
xarray.sel(comp="z")
<xarray.DataArray 'field' (x: 2, y: 2, z: 2)> array([[[1., 1.], [1., 1.]], [[1., 1.], [1., 1.]]]) Coordinates: * x (x) float64 2.5e-09 7.5e-09 * y (y) float64 2.5e-09 7.5e-09 * z (z) float64 2.5e-09 7.5e-09 comp <U1 'z' Attributes: units: None cell: (5e-09, 5e-09, 5e-09) p1: (0, 0, 0) p2: (1e-08, 1e-08, 1e-08)
One can take this a step forward and obtain a new DataArray
representing z component of the field on an x-axis cross-section of the geometry.
xarray.sel(comp="z").sel(x=2e-9, method="nearest")
<xarray.DataArray 'field' (y: 2, z: 2)> array([[1., 1.], [1., 1.]]) Coordinates: x float64 2.5e-09 * y (y) float64 2.5e-09 7.5e-09 * z (z) float64 2.5e-09 7.5e-09 comp <U1 'z' Attributes: units: None cell: (5e-09, 5e-09, 5e-09) p1: (0, 0, 0) p2: (1e-08, 1e-08, 1e-08)
One can also plot the above DataArray
as:
xarray.sel(comp="z").sel(x=2e-9, method="nearest").plot();
Notice that proper units and labels on the plot are displayed by default. The full list of xarray.DataArray
methods can be found in the API reference.
to_xarray
exceptions¶The to_xarray
method raises TypeError
if either name
or units
arguments are not str
.
discretisedfield.Field
from xarray.DataArray
¶It is possible to create a discretisedfield.Field
from an xarray.DataArray
with the help of class method from_xarray
. As a first step, we convert the xarray.DataArray
created by to_xarray
method to a new discretisedfield.Field
.
field_new = df.Field.from_xarray(xarray)
field_new
One can also first define an xarray.DataArray
and then convert it to discretisedfield.Field
.
import numpy as np
import xarray as xr
xdr = xr.DataArray(
np.ones((2, 2, 2, 3), dtype=float),
dims=["x", "y", "z", "comp"],
coords=dict(
x=np.arange(2.5e-9, 1e-8, 5e-9),
y=np.arange(2.5e-9, 1e-8, 5e-9),
z=np.arange(2.5e-9, 1e-8, 5e-9),
comp=["x", "y", "z"],
),
name="mag",
attrs=dict(cell=[5e-9, 5e-9, 5e-9], p1=[0.0, 0.0, 0.0], p2=[1e-8, 1e-8, 1e-8]),
)
xdr
<xarray.DataArray 'mag' (x: 2, y: 2, z: 2, comp: 3)> array([[[[1., 1., 1.], [1., 1., 1.]], [[1., 1., 1.], [1., 1., 1.]]], [[[1., 1., 1.], [1., 1., 1.]], [[1., 1., 1.], [1., 1., 1.]]]]) Coordinates: * x (x) float64 2.5e-09 7.5e-09 * y (y) float64 2.5e-09 7.5e-09 * z (z) float64 2.5e-09 7.5e-09 * comp (comp) <U1 'x' 'y' 'z' Attributes: cell: [5e-09, 5e-09, 5e-09] p1: [0.0, 0.0, 0.0] p2: [1e-08, 1e-08, 1e-08]
field_xdr = df.Field.from_xarray(xdr)
field_xdr
from_xarray
exceptions and expected properties of input xarray.DataArray
¶The input argument of from_xarray
must be an xarray.DataArray
, otherwise a TypeError
is raised. Further, the input xarray.DataArray
must have following properties in order to obtain a discretisedfield.Field
:
ValueError
is raised.x
, y
, z
, and comp
(if a vector field). If not, a ValueError
is raised.x
, y
, and z
dimensions must be equally spaced. If not, a ValueError
is raised.Lastly, it is strongly advised that the input xarray.DataArray
should have p1
, p2
, and cell
attributes required for the reconstruction of the mesh. The xarray.DataArray
in the above example has all these properties.