Additional functionality available with mammos_mumag#

Note: The functionality of mammos_mumag mentioned here is largly irrelevant to the average user.

At the moment, mammos_mumag primarily supports calculation of hysteresis loop and extrinsic magnetic properties (coercive field, remanent magnetisation, and maximum energy product). However, there are a few additional functionalities available with the package which are shown here using an example.

[1]:
import mammos_entity as me
import mammos_units as u
import pandas as pd
from mammos_mumag.materials import Materials
from mammos_mumag.mesh import Mesh
from mammos_mumag.parameters import Parameters
from mammos_mumag.simulation import Simulation

u.set_enabled_equivalencies(u.magnetic_flux_field());

Inputs#

Geometry and Mesh#

We load one of the meshes available with the package, cube20_singlegrain_msize2, defining a cube of length 20 with no grains and mesh size 2.

The magnetostatic potential is set to zero at infinity. To treat the open boundary condition numerically, we can apply the spherical shell transformation. A spherical shell surrounding the magnet is transformed to fill the space between the magnet and infinity. See: Imhoff, J. F., et al. “An original solution for unbounded electromagnetic 2D-and 3D-problems throughout the finite element method.” IEEE Transactions on Magnetics 26.5 (1990): 1659-1661.

So the mesh will effectively be a cube enclosed with a spherical shell. The region inside the cube corresponds to the magnetic material, while the one between the cube and the inner sphere of the shell corresponds to non-magnetic material.

This mesh was generated by Salome in the unv mesh format, and subsequently converted to the fly format.

[2]:
mesh = Mesh("cube20_singlegrain_msize2")

Material parameters#

We define the magnetic material parameters separately for the 3 domains: the magnetic material, the non-magnetic material and the shell. Therefore the parameters strictly related to the magnetic material will be only defined on the first domain. The other domain, representing vacuum, will have magnetic properties equal to zero.

[3]:
mat = Materials(
    domains=[
        {  # cube
            "theta": 0.0,
            "phi": 0.0,
            "K1": me.Ku(2805163.0590566983, unit=u.J / u.m**3),
            "K2": me.Ku(0.0, unit=u.J / u.m**3),
            "Ms": me.Ms(1189556.9164459957, unit=u.A / u.m),
            "A": me.A(6.015043651663635e-12, unit=u.J / u.m),
        },
        {  # sphere
            "theta": 0.0,
            "phi": 0.0,
            "K1": me.Ku(0.0, unit=u.J / u.m**3),
            "K2": me.Ku(0.0, unit=u.J / u.m**3),
            "Ms": me.Ms(0.0, unit=u.A / u.m),
            "A": me.A(0.0, unit=u.J / u.m),
        },
        {  # shell
            "theta": 0.0,
            "phi": 0.0,
            "K1": me.Ku(0.0, unit=u.J / u.m**3),
            "K2": me.Ku(0.0, unit=u.J / u.m**3),
            "Ms": me.Ms(0.0, unit=u.A / u.m),
            "A": me.A(0.0, unit=u.J / u.m),
        },
    ],
)

Simulation parameters#

Lastly, we define the simulation parameters. For an exhaustive explanation of these parameters, check the documentation.

[4]:
par = Parameters(
    size=1.0e-9,
    scale=0,
    h_start=(1 * u.T).to("A/m"),
    h_final=(-1 * u.T).to("A/m"),
    h_step=(-0.2 * u.T).to("A/m"),
    h_vect=[0.01745, 0, 0.99984],
    m_step=(0.4 * u.T).to("A/m"),
    m_final=(-1.2 * u.T).to("A/m"),
    m_vect=[0, 0, 1],
    tol_fun=1e-10,
    tol_h_mag_factor=1,
    precond_iter=10,
)

Putting all inputs in a Simulation object#

To define a Simulation object, we need to define a mesh, a Materials object, and a Parameters object.

[5]:
sim = Simulation(
    mesh=mesh,
    materials=mat,
    parameters=par,
)

Note that all of this could also have been defined using file paths:

sim = Simulation(
    mesh=Mesh("..."),
    materials_filepath="...",
    parameters_filepath="...",
)

Available methods beyond hysteresis loop calculation#

For all methods, we can specify the optional arguments outdir and name. While the first identifies the output directory where the input and output files will be stored (and where the script is executed), the name defines the output file names.

The outdir argument can be used to give an ID to the simulation.

Save the mesh and the materials#

To create a vtu file for the visualise the distribution of material properties within the geometry, we can use the run_materials method.

[6]:
sim.run_materials(outdir="out/materials", name="cube")

This creates a discretised representation of material properties as scalar functions and fields in the mesh and stores it as out/materials/cube_mat.vtu.

Compute the magnetostatic field#

To create the vtu file for visualisation of the magnetic scalar potential and the demagnetisation field, we use run_hmag method. With linear basis function for the magnetic scalar potential \(u\), the magnetostatic field \(h = - \nabla u\) is constant within a given finite-element cell. Further, the magnetostatic field is projected on a function space with linear basis function (the same as \(u\)) in order to obtain its values at the nodes. This projected function is defined as h_at_nodes in the vtu file.

[7]:
sim.run_hmag(outdir="out/hmag", name="cube")

The output cube_hmag.vtu will be stored in the output directory out/hmag.

The scripts creates two file: the magnetostatic field, as seen above, will be stored in cube_hmag.vtu. At the same time the software also gives the magnetostatic energy density computed with finite elements and compares it with the analytic solution:

  • from field: \begin{equation} E_{\mathsf{dmg}} = \frac{\mu_0}{2} \int_\Omega \frac{\mathbf{h} \cdot M_s \mathbf{m}}{V} \ \mathrm{d}x, \end{equation} where \(\Omega\) is the domain, \(\mathbf{h}\) is the demagnetization field, \(M_s\) is the spontaneous magnetisation, \(\mathbf{m}\) is the magnetisation vector field, and \(V\) is the volume of the domain.

  • from_gradient: \begin{equation} \frac{1}{2} \sum_i \mathbf{m}_i \cdot \mathbf{g}_i, \end{equation} where \(\mathbf{m}_i\) and \(\mathbf{g}_i\) are the unit vector of the magnetization and the gradient of the energy normalized by the volume of the energy with respect to \(\mathbf{m}_i\) at the nodes of the finite element mesh.

  • analytic: \(J_s^2 / (6 \mu_0)\)

This information is saved in cube.csv:

[8]:
pd.read_csv("out/hmag/cube.csv", skiprows=1)
[8]:
name value explanation
0 E_field 265453.504537 Energy density evaluated from field (J/m^3).
1 E_gradient 265453.504537 Energy density evaluated from gradient (J/m^3).
2 E_analytic 296366.469827 Energy density evaluated analytically (J/m^3).

Exchange and anisotropy energy#

To test the computation of the exchange and anisotropy energy density we can use the run_exani method.

This gives the exchange energy density of a vortex in the \(xy\)-plane and the anistropy energy density in the uniformly magnetized state. Here we have placed the anistropy direction paralle to to the \(z\)-axis. The anisotropy energy density is calculated as \(-K (\mathbf{m} \cdot \mathbf{k})^2\) where \(\mathbf{m}\) is the unit vector of magnetization and \(\mathbf{k}\) is the anisotropy direction. \(K\) is the magnetocrystalline anisotropy constant.

[9]:
sim.run_exani(outdir="out/exani", name="cube")

The accuracy of this script is then analyzed for two different magnetization, a vortex and a uniform vector.

[10]:
pd.read_csv("out/exani/cube_vortex.csv", skiprows=1)
[10]:
name value explanation
0 E_gradient 74103.192040 Energy evaluated from gradient (J/m^3).
1 E_analytic 74207.626622 Energy evaluated analytically (J/m^3).
[11]:
pd.read_csv("out/exani/cube_uniform.csv", skiprows=1)
[11]:
name value explanation
0 E_gradient -2.805163e+06 Energy evaluated from gradient (J/m^3).
1 E_analytic -2.805163e+06 Energy evaluated analytically (J/m^3).

Zeeman energy#

The method run_external calculates the Zeeman energy density for an external field of \(\mu_0 H_{\mathsf{ext}} = 1.2 \ T\) using finite elements and analytically.

[12]:
sim.run_external(outdir="out/external", name="cube")

The energy densities are written in the generated file cube.csv:

[13]:
pd.read_csv("out/external/cube.csv", skiprows=1)
[13]:
name value explanation
0 E_gradient -1.427251e+06 Energy evaluated from gradient (J/m^3).
1 E_analytic -1.427251e+06 Energy evaluated analytically (J/m^3).

jax implementation#

The above tools checked the energy calculation with the finite element backend. From the finite element backend system matrices are generated for micromagnetic simulations. The method run_mapping is used to test the energy calculations with matrices.

The mammos_mumag software uses sparse matrix methods from jax.

[14]:
sim.run_mapping(outdir="out/mapping", name="cube")

Information about the calculation are stored in different files. Among them, cube_energy.csv stores the total energy density for the uniformly magnetized state:

[15]:
pd.read_csv("out/mapping/cube_energy.csv", skiprows=1)
[15]:
name value explanation
0 E_jax -3966956.5 Energy evaluated with jax backend (J/m^3).
1 E_analytic -3966956.5 Energy evaluated analytically (J/m^3).

The file cube_stats.txt, on the other hand, shows information about memory and runtime.

[16]:
with open("out/mapping/cube_stats.txt") as file:
    print(file.read())
MAP FINITE ELEMENT BACKEND (esys-escript) TO JAX.
Memory before escript2jax: 403.2734375 MB.
Memory after  escript2jax: 609.01953125 MB.
Memory after garbage collection: 609.01953125 MB.
Timing and statistics.
elapsed time: 0.26778578758239746
function_calls: 1

Storing sparse matrices#

The sparse matrices used for computation can be stored and reused for simulations with the same finite element mesh. To store the matrices use the method run_store.

[17]:
sim.run_store(outdir="out/store", name="cube")