{ "cells": [ { "cell_type": "code", "execution_count": 21, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "0.14.3\n", "1.24.4\n" ] } ], "source": [ "from sys import getsizeof\n", "from timeit import default_timer as timer\n", "\n", "import sisl\n", "import sisl.viz\n", "from src.grogu_magn import *\n", "from mpi4py import MPI\n", "import warnings\n", "\n", "# runtime information\n", "times = dict()\n", "times[\"start_time\"] = timer()\n", "########################\n", "# it works if data is in downloads folder\n", "########################\n", "sisl.__version__\n", "\n", "try:\n", " print(sisl.__version__)\n", "except:\n", " print(\"sisl version unknown.\")\n", "\n", "try:\n", " print(np.__version__)\n", "except:\n", " print(\"numpy version unknown.\")" ] }, { "cell_type": "code", "execution_count": 22, "metadata": {}, "outputs": [], "source": [ "sisl.io.Sile(\"input.fdf\").read(\"scfOrientation\")" ] }, { "cell_type": "code", "execution_count": 23, "metadata": {}, "outputs": [], "source": [ "################################################################################\n", "#################################### INPUT #####################################\n", "################################################################################\n", "path = (\n", " \"/Users/danielpozsar/Downloads/nojij/Fe3GeTe2/monolayer/soc/lat3_791/Fe3GeTe2.fdf\"\n", ")\n", "outfile = \"./Fe3GeTe2_notebook\"\n", "\n", "# this information needs to be given at the input!!\n", "scf_xcf_orientation = np.array([0, 0, 1]) # z\n", "# list of reference directions for around which we calculate the derivatives\n", "# o is the quantization axis, v and w are two axes perpendicular to it\n", "# at this moment the user has to supply o,v,w on the input.\n", "# we can have some default for this\n", "ref_xcf_orientations = [\n", " dict(o=np.array([1, 0, 0]), vw=[np.array([0, 1, 0]), np.array([0, 0, 1])]),\n", " dict(o=np.array([0, 1, 0]), vw=[np.array([1, 0, 0]), np.array([0, 0, 1])]),\n", " dict(o=np.array([0, 0, 1]), vw=[np.array([1, 0, 0]), np.array([0, 1, 0])]),\n", "]\n", "magnetic_entities = [\n", " dict(atom=3, l=2),\n", " dict(atom=4, l=2),\n", " dict(atom=5, l=2),\n", "]\n", "pairs = [\n", " dict(ai=0, aj=1, Ruc=np.array([0, 0, 0])),\n", " dict(ai=0, aj=2, Ruc=np.array([0, 0, 0])),\n", " dict(ai=1, aj=2, Ruc=np.array([0, 0, 0])),\n", " dict(ai=0, aj=2, Ruc=np.array([-1, -1, 0])),\n", " dict(ai=1, aj=2, Ruc=np.array([-1, -1, 0])),\n", " dict(ai=0, aj=2, Ruc=np.array([-1, 0, 0])),\n", " dict(ai=1, aj=2, Ruc=np.array([-1, 0, 0])),\n", " dict(ai=1, aj=2, Ruc=np.array([-2, 0, 0])),\n", " dict(ai=1, aj=2, Ruc=np.array([-3, 0, 0])),\n", "]\n", "\n", "# Brilloun zone sampling and Green function contour integral\n", "kset = 3\n", "kdirs = \"xy\"\n", "ebot = -13\n", "eset = 300\n", "esetp = 1000\n", "################################################################################\n", "#################################### INPUT #####################################\n", "################################################################################" ] }, { "cell_type": "code", "execution_count": 24, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "================================================================================================================================================================\n", "Input file: \n", "/Users/danielpozsar/Downloads/nojij/Fe3GeTe2/monolayer/soc/lat3_791/Fe3GeTe2.fdf\n", "Output file: \n", "./Fe3GeTe2_notebook.pickle\n", "Number of nodes in the parallel cluster: 1\n", "================================================================================================================================================================\n", "Cell [Ang]: \n", "[[ 3.79100000e+00 0.00000000e+00 0.00000000e+00]\n", " [-1.89550000e+00 3.28310231e+00 0.00000000e+00]\n", " [ 1.25954923e-15 2.18160327e-15 2.05700000e+01]]\n", "================================================================================================================================================================\n", "DFT axis: \n", "[0 0 1]\n", "Quantization axis and perpendicular rotation directions:\n", "[1 0 0] --» [array([0, 1, 0]), array([0, 0, 1])]\n", "[0 1 0] --» [array([1, 0, 0]), array([0, 0, 1])]\n", "[0 0 1] --» [array([1, 0, 0]), array([0, 1, 0])]\n", "================================================================================================================================================================\n", "Parameters for the contour integral:\n", "Number of k points: 3\n", "k point directions: xy\n", "Ebot: -13\n", "Eset: 300\n", "Esetp: 1000\n", "================================================================================================================================================================\n", "Setup done. Elapsed time: 2435.031293583 s\n", "================================================================================================================================================================\n" ] } ], "source": [ "# MPI parameters\n", "comm = MPI.COMM_WORLD\n", "size = comm.Get_size()\n", "rank = comm.Get_rank()\n", "root_node = 0\n", "\n", "# rename outfile\n", "if not outfile.endswith(\".pickle\"):\n", " outfile += \".pickle\"\n", "\n", "simulation_parameters = dict(\n", " infile=path,\n", " outfile=outfile,\n", " scf_xcf_orientation=scf_xcf_orientation,\n", " ref_xcf_orientations=ref_xcf_orientations,\n", " kset=kset,\n", " kdirs=kdirs,\n", " ebot=ebot,\n", " eset=eset,\n", " esetp=esetp,\n", " parallel_size=size,\n", ")\n", "\n", "# if ebot is not given put it 0.1 eV under the smallest energy\n", "if simulation_parameters[\"ebot\"] is None:\n", " try:\n", " eigfile = simulation_parameters[\"infile\"][:-3] + \"EIG\"\n", " simulation_parameters[\"ebot\"] = read_siesta_emin(eigfile) - 0.1\n", " except:\n", " print(\"Could not determine ebot.\")\n", " print(\"Parameter was not given and .EIG file was not found.\")\n", "# digestion of the input\n", "# read sile\n", "fdf = sisl.get_sile(simulation_parameters[\"infile\"])\n", "# read in hamiltonian\n", "dh = fdf.read_hamiltonian()\n", "simulation_parameters[\"cell\"] = fdf.read_geometry().cell\n", "\n", "# unit cell index\n", "uc_in_sc_idx = dh.lattice.sc_index([0, 0, 0])\n", "\n", "if rank == root_node:\n", " print_parameters(simulation_parameters)\n", " times[\"setup_time\"] = timer()\n", " print(f\"Setup done. Elapsed time: {times['setup_time']} s\")\n", " print(\n", " \"================================================================================================================================================================\"\n", " )" ] }, { "cell_type": "code", "execution_count": 25, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Hamiltonian and exchange field rotated. Elapsed time: 2435.900105791 s\n", "================================================================================================================================================================\n" ] } ], "source": [ "hh, ss, NO = build_hh_ss(dh)\n", "\n", "\n", "# symmetrizing Hamiltonian and overlap matrix to make them hermitian\n", "for i in range(dh.lattice.sc_off.shape[0]):\n", " j = dh.lattice.sc_index(-dh.lattice.sc_off[i])\n", " h1, h1d = hh[i], hh[j]\n", " hh[i], hh[j] = (h1 + h1d.T.conj()) / 2, (h1d + h1.T.conj()) / 2\n", " s1, s1d = ss[i], ss[j]\n", " ss[i], ss[j] = (s1 + s1d.T.conj()) / 2, (s1d + s1.T.conj()) / 2\n", "\n", "\n", "# identifying TRS and TRB parts of the Hamiltonian\n", "TAUY = np.kron(np.eye(NO), tau_y)\n", "hTR = np.array([TAUY @ hh[i].conj() @ TAUY for i in range(dh.lattice.nsc.prod())])\n", "hTRS = (hh + hTR) / 2\n", "hTRB = (hh - hTR) / 2\n", "\n", "# extracting the exchange field\n", "traced = [spin_tracer(hTRB[i]) for i in range(dh.lattice.nsc.prod())] # equation 77\n", "XCF = np.array(\n", " [\n", " np.array([f[\"x\"] / 2 for f in traced]),\n", " np.array([f[\"y\"] / 2 for f in traced]),\n", " np.array([f[\"z\"] / 2 for f in traced]),\n", " ]\n", ") # equation 77\n", "\n", "\n", "# Check if exchange field has scalar part\n", "max_xcfs = abs(np.array(np.array([f[\"c\"] / 2 for f in traced]))).max()\n", "if max_xcfs > 1e-12:\n", " warnings.warn(\n", " f\"Exchange field has non negligible scalar part. Largest value is {max_xcfs}\"\n", " )\n", "\n", "if rank == root_node:\n", " times[\"H_and_XCF_time\"] = timer()\n", " print(\n", " f\"Hamiltonian and exchange field rotated. Elapsed time: {times['H_and_XCF_time']} s\"\n", " )\n", " print(\n", " \"================================================================================================================================================================\"\n", " )" ] }, { "cell_type": "code", "execution_count": 26, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Site and pair dictionaries created. Elapsed time: 2438.737295 s\n", "================================================================================================================================================================\n" ] } ], "source": [ "pairs, magnetic_entities = setup_pairs_and_magnetic_entities(\n", " magnetic_entities, pairs, dh, simulation_parameters\n", ")\n", "\n", "if rank == root_node:\n", " times[\"site_and_pair_dictionaries_time\"] = timer()\n", " print(\n", " f\"Site and pair dictionaries created. Elapsed time: {times['site_and_pair_dictionaries_time']} s\"\n", " )\n", " print(\n", " \"================================================================================================================================================================\"\n", " )" ] }, { "cell_type": "code", "execution_count": 27, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "k loop: 0%| | 0/9 [00:00iklm\", R, XCF)\n", " rot_H_XCF = sum(\n", " [np.kron(rot_XCF[i], tau) for i, tau in enumerate([tau_x, tau_y, tau_z])]\n", " )\n", " rot_H_XCF_uc = rot_H_XCF[uc_in_sc_idx]\n", "\n", " # obtain total Hamiltonian with the rotated exchange field\n", " rot_H = hTRS + rot_H_XCF # equation 76\n", "\n", " hamiltonians.append(\n", " dict(\n", " orient=orient[\"o\"],\n", " H=rot_H,\n", " GS=np.zeros(\n", " (simulation_parameters[\"eset\"], rot_H.shape[1], rot_H.shape[2]),\n", " dtype=\"complex128\",\n", " ),\n", " GS_tmp=np.zeros(\n", " (simulation_parameters[\"eset\"], rot_H.shape[1], rot_H.shape[2]),\n", " dtype=\"complex128\",\n", " ),\n", " )\n", " ) # store orientation and rotated Hamiltonian\n", "\n", " # these are the rotations (for now) perpendicular to the quantization axis\n", " for u in orient[\"vw\"]:\n", " Tu = np.kron(np.eye(NO, dtype=int), tau_u(u)) # section 2.H\n", "\n", " Vu1, Vu2 = calc_Vu(rot_H_XCF_uc, Tu)\n", "\n", " for mag_ent in magnetic_entities:\n", " idx = mag_ent[\"spin_box_indeces\"]\n", " # fill up the perturbed potentials (for now) based on the on-site projections\n", " mag_ent[\"Vu1\"][i].append(Vu1[:, idx][idx, :])\n", " mag_ent[\"Vu2\"][i].append(Vu2[:, idx][idx, :])\n", "\n", "if rank == root_node:\n", " times[\"reference_rotations_time\"] = timer()\n", " print(\n", " f\"Rotations done perpendicular to quantization axis. Elapsed time: {times['reference_rotations_time']} s\"\n", " )\n", " print(\n", " \"================================================================================================================================================================\"\n", " )" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Starting matrix inversions.\n", "Total number of k points: 9\n", "Number of energy samples per k point: 300\n", "Total number of directions: 3\n", "Total number of matrix inversions: 8100\n", "The shape of the Hamiltonian and the Greens function is 84x84=7056\n", "Memory taken by a single Hamiltonian is: 0.015625 KB\n", "Expected memory usage per matrix inversion: 0.5 KB\n", "Expected memory usage per k point for parallel inversion: 450.0 KB\n", "Expected memory usage on root node: 3.955078125 MB\n", "================================================================================================================================================================\n" ] }, { "name": "stderr", "output_type": "stream", "text": [ "k loop: 100%|██████████| 9/9 [00:31<00:00, 3.50s/it]" ] }, { "name": "stdout", "output_type": "stream", "text": [ "Calculated Greens functions. Elapsed time: 34.685783333 s\n", "================================================================================================================================================================\n" ] }, { "name": "stderr", "output_type": "stream", "text": [ "\n" ] } ], "source": [ "if rank == root_node:\n", " print(\"Starting matrix inversions.\")\n", " print(f\"Total number of k points: {kset.shape[0]}\")\n", " print(f\"Number of energy samples per k point: {simulation_parameters['eset']}\")\n", " print(f\"Total number of directions: {len(hamiltonians)}\")\n", " print(\n", " f\"Total number of matrix inversions: {kset.shape[0] * len(hamiltonians) * simulation_parameters['eset']}\"\n", " )\n", " print(f\"The shape of the Hamiltonian and the Greens function is {NO}x{NO}={NO*NO}\")\n", " # https://stackoverflow.com/questions/70746660/how-to-predict-memory-requirement-for-np-linalg-inv\n", " # memory is O(64 n**2) for complex matrices\n", " memory_size = getsizeof(hamiltonians[0][\"H\"].base) / 1024\n", " print(\n", " f\"Memory taken by a single Hamiltonian is: {getsizeof(hamiltonians[0]['H'].base) / 1024} KB\"\n", " )\n", " print(f\"Expected memory usage per matrix inversion: {memory_size * 32} KB\")\n", " print(\n", " f\"Expected memory usage per k point for parallel inversion: {memory_size * len(hamiltonians) * simulation_parameters['eset'] * 32} KB\"\n", " )\n", " print(\n", " f\"Expected memory usage on root node: {len(np.array_split(kset, size)[0]) * memory_size * len(hamiltonians) * simulation_parameters['eset'] * 32 / 1024} MB\"\n", " )\n", " print(\n", " \"================================================================================================================================================================\"\n", " )\n", "\n", "comm.Barrier()\n", "# ----------------------------------------------------------------------\n", "\n", "# make energy contour\n", "# we are working in eV now !\n", "# and sisl shifts E_F to 0 !\n", "cont = make_contour(\n", " emin=simulation_parameters[\"ebot\"],\n", " enum=simulation_parameters[\"eset\"],\n", " p=simulation_parameters[\"esetp\"],\n", ")\n", "eran = cont.ze\n", "\n", "# ----------------------------------------------------------------------\n", "# sampling the integrand on the contour and the BZ\n", "for k in kpcs[rank]:\n", " wk = wkset[rank] # weight of k point in BZ integral\n", " # iterate over reference directions\n", " for i, hamiltonian_orientation in enumerate(hamiltonians):\n", " # calculate Greens function\n", " H = hamiltonian_orientation[\"H\"]\n", " HK, SK = hsk(H, ss, dh.sc_off, k)\n", "\n", " # solve Greens function sequentially for the energies, because of memory bound\n", " Gk = sequential_GK(HK, SK, eran, simulation_parameters[\"eset\"])\n", "\n", " # saving this for total charge\n", " hamiltonian_orientation[\"GS_tmp\"] += Gk @ SK * wk\n", "\n", " # store the Greens function slice of the magnetic entities (for now) based on the on-site projections\n", " for mag_ent in magnetic_entities:\n", " mag_ent[\"Gii_tmp\"][i] += (\n", " Gk[:, mag_ent[\"spin_box_indeces\"], :][:, :, mag_ent[\"spin_box_indeces\"]]\n", " * wk\n", " )\n", "\n", " for pair in pairs:\n", " # add phase shift based on the cell difference\n", " phase = np.exp(1j * 2 * np.pi * k @ pair[\"Ruc\"].T)\n", "\n", " # get the pair orbital sizes from the magnetic entities\n", " ai = magnetic_entities[pair[\"ai\"]][\"spin_box_indeces\"]\n", " aj = magnetic_entities[pair[\"aj\"]][\"spin_box_indeces\"]\n", "\n", " # store the Greens function slice of the magnetic entities (for now) based on the on-site projections\n", " pair[\"Gij_tmp\"][i] += Gk[:, ai][..., aj] * phase * wk\n", " pair[\"Gji_tmp\"][i] += Gk[:, aj][..., ai] / phase * wk\n", "\n", "# summ reduce partial results of mpi nodes\n", "for i in range(len(hamiltonians)):\n", " # for total charge\n", " comm.Reduce(hamiltonians[i][\"GS_tmp\"], hamiltonians[i][\"GS\"], root=root_node)\n", "\n", " for mag_ent in magnetic_entities:\n", " comm.Reduce(mag_ent[\"Gii_tmp\"][i], mag_ent[\"Gii\"][i], root=root_node)\n", "\n", " for pair in pairs:\n", " comm.Reduce(pair[\"Gij_tmp\"][i], pair[\"Gij\"][i], root=root_node)\n", " comm.Reduce(pair[\"Gji_tmp\"][i], pair[\"Gji\"][i], root=root_node)\n", "\n", "if rank == root_node:\n", " times[\"green_function_inversion_time\"] = timer()\n", " print(\n", " f\"Calculated Greens functions. Elapsed time: {times['green_function_inversion_time']} s\"\n", " )\n", " print(\n", " \"================================================================================================================================================================\"\n", " )" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Total charge: 39.98745949985201\n", "Total charge: 39.987459508326964\n", "Total charge: 40.098563331816784\n", "Magnetic entities integrated.\n", "Pairs integrated.\n", "Magnetic parameters calculated.\n", "##################################################################### GROGU OUTPUT #############################################################################\n", "================================================================================================================================================================\n", "Input file: \n", "/Users/danielpozsar/Downloads/nojij/Fe3GeTe2/monolayer/soc/lat3_791/Fe3GeTe2.fdf\n", "Output file: \n", "./Fe3GeTe2_notebook.pickle\n", "Number of nodes in the parallel cluster: 1\n", "================================================================================================================================================================\n", "Cell [Ang]: \n", "[[ 3.79100000e+00 0.00000000e+00 0.00000000e+00]\n", " [-1.89550000e+00 3.28310231e+00 0.00000000e+00]\n", " [ 1.25954923e-15 2.18160327e-15 2.05700000e+01]]\n", "================================================================================================================================================================\n", "DFT axis: \n", "[0 0 1]\n", "Quantization axis and perpendicular rotation directions:\n", "[1 0 0] --» [array([0, 1, 0]), array([0, 0, 1])]\n", "[0 1 0] --» [array([1, 0, 0]), array([0, 0, 1])]\n", "[0 0 1] --» [array([1, 0, 0]), array([0, 1, 0])]\n", "================================================================================================================================================================\n", "Parameters for the contour integral:\n", "Number of k points: 3\n", "k point directions: xy\n", "Ebot: -13\n", "Eset: 300\n", "Esetp: 1000\n", "================================================================================================================================================================\n", "Atomic information: \n", "----------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "[atom index]Element(orbitals) x [Ang] y [Ang] z [Ang] Sx Sy Sz Q Lx Ly Lz Jx Jy Jz\n", "----------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "[3]Fe(2) -7.339158738013707e-06 4.149278510690423e-06 11.657585837928032\n", "\n", "[4]Fe(2) -7.326987662162937e-06 4.158274523275774e-06 8.912422537596708\n", "\n", "[5]Fe(2) 1.8954667088117545 1.0943913231921656 10.285002698393109\n", "\n", "================================================================================================================================================================\n", "Anisotropy [meV]\n", "----------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "Magnetic entity x [Ang] y [Ang] z [Ang]\n", "----------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "[3]Fe(2) -7.339158738013707e-06 4.149278510690423e-06 11.657585837928032\n", "Consistency check: 2.3860081328264116e-05\n", "Anisotropy diag: [0.02546754 0.00208022 0. ]\n", "\n", "[4]Fe(2) -7.326987662162937e-06 4.158274523275774e-06 8.912422537596708\n", "Consistency check: 7.768349833381372e-05\n", "Anisotropy diag: [-0.06673164 0.01049614 0. ]\n", "\n", "[5]Fe(2) 1.8954667088117545 1.0943913231921656 10.285002698393109\n", "Consistency check: 1.2052892692793193e-06\n", "Anisotropy diag: [-0.01943729 -0.0206467 0. ]\n", "\n", "================================================================================================================================================================\n", "Exchange [meV]\n", "----------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "Magnetic entity1 Magnetic entity2 [i j k] d [Ang]\n", "----------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "[3]Fe(2) [4]Fe(2) [0 0 0] d [Ang] 2.745163300331324\n", "Isotropic: -39.96748066866786\n", "DMI: [-4.16373454e-04 -8.52759845e-04 2.22666054e-07]\n", "Symmetric-anisotropy: [ 7.85170293e-01 -2.01786552e-05 -2.73439476e-07 -2.01786552e-05\n", " 7.85623794e-01 -9.79265975e-07 -2.73439476e-07 -9.79265975e-07\n", " -1.57079409e+00]\n", "J: [-3.91823104e+01 -2.01786552e-05 -2.73439476e-07 -2.01786552e-05\n", " -3.91818569e+01 -9.79265975e-07 -2.73439476e-07 -9.79265975e-07\n", " -4.15382748e+01]\n", "Energies for debugging: \n", "[[-4.15381735e-02 -4.15394188e-07 4.17352720e-07 -4.17291610e-02]\n", " [-4.15383760e-02 8.53033284e-07 -8.52486405e-07 -4.17300740e-02]\n", " [-3.66345528e-02 2.04013212e-08 1.99559891e-08 -3.66345467e-02]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.04173007 -0.03663455 -0.04153817]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.04173007401441678 -0.036634546737597175\n", "\n", "[3]Fe(2) [5]Fe(2) [0 0 0] d [Ang] 2.5835033632437767\n", "Isotropic: -54.8678926276239\n", "DMI: [ 8.10240490e-01 -1.36608099e+00 3.25215454e-05]\n", "Symmetric-anisotropy: [ 0.12903943 0.14426756 -0.06549931 0.14426756 -0.20896588 -0.05975815\n", " -0.06549931 -0.05975815 0.07992645]\n", "J: [-54.7388532 0.14426756 -0.06549931 0.14426756 -55.07685851\n", " -0.05975815 -0.06549931 -0.05975815 -54.78796617]\n", "Energies for debugging: \n", "[[-0.05497624 0.00087 -0.00075048 -0.0552075 ]\n", " [-0.0545997 0.00143158 -0.00130058 -0.05469815]\n", " [-0.05494622 -0.00014424 -0.0001443 -0.05477956]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.05469815 -0.05494622 -0.05497624]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.05469814671520505 -0.05477955968564471\n", "\n", "[4]Fe(2) [5]Fe(2) [0 0 0] d [Ang] 2.583501767937866\n", "Isotropic: -54.86287528016134\n", "DMI: [-8.08522805e-01 1.36388593e+00 3.14077512e-05]\n", "Symmetric-anisotropy: [ 0.12501449 0.14426595 0.05678873 0.14426595 -0.20956489 0.04990479\n", " 0.05678873 0.04990479 0.08455041]\n", "J: [-5.47378608e+01 1.44265952e-01 5.67887254e-02 1.44265952e-01\n", " -5.50724402e+01 4.99047869e-02 5.67887254e-02 4.99047869e-02\n", " -5.47783249e+01]\n", "Energies for debugging: \n", "[[-0.05496666 -0.00085843 0.00075862 -0.05519787]\n", " [-0.05458999 -0.00142067 0.0013071 -0.05469537]\n", " [-0.05494701 -0.00014423 -0.0001443 -0.05478035]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.05469537 -0.05494701 -0.05496666]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.054695368426152975 -0.05478035315936865\n", "\n", "[3]Fe(2) [5]Fe(2) [-1 -1 0] d [Ang] 2.5834973202859075\n", "Isotropic: -54.86725213337313\n", "DMI: [-1.55587908e+00 1.48356799e-05 -4.55339161e-05]\n", "Symmetric-anisotropy: [-3.76175058e-01 -2.58204985e-05 -2.58761719e-05 -2.58204985e-05\n", " 2.92919450e-01 4.97693839e-02 -2.58761719e-05 4.97693839e-02\n", " 8.32556083e-02]\n", "J: [-5.52434272e+01 -2.58204985e-05 -2.58761719e-05 -2.58204985e-05\n", " -5.45743327e+01 4.97693839e-02 -2.58761719e-05 4.97693839e-02\n", " -5.47839965e+01]\n", "Energies for debugging: \n", "[[-5.44093670e-02 -1.60564846e-03 1.50610969e-03 -5.44511949e-02]\n", " [-5.51586260e-02 1.10404920e-08 4.07118518e-08 -5.54561691e-02]\n", " [-5.46974705e-02 -1.97134176e-08 7.13544147e-08 -5.50306852e-02]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.05545617 -0.05469747 -0.05440937]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.0554561691350627 -0.05503068524779923\n", "\n", "[4]Fe(2) [5]Fe(2) [-1 -1 0] d [Ang] 2.583495745338251\n", "Isotropic: -54.87173513590157\n", "DMI: [ 1.55588615e+00 -3.50135692e-04 -4.46128995e-05]\n", "Symmetric-anisotropy: [-3.79737241e-01 -2.60848380e-05 -4.45024497e-05 -2.60848380e-05\n", " 2.96608300e-01 -4.97778119e-02 -4.45024497e-05 -4.97778119e-02\n", " 8.31289413e-02]\n", "J: [-5.52514724e+01 -2.60848380e-05 -4.45024497e-05 -2.60848380e-05\n", " -5.45751268e+01 -4.97778119e-02 -4.45024497e-05 -4.97778119e-02\n", " -5.47886062e+01]\n", "Energies for debugging: \n", "[[-5.44101654e-02 1.60566396e-03 -1.50610834e-03 -5.44519935e-02]\n", " [-5.51670469e-02 3.94638142e-07 -3.05633242e-07 -5.54714740e-02]\n", " [-5.46982601e-02 -1.85280615e-08 7.06977375e-08 -5.50314708e-02]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.05547147 -0.05469826 -0.05441017]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.05547147397315911 -0.05503147078078048\n", "\n", "[3]Fe(2) [5]Fe(2) [-1 0 0] d [Ang] 2.583541444641373\n", "Isotropic: -54.8543539372674\n", "DMI: [8.08531830e-01 1.36653494e+00 1.23153243e-05]\n", "Symmetric-anisotropy: [ 0.12565256 -0.14422424 0.06554219 -0.14422424 -0.2071834 -0.04991746\n", " 0.06554219 -0.04991746 0.08153084]\n", "J: [-5.47287014e+01 -1.44224241e-01 6.55421872e-02 -1.44224241e-01\n", " -5.50615373e+01 -4.99174557e-02 6.55421872e-02 -4.99174557e-02\n", " -5.47728231e+01]\n", "Energies for debugging: \n", "[[-0.0549559 0.00085845 -0.00075861 -0.05518723]\n", " [-0.05458975 -0.00143208 0.00130099 -0.05468811]\n", " [-0.05493584 0.00014424 0.00014421 -0.05476929]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.05468811 -0.05493584 -0.0549559 ]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.05468810820091341 -0.054769294557366933\n", "\n", "[4]Fe(2) [5]Fe(2) [-1 0 0] d [Ang] 2.5835398672184064\n", "Isotropic: -54.85632586904892\n", "DMI: [-8.10252027e-01 -1.36402738e+00 1.40111262e-05]\n", "Symmetric-anisotropy: [ 0.12852618 -0.14422269 -0.05676095 -0.14422269 -0.21124941 0.05978565\n", " -0.05676095 0.05978565 0.08272322]\n", "J: [-54.72779968 -0.14422269 -0.05676095 -0.14422269 -55.06757527\n", " 0.05978565 -0.05676095 0.05978565 -54.77360265]\n", "Energies for debugging: \n", "[[-0.05496711 -0.00087004 0.00075047 -0.0551985 ]\n", " [-0.05458009 0.00142079 -0.00130727 -0.05468549]\n", " [-0.05493665 0.00014424 0.00014421 -0.05477011]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.05468549 -0.05493665 -0.05496711]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.05468549220244033 -0.05477010716669793\n", "\n", "[4]Fe(2) [5]Fe(2) [-2 0 0] d [Ang] 5.951322298958084\n", "Isotropic: 0.7670124122315716\n", "DMI: [-2.53243316e-01 2.68939617e-04 -1.43892049e-04]\n", "Symmetric-anisotropy: [ 8.10949387e-02 -1.65488077e-05 -4.52585019e-05 -1.65488077e-05\n", " -3.13833695e-01 -2.50031263e-01 -4.52585019e-05 -2.50031263e-01\n", " 2.32738756e-01]\n", "J: [ 8.48107351e-01 -1.65488077e-05 -4.52585019e-05 -1.65488077e-05\n", " 4.53178717e-01 -2.50031263e-01 -4.52585019e-05 -2.50031263e-01\n", " 9.99751168e-01]\n", "Energies for debugging: \n", "[[ 1.27357111e-03 -3.21205279e-06 5.03274580e-04 6.56773632e-04]\n", " [ 7.25931228e-04 -2.23681115e-07 3.14198119e-07 7.58339032e-04]\n", " [ 2.49583803e-04 -1.27343241e-07 1.60440856e-07 9.37875670e-04]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[0.00075834 0.00024958 0.00127357]\n", "Test J_xx = E(y,z) = E(z,y)\n", "0.0007583390322003737 0.0009378756696001232\n", "\n", "[4]Fe(2) [5]Fe(2) [-3 0 0] d [Ang] 9.638732176310562\n", "Isotropic: -54.86287528016134\n", "DMI: [-8.08522805e-01 1.36388593e+00 3.14077512e-05]\n", "Symmetric-anisotropy: [ 0.12501449 0.14426595 0.05678873 0.14426595 -0.20956489 0.04990479\n", " 0.05678873 0.04990479 0.08455041]\n", "J: [-5.47378608e+01 1.44265952e-01 5.67887254e-02 1.44265952e-01\n", " -5.50724402e+01 4.99047869e-02 5.67887254e-02 4.99047869e-02\n", " -5.47783249e+01]\n", "Energies for debugging: \n", "[[-0.05496666 -0.00085843 0.00075862 -0.05519787]\n", " [-0.05458999 -0.00142067 0.0013071 -0.05469537]\n", " [-0.05494701 -0.00014423 -0.0001443 -0.05478035]]\n", "J_ii for debugging: (check if this is the same as in calculate_exchange_tensor)\n", "[-0.05469537 -0.05494701 -0.05496666]\n", "Test J_xx = E(y,z) = E(z,y)\n", "-0.05469536842615297 -0.05478035315936865\n", "\n", "================================================================================================================================================================\n", "Runtime information: \n", "Total runtime: 32.859570917000006 s\n", "----------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "Initial setup: 0.10918341700000012 s\n", "Hamiltonian conversion and XC field extraction: 0.612 s\n", "Pair and site datastructure creatrions: 0.030 s\n", "k set cration and distribution: 0.030 s\n", "Rotating XC potential: 0.278 s\n", "Greens function inversion: 31.432 s\n", "Calculate energies and magnetic components: 0.368 s\n" ] } ], "source": [ "if rank == root_node:\n", " # Calculate total charge\n", " for hamiltonian in hamiltonians:\n", " GS = hamiltonian[\"GS\"]\n", " traced = np.trace((GS), axis1=1, axis2=2)\n", " print(\"Total charge: \", int_de_ke(traced, cont.we))\n", "\n", " # iterate over the magnetic entities\n", " for tracker, mag_ent in enumerate(magnetic_entities):\n", " # iterate over the quantization axes\n", " for i, Gii in enumerate(mag_ent[\"Gii\"]):\n", " storage = []\n", " # iterate over the first and second order local perturbations\n", " for Vu1, Vu2 in zip(mag_ent[\"Vu1\"][i], mag_ent[\"Vu2\"][i]):\n", " # The Szunyogh-Lichtenstein formula\n", " traced = np.trace((Vu2 @ Gii + 0.5 * Gii @ Vu1 @ Gii), axis1=1, axis2=2)\n", " # evaluation of the contour integral\n", " storage.append(int_de_ke(traced, cont.we))\n", "\n", " # fill up the magnetic entities dictionary with the energies\n", " magnetic_entities[tracker][\"energies\"].append(storage)\n", " # convert to np array\n", " magnetic_entities[tracker][\"energies\"] = np.array(\n", " magnetic_entities[tracker][\"energies\"]\n", " )\n", " print(\"Magnetic entities integrated.\")\n", "\n", " # iterate over the pairs\n", " for tracker, pair in enumerate(pairs):\n", " # iterate over the quantization axes\n", " for i, (Gij, Gji) in enumerate(zip(pair[\"Gij\"], pair[\"Gji\"])):\n", " site_i = magnetic_entities[pair[\"ai\"]]\n", " site_j = magnetic_entities[pair[\"aj\"]]\n", "\n", " storage = []\n", " # iterate over the first order local perturbations in all possible orientations for the two sites\n", " for Vui in site_i[\"Vu1\"][i]:\n", " for Vuj in site_j[\"Vu1\"][i]:\n", " # The Szunyogh-Lichtenstein formula\n", " traced = np.trace((Vui @ Gij @ Vuj @ Gji), axis1=1, axis2=2)\n", " # evaluation of the contour integral\n", " storage.append(int_de_ke(traced, cont.we))\n", " # fill up the pairs dictionary with the energies\n", " pairs[tracker][\"energies\"].append(storage)\n", " # convert to np array\n", " pairs[tracker][\"energies\"] = np.array(pairs[tracker][\"energies\"])\n", "\n", " print(\"Pairs integrated.\")\n", "\n", " # calculate magnetic parameters\n", " for mag_ent in magnetic_entities:\n", " Kxx, Kyy, Kzz, consistency = calculate_anisotropy_tensor(mag_ent)\n", " mag_ent[\"K\"] = np.array([Kxx, Kyy, Kzz]) * sisl.unit_convert(\"eV\", \"meV\")\n", " mag_ent[\"K_consistency\"] = consistency\n", "\n", " for pair in pairs:\n", " J_iso, J_S, D, J = calculate_exchange_tensor(pair)\n", " pair[\"J_iso\"] = J_iso * sisl.unit_convert(\"eV\", \"meV\")\n", " pair[\"J_S\"] = J_S * sisl.unit_convert(\"eV\", \"meV\")\n", " pair[\"D\"] = D * sisl.unit_convert(\"eV\", \"meV\")\n", " pair[\"J\"] = J * sisl.unit_convert(\"eV\", \"meV\")\n", "\n", " print(\"Magnetic parameters calculated.\")\n", "\n", " times[\"end_time\"] = timer()\n", " print(\n", " \"##################################################################### GROGU OUTPUT #############################################################################\"\n", " )\n", "\n", " print_parameters(simulation_parameters)\n", " print_atoms_and_pairs(magnetic_entities, pairs)\n", " print_runtime_information(times)\n", "\n", " pairs, magnetic_entities = remove_clutter_for_save(pairs, magnetic_entities)\n", " # create output dictionary with all the relevant data\n", " results = dict(\n", " parameters=simulation_parameters,\n", " magnetic_entities=magnetic_entities,\n", " pairs=pairs,\n", " runtime=times,\n", " )\n", "\n", " save_pickle(simulation_parameters[\"outfile\"], results)" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "ename": "SyntaxError", "evalue": "invalid syntax (3105939143.py, line 1)", "output_type": "error", "traceback": [ "\u001b[0;36m Cell \u001b[0;32mIn[9], line 1\u001b[0;36m\u001b[0m\n\u001b[0;31m ========================================\u001b[0m\n\u001b[0m ^\u001b[0m\n\u001b[0;31mSyntaxError\u001b[0m\u001b[0;31m:\u001b[0m invalid syntax\n" ] } ], "source": [ "========================================\n", " \n", "Atom Angstrom\n", "# Label, x y z Sx Sy Sz #Q Lx Ly Lz Jx Jy Jz\n", "--------------------------------------------------------------------------------------------------------------------------------------------------------------------------\n", "Te1 1.8955 1.0943 13.1698 -0.0000 0.0000 -0.1543 # 5.9345 -0.0000 0.0000 -0.0537 -0.0000 0.0000 -0.2080 \n", "Te2 1.8955 1.0943 7.4002 0.0000 -0.0000 -0.1543 # 5.9345 0.0000 -0.0000 -0.0537 0.0000 -0.0000 -0.2080 \n", "Ge3 -0.0000 2.1887 10.2850 0.0000 0.0000 -0.1605 # 3.1927 -0.0000 0.0000 0.0012 0.0000 0.0000 -0.1593 \n", "Fe4 -0.0000 0.0000 11.6576 0.0001 -0.0001 2.0466 # 8.3044 0.0000 -0.0000 0.1606 0.0001 -0.0001 2.2072 \n", "Fe5 -0.0000 0.0000 8.9124 -0.0001 0.0001 2.0466 # 8.3044 -0.0000 0.0000 0.1606 -0.0001 0.0001 2.2072 \n", "Fe6 1.8955 1.0944 10.2850 0.0000 0.0000 1.5824 # 8.3296 -0.0000 -0.0000 0.0520 -0.0000 0.0000 1.6344 \n", "==================================================================================================================================\n", " \n", "Exchange meV\n", "--------------------------------------------------------------------------------\n", "# at1 at2 i j k # d (Ang)\n", "--------------------------------------------------------------------------------\n", "Fe4 Fe5 0 0 0 # 2.7452\n", "Isotropic -82.0854\n", "DMI 0.12557 -0.00082199 6.9668e-08\n", "Symmetric-anisotropy -0.60237 -0.83842 -0.00032278 -1.2166e-05 -3.3923e-05\n", "--------------------------------------------------------------------------------\n", "Fe4 Fe6 0 0 0 # 2.5835\n", "Isotropic -41.9627\n", "DMI 1.1205 -1.9532 0.0018386\n", "Symmetric-anisotropy 0.26007 -0.00013243 0.12977 -0.069979 -0.042066\n", "--------------------------------------------------------------------------------\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "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.9.6" } }, "nbformat": 4, "nbformat_minor": 2 }