Field Terms#
Effective-field contributions are subclasses of FieldTerm. Calling
term.register(state) adds two dynamic attributes to the state:
state.h_<name> for the effective field and state.E_<name> for the
energy. Multiple terms can be combined through TotalField.
See also
Discretization — nodal finite-difference and FIC schemes that underlie every field term.
Form Compiler — how
e_expris turned into vectorized backend code.Dynamic Attributes, resolve and remap — what registering a field term actually does to the
State.
- class neuralmag.FieldTerm(n_gauss=None)#
Base class of all effective field contributions.
In simple cases, a subclass is just required to implement the energy functional of a field contribution and a default_name which is used when registering the field term with a state. The form compiler of NeuralMag is then used to generate efficient code for the computation of both the energy and effective field by means of a just-in-time compiler.
- Parameters:
n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Examples
>>> import neuralmag as nm >>> from neuralmag.common import engine as en >>> >>> # Example subclass implementing a uniaxial anisotropy >>> class UniaxialAnisotropyField(nm.FieldTerm): ... default_name = "uaniso" ... ... @staticmethod ... def e_expr(m, dim, _options): ... K = en.Variable("material__Ku", "c" * dim) ... axis = en.Variable("material__Ku_axis", "c" * dim, (3,)) ... return -K * m.dot(axis) ** 2 * en.dV(dim) >>> >>> # Use instance of class to register dynamic attributes in state >>> state = nm.State(nm.Mesh((10, 10, 10), (1e-9, 1e-9, 1e-9))) >>> UniaxialAnisotropyField().register(state) >>> >>> # compute field and energy >>> h = state.h_uaniso >>> E = state.E_uaniso
- attr_name(attr, name=None)#
Return the attribute name for a given attribute.
Falls back to this instance’s
default_name. It must be an instance method, not a classmethod: subclasses such asLIFieldsetdefault_nameper instance (one name per Lifshitz invariant), and a classmethod readingcls.default_namewould ignore that and make sibling instances collide on one attribute.- Parameters:
attr (str) – The attribute name, e.g. “E” or “h”.
name (str, optional) – Overrides the registration name; defaults to
default_name.
- register(state, name=None)#
Register dynamic attributes for the effective field and energy.
Registers dynamic attributes for the computation of the effective field and energy with the given
Stateobject. By naming convention, these methods are registered asstate.h_{name}andstate.E_{name}orstate.handstate.Ein case that ofnamebeing an empty string.- Parameters:
state (
State) – The state.name (str, optional) – The name used for the registration, falls back to
default_nameattribute of the class.
- class neuralmag.TotalField(*field_names)#
Combine multiple field terms into a single total field term.
This class adds up the effective fields and energies of the given field contributions.
- Parameters:
*field_names – The names of the effective field contributions.
Examples
>>> state = nm.State(nm.Mesh((10, 10, 10), (1e-9, 1e-9, 1e-9))) >>> nm.ExchangeField().register(state, "exchange") >>> nm.DemagField().register(state, "demag") >>> nm.ExternalField(h_ext).register(state, "external") >>> nm.TotalField("exchange", "demag", "external").register(state) >>> # Compute total field and energy >>> h = state.h >>> E = state.E
- register(state, name=None)#
Register dynamic attributes for the effective field and energy.
Registers dynamic attributes for the computation of the effective field and energy with the given
Stateobject. By naming convention, these methods are registered asstate.h_{name}andstate.E_{name}orstate.handstate.Ein case that ofnamebeing an empty string.- Parameters:
state (
State) – The state.name (str, optional) – The name used for the registration, falls back to
default_nameattribute of the class.
Predefined Field Terms#
- class neuralmag.BulkDMIField(**kwargs)#
Effective field contribution for the micromagnetic bulk-DMI energy.
\[E = \int_\Omega D \vec{m} \cdot (\nabla \times \vec{m}) \dx\]with the DMI constant \(D\) given in units of \(\text{J/m}^2\).
- Parameters:
n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.Db– the DMI constant, in J/m^2 (cell scalar field)state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
- class neuralmag.DemagField(p=20, cache_dir=None, gpu_setup=None, **kwargs)#
Effective field contribution corresponding to the demagnetization field.
The demagnetization field (also referred to as magnetostatic field or stray field) is computed from the scalar potential \(u\) as
\[\vec{H}_\text{demag} = - \nabla u\]with \(u\) being calculated by the Poisson equation
\[\Delta u = \nabla \cdot (M_s \vec{m})\]with open boundary conditions by default. Periodic boundary conditions can be enabled on the
Meshviapbc: true PBC (requires a 3D mesh with all three directions set toTrue/inf) switches the convolution to the periodic kernel of Bruckner et al., Sci. Rep. 11, 9202 (2021); pseudo-PBC reuses the open-boundary kernel with image copies. Partial true PBC (e.g.pbc=(True, True, 0)) is not supported and raises aValueError— use pseudo-PBC (positive integer) for the periodic directions instead. See the PBC user guide for details.- Parameters:
p (int) – Distance threshold at which the demag tensor is approximated by a dipole field given in numbers of cells. Defaults to 20.
cache_dir (str or Path or None) – Directory for caching the demag tensor. If provided and a matching tensor exists (same n, dx, pbc, p), it is loaded instead of recomputed. The directory is created automatically if it does not exist. Caching is disabled by default (
None).gpu_setup (bool or None) – jax backend only: build the demag tensor on the GPU via PyTorch instead of the default NumPy CPU build (much faster for large meshes).
None(default) defers to theNM_JAX_GPU_SETUPenv var;True/Falseforce it. Requires PyTorch and a CUDA device; ignored by the torch backend.n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
- class neuralmag.CubicAnisotropyField(n_gauss=None)#
Effective field contribution corresponding to the cubic anisotropy energy.
\[E = - \int_\Omega K \big( (\vec{m} \cdot \vec{e}_1)^2(\vec{m} \cdot \vec{e}_2)^2 + (\vec{m} \cdot \vec{e}_2)^2(\vec{m} \cdot \vec{e}_3)^2 + (\vec{m} \cdot \vec{e}_1)^2(\vec{m} \cdot \vec{e}_3)^2 \big) \dx\]with the anisotropy constant \(K\) given in units of \(\text{J/m}^3\).
state.material.Kc_axis1,state.material.Kc_axis2, andstate.material.Kc_axis3should be orthogonal unit vectors.- Parameters:
n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.Kc– the anisotropy constant, in J/m^3 (cell scalar field)state.material.Kc_axis1– the first anisotropy axis as unit vector field (cell vector field)state.material.Kc_axis2– the second anisotropy axis as unit vector field (cell vector field)state.material.Kc_axis3– the third anisotropy axis as unit vector field (cell vector field)state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
- class neuralmag.ExchangeField(**kwargs)#
Effective field contribution corresponding to the micromagnetic exchange energy.
\[E = \int_\Omega A \big( \nabla m_x^2 + \nabla m_y^2 + \nabla m_z^2 \big) \dx\]with the exchange constant \(A\) given in units of \(\text{J/m}\).
- Parameters:
n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.A– the exchange constant, in J/m (cell scalar field)state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
- class neuralmag.ExternalField(h, **kwargs)#
Effective field contribution corresponding to the external field.
\[E = - \int_\Omega \mu_0 M_s \vec{m} \cdot \vec{h} \dx\]- Parameters:
h (config.backend.Tensor) – The field either given as a
config.backend.Tensoror callable in case of a field that depends e.g. onstate.t. The shape must be either (nx, ny, nz, 3) or (3,) in which case the field expanded to full size.n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
Examples
>>> state = nm.State(nm.Mesh((10, 10, 10), (1e-9, 1e-9, 1e-9))) >>> # define constant external field from expanded function >>> h_ext = nm.VectorFunction(state).fill((0, 0, 0), expand=True) >>> external = nm.ExternalField(h_ext) >>> # define external field in y-direction linearly increasing with time >>> external = nm.ExternalField(lambda t: t * state.tensor([0, 8e5 / 10e-9, 0]))
- class neuralmag.InterfaceDMIField(**kwargs)#
Effective field contribution corresponding to the micromagnetic interface-DMI energy.
\[E = \int_\Omega D \Big[ \vec{m} \cdot \nabla (\vec{e}_D \cdot \vec{m}) - (\nabla \cdot \vec{m}) (\vec{e}_D \cdot \vec{m}) \Big] \dx\]with the DMI constant \(D\) given in units of \(\text{J/m}^2\).
- Parameters:
n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.Di– the DMI constant, in J/m^2 (cell scalar field)state.material.Di_axis– the DMI surface normal as unit vector field (cell vector field)state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
- class neuralmag.InterlayerExchangeField(idx1, idx2, **kwargs)#
Effective field contribution for interface exchange up to biquadratic order.
\[E = - 0.5 \int_\Gamma A \vec{m} \cdot \vec{m}_\text{other} \ds\]where \(\Gamma\) denotes two coupled interfaces and \(\vec{m}_\text{other}\) denotes the magnetization at the nearest point of the other interface. This expression is equivalent to the more common expression
\[E = - \int_\Gamma A \vec{m}_1 \cdot \vec{m}_2 \ds\]where \(\Gamma\) denotes a shared interface and \(\vec{m}_1\) and \(\vec{m}_2\) denote the magnetization on the respective sides. Currently, only surfaces in the xy-plane are supported.
- Parameters:
idx1 (int) – z-index of the first interfaces.
idx2 (int) – z-index of the second interfaces.
n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.iA– the interface coupling constant, in J/m^2 (CCN scalar field)state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
- class neuralmag.LIField(LI: str, **kwargs)#
Effective field contribution for a single Lifshitz invariant (LI) term.
These terms can be combined to form the DMI energy for different crystal symmetries. Lifshitz invariants are of the form
\[L_{ij}^k(\vec m) = m_i \partial_k m_j - m_j \partial_k m_i,\]where \(i,j,k \in \{x,y,z\}\), \(i \neq j\), and \(m_i\) is the \(i\)-th Cartesian component of the unit magnetization vector \(\vec m\). This class constructs the energy density by multiplying the Lifshitz invariant by a DMI constant, giving the form
\[E = - \int_\Omega D_{ijk} L_{ij}^k(\vec m) \dx,\]with material constant \(D_{ijk}\) provided as a scalar cell field whose name encodes the chosen LI directions.
The specific invariant is selected by the string
LIpassed to the constructor. This must be a three-character string composed only of the lettersx,yandz. The first two characters specify the pair \((i,j)\) in order whilst the third character specifies \(k\) (the derivative direction).The corresponding material constant must be supplied in the simulation state as
state.material.D{LI}(e.g.state.material.Dxyz).- Parameters:
LI (str) – Three-character string selecting the Lifshitz invariant indices (e.g.
'xyz').n_gauss (int, optional) – Degree of Gauss quadrature used in the form compiler (inherited via kwargs).
- Raises:
TypeError – If
LIis not a string.ValueError – If
LIlength is not 3 or contains characters other thanx,y,z.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.D{LI}– the DMI constant associated with this invariant (cell scalar field)state.material.Ms– the saturation magnetization, in A/m (cell scalar field)
Examples
>>> state = nm.State(nm.Mesh((10, 10, 10), (1e-9, 1e-9, 1e-9))) >>> state.material.Dxyz = 1e-3 >>> nm.LIField("xyz").register(state, "LI_xyz")
- class neuralmag.OerstedField(p=20, cache_dir=None, region=None, gpu_setup=None, **kwargs)#
Effective field contribution corresponding to the Oersted field.
The Oersted field is created by a current density \(\vec{j}\) (state attribute
state.j) via the Biot–Savart law\[\vec{H}_\text{oe}(\vec{x}) = \frac{1}{4 \pi} \int \vec{j}(\vec{x}') \times \frac{\vec{x}-\vec{x}'}{\vert \vec{x}-\vec{x}'\vert^3} \, d\vec{x}'.\]This has the same convolution structure as the demagnetization field and reuses the same FFT infrastructure (
ConvolutionField), but the convolution tensor is the antisymmetric Krüger kernel and the source is the current density (in \(\text{A}/\text{m}^2\), no \(M_s\) and no \(\mu_0\)). The current is interpolated to the cells, convolved, and the resulting field projected back to the magnetization’s function space; since \(\vec{H}_\text{oe}\) is sourced by the (smooth) current and does not jump at the magnet boundary, the back-projection is a plain (unweighted) mass-lumped projection.state.jmay be a staticVectorFunction/VectorCellFunctionor a dynamic attribute (recomputed each step, e.g. solved from a Poisson problem); the kernel is precomputed once.Only open boundary conditions are supported; any periodic boundary raises a
NotImplementedError.- Parameters:
p (int) – Distance threshold (in numbers of cells) at which the Oersted tensor is approximated by a dipole field. Defaults to 20.
cache_dir (str or Path or None) – Directory for caching the Oersted tensor (same semantics as
DemagField). Disabled by default (None).region (str or None) – Optional name of a cell-scalar
CellFunctionattribute that masks the current to a conductor region (1inside,0outside), removing the one-cell leakage of the node→cell projection. This is not the magnetic regionrho(a conductor need not be magnetic). Defaults toNone(no mask).gpu_setup (bool or None) – jax backend only: build the Oersted tensor on the GPU via PyTorch instead of the default NumPy CPU build.
Nonedefers toNM_JAX_GPU_SETUP; requires PyTorch + CUDA; ignored by the torch backend.n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.j– the current density, in A/m² (cell or nodal vector field)state.material.Ms– the saturation magnetization (energy only) (cell scalar field)
- class neuralmag.UniaxialAnisotropyField(**kwargs)#
Effective field contribution corresponding to the quadratic uniaxial anisotropy energy.
\[E = - \int_\Omega K \big( \vec{m} \cdot \vec{e}_k \big)^2 \dx\]with the anisotropy constant \(K\) given in units of \(\text{J/m}^3\). This term covers the quadratic (second-order) uniaxial anisotropy only.
- Parameters:
n_gauss (int) – Degree of Gauss quadrature used in the form compiler.
Required state attributes
Assuming the default names (i.e. not renamed):
state.material.Ku– the anisotropy constant, in J/m^3 (cell scalar field)state.material.Ku_axis– the anisotropy axis as unit vector field (cell vector field)state.material.Ms– the saturation magnetization, in A/m (cell scalar field)