Importing/Exporting field from/to xarray.DataArray#

`xarray <https://docs.xarray.dev/en/stable/>`__ 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 and holoviews 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 rendered plots).

Exporting discretisedfield.Field to xarray.DataArray#

to_xarray method of the discretisedfield.Field object is utilised to export the field data to a DataArray. As an example, we use a three-dimensional field.

[1]:
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, nvdim=3, value=(0, 0, 1))
field
[1]:
Field
  • Mesh
    • Region
      • pmin = [0.0, 0.0, 0.0]
      • pmax = [1e-08, 1e-08, 1e-08]
      • dims = ['x', 'y', 'z']
      • units = ['m', 'm', 'm']
    • n = [2, 2, 2]
  • nvdim = 3
  • vdims:
    • x
    • y
    • z
[2]:
xarray = field.to_xarray()
xarray
[2]:
<xarray.DataArray 'field' (x: 2, y: 2, z: 2, vdims: 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
  * vdims    (vdims) <U1 'x' 'y' 'z'
Attributes:
    units:             None
    cell:              [5.e-09 5.e-09 5.e-09]
    pmin:              [0. 0. 0.]
    pmax:              [1.e-08 1.e-08 1.e-08]
    nvdim:             3
    tolerance_factor:  1e-12

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 vdims represents the components of the field. There will be no vdims dimension if the stored field is a scalar field.

By default, the DataArray has four attributes, namely units, cell, pmin, and pmax. The last three attributes store the information about the mesh, and they prove useful for mesh reconstruction from a ``DataArray` <#from_xarray>`__ 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 unit arguments.

[3]:
xarray_2 = field.to_xarray(name="m", unit="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.

Few useful xarray.DataArray methods#

The underlying array can be accessed through:

[4]:
xarray.values
[4]:
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:

[5]:
xarray.sel(x=2e-9, method="nearest")
[5]:
<xarray.DataArray 'field' (y: 2, z: 2, vdims: 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
  * vdims    (vdims) <U1 'x' 'y' 'z'
Attributes:
    units:             None
    cell:              [5.e-09 5.e-09 5.e-09]
    pmin:              [0. 0. 0.]
    pmax:              [1.e-08 1.e-08 1.e-08]
    nvdim:             3
    tolerance_factor:  1e-12

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:

[6]:
xarray.sel(vdims="z")
[6]:
<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
    vdims    <U1 'z'
Attributes:
    units:             None
    cell:              [5.e-09 5.e-09 5.e-09]
    pmin:              [0. 0. 0.]
    pmax:              [1.e-08 1.e-08 1.e-08]
    nvdim:             3
    tolerance_factor:  1e-12

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.

[7]:
xarray.sel(vdims="z").sel(x=2e-9, method="nearest")
[7]:
<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
    vdims    <U1 'z'
Attributes:
    units:             None
    cell:              [5.e-09 5.e-09 5.e-09]
    pmin:              [0. 0. 0.]
    pmax:              [1.e-08 1.e-08 1.e-08]
    nvdim:             3
    tolerance_factor:  1e-12

One can also plot the above DataArray as:

[8]:
xarray.sel(vdims="z").sel(x=2e-9, method="nearest").plot();
../../../_images/documentation_notebooks_discretisedfield_xarray-usage_14_0.png

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 unit arguments are not str.

## Importing 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 primary example, we convert the xarray.DataArray created by to_xarray method to a new discretisedfield.Field.

[9]:
field_new = df.Field.from_xarray(xarray)
field_new
[9]:
Field
  • Mesh
    • Region
      • pmin = [0.0, 0.0, 0.0]
      • pmax = [1e-08, 1e-08, 1e-08]
      • dims = ['x', 'y', 'z']
      • units = ['m', 'm', 'm']
    • n = [2, 2, 2]
  • nvdim = 3
  • vdims:
    • x
    • y
    • z

One can also first define an xarray.DataArray and then convert it to discretisedfield.Field.

[10]:
import numpy as np
import xarray as xr

xdr = xr.DataArray(
    np.ones((2, 2, 2, 3), dtype=float),
    dims=["x", "y", "z", "vdims"],
    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),
        vdims=["x", "y", "z"],
    ),
    name="mag",
    attrs=dict(
        nvdim=3, cell=[5e-9, 5e-9, 5e-9], pmin=[0.0, 0.0, 0.0], pmax=[1e-8, 1e-8, 1e-8]
    ),
)
xdr
[10]:
<xarray.DataArray 'mag' (x: 2, y: 2, z: 2, vdims: 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
  * vdims    (vdims) <U1 'x' 'y' 'z'
Attributes:
    nvdim:    3
    cell:     [5e-09, 5e-09, 5e-09]
    pmin:     [0.0, 0.0, 0.0]
    pmax:     [1e-08, 1e-08, 1e-08]
[11]:
field_xdr = df.Field.from_xarray(xdr)
field_xdr
[11]:
Field
  • Mesh
    • Region
      • pmin = [0.0, 0.0, 0.0]
      • pmax = [1e-08, 1e-08, 1e-08]
      • dims = ['x', 'y', 'z']
      • units = ['m', 'm', 'm']
    • n = [2, 2, 2]
  • nvdim = 3
  • vdims:
    • x
    • y
    • z

from_xarray exceptions and expected properties of input xarray.DataArray#

The input argument of from_xarray must be an xarray.DataArray. Further, the input xarray.DataArray must have following properties in order to obtain a discretisedfield.Field: 1. The attributes must have an nvdim entry equal to one for scalar fields and greater than one for vector fields. Further the value of nvdim must be an integer. 2. There must be a dimension vdims if nvdim is greater than one. 3. The coordinates of the geometric dimensions must be equally spaced.

Lastly, it is strongly advised that the input ``xarray.DataArray`` should have ``pmin``, ``pmax``, and ``cell`` attributes required for the reconstruction of the mesh. The xarray.DataArray in the above example has all these properties.