{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Calculate magnetic parameters" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "grogupy is specifically written in a modular way, to fully take advantage of python classes. In this simple example we show how to use grogupy to simulate the magnetic interaction between two atoms from a density functional theory (DFT) calculation and create an output file for [magnopy](https://docs.magnopy.org/en/latest/user-guide/start/about.html). In this example we will use a non-collinear Siesta calculation of the Fe3GeTe2 system with spin-orbit coupling." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "First import grogupy, which will import the most important classes, functions and variables in its namespace." ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "/Users/danielpozsar/Documents/studies/elte/phd/grogu/.venv/lib/python3.12/site-packages/grogupy/_tqdm.py:32: TqdmExperimentalWarning: Using `tqdm.autonotebook.tqdm` in notebook mode. Use `tqdm.tqdm` instead to force console mode (e.g. in jupyter console)\n", " from tqdm.autonotebook import tqdm\n" ] } ], "source": [ "import grogupy" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we start to create the basic classes based on the desired parameters. First create the ``Kspace`` class which contains the parameters fot the Brillouin zone integration. Because in this case we only take a single layer of Fe3GeTe2 and there is a large vacuum in the perpendicular direction it is sufficient to only \n", "integrate in the plane of the material. Furthermore for a fast initial calculation a 10x10 grid of k-points should be enough, but this is nowhere near converged." ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 2, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Fe3GeTe2_kspace = grogupy.Kspace(kset=[10, 10, 1])\n", "Fe3GeTe2_kspace" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we would like to set the parameters for the complex integral for the Green's \n", "function calculation. This can be done through the ``Contour`` class. Again, for an \n", "initial calculation 300 sample points should be enough. The energy minimum should be \n", "set below the smallest eigenvalue from the Siesta Hamiltonian.\n" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Fe3GeTe2_contour = grogupy.Contour(emin=-15, eset=300, esetp=10000)\n", "Fe3GeTe2_contour" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we can create the ``Hamiltonian`` class, which extracts and stores the Hamiltonian \n", "and geometrical data from the Siesta files using the ``sisl`` library. Furthermore we must \n", "provide the exchange field orientation in the DFT calculation, which is usually the \n", "perpendicular direction from the 2D material.\n" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "application/vnd.jupyter.widget-view+json": { "model_id": "94512d084447492b9d5c730842ba066b", "version_major": 2, "version_minor": 0 }, "text/plain": [ "Setting up Hamiltonian: 0%| | 0/243 [00:00" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "Fe3GeTe2_hamiltonian = grogupy.Hamiltonian(\n", " infile=\"./../../../benchmarks/Fe3GeTe2/Fe3GeTe2.fdf\",\n", " scf_xcf_orientation=[0, 0, 1],\n", ")\n", "Fe3GeTe2_hamiltonian" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "After these initial steps we can create the ``Builder`` class and initialize with the \n", "above information. We have to provide the directions of the rotated exchange field and the \n", "two perpendicular directions are calculated automatically.\n" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "orientations = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]\n", "\n", "Fe3GeTe2 = grogupy.Builder(ref_xcf_orientations=orientations)\n", "Fe3GeTe2.add_kspace(Fe3GeTe2_kspace)\n", "Fe3GeTe2.add_contour(Fe3GeTe2_contour)\n", "Fe3GeTe2.add_hamiltonian(Fe3GeTe2_hamiltonian)\n", "Fe3GeTe2" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we have to define the magnetic entities, which can be a list of orbitals from the \n", "Siesta Hamiltonian. We use ``sisl`` functions to extract this information based on a more \n", "human like definition. For example take two shells from two Fe atoms.\n" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "magnetic_entities = [dict(atom=3, l=2), dict(atom=4, l=2)]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "With this information we can create two ``MagneticEntity`` instances.\n" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "application/vnd.jupyter.widget-view+json": { "model_id": "0a7245040d864ecca1c35591a22804a2", "version_major": 2, "version_minor": 0 }, "text/plain": [ "Add magnetic entities: 0%| | 0/2 [00:00" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "pairs = [\n", " dict(ai=0, aj=1, Ruc=[0, 0, 0]),\n", " dict(ai=0, aj=1, Ruc=[1, 0, 0]),\n", "]\n", "Fe3GeTe2.add_pairs(pairs)\n", "Fe3GeTe2" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now every information is contained in the Fe3GeTe2 instance to run the simulation. \n", "There are multiple parameters that can be changed to tune the runtime and the precision \n", "of the simulation, but for now let us use the default parameters." ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "application/vnd.jupyter.widget-view+json": { "model_id": "175a5730dfe943328a78bd423b92ae55", "version_major": 2, "version_minor": 0 }, "text/plain": [ "Extracting exchange field: 0%| | 0/165 [00:00 [[ 0. 0. -1. ]\n", " [ 0. 1. 0. ]\n", " [ 0. 0.70710678 -0.70710678]]\n", "[0 1 0] --> [[ 1. 0. 0. ]\n", " [ 0. 0. -1. ]\n", " [ 0.70710678 0. -0.70710678]]\n", "[0 0 1] --> [[1. 0. 0. ]\n", " [0. 1. 0. ]\n", " [0.70710678 0.70710678 0. ]]\n", "========================================\n", "Parameters for the Brillouin zone sampling:\n", "Number of k points: 100\n", "K points in each directions: [10 10 1]\n", "Parameters for the contour integral:\n", "Eset: 300\n", "Esetp: 10000\n", "Ebot: -20\n", "Etop: 0\n", "========================================\n", "\n" ] } ], "source": [ "print(Fe3GeTe2)" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "================================================================================\n", "GROGU INFORMATION\n", "================================================================================\n", "Hamiltonian convention\n", "Double counting true\n", "Normalized spins true\n", "Intra-atomic factor +1\n", "Exchange factor +0.5\n", "================================================================================\n", "Cell (Ang)\n", "3.791001511088653242e+00 0.000000000000000000e+00 0.000000000000000000e+00\n", "-1.895500755544326621e+00 3.283103614407953064e+00 0.000000000000000000e+00\n", "0.000000000000000000e+00 0.000000000000000000e+00 2.057000819825037041e+01\n", "================================================================================\n", "Magnetic sites\n", "Number of sites 2\n", "Name x (Ang) y (Ang) z (Ang) s sx sy sz\n", "3Fe(l:2) 1.981762848288204e-06 -9.134332829322269e-08 11.653315176909826 2.011394466844794 3.205365646820619e-05 -0.00025700339056116973 0.9999999664609096\n", "4Fe(l:2) 1.990383996147767e-06 -1.089291427685217e-07 8.916695475519404 2.0114468870381006 -3.0682257099277274e-05 0.00025456982487434966 0.999999967126401\n", "================================================================================\n", "Intra-atomic anisotropy tensor (meV)\n", "--------------------------------------------------------------------------------\n", "3Fe(l:2)\n", "Matrix\n", " 0.19805059848871515 -2.5314299080536617e-05 0.00039210954835139367\n", " -2.5314299080536617e-05 1.0811832564281654 -0.1178227120808917\n", " 0.00039210954835139367 -0.1178227120808917 0.0\n", "--------------------------------------------------------------------------------\n", "4Fe(l:2)\n", "Matrix\n", " -0.6983450912837961 -2.545067293718839e-05 0.0003488429527985503\n", " -2.545067293718839e-05 1.0856995486172 0.11931122379468805\n", " 0.0003488429527985503 0.11931122379468805 0.0\n", "--------------------------------------------------------------------------------\n", "================================================================================\n", "Exchange tensor (meV)\n", "Number of pairs 2\n", "--------------------------------------------------------------------------------\n", "Name1 Name2 i j k d (Ang)\n", "--------------------------------------------------------------------------------\n", "3Fe(l:2) 4Fe(l:2) 0 0 0 2.7366197013904223\n", "Matrix\n", " -24.798782054677186 7.711382752243433e-05 -0.015536801936793711\n", " 7.702738396723201e-05 -80.07760111320495 2.4364059494482864\n", " 0.012837383630434405 -2.4372750205623954 -73.68415956008036\n", "--------------------------------------------------------------------------------\n", "3Fe(l:2) 4Fe(l:2) 1 0 0 4.675551295032514\n", "Matrix\n", " -5.912430191790732 5.230909275184269 -0.3293938893776426\n", " -5.230901487243636 0.014705197513052002 2.6937771914387367\n", " 1.9676391897287149 -2.6923809965197587 -10.339832228712105\n", "--------------------------------------------------------------------------------\n", "================================================================================\n", "\n" ] } ], "source": [ "print(Fe3GeTe2.to_magnopy(comments=False))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "More information can be found in the Tutorials." ] } ], "metadata": { "kernelspec": { "display_name": ".venv", "language": "python", "name": "python3" }, "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.12.10" } }, "nbformat": 4, "nbformat_minor": 2 }