References

Meshes

JuMag.FDMeshType
FDMesh(;dx=1e-9, dy=1e-9, dz=1e-9, nx=1, ny=1, nz=1, pbc="open")

Create a FDMesh for given parameters. pbc could be any combination of "x", "y" and "z".

source
JuMag.FDMeshGPUType
FDMeshGPU(;dx=1e-9, dy=1e-9, dz=1e-9, nx=1, ny=1, nz=1, pbc="open")

The GPU version of the FDMesh, which is needed for GPU compuation.

source
JuMag.TriangularMeshGPUType

Create a 3d triangular mesh.

The nearest neighbours are indexed as counterclockwise of the given spin:

| 1 2 3 4 5 6 7 8 | |right topright topleft left bottomleft bottomright below above |

and the next-nearest neighbours are

| 1 2 3 4 5 6 7 8 | |topright top topleft bottomleft bottom bottomright below above |

source

DataTypes

JuMag.NumberOrArrayOrFunctionConstant
NumberOrArrayOrFunction

In Micromagnetics, typical parameters such as saturation magnetization Ms and exchange stiffness constant A are constant. However, there are many cases that a spatial Ms is needed. For instance, if the simulated system is a circular disk the Ms in the corners should be set to zero. If the simulated system contains mutiple materials, the exchange constant A should be spatial as well. The Union NumberOrArrayOrFunction is designed to deal with such situations. As indicated from its name, it means that the input parameter could be a number or an array or a function:

  • Number: should be Real.

  • Array: the length of the array should be N where N is the total spin number of the system.

  • Function: the parameter of the function should be (i,j,k,dx,dy,dz) where i,j,k is the cell index and dx,dy,dz is the cellsize. The return value of the function should be a real number. For example,

    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
source
JuMag.NumberOrArrayConstant
NumberOrArray

Similar to Union NumberOrArrayOrFunction, the Union NumberOrArray is designed to deal with cases that a number or an array is needed.

source
JuMag.TupleOrArrayOrFunctionConstant
TupleOrArrayOrFunction

Similar to NumberOrArrayOrFunction, TupleOrArrayOrFunction means that the input parameter could be a tuple or an array or a function:

  • Tuple: should be Real with length 3. For example, (0,0,1e5).
  • Array: the length of the array should be 3N where N is the total spin number of the system.
  • Function: the parameter of the function should be (i,j,k,dx,dy,dz) and the return value should be a tuple with length 3. For example,
    function uniform_m0(i,j,k,dx,dy,dz)
        return (0,0,1)
    end
source

Interfaces

JuMag.SimFunction
Sim(mesh::Mesh; driver="LLG", name="dyn", integrator="DormandPrince")

Create a simulation instance for given mesh.

source
Sim(mesh::MeshGPU; driver="LLG", name="dyn", integrator="DormandPrince")

Create a simulation instance for GPU mesh.

source
JuMag.set_MsFunction
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)
source
set_Ms(sim::MicroSimGPU, Ms::NumberOrArrayOrFunction)
source
JuMag.set_Ms_cylindricalFunction
set_Ms_cylindrical(sim::MicroSim, Ms::Number, axis=ex)

Set the saturation magnetization Ms of the studied system to a cylindrical shape, axis could be JuMag.ex, JuMag.ey or JuMag.ez.

source
set_Ms_cylindrical(sim::MicroSimGPU, Ms::Number, axis=ex)
source
JuMag.set_mu_sFunction
set_mu_s(sim::AtomicSimGPU, Ms::NumberOrArrayOrFunction)

Set magnetic moment mu_s of the studied system. For example,

   set_mu_s(sim, 8.6e5)

or

function circular_shape(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_mu_s(sim, circular_shape)
source
JuMag.init_m0Function
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)
source
init_m0(sim::AbstractSimGPU, m0::TupleOrArrayOrFunction; norm=true)

Set the initial magnetization of the system. If norm=false the magnetization array will be not normalised.

source
JuMag.init_m0_skyrmionFunction
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 –> clockwise; -1 –> counter-clockwise

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")
source
JuMag.add_exchFunction
add_exch(sim::AbstractSim, A::NumberOrArrayOrFunction; name="exch")

Add exchange energy to the system.

source
add_exch(sim::AtomicSimGPU, J::Array; name="exch")

Add exchange energy to the system. The length of J should be equal to the length of neigbours.

source
add_exch(sim::AtomicSimGPU, J::Number; name="exch")

Add exchange energy to the system.

source
JuMag.add_anisFunction
add_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} (\vec{m} \cdot \hat{u})^2\]
source
JuMag.add_dmiFunction
add_dmi(sim::AbstractSim, D::Tuple{Real, Real, Real}; name="dmi")

Add DMI to the system. Example:

   add_dmi(sim, (1e-3, 1e-3, 0))
source
add_dmi(sim::AbstractSim, D::Real; name="dmi", type="bulk")

Add DMI to the system. type could be "bulk" or "interfacial" Examples:

   add_dmi(sim, 1e-3, type="interfacial")

or

   add_dmi(sim, 1e-3, type="bulk")
source
JuMag.add_demagFunction
add_demag(sim::MicroSim; name="demag", Nx=0, Ny=0, Nz=0)

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.

source
JuMag.add_zeemanFunction
add_zeeman(sim::AbstractSim, H0::TupleOrArrayOrFunction; name="zeeman")

Add a static Zeeman energy to the simulation.

source
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)
source
add_zeeman(sim::AbstractSimGPU, H0::TupleOrArrayOrFunction; name="zeeman")

Add a static Zeeman energy to the simulation.

source
JuMag.add_exch_vectorFunction
add_exch_vector(sim::AbstractSim, A::TupleOrArrayOrFunction; name="exch")

Add a vector form exchange energy to the system. The exchange constant of 3 directions can be different. For example:

add_exc_vector(sim, (2e-12,5e-12,0))
source
JuMag.add_exch_kagomeFunction
add_exch_kagome(sim::AtomicSimGPU, Jxy::Number, Jz::Number; name="exch")

Add exchange energy to the system.

source
JuMag.add_anis_kagomeFunction
add_anis_kagome(sim::AtomicSimGPU, 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\]
source
add_anis_kagome(sim::MonteCarlo, Ku::Float64)

Add Anisotropy for kagome system, where the energy density is given by

\[ E_\mathrm{anis} = - K_{u} (\vec{m} \cdot \hat{u})^2\]

and u is one of ax1=(-0.5,-sqrt(3)/2,0), ax2=(1,0,0) and ax3=(-0.5,sqrt(3)/2,0).

source
JuMag.update_zeemanFunction
update_zeeman(sim::AbstractSim, H0::Tuple; 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
source
JuMag.update_anisFunction
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 K2
source
JuMag.relaxFunction
relax(sim::AbstractSim; maxsteps=10000, stopping_dmdt=0.01, save_m_every = 10, save_ovf_every=-1, ovf_format = "binary", ovf_folder="ovfs", save_vtk_every=-1, vtk_folder="vtks",fields::Array{String, 1} = String[])

Relax the system using LLG or SD driver. The stop condition is determined by stopping_dmdt(both for LLG and SD drivers). Spins can be stored in ovfs or vtks. ovf format can be chosen in "binary"(float64),"binary8"(float64), "binary4"(float32), "text"

Fields can be stored in vtks by:

relax(sim, save_vtk_every = 10, fields = ["demag", "exch", "anis"])
source

DataSaving

JuMag.save_mFunction
save_m(sim::AbstractSim, fname::String; vtk::Bool = false, npy::Bool = false, vtk_folder = "vtks")

If vtk = true, save spins to dir ./vtks/ in vtk format;
If npy = true, save spins to current directory in npy format;

Example:

```julia
    save_m(sim, sim.name, vtk = true, npy= true)
```
source
JuMag.save_vtkFunction
save_vtk(sim::AbstractSim, fname::String; fields::Array{String, 1} = String[])

Save magnetization or other fields to vtk.

    save_vtk(sim, "m")
source
JuMag.save_ovfFunction

save_ovf(sim::AbstractSim, fname::String; dataformat::String = "text")

Save spins in format of ovf, which can be viewed by Muview. Dataformat could be "text", "binary4" or "binary8".

For example:

```julia
    save_ovf(sim, "m0")
```
source
JuMag.read_ovfFunction
read_ovf(sim, fname)

Initialize sim with an ovf file named of "fname.ovf".
source
read_ovf(fname)

Load ovf file as OVF2 Type, where spin is stored in OVF2.data
source

Tools

JuMag.ovf2vtkFunction
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)
source