-
Notifications
You must be signed in to change notification settings - Fork 5
/
examples.html
1103 lines (886 loc) · 31.7 KB
/
examples.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
<!DOCTYPE html>
<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
<title>mumax3</title>
<link rel="stylesheet" href="https://unpkg.com/[email protected]/build/pure-min.css">
<link rel="stylesheet" type="text/css" href="style.css">
<script>
(function(i,s,o,g,r,a,m){i['GoogleAnalyticsObject']=r;i[r]=i[r]||function(){
(i[r].q=i[r].q||[]).push(arguments)},i[r].l=1*new Date();a=s.createElement(o),
m=s.getElementsByTagName(o)[0];a.async=1;a.src=g;m.parentNode.insertBefore(a,m)
})(window,document,'script','//www.google-analytics.com/analytics.js','ga');
ga('create', 'UA-50169601-1', 'mumax.github.io');
ga('send', 'pageview');
</script>
</head>
<body>
<div style="float:left">
<img src="nimble-cubes128.png"/>
</div>
<div style="float:left">
<span style="font-size:48px"><b>mumax</b><sup>3</sup><br/></span>
GPU-accelerated micromagnetism <br/><br/>
<a class="pure-button pure-button-primary" href="index.html">Home</a>
<a class="pure-button pure-button-primary" href="download.html">Download</a>
<a class="pure-button pure-button-primary" href="examples.html">Examples</a>
<a class="pure-button pure-button-primary" href="api.html">API</a>
<a class="pure-button pure-button-primary" href="https://groups.google.com/forum/#!forum/mumax2">Forum</a>
</div>
<br style="clear:both"/>
<hr/>
<h1>mumax 3.10 examples</h1>
<p>These are example input scripts, the full API can be found <a href="http://mumax.github.io/api.html">here</a>.</p>
mumax<sup>3</sup> input files are run with the command
<pre>mumax3 myfile.mx3</pre>
Output is automatically stored in the "myfile.out" directory. Additionally, a web interface provides live output. Default is <code><a href="http://127.0.0.1:35367">http://localhost:35367</a></code>.<br/>
For more details, run <code>mumax3 -help</code> which will show the available command-line flags (e.g. to select a certain GPU).
<hr/><h2>Getting started with Standard Problem #4</h2>
Let's start with the classic mumag standard problem 4, as defined <a href="http://www.ctcms.nist.gov/~rdm/mumag.org.html">here</a>.
<a id=example1></a><pre>SetGridsize(128, 32, 1)
SetCellsize(500e-9/128, 125e-9/32, 3e-9)
Msat = 800e3
Aex = 13e-12
alpha = 0.02
m = uniform(1, .1, 0)
relax()
save(m) // relaxed state
autosave(m, 200e-12)
tableautosave(10e-12)
B_ext = vector(-24.6E-3, 4.3E-3, 0)
run(1e-9)</pre>
<p>This example should be pretty straight-forward to follow. Space-dependent output is stored in OVF format, which is compatible with OOMMF and can be converted with <a href="http://godoc.org/github.com/mumax/3/cmd/mumax3-convert"><code>mumax3-convert</code></a>. Below is the output converted to PNG.</p>
<p>The data table is stored in a simple text format compatible with <a href=http://www.gnuplot.info>gnuplot</a>, like used for the plot below.</p>
<h3>output</h3>
<figure style="float:left">
<img src="example1.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example1.out/m000001.png"/>
<figcaption> m000001 </figcaption>
</figure>
<figure style="float:left">
<img src="example1.out/m000002.png"/>
<figcaption> m000002 </figcaption>
</figure>
<figure style="float:left">
<img src="example1.out/m000003.png"/>
<figcaption> m000003 </figcaption>
</figure>
<figure style="float:left">
<img src="example1.out/m000004.png"/>
<figcaption> m000004 </figcaption>
</figure>
<figure style="float:left">
<img src="example1.out/m000005.png"/>
<figcaption> m000005 </figcaption>
</figure>
<figure style="float:left">
<img src="example1.out/m000006.png"/>
<figcaption> m000006 </figcaption>
</figure><br style="clear:both"/>
<figure>
<img src="example1.out/m.svg"/>
<figcaption> m.svg </figcaption>
</figure>
<hr/><h2>Standard Problem #2</h2>
Using the scripting language explained above, relatively complex input files can be easily defined. E.g. <a href="http://www.ctcms.nist.gov/~rdm/std2/spec2.html">micromagnetic standard problem #2</a> specifies the simulation size in exchange lengths. The script below calculates the exchange length and chooses cells not larger than 0.75 exchange lengths so that the number of cells is a power of two (for best performance).
<a id=example2></a><pre>Msat = 1000e3
Aex = 10e-12
// define exchange length
lex := sqrt(10e-12 / (0.5 * mu0 * pow(1000e3 ,2)))
d := 30 * lex // we test for d/lex = 30
Sizex := 5*d // magnet size x
Sizey := 1*d
Sizez := 0.1*d
nx := pow(2, ilogb(Sizex / (0.75*lex))) // power-of-two number of cells
ny := pow(2, ilogb(Sizey / (0.75*lex))) // not larger than 0.75 exchange lengths
SetGridSize(nx, ny, 1)
SetCellSize(Sizex/nx, Sizey/ny, Sizez)
m = Uniform(1, 0.1, 0) // initial mag
relax()
save(m) // remanent magnetization
print("<m> for d/lex=30: ", m.average())</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example2.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure><br style="clear:both"/>
This example saves and prints the remanent magnetization state so we can verify it against known values.
<hr/><h2>Hysteresis</h2>
Below is an example of a hysteresis loop where we step the applied field in small increments and find the magnetization ground state after each step. Minimize() finds the ground state using the conjugate gradient method, which is very fast. However, this method might fail on very high energy initial states like a random magnetization. In that case, Relax() is more robust (albeit much slower).
<a id=example3></a><pre>SetGridsize(128, 32, 1)
SetCellsize(4e-9, 4e-9, 30e-9)
Msat = 800e3
Aex = 13e-12
m = randomMag()
relax() // high-energy states best minimized by relax()
Bmax := 100.0e-3
Bstep := 1.0e-3
MinimizerStop = 1e-6
TableAdd(B_ext)
for B:=0.0; B<=Bmax; B+=Bstep{
B_ext = vector(B, 0, 0)
minimize() // small changes best minimized by minimize()
tablesave()
}
for B:=Bmax; B>=-Bmax; B-=Bstep{
B_ext = vector(B, 0, 0)
minimize() // small changes best minimized by minimize()
tablesave()
}
for B:=-Bmax; B<=Bmax; B+=Bstep{
B_ext = vector(B, 0, 0)
minimize() // small changes best minimized by minimize()
tablesave()
}</pre>
<h3>output</h3>
<figure>
<img src="example3.out/hysteresis.svg"/>
</figure>
<hr/><h2>Geometry</h2>
mumax3 has powerful API to programatically define geometries. A number of primitive shapes are defined, like ellipses, rectangles, etc. They can be transformed (rotated, translated) and combined using boolean logic (add, sub, inverse). All positions are specified in meters and the origin lies in the center of the simulation box. See the full <a href="http://mumax.github.io/api.html">API</a>.
Edges can be smoothed to reduce staircase effects. <code>EdgeSmooth=n</code> means <code>n³</code> samples per cell are used to determine its volume. <code>EdgeSmooth=0</code> implies a staircase approximation, while <code>EdgeSmooth=8</code> results in quite accurately resolved edges.
<a id=example4></a><pre>SetGridsize(100, 100, 50)
SetCellsize(1e-6/100, 1e-6/100, 1e-6/50)
EdgeSmooth = 8
setgeom( rect(800e-9, 500e-9) )
saveas(geom, "rect")
setgeom( cylinder(800e-9, inf) )
saveas(geom, "cylinder")
setgeom( circle(200e-9).repeat(300e-9, 400e-9, 0) )
saveas(geom, "circle_repeat")
setgeom( cylinder(800e-9, inf).inverse() )
saveas(geom, "cylinder_inverse")
setgeom( cylinder(800e-9, 600e-9).transl(200e-9, 100e-9, 0) )
saveas(geom, "cylinder_transl")
setgeom( ellipsoid(800e-9, 600e-9, 500e-9) )
saveas(geom, "ellipsoid")
setgeom( cuboid(800e-9, 600e-9, 500e-9) )
saveas(geom, "cuboid")
setgeom( cuboid(800e-9, 600e-9, 500e-9).rotz(-10*pi/180) )
saveas(geom, "cuboid_rotZ")
setgeom( layers(0, 25) )
saveas(geom, "layers")
setgeom( cell(50, 20, 0) )
saveas(geom, "cell")
setgeom( xrange(0, inf) )
saveas(geom, "xrange")
a := cylinder(600e-9, 600e-9).transl(-150e-9, 50e-9, 0 )
b := rect(600e-9, 600e-9).transl(150e-9, -50e-9, 0)
setgeom( a.add(b) )
saveas(geom, "logicAdd")
setgeom( a.sub(b) )
saveas(geom, "logicSub")
setgeom( a.intersect(b) )
saveas(geom, "logicAnd")
setgeom( a.xor(b) )
saveas(geom, "logicXor")
setgeom( imageShape("mask.png") )
saveas(geom, "imageShape")</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example4.out/cell.png"/>
<figcaption> cell </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/circle_repeat.png"/>
<figcaption> circle_repeat </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/cuboid.png"/>
<figcaption> cuboid </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/cuboid_rotZ.png"/>
<figcaption> cuboid_rotZ </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/cylinder.png"/>
<figcaption> cylinder </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/cylinder_inverse.png"/>
<figcaption> cylinder_inverse </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/cylinder_transl.png"/>
<figcaption> cylinder_transl </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/ellipsoid.png"/>
<figcaption> ellipsoid </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/imageShape.png"/>
<figcaption> imageShape </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/layers.png"/>
<figcaption> layers </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/logicAdd.png"/>
<figcaption> logicAdd </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/logicAnd.png"/>
<figcaption> logicAnd </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/logicSub.png"/>
<figcaption> logicSub </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/logicXor.png"/>
<figcaption> logicXor </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/rect.png"/>
<figcaption> rect </figcaption>
</figure>
<figure style="float:left">
<img src="example4.out/xrange.png"/>
<figcaption> xrange </figcaption>
</figure><br style="clear:both"/>
Note: these are 3D geometries seen from above. The displayed cell filling is averaged along the thickness (notable in ellipse and layers example). Black means empty space, white is filled.
<hr/><h2>Initial Magnetization</h2>
Some initial magnetization functions are provided, as well as transformations similar to those on Shapes. See the Config <a href="http://mumax.github.io/api.html">API</a>.
<a id=example5></a><pre>setgridsize(256, 128, 1)
setcellsize(5e-9, 5e-9, 5e-9)
m = Uniform(1, 1, 0) // no need to normalize length
saveas(m, "uniform")
m = Vortex(1, -1) // circulation, polarization
saveas(m, "vortex")
m = TwoDomain(1,0,0, 0,1,0, -1,0,0) // Néel wall
saveas(m, "twodomain")
m = RandomMag()
saveas(m, "randommag")
m = TwoDomain(1,0,0, 0,1,0, -1,0,0).rotz(-pi/4)
saveas(m, "twodomain_rot")
m = VortexWall(1, -1, 1, 1)
saveas(m, "vortexwall")
m = VortexWall(1, -1, 1, 1).scale(1/2, 1, 1)
saveas(m, "vortexwall_scale")
m = Vortex(1,-1).transl(100e-9, 50e-9, 0)
saveas(m, "vortex_transl")
m = Vortex(1,-1).Add(0.1, randomMag())
saveas(m, "vortex_add_random")
m = BlochSkyrmion(1, -1).scale(3,3,1)
saveas(m, "Bloch_skyrmion")
m = NeelSkyrmion(1,-1).scale(3,3,1)
saveas(m, "Néel_skyrmion")
// set m in only a part of space, or a single cell:
m = uniform(1, 1, 1)
m.setInShape(cylinder(400e-9, 100e-9), vortex(1, -1))
m.setCell(20, 10, 0, vector(0.1, 0.1, -0.9)) // set in cell index [20,10,0]
saveas(m, "setInShape_setCell")
//Read m form .ovf file.
m.loadfile("myfile.ovf")
saveas(m, "loadfile")</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example5.out/Bloch_skyrmion.png"/>
<figcaption> Bloch_skyrmion </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/Néel_skyrmion.png"/>
<figcaption> Néel_skyrmion </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/loadfile.png"/>
<figcaption> loadfile </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/randommag.png"/>
<figcaption> randommag </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/setInShape_setCell.png"/>
<figcaption> setInShape_setCell </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/twodomain.png"/>
<figcaption> twodomain </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/twodomain_rot.png"/>
<figcaption> twodomain_rot </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/uniform.png"/>
<figcaption> uniform </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/vortex.png"/>
<figcaption> vortex </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/vortex_add_random.png"/>
<figcaption> vortex_add_random </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/vortex_transl.png"/>
<figcaption> vortex_transl </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/vortexwall.png"/>
<figcaption> vortexwall </figcaption>
</figure>
<figure style="float:left">
<img src="example5.out/vortexwall_scale.png"/>
<figcaption> vortexwall_scale </figcaption>
</figure><br style="clear:both"/>
These initial states are approximate, after setting them it is a good idea to relax the magnetization to the actual ground state.
The magnetization can also be set in separate regions, see below.
<hr/><h2>Interlude: Rotating Cheese</h2>
In this example we define a geometry that looks like a slice of cheese and have it rotate in time.
<a id=example6></a><pre>setgridsize(128, 128, 1)
setcellsize(2e-9, 2e-9, 2e-9)
d := 200e-9
sq := rect(d, d) // square with side d
h := 50e-9
hole := cylinder(h, h) // circle with diameter h
hole1 := hole.transl(100e-9, 0, 0) // translated circle #1
hole2 := hole.transl(0, -50e-9, 0) // translated cricle #2
cheese:= sq.sub(hole1).sub(hole2)// subtract the circles from the square (makes holes).
setgeom(cheese)
msat = 600e3
aex = 12e-13
alpha = 3
// rotate the cheese.
for i:=0; i<=90; i=i+30{
angle := i*pi/180
setgeom(cheese.rotz(angle))
m = uniform(cos(angle), sin(angle), 0)
minimize()
save(m)
}</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example6.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example6.out/m000001.png"/>
<figcaption> m000001 </figcaption>
</figure>
<figure style="float:left">
<img src="example6.out/m000002.png"/>
<figcaption> m000002 </figcaption>
</figure>
<figure style="float:left">
<img src="example6.out/m000003.png"/>
<figcaption> m000003 </figcaption>
</figure><br style="clear:both"/>
<hr/><h2>Regions: Space-dependent Parameters</h2>
<p>Space-dependent parameters are defined using material <i>regions</i>. Regions are numbered 0-255 and represent different materials. Each cell can belong to only one region. At the start of a simulation all cells have region number 0.</p>
<p>Regions are defined with <code>defregion(number, shape)</code>, where <code>shape</code> is explained in the geometry example.</p>
<p>When you're not using regions, like in the above examples, you'll probably set parameters with a simple assign:
<pre>Aex = 12e-13</pre>
Behind the screens, this sets Aex in <i>all</i> regions.
</p>
<p>It's always a good idea to output the <code>regions</code> quantity, as well as all your material parameters. </p>
<a id=example7></a><pre>N := 128
setgridsize(N, N, 1)
c := 4e-9
setcellsize(c, c, c)
// disk with different anisotropy in left and right half
setgeom(circle(N*c))
defregion(1, xrange(0, inf)) // left half
defregion(2, xrange(-inf, 0)) // right half
save(regions)
Ku1.setregion(1, .1e6)
anisU.setRegion(1, vector(1, 0, 0))
Ku1.setregion(2, .2e6)
anisU.setRegion(2, vector(0, 1, 0))
save(Ku1)
save(anisU)
Msat = 800e3 // sets it everywhere
save(Msat)
Aex = 12e-13
alpha = 1
m.setRegion(1, uniform(1, 1, 0))
m.setRegion(2, uniform(-1, 1, 0))
saveas(m, "m_inital")
run(.1e-9)
saveas(m, "m_final")</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example7.out/Ku1000000.png"/>
<figcaption> Ku1000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example7.out/Msat000000.png"/>
<figcaption> Msat000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example7.out/anisU000000.png"/>
<figcaption> anisU000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example7.out/m_final.png"/>
<figcaption> m_final </figcaption>
</figure>
<figure style="float:left">
<img src="example7.out/m_inital.png"/>
<figcaption> m_inital </figcaption>
</figure>
<figure style="float:left">
<img src="example7.out/regions000000.png"/>
<figcaption> regions000000 </figcaption>
</figure><br style="clear:both"/>
<hr/><h2>Slicing and dicing output</h2>
The example below illustrates how to save only the part of the output you're interested in.
<a id=example8></a><pre>Nx := 256
Ny := 256
Nz := 1
setgridsize(Ny, Nx, Nz)
c := 4e-9
setcellsize(c, c, c)
setgeom(circle(Nx*c))
Msat = 800e3
Aex = 12e-13
alpha = 1
m = vortex(1, 1)
save(m)
save(m.Comp(0))
save(Crop(m, 0, Nx/2, 0, Ny/2, 0, Nz))
mx := m.Comp(0)
mx_center := CropY(mx, Ny/4, 3*Ny/4)
save(mx_center)</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example8.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example8.out/m_x000000.png"/>
<figcaption> m_x000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example8.out/m_x_yrange64-192_000000.png"/>
<figcaption> m_x_yrange64-192_000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example8.out/m_xrange0-128_yrange0-128_000000.png"/>
<figcaption> m_xrange0-128_yrange0-128_000000 </figcaption>
</figure><br style="clear:both"/>
<a id=MFM><hr/><h2>Magnetic Force Microscopy</h2></a>
<p>Mumax3 has built-in generation of MFM images from the magnetization. The MFM tip lift can be freely chosen. By default the tip magnetization is modeled as a point monopole at the apex. This is sufficient for most situations. Nevertheless, it is also possible to model partially magnetized tips by setting MFMDipole to the magnetized portion of the tip, in meters. E.g., if only the first 20nm of the tip is (vertically) magnetized, set MFMDipole=20e-9.</p>
<a id=example9></a><pre>setgridsize(256, 256, 1)
setcellsize(2e-9, 2e-9, 1e-9)
setgeom(rect(400e-9, 400e-9))
msat = 600e3
aex = 10e-12
m = vortex(1, 1)
relax()
save(m)
MFMLift = 10e-9
saveas(MFM, "lift_10nm")
MFMLift = 40e-9
saveas(MFM, "lift_40nm")
MFMLift = 90e-9
saveas(MFM, "lift_90nm")</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example9.out/lift_10nm.png"/>
<figcaption> lift_10nm </figcaption>
</figure>
<figure style="float:left">
<img src="example9.out/lift_40nm.png"/>
<figcaption> lift_40nm </figcaption>
</figure>
<figure style="float:left">
<img src="example9.out/lift_90nm.png"/>
<figcaption> lift_90nm </figcaption>
</figure>
<figure style="float:left">
<img src="example9.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure><br style="clear:both"/>
<hr/><h2>PMA Racetrack</h2>
In this example we drive a domain wall in PMA material by spin-transfer torque. We set up a post-step function that makes the simulation box "follow" the domain wall. Like this, only a small number of cells is needed to simulate an infinitely long magnetic wire.
<a id=example10></a><pre>setGridSize(128, 128, 1)
setCellSize(2e-9, 2e-9, 1e-9)
Msat = 600e3
Aex = 10e-12
anisU = vector(0, 0, 1)
Ku1 = 0.59e6
alpha = 0.02
Xi = 0.2
m = twoDomain(0, 0, 1, 1, 1, 0, 0, 0, -1) // up-down domains with wall between Bloch and Néél type
relax()
// Set post-step function that centers simulation window on domain wall.
ext_centerWall(2) // keep m[2] (= m_z) close to zero
// Schedule output
autosave(m, 100e-12)
// Run for 1ns with current through the sample
j = vector(1.5e13, 0, 0)
pol = 1
run(.5e-9)</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example10.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example10.out/m000001.png"/>
<figcaption> m000001 </figcaption>
</figure>
<figure style="float:left">
<img src="example10.out/m000002.png"/>
<figcaption> m000002 </figcaption>
</figure>
<figure style="float:left">
<img src="example10.out/m000003.png"/>
<figcaption> m000003 </figcaption>
</figure>
<figure style="float:left">
<img src="example10.out/m000004.png"/>
<figcaption> m000004 </figcaption>
</figure>
<figure style="float:left">
<img src="example10.out/m000005.png"/>
<figcaption> m000005 </figcaption>
</figure><br style="clear:both"/>
Since we center on the domain wall we can not see that it is actually moving, but the domain wall breakdown is visible.
<hr/><h2>Py Racetrack</h2>
In this example we drive a vortex wall in Permalloy by spin-transfer torque. The simulation box "follows" the domain wall. By removing surface charges at the left and right ends, we mimic an infintely long wire.
<a id=example11></a><pre>setGridSize(256, 64, 1)
setCellSize(3e-9, 3e-9, 10e-9)
Msat = 860e3
Aex = 13e-12
Xi = 0.1
alpha = 0.02
m = twodomain(1,0,0, 0,1,0, -1,0,0)
notches := rect(15e-9, 15e-9).RotZ(45*pi/180).Repeat(200e-9, 64*3e-9, 0).Transl(0, 32*3e-9, 0)
setGeom(notches.inverse())
// Remove surface charges from left (mx=1) and right (mx=-1) sides to mimic infinitely long wire. We have to specify the region (0) at the boundaries.
BoundaryRegion := 0
MagLeft := 1
MagRight := -1
ext_rmSurfaceCharge(BoundaryRegion, MagLeft, MagRight)
relax()
ext_centerWall(0) // keep m[0] (m_x) close to zero
// Schedule output
autosave(m, 50e-12)
tableadd(ext_dwpos) // domain wall position
tableautosave(10e-12)
// Run the simulation with current through the sample
pol = 0.56
J = vector(-10e12, 0, 0)
Run(0.5e-9)</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example11.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000001.png"/>
<figcaption> m000001 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000002.png"/>
<figcaption> m000002 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000003.png"/>
<figcaption> m000003 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000004.png"/>
<figcaption> m000004 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000005.png"/>
<figcaption> m000005 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000006.png"/>
<figcaption> m000006 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000007.png"/>
<figcaption> m000007 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000008.png"/>
<figcaption> m000008 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000009.png"/>
<figcaption> m000009 </figcaption>
</figure>
<figure style="float:left">
<img src="example11.out/m000010.png"/>
<figcaption> m000010 </figcaption>
</figure><br style="clear:both"/>
<figure>
<img src="example11.out/ext_dwpos.svg"/>
<figcaption> ext_dwpos.svg </figcaption>
</figure>
<figure>
<img src="example11.out/m.svg"/>
<figcaption> m.svg </figcaption>
</figure>
Since we center on the domain wall we can not really see the motion, despite the vortex wall moving pretty fast. Note the absence of closure domains at the edges due to the surface charges being removed there.
<hr/><h2>Voronoi tessellation</h2>
In this example we use regions to specify grains in a material. The built-in extension <code>ext_makegrains</code> is used to define grain-like regions using Voronoi tessellation. We vary the material parameters in each grain.
<a id=example12></a><pre>N := 256
c := 4e-9
d := 40e-9
setgridsize(N, N, 1)
setcellsize(c, c, d)
setGeom(circle(N*c))
// define grains with region number 0-255
grainSize := 40e-9 // m
randomSeed := 1234567
maxRegion := 255
ext_makegrains(grainSize, maxRegion, randomSeed)
defregion(256, circle(N*c).inverse()) // region 256 is outside, not really needed
alpha = 3
Kc1 = 1000
Aex = 13e-12
Msat = 860e3
// set random parameters per region
for i:=0; i<maxRegion; i++{
// random cubic anisotropy direction
axis1 := vector(randNorm(), randNorm(), randNorm())
helper := vector(randNorm(), randNorm(), randNorm())
axis2 := axis1.cross(helper) // perpendicular to axis1
AnisC1.SetRegion(i, axis1) // axes need not be normalized
AnisC2.SetRegion(i, axis2)
// random 10% anisotropy variation
K := 1e5
Kc1.SetRegion(i, K + randNorm() * 0.1 * K)
}
// reduce exchange coupling between grains by 10%
for i:=0; i<maxRegion; i++{
for j:=i+1; j<maxRegion; j++{
ext_ScaleExchange(i, j, 0.9)
}
}
m = vortex(1, 1)
run(.1e-9)
save(regions)
save(Kc1)
save(AnisC1)
save(AnisC2)
save(m)
save(exchCoupling)</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example12.out/ExchCoupling000000.png"/>
<figcaption> ExchCoupling000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example12.out/Kc1000000.png"/>
<figcaption> Kc1000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example12.out/anisC1000000.png"/>
<figcaption> anisC1000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example12.out/anisC2000000.png"/>
<figcaption> anisC2000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example12.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example12.out/regions000000.png"/>
<figcaption> regions000000 </figcaption>
</figure><br style="clear:both"/>
<hr/><h2>RKKY</h2>
Scaling the exchange coupling between regions can be used to obtain antiferromagnetic coupling like the RKKY interaction. In that case we only model the magnetic layers and do not explicitly add a spacer layer (which is negligibly thin). We scale the exchange coupling to get the desired RKKY strength: <code>scale = (RKKY * cellsize_z) / (2 * Aex)</code>.
<a id=example13></a><pre>N := 10
setgridsize(N, N, 2)
c := 1e-9
setcellsize(c, c, c)
defRegion(0, layer(0))
defRegion(1, layer(1))
Msat = 1e6
Aex = 10e-12
RKKY := -1e-3 // 1mJ/m2
scale := (RKKY * c) / (2 * Aex.Average())
ext_scaleExchange(0, 1, scale)
tableAdd(E_total)
m.setRegion(0, uniform(1, 0, 0))
for ang:=0; ang<360; ang++{
m.setRegion(1, uniform(cos(ang*pi/180), sin(ang*pi/180), 0))
t = ang * 1e-9 // output "time" is really angle
tablesave()
}</pre>
<h3>output</h3> <br style="clear:both"/>
<figure>
<img src="example13.out/E_total.svg"/>
<figcaption> E_total.svg </figcaption>
</figure>
<figure>
<img src="example13.out/m.svg"/>
<figcaption> m.svg </figcaption>
</figure>
<hr/><h2>Slonczewski STT</h2>
Example of a spin-torque MRAM stack consisting of a fixed layer, spacer and free layer. Only the free layer magnetization is explicitly modeled, so we use a 2D grid. The fixed layer polarization is set with <code>FixedLayer = ...</code>, which can be space-dependent. The spacer layer properties are modeled by setting the parameters <code>Lambda</code> and <code>EpsilonPrime</code>. Finally <code>Pol</code> sets the current polarization and <code>J</code> the current density, which should be along z in this case.
Below we switch an MRAM bit.
<a id=example14></a><pre>// geometry
sizeX := 160e-9
sizeY := 80e-9
sizeZ := 5e-9
Nx := 64
Ny := 32
setgridsize(Nx, Ny, 1)
setcellsize(sizeX/Nx, sizeY/Ny, sizeZ)
setGeom(ellipse(sizeX, sizeY))
// set up free layer
Msat = 800e3
Aex = 13e-12
alpha = 0.01
m = uniform(1, 0, 0)
// set up spacer layer parameters
lambda = 1
Pol = 0.5669
epsilonprime = 0
// set up fixed layer polarization
angle := 20
px := cos(angle * pi/180)
py := sin(angle * pi/180)
fixedlayer = vector(px, py, 0)
// send current
Jtot := -0.008 // total current in A
area := sizeX*sizeY*pi/4
jc := Jtot / area // current density in A/m2
J = vector(0, 0, jc)
// schedule output & run
autosave(m, 100e-12)
tableautosave(10e-12)
run(1e-9)</pre>
<h3>output</h3>
<figure style="float:left">
<img src="example14.out/m000000.png"/>
<figcaption> m000000 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000001.png"/>
<figcaption> m000001 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000002.png"/>
<figcaption> m000002 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000003.png"/>
<figcaption> m000003 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000004.png"/>
<figcaption> m000004 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000005.png"/>
<figcaption> m000005 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000006.png"/>
<figcaption> m000006 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000007.png"/>
<figcaption> m000007 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000008.png"/>
<figcaption> m000008 </figcaption>
</figure>
<figure style="float:left">
<img src="example14.out/m000009.png"/>