{ "cells": [ { "cell_type": "markdown", "source": [ "# AtomsBase integration" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "[AtomsBase.jl](https://github.com/JuliaMolSim/AtomsBase.jl) is a common interface\n", "for representing atomic structures in Julia. DFTK directly supports using such\n", "structures to run a calculation as is demonstrated here." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using DFTK\n", "using AtomsBuilder" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "## Feeding an AtomsBase AbstractSystem to DFTK\n", "\n", "In this example we construct a bulk silicon system using the `bulk` function\n", "from [AtomsBuilder](https://github.com/JuliaMolSim/AtomsBuilder.jl). This function\n", "uses tabulated data to set up a reasonable starting geometry and lattice for bulk silicon." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "FlexibleSystem(Si₂, periodicity = TTT):\n cell_vectors : [ 0 2.715 2.715;\n 2.715 0 2.715;\n 2.715 2.715 0]u\"Å\"\n\n Atom(Si, [ 0, 0, 0]u\"Å\")\n Atom(Si, [ 1.3575, 1.3575, 1.3575]u\"Å\")\n" }, "metadata": {}, "execution_count": 2 } ], "cell_type": "code", "source": [ "system = bulk(:Si)" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "By default the atoms of an `AbstractSystem` employ the bare Coulomb potential.\n", "To employ pseudpotential models (which is almost always advisable for\n", "plane-wave DFT) one employs the `pseudopotential` keyword argument in\n", "model constructors such as `model_DFT`.\n", "For example we can employ a `PseudoFamily` object\n", "from the [PseudoPotentialData](https://github.com/JuliaMolSim/PseudoPotentialData.jl)\n", "package. See its documentation for more information on the available\n", "pseudopotential families and how to select them." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Model(lda_x+lda_c_pw, 3D):\n lattice (in Bohr) : [0 , 5.13061 , 5.13061 ]\n [5.13061 , 0 , 5.13061 ]\n [5.13061 , 5.13061 , 0 ]\n unit cell volume : 270.11 Bohr³\n\n atoms : Si₂\n atom potentials : ElementPsp(Si, \"/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Si.upf\")\n ElementPsp(Si, \"/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Si.upf\")\n\n num. electrons : 8\n spin polarization : none\n temperature : 0.001 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": 3 } ], "cell_type": "code", "source": [ "using PseudoPotentialData # defines PseudoFamily\n", "\n", "pd_lda_family = PseudoFamily(\"dojo.nc.sr.lda.v0_4_1.standard.upf\")\n", "model = model_DFT(system; functionals=LDA(), temperature=1e-3,\n", " pseudopotentials=pd_lda_family)" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "Alternatively the `pseudopotentials` object also accepts a `Dict{Symbol,String}`,\n", "which provides for each element symbol the filename or identifier\n", "of the pseudopotential to be employed, e.g." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "Model(lda_x+lda_c_pw, 3D):\n lattice (in Bohr) : [0 , 5.13061 , 5.13061 ]\n [5.13061 , 0 , 5.13061 ]\n [5.13061 , 5.13061 , 0 ]\n unit cell volume : 270.11 Bohr³\n\n atoms : Si₂\n atom potentials : ElementPsp(Si, \"/home/runner/.julia/artifacts/966fd9cdcd7dbaba6dc2bf43ee50dd81e63e8837/Si.gth\")\n ElementPsp(Si, \"/home/runner/.julia/artifacts/966fd9cdcd7dbaba6dc2bf43ee50dd81e63e8837/Si.gth\")\n\n num. electrons : 8\n spin polarization : none\n temperature : 0.001 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": 4 } ], "cell_type": "code", "source": [ "path_to_pspfile = PseudoFamily(\"cp2k.nc.sr.lda.v0_1.semicore.gth\")[:Si]\n", "model = model_DFT(system; functionals=LDA(), temperature=1e-3,\n", " pseudopotentials=Dict(:Si => path_to_pspfile))" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "We can then discretise such a model and solve:" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) Diag Δtime\n", "--- --------------- --------- --------- ---- ------\n", " 1 -7.921681503863 -0.69 5.5 839ms\n", " 2 -7.926133575373 -2.35 -1.22 1.0 770ms\n", " 3 -7.926833524972 -3.15 -2.37 2.0 152ms\n", " 4 -7.926861274900 -4.56 -2.99 2.9 199ms\n", " 5 -7.926861645407 -6.43 -3.38 2.2 162ms\n", " 6 -7.926861670500 -7.60 -3.84 1.4 141ms\n", " 7 -7.926861679389 -8.05 -4.17 1.8 146ms\n", " 8 -7.926861681752 -8.63 -4.90 1.6 148ms\n", " 9 -7.926861681834 -10.08 -5.00 2.2 160ms\n", " 10 -7.926861681869 -10.46 -5.69 1.1 137ms\n", " 11 -7.926861681872 -11.48 -6.12 2.2 162ms\n", " 12 -7.926861681873 -12.59 -7.10 1.8 152ms\n", " 13 -7.926861681873 -13.52 -7.41 2.9 202ms\n", " 14 -7.926861681873 -14.75 -7.75 1.1 141ms\n", " 15 -7.926861681873 + -Inf -8.16 1.6 196ms\n" ] } ], "cell_type": "code", "source": [ "basis = PlaneWaveBasis(model; Ecut=15, kgrid=[4, 4, 4])\n", "scfres = self_consistent_field(basis, tol=1e-8);" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "If we did not want to use\n", "[AtomsBuilder](https://github.com/JuliaMolSim/AtomsBuilder.jl)\n", "we could of course use any other package\n", "which yields an AbstractSystem object. This includes:" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "### Reading a system using AtomsIO\n", "\n", "Read a file using [AtomsIO](https://github.com/mfherbst/AtomsIO.jl),\n", "which directly yields an AbstractSystem." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using AtomsIO\n", "system = load_system(\"Si.extxyz\");" ], "metadata": {}, "execution_count": 6 }, { "cell_type": "markdown", "source": [ "Run the LDA calculation:" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) Diag Δtime\n", "--- --------------- --------- --------- ---- ------\n", " 1 -7.921730371012 -0.69 5.8 187ms\n", " 2 -7.926131411206 -2.36 -1.22 1.0 133ms\n", " 3 -7.926833365786 -3.15 -2.37 2.1 154ms\n", " 4 -7.926861259195 -4.55 -3.00 3.1 241ms\n", " 5 -7.926861647203 -6.41 -3.39 2.1 159ms\n", " 6 -7.926861671123 -7.62 -3.84 1.4 139ms\n", " 7 -7.926861679640 -8.07 -4.18 1.8 172ms\n", " 8 -7.926861681772 -8.67 -4.90 1.4 140ms\n", " 9 -7.926861681856 -10.08 -5.18 1.8 150ms\n", " 10 -7.926861681871 -10.82 -5.84 1.8 148ms\n", " 11 -7.926861681873 -11.74 -6.41 2.0 155ms\n", " 12 -7.926861681873 -13.00 -7.52 1.8 155ms\n", " 13 -7.926861681873 -14.15 -8.00 3.5 197ms\n", " 14 -7.926861681873 + -14.75 -8.78 2.0 153ms\n" ] } ], "cell_type": "code", "source": [ "pseudopotentials = PseudoFamily(\"cp2k.nc.sr.lda.v0_1.semicore.gth\")\n", "model = model_DFT(system; pseudopotentials, functionals=LDA(), temperature=1e-3)\n", "basis = PlaneWaveBasis(model; Ecut=15, kgrid=[4, 4, 4])\n", "scfres = self_consistent_field(basis, tol=1e-8);" ], "metadata": {}, "execution_count": 7 }, { "cell_type": "markdown", "source": [ "The same could be achieved using [ExtXYZ](https://github.com/libAtoms/ExtXYZ.jl)\n", "by `system = Atoms(read_frame(\"Si.extxyz\"))`,\n", "since the `ExtXYZ.Atoms` object is directly AtomsBase-compatible." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "### Directly setting up a system in AtomsBase" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "n Energy log10(ΔE) log10(Δρ) Diag Δtime\n", "--- --------------- --------- --------- ---- ------\n", " 1 -7.921728831917 -0.69 5.6 218ms\n", " 2 -7.926138025114 -2.36 -1.22 1.0 143ms\n", " 3 -7.926836807033 -3.16 -2.37 2.0 147ms\n", " 4 -7.926864499287 -4.56 -2.99 3.0 195ms\n", " 5 -7.926865060810 -6.25 -3.37 2.1 156ms\n", " 6 -7.926865082027 -7.67 -3.84 1.2 128ms\n", " 7 -7.926865089458 -8.13 -4.06 1.9 141ms\n" ] } ], "cell_type": "code", "source": [ "using AtomsBase\n", "using Unitful\n", "using UnitfulAtomic\n", "\n", "# Construct a system in the AtomsBase world\n", "a = 10.26u\"bohr\" # Silicon lattice constant\n", "lattice = a / 2 * [[0, 1, 1.], # Lattice as vector of vectors\n", " [1, 0, 1.],\n", " [1, 1, 0.]]\n", "atoms = [:Si => ones(3)/8, :Si => -ones(3)/8]\n", "system = periodic_system(atoms, lattice; fractional=true)\n", "\n", "# Now run the LDA calculation:\n", "pseudopotentials = PseudoFamily(\"cp2k.nc.sr.lda.v0_1.semicore.gth\")\n", "model = model_DFT(system; pseudopotentials, functionals=LDA(), temperature=1e-3)\n", "basis = PlaneWaveBasis(model; Ecut=15, kgrid=[4, 4, 4])\n", "scfres = self_consistent_field(basis, tol=1e-4);" ], "metadata": {}, "execution_count": 8 }, { "cell_type": "markdown", "source": [ "## Obtaining an AbstractSystem from DFTK data" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "At any point we can also get back the DFTK model as an\n", "AtomsBase-compatible `AbstractSystem`:" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "FlexibleSystem(Si₂, periodicity = TTT):\n cell_vectors : [ 0 5.13 5.13;\n 5.13 0 5.13;\n 5.13 5.13 0]u\"a₀\"\n\n Atom(Si, [ 1.2825, 1.2825, 1.2825]u\"a₀\")\n Atom(Si, [ -1.2825, -1.2825, -1.2825]u\"a₀\")\n" }, "metadata": {}, "execution_count": 9 } ], "cell_type": "code", "source": [ "second_system = atomic_system(model)" ], "metadata": {}, "execution_count": 9 }, { "cell_type": "markdown", "source": [ "Similarly DFTK offers a method to the `atomic_system` and `periodic_system` functions\n", "(from AtomsBase), which enable a seamless conversion of the usual data structures for\n", "setting up DFTK calculations into an `AbstractSystem`:" ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "FlexibleSystem(Si₂, periodicity = TTT):\n cell_vectors : [ 0 5.13155 5.13155;\n 5.13155 0 5.13155;\n 5.13155 5.13155 0]u\"a₀\"\n\n Atom(Si, [ 1.28289, 1.28289, 1.28289]u\"a₀\")\n Atom(Si, [-1.28289, -1.28289, -1.28289]u\"a₀\")\n" }, "metadata": {}, "execution_count": 10 } ], "cell_type": "code", "source": [ "lattice = 5.431u\"Å\" / 2 * [[0 1 1.];\n", " [1 0 1.];\n", " [1 1 0.]];\n", "Si = ElementPsp(:Si, pseudopotentials)\n", "atoms = [Si, Si]\n", "positions = [ones(3)/8, -ones(3)/8]\n", "\n", "third_system = atomic_system(lattice, atoms, positions)" ], "metadata": {}, "execution_count": 10 } ], "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 }