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

dockonsurf / modules / screening.py @ cdc1edbe

Historique | Voir | Annoter | Télécharger (50,42 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
        if code == 'cp2k':
43
            conf_enrgs = collect_energies(calc_dirs, code, 'isolated')
44
        elif code == 'vasp':
45
            conf_enrgs = np.array([conf.get_total_energy()
46
                                   for conf in orig_conf_list])
47
    if 'moi' in magns:
48
        mois = np.array([conf.get_moments_of_inertia() for conf in conf_list])
49

    
50
    # Assign values
51
    for i, conf in enumerate(conf_list):
52
        assign_prop(conf, 'idx', i)
53
        assign_prop(conf, 'iso', calc_dirs[i])
54
        if 'energy' in magns:
55
            assign_prop(conf, 'energy', conf_enrgs[i])
56
        if 'moi' in magns:
57
            assign_prop(conf, 'moi', mois[i, 2])
58

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

    
70
    logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
71
    return [conf_list[idx] for idx in selected_ids]
72

    
73

    
74
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True):
75
    """Computes the angle between two vectors.
76

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

    
94
    return angle if not degrees else angle * 180 / np.pi
95

    
96

    
97
def vect_avg(vects):
98
    """Computes the element-wise mean of a set of vectors.
99

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

    
114

    
115
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None):
116
    """Gets the coordinates of the specified indices from a ase.Atoms object.
117

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

    
149

    
150
def compute_norm_vect(atoms, idxs, cell):
151
    """Computes the local normal vector of a surface at a given site.
152

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

    
183

    
184
def align_molec(orig_molec, ctr_coord, ref_vect):
185
    """Align a molecule to a vector by a center.
186

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

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

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

    
265

    
266
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
267
                  norm_vect=(0, 0, 1)):
268
    """Add an adsorbate to a surface.
269

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

    
320

    
321
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
322
                    nn_molec=0, coll_coeff=1.2, exclude_atom=False):
323
    """Checks whether a slab and a molecule collide or not.
324

325
    @param slab_molec: The system of adsorbate-slab for which to detect if there
326
        are collisions.
327
    @param nn_slab: Number of neigbors in the surface.
328
    @param nn_molec: Number of neighbors in the molecule.
329
    @param coll_coeff: The coefficient that multiplies the covalent radius of
330
        atoms resulting in a distance that two atoms being closer to that is
331
        considered as atomic collision.
332
    @param slab_num_atoms: Number of atoms of the bare slab.
333
    @param min_height: The minimum height atoms can have to not be considered as
334
        colliding.
335
    @param vect: The vector perpendicular to the slab.
336
    @param exclude_atom: Whether to exclude the adsorption center in the
337
        molecule in the collision detection.
338
    @return: bool, whether the surface and the molecule collide.
339
    """
340
    from copy import deepcopy
341
    from ase.neighborlist import natural_cutoffs, neighbor_list
342

    
343
    # Check structure overlap by height
344
    if min_height is not False:
345
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
346
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
347
        if vect.tolist() not in cart_axes:
348
            err_msg = "'min_coll_height' option is only implemented for " \
349
                      "'surf_norm_vect' to be one of the x, y or z axes. "
350
            logger.error(err_msg)
351
            raise ValueError(err_msg)
352
        for atom in slab_molec[slab_num_atoms:]:
353
            if exclude_atom is not False \
354
                    and atom.index == exclude_atom:
355
                continue
356
            for i, coord in enumerate(vect):
357
                if coord == 0:
358
                    continue
359
                if atom.position[i] * coord < min_height * coord:
360
                    return True
361

    
362
    # Check structure overlap by sphere collision
363
    if coll_coeff is not False:
364
        if exclude_atom is not False:
365
            slab_molec_wo_ctr = deepcopy(slab_molec)
366
            del slab_molec_wo_ctr[exclude_atom + slab_num_atoms]
367
            slab_molec_cutoffs = natural_cutoffs(slab_molec_wo_ctr,
368
                                                 mult=coll_coeff)
369
            slab_molec_nghbs = len(neighbor_list("i", slab_molec_wo_ctr,
370
                                                 slab_molec_cutoffs))
371
        else:
372
            slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
373
            slab_molec_nghbs = len(neighbor_list("i", slab_molec,
374
                                                 slab_molec_cutoffs))
375
        if slab_molec_nghbs > nn_slab + nn_molec:
376
            return True
377

    
378
    return False
379

    
380

    
381
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
382
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
383
                 coll_coeff, height=2.5, excl_atom=False):
384
    # TODO Rename this function
385
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
386
    small rotations.
387

388
    @param molec: ase.Atoms object of the molecule to adsorb
389
    @param slab: ase.Atoms object of the surface on which to adsorb the
390
        molecule
391
    @param ctr_coord: The coordinates of the molecule to use as adsorption
392
        center.
393
    @param site_coord: The coordinates of the surface on which to adsorb the
394
        molecule
395
    @param num_pts: Number on which to sample Euler angles.
396
    @param min_coll_height: The lowermost height for which to detect a collision
397
    @param norm_vect: The vector perpendicular to the surface.
398
    @param slab_nghbs: Number of neigbors in the surface.
399
    @param molec_nghbs: Number of neighbors in the molecule.
400
    @param coll_coeff: The coefficient that multiplies the covalent radius of
401
        atoms resulting in a distance that two atoms being closer to that is
402
        considered as atomic collision.
403
    @param height: Height on which to try adsorption.
404
    @param excl_atom: Whether to exclude the adsorption center in the
405
        molecule in the collision detection.
406
    @return collision: bool, whether the structure generated has collisions
407
        between slab and adsorbate.
408
    """
409
    from copy import deepcopy
410
    slab_num_atoms = len(slab)
411
    slab_molec = []
412
    collision = True
413
    max_corr = 6  # Should be an even number
414
    d_angle = 180 / ((max_corr / 2.0) * num_pts)
415
    num_corr = 0
416
    while collision and num_corr <= max_corr:
417
        k = num_corr * (-1) ** num_corr
418
        slab_molec = deepcopy(slab)
419
        molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
420
                           center=ctr_coord)
421
        add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height,
422
                      norm_vect=norm_vect)
423
        collision = check_collision(slab_molec, slab_num_atoms, min_coll_height,
424
                                    norm_vect, slab_nghbs, molec_nghbs,
425
                                    coll_coeff, excl_atom)
426
        num_corr += 1
427
    return slab_molec, collision
428

    
429

    
430
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
431
                 neigh_cutoff=1):
432
    # TODO rethink
433
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
434

435
    Tries to dissociate a H atom from the molecule and adsorb in on top of the
436
    surface if the distance is shorter than two times the neigh_cutoff value.
437
    @param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
438
    @param h_idx: The index of the hydrogen atom to carry out adsorption of.
439
    @param num_atoms_slab: The number of atoms of the slab without adsorbate.
440
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
441
    @param neigh_cutoff: half the maximum distance between the surface and the
442
        H for it to carry out dissociation.
443
    @return: An ase.Atoms object of the system adsorbate-surface with H
444
    """
445
    from copy import deepcopy
446
    from ase.neighborlist import NeighborList
447
    slab_molec = deepcopy(slab_molec_orig)
448
    cutoffs = len(slab_molec) * [neigh_cutoff]
449
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True, skin=0.0)
450
    nl.update(slab_molec)
451
    surf_h_vect = np.array([np.infty] * 3)
452
    if h_acceptor == 'all':
453
        h_acceptor = list(range(num_atoms_slab))
454
    for neigh_idx in nl.get_neighbors(h_idx)[0]:
455
        for elm in h_acceptor:
456
            if isinstance(elm, int):
457
                if neigh_idx == elm and neigh_idx < num_atoms_slab:
458
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
459
                                          slab_molec[h_idx].position)
460
                    if dist < np.linalg.norm(surf_h_vect):
461
                        surf_h_vect = slab_molec[neigh_idx].position \
462
                                      - slab_molec[h_idx].position
463
            else:
464
                if slab_molec[neigh_idx].symbol == elm \
465
                        and neigh_idx < num_atoms_slab:
466
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
467
                                          slab_molec[h_idx].position)
468
                    if dist < np.linalg.norm(surf_h_vect):
469
                        surf_h_vect = slab_molec[neigh_idx].position \
470
                                      - slab_molec[h_idx].position
471

    
472
    if np.linalg.norm(surf_h_vect) != np.infty:
473
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
474
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
475
        return slab_molec
476

    
477

    
478
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
479
    # TODO multiple dissociation
480
    """Decides which H atoms to dissociate according to a list of atoms.
481

482
    Given a list of chemical symbols or atom indices it checks for every atom
483
    or any of its neighbor if it's a H and calls dissociate_h to try to carry
484
    out dissociation of that H. For atom indices, it checks both whether
485
    the atom index or its neighbors are H, for chemical symbols, it only checks
486
    if there is a neighbor H.
487
    @param slab_molec: The ase.Atoms object of the system adsorbate-slab.
488
    @param h_donor: List of atom types or atom numbers that are H-donors.
489
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
490
    @param num_atoms_slab: Number of atoms of the bare slab.
491
    @return:
492
    """
493
    from ase.neighborlist import natural_cutoffs, NeighborList
494
    molec = slab_molec[num_atoms_slab:]
495
    cutoffs = natural_cutoffs(molec)
496
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True)
497
    nl.update(molec)
498
    disso_structs = []
499
    for el in h_donor:
500
        if isinstance(el, int):
501
            if molec[el].symbol == 'H':
502
                disso_struct = dissociate_h(slab_molec, el + num_atoms_slab,
503
                                            num_atoms_slab, h_acceptor)
504
                if disso_struct is not None:
505
                    disso_structs.append(disso_struct)
506
            else:
507
                for neigh_idx in nl.get_neighbors(el)[0]:
508
                    if molec[neigh_idx].symbol == 'H':
509
                        disso_struct = dissociate_h(slab_molec, neigh_idx +
510
                                                    num_atoms_slab,
511
                                                    num_atoms_slab, h_acceptor)
512
                        if disso_struct is not None:
513
                            disso_structs.append(disso_struct)
514
        else:
515
            for atom in molec:
516
                if atom.symbol.lower() == el.lower():
517
                    for neigh_idx in nl.get_neighbors(atom.index)[0]:
518
                        if molec[neigh_idx].symbol == 'H':
519
                            disso_struct = dissociate_h(slab_molec, neigh_idx
520
                                                        + num_atoms_slab,
521
                                                        num_atoms_slab,
522
                                                        h_acceptor)
523
                            if disso_struct is not None:
524
                                disso_structs.append(disso_struct)
525
    return disso_structs
526

    
527

    
528
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
529
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
530
              h_donor, h_acceptor, height, excl_atom):
531
    """Generates adsorbate-surface structures by sampling over Euler angles.
532

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

    
587
    return slab_ads_list
588

    
589

    
590
def internal_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
591
                    ctr2_surf, bond_vector, bond_angle_target,
592
                    dihedral_angle_target=None, mol_dihedral_angle_target=None):
593
    """Performs translation and rotation of an adsorbate defined by an
594
    adsorption bond length, direction, angle and dihedral angle
595

596
    Carles modification of chemcat's transform_adsorbate to work with
597
    coordinates instead of ase.Atom
598
    Parameters:
599
        molecule (ase.Atoms): The molecule to adsorb.
600

601
        surf (ase.Atoms): The surface ontop of which to adsorb.
602

603
        ctr1_mol (int/list): The position of the adsorption center in the
604
        molecule that will be bound to the surface.
605

606
        ctr2_mol (int/list): The position of a second center of the
607
        adsorbate used to define the adsorption bond angle, and the dihedral
608
        adsorption angle.
609

610
        ctr3_mol (int/list): The position of a third center in the
611
        adsorbate used to define the adsorbate dihedral angle.
612

613
        ctr1_surf (int/list): The position of the site on the surface that
614
        will be bound to the molecule.
615

616
        ctr2_surf (int/list): The position of a second center of the
617
        surface used to define the dihedral adsorption angle.
618

619
        bond_vector (numpy.ndarray): The adsorption bond desired.
620
            Details: offset = vect(atom1_surf, atom1_mol)
621

622
        bond_angle_target (float or int): The adsorption bond angle desired (in
623
            degrees).
624
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
625

626
        dihedral_angle_target (float or int): The dihedral adsorption angle
627
            desired (in degrees).
628
            Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
629
            atom1_mol, atom2_mol)
630
                The rotation vector is facing the adsorbate from the surface
631
                (i.e. counterclockwise rotation when facing the surface (i.e.
632
                view from top))
633

634
        mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
635
            desired (in degrees).
636
            Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
637
            atom2_mol, atom3_mol)
638
                The rotation vector is facing atom2_mol from atom1_mol
639

640
    Returns:
641
        None (the `molecule` object is modified in-place)
642
    """
643
    vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf,
644
                                                                   ctr1_surf)
645
    vect_inter = get_atom_coords(molecule, ctr1_mol) \
646
        - get_atom_coords(surf, ctr1_surf)
647
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
648
                                                                     ctr1_mol)
649
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
650
                                                                      ctr2_mol)
651

    
652
    # Check if dihedral angles can be defined
653
    do_dihedral = dihedral_angle_target is not None
654
    do_mol_dihedral = mol_dihedral_angle_target is not None
655
    dihedral_use_mol2 = False
656
    # Check if vect_surf and bond_vector are aligned
657
    if np.allclose(np.cross(vect_surf, bond_vector), 0):
658
        do_dihedral = False
659
    # Check if requested bond angle is not flat
660
    if np.isclose((bond_angle_target + 90) % 180 - 90, 0):
661
        do_mol_dihedral = False
662
        dihedral_use_mol2 = True
663
    # Check if vect_mol and vect2_mol are not aligned
664
    if np.allclose(np.cross(vect_mol, vect2_mol), 0):
665
        do_mol_dihedral = False
666

    
667
    ###########################
668
    #       Translation       #
669
    ###########################
670

    
671
    # Compute and apply translation of adsorbate
672
    translation = bond_vector - vect_inter
673
    molecule.translate(translation)
674

    
675
    # Update adsorption bond
676
    vect_inter = get_atom_coords(molecule, ctr1_mol) - \
677
        get_atom_coords(surf, ctr1_surf)
678

    
679
    # Check if translation was successful
680
    if np.allclose(vect_inter, bond_vector):
681
        pass  # print("Translation successfully applied (error: ~ {:.5g} unit "
682
        # "length)".format(np.linalg.norm(vect_inter - bond_vector)))
683
    else:
684
        err = 'An unknown error occured during the translation'
685
        logger.error(err)
686
        raise AssertionError(err)
687

    
688
    ###########################
689
    #   Bond angle rotation   #
690
    ###########################
691

    
692
    # Compute rotation vector
693
    rotation_vector = np.cross(-vect_inter, vect_mol)
694
    if np.allclose(rotation_vector, 0, atol=1e-5):
695
        # If molecular bonds are aligned, any vector orthogonal to vect_inter
696
        # can be used Such vector can be found as the orthogonal rejection of
697
        # either X-axis, Y-axis or Z-axis onto vect_inter (since they cannot
698
        # be all aligned)
699
        non_aligned_vector = np.zeros(3)
700
        # Select the most orthogonal axis (lowest dot product):
701
        non_aligned_vector[np.argmin(np.abs(vect_inter))] = 1
702
        rotation_vector = non_aligned_vector - np.dot(non_aligned_vector,
703
                                                      vect_inter) / np.dot(
704
            vect_inter, vect_inter) * vect_inter
705

    
706
    # Retrieve current bond angle
707
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
708

    
709
    # Apply rotation to reach desired bond_angle
710
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
711
                    center=get_atom_coords(molecule, ctr1_mol))
712

    
713
    # Update molecular bonds
714
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
715
                                                                     ctr1_mol)
716
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
717
                                                                      ctr2_mol)
718

    
719
    # Check if rotation was successful
720
    bond_angle = get_vect_angle(-vect_inter, vect_mol)
721
    if np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90, 0,
722
                  atol=1e-3) and np.allclose(get_atom_coords(molecule, ctr1_mol)
723
                                             - get_atom_coords(surf,
724
                                                               ctr1_surf),
725
                                             vect_inter):
726
        pass  # print("Rotation successfully applied (error: {:.5f}°)".format(
727
        # (bond_angle - bond_angle_target + 90) % 180 - 90))
728
    else:
729
        err = 'An unknown error occured during the rotation'
730
        logger.error(err)
731
        raise AssertionError(err)
732

    
733
    ###########################
734
    # Dihedral angle rotation #
735
    ###########################
736

    
737
    # Perform dihedral rotation if possible
738
    if do_dihedral:
739
        # Retrieve current dihedral angle (by computing the angle between the
740
        # vector rejection of vect_surf and vect_mol onto vect_inter)
741
        vect_inter_inner = np.dot(vect_inter, vect_inter)
742
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
743
            vect_inter_inner * vect_inter
744
        if dihedral_use_mol2:
745
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
746
                              vect_inter_inner * vect_inter
747
        else:
748
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
749
                              vect_inter_inner * vect_inter
750
        dihedral_angle_ini = get_vect_angle(vect_surf_reject, vect_mol_reject,
751
                                            vect_inter)
752

    
753
        # Apply dihedral rotation along vect_inter
754
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
755
                        v=vect_inter, center=get_atom_coords(molecule,
756
                                                             ctr1_mol))
757

    
758
        # Update molecular bonds
759
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
760
            get_atom_coords(molecule, ctr1_mol)
761
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
762
            get_atom_coords(molecule, ctr2_mol)
763

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

    
791
    #####################################
792
    # Adsorbate dihedral angle rotation #
793
    #####################################
794

    
795
    # Perform adsorbate dihedral rotation if possible
796
    if do_mol_dihedral:
797
        # Retrieve current adsorbate dihedral angle (by computing the angle
798
        # between the orthogonal rejection of vect_inter and vect2_mol onto
799
        # vect_mol)
800
        vect_mol_inner = np.dot(vect_mol, vect_mol)
801
        bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \
802
            vect_mol_inner * vect_mol
803
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
804
            vect_mol_inner * vect_mol
805
        dihedral_angle_ini = get_vect_angle(bond_inter_reject,
806
                                            bond2_mol_reject, vect_mol)
807

    
808
        # Apply dihedral rotation along vect_mol
809
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
810
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
811

    
812
        # Update molecular bonds
813
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
814
            - get_atom_coords(molecule, ctr1_mol)
815
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
816
            - get_atom_coords(molecule, ctr2_mol)
817

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

    
855

    
856
def ads_internal(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
857
                 ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
858
                 slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel, height,
859
                 excl_atom):
860
    """Generates adsorbate-surface structures by sampling over internal angles.
861

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

    
928
    return slab_ads_list
929

    
930

    
931
def adsorb_confs(conf_list, surf, inp_vars):
932
    """Generates a number of adsorbate-surface structure coordinates.
933

934
    Given a list of conformers, a surface, a list of atom indices (or list of
935
    list of indices) of both the surface and the adsorbate, it generates a
936
    number of adsorbate-surface structures for every possible combination of
937
    them at different orientations.
938
    @param conf_list: list of ase.Atoms of the different conformers
939
    @param surf: the ase.Atoms object of the surface
940
    @param inp_vars: Calculation parameters from input file.
941
    @return: list of ase.Atoms for the adsorbate-surface structures
942
    """
943
    from copy import deepcopy
944
    from ase.neighborlist import natural_cutoffs, neighbor_list
945
    molec_ctrs = inp_vars['molec_ctrs']
946
    sites = inp_vars['sites']
947
    angles = inp_vars['set_angles']
948
    num_pts = inp_vars['sample_points_per_angle']
949
    inp_norm_vect = inp_vars['surf_norm_vect']
950
    min_coll_height = inp_vars['min_coll_height']
951
    coll_coeff = inp_vars['collision_threshold']
952
    exclude_ads_ctr = inp_vars['exclude_ads_ctr']
953
    h_donor = inp_vars['h_donor']
954
    h_acceptor = inp_vars['h_acceptor']
955
    height = inp_vars['adsorption_height']
956

    
957
    if inp_vars['pbc_cell'] is not False:
958
        surf.set_pbc(True)
959
        surf.set_cell(inp_vars['pbc_cell'])
960

    
961
    surf_ads_list = []
962
    sites_coords = get_atom_coords(surf, sites)
963
    if coll_coeff is not False:
964
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
965
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
966
    else:
967
        surf_nghbs = 0
968
    for i, conf in enumerate(conf_list):
969
        molec_ctr_coords = get_atom_coords(conf, molec_ctrs)
970
        if inp_vars['pbc_cell'] is not False:
971
            conf.set_pbc(True)
972
            conf.set_cell(inp_vars['pbc_cell'])
973
        for s, site in enumerate(sites_coords):
974
            if isinstance(inp_norm_vect, str) and inp_norm_vect == 'auto':
975
                norm_vect = compute_norm_vect(surf, sites[s],
976
                                              inp_vars['pbc_cell'])
977
            else:
978
                norm_vect = inp_norm_vect
979
            for c, ctr in enumerate(molec_ctr_coords):
980
                if exclude_ads_ctr and isinstance(molec_ctrs[c], int):
981
                    exclude_atom = molec_ctrs[c]
982
                else:
983
                    exclude_atom = False
984
                    if exclude_ads_ctr and not isinstance(molec_ctrs[c], int):
985
                        logger.warning("'exclude_ads_ctr' only works for atomic"
986
                                       "centers and not for many-atoms "
987
                                       f"barycenters. {molec_ctrs[c]} are not "
988
                                       f"going to be excluded from collison.")
989
                if coll_coeff and exclude_atom is not False:
990
                    conf_wo_ctr = deepcopy(conf)
991
                    del conf_wo_ctr[exclude_atom]
992
                    conf_cutoffs = natural_cutoffs(conf_wo_ctr, mult=coll_coeff)
993
                    molec_nghbs = len(neighbor_list("i", conf_wo_ctr,
994
                                                    conf_cutoffs))
995
                elif coll_coeff and exclude_atom is False:
996
                    conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
997
                    molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
998
                else:
999
                    molec_nghbs = 0
1000
                if angles == 'euler':
1001
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
1002
                                                   num_pts, min_coll_height,
1003
                                                   coll_coeff, norm_vect,
1004
                                                   surf_nghbs, molec_nghbs,
1005
                                                   h_donor, h_acceptor, height,
1006
                                                   exclude_atom))
1007
                elif angles == 'internal':
1008
                    mol_ctr1 = molec_ctrs[c]
1009
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
1010
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
1011
                    surf_ctr1 = sites[s]
1012
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
1013
                    max_h = inp_vars["max_helic_angle"]
1014
                    surf_ads_list.extend(ads_internal(conf, surf, mol_ctr1,
1015
                                                      mol_ctr2, mol_ctr3,
1016
                                                      surf_ctr1, surf_ctr2,
1017
                                                      num_pts, min_coll_height,
1018
                                                      coll_coeff, norm_vect,
1019
                                                      surf_nghbs, molec_nghbs,
1020
                                                      h_donor, h_acceptor,
1021
                                                      max_h, height,
1022
                                                      exclude_atom))
1023
    return surf_ads_list
1024

    
1025

    
1026
def run_screening(inp_vars):
1027
    """Carries out the screening of adsorbate structures on a surface.
1028

1029
    @param inp_vars: Calculation parameters from input file.
1030
    """
1031
    import os
1032
    import random
1033
    from modules.formats import collect_coords, adapt_format
1034
    from modules.calculation import run_calc, check_finished_calcs
1035

    
1036
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
1037
                ' structures.')
1038
    if inp_vars['use_molec_file']:
1039
        selected_confs = [adapt_format('ase', inp_vars['use_molec_file'])]
1040
        logger.info(f"Using '{inp_vars['use_molec_file']}' as only conformer.")
1041
    else:
1042
        if not os.path.isdir("isolated"):
1043
            err = "'isolated' directory not found. It is needed in order to " \
1044
                  "carry out the screening of structures to be adsorbed"
1045
            logger.error(err)
1046
            raise FileNotFoundError(err)
1047

    
1048
        correct_calcs, failed_calcs = check_finished_calcs('isolated',
1049
                                                           inp_vars['code'])
1050
        if not correct_calcs:
1051
            err_msg = "No calculations on 'isolated' finished normally."
1052
            logger.error(err_msg)
1053
            raise FileNotFoundError(err_msg)
1054

    
1055
        logger.info(f"Found {len(correct_calcs)} structures of isolated "
1056
                    f"conformers whose calculation finished normally.")
1057
        if len(failed_calcs) != 0:
1058
            logger.warning(
1059
                f"Found {len(failed_calcs)} calculations more that "
1060
                f"did not finish normally: {failed_calcs}. \n"
1061
                f"Using only the ones that finished normally: "
1062
                f"{correct_calcs}.")
1063

    
1064
        conformer_atoms_list = collect_coords(correct_calcs, inp_vars['code'],
1065
                                              'isolated',
1066
                                              inp_vars['special_atoms'])
1067
        selected_confs = select_confs(conformer_atoms_list, correct_calcs,
1068
                                      inp_vars['select_magns'],
1069
                                      inp_vars['confs_per_magn'],
1070
                                      inp_vars['code'])
1071
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
1072
    surf.info = {}
1073
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
1074
    if len(surf_ads_list) > inp_vars['max_structures']:
1075
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
1076
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
1077
                f'configurations to carry out a calculation of.')
1078

    
1079
    run_calc('screening', inp_vars, surf_ads_list)
1080
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1081
                ' structures section.')