mumax3
GPU-accelerated micromagnetism

Home Download Examples Tutorial API Forum


mumax3.10 API

Basics

Advanced Features

Go to mumax3.9c API

  Back to top  

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( )  

examples: [3] [4] [5] [11] [13] [15]

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 )  

examples: [1] [6] [7] [8] [10] [11] [12] [14] [15]

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 )  

examples: [7] [10] [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 )  

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).

examples: [1] [10] [11] [14] [15]

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 )  

examples: [1] [3] [15]

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 )  

examples: [4] [6] [7] [8] [12]

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( )  

cos(float64) float64

Cos returns the cosine of the radian argument x.

examples: [6] [13] [14]

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 )  

examples: [4] [5] [6]

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

examples: [7] [12] [13]

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_total

total energy (J)

methods: Average( )   EvalTo( Slice )   Get( )  

examples: [13]

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

exp(float64) float64

Exp returns e**x, the base-e exponential of x.

examples: [15]

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

examples: [10] [11]

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_makegrains(float64, int, int)

Voronoi tesselation (grain size, num regions)

examples: [12] [15]

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.

examples: [12] [13] [15]

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( )  

examples: [4] [6] [7] [8] [9] [11] [12] [14]

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.

ilogb(float64) int

Ilogb returns the binary exponent of x as an integer.

examples: [2]

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.

examples: [4] [7] [11]

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 )  

examples: [10] [11] [12] [13] [14] [15]

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 )  

examples: [7] [10] [15]

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 )  

examples: [4] [13] [14]

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]

log(float64) float64

Log returns the natural logarithm of x.

examples: [2] [4]

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.

logb(float64) float64

Logb returns the binary exponent of x.

examples: [2]

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( )  

max(float64, float64) float64

Max returns the larger of x or y.

examples: [3] [12]

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( )  

MFM

MFM image (arb.)

methods: Average( )   EvalTo( Slice )   Region( int )  

examples: [9]

MFMDipole

Height of vertically magnetized part of MFM tip

MFMLift

MFM lift height

examples: [9]

min(float64, float64) float64

Min returns the smaller of x or y.

examples: [3] [6]

MinDt

Minimum time step the solver can take (s)

Minimize()

Use steepest conjugate gradient method to minimize the total energy

examples: [3] [6]

MinimizerSamples

Number of max dM to collect for Minimize convergence check.

MinimizerStop

Stopping max dM for Minimize

examples: [3]

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]

Mu0

Vacuum permeability (Tm/A)

examples: [2]

Mul(Quantity, Quantity) Quantity

Point-wise product of two quantities

methods: EvalTo( )  

examples: [10] [11]

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 )  

norm(float64) float64

Standard normal distribution

examples: [5] [12]

Normalized(Quantity) Quantity

Normalize quantity

methods: EvalTo( )  

examples: [12]

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( )  

pi

examples: [4] [5] [6] [11] [13] [14] [15]

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 )  

examples: [5] [10] [11] [14]

pow(float64, float64) float64

Pow returns x**y, the base-x exponential of y.

examples: [2] [15]

pow10(int) float64

Pow10 returns 10**n, the base-10 exponential of n.

Print(...interface {})

Print to standard output

examples: [2]

rand() float64

Random number between 0 and 1

examples: [3] [5] [12] [15]

randExp() float64

Exponentially distributed random number between 0 and +inf, mean=1

randInt(int) int

Random non-negative integer

randNorm() float64

Standard normal random number

examples: [12]

RandomMag() Config

Random magnetization

methods: Add( float64 Config )   RotZ( float64 )   Scale( float64 float64 float64 )   Transl( float64 float64 float64 )  

examples: [3] [5]

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 )  

examples: [4] [6] [9] [11] [12] [15]

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 )  

examples: [7] [12]

Relax()

Try to minimize the total energy

examples: [1] [2] [3] [9] [10] [11]

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

Run(float64)

Run the simulation for a time in seconds

examples: [1] [7] [10] [11] [12] [14] [15]

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

examples: [4] [5] [7] [9]

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]

SetGeom(Shape)

Sets the geometry to a given shape

examples: [4] [6] [7] [8] [9] [11] [12] [14]

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

Shift(int)

Shifts the simulation by +1/-1 cells along X

examples: [15]

Shifted(Quantity, int, int, int) Quantity

Shifted quantity

methods: EvalTo( )  

ShiftGeom

Whether Shift() acts on geometry

ShiftM

Whether Shift() acts on magnetization

examples: [15]

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.

examples: [5] [6] [13] [14] [15]

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.

sqrt(float64) float64

Sqrt returns the square root of x.

examples: [2]

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]

step

Total number of time steps taken

examples: [3] [10]

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]

TableAdd(Quantity)

Add quantity as a column to the data table.

examples: [3] [11] [13]

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.

examples: [1] [11] [14]

TablePrint(...interface {})

Print anyting in the data table

TableSave()

Save the data table right now (appends one line).

examples: [3] [13]

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 )  

examples: [5] [10] [11]

Uniform(float64, float64, float64) Config

Uniform magnetization in given direction

methods: Add( float64 Config )   RotZ( float64 )   Scale( float64 float64 float64 )   Transl( float64 float64 float64 )  

examples: [1] [2] [5] [6] [7] [13] [14] [15]

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( )  

examples: [1] [3] [5] [7] [10] [11] [12] [14] [15]

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 )  

examples: [5] [8] [9] [12]

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 )  

examples: [10] [11] [12]

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 )  

examples: [4] [7]

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 )