{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n\n# Surface adsorption study using the ASE database\n\nIn this tutorial we will adsorb C, N and O on 7 different FCC(111) surfaces\nwith 1, 2 and 3 layers and we will use database files to store the results.\n\n.. seealso::\n\n The :mod:`ase.db` module documentation.\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. centered:: |Cu1O| |Cu2O| |Cu3O|\n\n.. |Cu1O| image:: ../../../examples/tutorials/Cu1O.png\n.. |Cu2O| image:: ../../../examples/tutorials/Cu2O.png\n.. |Cu3O| image:: ../../../examples/tutorials/Cu3O.png\n\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Bulk\n\nFirst, we calculate the equilibrium bulk FCC lattice constants for the seven\nelements where the :mod:`EMT ` potential offers a\ncomputational affordable calculator. We store the resulting configurations\ntogether with its bulk modulus data in the database ``bulk.db``:\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from pathlib import Path\n\nfrom ase import Atoms\nfrom ase.build import add_adsorbate, bulk, fcc111\nfrom ase.calculators.emt import EMT\nfrom ase.constraints import FixAtoms\nfrom ase.db import connect\nfrom ase.eos import calculate_eos\nfrom ase.io import write\nfrom ase.optimize import BFGS\n\nbulk_syms = ['Al', 'Ni', 'Cu', 'Pd', 'Ag', 'Pt', 'Au']\npath_to_bulk_db = Path('bulk.db')\n\n# clean up if there is an old file\npath_to_bulk_db.unlink(missing_ok=True)\n# use context manage to ensure that connection is closed afterwards\nwith connect(path_to_bulk_db) as bulk_db:\n for symb in bulk_syms:\n atoms = bulk(symb, 'fcc')\n\n atoms.calc = EMT()\n eos = calculate_eos(atoms)\n v, e, B = eos.fit() # find minimum\n # do one more calculation at the minimum and write to database:\n atoms.cell *= (v / atoms.get_volume()) ** (1 / 3)\n atoms.get_potential_energy()\n bulk_db.write(atoms, bm=B)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. highlight:: bash\n\nHere is how we can inspect the results of the previous step::\n\n $ ase db bulk.db -c +bm # show also the bulk-modulus column\n id|age|formula|calculator|energy| fmax|pbc|volume|charge| mass| bm\n 1|10s|Al |emt |-0.005|0.000|TTT|15.932| 0.000| 26.982|0.249\n 2| 9s|Ni |emt |-0.013|0.000|TTT|10.601| 0.000| 58.693|1.105\n 3| 9s|Cu |emt |-0.007|0.000|TTT|11.565| 0.000| 63.546|0.839\n 4| 9s|Pd |emt |-0.000|0.000|TTT|14.588| 0.000|106.420|1.118\n 5| 9s|Ag |emt |-0.000|0.000|TTT|16.775| 0.000|107.868|0.625\n 6| 9s|Pt |emt |-0.000|0.000|TTT|15.080| 0.000|195.084|1.736\n 7| 9s|Au |emt |-0.000|0.000|TTT|16.684| 0.000|196.967|1.085\n Rows: 7\n Keys: bm\n $ ase gui bulk.db\n\nThe :file:`bulk.db` is an SQLite3_ database in a single file::\n\n $ file bulk.db\n bulk.db: SQLite 3.x database\n\n\nIf you want to see what's inside you can convert the database file to a json\nfile and open that in your text editor::\n\n $ ase db bulk.db --insert-into bulk.json\n Added 0 key-value pairs (0 pairs updated)\n Inserted 7 rows\n\nor, you can look at a single row like this::\n\n $ ase db bulk.db Cu -j\n {\"1\": {\n \"calculator\": \"emt\",\n \"energy\": -0.007036492048371201,\n \"forces\": [[0.0, 0.0, 0.0]],\n \"key_value_pairs\": {\"bm\": 0.8392875566787444},\n ...\n ...\n }\n\nThe json file format is human readable, but much less efficient to work with\ncompared to a SQLite3 file.\n\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Adsorbates\n\nNow we do the adsorption calculations and store the results\nin the ``ads.db`` database:\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "path_to_ads_db = Path('ads.db')\nads_syms = ['C', 'N', 'O']\nnlayers = [1, 2, 3]\n\n\ndef optimize_adsorbate(symb, alat, nlayer, ads):\n atoms = fcc111(symb, (1, 1, nlayer), a=alat)\n add_adsorbate(atoms, ads, height=1.0, position='fcc')\n\n # Constrain all atoms except the adsorbate:\n fixed = list(range(len(atoms) - 1))\n atoms.constraints = [FixAtoms(indices=fixed)]\n\n atoms.calc = EMT()\n opt = BFGS(atoms, logfile=None)\n opt.run(fmax=0.01)\n return atoms\n\n\n# clean up if there is an old file\npath_to_ads_db.unlink(missing_ok=True)\n# again use context-manager\nwith connect(path_to_ads_db) as ads_db:\n for row in bulk_db.select():\n alat = row.cell[0, 1] * 2 # lattice constant\n symb = row.symbols[0]\n for nlayer in nlayers:\n for ads in ads_syms:\n atoms = optimize_adsorbate(symb, alat, nlayer, ads)\n myid = ads_db.reserve(layers=nlayer, surf=symb, ads=ads)\n if myid is not None:\n ads_db.write(\n atoms, id=myid, layers=nlayer, surf=symb, ads=ads\n )" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We now have a new database file with 63 rows::\n\n $ ase db ads.db -n\n 63 rows\n\nThese 63 calculations only take a few seconds with EMT. Suppose you want to\nuse DFT and send the calculations to a supercomputer. In that case you may\nwant to run several calculations in different jobs on the computer. In\naddition, some of the jobs could time out and not finish. It's a good idea\nto modify the script a bit for this scenario. Therefore we added a couple of\nlines to the inner loop.\n\nThe :meth:`~ase.db.core.Database.reserve` method will check if there is a row\nwith the keys ``layers=n``, ``surf=symb`` and ``ads=ads``. If there is, then\nthe calculation will be skipped. If there is not, then an empty row with\nthose keys-values will be written and the calculation will start. When done,\nthe real row will be written into to the reserved row with ``id=myid``.\nThis modified script can run in several jobs all running in parallel and no\ncalculation will be done twice.\n\n.. highlight:: bash\n\nIn case a calculation crashes, you will see empty rows left in the database::\n\n $ ase db ads.db natoms=0 -c ++\n id|age|user |formula|pbc|charge| mass|ads|layers|surf\n 17|31s|jensj| |FFF| 0.000|0.000| N| 1| Cu\n Rows: 1\n Keys: ads, layers, surf\n\nDelete them, fix the problem and run the script again::\n\n $ ase db ads.db natoms=0 --delete\n Delete 1 row? (yes/No): yes\n Deleted 1 row\n $ python ads.py # or sbatch ...\n $ ase db ads.db natoms=0\n Rows: 0\n\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Reference energies\n\nLet's also calculate the energy of the clean surfaces and the isolated\nadsorbates and store them in the ``refs.db`` database:\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "path_to_refs_db = Path('refs.db')\n\n\ndef clean_surface_reference(symb, alat, nlayer):\n atoms = fcc111(symb, (1, 1, nlayer), a=alat)\n atoms.calc = EMT()\n atoms.get_forces()\n return atoms\n\n\n# clean slabs\npath_to_refs_db.unlink(missing_ok=True)\n# connect to refs_db for writing\nwith connect(path_to_refs_db) as refs_db:\n # connect bulk_db for reading\n with connect(path_to_bulk_db) as bulk_db:\n for row in bulk_db.select():\n alat = row.cell[0, 1] * 2 # lattice constant\n symb = row.symbols[0]\n for nlayer in nlayers:\n myid = refs_db.reserve(layers=nlayer, surf=symb, ads='clean')\n if myid is not None:\n atoms = clean_surface_reference(symb, alat, nlayer)\n refs_db.write(\n atoms, id=myid, layers=nlayer, surf=symb, ads='clean'\n )\n\n # isolated atoms:\n for ads in ads_syms:\n atoms = Atoms(ads)\n atoms.calc = EMT()\n atoms.get_potential_energy()\n myid = refs_db.reserve(ads=ads)\n if myid is not None:\n refs_db.write(atoms, id=myid, ads=ads)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The previous code snippet saves those 24 reference energies (clean surfaces\nand isolated adsorbates) in a :file:`refs.db` file.\nSuppose we want to extract the clean surfaces from\nthe ``refs.db`` and store them in ``clean_surfaces.db``.\nWe could change the script and run the calculations again,\nbut we can also manipulate the files using the ``ase db`` tool. First, we\nmove over the clean surfaces::\n\n $ ase db refs.db ads=clean --insert-into clean_surfaces.db\n Added 0 key-value pairs (0 pairs updated)\n Inserted 21 rows\n\nIf we want we can delete those from ``refs.db``::\n\n $ ase db refs.db ads=clean --delete --yes\n Deleted 21 rows\n\nSince this is only to demonstrate usage of the command line we will keep them.\n\nNow we might want to extract the three atoms (``pbc=FFF``, no periodicity)::\n\n $ ase db refs.db pbc=FFF --insert-into isolated_atoms.db\n Added 0 key-value pairs (0 pairs updated)\n Inserted 3 rows\n\nFinally, as we have not deleted anything these are the databases\nwhich we should have build so far::\n\n $ ase db ads.db -n\n 63 rows\n $ ase db refs.db -n\n 24 rows\n $ ase db clean_surfaces.db -n\n 21 rows\n $ ase db isolated_atoms.db -n\n 3 rows\n\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Analysis\n\nNow we have what we need to calculate the adsorption energies and heights:\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# connect to ads_db for updating\nwith connect(path_to_ads_db) as ads_db:\n # connect to refs_db for reading\n with connect(path_to_refs_db) as refs_db:\n for row in ads_db.select():\n # atoms\n e_ads = refs_db.get(formula=row.ads).energy\n # clean surface\n e_clean = refs_db.get(\n surf=row.surf, layers=row.layers, ads='clean'\n ).energy\n ea = row.energy - e_ads - e_clean\n h = row.positions[-1, 2] - row.positions[-2, 2]\n ads_db.update(row.id, ea=ea, height=h)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Once run we can inspect the results for three layers of Pt::\n\n $ ase db ads.db Pt,layers=3 -c formula,ea,height\n formula| ea|height\n Pt3C |-3.715| 1.504\n Pt3N |-5.419| 1.534\n Pt3O |-4.724| 1.706\n Rows: 3\n Keys: ads, ea, height, layers, surf\n\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Finally, we can load specific structures of the database and create figures:\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "with connect(path_to_ads_db) as ads_db:\n for nlayer in nlayers:\n atoms = ads_db.get(surf='Cu', ads='O', layers=nlayer).toatoms()\n atoms_sc = atoms * (4, 4, 1)\n renderer = write(f'Cu{nlayer}O.pov', atoms_sc, rotation='-80x')\n renderer.render()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Note

While the EMT description of Ni, Cu, Pd, Ag, Pt, Au and Al is alright, the\n parameters for C, N and O are not intended for real work!

\n\n" ] } ], "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.13.7" } }, "nbformat": 4, "nbformat_minor": 0 }