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

dockonsurf / modules / screening.py @ 365d5b9a

Historique | Voir | Annoter | Télécharger (47,49 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 1e9e784d Carles Marti
        conf_enrgs = collect_energies(calc_dirs, code, 'isolated')
43 bfe93f0d Carles Marti
    if 'moi' in magns:
44 bfe93f0d Carles Marti
        mois = np.array([conf.get_moments_of_inertia() for conf in conf_list])
45 bfe93f0d Carles Marti
46 bfe93f0d Carles Marti
    # Assign values
47 bfe93f0d Carles Marti
    for i, conf in enumerate(conf_list):
48 bfe93f0d Carles Marti
        assign_prop(conf, 'idx', i)
49 bfe93f0d Carles Marti
        if 'energy' in magns:
50 bfe93f0d Carles Marti
            assign_prop(conf, 'energy', conf_enrgs[i])
51 bfe93f0d Carles Marti
        if 'moi' in magns:
52 bfe93f0d Carles Marti
            assign_prop(conf, 'moi', mois[i, 2])
53 bfe93f0d Carles Marti
54 bfe93f0d Carles Marti
    # pick ids
55 bfe93f0d Carles Marti
    for magn in magns:
56 bfe93f0d Carles Marti
        sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn]))
57 bfe93f0d Carles Marti
        if sorted_list[-1].info['idx'] not in selected_ids:
58 bfe93f0d Carles Marti
            selected_ids.append(sorted_list[-1].info['idx'])
59 bfe93f0d Carles Marti
        if num_sel > 1:
60 bfe93f0d Carles Marti
            for i in range(0, len(sorted_list) - 1,
61 bfe93f0d Carles Marti
                           len(conf_list) // (num_sel - 1)):
62 bfe93f0d Carles Marti
                if sorted_list[i].info['idx'] not in selected_ids:
63 bfe93f0d Carles Marti
                    selected_ids.append(sorted_list[i].info['idx'])
64 bfe93f0d Carles Marti
65 695dcff8 Carles Marti
    logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
66 bfe93f0d Carles Marti
    return [conf_list[idx] for idx in selected_ids]
67 bfe93f0d Carles Marti
68 bfe93f0d Carles Marti
69 7dd94df7 Carles Marti
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True):
70 b4aef3d7 Carles Marti
    """Computes the angle between two vectors.
71 b4aef3d7 Carles Marti

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

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

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

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

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

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

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

369 e8bebcca Carles Marti
    @param molec: ase.Atoms object of the molecule to adsorb
370 e8bebcca Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the
371 e8bebcca Carles Marti
        molecule
372 e8bebcca Carles Marti
    @param ctr_coord: The coordinates of the molecule to use as adsorption
373 e8bebcca Carles Marti
        center.
374 e8bebcca Carles Marti
    @param site_coord: The coordinates of the surface on which to adsorb the
375 e8bebcca Carles Marti
        molecule
376 e8bebcca Carles Marti
    @param num_pts: Number on which to sample Euler angles.
377 e8bebcca Carles Marti
    @param min_coll_height: The lowermost height for which to detect a collision
378 e8bebcca Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
379 e8bebcca Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
380 e8bebcca Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
381 e8bebcca Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
382 e8bebcca Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
383 e8bebcca Carles Marti
        considered as atomic collision.
384 e8bebcca Carles Marti
    @param height: Height on which to try adsorption
385 e8bebcca Carles Marti
    @return collision: bool, whether the structure generated has collisions
386 e8bebcca Carles Marti
        between slab and adsorbate.
387 e8bebcca Carles Marti
    """
388 e8bebcca Carles Marti
    from copy import deepcopy
389 e8bebcca Carles Marti
    slab_num_atoms = len(slab)
390 e8bebcca Carles Marti
    slab_molec = []
391 e8bebcca Carles Marti
    collision = True
392 e8bebcca Carles Marti
    max_corr = 6  # Should be an even number
393 e8bebcca Carles Marti
    d_angle = 180 / ((max_corr / 2.0) * num_pts)
394 e8bebcca Carles Marti
    num_corr = 0
395 e8bebcca Carles Marti
    while collision and num_corr <= max_corr:
396 e8bebcca Carles Marti
        k = num_corr * (-1) ** num_corr
397 e8bebcca Carles Marti
        slab_molec = deepcopy(slab)
398 e8bebcca Carles Marti
        molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
399 e8bebcca Carles Marti
                           center=ctr_coord)
400 e8bebcca Carles Marti
        add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height,
401 e8bebcca Carles Marti
                      norm_vect=norm_vect)
402 e8bebcca Carles Marti
        collision = check_collision(slab_molec, slab_num_atoms, min_coll_height,
403 e8bebcca Carles Marti
                                    norm_vect, slab_nghbs, molec_nghbs,
404 e8bebcca Carles Marti
                                    coll_coeff)
405 e8bebcca Carles Marti
        num_corr += 1
406 e8bebcca Carles Marti
    return slab_molec, collision
407 5f3f4b69 Carles Marti
408 5f3f4b69 Carles Marti
409 c25aa299 Carles Marti
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
410 c25aa299 Carles Marti
                 neigh_cutoff=1):
411 b4b2f307 Carles Marti
    # TODO rethink
412 91ae8d86 Carles Marti
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
413 b4b2f307 Carles Marti

414 91ae8d86 Carles Marti
    Tries to dissociate a H atom from the molecule and adsorb in on top of the
415 91ae8d86 Carles Marti
    surface if the distance is shorter than two times the neigh_cutoff value.
416 b4b2f307 Carles Marti
    @param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
417 b4b2f307 Carles Marti
    @param h_idx: The index of the hydrogen atom to carry out adsorption of.
418 b4b2f307 Carles Marti
    @param num_atoms_slab: The number of atoms of the slab without adsorbate.
419 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
420 b4b2f307 Carles Marti
    @param neigh_cutoff: half the maximum distance between the surface and the
421 b4b2f307 Carles Marti
        H for it to carry out dissociation.
422 b4b2f307 Carles Marti
    @return: An ase.Atoms object of the system adsorbate-surface with H
423 b4b2f307 Carles Marti
    """
424 b4b2f307 Carles Marti
    from copy import deepcopy
425 b4b2f307 Carles Marti
    from ase.neighborlist import NeighborList
426 b4b2f307 Carles Marti
    slab_molec = deepcopy(slab_molec_orig)
427 b4b2f307 Carles Marti
    cutoffs = len(slab_molec) * [neigh_cutoff]
428 c25aa299 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True, skin=0.0)
429 b4b2f307 Carles Marti
    nl.update(slab_molec)
430 b4b2f307 Carles Marti
    surf_h_vect = np.array([np.infty] * 3)
431 c25aa299 Carles Marti
    if h_acceptor == 'all':
432 c25aa299 Carles Marti
        h_acceptor = list(range(num_atoms_slab))
433 b4b2f307 Carles Marti
    for neigh_idx in nl.get_neighbors(h_idx)[0]:
434 c25aa299 Carles Marti
        for elm in h_acceptor:
435 c25aa299 Carles Marti
            if isinstance(elm, int):
436 c25aa299 Carles Marti
                if neigh_idx == elm and neigh_idx < num_atoms_slab:
437 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
438 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
439 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
440 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
441 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
442 c25aa299 Carles Marti
            else:
443 c25aa299 Carles Marti
                if slab_molec[neigh_idx].symbol == elm \
444 c25aa299 Carles Marti
                        and neigh_idx < num_atoms_slab:
445 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
446 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
447 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
448 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
449 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
450 c25aa299 Carles Marti
451 b4b2f307 Carles Marti
    if np.linalg.norm(surf_h_vect) != np.infty:
452 b4b2f307 Carles Marti
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
453 b4b2f307 Carles Marti
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
454 b4b2f307 Carles Marti
        return slab_molec
455 b4b2f307 Carles Marti
456 b4b2f307 Carles Marti
457 c25aa299 Carles Marti
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
458 b4b2f307 Carles Marti
    # TODO multiple dissociation
459 b4b2f307 Carles Marti
    """Decides which H atoms to dissociate according to a list of atoms.
460 b4b2f307 Carles Marti

461 b4b2f307 Carles Marti
    Given a list of chemical symbols or atom indices it checks for every atom
462 b4b2f307 Carles Marti
    or any of its neighbor if it's a H and calls dissociate_h to try to carry
463 b4b2f307 Carles Marti
    out dissociation of that H. For atom indices, it checks both whether
464 b4b2f307 Carles Marti
    the atom index or its neighbors are H, for chemical symbols, it only checks
465 b4b2f307 Carles Marti
    if there is a neighbor H.
466 b4b2f307 Carles Marti
    @param slab_molec: The ase.Atoms object of the system adsorbate-slab.
467 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
468 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
469 c25aa299 Carles Marti
    @param num_atoms_slab: Number of atoms of the bare slab.
470 b4b2f307 Carles Marti
    @return:
471 b4b2f307 Carles Marti
    """
472 b4b2f307 Carles Marti
    from ase.neighborlist import natural_cutoffs, NeighborList
473 b4b2f307 Carles Marti
    molec = slab_molec[num_atoms_slab:]
474 b4b2f307 Carles Marti
    cutoffs = natural_cutoffs(molec)
475 b4b2f307 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True)
476 b4b2f307 Carles Marti
    nl.update(molec)
477 b4b2f307 Carles Marti
    disso_structs = []
478 c25aa299 Carles Marti
    for el in h_donor:
479 b4b2f307 Carles Marti
        if isinstance(el, int):
480 b4b2f307 Carles Marti
            if molec[el].symbol == 'H':
481 b4b2f307 Carles Marti
                disso_struct = dissociate_h(slab_molec, el + num_atoms_slab,
482 c25aa299 Carles Marti
                                            num_atoms_slab, h_acceptor)
483 b4b2f307 Carles Marti
                if disso_struct is not None:
484 b4b2f307 Carles Marti
                    disso_structs.append(disso_struct)
485 b4b2f307 Carles Marti
            else:
486 b4b2f307 Carles Marti
                for neigh_idx in nl.get_neighbors(el)[0]:
487 b4b2f307 Carles Marti
                    if molec[neigh_idx].symbol == 'H':
488 b4b2f307 Carles Marti
                        disso_struct = dissociate_h(slab_molec, neigh_idx +
489 b4b2f307 Carles Marti
                                                    num_atoms_slab,
490 c25aa299 Carles Marti
                                                    num_atoms_slab, h_acceptor)
491 b4b2f307 Carles Marti
                        if disso_struct is not None:
492 b4b2f307 Carles Marti
                            disso_structs.append(disso_struct)
493 b4b2f307 Carles Marti
        else:
494 b4b2f307 Carles Marti
            for atom in molec:
495 b4b2f307 Carles Marti
                if atom.symbol.lower() == el.lower():
496 b4b2f307 Carles Marti
                    for neigh_idx in nl.get_neighbors(atom.index)[0]:
497 b4b2f307 Carles Marti
                        if molec[neigh_idx].symbol == 'H':
498 5261a07f Carles Marti
                            disso_struct = dissociate_h(slab_molec, neigh_idx
499 b4b2f307 Carles Marti
                                                        + num_atoms_slab,
500 c25aa299 Carles Marti
                                                        num_atoms_slab,
501 c25aa299 Carles Marti
                                                        h_acceptor)
502 b4b2f307 Carles Marti
                            if disso_struct is not None:
503 b4b2f307 Carles Marti
                                disso_structs.append(disso_struct)
504 b4b2f307 Carles Marti
    return disso_structs
505 b4b2f307 Carles Marti
506 b4b2f307 Carles Marti
507 3ab0865c Carles Marti
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
508 b4b2f307 Carles Marti
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
509 c25aa299 Carles Marti
              h_donor, h_acceptor):
510 3ab0865c Carles Marti
    """Generates adsorbate-surface structures by sampling over Euler angles.
511 3ab0865c Carles Marti

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

572 7dd94df7 Carles Marti
    Carles modification of chemcat's transform_adsorbate to work with
573 7dd94df7 Carles Marti
    coordinates instead of ase.Atom
574 7dd94df7 Carles Marti
    Parameters:
575 7dd94df7 Carles Marti
        molecule (ase.Atoms): The molecule to adsorb.
576 7dd94df7 Carles Marti

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

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

582 7dd94df7 Carles Marti
        ctr2_mol (int/list): The position of a second center of the
583 7dd94df7 Carles Marti
        adsorbate used to define the adsorption bond angle, and the dihedral
584 7dd94df7 Carles Marti
        adsorption angle.
585 7dd94df7 Carles Marti

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

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

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

595 7dd94df7 Carles Marti
        bond_vector (numpy.ndarray): The adsorption bond desired.
596 7dd94df7 Carles Marti
            Details: offset = vect(atom1_surf, atom1_mol)
597 7dd94df7 Carles Marti

598 7dd94df7 Carles Marti
        bond_angle_target (float or int): The adsorption bond angle desired (in
599 7dd94df7 Carles Marti
            degrees).
600 7dd94df7 Carles Marti
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
601 7dd94df7 Carles Marti

602 7dd94df7 Carles Marti
        dihedral_angle_target (float or int): The dihedral adsorption angle
603 7dd94df7 Carles Marti
            desired (in degrees).
604 7dd94df7 Carles Marti
            Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
605 7dd94df7 Carles Marti
            atom1_mol, atom2_mol)
606 7dd94df7 Carles Marti
                The rotation vector is facing the adsorbate from the surface
607 7dd94df7 Carles Marti
                (i.e. counterclockwise rotation when facing the surface (i.e.
608 7dd94df7 Carles Marti
                view from top))
609 7dd94df7 Carles Marti

610 7dd94df7 Carles Marti
        mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
611 7dd94df7 Carles Marti
            desired (in degrees).
612 7dd94df7 Carles Marti
            Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
613 7dd94df7 Carles Marti
            atom2_mol, atom3_mol)
614 7dd94df7 Carles Marti
                The rotation vector is facing atom2_mol from atom1_mol
615 7dd94df7 Carles Marti

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

837 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb (adsorbate).
838 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the molecule
839 5261a07f Carles Marti
    @param ctr1_mol: The index/es of the center in the adsorbate to use as
840 5261a07f Carles Marti
        adsorption center.
841 5261a07f Carles Marti
    @param ctr2_mol: The index/es of the center in the adsorbate to use for the
842 5261a07f Carles Marti
        definition of the surf-adsorbate angle, surf-adsorbate dihedral angle
843 5261a07f Carles Marti
        and adsorbate dihedral angle.
844 5261a07f Carles Marti
    @param ctr3_mol: The index/es of the center in the adsorbate to use for the
845 5261a07f Carles Marti
        definition of the adsorbate dihedral angle.
846 5261a07f Carles Marti
    @param ctr1_surf: The index/es of the center in the surface to use as
847 5261a07f Carles Marti
        adsorption center.
848 5261a07f Carles Marti
    @param ctr2_surf: The index/es of the center in the surface to use for the
849 5261a07f Carles Marti
        definition of the surf-adsorbate dihedral angle.
850 5261a07f Carles Marti
    @param num_pts: Number on which to sample Euler angles.
851 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision
852 5261a07f Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
853 5261a07f Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
854 5261a07f Carles Marti
        considered as atomic collision.
855 5261a07f Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
856 5261a07f Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
857 5261a07f Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
858 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
859 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
860 c25aa299 Carles Marti
    @param max_hel: Maximum value for sampling the helicopter
861 5261a07f Carles Marti
        (surf-adsorbate dihedral) angle.
862 5261a07f Carles Marti
    @return: list of ase.Atoms object conatining all the orientations of a given
863 5261a07f Carles Marti
        conformer.
864 5261a07f Carles Marti
    """
865 7dd94df7 Carles Marti
    from copy import deepcopy
866 cf40df1b Carles Marti
    slab_ads_list = []
867 9dca524b Carles Marti
    # Rotation over bond angle
868 9dca524b Carles Marti
    for alpha in np.arange(90, 180+1, 90 / max(1, num_pts-1))[:num_pts]:
869 7dd94df7 Carles Marti
        # Rotation over surf-adsorbate dihedral angle
870 c25aa299 Carles Marti
        for beta in np.arange(0, max_hel, max_hel / num_pts):
871 7dd94df7 Carles Marti
            # Rotation over adsorbate bond dihedral angle
872 9dca524b Carles Marti
            for gamma in np.arange(90, 270+1, 180/max(1, num_pts-1))[:num_pts]:
873 9dca524b Carles Marti
                # Avoid duplicates as gamma rotation has no effect on plane
874 9dca524b Carles Marti
                # angles.
875 9dca524b Carles Marti
                if alpha == 180 and gamma > 90:
876 9dca524b Carles Marti
                    continue
877 7dd94df7 Carles Marti
                new_molec = deepcopy(orig_molec)
878 19727725 Carles Marti
                internal_rotate(new_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol,
879 19727725 Carles Marti
                                ctr1_surf, ctr2_surf, norm_vect, alpha,
880 19727725 Carles Marti
                                beta, gamma)
881 7dd94df7 Carles Marti
                site_coords = get_atom_coords(slab, ctr1_surf)
882 7dd94df7 Carles Marti
                ctr_coords = get_atom_coords(new_molec, ctr1_mol)
883 7dd94df7 Carles Marti
                slab_molec, collision = correct_coll(new_molec, slab,
884 7dd94df7 Carles Marti
                                                     ctr_coords, site_coords,
885 7dd94df7 Carles Marti
                                                     num_pts, min_coll_height,
886 7dd94df7 Carles Marti
                                                     norm_vect, slab_nghbs,
887 e8bebcca Carles Marti
                                                     molec_nghbs, coll_coeff)
888 cf40df1b Carles Marti
                if not collision and \
889 cf40df1b Carles Marti
                        not any([np.allclose(slab_molec.positions,
890 cf40df1b Carles Marti
                                             conf.positions)
891 cf40df1b Carles Marti
                                 for conf in slab_ads_list]):
892 cf40df1b Carles Marti
                    slab_ads_list.append(slab_molec)
893 cf40df1b Carles Marti
                    if h_donor is not False:
894 cf40df1b Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
895 cf40df1b Carles Marti
                                                          h_acceptor,
896 cf40df1b Carles Marti
                                                          len(slab)))
897 cf40df1b Carles Marti
898 cf40df1b Carles Marti
    return slab_ads_list
899 f3d1e601 Carles Marti
900 f3d1e601 Carles Marti
901 7dd94df7 Carles Marti
def adsorb_confs(conf_list, surf, inp_vars):
902 a5cc42ff Carles Marti
    """Generates a number of adsorbate-surface structure coordinates.
903 a5cc42ff Carles Marti

904 a5cc42ff Carles Marti
    Given a list of conformers, a surface, a list of atom indices (or list of
905 a5cc42ff Carles Marti
    list of indices) of both the surface and the adsorbate, it generates a
906 a5cc42ff Carles Marti
    number of adsorbate-surface structures for every possible combination of
907 a5cc42ff Carles Marti
    them at different orientations.
908 a5cc42ff Carles Marti
    @param conf_list: list of ase.Atoms of the different conformers
909 a5cc42ff Carles Marti
    @param surf: the ase.Atoms object of the surface
910 7dd94df7 Carles Marti
    @param inp_vars: Calculation parameters from input file.
911 a5cc42ff Carles Marti
    @return: list of ase.Atoms for the adsorbate-surface structures
912 a5cc42ff Carles Marti
    """
913 bb55f47c Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
914 7dd94df7 Carles Marti
    molec_ctrs = inp_vars['molec_ctrs']
915 7dd94df7 Carles Marti
    sites = inp_vars['sites']
916 7dd94df7 Carles Marti
    angles = inp_vars['set_angles']
917 7dd94df7 Carles Marti
    num_pts = inp_vars['sample_points_per_angle']
918 39df9c43 Carles Marti
    inp_norm_vect = inp_vars['surf_norm_vect']
919 7dd94df7 Carles Marti
    min_coll_height = inp_vars['min_coll_height']
920 7dd94df7 Carles Marti
    coll_coeff = inp_vars['collision_threshold']
921 c25aa299 Carles Marti
    h_donor = inp_vars['h_donor']
922 c25aa299 Carles Marti
    h_acceptor = inp_vars['h_acceptor']
923 7dd94df7 Carles Marti
924 bc703cab Carles Marti
    if inp_vars['pbc_cell'] is not False:
925 bc703cab Carles Marti
        surf.set_pbc(True)
926 bc703cab Carles Marti
        surf.set_cell(inp_vars['pbc_cell'])
927 bc703cab Carles Marti
928 cf40df1b Carles Marti
    surf_ads_list = []
929 f3d1e601 Carles Marti
    sites_coords = get_atom_coords(surf, sites)
930 e8bebcca Carles Marti
    if coll_coeff is not False:
931 5fb01677 Carles Marti
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
932 5fb01677 Carles Marti
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
933 5fb01677 Carles Marti
    else:
934 5fb01677 Carles Marti
        surf_nghbs = 0
935 7dd94df7 Carles Marti
    for i, conf in enumerate(conf_list):
936 bb55f47c Carles Marti
        molec_ctr_coords = get_atom_coords(conf, molec_ctrs)
937 bc703cab Carles Marti
        if inp_vars['pbc_cell'] is not False:
938 bc703cab Carles Marti
            conf.set_pbc(True)
939 bc703cab Carles Marti
            conf.set_cell(inp_vars['pbc_cell'])
940 e8bebcca Carles Marti
        if coll_coeff is not False:
941 5fb01677 Carles Marti
            conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
942 5fb01677 Carles Marti
            molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
943 5fb01677 Carles Marti
        else:
944 5fb01677 Carles Marti
            molec_nghbs = 0
945 7dd94df7 Carles Marti
        for s, site in enumerate(sites_coords):
946 8279a51d Carles Marti
            if isinstance(inp_norm_vect, str) and inp_norm_vect == 'auto':
947 d6da8693 Carles Marti
                norm_vect = compute_norm_vect(surf, sites[s],
948 d6da8693 Carles Marti
                                              inp_vars['pbc_cell'])
949 39df9c43 Carles Marti
            else:
950 39df9c43 Carles Marti
                norm_vect = inp_norm_vect
951 7dd94df7 Carles Marti
            for c, ctr in enumerate(molec_ctr_coords):
952 7dd94df7 Carles Marti
                if angles == 'euler':
953 cf40df1b Carles Marti
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
954 cf40df1b Carles Marti
                                                   num_pts, min_coll_height,
955 cf40df1b Carles Marti
                                                   coll_coeff, norm_vect,
956 cf40df1b Carles Marti
                                                   surf_nghbs, molec_nghbs,
957 cf40df1b Carles Marti
                                                   h_donor, h_acceptor))
958 609104e3 Carles Marti
                elif angles == 'internal':
959 7dd94df7 Carles Marti
                    mol_ctr1 = molec_ctrs[c]
960 7dd94df7 Carles Marti
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
961 7dd94df7 Carles Marti
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
962 7dd94df7 Carles Marti
                    surf_ctr1 = sites[s]
963 7dd94df7 Carles Marti
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
964 a98d4172 Carles Marti
                    max_h = inp_vars["max_helic_angle"]
965 d68dd4ad Carles Marti
                    surf_ads_list.extend(ads_internal(conf, surf, mol_ctr1,
966 d68dd4ad Carles Marti
                                                      mol_ctr2, mol_ctr3,
967 d68dd4ad Carles Marti
                                                      surf_ctr1, surf_ctr2,
968 d68dd4ad Carles Marti
                                                      num_pts, min_coll_height,
969 d68dd4ad Carles Marti
                                                      coll_coeff, norm_vect,
970 d68dd4ad Carles Marti
                                                      surf_nghbs, molec_nghbs,
971 d68dd4ad Carles Marti
                                                      h_donor, h_acceptor,
972 d68dd4ad Carles Marti
                                                      max_h))
973 cf40df1b Carles Marti
    return surf_ads_list
974 f3d1e601 Carles Marti
975 f3d1e601 Carles Marti
976 4614bb6a Carles
def run_screening(inp_vars):
977 91ae8d86 Carles Marti
    """Carries out the screening of adsorbate structures on a surface.
978 e07c09eb Carles

979 e07c09eb Carles
    @param inp_vars: Calculation parameters from input file.
980 e07c09eb Carles
    """
981 e07c09eb Carles
    import os
982 57e3a8c7 Carles Marti
    import random
983 fd2384fc Carles Marti
    from modules.formats import collect_coords, adapt_format
984 cf8fe0e3 Carles Marti
    from modules.calculation import run_calc, check_finished_calcs
985 e07c09eb Carles
986 76f4ac19 Carles Marti
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
987 76f4ac19 Carles Marti
                ' structures.')
988 b75bf97d Carles Marti
    if inp_vars['use_molec_file']:
989 b75bf97d Carles Marti
        selected_confs = [adapt_format('ase', inp_vars['use_molec_file'])]
990 b75bf97d Carles Marti
        logger.info(f"Using '{inp_vars['use_molec_file']}' as only conformer.")
991 b75bf97d Carles Marti
    else:
992 b75bf97d Carles Marti
        if not os.path.isdir("isolated"):
993 b75bf97d Carles Marti
            err = "'isolated' directory not found. It is needed in order to " \
994 b75bf97d Carles Marti
                  "carry out the screening of structures to be adsorbed"
995 b75bf97d Carles Marti
            logger.error(err)
996 b75bf97d Carles Marti
            raise FileNotFoundError(err)
997 e07c09eb Carles
998 b75bf97d Carles Marti
        correct_calcs, failed_calcs = check_finished_calcs('isolated',
999 b75bf97d Carles Marti
                                                           inp_vars['code'])
1000 b75bf97d Carles Marti
        if not correct_calcs:
1001 b75bf97d Carles Marti
            err_msg = "No calculations on 'isolated' finished normally."
1002 b75bf97d Carles Marti
            logger.error(err_msg)
1003 b75bf97d Carles Marti
            raise FileNotFoundError(err_msg)
1004 b75bf97d Carles Marti
1005 b75bf97d Carles Marti
        logger.info(f"Found {len(correct_calcs)} structures of isolated "
1006 b75bf97d Carles Marti
                    f"conformers whose calculation finished normally.")
1007 b75bf97d Carles Marti
        if len(failed_calcs) != 0:
1008 b75bf97d Carles Marti
            logger.warning(
1009 b75bf97d Carles Marti
                f"Found {len(failed_calcs)} calculations more that "
1010 b75bf97d Carles Marti
                f"did not finish normally: {failed_calcs}. \n"
1011 b75bf97d Carles Marti
                f"Using only the ones that finished normally: "
1012 b75bf97d Carles Marti
                f"{correct_calcs}.")
1013 b75bf97d Carles Marti
1014 b75bf97d Carles Marti
        conformer_atoms_list = collect_coords(correct_calcs, inp_vars['code'],
1015 b75bf97d Carles Marti
                                              'isolated',
1016 b75bf97d Carles Marti
                                              inp_vars['special_atoms'])
1017 b75bf97d Carles Marti
        selected_confs = select_confs(conformer_atoms_list, correct_calcs,
1018 b75bf97d Carles Marti
                                      inp_vars['select_magns'],
1019 b75bf97d Carles Marti
                                      inp_vars['confs_per_magn'],
1020 b75bf97d Carles Marti
                                      inp_vars['code'])
1021 90819cc3 Carles Marti
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
1022 cf40df1b Carles Marti
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
1023 7d97341d Carles Marti
    if len(surf_ads_list) > inp_vars['max_structures']:
1024 57e3a8c7 Carles Marti
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
1025 cf40df1b Carles Marti
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
1026 cf40df1b Carles Marti
                f'configurations to carry out a calculation of.')
1027 d9167fea Carles Marti
1028 f3d1e601 Carles Marti
    run_calc('screening', inp_vars, surf_ads_list)
1029 14f39d2a Carles Marti
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1030 14f39d2a Carles Marti
                ' structures section.')