API for users
Meshes – MicroMagnetic
MicroMagnetic.FDMesh Method
FDMesh(; 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.
sourceMeshes – Atomistic
MicroMagnetic.CubicMesh Method
CubicMesh(;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.
sourceMicroMagnetic.TriangularMesh Method
TriangularMesh(;dx=1e-9, dz=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 Method
CylindricalTubeMesh(;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 |
sourceShapes
MicroMagnetic.Cylinder Type
Cylinder(;center = (0.0, 0.0, 0.0), radius = 1.0, height= Inf, normal = (0.0, 0.0, 1.0))Create a Cylinder shape.
sourceMicroMagnetic.Box Type
Box(;center::Tuple=(0.0, 0.0, 0.0), sides::Tuple=(1.0,1.0,1.0), theta=0.0)Create a Box shape.
sourceMicroMagnetic.Torus Type
Torus(; center::Tuple=(0.0, 0.0, 0.0), R=1.0, r=0.2)Create a Torus shape.
sourceInterfaces
MicroMagnetic.sim_with Function
sim_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.AorJ: 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 (trueorfalse). 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 theDynamicstask.dt: The time interval of each step, so the total simulation time issteps * dtfor theDynamicstask.max_steps::Int: Maximum number of steps to run the simulation for theRelaxtask. Default is10000.saver_item: ASaverIteminstance or a list ofSaverIteminstances. 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 bothLLGandSDdrivers. For standard micromagnetic simulations, typical values range from0.01to1. InSDdriver mode, where time is not strictly defined, a factor ofγis applied to make it comparable to theLLGdriver. For atomistic models using dimensionless units, setusing_time_factortofalseto disable this factor.relax_data_interval::Int: Interval for saving overall data such as energies and average magnetization during aRelaxtask. A negative value disables data saving (e.g.,relax_data_interval = -1saves 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 theDynamicstask. Set totrueto enable, orfalseto disable.relax_m_interval::Int: Interval for saving magnetization data during aRelaxtask. A negative value disables magnetization saving.dynamic_m_interval::Int: Interval for saving magnetization data during aDynamicstask. A negative value disables magnetization saving.using_time_factor::Bool: Boolean flag to apply a time factor inSDmode for comparison withLLGmode. 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, anddrivercan be defined as arrays using the_sor_sweepsuffix. 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
argsparameter can be either aNamedTupleor aDict, providing flexibility in how you organize and pass the simulation parameters.Driver Selection: The
driverparameter (ordriver_sfor 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_dmdtparameter 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.01to1. For atomistic models, theusing_time_factorflag can be set tofalseto disable time scaling.Data Saving: The
relax_m_intervalanddynamic_m_intervalparameters control how frequently magnetization data is saved duringRelaxandDynamicstasks, respectively. Use negative values to disable data saving for these tasks.
MicroMagnetic.set_backend Function
set_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 CUDATo set the backend to use the CPU/CUDA:
set_backend("cpu")
set_backend("cuda")MicroMagnetic.set_precision Function
set_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 Function
set_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 loggingMicroMagnetic.Sim Function
Sim(mesh::Mesh; driver="LLG", name="dyn", integrator="DormandPrince", save_data=true)Create a simulation instance for the given mesh with specified driver and integrator.
Arguments
mesh::Mesh: The computational mesh defining the simulation geometry. Can be:FDMesh: Finite difference mesh for micromagneticsFEMesh: Finite element mesh for micromagneticsAtomisticMesh: Creates atomistic simulation
Keyword Arguments
driver::String="LLG": The simulation driver type. Options:"None": No driver (static simulation)"SD": Energy minimization (Steepest Descent)"LLG": Landau-Lifshitz-Gilbert equation"InertialLLG": Inertial LLG Equation"SpatialLLG": Spatial LLG equation allowing spatial damping constant."LLG_STT": LLG with spin transfer torque"LLG_CPP": LLG with CPP spin transfer torque
name::String="dyn": Name identifier for the simulationintegrator::String="DormandPrince": Time integration method. Options:Fixed step methods:
"Heun","RungeKutta","RungeKuttaCayley"Adaptive step methods:
"DormandPrince"(DOPRI54),"BS23","CashKarp54","Fehlberg54"(RKF54)Cayley-transform methods:
"DormandPrinceCayley","RungeKuttaCayley"
save_data::Bool=true: Whether to enable data saving during simulation
Returns
MicroSim{Float}: For FDMesh simulationsMicroSimFE{Float}: For FEMesh simulationsAtomisticSim{Float}: For AtomisticMesh types
Examples
# Create LLG simulation with Dormand-Prince integrator
mesh = FDMesh(nx=100, ny=100, nz=1)
sim = Sim(mesh, driver="LLG", integrator="DormandPrince")
# Create energy minimization simulation
sim_sd = Sim(mesh, driver="SD")
# Create STT simulation with adaptive stepping
sim_stt = Sim(mesh, driver="LLG_STT", integrator="BS23")MicroMagnetic.NEB Type
NEB(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 Function
set_driver(sim::AbstractSim; driver="LLG", integrator="DormandPrince", args...)Set the driver of the simulation, can be used to switch the driver.
sourceMicroMagnetic.set_alpha Method
set_alpha(sim::AbstractSim, alpha::ArrayOrFunction)Set the damping parameter α for a simulation with SpatialLLG driver.
This function updates the spatially varying damping coefficient array in the simulation's driver. The α parameter can be specified either as a uniform value (via a function) or as a spatially varying array.
Arguments
sim::AbstractSim: The simulation object containing the driveralpha::ArrayOrFunction: The damping parameter specification, which can be:A function
f(i, j, k, dx, dy, dz) -> Float64that returns the α value at each spatial coordinateAn array of Float64 values with length equal to
sim.n_total
Examples
# Set uniform α value of 0.1 using a function
set_alpha(sim, (i, j, k, dx, dy, dz) -> 0.1)
# Set spatially varying α using an array
alpha_array = rand(sim.n_total) * 0.2 + 0.01 # Random values between 0.01 and 0.21
set_alpha(sim, alpha_array)
# Set spatial α with a function
function spatial_alpha(i, j, k, dx, dy, dz)
if i < 10
return 0.01
else
return 0.99
end
end
set_alpha(sim, spatial_alpha)MicroMagnetic.init_m0 Function
init_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 Function
init_m0_random(sim::MicroSim)Set the initial magnetization with random direction.
sourceMicroMagnetic.init_m0_skyrmion Function
init_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 Function
add_zeeman(sim::AbstractSim, H0::TupleOrArrayOrFunction; name="zeeman")Add a static Zeeman energy to the simulation.
sourceadd_zeeman(sim::AbstractSim, H0::TupleOrArrayOrFunction, ft::Function; name="timezeeman")Add a time-varying Zeeman interaction to the simulation system.
Arguments
sim::AbstractSim: The simulation object to which the Zeeman interaction will be added.H0::TupleOrArrayOrFunction: The spatial field configuration, which can be:A tuple of 3 numbers representing a uniform field in x, y, z directions
An array representing a spatially varying field
A function that returns the field at specific spatial coordinates
ft::Function: The time-dependent function that modulates the field amplitude.name::String: (Optional) A name identifier for this interaction. Defaults to "timezeeman".
Time Function Specification
The time function ft should accept a time argument t and return either:
Scalar return: A single value that will be applied uniformly to all three field components
Tuple return: A 3-element tuple
(fx, fy, fz)specifying individual scaling factors for x, y, z components
The function is evaluated at each time step to modulate the applied field.
Examples
Example 1: Uniform rotating field with vector time function
function time_fun(t)
w = 2*pi*2.0e9 # 2 GHz frequency
return (sin(w*t), cos(w*t), 0) # Rotating in xy-plane
end
function spatial_H(i, j, k, dx, dy, dz)
H = 1e3 # 1000 A/m base field strength
if i <= 2
return (H, H, 0) # Apply field only in first two cells
end
return (0, 0, 0) # Zero field elsewhere
end
add_zeeman(sim, spatial_H, time_fun)Example 2: Uniformly oscillating field with scalar time function
function time_fun_scalar(t)
w = 2*pi*1.0e9 # 1 GHz frequency
return sin(w*t) # Same modulation for all components
end
# Apply uniform field in z-direction
add_zeeman(sim, (0, 0, 1e3), time_fun_scalar)Example 3: Complex modulated field
function complex_time_fun(t)
w1 = 2*pi*1.0e9
w2 = 2*pi*0.5e9
# Different modulation for each component
fx = sin(w1*t) * cos(w2*t)
fy = 0.5 * (1 + sin(w1*t))
fz = exp(-t/1e-9) # Exponential decay in z-component
return (fx, fy, fz)
end
# Spatially uniform initial field
add_zeeman(sim, (1e3, 2e3, 0.5e3), complex_time_fun, name="modulated_zeeman")Field Initialization
The H0 parameter defines the spatial distribution of the field:
For uniform fields: Use a 3-element tuple
(Hx, Hy, Hz)For spatially varying fields: Provide a function with signature
(i, j, k, dx, dy, dz) -> (Hx, Hy, Hz)wherei, j, kare grid indices anddx, dy, dzare spatial steps
Energy Calculation
The Zeeman energy density is calculated as: E = -μ₀ * M ⋅ H where M is the magnetization and H is the applied field.
sourceMicroMagnetic.add_anis Function
add_anis(sim::AbstractSim, Ku::NumberOrArrayOrFunction; axis::TupleOrArrayOrFunction, name="anis")Add Anisotropy to the system, where the energy density is given by
add_anis(sim::MicroSimFE, Ku::NumberOrArrayOrFunction; axis=(0,0,1), name="anis")Add Anisotropy to the system, where the energy density is given by
MicroMagnetic.add_cubic_anis Function
add_cubic_anis(sim::AbstractSim, Kc::NumberOrArrayOrFunction; 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.
Example
add_cubic_anis(sim, 1e3, (1, 1, 0), (1, -1, 0))MicroMagnetic.add_hex_anis Function
add_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:
Example
add_hex_anis(sim, K1=1e3, K2=0, K3=1e2)MicroMagnetic.update_zeeman Function
update_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/mMicroMagnetic.update_anis Function
update_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 K2MicroMagnetic.add_stt Function
add_stt(sim::AbstractSim; model::Symbol=:slonczewski, name="stt", kwargs...)Add spin transfer torque (STT) to the Landau-Lifshitz-Gilbert equation using various theoretical models.
Spin transfer torque is a fundamental phenomenon in spintronics where spin-polarized electric currents exert torques on magnetic moments, enabling precise manipulation of magnetization dynamics in nanoscale devices.
Supported Models
Zhang-Li Model (model=:zhang_li)
The Zhang-Li model describes spin transfer torque for in-plane current flows in continuous magnetic materials. This model is particularly suitable for describing current-induced domain wall motion, magnetization dynamics in nanowires, and spatially non-uniform current distributions.
The modified LLG equation becomes:
where
Parameters
Required:
J::TupleOrArrayOrFunction: Current density specification (A/m²)Constant vector:
TupleorVectorspecifying uniform current flowSpatially varying:
Functionthat returns current density at position(x,y,z)
Coefficient specification (choose one approach):
Physical parameters approach:
P::Real: Spin polarization (0 ≤ P ≤ 1)Ms::Real: Saturation magnetization (A/m)
Direct coefficient approach:
b::Real: Pre-calculated Zhang-Li torque coefficient
Optional:
xi::Real=0.0: Non-adiabatic parameter (default: 0.0)ft::Function: Time-dependent function modulating current density amplitude
Examples
# Basic usage with uniform current vector
add_stt(sim, model=:zhang_li, P=0.5, Ms=8e5, xi=0.05, J=(1e12, 0, 0))
# Using pre-calculated coefficient
add_stt(sim, model=:zhang_li, b=72.5, xi=0.1, J=(1, 0, 0))
# Time-dependent current pulse with spatial variation
function current_pulse(t)
return t < 1e-9 ? 1.0 : 0.0 # 1 ns pulse
end
function spatial_current(i, j, k, dx, dy, dz)
current = k > 5 ? 1e12 : 2e11
return (current, 0.0, 0.0)
end
add_stt(sim, model=:zhang_li, P=0.5, Ms=8e5, J=spatial_current, ft=current_pulse)Slonczewski Model (model=:slonczewski)
The Slonczewski model describes spin transfer torque for perpendicular current flows in magnetic tunnel junctions (MTJs) and spin valves. The extended LLG equation is:
where
where
Parameters
Required:
J::NumberOrArrayOrFunction: Current density magnitude (A/m²)tf::Real: Free layer thickness (m)Ms::Real: Saturation magnetization (A/m)p::Union{Tuple, Vector}: Polarization direction vectorP::Real: Spin polarization (0 ≤ P ≤ 1)
Optional:
xi::Real=0.0: Secondary spin-torque parameter (default: 0.0)Lambda::Real=2.0: Slonczewski parameter (default: 2.0)ft::Function: Time-dependent function modulating current density amplitude
Examples
# Const current density
add_stt(sim, model=:slonczewski, J=1e12, tf=2e-9, Ms=8e5, p=(0,0,1), P=0.7)
# With secondary spin-torque parameter
add_stt(sim, model=:slonczewski, J=8e11, tf=1.5e-9, Ms=6e5, p=(1,0,0), P=0.5, xi=0.05)
# Time-dependent current pulse with Gaussian current profile
function current_pulse(t)
return t < 1e-9 ? 1.0 : 0.0 # 1 ns pulse
end
function gaussian_current(i, j, k, dx, dy, dz)
x, y = (i-100)*dx, (j-100)*dy
sigma = 50e-9 # 50 nm width
r_sq = x^2 + y^2
return 1e12 * exp(-r_sq/(2*sigma^2))
end
add_stt(sim, model=:slonczewski, J=gaussian_current, tf=2e-9, Ms=8e5, p=(0,0,1), P=0.4, ft=current_pulse)MicroMagnetic.add_sot Function
add_sot(sim::AbstractSim, aj::NumberOrArrayOrFunction, bj::Number, p::Tuple{Real,Real,Real}; name="sot")Add damping-like and field-like torque to LLG equation, which is useful to model spin-orbit torques (SOT) or Slonczewski torque with constant
The equivalent effective field is
MicroMagnetic.add_torque Function
add_torque(sim::AbstractSim, fun::Function; name="torque")Add a torque
The modified LLG equation becomes:
where
Arguments
sim::AbstractSim: The simulation object to which the torque will be added.fun::Function: A function that defines the torque term. The function should have the signaturefun(y, m, t)where:y: The torquem: The magnetization vector mt: The current simulation time
Examples
function torque_fun(y, m, t)
y .= 0
y[1:3:end] .= m[1:3:end] .* sin(t) # x-component
y[2:3:end] .= 0.2 * cos(t) # y-component
end
add_torque(sim, torque_fun, name="test")MicroMagnetic.add_sahe_torque Function
add_sahe_torque(sim::AbstractSim, sigma_s::TupleOrArrayOrFunction, sigma_sa::TupleOrArrayOrFunction; a1=(1,0,0), a2=(0,1,0), beta=0, name="she")Add an effective field to represent the torque due to spin anomalous Hall effect, which is given by
The equivalent effective field is
where
Example
add_sahe_torque(sim, (0.1, 0.2, 0.3), (0.3, 0.4, 0.5), beta=0.1)MicroMagnetic.add_thermal_noise Function
add_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
and
For the atomistic model, the thermal noise is defined as
where
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.create_sim Function
create_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,uyoruz: 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, mu_s=2*mu_BAorJ: 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 Function
relax(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 bothLLGandSDdrivers. For standard micromagnetic simulations, typical values range from0.01to1. InSDdriver mode, where time is not strictly defined, a factor ofγis applied to make it comparable to theLLGdriver. For atomistic models using dimensionless units, setusing_time_factortofalseto 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=-1saves 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 inSDmode for comparison withLLGmode. Default istrue.
Examples:
relax(sim, max_steps=10000, stopping_dmdt=0.1)MicroMagnetic.run_sim Function
run_sim(sim::AbstractSim; steps=10, dt=1e-12, save_data=true, save_m_every=1, saver_item=nothing, call_back=nothing)Run the simulation for a total duration of steps * dt seconds.
steps: The total number of simulation steps. The total simulation time will besteps * dt.dt: The time step for each simulation step (in seconds). It controls the time interval between successive steps.save_data: A boolean flag to control whether the overall simulation data (such as total energy, exchange energy, and average magnetization) is saved at each step. The default istrue.save_m_every: Specifies how often to save the magnetization configuration. For example, ifsave_m_every = 1, the magnetization is saved at every step. A negative value will disable magnetization saving entirely.saver_item: ASaverIteminstance or a list ofSaverIteminstances. 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.
Example Usage
To compute a custom quantity, such as the guiding center of the magnetization, and save it in the output table, you can define a SaverItem as shown below:
run_sim(sim, saver_item=SaverItem(("Rx", "Ry"), ("<m>", "<m>"), compute_guiding_center))In this example:
("Rx", "Ry")are the labels for the data in the output.("<m>", "<m>")are the units of guiding center.compute_guiding_centeris the function that computes the guiding center at each step.
This setup will save the computed guiding center to the simulation output, in addition to the default data like energies and average magnetization.
sourceInterfaces – MicroMagnetic
MicroMagnetic.set_Ms Method
set_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 Method
set_Ms(sim::AbstractSim, shape::Union{CSGNode,Shape}, Ms::Number)Set the saturation magnetization Ms within the Shape.
sourceMicroMagnetic.add_exch Method
add_exch(sim::MicroSim, A::NumberOrTupleOrArrayOrFunction; name="exch")Add exchange energy to the system. The exchange energy is definded as
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 Method
add_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 Method
add_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 and2Nz+1times inx,yandz` direction, respectively.
MicroMagnetic.add_exch_int Method
add_exch_int(sim::AbstractSim, J::NumberOrArrayOrFunction; k1=1, k2=-1, name="rkky")Add an RKKY-type exchange for interlayers. The energy of RKKY-type exchange is defined as
where
The effective field is given then as
MicroMagnetic.add_dmi_int Method
add_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
where
The effective field is given
MicroMagnetic.voronoi Method
voronoi(mesh; min_dist=20, seed=123456, threshold=nothing) -> grain_ids, gb_mask, pointsGenerate a Voronoi tessellation on a 2D grid with grain boundaries detection using Poisson disk sampling.
Arguments
mesh:FDMeshobject containing grid informationmin_dist: Minimum distance between Voronoi seeds in nanometers (default: 20)seed: Random number generator seed for reproducible results (default: 123456)threshold: Optional threshold for continuous boundary detection. If provided, uses distance-based boundary detection instead of discrete neighbor comparison.
Returns
grain_ids: Matrix of integers where each element represents the ID of the Voronoi cell it belongs togb_mask: Boolean matrix marking grain boundary locations (true = boundary)points: Array of Voronoi seed points used for the tessellation
Example
mesh = FDMesh(; dx=2e-9, dy=2e-9, dz=2e-9, nx=100, ny=100, nz=2)
# Discrete boundary detection (default)
grain_ids, boundaries, seeds = voronoi(mesh; min_dist=20, seed=1000)
# Continuous boundary detection with custom threshold
grain_ids, boundaries, seeds = voronoi(mesh; min_dist=20, seed=1000, threshold=0.2)Interfaces – Atomistic
MicroMagnetic.set_mu_s Method
set_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 Method
add_exch(sim::AtomisticSim, J1::NumberOrArray; name="exch", J2=0, J3=0, J4=0)Add exchange interactions to the atomistic simulation sim.
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.
sourceMicroMagnetic.add_exch Method
add_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:
Jxrepresents the exchange interaction between the site at(i, j, k)and the site at(i+1, j, k),Jyrepresents the exchange interaction between the site at(i, j, k)and the site at(i, j+1, k),Jzrepresents 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. Ifiis less than 10, it is1meV, otherwise it is2meV.The interactions in the y- and z-directions (
JyandJz) are fixed at1meV.
MicroMagnetic.add_exch_bq Method
add_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
MicroMagnetic.add_dmi Method
add_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
where
Examples:
J = 0.1 * meV
add_dmi(sim, 0.01*J, type="bulk")MicroMagnetic.add_dmi Method
add_dmi(sim::AtomisticSim, Dij::Array{<:Real, 2}; name="dmi")Add Dzyaloshinskii-Moriya Interaction (DMI) to the system, defined as:
where 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:
D1is the DMI vector between the site at(i, j, k)and the site at(i-1, j, k),D2is the DMI vector between the site at(i, j, k)and the site at(i+1, j, k),D3is the DMI vector between the site at(i, j, k)and the site at(i, j-1, k),D4is the DMI vector between the site at(i, j, k)and the site at(i, j+1, k),D5is the DMI vector between the site at(i, j, k)and the site at(i, j, k-1),D6is 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 Method
add_dmi(sim::AtomisticSim, Dfun::Function; name="dmi")Add spatial DMI to the system. The DMI is defined as
where (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 Method
add_demag(sim::AtomisticSim; name="demag", Nx=0, Ny=0, Nz=0 )add dipolar interaction into the system.
MicroMagnetic.add_exch_kagome Method
add_exch_kagome(sim::AtomisticSim, Jxy::Number, Jz::Number; name="exch")Add exchange energy to the system.
sourceMicroMagnetic.add_anis_kagome Method
add_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
MicroMagnetic.add_anis_tube Method
add_anis_tube(sim::AtomisticSim, Ku::Float64; name="anis")add anisotropy to the system when the tube mesh is used. The anisotropy axis
DataSaving
MicroMagnetic.save_vtk Function
save_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")save_vtk(sim::MicroSimFEM, fname::String; fields::Array{String, 1} = String[])Save magnetization or other fields to vtk.
save_vtk(sim, "m")MicroMagnetic.save_ovf Function
save_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 Function
read_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.dataTools/Visualization
MicroMagnetic.voronoi Function
voronoi(mesh; min_dist=20, seed=123456, threshold=nothing) -> grain_ids, gb_mask, pointsGenerate a Voronoi tessellation on a 2D grid with grain boundaries detection using Poisson disk sampling.
Arguments
mesh:FDMeshobject containing grid informationmin_dist: Minimum distance between Voronoi seeds in nanometers (default: 20)seed: Random number generator seed for reproducible results (default: 123456)threshold: Optional threshold for continuous boundary detection. If provided, uses distance-based boundary detection instead of discrete neighbor comparison.
Returns
grain_ids: Matrix of integers where each element represents the ID of the Voronoi cell it belongs togb_mask: Boolean matrix marking grain boundary locations (true = boundary)points: Array of Voronoi seed points used for the tessellation
Example
mesh = FDMesh(; dx=2e-9, dy=2e-9, dz=2e-9, nx=100, ny=100, nz=2)
# Discrete boundary detection (default)
grain_ids, boundaries, seeds = voronoi(mesh; min_dist=20, seed=1000)
# Continuous boundary detection with custom threshold
grain_ids, boundaries, seeds = voronoi(mesh; min_dist=20, seed=1000, threshold=0.2)MicroMagnetic.ovf2vtk Function
ovf2vtk(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)MicroMagnetic.plot_ts Function
plot_ts(filename::String, keys::Vector{String};
x_key::String="time", x_unit::Real=1e9, y_units::Vector=[],
size=(400, 280), title="", xlabel="", ylabel="",
plot_type=:scatterlines, markersize=6, legend_position=:rt,
colors=nothing, linestyles=nothing, transparency=false)General time series plotting function
Arguments
filename: Data file namekeys: Data column keys to plotx_key: Time column key, default is "time"x_unit: Time unit conversion factor, default is 1e9 (seconds to nanoseconds)y_units: Y-axis data unit conversion factors for each key, same order askeysplot_type: Plot type,:lines,:scatter,:scatterlines,:linesmarkerslegend_position: Legend position,:rt(top right),:lt(top left),:rb(bottom right),:lb(bottom left)colors: Custom color sequencelinestyles: Custom line style sequencetransparency: Whether to use transparent background
Examples
# Plot magnetization components
plot_ts("std4_llg.txt", ["m_x", "m_y", "m_z"];
xlabel="Time (ns)", ylabel="m", plot_type=:scatter)
# Plot vortex center with unit conversion
plot_ts("std5_llg.txt", ["cx", "cy"];
xlabel="Time (ns)", ylabel="Vortex center (nm)",
y_units=[1e9, 1e9], plot_type=:scatterlines)
# Custom colors and line styles
plot_ts("data.txt", ["A", "B", "C"];
colors=[:red, :blue, :green],
linestyles=[:solid, :dash, :dot])MicroMagnetic.plot_m Function
MicroMagnetic.plot_m(spin; dx=1.0, dy=1.0, k=1, component='z', arrows=(-1, -1), figsize=(500, -1), fig=nothing, ax=nothing, kwargs...)Create a plot for the given magnetization.
Arguments
spin::Array{Float64, 4}: An array with dimensions (3, nx, ny, nz) representing the magnetization.
Keyword Arguments
dx::Float64: The spacing in the x-direction (default: 1.0).dy::Float64: The spacing in the y-direction (default: 1.0).k::Int: The layer index to plot (starting from 1) (default: 1).component::Char: The magnetization component to plot ('x', 'y', or 'z') (default: 'z').arrows::Tuple{Int, Int}: The number of arrows to plot, specified as a tuple. By default, arrows=(-1, -1), which auto-scales the number of arrows. If arrows is set tonothing, no arrows will be shown.figsize::Tuple{Int, Int}: The size of the figure, specified as a tuple (width, height). For example, figsize=(500, 400) or figsize=(500, -1) where -1 means auto-scaled height (default: (500, -1)).fig: An existing figure to plot on. If nothing, a new figure is created (default: nothing).ax: An existing axis to plot on. If nothing, a new axis is created (default: nothing).kwargs...: Additional keyword arguments that are passed toheatmap!.
Supported heatmap! Arguments
alpha: Transparency level of the heatmap.colormap: A list or function that defines the color mapping for the heatmap Colors.colorrange: The color range of the heatmap.interpolate: Whether to interpolate between values.colorscale: A scale factor for the color.transparency: Whether to apply transparency to the heatmap.
For more details on the supported arguments, please refer to the Makie documentation: Makie Heatmap Reference.
Returns
fig: The figure containing the plot.
Examples
spin = randn(3, 10, 10, 5) # Example spin data
# Creates a plot with default settings
MicroMagnetic.plot_m(spin, colorrange=[-1, 1])
# Creates a plot for the x-component of the second layer with custom settings
MicroMagnetic.plot_m(spin, dx=0.5, dy=0.5, k=2, component='x', arrows=(5, 5), figsize=(600, 400))MicroMagnetic.plot_m(sim::MicroMagnetic.AbstractSim; kwargs...)Create a plot for the given magnetization in a simulation object.
Arguments
sim::MicroMagnetic.AbstractSim: A simulation object containing the magnetization data and mesh information.
Keyword Arguments
This function forwards all keyword arguments to MicroMagnetic.plot_m. Refer to MicroMagnetic.plot_m for detailed descriptions of the keyword arguments.
Returns
fig: The figure containing the plot.
Examples
plot_m(sim, k=2, component='x', arrows=(20, 20), figsize=(600, 400))MicroMagnetic.ovf2png Function
ovf2png(ovf_name, output=nothing; k=1, arrows=(-1, -1), figsize=(500, -1))Create a png from the given ovf file.
kindicates the layer index (starting from 1)arrowsis the number of arrows, should be a Tuple of integers. By default, arrows=(-1, -1).figsizeshould be a Tuple of integers, for example, figsize=(500, 400) or figsize=(500, -1).
MicroMagnetic.ovf2movie Function
ovf2movie(folder; framerate=12, output=nothing, kwargs...)Create a moive from the given folder.
Keyword Arguments
This function forwards all keyword arguments to MicroMagnetic.plot_m. Refer to MicroMagnetic.plot_m for detailed descriptions of the keyword arguments.
output is the filename of the video and the support formats are 'mp4', 'avi' and 'gif'.
Examples
ovf2movie("std4_LLG"; output="std4.gif", component='x');