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

dockonsurf / modules / screening.py @ 07edc24f

Historique | Voir | Annoter | Télécharger (47,49 ko)

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

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

    
7

    
8
def assign_prop(atoms: ase.Atoms, prop_name: str, prop_val):  # TODO Needed?
9
    atoms.info[prop_name] = prop_val
10

    
11

    
12
def select_confs(orig_conf_list: list, calc_dirs: list, magns: list,
13
                 num_sel: int, code: str):
14
    """Takes a list ase.Atoms and selects the most different magnitude-wise.
15

16
    Given a list of ase.Atoms objects and a list of magnitudes, it selects a
17
    number of the most different conformers according to every magnitude
18
    specified.
19
     
20
    @param orig_conf_list: list of ase.Atoms objects to select among.
21
    @param calc_dirs: List of directories where to read the energies from.
22
    @param magns: list of str with the names of the magnitudes to use for the
23
        conformer selection. Supported magnitudes: 'energy', 'moi'.
24
    @param num_sel: number of conformers to select for every of the magnitudes.
25
    @param code: The code that generated the magnitude information.
26
         Supported codes: See formats.py
27
    @return: list of the selected ase.Atoms objects.
28
    """
29
    from copy import deepcopy
30
    from modules.formats import collect_energies
31

    
32
    conf_list = deepcopy(orig_conf_list)
33
    conf_enrgs, mois, selected_ids = [], [], []
34
    if num_sel >= len(conf_list):
35
        logger.warning('Number of conformers per magnitude is equal or larger '
36
                       'than the total number of conformers. Using all '
37
                       f'available conformers: {len(conf_list)}.')
38
        return conf_list
39

    
40
    # Read properties
41
    if 'energy' in magns:
42
        conf_enrgs = collect_energies(calc_dirs, code, 'isolated')
43
    if 'moi' in magns:
44
        mois = np.array([conf.get_moments_of_inertia() for conf in conf_list])
45

    
46
    # Assign values
47
    for i, conf in enumerate(conf_list):
48
        assign_prop(conf, 'idx', i)
49
        if 'energy' in magns:
50
            assign_prop(conf, 'energy', conf_enrgs[i])
51
        if 'moi' in magns:
52
            assign_prop(conf, 'moi', mois[i, 2])
53

    
54
    # pick ids
55
    for magn in magns:
56
        sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn]))
57
        if sorted_list[-1].info['idx'] not in selected_ids:
58
            selected_ids.append(sorted_list[-1].info['idx'])
59
        if num_sel > 1:
60
            for i in range(0, len(sorted_list) - 1,
61
                           len(conf_list) // (num_sel - 1)):
62
                if sorted_list[i].info['idx'] not in selected_ids:
63
                    selected_ids.append(sorted_list[i].info['idx'])
64

    
65
    logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
66
    return [conf_list[idx] for idx in selected_ids]
67

    
68

    
69
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True):
70
    """Computes the angle between two vectors.
71

72
    @param v1: The first vector.
73
    @param v2: The second vector.
74
    @param ref: Orthogonal vector to both v1 and v2,
75
        along which the sign of the rotation is defined (i.e. positive if
76
        counterclockwise angle when facing ref)
77
    @param degrees: Whether the result should be in radians (True) or in
78
        degrees (False).
79
    @return: The angle in radians if degrees = False, or in degrees if
80
        degrees =True
81
    """
82
    v1_u = v1 / np.linalg.norm(v1)
83
    v2_u = v2 / np.linalg.norm(v2)
84
    angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0))
85
    if ref is not None:
86
        # Give sign according to ref direction
87
        angle *= (1 if np.dot(np.cross(v1, v2), ref) >= 0 else -1)
88

    
89
    return angle if not degrees else angle * 180 / np.pi
90

    
91

    
92
def vect_avg(vects):
93
    """Computes the element-wise mean of a set of vectors.
94

95
    @param vects: list of lists-like: containing the vectors (num_vectors,
96
        length_vector).
97
    @return: vector average computed doing the element-wise mean.
98
    """
99
    from modules.utilities import try_command
100
    err = "vect_avg parameter vects must be a list-like, able to be converted" \
101
          " np.array"
102
    array = try_command(np.array, [(ValueError, err)], vects)
103
    if len(array.shape) == 1:
104
        return array
105
    else:
106
        num_vects = array.shape[1]
107
        return np.array([np.average(array[:, i]) for i in range(num_vects)])
108

    
109

    
110
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None):
111
    """Gets the coordinates of the specified indices from a ase.Atoms object.
112

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

    
144

    
145
def compute_norm_vect(atoms, idxs, cell):
146
    """Computes the local normal vector of a surface at a given site.
147

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

    
178

    
179
def align_molec(orig_molec, ctr_coord, ref_vect):
180
    """Align a molecule to a vector by a center.
181

182
    Given a reference vector to be aligned to and some coordinates acting as
183
    alignment center, it first averages the vectors pointing to neighboring
184
    atoms and then tries to align this average vector to the target. If the
185
    average vector is 0 it takes the vector to the nearest neighbor.
186
    @param orig_molec: The molecule to align.
187
    @param ctr_coord: The coordinates to use ase alignment center.
188
    @param ref_vect: The vector to be aligned with.
189
    @return: ase.Atoms of the aligned molecule.
190
    """
191
    from copy import deepcopy
192
    from ase import Atom
193
    from ase.neighborlist import natural_cutoffs, neighbor_list
194

    
195
    molec = deepcopy(orig_molec)
196
    if len(molec) == 1:
197
        err_msg = "Cannot align a single atom"
198
        logger.error(err_msg)
199
        ValueError(err_msg)
200
    cutoffs = natural_cutoffs(molec, mult=1.2)
201
    # Check if ctr_coord are the coordinates of an atom and if not creates a
202
    # dummy one to extract the neighboring atoms.
203
    ctr_idx = None
204
    dummy_atom = False
205
    for atom in molec:
206
        if np.allclose(ctr_coord, atom.position, rtol=1e-2):
207
            ctr_idx = atom.index
208
            break
209
    if ctr_idx is None:
210
        molec.append(Atom('X', position=ctr_coord))
211
        cutoffs.append(0.2)
212
        ctr_idx = len(molec) - 1
213
        dummy_atom = True
214
    # Builds the neighbors and computes the average vector
215
    refs, vects = neighbor_list("iD", molec, cutoffs, self_interaction=False)
216
    neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
217
    # If no neighbors are present, the cutoff of the alignment center is
218
    # set to a value where at least one atom is a neighbor and neighbors are
219
    # recalculated.
220
    if len(neigh_vects) == 0:
221
        min_dist, min_idx = (np.inf, np.inf)
222
        for atom in molec:
223
            if atom.index == ctr_idx:
224
                continue
225
            if molec.get_distance(ctr_idx, atom.index) < min_dist:
226
                min_dist = molec.get_distance(ctr_idx, atom.index)
227
                min_idx = atom.index
228
        cutoffs[ctr_idx] = min_dist - cutoffs[min_idx] + 0.05
229
        refs, vects = neighbor_list("iD", molec, cutoffs,
230
                                    self_interaction=False)
231
        neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
232
    target_vect = vect_avg(neigh_vects)
233
    # If the target vector is 0 (the center is at the baricenter of its
234
    # neighbors). Assuming the adsorption center is coplanar or colinear to its
235
    # neighbors (it would not make a lot of sense to chose a center which is
236
    # the baricenter of neighbors distributed in 3D), the target_vector is
237
    # chosen perpendicular to the nearest neighbor.
238
    if np.allclose(target_vect, 0, rtol=1e-3):
239
        nn_vect = np.array([np.inf] * 3)
240
        for vect in neigh_vects:
241
            if np.linalg.norm(vect) < np.linalg.norm(nn_vect):
242
                nn_vect = vect
243
        cart_axes = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]
244
        axis = cart_axes[int(np.argmax([np.linalg.norm(np.cross(ax, nn_vect))
245
                                        for ax in cart_axes]))]
246
        target_vect = np.cross(axis, nn_vect)
247

    
248
    rot_vect = np.cross(target_vect, ref_vect)
249
    if np.allclose(rot_vect, 0):
250
        cart_axes = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]
251
        axis = cart_axes[int(np.argmax([np.linalg.norm(np.cross(ax, ref_vect))
252
                                        for ax in cart_axes]))]
253
        rot_vect = np.cross(ref_vect, axis)
254
    rot_angle = -get_vect_angle(ref_vect, target_vect, rot_vect)
255
    molec.rotate(rot_angle, rot_vect, ctr_coord)
256
    if dummy_atom:
257
        del molec[-1]
258
    return molec
259

    
260

    
261
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
262
                  norm_vect=(0, 0, 1)):
263
    """Add an adsorbate to a surface.
264

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

    
315

    
316
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
317
                    nn_molec=0, coll_coeff=1.2):
318
    """Checks whether a slab and a molecule collide or not.
319

320
    @param slab_molec: The system of adsorbate-slab for which to detect if there
321
        are collisions.
322
    @param nn_slab: Number of neigbors in the surface.
323
    @param nn_molec: Number of neighbors in the molecule.
324
    @param coll_coeff: The coefficient that multiplies the covalent radius of
325
        atoms resulting in a distance that two atoms being closer to that is
326
        considered as atomic collision.
327
    @param slab_num_atoms: Number of atoms of the bare slab.
328
    @param min_height: The minimum height atoms can have to not be considered as
329
        colliding.
330
    @param vect: The vector perpendicular to the slab.
331
    @return: bool, whether the surface and the molecule collide.
332
    """
333
    from ase.neighborlist import natural_cutoffs, neighbor_list
334

    
335
    # Check structure overlap by height
336
    if min_height is not False:
337
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
338
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
339
        if vect.tolist() not in cart_axes:
340
            err_msg = "'min_coll_height' option is only implemented for " \
341
                      "'surf_norm_vect' to be one of the x, y or z axes. "
342
            logger.error(err_msg)
343
            raise ValueError(err_msg)
344
        for atom in slab_molec[slab_num_atoms:]:
345
            for i, coord in enumerate(vect):
346
                if coord == 0:
347
                    continue
348
                if atom.position[i] * coord < min_height * coord:
349
                    return True
350

    
351
    # Check structure overlap by sphere collision
352
    if coll_coeff is not False:
353
        slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
354
        slab_molec_nghbs = len(
355
            neighbor_list("i", slab_molec, 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):
365
    # TODO Rethink 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
    @return collision: bool, whether the structure generated has collisions
386
        between slab and adsorbate.
387
    """
388
    from copy import deepcopy
389
    slab_num_atoms = len(slab)
390
    slab_molec = []
391
    collision = True
392
    max_corr = 6  # Should be an even number
393
    d_angle = 180 / ((max_corr / 2.0) * num_pts)
394
    num_corr = 0
395
    while collision and num_corr <= max_corr:
396
        k = num_corr * (-1) ** num_corr
397
        slab_molec = deepcopy(slab)
398
        molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
399
                           center=ctr_coord)
400
        add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height,
401
                      norm_vect=norm_vect)
402
        collision = check_collision(slab_molec, slab_num_atoms, min_coll_height,
403
                                    norm_vect, slab_nghbs, molec_nghbs,
404
                                    coll_coeff)
405
        num_corr += 1
406
    return slab_molec, collision
407

    
408

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

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

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

    
456

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

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

    
506

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

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

    
563
    return slab_ads_list
564

    
565

    
566
def internal_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
567
                    ctr2_surf, bond_vector, bond_angle_target,
568
                    dihedral_angle_target=None, mol_dihedral_angle_target=None):
569
    """Performs translation and rotation of an adsorbate defined by an
570
    adsorption bond length, direction, angle and dihedral angle
571

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

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

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

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

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

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

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

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

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

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

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

616
    Returns:
617
        None (the `molecule` object is modified in-place)
618
    """
619
    vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf,
620
                                                                   ctr1_surf)
621
    vect_inter = get_atom_coords(molecule, ctr1_mol) \
622
        - get_atom_coords(surf, ctr1_surf)
623
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
624
                                                                     ctr1_mol)
625
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
626
                                                                      ctr2_mol)
627

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

    
643
    ###########################
644
    #       Translation       #
645
    ###########################
646

    
647
    # Compute and apply translation of adsorbate
648
    translation = bond_vector - vect_inter
649
    molecule.translate(translation)
650

    
651
    # Update adsorption bond
652
    vect_inter = get_atom_coords(molecule, ctr1_mol) - \
653
        get_atom_coords(surf, ctr1_surf)
654

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

    
664
    ###########################
665
    #   Bond angle rotation   #
666
    ###########################
667

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

    
682
    # Retrieve current bond angle
683
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
684

    
685
    # Apply rotation to reach desired bond_angle
686
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
687
                    center=get_atom_coords(molecule, ctr1_mol))
688

    
689
    # Update molecular bonds
690
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
691
                                                                     ctr1_mol)
692
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
693
                                                                      ctr2_mol)
694

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

    
709
    ###########################
710
    # Dihedral angle rotation #
711
    ###########################
712

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

    
729
        # Apply dihedral rotation along vect_inter
730
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
731
                        v=vect_inter, center=get_atom_coords(molecule,
732
                                                             ctr1_mol))
733

    
734
        # Update molecular bonds
735
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
736
            get_atom_coords(molecule, ctr1_mol)
737
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
738
            get_atom_coords(molecule, ctr2_mol)
739

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

    
767
    #####################################
768
    # Adsorbate dihedral angle rotation #
769
    #####################################
770

    
771
    # Perform adsorbate dihedral rotation if possible
772
    if do_mol_dihedral:
773
        # Retrieve current adsorbate dihedral angle (by computing the angle
774
        # between the orthogonal rejection of vect_inter and vect2_mol onto
775
        # vect_mol)
776
        vect_mol_inner = np.dot(vect_mol, vect_mol)
777
        bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \
778
            vect_mol_inner * vect_mol
779
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
780
            vect_mol_inner * vect_mol
781
        dihedral_angle_ini = get_vect_angle(bond_inter_reject,
782
                                            bond2_mol_reject, vect_mol)
783

    
784
        # Apply dihedral rotation along vect_mol
785
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
786
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
787

    
788
        # Update molecular bonds
789
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
790
            - get_atom_coords(molecule, ctr1_mol)
791
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
792
            - get_atom_coords(molecule, ctr2_mol)
793

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

    
831

    
832
def ads_internal(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
833
                 ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
834
                 slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel):
835
    """Generates adsorbate-surface structures by sampling over internal angles.
836

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

    
898
    return slab_ads_list
899

    
900

    
901
def adsorb_confs(conf_list, surf, inp_vars):
902
    """Generates a number of adsorbate-surface structure coordinates.
903

904
    Given a list of conformers, a surface, a list of atom indices (or list of
905
    list of indices) of both the surface and the adsorbate, it generates a
906
    number of adsorbate-surface structures for every possible combination of
907
    them at different orientations.
908
    @param conf_list: list of ase.Atoms of the different conformers
909
    @param surf: the ase.Atoms object of the surface
910
    @param inp_vars: Calculation parameters from input file.
911
    @return: list of ase.Atoms for the adsorbate-surface structures
912
    """
913
    from ase.neighborlist import natural_cutoffs, neighbor_list
914
    molec_ctrs = inp_vars['molec_ctrs']
915
    sites = inp_vars['sites']
916
    angles = inp_vars['set_angles']
917
    num_pts = inp_vars['sample_points_per_angle']
918
    inp_norm_vect = inp_vars['surf_norm_vect']
919
    min_coll_height = inp_vars['min_coll_height']
920
    coll_coeff = inp_vars['collision_threshold']
921
    h_donor = inp_vars['h_donor']
922
    h_acceptor = inp_vars['h_acceptor']
923

    
924
    if inp_vars['pbc_cell'] is not False:
925
        surf.set_pbc(True)
926
        surf.set_cell(inp_vars['pbc_cell'])
927

    
928
    surf_ads_list = []
929
    sites_coords = get_atom_coords(surf, sites)
930
    if coll_coeff is not False:
931
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
932
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
933
    else:
934
        surf_nghbs = 0
935
    for i, conf in enumerate(conf_list):
936
        molec_ctr_coords = get_atom_coords(conf, molec_ctrs)
937
        if inp_vars['pbc_cell'] is not False:
938
            conf.set_pbc(True)
939
            conf.set_cell(inp_vars['pbc_cell'])
940
        if coll_coeff is not False:
941
            conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
942
            molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
943
        else:
944
            molec_nghbs = 0
945
        for s, site in enumerate(sites_coords):
946
            if isinstance(inp_norm_vect, str) and inp_norm_vect == 'auto':
947
                norm_vect = compute_norm_vect(surf, sites[s],
948
                                              inp_vars['pbc_cell'])
949
            else:
950
                norm_vect = inp_norm_vect
951
            for c, ctr in enumerate(molec_ctr_coords):
952
                if angles == 'euler':
953
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
954
                                                   num_pts, min_coll_height,
955
                                                   coll_coeff, norm_vect,
956
                                                   surf_nghbs, molec_nghbs,
957
                                                   h_donor, h_acceptor))
958
                elif angles == 'internal':
959
                    mol_ctr1 = molec_ctrs[c]
960
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
961
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
962
                    surf_ctr1 = sites[s]
963
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
964
                    max_h = inp_vars["max_helic_angle"]
965
                    surf_ads_list.extend(ads_internal(conf, surf, mol_ctr1,
966
                                                      mol_ctr2, mol_ctr3,
967
                                                      surf_ctr1, surf_ctr2,
968
                                                      num_pts, min_coll_height,
969
                                                      coll_coeff, norm_vect,
970
                                                      surf_nghbs, molec_nghbs,
971
                                                      h_donor, h_acceptor,
972
                                                      max_h))
973
    return surf_ads_list
974

    
975

    
976
def run_screening(inp_vars):
977
    """Carries out the screening of adsorbate structures on a surface.
978

979
    @param inp_vars: Calculation parameters from input file.
980
    """
981
    import os
982
    import random
983
    from modules.formats import collect_coords, adapt_format
984
    from modules.calculation import run_calc, check_finished_calcs
985

    
986
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
987
                ' structures.')
988
    if inp_vars['use_molec_file']:
989
        selected_confs = [adapt_format('ase', inp_vars['use_molec_file'])]
990
        logger.info(f"Using '{inp_vars['use_molec_file']}' as only conformer.")
991
    else:
992
        if not os.path.isdir("isolated"):
993
            err = "'isolated' directory not found. It is needed in order to " \
994
                  "carry out the screening of structures to be adsorbed"
995
            logger.error(err)
996
            raise FileNotFoundError(err)
997

    
998
        correct_calcs, failed_calcs = check_finished_calcs('isolated',
999
                                                           inp_vars['code'])
1000
        if not correct_calcs:
1001
            err_msg = "No calculations on 'isolated' finished normally."
1002
            logger.error(err_msg)
1003
            raise FileNotFoundError(err_msg)
1004

    
1005
        logger.info(f"Found {len(correct_calcs)} structures of isolated "
1006
                    f"conformers whose calculation finished normally.")
1007
        if len(failed_calcs) != 0:
1008
            logger.warning(
1009
                f"Found {len(failed_calcs)} calculations more that "
1010
                f"did not finish normally: {failed_calcs}. \n"
1011
                f"Using only the ones that finished normally: "
1012
                f"{correct_calcs}.")
1013

    
1014
        conformer_atoms_list = collect_coords(correct_calcs, inp_vars['code'],
1015
                                              'isolated',
1016
                                              inp_vars['special_atoms'])
1017
        selected_confs = select_confs(conformer_atoms_list, correct_calcs,
1018
                                      inp_vars['select_magns'],
1019
                                      inp_vars['confs_per_magn'],
1020
                                      inp_vars['code'])
1021
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
1022
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
1023
    if len(surf_ads_list) > inp_vars['max_structures']:
1024
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
1025
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
1026
                f'configurations to carry out a calculation of.')
1027

    
1028
    run_calc('screening', inp_vars, surf_ads_list)
1029
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1030
                ' structures section.')