{ "cells": [ { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n# Flow for phonon-limited mobilities in semiconductors\n\nThis flow computes the phonon-limited mobility in AlAs.\nusing different dense k/q meshes\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import sys\nimport os\nimport abipy.data as abidata\nimport abipy.abilab as abilab\nimport abipy.flowtk as flowtk\nimport abipy.core.abinit_units as abu\n\n\ndef build_flow(options):\n\n # Working directory (default is the name of the script with '.py' removed and \"run_\" replaced by \"flow_\")\n if not options.workdir:\n options.workdir = os.path.basename(sys.argv[0]).replace(\".py\", \"\").replace(\"run_\", \"flow_\")\n\n # Initialize the flow\n flow = flowtk.Flow(workdir=options.workdir, manager=options.manager)\n\n # Initialize the tmesh, sigma_kerange, sigma_erange and dense meshes for the eph integrations\n tmesh = [300, 300, 1]\n sigma_kerange = [0, 0.5 * abu.eV_Ha] # 0.5 eV range for the WFK of electrons\n sigma_erange = [0, 0.25 * abu.eV_Ha] # 0.25 eV range for the EPH computation for electrons\n dense_meshes = [[30, 30, 30],\n [40, 40, 40]]\n\n # Initialize the structure and pseudos\n structure = abidata.structure_from_ucell(\"AlAs\")\n pseudos = abidata.pseudos(\"13al.981214.fhi\", \"33as.pspnc\")\n\n # Ground-state computation for 1) the phonons and 2) the WFK generation\n scf_input = abilab.AbinitInput(structure, pseudos=pseudos)\n\n scf_input.set_vars(\n nband=8,\n ecut=2.0,\n ngkpt=[4, 4, 4],\n shiftk=[0, 0, 0],\n tolvrs=1.0e-10,\n diemac=9.0,\n prtden=1,\n #iomode=3,\n )\n\n # Initialize the first work and add the ground-state task\n work0 = flowtk.Work()\n work0.register_scf_task(scf_input)\n\n # Band structure calculation to make sure everything is OK\n # Also allows to compare the results obtained with abitk to\n # check the SKW interpolation works as needed\n bs_input = scf_input.make_ebands_input(tolwfr=1e-12, ndivsm=10, nb_extra=5)\n bs_input.set_vars(nstep=100, nbdbuf=1)\n work0.register_nscf_task(bs_input, deps={work0[0]: \"DEN\"})\n\n # NSCF input for the WFK needed to interpolate with kerange\n nscf_input = abilab.AbinitInput(structure, pseudos)\n nscf_input.set_vars(\n ecut=2,\n nband=8,\n iscf=-2,\n tolwfr=1e-20,\n prtwf=1,\n ngkpt=[16, 16, 16], # Should be dense enough so that the kerange interpolation works\n shiftk=[0.0, 0.0, 0.0],\n )\n\n work0.register_nscf_task(nscf_input, deps={work0[0]: \"DEN\"})\n flow.register_work(work0)\n\n # Add the phonon work to the flow\n ddb_ngqpt = [4, 4, 4]\n ph_work = flowtk.PhononWork.from_scf_task(work0[0], qpoints=ddb_ngqpt,\n is_ngqpt=True, with_becs=True, with_quad=False)\n flow.register_work(ph_work)\n\n # We loop over the dense meshes\n for i, sigma_ngkpt in enumerate(dense_meshes):\n # Use the kerange trick to generate a WFK file\n multi = nscf_input.make_wfk_kerange_inputs(sigma_kerange=sigma_kerange,\n sigma_ngkpt=sigma_ngkpt)\n kerange_input, wfk_input = multi.split_datasets()\n\n work_eph = flowtk.Work()\n work_eph.register_kerange_task(kerange_input, deps={work0[2]: \"WFK\"})\n work_eph.register_nscf_task(wfk_input,\n deps={work0[0]: \"DEN\", work_eph[0]: \"KERANGE.nc\"})\n\n # Generate the input file for the transport calculation.\n # Use ibte_prep = 1 to activate the iterative BTE.\n eph_input = wfk_input.make_eph_transport_input(ddb_ngqpt=ddb_ngqpt,\n sigma_erange=sigma_erange,\n tmesh=tmesh,\n eph_ngqpt_fine=sigma_ngkpt,\n ibte_prep=1)\n\n # We compute the phonon dispersion in the EPH code to be able to check they are ok.\n if i == 0:\n eph_input.set_qpath(20)\n\n work_eph.register_eph_task(eph_input,\n deps={work_eph[1]: \"WFK\", ph_work: [\"DDB\", \"DVDB\"]})\n\n flow.register_work(work_eph)\n\n flow.allocate(use_smartio=True)\n\n return flow\n\n\n# This block generates the thumbnails in the AbiPy gallery.\n# You can safely REMOVE this part if you are using this script for production runs.\nif os.getenv(\"READTHEDOCS\", False):\n __name__ = None\n import tempfile\n options = flowtk.build_flow_main_parser().parse_args([\"-w\", tempfile.mkdtemp()])\n build_flow(options).graphviz_imshow()\n\n\n@flowtk.flow_main\ndef main(options):\n \"\"\"\n This is our main function that will be invoked by the script.\n flow_main is a decorator implementing the command line interface.\n Command line args are stored in `options`.\n \"\"\"\n\n return build_flow(options)\n\n\nif __name__ == \"__main__\":\n sys.exit(main())" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "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.8.0" } }, "nbformat": 4, "nbformat_minor": 0 }