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

dockonsurf / modules / screening.py @ 9cd032cf

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

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

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

    
7

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

    
11

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

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

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

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

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

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

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

    
68

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

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

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

    
91

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

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

    
109

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

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

    
144

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

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

    
178

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

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

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

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

    
260

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

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

    
315

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

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

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

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

    
373
    return False
374

    
375

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

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

    
424

    
425
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
426
                 neigh_cutoff=1):
427
    # TODO rethink
428
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
429

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

    
467
    if np.linalg.norm(surf_h_vect) != np.infty:
468
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
469
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
470
        return slab_molec
471

    
472

    
473
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
474
    # TODO multiple dissociation
475
    """Decides which H atoms to dissociate according to a list of atoms.
476

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

    
522

    
523
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
524
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
525
              h_donor, h_acceptor, height, excl_atom):
526
    """Generates adsorbate-surface structures by sampling over Euler angles.
527

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

    
582
    return slab_ads_list
583

    
584

    
585
def internal_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
586
                    ctr2_surf, bond_vector, bond_angle_target,
587
                    dihedral_angle_target=None, mol_dihedral_angle_target=None):
588
    """Performs translation and rotation of an adsorbate defined by an
589
    adsorption bond length, direction, angle and dihedral angle
590

591
    Carles modification of chemcat's transform_adsorbate to work with
592
    coordinates instead of ase.Atom
593
    Parameters:
594
        molecule (ase.Atoms): The molecule to adsorb.
595

596
        surf (ase.Atoms): The surface ontop of which to adsorb.
597

598
        ctr1_mol (int/list): The position of the adsorption center in the
599
        molecule that will be bound to the surface.
600

601
        ctr2_mol (int/list): The position of a second center of the
602
        adsorbate used to define the adsorption bond angle, and the dihedral
603
        adsorption angle.
604

605
        ctr3_mol (int/list): The position of a third center in the
606
        adsorbate used to define the adsorbate dihedral angle.
607

608
        ctr1_surf (int/list): The position of the site on the surface that
609
        will be bound to the molecule.
610

611
        ctr2_surf (int/list): The position of a second center of the
612
        surface used to define the dihedral adsorption angle.
613

614
        bond_vector (numpy.ndarray): The adsorption bond desired.
615
            Details: offset = vect(atom1_surf, atom1_mol)
616

617
        bond_angle_target (float or int): The adsorption bond angle desired (in
618
            degrees).
619
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
620

621
        dihedral_angle_target (float or int): The dihedral adsorption angle
622
            desired (in degrees).
623
            Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
624
            atom1_mol, atom2_mol)
625
                The rotation vector is facing the adsorbate from the surface
626
                (i.e. counterclockwise rotation when facing the surface (i.e.
627
                view from top))
628

629
        mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
630
            desired (in degrees).
631
            Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
632
            atom2_mol, atom3_mol)
633
                The rotation vector is facing atom2_mol from atom1_mol
634

635
    Returns:
636
        None (the `molecule` object is modified in-place)
637
    """
638
    vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf,
639
                                                                   ctr1_surf)
640
    vect_inter = get_atom_coords(molecule, ctr1_mol) \
641
        - get_atom_coords(surf, ctr1_surf)
642
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
643
                                                                     ctr1_mol)
644
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
645
                                                                      ctr2_mol)
646

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

    
662
    ###########################
663
    #       Translation       #
664
    ###########################
665

    
666
    # Compute and apply translation of adsorbate
667
    translation = bond_vector - vect_inter
668
    molecule.translate(translation)
669

    
670
    # Update adsorption bond
671
    vect_inter = get_atom_coords(molecule, ctr1_mol) - \
672
        get_atom_coords(surf, ctr1_surf)
673

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

    
683
    ###########################
684
    #   Bond angle rotation   #
685
    ###########################
686

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

    
701
    # Retrieve current bond angle
702
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
703

    
704
    # Apply rotation to reach desired bond_angle
705
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
706
                    center=get_atom_coords(molecule, ctr1_mol))
707

    
708
    # Update molecular bonds
709
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
710
                                                                     ctr1_mol)
711
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
712
                                                                      ctr2_mol)
713

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

    
728
    ###########################
729
    # Dihedral angle rotation #
730
    ###########################
731

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

    
748
        # Apply dihedral rotation along vect_inter
749
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
750
                        v=vect_inter, center=get_atom_coords(molecule,
751
                                                             ctr1_mol))
752

    
753
        # Update molecular bonds
754
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
755
            get_atom_coords(molecule, ctr1_mol)
756
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
757
            get_atom_coords(molecule, ctr2_mol)
758

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

    
786
    #####################################
787
    # Adsorbate dihedral angle rotation #
788
    #####################################
789

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

    
803
        # Apply dihedral rotation along vect_mol
804
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
805
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
806

    
807
        # Update molecular bonds
808
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
809
            - get_atom_coords(molecule, ctr1_mol)
810
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
811
            - get_atom_coords(molecule, ctr2_mol)
812

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

    
850

    
851
def ads_internal(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
852
                 ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
853
                 slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel, height,
854
                 excl_atom):
855
    """Generates adsorbate-surface structures by sampling over internal angles.
856

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

    
922
    return slab_ads_list
923

    
924

    
925
def adsorb_confs(conf_list, surf, inp_vars):
926
    """Generates a number of adsorbate-surface structure coordinates.
927

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

    
951
    if inp_vars['pbc_cell'] is not False:
952
        surf.set_pbc(True)
953
        surf.set_cell(inp_vars['pbc_cell'])
954

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

    
1019

    
1020
def run_screening(inp_vars):
1021
    """Carries out the screening of adsorbate structures on a surface.
1022

1023
    @param inp_vars: Calculation parameters from input file.
1024
    """
1025
    import os
1026
    import random
1027
    from modules.formats import collect_coords, adapt_format
1028
    from modules.calculation import run_calc, check_finished_calcs
1029

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

    
1042
        correct_calcs, failed_calcs = check_finished_calcs('isolated',
1043
                                                           inp_vars['code'])
1044
        if not correct_calcs:
1045
            err_msg = "No calculations on 'isolated' finished normally."
1046
            logger.error(err_msg)
1047
            raise FileNotFoundError(err_msg)
1048

    
1049
        logger.info(f"Found {len(correct_calcs)} structures of isolated "
1050
                    f"conformers whose calculation finished normally.")
1051
        if len(failed_calcs) != 0:
1052
            logger.warning(
1053
                f"Found {len(failed_calcs)} calculations more that "
1054
                f"did not finish normally: {failed_calcs}. \n"
1055
                f"Using only the ones that finished normally: "
1056
                f"{correct_calcs}.")
1057

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

    
1072
    run_calc('screening', inp_vars, surf_ads_list)
1073
    logger.info('Finished the procedures for the screening of adsorbate-surface'
1074
                ' structures section.')