def id_adsorbates(self, classifier='radial', return_atoms=False): """Return a list of Gratoms objects for each adsorbate classified on a surface. Requires classification of adsorbate atoms. Parameters ---------- classifier : str Classification technique to identify individual adsorbates. 'radial': Use standard cutoff distances to identify neighboring atoms. return_atoms : bool Return Gratoms objects instead of adsorbate indices. Returns ------- adsorabtes : list (n,) Adsorbate indices of adsorbates in unit cell. """ atoms = self.atoms.copy() # Remove the slab atoms ads_atoms = self.ads_atoms if ads_atoms is None: ads_atoms = self.id_adsorbate_atoms() if classifier == 'radial': con = utils.get_cutoff_neighbors(atoms) ads_con = con[ads_atoms][:, ads_atoms] G = nx.Graph() G.add_nodes_from(ads_atoms) edges = utils.connectivity_to_edges(ads_con, indices=ads_atoms) G.add_weighted_edges_from(edges, weight='bonds') SG = nx.connected_component_subgraphs(G) adsorabtes = [] for sg in SG: nodes = list(sg.nodes) if return_atoms: edges = list(sg.edges) ads = Gratoms(numbers=atoms.numbers[nodes], positions=atoms.positions[nodes], edges=edges) ads.center(vacuum=5) else: ads = nodes adsorabtes += [ads] return adsorabtes
def get_slab(self, size=(1, 1), root=None, iterm=None, primitive=False): """Generate a slab object with a certain number of layers. Parameters ---------- size : tuple (2,) Repeat the x and y lattice vectors by the indicated dimensions root : int Produce a slab with a primitive a1 basis vector multiplied by the square root of a provided value. Uses primitive unit cell. iterm : int A termination index in reference to the list of possible terminations. primitive : bool Whether to reduce the unit cell to its primitive form. Returns ------- slab : atoms object The modified basis slab produced based on the layer specifications given. """ slab = self._basis.copy() if iterm: if self.unique_terminations is None: terminations = self.get_unique_terminations() else: terminations = self.unique_terminations zshift = terminations[iterm] slab.translate([0, 0, -zshift]) slab.wrap(pbc=True) # Get the minimum number of layers needed zlayers = utils.get_unique_coordinates(slab, direct=False, tol=self.tol) if self.min_width: width = slab.cell[2][2] z_repetitions = np.ceil(width / len(zlayers) * self.min_width) else: z_repetitions = np.ceil(self.layers / len(zlayers)) slab *= (1, 1, int(z_repetitions)) if primitive or root: if self.vacuum: slab.center(vacuum=self.vacuum, axis=2) else: raise ValueError('Primitive slab generation requires vacuum') nslab = utils.get_primitive_cell(slab) if nslab is not None: slab = nslab # spglib occasionally returns a split slab zpos = slab.get_scaled_positions() if zpos[:, 2].max() > 0.9 or zpos[:, 2].min() < 0.1: zpos[:, 2] -= 0.5 zpos[:, 2] %= 1 slab.set_scaled_positions(zpos) slab.center(vacuum=self.vacuum, axis=2) # For hcp(1, 1, 0), primitive alters z-axis d = norm(slab.cell, axis=0) maxd = np.argwhere(d == d.max())[0][0] if maxd != 2: slab.rotate(slab.cell[maxd], 'z', rotate_cell=True) slab.cell[[maxd, 2]] = slab.cell[[2, maxd]] slab.cell[maxd] = -slab.cell[maxd] slab.wrap(pbc=True) slab.rotate(slab.cell[0], 'x', rotate_cell=True) # Orthogonalize the z-coordinate # Warning: bulk symmetry is lost at this point a, b, c = slab.cell nab = np.cross(a, b) c = (nab * np.dot(c, nab) / norm(nab)**2) slab.cell[2] = c # Align the longest remaining basis vector with x vdist = norm(slab.cell[:2], axis=1) if vdist[1] > vdist[0]: slab.rotate(slab.cell[1], 'x', rotate_cell=True) slab.cell[0] *= -1 slab.cell[[0, 1]] = slab.cell[[1, 0]] # Enforce that the angle between basis vectors is acute. if slab.cell[1][0] < 0: slab.rotate(slab.cell[2], '-z') slab.cell *= [[1, 0, 0], [-1, 1, 0], [0, 0, 1]] if root: roots, vectors = root_surface_analysis(slab, return_vectors=True) if root not in roots: raise ValueError( 'Requested root structure unavailable for this system.' 'Try: {}'.format(roots)) vect = vectors[np.where(root == roots)][0] slab = self.root_surface(slab, root, vect) # Get the direct z-coordinate of the requested layer zlayers = utils.get_unique_coordinates(slab, direct=False, tag=True, tol=self.tol) if not self.fix_stoichiometry: reverse_sort = np.sort(zlayers)[::-1] if self.min_width: n = np.where(zlayers < self.min_width, 1, 0).sum() ncut = reverse_sort[n] else: ncut = reverse_sort[:self.layers][-1] zpos = slab.positions[:, 2] index = np.arange(len(slab)) del slab[index[zpos - ncut < -self.tol]] slab.cell[2][2] -= ncut slab.translate([0, 0, -ncut]) slab *= (size[0], size[1], 1) tags = slab.get_tags() m = np.where(tags == 1)[0][0] translation = slab[m].position.copy() translation[2] = 0 slab.translate(-translation) slab.wrap() ind = np.lexsort( (slab.positions[:, 0], slab.positions[:, 1], slab.positions[:, 2])) slab = Gratoms(positions=slab.positions[ind], numbers=slab.numbers[ind], cell=slab.cell, pbc=[1, 1, 0], tags=tags[ind]) fix = tags.max() - self.fixed constraints = FixAtoms(indices=[a.index for a in slab if a.tag > fix]) slab.set_constraint(constraints) if self.vacuum: slab.center(vacuum=self.vacuum, axis=2) self.slab = slab return slab