{ "cells": [ { "cell_type": "markdown", "source": [ "# Input and output formats" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "This section provides an overview of the input and output formats\n", "supported by DFTK, usually via integration with a third-party library.\n", "\n", "## Reading / writing files supported by AtomsIO\n", "[AtomsIO](https://github.com/mfherbst/AtomsIO.jl) is a Julia package which supports\n", "reading / writing atomistic structures from / to a large range of file formats.\n", "Supported formats include Crystallographic Information Framework (CIF),\n", "XYZ and extxyz files, ASE / Gromacs / LAMMPS / Amber trajectory files\n", "or input files of various other codes (e.g. Quantum Espresso, VASP, ABINIT, CASTEP, …).\n", "The full list of formats is is available in the\n", "[AtomsIO documentation](https://mfherbst.github.io/AtomsIO.jl/stable).\n", "\n", "The AtomsIO functionality is split into two packages. The main package, `AtomsIO` itself,\n", "only depends on packages, which are registered in the Julia General package registry.\n", "In contrast `AtomsIOPython` extends `AtomsIO` by parsers depending on python packages,\n", "which are automatically managed via `PythonCall`. While it thus provides the full set of\n", "supported IO formats, this also adds additional practical complications, so some users\n", "may choose not to use `AtomsIOPython`.\n", "\n", "As an example we start the calculation of a simple antiferromagnetic iron crystal\n", "using a Quantum-Espresso input file, [Fe_afm.pwi](Fe_afm.pwi).\n", "For more details about calculations on magnetic systems\n", "using collinear spin, see Collinear spin and magnetic systems.\n", "\n", "First we parse the Quantum Espresso input file using AtomsIO,\n", "which reads the lattice, atomic positions and initial magnetisation\n", "from the input file and returns it as an\n", "[AtomsBase](https://github.com/JuliaMolSim/AtomsBase.jl) `AbstractSystem`,\n", "the JuliaMolSim community standard for representing atomic systems." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "FlexibleSystem(Fe₂, periodicity = TTT):\n cell_vectors : [ 2.86814 0 0;\n 0 2.86814 0;\n 0 0 2.86814]u\"Å\"\n\n Atom(Fe, [ 0, 0, 0]u\"Å\")\n Atom(Fe, [ 1.43407, 1.43407, 1.43407]u\"Å\")\n" }, "metadata": {}, "execution_count": 1 } ], "cell_type": "code", "source": [ "using AtomsIO # Use Julia-only IO parsers\n", "using AtomsIOPython # Use python-based IO parsers (e.g. ASE)\n", "system = load_system(\"Fe_afm.pwi\")" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "Next we build a model making use of the AtomsBase integration\n", "of DFTK and supplying a pseudpotential family\n", "from the [PseudoPotentialData](https://github.com/JuliaMolSim/PseudoPotentialData.jl)" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Model(lda_x+lda_c_pw, 3D):\n lattice (in Bohr) : [5.42 , 0 , 0 ]\n [0 , 5.42 , 0 ]\n [0 , 0 , 5.42 ]\n unit cell volume : 159.22 Bohr³\n\n atoms : Fe₂\n atom potentials : ElementPsp(Fe, \"/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Fe.upf\")\n ElementPsp(Fe, \"/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Fe.upf\")\n\n num. electrons : 32\n spin polarization : collinear\n temperature : 0.01 Ha\n smearing : DFTK.Smearing.FermiDirac()\n\n terms : Kinetic()\n AtomicLocal()\n AtomicNonlocal()\n Ewald(nothing)\n PspCorrection()\n Hartree()\n Xc(lda_x, lda_c_pw)\n Entropy()" }, "metadata": {}, "execution_count": 2 } ], "cell_type": "code", "source": [ "using DFTK\n", "using PseudoPotentialData\n", "pd_lda_family = PseudoFamily(\"dojo.nc.sr.lda.v0_4_1.standard.upf\")\n", "model = model_DFT(system;\n", " functionals=LDA(),\n", " temperature=0.01,\n", " pseudopotentials=pd_lda_family)" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "Finally we run the calculation:" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "┌ Warning: Negative ρcore detected: -0.010160696831625975\n", "└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/terms/xc.jl:39\n", "n Energy log10(ΔE) log10(Δρ) Magnet Diag Δtime\n", "--- --------------- --------- --------- ------ ---- ------\n", " 1 -230.6877319769 0.31 -0.519 5.2 166ms\n", " 2 -234.1255799566 0.54 -0.32 -0.461 2.5 192ms\n", " 3 -234.2151242682 -1.05 -1.06 -0.449 3.4 162ms\n", " 4 -234.2161692706 -2.98 -1.94 -0.452 1.2 130ms\n", " 5 -234.2162035304 -4.47 -2.79 -0.455 1.9 111ms\n", " 6 -234.2162059315 -5.62 -2.96 -0.456 3.8 138ms\n", " 7 -234.2162066780 -6.13 -3.09 -0.458 1.2 112ms\n", " 8 -234.2162069457 -6.57 -3.18 -0.459 1.0 112ms\n", " 9 -234.2162070616 -6.94 -3.25 -0.460 1.0 107ms\n", " 10 -234.2162074006 -6.47 -3.81 -0.466 1.4 114ms\n", " 11 -234.2162074505 -7.30 -3.78 -0.470 1.8 122ms\n", " 12 -234.2162075074 -7.24 -3.78 -0.475 1.0 102ms\n", " 13 -234.2162077177 -6.68 -3.86 -0.496 2.1 132ms\n", " 14 -234.2162079115 -6.71 -4.13 -0.527 2.9 150ms\n", " 15 -234.2162079637 -7.28 -4.59 -0.550 3.0 155ms\n", " 16 -234.2162079648 -8.97 -4.98 -0.550 1.6 137ms\n", " 17 -234.2162079651 -9.53 -5.41 -0.552 2.6 144ms\n", " 18 -234.2162079652 -10.21 -5.66 -0.552 2.0 129ms\n", " 19 -234.2162079652 -10.56 -5.76 -0.552 2.5 142ms\n", " 20 -234.2162079652 -10.55 -5.80 -0.552 1.9 121ms\n", " 21 -234.2162079653 -10.48 -5.85 -0.552 1.0 107ms\n", " 22 -234.2162079653 -10.54 -5.88 -0.552 1.0 106ms\n", " 23 -234.2162079654 -10.20 -5.99 -0.552 1.8 112ms\n", " 24 -234.2162079654 -10.54 -6.06 -0.552 1.2 129ms\n" ] } ], "cell_type": "code", "source": [ "basis = PlaneWaveBasis(model; Ecut=10, kgrid=(2, 2, 2))\n", "ρ0 = guess_density(basis, system)\n", "scfres = self_consistent_field(basis, ρ=ρ0);" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "> **DFTK data formats are not yet fully matured**\n", ">\n", "> The data format in which DFTK saves data as well as the general interface\n", "> of the `load_scfres` and `save_scfres` pair of functions\n", "> are not yet fully matured. They have become fairly stable, but we reserve\n", "> the right to change them incompatibly between DFTK versions (including\n", "> patch versions) at this point." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Writing VTK files for visualization\n", "For visualizing the density or the Kohn-Sham orbitals DFTK supports storing\n", "the result of an SCF calculations in the form of VTK files.\n", "These can afterwards be visualized using tools such\n", "as [paraview](https://www.paraview.org/).\n", "Using this feature requires\n", "the [WriteVTK.jl](https://github.com/jipolanco/WriteVTK.jl/) Julia package." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using WriteVTK\n", "save_scfres(\"iron_afm.vts\", scfres; save_ψ=true);" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "This will save the iron calculation above into the file `iron_afm.vts`,\n", "using `save_ψ=true` to also include the KS orbitals." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Parsable data-export using json\n", "Many structures in DFTK support the (unexported) `todict` function,\n", "which returns a simplified dictionary representation of the data." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Dict{String, Float64} with 9 entries:\n \"AtomicNonlocal\" => 0.657977\n \"PspCorrection\" => 6.08845\n \"Ewald\" => -171.89\n \"total\" => -234.216\n \"Entropy\" => -0.0383957\n \"Kinetic\" => 79.7543\n \"AtomicLocal\" => -154.423\n \"Hartree\" => 36.1757\n \"Xc\" => -30.5411" }, "metadata": {}, "execution_count": 5 } ], "cell_type": "code", "source": [ "DFTK.todict(scfres.energies)" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "This in turn can be easily written to disk using a JSON library.\n", "Currently we integrate most closely with `JSON3`,\n", "which is thus recommended." ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "{\n", " \"AtomicNonlocal\": 0.6579767812496529,\n", " \"PspCorrection\": 6.088454244905408,\n", " \"Ewald\": -171.88999401525265,\n", " \"total\": -234.2162079653813,\n", " \"Entropy\": -0.0383957450517167,\n", " \"Kinetic\": 79.7543109320243,\n", " \"AtomicLocal\": -154.4232165863035,\n", " \"Hartree\": 36.17571837639441,\n", " \"Xc\": -30.541061953347185\n", "}\n" ] } ], "cell_type": "code", "source": [ "using JSON3\n", "open(\"iron_afm_energies.json\", \"w\") do io\n", " JSON3.pretty(io, DFTK.todict(scfres.energies))\n", "end\n", "println(read(\"iron_afm_energies.json\", String))" ], "metadata": {}, "execution_count": 6 }, { "cell_type": "markdown", "source": [ "Once JSON3 is loaded, additionally a convenience function for saving\n", "a summary of `scfres` objects using `save_scfres` is available:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using JSON3\n", "save_scfres(\"iron_afm.json\", scfres)" ], "metadata": {}, "execution_count": 7 }, { "cell_type": "markdown", "source": [ "Similarly a summary of the band data (occupations, eigenvalues, εF, etc.)\n", "for post-processing can be dumped using `save_bands`:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "save_bands(\"iron_afm_scfres.json\", scfres)" ], "metadata": {}, "execution_count": 8 }, { "cell_type": "markdown", "source": [ "Notably this function works both for the results obtained\n", "by `self_consistent_field` as well as `compute_bands`:" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "┌ Warning: The provided cell is a supercell: the returned k-path is the standard k-path of the associated primitive cell in the basis of the supercell reciprocal lattice.\n", "│ cell =\n", "│ Spglib.SpglibCell{Float64, Float64, Int64, Float64}\n", "│ lattice:\n", "│ 5.41999997388764 0.0 0.0\n", "│ 0.0 5.41999997388764 0.0\n", "│ 0.0 0.0 5.41999997388764\n", "│ 2 atomic positions:\n", "│ 0.0 0.0 0.0\n", "│ 0.5 0.5 0.5\n", "│ 2 atoms:\n", "│ 1 1\n", "│ 2 magmoms:\n", "│ 0.0 0.0\n", "│ \n", "└ @ BrillouinSpglibExt ~/.julia/packages/Brillouin/FwrEN/ext/BrillouinSpglibExt.jl:28\n", "┌ Warning: Negative ρcore detected: -0.010160696831625975\n", "└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/terms/xc.jl:39\n" ] } ], "cell_type": "code", "source": [ "bands = compute_bands(scfres, kline_density=10)\n", "save_bands(\"iron_afm_bands.json\", bands)" ], "metadata": {}, "execution_count": 9 }, { "cell_type": "markdown", "source": [ "## Writing and reading JLD2 files\n", "The full state of a DFTK self-consistent field calculation can be\n", "stored on disk in form of an [JLD2.jl](https://github.com/JuliaIO/JLD2.jl) file.\n", "This file can be read from other Julia scripts\n", "as well as other external codes supporting the HDF5 file format\n", "(since the JLD2 format is based on HDF5). This includes notably `h5py`\n", "to read DFTK output from python." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using JLD2\n", "save_scfres(\"iron_afm.jld2\", scfres)" ], "metadata": {}, "execution_count": 10 }, { "cell_type": "markdown", "source": [ "Saving such JLD2 files supports some options, such as `save_ψ=false`, which avoids saving\n", "the Bloch waves (much faster and smaller files). Notice that JLD2 files can also be used\n", "with `save_bands`." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "Since such JLD2 can also be read by DFTK to start or continue a calculation,\n", "these can also be used for checkpointing or for transferring results\n", "to a different computer.\n", "See Saving SCF results on disk and SCF checkpoints for details." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "(Cleanup files generated by this notebook.)" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "rm.([\"iron_afm.vts\", \"iron_afm.jld2\",\n", " \"iron_afm.json\", \"iron_afm_energies.json\", \"iron_afm_scfres.json\",\n", " \"iron_afm_bands.json\"]);" ], "metadata": {}, "execution_count": 11 } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.11.4" }, "kernelspec": { "name": "julia-1.11", "display_name": "Julia 1.11.4", "language": "julia" } }, "nbformat": 4 }