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

dockonsurf / modules / screening.py @ 58ede1f9

Historique | Voir | Annoter | Télécharger (50,51 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, center=None):
116
    """Gets the coordinates of the specified center for an ase.Atoms object.
117

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

    
145

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

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

    
179

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

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

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

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

    
261

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

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

    
316

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

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

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

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

    
374
    return False
375

    
376

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

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

    
425

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

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

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

    
473

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

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

    
523

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

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

    
583
    return slab_ads_list
584

    
585

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

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

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

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

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

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

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

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

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

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

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

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

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

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

    
663
    ###########################
664
    #       Translation       #
665
    ###########################
666

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

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

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

    
684
    ###########################
685
    #   Bond angle rotation   #
686
    ###########################
687

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

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

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

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

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

    
729
    ###########################
730
    # Dihedral angle rotation #
731
    ###########################
732

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

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

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

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

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

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

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

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

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

    
851

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

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

    
924
    return slab_ads_list
925

    
926

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

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

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

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

    
1021

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

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

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

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

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

    
1061
        conformer_atoms_list = collect_coords(correct_calcs, inp_vars['code'],
1062
                                              'isolated',
1063
                                              inp_vars['special_atoms'])
1064
        selected_confs = select_confs(conformer_atoms_list, correct_calcs,
1065
                                      inp_vars['select_magns'],
1066
                                      inp_vars['confs_per_magn'],
1067
                                      inp_vars['code'])
1068
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
1069
    surf.info = {}
1070
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
1071
    if len(surf_ads_list) > inp_vars['max_structures']:
1072
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
1073
    elif len(surf_ads_list) == 0:
1074
        err_msg = "No configurations were generated: Check the parameters in" \
1075
                  "dockonsurf.inp"
1076
        logger.error(err_msg)
1077
        raise ValueError(err_msg)
1078
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
1079
                f'configurations to carry out a calculation of.')
1080

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