engine

package
v0.0.0-...-f85407f Latest Latest
Warning

This package is not in the latest version of its module.

Go to latest
Published: Dec 11, 2023 License: CC-BY-3.0, Freetype, GPL-3.0-or-later Imports: 38 Imported by: 0

Documentation

Overview

engine does the simulation bookkeeping, I/O and GUI.

space-dependence: value: space-independent param: region-dependent parameter (always input) field: fully space-dependent field

TODO: godoc everything

Index

Constants

View Source
const (
	TILE   = 2           // tile size in grains
	LAMBDA = TILE * TILE // expected grains per tile
)
View Source
const (
	BACKWARD_EULER = -1
	EULER          = 1
	HEUN           = 2
	BOGAKISHAMPINE = 3
	RUNGEKUTTA     = 4
	DORMANDPRINCE  = 5
	FEHLBERG       = 6
)

Arguments for SetSolver

View Source
const (
	X = 0
	Y = 1
	Z = 2
)
View Source
const (
	SCALAR = 1
	VECTOR = 3
)
View Source
const CSS = `` /* 756-byte string literal not displayed */
View Source
const MAX_GNUPLOTS = 5 // maximum allowed number of gnuplot processes
View Source
const NREGION = 256 // maximum number of regions, limited by size of byte.
View Source
const TableAutoflushRate = 5 // auto-flush table every X seconds
View Source
const VERSION = "mumax 3.10"

Variables

View Source
var (
	Ku1        = NewScalarParam("Ku1", "J/m3", "1st order uniaxial anisotropy constant")
	Ku2        = NewScalarParam("Ku2", "J/m3", "2nd order uniaxial anisotropy constant")
	Kc1        = NewScalarParam("Kc1", "J/m3", "1st order cubic anisotropy constant")
	Kc2        = NewScalarParam("Kc2", "J/m3", "2nd order cubic anisotropy constant")
	Kc3        = NewScalarParam("Kc3", "J/m3", "3rd order cubic anisotropy constant")
	AnisU      = NewVectorParam("anisU", "", "Uniaxial anisotropy direction")
	AnisC1     = NewVectorParam("anisC1", "", "Cubic anisotropy direction #1")
	AnisC2     = NewVectorParam("anisC2", "", "Cubic anisotorpy directon #2")
	B_anis     = NewVectorField("B_anis", "T", "Anisotropy field", AddAnisotropyField)
	Edens_anis = NewScalarField("Edens_anis", "J/m3", "Anisotropy energy density", AddAnisotropyEnergyDensity)
	E_anis     = NewScalarValue("E_anis", "J", "total anisotropy energy", GetAnisotropyEnergy)
)

Anisotropy variables

View Source
var (
	B_custom     = NewVectorField("B_custom", "T", "User-defined field", AddCustomField)
	Edens_custom = NewScalarField("Edens_custom", "J/m3", "Energy density of user-defined field.", AddCustomEnergyDensity)
	E_custom     = NewScalarValue("E_custom", "J", "total energy of user-defined field", GetCustomEnergy)
)
View Source
var (
	Msat        = NewScalarParam("Msat", "A/m", "Saturation magnetization")
	M_full      = NewVectorField("m_full", "A/m", "Unnormalized magnetization", SetMFull)
	B_demag     = NewVectorField("B_demag", "T", "Magnetostatic field", SetDemagField)
	Edens_demag = NewScalarField("Edens_demag", "J/m3", "Magnetostatic energy density", AddEdens_demag)
	E_demag     = NewScalarValue("E_demag", "J", "Magnetostatic energy", GetDemagEnergy)

	EnableDemag  = true // enable/disable global demag field
	NoDemagSpins = NewScalarParam("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.")

	DemagAccuracy = 6.0 // Demag accuracy (divide cubes in at most N^3 points)
)

Demag variables

View Source
var (
	Edens_total = NewScalarField("Edens_total", "J/m3", "Total energy density", SetTotalEdens)
	E_total     = NewScalarValue("E_total", "J", "total energy", GetTotalEnergy)
)
View Source
var (
	Aex   = NewScalarParam("Aex", "J/m", "Exchange stiffness", &lex2)
	Dind  = NewScalarParam("Dind", "J/m2", "Interfacial Dzyaloshinskii-Moriya strength", &din2)
	Dbulk = NewScalarParam("Dbulk", "J/m2", "Bulk Dzyaloshinskii-Moriya strength", &dbulk2)

	B_exch     = NewVectorField("B_exch", "T", "Exchange field", AddExchangeField)
	E_exch     = NewScalarValue("E_exch", "J", "Total exchange energy (including the DMI energy)", GetExchangeEnergy)
	Edens_exch = NewScalarField("Edens_exch", "J/m3", "Total exchange energy density (including the DMI energy density)", AddExchangeEnergyDensity)

	// Average exchange coupling with neighbors. Useful to debug inter-region exchange
	ExchCoupling = NewScalarField("ExchCoupling", "arb.", "Average exchange coupling with neighbors", exchangeDecode)
	DindCoupling = NewScalarField("DindCoupling", "arb.", "Average DMI coupling with neighbors", dindDecode)

	OpenBC = false
)
View Source
var (
	BubblePos      = NewVectorValue("ext_bubblepos", "m", "Bubble core position", bubblePos)
	BubbleDist     = NewScalarValue("ext_bubbledist", "m", "Bubble traveled distance", bubbleDist)
	BubbleSpeed    = NewScalarValue("ext_bubblespeed", "m/s", "Bubble velocity", bubbleSpeed)
	BubbleMz       = 1.0
	BackGroundTilt = 0.25
)
View Source
var (
	DWPos   = NewScalarValue("ext_dwpos", "m", "Position of the simulation window while following a domain wall", GetShiftPos) // TODO: make more accurate
	DWxPos  = NewScalarValue("ext_dwxpos", "m", "Position of the simulation window while following a domain wall", GetDWxPos)
	DWSpeed = NewScalarValue("ext_dwspeed", "m/s", "Speed of the simulation window while following a domain wall", getShiftSpeed)
)
View Source
var (
	B1 = NewScalarParam("B1", "J/m3", "First magneto-elastic coupling constant")
	B2 = NewScalarParam("B2", "J/m3", "Second magneto-elastic coupling constant")

	B_mel     = NewVectorField("B_mel", "T", "Magneto-elastic filed", AddMagnetoelasticField)
	F_mel     = NewVectorField("F_mel", "N/m3", "Magneto-elastic force density", GetMagnetoelasticForceDensity)
	Edens_mel = NewScalarField("Edens_mel", "J/m3", "Magneto-elastic energy density", AddMagnetoelasticEnergyDensity)
	E_mel     = NewScalarValue("E_mel", "J", "Magneto-elastic energy", GetMagnetoelasticEnergy)
)
View Source
var (
	Ext_TopologicalCharge        = NewScalarValue("ext_topologicalcharge", "", "2D topological charge", GetTopologicalCharge)
	Ext_TopologicalChargeDensity = NewScalarField("ext_topologicalchargedensity", "1/m2",
		"2D topological charge density m·(∂m/∂x ✕ ∂m/∂y)", SetTopologicalChargeDensity)
)
View Source
var (
	Ext_TopologicalChargeLattice        = NewScalarValue("ext_topologicalchargelattice", "", "2D topological charge according to Berg and Lüscher", GetTopologicalChargeLattice)
	Ext_TopologicalChargeDensityLattice = NewScalarField("ext_topologicalchargedensitylattice", "1/m2",
		"2D topological charge density according to Berg and Lüscher", SetTopologicalChargeDensityLattice)
)
View Source
var (
	// These flags are shared between cmd/mumax3 and Go input files.
	Flag_cachedir    = flag.String("cache", os.TempDir(), "Kernel cache directory (empty disables caching)")
	Flag_gpu         = flag.Int("gpu", 0, "Specify GPU")
	Flag_interactive = flag.Bool("i", false, "Open interactive browser session")
	Flag_od          = flag.String("o", "", "Override output directory")
	Flag_port        = flag.String("http", ":35367", "Port to serve web gui")
	Flag_selftest    = flag.Bool("paranoid", false, "Enable convolution self-test for cuFFT sanity.")
	Flag_silent      = flag.Bool("s", false, "Silent") // provided for backwards compatibility
	Flag_sync        = flag.Bool("sync", false, "Synchronize all CUDA calls (debug)")
	Flag_forceclean  = flag.Bool("f", false, "Force start, clean existing output directory")
)
View Source
var (
	MaxAngle  = NewScalarValue("MaxAngle", "rad", "maximum angle between neighboring spins", GetMaxAngle)
	SpinAngle = NewScalarField("spinAngle", "rad", "Angle between neighboring spins", SetSpinAngle)
)
View Source
var (
	MFM        = NewScalarField("MFM", "arb.", "MFM image", SetMFM)
	MFMLift    inputValue
	MFMTipSize inputValue
)
View Source
var (
	DmSamples int     = 10   // number of dm to keep for convergence check
	StopMaxDm float64 = 1e-6 // stop minimizer if sampled dm is smaller than this
)
View Source
var (
	Time float64 // time in seconds

	Inject                          = make(chan func()) // injects code in between time steps. Used by web interface.
	Dt_si                   float64 = 1e-15             // time step = dt_si (seconds) *dt_mul, which should be nice float32
	MinDt, MaxDt            float64                     // minimum and maximum time step
	MaxErr                  float64 = 1e-5              // maximum error/step
	Headroom                float64 = 0.8               // solver headroom, (Gustafsson, 1992, Control of Error and Convergence in ODE Solvers)
	LastErr, PeakErr        float64                     // error of last step, highest error ever
	LastTorque              float64                     // maxTorque of last time step
	NSteps, NUndone, NEvals int                         // number of good steps, undone steps
	FixDt                   float64                     // fixed time step?

)

Solver globals

View Source
var (
	FilenameFormat = "%s%06d" // formatting string for auto filenames.
	SnapshotFormat = "jpg"    // user-settable snapshot format

)
View Source
var (
	TotalShift, TotalYShift                    float64                        // accumulated window shift (X and Y) in meter
	ShiftMagL, ShiftMagR, ShiftMagU, ShiftMagD data.Vector                    // when shifting m, put these value at the left/right edge.
	ShiftM, ShiftGeom, ShiftRegions            bool        = true, true, true // should shift act on magnetization, geometry, regions?
)
View Source
var (
	Temp        = NewScalarParam("Temp", "K", "Temperature")
	E_therm     = NewScalarValue("E_therm", "J", "Thermal energy", GetThermalEnergy)
	Edens_therm = NewScalarField("Edens_therm", "J/m3", "Thermal energy density", AddThermalEnergyDensity)
	B_therm     thermField // Thermal effective field (T)
)
View Source
var (
	Alpha                            = NewScalarParam("alpha", "", "Landau-Lifshitz damping constant")
	Xi                               = NewScalarParam("xi", "", "Non-adiabaticity of spin-transfer-torque")
	Pol                              = NewScalarParam("Pol", "", "Electrical current polarization")
	Lambda                           = NewScalarParam("Lambda", "", "Slonczewski Λ parameter")
	EpsilonPrime                     = NewScalarParam("EpsilonPrime", "", "Slonczewski secondairy STT term ε'")
	FrozenSpins                      = NewScalarParam("frozenspins", "", "Defines spins that should be fixed") // 1 - frozen, 0 - free. TODO: check if it only contains 0/1 values
	FreeLayerThickness               = NewScalarParam("FreeLayerThickness", "m", "Slonczewski free layer thickness (if set to zero (default), then the thickness will be deduced from the mesh size)")
	FixedLayer                       = NewExcitation("FixedLayer", "", "Slonczewski fixed layer polarization")
	Torque                           = NewVectorField("torque", "T", "Total torque/γ0", SetTorque)
	LLTorque                         = NewVectorField("LLtorque", "T", "Landau-Lifshitz torque/γ0", SetLLTorque)
	STTorque                         = NewVectorField("STTorque", "T", "Spin-transfer torque/γ0", AddSTTorque)
	J                                = NewExcitation("J", "A/m2", "Electrical current density")
	MaxTorque                        = NewScalarValue("maxTorque", "T", "Maximum torque/γ0, over all cells", GetMaxTorque)
	GammaLL                  float64 = 1.7595e11 // Gyromagnetic ratio of spins, in rad/Ts
	Precess                          = true
	DisableZhangLiTorque             = false
	DisableSlonczewskiTorque         = false
)
View Source
var (
	B_ext        = NewExcitation("B_ext", "T", "Externally applied field")
	Edens_zeeman = NewScalarField("Edens_Zeeman", "J/m3", "Zeeman energy density", AddEdens_zeeman)
	E_Zeeman     = NewScalarValue("E_Zeeman", "J", "Zeeman energy", GetZeemanEnergy)
)
View Source
var AddEdens_demag = makeEdensAdder(&B_demag, -0.5)
View Source
var AddEdens_zeeman = makeEdensAdder(B_ext, -1)
View Source
var AddExchangeEnergyDensity = makeEdensAdder(&B_exch, -0.5) // TODO: normal func
View Source
var AddThermalEnergyDensity = makeEdensAdder(&B_therm, -1)
View Source
var B_eff = NewVectorField("B_eff", "T", "Effective field", SetEffectiveField)
View Source
var CorePos = NewVectorValue("ext_corepos", "m", "Vortex core position (x,y) + polarization (z)", corePos)
View Source
var (
	CurrentSignFromFixedLayerPosition = map[FixedLayerPosition]float64{
		FIXEDLAYER_TOP:    1.0,
		FIXEDLAYER_BOTTOM: -1.0,
	}
)
View Source
var DWTiltPMA = NewScalarValue("ext_dwtilt", "rad", "PMA domain wall tilt", dwTiltPMA)

PMA domain wall tilt assuming straight wall.

View Source
var (
	InputFile string
)
View Source
var M magnetization // reduced magnetization (unit length)
View Source
var RelaxTorqueThreshold float64 = -1.

Stopping relax Maxtorque in T. The user can check MaxTorque for sane values (e.g. 1e-3). If set to 0, relax() will stop when the average torque is steady or increasing.

View Source
var StartTime = time.Now()
View Source
var (
	StringFromOutputFormat = map[OutputFormat]string{
		OVF1_TEXT:   "ovf",
		OVF1_BINARY: "ovf",
		OVF2_TEXT:   "ovf",
		OVF2_BINARY: "ovf",
		DUMP:        "dump"}
)
View Source
var Table = *newTable("table") // output handle for tabular data (average magnetization etc.)
View Source
var (
	Timeout = 3 * time.Second // exit finished simulation this long after browser was closed
)

global GUI state stores what is currently shown in the web page.

View Source
var UNAME = fmt.Sprintf("%s [%s_%s %s(%s) CUDA-%d.%d]",
	VERSION, runtime.GOOS, runtime.GOARCH, runtime.Version(), runtime.Compiler,
	cu.CUDA_VERSION/1000, (cu.CUDA_VERSION%1000)/10)
View Source
var World = script.NewWorld()

holds the script state (variables etc)

Functions

func AddAnisotropyEnergyDensity

func AddAnisotropyEnergyDensity(dst *data.Slice)

Add the anisotropy energy density to dst

func AddAnisotropyField

func AddAnisotropyField(dst *data.Slice)

Add the anisotropy field to dst

func AddCustomEnergyDensity

func AddCustomEnergyDensity(dst *data.Slice)

Adds the custom energy densities (defined with AddCustomE

func AddCustomField

func AddCustomField(dst *data.Slice)

AddCustomField evaluates the user-defined custom field terms and adds the result to dst.

func AddEdensTerm

func AddEdensTerm(e Quantity)

AddEnergyTerm adds an energy density function (returning Joules/m³) to Edens_total. Needed when AddFieldTerm was used and a correct energy is needed (e.g. for Relax, Minimize, ...).

func AddExchangeField

func AddExchangeField(dst *data.Slice)

Adds the current exchange field to dst

func AddFieldTerm

func AddFieldTerm(b Quantity)

AddFieldTerm adds an effective field function (returning Teslas) to B_eff. Be sure to also add the corresponding energy term using AddEnergyTerm.

func AddMagnetoelasticEnergyDensity

func AddMagnetoelasticEnergyDensity(dst *data.Slice)

func AddMagnetoelasticField

func AddMagnetoelasticField(dst *data.Slice)

func AddSTTorque

func AddSTTorque(dst *data.Slice)

Adds the current spin transfer torque to dst

func AutoSave

func AutoSave(q Quantity, period float64)

Register quant to be auto-saved every period. period == 0 stops autosaving.

func AutoSnapshot

func AutoSnapshot(q Quantity, period float64)

Register quant to be auto-saved as image, every period.

func AverageOf

func AverageOf(q Quantity) []float64

func Break

func Break()

func CenterBubble

func CenterBubble()

This post-step function centers the simulation window on a bubble

func CenterWall

func CenterWall(magComp int)

This post-step function centers the simulation window on a domain wall between up-down (or down-up) domains (like in perpendicular media). E.g.:

PostStep(CenterPMAWall)

func CheckRecoverable

func CheckRecoverable(err error)

func Close

func Close()

Cleanly exits the simulation, assuring all output is flushed.

func CompileFile

func CompileFile(fname string) (*script.BlockStmt, error)

func Crop

func Crop(parent Quantity, x1, x2, y1, y2, z1, z2 int) *cropped

func CropLayer

func CropLayer(parent Quantity, layer int) *cropped

func CropRegion

func CropRegion(parent Quantity, region int) *cropped

Crop quantity to a box enclosing the given region. Used to output a region of interest, even if the region is non-rectangular.

func CropX

func CropX(parent Quantity, x1, x2 int) *cropped

func CropY

func CropY(parent Quantity, y1, y2 int) *cropped

func CropZ

func CropZ(parent Quantity, z1, z2 int) *cropped

func DeclConst

func DeclConst(name string, value float64, doc string)

Add a constant to the script world

func DeclFunc

func DeclFunc(name string, f interface{}, doc string)

Add a function to the script world

func DeclLValue

func DeclLValue(name string, value LValue, doc string)

Add an LValue to the script world. Assign to LValue invokes SetValue()

func DeclROnly

func DeclROnly(name string, value interface{}, doc string)

Add a read-only variable to the script world. It can be changed, but not by the user.

func DeclTVar

func DeclTVar(name string, value interface{}, doc string)

Hack for fixing the closure caveat: Defines "t", the time variable, handled specially by Fix()

func DeclVar

func DeclVar(name string, value interface{}, doc string)

Add a (pointer to) variable to the script world

func DefRegion

func DefRegion(id int, s Shape)

Define a region with id (0-255) to be inside the Shape.

func DefRegionCell

func DefRegionCell(id int, x, y, z int)

func DoOutput

func DoOutput()

Periodically called by run loop to save everything that's needed at this time.

func Download

func Download(q Quantity) *data.Slice

Download a quantity to host, or just return its data when already on host.

func EnableUnsafe

func EnableUnsafe()

func Eval

func Eval(code string)

func Eval1Line

func Eval1Line(code string) interface{}

func EvalFile

func EvalFile(code *script.BlockStmt)

evaluate code, exit on error (behavior for input files)

func EvalTo

func EvalTo(q interface {
	Slice() (*data.Slice, bool)
}, dst *data.Slice)

Temporary shim to fit Slice into EvalTo

func Exit

func Exit()

func Expect

func Expect(msg string, have, want, maxError float64)

Test if have lies within want +/- maxError, and print suited message.

func ExpectV

func ExpectV(msg string, have, want data.Vector, maxErr float64)

func Export

func Export(q interface {
	Name() string
	Unit() string
}, doc string)

func FifoRing

func FifoRing(length int) fifoRing

func Fprintln

func Fprintln(filename string, msg ...interface{})

Append msg to file. Used to write aggregated output of many simulations in one file.

func FreezeSpins

func FreezeSpins(dst *data.Slice)

func GUIAdd

func GUIAdd(name string, value interface{})

func GetAnisotropyEnergy

func GetAnisotropyEnergy() float64

Returns anisotropy energy in joules.

func GetBusy

func GetBusy() bool

func GetCustomEnergy

func GetCustomEnergy() float64

func GetDWxPos

func GetDWxPos() float64

func GetDemagEnergy

func GetDemagEnergy() float64

Returns the current demag energy in Joules.

func GetExchangeEnergy

func GetExchangeEnergy() float64

Returns the current exchange energy in Joules.

func GetMagnetoelasticEnergy

func GetMagnetoelasticEnergy() float64

Returns magneto-ell energy in joules.

func GetMagnetoelasticForceDensity

func GetMagnetoelasticForceDensity(dst *data.Slice)

func GetMaxAngle

func GetMaxAngle() float64

func GetMaxTorque

func GetMaxTorque() float64

func GetShiftPos

func GetShiftPos() float64

position of the window lab frame

func GetShiftYPos

func GetShiftYPos() float64

func GetThermalEnergy

func GetThermalEnergy() float64

func GetTopologicalCharge

func GetTopologicalCharge() float64

func GetTopologicalChargeLattice

func GetTopologicalChargeLattice() float64

func GetTotalEnergy

func GetTotalEnergy() float64

Returns the total energy in J.

func GetZeemanEnergy

func GetZeemanEnergy() float64

func GoServe

func GoServe(addr string) string

func Index2Coord

func Index2Coord(ix, iy, iz int) data.Vector

converts cell index to coordinate, internal coordinates

func InitAndClose

func InitAndClose() func()

Usage: in every Go input file, write:

func main(){
	defer InitAndClose()()
	// ...
}

This initialises the GPU, output directory, etc, and makes sure pending output will get flushed.

func InitIO

func InitIO(inputfile, od string, force bool)

SetOD sets the output directory where auto-saved files will be stored. The -o flag can also be used for this purpose.

func InjectAndWait

func InjectAndWait(task func())

inject code into engine and wait for it to complete.

func InterDind

func InterDind(region1, region2 int, value float64)

Sets the DMI interaction between region 1 and 2.

func InterExchange

func InterExchange(region1, region2 int, value float64)

Sets the exchange interaction between region 1 and 2.

func IsConst

func IsConst(e script.Expr) bool

checks if a script expression contains t (time)

func LoadFile

func LoadFile(fname string) *data.Slice

Read a magnetization state from .dump file.

func LogErr

func LogErr(msg ...interface{})

func LogIn

func LogIn(msg ...interface{})

func LogOut

func LogOut(msg ...interface{})

func LogUsedRefs

func LogUsedRefs()

func Madd

func Madd(a, b Quantity, fac1, fac2 float64) *mAddition

func Mesh

func Mesh() *data.Mesh

func MeshOf

func MeshOf(q Quantity) *data.Mesh

func MeshSize

func MeshSize() [3]int

func Minimize

func Minimize()

func NameOf

func NameOf(q Quantity) string

func NewScalarMask

func NewScalarMask(Nx, Ny, Nz int) *data.Slice

func NewSlice

func NewSlice(ncomp, Nx, Ny, Nz int) *data.Slice

Returns a new new slice (3D array) with given number of components and size.

func NewVectorMask

func NewVectorMask(Nx, Ny, Nz int) *data.Slice

func OD

func OD() string

func PostStep

func PostStep(f func())

Register function f to be called after every time step. Typically used, e.g., to manipulate the magnetization.

func Refer

func Refer(tag string)

func Relax

func Relax()

func RemoveCustomEnergies

func RemoveCustomEnergies()

Removes all customenergies

func RemoveCustomFields

func RemoveCustomFields()

Removes all customfields

func RemoveLRSurfaceCharge

func RemoveLRSurfaceCharge(region int, mxLeft, mxRight float64)

For a nanowire magnetized in-plane, with mx = mxLeft on the left end and mx = mxRight on the right end (both -1 or +1), add a B field needed to compensate for the surface charges on the left and right edges. This will mimic an infinitely long wire.

func Run

func Run(seconds float64)

Run the simulation for a number of seconds.

func RunInteractive

func RunInteractive()

Runs as long as browser is connected to gui.

func RunWhile

func RunWhile(condition func() bool)

Runs as long as condition returns true, saves output.

func SanityCheck

func SanityCheck()

func Save

func Save(q Quantity)

Save once, with auto file name

func SaveAs

func SaveAs(q Quantity, fname string)

Save under given file name (transparent async I/O).

func ScaleInterDind

func ScaleInterDind(region1, region2 int, scale float64)

Scales the DMI interaction between region 1 and 2.

func ScaleInterExchange

func ScaleInterExchange(region1, region2 int, scale float64)

Scales the heisenberg exchange interaction between region1 and 2. Scale = 1 means the harmonic mean over the regions of Aex.

func SetBusy

func SetBusy(b bool)

We set SetBusy(true) when the simulation is too busy too accept GUI input on Inject channel. E.g. during kernel init.

func SetCellSize

func SetCellSize(cx, cy, cz float64)

func SetDemagField

func SetDemagField(dst *data.Slice)

Sets dst to the current demag field

func SetEffectiveField

func SetEffectiveField(dst *data.Slice)

Sets dst to the current effective field, in Tesla. This is the sum of all effective field terms, like demag, exchange, ...

func SetGeom

func SetGeom(s Shape)

func SetGridSize

func SetGridSize(Nx, Ny, Nz int)

func SetLLTorque

func SetLLTorque(dst *data.Slice)

Sets dst to the current Landau-Lifshitz torque

func SetMFM

func SetMFM(dst *data.Slice)

func SetMFull

func SetMFull(dst *data.Slice)

Sets dst to the full (unnormalized) magnetization in A/m

func SetMesh

func SetMesh(Nx, Ny, Nz int, cellSizeX, cellSizeY, cellSizeZ float64, pbcx, pbcy, pbcz int)

Set the simulation mesh to Nx x Ny x Nz cells of given size. Can be set only once at the beginning of the simulation. TODO: dedup arguments from globals

func SetPBC

func SetPBC(nx, ny, nz int)

func SetPhi

func SetPhi(dst *data.Slice)

func SetSolver

func SetSolver(typ int)

func SetSpinAngle

func SetSpinAngle(dst *data.Slice)

func SetTheta

func SetTheta(dst *data.Slice)

func SetTopologicalChargeDensity

func SetTopologicalChargeDensity(dst *data.Slice)

func SetTopologicalChargeDensityLattice

func SetTopologicalChargeDensityLattice(dst *data.Slice)

func SetTorque

func SetTorque(dst *data.Slice)

Sets dst to the current total torque

func SetTotalEdens

func SetTotalEdens(dst *data.Slice)

Set dst to total energy density in J/m3

func Shift

func Shift(dx int)

shift the simulation window over dx cells in X direction

func SizeOf

func SizeOf(q Quantity) [3]int

func Snapshot

func Snapshot(q Quantity)

Save image once, with auto file name

func SnapshotAs

func SnapshotAs(q Quantity, fname string)

func Steps

func Steps(n int)

Run the simulation for a number of steps.

func TableAdd

func TableAdd(col Quantity)

func TableAddVariable

func TableAddVariable(x script.ScalarFunction, name, unit string)

func TableAutoSave

func TableAutoSave(period float64)

func TablePrint

func TablePrint(msg ...interface{})

func TableSave

func TableSave()

func ThermSeed

func ThermSeed(seed int)

Seeds the thermal noise generator

func UnitOf

func UnitOf(q Quantity) string

func ValueOf

func ValueOf(q Quantity) *data.Slice

func Vector

func Vector(x, y, z float64) data.Vector

Constructs a vector

func Voronoi

func Voronoi(grainsize float64, numRegions, seed int)

func Voronoi3d

func Voronoi3d(grainsize float64, startRegion int, numRegions int, inputShape Shape, seed int)

func YShift

func YShift(dy int)

shift the simulation window over dy cells in Y direction

Types

type BackwardEuler

type BackwardEuler struct {
	// contains filtered or unexported fields
}

Implicit midpoint solver.

func (*BackwardEuler) Free

func (s *BackwardEuler) Free()

func (*BackwardEuler) Step

func (s *BackwardEuler) Step()

Euler method, can be used as solver.Step.

type Config

type Config func(x, y, z float64) data.Vector

Magnetic configuration returns m vector for position (x,y,z)

func AntiVortex

func AntiVortex(circ, pol int) Config

func BlochSkyrmion

func BlochSkyrmion(charge, pol int) Config

func Conical

func Conical(q, coneDirection data.Vector, coneAngle float64) Config

Conical magnetization configuration. The magnetization rotates on a cone defined by coneAngle and coneDirection. q is the wave vector of the conical magnetization configuration. The magnetization is

m = u*cos(coneAngle) + sin(coneAngle)*( ua*cos(q*r) + ub*sin(q*r) )

with ua and ub unit vectors perpendicular to u (normalized coneDirection)

func Helical

func Helical(q data.Vector) Config

func NeelSkyrmion

func NeelSkyrmion(charge, pol int) Config

func RandomMag

func RandomMag() Config

Random initial magnetization.

func RandomMagSeed

func RandomMagSeed(seed int) Config

Random initial magnetization, generated from random seed.

func TwoDomain

func TwoDomain(mx1, my1, mz1, mxwall, mywall, mzwall, mx2, my2, mz2 float64) Config

Make a 2-domain configuration with domain wall. (mx1, my1, mz1) and (mx2, my2, mz2) are the magnetizations in the left and right domain, respectively. (mxwall, mywall, mzwall) is the magnetization in the wall. The wall is smoothed over a few cells so it will easily relax to its ground state. E.g.:

TwoDomain(1,0,0,  0,1,0,  -1,0,0) // head-to-head domains with transverse (Néel) wall
TwoDomain(1,0,0,  0,0,1,  -1,0,0) // head-to-head domains with perpendicular (Bloch) wall
TwoDomain(0,0,1,  1,0,0,   0,0,-1)// up-down domains with Bloch wall

func Uniform

func Uniform(mx, my, mz float64) Config

Returns a uniform magnetization state. E.g.:

M = Uniform(1, 0, 0)) // saturated along X

func Vortex

func Vortex(circ, pol int) Config

Make a vortex magnetization with given circulation and core polarization (+1 or -1). The core is smoothed over a few exchange lengths and should easily relax to its ground state.

func VortexWall

func VortexWall(mleft, mright float64, circ, pol int) Config

Make a vortex wall configuration.

func (Config) Add

func (c Config) Add(weight float64, other Config) Config

Returns a new magnetization equal to c + weight * other. E.g.:

Uniform(1, 0, 0).Add(0.2, RandomMag())

for a uniform state with 20% random distortion.

func (Config) RotZ

func (c Config) RotZ(θ float64) Config

Rotates the configuration around the Z-axis, over θ radians.

func (Config) Scale

func (c Config) Scale(sx, sy, sz float64) Config

Scale returns a scaled copy of configuration c.

func (Config) Transl

func (c Config) Transl(dx, dy, dz float64) Config

Transl returns a translated copy of configuration c. E.g.:

M = Vortex(1, 1).Transl(100e-9, 0, 0)  // vortex with center at x=100nm

type DataTable

type DataTable struct {
	// contains filtered or unexported fields
}

func (*DataTable) Add

func (t *DataTable) Add(output Quantity)

func (*DataTable) AddVariable

func (t *DataTable) AddVariable(x script.ScalarFunction, name, unit string)

func (*DataTable) Flush

func (t *DataTable) Flush() error

func (*DataTable) NComp

func (i *DataTable) NComp() int

func (*DataTable) Name

func (i *DataTable) Name() string

func (*DataTable) Println

func (t *DataTable) Println(msg ...interface{})

func (*DataTable) Save

func (t *DataTable) Save()

func (*DataTable) Unit

func (i *DataTable) Unit() string

func (*DataTable) Write

func (t *DataTable) Write(p []byte) (int, error)

type DerivedParam

type DerivedParam struct {
	// contains filtered or unexported fields
}

parameter derived from others (not directly settable). E.g.: Bsat derived from Msat

func NewDerivedParam

func NewDerivedParam(nComp int, parents []parent, updater func(*DerivedParam)) *DerivedParam

func (*DerivedParam) EvalTo

func (p *DerivedParam) EvalTo(dst *data.Slice)

uncompress the table to a full array in the dst Slice with parameter values per cell.

func (*DerivedParam) GetRegion

func (p *DerivedParam) GetRegion(r int) []float64

Get value in region r.

func (*DerivedParam) NComp

func (b *DerivedParam) NComp() int

func (*DerivedParam) Slice

func (p *DerivedParam) Slice() (*data.Slice, bool)

uncompress the table to a full array with parameter values per cell.

type Euler

type Euler struct{}

func (*Euler) Free

func (_ *Euler) Free()

func (*Euler) Step

func (_ *Euler) Step()

Euler method, can be used as solver.Step.

type Excitation

type Excitation struct {
	// contains filtered or unexported fields
}

An excitation, typically field or current, can be defined region-wise plus extra mask*multiplier terms.

func NewExcitation

func NewExcitation(name, unit, desc string) *Excitation

func (*Excitation) Add

func (e *Excitation) Add(mask *data.Slice, f script.ScalarFunction)

Add an extra mask*multiplier term to the excitation.

func (*Excitation) AddGo

func (e *Excitation) AddGo(mask *data.Slice, mul func() float64)

An Add(mask, f) equivalent for Go use

func (*Excitation) AddTo

func (e *Excitation) AddTo(dst *data.Slice)

func (*Excitation) Average

func (e *Excitation) Average() data.Vector

func (*Excitation) Comp

func (e *Excitation) Comp(c int) ScalarField

func (*Excitation) Eval

func (e *Excitation) Eval() interface{}

func (*Excitation) EvalTo

func (e *Excitation) EvalTo(dst *data.Slice)

func (*Excitation) InputType

func (e *Excitation) InputType() reflect.Type

func (*Excitation) IsUniform

func (e *Excitation) IsUniform() bool

func (*Excitation) MSlice

func (p *Excitation) MSlice() cuda.MSlice

func (*Excitation) Mesh

func (e *Excitation) Mesh() *data.Mesh

func (*Excitation) NComp

func (e *Excitation) NComp() int

func (*Excitation) Name

func (e *Excitation) Name() string

func (*Excitation) Region

func (e *Excitation) Region(r int) *vOneReg

func (*Excitation) RemoveExtraTerms

func (e *Excitation) RemoveExtraTerms()

After resizing the mesh, the extra terms don't fit the grid anymore and there is no reasonable way to resize them. So remove them and have the user re-add them.

func (*Excitation) Set

func (e *Excitation) Set(v data.Vector)

func (*Excitation) SetRegion

func (e *Excitation) SetRegion(region int, f script.VectorFunction)

func (*Excitation) SetRegionFn

func (e *Excitation) SetRegionFn(region int, f func() [3]float64)

func (*Excitation) SetValue

func (e *Excitation) SetValue(v interface{})

func (*Excitation) Slice

func (e *Excitation) Slice() (*data.Slice, bool)

func (*Excitation) Type

func (e *Excitation) Type() reflect.Type

func (*Excitation) Unit

func (e *Excitation) Unit() string

type FixedLayerPosition

type FixedLayerPosition int
const (
	FIXEDLAYER_TOP FixedLayerPosition = iota + 1
	FIXEDLAYER_BOTTOM
)

type Heun

type Heun struct{}

Adaptive Heun solver.

func (*Heun) Free

func (_ *Heun) Free()

func (*Heun) Step

func (_ *Heun) Step()

Adaptive Heun method, can be used as solver.Step

type Info

type Info interface {
	Name() string // number of components (scalar, vector, ...)
	Unit() string // name used for output file (e.g. "m")
	NComp() int   // unit, e.g. "A/m"
}

The Info interface defines the bare minimum methods a quantity must implement to be accessible for scripting and I/O.

type LValue

type LValue interface {
	SetValue(interface{}) // assigns a new value
	Eval() interface{}    // evaluate and return result (nil for void)
	Type() reflect.Type   // type that can be assigned and will be returned by Eval
}

LValue is settable

type Minimizer

type Minimizer struct {
	// contains filtered or unexported fields
}

func (*Minimizer) Free

func (mini *Minimizer) Free()

func (*Minimizer) Step

func (mini *Minimizer) Step()

type OutputFormat

type OutputFormat int
const (
	OVF1_TEXT OutputFormat = iota + 1
	OVF1_BINARY
	OVF2_TEXT
	OVF2_BINARY
	DUMP
)

type Param

type Param interface {
	NComp() int
	Name() string
	Unit() string

	IsUniform() bool
	// contains filtered or unexported methods
}

displayable quantity in GUI Parameters section

type Quantity

type Quantity interface {
	NComp() int
	EvalTo(dst *data.Slice)
}

Arbitrary physical quantity.

func Add

func Add(a, b Quantity) Quantity

func Const

func Const(v float64) Quantity

Const returns a constant (uniform) scalar quantity, that can be used to construct custom field terms.

func ConstVector

func ConstVector(x, y, z float64) Quantity

ConstVector returns a constant (uniform) vector quantity, that can be used to construct custom field terms.

func Cross

func Cross(a, b Quantity) Quantity

CrossProduct creates a new quantity that is the cross product of quantities a and b. E.g.:

CrossProct(&M, &B_ext)

func Div

func Div(a, b Quantity) Quantity

Div returns a new quantity that evaluates to the pointwise product a and b.

func Dot

func Dot(a, b Quantity) Quantity

DotProduct creates a new quantity that is the dot product of quantities a and b. E.g.:

DotProct(&M, &B_ext)

func Masked

func Masked(q Quantity, shape Shape) Quantity

Masks a quantity with a shape The shape will be only evaluated once on the mesh, and will be re-evaluated after mesh change, because otherwise too slow

func Mul

func Mul(a, b Quantity) Quantity

Mul returns a new quantity that evaluates to the pointwise product a and b.

func MulMV

func MulMV(Ax, Ay, Az, b Quantity) Quantity

MulMV returns a new Quantity that evaluates to the matrix-vector product (Ax·b, Ay·b, Az·b).

func Normalized

func Normalized(q Quantity) Quantity

Normalized returns a quantity that evaluates to the unit vector of q

func Shifted

func Shifted(q Quantity, dx, dy, dz int) Quantity

Shifted returns a new Quantity that evaluates to the original, shifted over dx, dy, dz cells.

type RK23

type RK23 struct {
	// contains filtered or unexported fields
}

Bogacki-Shampine solver. 3rd order, 3 evaluations per step, adaptive step.

http://en.wikipedia.org/wiki/Bogacki-Shampine_method

k1 = f(tn, yn)
k2 = f(tn + 1/2 h, yn + 1/2 h k1)
k3 = f(tn + 3/4 h, yn + 3/4 h k2)
y{n+1}  = yn + 2/9 h k1 + 1/3 h k2 + 4/9 h k3            // 3rd order
k4 = f(tn + h, y{n+1})
z{n+1} = yn + 7/24 h k1 + 1/4 h k2 + 1/3 h k3 + 1/8 h k4 // 2nd order

func (*RK23) Free

func (rk *RK23) Free()

func (*RK23) Step

func (rk *RK23) Step()

type RK4

type RK4 struct {
}

Classical 4th order RK solver.

func (*RK4) Free

func (_ *RK4) Free()

func (*RK4) Step

func (rk *RK4) Step()

type RK45DP

type RK45DP struct {
	// contains filtered or unexported fields
}

func (*RK45DP) Free

func (rk *RK45DP) Free()

func (*RK45DP) Step

func (rk *RK45DP) Step()

type RK56

type RK56 struct {
}

func (*RK56) Free

func (rk *RK56) Free()

func (*RK56) Step

func (rk *RK56) Step()

type Regions

type Regions struct {
	// contains filtered or unexported fields
}

stores the region index for each cell

func (*Regions) Average

func (r *Regions) Average() float64

func (*Regions) EvalTo

func (r *Regions) EvalTo(dst *data.Slice)

func (*Regions) GetCell

func (r *Regions) GetCell(ix, iy, iz int) int

func (*Regions) Gpu

func (r *Regions) Gpu() *cuda.Bytes

Get the region data on GPU

func (*Regions) HostArray

func (r *Regions) HostArray() [][][]byte

func (*Regions) HostList

func (r *Regions) HostList() []byte

func (*Regions) LoadFile

func (r *Regions) LoadFile(fname string)

Load regions from ovf file, use first component. Regions should be between 0 and 256

func (*Regions) Mesh

func (r *Regions) Mesh() *data.Mesh

func (*Regions) NComp

func (i *Regions) NComp() int

func (*Regions) Name

func (i *Regions) Name() string

func (*Regions) SetCell

func (r *Regions) SetCell(ix, iy, iz int, region int)

Set the region of one cell

func (*Regions) Slice

func (r *Regions) Slice() (*data.Slice, bool)

Get returns the regions as a slice of floats, so it can be output.

func (*Regions) Unit

func (i *Regions) Unit() string

type RegionwiseScalar

type RegionwiseScalar struct {
	// contains filtered or unexported fields
}

specialized param with 1 component

func NewScalarParam

func NewScalarParam(name, unit, desc string, children ...derived) *RegionwiseScalar

TODO: auto derived

func (*RegionwiseScalar) Average

func (p *RegionwiseScalar) Average() float64

func (*RegionwiseScalar) Eval

func (p *RegionwiseScalar) Eval() interface{}

func (*RegionwiseScalar) GetRegion

func (p *RegionwiseScalar) GetRegion(region int) float64

func (*RegionwiseScalar) InputType

func (p *RegionwiseScalar) InputType() reflect.Type

func (*RegionwiseScalar) IsUniform

func (p *RegionwiseScalar) IsUniform() bool

func (*RegionwiseScalar) MSlice

func (p *RegionwiseScalar) MSlice() cuda.MSlice

func (*RegionwiseScalar) Mesh

func (p *RegionwiseScalar) Mesh() *data.Mesh

func (*RegionwiseScalar) Name

func (p *RegionwiseScalar) Name() string

func (*RegionwiseScalar) Region

func (p *RegionwiseScalar) Region(r int) *sOneReg

func (*RegionwiseScalar) Set

func (p *RegionwiseScalar) Set(v float64)

func (*RegionwiseScalar) SetRegion

func (p *RegionwiseScalar) SetRegion(region int, f script.ScalarFunction)

func (*RegionwiseScalar) SetRegionFuncGo

func (p *RegionwiseScalar) SetRegionFuncGo(region int, f func() float64)

func (*RegionwiseScalar) SetRegionValueGo

func (p *RegionwiseScalar) SetRegionValueGo(region int, v float64)

func (*RegionwiseScalar) SetValue

func (p *RegionwiseScalar) SetValue(v interface{})

func (*RegionwiseScalar) Type

func (p *RegionwiseScalar) Type() reflect.Type

func (*RegionwiseScalar) Unit

func (p *RegionwiseScalar) Unit() string

type RegionwiseVector

type RegionwiseVector struct {
	// contains filtered or unexported fields
}

vector input parameter, settable by user

func NewVectorParam

func NewVectorParam(name, unit, desc string) *RegionwiseVector

func (*RegionwiseVector) Average

func (p *RegionwiseVector) Average() data.Vector

func (*RegionwiseVector) Comp

func (p *RegionwiseVector) Comp(c int) ScalarField

func (*RegionwiseVector) Eval

func (p *RegionwiseVector) Eval() interface{}

func (*RegionwiseVector) GetRegion

func (p *RegionwiseVector) GetRegion(region int) [3]float64

func (*RegionwiseVector) InputType

func (p *RegionwiseVector) InputType() reflect.Type

func (*RegionwiseVector) IsUniform

func (p *RegionwiseVector) IsUniform() bool

func (*RegionwiseVector) MSlice

func (p *RegionwiseVector) MSlice() cuda.MSlice

func (*RegionwiseVector) Mesh

func (p *RegionwiseVector) Mesh() *data.Mesh

func (*RegionwiseVector) Name

func (p *RegionwiseVector) Name() string

func (*RegionwiseVector) Region

func (p *RegionwiseVector) Region(r int) *vOneReg

func (*RegionwiseVector) SetRegion

func (p *RegionwiseVector) SetRegion(region int, f script.VectorFunction)

func (*RegionwiseVector) SetRegionFn

func (p *RegionwiseVector) SetRegionFn(region int, f func() [3]float64)

func (*RegionwiseVector) SetValue

func (p *RegionwiseVector) SetValue(v interface{})

func (*RegionwiseVector) Type

func (p *RegionwiseVector) Type() reflect.Type

func (*RegionwiseVector) Unit

func (p *RegionwiseVector) Unit() string

type ScalarExcitation

type ScalarExcitation struct {
	// contains filtered or unexported fields
}

An excitation, typically field or current, can be defined region-wise plus extra mask*multiplier terms.

func NewScalarExcitation

func NewScalarExcitation(name, unit, desc string) *ScalarExcitation

func (*ScalarExcitation) Add

func (e *ScalarExcitation) Add(mask *data.Slice, f script.ScalarFunction)

Add an extra mask*multiplier term to the excitation.

func (*ScalarExcitation) AddGo

func (e *ScalarExcitation) AddGo(mask *data.Slice, mul func() float64)

An Add(mask, f) equivalent for Go use

func (*ScalarExcitation) AddTo

func (e *ScalarExcitation) AddTo(dst *data.Slice)

func (*ScalarExcitation) Average

func (e *ScalarExcitation) Average() float64

func (*ScalarExcitation) Comp

func (e *ScalarExcitation) Comp(c int) ScalarField

func (*ScalarExcitation) Eval

func (e *ScalarExcitation) Eval() interface{}

func (*ScalarExcitation) EvalTo

func (e *ScalarExcitation) EvalTo(dst *data.Slice)

func (*ScalarExcitation) InputType

func (e *ScalarExcitation) InputType() reflect.Type

func (*ScalarExcitation) IsUniform

func (e *ScalarExcitation) IsUniform() bool

func (*ScalarExcitation) MSlice

func (p *ScalarExcitation) MSlice() cuda.MSlice

func (*ScalarExcitation) Mesh

func (e *ScalarExcitation) Mesh() *data.Mesh

func (*ScalarExcitation) NComp

func (e *ScalarExcitation) NComp() int

func (*ScalarExcitation) Name

func (e *ScalarExcitation) Name() string

func (*ScalarExcitation) Region

func (e *ScalarExcitation) Region(r int) *vOneReg

func (*ScalarExcitation) RemoveExtraTerms

func (e *ScalarExcitation) RemoveExtraTerms()

After resizing the mesh, the extra terms don't fit the grid anymore and there is no reasonable way to resize them. So remove them and have the user re-add them.

func (*ScalarExcitation) Set

func (e *ScalarExcitation) Set(v float64)

func (*ScalarExcitation) SetRegion

func (e *ScalarExcitation) SetRegion(region int, f script.ScalarFunction)

func (*ScalarExcitation) SetRegionFn

func (e *ScalarExcitation) SetRegionFn(region int, f func() [3]float64)

func (*ScalarExcitation) SetValue

func (e *ScalarExcitation) SetValue(v interface{})

func (*ScalarExcitation) Slice

func (e *ScalarExcitation) Slice() (*data.Slice, bool)

func (*ScalarExcitation) Type

func (e *ScalarExcitation) Type() reflect.Type

func (*ScalarExcitation) Unit

func (e *ScalarExcitation) Unit() string

type ScalarField

type ScalarField struct {
	Quantity
}

ScalarField enhances an outputField with methods specific to a space-dependent scalar quantity.

func AsScalarField

func AsScalarField(q Quantity) ScalarField

AsScalarField promotes a quantity to a ScalarField, enabling convenience methods particular to scalars.

func Comp

func Comp(parent Quantity, c int) ScalarField

Comp returns vector component c of the parent Quantity

func NewScalarField

func NewScalarField(name, unit, desc string, f func(dst *data.Slice)) ScalarField

NewVectorField constructs an outputable space-dependent scalar quantity whose value is provided by function f.

func (ScalarField) Average

func (s ScalarField) Average() float64

func (ScalarField) Name

func (s ScalarField) Name() string

func (ScalarField) Region

func (s ScalarField) Region(r int) ScalarField

func (ScalarField) Unit

func (s ScalarField) Unit() string

type ScalarValue

type ScalarValue struct {
	// contains filtered or unexported fields
}

ScalarValue enhances an outputValue with methods specific to a space-independent scalar quantity (e.g. total energy).

func NewScalarValue

func NewScalarValue(name, unit, desc string, f func() float64) *ScalarValue

NewScalarValue constructs an outputable space-independent scalar quantity whose value is provided by function f.

func (ScalarValue) Average

func (s ScalarValue) Average() float64

func (ScalarValue) EvalTo

func (g ScalarValue) EvalTo(dst *data.Slice)

func (ScalarValue) Get

func (s ScalarValue) Get() float64

type Shape

type Shape func(x, y, z float64) bool

geometrical shape for setting sample geometry

func Cell

func Cell(ix, iy, iz int) Shape

Single cell with given index

func Circle

func Circle(diam float64) Shape

func Cone

func Cone(diam, height float64) Shape

3D Cone with the vertex down.

func Cuboid

func Cuboid(sidex, sidey, sidez float64) Shape

3D Rectangular slab with given sides.

func Cylinder

func Cylinder(diam, height float64) Shape

cylinder along z.

func Ellipse

func Ellipse(diamx, diamy float64) Shape

func Ellipsoid

func Ellipsoid(diamx, diamy, diamz float64) Shape

Ellipsoid with given diameters

func GrainRoughness

func GrainRoughness(grainsize, zmin, zmax float64, seed int) Shape

func ImageShape

func ImageShape(fname string) Shape

func Layer

func Layer(index int) Shape

func Layers

func Layers(a, b int) Shape

Cell layers #a (inclusive) up to #b (exclusive).

func Rect

func Rect(sidex, sidey float64) Shape

2D Rectangle with given sides.

func Square

func Square(side float64) Shape

2D square with given side.

func Universe

func Universe() Shape

func XRange

func XRange(a, b float64) Shape

All cells with x-coordinate between a and b

func YRange

func YRange(a, b float64) Shape

All cells with y-coordinate between a and b

func ZRange

func ZRange(a, b float64) Shape

All cells with z-coordinate between a and b

func (Shape) Add

func (a Shape) Add(b Shape) Shape

Union of shapes a and b (logical OR).

func (Shape) Intersect

func (a Shape) Intersect(b Shape) Shape

Intersection of shapes a and b (logical AND).

func (Shape) Inverse

func (s Shape) Inverse() Shape

Inverse (outside) of shape (logical NOT).

func (Shape) Repeat

func (s Shape) Repeat(periodX, periodY, periodZ float64) Shape

Infinitely repeats the shape with given period in x, y, z. A period of 0 or infinity means no repetition.

func (Shape) RotX

func (s Shape) RotX(θ float64) Shape

Rotates the shape around the X-axis, over θ radians.

func (Shape) RotY

func (s Shape) RotY(θ float64) Shape

Rotates the shape around the Y-axis, over θ radians.

func (Shape) RotZ

func (s Shape) RotZ(θ float64) Shape

Rotates the shape around the Z-axis, over θ radians.

func (Shape) Scale

func (s Shape) Scale(sx, sy, sz float64) Shape

Scale returns a scaled copy of the shape.

func (Shape) Sub

func (a Shape) Sub(b Shape) Shape

Removes b from a (logical a AND NOT b)

func (Shape) Transl

func (s Shape) Transl(dx, dy, dz float64) Shape

Transl returns a translated copy of the shape.

func (Shape) Xor

func (a Shape) Xor(b Shape) Shape

Logical XOR of shapes a and b

type Stepper

type Stepper interface {
	Step() // take time step using solver globals
	Free() // free resources, if any (e.g.: RK23 previous torque)
}

Time stepper like Euler, Heun, RK23

type UserErr

type UserErr string

Special error that is not fatal when paniced on and called from GUI E.g.: try to set bad grid size: panic on UserErr, recover, print error, carry on.

func (UserErr) Error

func (e UserErr) Error() string

type VectorField

type VectorField struct {
	Quantity
}

VectorField enhances an outputField with methods specific to a space-dependent vector quantity.

func AsVectorField

func AsVectorField(q Quantity) VectorField

AsVectorField promotes a quantity to a VectorField, enabling convenience methods particular to vectors.

func NewVectorField

func NewVectorField(name, unit, desc string, f func(dst *data.Slice)) VectorField

NewVectorField constructs an outputable space-dependent vector quantity whose value is provided by function f.

func (VectorField) Average

func (v VectorField) Average() data.Vector

func (VectorField) Comp

func (v VectorField) Comp(c int) ScalarField

func (VectorField) HostCopy

func (v VectorField) HostCopy() *data.Slice

func (VectorField) Mesh

func (v VectorField) Mesh() *data.Mesh

func (VectorField) Name

func (v VectorField) Name() string

func (VectorField) Region

func (v VectorField) Region(r int) VectorField

func (VectorField) Unit

func (v VectorField) Unit() string

type VectorValue

type VectorValue struct {
	// contains filtered or unexported fields
}

VectorValue enhances an outputValue with methods specific to a space-independent vector quantity (e.g. averaged magnetization).

func NewVectorValue

func NewVectorValue(name, unit, desc string, f func() []float64) *VectorValue

NewVectorValue constructs an outputable space-independent vector quantity whose value is provided by function f.

func (*VectorValue) Average

func (v *VectorValue) Average() data.Vector

func (VectorValue) EvalTo

func (g VectorValue) EvalTo(dst *data.Slice)

func (*VectorValue) Get

func (v *VectorValue) Get() data.Vector

Jump to

Keyboard shortcuts

? : This menu
/ : Search site
f or F : Jump to
y or Y : Canonical URL