3. Creating and defining a sample#

3.1. The sample object#

3.1.1. Initialization#

The central object in a TetraX workflow is the sample, represented by an instance of the AbstractSample class. Much like physical samples, they are characterized by a geometry, a set of material parameters, a magnetization state and the involved magnetic interactions. In numerical experiments, the sample is additionally characterized by an underlying mesh, a scale, a set of differential operators (gradient, divergence, etc.) tied to the mesh and a number of other attributes necessary for computation. We can create a sample using

>>> sample = tx.create_sample("name of my sample")

This will initialize a sample (without a geometry yet) with a set of default material parameters. These parameters can be accessed and set in the usual pythonic way

>>> sample.Msat = 800e3 # A/m
>>> sample.alpha = 0.08
>>> print(sample.alpha)
0.08

>>> print(sample.name)
"name of my sample"

The name of the sample is set purely for the purpose of storing data in a meaningful way, as all data obtained in numerical experiments is saved in a directory named after the sample. Note that, as of version 1.3.0, several basic material parameters such as the saturation magnetization sample.Msat or the exchange-stiffness constant sample.Aex can be inhomogeneous within the sample (more details below).

3.1.2. Setting a geometry#

Before proceeding with anything further, a specific geometry for the sample needs to be set. This is done by supplying a mesh (in meshio format) to the sample using

>>> sample.set_geom(mesh)

In the tetrax.geometries submodule, we supply a couple of templates to generate common meshes, such as cuboids, prisms, tubes, wires, and so forth. However, desired meshes can be generated also using pygmsh or even any external software, as long as the generated mesh files are readable by meshio. Mesh files can be read in using

>>> sample.read_mesh(filename)

When setting a mesh/geometry, under the hood, certain preprocessing necessary for later computations is performed. For example, the differential operators on the mesh or the normal vector field of the sample are calculated. After these calculations have been performed, important numerical parameters of the mesh can be accessed, including

  • mesh, the mesh itself in meshio format,

  • nx, the total number of nodes in the mesh,

  • nb, the number of boundary nodes in the mesh, and

  • xyz, the coordinates at each node as a MeshVector. These coordinates are expressed in units of

  • scale, the characteristic length of the sample. The default is 1e-9 which means that the coordinates on the mesh are given in nanometers.

Note

The mesh of the sample can visualized using the show method on the sample, to allow for checking if the mesh was created/obtained according to the request or mesh in mind:

>>> sample.show()

3.1.3. Magnetization vector field#

Once a mesh has been assigned to the sample, the magnetization vector field of the sample (mag) can be set for example as

>>> sample.mag = [0, 0, 1]

Note

A magnetization vector field can only be set for samples which already have a mesh.

In antiferromagnetic cases, the sample has a magnetization vector field for each sublattice, mag1 and mag2. To initialize an inhomogeneous magnetization, an array of shape (nx, 3) needs to supplied which can be seen as a list of triplets with each triplet representing the vector at a specific node. For convenience, a couple of template vector fields (such as domain walls, vortices, etc.) can be found in the tetrax.vectorfields submodule. These fields can be called by supplying xyz, the coordinates on the mesh.

>>> sample.mag = tx.vectorfields.bloch_wall(sample.xyz, ...)

When setting a magnetization field, the input field is automatically normalized and internally converted to a MeshVector, which is a data structure that allows to easily access the different components of the magnetization.

>>> sample.mag.x
MeshScalar([0.0, 0.0, 0.12, 0.15, ...])

Here, we see that the individual components of a mesh vector are MeshScalars, another basic data structure in TetraX. We can easily calculate the averages of these objects using the sample_average function from the tetrax.helpers.math module, which takes some vector- or scalar field and the sample itself as an input.

>>> m = tx.sample_average(sample.mag, sample)
>>> m
[0.98, 0.02, 0.0]

>>> mx = tx.sample_average(sample.mag.x, sample)
>>> mx
0.98

Note

The reason, why sample_average is not called volume_average is because not all samples have to be three-dimensional. For example, in waveguide samples, sample_average calculates the average in a cross section.

The mesh and the current magnetization of a sample (if available) can always be visualized using the show method.

>>> sample.show()

For further details, see Visualization and evaluation.

3.2. Different types of samples#

Depending on the specific type of geometry and magnetic order, concrete subclasses of AbstractSample can be generated which are distinguished by their set of possible material parameters (ferromagnetic or antiferromagnetic), dimensionality, or numerical solvers needed to calculate magnetization dynamics or statics. The type of sample is chosen when creating the sample. At the beginning of this section, a sample was created by calling tetrax.create_sample(). Per default, this function creates a ferromagnetic waveguide sample. However, different options are possible.

tetrax.core.sample.create_sample(name='my_sample', geometry='waveguide', magnetic_order='FM')#

Create a sample with certain geometry and magnetic order.

Parameters
geometrystr, {“confined”, “waveguide”, “layer”}

Geometry of the sample (default is “waveguide”).

magnetic_orderstr {“FM”, “AFM”}, optional

Magnetic order of the material. Can be either ferromagnetic (“FM”) or antiferromagnetic (“AFM”). The default is “FM”.

Returns
sampleAbstractSample

We can see from the geometry parameter, that different types of samples can be selected:

  • confined samples such as disks, cuboids, spheres or Möbius ribbons. These types of samples are often the subject of standard three-dimensional micromagnetic simulations. Only in such samples, an LLG solver can be used. Using a confined-wave dynamic-matrix approach, the normal modes of confined samples can be calculated.

  • waveguide samples are infinitely long in one direction (here, the \(z\) direction). Moreover, their (equilibrium) magnetization is taken to be translationally invariant in this direction. Such samples are modelled by discretizing only a single cross section of the waveguide. A propagating-wave dynamic-matrix approach can be used to calculate the normal modes as a function of the wave vector \(k\) along the \(z\) direction.

  • layer samples consist of one or multiple infinitely extended layers which are modelled by discretizing a line trace along the thickness of the layer(s). The (equilibrium) magnetization is taken to be homogeneous within the layer planes. A propagating-wave dynamic-matrix approach can be used to calculate the normal modes as a function of the wave vector \(k\) along the \(z\) direction.

For all three different types of samples, energy minimizers are implemented to calculated the magnetic ground states.

../_images/sample_types.png

Note

As of version 1.3.0, only ferromagnetic waveguides and multilayers are fully supported, while antiferromagnets and confined samples (without dipolar interaction) are included as experimental features.

3.3. Specific geometries#

Template geometries (mesh generators) are all available in the tetrax.geometries module. Even though they are implemented in submodules, all of them can be accessed in the namespace of tetrax.geometries. For examples

sample.set_geom(tx.geometries.tube_cross_section(r=20, R=30, lc=3))

Note

For a given sample, only such meshes can be set which fit the type of the sample selected. For example, it is not possible to set a three-dimensional cuboid mesh for a sample which represents the cross section of an infinitely long waveguide.

3.3.1. Geometries for confined samples#

The available geometries for three-dimensional confined samples are listed here.

Warning

Three-dimensional confined samples are not fully supported yet. While antiferromagnets are, in principle, fully supported accross all types of geometries, their magnetic tensors are not properly tested yet. Ferromagnetic confined samples do not include dipolar interactions.

tetrax.geometries.geometries3D.disk(R, t, lc=5, nl=2, x0=0, y0=0)#

Creates a mesh for a round disk.

Parameters
Rfloat

Radius of the round wire.

tfloat

thickness of the disk.

lcfloat

Characteristic length of mesh in the (xy) plane (default is 5).

nlinteger

The number of layers along the thickness (default is 3).

x0float

x coordinate of the origin (Default is 0).

y0float

y coordinate of the origin (Default is 0).

Returns
meshmeshio.Mesh
tetrax.geometries.geometries3D.prism(a, b, c, lca=5, lcb=5, lcc=5)#

Creates a regular mesh of a prism.

usage/prism.png
Parameters
afloat

Length of the prism.

bfloat

Width of the prism.

cfloat

Thickness of the prism.

lcafloat

characteristic mesh length along the length (Default is 5). The number of FEM layers along the length will be a/lca + 1.

lcbfloat

characteristic mesh length along the width (Default is 5). The number of FEM layers along the width will be b/lcb + 1.

lccfloat

characteristic mesh length along the thickness (Default is 5). The number of FEM layers along the thickness will be c/lcc + 1.

Returns
——-
meshmeshio.Mesh

See also

round_wire_cross_section, polygonal_cross_section
tetrax.geometries.geometries3D.tube(r, R, t, lc=5, nl=2, x0=0, y0=0)#

Creates a mesh with the cross section of a tube.

../_images/tube_cross_section.png
Parameters
rfloat

Inner radius of the tube.

Rfloat

Outer radius of the tube.

lcfloat

Characteristic length of mesh (default is 5).

x0float

x coordinate of the origin (Default is 0).

y0float

y coordinate of the origin (Default is 0).

Returns
meshmeshio.Mesh

See also

tube_segment_cross_section, round_wire_cross_section, polygonal_cross_section

3.3.2. Geometries for waveguide samples#

tetrax.geometries.geometries2D.polygonal_cross_section(points, lc=5)#

Creates a mesh with the cross section of a polygonal wire.

../_images/polygonal_cross_section.png
Parameters
pointslist (Point)

List of points defining the corners of the polygon.

lcfloat

Characteristic length of mesh (default is 5).

Returns
meshmeshio.Mesh
tetrax.geometries.geometries2D.rectangle_cross_section(a, b, lca=5, lcb=5)#

Creates a regular mesh with rectangular cross section.

../_images/rectangle_cross_section.png
Parameters
afloat

Width of the rectangle.

bfloat

Thickness of the rectangle.

lcafloat

characteristic mesh length along the width (Default is 5). The number of FEM layers along the width will be a/lca + 1.

lcbfloat

characteristic mesh length along the thickness (Default is 5). The number of FEM layers along the thickness will be b/lcb + 1.

Returns
meshmeshio.Mesh
tetrax.geometries.geometries2D.round_wire_cross_section(R, lc=5, x0=0, y0=0)#

Creates a mesh with the cross section of a round wire.

../_images/round_wire_cross_section.png
Parameters
Rfloat

Radius of the round wire.

lcfloat

Characteristic length of mesh (default is 5).

x0float

x coordinate of the origin (Default is 0).

y0float

y coordinate of the origin (Default is 0).

Returns
meshmeshio.Mesh
tetrax.geometries.geometries2D.round_wire_cross_section_refined(R, func, x0=0, y0=0)#

Creates a refined mesh with the cross section of a round wire.

Instead of a fixed characteristic mesh length lc, a scalar field func(x,y,z) is supplied to define lc at each position of the mesh. This mesh is preferable, for example, for round wires housing a vortex string.

../_images/round_wire_cross_section.png
Parameters
Rfloat

Radius of the round wire.

funcfunc

Scalar field of signture func(x,y,z) which returns the characteristic mesh length at each position.

x0float

x coordinate of the origin (Default is 0).

y0float

y coordinate of the origin (Default is 0).

Returns
meshmeshio.Mesh
tetrax.geometries.geometries2D.tube_cross_section(r, R, lc=5, x0=0, y0=0)#

Creates a mesh with the cross section of a tube.

../_images/tube_cross_section.png
Parameters
rfloat

Inner radius of the tube.

Rfloat

Outer radius of the tube.

lcfloat

Characteristic length of mesh (default is 5).

x0float

x coordinate of the origin (Default is 0).

y0float

y coordinate of the origin (Default is 0).

Returns
meshmeshio.Mesh
tetrax.geometries.geometries2D.tube_segment_cross_section(arc_length, thickness, angle, lc=5, shifted_to_origin=True)#

Creates a mesh with the cross section of a tube segment.

../_images/tube_segment_cross_section.png
Parameters
arc_lengthfloat

Average arc length of the tube segment (at the center surface).

thicknessfloat

Thickness of the tube segment.

anglefloat

Opening angle of the tube segment (degrees).

lcfloat

Characteristic length of mesh (default is 5).

shifted_to_originbool

If True, shift the tube segment to the coordinate origin (default is True). Otherwise, the origin will be the symmetry axis of the tube.

Returns
meshmeshio.Mesh

3.3.3. Geometries for layer samples#

The available geometries for one-dimensional layer samples are listed here.

tetrax.geometries.geometries1D.bilayer_line_trace(thickness1, thickness2, spacing, lcs, y0=0)#

Creates a mesh for the line trace along the normal direction of a magnetic bilayer.

../_images/bilayer_line_trace.png
Parameters
thickness1float

Thickness \(d_1\) of the bottom magnetic layer.

thickness2float

Thickness \(d_2\) of the top magnetic layer.

spacingfloat

Thickness \(s\) of the interlayer.

lcsfloat or list(float)

Characteristic lengths for the different magnetic layers. If list, has to contain two elements. If not a list it will be same for both layers.

y0float

y coordinate of the center (Default is 0).

tetrax.geometries.geometries1D.monolayer_line_trace(thickness, lc=2, y0=0)#

Creates a mesh for the line trace along the normal direction of a magnetic monolayer.

../_images/monolayer_line_trace.png
Parameters
thicknessfloat

Thickness \(d\) of the monolayer.

lcfloat

Characteristic mesh length.

y0float

y coordinate of the center (Default is 0).

tetrax.geometries.geometries1D.multilayer_line_trace(thicknesses, spacings, lcs, y0=0)#

Creates a mesh for the line trace along the normal direction of a magnetic multi-layer.

../_images/multilayer_line_trace.png
Parameters
thicknesseslist(float)

Thicknesses \([d_1, d_2, d_3, ...]\) of the magnetic layers (ordered with increasing \(y\) coordinate). List must be one element longer than spacings.

spacingslist(float)

Spacings \([s_1, s_2, s_3, ...]\) between the magnetic layers, which is are the thicknesses of the non-magnetic interlayers (ordered with increasing \(y\) coordinate). List must be one element shorter than thicknesses.

lcslist(float) or float

Characteristic lengths for the different magnetic layers. If list, has to be same length as thicknesses. If not a list it will be same for all layers.

y0float

y coordinate of the center (Default is 0).

Returns
meshmeshio.Mesh

3.4. Material parameters#

Once a sample object is created, being ferro- or antiferromagnetic, it will be initialized with the material parameters of Ni80Fe20(permalloy). However, the material parameters can be set individually to match the system in study.

3.4.1. Ferromagnetic samples#

Typing the following command one can easily check the current material parameters as well as some of the underlaying mesh parameters.

>>> sample

Sample: Material parameters

Material parameters

Attribute name

(Average) value

Saturation magnetization

Msat

796000.0 A/m

Exch. stiffness constant

Aex

1.3e-11 J/m

Gyromagnetic ratio

gamma

176085964400.0 rad Hz/T

Gilbert damping parameter

alpha

0.01

Uniaxial anisotropy constant

Ku1

0.0 J/m3

Uniaxial anisotropy direction

e_u

array([0,0,1])

Cubic anisotropy constant

Kc1

0.0 J/m3

Cubic anisotropy axes

v1Kc, v2Kc

(1,0,0), (0,1,0)

Bulk DMI constant

Dbulk

0.0 J/m2

Interfacial DMI constants

Didmi

0.0 J/m2

Interlayer exchange constant

J1

0.0 J/m2

Sample: Mesh parameters

Mesh parameters

Attribute name

Value

Mesh scaling

scale

1e-9

Number of nodes

nx

283

Number of boundary nodes

nb

83

Any of the listed parameters can easily be changed using:

sample.Msat = 796e3
sample.Aex = 13e-12
sample.alpha = 0.007
sample.scale = 10e-9

The uniaxial anisotropy can be defined to be homogeneous over the whole sample or to be spatially varying, both in magnitude and direction. The direction doesn’t need to be a unit vector and will not be normalized in the code. Therefore, the implementation allows to define directions in units of the anisotropy and set the anisotropy constant to 1. For details please check the next sections.

3.4.1.1. Spatial-dependent saturation and exchange stiffness#

Eventhough the attributes sample.Msat and sample.Aex can be set by specifying a single number, they are internally converted to MeshScalars of length sample.nx. One can check this by

>>> sample.Msat = 796e3
>>> sample.Msat
MeshScalar([796000., 796000., 796000., 796000., 796000., 796000., 796000.,
         .....
         796000., 796000., 796000., 796000., 796000., 796000.])

or, similarly,

>>> sample.Aex
MeshScalar([1.3e-11, 1.3e-11, 1.3e-11, 1.3e-11, 1.3e-11, 1.3e-11, 1.3e-11,
            ....
            1.3e-11, 1.3e-11, 1.3e-11, 1.3e-11, 1.3e-11, 1.3e-11])

The average of these scalar fields can always be obtained using

>>> sample.average(sample.Aex)
MeshScalar(1.3e-11)

>>> sample.average(sample.Msat)
MeshScalar(796000.)

For the saturation, we also supply the shorthand sample.Msat_avrg. Of course, by supplying the saturation or the exchange stiffness as an array, list or MeshScalars one can specify inhomogeneous material parameters and, thereby, different materials. For example, we can create a tube with a saturation magnetization that is linearly varying along the \(y\) direction.

Mavrg = 1200e3
dMdy = 20e3
sample.mag = [0, 0, 1]
sample.Msat = Mavrg + dMdy * sample.xyz.y
sample.show(show_extrusion=False)

Notice that the arrows in the 3D plot are now scaled according to their local value of sample.Msat (can be deactivated using show_scaled_mag=False). Importantly, the attribute sample.mag is still a unit vector (\(\mathbf{m}\)) and denotes the local direction of magnetization!

>>> sample.mag
MeshVector([[0., 0., 1.],
            [0., 0., 1.],
            [0., 0., 1.],
            ...

As result, when initializing a magnetization state, we supply only the unit vector (field) for the direction, as above. The full magnetization \(\mathbf{M}=M_\mathrm{s}\mathbf{m}\) (in units of A/m) can be retrieved using

>>> sample.mag_full
MeshVector([[0., 0., 0.12],
            [0., 0., 0.12],
            [0., 0., 0.125],
            ...

This quantity cannot be set directly. For examples including inhomogeneous exchange and/or saturation magnetization, see ../examples/Py25_CoFeB25_double_layer_Grassi and Magnetization-graded waveguides.

3.4.1.2. Homogeneous uniaxal anisotropy#

A homogeneous uniaxial anisotropy can be set using a triplet for the direction and a single constant. For a single crystal ferromagnetic sample would be like:

sample.Ku = 48e3 # in J/m^3
sample.e_u = [1, 0, 0] # along the x axis

Easy plane anisotropies can also be defined, by setting the anisotropy constant to negative value. An example for an xy easy plane would be:

sample.Ku = -20e3 # in J/m^3
sample.e_u = [0, 0, 1] # along the z axis

3.4.1.3. Spatial-dependent uniaxal anisotropy#

A spatial dependent uniaxial anisotropy can be defined by a list of triplets defining the orientation at each node position and by a list of the constants. The template vectorfields in tetrax.vectorfields can be used to define spatial dependent vectorfields.

Here is an example for a radial uniaxial anisotropy with a homogeneous anisotropy constant:

sample.Ku1 = 2e4 #J/m^3
sample.e_u = tx.vectorfields.radial(sample.xyz,1)

The defined vectorfield of the uniaxial anisotropy can be visualized to check if the aimed spatial dependence has been set or not. This can be done by the plot() method. For further details please check the provided example: Spatially dependent uniaxial anisotropy.

3.4.1.4. Cubic anisotropy#

A homogeneous cubic anisotropy can be defined by an anisotropy constant Kc1 and two triplets defining the orientation of the easy axes (the third axis will be automatically calculated as the cross product of the other two). The higher order term is not yet included.

Here is an example for a cubic anisotropy with axes along the (100) directions, typical for iron (Fe):

sample.Kc1 = 48e3 #J/m^3
sample.v1Kc = (1,0,0)
sample.v2Kc = (0,1,0)

For further details please check the provided example: Cubic anisotropy.

3.4.1.5. Interlayer exchange#

Bilinear interlayer-exchange interaction can be used, implemented for layered systems only. The bilinear term can be set using the J1 constant (in J/m^). The constant can also have different values between different layers in a multilayered sample, by supplying a list of values.

Example for bilayers:

sample = tx.create_sample(geometry="layer", name="bilayer_iec")
sample.J1 = -3e-4
mesh = tx.geometries.bilayer_line_trace(20,20,2,1)

Example for multilayers consisting of 3 magnetic layers and 2 non-magnetic spacers: - antiferromagnetic coupling between the first two layers and ferromagnetic coupling between the second and the third layer

sample = tx.create_sample(geometry="layer", name="trilayer_iec")
sample.J1 = [-3e-4,1e-3]
mesh = tx.geometries.multilayer_line_trace([50,20,10],[2,1],[5,5,5])

For further details please check the provided example: ../examples/AFM_bilayers.

3.4.2. Antiferromagnetic samples#

Warning

As of version 1.3.0 antiferromagnetic samples are implemented, but not fully tested. Therefore, we do not supply any documentation yet. Use on your own risk (and tell us your experiences).