Channeling mixed Néel-Bloch domain walls with iDMI#

In this example, we reproduce the dispersion of magnons propagating along a waveguide housing a mixed Néel-Bloch domain wall, stabilized by interfacial Dzyaloshinskii-Moryia interaction (iDMI), as reported by Garcia-Sanchez et al. in Ref. 1. This example was adapted from a notebook kindly provided by Joo-Von Kim.

Creating the sample#

First, we create the magnetic sample. To model spin waves propagating along a domain wall, we can use the waveguide geometry.

[1]:
# Define sample
import tetrax as tx
import numpy as np

width = 250 # nm
thickness = 1 # nm

sample = tx.Sample(
    tx.geometries.waveguide.rectangular(width,
                                        thickness,
                                        cell_size_width=1.0,
                                        cell_size_thickness=.5
                                       ),
    name="DW_waveguide_iDMI"
)

# Micromagnetic parameters taken from Ref. 1
sample.material["Msat"] = 1e6
sample.material["Aex"] = 15e-12
sample.material["Ku1"] = 1e6
sample.material["e_u"] = (0, 1, 0)
sample.material["Didmi"] = 1.5e-3
sample.material["gamma"] = 1.760859644e11

Preparing the equilibrium state#

We initialize a left-handed domain wall at the center of the waveguide in the following way, using the theoretical domain-wall width given by

\[\Delta = \sqrt{\frac{A_\mathrm{ex}}{K_\mathrm{eff}}} \simeq 6.35 \; \mathrm{nm},\]

where

\[K_\mathrm{eff} = K_\mathrm{u1} - \frac{1}{2}\mu_0 M_s^2.\]

A small magnetic field is applied to help the calculations.

[2]:
def mixed_wall(xyz, thickness=10, x0=0, phi=1.5708):
    x = xyz.x
    theta = 2 * np.arctan(np.exp((x - x0) / thickness))
    mz = np.sin(theta) * np.cos(phi)
    mx = np.sin(theta) * np.sin(phi)
    my = np.cos(theta)

    return np.array([mx, my, mz]).T

sample.mag = mixed_wall(sample.xyz, thickness=12.7, x0=0.0, phi=1.5*np.pi)
sample.Bext = (0,0,1e-5)

Systems with DMI require special attention when equilibrating them. Simple energy minimization using the relax() experiment is insufficient and can lead to shifts in the spin-wave dispersion. A more precise equilibrium state can be calculated by evolving the system according to the over-damped Landau-Lifshitz-Gilbert equation, minimizing the overall torque. This is implemented in the relax_dynamic() experiment. Torque minimization usually takes longer but can be accelerated by precursing it with an energy minimization.

[3]:
relaxation_energy = tx.experiments.relax(sample, tolerance=1e-10)
relaxation_torque = tx.experiments.relax_dynamic(sample, tolerance=1e-10)

relaxation_torque.plot(renderer="notebook")
Minimizing energy in using 'L-BFGS-B' (tolerance 1e-10) ...
energy length density: -9.094820983684039e-11 J/m, <mag.x> = -0.07, <mag.y> = -0.00, <mag.z> = 0.00
Success!

Minimizing torque with over-damped LLG (alpha=1, tolerance=1e-10) ...
t = 1.5477691204528555e-11 (s), dm/dt = 114444063.51931995 [Hz], <mag.x> = -0.07, <mag.y> = -0.00, <mag.z> = 0.00
Success!

The torque still significantly decreases, although the average magnetization components barely change. We have a look at the equilibrium by taking a line scan along the transversal (\(x\)) direction of the waveguide and comparing it with MuMax3 (also using torque minimization). The MuMax3 data will be loaded using the python package `pooch <>`__, which allows to read csv files from

[4]:
mumax_eq = tx.fetch_reference_data("example_equilibrium_mumax3_neel_wall_iDMI")

scan_tx, curve_tx = sample.scan_along_curve(sample.mag, ((-width/2,0,0),(width/2,0,0)),num_points=500, return_curve=True)

from plotly import graph_objects as go

fig = go.Figure()

for i, name in enumerate(["x", "y", "z"]):
    fig.add_trace(
                go.Scatter(
                    x=mumax_eq["Points:1"]/255*width-width/2,
                    y=mumax_eq[f'm/1:{i}'],
                    mode="lines",
                    line = dict(color="black",width=1),
                    legendgroup="0",
                    legendgrouptitle_text="MuMax",
                    name=f"<i>m<sub>{name}</sub><i>"
                )
    )

for comp, name in zip(scan_tx.components, scan_tx.component_names):
    fig.add_trace(
                    go.Scatter(
                        x=curve_tx.x,
                        y=comp,
                        mode="lines",
                        line = dict(color="black",width=2, dash="dot"),
                        legendgroup="2",
                        legendgrouptitle_text="TetraX",
                        name=f"<i>m<sub>{name}</sub><i>"
                    )
        )



fig.update_xaxes(title_text="<i>y</i> (nm)", exponentformat="power")
fig.update_yaxes(title_text="Equilibrium magnetization <i><b>m</b></i><sub>0</sub>", exponentformat="power")

fig.update_layout(
        template="simple_white",
        height=600,
        width=600,
        hoverlabel={"bordercolor": "rgba(255,255,255,1)"},
    )

fig.show(renderer="notebook")