Field Terms
The effective field can be comprised of many different fields such as the quantum mechanical exchange interaction or an anisotropy field due to spin orbit coupling. FieldTerm
is a class representing the possible field terms. magnum.np comes with a number of predefined field term classes that cover a variety of use cases. In most applications a list of FieldTerm
objects is usually passed to the solver in order to solve micromagnetic problems under the influence of different effects. Furthermore, each FieldTerm
object has methods to compute the corresponding effective field and energy for a given State object.
Here is an examples on how to use the field term classes:
# initialize state with uniform magnetization configuration
state = State(mesh, m = Constant((1.0, 0.0, 0.0)))
# initialize demagnetization field object
demag = DemagField()
# compute field and save to file
File("h.pvd") << demag.h(state)
# compute energy and print on screen
print demag.E(state)
Both, the field method h
and the energy method E
can also be used as virtual attribues:
state.h_demag = demag.h()
state.E_demag = demag.E()
The computation of the field is carried out in the default region which usually is the ‘magnetic’ region. However, if the user wishes to restrict the computation of the field to a certain region, this default can be overriden in the initialization, like so:
exchange = ExchangeField(region = 'ferromagnetic')
where the region called ‘ferromagnetic’ would have to be defined by the user beforehand. How this can be done is explained under Domains.
Some field term classes require certain material parameters to be defined. The following documentatin offers more detailed information about each class.
LinearFieldTerm
UniaxialAnisotropyField
- class magnumnp.UniaxialAnisotropyField(**kwargs)[source]
Bases:
LinearFieldTerm
Uniaxial Anisotropy Field:
\[\vec{h}^\text{u} = \frac{2 K_\text{u}}{\mu_0 \, M_s} \; \vec{e}_\text{u} \; (\vec{e}_\text{u} \cdot \vec{m}),\]with the anisotropy constant \(K_\text{u}\) given in units of \(\text{J/m}^3\).
- Parameters:
Ku (str, optional) – Name of the material parameter for the anisotropy constant \(K_\text{u}\), defaults to “Ku”
Ku_axis (str, optional) – Name of the material parameter for the anisotropy axis \(\vec{e}_\text{u}\), defaults to “Ku_axis”
CubicAnisotropyField
- class magnumnp.CubicAnisotropyField(**kwargs)[source]
Bases:
FieldTerm
Cubic Anisotropy Field:
\[\begin{split}\vec{h}^\text{c} = -\frac{2 K_\text{c1}}{\mu_0 \, M_s} \; \begin{pmatrix} m_1 \, m_2^2 + m_1 \, m_3^2 \\ m_2 \, m_3^2 + m_2 \, m_1^2 \\ m_3 \, m_1^2 + m_3 \, m_2^2\end{pmatrix} -\frac{2 K_\text{c2}}{\mu_0 \, M_s} \; \begin{pmatrix} m_1 \, m_2^2 \, m_3^2 \\ m_1^2 \, m_2 \, m_3^2 \\ m_1^2 \, m_2^2 \, m_3\end{pmatrix},\end{split}\]with the anisotropy constants \(K_{c1}\) and \(K_{c2}\) given in units of \(\text{J/m}^3\). If Euler angles \(\alpha\), \(\beta\) and \(\gamma\) are provided, the cubic anisotropy axes are rotated such that the effective magnetization components \(m_i\) with \(i \in \{x,y,z\}\) read
\[m_i = (\mat{A} \vec{e}_i) \cdot \vec{m}\]with
\[\begin{split}\small \mat{A} = \begin{pmatrix} \cos(\alpha) \cos(\gamma) - \cos(\beta) \sin(\alpha) \sin(\gamma) & -\cos(\beta) \cos(\gamma) \sin(\alpha) - \cos(\alpha) \sin(\gamma) & \sin(\alpha) \sin(\beta)\\ \cos(\gamma) \sin(\alpha) + \cos(\alpha) \cos(\beta) \sin(\gamma) & \cos(\alpha) \cos(\beta) \cos(\gamma) - \sin(\alpha) \sin(\gamma) & -\cos(\alpha) \sin(\beta)\\ \sin(\beta) \sin(\gamma) & \cos(\gamma) \sin(\beta) & \cos(\beta) \end{pmatrix}\end{split}\]- Parameters:
Kc1 (str, optional) – Name of the material parameter for the anisotropy constant \(K_\text{c1}\), defaults to “Kc1”
Kc2 (str, optional) – Name of the material parameter for the anisotropy constant \(K_\text{c2}\), defaults to “Kc2”
Kc_alpha (str, optional) – Euler angle \(\alpha\) for the rotation of the anisotropy axes
Kc_beta (str, optional) – Euler angle \(\beta\) for the rotation of the anisotropy axes
Kc_gamma (str, optional) – Euler angle \(\gamma\) for the rotation of the anisotropy axes
CubicAnisotropyField2
- class magnumnp.CubicAnisotropyField2(**kwargs)[source]
Bases:
FieldTerm
Cubic Anisotropy Field (alternative definition):
\[\begin{split}\vec{h}^\text{c} = - \frac{2 K_\text{c1}}{\mu_0 \, M_s} \, (&((\vec{c_2} \cdot \vec{m})^2 + (\vec{c_3} \cdot \vec{m})^2) \; (\vec{c_1} \cdot \vec{m}) \; \vec{c_1} + \\ &((\vec{c_1} \cdot \vec{m})^2 + (\vec{c_3} \cdot \vec{m})^2) \; (\vec{c_2} \cdot \vec{m}) \; \vec{c_2} + \\ &((\vec{c_1} \cdot \vec{m})^2 + (\vec{c_2} \cdot \vec{m})^2) \; (\vec{c_3} \cdot \vec{m}) \; \vec{c_3}) \\ - \frac{2 K_\text{c2}}{\mu_0 \, M_s} \, (&(\vec{c_2} \cdot \vec{m})^2 \; (\vec{c_3} \cdot \vec{m})^2 \; (\vec{c_1} \cdot \vec{m}) \; \vec{c_1} + \\ &(\vec{c_1} \cdot \vec{m})^2 \; (\vec{c_3} \cdot \vec{m})^2 \; (\vec{c_2} \cdot \vec{m}) \; \vec{c_2} + \\ &(\vec{c_1} \cdot \vec{m})^2 \; (\vec{c_2} \cdot \vec{m})^2 \; (\vec{c_3} \cdot \vec{m}) \; \vec{c_3}) \\ - \frac{4 K_\text{c3}}{\mu_0 \, M_s} \, (&((\vec{c_2} \cdot \vec{m})^4 + (\vec{c_3} \cdot \vec{m})^4) \; (\vec{c_1} \cdot \vec{m})^3 \; \vec{c_1} + \\ &((\vec{c_1} \cdot \vec{m})^4 + (\vec{c_3} \cdot \vec{m})^4) \; (\vec{c_2} \cdot \vec{m})^3 \; \vec{c_2} + \\ &((\vec{c_1} \cdot \vec{m})^4 + (\vec{c_2} \cdot \vec{m})^4) \; (\vec{c_3} \cdot \vec{m})^3 \; \vec{c_3})\end{split}\]with the anisotropy constants \(K_{c1}\) and \(K_{c2}\) given in units of \(\text{J/m}^3\).
NOTE: Orthonormal axes need to be provided by the user!
- Parameters:
Kc1 (str, optional) – Name of the material parameter for the anisotropy constant \(K_\text{c1}\), defaults to “Kc1”
Kc2 (str, optional) – Name of the material parameter for the anisotropy constant \(K_\text{c2}\), defaults to “Kc2”
Kc_axis1 (str, optional) – Name of the material parameter for the 1. axis
Kc_axis2 (str, optional) – Name of the material parameter for the 2. axis
DemagField
- class magnumnp.DemagField(p=20, cache_dir=None)[source]
Bases:
LinearFieldTerm
Demagnetization Field:
The dipole-dipole interaction gives rise to a long-range interaction. The integral formulation of the corresponding Maxwell equations can be represented as convolution of the magnetization \(\vec{M} = M_s \; \vec{m}\) with a proper demagnetization kernel \(\vec{N}\)
\[\vec{h}^\text{dem}_{\vec{i}} = \sum\limits_{\vec{j}} \vec{N}_{\vec{i} - \vec{j}} \, \vec{M}_{\vec{j}},\]The convolution can be evaluated efficiently using an FFT method.
- Parameters:
p (int, optional) – number of next neighbors for near field via Newell’s equation (default = 20)
DemagFieldPBC
- class magnumnp.DemagFieldPBC(**kwargs)[source]
Bases:
LinearFieldTerm
DemagFieldNonEquidistant
- class magnumnp.DemagFieldNonEquidistant(p=20)[source]
Bases:
LinearFieldTerm
Demagnetization Field:
The dipole-dipole interaction gives rise to a long-range interaction. The integral formulation of the corresponding Maxwell equations can be represented as convolution of the magneti_dstation \(\vec{M} = M_s \; \vec{m}\) with a proper demagneti_dstation kernel \(\vec{N}\)
\[\vec{h}^\text{dem}_{\vec{i}} = \sum\limits_{\vec{j}} \vec{N}_{\vec{i} - \vec{j}} \, \vec{M}_{\vec{j}},\]The convolution can be evaluated efficiently using an FFT method.
- Parameters:
p (int, optional) – number of next neighbors for near field via Newell’s equation (default = 20)
ExchangeField
- class magnumnp.ExchangeField(domain=None, **kwargs)[source]
Bases:
LinearFieldTerm
Exchange Field
\[\vec{h}^\text{ex}_i = \frac{2}{\mu_0 \, M_{s,i}} \; \sum_{k=\pm x, \pm y,\pm z} \frac{2}{\Delta_k} \frac{A_{i+\vec{e}_k} \; A_i}{A_{i+\vec{e}_k} + A_i} \; \left( \vec{m}_{i+\vec{e}_k} - \vec{m}_i \right),\]with the vacuum permeability \(\mu_0\), the saturation magnetization \(M_s\), and the exchange constant \(A\). \(\Delta_k\) and \(\vec{e}_k\) represent the grid spacing and the unit vector in direction \(k\), respectively.
- Parameters:
A (str, optional) – Name of the material parameter for the exchange constant \(A\), defaults to “A”
ExternalField
- class magnumnp.ExternalField(h=None)[source]
Bases:
object
External Field
- Parameters:
h (list or tuple or
Tensor
or function) – External Field- Examples:
# homogenious, constant field external = ExternalField([Hx, 0, 0]) # homogenious, time-dependent field external = ExternalField(lambda state: [Hx*state.t, 0, 0]) # inhomogenious, constant field x, y, z = SpatialCoordinate(state) external = ExternalField(Expression([x,y,z]))
DMIField
- class magnumnp.DMIField(dmi_vector, **kwargs)[source]
Bases:
LinearFieldTerm
General Dzyaloshinskii-Moriya interaction
The general expression for the DMI field can be expressed as
\[\vec{h}^\text{dmi}(\vec{x}) = \frac{2 \, D}{\mu_0 \, M_s} \; \sum_{k=x,y,z} \vec{e}^\text{dmi}_k \times \frac{\partial \vec{m}}{\partial \vec{e}_k},\]with the DMI strength \(D\) and the DMI vectors \(\vec{e}^\text{dmi}_k\), which describe which components of the gradient of \(\vec{m}\) contribute to which component of the corresponding field. It is assumed that \(\vec{e}^\text{dmi}_{-k} = -\vec{e}^\text{dmi}_k\).
The occuring gradient is discretized using central differences which finally yields
\[\vec{h}^\text{dmi}_i = \frac{2 \, D_i}{\mu_0 \, M_{s,i}} \; \sum_{k=\pm x, \pm y,\pm z} \frac{\vec{e}^\text{dmi}_k \times \vec{m}_{i+\vec{e}_k}}{2 \, \Delta_k}.\]- Parameters:
Ku (str, optional) – Name of the material parameter for the anisotropy constant \(K_\text{u}\), defaults to “Ku”
Ku_axis (str, optional) – Name of the material parameter for the anisotropy axis \(\vec{e}_\text{u}\), defaults to “Ku_axis”
BulkDMIField
- class magnumnp.BulkDMIField(Db='Db')[source]
Bases:
DMIField
Bulk Dzyaloshinskii-Moriya interaction
\[\vec{h}^\text{dmib}(\vec{x}) = -\frac{2 \, D_b}{\mu_0 \, M_s} \; \nabla \times \vec{m},\]with the DMI strength \(D_b\). The corresponding DMI vectors are \(\vec{e}^\text{dmi}_x = [1, 0, 0]\), \(\vec{e}^\text{dmi}_y = [0, 1, 0]\), and \(\vec{e}^\text{dmi}_z = [0, 0, 1]\).
- Parameters:
Db (str, optional) – Name of the material parameter for the anisotropy constant \(D_i\), defaults to “Di”
InterfaceDMIField
- class magnumnp.InterfaceDMIField(Di='Di')[source]
Bases:
DMIField
Interface Dzyaloshinskii-Moriya interaction
\[\vec{h}^\text{dmii}(\vec{x}) = -\frac{2 \, D_i}{\mu_0 \, M_s} \; \left[ \nabla \left(\vec{e}_z \cdot \vec{m} \right) - \left(\nabla \cdot \vec{m} \right) \, \vec{e}_z\right],\]with the DMI strength \(D_i\) and an interface normal \(\vec{e}_z\) in z-direction. The corresponding DMI vectors are \(\vec{e}^\text{dmi}_x = [ 0, 1, 0]\), \(\vec{e}^\text{dmi}_y = [-1, 0, 0]\), and \(\vec{e}^\text{dmi}_z = [0, 0, 0]\).
- Parameters:
Di (str, optional) – Name of the material parameter for the anisotropy constant \(D_i\), defaults to “Di”
D2dDMIField
- class magnumnp.D2dDMIField(DD2d='DD2d')[source]
Bases:
DMIField
D2d Dzyaloshinskii-Moriya interaction
\[\vec{h}^\text{dmib}(\vec{x}) = -\frac{2 \, D_b}{\mu_0 \, M_s} \; \nabla \times \vec{m},\]with the DMI strength \(D_{D2d}\). The corresponding DMI vectors are \(\vec{e}^\text{dmi}_x = [-1, 0, 0]\), \(\vec{e}^\text{dmi}_y = [0, 1, 0]\), and \(\vec{e}^\text{dmi}_z = [0, 0, 0]\).
- Parameters:
DD2d (str, optional) – Name of the material parameter for the anisotropy constant \(D_{D2d}\), defaults to “DD2d”
RKKYField
- class magnumnp.RKKYField(J_rkky, dir, id1, id2, order=0)[source]
Interlayer-Exchange interaction between two layers gives rise to the following energy contribution:
\[E^\text{rkky} = -\int\limits_\Gamma J_\text{rkky} \, \vec{m}_i \cdot \vec{m}_j \, d\vec{A},\]where \(\Gamma\) is the interface between two layers \(i\) and \(j\) with magnetizations \(\vec{m}_i\) and \(\vec{m}_j\), respectively.
Special care has to be taken, if domain walls or partial domain walls are formed across the RKKY interface. In this case higher order approximations of the magnetization needs to be used near the interface in order to accurately describe e.g. the equilibrium magnetization or critical switching fields. (see Suess et al. “Accurate finite difference micromagnetic of magnetis including RKKY interaction – analytical solution and comparision to standard micromagnetic codes”)
- Parameters:
J_rkky (float) – Interlayer-Exchange constant \(J_\text{rkky}\)
dir – normal direction of the interface (currently “z” is hard-coded”)
id1 (int) – Index of the first layer
id2 (int) – Index of the second layer
order (int, optional) – appoximation order of the magnetization near the interface (default = 0)
- Example:
# create state with named domains from mesh state = State(mesh) # create domains as bool arrays, e.g: domain1 = torch.zeros(n, dtype=torch.bool) domain1[n[0]//2:,:,:] = True domain2 = torch.zeros(n, dtype=torch.bool) domain2[:-n[0]//2:,:,:] = True # rotate magnetization within one subdomain state.m[domain1] = torch.tensor([np.cos(phi), np.sin(phi), 0]) # without interface layer, two seperate exchange fields need to be defined exchange1 = ExchangeField(Aex1, domain1) exchange2 = ExchangeField(Aex2, domain2) rkky = RKKYField(J_rkky, "z", id1, id2)
BiquadraticRKKYField
- class magnumnp.BiquadraticRKKYField(J_rkky_BQ, dir, id1, id2)[source]
Biquadratic surface exchange couplong between two layers gives rise to the following energy contribution:
\[E^\text{biquadratic} = -\int\limits_\Gamma J_\text{biquadratic} \, (\vec{m}_i \cdot \vec{m}_j)^2 \, d\vec{A},\]where \(\Gamma\) is the interface between two layers \(i\) and \(j\) with magnetizations \(\vec{m}_i\) and \(\vec{m}_j\), respectively.
The effective field is given by:
\[\vec{h}^\text{biquadratic}_i = \frac{2 J_\text{biquadratic}} {M_s \Delta z \mu_0} \, (\vec{m}_i \cdot \vec{m}_j) \, \vec{m}_j,\]with the interlayer exchange constant \(J_\text{biquadratic}\).
SpinOrbitTorque
- class magnumnp.SpinOrbitTorque[source]
General spin torque contributions can be described by the following field
\[\vec{h}^\text{sot} = -\frac{j_e \hbar}{2 e \mu_0 M_s} \left[\eta_\text{damp} \, \vec{m} \times \vec{p} + \eta_\text{field} \, \vec{p} \right],\]with the current density \(j_e\), the reduced Planck constant \(\hbar\), the elementary charge \(e\), and the polarization of the electrons \(\vec{p}\). \(\eta_\text{damp}\) and \(\eta_\text{field}\) are material parameters which describe the amplitude of damping- and field-like torque.
In case of Spin-Orbit-Torqe (SOT) \(\eta_\text{field}\) and \(\eta_\text{damp}\) are constant material parameters.
SpinTorqueZhangLi
- class magnumnp.SpinTorqueZhangLi[source]
Zhang Lie spin torque contributions can be described by the following field
\[\vec{h}^\text{stt,zl} = \frac{b}{\gamma} \left[\vec{m} \times (\vec{j}_e \cdot \nabla) \vec{m} + \xi \; (\vec{j}_e \cdot \nabla) \vec{m} \right],\]with the reduced gyromagnetic ratio \(\gamma\), the degree of nonadiabacity \(\xi\). \(b\) is the polarization rate of the conducting electrons and can be written as
\[b = \frac{\beta \mu_B}{e M_s (1+\xi^2)},\]with the Bohr magneton \(\mu_B\), and the dimensionless polarization rate \(\beta\).
SpinTorqueSlonczewski
- class magnumnp.SpinTorqueSlonczewski[source]
Slonczewski spin torque contributions can be described by the following field:
\[\vec{h}^\text{stt,slonczewski} = \beta \left[ \epsilon (\vec{m} \times \vec{p}_\text{p}) + \epsilon' \vec{p}_\text{p} \right],\]with the polarization vector of the conducting electrons \(\vec{p}_\text{p}\) and the secondary spin torque efficiency \(\epsilon'\).
The \(\beta\) is given by current density \(J\), the thickness of the fixed layer \(d\) (default using the thickness of the mesh), the reduced Planck constant \(\hbar\), the elementary charge \(e\), the permeability of free space \(\mu_0\),
\[\beta = \frac{\hbar J} {e M_s d \mu_0},\]The \(\epsilon\) is computed by assuming both the fixed layer and free layer have the same spin polarization \(P\) and the spin diffusion length \(\Lambda\): across the spacer layer
\[\epsilon = \frac{P \Lambda^2} {(\Lambda^2 + 1) + (\Lambda^2 - 1) \vec{m} \cdot \vec{p}_\text{p}}.\]
OerstedField
- class magnumnp.OerstedField(p=20, cache_dir=None)[source]
The Oersted field created by some current density \(\vec{j}\) can be calculated by means of the Biot-Savart law
\[\vec{h}^\text{oersted}(\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}'.\]The occuring equations [krueger] look very similar to those of the demag field, and the occuring convolution can be efficiently calculated by means of an FFT method.
- Parameters:
p (int, optional) – number of next neighbors for near field via Krueger’s equations (default = 20)