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

dockonsurf / modules / screening.py @ f16ebc80

Historique | Voir | Annoter | Télécharger (47,38 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 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 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(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 molec: The molecule to alugn
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 ase import Atom
192
    from ase.neighborlist import natural_cutoffs, neighbor_list
193

    
194
    if len(molec) == 1:
195
        err_msg = "Cannot align a single atom"
196
        logger.error(err_msg)
197
        ValueError(err_msg)
198
    cutoffs = natural_cutoffs(molec, mult=1.2)
199
    # Check if ctr_coord are the coordinates of an atom and if not creates a
200
    # dummy one to extract the neighboring atoms.
201
    ctr_idx = None
202
    dummy_atom = False
203
    for atom in molec:
204
        if np.allclose(ctr_coord, atom.position, rtol=1e-2):
205
            ctr_idx = atom.index
206
            break
207
    if ctr_idx is None:
208
        molec.append(Atom('X', position=ctr_coord))
209
        cutoffs.append(0.2)
210
        ctr_idx = len(molec) - 1
211
        dummy_atom = True
212
    # Builds the neighbors and computes the average vector
213
    refs, vects = neighbor_list("iD", molec, cutoffs, self_interaction=False)
214
    neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
215
    # If no neighbors are present, the cutoff of the alignment center is
216
    # set to a value where at least one atom is a neighbor and neighbors are
217
    # recalculated.
218
    if len(neigh_vects) == 0:
219
        min_dist, min_idx = (np.inf, np.inf)
220
        for atom in molec:
221
            if atom.index == ctr_idx:
222
                continue
223
            if molec.get_distance(ctr_idx, atom) < min_dist:
224
                min_dist = molec.get_distance(ctr_idx, atom.index)
225
                min_idx = atom.index
226
        cutoffs[ctr_idx] = min_dist - cutoffs[min_idx] + 0.05
227
        refs, vects = neighbor_list("iD", molec, cutoffs,
228
                                    self_interaction=False)
229
        neigh_vects = [vects[i] for i, atm in enumerate(refs) if atm == ctr_idx]
230
    target_vect = vect_avg(neigh_vects)
231
    # If the target vector is 0 it gets reassigned to the vector pointing
232
    # towards the nearest atom.
233
    if target_vect.tolist() == [0, 0, 0]:
234
        target_vect = np.array([np.inf] * 3)
235
        for vect in neigh_vects:
236
            if np.linalg.norm(vect) < np.linalg.norm(target_vect):
237
                target_vect = vect
238

    
239
    rot_vect = np.cross(target_vect, ref_vect)
240
    if np.allclose(rot_vect, 0):
241
        cart_axes = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]
242
        axis = cart_axes[int(np.argmax([np.linalg.norm(np.cross(ax, ref_vect))
243
                                        for ax in cart_axes]))]
244
        rot_vect = np.cross(ref_vect, axis)
245
    rot_angle = -get_vect_angle(ref_vect, target_vect, rot_vect)
246
    molec.rotate(rot_angle, rot_vect, ctr_coord)
247
    if dummy_atom:
248
        del molec[-1]
249
    return molec
250

    
251

    
252
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
253
                  norm_vect=(0, 0, 1)):
254
    """Add an adsorbate to a surface.
255

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

    
306

    
307
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
308
                    nn_molec=0, coll_coeff=1.2):
309
    """Checks whether a slab and a molecule collide or not.
310

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

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

    
342
    # Check structure overlap by sphere collision
343
    if coll_coeff is not False:
344
        slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
345
        slab_molec_nghbs = len(
346
            neighbor_list("i", slab_molec, slab_molec_cutoffs))
347
        if slab_molec_nghbs > nn_slab + nn_molec:
348
            return True
349

    
350
    return False
351

    
352

    
353
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
354
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
355
                 coll_coeff, height=2.5):
356
    # TODO Rethink this function
357
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
358
    small rotations.
359

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

    
399

    
400
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
401
                 neigh_cutoff=1):
402
    # TODO rethink
403
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
404

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

    
442
    if np.linalg.norm(surf_h_vect) != np.infty:
443
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
444
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
445
        return slab_molec
446

    
447

    
448
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
449
    # TODO multiple dissociation
450
    """Decides which H atoms to dissociate according to a list of atoms.
451

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

    
497

    
498
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
499
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
500
              h_donor, h_acceptor):
501
    """Generates adsorbate-surface structures by sampling over Euler angles.
502

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

    
560
    return slab_ads_list
561

    
562

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

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

574
        surf (ase.Atoms): The surface ontop of which to adsorb.
575

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

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

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

586
        ctr1_surf (int/list): The position of the site on the surface that
587
        will be bound to the molecule.
588

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

592
        bond_vector (numpy.ndarray): The adsorption bond desired.
593
            Details: offset = vect(atom1_surf, atom1_mol)
594

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

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

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

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

    
625
    # Check if dihedral angles can be defined
626
    do_dihedral = dihedral_angle_target is not None
627
    do_mol_dihedral = mol_dihedral_angle_target is not None
628
    dihedral_use_mol2 = False
629
    # Check if vect_surf and vect_inter are not aligned
630
    if np.allclose(np.cross(vect_surf, vect_inter), 0):
631
        logger.warning(
632
            "Surface atoms are incompatible with adsorption "
633
            "direction/bond. An adsorption dihedral angle cannot be defined.")
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
        logger.warning("Requested bond angle is flat. Only a single dihedral "
638
                       "angle can be defined (ctr2_surf, ctr1_surf, ctr2_mol, "
639
                       "ctr3_mol).")
640
        do_mol_dihedral = False
641
        dihedral_use_mol2 = True
642
        logger.warning("Dihedral adsorption angle rotation will be perfomed "
643
                       "with (ctr2_surf, ctr1_surf, ctr2_mol, ctr3_mol).")
644
    # Check if vect_mol and vect2_mol are not aligned
645
    if np.allclose(np.cross(vect_mol, vect2_mol), 0):
646
        logger.warning("Adsorbates atoms are aligned. An adsorbate dihedral "
647
                       "angle cannot be defined.")
648
        do_mol_dihedral = False
649
    if not do_dihedral:
650
        logger.warning("Dihedral adsorption angle rotation will not be "
651
                       "performed.")
652
    if not do_mol_dihedral:
653
        logger.warning("Adsorbate dihedral angle rotation will not be "
654
                       "performed.")
655

    
656
    ###########################
657
    #       Translation       #
658
    ###########################
659

    
660
    # Compute and apply translation of adsorbate
661
    translation = bond_vector - vect_inter
662
    molecule.translate(translation)
663

    
664
    # Update adsorption bond
665
    vect_inter = get_atom_coords(molecule, ctr1_mol) - \
666
        get_atom_coords(surf, ctr1_surf)
667

    
668
    # Check if translation was successful
669
    if np.allclose(vect_inter, bond_vector):
670
        pass  # print("Translation successfully applied (error: ~ {:.5g} unit "
671
        # "length)".format(np.linalg.norm(vect_inter - bond_vector)))
672
    else:
673
        err = 'An unknown error occured during the translation'
674
        logger.error(err)
675
        raise AssertionError(err)
676

    
677
    ###########################
678
    #   Bond angle rotation   #
679
    ###########################
680

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

    
695
    # Retrieve current bond angle
696
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
697

    
698
    # Apply rotation to reach desired bond_angle
699
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
700
                    center=get_atom_coords(molecule, ctr1_mol))
701

    
702
    # Update molecular bonds
703
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
704
                                                                     ctr1_mol)
705
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
706
                                                                      ctr2_mol)
707

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

    
722
    ###########################
723
    # Dihedral angle rotation #
724
    ###########################
725

    
726
    # Perform dihedral rotation if possible
727
    if do_dihedral:
728
        # Retrieve current dihedral angle (by computing the angle between the
729
        # vector rejection of vect_surf and vect_mol onto vect_inter)
730
        vect_inter_inner = np.dot(vect_inter, vect_inter)
731
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
732
            vect_inter_inner * vect_inter
733
        if dihedral_use_mol2:
734
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
735
                              vect_inter_inner * vect_inter
736
        else:
737
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
738
                              vect_inter_inner * vect_inter
739
        dihedral_angle_ini = get_vect_angle(vect_surf_reject, vect_mol_reject,
740
                                            vect_inter)
741

    
742
        # Apply dihedral rotation along vect_inter
743
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
744
                        v=vect_inter, center=get_atom_coords(molecule,
745
                                                             ctr1_mol))
746

    
747
        # Update molecular bonds
748
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
749
            get_atom_coords(molecule, ctr1_mol)
750
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
751
            get_atom_coords(molecule, ctr2_mol)
752

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

    
780
    #####################################
781
    # Adsorbate dihedral angle rotation #
782
    #####################################
783

    
784
    # Perform adsorbate dihedral rotation if possible
785
    if do_mol_dihedral:
786
        # Retrieve current adsorbate dihedral angle (by computing the angle
787
        # between the orthogonal rejection of vect_inter and vect2_mol onto
788
        # vect_mol)
789
        vect_mol_inner = np.dot(vect_mol, vect_mol)
790
        bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \
791
            vect_mol_inner * vect_mol
792
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
793
            vect_mol_inner * vect_mol
794
        dihedral_angle_ini = get_vect_angle(bond_inter_reject,
795
                                            bond2_mol_reject, vect_mol)
796

    
797
        # Apply dihedral rotation along vect_mol
798
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
799
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
800

    
801
        # Update molecular bonds
802
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
803
            - get_atom_coords(molecule, ctr1_mol)
804
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
805
            - get_atom_coords(molecule, ctr2_mol)
806

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

    
844

    
845
def ads_chemcat(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
846
                ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
847
                slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel):
848
    """Generates adsorbate-surface structures by sampling over chemcat angles.
849

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

    
907
    return slab_ads_list
908

    
909

    
910
def adsorb_confs(conf_list, surf, inp_vars):
911
    """Generates a number of adsorbate-surface structure coordinates.
912

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

    
933
    if inp_vars['pbc_cell'] is not False:
934
        surf.set_pbc(True)
935
        surf.set_cell(inp_vars['pbc_cell'])
936

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

    
982

    
983
def run_screening(inp_vars):
984
    """Carries out the screening of adsorbate structures on a surface.
985

986
    @param inp_vars: Calculation parameters from input file.
987
    """
988
    import os
989
    import random
990
    from modules.formats import collect_coords, adapt_format
991
    from modules.calculation import run_calc, check_finished_calcs
992

    
993
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
994
                ' structures.')
995
    if not os.path.isdir("isolated"):
996
        err = "'isolated' directory not found. It is needed in order to carry "
997
        "out the screening of structures to be adsorbed"
998
        logger.error(err)
999
        raise FileNotFoundError(err)
1000

    
1001
    finished_calcs, unfinished_calcs = check_finished_calcs('isolated',
1002
                                                            inp_vars['code'])
1003
    logger.info(f"Found {len(finished_calcs)} structures of isolated "
1004
                f"conformers whose calculation finished normally.")
1005
    if len(unfinished_calcs) != 0:
1006
        logger.warning(f"Found {len(unfinished_calcs)} calculations more that "
1007
                       f"did not finish normally: {unfinished_calcs}. \n"
1008
                       f"Using only the ones that finished normally: "
1009
                       f"{finished_calcs}.")
1010

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

    
1024
    run_calc('screening', inp_vars, surf_ads_list)
1025
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1026
                ' structures section.')