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

dockonsurf / modules / screening.py @ 587dca22

Historique | Voir | Annoter | Télécharger (43,79 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 12a12bdd Carles Marti
    from 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 d6da8693 Carles Marti
    from 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 dadc6016 Carles Marti
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
180 1d22a086 Carles Marti
                  norm_vect=(0, 0, 1)):
181 1d22a086 Carles Marti
    """Add an adsorbate to a surface.
182 1d22a086 Carles Marti

183 1d22a086 Carles Marti
    This function extends the functionality of ase.build.add_adsorbate
184 1d22a086 Carles Marti
    (https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
185 1d22a086 Carles Marti
    by enabling to change the z coordinate and the axis perpendicular to the
186 1d22a086 Carles Marti
    surface.
187 1d22a086 Carles Marti
    @param slab: ase.Atoms object containing the coordinates of the surface
188 1d22a086 Carles Marti
    @param adsorbate: ase.Atoms object containing the coordinates of the
189 1d22a086 Carles Marti
        adsorbate.
190 dadc6016 Carles Marti
    @param site_coord: The coordinates of the adsorption site on the surface.
191 dadc6016 Carles Marti
    @param ctr_coord: The coordinates of the adsorption center in the molecule.
192 dadc6016 Carles Marti
    @param height: The height above the surface where to adsorb.
193 1d22a086 Carles Marti
    @param offset: Offsets the adsorbate by a number of unit cells. Mostly
194 1d22a086 Carles Marti
        useful when adding more than one adsorbate.
195 1d22a086 Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
196 1d22a086 Carles Marti
    """
197 36d92f4f Carles Marti
    from copy import deepcopy
198 1d22a086 Carles Marti
    info = slab.info.get('adsorbate_info', {})
199 dadc6016 Carles Marti
    pos = np.array([0.0, 0.0, 0.0])  # part of absolute coordinates
200 1d22a086 Carles Marti
    spos = np.array([0.0, 0.0, 0.0])  # part relative to unit cell
201 dadc6016 Carles Marti
    norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect)
202 1d22a086 Carles Marti
    if offset is not None:
203 1d22a086 Carles Marti
        spos += np.asarray(offset, float)
204 dadc6016 Carles Marti
    if isinstance(site_coord, str):
205 1d22a086 Carles Marti
        # A site-name:
206 1d22a086 Carles Marti
        if 'sites' not in info:
207 1d22a086 Carles Marti
            raise TypeError('If the atoms are not made by an ase.build '
208 1d22a086 Carles Marti
                            'function, position cannot be a name.')
209 dadc6016 Carles Marti
        if site_coord not in info['sites']:
210 dadc6016 Carles Marti
            raise TypeError('Adsorption site %s not supported.' % site_coord)
211 dadc6016 Carles Marti
        spos += info['sites'][site_coord]
212 1d22a086 Carles Marti
    else:
213 dadc6016 Carles Marti
        pos += site_coord
214 1d22a086 Carles Marti
    if 'cell' in info:
215 1d22a086 Carles Marti
        cell = info['cell']
216 1d22a086 Carles Marti
    else:
217 1d22a086 Carles Marti
        cell = slab.get_cell()
218 1d22a086 Carles Marti
    pos += np.dot(spos, cell)
219 1d22a086 Carles Marti
    # Convert the adsorbate to an Atoms object
220 1d22a086 Carles Marti
    if isinstance(adsorbate, ase.Atoms):
221 36d92f4f Carles Marti
        ads = deepcopy(adsorbate)
222 1d22a086 Carles Marti
    elif isinstance(adsorbate, ase.Atom):
223 1d22a086 Carles Marti
        ads = ase.Atoms([adsorbate])
224 1d22a086 Carles Marti
    else:
225 1d22a086 Carles Marti
        # Assume it is a string representing a single Atom
226 1d22a086 Carles Marti
        ads = ase.Atoms([ase.Atom(adsorbate)])
227 dadc6016 Carles Marti
    pos += height * norm_vect_u
228 1d22a086 Carles Marti
    # Move adsorbate into position
229 dadc6016 Carles Marti
    ads.translate(pos - ctr_coord)
230 1d22a086 Carles Marti
    # Attach the adsorbate
231 1d22a086 Carles Marti
    slab.extend(ads)
232 1d22a086 Carles Marti
233 1d22a086 Carles Marti
234 a4b57124 Carles Marti
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
235 e8bebcca Carles Marti
                    nn_molec=0, coll_coeff=1.2):
236 5f3f4b69 Carles Marti
    """Checks whether a slab and a molecule collide or not.
237 5f3f4b69 Carles Marti

238 5f3f4b69 Carles Marti
    @param slab_molec: The system of adsorbate-slab for which to detect if there
239 5f3f4b69 Carles Marti
        are collisions.
240 5f3f4b69 Carles Marti
    @param nn_slab: Number of neigbors in the surface.
241 5f3f4b69 Carles Marti
    @param nn_molec: Number of neighbors in the molecule.
242 5f3f4b69 Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
243 5f3f4b69 Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
244 5f3f4b69 Carles Marti
        considered as atomic collision.
245 5f3f4b69 Carles Marti
    @param slab_num_atoms: Number of atoms of the bare slab.
246 5f3f4b69 Carles Marti
    @param min_height: The minimum height atoms can have to not be considered as
247 5f3f4b69 Carles Marti
        colliding.
248 5f3f4b69 Carles Marti
    @param vect: The vector perpendicular to the slab.
249 5f3f4b69 Carles Marti
    @return: bool, whether the surface and the molecule collide.
250 e8bebcca Carles Marti
    """
251 5f3f4b69 Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
252 e8bebcca Carles Marti
253 e8bebcca Carles Marti
    # Check structure overlap by height
254 5fb01677 Carles Marti
    if min_height is not False:
255 a4b57124 Carles Marti
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
256 a4b57124 Carles Marti
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
257 e8bebcca Carles Marti
        if vect.tolist() not in cart_axes:
258 e8bebcca Carles Marti
            err_msg = "'min_coll_height' option is only implemented for " \
259 e8bebcca Carles Marti
                      "'surf_norm_vect' to be one of the x, y or z axes. "
260 e8bebcca Carles Marti
            logger.error(err_msg)
261 e8bebcca Carles Marti
            raise ValueError(err_msg)
262 e8bebcca Carles Marti
        for atom in slab_molec[slab_num_atoms:]:
263 e8bebcca Carles Marti
            for i, coord in enumerate(vect):
264 e8bebcca Carles Marti
                if coord == 0:
265 e8bebcca Carles Marti
                    continue
266 e8bebcca Carles Marti
                if atom.position[i] * coord < min_height * coord:
267 e8bebcca Carles Marti
                    return True
268 e8bebcca Carles Marti
269 e8bebcca Carles Marti
    # Check structure overlap by sphere collision
270 e8bebcca Carles Marti
    if coll_coeff is not False:
271 5fb01677 Carles Marti
        slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
272 b4b2f307 Carles Marti
        slab_molec_nghbs = len(
273 b4b2f307 Carles Marti
            neighbor_list("i", slab_molec, slab_molec_cutoffs))
274 5fb01677 Carles Marti
        if slab_molec_nghbs > nn_slab + nn_molec:
275 5fb01677 Carles Marti
            return True
276 e8bebcca Carles Marti
277 e8bebcca Carles Marti
    return False
278 e8bebcca Carles Marti
279 e8bebcca Carles Marti
280 e8bebcca Carles Marti
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
281 e8bebcca Carles Marti
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
282 e8bebcca Carles Marti
                 coll_coeff, height=2.5):
283 e8bebcca Carles Marti
    # TODO Rethink this function
284 e8bebcca Carles Marti
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
285 e8bebcca Carles Marti
    small rotations.
286 e8bebcca Carles Marti

287 e8bebcca Carles Marti
    @param molec: ase.Atoms object of the molecule to adsorb
288 e8bebcca Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the
289 e8bebcca Carles Marti
        molecule
290 e8bebcca Carles Marti
    @param ctr_coord: The coordinates of the molecule to use as adsorption
291 e8bebcca Carles Marti
        center.
292 e8bebcca Carles Marti
    @param site_coord: The coordinates of the surface on which to adsorb the
293 e8bebcca Carles Marti
        molecule
294 e8bebcca Carles Marti
    @param num_pts: Number on which to sample Euler angles.
295 e8bebcca Carles Marti
    @param min_coll_height: The lowermost height for which to detect a collision
296 e8bebcca Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
297 e8bebcca Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
298 e8bebcca Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
299 e8bebcca Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
300 e8bebcca Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
301 e8bebcca Carles Marti
        considered as atomic collision.
302 e8bebcca Carles Marti
    @param height: Height on which to try adsorption
303 e8bebcca Carles Marti
    @return collision: bool, whether the structure generated has collisions
304 e8bebcca Carles Marti
        between slab and adsorbate.
305 e8bebcca Carles Marti
    """
306 e8bebcca Carles Marti
    from copy import deepcopy
307 e8bebcca Carles Marti
    slab_num_atoms = len(slab)
308 e8bebcca Carles Marti
    slab_molec = []
309 e8bebcca Carles Marti
    collision = True
310 e8bebcca Carles Marti
    max_corr = 6  # Should be an even number
311 e8bebcca Carles Marti
    d_angle = 180 / ((max_corr / 2.0) * num_pts)
312 e8bebcca Carles Marti
    num_corr = 0
313 e8bebcca Carles Marti
    while collision and num_corr <= max_corr:
314 e8bebcca Carles Marti
        k = num_corr * (-1) ** num_corr
315 e8bebcca Carles Marti
        slab_molec = deepcopy(slab)
316 e8bebcca Carles Marti
        molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
317 e8bebcca Carles Marti
                           center=ctr_coord)
318 e8bebcca Carles Marti
        add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height,
319 e8bebcca Carles Marti
                      norm_vect=norm_vect)
320 e8bebcca Carles Marti
        collision = check_collision(slab_molec, slab_num_atoms, min_coll_height,
321 e8bebcca Carles Marti
                                    norm_vect, slab_nghbs, molec_nghbs,
322 e8bebcca Carles Marti
                                    coll_coeff)
323 e8bebcca Carles Marti
        num_corr += 1
324 e8bebcca Carles Marti
    return slab_molec, collision
325 5f3f4b69 Carles Marti
326 5f3f4b69 Carles Marti
327 c25aa299 Carles Marti
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
328 c25aa299 Carles Marti
                 neigh_cutoff=1):
329 b4b2f307 Carles Marti
    # TODO rethink
330 91ae8d86 Carles Marti
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
331 b4b2f307 Carles Marti

332 91ae8d86 Carles Marti
    Tries to dissociate a H atom from the molecule and adsorb in on top of the
333 91ae8d86 Carles Marti
    surface if the distance is shorter than two times the neigh_cutoff value.
334 b4b2f307 Carles Marti
    @param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
335 b4b2f307 Carles Marti
    @param h_idx: The index of the hydrogen atom to carry out adsorption of.
336 b4b2f307 Carles Marti
    @param num_atoms_slab: The number of atoms of the slab without adsorbate.
337 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
338 b4b2f307 Carles Marti
    @param neigh_cutoff: half the maximum distance between the surface and the
339 b4b2f307 Carles Marti
        H for it to carry out dissociation.
340 b4b2f307 Carles Marti
    @return: An ase.Atoms object of the system adsorbate-surface with H
341 b4b2f307 Carles Marti
    """
342 b4b2f307 Carles Marti
    from copy import deepcopy
343 b4b2f307 Carles Marti
    from ase.neighborlist import NeighborList
344 b4b2f307 Carles Marti
    slab_molec = deepcopy(slab_molec_orig)
345 b4b2f307 Carles Marti
    cutoffs = len(slab_molec) * [neigh_cutoff]
346 c25aa299 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True, skin=0.0)
347 b4b2f307 Carles Marti
    nl.update(slab_molec)
348 b4b2f307 Carles Marti
    surf_h_vect = np.array([np.infty] * 3)
349 c25aa299 Carles Marti
    if h_acceptor == 'all':
350 c25aa299 Carles Marti
        h_acceptor = list(range(num_atoms_slab))
351 b4b2f307 Carles Marti
    for neigh_idx in nl.get_neighbors(h_idx)[0]:
352 c25aa299 Carles Marti
        for elm in h_acceptor:
353 c25aa299 Carles Marti
            if isinstance(elm, int):
354 c25aa299 Carles Marti
                if neigh_idx == elm and neigh_idx < num_atoms_slab:
355 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
356 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
357 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
358 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
359 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
360 c25aa299 Carles Marti
            else:
361 c25aa299 Carles Marti
                if slab_molec[neigh_idx].symbol == elm \
362 c25aa299 Carles Marti
                        and neigh_idx < num_atoms_slab:
363 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
364 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
365 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
366 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
367 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
368 c25aa299 Carles Marti
369 b4b2f307 Carles Marti
    if np.linalg.norm(surf_h_vect) != np.infty:
370 b4b2f307 Carles Marti
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
371 b4b2f307 Carles Marti
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
372 b4b2f307 Carles Marti
        return slab_molec
373 b4b2f307 Carles Marti
374 b4b2f307 Carles Marti
375 c25aa299 Carles Marti
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
376 b4b2f307 Carles Marti
    # TODO multiple dissociation
377 b4b2f307 Carles Marti
    """Decides which H atoms to dissociate according to a list of atoms.
378 b4b2f307 Carles Marti

379 b4b2f307 Carles Marti
    Given a list of chemical symbols or atom indices it checks for every atom
380 b4b2f307 Carles Marti
    or any of its neighbor if it's a H and calls dissociate_h to try to carry
381 b4b2f307 Carles Marti
    out dissociation of that H. For atom indices, it checks both whether
382 b4b2f307 Carles Marti
    the atom index or its neighbors are H, for chemical symbols, it only checks
383 b4b2f307 Carles Marti
    if there is a neighbor H.
384 b4b2f307 Carles Marti
    @param slab_molec: The ase.Atoms object of the system adsorbate-slab.
385 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
386 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
387 c25aa299 Carles Marti
    @param num_atoms_slab: Number of atoms of the bare slab.
388 b4b2f307 Carles Marti
    @return:
389 b4b2f307 Carles Marti
    """
390 b4b2f307 Carles Marti
    from ase.neighborlist import natural_cutoffs, NeighborList
391 b4b2f307 Carles Marti
    molec = slab_molec[num_atoms_slab:]
392 b4b2f307 Carles Marti
    cutoffs = natural_cutoffs(molec)
393 b4b2f307 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True)
394 b4b2f307 Carles Marti
    nl.update(molec)
395 b4b2f307 Carles Marti
    disso_structs = []
396 c25aa299 Carles Marti
    for el in h_donor:
397 b4b2f307 Carles Marti
        if isinstance(el, int):
398 b4b2f307 Carles Marti
            if molec[el].symbol == 'H':
399 b4b2f307 Carles Marti
                disso_struct = dissociate_h(slab_molec, el + num_atoms_slab,
400 c25aa299 Carles Marti
                                            num_atoms_slab, h_acceptor)
401 b4b2f307 Carles Marti
                if disso_struct is not None:
402 b4b2f307 Carles Marti
                    disso_structs.append(disso_struct)
403 b4b2f307 Carles Marti
            else:
404 b4b2f307 Carles Marti
                for neigh_idx in nl.get_neighbors(el)[0]:
405 b4b2f307 Carles Marti
                    if molec[neigh_idx].symbol == 'H':
406 b4b2f307 Carles Marti
                        disso_struct = dissociate_h(slab_molec, neigh_idx +
407 b4b2f307 Carles Marti
                                                    num_atoms_slab,
408 c25aa299 Carles Marti
                                                    num_atoms_slab, h_acceptor)
409 b4b2f307 Carles Marti
                        if disso_struct is not None:
410 b4b2f307 Carles Marti
                            disso_structs.append(disso_struct)
411 b4b2f307 Carles Marti
        else:
412 b4b2f307 Carles Marti
            for atom in molec:
413 b4b2f307 Carles Marti
                if atom.symbol.lower() == el.lower():
414 b4b2f307 Carles Marti
                    for neigh_idx in nl.get_neighbors(atom.index)[0]:
415 b4b2f307 Carles Marti
                        if molec[neigh_idx].symbol == 'H':
416 5261a07f Carles Marti
                            disso_struct = dissociate_h(slab_molec, neigh_idx
417 b4b2f307 Carles Marti
                                                        + num_atoms_slab,
418 c25aa299 Carles Marti
                                                        num_atoms_slab,
419 c25aa299 Carles Marti
                                                        h_acceptor)
420 b4b2f307 Carles Marti
                            if disso_struct is not None:
421 b4b2f307 Carles Marti
                                disso_structs.append(disso_struct)
422 b4b2f307 Carles Marti
    return disso_structs
423 b4b2f307 Carles Marti
424 b4b2f307 Carles Marti
425 3ab0865c Carles Marti
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
426 b4b2f307 Carles Marti
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
427 c25aa299 Carles Marti
              h_donor, h_acceptor):
428 3ab0865c Carles Marti
    """Generates adsorbate-surface structures by sampling over Euler angles.
429 3ab0865c Carles Marti

430 3ab0865c Carles Marti
    This function generates a number of adsorbate-surface structures at
431 3ab0865c Carles Marti
    different orientations of the adsorbate sampled at multiple Euler (zxz)
432 3ab0865c Carles Marti
    angles.
433 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb.
434 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the
435 5261a07f Carles Marti
        molecule.
436 3ab0865c Carles Marti
    @param ctr_coord: The coordinates of the molecule to use as adsorption
437 3ab0865c Carles Marti
        center.
438 3ab0865c Carles Marti
    @param site_coord: The coordinates of the surface on which to adsorb the
439 5261a07f Carles Marti
        molecule.
440 3ab0865c Carles Marti
    @param num_pts: Number on which to sample Euler angles.
441 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision.
442 3ab0865c Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
443 3ab0865c Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
444 3ab0865c Carles Marti
        considered as atomic collision.
445 3ab0865c Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
446 3ab0865c Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
447 3ab0865c Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
448 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
449 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
450 3ab0865c Carles Marti
    @return: list of ase.Atoms object conatining all the orientations of a given
451 5261a07f Carles Marti
        conformer.
452 3ab0865c Carles Marti
    """
453 3ab0865c Carles Marti
    from copy import deepcopy
454 3ab0865c Carles Marti
    slab_ads_list = []
455 3ab0865c Carles Marti
    # rotation around z
456 3ab0865c Carles Marti
    for alpha in np.arange(0, 360, 360 / num_pts):
457 3ab0865c Carles Marti
        # rotation around x'
458 3ab0865c Carles Marti
        for beta in np.arange(0, 180, 180 / num_pts):
459 3ab0865c Carles Marti
            # rotation around z'
460 3ab0865c Carles Marti
            for gamma in np.arange(0, 360, 360 / num_pts):
461 3ab0865c Carles Marti
                molec = deepcopy(orig_molec)
462 3ab0865c Carles Marti
                molec.euler_rotate(alpha, beta, gamma, center=ctr_coord)
463 3ab0865c Carles Marti
                slab_molec, collision = correct_coll(molec, slab,
464 5fb01677 Carles Marti
                                                     ctr_coord, site_coord,
465 5fb01677 Carles Marti
                                                     num_pts, min_coll_height,
466 5fb01677 Carles Marti
                                                     norm_vect,
467 5fb01677 Carles Marti
                                                     slab_nghbs, molec_nghbs,
468 5fb01677 Carles Marti
                                                     coll_coeff)
469 5261a07f Carles Marti
                if not collision and not any([np.allclose(slab_molec.positions,
470 5261a07f Carles Marti
                                                          conf.positions)
471 5261a07f Carles Marti
                                              for conf in slab_ads_list]):
472 3ab0865c Carles Marti
                    slab_ads_list.append(slab_molec)
473 c25aa299 Carles Marti
                    if h_donor is not False:
474 c25aa299 Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
475 c25aa299 Carles Marti
                                                          h_acceptor,
476 c25aa299 Carles Marti
                                                          len(slab)))
477 3ab0865c Carles Marti
478 3ab0865c Carles Marti
    return slab_ads_list
479 f3d1e601 Carles Marti
480 f3d1e601 Carles Marti
481 7dd94df7 Carles Marti
def chemcat_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
482 7dd94df7 Carles Marti
                   ctr2_surf, bond_vector, bond_angle_target,
483 7dd94df7 Carles Marti
                   dihedral_angle_target=None, mol_dihedral_angle_target=None):
484 7dd94df7 Carles Marti
    """Performs translation and rotation of an adsorbate defined by an
485 7dd94df7 Carles Marti
    adsorption bond length, direction, angle and dihedral angle
486 7dd94df7 Carles Marti

487 7dd94df7 Carles Marti
    Carles modification of chemcat's transform_adsorbate to work with
488 7dd94df7 Carles Marti
    coordinates instead of ase.Atom
489 7dd94df7 Carles Marti
    Parameters:
490 7dd94df7 Carles Marti
        molecule (ase.Atoms): The molecule to adsorb.
491 7dd94df7 Carles Marti

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

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

497 7dd94df7 Carles Marti
        ctr2_mol (int/list): The position of a second center of the
498 7dd94df7 Carles Marti
        adsorbate used to define the adsorption bond angle, and the dihedral
499 7dd94df7 Carles Marti
        adsorption angle.
500 7dd94df7 Carles Marti

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

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

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

510 7dd94df7 Carles Marti
        bond_vector (numpy.ndarray): The adsorption bond desired.
511 7dd94df7 Carles Marti
            Details: offset = vect(atom1_surf, atom1_mol)
512 7dd94df7 Carles Marti

513 7dd94df7 Carles Marti
        bond_angle_target (float or int): The adsorption bond angle desired (in
514 7dd94df7 Carles Marti
            degrees).
515 7dd94df7 Carles Marti
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
516 7dd94df7 Carles Marti

517 7dd94df7 Carles Marti
        dihedral_angle_target (float or int): The dihedral adsorption angle
518 7dd94df7 Carles Marti
            desired (in degrees).
519 7dd94df7 Carles Marti
            Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
520 7dd94df7 Carles Marti
            atom1_mol, atom2_mol)
521 7dd94df7 Carles Marti
                The rotation vector is facing the adsorbate from the surface
522 7dd94df7 Carles Marti
                (i.e. counterclockwise rotation when facing the surface (i.e.
523 7dd94df7 Carles Marti
                view from top))
524 7dd94df7 Carles Marti

525 7dd94df7 Carles Marti
        mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
526 7dd94df7 Carles Marti
            desired (in degrees).
527 7dd94df7 Carles Marti
            Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
528 7dd94df7 Carles Marti
            atom2_mol, atom3_mol)
529 7dd94df7 Carles Marti
                The rotation vector is facing atom2_mol from atom1_mol
530 7dd94df7 Carles Marti

531 7dd94df7 Carles Marti
    Returns:
532 7dd94df7 Carles Marti
        None (the `molecule` object is modified in-place)
533 7dd94df7 Carles Marti
    """
534 7dd94df7 Carles Marti
    vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf,
535 7dd94df7 Carles Marti
                                                                   ctr1_surf)
536 d6da8693 Carles Marti
    vect_inter = get_atom_coords(molecule, ctr1_mol) \
537 d6da8693 Carles Marti
        - get_atom_coords(surf, ctr1_surf)
538 7dd94df7 Carles Marti
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
539 7dd94df7 Carles Marti
                                                                     ctr1_mol)
540 7dd94df7 Carles Marti
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
541 7dd94df7 Carles Marti
                                                                      ctr2_mol)
542 7dd94df7 Carles Marti
543 7dd94df7 Carles Marti
    # Check if dihedral angles can be defined
544 7dd94df7 Carles Marti
    do_dihedral = dihedral_angle_target is not None
545 7dd94df7 Carles Marti
    do_mol_dihedral = mol_dihedral_angle_target is not None
546 7dd94df7 Carles Marti
    dihedral_use_mol2 = False
547 7dd94df7 Carles Marti
    # Check if vect_surf and vect_inter are not aligned
548 7dd94df7 Carles Marti
    if np.allclose(np.cross(vect_surf, vect_inter), 0):
549 7dd94df7 Carles Marti
        logger.warning(
550 7dd94df7 Carles Marti
            "Surface atoms are incompatible with adsorption "
551 7dd94df7 Carles Marti
            "direction/bond. An adsorption dihedral angle cannot be defined.")
552 7dd94df7 Carles Marti
        do_dihedral = False
553 7dd94df7 Carles Marti
    # Check if requested bond angle is not flat
554 7dd94df7 Carles Marti
    if np.isclose((bond_angle_target + 90) % 180 - 90, 0):
555 7dd94df7 Carles Marti
        logger.warning("Requested bond angle is flat. Only a single dihedral "
556 7dd94df7 Carles Marti
                       "angle can be defined (ctr2_surf, ctr1_surf, ctr2_mol, "
557 7dd94df7 Carles Marti
                       "ctr3_mol).")
558 7dd94df7 Carles Marti
        do_mol_dihedral = False
559 7dd94df7 Carles Marti
        dihedral_use_mol2 = True
560 7dd94df7 Carles Marti
        logger.warning("Dihedral adsorption angle rotation will be perfomed "
561 7dd94df7 Carles Marti
                       "with (ctr2_surf, ctr1_surf, ctr2_mol, ctr3_mol).")
562 7dd94df7 Carles Marti
    # Check if vect_mol and vect2_mol are not aligned
563 7dd94df7 Carles Marti
    if np.allclose(np.cross(vect_mol, vect2_mol), 0):
564 7dd94df7 Carles Marti
        logger.warning("Adsorbates atoms are aligned. An adsorbate dihedral "
565 7dd94df7 Carles Marti
                       "angle cannot be defined.")
566 7dd94df7 Carles Marti
        do_mol_dihedral = False
567 7dd94df7 Carles Marti
    if not do_dihedral:
568 7dd94df7 Carles Marti
        logger.warning("Dihedral adsorption angle rotation will not be "
569 7dd94df7 Carles Marti
                       "performed.")
570 7dd94df7 Carles Marti
    if not do_mol_dihedral:
571 7dd94df7 Carles Marti
        logger.warning("Adsorbate dihedral angle rotation will not be "
572 7dd94df7 Carles Marti
                       "performed.")
573 7dd94df7 Carles Marti
574 7dd94df7 Carles Marti
    ###########################
575 7dd94df7 Carles Marti
    #       Translation       #
576 7dd94df7 Carles Marti
    ###########################
577 7dd94df7 Carles Marti
578 7dd94df7 Carles Marti
    # Compute and apply translation of adsorbate
579 7dd94df7 Carles Marti
    translation = bond_vector - vect_inter
580 7dd94df7 Carles Marti
    molecule.translate(translation)
581 7dd94df7 Carles Marti
582 7dd94df7 Carles Marti
    # Update adsorption bond
583 d6da8693 Carles Marti
    vect_inter = get_atom_coords(molecule, ctr1_mol) - \
584 d6da8693 Carles Marti
        get_atom_coords(surf, ctr1_surf)
585 7dd94df7 Carles Marti
586 7dd94df7 Carles Marti
    # Check if translation was successful
587 7dd94df7 Carles Marti
    if np.allclose(vect_inter, bond_vector):
588 7dd94df7 Carles Marti
        pass  # print("Translation successfully applied (error: ~ {:.5g} unit "
589 7dd94df7 Carles Marti
        # "length)".format(np.linalg.norm(vect_inter - bond_vector)))
590 7dd94df7 Carles Marti
    else:
591 7dd94df7 Carles Marti
        err = 'An unknown error occured during the translation'
592 7dd94df7 Carles Marti
        logger.error(err)
593 7dd94df7 Carles Marti
        raise AssertionError(err)
594 7dd94df7 Carles Marti
595 7dd94df7 Carles Marti
    ###########################
596 7dd94df7 Carles Marti
    #   Bond angle rotation   #
597 7dd94df7 Carles Marti
    ###########################
598 7dd94df7 Carles Marti
599 7dd94df7 Carles Marti
    # Compute rotation vector
600 7dd94df7 Carles Marti
    rotation_vector = np.cross(-vect_inter, vect_mol)
601 7dd94df7 Carles Marti
    if np.allclose(rotation_vector, 0, atol=1e-5):
602 7dd94df7 Carles Marti
        # If molecular bonds are aligned, any vector orthogonal to vect_inter
603 7dd94df7 Carles Marti
        # can be used Such vector can be found as the orthogonal rejection of
604 7dd94df7 Carles Marti
        # either X-axis, Y-axis or Z-axis onto vect_inter (since they cannot
605 7dd94df7 Carles Marti
        # be all aligned)
606 7dd94df7 Carles Marti
        non_aligned_vector = np.zeros(3)
607 7dd94df7 Carles Marti
        # Select the most orthogonal axis (lowest dot product):
608 7dd94df7 Carles Marti
        non_aligned_vector[np.argmin(np.abs(vect_inter))] = 1
609 7dd94df7 Carles Marti
        rotation_vector = non_aligned_vector - np.dot(non_aligned_vector,
610 7dd94df7 Carles Marti
                                                      vect_inter) / np.dot(
611 7dd94df7 Carles Marti
            vect_inter, vect_inter) * vect_inter
612 7dd94df7 Carles Marti
613 7dd94df7 Carles Marti
    # Retrieve current bond angle
614 7dd94df7 Carles Marti
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
615 7dd94df7 Carles Marti
616 7dd94df7 Carles Marti
    # Apply rotation to reach desired bond_angle
617 7dd94df7 Carles Marti
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
618 7dd94df7 Carles Marti
                    center=get_atom_coords(molecule, ctr1_mol))
619 7dd94df7 Carles Marti
620 7dd94df7 Carles Marti
    # Update molecular bonds
621 7dd94df7 Carles Marti
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
622 7dd94df7 Carles Marti
                                                                     ctr1_mol)
623 7dd94df7 Carles Marti
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
624 7dd94df7 Carles Marti
                                                                      ctr2_mol)
625 7dd94df7 Carles Marti
626 7dd94df7 Carles Marti
    # Check if rotation was successful
627 7dd94df7 Carles Marti
    bond_angle = get_vect_angle(-vect_inter, vect_mol)
628 7dd94df7 Carles Marti
    if np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90, 0,
629 7dd94df7 Carles Marti
                  atol=1e-3) and np.allclose(get_atom_coords(molecule, ctr1_mol)
630 7dd94df7 Carles Marti
                                             - get_atom_coords(surf,
631 7dd94df7 Carles Marti
                                                               ctr1_surf),
632 7dd94df7 Carles Marti
                                             vect_inter):
633 7dd94df7 Carles Marti
        pass  # print("Rotation successfully applied (error: {:.5f}°)".format(
634 7dd94df7 Carles Marti
        # (bond_angle - bond_angle_target + 90) % 180 - 90))
635 7dd94df7 Carles Marti
    else:
636 7dd94df7 Carles Marti
        err = 'An unknown error occured during the rotation'
637 7dd94df7 Carles Marti
        logger.error(err)
638 7dd94df7 Carles Marti
        raise AssertionError(err)
639 7dd94df7 Carles Marti
640 7dd94df7 Carles Marti
    ###########################
641 7dd94df7 Carles Marti
    # Dihedral angle rotation #
642 7dd94df7 Carles Marti
    ###########################
643 7dd94df7 Carles Marti
644 7dd94df7 Carles Marti
    # Perform dihedral rotation if possible
645 7dd94df7 Carles Marti
    if do_dihedral:
646 7dd94df7 Carles Marti
        # Retrieve current dihedral angle (by computing the angle between the
647 7dd94df7 Carles Marti
        # vector rejection of vect_surf and vect_mol onto vect_inter)
648 7dd94df7 Carles Marti
        vect_inter_inner = np.dot(vect_inter, vect_inter)
649 7dd94df7 Carles Marti
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
650 d6da8693 Carles Marti
            vect_inter_inner * vect_inter
651 7dd94df7 Carles Marti
        if dihedral_use_mol2:
652 7dd94df7 Carles Marti
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
653 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
654 7dd94df7 Carles Marti
        else:
655 7dd94df7 Carles Marti
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
656 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
657 7dd94df7 Carles Marti
        dihedral_angle_ini = get_vect_angle(vect_surf_reject, vect_mol_reject,
658 7dd94df7 Carles Marti
                                            vect_inter)
659 7dd94df7 Carles Marti
660 7dd94df7 Carles Marti
        # Apply dihedral rotation along vect_inter
661 7dd94df7 Carles Marti
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
662 7dd94df7 Carles Marti
                        v=vect_inter, center=get_atom_coords(molecule,
663 7dd94df7 Carles Marti
                                                             ctr1_mol))
664 7dd94df7 Carles Marti
665 7dd94df7 Carles Marti
        # Update molecular bonds
666 7dd94df7 Carles Marti
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
667 d6da8693 Carles Marti
            get_atom_coords(molecule, ctr1_mol)
668 7dd94df7 Carles Marti
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
669 d6da8693 Carles Marti
            get_atom_coords(molecule, ctr2_mol)
670 7dd94df7 Carles Marti
671 7dd94df7 Carles Marti
        # Check if rotation was successful
672 7dd94df7 Carles Marti
        # Check dihedral rotation
673 7dd94df7 Carles Marti
        if dihedral_use_mol2:
674 7dd94df7 Carles Marti
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
675 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
676 7dd94df7 Carles Marti
        else:
677 7dd94df7 Carles Marti
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
678 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
679 7dd94df7 Carles Marti
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
680 7dd94df7 Carles Marti
                                        vect_inter)
681 7dd94df7 Carles Marti
        # Check bond rotation is unmodified
682 7dd94df7 Carles Marti
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
683 7dd94df7 Carles Marti
        if np.isclose((dihedral_angle - dihedral_angle_target + 90) % 180 - 90,
684 5261a07f Carles Marti
                      0, atol=1e-3) \
685 5261a07f Carles Marti
                and np.isclose((bond_angle - bond_angle_target + 90) %
686 5261a07f Carles Marti
                               180 - 90, 0, atol=1e-5) \
687 c25aa299 Carles Marti
                and np.allclose(get_atom_coords(molecule, ctr1_mol)
688 c25aa299 Carles Marti
                                - get_atom_coords(surf, ctr1_surf),
689 c25aa299 Carles Marti
                                vect_inter):
690 7dd94df7 Carles Marti
            pass  # print( "Dihedral rotation successfully applied (error: {
691 7dd94df7 Carles Marti
            # :.5f}°)".format((dihedral_angle - dihedral_angle_target + 90) %
692 7dd94df7 Carles Marti
            # 180 - 90))
693 7dd94df7 Carles Marti
        else:
694 7dd94df7 Carles Marti
            err = 'An unknown error occured during the dihedral rotation'
695 7dd94df7 Carles Marti
            logger.error(err)
696 7dd94df7 Carles Marti
            raise AssertionError(err)
697 7dd94df7 Carles Marti
698 7dd94df7 Carles Marti
    #####################################
699 7dd94df7 Carles Marti
    # Adsorbate dihedral angle rotation #
700 7dd94df7 Carles Marti
    #####################################
701 7dd94df7 Carles Marti
702 7dd94df7 Carles Marti
    # Perform adsorbate dihedral rotation if possible
703 7dd94df7 Carles Marti
    if do_mol_dihedral:
704 7dd94df7 Carles Marti
        # Retrieve current adsorbate dihedral angle (by computing the angle
705 7dd94df7 Carles Marti
        # between the orthogonal rejection of vect_inter and vect2_mol onto
706 7dd94df7 Carles Marti
        # vect_mol)
707 7dd94df7 Carles Marti
        vect_mol_inner = np.dot(vect_mol, vect_mol)
708 7dd94df7 Carles Marti
        bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \
709 5261a07f Carles Marti
            vect_mol_inner * vect_mol
710 7dd94df7 Carles Marti
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
711 5261a07f Carles Marti
            vect_mol_inner * vect_mol
712 7dd94df7 Carles Marti
        dihedral_angle_ini = get_vect_angle(bond_inter_reject,
713 7dd94df7 Carles Marti
                                            bond2_mol_reject, vect_mol)
714 7dd94df7 Carles Marti
715 7dd94df7 Carles Marti
        # Apply dihedral rotation along vect_mol
716 7dd94df7 Carles Marti
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
717 7dd94df7 Carles Marti
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
718 7dd94df7 Carles Marti
719 7dd94df7 Carles Marti
        # Update molecular bonds
720 7dd94df7 Carles Marti
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
721 5261a07f Carles Marti
            - get_atom_coords(molecule, ctr1_mol)
722 7dd94df7 Carles Marti
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
723 5261a07f Carles Marti
            - get_atom_coords(molecule, ctr2_mol)
724 7dd94df7 Carles Marti
725 7dd94df7 Carles Marti
        # Check if rotation was successful
726 7dd94df7 Carles Marti
        # Check adsorbate dihedral rotation
727 7dd94df7 Carles Marti
        vect_mol_inner = np.dot(vect_mol, vect_mol)
728 7dd94df7 Carles Marti
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
729 5261a07f Carles Marti
            vect_mol_inner * vect_mol
730 7dd94df7 Carles Marti
        mol_dihedral_angle = get_vect_angle(bond_inter_reject,
731 7dd94df7 Carles Marti
                                            bond2_mol_reject, vect_mol)
732 7dd94df7 Carles Marti
        # Check dihedral rotation
733 7dd94df7 Carles Marti
        vect_inter_inner = np.dot(vect_inter, vect_inter)
734 7dd94df7 Carles Marti
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
735 5261a07f Carles Marti
            vect_inter_inner * vect_inter
736 7dd94df7 Carles Marti
        vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
737 5261a07f Carles Marti
            vect_inter_inner * vect_inter
738 7dd94df7 Carles Marti
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
739 7dd94df7 Carles Marti
                                        vect_inter)
740 7dd94df7 Carles Marti
        # Check bond rotation is unmodified
741 7dd94df7 Carles Marti
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
742 7dd94df7 Carles Marti
        if np.isclose((mol_dihedral_angle - mol_dihedral_angle_target + 90) %
743 7dd94df7 Carles Marti
                      180 - 90, 0, atol=1e-3) \
744 7dd94df7 Carles Marti
                and np.isclose((dihedral_angle -
745 7dd94df7 Carles Marti
                                dihedral_angle_target + 90) % 180 - 90, 0,
746 7dd94df7 Carles Marti
                               atol=1e-5) \
747 7dd94df7 Carles Marti
                and np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90,
748 7dd94df7 Carles Marti
                               0, atol=1e-5) \
749 7dd94df7 Carles Marti
                and np.allclose(get_atom_coords(molecule, ctr1_mol) -
750 7dd94df7 Carles Marti
                                get_atom_coords(surf, ctr1_surf),
751 7dd94df7 Carles Marti
                                vect_inter):
752 7dd94df7 Carles Marti
            pass  # print(
753 7dd94df7 Carles Marti
            # "Adsorbate dihedral rotation successfully applied (error:
754 7dd94df7 Carles Marti
            # {:.5f}°)".format((mol_dihedral_angle - mol_dihedral_angle_target
755 7dd94df7 Carles Marti
            # + 90) % 180 - 90))
756 7dd94df7 Carles Marti
        else:
757 7dd94df7 Carles Marti
            err = 'An unknown error occured during the adsorbate dihedral ' \
758 7dd94df7 Carles Marti
                  'rotation'
759 7dd94df7 Carles Marti
            logger.error(err)
760 7dd94df7 Carles Marti
            raise AssertionError(err)
761 7dd94df7 Carles Marti
762 7dd94df7 Carles Marti
763 7dd94df7 Carles Marti
def ads_chemcat(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
764 7dd94df7 Carles Marti
                ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
765 c25aa299 Carles Marti
                slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel):
766 5261a07f Carles Marti
    """Generates adsorbate-surface structures by sampling over chemcat angles.
767 5261a07f Carles Marti

768 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb (adsorbate).
769 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the molecule
770 5261a07f Carles Marti
    @param ctr1_mol: The index/es of the center in the adsorbate to use as
771 5261a07f Carles Marti
        adsorption center.
772 5261a07f Carles Marti
    @param ctr2_mol: The index/es of the center in the adsorbate to use for the
773 5261a07f Carles Marti
        definition of the surf-adsorbate angle, surf-adsorbate dihedral angle
774 5261a07f Carles Marti
        and adsorbate dihedral angle.
775 5261a07f Carles Marti
    @param ctr3_mol: The index/es of the center in the adsorbate to use for the
776 5261a07f Carles Marti
        definition of the adsorbate dihedral angle.
777 5261a07f Carles Marti
    @param ctr1_surf: The index/es of the center in the surface to use as
778 5261a07f Carles Marti
        adsorption center.
779 5261a07f Carles Marti
    @param ctr2_surf: The index/es of the center in the surface to use for the
780 5261a07f Carles Marti
        definition of the surf-adsorbate dihedral angle.
781 5261a07f Carles Marti
    @param num_pts: Number on which to sample Euler angles.
782 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision
783 5261a07f Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
784 5261a07f Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
785 5261a07f Carles Marti
        considered as atomic collision.
786 5261a07f Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
787 5261a07f Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
788 5261a07f Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
789 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
790 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
791 c25aa299 Carles Marti
    @param max_hel: Maximum value for sampling the helicopter
792 5261a07f Carles Marti
        (surf-adsorbate dihedral) angle.
793 5261a07f Carles Marti
    @return: list of ase.Atoms object conatining all the orientations of a given
794 5261a07f Carles Marti
        conformer.
795 5261a07f Carles Marti
    """
796 7dd94df7 Carles Marti
    from copy import deepcopy
797 7dd94df7 Carles Marti
    slab_ads_list = []
798 70d3cf53 Carles Marti
    # Rotation over bond angle # TODO Check sampling
799 7dd94df7 Carles Marti
    for alpha in np.arange(90, 180, 90 / min(num_pts, 2)):
800 7dd94df7 Carles Marti
        # Rotation over surf-adsorbate dihedral angle
801 c25aa299 Carles Marti
        for beta in np.arange(0, max_hel, max_hel / num_pts):
802 7dd94df7 Carles Marti
            # Rotation over adsorbate bond dihedral angle
803 7dd94df7 Carles Marti
            for gamma in np.arange(0, 360, 360 / num_pts):
804 7dd94df7 Carles Marti
                new_molec = deepcopy(orig_molec)
805 7dd94df7 Carles Marti
                chemcat_rotate(new_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol,
806 7dd94df7 Carles Marti
                               ctr1_surf, ctr2_surf, norm_vect, alpha,
807 7dd94df7 Carles Marti
                               beta, gamma)
808 7dd94df7 Carles Marti
                site_coords = get_atom_coords(slab, ctr1_surf)
809 7dd94df7 Carles Marti
                ctr_coords = get_atom_coords(new_molec, ctr1_mol)
810 7dd94df7 Carles Marti
                slab_molec, collision = correct_coll(new_molec, slab,
811 7dd94df7 Carles Marti
                                                     ctr_coords, site_coords,
812 7dd94df7 Carles Marti
                                                     num_pts, min_coll_height,
813 7dd94df7 Carles Marti
                                                     norm_vect, slab_nghbs,
814 e8bebcca Carles Marti
                                                     molec_nghbs, coll_coeff)
815 7dd94df7 Carles Marti
                if not collision and \
816 5261a07f Carles Marti
                        not any([np.allclose(slab_molec.positions,
817 5261a07f Carles Marti
                                             conf.positions)
818 7dd94df7 Carles Marti
                                 for conf in slab_ads_list]):
819 7dd94df7 Carles Marti
                    slab_ads_list.append(slab_molec)
820 c25aa299 Carles Marti
                    if h_donor is not False:
821 c25aa299 Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
822 c25aa299 Carles Marti
                                                          h_acceptor,
823 c25aa299 Carles Marti
                                                          len(slab)))
824 7dd94df7 Carles Marti
825 7dd94df7 Carles Marti
    return slab_ads_list
826 f3d1e601 Carles Marti
827 f3d1e601 Carles Marti
828 7dd94df7 Carles Marti
def adsorb_confs(conf_list, surf, inp_vars):
829 a5cc42ff Carles Marti
    """Generates a number of adsorbate-surface structure coordinates.
830 a5cc42ff Carles Marti

831 a5cc42ff Carles Marti
    Given a list of conformers, a surface, a list of atom indices (or list of
832 a5cc42ff Carles Marti
    list of indices) of both the surface and the adsorbate, it generates a
833 a5cc42ff Carles Marti
    number of adsorbate-surface structures for every possible combination of
834 a5cc42ff Carles Marti
    them at different orientations.
835 a5cc42ff Carles Marti
    @param conf_list: list of ase.Atoms of the different conformers
836 a5cc42ff Carles Marti
    @param surf: the ase.Atoms object of the surface
837 7dd94df7 Carles Marti
    @param inp_vars: Calculation parameters from input file.
838 a5cc42ff Carles Marti
    @return: list of ase.Atoms for the adsorbate-surface structures
839 a5cc42ff Carles Marti
    """
840 bb55f47c Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
841 7dd94df7 Carles Marti
    molec_ctrs = inp_vars['molec_ctrs']
842 7dd94df7 Carles Marti
    sites = inp_vars['sites']
843 7dd94df7 Carles Marti
    angles = inp_vars['set_angles']
844 7dd94df7 Carles Marti
    num_pts = inp_vars['sample_points_per_angle']
845 7dd94df7 Carles Marti
    norm_vect = inp_vars['surf_norm_vect']
846 7dd94df7 Carles Marti
    min_coll_height = inp_vars['min_coll_height']
847 7dd94df7 Carles Marti
    coll_coeff = inp_vars['collision_threshold']
848 c25aa299 Carles Marti
    h_donor = inp_vars['h_donor']
849 c25aa299 Carles Marti
    h_acceptor = inp_vars['h_acceptor']
850 7dd94df7 Carles Marti
851 bc703cab Carles Marti
    if inp_vars['pbc_cell'] is not False:
852 bc703cab Carles Marti
        surf.set_pbc(True)
853 bc703cab Carles Marti
        surf.set_cell(inp_vars['pbc_cell'])
854 bc703cab Carles Marti
855 f3d1e601 Carles Marti
    surf_ads_list = []
856 f3d1e601 Carles Marti
    sites_coords = get_atom_coords(surf, sites)
857 e8bebcca Carles Marti
    if coll_coeff is not False:
858 5fb01677 Carles Marti
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
859 5fb01677 Carles Marti
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
860 5fb01677 Carles Marti
    else:
861 5fb01677 Carles Marti
        surf_nghbs = 0
862 7dd94df7 Carles Marti
    for i, conf in enumerate(conf_list):
863 bb55f47c Carles Marti
        molec_ctr_coords = get_atom_coords(conf, molec_ctrs)
864 bc703cab Carles Marti
        if inp_vars['pbc_cell'] is not False:
865 bc703cab Carles Marti
            conf.set_pbc(True)
866 bc703cab Carles Marti
            conf.set_cell(inp_vars['pbc_cell'])
867 e8bebcca Carles Marti
        if coll_coeff is not False:
868 5fb01677 Carles Marti
            conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
869 5fb01677 Carles Marti
            molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
870 5fb01677 Carles Marti
        else:
871 5fb01677 Carles Marti
            molec_nghbs = 0
872 7dd94df7 Carles Marti
        for s, site in enumerate(sites_coords):
873 d6da8693 Carles Marti
            if norm_vect == 'auto':
874 d6da8693 Carles Marti
                norm_vect = compute_norm_vect(surf, sites[s],
875 d6da8693 Carles Marti
                                              inp_vars['pbc_cell'])
876 7dd94df7 Carles Marti
            for c, ctr in enumerate(molec_ctr_coords):
877 7dd94df7 Carles Marti
                if angles == 'euler':
878 bb55f47c Carles Marti
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
879 5fb01677 Carles Marti
                                                   num_pts, min_coll_height,
880 bb55f47c Carles Marti
                                                   coll_coeff, norm_vect,
881 b4b2f307 Carles Marti
                                                   surf_nghbs, molec_nghbs,
882 c25aa299 Carles Marti
                                                   h_donor, h_acceptor))
883 7dd94df7 Carles Marti
                elif angles == 'chemcat':
884 7dd94df7 Carles Marti
                    mol_ctr1 = molec_ctrs[c]
885 7dd94df7 Carles Marti
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
886 7dd94df7 Carles Marti
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
887 7dd94df7 Carles Marti
                    surf_ctr1 = sites[s]
888 7dd94df7 Carles Marti
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
889 a98d4172 Carles Marti
                    max_h = inp_vars["max_helic_angle"]
890 7dd94df7 Carles Marti
                    surf_ads_list.extend(ads_chemcat(conf, surf, mol_ctr1,
891 7dd94df7 Carles Marti
                                                     mol_ctr2, mol_ctr3,
892 7dd94df7 Carles Marti
                                                     surf_ctr1, surf_ctr2,
893 7dd94df7 Carles Marti
                                                     num_pts, min_coll_height,
894 7dd94df7 Carles Marti
                                                     coll_coeff, norm_vect,
895 7dd94df7 Carles Marti
                                                     surf_nghbs, molec_nghbs,
896 c25aa299 Carles Marti
                                                     h_donor, h_acceptor,
897 c25aa299 Carles Marti
                                                     max_h))
898 f3d1e601 Carles Marti
    return surf_ads_list
899 f3d1e601 Carles Marti
900 f3d1e601 Carles Marti
901 4614bb6a Carles
def run_screening(inp_vars):
902 91ae8d86 Carles Marti
    """Carries out the screening of adsorbate structures on a surface.
903 e07c09eb Carles

904 e07c09eb Carles
    @param inp_vars: Calculation parameters from input file.
905 e07c09eb Carles
    """
906 e07c09eb Carles
    import os
907 57e3a8c7 Carles Marti
    import random
908 fd2384fc Carles Marti
    from modules.formats import collect_coords, adapt_format
909 cf8fe0e3 Carles Marti
    from modules.calculation import run_calc, check_finished_calcs
910 e07c09eb Carles
911 76f4ac19 Carles Marti
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
912 76f4ac19 Carles Marti
                ' structures.')
913 e07c09eb Carles
    if not os.path.isdir("isolated"):
914 e07c09eb Carles
        err = "'isolated' directory not found. It is needed in order to carry "
915 e07c09eb Carles
        "out the screening of structures to be adsorbed"
916 e07c09eb Carles
        logger.error(err)
917 cf8fe0e3 Carles Marti
        raise FileNotFoundError(err)
918 e07c09eb Carles
919 cf8fe0e3 Carles Marti
    finished_calcs, unfinished_calcs = check_finished_calcs('isolated',
920 cf8fe0e3 Carles Marti
                                                            inp_vars['code'])
921 1a1164e0 Carles Marti
    logger.info(f"Found {len(finished_calcs)} structures of isolated "
922 1a1164e0 Carles Marti
                f"conformers whose calculation finished normally.")
923 1a1164e0 Carles Marti
    if len(unfinished_calcs) != 0:
924 cf8fe0e3 Carles Marti
        logger.warning(f"Found {len(unfinished_calcs)} calculations more that "
925 76f4ac19 Carles Marti
                       f"did not finish normally: {unfinished_calcs}. \n"
926 76f4ac19 Carles Marti
                       f"Using only the ones that finished normally: "
927 76f4ac19 Carles Marti
                       f"{finished_calcs}.")
928 cf8fe0e3 Carles Marti
929 fd2384fc Carles Marti
    conformer_atoms_list = collect_coords(finished_calcs, inp_vars['code'],
930 fd2384fc Carles Marti
                                          'isolated', inp_vars['special_atoms'])
931 1e9e784d Carles Marti
    selected_confs = select_confs(conformer_atoms_list, finished_calcs,
932 fd2384fc Carles Marti
                                  inp_vars['select_magns'],
933 bfe93f0d Carles Marti
                                  inp_vars['confs_per_magn'],
934 bfe93f0d Carles Marti
                                  inp_vars['code'])
935 90819cc3 Carles Marti
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
936 7dd94df7 Carles Marti
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
937 7d97341d Carles Marti
    if len(surf_ads_list) > inp_vars['max_structures']:
938 57e3a8c7 Carles Marti
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
939 bfe93f0d Carles Marti
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
940 d9167fea Carles Marti
                f'configurations to carry out a calculation of.')
941 d9167fea Carles Marti
942 f3d1e601 Carles Marti
    run_calc('screening', inp_vars, surf_ads_list)
943 14f39d2a Carles Marti
    logger.info('Finished the procedures for the screening of adsorbate-surface'
944 14f39d2a Carles Marti
                ' structures section.')