Python API Reference

Field Classes#

All field classes follow the same interface:

  • compute(magnetization) – returns the effective field as an (n, 3) array
  • energy(magnetization) – returns the energy as a scalar (Joules)

Magnetization m is always an (n, 3) array of unit vectors.


DemagFieldFK#

Demagnetizing field via the Fredkin-Koehler scalar potential method. This is the primary solver and supports multiple BEM backends.

DemagFieldFK(coordinates, elements, Ms, *,
           scale=1.0, tolerance=1e-8,
           solver=SolverType.AMG_CG,
           h2_cache="", h2params=None,
           fmm_params=None, backend=None,
           tree_type=None)
ParameterTypeDescription
coordinates(n, 3)Node positions
elements(ntet, 4)Tetrahedral connectivity, 0-indexed
Ms(n,)Saturation magnetization per node (A/m)
scalefloatLength unit: 1.0 = m, 1e-9 = nm
tolerancefloatIterative solver tolerance
solverSolverTypeGPU solver preset
h2_cachestrH2 matrix cache file path (empty = no caching)
h2paramsParamsFKH2LibH2-matrix compression parameters
fmm_paramsParamsFKFMMGPU FMM parameters
backendBemBackendBEM backend (default: FMM_GPU)
tree_typeTreeType or strFMM tree strategy

Additional methods:

  • energy_weak_form(m) – energy via element-level FEM weak form, more stable on meshes with sharp potential transitions

Convenience subclasses select a specific backend:

  • DemagFieldFKH2Lib – CPU H2-matrix BEM
  • DemagFieldFKFMM – GPU FMM (same as default)
  • DemagFieldFKDense – dense BEM, $O(n^2)$, exact reference
  • DemagFieldFKGalerkin – Galerkin BEM, $O(n^2)$, smoothed

DemagFieldShell#

Shell-transform demagnetizing field solver. Replaces the BEM boundary integral with a mapped shell region. Fully GPU – no CPU round-trip.

DemagFieldShell(coordinates, elements, element_zones, Ms, *,
                scale=1.0, tolerance=1e-8,
                solver=SolverType.AMG_CG,
                transform=ShellTransform.PARALLELEPIPED)
ParameterTypeDescription
element_zones(ntet,)Zone tag: >0 = magnetic, 0 = air, -1 = shell
transformShellTransformShell mapping type

DemagFieldVolFMM#

Volume FMM demagnetizing field solver. Bypasses the Fredkin-Koehler decomposition entirely – computes the demagnetizing field directly via volumetric FMM. Fully GPU.

DemagFieldVolFMM(coordinates, elements, Ms, *,
                 scale=1.0, fmm_vol_params=None)
ParameterTypeDescription
fmm_vol_paramsParamsVolFMMVolume FMM parameters

ExchangeField#

Exchange field via the FEM stiffness matrix with GPU SpMV.

ExchangeField(coordinates, elements, A, *, scale=1.0)
ParameterTypeDescription
A(n,)Exchange stiffness constant per node (J/m)

AnisotropyField#

Uniaxial magnetocrystalline anisotropy field.

AnisotropyField(coordinates, elements, Ku, axis, *, scale=1.0)
ParameterTypeDescription
Ku(n,)Anisotropy constant per node (J/m$^3$)
axis(n, 3)Easy axis unit vectors per node

ZeemanField#

External (Zeeman) field.

ZeemanField(coordinates, elements, Ms, field, *, scale=1.0)
ParameterTypeDescription
Ms(n,)Saturation magnetization per node (A/m)
field(n, 3)External field per node (A/m)

compute() returns $\mu_0 M_s \mathbf{H}_\text{ext}$.

Additional methods:

  • set_field(field) – update the external field without rebuilding

Composable Effective Fields#

EffectiveField#

Low-level compositor. You construct individual field objects and add them by name.

eff = EffectiveField(coordinates, elements, scale=1.0)
eff.add("demag", DemagFieldFK(coords, elems, Ms, scale=scale))
eff.add("exchange", ExchangeField(coords, elems, A, scale=scale))
Method / PropertyDescription
add(name, field)Add a named field component (ownership is moved)
compute(m)Total effective field (n, 3)
energy(m)Total energy (Joules)
compute_component(name, m)Single component’s field
energy_component(name, m)Single component’s energy
component_namesList of registered component names
node_volumePer-node volumes (n,)
num_nodesNumber of mesh nodes

EffectiveFieldFK#

High-level compositor for BEM-based models. Owns the mesh and $M_s$, and provides add_*() / set_*() methods.

eff = EffectiveFieldFK(coordinates, elements, Ms, scale=1.0)
MethodDescription
add_demag(...)Add demagnetizing field (see DemagFieldFK for parameters)
add_exchange(A)Add exchange field
add_anisotropy(Ku, axis)Add uniaxial anisotropy
add_anisotropy_cubic(Kc1, Kc2, cubic_axes)Add cubic anisotropy
add_zeeman(Hext)Add Zeeman field
set_Ms(Ms)Update $M_s$ (call recompute() to apply)
set_exchange(A)Update exchange stiffness (call recompute())
set_anisotropy(Ku, axis)Update uniaxial anisotropy (call recompute())
set_zeeman(Hext)Update external field (no recompute needed)
recompute()Rebuild dirty components

Also has compute(), energy(), compute_component(), energy_component(), component_names, node_volume, num_nodes – same as EffectiveField.


EffectiveFieldShell#

Same interface as EffectiveFieldFK, but for shell-transform meshes. Exchange, anisotropy, and Zeeman operate on body-only elements.

eff = EffectiveFieldShell(coordinates, elements, element_zones, Ms, scale=1.0)

The element_zones array tags each tetrahedron: >0 = magnetic body, 0 = air, -1 = shell.


EffectiveFieldVolFMM#

Same interface as EffectiveFieldFK, but uses DemagFieldVolFMM (direct volumetric FMM, no BEM or FK decomposition) for the demagnetizing field.

eff = EffectiveFieldVolFMM(coordinates, elements, Ms, scale=1.0)
MethodDescription
add_demag(fmm_vol_params=None)Add demagnetizing field (volume FMM)

All other methods (add_exchange, add_anisotropy, add_zeeman, set_*, recompute, compute, energy, etc.) are the same as EffectiveFieldFK.


Enums#

SolverType#

GPU linear solver presets (AMGCL with CUDA backend).

ValuePreconditionerKrylovNotes
AMG_GMRESAMG + ILU0GMRES
AMG_CGAMG + ILU0CGDefault. Faster for SPD systems.
AMG_JACOBI_CGAMG + JacobiCGAvoids cuSPARSE ILU0.

BemBackend#

ValueDescription
FMM_GPUGPU Fast Multipole Method (default)
H2LIBCPU H2-matrix
DENSE_BEMDense $O(n^2)$ reference
GALERKIN_BEMGalerkin DLP BEM

TreeType#

FMM tree construction strategy (only affects FMM_GPU backend).

ValueDescription
OCTREEUniform cubic octree (default)
ADAPTIVE_OCTREEAxis-aligned bounding box
KD_TREEBinary k-d tree (best for high aspect ratio)

String shortcuts: "octree", "adaptive", "kdtree".

ShellTransform#

ValueDescription
PARALLELEPIPEDParallelepiped mapping
SPHERICALSpherical mapping

Parameter Structs#

ParamsFKH2Lib#

H2-matrix compression tuning (only affects H2LIB backend).

FieldDefaultDescription
q_reg3Quadrature order for regular integrals (singular = q_reg + 2)
clf32Cluster leaf size
eta2.0Admissibility parameter
eps_aca1e-8ACA approximation tolerance
eps_recomp1e-5H2 recompression tolerance
m_hca5HCA interpolation order

ParamsFKFMM#

GPU FMM accuracy and performance.

FieldDefaultDescription
expansion_order10Solid harmonic truncation order $p$; $(p+1)^2$ coefficients
leaf_size64Maximum boundary nodes per leaf cell
eta0.5Admissibility: cells are “far” if diam $< \eta \times$ dist
quad_order3Triangle quadrature order
p2p_onlyFalseForce P2P-only mode (no far-field)
tree_typeOCTREETree construction strategy

ParamsVolFMM#

Volume FMM parameters.

FieldDefaultDescription
expansion_order8Spherical harmonic truncation $p$
leaf_size64Octree leaf size