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

dockonsurf / modules / screening.py @ cf980c86

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

1
import logging
2
import numpy as np
3
import ase
4

    
5
logger = logging.getLogger('DockOnSurf')
6

    
7

    
8
def select_confs(conf_list: list, magns: list, num_sel: int):
9
    """Takes a list ase.Atoms and selects the most different magnitude-wise.
10

11
    Given a list of ase.Atoms objects and a list of magnitudes, it selects a
12
    number of the most different conformers according to every magnitude
13
    specified.
14
     
15
    @param conf_list: list of ase.Atoms objects to select among.
16
    @param magns: list of str with the names of the magnitudes to use for the
17
        conformer selection. Supported magnitudes: 'energy', 'moi'.
18
    @param num_sel: number of conformers to select for every of the magnitudes.
19
    @return: list of the selected ase.Atoms objects.
20
    """
21
    selected_ids = []
22
    if num_sel >= len(conf_list):
23
        logger.warning('Number of conformers per magnitude is equal or larger '
24
                       'than the total number of conformers. Using all '
25
                       f'available conformers: {len(conf_list)}.')
26
        return conf_list
27

    
28
    # Assign mois
29
    if 'moi' in magns:
30
        for conf in conf_list:
31
            conf.info["moi"] = conf.get_moments_of_inertia().sum()
32

    
33
    # pick ids
34
    for magn in magns:
35
        sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn]))
36
        if sorted_list[-1].info['iso'] not in selected_ids:
37
            selected_ids.append(sorted_list[-1].info['iso'])
38
        if num_sel > 1:
39
            for i in range(0, len(sorted_list) - 1,
40
                           len(conf_list) // (num_sel - 1)):
41
                if sorted_list[i].info['iso'] not in selected_ids:
42
                    selected_ids.append(sorted_list[i].info['iso'])
43

    
44
    logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
45
    return [conf for conf in conf_list if conf.info["iso"] in selected_ids]
46

    
47

    
48
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True):
49
    """Computes the angle between two vectors.
50

51
    @param v1: The first vector.
52
    @param v2: The second vector.
53
    @param ref: Orthogonal vector to both v1 and v2,
54
        along which the sign of the rotation is defined (i.e. positive if
55
        counterclockwise angle when facing ref)
56
    @param degrees: Whether the result should be in radians (True) or in
57
        degrees (False).
58
    @return: The angle in radians if degrees = False, or in degrees if
59
        degrees =True
60
    """
61
    v1_u = v1 / np.linalg.norm(v1)
62
    v2_u = v2 / np.linalg.norm(v2)
63
    angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0))
64
    if ref is not None:
65
        # Give sign according to ref direction
66
        angle *= (1 if np.dot(np.cross(v1, v2), ref) >= 0 else -1)
67

    
68
    return angle if not degrees else angle * 180 / np.pi
69

    
70

    
71
def vect_avg(vects):
72
    """Computes the element-wise mean of a set of vectors.
73

74
    @param vects: list of lists-like: containing the vectors (num_vectors,
75
        length_vector).
76
    @return: vector average computed doing the element-wise mean.
77
    """
78
    from modules.utilities import try_command
79
    err = "vect_avg parameter vects must be a list-like, able to be converted" \
80
          " np.array"
81
    array = try_command(np.array, [(ValueError, err)], vects)
82
    if len(array.shape) == 1:
83
        return array
84
    else:
85
        num_vects = array.shape[1]
86
        return np.array([np.average(array[:, i]) for i in range(num_vects)])
87

    
88

    
89
def get_atom_coords(atoms: ase.Atoms, center=None):
90
    """Gets the coordinates of the specified center for an ase.Atoms object.
91

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

    
119

    
120
def compute_norm_vect(atoms, idxs, cell):
121
    """Computes the local normal vector of a surface at a given site.
122

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

    
154

    
155
def align_molec(orig_molec, ctr_coord, ref_vect):
156
    """Align a molecule to a vector by a center.
157

158
    Given a reference vector to be aligned to and some coordinates acting as
159
    alignment center, it first averages the vectors pointing to neighboring
160
    atoms and then tries to align this average vector to the target. If the
161
    average vector is 0 it takes the vector to the nearest neighbor.
162
    @param orig_molec: The molecule to align.
163
    @param ctr_coord: The coordinates to use ase alignment center.
164
    @param ref_vect: The vector to be aligned with.
165
    @return: ase.Atoms of the aligned molecule.
166
    """
167
    from copy import deepcopy
168
    from ase import Atom
169
    from ase.neighborlist import natural_cutoffs, neighbor_list
170

    
171
    molec = deepcopy(orig_molec)
172
    if len(molec) == 1:
173
        err_msg = "Cannot align a single atom"
174
        logger.error(err_msg)
175
        ValueError(err_msg)
176
    cutoffs = natural_cutoffs(molec, mult=1.2)
177
    # Check if ctr_coord are the coordinates of an atom and if not creates a
178
    # dummy one to extract the neighboring atoms.
179
    ctr_idx = None
180
    dummy_atom = False
181
    for atom in molec:
182
        if np.allclose(ctr_coord, atom.position, rtol=1e-2):
183
            ctr_idx = atom.index
184
            break
185
    if ctr_idx is None:
186
        molec.append(Atom('X', position=ctr_coord))
187
        cutoffs.append(0.2)
188
        ctr_idx = len(molec) - 1
189
        dummy_atom = True
190
    # Builds the neighbors and computes the average vector
191
    refs, vects = neighbor_list("iD", molec, cutoffs, self_interaction=False)
192
    neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
193
    # If no neighbors are present, the cutoff of the alignment center is
194
    # set to a value where at least one atom is a neighbor and neighbors are
195
    # recalculated.
196
    if len(neigh_vects) == 0:
197
        min_dist, min_idx = (np.inf, np.inf)
198
        for atom in molec:
199
            if atom.index == ctr_idx:
200
                continue
201
            if molec.get_distance(ctr_idx, atom.index) < min_dist:
202
                min_dist = molec.get_distance(ctr_idx, atom.index)
203
                min_idx = atom.index
204
        cutoffs[ctr_idx] = min_dist - cutoffs[min_idx] + 0.05
205
        refs, vects = neighbor_list("iD", molec, cutoffs,
206
                                    self_interaction=False)
207
        neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
208
    target_vect = vect_avg(neigh_vects)
209
    # If the target vector is 0 (the center is at the baricenter of its
210
    # neighbors). Assuming the adsorption center is coplanar or colinear to its
211
    # neighbors (it would not make a lot of sense to chose a center which is
212
    # the baricenter of neighbors distributed in 3D), the target_vector is
213
    # chosen perpendicular to the nearest neighbor.
214
    if np.allclose(target_vect, 0, rtol=1e-3):
215
        nn_vect = np.array([np.inf] * 3)
216
        for vect in neigh_vects:
217
            if np.linalg.norm(vect) < np.linalg.norm(nn_vect):
218
                nn_vect = vect
219
        cart_axes = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]
220
        axis = cart_axes[int(np.argmax([np.linalg.norm(np.cross(ax, nn_vect))
221
                                        for ax in cart_axes]))]
222
        target_vect = np.cross(axis, nn_vect)
223

    
224
    rot_vect = np.cross(target_vect, ref_vect)
225
    if np.allclose(rot_vect, 0):
226
        cart_axes = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]
227
        axis = cart_axes[int(np.argmax([np.linalg.norm(np.cross(ax, ref_vect))
228
                                        for ax in cart_axes]))]
229
        rot_vect = np.cross(ref_vect, axis)
230
    rot_angle = -get_vect_angle(ref_vect, target_vect, rot_vect)
231
    molec.rotate(rot_angle, rot_vect, ctr_coord)
232
    if dummy_atom:
233
        del molec[-1]
234
    return molec
235

    
236

    
237
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
238
                  norm_vect=(0, 0, 1)):
239
    """Add an adsorbate to a surface.
240

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

    
291

    
292
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
293
                    nn_molec=0, coll_coeff=1.2, exclude_atom=False):
294
    """Checks whether a slab and a molecule collide or not.
295

296
    @param slab_molec: The system of adsorbate-slab for which to detect if there
297
        are collisions.
298
    @param nn_slab: Number of neigbors in the surface.
299
    @param nn_molec: Number of neighbors in the molecule.
300
    @param coll_coeff: The coefficient that multiplies the covalent radius of
301
        atoms resulting in a distance that two atoms being closer to that is
302
        considered as atomic collision.
303
    @param slab_num_atoms: Number of atoms of the bare slab.
304
    @param min_height: The minimum height atoms can have to not be considered as
305
        colliding.
306
    @param vect: The vector perpendicular to the slab.
307
    @param exclude_atom: Whether to exclude the adsorption center in the
308
        molecule in the collision detection.
309
    @return: bool, whether the surface and the molecule collide.
310
    """
311
    from copy import deepcopy
312
    from ase.neighborlist import natural_cutoffs, neighbor_list
313

    
314
    normal = 0
315
    for i, coord in enumerate(vect):
316
        if coord == 0:
317
            continue
318
        normal = i
319
    if vect[normal] > 0:
320
        surf_height = max(slab_molec[:slab_num_atoms].positions[:, normal])
321
    else:
322
        surf_height = min(slab_molec[:slab_num_atoms].positions[:, normal])
323

    
324
    # Check structure overlap by height
325
    if min_height is not False:
326
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
327
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
328
        if vect.tolist() not in cart_axes:
329
            err_msg = "'min_coll_height' option is only implemented for " \
330
                      "'surf_norm_vect' to be one of the x, y or z axes. "
331
            logger.error(err_msg)
332
            raise ValueError(err_msg)
333
        for atom in slab_molec[slab_num_atoms:]:
334
            if exclude_atom is not False \
335
                    and atom.index == exclude_atom:
336
                continue
337
            for i, coord in enumerate(vect):
338
                if coord == 0:
339
                    continue
340
                if (atom.position[i] - surf_height) * coord < min_height:
341
                    return True
342

    
343
    # Check structure overlap by sphere collision
344
    if coll_coeff is not False:
345
        if exclude_atom is not False:
346
            slab_molec_wo_ctr = deepcopy(slab_molec)
347
            del slab_molec_wo_ctr[exclude_atom + slab_num_atoms]
348
            slab_molec_cutoffs = natural_cutoffs(slab_molec_wo_ctr,
349
                                                 mult=coll_coeff)
350
            slab_molec_nghbs = len(neighbor_list("i", slab_molec_wo_ctr,
351
                                                 slab_molec_cutoffs))
352
        else:
353
            slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
354
            slab_molec_nghbs = len(neighbor_list("i", slab_molec,
355
                                                 slab_molec_cutoffs))
356
        if slab_molec_nghbs > nn_slab + nn_molec:
357
            return True
358

    
359
    return False
360

    
361

    
362
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
363
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
364
                 coll_coeff, height=2.5, excl_atom=False):
365
    # TODO Rename this function
366
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
367
    small rotations.
368

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

    
410

    
411
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
412
                 neigh_cutoff=1):
413
    # TODO rethink
414
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
415

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

    
453
    if np.linalg.norm(surf_h_vect) != np.infty:
454
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
455
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
456
        return slab_molec
457

    
458

    
459
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
460
    # TODO multiple dissociation
461
    """Decides which H atoms to dissociate according to a list of atoms.
462

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

    
508

    
509
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
510
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
511
              h_donor, h_acceptor, height, excl_atom):
512
    """Generates adsorbate-surface structures by sampling over Euler angles.
513

514
    This function generates a number of adsorbate-surface structures at
515
    different orientations of the adsorbate sampled at multiple Euler (zxz)
516
    angles.
517
    @param orig_molec: ase.Atoms object of the molecule to adsorb.
518
    @param slab: ase.Atoms object of the surface on which to adsorb the
519
        molecule.
520
    @param ctr_coord: The coordinates of the molecule to use as adsorption
521
        center.
522
    @param site_coord: The coordinates of the surface on which to adsorb the
523
        molecule.
524
    @param num_pts: Number on which to sample Euler angles.
525
    @param min_coll_height: The lowest height for which to detect a collision.
526
    @param coll_coeff: The coefficient that multiplies the covalent radius of
527
        atoms resulting in a distance that two atoms being closer to that is
528
        considered as atomic collision.
529
    @param norm_vect: The vector perpendicular to the surface.
530
    @param slab_nghbs: Number of neigbors in the surface.
531
    @param molec_nghbs: Number of neighbors in the molecule.
532
    @param h_donor: List of atom types or atom numbers that are H-donors.
533
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
534
    @param height: Height on which to try adsorption.
535
    @param excl_atom: Whether to exclude the adsorption center in the
536
        molecule in the collision detection.
537
    @return: list of ase.Atoms object conatining all the orientations of a given
538
        conformer.
539
    """
540
    from copy import deepcopy
541
    slab_ads_list = []
542
    prealigned_molec = align_molec(orig_molec, ctr_coord, norm_vect)
543
    # rotation around z
544
    for alpha in np.arange(0, 360, 360 / num_pts):
545
        # rotation around x'
546
        for beta in np.arange(0, 180, 180 / num_pts):
547
            # rotation around z'
548
            for gamma in np.arange(0, 360, 360 / num_pts):
549
                if beta == 0 and gamma > 0:
550
                    continue
551
                molec = deepcopy(prealigned_molec)
552
                molec.euler_rotate(alpha, beta, gamma, center=ctr_coord)
553
                slab_molec, collision = correct_coll(molec, slab, ctr_coord,
554
                                                     site_coord, num_pts,
555
                                                     min_coll_height, norm_vect,
556
                                                     slab_nghbs, molec_nghbs,
557
                                                     coll_coeff, height,
558
                                                     excl_atom)
559
                if not collision and not any([np.allclose(slab_molec.positions,
560
                                                          conf.positions)
561
                                              for conf in slab_ads_list]):
562
                    slab_molec.info = {**slab_molec.info, **molec.info}
563
                    slab_ads_list.append(slab_molec)
564
                    if h_donor is not False:
565
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
566
                                                          h_acceptor,
567
                                                          len(slab)))
568

    
569
    return slab_ads_list
570

    
571

    
572
def internal_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
573
                    ctr2_surf, bond_vector, bond_angle_target,
574
                    dihedral_angle_target=None, mol_dihedral_angle_target=None):
575
    """Performs translation and rotation of an adsorbate defined by an
576
    adsorption bond length, direction, angle and dihedral angle
577

578
    Carles modification of chemcat's transform_adsorbate to work with
579
    coordinates instead of ase.Atom
580
    Parameters:
581
        molecule (ase.Atoms): The molecule to adsorb.
582

583
        surf (ase.Atoms): The surface ontop of which to adsorb.
584

585
        ctr1_mol (int/list): The position of the adsorption center in the
586
        molecule that will be bound to the surface.
587

588
        ctr2_mol (int/list): The position of a second center of the
589
        adsorbate used to define the adsorption bond angle, and the dihedral
590
        adsorption angle.
591

592
        ctr3_mol (int/list): The position of a third center in the
593
        adsorbate used to define the adsorbate dihedral angle.
594

595
        ctr1_surf (int/list): The position of the site on the surface that
596
        will be bound to the molecule.
597

598
        ctr2_surf (int/list): The position of a second center of the
599
        surface used to define the dihedral adsorption angle.
600

601
        bond_vector (numpy.ndarray): The adsorption bond desired.
602
            Details: offset = vect(atom1_surf, atom1_mol)
603

604
        bond_angle_target (float or int): The adsorption bond angle desired (in
605
            degrees).
606
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
607

608
        dihedral_angle_target (float or int): The dihedral adsorption angle
609
            desired (in degrees).
610
            Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
611
            atom1_mol, atom2_mol)
612
                The rotation vector is facing the adsorbate from the surface
613
                (i.e. counterclockwise rotation when facing the surface (i.e.
614
                view from top))
615

616
        mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
617
            desired (in degrees).
618
            Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
619
            atom2_mol, atom3_mol)
620
                The rotation vector is facing atom2_mol from atom1_mol
621

622
    Returns:
623
        None (the `molecule` object is modified in-place)
624
    """
625
    vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf,
626
                                                                   ctr1_surf)
627
    vect_inter = get_atom_coords(molecule, ctr1_mol) \
628
        - get_atom_coords(surf, ctr1_surf)
629
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
630
                                                                     ctr1_mol)
631
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
632
                                                                      ctr2_mol)
633

    
634
    # Check if dihedral angles can be defined
635
    do_dihedral = dihedral_angle_target is not None
636
    do_mol_dihedral = mol_dihedral_angle_target is not None
637
    dihedral_use_mol2 = False
638
    # Check if vect_surf and bond_vector are aligned
639
    if np.allclose(np.cross(vect_surf, bond_vector), 0):
640
        do_dihedral = False
641
    # Check if requested bond angle is not flat
642
    if np.isclose((bond_angle_target + 90) % 180 - 90, 0):
643
        do_mol_dihedral = False
644
        dihedral_use_mol2 = True
645
    # Check if vect_mol and vect2_mol are not aligned
646
    if np.allclose(np.cross(vect_mol, vect2_mol), 0):
647
        do_mol_dihedral = False
648

    
649
    ###########################
650
    #       Translation       #
651
    ###########################
652

    
653
    # Compute and apply translation of adsorbate
654
    translation = bond_vector - vect_inter
655
    molecule.translate(translation)
656

    
657
    # Update adsorption bond
658
    vect_inter = get_atom_coords(molecule, ctr1_mol) - \
659
        get_atom_coords(surf, ctr1_surf)
660

    
661
    # Check if translation was successful
662
    if np.allclose(vect_inter, bond_vector):
663
        pass  # print("Translation successfully applied (error: ~ {:.5g} unit "
664
        # "length)".format(np.linalg.norm(vect_inter - bond_vector)))
665
    else:
666
        err = 'An unknown error occured during the translation'
667
        logger.error(err)
668
        raise AssertionError(err)
669

    
670
    ###########################
671
    #   Bond angle rotation   #
672
    ###########################
673

    
674
    # Compute rotation vector
675
    rotation_vector = np.cross(-vect_inter, vect_mol)
676
    if np.allclose(rotation_vector, 0, atol=1e-5):
677
        # If molecular bonds are aligned, any vector orthogonal to vect_inter
678
        # can be used Such vector can be found as the orthogonal rejection of
679
        # either X-axis, Y-axis or Z-axis onto vect_inter (since they cannot
680
        # be all aligned)
681
        non_aligned_vector = np.zeros(3)
682
        # Select the most orthogonal axis (lowest dot product):
683
        non_aligned_vector[np.argmin(np.abs(vect_inter))] = 1
684
        rotation_vector = non_aligned_vector - np.dot(non_aligned_vector,
685
                                                      vect_inter) / np.dot(
686
            vect_inter, vect_inter) * vect_inter
687

    
688
    # Retrieve current bond angle
689
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
690

    
691
    # Apply rotation to reach desired bond_angle
692
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
693
                    center=get_atom_coords(molecule, ctr1_mol))
694

    
695
    # Update molecular bonds
696
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
697
                                                                     ctr1_mol)
698
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
699
                                                                      ctr2_mol)
700

    
701
    # Check if rotation was successful
702
    bond_angle = get_vect_angle(-vect_inter, vect_mol)
703
    if np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90, 0,
704
                  atol=1e-3) and np.allclose(get_atom_coords(molecule, ctr1_mol)
705
                                             - get_atom_coords(surf,
706
                                                               ctr1_surf),
707
                                             vect_inter):
708
        pass  # print("Rotation successfully applied (error: {:.5f}°)".format(
709
        # (bond_angle - bond_angle_target + 90) % 180 - 90))
710
    else:
711
        err = 'An unknown error occured during the rotation'
712
        logger.error(err)
713
        raise AssertionError(err)
714

    
715
    ###########################
716
    # Dihedral angle rotation #
717
    ###########################
718

    
719
    # Perform dihedral rotation if possible
720
    if do_dihedral:
721
        # Retrieve current dihedral angle (by computing the angle between the
722
        # vector rejection of vect_surf and vect_mol onto vect_inter)
723
        vect_inter_inner = np.dot(vect_inter, vect_inter)
724
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
725
            vect_inter_inner * vect_inter
726
        if dihedral_use_mol2:
727
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
728
                              vect_inter_inner * vect_inter
729
        else:
730
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
731
                              vect_inter_inner * vect_inter
732
        dihedral_angle_ini = get_vect_angle(vect_surf_reject, vect_mol_reject,
733
                                            vect_inter)
734

    
735
        # Apply dihedral rotation along vect_inter
736
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
737
                        v=vect_inter, center=get_atom_coords(molecule,
738
                                                             ctr1_mol))
739

    
740
        # Update molecular bonds
741
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
742
            get_atom_coords(molecule, ctr1_mol)
743
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
744
            get_atom_coords(molecule, ctr2_mol)
745

    
746
        # Check if rotation was successful
747
        # Check dihedral rotation
748
        if dihedral_use_mol2:
749
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
750
                              vect_inter_inner * vect_inter
751
        else:
752
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
753
                              vect_inter_inner * vect_inter
754
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
755
                                        vect_inter)
756
        # Check bond rotation is unmodified
757
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
758
        if np.isclose((dihedral_angle - dihedral_angle_target + 90) % 180 - 90,
759
                      0, atol=1e-3) \
760
                and np.isclose((bond_angle - bond_angle_target + 90) %
761
                               180 - 90, 0, atol=1e-5) \
762
                and np.allclose(get_atom_coords(molecule, ctr1_mol)
763
                                - get_atom_coords(surf, ctr1_surf),
764
                                vect_inter):
765
            pass  # print( "Dihedral rotation successfully applied (error: {
766
            # :.5f}°)".format((dihedral_angle - dihedral_angle_target + 90) %
767
            # 180 - 90))
768
        else:
769
            err = 'An unknown error occured during the dihedral rotation'
770
            logger.error(err)
771
            raise AssertionError(err)
772

    
773
    #####################################
774
    # Adsorbate dihedral angle rotation #
775
    #####################################
776

    
777
    # Perform adsorbate dihedral rotation if possible
778
    if do_mol_dihedral:
779
        # Retrieve current adsorbate dihedral angle (by computing the angle
780
        # between the orthogonal rejection of vect_inter and vect2_mol onto
781
        # vect_mol)
782
        vect_mol_inner = np.dot(vect_mol, vect_mol)
783
        bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \
784
            vect_mol_inner * vect_mol
785
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
786
            vect_mol_inner * vect_mol
787
        dihedral_angle_ini = get_vect_angle(bond_inter_reject,
788
                                            bond2_mol_reject, vect_mol)
789

    
790
        # Apply dihedral rotation along vect_mol
791
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
792
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
793

    
794
        # Update molecular bonds
795
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
796
            - get_atom_coords(molecule, ctr1_mol)
797
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
798
            - get_atom_coords(molecule, ctr2_mol)
799

    
800
        # Check if rotation was successful
801
        # Check adsorbate dihedral rotation
802
        vect_mol_inner = np.dot(vect_mol, vect_mol)
803
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
804
            vect_mol_inner * vect_mol
805
        mol_dihedral_angle = get_vect_angle(bond_inter_reject,
806
                                            bond2_mol_reject, vect_mol)
807
        # Check dihedral rotation
808
        vect_inter_inner = np.dot(vect_inter, vect_inter)
809
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
810
            vect_inter_inner * vect_inter
811
        vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
812
            vect_inter_inner * vect_inter
813
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
814
                                        vect_inter)
815
        # Check bond rotation is unmodified
816
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
817
        if np.isclose((mol_dihedral_angle - mol_dihedral_angle_target + 90) %
818
                      180 - 90, 0, atol=1e-3) \
819
                and np.isclose((dihedral_angle -
820
                                dihedral_angle_target + 90) % 180 - 90, 0,
821
                               atol=1e-5) \
822
                and np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90,
823
                               0, atol=1e-5) \
824
                and np.allclose(get_atom_coords(molecule, ctr1_mol) -
825
                                get_atom_coords(surf, ctr1_surf),
826
                                vect_inter):
827
            pass  # print(
828
            # "Adsorbate dihedral rotation successfully applied (error:
829
            # {:.5f}°)".format((mol_dihedral_angle - mol_dihedral_angle_target
830
            # + 90) % 180 - 90))
831
        else:
832
            err = 'An unknown error occured during the adsorbate dihedral ' \
833
                  'rotation'
834
            logger.error(err)
835
            raise AssertionError(err)
836

    
837

    
838
def ads_internal(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
839
                 ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
840
                 slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel, height,
841
                 excl_atom):
842
    """Generates adsorbate-surface structures by sampling over internal angles.
843

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

    
910
    return slab_ads_list
911

    
912

    
913
def adsorb_confs(conf_list, surf, inp_vars):
914
    """Generates a number of adsorbate-surface structure coordinates.
915

916
    Given a list of conformers, a surface, a list of atom indices (or list of
917
    list of indices) of both the surface and the adsorbate, it generates a
918
    number of adsorbate-surface structures for every possible combination of
919
    them at different orientations.
920
    @param conf_list: list of ase.Atoms of the different conformers
921
    @param surf: the ase.Atoms object of the surface
922
    @param inp_vars: Calculation parameters from input file.
923
    @return: list of ase.Atoms for the adsorbate-surface structures
924
    """
925
    from copy import deepcopy
926
    from ase.neighborlist import natural_cutoffs, neighbor_list
927
    molec_ctrs = inp_vars['molec_ctrs']
928
    sites = inp_vars['sites']
929
    angles = inp_vars['set_angles']
930
    num_pts = inp_vars['sample_points_per_angle']
931
    inp_norm_vect = inp_vars['surf_norm_vect']
932
    min_coll_height = inp_vars['min_coll_height']
933
    coll_coeff = inp_vars['collision_threshold']
934
    exclude_ads_ctr = inp_vars['exclude_ads_ctr']
935
    h_donor = inp_vars['h_donor']
936
    h_acceptor = inp_vars['h_acceptor']
937
    height = inp_vars['adsorption_height']
938

    
939
    if inp_vars['pbc_cell'] is not False:
940
        surf.set_pbc(True)
941
        surf.set_cell(inp_vars['pbc_cell'])
942

    
943
    surf_ads_list = []
944
    sites_coords = np.array([get_atom_coords(surf, site) for site in sites])
945
    if coll_coeff is not False:
946
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
947
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
948
    else:
949
        surf_nghbs = 0
950
    for i, conf in enumerate(conf_list):
951
        molec_ctr_coords = np.array([get_atom_coords(conf, ctr)
952
                                     for ctr in molec_ctrs])
953
        if inp_vars['pbc_cell'] is not False:
954
            conf.set_pbc(True)
955
            conf.set_cell(inp_vars['pbc_cell'])
956
        for s, site in enumerate(sites_coords):
957
            if isinstance(inp_norm_vect, str) and inp_norm_vect == 'auto':
958
                norm_vect = compute_norm_vect(surf, sites[s],
959
                                              inp_vars['pbc_cell'])
960
                if np.isnan(norm_vect).any():
961
                    logger.warning(f"Could not compute the normal vector to "
962
                                   f"site '{sites[s]}'. Skipping site.")
963
                    continue
964
            else:
965
                norm_vect = inp_norm_vect
966
            for c, ctr in enumerate(molec_ctr_coords):
967
                if exclude_ads_ctr and isinstance(molec_ctrs[c], int):
968
                    exclude_atom = molec_ctrs[c]
969
                else:
970
                    exclude_atom = False
971
                    if exclude_ads_ctr and not isinstance(molec_ctrs[c], int):
972
                        logger.warning("'exclude_ads_ctr' only works for atomic"
973
                                       "centers and not for many-atoms "
974
                                       f"barycenters. {molec_ctrs[c]} are not "
975
                                       f"going to be excluded from collison.")
976
                if coll_coeff and exclude_atom is not False:
977
                    conf_wo_ctr = deepcopy(conf)
978
                    del conf_wo_ctr[exclude_atom]
979
                    conf_cutoffs = natural_cutoffs(conf_wo_ctr, mult=coll_coeff)
980
                    molec_nghbs = len(neighbor_list("i", conf_wo_ctr,
981
                                                    conf_cutoffs))
982
                elif coll_coeff and exclude_atom is False:
983
                    conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
984
                    molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
985
                else:
986
                    molec_nghbs = 0
987
                if angles == 'euler':
988
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
989
                                                   num_pts, min_coll_height,
990
                                                   coll_coeff, norm_vect,
991
                                                   surf_nghbs, molec_nghbs,
992
                                                   h_donor, h_acceptor, height,
993
                                                   exclude_atom))
994
                elif angles == 'internal':
995
                    mol_ctr1 = molec_ctrs[c]
996
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
997
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
998
                    surf_ctr1 = sites[s]
999
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
1000
                    max_h = inp_vars["max_helic_angle"]
1001
                    surf_ads_list.extend(ads_internal(conf, surf, mol_ctr1,
1002
                                                      mol_ctr2, mol_ctr3,
1003
                                                      surf_ctr1, surf_ctr2,
1004
                                                      num_pts, min_coll_height,
1005
                                                      coll_coeff, norm_vect,
1006
                                                      surf_nghbs, molec_nghbs,
1007
                                                      h_donor, h_acceptor,
1008
                                                      max_h, height,
1009
                                                      exclude_atom))
1010
    return surf_ads_list
1011

    
1012

    
1013
def run_screening(inp_vars):
1014
    """Carries out the screening of adsorbate structures on a surface.
1015

1016
    @param inp_vars: Calculation parameters from input file.
1017
    """
1018
    import os
1019
    import random
1020
    from modules.formats import collect_confs, adapt_format
1021
    from modules.calculation import run_calc, check_finished_calcs
1022

    
1023
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
1024
                ' structures.')
1025
    if inp_vars['use_molec_file']:
1026
        selected_confs = [adapt_format('ase', inp_vars['use_molec_file'],
1027
                                       inp_vars['special_atoms'])]
1028
        logger.info(f"Using '{inp_vars['use_molec_file']}' as only conformer.")
1029
    else:
1030
        if not os.path.isdir("isolated"):
1031
            err = "'isolated' directory not found. It is needed in order to " \
1032
                  "carry out the screening of structures to be adsorbed"
1033
            logger.error(err)
1034
            raise FileNotFoundError(err)
1035

    
1036
        finished_calcs, failed_calcs = check_finished_calcs('isolated',
1037
                                                            inp_vars['code'])
1038
        if not finished_calcs:
1039
            err_msg = "No calculations on 'isolated' finished normally."
1040
            logger.error(err_msg)
1041
            raise FileNotFoundError(err_msg)
1042

    
1043
        logger.info(f"Found {len(finished_calcs)} structures of isolated "
1044
                    f"conformers whose calculation finished normally.")
1045
        if len(failed_calcs) != 0:
1046
            logger.warning(
1047
                f"Found {len(failed_calcs)} calculations more that "
1048
                f"did not finish normally: {failed_calcs}. \n"
1049
                f"Using only the ones that finished normally: "
1050
                f"{finished_calcs}.")
1051

    
1052
        conf_list = collect_confs(finished_calcs, inp_vars['code'], 'isolated',
1053
                                  inp_vars['special_atoms'])
1054
        selected_confs = select_confs(conf_list, inp_vars['select_magns'],
1055
                                      inp_vars['confs_per_magn'])
1056
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
1057
    surf.info = {}
1058
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
1059
    if len(surf_ads_list) > inp_vars['max_structures']:
1060
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
1061
    elif len(surf_ads_list) == 0:
1062
        err_msg = "No configurations were generated: Check the parameters in" \
1063
                  "dockonsurf.inp"
1064
        logger.error(err_msg)
1065
        raise ValueError(err_msg)
1066
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
1067
                f'configurations to carry out a calculation of.')
1068

    
1069
    run_calc('screening', inp_vars, surf_ads_list)
1070
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1071
                ' structures section.')