Statistiques
| Branche: | Tag: | Révision :

dockonsurf / modules / screening.py @ a44ad3c2

Historique | Voir | Annoter | Télécharger (50,63 ko)

1 e07c09eb Carles
import logging
2 f3d1e601 Carles Marti
import numpy as np
3 f3d1e601 Carles Marti
import ase
4 e07c09eb Carles
5 e07c09eb Carles
logger = logging.getLogger('DockOnSurf')
6 e07c09eb Carles
7 e07c09eb Carles
8 bfe93f0d Carles Marti
def assign_prop(atoms: ase.Atoms, prop_name: str, prop_val):  # TODO Needed?
9 bfe93f0d Carles Marti
    atoms.info[prop_name] = prop_val
10 bfe93f0d Carles Marti
11 bfe93f0d Carles Marti
12 e1c5f171 Carles Marti
def select_confs(orig_conf_list: list, calc_dirs: list, magns: list,
13 e1c5f171 Carles Marti
                 num_sel: int, code: str):
14 bfe93f0d Carles Marti
    """Takes a list ase.Atoms and selects the most different magnitude-wise.
15 bfe93f0d Carles Marti

16 bfe93f0d Carles Marti
    Given a list of ase.Atoms objects and a list of magnitudes, it selects a
17 bfe93f0d Carles Marti
    number of the most different conformers according to every magnitude
18 bfe93f0d Carles Marti
    specified.
19 1e9e784d Carles Marti
     
20 bfe93f0d Carles Marti
    @param orig_conf_list: list of ase.Atoms objects to select among.
21 1e9e784d Carles Marti
    @param calc_dirs: List of directories where to read the energies from.
22 bfe93f0d Carles Marti
    @param magns: list of str with the names of the magnitudes to use for the
23 1e9e784d Carles Marti
        conformer selection. Supported magnitudes: 'energy', 'moi'.
24 bfe93f0d Carles Marti
    @param num_sel: number of conformers to select for every of the magnitudes.
25 bfe93f0d Carles Marti
    @param code: The code that generated the magnitude information.
26 bfe93f0d Carles Marti
         Supported codes: See formats.py
27 bfe93f0d Carles Marti
    @return: list of the selected ase.Atoms objects.
28 bfe93f0d Carles Marti
    """
29 bfe93f0d Carles Marti
    from copy import deepcopy
30 fd2384fc Carles Marti
    from modules.formats import collect_energies
31 bfe93f0d Carles Marti
32 bfe93f0d Carles Marti
    conf_list = deepcopy(orig_conf_list)
33 e1c5f171 Carles Marti
    conf_enrgs, mois, selected_ids = [], [], []
34 bfe93f0d Carles Marti
    if num_sel >= len(conf_list):
35 bfe93f0d Carles Marti
        logger.warning('Number of conformers per magnitude is equal or larger '
36 bfe93f0d Carles Marti
                       'than the total number of conformers. Using all '
37 695dcff8 Carles Marti
                       f'available conformers: {len(conf_list)}.')
38 bfe93f0d Carles Marti
        return conf_list
39 bfe93f0d Carles Marti
40 bfe93f0d Carles Marti
    # Read properties
41 bfe93f0d Carles Marti
    if 'energy' in magns:
42 cdc1edbe Carles Martí
        if code == 'cp2k':
43 cdc1edbe Carles Martí
            conf_enrgs = collect_energies(calc_dirs, code, 'isolated')
44 cdc1edbe Carles Martí
        elif code == 'vasp':
45 cdc1edbe Carles Martí
            conf_enrgs = np.array([conf.get_total_energy()
46 cdc1edbe Carles Martí
                                   for conf in orig_conf_list])
47 bfe93f0d Carles Marti
    if 'moi' in magns:
48 bfe93f0d Carles Marti
        mois = np.array([conf.get_moments_of_inertia() for conf in conf_list])
49 bfe93f0d Carles Marti
50 bfe93f0d Carles Marti
    # Assign values
51 bfe93f0d Carles Marti
    for i, conf in enumerate(conf_list):
52 bfe93f0d Carles Marti
        assign_prop(conf, 'idx', i)
53 bb387578 Carles Marti
        assign_prop(conf, 'iso', calc_dirs[i])
54 bfe93f0d Carles Marti
        if 'energy' in magns:
55 bfe93f0d Carles Marti
            assign_prop(conf, 'energy', conf_enrgs[i])
56 bfe93f0d Carles Marti
        if 'moi' in magns:
57 bfe93f0d Carles Marti
            assign_prop(conf, 'moi', mois[i, 2])
58 bfe93f0d Carles Marti
59 bfe93f0d Carles Marti
    # pick ids
60 bfe93f0d Carles Marti
    for magn in magns:
61 bfe93f0d Carles Marti
        sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn]))
62 bfe93f0d Carles Marti
        if sorted_list[-1].info['idx'] not in selected_ids:
63 bfe93f0d Carles Marti
            selected_ids.append(sorted_list[-1].info['idx'])
64 bfe93f0d Carles Marti
        if num_sel > 1:
65 bfe93f0d Carles Marti
            for i in range(0, len(sorted_list) - 1,
66 bfe93f0d Carles Marti
                           len(conf_list) // (num_sel - 1)):
67 bfe93f0d Carles Marti
                if sorted_list[i].info['idx'] not in selected_ids:
68 bfe93f0d Carles Marti
                    selected_ids.append(sorted_list[i].info['idx'])
69 bfe93f0d Carles Marti
70 695dcff8 Carles Marti
    logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
71 bfe93f0d Carles Marti
    return [conf_list[idx] for idx in selected_ids]
72 bfe93f0d Carles Marti
73 bfe93f0d Carles Marti
74 7dd94df7 Carles Marti
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True):
75 b4aef3d7 Carles Marti
    """Computes the angle between two vectors.
76 b4aef3d7 Carles Marti

77 b4aef3d7 Carles Marti
    @param v1: The first vector.
78 b4aef3d7 Carles Marti
    @param v2: The second vector.
79 b516a1d6 Carles Marti
    @param ref: Orthogonal vector to both v1 and v2,
80 b516a1d6 Carles Marti
        along which the sign of the rotation is defined (i.e. positive if
81 b516a1d6 Carles Marti
        counterclockwise angle when facing ref)
82 b4aef3d7 Carles Marti
    @param degrees: Whether the result should be in radians (True) or in
83 b4aef3d7 Carles Marti
        degrees (False).
84 b4aef3d7 Carles Marti
    @return: The angle in radians if degrees = False, or in degrees if
85 b4aef3d7 Carles Marti
        degrees =True
86 b4aef3d7 Carles Marti
    """
87 b4aef3d7 Carles Marti
    v1_u = v1 / np.linalg.norm(v1)
88 b4aef3d7 Carles Marti
    v2_u = v2 / np.linalg.norm(v2)
89 b4aef3d7 Carles Marti
    angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0))
90 b516a1d6 Carles Marti
    if ref is not None:
91 b516a1d6 Carles Marti
        # Give sign according to ref direction
92 b516a1d6 Carles Marti
        angle *= (1 if np.dot(np.cross(v1, v2), ref) >= 0 else -1)
93 b516a1d6 Carles Marti
94 b4aef3d7 Carles Marti
    return angle if not degrees else angle * 180 / np.pi
95 b4aef3d7 Carles Marti
96 b4aef3d7 Carles Marti
97 12a12bdd Carles Marti
def vect_avg(vects):
98 a5cc42ff Carles Marti
    """Computes the element-wise mean of a set of vectors.
99 12a12bdd Carles Marti

100 dadc6016 Carles Marti
    @param vects: list of lists-like: containing the vectors (num_vectors,
101 dadc6016 Carles Marti
        length_vector).
102 a5cc42ff Carles Marti
    @return: vector average computed doing the element-wise mean.
103 12a12bdd Carles Marti
    """
104 d26885d3 Carles Marti
    from modules.utilities import try_command
105 12a12bdd Carles Marti
    err = "vect_avg parameter vects must be a list-like, able to be converted" \
106 12a12bdd Carles Marti
          " np.array"
107 dadc6016 Carles Marti
    array = try_command(np.array, [(ValueError, err)], vects)
108 dadc6016 Carles Marti
    if len(array.shape) == 1:
109 dadc6016 Carles Marti
        return array
110 12a12bdd Carles Marti
    else:
111 dadc6016 Carles Marti
        num_vects = array.shape[1]
112 dadc6016 Carles Marti
        return np.array([np.average(array[:, i]) for i in range(num_vects)])
113 12a12bdd Carles Marti
114 12a12bdd Carles Marti
115 8fdff302 Carles Marti
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None):
116 8fdff302 Carles Marti
    """Gets the coordinates of the specified indices from a ase.Atoms object.
117 a5cc42ff Carles Marti

118 a5cc42ff Carles Marti
    Given an ase.Atoms object and a list of atom indices specified in ctrs_list
119 a5cc42ff Carles Marti
    it gets the coordinates of the specified atoms. If the element in the
120 a5cc42ff Carles Marti
    ctrs_list is not an index but yet a list of indices, it computes the
121 a5cc42ff Carles Marti
    element-wise mean of the coordinates of the atoms specified in the inner
122 a5cc42ff Carles Marti
    list.
123 a5cc42ff Carles Marti
    @param atoms: ase.Atoms object for which to obtain the coordinates of.
124 a5cc42ff Carles Marti
    @param ctrs_list: list of (indices/list of indices) of the atoms for which
125 a5cc42ff Carles Marti
                      the coordinates should be extracted.
126 a5cc42ff Carles Marti
    @return: np.ndarray of atomic coordinates.
127 a5cc42ff Carles Marti
    """
128 12a12bdd Carles Marti
    coords = []
129 8fdff302 Carles Marti
    err = "'ctrs_list' argument must be an integer, a list of integers or a " \
130 8fdff302 Carles Marti
          "list of lists of integers. Every integer must be in the range " \
131 8fdff302 Carles Marti
          "[0, num_atoms)"
132 8fdff302 Carles Marti
    if ctrs_list is None:
133 8fdff302 Carles Marti
        ctrs_list = range(len(atoms))
134 8fdff302 Carles Marti
    elif isinstance(ctrs_list, int):
135 8fdff302 Carles Marti
        if ctrs_list not in range(len(atoms)):
136 8fdff302 Carles Marti
            logger.error(err)
137 8fdff302 Carles Marti
            raise ValueError(err)
138 8fdff302 Carles Marti
        return atoms[ctrs_list].position
139 8fdff302 Carles Marti
    for elem in ctrs_list:
140 12a12bdd Carles Marti
        if isinstance(elem, list):
141 8fdff302 Carles Marti
            coords.append(vect_avg([atoms[c].position for c in elem]))
142 12a12bdd Carles Marti
        elif isinstance(elem, int):
143 12a12bdd Carles Marti
            coords.append(atoms[elem].position)
144 12a12bdd Carles Marti
        else:
145 12a12bdd Carles Marti
            logger.error(err)
146 12a12bdd Carles Marti
            raise ValueError
147 12a12bdd Carles Marti
    return np.array(coords)
148 f3d1e601 Carles Marti
149 f3d1e601 Carles Marti
150 d6da8693 Carles Marti
def compute_norm_vect(atoms, idxs, cell):
151 d6da8693 Carles Marti
    """Computes the local normal vector of a surface at a given site.
152 d6da8693 Carles Marti

153 d6da8693 Carles Marti
    Given an ase.Atoms object and a site defined as a linear combination of
154 d6da8693 Carles Marti
    atoms it computes the vector perpendicular to the surface, considering the
155 d6da8693 Carles Marti
    local environment of the site.
156 d6da8693 Carles Marti
    @param atoms: ase.Atoms object of the surface.
157 d6da8693 Carles Marti
    @param idxs: list or int: Index or list of indices of the atom/s that define
158 d6da8693 Carles Marti
        the site
159 d6da8693 Carles Marti
    @param cell: Unit cell. A 3x3 matrix (the three unit cell vectors)
160 d6da8693 Carles Marti
    @return: numpy.ndarray of the coordinates of the vector locally
161 d6da8693 Carles Marti
    perpendicular to the surface.
162 d6da8693 Carles Marti
    """
163 6d1bd296 Carles Marti
    from modules.ASANN import coordination_numbers as coord_nums
164 d6da8693 Carles Marti
    if isinstance(idxs, list):
165 d6da8693 Carles Marti
        atm_vect = [-np.round(coord_nums(atoms.get_scaled_positions(),
166 d6da8693 Carles Marti
                                         pbc=np.any(cell),
167 d6da8693 Carles Marti
                                         cell_vectors=cell)[3][i], 2)
168 d6da8693 Carles Marti
                    for i in idxs]
169 d6da8693 Carles Marti
        norm_vec = vect_avg([vect / np.linalg.norm(vect) for vect in atm_vect])
170 d6da8693 Carles Marti
    elif isinstance(idxs, int):
171 d6da8693 Carles Marti
        norm_vec = -coord_nums(atoms.get_scaled_positions(),
172 d6da8693 Carles Marti
                               pbc=np.any(cell),
173 d6da8693 Carles Marti
                               cell_vectors=cell)[3][idxs]
174 d6da8693 Carles Marti
    else:
175 d6da8693 Carles Marti
        err = "'idxs' must be either an int or a list"
176 d6da8693 Carles Marti
        logger.error(err)
177 d6da8693 Carles Marti
        raise ValueError(err)
178 587dca22 Carles Marti
    norm_vec = np.round(norm_vec, 2) / np.linalg.norm(np.round(norm_vec, 2))
179 587dca22 Carles Marti
    logger.info(f"The perpendicular vector to the surface at site '{idxs}' is "
180 587dca22 Carles Marti
                f"{norm_vec}")
181 587dca22 Carles Marti
    return norm_vec
182 d6da8693 Carles Marti
183 d6da8693 Carles Marti
184 d8d92cfb Carles Marti
def align_molec(orig_molec, ctr_coord, ref_vect):
185 4918e2ad Carles Marti
    """Align a molecule to a vector by a center.
186 4918e2ad Carles Marti

187 4918e2ad Carles Marti
    Given a reference vector to be aligned to and some coordinates acting as
188 4918e2ad Carles Marti
    alignment center, it first averages the vectors pointing to neighboring
189 4918e2ad Carles Marti
    atoms and then tries to align this average vector to the target. If the
190 4918e2ad Carles Marti
    average vector is 0 it takes the vector to the nearest neighbor.
191 5412d6ef Carles Marti
    @param orig_molec: The molecule to align.
192 4918e2ad Carles Marti
    @param ctr_coord: The coordinates to use ase alignment center.
193 4918e2ad Carles Marti
    @param ref_vect: The vector to be aligned with.
194 5412d6ef Carles Marti
    @return: ase.Atoms of the aligned molecule.
195 4918e2ad Carles Marti
    """
196 d8d92cfb Carles Marti
    from copy import deepcopy
197 4918e2ad Carles Marti
    from ase import Atom
198 4918e2ad Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
199 4918e2ad Carles Marti
200 d8d92cfb Carles Marti
    molec = deepcopy(orig_molec)
201 4918e2ad Carles Marti
    if len(molec) == 1:
202 4918e2ad Carles Marti
        err_msg = "Cannot align a single atom"
203 4918e2ad Carles Marti
        logger.error(err_msg)
204 4918e2ad Carles Marti
        ValueError(err_msg)
205 4918e2ad Carles Marti
    cutoffs = natural_cutoffs(molec, mult=1.2)
206 4918e2ad Carles Marti
    # Check if ctr_coord are the coordinates of an atom and if not creates a
207 4918e2ad Carles Marti
    # dummy one to extract the neighboring atoms.
208 4918e2ad Carles Marti
    ctr_idx = None
209 f16ebc80 Carles Marti
    dummy_atom = False
210 4918e2ad Carles Marti
    for atom in molec:
211 4918e2ad Carles Marti
        if np.allclose(ctr_coord, atom.position, rtol=1e-2):
212 4918e2ad Carles Marti
            ctr_idx = atom.index
213 4918e2ad Carles Marti
            break
214 4918e2ad Carles Marti
    if ctr_idx is None:
215 4918e2ad Carles Marti
        molec.append(Atom('X', position=ctr_coord))
216 4918e2ad Carles Marti
        cutoffs.append(0.2)
217 4918e2ad Carles Marti
        ctr_idx = len(molec) - 1
218 f16ebc80 Carles Marti
        dummy_atom = True
219 4918e2ad Carles Marti
    # Builds the neighbors and computes the average vector
220 4918e2ad Carles Marti
    refs, vects = neighbor_list("iD", molec, cutoffs, self_interaction=False)
221 4918e2ad Carles Marti
    neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
222 4918e2ad Carles Marti
    # If no neighbors are present, the cutoff of the alignment center is
223 4918e2ad Carles Marti
    # set to a value where at least one atom is a neighbor and neighbors are
224 4918e2ad Carles Marti
    # recalculated.
225 4918e2ad Carles Marti
    if len(neigh_vects) == 0:
226 4918e2ad Carles Marti
        min_dist, min_idx = (np.inf, np.inf)
227 4918e2ad Carles Marti
        for atom in molec:
228 4918e2ad Carles Marti
            if atom.index == ctr_idx:
229 4918e2ad Carles Marti
                continue
230 e5faf87a Carles Marti
            if molec.get_distance(ctr_idx, atom.index) < min_dist:
231 4918e2ad Carles Marti
                min_dist = molec.get_distance(ctr_idx, atom.index)
232 4918e2ad Carles Marti
                min_idx = atom.index
233 4918e2ad Carles Marti
        cutoffs[ctr_idx] = min_dist - cutoffs[min_idx] + 0.05
234 4918e2ad Carles Marti
        refs, vects = neighbor_list("iD", molec, cutoffs,
235 4918e2ad Carles Marti
                                    self_interaction=False)
236 4918e2ad Carles Marti
        neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
237 4918e2ad Carles Marti
    target_vect = vect_avg(neigh_vects)
238 e5faf87a Carles Marti
    # If the target vector is 0 (the center is at the baricenter of its
239 e5faf87a Carles Marti
    # neighbors). Assuming the adsorption center is coplanar or colinear to its
240 e5faf87a Carles Marti
    # neighbors (it would not make a lot of sense to chose a center which is
241 e5faf87a Carles Marti
    # the baricenter of neighbors distributed in 3D), the target_vector is
242 e5faf87a Carles Marti
    # chosen perpendicular to the nearest neighbor.
243 e5faf87a Carles Marti
    if np.allclose(target_vect, 0, rtol=1e-3):
244 e5faf87a Carles Marti
        nn_vect = np.array([np.inf] * 3)
245 4918e2ad Carles Marti
        for vect in neigh_vects:
246 e5faf87a Carles Marti
            if np.linalg.norm(vect) < np.linalg.norm(nn_vect):
247 e5faf87a Carles Marti
                nn_vect = vect
248 e5faf87a Carles Marti
        cart_axes = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]
249 e5faf87a Carles Marti
        axis = cart_axes[int(np.argmax([np.linalg.norm(np.cross(ax, nn_vect))
250 e5faf87a Carles Marti
                                        for ax in cart_axes]))]
251 e5faf87a Carles Marti
        target_vect = np.cross(axis, nn_vect)
252 4918e2ad Carles Marti
253 4918e2ad Carles Marti
    rot_vect = np.cross(target_vect, ref_vect)
254 4918e2ad Carles Marti
    if np.allclose(rot_vect, 0):
255 4918e2ad Carles Marti
        cart_axes = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]
256 4918e2ad Carles Marti
        axis = cart_axes[int(np.argmax([np.linalg.norm(np.cross(ax, ref_vect))
257 4918e2ad Carles Marti
                                        for ax in cart_axes]))]
258 4918e2ad Carles Marti
        rot_vect = np.cross(ref_vect, axis)
259 4918e2ad Carles Marti
    rot_angle = -get_vect_angle(ref_vect, target_vect, rot_vect)
260 4918e2ad Carles Marti
    molec.rotate(rot_angle, rot_vect, ctr_coord)
261 f16ebc80 Carles Marti
    if dummy_atom:
262 f16ebc80 Carles Marti
        del molec[-1]
263 4918e2ad Carles Marti
    return molec
264 4918e2ad Carles Marti
265 4918e2ad Carles Marti
266 dadc6016 Carles Marti
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
267 1d22a086 Carles Marti
                  norm_vect=(0, 0, 1)):
268 1d22a086 Carles Marti
    """Add an adsorbate to a surface.
269 1d22a086 Carles Marti

270 1d22a086 Carles Marti
    This function extends the functionality of ase.build.add_adsorbate
271 1d22a086 Carles Marti
    (https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
272 1d22a086 Carles Marti
    by enabling to change the z coordinate and the axis perpendicular to the
273 1d22a086 Carles Marti
    surface.
274 1d22a086 Carles Marti
    @param slab: ase.Atoms object containing the coordinates of the surface
275 1d22a086 Carles Marti
    @param adsorbate: ase.Atoms object containing the coordinates of the
276 1d22a086 Carles Marti
        adsorbate.
277 dadc6016 Carles Marti
    @param site_coord: The coordinates of the adsorption site on the surface.
278 dadc6016 Carles Marti
    @param ctr_coord: The coordinates of the adsorption center in the molecule.
279 dadc6016 Carles Marti
    @param height: The height above the surface where to adsorb.
280 1d22a086 Carles Marti
    @param offset: Offsets the adsorbate by a number of unit cells. Mostly
281 1d22a086 Carles Marti
        useful when adding more than one adsorbate.
282 1d22a086 Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
283 1d22a086 Carles Marti
    """
284 36d92f4f Carles Marti
    from copy import deepcopy
285 1d22a086 Carles Marti
    info = slab.info.get('adsorbate_info', {})
286 dadc6016 Carles Marti
    pos = np.array([0.0, 0.0, 0.0])  # part of absolute coordinates
287 1d22a086 Carles Marti
    spos = np.array([0.0, 0.0, 0.0])  # part relative to unit cell
288 dadc6016 Carles Marti
    norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect)
289 1d22a086 Carles Marti
    if offset is not None:
290 1d22a086 Carles Marti
        spos += np.asarray(offset, float)
291 dadc6016 Carles Marti
    if isinstance(site_coord, str):
292 1d22a086 Carles Marti
        # A site-name:
293 1d22a086 Carles Marti
        if 'sites' not in info:
294 1d22a086 Carles Marti
            raise TypeError('If the atoms are not made by an ase.build '
295 1d22a086 Carles Marti
                            'function, position cannot be a name.')
296 dadc6016 Carles Marti
        if site_coord not in info['sites']:
297 dadc6016 Carles Marti
            raise TypeError('Adsorption site %s not supported.' % site_coord)
298 dadc6016 Carles Marti
        spos += info['sites'][site_coord]
299 1d22a086 Carles Marti
    else:
300 dadc6016 Carles Marti
        pos += site_coord
301 1d22a086 Carles Marti
    if 'cell' in info:
302 1d22a086 Carles Marti
        cell = info['cell']
303 1d22a086 Carles Marti
    else:
304 1d22a086 Carles Marti
        cell = slab.get_cell()
305 1d22a086 Carles Marti
    pos += np.dot(spos, cell)
306 1d22a086 Carles Marti
    # Convert the adsorbate to an Atoms object
307 1d22a086 Carles Marti
    if isinstance(adsorbate, ase.Atoms):
308 36d92f4f Carles Marti
        ads = deepcopy(adsorbate)
309 1d22a086 Carles Marti
    elif isinstance(adsorbate, ase.Atom):
310 1d22a086 Carles Marti
        ads = ase.Atoms([adsorbate])
311 1d22a086 Carles Marti
    else:
312 1d22a086 Carles Marti
        # Assume it is a string representing a single Atom
313 1d22a086 Carles Marti
        ads = ase.Atoms([ase.Atom(adsorbate)])
314 dadc6016 Carles Marti
    pos += height * norm_vect_u
315 1d22a086 Carles Marti
    # Move adsorbate into position
316 dadc6016 Carles Marti
    ads.translate(pos - ctr_coord)
317 1d22a086 Carles Marti
    # Attach the adsorbate
318 1d22a086 Carles Marti
    slab.extend(ads)
319 1d22a086 Carles Marti
320 1d22a086 Carles Marti
321 a4b57124 Carles Marti
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
322 9cd032cf Carles Marti
                    nn_molec=0, coll_coeff=1.2, exclude_atom=False):
323 5f3f4b69 Carles Marti
    """Checks whether a slab and a molecule collide or not.
324 5f3f4b69 Carles Marti

325 5f3f4b69 Carles Marti
    @param slab_molec: The system of adsorbate-slab for which to detect if there
326 5f3f4b69 Carles Marti
        are collisions.
327 5f3f4b69 Carles Marti
    @param nn_slab: Number of neigbors in the surface.
328 5f3f4b69 Carles Marti
    @param nn_molec: Number of neighbors in the molecule.
329 5f3f4b69 Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
330 5f3f4b69 Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
331 5f3f4b69 Carles Marti
        considered as atomic collision.
332 5f3f4b69 Carles Marti
    @param slab_num_atoms: Number of atoms of the bare slab.
333 5f3f4b69 Carles Marti
    @param min_height: The minimum height atoms can have to not be considered as
334 5f3f4b69 Carles Marti
        colliding.
335 5f3f4b69 Carles Marti
    @param vect: The vector perpendicular to the slab.
336 9cd032cf Carles Marti
    @param exclude_atom: Whether to exclude the adsorption center in the
337 9cd032cf Carles Marti
        molecule in the collision detection.
338 5f3f4b69 Carles Marti
    @return: bool, whether the surface and the molecule collide.
339 e8bebcca Carles Marti
    """
340 9cd032cf Carles Marti
    from copy import deepcopy
341 5f3f4b69 Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
342 e8bebcca Carles Marti
343 e8bebcca Carles Marti
    # Check structure overlap by height
344 5fb01677 Carles Marti
    if min_height is not False:
345 a4b57124 Carles Marti
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
346 a4b57124 Carles Marti
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
347 e8bebcca Carles Marti
        if vect.tolist() not in cart_axes:
348 e8bebcca Carles Marti
            err_msg = "'min_coll_height' option is only implemented for " \
349 e8bebcca Carles Marti
                      "'surf_norm_vect' to be one of the x, y or z axes. "
350 e8bebcca Carles Marti
            logger.error(err_msg)
351 e8bebcca Carles Marti
            raise ValueError(err_msg)
352 e8bebcca Carles Marti
        for atom in slab_molec[slab_num_atoms:]:
353 9cd032cf Carles Marti
            if exclude_atom is not False \
354 9cd032cf Carles Marti
                    and atom.index == exclude_atom:
355 9cd032cf Carles Marti
                continue
356 e8bebcca Carles Marti
            for i, coord in enumerate(vect):
357 e8bebcca Carles Marti
                if coord == 0:
358 e8bebcca Carles Marti
                    continue
359 e8bebcca Carles Marti
                if atom.position[i] * coord < min_height * coord:
360 e8bebcca Carles Marti
                    return True
361 e8bebcca Carles Marti
362 e8bebcca Carles Marti
    # Check structure overlap by sphere collision
363 e8bebcca Carles Marti
    if coll_coeff is not False:
364 9cd032cf Carles Marti
        if exclude_atom is not False:
365 9cd032cf Carles Marti
            slab_molec_wo_ctr = deepcopy(slab_molec)
366 9cd032cf Carles Marti
            del slab_molec_wo_ctr[exclude_atom + slab_num_atoms]
367 9cd032cf Carles Marti
            slab_molec_cutoffs = natural_cutoffs(slab_molec_wo_ctr,
368 9cd032cf Carles Marti
                                                 mult=coll_coeff)
369 9cd032cf Carles Marti
            slab_molec_nghbs = len(neighbor_list("i", slab_molec_wo_ctr,
370 9cd032cf Carles Marti
                                                 slab_molec_cutoffs))
371 9cd032cf Carles Marti
        else:
372 9cd032cf Carles Marti
            slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
373 9cd032cf Carles Marti
            slab_molec_nghbs = len(neighbor_list("i", slab_molec,
374 9cd032cf Carles Marti
                                                 slab_molec_cutoffs))
375 5fb01677 Carles Marti
        if slab_molec_nghbs > nn_slab + nn_molec:
376 5fb01677 Carles Marti
            return True
377 e8bebcca Carles Marti
378 e8bebcca Carles Marti
    return False
379 e8bebcca Carles Marti
380 e8bebcca Carles Marti
381 e8bebcca Carles Marti
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
382 e8bebcca Carles Marti
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
383 9cd032cf Carles Marti
                 coll_coeff, height=2.5, excl_atom=False):
384 9cd032cf Carles Marti
    # TODO Rename this function
385 e8bebcca Carles Marti
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
386 e8bebcca Carles Marti
    small rotations.
387 e8bebcca Carles Marti

388 e8bebcca Carles Marti
    @param molec: ase.Atoms object of the molecule to adsorb
389 e8bebcca Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the
390 e8bebcca Carles Marti
        molecule
391 e8bebcca Carles Marti
    @param ctr_coord: The coordinates of the molecule to use as adsorption
392 e8bebcca Carles Marti
        center.
393 e8bebcca Carles Marti
    @param site_coord: The coordinates of the surface on which to adsorb the
394 e8bebcca Carles Marti
        molecule
395 e8bebcca Carles Marti
    @param num_pts: Number on which to sample Euler angles.
396 e8bebcca Carles Marti
    @param min_coll_height: The lowermost height for which to detect a collision
397 e8bebcca Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
398 e8bebcca Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
399 e8bebcca Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
400 e8bebcca Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
401 e8bebcca Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
402 e8bebcca Carles Marti
        considered as atomic collision.
403 fe91ddb2 Carles Marti
    @param height: Height on which to try adsorption.
404 9cd032cf Carles Marti
    @param excl_atom: Whether to exclude the adsorption center in the
405 9cd032cf Carles Marti
        molecule in the collision detection.
406 e8bebcca Carles Marti
    @return collision: bool, whether the structure generated has collisions
407 e8bebcca Carles Marti
        between slab and adsorbate.
408 e8bebcca Carles Marti
    """
409 e8bebcca Carles Marti
    from copy import deepcopy
410 e8bebcca Carles Marti
    slab_num_atoms = len(slab)
411 e8bebcca Carles Marti
    slab_molec = []
412 e8bebcca Carles Marti
    collision = True
413 e8bebcca Carles Marti
    max_corr = 6  # Should be an even number
414 e8bebcca Carles Marti
    d_angle = 180 / ((max_corr / 2.0) * num_pts)
415 e8bebcca Carles Marti
    num_corr = 0
416 e8bebcca Carles Marti
    while collision and num_corr <= max_corr:
417 e8bebcca Carles Marti
        k = num_corr * (-1) ** num_corr
418 e8bebcca Carles Marti
        slab_molec = deepcopy(slab)
419 e8bebcca Carles Marti
        molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
420 e8bebcca Carles Marti
                           center=ctr_coord)
421 e8bebcca Carles Marti
        add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height,
422 e8bebcca Carles Marti
                      norm_vect=norm_vect)
423 e8bebcca Carles Marti
        collision = check_collision(slab_molec, slab_num_atoms, min_coll_height,
424 e8bebcca Carles Marti
                                    norm_vect, slab_nghbs, molec_nghbs,
425 9cd032cf Carles Marti
                                    coll_coeff, excl_atom)
426 e8bebcca Carles Marti
        num_corr += 1
427 e8bebcca Carles Marti
    return slab_molec, collision
428 5f3f4b69 Carles Marti
429 5f3f4b69 Carles Marti
430 c25aa299 Carles Marti
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
431 c25aa299 Carles Marti
                 neigh_cutoff=1):
432 b4b2f307 Carles Marti
    # TODO rethink
433 91ae8d86 Carles Marti
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
434 b4b2f307 Carles Marti

435 91ae8d86 Carles Marti
    Tries to dissociate a H atom from the molecule and adsorb in on top of the
436 91ae8d86 Carles Marti
    surface if the distance is shorter than two times the neigh_cutoff value.
437 b4b2f307 Carles Marti
    @param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
438 b4b2f307 Carles Marti
    @param h_idx: The index of the hydrogen atom to carry out adsorption of.
439 b4b2f307 Carles Marti
    @param num_atoms_slab: The number of atoms of the slab without adsorbate.
440 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
441 b4b2f307 Carles Marti
    @param neigh_cutoff: half the maximum distance between the surface and the
442 b4b2f307 Carles Marti
        H for it to carry out dissociation.
443 b4b2f307 Carles Marti
    @return: An ase.Atoms object of the system adsorbate-surface with H
444 b4b2f307 Carles Marti
    """
445 b4b2f307 Carles Marti
    from copy import deepcopy
446 b4b2f307 Carles Marti
    from ase.neighborlist import NeighborList
447 b4b2f307 Carles Marti
    slab_molec = deepcopy(slab_molec_orig)
448 b4b2f307 Carles Marti
    cutoffs = len(slab_molec) * [neigh_cutoff]
449 c25aa299 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True, skin=0.0)
450 b4b2f307 Carles Marti
    nl.update(slab_molec)
451 b4b2f307 Carles Marti
    surf_h_vect = np.array([np.infty] * 3)
452 c25aa299 Carles Marti
    if h_acceptor == 'all':
453 c25aa299 Carles Marti
        h_acceptor = list(range(num_atoms_slab))
454 b4b2f307 Carles Marti
    for neigh_idx in nl.get_neighbors(h_idx)[0]:
455 c25aa299 Carles Marti
        for elm in h_acceptor:
456 c25aa299 Carles Marti
            if isinstance(elm, int):
457 c25aa299 Carles Marti
                if neigh_idx == elm and neigh_idx < num_atoms_slab:
458 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
459 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
460 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
461 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
462 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
463 c25aa299 Carles Marti
            else:
464 c25aa299 Carles Marti
                if slab_molec[neigh_idx].symbol == elm \
465 c25aa299 Carles Marti
                        and neigh_idx < num_atoms_slab:
466 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
467 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
468 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
469 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
470 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
471 c25aa299 Carles Marti
472 b4b2f307 Carles Marti
    if np.linalg.norm(surf_h_vect) != np.infty:
473 b4b2f307 Carles Marti
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
474 b4b2f307 Carles Marti
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
475 b4b2f307 Carles Marti
        return slab_molec
476 b4b2f307 Carles Marti
477 b4b2f307 Carles Marti
478 c25aa299 Carles Marti
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
479 b4b2f307 Carles Marti
    # TODO multiple dissociation
480 b4b2f307 Carles Marti
    """Decides which H atoms to dissociate according to a list of atoms.
481 b4b2f307 Carles Marti

482 b4b2f307 Carles Marti
    Given a list of chemical symbols or atom indices it checks for every atom
483 b4b2f307 Carles Marti
    or any of its neighbor if it's a H and calls dissociate_h to try to carry
484 b4b2f307 Carles Marti
    out dissociation of that H. For atom indices, it checks both whether
485 b4b2f307 Carles Marti
    the atom index or its neighbors are H, for chemical symbols, it only checks
486 b4b2f307 Carles Marti
    if there is a neighbor H.
487 b4b2f307 Carles Marti
    @param slab_molec: The ase.Atoms object of the system adsorbate-slab.
488 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
489 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
490 c25aa299 Carles Marti
    @param num_atoms_slab: Number of atoms of the bare slab.
491 b4b2f307 Carles Marti
    @return:
492 b4b2f307 Carles Marti
    """
493 b4b2f307 Carles Marti
    from ase.neighborlist import natural_cutoffs, NeighborList
494 b4b2f307 Carles Marti
    molec = slab_molec[num_atoms_slab:]
495 b4b2f307 Carles Marti
    cutoffs = natural_cutoffs(molec)
496 b4b2f307 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True)
497 b4b2f307 Carles Marti
    nl.update(molec)
498 b4b2f307 Carles Marti
    disso_structs = []
499 c25aa299 Carles Marti
    for el in h_donor:
500 b4b2f307 Carles Marti
        if isinstance(el, int):
501 b4b2f307 Carles Marti
            if molec[el].symbol == 'H':
502 b4b2f307 Carles Marti
                disso_struct = dissociate_h(slab_molec, el + num_atoms_slab,
503 c25aa299 Carles Marti
                                            num_atoms_slab, h_acceptor)
504 b4b2f307 Carles Marti
                if disso_struct is not None:
505 b4b2f307 Carles Marti
                    disso_structs.append(disso_struct)
506 b4b2f307 Carles Marti
            else:
507 b4b2f307 Carles Marti
                for neigh_idx in nl.get_neighbors(el)[0]:
508 b4b2f307 Carles Marti
                    if molec[neigh_idx].symbol == 'H':
509 b4b2f307 Carles Marti
                        disso_struct = dissociate_h(slab_molec, neigh_idx +
510 b4b2f307 Carles Marti
                                                    num_atoms_slab,
511 c25aa299 Carles Marti
                                                    num_atoms_slab, h_acceptor)
512 b4b2f307 Carles Marti
                        if disso_struct is not None:
513 b4b2f307 Carles Marti
                            disso_structs.append(disso_struct)
514 b4b2f307 Carles Marti
        else:
515 b4b2f307 Carles Marti
            for atom in molec:
516 b4b2f307 Carles Marti
                if atom.symbol.lower() == el.lower():
517 b4b2f307 Carles Marti
                    for neigh_idx in nl.get_neighbors(atom.index)[0]:
518 b4b2f307 Carles Marti
                        if molec[neigh_idx].symbol == 'H':
519 5261a07f Carles Marti
                            disso_struct = dissociate_h(slab_molec, neigh_idx
520 b4b2f307 Carles Marti
                                                        + num_atoms_slab,
521 c25aa299 Carles Marti
                                                        num_atoms_slab,
522 c25aa299 Carles Marti
                                                        h_acceptor)
523 b4b2f307 Carles Marti
                            if disso_struct is not None:
524 b4b2f307 Carles Marti
                                disso_structs.append(disso_struct)
525 b4b2f307 Carles Marti
    return disso_structs
526 b4b2f307 Carles Marti
527 b4b2f307 Carles Marti
528 3ab0865c Carles Marti
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
529 b4b2f307 Carles Marti
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
530 9cd032cf Carles Marti
              h_donor, h_acceptor, height, excl_atom):
531 3ab0865c Carles Marti
    """Generates adsorbate-surface structures by sampling over Euler angles.
532 3ab0865c Carles Marti

533 3ab0865c Carles Marti
    This function generates a number of adsorbate-surface structures at
534 3ab0865c Carles Marti
    different orientations of the adsorbate sampled at multiple Euler (zxz)
535 3ab0865c Carles Marti
    angles.
536 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb.
537 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the
538 5261a07f Carles Marti
        molecule.
539 3ab0865c Carles Marti
    @param ctr_coord: The coordinates of the molecule to use as adsorption
540 3ab0865c Carles Marti
        center.
541 3ab0865c Carles Marti
    @param site_coord: The coordinates of the surface on which to adsorb the
542 5261a07f Carles Marti
        molecule.
543 3ab0865c Carles Marti
    @param num_pts: Number on which to sample Euler angles.
544 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision.
545 3ab0865c Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
546 3ab0865c Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
547 3ab0865c Carles Marti
        considered as atomic collision.
548 3ab0865c Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
549 3ab0865c Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
550 3ab0865c Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
551 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
552 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
553 fe91ddb2 Carles Marti
    @param height: Height on which to try adsorption.
554 9cd032cf Carles Marti
    @param excl_atom: Whether to exclude the adsorption center in the
555 9cd032cf Carles Marti
        molecule in the collision detection.
556 3ab0865c Carles Marti
    @return: list of ase.Atoms object conatining all the orientations of a given
557 5261a07f Carles Marti
        conformer.
558 3ab0865c Carles Marti
    """
559 3ab0865c Carles Marti
    from copy import deepcopy
560 cf40df1b Carles Marti
    slab_ads_list = []
561 d8d92cfb Carles Marti
    prealigned_molec = align_molec(orig_molec, ctr_coord, norm_vect)
562 3ab0865c Carles Marti
    # rotation around z
563 3ab0865c Carles Marti
    for alpha in np.arange(0, 360, 360 / num_pts):
564 3ab0865c Carles Marti
        # rotation around x'
565 3ab0865c Carles Marti
        for beta in np.arange(0, 180, 180 / num_pts):
566 3ab0865c Carles Marti
            # rotation around z'
567 3ab0865c Carles Marti
            for gamma in np.arange(0, 360, 360 / num_pts):
568 5864c86e Carles Marti
                if beta == 0 and gamma > 0:
569 5864c86e Carles Marti
                    continue
570 d8d92cfb Carles Marti
                molec = deepcopy(prealigned_molec)
571 3ab0865c Carles Marti
                molec.euler_rotate(alpha, beta, gamma, center=ctr_coord)
572 9cd032cf Carles Marti
                slab_molec, collision = correct_coll(molec, slab, ctr_coord,
573 9cd032cf Carles Marti
                                                     site_coord, num_pts,
574 9cd032cf Carles Marti
                                                     min_coll_height, norm_vect,
575 5fb01677 Carles Marti
                                                     slab_nghbs, molec_nghbs,
576 9cd032cf Carles Marti
                                                     coll_coeff, height,
577 9cd032cf Carles Marti
                                                     excl_atom)
578 cf40df1b Carles Marti
                if not collision and not any([np.allclose(slab_molec.positions,
579 cf40df1b Carles Marti
                                                          conf.positions)
580 cf40df1b Carles Marti
                                              for conf in slab_ads_list]):
581 cf40df1b Carles Marti
                    slab_ads_list.append(slab_molec)
582 cf40df1b Carles Marti
                    if h_donor is not False:
583 cf40df1b Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
584 cf40df1b Carles Marti
                                                          h_acceptor,
585 cf40df1b Carles Marti
                                                          len(slab)))
586 3ab0865c Carles Marti
587 cf40df1b Carles Marti
    return slab_ads_list
588 f3d1e601 Carles Marti
589 d68dd4ad Carles Marti
590 19727725 Carles Marti
def internal_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
591 19727725 Carles Marti
                    ctr2_surf, bond_vector, bond_angle_target,
592 19727725 Carles Marti
                    dihedral_angle_target=None, mol_dihedral_angle_target=None):
593 7dd94df7 Carles Marti
    """Performs translation and rotation of an adsorbate defined by an
594 7dd94df7 Carles Marti
    adsorption bond length, direction, angle and dihedral angle
595 7dd94df7 Carles Marti

596 7dd94df7 Carles Marti
    Carles modification of chemcat's transform_adsorbate to work with
597 7dd94df7 Carles Marti
    coordinates instead of ase.Atom
598 7dd94df7 Carles Marti
    Parameters:
599 7dd94df7 Carles Marti
        molecule (ase.Atoms): The molecule to adsorb.
600 7dd94df7 Carles Marti

601 7dd94df7 Carles Marti
        surf (ase.Atoms): The surface ontop of which to adsorb.
602 7dd94df7 Carles Marti

603 7dd94df7 Carles Marti
        ctr1_mol (int/list): The position of the adsorption center in the
604 7dd94df7 Carles Marti
        molecule that will be bound to the surface.
605 7dd94df7 Carles Marti

606 7dd94df7 Carles Marti
        ctr2_mol (int/list): The position of a second center of the
607 7dd94df7 Carles Marti
        adsorbate used to define the adsorption bond angle, and the dihedral
608 7dd94df7 Carles Marti
        adsorption angle.
609 7dd94df7 Carles Marti

610 7dd94df7 Carles Marti
        ctr3_mol (int/list): The position of a third center in the
611 7dd94df7 Carles Marti
        adsorbate used to define the adsorbate dihedral angle.
612 7dd94df7 Carles Marti

613 7dd94df7 Carles Marti
        ctr1_surf (int/list): The position of the site on the surface that
614 7dd94df7 Carles Marti
        will be bound to the molecule.
615 7dd94df7 Carles Marti

616 7dd94df7 Carles Marti
        ctr2_surf (int/list): The position of a second center of the
617 7dd94df7 Carles Marti
        surface used to define the dihedral adsorption angle.
618 7dd94df7 Carles Marti

619 7dd94df7 Carles Marti
        bond_vector (numpy.ndarray): The adsorption bond desired.
620 7dd94df7 Carles Marti
            Details: offset = vect(atom1_surf, atom1_mol)
621 7dd94df7 Carles Marti

622 7dd94df7 Carles Marti
        bond_angle_target (float or int): The adsorption bond angle desired (in
623 7dd94df7 Carles Marti
            degrees).
624 7dd94df7 Carles Marti
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
625 7dd94df7 Carles Marti

626 7dd94df7 Carles Marti
        dihedral_angle_target (float or int): The dihedral adsorption angle
627 7dd94df7 Carles Marti
            desired (in degrees).
628 7dd94df7 Carles Marti
            Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
629 7dd94df7 Carles Marti
            atom1_mol, atom2_mol)
630 7dd94df7 Carles Marti
                The rotation vector is facing the adsorbate from the surface
631 7dd94df7 Carles Marti
                (i.e. counterclockwise rotation when facing the surface (i.e.
632 7dd94df7 Carles Marti
                view from top))
633 7dd94df7 Carles Marti

634 7dd94df7 Carles Marti
        mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
635 7dd94df7 Carles Marti
            desired (in degrees).
636 7dd94df7 Carles Marti
            Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
637 7dd94df7 Carles Marti
            atom2_mol, atom3_mol)
638 7dd94df7 Carles Marti
                The rotation vector is facing atom2_mol from atom1_mol
639 7dd94df7 Carles Marti

640 7dd94df7 Carles Marti
    Returns:
641 7dd94df7 Carles Marti
        None (the `molecule` object is modified in-place)
642 7dd94df7 Carles Marti
    """
643 7dd94df7 Carles Marti
    vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf,
644 7dd94df7 Carles Marti
                                                                   ctr1_surf)
645 d6da8693 Carles Marti
    vect_inter = get_atom_coords(molecule, ctr1_mol) \
646 d6da8693 Carles Marti
        - get_atom_coords(surf, ctr1_surf)
647 7dd94df7 Carles Marti
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
648 7dd94df7 Carles Marti
                                                                     ctr1_mol)
649 7dd94df7 Carles Marti
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
650 7dd94df7 Carles Marti
                                                                      ctr2_mol)
651 7dd94df7 Carles Marti
652 7dd94df7 Carles Marti
    # Check if dihedral angles can be defined
653 7dd94df7 Carles Marti
    do_dihedral = dihedral_angle_target is not None
654 7dd94df7 Carles Marti
    do_mol_dihedral = mol_dihedral_angle_target is not None
655 7dd94df7 Carles Marti
    dihedral_use_mol2 = False
656 8c519241 Carles Marti
    # Check if vect_surf and bond_vector are aligned
657 8c519241 Carles Marti
    if np.allclose(np.cross(vect_surf, bond_vector), 0):
658 7dd94df7 Carles Marti
        do_dihedral = False
659 7dd94df7 Carles Marti
    # Check if requested bond angle is not flat
660 7dd94df7 Carles Marti
    if np.isclose((bond_angle_target + 90) % 180 - 90, 0):
661 7dd94df7 Carles Marti
        do_mol_dihedral = False
662 7dd94df7 Carles Marti
        dihedral_use_mol2 = True
663 7dd94df7 Carles Marti
    # Check if vect_mol and vect2_mol are not aligned
664 7dd94df7 Carles Marti
    if np.allclose(np.cross(vect_mol, vect2_mol), 0):
665 7dd94df7 Carles Marti
        do_mol_dihedral = False
666 7dd94df7 Carles Marti
667 7dd94df7 Carles Marti
    ###########################
668 7dd94df7 Carles Marti
    #       Translation       #
669 7dd94df7 Carles Marti
    ###########################
670 7dd94df7 Carles Marti
671 7dd94df7 Carles Marti
    # Compute and apply translation of adsorbate
672 7dd94df7 Carles Marti
    translation = bond_vector - vect_inter
673 7dd94df7 Carles Marti
    molecule.translate(translation)
674 7dd94df7 Carles Marti
675 7dd94df7 Carles Marti
    # Update adsorption bond
676 d6da8693 Carles Marti
    vect_inter = get_atom_coords(molecule, ctr1_mol) - \
677 d6da8693 Carles Marti
        get_atom_coords(surf, ctr1_surf)
678 7dd94df7 Carles Marti
679 7dd94df7 Carles Marti
    # Check if translation was successful
680 7dd94df7 Carles Marti
    if np.allclose(vect_inter, bond_vector):
681 7dd94df7 Carles Marti
        pass  # print("Translation successfully applied (error: ~ {:.5g} unit "
682 7dd94df7 Carles Marti
        # "length)".format(np.linalg.norm(vect_inter - bond_vector)))
683 7dd94df7 Carles Marti
    else:
684 7dd94df7 Carles Marti
        err = 'An unknown error occured during the translation'
685 7dd94df7 Carles Marti
        logger.error(err)
686 7dd94df7 Carles Marti
        raise AssertionError(err)
687 7dd94df7 Carles Marti
688 7dd94df7 Carles Marti
    ###########################
689 7dd94df7 Carles Marti
    #   Bond angle rotation   #
690 7dd94df7 Carles Marti
    ###########################
691 7dd94df7 Carles Marti
692 7dd94df7 Carles Marti
    # Compute rotation vector
693 7dd94df7 Carles Marti
    rotation_vector = np.cross(-vect_inter, vect_mol)
694 7dd94df7 Carles Marti
    if np.allclose(rotation_vector, 0, atol=1e-5):
695 7dd94df7 Carles Marti
        # If molecular bonds are aligned, any vector orthogonal to vect_inter
696 7dd94df7 Carles Marti
        # can be used Such vector can be found as the orthogonal rejection of
697 7dd94df7 Carles Marti
        # either X-axis, Y-axis or Z-axis onto vect_inter (since they cannot
698 7dd94df7 Carles Marti
        # be all aligned)
699 7dd94df7 Carles Marti
        non_aligned_vector = np.zeros(3)
700 7dd94df7 Carles Marti
        # Select the most orthogonal axis (lowest dot product):
701 7dd94df7 Carles Marti
        non_aligned_vector[np.argmin(np.abs(vect_inter))] = 1
702 7dd94df7 Carles Marti
        rotation_vector = non_aligned_vector - np.dot(non_aligned_vector,
703 7dd94df7 Carles Marti
                                                      vect_inter) / np.dot(
704 7dd94df7 Carles Marti
            vect_inter, vect_inter) * vect_inter
705 7dd94df7 Carles Marti
706 7dd94df7 Carles Marti
    # Retrieve current bond angle
707 7dd94df7 Carles Marti
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
708 7dd94df7 Carles Marti
709 7dd94df7 Carles Marti
    # Apply rotation to reach desired bond_angle
710 7dd94df7 Carles Marti
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
711 7dd94df7 Carles Marti
                    center=get_atom_coords(molecule, ctr1_mol))
712 7dd94df7 Carles Marti
713 7dd94df7 Carles Marti
    # Update molecular bonds
714 7dd94df7 Carles Marti
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
715 7dd94df7 Carles Marti
                                                                     ctr1_mol)
716 7dd94df7 Carles Marti
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
717 7dd94df7 Carles Marti
                                                                      ctr2_mol)
718 7dd94df7 Carles Marti
719 7dd94df7 Carles Marti
    # Check if rotation was successful
720 7dd94df7 Carles Marti
    bond_angle = get_vect_angle(-vect_inter, vect_mol)
721 7dd94df7 Carles Marti
    if np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90, 0,
722 7dd94df7 Carles Marti
                  atol=1e-3) and np.allclose(get_atom_coords(molecule, ctr1_mol)
723 7dd94df7 Carles Marti
                                             - get_atom_coords(surf,
724 7dd94df7 Carles Marti
                                                               ctr1_surf),
725 7dd94df7 Carles Marti
                                             vect_inter):
726 7dd94df7 Carles Marti
        pass  # print("Rotation successfully applied (error: {:.5f}°)".format(
727 7dd94df7 Carles Marti
        # (bond_angle - bond_angle_target + 90) % 180 - 90))
728 7dd94df7 Carles Marti
    else:
729 7dd94df7 Carles Marti
        err = 'An unknown error occured during the rotation'
730 7dd94df7 Carles Marti
        logger.error(err)
731 7dd94df7 Carles Marti
        raise AssertionError(err)
732 7dd94df7 Carles Marti
733 7dd94df7 Carles Marti
    ###########################
734 7dd94df7 Carles Marti
    # Dihedral angle rotation #
735 7dd94df7 Carles Marti
    ###########################
736 7dd94df7 Carles Marti
737 7dd94df7 Carles Marti
    # Perform dihedral rotation if possible
738 7dd94df7 Carles Marti
    if do_dihedral:
739 7dd94df7 Carles Marti
        # Retrieve current dihedral angle (by computing the angle between the
740 7dd94df7 Carles Marti
        # vector rejection of vect_surf and vect_mol onto vect_inter)
741 7dd94df7 Carles Marti
        vect_inter_inner = np.dot(vect_inter, vect_inter)
742 7dd94df7 Carles Marti
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
743 d6da8693 Carles Marti
            vect_inter_inner * vect_inter
744 7dd94df7 Carles Marti
        if dihedral_use_mol2:
745 7dd94df7 Carles Marti
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
746 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
747 7dd94df7 Carles Marti
        else:
748 7dd94df7 Carles Marti
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
749 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
750 7dd94df7 Carles Marti
        dihedral_angle_ini = get_vect_angle(vect_surf_reject, vect_mol_reject,
751 7dd94df7 Carles Marti
                                            vect_inter)
752 7dd94df7 Carles Marti
753 7dd94df7 Carles Marti
        # Apply dihedral rotation along vect_inter
754 7dd94df7 Carles Marti
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
755 7dd94df7 Carles Marti
                        v=vect_inter, center=get_atom_coords(molecule,
756 7dd94df7 Carles Marti
                                                             ctr1_mol))
757 7dd94df7 Carles Marti
758 7dd94df7 Carles Marti
        # Update molecular bonds
759 7dd94df7 Carles Marti
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
760 d6da8693 Carles Marti
            get_atom_coords(molecule, ctr1_mol)
761 7dd94df7 Carles Marti
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
762 d6da8693 Carles Marti
            get_atom_coords(molecule, ctr2_mol)
763 7dd94df7 Carles Marti
764 7dd94df7 Carles Marti
        # Check if rotation was successful
765 7dd94df7 Carles Marti
        # Check dihedral rotation
766 7dd94df7 Carles Marti
        if dihedral_use_mol2:
767 7dd94df7 Carles Marti
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
768 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
769 7dd94df7 Carles Marti
        else:
770 7dd94df7 Carles Marti
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
771 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
772 7dd94df7 Carles Marti
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
773 7dd94df7 Carles Marti
                                        vect_inter)
774 7dd94df7 Carles Marti
        # Check bond rotation is unmodified
775 7dd94df7 Carles Marti
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
776 7dd94df7 Carles Marti
        if np.isclose((dihedral_angle - dihedral_angle_target + 90) % 180 - 90,
777 5261a07f Carles Marti
                      0, atol=1e-3) \
778 5261a07f Carles Marti
                and np.isclose((bond_angle - bond_angle_target + 90) %
779 5261a07f Carles Marti
                               180 - 90, 0, atol=1e-5) \
780 c25aa299 Carles Marti
                and np.allclose(get_atom_coords(molecule, ctr1_mol)
781 c25aa299 Carles Marti
                                - get_atom_coords(surf, ctr1_surf),
782 c25aa299 Carles Marti
                                vect_inter):
783 7dd94df7 Carles Marti
            pass  # print( "Dihedral rotation successfully applied (error: {
784 7dd94df7 Carles Marti
            # :.5f}°)".format((dihedral_angle - dihedral_angle_target + 90) %
785 7dd94df7 Carles Marti
            # 180 - 90))
786 7dd94df7 Carles Marti
        else:
787 7dd94df7 Carles Marti
            err = 'An unknown error occured during the dihedral rotation'
788 7dd94df7 Carles Marti
            logger.error(err)
789 7dd94df7 Carles Marti
            raise AssertionError(err)
790 7dd94df7 Carles Marti
791 7dd94df7 Carles Marti
    #####################################
792 7dd94df7 Carles Marti
    # Adsorbate dihedral angle rotation #
793 7dd94df7 Carles Marti
    #####################################
794 7dd94df7 Carles Marti
795 7dd94df7 Carles Marti
    # Perform adsorbate dihedral rotation if possible
796 7dd94df7 Carles Marti
    if do_mol_dihedral:
797 7dd94df7 Carles Marti
        # Retrieve current adsorbate dihedral angle (by computing the angle
798 7dd94df7 Carles Marti
        # between the orthogonal rejection of vect_inter and vect2_mol onto
799 7dd94df7 Carles Marti
        # vect_mol)
800 7dd94df7 Carles Marti
        vect_mol_inner = np.dot(vect_mol, vect_mol)
801 7dd94df7 Carles Marti
        bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \
802 5261a07f Carles Marti
            vect_mol_inner * vect_mol
803 7dd94df7 Carles Marti
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
804 5261a07f Carles Marti
            vect_mol_inner * vect_mol
805 7dd94df7 Carles Marti
        dihedral_angle_ini = get_vect_angle(bond_inter_reject,
806 7dd94df7 Carles Marti
                                            bond2_mol_reject, vect_mol)
807 7dd94df7 Carles Marti
808 7dd94df7 Carles Marti
        # Apply dihedral rotation along vect_mol
809 7dd94df7 Carles Marti
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
810 7dd94df7 Carles Marti
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
811 7dd94df7 Carles Marti
812 7dd94df7 Carles Marti
        # Update molecular bonds
813 7dd94df7 Carles Marti
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
814 5261a07f Carles Marti
            - get_atom_coords(molecule, ctr1_mol)
815 7dd94df7 Carles Marti
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
816 5261a07f Carles Marti
            - get_atom_coords(molecule, ctr2_mol)
817 7dd94df7 Carles Marti
818 7dd94df7 Carles Marti
        # Check if rotation was successful
819 7dd94df7 Carles Marti
        # Check adsorbate dihedral rotation
820 7dd94df7 Carles Marti
        vect_mol_inner = np.dot(vect_mol, vect_mol)
821 7dd94df7 Carles Marti
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
822 5261a07f Carles Marti
            vect_mol_inner * vect_mol
823 7dd94df7 Carles Marti
        mol_dihedral_angle = get_vect_angle(bond_inter_reject,
824 7dd94df7 Carles Marti
                                            bond2_mol_reject, vect_mol)
825 7dd94df7 Carles Marti
        # Check dihedral rotation
826 7dd94df7 Carles Marti
        vect_inter_inner = np.dot(vect_inter, vect_inter)
827 7dd94df7 Carles Marti
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
828 5261a07f Carles Marti
            vect_inter_inner * vect_inter
829 7dd94df7 Carles Marti
        vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
830 5261a07f Carles Marti
            vect_inter_inner * vect_inter
831 7dd94df7 Carles Marti
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
832 7dd94df7 Carles Marti
                                        vect_inter)
833 7dd94df7 Carles Marti
        # Check bond rotation is unmodified
834 7dd94df7 Carles Marti
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
835 7dd94df7 Carles Marti
        if np.isclose((mol_dihedral_angle - mol_dihedral_angle_target + 90) %
836 7dd94df7 Carles Marti
                      180 - 90, 0, atol=1e-3) \
837 7dd94df7 Carles Marti
                and np.isclose((dihedral_angle -
838 7dd94df7 Carles Marti
                                dihedral_angle_target + 90) % 180 - 90, 0,
839 7dd94df7 Carles Marti
                               atol=1e-5) \
840 7dd94df7 Carles Marti
                and np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90,
841 7dd94df7 Carles Marti
                               0, atol=1e-5) \
842 7dd94df7 Carles Marti
                and np.allclose(get_atom_coords(molecule, ctr1_mol) -
843 7dd94df7 Carles Marti
                                get_atom_coords(surf, ctr1_surf),
844 7dd94df7 Carles Marti
                                vect_inter):
845 7dd94df7 Carles Marti
            pass  # print(
846 7dd94df7 Carles Marti
            # "Adsorbate dihedral rotation successfully applied (error:
847 7dd94df7 Carles Marti
            # {:.5f}°)".format((mol_dihedral_angle - mol_dihedral_angle_target
848 7dd94df7 Carles Marti
            # + 90) % 180 - 90))
849 7dd94df7 Carles Marti
        else:
850 7dd94df7 Carles Marti
            err = 'An unknown error occured during the adsorbate dihedral ' \
851 7dd94df7 Carles Marti
                  'rotation'
852 7dd94df7 Carles Marti
            logger.error(err)
853 7dd94df7 Carles Marti
            raise AssertionError(err)
854 7dd94df7 Carles Marti
855 7dd94df7 Carles Marti
856 19727725 Carles Marti
def ads_internal(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
857 19727725 Carles Marti
                 ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
858 9cd032cf Carles Marti
                 slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel, height,
859 9cd032cf Carles Marti
                 excl_atom):
860 609104e3 Carles Marti
    """Generates adsorbate-surface structures by sampling over internal angles.
861 5261a07f Carles Marti

862 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb (adsorbate).
863 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the molecule
864 5261a07f Carles Marti
    @param ctr1_mol: The index/es of the center in the adsorbate to use as
865 5261a07f Carles Marti
        adsorption center.
866 5261a07f Carles Marti
    @param ctr2_mol: The index/es of the center in the adsorbate to use for the
867 5261a07f Carles Marti
        definition of the surf-adsorbate angle, surf-adsorbate dihedral angle
868 5261a07f Carles Marti
        and adsorbate dihedral angle.
869 5261a07f Carles Marti
    @param ctr3_mol: The index/es of the center in the adsorbate to use for the
870 5261a07f Carles Marti
        definition of the adsorbate dihedral angle.
871 5261a07f Carles Marti
    @param ctr1_surf: The index/es of the center in the surface to use as
872 5261a07f Carles Marti
        adsorption center.
873 5261a07f Carles Marti
    @param ctr2_surf: The index/es of the center in the surface to use for the
874 5261a07f Carles Marti
        definition of the surf-adsorbate dihedral angle.
875 5261a07f Carles Marti
    @param num_pts: Number on which to sample Euler angles.
876 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision
877 5261a07f Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
878 5261a07f Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
879 5261a07f Carles Marti
        considered as atomic collision.
880 5261a07f Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
881 5261a07f Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
882 5261a07f Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
883 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
884 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
885 c25aa299 Carles Marti
    @param max_hel: Maximum value for sampling the helicopter
886 5261a07f Carles Marti
        (surf-adsorbate dihedral) angle.
887 fe91ddb2 Carles Marti
    @param height: Height on which to try adsorption.
888 9cd032cf Carles Marti
    @param excl_atom: Whether to exclude the adsorption center in the
889 9cd032cf Carles Marti
        molecule in the collision detection.
890 5261a07f Carles Marti
    @return: list of ase.Atoms object conatining all the orientations of a given
891 5261a07f Carles Marti
        conformer.
892 5261a07f Carles Marti
    """
893 7dd94df7 Carles Marti
    from copy import deepcopy
894 cf40df1b Carles Marti
    slab_ads_list = []
895 9dca524b Carles Marti
    # Rotation over bond angle
896 9dca524b Carles Marti
    for alpha in np.arange(90, 180+1, 90 / max(1, num_pts-1))[:num_pts]:
897 7dd94df7 Carles Marti
        # Rotation over surf-adsorbate dihedral angle
898 c25aa299 Carles Marti
        for beta in np.arange(0, max_hel, max_hel / num_pts):
899 7dd94df7 Carles Marti
            # Rotation over adsorbate bond dihedral angle
900 9dca524b Carles Marti
            for gamma in np.arange(90, 270+1, 180/max(1, num_pts-1))[:num_pts]:
901 9dca524b Carles Marti
                # Avoid duplicates as gamma rotation has no effect on plane
902 9dca524b Carles Marti
                # angles.
903 9dca524b Carles Marti
                if alpha == 180 and gamma > 90:
904 9dca524b Carles Marti
                    continue
905 7dd94df7 Carles Marti
                new_molec = deepcopy(orig_molec)
906 19727725 Carles Marti
                internal_rotate(new_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol,
907 19727725 Carles Marti
                                ctr1_surf, ctr2_surf, norm_vect, alpha,
908 19727725 Carles Marti
                                beta, gamma)
909 7dd94df7 Carles Marti
                site_coords = get_atom_coords(slab, ctr1_surf)
910 7dd94df7 Carles Marti
                ctr_coords = get_atom_coords(new_molec, ctr1_mol)
911 7dd94df7 Carles Marti
                slab_molec, collision = correct_coll(new_molec, slab,
912 7dd94df7 Carles Marti
                                                     ctr_coords, site_coords,
913 7dd94df7 Carles Marti
                                                     num_pts, min_coll_height,
914 7dd94df7 Carles Marti
                                                     norm_vect, slab_nghbs,
915 fe91ddb2 Carles Marti
                                                     molec_nghbs, coll_coeff,
916 9cd032cf Carles Marti
                                                     height, excl_atom)
917 bb387578 Carles Marti
                slab_molec.info = {**slab_molec.info, **new_molec.info}
918 cf40df1b Carles Marti
                if not collision and \
919 cf40df1b Carles Marti
                        not any([np.allclose(slab_molec.positions,
920 cf40df1b Carles Marti
                                             conf.positions)
921 cf40df1b Carles Marti
                                 for conf in slab_ads_list]):
922 cf40df1b Carles Marti
                    slab_ads_list.append(slab_molec)
923 cf40df1b Carles Marti
                    if h_donor is not False:
924 cf40df1b Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
925 cf40df1b Carles Marti
                                                          h_acceptor,
926 cf40df1b Carles Marti
                                                          len(slab)))
927 cf40df1b Carles Marti
928 cf40df1b Carles Marti
    return slab_ads_list
929 f3d1e601 Carles Marti
930 f3d1e601 Carles Marti
931 7dd94df7 Carles Marti
def adsorb_confs(conf_list, surf, inp_vars):
932 a5cc42ff Carles Marti
    """Generates a number of adsorbate-surface structure coordinates.
933 a5cc42ff Carles Marti

934 a5cc42ff Carles Marti
    Given a list of conformers, a surface, a list of atom indices (or list of
935 a5cc42ff Carles Marti
    list of indices) of both the surface and the adsorbate, it generates a
936 a5cc42ff Carles Marti
    number of adsorbate-surface structures for every possible combination of
937 a5cc42ff Carles Marti
    them at different orientations.
938 a5cc42ff Carles Marti
    @param conf_list: list of ase.Atoms of the different conformers
939 a5cc42ff Carles Marti
    @param surf: the ase.Atoms object of the surface
940 7dd94df7 Carles Marti
    @param inp_vars: Calculation parameters from input file.
941 a5cc42ff Carles Marti
    @return: list of ase.Atoms for the adsorbate-surface structures
942 a5cc42ff Carles Marti
    """
943 9cd032cf Carles Marti
    from copy import deepcopy
944 bb55f47c Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
945 7dd94df7 Carles Marti
    molec_ctrs = inp_vars['molec_ctrs']
946 7dd94df7 Carles Marti
    sites = inp_vars['sites']
947 7dd94df7 Carles Marti
    angles = inp_vars['set_angles']
948 7dd94df7 Carles Marti
    num_pts = inp_vars['sample_points_per_angle']
949 39df9c43 Carles Marti
    inp_norm_vect = inp_vars['surf_norm_vect']
950 7dd94df7 Carles Marti
    min_coll_height = inp_vars['min_coll_height']
951 7dd94df7 Carles Marti
    coll_coeff = inp_vars['collision_threshold']
952 9cd032cf Carles Marti
    exclude_ads_ctr = inp_vars['exclude_ads_ctr']
953 c25aa299 Carles Marti
    h_donor = inp_vars['h_donor']
954 c25aa299 Carles Marti
    h_acceptor = inp_vars['h_acceptor']
955 fe91ddb2 Carles Marti
    height = inp_vars['adsorption_height']
956 7dd94df7 Carles Marti
957 bc703cab Carles Marti
    if inp_vars['pbc_cell'] is not False:
958 bc703cab Carles Marti
        surf.set_pbc(True)
959 bc703cab Carles Marti
        surf.set_cell(inp_vars['pbc_cell'])
960 bc703cab Carles Marti
961 cf40df1b Carles Marti
    surf_ads_list = []
962 f3d1e601 Carles Marti
    sites_coords = get_atom_coords(surf, sites)
963 e8bebcca Carles Marti
    if coll_coeff is not False:
964 5fb01677 Carles Marti
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
965 5fb01677 Carles Marti
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
966 5fb01677 Carles Marti
    else:
967 5fb01677 Carles Marti
        surf_nghbs = 0
968 7dd94df7 Carles Marti
    for i, conf in enumerate(conf_list):
969 bb55f47c Carles Marti
        molec_ctr_coords = get_atom_coords(conf, molec_ctrs)
970 bc703cab Carles Marti
        if inp_vars['pbc_cell'] is not False:
971 bc703cab Carles Marti
            conf.set_pbc(True)
972 bc703cab Carles Marti
            conf.set_cell(inp_vars['pbc_cell'])
973 7dd94df7 Carles Marti
        for s, site in enumerate(sites_coords):
974 8279a51d Carles Marti
            if isinstance(inp_norm_vect, str) and inp_norm_vect == 'auto':
975 d6da8693 Carles Marti
                norm_vect = compute_norm_vect(surf, sites[s],
976 d6da8693 Carles Marti
                                              inp_vars['pbc_cell'])
977 39df9c43 Carles Marti
            else:
978 39df9c43 Carles Marti
                norm_vect = inp_norm_vect
979 7dd94df7 Carles Marti
            for c, ctr in enumerate(molec_ctr_coords):
980 9cd032cf Carles Marti
                if exclude_ads_ctr and isinstance(molec_ctrs[c], int):
981 9cd032cf Carles Marti
                    exclude_atom = molec_ctrs[c]
982 9cd032cf Carles Marti
                else:
983 9cd032cf Carles Marti
                    exclude_atom = False
984 9cd032cf Carles Marti
                    if exclude_ads_ctr and not isinstance(molec_ctrs[c], int):
985 9cd032cf Carles Marti
                        logger.warning("'exclude_ads_ctr' only works for atomic"
986 9cd032cf Carles Marti
                                       "centers and not for many-atoms "
987 9cd032cf Carles Marti
                                       f"barycenters. {molec_ctrs[c]} are not "
988 9cd032cf Carles Marti
                                       f"going to be excluded from collison.")
989 9cd032cf Carles Marti
                if coll_coeff and exclude_atom is not False:
990 9cd032cf Carles Marti
                    conf_wo_ctr = deepcopy(conf)
991 9cd032cf Carles Marti
                    del conf_wo_ctr[exclude_atom]
992 9cd032cf Carles Marti
                    conf_cutoffs = natural_cutoffs(conf_wo_ctr, mult=coll_coeff)
993 9cd032cf Carles Marti
                    molec_nghbs = len(neighbor_list("i", conf_wo_ctr,
994 9cd032cf Carles Marti
                                                    conf_cutoffs))
995 9cd032cf Carles Marti
                elif coll_coeff and exclude_atom is False:
996 9cd032cf Carles Marti
                    conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
997 9cd032cf Carles Marti
                    molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
998 9cd032cf Carles Marti
                else:
999 9cd032cf Carles Marti
                    molec_nghbs = 0
1000 7dd94df7 Carles Marti
                if angles == 'euler':
1001 cf40df1b Carles Marti
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
1002 cf40df1b Carles Marti
                                                   num_pts, min_coll_height,
1003 cf40df1b Carles Marti
                                                   coll_coeff, norm_vect,
1004 cf40df1b Carles Marti
                                                   surf_nghbs, molec_nghbs,
1005 9cd032cf Carles Marti
                                                   h_donor, h_acceptor, height,
1006 9cd032cf Carles Marti
                                                   exclude_atom))
1007 609104e3 Carles Marti
                elif angles == 'internal':
1008 7dd94df7 Carles Marti
                    mol_ctr1 = molec_ctrs[c]
1009 7dd94df7 Carles Marti
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
1010 7dd94df7 Carles Marti
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
1011 7dd94df7 Carles Marti
                    surf_ctr1 = sites[s]
1012 7dd94df7 Carles Marti
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
1013 a98d4172 Carles Marti
                    max_h = inp_vars["max_helic_angle"]
1014 d68dd4ad Carles Marti
                    surf_ads_list.extend(ads_internal(conf, surf, mol_ctr1,
1015 d68dd4ad Carles Marti
                                                      mol_ctr2, mol_ctr3,
1016 d68dd4ad Carles Marti
                                                      surf_ctr1, surf_ctr2,
1017 d68dd4ad Carles Marti
                                                      num_pts, min_coll_height,
1018 d68dd4ad Carles Marti
                                                      coll_coeff, norm_vect,
1019 d68dd4ad Carles Marti
                                                      surf_nghbs, molec_nghbs,
1020 d68dd4ad Carles Marti
                                                      h_donor, h_acceptor,
1021 9cd032cf Carles Marti
                                                      max_h, height,
1022 9cd032cf Carles Marti
                                                      exclude_atom))
1023 cf40df1b Carles Marti
    return surf_ads_list
1024 f3d1e601 Carles Marti
1025 f3d1e601 Carles Marti
1026 4614bb6a Carles
def run_screening(inp_vars):
1027 91ae8d86 Carles Marti
    """Carries out the screening of adsorbate structures on a surface.
1028 e07c09eb Carles

1029 e07c09eb Carles
    @param inp_vars: Calculation parameters from input file.
1030 e07c09eb Carles
    """
1031 e07c09eb Carles
    import os
1032 57e3a8c7 Carles Marti
    import random
1033 fd2384fc Carles Marti
    from modules.formats import collect_coords, adapt_format
1034 cf8fe0e3 Carles Marti
    from modules.calculation import run_calc, check_finished_calcs
1035 e07c09eb Carles
1036 76f4ac19 Carles Marti
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
1037 76f4ac19 Carles Marti
                ' structures.')
1038 b75bf97d Carles Marti
    if inp_vars['use_molec_file']:
1039 b75bf97d Carles Marti
        selected_confs = [adapt_format('ase', inp_vars['use_molec_file'])]
1040 b75bf97d Carles Marti
        logger.info(f"Using '{inp_vars['use_molec_file']}' as only conformer.")
1041 b75bf97d Carles Marti
    else:
1042 b75bf97d Carles Marti
        if not os.path.isdir("isolated"):
1043 b75bf97d Carles Marti
            err = "'isolated' directory not found. It is needed in order to " \
1044 b75bf97d Carles Marti
                  "carry out the screening of structures to be adsorbed"
1045 b75bf97d Carles Marti
            logger.error(err)
1046 b75bf97d Carles Marti
            raise FileNotFoundError(err)
1047 e07c09eb Carles
1048 b75bf97d Carles Marti
        correct_calcs, failed_calcs = check_finished_calcs('isolated',
1049 b75bf97d Carles Marti
                                                           inp_vars['code'])
1050 b75bf97d Carles Marti
        if not correct_calcs:
1051 b75bf97d Carles Marti
            err_msg = "No calculations on 'isolated' finished normally."
1052 b75bf97d Carles Marti
            logger.error(err_msg)
1053 b75bf97d Carles Marti
            raise FileNotFoundError(err_msg)
1054 b75bf97d Carles Marti
1055 b75bf97d Carles Marti
        logger.info(f"Found {len(correct_calcs)} structures of isolated "
1056 b75bf97d Carles Marti
                    f"conformers whose calculation finished normally.")
1057 b75bf97d Carles Marti
        if len(failed_calcs) != 0:
1058 b75bf97d Carles Marti
            logger.warning(
1059 b75bf97d Carles Marti
                f"Found {len(failed_calcs)} calculations more that "
1060 b75bf97d Carles Marti
                f"did not finish normally: {failed_calcs}. \n"
1061 b75bf97d Carles Marti
                f"Using only the ones that finished normally: "
1062 b75bf97d Carles Marti
                f"{correct_calcs}.")
1063 b75bf97d Carles Marti
1064 b75bf97d Carles Marti
        conformer_atoms_list = collect_coords(correct_calcs, inp_vars['code'],
1065 b75bf97d Carles Marti
                                              'isolated',
1066 b75bf97d Carles Marti
                                              inp_vars['special_atoms'])
1067 b75bf97d Carles Marti
        selected_confs = select_confs(conformer_atoms_list, correct_calcs,
1068 b75bf97d Carles Marti
                                      inp_vars['select_magns'],
1069 b75bf97d Carles Marti
                                      inp_vars['confs_per_magn'],
1070 b75bf97d Carles Marti
                                      inp_vars['code'])
1071 90819cc3 Carles Marti
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
1072 bb387578 Carles Marti
    surf.info = {}
1073 cf40df1b Carles Marti
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
1074 7d97341d Carles Marti
    if len(surf_ads_list) > inp_vars['max_structures']:
1075 57e3a8c7 Carles Marti
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
1076 a44ad3c2 Carles Martí
    elif len(surf_ads_list) == 0:
1077 a44ad3c2 Carles Martí
        err_msg = "No configurations were generated: Check the parameters in" \
1078 a44ad3c2 Carles Martí
                  "dockonsurf.inp"
1079 a44ad3c2 Carles Martí
        logger.error(err_msg)
1080 a44ad3c2 Carles Martí
        raise ValueError(err_msg)
1081 cf40df1b Carles Marti
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
1082 cf40df1b Carles Marti
                f'configurations to carry out a calculation of.')
1083 d9167fea Carles Marti
1084 f3d1e601 Carles Marti
    run_calc('screening', inp_vars, surf_ads_list)
1085 14f39d2a Carles Marti
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1086 14f39d2a Carles Marti
                ' structures section.')