The mumax3 input syntax is a subset of Go's syntax, somewhat similar to C. It is case-independent however, so msat is the same as Msat or MSAT.
:=
. Variables have a fixed type, inferred from the declaration's right-hand-side. Assigning to existing variables is done using =
. E.g.:
i := 7 // defines a new variable i, type automatically detected to be int
print(i) // now we can use i
i = 5 // assign new value, don't use ':=' (attempt to re-declare)
str := "hello" // defines str, type automatically is string
//str = 1 // would fail, cannot assign int to string
x := pi*(3+4)/5
x = pow(x, 3)
x++
y := abs(cbrt(cosh(erf(erfc(gamma(J0(Y0(2))))))))
for i:=0; i<10; i++{
print(i)
}
RunWhile(func()bool)
, or all input parameters). In that case, and only in this case, the argument is implicitly converted to a function, which is re-evaluated each time it's needed. E.g.:
value := sin(pi*t) // value is a float64, RHS evaluated only once
Msat = value // time-independent Msat
versus:
Msat = sin(pi*t) // RHS converted to function, re-evaluted every time
Msat
has the method SetRegion(int, float)
to set the value of the material parameter in a certain region:
Msat.SetRegion(1, 800e3) // Set Msat=520e3 in region 1
Nx := 128
Ny := 64
Nz := 2
sizeX := 500e-9
sizeY := 250e-9
sizeZ := 10e-9
SetGridSize(Nx, Ny, Nz)
SetCellSize(sizeX/Nx, sizeY/Ny, sizeZ/Nz)
SetPBC(5, 0, 0) // 5 extra images on left and right sides.
SetGridSize(128, 64, 1)
SetCellSize(5e-9, 5e-9, 5e-9)
Setting a nonzero PBC value in a direction enables wrap-around in that direction. The precise value passed determines how many repetitions are seen by the demag field. E.g., in the above example the demag field behaves as if 5 repetitions are present to the left and to the right side. Choosing a large number may cause long initialization time.
geometryShape := cylinder(400e-9, 20e-9).RotX(45*pi/180).Transl(1e-6,0,0)
SetGeom(geometryShape)
SetCellSize(float64, float64, float64)
Sets the X,Y,Z cell size in meters
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
SetGridSize(int, int, int)
Sets the number of cells for X,Y,Z
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
SetMesh(int, int, int, float64, float64, float64, int, int, int)
Sets GridSize, CellSize and PBC at the same time
SetPBC(int, int, int)
Sets the number of repetitions in X,Y,Z to create periodic boundary conditions. The number of repetitions determines the cutoff range for the demagnetization.
EdgeSmooth
Geometry edge smoothing with edgeSmooth^3 samples per cell, 0=staircase, ~8=very smooth
examples: [4]
myShape := cylinder(400e-9, 20e-9).RotX(45*pi/180).Transl(1e-6,0,0))
anotherShape := Circle(400e-9).sub(Circle(200e-9))
Cell(int, int, int) Shape
Single cell with given integer index (i, j, k)
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
Circle(float64) Shape
2D Circle with diameter in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Cone(float64, float64) Shape
3D Cone with diameter and height in meter. The top of the cone points in the +z direction.
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Cuboid(float64, float64, float64) Shape
Cuboid with sides in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
Cylinder(float64, float64) Shape
3D Cylinder with diameter and height in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Ellipse(float64, float64) Shape
2D Ellipse with axes in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [14]
Ellipsoid(float64, float64, float64) Shape
3D Ellipsoid with axes in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
GrainRoughness(float64, float64, float64, int) Shape
Grainy surface with different heights per grain with a typical grain size (first argument), minimal height (second argument), and maximal height (third argument). The last argument is a seed for the random number generator.
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
ImageShape(string) Shape
Use black/white image as shape
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
Layer(int) Shape
Single layer (along z), by integer index starting from 0
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Layers(int, int) Shape
Part of space between cell layer1 (inclusive) and layer2 (exclusive), in integer indices
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
Rect(float64, float64) Shape
2D rectangle with size in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Square(float64) Shape
2D square with size in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [6]
Universe() Shape
Entire space
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
XRange(float64, float64) Shape
Part of space between x1 (inclusive) and x2 (exclusive), in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
YRange(float64, float64) Shape
Part of space between y1 (inclusive) and y2 (exclusive), in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
ZRange(float64, float64) Shape
Part of space between z1 (inclusive) and z2 (exclusive), in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
regions
to verify whether each cell is assigned to the intended region. Each region can have its own material parameters, and we can output averages over each region. E.g.:
DefRegion(1, circle(1e-6))
DefRegion(0, circle(1e-6).Inverse()) // redundant
save(regions)
Msat.SetRegion(1, 800e6)
tableAdd(m.Region(1)) // add average m over region 1 to table
DefRegion(int, Shape)
Define a material region with given index (0-255) and shape
DefRegionCell(int, int, int, int)
Set a material region (first argument) in one cell by the index of the cell (last three arguments)
Config
to m, setting it in separate regions, or by loading a file directly.
m = uniform(1, 0, 0)
m.SetRegion(1, vortex(1, 1))
m.LoadFile("config.ovf")
m.SetInShape(circle(50e-9), uniform(0,0,1))
m
Reduced magnetization (unit length)
methods: Average( ) Buffer( ) Comp( int ) EvalTo( Slice ) GetCell( int int int ) LoadFile( string ) Quantity( ) Region( int ) Set( Config ) SetArray( Slice ) SetCell( int int int data.Vector ) SetInShape( Shape Config ) SetRegion( int Config )
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
Antivortex(int, int) Config
Antivortex magnetization with given circulation and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
BlochSkyrmion(int, int) Config
Bloch skyrmion magnetization with given chirality and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
examples: [5]
Conical(data.Vector, data.Vector, float64) Config
Conical state for given wave vector, cone direction, and cone angle
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
Helical(data.Vector) Config
Helical state for given wave vector
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
NeelSkyrmion(int, int) Config
Néél skyrmion magnetization with given charge and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
examples: [5]
RandomMag() Config
Random magnetization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
RandomMagSeed(int) Config
Random magnetization with given seed
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
TwoDomain(float64, float64, float64, float64, float64, float64, float64, float64, float64) Config
Twodomain magnetization with with given magnetization in left domain, wall, and right domain
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
Uniform(float64, float64, float64) Config
Uniform magnetization in given direction
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
Vortex(int, int) Config
Vortex magnetization with given circulation and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
VortexWall(float64, float64, int, int) Config
Vortex wall magnetization with given mx in left and right domain and core circulation and polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
examples: [5]
Msat = 800e3
AnisU = vector(1, 0, 0)
When regions are defined, they can also be set region-wise:
Msat.SetRegion(0, 800e3)
Msat.SetRegion(1, 540e3)
Material parameters can be functions of time as well. E.g.:
f := 500e6
Ku1 = 500 * sin(2*pi*f*t)
Aex
Exchange stiffness (J/m)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [1] [2] [3] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
alpha
Landau-Lifshitz damping constant
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
anisC1
Cubic anisotropy direction #1
methods: Average( ) Comp( int ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
examples: [12]
anisC2
Cubic anisotorpy directon #2
methods: Average( ) Comp( int ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
examples: [12]
anisU
Uniaxial anisotropy direction
methods: Average( ) Comp( int ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
B1
First magneto-elastic coupling constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
B2
Second magneto-elastic coupling constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Dbulk
Bulk Dzyaloshinskii-Moriya strength (J/m2)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Dind
Interfacial Dzyaloshinskii-Moriya strength (J/m2)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
EpsilonPrime
Slonczewski secondairy STT term ε'
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [14]
FreeLayerThickness
Slonczewski free layer thickness (if set to zero (default), then the thickness will be deduced from the mesh size) (m)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
frozenspins
Defines spins that should be fixed
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Kc1
1st order cubic anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [12]
Kc2
2nd order cubic anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Kc3
3rd order cubic anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Ku1
1st order uniaxial anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Ku2
2nd order uniaxial anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Lambda
Slonczewski Λ parameter
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [14]
Msat
Saturation magnetization (A/m)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [1] [2] [3] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
NoDemagSpins
Disable magnetostatic interaction per region (default=0, set to 1 to disable). E.g.: NoDemagSpins.SetRegion(5, 1) disables the magnetostatic interaction in region 5.
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Pol
Electrical current polarization
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Temp
Temperature (K)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
B_ext = vector(0.01, 1e-6*sin(2*pi*f*t), 0)
B_ext.SetRegion(1, vector(0, 0, 0.1))
Additionally, an arbitrary number of time- and space-dependent vector fields of the form g(x,y,z) * f(t)
may be added. (E.g., to simulate the field of an antenna or an arbitrary current running through the magnet)
B_ext.Add(LoadFile("antenna.ovf"), sin(2*pi*f*t))
J.Add(LoadFile("current.ovf"), 1)
B_ext
Externally applied field (T)
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
FixedLayer
Slonczewski fixed layer polarization
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
examples: [14]
J
Electrical current density (A/m2)
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
J
and the polarization Pol
.
xi
:
J = vector(1e12, 0, 0)
Pol = 1
xi = 0.1
Lambda
and EpsilonPrime
respectively:
DisableZhangLiTorque = true
J = vector(1e12, 0, 0)
Pol = 0.6
FixedLayer = vector(1,0,0)
FixedLayerPosition = FIXEDLAYER_TOP
EpsilonPrime = 0.02
Lambda = 1
DisableSlonczewskiTorque
Disables Slonczewski torque (default=false)
DisableZhangLiTorque
Disables Zhang-Li torque (default=false)
EpsilonPrime
Slonczewski secondairy STT term ε'
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [14]
FixedLayer
Slonczewski fixed layer polarization
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
examples: [14]
FIXEDLAYER_BOTTOM
FixedLayerPosition = FIXEDLAYER_BOTTOM instructs mumax3 that fixed layer is underneath of the free layer
FIXEDLAYER_TOP
FixedLayerPosition = FIXEDLAYER_TOP instructs mumax3 that fixed layer is on top of the free layer
FixedLayerPosition
Position of the fixed layer: FIXEDLAYER_TOP, FIXEDLAYER_BOTTOM (default=FIXEDLAYER_TOP)
FreeLayerThickness
Slonczewski free layer thickness (if set to zero (default), then the thickness will be deduced from the mesh size) (m)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
J
Electrical current density (A/m2)
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
Lambda
Slonczewski Λ parameter
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [14]
Pol
Electrical current polarization
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Mumax3 has built-in generation of MFM images from a 2D magnetization. The MFM tip lift can be freely chosen. By default the tip magnetization is modeled as a point monopole at the apex. This is sufficient for most situations. Nevertheless, it is also possible to model partially magnetized tips by setting MFMDipole to the magnetized portion of the tip, in meters. E.g., if only the first 20nm of the tip is (vertically) magnetized, set MFMDipole=20e-9.
MFMDipole
Height of vertically magnetized part of MFM tip
m // magnetization quantity
m.Comp(0) // x-component
m.Region(1) // magnetization in region 1 (0 elsewhere)
B_anis
Anisotropy field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_custom
User-defined field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_demag
Magnetostatic field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_eff
Effective field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_exch
Exchange field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_mel
Magneto-elastic filed (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_therm
Thermal field (T)
methods: AddTo( Slice ) EvalTo( Slice )
DindCoupling
Average DMI coupling with neighbors (arb.)
methods: Average( ) EvalTo( Slice ) Region( int )
dt
Time Step (s)
methods: Average( ) EvalTo( Slice ) Get( )
E_anis
total anisotropy energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_custom
total energy of user-defined field (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_demag
Magnetostatic energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_exch
Total exchange energy (including the DMI energy) (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_mel
Magneto-elastic energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_therm
Thermal energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_Zeeman
Zeeman energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
Edens_anis
Anisotropy energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_custom
Energy density of user-defined field. (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_demag
Magnetostatic energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_exch
Total exchange energy density (including the DMI energy density) (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_mel
Magneto-elastic energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_therm
Thermal energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_total
Total energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_Zeeman
Zeeman energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
ExchCoupling
Average exchange coupling with neighbors (arb.)
methods: Average( ) EvalTo( Slice ) Region( int )
examples: [12]
F_mel
Magneto-elastic force density (N/m3)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
geom
Cell fill fraction (0..1)
methods: Average( ) EvalTo( Slice ) Gpu( )
LastErr
Error of last step
methods: Average( ) EvalTo( Slice ) Get( )
LLtorque
Landau-Lifshitz torque/γ0 (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
m_full
Unnormalized magnetization (A/m)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
MaxAngle
maximum angle between neighboring spins (rad)
methods: Average( ) EvalTo( Slice ) Get( )
maxTorque
Maximum torque/γ0, over all cells (T)
methods: Average( ) EvalTo( Slice ) Get( )
NEval
Total number of torque evaluations
methods: Average( ) EvalTo( Slice ) Get( )
PeakErr
Overall maxium error per step
methods: Average( ) EvalTo( Slice ) Get( )
spinAngle
Angle between neighboring spins (rad)
methods: Average( ) EvalTo( Slice ) Region( int )
STTorque
Spin-transfer torque/γ0 (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
torque
Total torque/γ0 (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
save(m) // save full magnetization
save(m.Comp(0)) // save only x-component
save(CropLayer(m, 13)) // save only layer 13
save(CropLayer(m.Comp(0), 13)) // save only x-component of layer 13
Or even:
mx := m.Comp(0)
mx13 := CropLayer(mx, 13)
save(mx13)
tableAdd(mx13)
Crop(Quantity, int, int, int, int, int, int) *cropped
Crops a quantity to cell ranges [x1,x2[, [y1,y2[, [z1,z2[
methods: Average( ) EvalTo( Slice )
examples: [8]
CropLayer(Quantity, int) *cropped
Crops a quantity to a single layer
methods: Average( ) EvalTo( Slice )
CropRegion(Quantity, int) *cropped
Crops a quantity to a region
methods: Average( ) EvalTo( Slice )
CropX(Quantity, int, int) *cropped
Crops a quantity to cell ranges [x1,x2[
methods: Average( ) EvalTo( Slice )
CropY(Quantity, int, int) *cropped
Crops a quantity to cell ranges [y1,y2[
methods: Average( ) EvalTo( Slice )
examples: [8]
CropZ(Quantity, int, int) *cropped
Crops a quantity to cell ranges [z1,z2[
methods: Average( ) EvalTo( Slice )
TableAdd()
.
save(B_ext)
tableadd(B_ext)
tablesave()
Optionally, the output/averaging can be done over a single region:
save(m.Region(1))
TableAdd(m.Region(1))
User-defined variables can be added to the table with TableAddVar()
.
myField := 0.42
TableAddVar(myField, "B_extra", "T")
myField = ...
AutoSave(Quantity, float64)
Auto save space-dependent quantity every period (s).
AutoSnapshot(Quantity, float64)
Auto save image of quantity every period (s).
DUMP
OutputFormat = DUMP sets text DUMP output
FilenameFormat
printf formatting string for output filenames.
Flush()
Flush all pending output to disk.
Fprintln(string, ...interface {})
Print to file
OutputFormat
Format for data files: OVF1_TEXT, OVF1_BINARY, OVF2_TEXT or OVF2_BINARY
OVF1_BINARY
OutputFormat = OVF1_BINARY sets binary OVF1 output
OVF1_TEXT
OutputFormat = OVF1_TEXT sets text OVF1 output
OVF2_BINARY
OutputFormat = OVF2_BINARY sets binary OVF2 output
OVF2_TEXT
OutputFormat = OVF2_TEXT sets text OVF2 output
Save(Quantity)
Save space-dependent quantity once, with auto filename
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
SaveAs(Quantity, string)
Save space-dependent quantity with custom filename
Snapshot(Quantity)
Save image of quantity
SnapshotAs(Quantity, string)
Save image of quantity with custom filename
SnapshotFormat
Image format for snapshots: jpg, png or gif.
sprint(...interface {}) string
Print all arguments to string with automatic formatting
sprintf(string, ...interface {}) string
Print to string with C-style formatting.
TableAddVar(ScalarFunction, string, string)
Add user-defined variable + name + unit to data table.
TableAutoSave(float64)
Auto-save the data table every period (s). Zero disables save.
TablePrint(...interface {})
Print anyting in the data table
Run(time)
runs the simulation for a given time in seconds, using sensible error settings.
Run(1e-9)
More fine-grained control is provided by RunWhile(condition)
, which runs as long as an arbitrary condition is met. E.g.:
mx := m.comp(0)
RunWhile(mx.average() < 0) // search for switching field during reversal
Optionally, the solver accuracy may be fine-tuned. E.g.:
MaxDt = 1e-12
MinDt = 1e-15
MaxErr = 1e-6
Optionally, a different solver may be chosen (at any point) with SetSolver(int)
. Currently available solver types:
6
: RK56 (Fehlberg) solver. This is the highest order solver available, but which is typically not faster than the RK45 solver.5
: RK45 (Dormand-Prince) solver (the default). An accurate solver, very fast for magnetization dynamics at the cost of some memory usage. 4
: Classical 4th-order Runge-Kutta method. Intended for simulations where a fixed, relatively large time step is desired.3
: RK23 (Bogacki-Shampine) solver. A robust and reasonably fast solver with low memory requirements. Typically outperforms RK45 when relaxing the magnetization with little dynamics, so it used internally by Relax()
. 2
: Adaptive Heun solver. Robust and uses very little memory but takes smaller time steps than the higher-order solvers. Also suited when a fixed, relatively small time step is desired. 1
: Euler solver (requires FixDt = ...
, ignores other settings). Only useful in exceptional situations or for debugging. SetSolver(2) // Heun
FixDt = 1e-15
Relax()
tries to evolve the magnetization as closely as possible to the minimum energy state. This function assumes all excitations have been turned off (temperature, electrical current, time-dependent magnetic fields). During relax precession is disabled and the time t
does not increase. There is no need to set high damping.
In general it is difficult to be sure the minimum energy state has been truly reached. Hence, relax may occasionally return after the energy has reached a local minimum, a saddle point, or a rather flat valley in the energy landscape.
Minimize()
is like Relax, but uses the conjugate gradient method to find the energy minimum. It is usually much faster than Relax, but is a bit less robust against divergence. E.g., a random starting configuration can be Relaxed, but may fail with Minimize. Minimize is very well suited for hysteresis calculations, where we are never far away from the ground state.
RunWhile(func() bool)
Run while condition function is true
Steps(int)
Run the simulation for a number of time steps
dt
Time Step (s)
methods: Average( ) EvalTo( Slice ) Get( )
FixDt
Set a fixed time step, 0 disables fixed step (which is the default)
Headroom
Solver headroom (default = 0.8)
LastErr
Error of last step
methods: Average( ) EvalTo( Slice ) Get( )
MaxDt
Maximum time step the solver can take (s)
MaxErr
Maximum error per step the solver can tolerate (default = 1e-5)
MinDt
Minimum time step the solver can take (s)
MinimizerSamples
Number of max dM to collect for Minimize convergence check.
NEval
Total number of torque evaluations
methods: Average( ) EvalTo( Slice ) Get( )
PeakErr
Overall maxium error per step
methods: Average( ) EvalTo( Slice ) Get( )
RelaxTorqueThreshold
MaxTorque threshold for relax(). If set to -1 (default), relax() will stop when the average torque is steady or increasing.
t
Total simulated time (s)
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
SetSolver(int)
Set solver type. 1:Euler, 2:Heun, 3:Bogaki-Shampine, 4: Runge-Kutta (RK45), 5: Dormand-Prince, 6: Fehlberg, -1: Backward Euler
ext_centerwall(0)
ext_rmSurfaceCharge(0, -1, 1)
TableAdd(TotalShift)
will try to keep mx
(component 0, counting from 0) close to zero. If desired, one can override which "new" magnetization is inserted from the sides by setting ShiftMagL
and ShiftMagR
, though the default behaviour is usually OK.
ShiftGeom
Whether Shift() acts on geometry
ShiftMagD
Upon shift, insert this magnetization from the bottom
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
ShiftMagL
Upon shift, insert this magnetization from the left
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
ShiftMagR
Upon shift, insert this magnetization from the right
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
examples: [15]
ShiftMagU
Upon shift, insert this magnetization from the top
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
ShiftRegions
Whether Shift() acts on regions
TotalShift
Amount by which the simulation has been shifted (m).
ext_bubbledist
Bubble traveled distance (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_BubbleMz
Center magnetization 1.0 or -1.0 (default = 1.0)
ext_bubblepos
Bubble core position (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_bubblespeed
Bubble velocity (m/s)
methods: Average( ) EvalTo( Slice ) Get( )
ext_centerBubble()
centerBubble shifts m after each step to keep the bubble position close to the center of the window
ext_centerWall(int)
centerWall(c) shifts m after each step to keep m_c close to zero
ext_corepos
Vortex core position (x,y) + polarization (z) (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_dwpos
Position of the simulation window while following a domain wall (m)
methods: Average( ) EvalTo( Slice ) Get( )
examples: [11]
ext_dwspeed
Speed of the simulation window while following a domain wall (m/s)
methods: Average( ) EvalTo( Slice ) Get( )
ext_dwtilt
PMA domain wall tilt (rad)
methods: Average( ) EvalTo( Slice ) Get( )
ext_dwxpos
Position of the simulation window while following a domain wall (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_EnableUnsafe()
Deprecated. Only here to ensure maximal backwards compatibility with mumax3.9c.
ext_InterDind(int, int, float64)
Sets Dind coupling between two regions.
ext_InterExchange(int, int, float64)
Sets exchange coupling between two regions.
ext_make3dgrains(float64, int, int, Shape, int)
3D Voronoi tesselation over shape (grain size, starting region number, num regions, shape, seed)
ext_phi
Azimuthal angle (rad)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_rmSurfaceCharge(int, float64, float64)
Compensate magnetic charges on the left and right sides of an in-plane magnetized wire. Arguments: region, mx on left and right side, resp.
examples: [11]
ext_ScaleDind(int, int, float64)
Re-scales Dind coupling between two regions.
ext_ScaleExchange(int, int, float64)
Re-scales exchange coupling between two regions.
ext_theta
Polar angle (rad)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_topologicalcharge
2D topological charge
methods: Average( ) EvalTo( Slice ) Get( )
ext_topologicalchargedensity
2D topological charge density m·(∂m/∂x ✕ ∂m/∂y) (1/m2)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_topologicalchargedensitylattice
2D topological charge density according to Berg and Lüscher (1/m2)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_topologicalchargelattice
2D topological charge according to Berg and Lüscher
methods: Average( ) EvalTo( Slice ) Get( )
ext_topologicalchargedensity
quantity,
it is possible to define this quantity yourselves inside an input script:
cs := 1e-9
setcellsize(cs,cs,cs)
setgridsize(64,64,1)
// Use central finite differences to approximate the spatial derivatives of m
mL := Shifted(m,-1,0,0) // shift left
mR := Shifted(m,1,0,0) // shift right
mD := Shifted(m,0,-1,0) // shift up
mU := Shifted(m,0,1,0) // shift down
dmdx := Mul( Const(1/(2*cs)), Madd(mR,mL,1,-1) )
dmdy := Mul( Const(1/(2*cs)), Madd(mU,mD,1,-1) )
// Define the topological charge density
chargeDensity := Mul( Const(1/(4*pi)), Dot(m, Cross(dmdx,dmdy)))
// Save the topological charge density of a skyrmion
m = neelskyrmion(1,-1)
saveas(chargeDensity, "chargeDensity.ovf")
Add(Quantity, Quantity) Quantity
Add two quantities
methods: EvalTo( )
Const(float64) Quantity
Constant, uniform number
methods: EvalTo( )
ConstVector(float64, float64, float64) Quantity
Constant, uniform vector
methods: EvalTo( )
Cross(Quantity, Quantity) Quantity
Cross product of two vector quantities
methods: EvalTo( )
examples: [12]
Div(Quantity, Quantity) Quantity
Point-wise division of two quantities
methods: EvalTo( )
Dot(Quantity, Quantity) Quantity
Dot product of two vector quantities
methods: EvalTo( )
Madd(Quantity, Quantity, float64, float64) *mAddition
Weighted addition: Madd(Q1,Q2,c1,c2) = c1*Q1 + c2*Q2
methods: EvalTo( Slice )
Masked(Quantity, Shape) Quantity
Mask quantity with shape
methods: EvalTo( )
Mul(Quantity, Quantity) Quantity
Point-wise product of two quantities
methods: EvalTo( )
MulMV(Quantity, Quantity, Quantity, Quantity) Quantity
Matrix-Vector product: MulMV(AX, AY, AZ, m) = (AX·m, AY·m, AZ·m). The arguments Ax, Ay, Az and m are quantities with 3 componets.
methods: EvalTo( )
Shifted(Quantity, int, int, int) Quantity
Shifted quantity
methods: EvalTo( )
Ms := 1100e3
K := 0.5e6
u := ConstVector(1, 0, 0)
anisField := Mul( Const(2*K/Ms) , Mul( Dot(u, m), u))
anisEdens := Mul( Const(-0.5*Ms) , Dot( anisField, m))
AddFieldTerm(anisField) // promote anisField to an effective field term
AddEdensTerm(anisEdens) // promote anisEdens to an energy density term
tableAdd(E_custom) // Add a column with the energy related to the custom field
AddEdensTerm(Quantity)
Add an expression to Edens.
AddFieldTerm(Quantity)
Add an expression to B_eff.
B_custom
User-defined field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
E_custom
total energy of user-defined field (J)
methods: Average( ) EvalTo( Slice ) Get( )
Edens_custom
Energy density of user-defined field. (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
RemoveCustomFields()
Removes all custom fields again
abs(float64) float64
Abs returns the absolute value of x.
acos(float64) float64
Acos returns the arccosine, in radians, of x.
acosh(float64) float64
Acosh returns the inverse hyperbolic cosine of x.
asin(float64) float64
Asin returns the arcsine, in radians, of x.
asinh(float64) float64
Asinh returns the inverse hyperbolic sine of x.
atan(float64) float64
Atan returns the arctangent, in radians, of x.
atan2(float64, float64) float64
Atan2 returns the arc tangent of y/x, using the signs of the two to determine the quadrant of the return value.
atanh(float64) float64
Atanh returns the inverse hyperbolic tangent of x.
cbrt(float64) float64
Cbrt returns the cube root of x.
ceil(float64) float64
Ceil returns the least integer value greater than or equal to x.
cosh(float64) float64
Cosh returns the hyperbolic cosine of x.
DemagAccuracy
Controls accuracy of demag kernel
DoPrecess
Enables LL precession (default=true)
EnableDemag
Enables/disables demag (default=true)
erf(float64) float64
Erf returns the error function of x.
erfc(float64) float64
Erfc returns the complementary error function of x.
Exit()
Exit from the program
exp2(float64) float64
Exp2 returns 2**x, the base-2 exponential of x.
Expect(string, float64, float64, float64)
Used for automated tests: checks if a value is close enough to the expected value
ExpectV(string, data.Vector, data.Vector, float64)
Used for automated tests: checks if a vector is close enough to the expected value
expm1(float64) float64
Expm1 returns e**x - 1, the base-e exponential of x minus 1. It is more accurate than Exp(x) - 1 when x is near zero.
exx
exx component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
exy
exy component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
exz
exz component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
eyy
eyy component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
eyz
eyz component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
ezz
ezz component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
false
floor(float64) float64
Floor returns the greatest integer value less than or equal to x.
gamma(float64) float64
Gamma returns the Gamma function of x.
GammaLL
Gyromagnetic ratio in rad/Ts
heaviside(float64) float64
Returns 1 if x>0, 0 if x<0, and 0.5 if x==0
hypot(float64, float64) float64
Hypot returns Sqrt(p*p + q*q), taking care to avoid unnecessary overflow and underflow.
Index2Coord(int, int, int) data.Vector
Convert cell index to x,y,z coordinate in meter
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
examples: [15]
inf
Inf returns positive infinity if sign >= 0, negative infinity if sign < 0.
isInf(float64, int) bool
IsInf reports whether f is an infinity, according to sign. If sign > 0, IsInf reports whether f is positive infinity. If sign < 0, IsInf reports
isNaN(float64) bool
IsNaN reports whether f is an IEEE 754 “not-a-number” value.
j0(float64) float64
J0 returns the order-zero Bessel function of the first kind.
j1(float64) float64
J1 returns the order-one Bessel function of the first kind.
jn(int, float64) float64
Jn returns the order-n Bessel function of the first kind.
ldexp(float64, int) float64
Ldexp is the inverse of Frexp. It returns frac × 2**exp.
LoadFile(string) Slice
Load a data file (ovf or dump)
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
examples: [5]
log10(float64) float64
Log10 returns the decimal logarithm of x. The special cases are the same as for Log.
log1p(float64) float64
Log1p returns the natural logarithm of 1 plus its argument x. It is more accurate than Log(1 + x) when x is near zero.
log2(float64) float64
Log2 returns the binary logarithm of x. The special cases are the same as for Log.
mod(float64, float64) float64
Mod returns the floating-point remainder of x/y. The magnitude of the result is less than y and its sign agrees with that of x.
NewScalarMask(int, int, int) Slice
Makes a 3D array of scalars
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
NewSlice(int, int, int, int) Slice
Makes a 4D array with a specified number of components (first argument) and a specified size nx,ny,nz (remaining arguments)
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
NewVectorMask(int, int, int) Slice
Makes a 3D array of vectors
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
examples: [15]
now() time.Time
Returns the current time
methods: Add( time.Duration ) AddDate( int int int ) After( time.Time ) AppendFormat( []uint8 string ) Before( time.Time ) Clock( ) Date( ) Day( ) Equal( time.Time ) Format( string ) GobEncode( ) Hour( ) ISOWeek( ) In( *time.Location ) IsZero( ) Local( ) Location( ) MarshalBinary( ) MarshalJSON( ) MarshalText( ) Minute( ) Month( ) Nanosecond( ) Round( time.Duration ) Second( ) Sub( time.Time ) Truncate( time.Duration ) UTC( ) Unix( ) UnixNano( ) Weekday( ) Year( ) YearDay( ) Zone( )
OpenBC
Use open boundary conditions (default=false)
pow10(int) float64
Pow10 returns 10**n, the base-10 exponential of n.
randExp() float64
Exponentially distributed random number between 0 and +inf, mean=1
randInt(int) int
Random non-negative integer
randSeed(int)
Sets the random number seed
remainder(float64, float64) float64
Remainder returns the IEEE 754 floating-point remainder of x/y.
Sign(float64) float64
Signum function
sin(float64) float64
Sin returns the sine of the radian argument x.
sinc(float64) float64
Sinc returns sin(x)/x. If x=0, then Sinc(x) returns 1.
since(time.Time) time.Duration
Returns the time elapsed since argument
methods: Hours( ) Microseconds( ) Milliseconds( ) Minutes( ) Nanoseconds( ) Round( time.Duration ) Seconds( ) Truncate( time.Duration )
sinh(float64) float64
Sinh returns the hyperbolic sine of x.
tan(float64) float64
Tan returns the tangent of the radian argument x.
tanh(float64) float64
Tanh returns the hyperbolic tangent of x.
ThermSeed(int)
Set a random seed for thermal noise
true
trunc(float64) float64
Trunc returns the integer value of x.
Vector(float64, float64, float64) data.Vector
Constructs a vector with given components
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
y0(float64) float64
Y0 returns the order-zero Bessel function of the second kind.
y1(float64) float64
Y1 returns the order-one Bessel function of the second kind.
yn(int, float64) float64
Yn returns the order-n Bessel function of the second kind.
abs(float64) float64
Abs returns the absolute value of x.
acos(float64) float64
Acos returns the arccosine, in radians, of x.
acosh(float64) float64
Acosh returns the inverse hyperbolic cosine of x.
Add(Quantity, Quantity) Quantity
Add two quantities
methods: EvalTo( )
AddEdensTerm(Quantity)
Add an expression to Edens.
AddFieldTerm(Quantity)
Add an expression to B_eff.
Aex
Exchange stiffness (J/m)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [1] [2] [3] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
alpha
Landau-Lifshitz damping constant
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
anisC1
Cubic anisotropy direction #1
methods: Average( ) Comp( int ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
examples: [12]
anisC2
Cubic anisotorpy directon #2
methods: Average( ) Comp( int ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
examples: [12]
anisU
Uniaxial anisotropy direction
methods: Average( ) Comp( int ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
Antivortex(int, int) Config
Antivortex magnetization with given circulation and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
asin(float64) float64
Asin returns the arcsine, in radians, of x.
asinh(float64) float64
Asinh returns the inverse hyperbolic sine of x.
atan(float64) float64
Atan returns the arctangent, in radians, of x.
atan2(float64, float64) float64
Atan2 returns the arc tangent of y/x, using the signs of the two to determine the quadrant of the return value.
atanh(float64) float64
Atanh returns the inverse hyperbolic tangent of x.
AutoSave(Quantity, float64)
Auto save space-dependent quantity every period (s).
AutoSnapshot(Quantity, float64)
Auto save image of quantity every period (s).
B1
First magneto-elastic coupling constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
B2
Second magneto-elastic coupling constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
B_anis
Anisotropy field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_custom
User-defined field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_demag
Magnetostatic field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_eff
Effective field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_exch
Exchange field (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_ext
Externally applied field (T)
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
B_mel
Magneto-elastic filed (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
B_therm
Thermal field (T)
methods: AddTo( Slice ) EvalTo( Slice )
BlochSkyrmion(int, int) Config
Bloch skyrmion magnetization with given chirality and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
examples: [5]
cbrt(float64) float64
Cbrt returns the cube root of x.
ceil(float64) float64
Ceil returns the least integer value greater than or equal to x.
Cell(int, int, int) Shape
Single cell with given integer index (i, j, k)
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
Circle(float64) Shape
2D Circle with diameter in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Cone(float64, float64) Shape
3D Cone with diameter and height in meter. The top of the cone points in the +z direction.
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Conical(data.Vector, data.Vector, float64) Config
Conical state for given wave vector, cone direction, and cone angle
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
Const(float64) Quantity
Constant, uniform number
methods: EvalTo( )
ConstVector(float64, float64, float64) Quantity
Constant, uniform vector
methods: EvalTo( )
cosh(float64) float64
Cosh returns the hyperbolic cosine of x.
Crop(Quantity, int, int, int, int, int, int) *cropped
Crops a quantity to cell ranges [x1,x2[, [y1,y2[, [z1,z2[
methods: Average( ) EvalTo( Slice )
examples: [8]
CropLayer(Quantity, int) *cropped
Crops a quantity to a single layer
methods: Average( ) EvalTo( Slice )
CropRegion(Quantity, int) *cropped
Crops a quantity to a region
methods: Average( ) EvalTo( Slice )
CropX(Quantity, int, int) *cropped
Crops a quantity to cell ranges [x1,x2[
methods: Average( ) EvalTo( Slice )
CropY(Quantity, int, int) *cropped
Crops a quantity to cell ranges [y1,y2[
methods: Average( ) EvalTo( Slice )
examples: [8]
CropZ(Quantity, int, int) *cropped
Crops a quantity to cell ranges [z1,z2[
methods: Average( ) EvalTo( Slice )
Cross(Quantity, Quantity) Quantity
Cross product of two vector quantities
methods: EvalTo( )
examples: [12]
Cuboid(float64, float64, float64) Shape
Cuboid with sides in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
Cylinder(float64, float64) Shape
3D Cylinder with diameter and height in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Dbulk
Bulk Dzyaloshinskii-Moriya strength (J/m2)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
DefRegion(int, Shape)
Define a material region with given index (0-255) and shape
DefRegionCell(int, int, int, int)
Set a material region (first argument) in one cell by the index of the cell (last three arguments)
DemagAccuracy
Controls accuracy of demag kernel
Dind
Interfacial Dzyaloshinskii-Moriya strength (J/m2)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
DindCoupling
Average DMI coupling with neighbors (arb.)
methods: Average( ) EvalTo( Slice ) Region( int )
DisableSlonczewskiTorque
Disables Slonczewski torque (default=false)
DisableZhangLiTorque
Disables Zhang-Li torque (default=false)
Div(Quantity, Quantity) Quantity
Point-wise division of two quantities
methods: EvalTo( )
DoPrecess
Enables LL precession (default=true)
Dot(Quantity, Quantity) Quantity
Dot product of two vector quantities
methods: EvalTo( )
dt
Time Step (s)
methods: Average( ) EvalTo( Slice ) Get( )
DUMP
OutputFormat = DUMP sets text DUMP output
E_anis
total anisotropy energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_custom
total energy of user-defined field (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_demag
Magnetostatic energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_exch
Total exchange energy (including the DMI energy) (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_mel
Magneto-elastic energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_therm
Thermal energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
E_Zeeman
Zeeman energy (J)
methods: Average( ) EvalTo( Slice ) Get( )
Edens_anis
Anisotropy energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_custom
Energy density of user-defined field. (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_demag
Magnetostatic energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_exch
Total exchange energy density (including the DMI energy density) (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_mel
Magneto-elastic energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_therm
Thermal energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_total
Total energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
Edens_Zeeman
Zeeman energy density (J/m3)
methods: Average( ) EvalTo( Slice ) Region( int )
EdgeSmooth
Geometry edge smoothing with edgeSmooth^3 samples per cell, 0=staircase, ~8=very smooth
examples: [4]
Ellipse(float64, float64) Shape
2D Ellipse with axes in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [14]
Ellipsoid(float64, float64, float64) Shape
3D Ellipsoid with axes in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
EnableDemag
Enables/disables demag (default=true)
EpsilonPrime
Slonczewski secondairy STT term ε'
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [14]
erf(float64) float64
Erf returns the error function of x.
erfc(float64) float64
Erfc returns the complementary error function of x.
ExchCoupling
Average exchange coupling with neighbors (arb.)
methods: Average( ) EvalTo( Slice ) Region( int )
examples: [12]
Exit()
Exit from the program
exp2(float64) float64
Exp2 returns 2**x, the base-2 exponential of x.
Expect(string, float64, float64, float64)
Used for automated tests: checks if a value is close enough to the expected value
ExpectV(string, data.Vector, data.Vector, float64)
Used for automated tests: checks if a vector is close enough to the expected value
expm1(float64) float64
Expm1 returns e**x - 1, the base-e exponential of x minus 1. It is more accurate than Exp(x) - 1 when x is near zero.
ext_bubbledist
Bubble traveled distance (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_BubbleMz
Center magnetization 1.0 or -1.0 (default = 1.0)
ext_bubblepos
Bubble core position (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_bubblespeed
Bubble velocity (m/s)
methods: Average( ) EvalTo( Slice ) Get( )
ext_centerBubble()
centerBubble shifts m after each step to keep the bubble position close to the center of the window
ext_centerWall(int)
centerWall(c) shifts m after each step to keep m_c close to zero
ext_corepos
Vortex core position (x,y) + polarization (z) (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_dwpos
Position of the simulation window while following a domain wall (m)
methods: Average( ) EvalTo( Slice ) Get( )
examples: [11]
ext_dwspeed
Speed of the simulation window while following a domain wall (m/s)
methods: Average( ) EvalTo( Slice ) Get( )
ext_dwtilt
PMA domain wall tilt (rad)
methods: Average( ) EvalTo( Slice ) Get( )
ext_dwxpos
Position of the simulation window while following a domain wall (m)
methods: Average( ) EvalTo( Slice ) Get( )
ext_EnableUnsafe()
Deprecated. Only here to ensure maximal backwards compatibility with mumax3.9c.
ext_InterDind(int, int, float64)
Sets Dind coupling between two regions.
ext_InterExchange(int, int, float64)
Sets exchange coupling between two regions.
ext_make3dgrains(float64, int, int, Shape, int)
3D Voronoi tesselation over shape (grain size, starting region number, num regions, shape, seed)
ext_phi
Azimuthal angle (rad)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_rmSurfaceCharge(int, float64, float64)
Compensate magnetic charges on the left and right sides of an in-plane magnetized wire. Arguments: region, mx on left and right side, resp.
examples: [11]
ext_ScaleDind(int, int, float64)
Re-scales Dind coupling between two regions.
ext_ScaleExchange(int, int, float64)
Re-scales exchange coupling between two regions.
ext_theta
Polar angle (rad)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_topologicalcharge
2D topological charge
methods: Average( ) EvalTo( Slice ) Get( )
ext_topologicalchargedensity
2D topological charge density m·(∂m/∂x ✕ ∂m/∂y) (1/m2)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_topologicalchargedensitylattice
2D topological charge density according to Berg and Lüscher (1/m2)
methods: Average( ) EvalTo( Slice ) Region( int )
ext_topologicalchargelattice
2D topological charge according to Berg and Lüscher
methods: Average( ) EvalTo( Slice ) Get( )
exx
exx component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
exy
exy component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
exz
exz component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
eyy
eyy component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
eyz
eyz component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
ezz
ezz component of the strain tensor
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFn( int func() [3]float64 )
F_mel
Magneto-elastic force density (N/m3)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
false
FilenameFormat
printf formatting string for output filenames.
FixDt
Set a fixed time step, 0 disables fixed step (which is the default)
FixedLayer
Slonczewski fixed layer polarization
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
examples: [14]
FIXEDLAYER_BOTTOM
FixedLayerPosition = FIXEDLAYER_BOTTOM instructs mumax3 that fixed layer is underneath of the free layer
FIXEDLAYER_TOP
FixedLayerPosition = FIXEDLAYER_TOP instructs mumax3 that fixed layer is on top of the free layer
FixedLayerPosition
Position of the fixed layer: FIXEDLAYER_TOP, FIXEDLAYER_BOTTOM (default=FIXEDLAYER_TOP)
floor(float64) float64
Floor returns the greatest integer value less than or equal to x.
Flush()
Flush all pending output to disk.
Fprintln(string, ...interface {})
Print to file
FreeLayerThickness
Slonczewski free layer thickness (if set to zero (default), then the thickness will be deduced from the mesh size) (m)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
frozenspins
Defines spins that should be fixed
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
gamma(float64) float64
Gamma returns the Gamma function of x.
GammaLL
Gyromagnetic ratio in rad/Ts
geom
Cell fill fraction (0..1)
methods: Average( ) EvalTo( Slice ) Gpu( )
GrainRoughness(float64, float64, float64, int) Shape
Grainy surface with different heights per grain with a typical grain size (first argument), minimal height (second argument), and maximal height (third argument). The last argument is a seed for the random number generator.
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Headroom
Solver headroom (default = 0.8)
heaviside(float64) float64
Returns 1 if x>0, 0 if x<0, and 0.5 if x==0
Helical(data.Vector) Config
Helical state for given wave vector
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
hypot(float64, float64) float64
Hypot returns Sqrt(p*p + q*q), taking care to avoid unnecessary overflow and underflow.
ImageShape(string) Shape
Use black/white image as shape
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
Index2Coord(int, int, int) data.Vector
Convert cell index to x,y,z coordinate in meter
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
examples: [15]
inf
Inf returns positive infinity if sign >= 0, negative infinity if sign < 0.
isInf(float64, int) bool
IsInf reports whether f is an infinity, according to sign. If sign > 0, IsInf reports whether f is positive infinity. If sign < 0, IsInf reports
isNaN(float64) bool
IsNaN reports whether f is an IEEE 754 “not-a-number” value.
J
Electrical current density (A/m2)
methods: Add( Slice ScalarFunction ) AddGo( Slice func() float64 ) AddTo( Slice ) Average( ) Comp( int ) EvalTo( Slice ) IsUniform( ) MSlice( ) Region( int ) RemoveExtraTerms( ) Set( data.Vector ) SetRegion( int VectorFunction ) SetRegionFn( int func() [3]float64 )
j0(float64) float64
J0 returns the order-zero Bessel function of the first kind.
j1(float64) float64
J1 returns the order-one Bessel function of the first kind.
jn(int, float64) float64
Jn returns the order-n Bessel function of the first kind.
Kc1
1st order cubic anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [12]
Kc2
2nd order cubic anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Kc3
3rd order cubic anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Ku1
1st order uniaxial anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Ku2
2nd order uniaxial anisotropy constant (J/m3)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
Lambda
Slonczewski Λ parameter
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [14]
LastErr
Error of last step
methods: Average( ) EvalTo( Slice ) Get( )
Layer(int) Shape
Single layer (along z), by integer index starting from 0
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Layers(int, int) Shape
Part of space between cell layer1 (inclusive) and layer2 (exclusive), in integer indices
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [4]
ldexp(float64, int) float64
Ldexp is the inverse of Frexp. It returns frac × 2**exp.
LLtorque
Landau-Lifshitz torque/γ0 (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
LoadFile(string) Slice
Load a data file (ovf or dump)
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
examples: [5]
log10(float64) float64
Log10 returns the decimal logarithm of x. The special cases are the same as for Log.
log1p(float64) float64
Log1p returns the natural logarithm of 1 plus its argument x. It is more accurate than Log(1 + x) when x is near zero.
log2(float64) float64
Log2 returns the binary logarithm of x. The special cases are the same as for Log.
m
Reduced magnetization (unit length)
methods: Average( ) Buffer( ) Comp( int ) EvalTo( Slice ) GetCell( int int int ) LoadFile( string ) Quantity( ) Region( int ) Set( Config ) SetArray( Slice ) SetCell( int int int data.Vector ) SetInShape( Shape Config ) SetRegion( int Config )
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
m_full
Unnormalized magnetization (A/m)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
Madd(Quantity, Quantity, float64, float64) *mAddition
Weighted addition: Madd(Q1,Q2,c1,c2) = c1*Q1 + c2*Q2
methods: EvalTo( Slice )
Masked(Quantity, Shape) Quantity
Mask quantity with shape
methods: EvalTo( )
MaxAngle
maximum angle between neighboring spins (rad)
methods: Average( ) EvalTo( Slice ) Get( )
MaxDt
Maximum time step the solver can take (s)
MaxErr
Maximum error per step the solver can tolerate (default = 1e-5)
maxTorque
Maximum torque/γ0, over all cells (T)
methods: Average( ) EvalTo( Slice ) Get( )
MFMDipole
Height of vertically magnetized part of MFM tip
MinDt
Minimum time step the solver can take (s)
MinimizerSamples
Number of max dM to collect for Minimize convergence check.
mod(float64, float64) float64
Mod returns the floating-point remainder of x/y. The magnitude of the result is less than y and its sign agrees with that of x.
Msat
Saturation magnetization (A/m)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
examples: [1] [2] [3] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
Mul(Quantity, Quantity) Quantity
Point-wise product of two quantities
methods: EvalTo( )
MulMV(Quantity, Quantity, Quantity, Quantity) Quantity
Matrix-Vector product: MulMV(AX, AY, AZ, m) = (AX·m, AY·m, AZ·m). The arguments Ax, Ay, Az and m are quantities with 3 componets.
methods: EvalTo( )
NeelSkyrmion(int, int) Config
Néél skyrmion magnetization with given charge and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
examples: [5]
NEval
Total number of torque evaluations
methods: Average( ) EvalTo( Slice ) Get( )
NewScalarMask(int, int, int) Slice
Makes a 3D array of scalars
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
NewSlice(int, int, int, int) Slice
Makes a 4D array with a specified number of components (first argument) and a specified size nx,ny,nz (remaining arguments)
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
NewVectorMask(int, int, int) Slice
Makes a 3D array of vectors
methods: CPUAccess( ) Comp( int ) DevPtr( int ) Disable( ) Free( ) GPUAccess( ) Get( int int int int ) Host( ) HostCopy( ) Index( int int int ) IsNil( ) Len( ) MemType( ) Scalars( ) Set( int int int int float64 ) SetScalar( int int int float64 ) SetVector( int int int data.Vector ) Size( ) Tensors( ) Vectors( )
examples: [15]
NoDemagSpins
Disable magnetostatic interaction per region (default=0, set to 1 to disable). E.g.: NoDemagSpins.SetRegion(5, 1) disables the magnetostatic interaction in region 5.
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
now() time.Time
Returns the current time
methods: Add( time.Duration ) AddDate( int int int ) After( time.Time ) AppendFormat( []uint8 string ) Before( time.Time ) Clock( ) Date( ) Day( ) Equal( time.Time ) Format( string ) GobEncode( ) Hour( ) ISOWeek( ) In( *time.Location ) IsZero( ) Local( ) Location( ) MarshalBinary( ) MarshalJSON( ) MarshalText( ) Minute( ) Month( ) Nanosecond( ) Round( time.Duration ) Second( ) Sub( time.Time ) Truncate( time.Duration ) UTC( ) Unix( ) UnixNano( ) Weekday( ) Year( ) YearDay( ) Zone( )
OpenBC
Use open boundary conditions (default=false)
OutputFormat
Format for data files: OVF1_TEXT, OVF1_BINARY, OVF2_TEXT or OVF2_BINARY
OVF1_BINARY
OutputFormat = OVF1_BINARY sets binary OVF1 output
OVF1_TEXT
OutputFormat = OVF1_TEXT sets text OVF1 output
OVF2_BINARY
OutputFormat = OVF2_BINARY sets binary OVF2 output
OVF2_TEXT
OutputFormat = OVF2_TEXT sets text OVF2 output
PeakErr
Overall maxium error per step
methods: Average( ) EvalTo( Slice ) Get( )
Pol
Electrical current polarization
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
pow10(int) float64
Pow10 returns 10**n, the base-10 exponential of n.
randExp() float64
Exponentially distributed random number between 0 and +inf, mean=1
randInt(int) int
Random non-negative integer
RandomMag() Config
Random magnetization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
RandomMagSeed(int) Config
Random magnetization with given seed
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
randSeed(int)
Sets the random number seed
Rect(float64, float64) Shape
2D rectangle with size in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
regions
Outputs the region index for each cell
methods: Average( ) EvalTo( Slice ) GetCell( int int int ) Gpu( ) HostArray( ) HostList( ) LoadFile( string ) SetCell( int int int int )
RelaxTorqueThreshold
MaxTorque threshold for relax(). If set to -1 (default), relax() will stop when the average torque is steady or increasing.
remainder(float64, float64) float64
Remainder returns the IEEE 754 floating-point remainder of x/y.
RemoveCustomFields()
Removes all custom fields again
RunWhile(func() bool)
Run while condition function is true
Save(Quantity)
Save space-dependent quantity once, with auto filename
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
SaveAs(Quantity, string)
Save space-dependent quantity with custom filename
SetCellSize(float64, float64, float64)
Sets the X,Y,Z cell size in meters
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
SetGridSize(int, int, int)
Sets the number of cells for X,Y,Z
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
SetMesh(int, int, int, float64, float64, float64, int, int, int)
Sets GridSize, CellSize and PBC at the same time
SetPBC(int, int, int)
Sets the number of repetitions in X,Y,Z to create periodic boundary conditions. The number of repetitions determines the cutoff range for the demagnetization.
SetSolver(int)
Set solver type. 1:Euler, 2:Heun, 3:Bogaki-Shampine, 4: Runge-Kutta (RK45), 5: Dormand-Prince, 6: Fehlberg, -1: Backward Euler
Shifted(Quantity, int, int, int) Quantity
Shifted quantity
methods: EvalTo( )
ShiftGeom
Whether Shift() acts on geometry
ShiftMagD
Upon shift, insert this magnetization from the bottom
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
ShiftMagL
Upon shift, insert this magnetization from the left
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
ShiftMagR
Upon shift, insert this magnetization from the right
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
examples: [15]
ShiftMagU
Upon shift, insert this magnetization from the top
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
ShiftRegions
Whether Shift() acts on regions
Sign(float64) float64
Signum function
sin(float64) float64
Sin returns the sine of the radian argument x.
sinc(float64) float64
Sinc returns sin(x)/x. If x=0, then Sinc(x) returns 1.
since(time.Time) time.Duration
Returns the time elapsed since argument
methods: Hours( ) Microseconds( ) Milliseconds( ) Minutes( ) Nanoseconds( ) Round( time.Duration ) Seconds( ) Truncate( time.Duration )
sinh(float64) float64
Sinh returns the hyperbolic sine of x.
Snapshot(Quantity)
Save image of quantity
SnapshotAs(Quantity, string)
Save image of quantity with custom filename
SnapshotFormat
Image format for snapshots: jpg, png or gif.
spinAngle
Angle between neighboring spins (rad)
methods: Average( ) EvalTo( Slice ) Region( int )
sprint(...interface {}) string
Print all arguments to string with automatic formatting
sprintf(string, ...interface {}) string
Print to string with C-style formatting.
Square(float64) Shape
2D square with size in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
examples: [6]
Steps(int)
Run the simulation for a number of time steps
STTorque
Spin-transfer torque/γ0 (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
t
Total simulated time (s)
examples: [1] [2] [3] [4] [5] [6] [7] [8] [9] [10] [11] [12] [13] [14] [15]
TableAddVar(ScalarFunction, string, string)
Add user-defined variable + name + unit to data table.
TableAutoSave(float64)
Auto-save the data table every period (s). Zero disables save.
TablePrint(...interface {})
Print anyting in the data table
tan(float64) float64
Tan returns the tangent of the radian argument x.
tanh(float64) float64
Tanh returns the hyperbolic tangent of x.
Temp
Temperature (K)
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
ThermSeed(int)
Set a random seed for thermal noise
torque
Total torque/γ0 (T)
methods: Average( ) Comp( int ) EvalTo( Slice ) HostCopy( ) Region( int )
TotalShift
Amount by which the simulation has been shifted (m).
true
trunc(float64) float64
Trunc returns the integer value of x.
TwoDomain(float64, float64, float64, float64, float64, float64, float64, float64, float64) Config
Twodomain magnetization with with given magnetization in left domain, wall, and right domain
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
Uniform(float64, float64, float64) Config
Uniform magnetization in given direction
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
Universe() Shape
Entire space
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
Vector(float64, float64, float64) data.Vector
Constructs a vector with given components
methods: Add( data.Vector ) Cross( data.Vector ) Div( float64 ) Dot( data.Vector ) Len( ) MAdd( float64 data.Vector ) Mul( float64 ) Sub( data.Vector ) X( ) Y( ) Z( )
Vortex(int, int) Config
Vortex magnetization with given circulation and core polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
VortexWall(float64, float64, int, int) Config
Vortex wall magnetization with given mx in left and right domain and core circulation and polarization
methods: Add( float64 Config ) RotZ( float64 ) Scale( float64 float64 float64 ) Transl( float64 float64 float64 )
examples: [5]
xi
Non-adiabaticity of spin-transfer-torque
methods: Average( ) EvalTo( Slice ) GetRegion( int ) IsUniform( ) MSlice( ) Region( int ) Set( float64 ) SetRegion( int ScalarFunction ) SetRegionFuncGo( int func() float64 ) SetRegionValueGo( int float64 )
XRange(float64, float64) Shape
Part of space between x1 (inclusive) and x2 (exclusive), in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
y0(float64) float64
Y0 returns the order-zero Bessel function of the second kind.
y1(float64) float64
Y1 returns the order-one Bessel function of the second kind.
yn(int, float64) float64
Yn returns the order-n Bessel function of the second kind.
YRange(float64, float64) Shape
Part of space between y1 (inclusive) and y2 (exclusive), in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )
ZRange(float64, float64) Shape
Part of space between z1 (inclusive) and z2 (exclusive), in meter
methods: Add( Shape ) Intersect( Shape ) Inverse( ) Repeat( float64 float64 float64 ) RotX( float64 ) RotY( float64 ) RotZ( float64 ) Scale( float64 float64 float64 ) Sub( Shape ) Transl( float64 float64 float64 ) Xor( Shape )