Steepest Descent minimiser

Notebook showing how to use the implementation of an optimised version of the Steepest Descent (SD) algorithm in Fidimag, for an atomistic system.

This SD algorithm is based on Exl et al. Journal of Applied Physics 115, 17D118 (2014) https://aip.scitation.org/doi/10.1063/1.4862839

The following links are also relevant:

The SD driver inherits from the minimiser base class in fidimag/common/ and can be specified in the driver option of the Simulation class

[1]:
import matplotlib.pyplot as plt
import fidimag
import fidimag.common.constant as C
import numpy as np
%matplotlib inline
[2]:
# import imp
# imp.reload(fidimag)

1D example

We start defining parameters for an atomistic simulation

[3]:
# System parameters
L = 100

# Some atomistic parameters
J = 5.88 * C.meV
D = 1.56 * C.meV
Ku = 0.41 * C.meV
mus = 3 * C.mu_B

# Lattice constants (in nm)
a = 0.2715
az = 0.408

# Magnetic field in Tesla
B = 2

# Free electron gyrom ratio
gamma = 1.76e11

Define the mesh of the system

[4]:
nx, ny, nz = 100, 1, 1
dx, dy, dz = a, a, az

mesh = fidimag.common.CuboidMesh(nx=nx, ny=ny, nz=nz, dx=dx, dy=dy, dz=dz,
                                periodicity=(False, False, False),
                                unit_length=1e-9)

Steepest Descent

Here we set the minimiser by specifying it in the driver argument in the Simulation class:

[5]:
sim = fidimag.atomistic.Sim(mesh, name='one_dim_SD', driver='steepest_descent')

# Define the magnetisation
sim.set_mu_s(mus)

# Add the magnetic interactions
sim.add(fidimag.atomistic.Exchange(J))
sim.add(fidimag.atomistic.Anisotropy(Ku, axis=(0, 0, 1)))
sim.add(fidimag.atomistic.DMI(D, dmi_type='interfacial'))
sim.add(fidimag.atomistic.Zeeman((0, 0, B)))

xs = mesh.coordinates[:, 0]
centre_x = (xs.max() + xs.min()) * 0.5 + xs.min()

def m_initial(r):
    x, y, z = r[0], r[1], r[2]
    if np.abs(x - centre_x) < 2:
        return (0, 0.1, -.9)
    else:
        return (0, 0.1, .9)

# sim.set_m((0.1, 0, 0.9))
sim.set_m(m_initial)

The initial configuration showing \(m_{z}\), we want a domain wall after relaxation:

[6]:
plt.plot(sim.mesh.coordinates[:, 0], sim.spin.reshape(-1, 3)[:, 2], 'o-')
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_12_0.png

Here we relax the system with the steepest descent. In this case, the relevant parameter to stop the minimisation is the dm difference of the magnetisation with the previous step:

[7]:
sim.driver.tmax = 1e-1
sim.driver.minimise(max_steps=20000, stopping_dm=1e-10, initial_t_step=1e-2)
#0    max_tau=0.01     max_dm=0.191
#711  max_tau=0.01     max_dm=9.75e-11

The performance of the algorithm can be tuned by modifying the tmax and tmin tolerances in the driver. By defalt, tmax is around 0.01 (MuMax3 uses this magnitude) and tmin must be significantly small. These values seem to work reasonably well for most atomistic simulations:

[8]:
# sim.driver.tmax =
# sim.driver.tmin =

We finally obtain the \(360^\circ\) domain wall:

[9]:
plt.plot(sim.mesh.coordinates[:, 0], sim.spin.reshape(-1, 3)[:, 0], 'o-')
plt.plot(sim.mesh.coordinates[:, 0], sim.spin.reshape(-1, 3)[:, 2], 'o-')
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_18_0.png

LLG

We can compare the previous result with the LLG driver, which is widely tested:

[10]:
sim = fidimag.atomistic.Sim(mesh, name='one_dim', driver='llg')

# Define the magnetisation
sim.set_mu_s(mus)

# Add the magnetic interactions
sim.add(fidimag.atomistic.Exchange(J))
sim.add(fidimag.atomistic.Anisotropy(Ku, axis=(0, 0, 1)))
sim.add(fidimag.atomistic.DMI(D, dmi_type='interfacial'))
sim.add(fidimag.atomistic.Zeeman((0, 0, B)))

xs = mesh.coordinates[:, 0]
centre_x = (xs.max() + xs.min()) * 0.5 + xs.min()

def m_initial(r):
    x, y, z = r[0], r[1], r[2]
    if np.abs(x - centre_x) < 2:
        return (0, 0.1, -.9)
    else:
        return (0, 0.1, .9)

# sim.set_m((0.1, 0, 0.9))
sim.set_m(m_initial)
[11]:
plt.plot(sim.mesh.coordinates[:, 0], sim.spin.reshape(-1, 3)[:, 2], 'o-')
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_22_0.png
[12]:
sim.driver.do_precession = False
sim.relax()
#1    t=1e-11    dt=1e-11 max_dmdt=1.77
#2    t=3.53e-10 dt=3.43e-10 max_dmdt=1.77
#3    t=3.79e-09 dt=3.43e-09 max_dmdt=1.77
#4    t=3.81e-08 dt=3.43e-08 max_dmdt=1.77
#5    t=3.81e-07 dt=3.43e-07 max_dmdt=1.77
#6    t=3.81e-06 dt=3.43e-06 max_dmdt=1.77
#7    t=2.51e-05 dt=2.13e-05 max_dmdt=1.77
#8    t=8.19e-05 dt=5.67e-05 max_dmdt=1.77
#9    t=0.000197 dt=0.000115 max_dmdt=1.77
#10   t=0.000417 dt=0.00022 max_dmdt=1.77
#11   t=0.000637 dt=0.00022 max_dmdt=1.77
#12   t=0.000856 dt=0.00022 max_dmdt=1.77
#13   t=0.00108  dt=0.00022 max_dmdt=1.77
#14   t=0.00161  dt=0.000532 max_dmdt=1.78
#15   t=0.00214  dt=0.000532 max_dmdt=1.78
#16   t=0.00267  dt=0.000532 max_dmdt=1.78
#17   t=0.0032   dt=0.000532 max_dmdt=1.78
#18   t=0.00374  dt=0.000532 max_dmdt=1.79
#19   t=0.00483  dt=0.00109 max_dmdt=1.79
#20   t=0.00591  dt=0.00109 max_dmdt=1.8
#21   t=0.007    dt=0.00109 max_dmdt=1.8
#22   t=0.00809  dt=0.00109 max_dmdt=1.81
#23   t=0.00918  dt=0.00109 max_dmdt=1.82
#24   t=0.0109   dt=0.0017 max_dmdt=1.82
#25   t=0.0126   dt=0.0017 max_dmdt=1.83
#26   t=0.0143   dt=0.0017 max_dmdt=1.84
#27   t=0.016    dt=0.0017 max_dmdt=1.85
#28   t=0.0177   dt=0.0017 max_dmdt=1.86
#29   t=0.0203   dt=0.00262 max_dmdt=1.87
#30   t=0.0229   dt=0.00262 max_dmdt=1.88
#31   t=0.027    dt=0.00405 max_dmdt=1.9
#32   t=0.031    dt=0.00405 max_dmdt=1.92
#33   t=0.0351   dt=0.00405 max_dmdt=1.95
#34   t=0.0391   dt=0.00405 max_dmdt=1.97
#35   t=0.0432   dt=0.00405 max_dmdt=1.99
#36   t=0.0472   dt=0.00405 max_dmdt=2.01
#37   t=0.0513   dt=0.00405 max_dmdt=2.03
#38   t=0.0553   dt=0.00405 max_dmdt=2.05
#39   t=0.0594   dt=0.00405 max_dmdt=2.07
#40   t=0.0658   dt=0.0064 max_dmdt=2.1
#41   t=0.0722   dt=0.0064 max_dmdt=2.13
#42   t=0.0786   dt=0.0064 max_dmdt=2.16
#43   t=0.085    dt=0.0064 max_dmdt=2.19
#44   t=0.0913   dt=0.0064 max_dmdt=2.22
#45   t=0.0977   dt=0.0064 max_dmdt=2.24
#46   t=0.104    dt=0.0064 max_dmdt=2.27
#47   t=0.111    dt=0.0064 max_dmdt=2.29
#48   t=0.117    dt=0.0064 max_dmdt=2.31
#49   t=0.123    dt=0.0064 max_dmdt=2.32
#50   t=0.13     dt=0.0064 max_dmdt=2.34
#51   t=0.136    dt=0.0064 max_dmdt=2.35
#52   t=0.143    dt=0.0064 max_dmdt=2.36
#53   t=0.149    dt=0.0064 max_dmdt=2.37
#54   t=0.155    dt=0.0064 max_dmdt=2.37
#55   t=0.162    dt=0.0064 max_dmdt=2.37
#56   t=0.168    dt=0.0064 max_dmdt=2.37
#57   t=0.174    dt=0.0064 max_dmdt=2.36
#58   t=0.181    dt=0.0064 max_dmdt=2.35
#59   t=0.187    dt=0.0064 max_dmdt=2.34
#60   t=0.194    dt=0.0064 max_dmdt=2.32
#61   t=0.2      dt=0.0064 max_dmdt=2.3
#62   t=0.206    dt=0.0064 max_dmdt=2.28
#63   t=0.213    dt=0.0064 max_dmdt=2.25
#64   t=0.219    dt=0.0064 max_dmdt=2.22
#65   t=0.226    dt=0.0064 max_dmdt=2.19
#66   t=0.232    dt=0.0064 max_dmdt=2.16
#67   t=0.238    dt=0.0064 max_dmdt=2.12
#68   t=0.245    dt=0.0064 max_dmdt=2.08
#69   t=0.251    dt=0.0064 max_dmdt=2.04
#70   t=0.258    dt=0.0064 max_dmdt=2
#71   t=0.264    dt=0.0064 max_dmdt=1.95
#72   t=0.27     dt=0.0064 max_dmdt=1.91
#73   t=0.277    dt=0.0064 max_dmdt=1.86
#74   t=0.283    dt=0.0064 max_dmdt=1.81
#75   t=0.29     dt=0.0064 max_dmdt=1.77
#76   t=0.296    dt=0.0064 max_dmdt=1.72
#77   t=0.302    dt=0.0064 max_dmdt=1.67
#78   t=0.309    dt=0.0064 max_dmdt=1.62
#79   t=0.315    dt=0.0064 max_dmdt=1.6
#80   t=0.322    dt=0.0064 max_dmdt=1.59
#81   t=0.328    dt=0.0064 max_dmdt=1.58
#82   t=0.334    dt=0.0064 max_dmdt=1.57
#83   t=0.341    dt=0.0064 max_dmdt=1.55
#84   t=0.347    dt=0.0064 max_dmdt=1.53
#85   t=0.354    dt=0.0064 max_dmdt=1.52
#86   t=0.36     dt=0.0064 max_dmdt=1.5
#87   t=0.366    dt=0.0064 max_dmdt=1.47
#88   t=0.376    dt=0.00964 max_dmdt=1.45
#89   t=0.386    dt=0.00964 max_dmdt=1.41
#90   t=0.395    dt=0.00964 max_dmdt=1.37
#91   t=0.405    dt=0.00964 max_dmdt=1.34
#92   t=0.415    dt=0.00964 max_dmdt=1.3
#93   t=0.424    dt=0.00964 max_dmdt=1.26
#94   t=0.434    dt=0.00964 max_dmdt=1.22
#95   t=0.443    dt=0.00964 max_dmdt=1.18
#96   t=0.453    dt=0.00964 max_dmdt=1.14
#97   t=0.463    dt=0.00964 max_dmdt=1.11
#98   t=0.472    dt=0.00964 max_dmdt=1.07
#99   t=0.482    dt=0.00964 max_dmdt=1.03
#100  t=0.492    dt=0.00964 max_dmdt=0.997
#101  t=0.501    dt=0.00964 max_dmdt=0.963
#102  t=0.511    dt=0.00964 max_dmdt=0.93
#103  t=0.521    dt=0.00964 max_dmdt=0.907
#104  t=0.53     dt=0.00964 max_dmdt=0.89
#105  t=0.54     dt=0.00964 max_dmdt=0.874
#106  t=0.549    dt=0.00964 max_dmdt=0.857
#107  t=0.559    dt=0.00964 max_dmdt=0.84
#108  t=0.569    dt=0.00964 max_dmdt=0.823
#109  t=0.578    dt=0.00964 max_dmdt=0.807
#110  t=0.588    dt=0.00964 max_dmdt=0.791
#111  t=0.598    dt=0.00964 max_dmdt=0.775
#112  t=0.607    dt=0.00964 max_dmdt=0.759
#113  t=0.617    dt=0.00964 max_dmdt=0.744
#114  t=0.627    dt=0.00964 max_dmdt=0.728
#115  t=0.636    dt=0.00964 max_dmdt=0.713
#116  t=0.651    dt=0.0145 max_dmdt=0.695
#117  t=0.665    dt=0.0145 max_dmdt=0.674
#118  t=0.68     dt=0.0145 max_dmdt=0.653
#119  t=0.694    dt=0.0145 max_dmdt=0.633
#120  t=0.708    dt=0.0145 max_dmdt=0.614
#121  t=0.723    dt=0.0145 max_dmdt=0.595
#122  t=0.737    dt=0.0145 max_dmdt=0.577
#123  t=0.752    dt=0.0145 max_dmdt=0.56
#124  t=0.766    dt=0.0145 max_dmdt=0.543
#125  t=0.781    dt=0.0145 max_dmdt=0.527
#126  t=0.795    dt=0.0145 max_dmdt=0.512
#127  t=0.81     dt=0.0145 max_dmdt=0.502
#128  t=0.824    dt=0.0145 max_dmdt=0.495
#129  t=0.839    dt=0.0145 max_dmdt=0.489
#130  t=0.853    dt=0.0145 max_dmdt=0.482
#131  t=0.867    dt=0.0145 max_dmdt=0.475
#132  t=0.882    dt=0.0145 max_dmdt=0.469
#133  t=0.896    dt=0.0145 max_dmdt=0.462
#134  t=0.911    dt=0.0145 max_dmdt=0.456
#135  t=0.925    dt=0.0145 max_dmdt=0.449
#136  t=0.94     dt=0.0145 max_dmdt=0.443
#137  t=0.954    dt=0.0145 max_dmdt=0.436
#138  t=0.969    dt=0.0145 max_dmdt=0.43
#139  t=0.983    dt=0.0145 max_dmdt=0.424
#140  t=1.01     dt=0.0219 max_dmdt=0.416
#141  t=1.03     dt=0.0219 max_dmdt=0.407
#142  t=1.05     dt=0.0219 max_dmdt=0.398
#143  t=1.07     dt=0.0219 max_dmdt=0.389
#144  t=1.09     dt=0.0219 max_dmdt=0.381
#145  t=1.11     dt=0.0219 max_dmdt=0.373
#146  t=1.14     dt=0.0219 max_dmdt=0.365
#147  t=1.16     dt=0.0219 max_dmdt=0.357
#148  t=1.18     dt=0.0219 max_dmdt=0.35
#149  t=1.2      dt=0.0219 max_dmdt=0.342
#150  t=1.22     dt=0.0219 max_dmdt=0.335
#151  t=1.25     dt=0.0219 max_dmdt=0.329
#152  t=1.27     dt=0.0219 max_dmdt=0.322
#153  t=1.29     dt=0.0219 max_dmdt=0.315
#154  t=1.31     dt=0.0219 max_dmdt=0.309
#155  t=1.33     dt=0.0219 max_dmdt=0.303
#156  t=1.36     dt=0.0219 max_dmdt=0.297
#157  t=1.38     dt=0.0219 max_dmdt=0.292
#158  t=1.4      dt=0.0219 max_dmdt=0.286
#159  t=1.42     dt=0.0219 max_dmdt=0.281
#160  t=1.44     dt=0.0219 max_dmdt=0.276
#161  t=1.46     dt=0.0219 max_dmdt=0.271
#162  t=1.49     dt=0.0219 max_dmdt=0.266
#163  t=1.51     dt=0.0219 max_dmdt=0.261
#164  t=1.54     dt=0.0331 max_dmdt=0.256
#165  t=1.57     dt=0.0331 max_dmdt=0.249
#166  t=1.61     dt=0.0331 max_dmdt=0.245
#167  t=1.64     dt=0.0331 max_dmdt=0.242
#168  t=1.67     dt=0.0331 max_dmdt=0.239
#169  t=1.71     dt=0.0331 max_dmdt=0.237
#170  t=1.74     dt=0.0331 max_dmdt=0.234
#171  t=1.77     dt=0.0331 max_dmdt=0.232
#172  t=1.81     dt=0.0331 max_dmdt=0.23
#173  t=1.84     dt=0.0331 max_dmdt=0.227
#174  t=1.87     dt=0.0331 max_dmdt=0.226
#175  t=1.91     dt=0.0331 max_dmdt=0.224
#176  t=1.94     dt=0.0331 max_dmdt=0.222
#177  t=1.97     dt=0.0331 max_dmdt=0.22
#178  t=2.01     dt=0.0331 max_dmdt=0.219
#179  t=2.04     dt=0.0331 max_dmdt=0.217
#180  t=2.07     dt=0.0331 max_dmdt=0.216
#181  t=2.1      dt=0.0331 max_dmdt=0.214
#182  t=2.14     dt=0.0331 max_dmdt=0.213
#183  t=2.17     dt=0.0331 max_dmdt=0.212
#184  t=2.2      dt=0.0331 max_dmdt=0.211
#185  t=2.25     dt=0.0502 max_dmdt=0.209
#186  t=2.3      dt=0.0502 max_dmdt=0.208
#187  t=2.35     dt=0.0502 max_dmdt=0.207
#188  t=2.4      dt=0.0502 max_dmdt=0.205
#189  t=2.45     dt=0.0502 max_dmdt=0.204
#190  t=2.5      dt=0.0502 max_dmdt=0.203
#191  t=2.56     dt=0.0502 max_dmdt=0.202
#192  t=2.61     dt=0.0502 max_dmdt=0.201
#193  t=2.66     dt=0.0502 max_dmdt=0.2
#194  t=2.71     dt=0.0502 max_dmdt=0.199
#195  t=2.76     dt=0.0502 max_dmdt=0.198
#196  t=2.81     dt=0.0502 max_dmdt=0.197
#197  t=2.86     dt=0.0502 max_dmdt=0.197
#198  t=2.91     dt=0.0502 max_dmdt=0.196
#199  t=2.96     dt=0.0502 max_dmdt=0.195
#200  t=3.01     dt=0.0502 max_dmdt=0.195
#201  t=3.06     dt=0.0502 max_dmdt=0.194
#202  t=3.11     dt=0.0502 max_dmdt=0.193
#203  t=3.16     dt=0.0502 max_dmdt=0.193
#204  t=3.21     dt=0.0502 max_dmdt=0.192
#205  t=3.26     dt=0.0502 max_dmdt=0.192
#206  t=3.31     dt=0.0502 max_dmdt=0.191
#207  t=3.36     dt=0.0502 max_dmdt=0.19
#208  t=3.41     dt=0.0502 max_dmdt=0.19
#209  t=3.46     dt=0.0502 max_dmdt=0.189
#210  t=3.51     dt=0.0502 max_dmdt=0.189
#211  t=3.56     dt=0.0502 max_dmdt=0.188
#212  t=3.63     dt=0.0763 max_dmdt=0.188
#213  t=3.71     dt=0.0763 max_dmdt=0.187
#214  t=3.79     dt=0.0763 max_dmdt=0.186
#215  t=3.86     dt=0.0763 max_dmdt=0.185
#216  t=3.94     dt=0.0763 max_dmdt=0.184
#217  t=4.02     dt=0.0763 max_dmdt=0.183
#218  t=4.09     dt=0.0763 max_dmdt=0.182
#219  t=4.17     dt=0.0763 max_dmdt=0.181
#220  t=4.24     dt=0.0763 max_dmdt=0.18
#221  t=4.32     dt=0.0763 max_dmdt=0.179
#222  t=4.4      dt=0.0763 max_dmdt=0.178
#223  t=4.47     dt=0.0763 max_dmdt=0.177
#224  t=4.55     dt=0.0763 max_dmdt=0.176
#225  t=4.63     dt=0.0763 max_dmdt=0.175
#226  t=4.7      dt=0.0763 max_dmdt=0.173
#227  t=4.78     dt=0.0763 max_dmdt=0.172
#228  t=4.86     dt=0.0763 max_dmdt=0.171
#229  t=4.93     dt=0.0763 max_dmdt=0.169
#230  t=5.01     dt=0.0763 max_dmdt=0.168
#231  t=5.08     dt=0.0763 max_dmdt=0.166
#232  t=5.16     dt=0.0763 max_dmdt=0.165
#233  t=5.24     dt=0.0763 max_dmdt=0.163
#234  t=5.31     dt=0.0763 max_dmdt=0.163
#235  t=5.39     dt=0.0763 max_dmdt=0.163
#236  t=5.47     dt=0.0763 max_dmdt=0.163
#237  t=5.54     dt=0.0763 max_dmdt=0.162
#238  t=5.66     dt=0.115 max_dmdt=0.162
#239  t=5.77     dt=0.115 max_dmdt=0.161
#240  t=5.89     dt=0.115 max_dmdt=0.16
#241  t=6        dt=0.115 max_dmdt=0.159
#242  t=6.12     dt=0.115 max_dmdt=0.158
#243  t=6.23     dt=0.115 max_dmdt=0.157
#244  t=6.35     dt=0.115 max_dmdt=0.156
#245  t=6.46     dt=0.115 max_dmdt=0.155
#246  t=6.58     dt=0.115 max_dmdt=0.153
#247  t=6.69     dt=0.115 max_dmdt=0.152
#248  t=6.81     dt=0.115 max_dmdt=0.15
#249  t=6.92     dt=0.115 max_dmdt=0.148
#250  t=7.04     dt=0.115 max_dmdt=0.146
#251  t=7.15     dt=0.115 max_dmdt=0.144
#252  t=7.27     dt=0.115 max_dmdt=0.142
#253  t=7.38     dt=0.115 max_dmdt=0.14
#254  t=7.5      dt=0.115 max_dmdt=0.138
#255  t=7.61     dt=0.115 max_dmdt=0.136
#256  t=7.73     dt=0.115 max_dmdt=0.133
#257  t=7.84     dt=0.115 max_dmdt=0.131
#258  t=7.96     dt=0.115 max_dmdt=0.129
#259  t=8.07     dt=0.115 max_dmdt=0.126
#260  t=8.19     dt=0.115 max_dmdt=0.124
#261  t=8.3      dt=0.115 max_dmdt=0.121
#262  t=8.42     dt=0.115 max_dmdt=0.12
#263  t=8.53     dt=0.115 max_dmdt=0.118
#264  t=8.65     dt=0.115 max_dmdt=0.117
#265  t=8.76     dt=0.115 max_dmdt=0.116
#266  t=8.88     dt=0.115 max_dmdt=0.114
#267  t=8.99     dt=0.115 max_dmdt=0.113
#268  t=9.11     dt=0.115 max_dmdt=0.111
#269  t=9.22     dt=0.115 max_dmdt=0.109
#270  t=9.34     dt=0.115 max_dmdt=0.108
#271  t=9.45     dt=0.115 max_dmdt=0.106
#272  t=9.57     dt=0.115 max_dmdt=0.104
#273  t=9.68     dt=0.115 max_dmdt=0.102
#274  t=9.8      dt=0.115 max_dmdt=0.101
#275  t=9.91     dt=0.115 max_dmdt=0.0988
#276  t=10       dt=0.115 max_dmdt=0.0969
#277  t=10.1     dt=0.115 max_dmdt=0.095
#278  t=10.3     dt=0.115 max_dmdt=0.0931
#279  t=10.4     dt=0.115 max_dmdt=0.0912
#280  t=10.5     dt=0.115 max_dmdt=0.0893
#281  t=10.6     dt=0.115 max_dmdt=0.0874
#282  t=10.7     dt=0.115 max_dmdt=0.0855
#283  t=10.9     dt=0.179 max_dmdt=0.083
#284  t=11.1     dt=0.179 max_dmdt=0.0801
#285  t=11.3     dt=0.179 max_dmdt=0.0771
#286  t=11.4     dt=0.179 max_dmdt=0.0742
#287  t=11.6     dt=0.179 max_dmdt=0.0714
#288  t=11.8     dt=0.179 max_dmdt=0.0685
#289  t=12       dt=0.179 max_dmdt=0.0658
#290  t=12.2     dt=0.179 max_dmdt=0.0631
#291  t=12.3     dt=0.179 max_dmdt=0.0605
#292  t=12.5     dt=0.179 max_dmdt=0.0579
#293  t=12.7     dt=0.179 max_dmdt=0.0554
#294  t=12.9     dt=0.179 max_dmdt=0.053
#295  t=13       dt=0.179 max_dmdt=0.0506
#296  t=13.2     dt=0.179 max_dmdt=0.0484
#297  t=13.4     dt=0.179 max_dmdt=0.0462
#298  t=13.6     dt=0.179 max_dmdt=0.0441
#299  t=13.8     dt=0.179 max_dmdt=0.042
#300  t=13.9     dt=0.179 max_dmdt=0.0401
#301  t=14.1     dt=0.179 max_dmdt=0.0382
#302  t=14.3     dt=0.179 max_dmdt=0.0364
#303  t=14.5     dt=0.179 max_dmdt=0.0347
#304  t=14.7     dt=0.179 max_dmdt=0.033
#305  t=14.8     dt=0.179 max_dmdt=0.0314
#306  t=15       dt=0.179 max_dmdt=0.0299
#307  t=15.2     dt=0.179 max_dmdt=0.0285
#308  t=15.4     dt=0.179 max_dmdt=0.0272
#309  t=15.6     dt=0.179 max_dmdt=0.026
#310  t=15.7     dt=0.179 max_dmdt=0.0248
#311  t=15.9     dt=0.179 max_dmdt=0.0237
#312  t=16.1     dt=0.179 max_dmdt=0.0226
#313  t=16.3     dt=0.179 max_dmdt=0.0215
#314  t=16.5     dt=0.179 max_dmdt=0.0205
#315  t=16.6     dt=0.179 max_dmdt=0.0195
#316  t=16.8     dt=0.179 max_dmdt=0.0186
#317  t=17       dt=0.179 max_dmdt=0.0177
#318  t=17.2     dt=0.179 max_dmdt=0.0168
#319  t=17.3     dt=0.179 max_dmdt=0.016
#320  t=17.5     dt=0.179 max_dmdt=0.0152
#321  t=17.7     dt=0.179 max_dmdt=0.0145
#322  t=18       dt=0.282 max_dmdt=0.0136
#323  t=18.3     dt=0.282 max_dmdt=0.0125
#324  t=18.6     dt=0.282 max_dmdt=0.0116
#325  t=18.8     dt=0.282 max_dmdt=0.0107
#326  t=19.1     dt=0.282 max_dmdt=0.00983
[13]:
plt.plot(sim.mesh.coordinates[:, 0], sim.spin.reshape(-1, 3)[:, 0], 'o-')
plt.plot(sim.mesh.coordinates[:, 0], sim.spin.reshape(-1, 3)[:, 2], 'o-')
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_24_0.png

2D skyrmion

A more complex example is a 2D skyrmion in a disk. We will use the minimiser to find the equilibrium state:

[14]:
nx, ny, nz = 100, 100, 1
dx, dy, dz = a, a, az

mesh = fidimag.common.CuboidMesh(nx=nx, ny=ny, nz=nz, dx=dx, dy=dy, dz=dz,
                                periodicity=(False, False, False),
                                unit_length=1e-9)
xs = mesh.coordinates[:, 0]
ys = mesh.coordinates[:, 1]
centre_x = (xs.max() + xs.min()) * 0.5 + xs.min()

LLG

We start finding the solution with the LLG driver:

[15]:
sim = fidimag.atomistic.Sim(mesh, name='two_dim', driver='llg')

# Define the magnetisation
def material(r):
    x, y = r[0] - centre_x, r[1] - centre_x

    if x ** 2 + y ** 2 < (np.max(xs) - centre_x) ** 2:
        return mus
    else:
        return 0

sim.set_mu_s(material)
# sim.set_mu_s(mus)

# Add the magnetic interactions
sim.add(fidimag.atomistic.Exchange(J))
sim.add(fidimag.atomistic.Anisotropy(Ku, axis=(0, 0, 1)))
sim.add(fidimag.atomistic.DMI(D, dmi_type='interfacial'))
sim.add(fidimag.atomistic.Zeeman((0, 0, B)))

def m_initial(r):
    x, y = r[0] - centre_x, r[1] - centre_x
    if x ** 2 + y ** 2 < 1 ** 2:
        return (0, 0, -1)
    else:
        return (0, 0, 1)

sim.set_m(m_initial)
# sim.set_m((0, 0, 1))

The initial configuration with a dot at the centre of the sample:

[16]:
plt.figure(figsize=(5, 5))
plt.scatter(xs, ys, c=sim.spin.reshape(-1, 3)[:, 2], s=20, marker='s')
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_32_0.png
[17]:
%%capture
sim.driver.do_precession = False
sim.driver.relax()

After relaxation we see a skyrmion and the spin tilting at the boundary, characteristic of the DMI

[18]:
plt.figure(figsize=(6, 6))
plt.scatter(xs, ys, c=sim.spin.reshape(-1, 3)[:, 2], vmin=-1, vmax=1)

plt.quiver(xs, ys, sim.spin.reshape(-1, 3)[:, 0], sim.spin.reshape(-1, 3)[:, 1])
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_35_0.png

Steepest Descent

We can get the same result using the minimiser:

[19]:
sim = fidimag.atomistic.Sim(mesh, name='two_dim_SD', driver='steepest_descent')

# Define the magnetisation
def material(r):
    x, y = r[0] - centre_x, r[1] - centre_x

    if x ** 2 + y ** 2 < (np.max(xs) - centre_x) ** 2:
        return mus
    else:
        return 0

sim.set_mu_s(material)
# sim.set_mu_s(mus)

# Add the magnetic interactions
sim.add(fidimag.atomistic.Exchange(J))
sim.add(fidimag.atomistic.Anisotropy(Ku, axis=(0, 0, 1)))
sim.add(fidimag.atomistic.DMI(D, dmi_type='interfacial'))
sim.add(fidimag.atomistic.Zeeman((0, 0, B)))

def m_initial(r):
    x, y = r[0] - centre_x, r[1] - centre_x
    if x ** 2 + y ** 2 < 1:
        return (0, 0, -1)
    else:
        return (0, 0, 1)

sim.set_m(m_initial)
# sim.set_m((0, 0, 1))
[20]:
plt.figure(figsize=(2, 2))
plt.scatter(xs, ys, c=sim.spin.reshape(-1, 3)[:, 2])
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_39_0.png
[21]:
sim.driver.tmax = 1e-1
sim.driver.minimise(max_steps=6000, stopping_dm=1e-5, initial_t_step=1e-4)
#0    max_tau=0.0001   max_dm=0.00254
#1000 max_tau=0.0001   max_dm=0.000413
#2000 max_tau=0.0001   max_dm=0.000109
#3000 max_tau=0.0001   max_dm=6.32e-05
#4000 max_tau=0.0001   max_dm=4.88e-05
#5000 max_tau=0.0001   max_dm=4.15e-05
Warning: minimise did not converge in 6000 steps - maxdm = 3.653240142795239e-05
[22]:
sim.driver.spin.shape
[22]:
(30000,)
[23]:
plt.figure(figsize=(6, 6))
plt.scatter(xs, ys, c=sim.spin.reshape(-1, 3)[:, 2], vmin=-1, vmax=1)

plt.quiver(xs, ys, sim.spin.reshape(-1, 3)[:, 0], sim.spin.reshape(-1, 3)[:, 1])
plt.show()
../../_images/user_guide_ipynb_steepest_descent_atomistic_42_0.png
[ ]: