API for users
Meshes – MicroMagnetic
MicroMagnetic.FDMesh
— MethodFDMesh(; dx=1e-9, dy=1e-9, dz=1e-9, nx=1, ny=1, nz=1, pbc="open")
Create a finite difference mesh (FDMesh
) for simulations. This method is designed for the Micromagnetic model.
This constructor initializes a finite difference mesh with the specified parameters. The pbc
(periodic boundary conditions) parameter can be any combination of the characters "x"
, "y"
, and "z"
, which determines whether the mesh wraps around the corresponding axis.
Parameters:
dx::Float64
: Grid spacing in the x-direction (default:1e-9
).dy::Float64
: Grid spacing in the y-direction (default:1e-9
).dz::Float64
: Grid spacing in the z-direction (default:1e-9
).nx::Int
: Number of grid points in the x-direction (default:1
).ny::Int
: Number of grid points in the y-direction (default:1
).nz::Int
: Number of grid points in the z-direction (default:1
).pbc::String
: Periodic boundary conditions, which could be any combination of"x"
,"y"
, and"z"
(default:"open"
).
Examples:
mesh = FDMesh(dx=1e-8, dy=1e-8, dz=1e-8, nx=10, ny=10, nz=10, pbc="xyz")
This creates a FDMesh with 10 grid points in each direction and periodic boundary conditions in all three axes.
Meshes – Atomistic
MicroMagnetic.CubicMesh
— MethodCubicMesh(;dx=1e-9, dy=1e-9, dz=1e-9, nx=1, ny=1, nz=1, pbc="open")
Create a simple cubic mesh, in which each cell is only connected to six nearest neighbors and forming a simple cube structure.
MicroMagnetic.TriangularMesh
— MethodTriangularMesh(; dx=1e-9, nx=1, ny=1, nz=1, pbc="open")
Create a triangular mesh. The index of the nearest neighbours and the next-nearest neighbours are given as follows:
nearest index | location | next-nearest index | location |
---|---|---|---|
1 | right | 1 | top-right |
2 | top-right | 2 | top |
3 | top-left | 3 | top-left |
4 | left | 4 | bottom-left |
5 | bottom-left | 5 | bottom |
6 | bottom-right | 6 | bottom-right |
7 | above | ||
8 | below |
MicroMagnetic.CylindricalTubeMesh
— MethodCylindricalTubeMesh(;dz=1e-9, R=20e-9, nz=1, nr=10, pbc="open")
Create a cylindrical tube mesh along the +z direction.
The spins are located on the cylindrical tube uniformly, and are indexed as follows:
id = index(i, 1, k, nr, 1, nz)
which means that the spins are labelled in a ring firstly, then nz
is the number of rings. The coordinates of the spins at each ring are given as (R cos(2*pi*(i-1)/nr), R sin(2*pi*(nr-1)/nr))
.
The nearest neighbours are indexed as follows:
| 1 2 3 4 |left right bottom top |
Shapes
MicroMagnetic.Cylinder
— TypeCylinder(;center = (0.0, 0.0, 0.0), radius = 1.0, height= Inf, normal = (0.0, 0.0, 1.0))
Create a Cylinder shape.
MicroMagnetic.Box
— TypeBox(;center::Tuple=(0.0, 0.0, 0.0), sides::Tuple=(1.0,1.0,1.0), theta=0.0)
Create a Box shape.
MicroMagnetic.Torus
— TypeTorus(; center::Tuple=(0.0, 0.0, 0.0), R=1.0, r=0.2)
Create a Torus shape.
Interfaces
MicroMagnetic.sim_with
— Functionsim_with(args::Union{NamedTuple, Dict})
High-Level Interface for starting a typical micromagnetic simulation. All parameters are set using args
, which can be either a NamedTuple
or a Dict
.
Keywords
mesh
: A mesh must be provided to start the simulation. The mesh could beFDMesh
,CubicMesh
, orTriangularMesh
.name
: The name of the simulation, provided as a string.task
: The type of simulation task, which can be"Relax"
or"Dynamics"
. The default is "Relax".driver
: The name of the driver, which should be "SD", "LLG", or "LLG_STT". The default is "SD".alpha
: The Gilbert damping parameter in the LLG equation, provided as a number.beta
: The nonadiabatic strength in the LLG equation with spin-transfer torques (Zhang-Li model), provided as a number.gamma
: The gyromagnetic ratio, with a default value of 2.21e5.ux
,uy
,uz
: The components of the spin-transfer torque strength.ufun
: A time-dependent function foru
.Ms
: The saturation magnetization, which should be aNumberOrArrayOrFunction
. The default isMs=8e5
.mu_s
: The magnetic moment, which should be aNumberOrArrayOrFunction
. The default ismu_s=2*mu_B
.A
orJ
: The exchange constant, which should be aNumberOrArrayOrFunction
.D
: The Dzyaloshinskii-Moriya interaction (DMI) constant, which should be aNumberOrArrayOrFunction
.dmi_type
: The type of DMI, either "bulk" or "interfacial".Ku
: The anisotropy constant, which should be aNumberOrArrayOrFunction
.axis
: The anisotropy axis, provided as a tuple, e.g.,(0, 0, 1)
.Kc
: the cubic anisotropy constant, should beNumberOrArrayOrFunction
.axis1
: the cubic anisotropy axis1, should be a tuple, such as (1,0,0)axis2
: the cubic anisotropy axis2, should be a tuple, such as (0,1,0)demag
: Whether to include demagnetization. This should be a boolean (true
orfalse
). The default isdemag=false
.H
: The external magnetic field, which should be a tuple or function, i.e.,TupleOrArrayOrFunction
.m0
: The initial magnetization, which should be a tuple or function, i.e.,TupleOrArrayOrFunction
.T
: The temperature, which should be aNumberOrArrayOrFunction
.shape
: The shape defines the geometry of the sample, where parameters are configured.steps
: The total number of simulation steps for theDynamics
task.dt
: The time interval of each step, so the total simulation time issteps * dt
for theDynamics
task.max_steps::Int
: Maximum number of steps to run the simulation for theRelax
task. Default is10000
.saver_item
: ASaverItem
instance or a list ofSaverItem
instances. These are custom data-saving utilities that can be used to store additional quantities during the simulation (e.g., guiding centers or other derived values). Ifnothing
, no additional data is saved beyond the default.call_back
: A user-defined function ornothing
. If provided, this function will be called at every step, allowing for real-time inspection or manipulation of the simulation state.stopping_dmdt::Float64
: Primary stopping condition for bothLLG
andSD
drivers. For standard micromagnetic simulations, typical values range from0.01
to1
. InSD
driver mode, where time is not strictly defined, a factor ofγ
is applied to make it comparable to theLLG
driver. For atomistic models using dimensionless units, setusing_time_factor
tofalse
to disable this factor.relax_data_interval::Int
: Interval for saving overall data such as energies and average magnetization during aRelax
task. A negative value disables data saving (e.g.,relax_data_interval = -1
saves data only at the end of the relaxation).dynamic_data_save::Bool
: Boolean flag to enable or disable saving overall data such as energies and average magnetization during theDynamics
task. Set totrue
to enable, orfalse
to disable.relax_m_interval::Int
: Interval for saving magnetization data during aRelax
task. A negative value disables magnetization saving.dynamic_m_interval::Int
: Interval for saving magnetization data during aDynamics
task. A negative value disables magnetization saving.using_time_factor::Bool
: Boolean flag to apply a time factor inSD
mode for comparison withLLG
mode. Default istrue
.save_vtk::Bool
: Boolean flag to save the magnetization to vtk files after finishing each task. Default isfalse
.
Example
See examples at High-Level Interface.
Notes
- Suffix Usage: Parameters like
H
,Ms
,Ku
,A
,D
,task
, anddriver
can be defined as arrays using the_s
or_sweep
suffix. When using these suffixes, ensure that the corresponding array lengths match. This allows the simulation to iterate over different values for these parameters. - Argument Types: The
args
parameter can be either aNamedTuple
or aDict
, providing flexibility in how you organize and pass the simulation parameters. - Driver Selection: The
driver
parameter (ordriver_s
for multiple drivers) specifies the simulation type. Options include"SD"
for the steepest-descent method,"LLG"
for the Landau-Lifshitz-Gilbert equation, and"LLG_STT"
for simulations involving spin-transfer torques. - Stopping Criterion: The
stopping_dmdt
parameter is critical for determining when to stop a simulation, particularly in relaxation tasks. It measures the rate of change in magnetization, with typical values ranging from0.01
to1
. For atomistic models, theusing_time_factor
flag can be set tofalse
to disable time scaling. - Data Saving: The
relax_m_interval
anddynamic_m_interval
parameters control how frequently magnetization data is saved duringRelax
andDynamics
tasks, respectively. Use negative values to disable data saving for these tasks.
MicroMagnetic.set_backend
— Functionset_backend(backend="cuda")
Set the backend of MicroMagnetic. This function allows you to specify the backend for MicroMagnetic simulations.
The available options and their corresponding hardware and backends are shown below:
Option | Hardware | Backend |
---|---|---|
"cpu" | CPU | KernelAbstractions.CPU() |
"cuda" or "nvidia" | NVIDIA GPU | CUDA.CUDABackend() |
"amd" or "roc" | AMD GPU | AMDGPU.ROCBackend() |
"oneAPI" or "intel" | Intel GPU | oneAPI.oneAPIBackend() |
"metal" or "apple" | Apple GPU | Metal.MetalBackend() |
Examples
To set the backend to use CUDA (NVIDIA GPU):
using MicroMagnetic
using CUDA
To set the backend to use the CPU/CUDA:
set_backend("cpu")
set_backend("cuda")
MicroMagnetic.set_precision
— Functionset_precision(x::Type{<:AbstractFloat}=Float64)
Set the precision for MicroMagnetic simulations.
This function allows you to specify the precision of floating-point numbers to be used in MicroMagnetic simulations. By default, it sets the precision to Float64
. If single-precision computation is required, you can specify Float32
.
Example
set_precision(Float32)
MicroMagnetic.set_verbose
— Functionset_verbose(verbose::Bool=false)
Sets the global verbosity level for logging.
Arguments
verbose::Bool
: If true, enables verbose logging; if false, disables it. Default is false.
Example
set_verbose(true) # Enables verbose logging
set_verbose(false) # Disables verbose logging
MicroMagnetic.Sim
— FunctionSim(mesh::Mesh; driver="LLG", name="dyn", integrator="DormandPrince")
Create a simulation instance for given mesh.
MicroMagnetic.NEB
— TypeNEB(sim::AbstractSim, init_images::TupleOrArray, frames_between_images::TupleOrArray;
name="NEB", spring_constant=1.0e5, driver="LLG", clib_image=-1)
Create a NEB instance.
Arguments
sim::AbstractSim
: Simulation instance that includes the interactions.init_images::TupleOrArray
: Tuple or array of initial and final state images. Intermediate state images can also be included.frames_between_images::TupleOrArray
: Tuple or array specifying the number of frames between each pair of images.name::String="NEB"
: Name for the NEB instance.spring_constant::Float64=1.0e5
: Spring constant used in the NEB simulation.driver::String="LLG"
: Driver for the NEB simulation. Options are"SD"
or"LLG"
.clib_image::Int=-1
: Optional parameter for specifying a particular image in the library (default is -1).
Returns
A NEB instance configured with the provided parameters.
Example
# define mesh and the corresponding parameters
mesh = FDMesh(nx=60, ny=60, nz=1, dx=2e-9, dy=2e-9, dz=2e-9, pbc="xy")
params = Dict(
:Ms => 3.84e5,
:A => 3.25e-12,
:D => 5.83e-4,
:H => (0, 0, 120 * mT)
)
# Define the initial and final state, stored in the init_images list.
# Any acceptable object, such as a function, a tuple, or an array, can be used.
# The init_images list can also contain the intermediate state if you have one.
init_images = [read_vtk("skx.vts"), (0, 0, 1)]
# Define the interpolation array to specify the number of images used in the NEB simulation.
# The length of the interpolation array should be the length of init_images minus one.
# For example, if init_images = [read_vtk("skx.vts"), read_vtk("skx2.vts"), (0, 0, 1)],
# the length of interpolation should be 2, i.e., something like interpolation = [5,5].
interpolation = [6]
# Use the create_sim method to create a Sim instance.
sim = create_sim(mesh; params...)
# Create the NEB instance and set the spring_constant. The driver can be "SD" or "LLG".
neb = NEB(sim, init_images, interpolation; name="skx_fm", driver="SD")
# neb.spring_constant = 1e7
# Relax the entire system.
relax(neb; stopping_dmdt=0.1, save_vtk_every=1000, max_steps=5000)
MicroMagnetic.set_driver
— Functionset_driver(sim::AbstractSim; driver="LLG", integrator="DormandPrince", args...)
Set the driver of the simulation, can be used to switch the driver.
MicroMagnetic.init_m0
— Functioninit_m0(sim::MicroSim, m0::TupleOrArrayOrFunction; norm=true)
Set the initial magnetization of the system. If norm=false
the magnetization array will be not normalised. Examples:
init_m0(sim, (1,1,1))
or
init_m0(sim, (1,1,1), norm=false)
or
function uniform_m0(i,j,k,dx,dy,dz)
return (0,0,1)
end
init_m0(sim, uniform_m0)
MicroMagnetic.init_m0_random
— Functioninit_m0_random(sim::MicroSim)
Set the initial magnetization with random direction.
MicroMagnetic.init_m0_skyrmion
— Functioninit_m0_skyrmion(sim::AbstractSim, center::Tuple, R::Float64; ratio=0.7, p=-1, c=1, type="B")
Set the magnetization with skyrmions. Note that this function can be called mulitple times to add more skyrmons.
center : the skyrmion center, should be a Tuple. For example, center = (50e-9,50e-9)
R : the skyrmion radius.
ratio : ratio=w/R where w is the width of domain wall. By default ratio = 0.7
p : polarity, +1 –> core up; -1 –> core down
c : chirality, +1 –> lefthand,for positive D; -1 –> righthand,for negative D
type : "B" or "N", representing Bloch or Neel skyrmions.
For example:
init_m0_skyrmion(sim, (50e-9,50e-9), 2e-8, ratio=0.5, p=-1, c=1, type="B")
MicroMagnetic.add_zeeman
— Functionadd_zeeman(sim::AbstractSim, H0::TupleOrArrayOrFunction; name="zeeman")
Add a static Zeeman energy to the simulation.
add_zeeman(sim::AbstractSim, H0::TupleOrArrayOrFunction, ft::Function; name="timezeeman")
Add a time varying zeeman to system.
The input ft
is a function of time t
and its return value should be a tuple with length 3.
Example:
function time_fun(t)
w = 2*pi*2.0e9
return (sin(w*t), cos(w*t), 0)
end
function spatial_H(i, j, k, dx, dy, dz)
H = 1e3
if i<=2
return (H, H, 0)
end
return (0, 0, 0)
end
add_zeeman(sim, spatial_H, time_fun)
MicroMagnetic.add_anis
— Functionadd_anis(sim::AbstractSim, Ku::NumberOrArrayOrFunction; axis=(0,0,1), name="anis")
Add Anisotropy to the system, where the energy density is given by
\[ E_\mathrm{anis} = K_{u} (1 - \vec{m} \cdot \hat{u})^2\]
MicroMagnetic.add_cubic_anis
— Functionadd_cubic_anis(sim::AbstractSim, Kc::Float64; axis1=(1,0,0), axis2=(0,1,0), name="cubic")
add a cubic anisotropy with default axis (1,0,0) , (0,1,0), and (0,0,1). The third axis is defined as axis3 = axis1 x axis2.
\[ E_\mathrm{cubic} = -\int_{V} K_c (m_x^4 + m_y^4 + m_z^4) \, dV\]
Example
add_cubic_anis(sim, 1e3, (1, 1, 0), (1, -1, 0))
MicroMagnetic.add_hex_anis
— Functionadd_hex_anis(sim::AbstractSim; K1=0, K2=0, K3=0, name="hex")
Add hexagonal anisotropy to a simulation. The energy density of the anisotropy is defined as:
\[E = K_1 \sin^2 \theta + K_2 \sin^4 \theta + K_3 \sin^6 \theta \cos 6\phi\]
Example
add_hex_anis(sim, K1=1e3, K2=0, K3=1e2)
MicroMagnetic.add_she_torque
— Functionadd_she_torque(sim::AbstractSim, sigma_s::Tuple{Real,Real,Real}, sigma_sa::Tuple{Real,Real,Real}; a1=(1,0,0), a2=(0,1,0), beta=0, name="she")
Add an effective field to represent the SHE torque given by
\[\frac{\partial \mathbf{m}}{\partial t} = - \gamma \mathbf{m} \times \mathbf{H} + \alpha \mathbf{m} \times \frac{\partial \mathbf{m}}{\partial t} - \mathbf{m} \times (\mathbf{m} \times \boldsymbol{\sigma}_{x y}^\mathrm{S H E}) - \beta\mathbf{m} \times \boldsymbol{\sigma}_{x y}^\mathrm{S H E}\]
The equivalent effective field is
\[\mathbf{H}_\mathrm{she} = (1/\gamma)(\beta \boldsymbol{\sigma}_{x y}^\mathrm{S H E} + \mathbf{m} \times \boldsymbol{\sigma}_{x y}^\mathrm{S H E})\]
where
\[\boldsymbol{\sigma}_{x y}^\mathrm{S H E}=\widetilde{\boldsymbol{\sigma}}_\mathrm{S} \times \hat{\mathbf{a}}_2-\left({\mathbf{m}} \cdot\left(\hat{\mathbf{a}}_2 \times \hat{\mathbf{a}}_1\right)\right)^2\left(\widetilde{\boldsymbol{\sigma}}_\mathrm{SA} \times \hat{\mathbf{a}}_2\right)\]
Example
add_she_torque(sim, (0.1, 0.2, 0.3), (0.3, 0.4, 0.5), beta=0.1)
MicroMagnetic.add_thermal_noise
— Functionadd_thermal_noise(sim::AbstractSim, Temp::NumberOrArrayOrFunction; name="thermal", scaling=t -> 1.0, k_B=k_B)
Adds thermal noise fields to the simulation. For micromagnetic model, the thermal noise is defined as
\[\mathbf{b}^u = \eta \sqrt \frac{2 \alpha k_B T}{\mu_0 M_s \gamma \Delta V dt}\]
and $\eta$ is a random number follows the normal distribution where $\gamma=2.211\times 10^5$ m/(A·s) is the gyromagnetic ratio.
For the atomistic model, the thermal noise is defined as
\[\mathbf{b}^u = \eta \sqrt \frac{2 \alpha k_B T}{\gamma \mu_s dt}.\]
where $\eta$ is a random number follows the normal distribution and $\gamma=1.76\times 10^{11}$ rad/(T·s).
Arguments
sim::AbstractSim
: The simulation object.Temp::NumberOrArrayOrFunction
: Temperature value (can be a constant, array, or function).name::String
: Name for the noise field (default:"thermal"
).scaling::Function
: A function to scale the noise over time (default:t -> 1.0
).k_B::Float64
: Boltzmann constant (default:k_B
).
Example
# Add thermal noise with a constant temperature of 100 K and a scaling function
add_thermal_noise(sim, 100.0, scaling = t -> exp(-t/1e-9))
MicroMagnetic.update_zeeman
— Functionupdate_zeeman(sim::AbstractSim, H0::TupleOrArrayOrFunction; name="zeeman")
Set the Zeeman field to H0 where H0 is TupleOrArrayOrFunction according to its name. For example,
add_zeeman(sim, (0,0,0), name="my_H") #create a zeeman energy with field (0,0,0) A/m
update_zeeman(sim, (0,0,1e5), name="my_H") #change the field to (0,0,1e5) A/m
MicroMagnetic.update_anis
— Functionupdate_anis(sim::MicroSim, Ku::NumberOrArrayOrFunction; name = "anis")
update anisotropy constant Ku according to its name.
Example:
mesh = FDMesh(nx=200, ny=200, nz=12, dx=5e-9, dy=5e-9, dz=5e-9)
sim = Sim(mesh)
add_anis(sim, 3e4, axis = (0,0,1), name="K1")
add_anis(sim, 1e5, axis = (1,0,0), name="K2")
update_anis(sim, 5e4, name="K2") #update anisotropy K2
MicroMagnetic.create_sim
— Functioncreate_sim(mesh; args...)
Create a micromagnetic simulation instance with given arguments.
mesh
: a mesh has to be provided to start the simulation. The mesh could beFDMesh
,CubicMesh
, orTriangularMesh
.
Arguments
name
: the simulation name, should be a string.driver
: the driver name, should be a string. By default, the driver is "SD".alpha
: the Gilbert damping in the LLG equation, should be a number.beta
: the nonadiabatic strength in the LLG equation with spin transfer torques (zhang-li model), should be a number.gamma
: the gyromagnetic ratio, default value = 2.21e5.ux
,uy
oruz
: the strengths of the spin transfer torque.ufun
: the time-dependent function foru
.Ms
: the saturation magnetization, should beNumberOrArrayOrFunction
. By default, Ms=8e5mu_s
: the magnetic moment, should beNumberOrArrayOrFunction
. By default, mus=2*muBA
orJ
: the exchange constant, should beNumberOrArrayOrFunction
.D
: the DMI constant, should beNumberOrArrayOrFunction
.dmi_type
: the type of DMI, could be "bulk", "interfacial" or "D2d".Ku
: the anisotropy constant, should beNumberOrArrayOrFunction
.axis
: the anisotropy axis, should be a tuple, such as (0,0, 1)Kc
: the cubic anisotropy constant, should beNumberOrArrayOrFunction
.axis1
: the cubic anisotropy axis1, should be a tuple, such as (1,0,0)axis2
: the cubic anisotropy axis2, should be a tuple, such as (0,1,0)demag
: include demagnetization or not, should be a boolean, i.e., true or false. By default, demag=false.H
: the external field, should be a tuple or function, i.e.,TupleOrArrayOrFunction
.m0
: the initial magnetization, should be a tuple or function, i.e.,TupleOrArrayOrFunction
.T
: the temperature, should be should beNumberOrArrayOrFunction
.shape
: the shape defines the geometry of the sample, where parameters are configured.
MicroMagnetic.relax
— Functionrelax(sim::AbstractSim; max_steps=10000, stopping_dmdt=0.01, save_data_every=100,
save_m_every=-1, using_time_factor=true)
Relaxes the system using either the LLG
or SD
driver. This function is compatible with both Micromagnetic model and Atomistic spin model.
Arguments:
max_steps::Int
: Maximum number of steps to run the simulation. Default is10000
.stopping_dmdt::Float64
: Primary stopping condition for bothLLG
andSD
drivers. For standard micromagnetic simulations, typical values range from0.01
to1
. InSD
driver mode, where time is not strictly defined, a factor ofγ
is applied to make it comparable to theLLG
driver. For atomistic models using dimensionless units, setusing_time_factor
tofalse
to disable this factor.save_data_every::Int
: Interval for saving overall data such as energies and average magnetization. A negative value disables data saving (e.g.,save_data_every=-1
saves data only at the end of the relaxation).save_m_every::Int
: Interval for saving magnetization data. A negative value disables magnetization saving.using_time_factor::Bool
: Boolean flag to apply a time factor inSD
mode for comparison withLLG
mode. Default istrue
.
Examples:
relax(sim, max_steps=10000, stopping_dmdt=0.1)
Missing docstring for `run_sim
Interfaces – MicroMagnetic
MicroMagnetic.set_Ms
— Methodset_Ms(sim::MicroSim, Ms::NumberOrArrayOrFunction)
Set the saturation magnetization Ms of the studied system. For example,
set_Ms(sim, 8.6e5)
or
function circular_Ms(i,j,k,dx,dy,dz)
if (i-50.5)^2 + (j-50.5)^2 <= 50^2
return 8.6e5
end
return 0.0
end
set_Ms(sim, circular_Ms)
MicroMagnetic.set_Ms
— Methodset_Ms(sim::AbstractSim, shape::Union{CSGNode,Shape}, Ms::Number)
Set the saturation magnetization Ms within the Shape.
MicroMagnetic.add_exch
— Methodadd_exch(sim::MicroSim, A::NumberOrTupleOrArrayOrFunction; name="exch")
Add exchange energy to the system. The exchange energy is definded as
\[ E_\mathrm{ex} = \int_{V} A (\nabla \mathbf{m})^2 \mathrm{d}V\]
Examples:
add_exch(sim, 1e-11)
or
add_exch(sim, (2e-12,5e-12,0))
or
function spatial_A(i,j,k,dx,dy,dz)
if i<10
return 1e-11
else
return 2e-11
end
end
add_exch(sim, spatial_A)
MicroMagnetic.add_dmi
— Methodadd_dmi(sim::MicroSim, D::NumberOrTupleOrArrayOrFunction; name="dmi", type="bulk")
Add DMI to the system. type
could be "bulk", "interfacial" or "D2d".
Examples:
add_dmi(sim, 1e-3, type="interfacial")
or
add_dmi(sim, 1e-3, type="D2d")
add_dmi(sim, (1e-3, 1e-3, 0), type="bulk")
MicroMagnetic.add_demag
— Methodadd_demag(sim::MicroSim; name="demag", Nx=0, Ny=0, Nz=0, fft=true)
Add Demag to the system. Nx
, Ny
and Nz
can be used to describe the macro boundary conditions which means that the given mesh is repeated 2Nx+1
, 2Ny+1 and
2Nz+1times in
x,
yand
z` direction, respectively.
MicroMagnetic.add_exch_int
— Methodadd_exch_int(sim::AbstractSim, J::Float64; k1=1, k2=-1, name="rkky")
Add an RKKY-type exchange for interlayers. The energy of RKKY-type exchange is defined as
\[E_\mathrm{rkky} = - \int_\Gamma J_\mathrm{rkky} \mathbf{m}_{i} \cdot \mathbf{m}_{j} dA\]
where $\Gamma$ is the interface between two layers with magnetizations $\mathbf{m}_{i}$ and $\mathbf{m}_{j}$, $J_\mathrm{rkky}$ is the coupling constant which is related to the spacer layer thickness.
The effective field is given then as
\[\mathbf{H}_i = \frac{1}{\mu_0 M_s} \frac{J_\mathrm{rkky}}{\Delta_z} \mathbf{m}_{j} \]
MicroMagnetic.add_dmi_int
— Methodadd_dmi_int(sim::MicroSimGPU, D::Tuple{Real, Real, Real}; k1=1, k2=-1, name="dmi_int")
Add an interlayer DMI to the system. The energy of interlayer DMI is defined as
\[E_\mathrm{dmi-int} = \int_\Gamma \mathbf{D} \cdot \left(\mathbf{m}_{i} \times \mathbf{m}_{j}\right) dA\]
where $\Gamma$ is the interface between two layers with magnetizations $\mathbf{m}_{i}$ and $\mathbf{m}_{j}$. $\mathbf{D}$ is the effective DMI vector.
The effective field is given
\[\mathbf{H}_i = \frac{1}{\mu_0 M_s \Delta_z} \mathbf{D} \times \mathbf{m}_{j} \]
Interfaces – Atomistic
MicroMagnetic.set_mu_s
— Methodset_mu_s(sim::AtomisticSim, init::NumberOrArrayOrFunction)
Set the magnetic moment mu_s
for the given atomistic system sim
.
Usage Examples
You can set a uniform magnetic moment for the entire system as follows:
set_mu_s(sim, 2 * mu_B)
Alternatively, you can use a function to define spatially varying magnetic moments. For example, to set a circular region with a specific magnetic moment:
function circular_shape(i, j, k, dx, dy, dz)
if (i - 50.5)^2 + (j - 50.5)^2 <= 50^2
return 2 * mu_B
end
return 0.0
end
set_mu_s(sim, circular_shape)
The circular_shape
function here assigns 2 * mu_B
to atoms within a circle of radius 50 (in lattice units), centered at (50.5, 50.5)
, while setting it to 0.0
outside this region.
MicroMagnetic.add_exch
— Methodadd_exch(sim::AtomisticSim, J1::NumberOrArray; name="exch", J2=0, J3=0, J4=0)
Add exchange interactions to the atomistic simulation sim
.
\[\mathcal{H}_\mathrm{ex} = -J \sum_{\langle i, j\rangle} \mathbf{m}_{i} \cdot \mathbf{m}_{j}\]
Parameters
J1::NumberOrArray
: The first nearest-neighbor exchange coupling constant. Can be a single number (applied uniformly) or an array with a length equal to the number of first nearest neighbors.J2, J3, J4
(optional): Exchange coupling constants for second, third, and fourth nearest neighbors, respectively. Each can be a single number or an array with lengths matching their corresponding neighbor counts.name::String
: Optional identifier for the exchange interaction (default is"exch"
).
Usage Examples
To set a uniform exchange coupling for the first nearest neighbors:
add_exch(sim, 0.1*meV)
To set different exchange couplings for multiple neighbor shells:
Js4 = [0.1*meV, 0.2*meV, 0.3*meV, 0.1*meV, 0.2*meV, 0.3*meV]
add_exch(sim, 0.1*meV, J2=0.2meV, J3=0.0, J4=Js4)
Notes
If an array is provided for J1, J2, J3, or J4, its length must match the number of corresponding neighbors defined in the system's mesh. Otherwise, an error will be thrown.
MicroMagnetic.add_exch
— Methodadd_exch(sim::AtomisticSim, Jfun::Function; name="exch")
Add spatial exchange interaction to the system by accepting a function that defines the exchange parameters at different spatial points. The function should accept the indices (i, j, k)
representing the position of the spin on the mesh and return an array. The length of this array should be half the number of neighbors.
For example, for a CubicMesh
, the function should return an array [Jx, Jy, Jz]
, where:
Jx
represents the exchange interaction between the site at(i, j, k)
and the site at(i+1, j, k)
,Jy
represents the exchange interaction between the site at(i, j, k)
and the site at(i, j+1, k)
,Jz
represents the exchange interaction between the site at(i, j, k)
and the site at(i, j, k+1)
.
The array length corresponds to half the number of neighbors, because exchange interactions are considered only in the positive direction (i.e., the interactions with neighboring sites in increasing index directions).
Examples:
function spatial_J(i, j, k)
Jx = i < 10 ? 1meV : 2meV
Jy = 1meV
Jz = 1meV
return [Jx, Jy, Jz] # Returns an array
end
add_exch(sim, spatial_J)
In this example:
- The exchange interaction in the x-direction (
Jx
) depends on the value ofi
. Ifi
is less than 10, it is1meV
, otherwise it is2meV
. - The interactions in the y- and z-directions (
Jy
andJz
) are fixed at1meV
.
MicroMagnetic.add_exch_bq
— Methodadd_exch_bq(sim::AtomisticSim, K::NumberOrArray; name="exch_bq")
Add biquadratic exchange interaction to the atomistic simulation sim
. The biquadratic exchange interaction is defined as
\[\mathcal{H}_\mathrm{ex} = - \sum_{\langle i, j\rangle} K_{ij} (\mathbf{m}_{i} \cdot \mathbf{m}_{j})^2\]
MicroMagnetic.add_dmi
— Methodadd_dmi(sim::AtomisticSim, D::Real; name="dmi", type="bulk")
Add DMI to the system. type
could be "bulk" or "interfacial". The DMI is defined as
\[\mathcal{H}_\mathrm{dmi} = \sum_{\langle i, j\rangle} \mathbf{D}_{i j} \cdot\left(\mathbf{m}_{i} \times \mathbf{m}_{j}\right)\]
where $\mathbf{D}_{i j}$ is the DM vector. For the bulk dmi $\mathbf{D}_{i j} = D \hat{r}_{ij}$ and for interfacial dmi $\mathbf{D}_{i j} = D \hat{r}_{ij} \times \hat{z}$.
Examples:
J = 0.1 * meV
add_dmi(sim, 0.01*J, type="bulk")
MicroMagnetic.add_dmi
— Methodadd_dmi(sim::AtomisticSim, Dij::Array{<:Real, 2}; name="dmi")
Add Dzyaloshinskii-Moriya Interaction (DMI) to the system, defined as:
\[\mathcal{H}_\mathrm{dmi} = \sum_{\langle i, j\rangle} \mathbf{D}_{ij} \cdot \left(\mathbf{m}_{i} \times \mathbf{m}_{j}\right)\]
where $\mathbf{D}_{ij}$ is the DMI vector. The Dij
array should have size (3, mesh.n_ngbs)
and represents the DMI vectors between each pair of neighboring sites. For example, for a CubicMesh
, the Dij
array should be an array of DMI vectors [D1, D2, D3, D4, D5, D6]
, where:
D1
is the DMI vector between the site at(i, j, k)
and the site at(i-1, j, k)
,D2
is the DMI vector between the site at(i, j, k)
and the site at(i+1, j, k)
,D3
is the DMI vector between the site at(i, j, k)
and the site at(i, j-1, k)
,D4
is the DMI vector between the site at(i, j, k)
and the site at(i, j+1, k)
,D5
is the DMI vector between the site at(i, j, k)
and the site at(i, j, k-1)
,D6
is the DMI vector between the site at(i, j, k)
and the site at(i, j, k+1)
.
Examples
# Define bulk DMI for a CubicMesh
D = 0.1 * meV
Dij = hcat([-D, 0, 0], # DMI vector for (i-1, j, k)
[D, 0, 0], # DMI vector for (i+1, j, k)
[0, -D, 0], # DMI vector for (i, j-1, k)
[0, D, 0], # DMI vector for (i, j+1, k)
[0, 0, -D], # DMI vector for (i, j, k-1)
[0, 0, D]) # DMI vector for (i, j, k+1)
add_dmi(sim, Dij)
MicroMagnetic.add_dmi
— Methodadd_dmi(sim::AtomisticSim, Dfun::Function; name="dmi")
Add spatial DMI to the system. The DMI is defined as
\[\mathcal{H}_\mathrm{dmi} = \sum_{\langle i, j\rangle} \mathbf{D}_{i j} \cdot\left(\mathbf{m}_{i} \times \mathbf{m}_{j}\right)\]
where $\mathbf{D}_{i j}$ is the DM vector. The function should accept the indices (i, j, k)
representing the position of the spin on the mesh and return an array. The length of this array should be half the number of neighbors.
For example, for a CubicMesh
, the function should return an array [(D1, D2, D3), (D4, D5, D6), (D4, D5, D6)]
, where:
(D1, D2, D3)
represents the DM vector between the site at(i, j, k)
and the site at(i+1, j, k)
,(D4, D5, D6)
represents the DM vector between the site at(i, j, k)
and the site at(i, j+1, k)
,(D7, D8, D9)
represents the DM vector between the site at(i, j, k)
and the site at(i, j, k+1)
.
The array length corresponds to half the number of neighbors, because exchange interactions are considered only in the positive direction (i.e., the interactions with neighboring sites in increasing index directions).
Examples:
function spatial_bulk_DMI(i, j, k)
Dx = i < 10 ? 0.1meV : 0.2meV
Dy = 0.1meV
Dz = 0.3meV
return [(Dx, 0, 0), (0, Dy, 0), (0, 0, Dz)] # Returns an array
end
add_dmi(sim, spatial_bulk_DMI)
MicroMagnetic.add_demag
— Methodadd_demag(sim::AtomisticSim; name="demag", Nx=0, Ny=0, Nz=0 )
add dipolar interaction into the system.
\[\mathcal{H}_{\mathrm{d}}=-\frac{\mu_0 \mu_s^2}{4 \pi} \sum_{i<j} \frac{3\left(\mathbf{m}_i \cdot \hat{\mathbf{r}}_{i j}\right)\left(\mathbf{m}_j \cdot \hat{\mathbf{r}}_{i j}\right)-\mathbf{m}_i \cdot \mathbf{m}_j}{r_{i j}^3}\]
MicroMagnetic.add_exch_kagome
— Methodadd_exch_kagome(sim::AtomisticSim, Jxy::Number, Jz::Number; name="exch")
Add exchange energy to the system.
MicroMagnetic.add_anis_kagome
— Methodadd_anis_kagome(sim::AtomisticSim, Ku::Float64; ax1=(-0.5,-sqrt(3)/2,0), ax2=(1,0,0), ax3=(-0.5,sqrt(3)/2,0), name="anis")
Add Anisotropy for kagome system, where the energy density is given by
\[E_\mathrm{anis} = - K_{u} (\vec{m} \cdot \hat{u})^2\]
MicroMagnetic.add_anis_tube
— Methodadd_anis_tube(sim::AtomisticSim, Ku::Float64; name="anis")
add anisotropy to the system when the tube mesh is used. The anisotropy axis $u$ is along with the radial direction.
\[E_\mathrm{anis} = - K_{u} (\vec{m} \cdot \hat{u})^2\]
DataSaving
MicroMagnetic.save_vtk
— Functionsave_vtk(mesh::Mesh, shape::Union{CSGNode,Shape}, fname::String)
Save the shape to vtk.
mesh = FDMesh(dx=2e-9, dy=2e-9, dz=2e-9, nx=100, ny=100, nz=50)
t1 = Torus(R = 60e-9, r=20e-9)
save_vtk(mesh, t1, "torus")
save_vtk(sim::AbstractSim, fname::String; fields::Array{String, 1} = String[])
Save magnetization or other fields to vtk.
save_vtk(sim, "m")
MicroMagnetic.save_ovf
— Functionsave_ovf(sim::AbstractSim, fname::String; type::DataType = Float64)
Save spins by ovf, which can be viewed by Muview.
Parameters:
Sim : Sim struct whose spin to be saved.
fname : Save file name.
Optional:
type : Data type of ovf2 file. Can be chosen from Float32, Float64 or String.
For example:
```julia
save_ovf(sim, "ovf_example")
```
Or to specify a certain data type:
```julia
save_ovf(sim, "ovf_example", type = String)
```
MicroMagnetic.read_ovf
— Functionread_ovf(sim, fname)
Initialize sim with an ovf file named of "fname.ovf".
read_ovf(fname)
Load ovf file as OVF2 Type, where spin is stored in OVF2.data
Tools
MicroMagnetic.ovf2vtk
— Functionovf2vtk(ovf_name, vtk_name=nothing; point_data=false, box=noting)
Convert ovf file to vtk format. The data will be saved to points if point_data == true otherwise the data will be saved to cells.
If box is not nothing, it should be a tuple. For instance, box = (nx1, nx2, ny1, ny2, nz1, nz2). In this case, the generated vtk only contains the spins inside the box (including the boundary).
ovf2vtk("my.ovf", "test.vts")
ovf2vtk("my.ovf", point_data=true)
Others
MicroMagnetic.MicroSim
— TypeMicroSim{T<:AbstractFloat} <: AbstractSim
MicroMagnetic.AtomisticSim
— TypeAtomisticSim{T<:AbstractFloat} <: AbstractSim