{ "cells": [ { "cell_type": "code", "metadata": {}, "source": [], "outputs": [], "execution_count": null }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "\n", "# Noncollinear magnetism - VI$_2$\n", "\n", "For this part of the project, we will move on to a magnetic monolayer in which the dominant exchange coupling is antiferromagnetic, namely VI$_2$. We will still use the localized spin Hamiltonian\n", "\n", "$$H = -\\frac{1}{2}\\sum_{i,j}J_{ij}\\mathbf{S}_i\\cdot \\mathbf{S}_j+A\\sum_i(S_i^z)^2$$\n", "\n", "but here in the antiferromagnetic case, $J<0$.\n", "\n", "## Optimizing the atomic structure\n", "\n", "Since we will need to do LDA calculations later on, we will start of this part of the project by relaxing the atomic structure of VI$_2$ using the LDA functional. Usually, the difference in crystal structure between different magnetically ordered states is small, so we will perform the relaxation in the ferromagnetic state, which has a smaller unit cell.\n", "\n", "1. First you should download the relaxed PBE crystal structure. Either, browse the C2DB at https://cmrdb.fysik.dtu.dk/c2db and download the `.xyz` file for VI$_2$ or dowload it directly from the summer school tutorial website [here](https://wiki.fysik.dtu.dk/gpaw/summerschools/summerschool22/magnetism/magnetism.html)\n", "2. Fill in the expected ionic value for the V spins `S` below and run the cell to relax the crystal structure. The calculation takes about 17 minutes. (Hint: V has the electronic configuration [Ar]3d$^3$4s$^2$)\n", "\n" ] }, { "cell_type": "code", "metadata": {}, "source": [ "from ase.io import read\n", "from ase.visualize import view\n", "from ase.optimize import BFGS\n", "from ase.filters import UnitCellFilter\n", "from gpaw import GPAW, PW\n", "\n", "S = ???\n", "m = S * 2\n", "layer = read('VI2.xyz')\n", "layer.set_initial_magnetic_moments([m, 0, 0])\n", "view(layer)\n", "\n", "calc = GPAW(mode=PW(400), xc='LDA', kpts=(4, 4, 1))\n", "layer.calc = calc\n", "\n", "uf = UnitCellFilter(layer, mask=[1, 1, 0, 0, 0, 1])\n", "opt = BFGS(uf)\n", "opt.run(fmax=0.1)\n", "\n", "calc.set(symmetry='off')\n", "calc.get_potential_energy()\n", "\n", "calc.write('VI2_relaxed.gpw')" ], "outputs": [], "execution_count": null }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "## Magnetic anisotropy\n", "\n", "In the cell above, we not only performed a structural optimization, but we also took the relaxed structure, switched off the $k$-point symmetries and did a final DFT calculation with all the $k$-points in the Brillouin zone. We did this in order to be able to evaluate the the magnetic anisotropy arising from the spin-orbit coupling.\n", "\n", "1. Adapt the code you used for CrI$_3$ to calculate the magnetic anisotropy in the cell below.\n", "2. Is the easy-axis in-plane or out-of-plane?\n", "3. Do you expect to the VI$_2$ to exhibit magnetic order at finite temperatures?\n", "\n" ] }, { "cell_type": "code", "metadata": {}, "source": [ "from gpaw.spinorbit import soc_eigenstates\n" ], "outputs": [], "execution_count": null }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "## DFT calculations in a repeated cell\n", "\n", "To realize that VI$_2$ is in fact an antiferromagnetic material, we need to do a calculation starting from an antiferromagnetic alignment of the spins. To do so, we need more than a single V atom in the unit cell. In the cell below, the atomic structure is repeated once to obtain two V atoms in the unit cell and it is shown how do a DFT calculation for the antiferromagnetic state.\n", "\n", "1. Fill in the `kpts` for the repeated unit cell.\n", "2. Replace the `...` with code to calculate the ferromagnetic state.\n", "3. When you have finalzed the code for the ferromagnetic state, run the cell. The calculation takes about 7 minutes.\n", "\n" ] }, { "cell_type": "code", "metadata": {}, "source": [ "layer_afm = layer.repeat((2, 1, 1))\n", "layer_afm.set_initial_magnetic_moments([m, 0, 0, -m, 0, 0])\n", "view(layer_afm)\n", "\n", "calc = GPAW(mode=PW(400),\n", " xc='LDA',\n", " kpts=???))\n", "layer_afm.calc = calc\n", "layer_afm.get_potential_energy()\n", "calc.write('V2I4_afm.gpw')\n", "\n", "...\n" ], "outputs": [], "execution_count": null }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "## Calculating J\n", "\n", "Finally, we are in a position to calculate the nearest neighbour Heisenberg exchange coupling $J$. Before doing so, please compare the output of the antiferromagnetic and ferromagnetic calculations from the cell above.\n", "\n", "1. Which state has the lowest energy?\n", "2. What will the sign of $J$ be?\n", "\n", "You should find that the antiferromagnetic state has a lower energy than the ferromagnetic one, but does this also mean that the calculated configuration is the correct magnetic ground state? Rather, it implies that the spins prefer to be antialigned, i.e. that the exchange coupling $J$ is negative.\n", "\n", "3. Draw the structural arrangement of the V atoms on a piece of paper. Which type of magnetic lattice do they form?\n", "4. Fill in the spin configuration of the magnetic lattice and convince yourself that all spins cannot be antialigned to their nearest neighbours.\n", "\n", "The latter finding means that the antiferromagnetic system is frustrated and the antiferromagnetic configuration we have computed will not be the true ground state of the system.\n", "\n", "Leaving the magnetic frustration aside for the moment, we will first calculate $J$, which we can still do even though that we have not found the ground state yet.\n", "\n", "5. Use the classical Heisenberg model with nearest neighbor interaction only (and $A=0$) to derive the energy per magnetic site of the ferromagnetic and antiferromagnetic configurations respectively.\n", "\n", "You should obtain the following:\n", "\n", "$$E_{\\mathrm{FM}} = E_0 - 3JS^2$$\n", "\n", "and\n", "\n", "$$E_{\\mathrm{AFM}} = E_0 + JS^2$$\n", "\n", "where $E_0$ is some reference energy.\n", "\n", "6. Use these expressions to eliminate $E_0$ and express $J$ in terms of the energy difference per magnetic site of the two configurations.\n", "7. Write code to extract `E_fm` and `E_afm` in the cell below.\n", "8. Fill in the formula for `J` and evaluate the cell to calculate it.\n", "\n", "You should get a value for $J$ around -1.4 meV.\n" ] }, { "cell_type": "code", "metadata": {}, "source": [ "E_fm = ???\n", "E_afm = ???\n", "\n", "J = ???\n", "print(f'J = {J * 1000:1.2f} meV')" ], "outputs": [], "execution_count": null }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "## Noncollinear configuration\n", "\n", "As it turn out, the optimal spin structure on a trigonal lattice with antiferromagntice exchange coupling is to place all spins at 120$^\\circ$ angles with respect its neighbors.\n", "\n", "1. Draw this structure and convince yourself that it is indeed possible to put every spin at a 120$^\\circ$ angle with respect to its neighbors.\n", "2. What is the minimal number of magnetic atoms required in the magnetic unit cell to represent such a state?\n", "3. Verrify that the Heisenberg model with classical spins gives a lower energy with this configuration than in the antialigned structure calculated above.\n", "\n", "You should obtain the following energy for the 120$^\\circ$ noncolinnear configuration:\n", "\n", "$$E_{\\mathrm{NC}}=E_0+\\frac{3}{2}JS^2.$$\n", "\n", "We will now check if we can verify this prediction within the LSDA. To do that we need to perform a noncollinear DFT calculation, which is done in the cell below.\n", "\n", "4. Read and try to understand the code to perform a noncollinear LSDA calculation in the 120$^\\circ$ noncolinnear configuration.\n", "5. Replace the `...` with code to make a noncolinnear LSDA calculation for the ferromagnetic state as well (we will need this for a late comparison).\n", "\n", "Run the cell and verify that the energy per magnetic atom is lower in the 120$^\\circ$ noncolinnear configuration compared to both of the previous calculated states.\n" ] }, { "cell_type": "code", "metadata": {}, "source": [ "import numpy as np\n", "from ase.io import read\n", "from ase.visualize import view\n", "from gpaw import GPAW, PW, MixerDif\n", "\n", "m = 3\n", "cell_cv = layer.get_cell()\n", "layer_nc = layer.repeat((3, 1, 1))\n", "new_cell_cv = [[3 * cell_cv[0, 0] / 2, 3**0.5 * cell_cv[0, 0] / 2, 0.0],\n", " [3 * cell_cv[0, 0] / 2, -3**0.5 * cell_cv[0, 0] / 2, 0.0],\n", " [0.0, 0.0, cell_cv[2, 2]]]\n", "layer_nc.set_cell(new_cell_cv)\n", "view(layer_nc)\n", "\n", "magmoms = np.zeros((len(layer_nc), 3), float)\n", "magmoms[0] = [m, 0, 0]\n", "magmoms[3] = [m * np.cos(2 * np.pi / 3), m * np.sin(2 * np.pi / 3), 0]\n", "magmoms[6] = [m * np.cos(2 * np.pi / 3), -m * np.sin(2 * np.pi / 3), 0]\n", "\n", "calc = GPAW(mode=PW(400),\n", " xc='LDA',\n", " mixer=MixerDif(),\n", " symmetry='off',\n", " experimental={'magmoms': magmoms, 'soc': False},\n", " parallel={'domain': 1, 'band': 1},\n", " kpts=(2, 2, 1))\n", "layer_nc.calc = calc\n", "layer_nc.get_potential_energy()\n", "calc.write('nc_nosoc.gpw')\n", "\n", "...\n" ], "outputs": [], "execution_count": null }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "## Anisotropy and exchange coupling from noncollinear DFT\n", "\n", "In the cell above we could have included spin-orbit coupling in a self-consistent way by setting 'soc'=True'. However, it is more convenient for us to exclude it such that we can explicitly compute the single-ion anisotropy parameter $A$ afterwards. If the localized spin Hamiltonian with nearest neighbor exchange interactions is a good model, we should be able to obtain both $J$ and $A$ from the noncollinear calculation as well.\n", "\n", "1. Fill in the cell below with code to compute the nearest neighbour exchange coupling $J$ based on the noncollinear calculations performed above.\n", "2. Fill in code that computes the spin-orbit coupling corrected energies of the ferromagnetic state from the noncollinear calculation with spins directed along the $x$, $y$ and $z$ directions (Hint: the anisotropy is calculated by rotating the entire initial spin configuration first by $\\theta$ and then by $\\varphi$).\n", "3. Repeat point 2, but for the 120$^\\circ$ noncolinnear configuration. In this case, you cannot align all spins to one direction, but takes as a reference the first V atom.\n", "\n", "With the above code in place, please evaluate the cell. You should obtain a nearest neighbour exchange coupling of about -1.9 meV.\n", "\n" ] }, { "cell_type": "code", "metadata": {}, "source": [ "\n" ], "outputs": [], "execution_count": null }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n", "## Critical temperature?\n", "\n", "Based on the noncollinear calculations above:\n", "\n", "1. What is the easy axis of the system?\n", "2. Does the easy axis agree with your initial findings for the simple ferromagnetic state?\n", "\n", "You should find that the energy of the 120$^\\circ$ noncolinnear configuration does not depend strongly on whether the first V atom is aligned to the $x$ or the $y$ direction.\n", "\n", "3. If we assume full in-plane isotropy is there any rotational freedom left in the noncollinear ground state?\n", "4. What implications does this have for the critical temperature of the monolayer?\n", "\n", "You might be able to convince yourself that some degree of in-plane anisotropy is required to obtain a finite critical temperature for the 120$^\\circ$ noncolinnear magnetic order. Again, bear in mind that all the calculations in the present notebook ought to be properly converged with respect to $k$-points, plane wave cutoff etc. to achieve an accurate estimate of e.g. the in-plane anisotropy.\n", "\n", "Clearly the noncollinear spin state of VI$_2$ is more difficult to describe than the ferromagnetic state in CrI$_3$ and we do not yet have a simple theoretical expression for the critical temperature as a function of anisotropy and exchange coupling constants. However, with the rapid development of experimental techniques to synthesize and characterize 2D materials it does seem plausible that such a noncollinear 2D material may be observed in the future.\n" ] } ], "metadata": { "kernelspec": { "display_name": "CAMD2022", "language": "python", "name": "camd2022" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.6.1" } }, "nbformat": 4, "nbformat_minor": 1 }