1. Spin wave simulations of CoRh₂O₄
This tutorial introduces Sunny through its features for performing conventional spin wave theory calculations. We consider the crystal CoRh₂O₄ and reproduce the calculations of Ge et al., Phys. Rev. B 96, 064413 (2017) .
Get Julia and Sunny
Sunny is implemented in Julia, which allows for interactive development (like Python or Matlab) while also providing high numerical efficiency (like C++ or Fortran). New Julia users should begin with our Getting Started guide. Sunny requires Julia 1.10 or later.
From the Julia prompt, load Sunny and also GLMakie for graphics.
If these packages are not yet installed, Julia will offer to install them. If executing this tutorial gives an error, you may need to update Sunny and GLMakie from the built-in package manager.
Units system
The
Units
object selects reference energy and length scales, and uses these to provide physical constants. For example,
units.K
returns one kelvin as 0.086 meV, where the Boltzmann constant is implicit.
Crystal cell
A crystallographic cell may be loaded from a
.cif
file, or specified from atom positions and types.
Start by defining the shape of the conventional chemical cell. CoRh₂O₄ has cubic spacegroup 227 (Fd-3m). Its lattice constants are 8.5 Å and the cell angles are 90°. With this information,
lattice_vectors
constructs a 3×3 matrix
latvecs
. Columns of
latvecs
define the lattice vectors
in the global Cartesian coordinate system. Conversely, columns of
inv(latvecs)
define the global Cartesian axes
in components of the lattice vectors.
3×3 StaticArraysCore.SMatrix{3, 3, Float64, 9} with indices SOneTo(3)×SOneTo(3): 8.5031 0.0 0.0 0.0 8.5031 0.0 0.0 0.0 8.5031
Construct a
Crystal
cell from spacegroup 227 in the ITA standard setting. Cobalt atoms belong to Wyckoff 8a, which is the diamond cubic lattice.
[4m[1mCrystal[22m[24m Spacegroup 'F d -3 m' (227) Lattice params a=8.503, b=8.503, c=8.503, α=90°, β=90°, γ=90° Cell volume 614.8 Type 'Co', Wyckoff 8a (site sym. '-43m'): 1. [1/8, 1/8, 1/8] 2. [5/8, 5/8, 1/8] 3. [7/8, 3/8, 3/8] 4. [3/8, 7/8, 3/8] 5. [5/8, 1/8, 5/8] 6. [1/8, 5/8, 5/8] 7. [3/8, 3/8, 7/8] 8. [7/8, 7/8, 7/8]
view_crystal
launches an interface for interactive inspection and symmetry analysis.
Spin system
A
System
will define the spin model. Each cobalt atom carries quantum spin
, with a
-factor of 2. Specify this
Moment
data for cobalt atom 1. By symmetry, the same moment data also applies to cobalt atoms 2, 3, ... 7. The option
:dipole
indicates a traditional model type, for which quantum spin is modeled as a dipole expectation value.
[4m[1mSystem [Dipole mode][22m[24m Supercell (1×1×1)×8 Energy per site -567/200
Ge et al. demonstrated that inelastic neutron scattering data for CoRh₂O₄ is well modeled by antiferromagnetic nearest neighbor exchange,
J = 0.63
meV. Call
set_exchange!
with the bond that connects atom 1 to atom 3, and has zero displacement between chemical cells. Consistent with the symmetries of spacegroup 227, this interaction will be propagated to all other nearest-neighbor bonds. Calling
view_crystal
with
sys
now shows the antiferromagnetic Heisenberg interactions as blue polkadot spheres.
Optimizing spins
To search for the ground state, call
randomize_spins!
and
minimize_energy!
in sequence. For this problem, optimization converges rapidly to the expected Néel order. See this with
plot_spins
, where spins are colored according to their global
-component.
The diamond lattice is bipartite, allowing each spin to perfectly anti-align with its 4 nearest-neighbors. Each of these 4 bonds contribute
to the total energy. Two sites participate in each bond, so the energy per site is
. Check this by calling
energy_per_site
.
Reshaping the magnetic cell
The most compact magnetic cell for this Néel order is the primitive unit cell. Columns of the
primitive_cell
matrix provide the primitive lattice vectors as multiples of the conventional cubic lattice vectors.
3×3 StaticArraysCore.SMatrix{3, 3, Float64, 9} with indices SOneTo(3)×SOneTo(3): 0.0 0.5 0.5 0.5 0.0 0.5 0.5 0.5 0.0
Reduce the magnetic cell size using
reshape_supercell
. Verify that the energy per site is unchanged after the reshaping the supercell.
Plotting
sys_prim
shows the two spins within the primitive cell, as well as the larger conventional cubic cell for context.
Spin wave theory
With this primitive cell, we will perform a
SpinWaveTheory
calculation of the structure factor
. The measurement
ssf_perp
indicates projection of the spin structure factor
perpendicular to the direction of momentum transfer, as appropriate for unpolarized neutron scattering. The isotropic
FormFactor
for Co²⁺ dampens intensities at large
.
[4m[1mSpinWaveTheory [Dipole mode][22m[24m 2 atoms
Select
lorentzian
broadening with a full-width at half-maximum (FWHM) of 0.8 meV.
Lorentzian kernel, fwhm=0.8
Define a
q_space_path
that connects high-symmetry points in reciprocal space. The
-points are given in reciprocal lattice units (RLU) for the
original cubic cell. For example,
[1/2, 1/2, 0]
denotes the sum of the first two reciprocal lattice vectors,
. A total of 500
-points will be sampled along the path.
[4m[1mQPath (500 samples)[22m[24m [0, 0, 0] → [1/2, 0, 0] → [1/2, 1/2, 0] → [0, 0, 0]
Calculate single-crystal scattering
intensities
along this path, for energies between 0 and 6 meV. Use
plot_intensities
to visualize the result.
Sometimes experimental data is only available as a powder average, i.e., as an average over all possible crystal orientations. Use
powder_average
to simulate these intensities. Each
-magnitude defines a spherical shell in reciprocal space. Consider 200 radii from 0 to 3 inverse angstroms and collect 2000 random samples per spherical shell. As configured, this calculation completes in about two seconds. Had we used the conventional cubic cell, the calculation would be an order of magnitude slower.
This result can be compared to experimental neutron scattering data from Fig. 5 of Ge et al.
<img width="95%" src="https://raw.githubusercontent.com/SunnySuite/Sunny.jl/main/docs/src/assets/CoRh2O4_intensity.jpg">
What's next?
-
For more spin wave calculations of this type, browse the SpinW tutorials ported to Sunny.
-
Spin wave theory neglects thermal fluctuations of the magnetic order. The next CoRh₂O₄ tutorial demonstrates how to sample spins in thermal equilibrium and measure correlations from the classical spin dynamics.
-
Sunny also offers features that go beyond the dipole approximation of a quantum spin via the theory of SU( N) coherent states. This can be especially useful for systems with strong single-ion anisotropy, as demonstrated in the FeI₂ tutorial.