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

dockonsurf / modules / screening.py @ 8d2efb97

Historique | Voir | Annoter | Télécharger (49,64 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 1d8c374e Carles Martí
def select_confs(conf_list: list, magns: list, num_sel: int):
9 bfe93f0d Carles Marti
    """Takes a list ase.Atoms and selects the most different magnitude-wise.
10 bfe93f0d Carles Marti

11 bfe93f0d Carles Marti
    Given a list of ase.Atoms objects and a list of magnitudes, it selects a
12 bfe93f0d Carles Marti
    number of the most different conformers according to every magnitude
13 bfe93f0d Carles Marti
    specified.
14 1e9e784d Carles Marti
     
15 1d8c374e Carles Martí
    @param conf_list: list of ase.Atoms objects to select among.
16 bfe93f0d Carles Marti
    @param magns: list of str with the names of the magnitudes to use for the
17 1e9e784d Carles Marti
        conformer selection. Supported magnitudes: 'energy', 'moi'.
18 bfe93f0d Carles Marti
    @param num_sel: number of conformers to select for every of the magnitudes.
19 bfe93f0d Carles Marti
    @return: list of the selected ase.Atoms objects.
20 bfe93f0d Carles Marti
    """
21 1d8c374e Carles Martí
    selected_ids = []
22 bfe93f0d Carles Marti
    if num_sel >= len(conf_list):
23 bfe93f0d Carles Marti
        logger.warning('Number of conformers per magnitude is equal or larger '
24 bfe93f0d Carles Marti
                       'than the total number of conformers. Using all '
25 695dcff8 Carles Marti
                       f'available conformers: {len(conf_list)}.')
26 bfe93f0d Carles Marti
        return conf_list
27 bfe93f0d Carles Marti
28 1d8c374e Carles Martí
    # Assign mois
29 bfe93f0d Carles Marti
    if 'moi' in magns:
30 1d8c374e Carles Martí
        for conf in conf_list:
31 1d8c374e Carles Martí
            conf.info["moi"] = conf.get_moments_of_inertia()[2]
32 bfe93f0d Carles Marti
33 bfe93f0d Carles Marti
    # pick ids
34 bfe93f0d Carles Marti
    for magn in magns:
35 bfe93f0d Carles Marti
        sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn]))
36 1d8c374e Carles Martí
        if sorted_list[-1].info['iso'] not in selected_ids:
37 1d8c374e Carles Martí
            selected_ids.append(sorted_list[-1].info['iso'])
38 bfe93f0d Carles Marti
        if num_sel > 1:
39 bfe93f0d Carles Marti
            for i in range(0, len(sorted_list) - 1,
40 bfe93f0d Carles Marti
                           len(conf_list) // (num_sel - 1)):
41 1d8c374e Carles Martí
                if sorted_list[i].info['iso'] not in selected_ids:
42 1d8c374e Carles Martí
                    selected_ids.append(sorted_list[i].info['iso'])
43 bfe93f0d Carles Marti
44 695dcff8 Carles Marti
    logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
45 1d8c374e Carles Martí
    return [conf for conf in conf_list if conf.info["iso"] in selected_ids]
46 bfe93f0d Carles Marti
47 bfe93f0d Carles Marti
48 7dd94df7 Carles Marti
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True):
49 b4aef3d7 Carles Marti
    """Computes the angle between two vectors.
50 b4aef3d7 Carles Marti

51 b4aef3d7 Carles Marti
    @param v1: The first vector.
52 b4aef3d7 Carles Marti
    @param v2: The second vector.
53 b516a1d6 Carles Marti
    @param ref: Orthogonal vector to both v1 and v2,
54 b516a1d6 Carles Marti
        along which the sign of the rotation is defined (i.e. positive if
55 b516a1d6 Carles Marti
        counterclockwise angle when facing ref)
56 b4aef3d7 Carles Marti
    @param degrees: Whether the result should be in radians (True) or in
57 b4aef3d7 Carles Marti
        degrees (False).
58 b4aef3d7 Carles Marti
    @return: The angle in radians if degrees = False, or in degrees if
59 b4aef3d7 Carles Marti
        degrees =True
60 b4aef3d7 Carles Marti
    """
61 b4aef3d7 Carles Marti
    v1_u = v1 / np.linalg.norm(v1)
62 b4aef3d7 Carles Marti
    v2_u = v2 / np.linalg.norm(v2)
63 b4aef3d7 Carles Marti
    angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0))
64 b516a1d6 Carles Marti
    if ref is not None:
65 b516a1d6 Carles Marti
        # Give sign according to ref direction
66 b516a1d6 Carles Marti
        angle *= (1 if np.dot(np.cross(v1, v2), ref) >= 0 else -1)
67 b516a1d6 Carles Marti
68 b4aef3d7 Carles Marti
    return angle if not degrees else angle * 180 / np.pi
69 b4aef3d7 Carles Marti
70 b4aef3d7 Carles Marti
71 12a12bdd Carles Marti
def vect_avg(vects):
72 a5cc42ff Carles Marti
    """Computes the element-wise mean of a set of vectors.
73 12a12bdd Carles Marti

74 dadc6016 Carles Marti
    @param vects: list of lists-like: containing the vectors (num_vectors,
75 dadc6016 Carles Marti
        length_vector).
76 a5cc42ff Carles Marti
    @return: vector average computed doing the element-wise mean.
77 12a12bdd Carles Marti
    """
78 d26885d3 Carles Marti
    from modules.utilities import try_command
79 12a12bdd Carles Marti
    err = "vect_avg parameter vects must be a list-like, able to be converted" \
80 12a12bdd Carles Marti
          " np.array"
81 dadc6016 Carles Marti
    array = try_command(np.array, [(ValueError, err)], vects)
82 dadc6016 Carles Marti
    if len(array.shape) == 1:
83 dadc6016 Carles Marti
        return array
84 12a12bdd Carles Marti
    else:
85 dadc6016 Carles Marti
        num_vects = array.shape[1]
86 dadc6016 Carles Marti
        return np.array([np.average(array[:, i]) for i in range(num_vects)])
87 12a12bdd Carles Marti
88 12a12bdd Carles Marti
89 ae097639 Carles Martí
def get_atom_coords(atoms: ase.Atoms, center=None):
90 ae097639 Carles Martí
    """Gets the coordinates of the specified center for an ase.Atoms object.
91 a5cc42ff Carles Marti

92 ae097639 Carles Martí
    If center is not an index but a list of indices, it computes the
93 a5cc42ff Carles Marti
    element-wise mean of the coordinates of the atoms specified in the inner
94 a5cc42ff Carles Marti
    list.
95 a5cc42ff Carles Marti
    @param atoms: ase.Atoms object for which to obtain the coordinates of.
96 ae097639 Carles Martí
    @param center: index/list of indices of the atoms for which the coordinates
97 ae097639 Carles Martí
                   should be extracted.
98 a5cc42ff Carles Marti
    @return: np.ndarray of atomic coordinates.
99 a5cc42ff Carles Marti
    """
100 ae097639 Carles Martí
    err_msg = "Argument 'ctr' must be an integer or a list of integers. "\
101 ae097639 Carles Martí
              "Every integer must be in the range [0, num_atoms)"
102 ae097639 Carles Martí
    if center is None:
103 ae097639 Carles Martí
        center = list(range(len(atoms)))
104 ae097639 Carles Martí
    if isinstance(center, int):
105 ae097639 Carles Martí
        if center not in list(range(len(atoms))):
106 ae097639 Carles Martí
            logger.error(err_msg)
107 ae097639 Carles Martí
            raise ValueError(err_msg)
108 ae097639 Carles Martí
        return atoms[center].position
109 ae097639 Carles Martí
    elif isinstance(center, list):
110 ae097639 Carles Martí
        for elm in center:
111 ae097639 Carles Martí
            if elm not in list(range(len(atoms))):
112 ae097639 Carles Martí
                logger.error(err_msg)
113 ae097639 Carles Martí
                raise ValueError(err_msg)
114 ae097639 Carles Martí
        return vect_avg([atoms[idx].position for idx in center])
115 ae097639 Carles Martí
    else:
116 ae097639 Carles Martí
        logger.error(err_msg)
117 ae097639 Carles Martí
        raise ValueError(err_msg)
118 f3d1e601 Carles Marti
119 f3d1e601 Carles Marti
120 d6da8693 Carles Marti
def compute_norm_vect(atoms, idxs, cell):
121 d6da8693 Carles Marti
    """Computes the local normal vector of a surface at a given site.
122 d6da8693 Carles Marti

123 d6da8693 Carles Marti
    Given an ase.Atoms object and a site defined as a linear combination of
124 d6da8693 Carles Marti
    atoms it computes the vector perpendicular to the surface, considering the
125 d6da8693 Carles Marti
    local environment of the site.
126 d6da8693 Carles Marti
    @param atoms: ase.Atoms object of the surface.
127 d6da8693 Carles Marti
    @param idxs: list or int: Index or list of indices of the atom/s that define
128 d6da8693 Carles Marti
        the site
129 d6da8693 Carles Marti
    @param cell: Unit cell. A 3x3 matrix (the three unit cell vectors)
130 d6da8693 Carles Marti
    @return: numpy.ndarray of the coordinates of the vector locally
131 d6da8693 Carles Marti
    perpendicular to the surface.
132 d6da8693 Carles Marti
    """
133 6d1bd296 Carles Marti
    from modules.ASANN import coordination_numbers as coord_nums
134 d6da8693 Carles Marti
    if isinstance(idxs, list):
135 d6da8693 Carles Marti
        atm_vect = [-np.round(coord_nums(atoms.get_scaled_positions(),
136 d6da8693 Carles Marti
                                         pbc=np.any(cell),
137 d6da8693 Carles Marti
                                         cell_vectors=cell)[3][i], 2)
138 d6da8693 Carles Marti
                    for i in idxs]
139 d6da8693 Carles Marti
        norm_vec = vect_avg([vect / np.linalg.norm(vect) for vect in atm_vect])
140 d6da8693 Carles Marti
    elif isinstance(idxs, int):
141 d6da8693 Carles Marti
        norm_vec = -coord_nums(atoms.get_scaled_positions(),
142 d6da8693 Carles Marti
                               pbc=np.any(cell),
143 d6da8693 Carles Marti
                               cell_vectors=cell)[3][idxs]
144 d6da8693 Carles Marti
    else:
145 d6da8693 Carles Marti
        err = "'idxs' must be either an int or a list"
146 d6da8693 Carles Marti
        logger.error(err)
147 d6da8693 Carles Marti
        raise ValueError(err)
148 587dca22 Carles Marti
    norm_vec = np.round(norm_vec, 2) / np.linalg.norm(np.round(norm_vec, 2))
149 77f8c47b Carles Martí
    if not np.isnan(norm_vec).any():
150 77f8c47b Carles Martí
        logger.info(f"The perpendicular vector to the surface at site '{idxs}' "
151 77f8c47b Carles Martí
                    f"is {norm_vec}")
152 587dca22 Carles Marti
    return norm_vec
153 d6da8693 Carles Marti
154 d6da8693 Carles Marti
155 d8d92cfb Carles Marti
def align_molec(orig_molec, ctr_coord, ref_vect):
156 4918e2ad Carles Marti
    """Align a molecule to a vector by a center.
157 4918e2ad Carles Marti

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

241 1d22a086 Carles Marti
    This function extends the functionality of ase.build.add_adsorbate
242 1d22a086 Carles Marti
    (https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
243 1d22a086 Carles Marti
    by enabling to change the z coordinate and the axis perpendicular to the
244 1d22a086 Carles Marti
    surface.
245 1d22a086 Carles Marti
    @param slab: ase.Atoms object containing the coordinates of the surface
246 1d22a086 Carles Marti
    @param adsorbate: ase.Atoms object containing the coordinates of the
247 1d22a086 Carles Marti
        adsorbate.
248 dadc6016 Carles Marti
    @param site_coord: The coordinates of the adsorption site on the surface.
249 dadc6016 Carles Marti
    @param ctr_coord: The coordinates of the adsorption center in the molecule.
250 dadc6016 Carles Marti
    @param height: The height above the surface where to adsorb.
251 1d22a086 Carles Marti
    @param offset: Offsets the adsorbate by a number of unit cells. Mostly
252 1d22a086 Carles Marti
        useful when adding more than one adsorbate.
253 1d22a086 Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
254 1d22a086 Carles Marti
    """
255 36d92f4f Carles Marti
    from copy import deepcopy
256 1d22a086 Carles Marti
    info = slab.info.get('adsorbate_info', {})
257 dadc6016 Carles Marti
    pos = np.array([0.0, 0.0, 0.0])  # part of absolute coordinates
258 1d22a086 Carles Marti
    spos = np.array([0.0, 0.0, 0.0])  # part relative to unit cell
259 dadc6016 Carles Marti
    norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect)
260 1d22a086 Carles Marti
    if offset is not None:
261 1d22a086 Carles Marti
        spos += np.asarray(offset, float)
262 dadc6016 Carles Marti
    if isinstance(site_coord, str):
263 1d22a086 Carles Marti
        # A site-name:
264 1d22a086 Carles Marti
        if 'sites' not in info:
265 1d22a086 Carles Marti
            raise TypeError('If the atoms are not made by an ase.build '
266 1d22a086 Carles Marti
                            'function, position cannot be a name.')
267 dadc6016 Carles Marti
        if site_coord not in info['sites']:
268 dadc6016 Carles Marti
            raise TypeError('Adsorption site %s not supported.' % site_coord)
269 dadc6016 Carles Marti
        spos += info['sites'][site_coord]
270 1d22a086 Carles Marti
    else:
271 dadc6016 Carles Marti
        pos += site_coord
272 1d22a086 Carles Marti
    if 'cell' in info:
273 1d22a086 Carles Marti
        cell = info['cell']
274 1d22a086 Carles Marti
    else:
275 1d22a086 Carles Marti
        cell = slab.get_cell()
276 1d22a086 Carles Marti
    pos += np.dot(spos, cell)
277 1d22a086 Carles Marti
    # Convert the adsorbate to an Atoms object
278 1d22a086 Carles Marti
    if isinstance(adsorbate, ase.Atoms):
279 36d92f4f Carles Marti
        ads = deepcopy(adsorbate)
280 1d22a086 Carles Marti
    elif isinstance(adsorbate, ase.Atom):
281 1d22a086 Carles Marti
        ads = ase.Atoms([adsorbate])
282 1d22a086 Carles Marti
    else:
283 1d22a086 Carles Marti
        # Assume it is a string representing a single Atom
284 1d22a086 Carles Marti
        ads = ase.Atoms([ase.Atom(adsorbate)])
285 dadc6016 Carles Marti
    pos += height * norm_vect_u
286 1d22a086 Carles Marti
    # Move adsorbate into position
287 dadc6016 Carles Marti
    ads.translate(pos - ctr_coord)
288 1d22a086 Carles Marti
    # Attach the adsorbate
289 1d22a086 Carles Marti
    slab.extend(ads)
290 1d22a086 Carles Marti
291 1d22a086 Carles Marti
292 a4b57124 Carles Marti
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
293 9cd032cf Carles Marti
                    nn_molec=0, coll_coeff=1.2, exclude_atom=False):
294 5f3f4b69 Carles Marti
    """Checks whether a slab and a molecule collide or not.
295 5f3f4b69 Carles Marti

296 5f3f4b69 Carles Marti
    @param slab_molec: The system of adsorbate-slab for which to detect if there
297 5f3f4b69 Carles Marti
        are collisions.
298 5f3f4b69 Carles Marti
    @param nn_slab: Number of neigbors in the surface.
299 5f3f4b69 Carles Marti
    @param nn_molec: Number of neighbors in the molecule.
300 5f3f4b69 Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
301 5f3f4b69 Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
302 5f3f4b69 Carles Marti
        considered as atomic collision.
303 5f3f4b69 Carles Marti
    @param slab_num_atoms: Number of atoms of the bare slab.
304 5f3f4b69 Carles Marti
    @param min_height: The minimum height atoms can have to not be considered as
305 5f3f4b69 Carles Marti
        colliding.
306 5f3f4b69 Carles Marti
    @param vect: The vector perpendicular to the slab.
307 9cd032cf Carles Marti
    @param exclude_atom: Whether to exclude the adsorption center in the
308 9cd032cf Carles Marti
        molecule in the collision detection.
309 5f3f4b69 Carles Marti
    @return: bool, whether the surface and the molecule collide.
310 e8bebcca Carles Marti
    """
311 9cd032cf Carles Marti
    from copy import deepcopy
312 5f3f4b69 Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
313 e8bebcca Carles Marti
314 e8bebcca Carles Marti
    # Check structure overlap by height
315 5fb01677 Carles Marti
    if min_height is not False:
316 a4b57124 Carles Marti
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
317 a4b57124 Carles Marti
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
318 e8bebcca Carles Marti
        if vect.tolist() not in cart_axes:
319 e8bebcca Carles Marti
            err_msg = "'min_coll_height' option is only implemented for " \
320 e8bebcca Carles Marti
                      "'surf_norm_vect' to be one of the x, y or z axes. "
321 e8bebcca Carles Marti
            logger.error(err_msg)
322 e8bebcca Carles Marti
            raise ValueError(err_msg)
323 e8bebcca Carles Marti
        for atom in slab_molec[slab_num_atoms:]:
324 9cd032cf Carles Marti
            if exclude_atom is not False \
325 9cd032cf Carles Marti
                    and atom.index == exclude_atom:
326 9cd032cf Carles Marti
                continue
327 e8bebcca Carles Marti
            for i, coord in enumerate(vect):
328 e8bebcca Carles Marti
                if coord == 0:
329 e8bebcca Carles Marti
                    continue
330 e8bebcca Carles Marti
                if atom.position[i] * coord < min_height * coord:
331 e8bebcca Carles Marti
                    return True
332 e8bebcca Carles Marti
333 e8bebcca Carles Marti
    # Check structure overlap by sphere collision
334 e8bebcca Carles Marti
    if coll_coeff is not False:
335 9cd032cf Carles Marti
        if exclude_atom is not False:
336 9cd032cf Carles Marti
            slab_molec_wo_ctr = deepcopy(slab_molec)
337 9cd032cf Carles Marti
            del slab_molec_wo_ctr[exclude_atom + slab_num_atoms]
338 9cd032cf Carles Marti
            slab_molec_cutoffs = natural_cutoffs(slab_molec_wo_ctr,
339 9cd032cf Carles Marti
                                                 mult=coll_coeff)
340 9cd032cf Carles Marti
            slab_molec_nghbs = len(neighbor_list("i", slab_molec_wo_ctr,
341 9cd032cf Carles Marti
                                                 slab_molec_cutoffs))
342 9cd032cf Carles Marti
        else:
343 9cd032cf Carles Marti
            slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
344 9cd032cf Carles Marti
            slab_molec_nghbs = len(neighbor_list("i", slab_molec,
345 9cd032cf Carles Marti
                                                 slab_molec_cutoffs))
346 5fb01677 Carles Marti
        if slab_molec_nghbs > nn_slab + nn_molec:
347 5fb01677 Carles Marti
            return True
348 e8bebcca Carles Marti
349 e8bebcca Carles Marti
    return False
350 e8bebcca Carles Marti
351 e8bebcca Carles Marti
352 e8bebcca Carles Marti
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
353 e8bebcca Carles Marti
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
354 9cd032cf Carles Marti
                 coll_coeff, height=2.5, excl_atom=False):
355 9cd032cf Carles Marti
    # TODO Rename this function
356 e8bebcca Carles Marti
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
357 e8bebcca Carles Marti
    small rotations.
358 e8bebcca Carles Marti

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

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

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

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

568 7dd94df7 Carles Marti
    Carles modification of chemcat's transform_adsorbate to work with
569 7dd94df7 Carles Marti
    coordinates instead of ase.Atom
570 7dd94df7 Carles Marti
    Parameters:
571 7dd94df7 Carles Marti
        molecule (ase.Atoms): The molecule to adsorb.
572 7dd94df7 Carles Marti

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

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

578 7dd94df7 Carles Marti
        ctr2_mol (int/list): The position of a second center of the
579 7dd94df7 Carles Marti
        adsorbate used to define the adsorption bond angle, and the dihedral
580 7dd94df7 Carles Marti
        adsorption angle.
581 7dd94df7 Carles Marti

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

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

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

591 7dd94df7 Carles Marti
        bond_vector (numpy.ndarray): The adsorption bond desired.
592 7dd94df7 Carles Marti
            Details: offset = vect(atom1_surf, atom1_mol)
593 7dd94df7 Carles Marti

594 7dd94df7 Carles Marti
        bond_angle_target (float or int): The adsorption bond angle desired (in
595 7dd94df7 Carles Marti
            degrees).
596 7dd94df7 Carles Marti
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
597 7dd94df7 Carles Marti

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

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

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

834 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb (adsorbate).
835 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the molecule
836 5261a07f Carles Marti
    @param ctr1_mol: The index/es of the center in the adsorbate to use as
837 5261a07f Carles Marti
        adsorption center.
838 5261a07f Carles Marti
    @param ctr2_mol: The index/es of the center in the adsorbate to use for the
839 5261a07f Carles Marti
        definition of the surf-adsorbate angle, surf-adsorbate dihedral angle
840 5261a07f Carles Marti
        and adsorbate dihedral angle.
841 5261a07f Carles Marti
    @param ctr3_mol: The index/es of the center in the adsorbate to use for the
842 5261a07f Carles Marti
        definition of the adsorbate dihedral angle.
843 5261a07f Carles Marti
    @param ctr1_surf: The index/es of the center in the surface to use as
844 5261a07f Carles Marti
        adsorption center.
845 5261a07f Carles Marti
    @param ctr2_surf: The index/es of the center in the surface to use for the
846 5261a07f Carles Marti
        definition of the surf-adsorbate dihedral angle.
847 5261a07f Carles Marti
    @param num_pts: Number on which to sample Euler angles.
848 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision
849 5261a07f Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
850 5261a07f Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
851 5261a07f Carles Marti
        considered as atomic collision.
852 5261a07f Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
853 5261a07f Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
854 5261a07f Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
855 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
856 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
857 c25aa299 Carles Marti
    @param max_hel: Maximum value for sampling the helicopter
858 5261a07f Carles Marti
        (surf-adsorbate dihedral) angle.
859 fe91ddb2 Carles Marti
    @param height: Height on which to try adsorption.
860 9cd032cf Carles Marti
    @param excl_atom: Whether to exclude the adsorption center in the
861 9cd032cf Carles Marti
        molecule in the collision detection.
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 fe91ddb2 Carles Marti
                                                     molec_nghbs, coll_coeff,
888 9cd032cf Carles Marti
                                                     height, excl_atom)
889 bb387578 Carles Marti
                slab_molec.info = {**slab_molec.info, **new_molec.info}
890 cf40df1b Carles Marti
                if not collision and \
891 cf40df1b Carles Marti
                        not any([np.allclose(slab_molec.positions,
892 cf40df1b Carles Marti
                                             conf.positions)
893 cf40df1b Carles Marti
                                 for conf in slab_ads_list]):
894 cf40df1b Carles Marti
                    slab_ads_list.append(slab_molec)
895 cf40df1b Carles Marti
                    if h_donor is not False:
896 cf40df1b Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
897 cf40df1b Carles Marti
                                                          h_acceptor,
898 cf40df1b Carles Marti
                                                          len(slab)))
899 cf40df1b Carles Marti
900 cf40df1b Carles Marti
    return slab_ads_list
901 f3d1e601 Carles Marti
902 f3d1e601 Carles Marti
903 7dd94df7 Carles Marti
def adsorb_confs(conf_list, surf, inp_vars):
904 a5cc42ff Carles Marti
    """Generates a number of adsorbate-surface structure coordinates.
905 a5cc42ff Carles Marti

906 a5cc42ff Carles Marti
    Given a list of conformers, a surface, a list of atom indices (or list of
907 a5cc42ff Carles Marti
    list of indices) of both the surface and the adsorbate, it generates a
908 a5cc42ff Carles Marti
    number of adsorbate-surface structures for every possible combination of
909 a5cc42ff Carles Marti
    them at different orientations.
910 a5cc42ff Carles Marti
    @param conf_list: list of ase.Atoms of the different conformers
911 a5cc42ff Carles Marti
    @param surf: the ase.Atoms object of the surface
912 7dd94df7 Carles Marti
    @param inp_vars: Calculation parameters from input file.
913 a5cc42ff Carles Marti
    @return: list of ase.Atoms for the adsorbate-surface structures
914 a5cc42ff Carles Marti
    """
915 9cd032cf Carles Marti
    from copy import deepcopy
916 bb55f47c Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
917 7dd94df7 Carles Marti
    molec_ctrs = inp_vars['molec_ctrs']
918 7dd94df7 Carles Marti
    sites = inp_vars['sites']
919 7dd94df7 Carles Marti
    angles = inp_vars['set_angles']
920 7dd94df7 Carles Marti
    num_pts = inp_vars['sample_points_per_angle']
921 39df9c43 Carles Marti
    inp_norm_vect = inp_vars['surf_norm_vect']
922 7dd94df7 Carles Marti
    min_coll_height = inp_vars['min_coll_height']
923 7dd94df7 Carles Marti
    coll_coeff = inp_vars['collision_threshold']
924 9cd032cf Carles Marti
    exclude_ads_ctr = inp_vars['exclude_ads_ctr']
925 c25aa299 Carles Marti
    h_donor = inp_vars['h_donor']
926 c25aa299 Carles Marti
    h_acceptor = inp_vars['h_acceptor']
927 fe91ddb2 Carles Marti
    height = inp_vars['adsorption_height']
928 7dd94df7 Carles Marti
929 bc703cab Carles Marti
    if inp_vars['pbc_cell'] is not False:
930 bc703cab Carles Marti
        surf.set_pbc(True)
931 bc703cab Carles Marti
        surf.set_cell(inp_vars['pbc_cell'])
932 bc703cab Carles Marti
933 cf40df1b Carles Marti
    surf_ads_list = []
934 ae097639 Carles Martí
    sites_coords = [get_atom_coords(surf, site) for site in sites]
935 e8bebcca Carles Marti
    if coll_coeff is not False:
936 5fb01677 Carles Marti
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
937 5fb01677 Carles Marti
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
938 5fb01677 Carles Marti
    else:
939 5fb01677 Carles Marti
        surf_nghbs = 0
940 7dd94df7 Carles Marti
    for i, conf in enumerate(conf_list):
941 ae097639 Carles Martí
        molec_ctr_coords = [get_atom_coords(conf, ctr) for ctr in molec_ctrs]
942 bc703cab Carles Marti
        if inp_vars['pbc_cell'] is not False:
943 bc703cab Carles Marti
            conf.set_pbc(True)
944 bc703cab Carles Marti
            conf.set_cell(inp_vars['pbc_cell'])
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 77f8c47b Carles Martí
                if np.isnan(norm_vect).any():
950 77f8c47b Carles Martí
                    logger.warning(f"Could not compute the normal vector to "
951 77f8c47b Carles Martí
                                   f"site '{sites[s]}'. Skipping site.")
952 77f8c47b Carles Martí
                    continue
953 39df9c43 Carles Marti
            else:
954 39df9c43 Carles Marti
                norm_vect = inp_norm_vect
955 7dd94df7 Carles Marti
            for c, ctr in enumerate(molec_ctr_coords):
956 9cd032cf Carles Marti
                if exclude_ads_ctr and isinstance(molec_ctrs[c], int):
957 9cd032cf Carles Marti
                    exclude_atom = molec_ctrs[c]
958 9cd032cf Carles Marti
                else:
959 9cd032cf Carles Marti
                    exclude_atom = False
960 9cd032cf Carles Marti
                    if exclude_ads_ctr and not isinstance(molec_ctrs[c], int):
961 9cd032cf Carles Marti
                        logger.warning("'exclude_ads_ctr' only works for atomic"
962 9cd032cf Carles Marti
                                       "centers and not for many-atoms "
963 9cd032cf Carles Marti
                                       f"barycenters. {molec_ctrs[c]} are not "
964 9cd032cf Carles Marti
                                       f"going to be excluded from collison.")
965 9cd032cf Carles Marti
                if coll_coeff and exclude_atom is not False:
966 9cd032cf Carles Marti
                    conf_wo_ctr = deepcopy(conf)
967 9cd032cf Carles Marti
                    del conf_wo_ctr[exclude_atom]
968 9cd032cf Carles Marti
                    conf_cutoffs = natural_cutoffs(conf_wo_ctr, mult=coll_coeff)
969 9cd032cf Carles Marti
                    molec_nghbs = len(neighbor_list("i", conf_wo_ctr,
970 9cd032cf Carles Marti
                                                    conf_cutoffs))
971 9cd032cf Carles Marti
                elif coll_coeff and exclude_atom is False:
972 9cd032cf Carles Marti
                    conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
973 9cd032cf Carles Marti
                    molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
974 9cd032cf Carles Marti
                else:
975 9cd032cf Carles Marti
                    molec_nghbs = 0
976 7dd94df7 Carles Marti
                if angles == 'euler':
977 cf40df1b Carles Marti
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
978 cf40df1b Carles Marti
                                                   num_pts, min_coll_height,
979 cf40df1b Carles Marti
                                                   coll_coeff, norm_vect,
980 cf40df1b Carles Marti
                                                   surf_nghbs, molec_nghbs,
981 9cd032cf Carles Marti
                                                   h_donor, h_acceptor, height,
982 9cd032cf Carles Marti
                                                   exclude_atom))
983 609104e3 Carles Marti
                elif angles == 'internal':
984 7dd94df7 Carles Marti
                    mol_ctr1 = molec_ctrs[c]
985 7dd94df7 Carles Marti
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
986 7dd94df7 Carles Marti
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
987 7dd94df7 Carles Marti
                    surf_ctr1 = sites[s]
988 7dd94df7 Carles Marti
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
989 a98d4172 Carles Marti
                    max_h = inp_vars["max_helic_angle"]
990 d68dd4ad Carles Marti
                    surf_ads_list.extend(ads_internal(conf, surf, mol_ctr1,
991 d68dd4ad Carles Marti
                                                      mol_ctr2, mol_ctr3,
992 d68dd4ad Carles Marti
                                                      surf_ctr1, surf_ctr2,
993 d68dd4ad Carles Marti
                                                      num_pts, min_coll_height,
994 d68dd4ad Carles Marti
                                                      coll_coeff, norm_vect,
995 d68dd4ad Carles Marti
                                                      surf_nghbs, molec_nghbs,
996 d68dd4ad Carles Marti
                                                      h_donor, h_acceptor,
997 9cd032cf Carles Marti
                                                      max_h, height,
998 9cd032cf Carles Marti
                                                      exclude_atom))
999 cf40df1b Carles Marti
    return surf_ads_list
1000 f3d1e601 Carles Marti
1001 f3d1e601 Carles Marti
1002 4614bb6a Carles
def run_screening(inp_vars):
1003 91ae8d86 Carles Marti
    """Carries out the screening of adsorbate structures on a surface.
1004 e07c09eb Carles

1005 e07c09eb Carles
    @param inp_vars: Calculation parameters from input file.
1006 e07c09eb Carles
    """
1007 e07c09eb Carles
    import os
1008 57e3a8c7 Carles Marti
    import random
1009 1d8c374e Carles Martí
    from modules.formats import collect_confs, adapt_format
1010 cf8fe0e3 Carles Marti
    from modules.calculation import run_calc, check_finished_calcs
1011 e07c09eb Carles
1012 76f4ac19 Carles Marti
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
1013 76f4ac19 Carles Marti
                ' structures.')
1014 b75bf97d Carles Marti
    if inp_vars['use_molec_file']:
1015 7f5fcd41 Carles Martí
        selected_confs = [adapt_format('ase', inp_vars['use_molec_file'],
1016 7f5fcd41 Carles Martí
                                       inp_vars['special_atoms'])]
1017 b75bf97d Carles Marti
        logger.info(f"Using '{inp_vars['use_molec_file']}' as only conformer.")
1018 b75bf97d Carles Marti
    else:
1019 b75bf97d Carles Marti
        if not os.path.isdir("isolated"):
1020 b75bf97d Carles Marti
            err = "'isolated' directory not found. It is needed in order to " \
1021 b75bf97d Carles Marti
                  "carry out the screening of structures to be adsorbed"
1022 b75bf97d Carles Marti
            logger.error(err)
1023 b75bf97d Carles Marti
            raise FileNotFoundError(err)
1024 e07c09eb Carles
1025 1d8c374e Carles Martí
        finished_calcs, failed_calcs = check_finished_calcs('isolated',
1026 1d8c374e Carles Martí
                                                            inp_vars['code'])
1027 1d8c374e Carles Martí
        if not finished_calcs:
1028 b75bf97d Carles Marti
            err_msg = "No calculations on 'isolated' finished normally."
1029 b75bf97d Carles Marti
            logger.error(err_msg)
1030 b75bf97d Carles Marti
            raise FileNotFoundError(err_msg)
1031 b75bf97d Carles Marti
1032 1d8c374e Carles Martí
        logger.info(f"Found {len(finished_calcs)} structures of isolated "
1033 b75bf97d Carles Marti
                    f"conformers whose calculation finished normally.")
1034 b75bf97d Carles Marti
        if len(failed_calcs) != 0:
1035 b75bf97d Carles Marti
            logger.warning(
1036 b75bf97d Carles Marti
                f"Found {len(failed_calcs)} calculations more that "
1037 b75bf97d Carles Marti
                f"did not finish normally: {failed_calcs}. \n"
1038 b75bf97d Carles Marti
                f"Using only the ones that finished normally: "
1039 1d8c374e Carles Martí
                f"{finished_calcs}.")
1040 1d8c374e Carles Martí
1041 1d8c374e Carles Martí
        conf_list = collect_confs(finished_calcs, inp_vars['code'], 'isolated',
1042 1d8c374e Carles Martí
                                  inp_vars['special_atoms'])
1043 1d8c374e Carles Martí
        selected_confs = select_confs(conf_list, inp_vars['select_magns'],
1044 1d8c374e Carles Martí
                                      inp_vars['confs_per_magn'])
1045 90819cc3 Carles Marti
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
1046 bb387578 Carles Marti
    surf.info = {}
1047 cf40df1b Carles Marti
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
1048 7d97341d Carles Marti
    if len(surf_ads_list) > inp_vars['max_structures']:
1049 57e3a8c7 Carles Marti
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
1050 a44ad3c2 Carles Martí
    elif len(surf_ads_list) == 0:
1051 a44ad3c2 Carles Martí
        err_msg = "No configurations were generated: Check the parameters in" \
1052 a44ad3c2 Carles Martí
                  "dockonsurf.inp"
1053 a44ad3c2 Carles Martí
        logger.error(err_msg)
1054 a44ad3c2 Carles Martí
        raise ValueError(err_msg)
1055 cf40df1b Carles Marti
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
1056 cf40df1b Carles Marti
                f'configurations to carry out a calculation of.')
1057 d9167fea Carles Marti
1058 f3d1e601 Carles Marti
    run_calc('screening', inp_vars, surf_ads_list)
1059 14f39d2a Carles Marti
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1060 14f39d2a Carles Marti
                ' structures section.')