{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Introduction\n", "\n", "This notebook demonstrates how to carry out an ordering of a disordered structure using pymatgen." ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Full Formula (Cu2 Au2)\n", "Reduced Formula: CuAu\n", "abc : 3.677000 3.677000 3.677000\n", "angles: 90.000000 90.000000 90.000000\n", "Sites (4)\n", " # SP a b c\n", "--- ---------------------- --- --- ---\n", " 0 Cu0+:0.500, Au0+:0.500 0 0 0\n", " 1 Cu0+:0.500, Au0+:0.500 0 0.5 0.5\n", " 2 Cu0+:0.500, Au0+:0.500 0.5 0 0.5\n", " 3 Cu0+:0.500, Au0+:0.500 0.5 0.5 0\n" ] } ], "source": [ "# Let us start by creating a disordered CuAu fcc structure.\n", "\n", "from pymatgen import Structure, Lattice\n", "\n", "specie = {\"Cu0+\": 0.5, \"Au0+\": 0.5}\n", "cuau = Structure.from_spacegroup(\"Fm-3m\", Lattice.cubic(3.677), [specie], [[0, 0, 0]])\n", "print cuau" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that each site is now 50% occupied by Cu and Au. Because the ordering algorithms uses an Ewald summation to rank the structures, you need to explicitly specify the oxidation state for each species, even if it is 0. Let us now perform ordering of these sites using two methods." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Method 1 - Using the OrderDisorderedStructureTransformation\n", "\n", "The first method is to use the OrderDisorderedStructureTransformation." ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "6\n", "{u'energy': 0.0, u'energy_above_minimum': 0.0, u'structure': Structure Summary\n", "Lattice\n", " abc : 3.677 3.677 3.677\n", " angles : 90.0 90.0 90.0\n", " volume : 49.714249733000003\n", " A : 3.677 0.0 0.0\n", " B : 0.0 3.677 0.0\n", " C : 0.0 0.0 3.677\n", "PeriodicSite: Cu0+ (0.0000, 0.0000, 0.0000) [0.0000, 0.0000, 0.0000]\n", "PeriodicSite: Cu0+ (0.0000, 1.8385, 1.8385) [0.0000, 0.5000, 0.5000]\n", "PeriodicSite: Au0+ (1.8385, 0.0000, 1.8385) [0.5000, 0.0000, 0.5000]\n", "PeriodicSite: Au0+ (1.8385, 1.8385, 0.0000) [0.5000, 0.5000, 0.0000]}\n" ] } ], "source": [ "from pymatgen.transformations.standard_transformations import OrderDisorderedStructureTransformation\n", "\n", "trans = OrderDisorderedStructureTransformation()\n", "\n", "ss = trans.apply_transformation(cuau, return_ranked_list=100)\n", "\n", "print(len(ss))\n", "\n", "print(ss[0])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that the OrderDisorderedTransformation (with a sufficiently large return_ranked_list parameter) returns all orderings, including duplicates without accounting for symmetry. A computed ewald energy is returned together with each structure. To eliminate duplicates, the best way is to use StructureMatcher's group_structures method, as demonstrated below." ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "1\n", "Full Formula (Cu2 Au2)\n", "Reduced Formula: CuAu\n", "abc : 3.677000 3.677000 3.677000\n", "angles: 90.000000 90.000000 90.000000\n", "Sites (4)\n", " # SP a b c\n", "--- ---- --- --- ---\n", " 0 Cu0+ 0 0 0\n", " 1 Cu0+ 0 0.5 0.5\n", " 2 Au0+ 0.5 0 0.5\n", " 3 Au0+ 0.5 0.5 0\n" ] } ], "source": [ "from pymatgen.analysis.structure_matcher import StructureMatcher\n", "\n", "matcher = StructureMatcher()\n", "groups = matcher.group_structures([d[\"structure\"] for d in ss])\n", "print(len(groups))\n", "print(groups[0][0])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Method 2 - Using the EnumerateStructureTransformation\n", "\n", "If you have enumlib installed, you can use the EnumerateStructureTransformation. This automatically takes care of symmetrically equivalent orderings and can enumerate supercells, but is much more prone to parameter sensitivity and cannot handle very large structures. The example below shows an enumerate of CuAu up to cell sizes of 4." ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from pymatgen.transformations.advanced_transformations import EnumerateStructureTransformation\n", "specie = {\"Cu\": 0.5, \"Au\": 0.5}\n", "cuau = Structure.from_spacegroup(\"Fm-3m\", Lattice.cubic(3.677), [specie], [[0, 0, 0]])\n", "\n", "trans = EnumerateStructureTransformation(max_cell_size=3)\n", "ss = trans.apply_transformation(cuau, return_ranked_list=1000)" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "94\n", "cell sizes are [4, 8, 8, 8, 8, 8, 8, 8, 8, 8, 8, 8, 8, 8, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12, 12]\n" ] } ], "source": [ "print(len(ss))\n", "print(\"cell sizes are %s\" % ([len(d[\"structure\"]) for d in ss]))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that structures with cell sizes ranging from 1-3x the unit cell size is generated." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Conclusion\n", "\n", "This notebook illustrates two basic ordering/enumeration approaches. In general, OrderDisorderedTransformation works better for large cells and is useful if you need just any quick plausible ordering. EnumerateStructureTransformation is more rigorous, but is prone to sensitivity errors and may require fiddling with various parameters." ] } ], "metadata": { "kernelspec": { "display_name": "Python 2", "language": "python", "name": "python2" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 2 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython2", "version": "2.7.11" } }, "nbformat": 4, "nbformat_minor": 0 }