// Finds the average exchange strength around each cell, for debugging. func ExchangeDecode(dst *data.Slice, Aex_red SymmLUT, regions *Bytes, mesh *data.Mesh) { c := mesh.CellSize() wx := float32(2 * 1e-18 / (c[X] * c[X])) wy := float32(2 * 1e-18 / (c[Y] * c[Y])) wz := float32(2 * 1e-18 / (c[Z] * c[Z])) N := mesh.Size() pbc := mesh.PBC_code() cfg := make3DConf(N) k_exchangedecode_async(dst.DevPtr(0), unsafe.Pointer(Aex_red), regions.Ptr, wx, wy, wz, N[X], N[Y], N[Z], pbc, cfg) }
// Add effective field of Dzyaloshinskii-Moriya interaction to Beff (Tesla). // According to Bagdanov and Röβler, PRL 87, 3, 2001. eq.8 (out-of-plane symmetry breaking). // See dmi.cu func AddDMI(Beff *data.Slice, m *data.Slice, Aex_red, Dex_red SymmLUT, regions *Bytes, mesh *data.Mesh) { cellsize := mesh.CellSize() N := Beff.Size() util.Argument(m.Size() == N) cfg := make3DConf(N) k_adddmi_async(Beff.DevPtr(X), Beff.DevPtr(Y), Beff.DevPtr(Z), m.DevPtr(X), m.DevPtr(Y), m.DevPtr(Z), unsafe.Pointer(Aex_red), unsafe.Pointer(Dex_red), regions.Ptr, float32(cellsize[X]*1e9), float32(cellsize[Y]*1e9), float32(cellsize[Z]*1e9), N[X], N[Y], N[Z], mesh.PBC_code(), cfg) }
// Set s to the toplogogical charge density s = m · (m/∂x ❌ ∂m/∂y) // See topologicalcharge.cu func SetTopologicalCharge(s *data.Slice, m *data.Slice, mesh *data.Mesh) { cellsize := mesh.CellSize() N := s.Size() util.Argument(m.Size() == N) cfg := make3DConf(N) icxcy := float32(1.0 / (cellsize[X] * cellsize[Y])) k_settopologicalcharge_async(s.DevPtr(X), m.DevPtr(X), m.DevPtr(Y), m.DevPtr(Z), icxcy, N[X], N[Y], N[Z], mesh.PBC_code(), cfg) }
// Add Zhang-Li ST torque (Tesla) to torque. // see zhangli.cu func AddZhangLiTorque(torque, m, J *data.Slice, bsat, alpha, xi, pol LUTPtr, regions *Bytes, mesh *data.Mesh) { c := mesh.CellSize() N := mesh.Size() cfg := make3DConf(N) k_addzhanglitorque_async(torque.DevPtr(X), torque.DevPtr(Y), torque.DevPtr(Z), m.DevPtr(X), m.DevPtr(Y), m.DevPtr(Z), J.DevPtr(X), J.DevPtr(Y), J.DevPtr(Z), float32(c[X]), float32(c[Y]), float32(c[Z]), unsafe.Pointer(bsat), unsafe.Pointer(alpha), unsafe.Pointer(xi), unsafe.Pointer(pol), regions.Ptr, N[X], N[Y], N[Z], mesh.PBC_code(), cfg) }
// Add exchange field to Beff. // m: normalized magnetization // B: effective field in Tesla // Aex_red: Aex / (Msat * 1e18 m2) // see exchange.cu func AddExchange(B, m *data.Slice, Aex_red SymmLUT, regions *Bytes, mesh *data.Mesh) { c := mesh.CellSize() wx := float32(2 * 1e-18 / (c[X] * c[X])) wy := float32(2 * 1e-18 / (c[Y] * c[Y])) wz := float32(2 * 1e-18 / (c[Z] * c[Z])) N := mesh.Size() pbc := mesh.PBC_code() cfg := make3DConf(N) k_addexchange_async(B.DevPtr(X), B.DevPtr(Y), B.DevPtr(Z), m.DevPtr(X), m.DevPtr(Y), m.DevPtr(Z), unsafe.Pointer(Aex_red), regions.Ptr, wx, wy, wz, N[X], N[Y], N[Z], pbc, cfg) }
// Add interlayer exchange field to Beff. // see interlayer.cu func AddInterlayerExchange(Beff, m *data.Slice, J1_red, J2_red, toplayer, bottomlayer LUTPtr, direc LUTPtrs, regions *Bytes, mesh *data.Mesh) { cellsize := mesh.CellSize() N := Beff.Size() util.Argument(m.Size() == N) cfg := make3DConf(N) k_addinterlayerexchange_async(Beff.DevPtr(X), Beff.DevPtr(Y), Beff.DevPtr(Z), m.DevPtr(X), m.DevPtr(Y), m.DevPtr(Z), unsafe.Pointer(J1_red), unsafe.Pointer(J2_red), unsafe.Pointer(toplayer), unsafe.Pointer(bottomlayer), direc[X], direc[Y], direc[Z], float32(cellsize[X])*1e9, float32(cellsize[Y])*1e9, float32(cellsize[Z])*1e9, N[X], N[Y], N[Z], regions.Ptr, cfg) }
func compensateLRSurfaceCharges(m *data.Mesh, mxLeft, mxRight float64, bsat float64) *data.Slice { h := data.NewSlice(3, m.Size()) H := h.Vectors() world := m.WorldSize() cell := m.CellSize() size := m.Size() q := cell[Z] * cell[Y] * bsat q1 := q * mxLeft q2 := q * (-mxRight) prog, maxProg := 0, (size[Z]+1)*(size[Y]+1) // surface loop (source) for I := 0; I < size[Z]; I++ { for J := 0; J < size[Y]; J++ { prog++ util.Progress(prog, maxProg, "removing surface charges") y := (float64(J) + 0.5) * cell[Y] z := (float64(I) + 0.5) * cell[Z] source1 := [3]float64{0, y, z} // left surface source source2 := [3]float64{world[X], y, z} // right surface source // volume loop (destination) for iz := range H[0] { for iy := range H[0][iz] { for ix := range H[0][iz][iy] { dst := [3]float64{ // destination coordinate (float64(ix) + 0.5) * cell[X], (float64(iy) + 0.5) * cell[Y], (float64(iz) + 0.5) * cell[Z]} h1 := hfield(q1, source1, dst) h2 := hfield(q2, source2, dst) // add this surface charges' field to grand total for c := 0; c < 3; c++ { H[c][iz][iy][ix] += float32(h1[c] + h2[c]) } } } } } } return h }
// Add Zhang-Li ST torque (Tesla) to torque. // see zhangli.cu func AddZhangLiTorque(torque, m *data.Slice, Msat, J, alpha, xi, pol MSlice, mesh *data.Mesh) { c := mesh.CellSize() N := mesh.Size() cfg := make3DConf(N) k_addzhanglitorque2_async( torque.DevPtr(X), torque.DevPtr(Y), torque.DevPtr(Z), m.DevPtr(X), m.DevPtr(Y), m.DevPtr(Z), Msat.DevPtr(0), Msat.Mul(0), J.DevPtr(X), J.Mul(X), J.DevPtr(Y), J.Mul(Y), J.DevPtr(Z), J.Mul(Z), alpha.DevPtr(0), alpha.Mul(0), xi.DevPtr(0), xi.Mul(0), pol.DevPtr(0), pol.Mul(0), float32(c[X]), float32(c[Y]), float32(c[Z]), N[X], N[Y], N[Z], mesh.PBC_code(), cfg) }
// Kernel for the vertical derivative of the force on an MFM tip due to mx, my, mz. // This is the 2nd derivative of the energy w.r.t. z. func MFMKernel(mesh *d.Mesh, lift, tipsize float64) (kernel [3]*d.Slice) { const TipCharge = 1 / Mu0 // tip charge const Δ = 1e-9 // tip oscillation, take 2nd derivative over this distance util.AssertMsg(lift > 0, "MFM tip crashed into sample, please lift the new one higher") { // Kernel mesh is 2x larger than input, instead in case of PBC pbc := mesh.PBC() sz := padSize(mesh.Size(), pbc) cs := mesh.CellSize() mesh = d.NewMesh(sz[X], sz[Y], sz[Z], cs[X], cs[Y], cs[Z], pbc[:]...) } // Shorthand size := mesh.Size() pbc := mesh.PBC() cellsize := mesh.CellSize() volume := cellsize[X] * cellsize[Y] * cellsize[Z] fmt.Println("calculating MFM kernel") // Sanity check { util.Assert(size[Z] >= 1 && size[Y] >= 2 && size[X] >= 2) util.Assert(cellsize[X] > 0 && cellsize[Y] > 0 && cellsize[Z] > 0) util.AssertMsg(size[X]%2 == 0 && size[Y]%2 == 0, "Even kernel size needed") if size[Z] > 1 { util.AssertMsg(size[Z]%2 == 0, "Even kernel size needed") } } // Allocate only upper diagonal part. The rest is symmetric due to reciprocity. var K [3][][][]float32 for i := 0; i < 3; i++ { kernel[i] = d.NewSlice(1, mesh.Size()) K[i] = kernel[i].Scalars() } r1, r2 := kernelRanges(size, pbc) progress, progmax := 0, (1+r2[Y]-r1[Y])*(1+r2[Z]-r1[Z]) for iz := r1[Z]; iz <= r2[Z]; iz++ { zw := wrap(iz, size[Z]) z := float64(iz) * cellsize[Z] for iy := r1[Y]; iy <= r2[Y]; iy++ { yw := wrap(iy, size[Y]) y := float64(iy) * cellsize[Y] progress++ util.Progress(progress, progmax, "Calculating MFM kernel") for ix := r1[X]; ix <= r2[X]; ix++ { x := float64(ix) * cellsize[X] xw := wrap(ix, size[X]) for s := 0; s < 3; s++ { // source index Ksxyz m := d.Vector{0, 0, 0} m[s] = 1 var E [3]float64 // 3 energies for 2nd derivative for i := -1; i <= 1; i++ { I := float64(i) R := d.Vector{-x, -y, z - (lift + (I * Δ))} r := R.Len() B := R.Mul(TipCharge / (4 * math.Pi * r * r * r)) R = d.Vector{-x, -y, z - (lift + tipsize + (I * Δ))} r = R.Len() B = B.Add(R.Mul(-TipCharge / (4 * math.Pi * r * r * r))) E[i+1] = B.Dot(m) * volume // i=-1 stored in E[0] } dFdz_tip := ((E[0] - E[1]) + (E[2] - E[1])) / (Δ * Δ) // dFz/dz = d2E/dz2 K[s][zw][yw][xw] += float32(dFdz_tip) // += needed in case of PBC } } } } return kernel }