Building things¶
Operate adsorbate¶
add_adsorbate function¶
acat.build.action.
add_adsorbate
(atoms, adsorbate, site=None, surface=None, morphology=None, indices=None, height=None, composition=None, orientation=None, tilt_angle=0.0, subsurf_element=None, both_sides=False, all_sites=None)[source]¶A general function for adding one adsorbate to the surface. Note that this function adds one adsorbate to a random site that meets the specified condition regardless of it is already occupied or not. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
atoms (ase.Atoms object) – Accept any ase.Atoms object. No need to be built-in.
adsorbate (str or ase.Atom object or ase.Atoms object) – The adsorbate species to be added onto the surface.
site (str, default None) – The site type that the adsorbate should be added to.
surface (str, default None) – The surface type (crystal structure + Miller indices) If the structure is a periodic surface slab, this is required. If the structure is a nanoparticle, the function enumerates only the sites on the specified surface.
morphology (str, default None) – The morphology type that the adsorbate should be added to. Only available for surface slabs.
indices (list or tuple) – The indices of the atoms that contribute to the site that you want to add adsorbate to. This has the highest priority.
height (float, default None) – The height of the added adsorbate from the surface. Use the default settings if not specified.
composition (str, default None) – The elemental of the site that should be added to.
orientation (list or numpy.array, default None) – The vector that the multidentate adsorbate is aligned to.
tilt_angle (float, default 0.) – Tilt the adsorbate with an angle (in degress) relative to the surface normal.
subsurf_element (str, default None) – The subsurface element of the hcp or 4fold hollow site that should be added to.
both_sides (bool, default False) – Whether to consider sites on both top and bottom sides of the slab. Only relevant for periodic surface slabs.
all_sites (list of dicts, default None) – The list of all sites. Provide this to make the function much faster. Useful when the function is called many times.
Example
To add a NO molecule to a bridge site consists of one Pt and one Ni on the fcc111 surface of a bimetallic truncated octahedron:
>>> from acat.build.action import add_adsorbate >>> from ase.cluster import Octahedron >>> from ase.visualize import view >>> atoms = Octahedron('Ni', length=7, cutoff=2) >>> for atom in atoms: ... if atom.index % 2 == 0: ... atom.symbol = 'Pt' >>> add_adsorbate(atoms, adsorbate='NO', site='bridge', ... surface='fcc111', composition='NiPt') >>> view(atoms)Output:
![]()
add_adsorbate_to_site function¶
acat.build.action.
add_adsorbate_to_site
(atoms, adsorbate, site, height=None, orientation=None, tilt_angle=0.0)[source]¶The base function for adding one adsorbate to a site. Site must include information of ‘normal’ and ‘position’. Useful for adding adsorbate to multiple sites or adding multidentate adsorbates.
- Parameters
atoms (ase.Atoms object) – Accept any ase.Atoms object. No need to be built-in.
adsorbate (str or ase.Atom object or ase.Atoms object) – The adsorbate species to be added onto the surface.
site (dict) – The site that the adsorbate should be added to. Must contain information of the position and the normal vector of the site.
height (float, default None) – The height of the added adsorbate from the surface. Use the default settings if not specified.
orientation (list or numpy.array, default None) – The vector that the multidentate adsorbate is aligned to.
tilt_angle (float, default None) – Tilt the adsorbate with an angle (in degress) relative to the surface normal.
Example
To add CO to all fcc sites of an icosahedral nanoparticle:
>>> from acat.adsorption_sites import ClusterAdsorptionSites >>> from acat.build.action import add_adsorbate_to_site >>> from ase.cluster import Icosahedron >>> from ase.visualize import view >>> atoms = Icosahedron('Pt', noshells=5) >>> atoms.center(vacuum=5.) >>> cas = ClusterAdsorptionSites(atoms) >>> fcc_sites = cas.get_sites(site='fcc') >>> for site in fcc_sites: ... add_adsorbate_to_site(atoms, adsorbate='CO', site=site) >>> view(atoms)Output:
![]()
To add a bidentate CH3OH to the (54, 57, 58) site on a Pt fcc111 surface slab and rotate the orientation to a neighbor site:
>>> from acat.adsorption_sites import SlabAdsorptionSites >>> from acat.adsorption_sites import get_adsorption_site >>> from acat.build.action import add_adsorbate_to_site >>> from acat.utilities import get_mic >>> from ase.build import fcc111 >>> from ase.visualize import view >>> atoms = fcc111('Pt', (4, 4, 4), vacuum=5.) >>> i, site = get_adsorption_site(atoms, indices=(54, 57, 58), ... surface='fcc111', ... return_index=True) >>> sas = SlabAdsorptionSites(atoms, surface='fcc111') >>> sites = sas.get_sites() >>> nbsites = sas.get_neighbor_site_list(neighbor_number=1) >>> nbsite = sites[nbsites[i][0]] # Choose the first neighbor site >>> ori = get_mic(site['position'], nbsite['position'], atoms.cell) >>> add_adsorbate_to_site(atoms, adsorbate='CH3OH', site=site, ... orientation=ori) >>> view(atoms)Output:
![]()
add_adsorbate_to_label function¶
acat.build.action.
add_adsorbate_to_label
(atoms, adsorbate, label, surface=None, height=None, orientation=None, tilt_angle=0.0, composition_effect=False, both_sides=False, all_sites=None)[source]¶Same as add_adsorbate function, except that the site type is represented by a numerical label. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
atoms (ase.Atoms object) – Accept any ase.Atoms object. No need to be built-in.
adsorbate (str or ase.Atom object or ase.Atoms object) – The adsorbate species to be added onto the surface.
label (int or str) – The label of the site that the adsorbate should be added to.
surface (str, default None) – The surface type (crystal structure + Miller indices) If the structure is a periodic surface slab, this is required. If the structure is a nanoparticle, the function enumerates only the sites on the specified surface.
height (float, default None) – The height of the added adsorbate from the surface. Use the default settings if not specified.
orientation (list or numpy.array, default None) – The vector that the multidentate adsorbate is aligned to.
tilt_angle (float, default 0.) – Tilt the adsorbate with an angle (in degress) relative to the surface normal.
composition_effect (bool, default False) – Whether the label is defined in bimetallic labels or not.
both_sides (bool, default False) – Whether to consider sites on both top and bottom sides of the slab. Only relevant for periodic surface slabs.
all_sites (list of dicts, default None) – The list of all sites. Provide this to make the function much faster. Useful when the function is called many times.
Example
To add a NH molecule to a site with bimetallic label 14 (an hcp CuCuAu site) on a bimetallic fcc110 surface slab:
>>> from acat.build.action import add_adsorbate_to_label >>> from ase.build import fcc110 >>> from ase.visualize import view >>> atoms = fcc110('Cu', (3, 3, 8), vacuum=5.) >>> for atom in atoms: ... if atom.index % 2 == 0: ... atom.symbol = 'Au' ... atoms.center() >>> add_adsorbate_to_label(atoms, adsorbate='NH', label=14, ... surface='fcc110', composition_effect=True) >>> view(atoms)Output:
![]()
remove_adsorbate_from_site function¶
acat.build.action.
remove_adsorbate_from_site
(atoms, site, remove_fragment=False)[source]¶The base function for removing one adsorbate from an occupied site. The site must include information of ‘adsorbate_indices’ or ‘fragment_indices’. Note that if you want to remove adsorbates from multiple sites, call this function multiple times will return the wrong result. Please use remove_adsorbates_from_sites instead.
- Parameters
atoms (ase.Atoms object) – Accept any ase.Atoms object. No need to be built-in.
site (dict) – The site that the adsorbate should be removed from. Must contain information of the adsorbate indices.
remove_fragment (bool, default False) – Remove the fragment of a multidentate adsorbate instead of the whole adsorbate.
Example
To remove a CO molecule from a fcc111 surface slab with one CO and one OH:
>>> from acat.adsorption_sites import SlabAdsorptionSites >>> from acat.adsorbate_coverage import SlabAdsorbateCoverage >>> from acat.build.action import add_adsorbate_to_site >>> from acat.build.action import remove_adsorbate_from_site >>> from ase.build import fcc111 >>> from ase.visualize import view >>> atoms = fcc111('Pt', (6, 6, 4), 4, vacuum=5.) >>> atoms.center() >>> sas = SlabAdsorptionSites(atoms, surface='fcc111') >>> sites = sas.get_sites() >>> add_adsorbate_to_site(atoms, adsorbate='OH', site=sites[0]) >>> add_adsorbate_to_site(atoms, adsorbate='CO', site=sites[-1]) >>> sac = SlabAdsorbateCoverage(atoms, sas) >>> occupied_sites = sac.get_sites(occupied_only=True) >>> CO_site = next((s for s in occupied_sites if s['adsorbate'] == 'CO')) >>> remove_adsorbate_from_site(atoms, site=CO_site) >>> view(atoms)![]()
remove_adsorbates_from_sites function¶
acat.build.action.
remove_adsorbates_from_sites
(atoms, sites, remove_fragments=False)[source]¶The base function for removing multiple adsorbates from an occupied site. The sites must include information of ‘adsorbate_indices’ or ‘fragment_indices’.
- Parameters
atoms (ase.Atoms object) – Accept any ase.Atoms object. No need to be built-in.
sites (list of dicts) – The site that the adsorbate should be removed from. Must contain information of the adsorbate indices.
remove_fragments (bool, default False) – Remove the fragment of a multidentate adsorbate instead of the whole adsorbate.
Example
To remove all CO species from a fcc111 surface slab covered with both CO and OH:
>>> from acat.adsorption_sites import SlabAdsorptionSites >>> from acat.adsorbate_coverage import SlabAdsorbateCoverage >>> from acat.build.overlayer import random_coverage_pattern >>> from acat.build.action import remove_adsorbates_from_sites >>> from ase.build import fcc111 >>> from ase.visualize import view >>> slab = fcc111('Pt', (6, 6, 4), 4, vacuum=5.) >>> slab.center() >>> atoms = random_coverage_pattern(slab, adsorbate_species=['OH','CO'], ... surface='fcc111', ... min_adsorbate_distance=5.) >>> sas = SlabAdsorptionSites(atoms, surface='fcc111') >>> sac = SlabAdsorbateCoverage(atoms, sas) >>> occupied_sites = sac.get_sites(occupied_only=True) >>> CO_sites = [s for s in occupied_sites if s['adsorbate'] == 'CO'] >>> remove_adsorbates_from_sites(atoms, sites=CO_sites) >>> view(atoms)Output:
![]()
remove_adsorbates_too_close function¶
acat.build.action.
remove_adsorbates_too_close
(atoms, adsorbate_coverage=None, surface=None, min_adsorbate_distance=0.5)[source]¶Find adsorbates that are too close, remove one set of them. The function is intended to remove atoms that are unphysically close. Please do not use a min_adsorbate_distace larger than 2. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
atoms (ase.Atoms object) – The nanoparticle or surface slab onto which the adsorbates are added. Accept any ase.Atoms object. No need to be built-in.
adsorbate_coverage (acat.adsorbate_coverage.ClusterAdsorbateCoverage object or acat.adsorbate_coverage.SlabAdsorbateCoverage object, default None) – The built-in adsorbate coverage class.
surface (str, default None) – The surface type (crystal structure + Miller indices). If the structure is a periodic surface slab, this is required. If the structure is a nanoparticle, the function only add adsorbates to the sites on the specified surface.
min_adsorbate_distance (float, default 0.) – The minimum distance between two atoms that is not considered to be to close. This distance has to be small.
Example
To remove unphysically close adsorbates on the edges of a Marks decahedron with 0.75 ML ordered CO coverage:
>>> from acat.build.overlayer import ordered_coverage_pattern >>> from acat.build.action import remove_adsorbates_too_close >>> from ase.cluster import Decahedron >>> from ase.visualize import view >>> atoms = Decahedron('Pt', p=4, q=3, r=1) >>> atoms.center(vacuum=5.) >>> pattern = ordered_coverage_pattern(atoms, adsorbate='CO', ... coverage=0.75) >>> remove_adsorbates_too_close(pattern, min_adsorbate_distance=1.) >>> view(pattern)Output:
![]()
Generate adsorbate overlayer patterns¶
StochasticPatternGenerator class¶
- class
acat.build.overlayer.
StochasticPatternGenerator
(images, adsorbate_species, image_probabilities=None, species_probabilities=None, min_adsorbate_distance=1.5, adsorption_sites=None, surface=None, heights={'3fold': 1.3, '4fold': 1.3, '5fold': 1.5, '6fold': 0.0, 'bridge': 1.5, 'fcc': 1.3, 'hcp': 1.3, 'longbridge': 1.5, 'ontop': 1.8, 'shortbridge': 1.5}, allow_6fold=False, composition_effect=True, both_sides=False, dmax=2.5, species_forbidden_sites=None, species_forbidden_labels=None, fragmentation=True, trajectory='patterns.traj', append_trajectory=False, logfile='patterns.log')[source]¶StochasticPatternGenerator is a class for generating adsorbate overlayer patterns stochastically. Graph isomorphism is implemented to identify identical overlayer patterns. 4 adsorbate actions are supported: add, remove, move, replace. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
images (ase.Atoms object or list of ase.Atoms objects) – The structure to perform the adsorbate actions on. If a list of structures is provided, perform one adsorbate action on one of the structures in each step. Accept any ase.Atoms object. No need to be built-in.
adsorbate_species (str or list of strs) – A list of adsorbate species to be randomly added to the surface.
image_probabilities (listt, default None) – A list of the probabilities of selecting each structure. Selecting structure with equal probability if not specified.
species_probabilities (dict, default None) – A dictionary that contains keys of each adsorbate species and values of their probabilities of adding onto the surface. Adding adsorbate species with equal probability if not specified.
min_adsorbate_distance (float, default 1.5) – The minimum distance constraint between two atoms that belongs to two adsorbates.
adsorption_sites (acat.adsorption_sites.ClusterAdsorptionSites object or acat.adsorption_sites.SlabAdsorptionSites object, default None) – Provide the built-in adsorption sites class to accelerate the pattern generation. Make sure all the structures have the same periodicity and atom indexing. If composition_effect=True, you should only provide adsorption_sites when the surface composition is fixed.
surface (str, default None) – The surface type (crystal structure + Miller indices). Only required if the structure is a periodic surface slab.
heights (dict, default acat.settings.site_heights) – A dictionary that contains the adsorbate height for each site type. Use the default height settings if the height for a site type is not specified.
allow_6fold (bool, default False) – Whether to allow the adsorption on 6-fold subsurf sites underneath fcc hollow sites.
composition_effect (bool, default False) – Whether to consider sites with different elemental compositions as different sites. It is recommended to set composition=False for monometallics.
both_sides (bool, default False) – Whether to add adsorbate to both top and bottom sides of the slab. Only works if adsorption_sites is not provided, otherwise consistent with what is specified in adsorption_sites. Only relvevant for periodic surface slabs.
dmax (float, default 2.5) – The maximum bond length (in Angstrom) between the site and the bonding atom that should be considered as an adsorbate.
species_forbidden_sites (dict, default None) – A dictionary that contains keys of each adsorbate species and values of the site (can be one or multiple site types) that the speices is not allowed to add to. All sites are availabe for a species if not specified. Note that this does not differentiate sites with different compositions.
species_forbidden_labels (dict, default None) – Same as species_forbidden_sites except that the adsorption sites are written as numerical labels according to acat.labels. Useful when you need to differentiate sites with different compositions.
fragmentation (bool, default True) – Whether to cut multidentate species into fragments. This ensures that multidentate species with different orientations are considered as different overlayer patterns.
trajectory (str, default 'patterns.traj') – The name of the output ase trajectory file.
append_trajectory (bool, default False) –
Whether to append structures to the existing trajectory. If only unique patterns are accepted, the code will also check graph isomorphism for the existing structures in the trajectory. This is also useful when you want to generate overlayer patterns stochastically but for all images systematically, e.g. generating 10 stochastic overlayer patterns for each image:
>>> from acat.build.overlayer import StochasticPatternGenerator as SPG >>> for atoms in images: ... spg = SPG(atoms, ..., append_trajectory=True) ... spg.run(num_gen = 10)logfile (str, default 'patterns.log') – The name of the log file.
run
(num_gen, action='add', num_act=1, action_probabilities=None, unique=True, subsurf_effect=False)[source]¶Run the pattern generator.
- Parameters
num_gen (int) – Number of patterns to generate.
action (str or list of strs, default 'add') – Action(s) to perform. If a list of actions is provided, select actions from the list randomly or with probabilities.
num_act (int, default 1) – Number of times performed for each action. Useful for operating more than one adsorbates at a time.
action_probabilities (dict, default None) – A dictionary that contains keys of each action and values of the corresponding probabilities. Select actions with equal probability if not specified.
unique (bool, default True) – Whether to discard duplicate patterns based on isomorphism.
subsurf_effect (bool, default False) – Take subsurface atoms into consideration when checking uniqueness. Could be important for surfaces like fcc100.
Example
The following example illustrates how to generate 100 stochastic adsorbate overlayer patterns with CO, OH, CH3 and CHO, based on 10 Pt fcc111 surface slabs with random C and O coverages, where CH3 is forbidden to be added to ontop and bridge sites:
>>> from acat.build.overlayer import StochasticPatternGenerator as SPG >>> from acat.build.overlayer import random_coverage_pattern >>> from ase.build import fcc111 >>> from ase.io import read >>> from ase.visualize import view >>> slab = fcc111('Pt', (6, 6, 4), 4, vacuum=5.) >>> slab.center() >>> images = [] >>> for _ in range(10): ... atoms = slab.copy() ... image = random_coverage_pattern(atoms, adsorbate_species=['C','O'], ... surface='fcc111', ... min_adsorbate_distance=5.) ... images.append(image) >>> spg = SPG(images, adsorbate_species=['CO','OH','CH3','CHO'], ... species_probabilities={'CO': 0.3, 'OH': 0.3, ... 'CH3': 0.2, 'CHO': 0.2}, ... min_adsorbate_distance=1.5, ... surface='fcc111', ... composition_effect=False, ... species_forbidden_sites={'CH3': ['ontop','bridge']}) >>> spg.run(num_gen=100, action='add') >>> images = read('patterns.traj', index=':') >>> view(images)Output:
![]()
The following example illustrates how to generate 20 unique coverage patterns, each adding 4 adsorbates (randomly chosen from H, OH and H2O) onto a fcc100 Ni2Cu surface slab on both top and bottom interfaces (
bulk water
in between) with probabilities of 0.25, 0.25, 0.5, respectively, and a minimum adsorbate distance of 2.5 Angstrom:>>> from acat.build.overlayer import StochasticPatternGenerator as SPG >>> from acat.adsorption_sites import SlabAdsorptionSites >>> from ase.io import read >>> from ase.build import fcc100 >>> from ase.visualize import view >>> water_bulk = read('water_bulk.xyz') >>> water_bulk.center(vacuum=11., axis=2) >>> slab = fcc100('Ni', (4, 4, 8), vacuum=9.5, periodic=True) >>> for atom in slab: ... if atom.index % 3 == 0: ... atom.symbol = 'Cu' >>> slab.translate(-slab.cell[2] / 2) >>> slab.wrap() >>> atoms = slab + water_bulk >>> sas = SlabAdsorptionSites(atoms, surface='fcc100', ... composition_effect=True, ... both_sides=True) >>> spg = SPG(atoms, adsorbate_species=['H','OH','OH2'], ... species_probabilities={'H': 0.25, 'OH': 0.25, 'OH2': 0.5}, ... min_adsorbate_distance=2.5, ... adsorption_sites=sas, ... surface='fcc100') >>> spg.run(num_gen=20, action='add', num_act=4) >>> images = read('patterns.traj', index=':') >>> view(images)Output:
![]()
SystematicPatternGenerator class¶
- class
acat.build.overlayer.
SystematicPatternGenerator
(images, adsorbate_species, min_adsorbate_distance=1.5, adsorption_sites=None, surface=None, heights={'3fold': 1.3, '4fold': 1.3, '5fold': 1.5, '6fold': 0.0, 'bridge': 1.5, 'fcc': 1.3, 'hcp': 1.3, 'longbridge': 1.5, 'ontop': 1.8, 'shortbridge': 1.5}, allow_6fold=False, composition_effect=True, both_sides=False, dmax=2.5, species_forbidden_sites=None, species_forbidden_labels=None, enumerate_orientations=True, trajectory='patterns.traj', append_trajectory=False, logfile='patterns.log')[source]¶SystematicPatternGenerator is a class for generating adsorbate overlayer patterns systematically. This is useful to enumerate all unique patterns at low coverage, but explodes at higher coverages. Graph isomorphism is implemented to identify identical overlayer patterns. 4 adsorbate actions are supported: add, remove, move, replace. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
images (ase.Atoms object or list of ase.Atoms objects) – The structure to perform the adsorbate actions on. If a list of structures is provided, perform one adsorbate action on one of the structures in each step. Accept any ase.Atoms object. No need to be built-in.
adsorbate_species (str or list of strs) – A list of adsorbate species to be randomly added to the surface.
min_adsorbate_distance (float, default 1.5) – The minimum distance constraint between two atoms that belongs to two adsorbates.
adsorption_sites (acat.adsorption_sites.ClusterAdsorptionSites object or acat.adsorption_sites.SlabAdsorptionSites object, default None) – Provide the built-in adsorption sites class to accelerate the pattern generation. Make sure all the structures have the same periodicity and atom indexing. If composition_effect=True, you should only provide adsorption_sites when the surface composition is fixed.
surface (str, default None) – The surface type (crystal structure + Miller indices). Only required if the structure is a periodic surface slab.
heights (dict, default acat.settings.site_heights) – A dictionary that contains the adsorbate height for each site type. Use the default height settings if the height for a site type is not specified.
allow_6fold (bool, default False) – Whether to allow the adsorption on 6-fold subsurf sites underneath fcc hollow sites.
composition_effect (bool, default False) – Whether to consider sites with different elemental compositions as different sites. It is recommended to set composition=False for monometallics.
both_sides (bool, default False) – Whether to add adsorbate to both top and bottom sides of the slab. Only works if adsorption_sites is not provided, otherwise consistent with what is specified in adsorption_sites. Only relevant for periodic surface slabs.
dmax (float, default 2.5) – The maximum bond length (in Angstrom) between the site and the bonding atom that should be considered as an adsorbate.
species_forbidden_sites (dict, default None) – A dictionary that contains keys of each adsorbate species and values of the site (can be one or multiple site types) that the speices is not allowed to add to. All sites are availabe for a species if not specified. Note that this does not differentiate sites with different compositions.
species_forbidden_labels (dict, default None) – Same as species_forbidden_sites except that the adsorption sites are written as numerical labels according to acat.labels. Useful when you need to differentiate sites with different compositions.
enumerate_orientations (bool, default True) – Whether to enumerate all orientations of multidentate species. This ensures that multidentate species with different orientations are all enumerated.
trajectory (str, default 'patterns.traj') – The name of the output ase trajectory file.
append_trajectory (bool, default False) – Whether to append structures to the existing trajectory. If only unique patterns are accepted, the code will also check graph isomorphism for the existing structures in the trajectory.
logfile (str, default 'patterns.log') – The name of the log file.
run
(max_gen_per_image=None, action='add', num_act=1, unique=True, subsurf_effect=False)[source]¶Run the pattern generator.
- Parameters
max_gen_per_image (int, default None) – Maximum number of patterns to generate for each image. Enumerate all possible patterns if not specified.
action (str, defualt 'add') – Action to perform.
num_act (int, default 1) – Number of times performed for the action. Useful for operating more than one adsorbates at a time.
unique (bool, default True) – Whether to discard duplicate patterns based on isomorphism.
subsurf_effect (bool, default False) – Take subsurface atoms into consideration when checking uniqueness. Could be important for surfaces like fcc100.
Example
The following example illustrates how to add CO to all unique sites on a cuboctahedral bimetallic nanoparticle with a minimum adsorbate distance of 2 Angstrom:
>>> from acat.adsorption_sites import ClusterAdsorptionSites >>> from acat.build.overlayer import SystematicPatternGenerator as SPG >>> from ase.cluster import Octahedron >>> from ase.io import read >>> from ase.visualize import view >>> atoms = Octahedron('Cu', length=7, cutoff=3) >>> for atom in atoms: ... if atom.index % 2 == 0: ... atom.symbol = 'Au' >>> atoms.center(vacuum=5.) >>> cas = ClusterAdsorptionSites(atoms, composition_effect=True) >>> spg = SPG(atoms, adsorbate_species='CO', ... min_adsorbate_distance=2., ... adsorption_sites=cas, ... composition_effect=True) >>> spg.run(action='add') >>> images = read('patterns.traj', index=':') >>> view(images)Output:
![]()
The following example illustrates how to enumerate all unique coverage paterns consists of 3 adsorbates choosen from C, N and O on a bimetallic bcc111 surface slab with a minimum adsorbate distance of 2 Angstrom (here only generate a maximum of 100 unique patterns):
>>> from acat.build.overlayer import SystematicPatternGenerator as SPG >>> from acat.adsorption_sites import SlabAdsorptionSites >>> from ase.io import read >>> from ase.build import bcc111 >>> from ase.visualize import view >>> atoms = bcc111('Fe', (2, 2, 12), vacuum=5.) >>> for atom in atoms: ... if atom.index % 2 == 0: ... atom.symbol = 'Mo' >>> atoms.center() >>> sas = SlabAdsorptionSites(atoms, surface='bcc111', ... composition_effect=True) >>> spg = SPG(atoms, adsorbate_species=['C','N','O'], ... min_adsorbate_distance=2., ... adsorption_sites=sas, ... surface='bcc111', ... composition_effect=True) >>> spg.run(max_gen_per_image=100, action='add', num_act=3) >>> images = read('patterns.traj', index=':') >>> view(images)Output:
![]()
ordered_coverage_pattern function¶
acat.build.overlayer.
ordered_coverage_pattern
(atoms, adsorbate, coverage=1.0, surface=None, height=None, min_adsorbate_distance=0.0, both_sides=False)[source]¶A function for generating representative ordered adsorbate overlayer patterns. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
atoms (ase.Atoms object) – The nanoparticle or surface slab onto which the adsorbates are added. Accept any ase.Atoms object. No need to be built-in.
adsorbate (str or ase.Atom object or ase.Atoms object) – The adsorbate species to be added onto the surface. For now only support adding one type of adsorbate species.
coverage (float, default 1.) – The coverage (ML) of the adsorbate (N_adsorbate / N_surf_atoms). Support 4 overlayer patterns (0.25 for p(2x2) pattern; 0.5 for c(2x2) pattern on fcc100 or honeycomb pattern on fcc111; 0.75 for (2x2) pattern on fcc100 or Kagome pattern on fcc111; 1 for p(1x1) pattern. Note that for small nanoparticles, the function might give results that do not correspond to the coverage. This is normal since the surface area can be too small to encompass the overlayer pattern properly. We expect this function to work well on large nanoparticles and surface slabs.
surface (str, default None) – The surface type (crystal structure + Miller indices). For now only support 2 common surfaces: fcc100 and fcc111. If the structure is a periodic surface slab, this is required. If the structure is a nanoparticle, the function only add adsorbates to the sites on the specified surface.
height (float, default None) – The height of the added adsorbate from the surface. Use the default settings if not specified.
min_adsorbate_distance (float, default 0.) – The minimum distance between two atoms that belongs to two adsorbates.
both_sides (bool, default False) – Whether to add adsorbate to both top and bottom sides of the slab. Only relvevant for periodic surface slabs.
Example
To add a 0.5 ML CO overlayer pattern on a cuboctahedron:
>>> from acat.build.overlayer import ordered_coverage_pattern >>> from ase.cluster import Octahedron >>> from ase.visualize import view >>> atoms = Octahedron('Au', length=9, cutoff=4) >>> atoms.center(vacuum=5.) >>> pattern = ordered_coverage_pattern(atoms, adsorbate='CO', ... coverage=0.5) >>> view(pattern)Output:
![]()
To add a 0.75 ML CO overlayer pattern on a fcc111 surface slab:
>>> from acat.build.overlayer import ordered_coverage_pattern >>> from ase.build import fcc111 >>> from ase.visualize import view >>> atoms = fcc111('Cu', (8, 8, 4), vacuum=5.) >>> atoms.center() >>> pattern = ordered_coverage_pattern(atoms, adsorbate='CO', ... coverage=0.5, ... surface='fcc111') >>> view(pattern)Output:
![]()
full_coverage_pattern function¶
acat.build.overlayer.
full_coverage_pattern
(atoms, adsorbate, site, surface=None, height=None, min_adsorbate_distance=0.0, both_sides=False)[source]¶A function for generating adsorbate overlayer patterns on sites that are of a same type. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
atoms (ase.Atoms object) – The nanoparticle or surface slab onto which the adsorbates are added. Accept any ase.Atoms object. No need to be built-in.
adsorbate (str or ase.Atom object or ase.Atoms object) – The adsorbate species to be added onto the surface. For now only support adding one type of adsorbate species.
site (str) – The site type that the adsorbates should be added to.
surface (str, default None) – The surface type (crystal structure + Miller indices). If the structure is a periodic surface slab, this is required. If the structure is a nanoparticle, the function only add adsorbates to the sites on the specified surface.
height (float, default None) – The height of the added adsorbate from the surface. Use the default settings if not specified.
min_adsorbate_distance (float, default 0.) – The minimum distance between two atoms that belongs to two adsorbates.
both_sides (bool, default False) – Whether to add adsorbate to both top and bottom sides of the slab. Only relvevant for periodic surface slabs.
Example
To add CO to all hcp sites on a icosahedron:
>>> from acat.build.overlayer import full_coverage_pattern >>> from ase.cluster import Icosahedron >>> from ase.visualize import view >>> atoms = Icosahedron('Au', noshells=5) >>> atoms.center(vacuum=5.) >>> pattern = full_coverage_pattern(atoms, adsorbate='CO', site='hcp') >>> view(pattern)Output:
![]()
To add CO to all 3fold sites on a bcc110 surface slab:
>>> from acat.build.overlayer import full_coverage_pattern >>> from ase.build import bcc110 >>> from ase.visualize import view >>> atoms = bcc110('Mo', (8, 8, 4), vacuum=5.) >>> atoms.center() >>> pattern = full_coverage_pattern(atoms, adsorbate='CO', ... surface='bcc110', site='3fold') >>> view(pattern)Output:
![]()
random_coverage_pattern function¶
acat.build.overlayer.
random_coverage_pattern
(atoms, adsorbate_species, species_probabilities=None, surface=None, min_adsorbate_distance=1.5, heights={'3fold': 1.3, '4fold': 1.3, '5fold': 1.5, '6fold': 0.0, 'bridge': 1.5, 'fcc': 1.3, 'hcp': 1.3, 'longbridge': 1.5, 'ontop': 1.8, 'shortbridge': 1.5}, allow_6fold=False, both_sides=False)[source]¶A function for generating random overlayer patterns with a minimum distance constraint. The function is generalized for both periodic and non-periodic systems (distinguished by atoms.pbc).
- Parameters
atoms (ase.Atoms object) – The nanoparticle or surface slab onto which the adsorbates are added. Accept any ase.Atoms object. No need to be built-in.
adsorbate_species (str or list of strs) – A list of adsorbate species to be randomly added to the surface.
species_probabilities (dict, default None) – A dictionary that contains keys of each adsorbate species and values of their probabilities of adding onto the surface.
surface (str, default None) – The surface type (crystal structure + Miller indices). If the structure is a periodic surface slab, this is required. If the structure is a nanoparticle, the function only add adsorbates to the sites on the specified surface.
heights (dict, default acat.settings.site_heights) – A dictionary that contains the adsorbate height for each site type. Use the default height settings if the height for a site type is not specified.
min_adsorbate_distance (float, default 1.5) – The minimum distance constraint between two atoms that belongs to two adsorbates.
allow_6fold (bool, default False) – Whether to allow the adsorption on 6-fold subsurf sites underneath fcc hollow sites.
both_sides (bool, default False) – Whether to add adsorbate to both top and bottom sides of the slab. Only relvevant for periodic surface slabs.
Example
To add CO randomly onto a cuboctahedron with a minimum adsorbate distance of 5 Angstrom:
>>> from acat.build.overlayer import random_coverage_pattern >>> from ase.cluster import Octahedron >>> from ase.visualize import view >>> atoms = Octahedron('Au', length=9, cutoff=4) >>> atoms.center(vacuum=5.) >>> pattern = random_coverage_pattern(atoms, adsorbate_species='CO', ... min_adsorbate_distance=5.) >>> view(pattern)Output:
![]()
To add C, N, O randomly onto a hcp0001 surface slab with probabilities of 0.25, 0.25, 0.5, respectively, and a minimum adsorbate distance of 2 Angstrom:
>>> from acat.build.overlayer import random_coverage_pattern >>> from ase.build import hcp0001 >>> from ase.visualize import view >>> atoms = hcp0001('Ru', (8, 8, 4), vacuum=5.) >>> atoms.center() >>> pattern = random_coverage_pattern(atoms, adsorbate_species=['C','N','O'], ... species_probabilities={'C': 0.25, ... 'N': 0.25, ... 'O': 0.5}, ... surface='hcp0001', ... min_adsorbate_distance=2.) >>> view(pattern)Output:
![]()
To add H, OH and H2O randomly onto a fcc100 Ni2Cu surface slab on both top and bottom interfaces (
bulk water
in between) with probabilities of 0.25, 0.25, 0.5, respectively, and a minimum adsorbate distance of 2 Angstrom:>>> from acat.build.overlayer import random_coverage_pattern >>> from ase.io import read >>> from ase.build import fcc100 >>> from ase.visualize import view >>> water_bulk = read('water_bulk.xyz') >>> water_bulk.center(vacuum=11., axis=2) >>> slab = fcc100('Ni', (4, 4, 8), vacuum=9.5, periodic=True) >>> for atom in slab: ... if atom.index % 3 == 0: ... atom.symbol = 'Cu' >>> slab.translate(-slab.cell[2] / 2) >>> slab.wrap() >>> atoms = slab + water_bulk >>> pattern = random_coverage_pattern(atoms, adsorbate_species=['H','OH','OH2'], ... species_probabilities={'H': 0.25, ... 'OH': 0.25, ... 'OH2': 0.5}, ... surface='fcc100', ... min_adsorbate_distance=2., ... both_sides=True) >>> view(pattern)Output:
![]()
Generate alloy chemical orderings¶
-
class
acat.build.ordering.
OrderedSlabOrderingGenerator
(atoms, elements, reducing_size=2, 2, composition=None, dtol=0.01, ztol=0.1, trajectory='orderings.traj', append_trajectory=False)[source]¶ Bases:
object
OrderedSlabOrderingGenerator is a class for generating ordered chemical orderings for a alloy surface slab. There is no limitation of the number of metal components.
- Parameters
atoms (ase.Atoms object) – The surface slab to use as a template to generate ordered chemical orderings. Accept any ase.Atoms object. No need to be built-in.
elements (list of strs) – The metal elements of the alloy catalyst.
reducing_size (list of ints or tuple of ints, default (2, 2)) – The multiples that groups the symmetry-equivalent atoms in the x and y directions. The x or y length of the cell must be this multiple of the distance between each pair of symmetry-equivalent atoms. Larger reducing size generates less structures.
composition (dict, default None) – Generate ordered orderings only at a certain composition. The dictionary contains the metal elements as keys and their concentrations as values. Generate orderings at all compositions if not specified. Note that the computational cost scales badly with the number of groups for a fixed-composition search.
dtol (float, default 0.01) – The distance tolerance (in Angstrom) when comparing with (cell length / multiple). Use a larger value if the structure is distorted.
ztol (float, default 0.1) – The tolerance (in Angstrom) when comparing z values. Use a larger ztol if the structure is distorted.
trajectory (str, default 'orderings.traj') – The name of the output ase trajectory file.
append_trajectory (bool, default False) – Whether to append structures to the existing trajectory.
-
get_groups
()[source]¶ Get the groups (a list of lists of atom indices) of all symmetry-equivalent atoms.
-
run
(max_gen=None, mode='systematic', verbose=False)[source]¶ Run the chemical ordering generator.
- Parameters
max_gen (int, default None) – Maximum number of chemical orderings to generate. Enumerate all symetric patterns if not specified.
mode (str, default 'systematic') –
‘systematic’: enumerate all possible unique chemical orderings. Recommended when there are not many groups. Switch to stochastic mode automatically if the number of groups is more than 20. This mode is the only option when the composition is fixed.
’stochastic’: sample chemical orderings stochastically. Duplicate structures can be generated. Recommended when there are many groups. Switch to systematic mode automatically if the composition is fixed.
verbose (bool, default False) – Whether to print out information about number of groups and number of generated structures.
SymmetricClusterOrderingGenerator class¶
- class
acat.build.ordering.
SymmetricClusterOrderingGenerator
(atoms, elements, symmetry='spherical', cutoff=1.0, secondary_symmetry=None, secondary_cutoff=1.0, composition=None, trajectory='orderings.traj', append_trajectory=False)[source]¶SymmetricClusterOrderingGenerator is a class for generating symmetric chemical orderings for a nanoalloy. There is no limitation of the number of metal components. Please always align the z direction to the symmetry axis of the nanocluster.
- Parameters
atoms (ase.Atoms object) – The nanoparticle to use as a template to generate symmetric chemical orderings. Accept any ase.Atoms object. No need to be built-in.
elements (list of strs) – The metal elements of the nanoalloy.
symmetry (str, default 'spherical') –
Support 9 symmetries:
’spherical’: centrosymmetry (groups defined by the distances to the geometric center);
’cylindrical’: cylindrical symmetry around z axis (groups defined by the distances to the z axis);
’planar’: planar symmetry around z axis (groups defined by the z coordinates), common for phase-separated nanoalloys;
’mirror_planar’: mirror planar symmetry around both z and xy plane (groups defined by the absolute z coordinate), high symmetry subset of ‘planar’; ‘circular’ = ring symmetry around z axis (groups defined by both z coordinate and distance to z axis);
’mirror_circular’: mirror ring symmetry around both z and xy plane (groups defined by both absolute z coordinate and distance to z axis);
’chemical’: symmetry w.r.t chemical environment (groups defined by the atomic energies given by a Gupta potential)
’geometrical’: symmetry w.r.t geometrical environment (groups defined by vertex / edge / fcc111 / fcc100 / bulk identified by CNA analysis);
’concentric’: conventional definition of the concentric shells (surface / subsurface, subsubsurface, …, core).
cutoff (float, default 1.0) – Maximum thickness (in Angstrom) of a single group. The thickness is calculated as the difference between the “distances” of the closest-to-center atoms in two neighbor groups. Note that the criterion of “distance” depends on the symmetry. This parameter works differently if the symmetry is ‘chemical’, ‘geometrical’ or ‘concentric’. For ‘chemical’ it is defined as the maximum atomic energy difference (in eV) of a single group predicted by a Gupta potential. For ‘geometrical’ and ‘concentric’ it is defined as the cutoff radius (in Angstrom) for CNA, and a reasonable cutoff based on the lattice constant of the material will automatically be used if cutoff <= 1. Use a larger cutoff if the structure is distorted.
secondary_symmetry (str, default None) – Add a secondary symmetry check to define groups hierarchically. For example, even if two atoms are classifed in the same group defined by the primary symmetry, they can still end up in different groups if they fall into two different groups defined by the secondary symmetry. Support 7 symmetries: ‘spherical’, ‘cylindrical’, ‘planar’, ‘mirror_planar’, ‘chemical’, ‘geometrical’ and ‘concentric’. Note that secondary symmetry has the same importance as the primary symmetry, so you can set either of the two symmetries of interest as the secondary symmetry. Useful for combining symmetries of different types (e.g. circular + chemical) or combining symmetries with different cutoffs.
secondary_cutoff (float, default 1.0) – Same as cutoff, except that it is for the secondary symmetry.
composition (dict, default None) – Generate symmetric orderings only at a certain composition. The dictionary contains the metal elements as keys and their concentrations as values. Generate orderings at all compositions if not specified. Note that the computational cost scales badly with the number of groups for a fixed-composition search.
trajectory (str, default 'orderings.traj') – The name of the output ase trajectory file.
append_trajectory (bool, default False) – Whether to append structures to the existing trajectory.
get_sorted_indices
(symmetry)[source]¶Returns the indices sorted by the metric that defines different groups, together with the corresponding vlues, given a specific symmetry. Returns the indices sorted by geometrical environment if symmetry=’geometrical’. Returns the indices sorted by surface, subsurface, subsubsurface, …, core if symmetry=’concentric’.
- Parameters
symmetry (str) – Support 7 symmetries: spherical, cylindrical, planar, mirror_planar, chemical, geometrical, concentric.
get_groups
()[source]¶Get the groups (a list of lists of atom indices) of all symmetry-equivalent atoms.
run
(max_gen=None, mode='systematic', verbose=False)[source]¶Run the chemical ordering generator.
- Parameters
max_gen (int, default None) – Maximum number of chemical orderings to generate. Enumerate all symetric patterns if not specified.
mode (str, default 'systematic') –
‘systematic’: enumerate all possible unique chemical orderings. Recommended when there are not many groups. Switch to stochastic mode automatically if the number of groups is more than 20. This mode is the only option when the composition is fixed.
’stochastic’: sample chemical orderings stochastically. Duplicate structures can be generated. Recommended when there are many groups. Switch to systematic mode automatically if the composition is fixed.
verbose (bool, default False) – Whether to print out information about number of groups and number of generated structures.
Example
To generate 100 symmetric chemical orderings for truncated octahedral NiPt nanoalloys with spherical symmetry:
>>> from acat.build.ordering import SymmetricClusterOrderingGenerator as SCOG >>> from ase.cluster import Octahedron >>> from ase.io import read >>> from ase.visualize import view >>> atoms = Octahedron('Ni', length=8, cutoff=3) >>> sog = SCOG(atoms, elements=['Ni', 'Pt'], symmetry='spherical') >>> sog.run(max_gen=100, verbose=True) >>> images = read('orderings.traj', index=':') >>> view(images)Output:
10 groups classified100 symmetric chemical orderings generated![]()
To generate 50 symmetric chemical orderings for quaternary truncated octahedral Ni0.4Cu0.3Pt0.2Au0.1 nanoalloys with mirror circular symmetry:
>>> from acat.build.ordering import SymmetricClusterOrderingGenerator as SCOG >>> from ase.cluster import Octahedron >>> from ase.io import read >>> from ase.visualize import view >>> atoms = Octahedron('Ni', 7, 2) >>> sog = SCOG(atoms, elements=['Ni', 'Cu', 'Pt', 'Au'], ... symmetry='mirror_circular', ... composition={'Ni': 0.4, 'Cu': 0.3, 'Pt': 0.2, 'Au': 0.1}) >>> sog.run(max_gen=50, verbose=True) >>> images = read('orderings.traj', index=':') >>> view(images)Output:
25 groups classified50 symmetric chemical orderings generated![]()
Sometimes it is also useful to get the structure of each group. For instance, to visualize the concentric shells of a truncated octahedral nanoparticle:
>>> from acat.build.ordering import SymmetricClusterOrderingGenerator as SCOG >>> from ase.cluster import Octahedron >>> from ase.visualize import view >>> atoms = Octahedron('Ni', 12, 5) >>> sog = SCOG(atoms, elements=['Ni', 'Pt'], symmetry='concentric') >>> groups = sog.get_groups() >>> images = [atoms[g] for g in groups] >>> view(images)Output:
![]()
OrderedSlabOrderingGenerator class¶
- class
acat.build.ordering.
OrderedSlabOrderingGenerator
(atoms, elements, reducing_size=2, 2, composition=None, dtol=0.01, ztol=0.1, trajectory='orderings.traj', append_trajectory=False)[source]¶OrderedSlabOrderingGenerator is a class for generating ordered chemical orderings for a alloy surface slab. There is no limitation of the number of metal components.
- Parameters
atoms (ase.Atoms object) – The surface slab to use as a template to generate ordered chemical orderings. Accept any ase.Atoms object. No need to be built-in.
elements (list of strs) – The metal elements of the alloy catalyst.
reducing_size (list of ints or tuple of ints, default (2, 2)) – The multiples that groups the symmetry-equivalent atoms in the x and y directions. The x or y length of the cell must be this multiple of the distance between each pair of symmetry-equivalent atoms. Larger reducing size generates less structures.
composition (dict, default None) – Generate ordered orderings only at a certain composition. The dictionary contains the metal elements as keys and their concentrations as values. Generate orderings at all compositions if not specified. Note that the computational cost scales badly with the number of groups for a fixed-composition search.
dtol (float, default 0.01) – The distance tolerance (in Angstrom) when comparing with (cell length / multiple). Use a larger value if the structure is distorted.
ztol (float, default 0.1) – The tolerance (in Angstrom) when comparing z values. Use a larger ztol if the structure is distorted.
trajectory (str, default 'orderings.traj') – The name of the output ase trajectory file.
append_trajectory (bool, default False) – Whether to append structures to the existing trajectory.
get_groups
()[source]¶Get the groups (a list of lists of atom indices) of all symmetry-equivalent atoms.
run
(max_gen=None, mode='systematic', verbose=False)[source]¶Run the chemical ordering generator.
- Parameters
max_gen (int, default None) – Maximum number of chemical orderings to generate. Enumerate all symetric patterns if not specified.
mode (str, default 'systematic') –
‘systematic’: enumerate all possible unique chemical orderings. Recommended when there are not many groups. Switch to stochastic mode automatically if the number of groups is more than 20. This mode is the only option when the composition is fixed.
’stochastic’: sample chemical orderings stochastically. Duplicate structures can be generated. Recommended when there are many groups. Switch to systematic mode automatically if the composition is fixed.
verbose (bool, default False) – Whether to print out information about number of groups and number of generated structures.
Example
To stochastically generate 50 ordered chemical orderings for ternary NixPtyAu1-x-y fcc111 surface slabs:
>>> from acat.build.ordering import OrderedSlabOrderingGenerator as OSOG >>> from ase.build import fcc111 >>> from ase.io import read >>> from ase.visualize import view >>> atoms = fcc111('Ni', (6, 6, 4), vacuum=5.) >>> atoms.center() >>> osog = OSOG(atoms, elements=['Ni', 'Pt', 'Au'], ... reducing_size=(3, 3)) >>> osog.run(max_gen=50, mode='stochastic', verbose=True) >>> images = read('orderings.traj', index=':') >>> view(images)Output:
16 groups classified50 ordered chemical orderings generated![]()
To generate 50 ordered chemical orderings for Ni0.75Pt0.25 fcc110 surface slabs:
>>> from acat.build.ordering import OrderedSlabOrderingGenerator as OSOG >>> from ase.build import fcc110 >>> from ase.io import read >>> from ase.visualize import view >>> atoms = fcc110('Ni', (4, 4, 4), vacuum=5.) >>> atoms.center() >>> osog = OSOG(atoms, elements=['Ni', 'Pt'], ... composition={'Ni': 0.75, 'Pt': 0.25}, ... reducing_size=(2, 2)) >>> osog.run(max_gen=50, verbose=True) >>> images = read('orderings.traj', index=':') >>> view(images)Output:
16 groups classified50 ordered chemical orderings generated![]()
RandomOrderingGenerator class¶
- class
acat.build.ordering.
RandomOrderingGenerator
(atoms, elements, composition=None, trajectory='orderings.traj', append_trajectory=False)[source]¶RandomOrderingGenerator is a class for generating random chemical orderings for an alloy catalyst. The function is generalized for both periodic and non-periodic systems, and there is no limitation of the number of metal components.
- Parameters
atoms (ase.Atoms object) – The nanoparticle or surface slab to use as a template to generate random chemical orderings. Accept any ase.Atoms object. No need to be built-in.
elements (list of strs) – The metal elements of the alloy catalyst.
composition (dict, default None) – Generate random orderings only at a certain composition. The dictionary contains the metal elements as keys and their concentrations as values. Generate orderings at all compositions if not specified.
trajectory (str, default 'orderings.traj') – The name of the output ase trajectory file.
append_trajectory (bool, default False) – Whether to append structures to the existing trajectory.
Example
To generate 50 random chemical orderings for icosahedral Ni0.5Pt0.2Au0.3 nanoalloys:
>>> from acat.build.ordering import RandomOrderingGenerator as ROG >>> from ase.cluster import Icosahedron >>> from ase.io import read >>> from ase.visualize import view >>> atoms = Icosahedron('Ni', noshells=4) >>> rog = ROG(atoms, elements=['Ni', 'Pt', 'Au'], ... composition={'Ni': 0.5, 'Pt': 0.2, 'Au': 0.3}) >>> rog.run(num_gen=50) >>> images = read('orderings.traj', index=':') >>> view(images)Output:
![]()
To generate 50 random chemical orderings for Pt0.5Au0.5 fcc111 surface slabs:
>>> from acat.build.ordering import RandomOrderingGenerator as ROG >>> from ase.build import fcc111 >>> from ase.io import read >>> from ase.visualize import view >>> atoms = fcc111('Ni', (4, 4, 4), vacuum=5.) >>> atoms.center() >>> rog = ROG(atoms, elements=['Pt', 'Au'], ... composition={'Pt': 0.5, 'Au': 0.5}) >>> rog.run(num_gen=50) >>> images = read('orderings.traj', index=':') >>> view(images)Output:
![]()