pymatgen.analysis.surface_analysis module

class NanoscaleStability(se_analyzers, symprec=1e-05)[source]

Bases: object

A class for analyzing the stability of nanoparticles of different

polymorphs with respect to size. The Wulff shape will be the model for the nanoparticle. Stability will be determined by an energetic competition between the weighted surface energy (surface energy of the Wulff shape) and the bulk energy. A future release will include a 2D phase diagram (e.g. wrt size vs chempot for adsorbed or nonstoichiometric surfaces). Based on the following work:

Kang, S., Mo, Y., Ong, S. P., & Ceder, G. (2014). Nanoscale
stabilization of sodium oxides: Implications for Na-O2 batteries. Nano Letters, 14(2), 1016–1020. https://doi.org/10.1021/nl404557w
se_analyzers
List of SurfaceEnergyPlotter objects. Each item corresponds to a
different polymorph.
symprec

See WulffShape.

Analyzes the nanoscale stability of different polymorphs.

bulk_gform(bulk_entry)[source]

Returns the formation energy of the bulk :param bulk_entry: Entry of the corresponding bulk. :type bulk_entry: ComputedStructureEntry

plot_all_stability_map(max_r, increments=50, u_dict={}, u_default=0, plt=None, labels=[], from_sphere_area=False)[source]
Returns the plot of the formation energy of a particles
of different polymorphs against its effect radius
Parameters:
  • max_r (float) – The maximum radius of the particle to plot up to.
  • increments (int) – Number of plot points
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
  • plt (pylab) – Plot
  • labels (list) – List of labels for each plot, corresponds to the list of se_analyzers
  • from_sphere_area (bool) – There are two ways to calculate the bulk formation energy. Either by treating the volume and thus surface area of the particle as a perfect sphere, or as a Wulff shape.
plot_one_stability_map(analyzer, max_r, u_dict={}, label='', increments=50, u_default=0, plt=None, from_sphere_area=False)[source]
Returns the plot of the formation energy of a particle against its
effect radius
Parameters:
  • analyzer (SurfaceEnergyPlotter) – Analyzer associated with the first polymorph
  • max_r (float) – The maximum radius of the particle to plot up to.
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • label (str) – Label of the plot for legend
  • increments (int) – Number of plot points
  • u_default (float) – Default value for all unset chemical potentials
  • plt (pylab) – Plot
  • from_sphere_area (bool) – There are two ways to calculate the bulk formation energy. Either by treating the volume and thus surface area of the particle as a perfect sphere, or as a Wulff shape.
scaled_wulff(wulffshape, r)[source]
Scales the Wulff shape with an effective radius r. Note that the resulting
Wulff does not neccesarily have the same effective radius as the one provided. The Wulff shape is scaled by its surface energies where first the surface energies are scale by the minimum surface energy and then multiplied by the given effective radius.
Parameters:
  • wulffshape (WulffShape) – Initial, unscaled WulffShape
  • r (float) – Arbitrary effective radius of the WulffShape
Returns:

WulffShape (scaled by r)

solve_equilibrium_point(analyzer1, analyzer2, u_dict={}, u_default=0)[source]
Gives the radial size of two particles where equilibrium is reached
between both particles. NOTE: the solution here is not the same as the solution visualized in the plot because solving for r requires that both the total surface area and volume of the particles are functions of r.
Parameters:
  • analyzer1 (SurfaceEnergyPlotter) – Analyzer associated with the first polymorph
  • analyzer2 (SurfaceEnergyPlotter) – Analyzer associated with the second polymorph
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
Returns:

Particle radius in nm

wulff_gform_and_r(wulffshape, bulk_entry, r, from_sphere_area=False)[source]

Calculates the formation energy of the particle with arbitrary radius r.

Parameters:
  • wulffshape (WulffShape) – Initial, unscaled WulffShape
  • bulk_entry (ComputedStructureEntry) – Entry of the corresponding bulk.
  • r (float (Ang)) – Arbitrary effective radius of the WulffShape
  • from_sphere_area (bool) – There are two ways to calculate the bulk formation energy. Either by treating the volume and thus surface area of the particle as a perfect sphere, or as a Wulff shape.
Returns:

particle formation energy (float in keV), effective radius (float in nm)

class SlabEntry(structure, energy, miller_index, correction=0.0, parameters=None, data=None, entry_id=None, label=None, adsorbates=[], clean_entry=None)[source]

Bases: pymatgen.entries.computed_entries.ComputedStructureEntry

A ComputedStructureEntry object encompassing all data relevant to a
slab for analyzing surface thermodynamics.
miller_index

Miller index of plane parallel to surface.

label

Brief description for this slab.

adsorbates

List of ComputedStructureEntry for the types of adsorbates

..attribute:: clean_entry

SlabEntry for the corresponding clean slab for an adsorbed slab

..attribute:: ads_entries_dict

Dictionary where the key is the reduced composition of the
adsorbate entry and value is the entry itself

Make a SlabEntry containing all relevant surface thermodynamics data.

Parameters:
  • structure (Slab) – The primary slab associated with this entry.
  • energy (float) – Energy from total energy calculation
  • miller_index (tuple(h, k, l)) – Miller index of plane parallel to surface
  • correction (float) – See ComputedSlabEntry
  • parameters (dict) – See ComputedSlabEntry
  • data (dict) – See ComputedSlabEntry
  • entry_id (str) – See ComputedSlabEntry
  • data – See ComputedSlabEntry
  • entry_id – See ComputedSlabEntry
  • label (str) – Any particular label for this slab, e.g. “Tasker 2”, “non-stoichiometric”, “reconstructed”
  • adsorbates ([ComputedStructureEntry]) – List of reference entries for the adsorbates on the slab, can be an isolated molecule (e.g. O2 for O or O2 adsorption), a bulk structure (eg. fcc Cu for Cu adsorption) or anything.
  • clean_entry (ComputedStructureEntry) – If the SlabEntry is for an adsorbed slab, this is the corresponding SlabEntry for the clean slab
Nads_in_slab

Returns the TOTAL number of adsorbates in the slab on BOTH sides

Nsurfs_ads_in_slab

Returns the TOTAL number of adsorbed surfaces in the slab

as_dict()[source]

Returns dict which contains Slab Entry data.

cleaned_up_slab

Returns a slab with the adsorbates removed

create_slab_label

Returns a label (str) for this particular slab based on composition, coverage and Miller index.

classmethod from_dict(d)[source]

Returns a SlabEntry by reading in an dictionary

get_monolayer

Returns the primitive unit surface area density of the adsorbate.

get_unit_primitive_area

Returns the surface area of the adsorbed system per unit area of the primitive slab system.

gibbs_binding_energy(eads=False)[source]
Returns the adsorption energy or Gibb’s binding energy
of an adsorbate on a surface
Parameters:eads (bool) – Whether to calculate the adsorption energy (True) or the binding energy (False) which is just adsorption energy normalized by number of adsorbates.
surface_area

Calculates the surface area of the slab

surface_energy(ucell_entry, ref_entries=[])[source]

Calculates the surface energy of this SlabEntry. :param ucell_entry: An entry object for the bulk :type ucell_entry: entry :param ref_entries (list: [entry]): A list of entries for each type

of element to be used as a reservoir for nonstoichiometric systems. The length of this list MUST be n-1 where n is the number of different elements in the bulk entry. The chempot of the element ref_entry that is not in the list will be treated as a variable.

Returns (Add (Sympy class)): Surface energy

class SurfaceEnergyPlotter(entry_dict, ucell_entry, ref_entries=[])[source]

Bases: object

A class used for generating plots to analyze the thermodynamics of surfaces
of a material. Produces stability maps of different slab configurations, phases diagrams of two parameters to determine stability of configurations (future release), and Wulff shapes.
entry_dict
Nested dictionary containing a list of entries for slab calculations as

items and the corresponding Miller index of the slab as the key. To account for adsorption, each value is a sub-dictionary with the entry of a clean slab calculation as the sub-key and a list of entries for adsorption calculations as the sub-value. The sub-value can contain different adsorption configurations such as a different site or a different coverage, however, ordinarily only the most stable configuration for a particular coverage will be considered as the function of the adsorbed surface energy has an intercept dependent on the adsorption energy (ie an adsorption site with a higher adsorption energy will always provide a higher surface energy than a site with a lower adsorption energy). An example parameter is provided: {(h1,k1,l1): {clean_entry1: [ads_entry1, ads_entry2, …],

clean_entry2: […], …}, (h2,k2,l2): {…}}

where clean_entry1 can be a pristine surface and clean_entry2 can be a reconstructed surface while ads_entry1 can be adsorption at site 1 with a 2x2 coverage while ads_entry2 can have a 3x3 coverage. If adsorption entries are present (i.e. if entry_dict[(h,k,l)][clean_entry1]), we consider adsorption in all plots and analysis for this particular facet.

..attribute:: color_dict

Dictionary of colors (r,g,b,a) when plotting surface energy stability. The
keys are individual surface entries where clean surfaces have a solid color while the corresponding adsorbed surface will be transparent.
ucell_entry

ComputedStructureEntry of the bulk reference for this particular material.

ref_entries

List of ComputedStructureEntries to be used for calculating chemical potential.

color_dict

Randomly generated dictionary of colors associated with each facet.

Object for plotting surface energy in different ways for clean and
adsorbed surfaces.
Parameters:
  • entry_dict (dict) – Dictionary containing a list of entries for slab calculations. See attributes.
  • ucell_entry (ComputedStructureEntry) – ComputedStructureEntry of the bulk reference for this particular material.
  • ref_entries ([ComputedStructureEntries]) – A list of entries for each type of element to be used as a reservoir for nonstoichiometric systems. The length of this list MUST be n-1 where n is the number of different elements in the bulk entry. The chempot of the element ref_entry that is not in the list will be treated as a variable. e.g. if your ucell_entry is for LiFePO4 than your ref_entries should have an entry for Li, Fe, and P if you want to use the chempot of O as the variable.
BE_vs_clean_SE(u_dict, u_default=0, plot_eads=False, annotate_monolayer=True, JPERM2=False)[source]
For each facet, plot the clean surface energy against the most
stable binding energy.
Parameters:
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
  • plot_eads (bool) – Option to plot the adsorption energy (binding energy multiplied by number of adsorbates) instead.
  • annotate_monolayer (bool) – Whether or not to label each data point with its monolayer (adsorbate density per unit primiitve area)
  • JPERM2 (bool) – Whether to plot surface energy in /m^2 (True) or eV/A^2 (False)
Returns:

Plot of clean surface energy vs binding energy for

all facets.

Return type:

(Plot)

area_frac_vs_chempot_plot(ref_delu, chempot_range, u_dict={}, u_default=0, increments=10)[source]

1D plot. Plots the change in the area contribution of each facet as a function of chemical potential.

Parameters:
  • ref_delu (sympy Symbol) – The free variable chempot with the format: Symbol(“delu_el”) where el is the name of the element.
  • chempot_range (list) – Min/max range of chemical potential to plot along
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
  • increments (int) – Number of data points between min/max or point of intersection. Defaults to 10 points.
Returns:

Plot of area frac on the Wulff shape

for each facet vs chemical potential.

Return type:

(Pylab)

chempot_plot_addons(plt, xrange, ref_el, axes, pad=2.4, rect=[-0.047, 0, 0.84, 1], ylim=[])[source]

Helper function to a chempot plot look nicer.

Parameters:
  • plt (Plot) –
  • xrange (list) – xlim parameter
  • ref_el (str) – Element of the referenced chempot.
  • axes (axes) –
  • pad (float) –
  • rect (list) – For tight layout
  • ylim (ylim parameter) –

return (Plot): Modified plot with addons.

chempot_vs_gamma(ref_delu, chempot_range, miller_index=(), u_dict={}, u_default=0, JPERM2=False, show_unstable=False, ylim=[], plt=None, no_clean=False, no_doped=False, use_entry_labels=False, no_label=False)[source]
Plots the surface energy as a function of chemical potential.
Each facet will be associated with its own distinct colors. Dashed lines will represent stoichiometries different from that of the mpid’s compound. Transparent lines indicates adsorption.
Parameters:
  • ref_delu (sympy Symbol) – The range stability of each slab is based on the chempot range of this chempot. Should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element
  • chempot_range ([max_chempot, min_chempot]) – Range to consider the stability of the slabs.
  • miller_index (list) – Miller index for a specific facet to get a dictionary for.
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
  • JPERM2 (bool) – Whether to plot surface energy in /m^2 (True) or eV/A^2 (False)
  • show_unstable (bool) – Whether or not to show parts of the surface energy plot outside the region of stability.
  • ylim ([ymax, ymin]) – Range of y axis
  • no_doped (bool) – Whether to plot for the clean slabs only.
  • no_clean (bool) – Whether to plot for the doped slabs only.
  • use_entry_labels (bool) – If True, will label each slab configuration according to their given label in the SlabEntry object.
  • no_label (bool) – Option to turn off labels.
Returns:

Plot of surface energy vs chempot for all entries.

Return type:

(Plot)

chempot_vs_gamma_plot_one(plt, entry, ref_delu, chempot_range, u_dict={}, u_default=0, label='', JPERM2=False)[source]

Helper function to help plot the surface energy of a single SlabEntry as a function of chemical potential.

Parameters:
  • plt (Plot) – A plot.
  • entry (SlabEntry) – Entry of the slab whose surface energy we want to plot
  • ref_delu (sympy Symbol) – The range stability of each slab is based on the chempot range of this chempot. Should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element
  • chempot_range ([max_chempot, min_chempot]) – Range to consider the stability of the slabs.
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
  • label (str) – Label of the slab for the legend.
  • JPERM2 (bool) – Whether to plot surface energy in /m^2 (True) or eV/A^2 (False)
Returns:

Plot of surface energy vs chemical potential for one entry.

Return type:

(Plot)

color_palette_dict(alpha=0.35)[source]

Helper function to assign each facet a unique color using a dictionary.

Parameters:alpha (float) – Degree of transparency
return (dict): Dictionary of colors (r,g,b,a) when plotting surface
energy stability. The keys are individual surface entries where clean surfaces have a solid color while the corresponding adsorbed surface will be transparent.
get_stable_entry_at_u(miller_index, u_dict={}, u_default=0, no_doped=False, no_clean=False)[source]
Returns the entry corresponding to the most stable slab for a particular
facet at a specific chempot. We assume that surface energy is constant so all free variables must be set with u_dict, otherwise they are assumed to be equal to u_default.
Parameters:
  • miller_index ((h,k,l)) – The facet to find the most stable slab in
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
  • no_doped (bool) – Consider stability of clean slabs only.
  • no_clean (bool) – Consider stability of doped slabs only.
Returns:

SlabEntry, surface_energy (float)

get_surface_equilibrium(slab_entries, u_dict={})[source]
Takes in a list of SlabEntries and calculates the chemical potentials
at which all slabs in the list coexists simultaneously. Useful for building surface phase diagrams. Note that to solve for x equations (x slab_entries), there must be x free variables (chemical potentials). Adjust u_dict as need be to get the correct number of free variables.
Parameters:
  • slab_entries (array) – The coefficients of the first equation
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
Returns:

Array containing a solution to x equations with x

variables (x-1 chemical potential and 1 surface energy)

Return type:

(array)

monolayer_vs_BE(plot_eads=False)[source]
Plots the binding energy energy as a function of monolayers (ML), i.e.
the fractional area adsorbate density for all facets. For each facet at a specific monlayer, only plot the lowest binding energy.
Parameters:plot_eads (bool) – Option to plot the adsorption energy (binding energy multiplied by number of adsorbates) instead.
Returns:Plot of binding energy vs monolayer for all facets.
Return type:(Plot)
set_all_variables(entry, u_dict, u_default)[source]
Sets all chemical potential values and returns a dictionary where
the key is a sympy Symbol and the value is a float (chempot).
Parameters:
  • entry (SlabEntry) – Computed structure entry of the slab
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
Returns:

Dictionary of set chemical potential values

stable_u_range_dict(chempot_range, ref_delu, no_doped=True, u_dict={}, miller_index=())[source]

Creates a dictionary where each entry is a key pointing to a chemical potential range where the surface of that entry is stable. Does so by enumerating through all possible solutions (intersect) for surface energies of a specific facet.

Parameters:
  • chempot_range ([max_chempot, min_chempot]) – Range to consider the stability of the slabs.
  • ref_delu (sympy Symbol) – The range stability of each slab is based on the chempot range of this chempot. Should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element
  • no_doped (bool) – Consider stability of clean slabs only.
  • no_clean (bool) – Consider stability of doped slabs only.
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • miller_index (list) – Miller index for a specific facet to get a dictionary for.
wulff_from_chempot(u_dict={}, u_default=0, symprec=1e-05, no_clean=False, no_doped=False)[source]

Method to get the Wulff shape at a specific chemical potential.

Parameters:
  • u_dict (Dict) – Dictionary of the chemical potentials to be set as constant. Note the key should be a sympy Symbol object of the format: Symbol(“delu_el”) where el is the name of the element.
  • u_default (float) – Default value for all unset chemical potentials
  • symprec (float) – See WulffShape.
  • no_doped (bool) – Consider stability of clean slabs only.
  • no_clean (bool) – Consider stability of doped slabs only.
Returns:

The WulffShape at u_ref and u_ads.

Return type:

(WulffShape)