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

dockonsurf / modules / screening.py @ b4b2f307

Historique | Voir | Annoter | Télécharger (21,84 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, magns: list, num_sel: int, code: str):
13
    """Takes a list ase.Atoms and selects the most different magnitude-wise.
14

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

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

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

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

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

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

    
66

    
67
def get_vect_angle(v1, v2, degrees=False):
68
    """Computes the angle between two vectors.
69

70
    @param v1: The first vector.
71
    @param v2: The second vector.
72
    @param degrees: Whether the result should be in radians (True) or in
73
        degrees (False).
74
    @return: The angle in radians if degrees = False, or in degrees if
75
        degrees =True
76
    """
77
    v1_u = v1 / np.linalg.norm(v1)
78
    v2_u = v2 / np.linalg.norm(v2)
79
    angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0))
80
    return angle if not degrees else angle * 180 / np.pi
81

    
82

    
83
def vect_avg(vects):
84
    """Computes the element-wise mean of a set of vectors.
85

86
    @param vects: list of lists-like: containing the vectors (num_vectors,
87
        length_vector).
88
    @return: vector average computed doing the element-wise mean.
89
    """
90
    from utilities import try_command
91
    err = "vect_avg parameter vects must be a list-like, able to be converted" \
92
          " np.array"
93
    array = try_command(np.array, [(ValueError, err)], vects)
94
    if len(array.shape) == 1:
95
        return array
96
    else:
97
        num_vects = array.shape[1]
98
        return np.array([np.average(array[:, i]) for i in range(num_vects)])
99

    
100

    
101
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None):
102
    """Gets the coordinates of the specified indices from a ase.Atoms object.
103

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

    
132
            logger.error(err)
133
            raise ValueError
134
    return np.array(coords)
135

    
136

    
137
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
138
                  norm_vect=(0, 0, 1)):
139
    """Add an adsorbate to a surface.
140

141
    This function extends the functionality of ase.build.add_adsorbate
142
    (https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
143
    by enabling to change the z coordinate and the axis perpendicular to the
144
    surface.
145
    @param slab: ase.Atoms object containing the coordinates of the surface
146
    @param adsorbate: ase.Atoms object containing the coordinates of the
147
        adsorbate.
148
    @param site_coord: The coordinates of the adsorption site on the surface.
149
    @param ctr_coord: The coordinates of the adsorption center in the molecule.
150
    @param height: The height above the surface where to adsorb.
151
    @param offset: Offsets the adsorbate by a number of unit cells. Mostly
152
        useful when adding more than one adsorbate.
153
    @param norm_vect: The vector perpendicular to the surface.
154
    """
155
    from copy import deepcopy
156
    info = slab.info.get('adsorbate_info', {})
157
    pos = np.array([0.0, 0.0, 0.0])  # part of absolute coordinates
158
    spos = np.array([0.0, 0.0, 0.0])  # part relative to unit cell
159
    norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect)
160
    if offset is not None:
161
        spos += np.asarray(offset, float)
162
    if isinstance(site_coord, str):
163
        # A site-name:
164
        if 'sites' not in info:
165
            raise TypeError('If the atoms are not made by an ase.build '
166
                            'function, position cannot be a name.')
167
        if site_coord not in info['sites']:
168
            raise TypeError('Adsorption site %s not supported.' % site_coord)
169
        spos += info['sites'][site_coord]
170
    else:
171
        pos += site_coord
172
    if 'cell' in info:
173
        cell = info['cell']
174
    else:
175
        cell = slab.get_cell()
176
    pos += np.dot(spos, cell)
177
    # Convert the adsorbate to an Atoms object
178
    if isinstance(adsorbate, ase.Atoms):
179
        ads = deepcopy(adsorbate)
180
    elif isinstance(adsorbate, ase.Atom):
181
        ads = ase.Atoms([adsorbate])
182
    else:
183
        # Assume it is a string representing a single Atom
184
        ads = ase.Atoms([ase.Atom(adsorbate)])
185
    pos += height * norm_vect_u
186
    # Move adsorbate into position
187
    ads.translate(pos - ctr_coord)
188
    # Attach the adsorbate
189
    slab.extend(ads)
190

    
191

    
192
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
193
                    nn_molec=0, coll_coeff=0.9):
194
    """Checks whether a slab and a molecule collide or not.
195

196
    @param slab_molec: The system of adsorbate-slab for which to detect if there
197
        are collisions.
198
    @param nn_slab: Number of neigbors in the surface.
199
    @param nn_molec: Number of neighbors in the molecule.
200
    @param coll_coeff: The coefficient that multiplies the covalent radius of
201
        atoms resulting in a distance that two atoms being closer to that is
202
        considered as atomic collision.
203
    @param slab_num_atoms: Number of atoms of the bare slab.
204
    @param min_height: The minimum height atoms can have to not be considered as
205
        colliding.
206
    @param vect: The vector perpendicular to the slab.
207
    @return: bool, whether the surface and the molecule collide.
208
    """
209
    from ase.neighborlist import natural_cutoffs, neighbor_list
210
    if min_height is not False:
211
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
212
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
213
        if vect.tolist() in cart_axes:
214
            for atom in slab_molec[slab_num_atoms:]:
215
                for i, coord in enumerate(vect):
216
                    if coord == 0:
217
                        continue
218
                    if atom.position[i] * coord < min_height:
219
                        return True
220
            return False
221
    else:
222
        slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
223
        slab_molec_nghbs = len(
224
            neighbor_list("i", slab_molec, slab_molec_cutoffs))
225
        if slab_molec_nghbs > nn_slab + nn_molec:
226
            return True
227
        else:
228
            return False
229

    
230

    
231
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, neigh_cutoff=1):
232
    # TODO rethink
233
    """Moves a H atom from a molecule to the surface if the distance is short
234
    enough.
235

236
    @param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
237
    @param h_idx: The index of the hydrogen atom to carry out adsorption of.
238
    @param num_atoms_slab: The number of atoms of the slab without adsorbate.
239
    @param neigh_cutoff: half the maximum distance between the surface and the
240
        H for it to carry out dissociation.
241
    @return: An ase.Atoms object of the system adsorbate-surface with H
242
    """
243
    from copy import deepcopy
244
    from ase.neighborlist import NeighborList
245
    slab_molec = deepcopy(slab_molec_orig)
246
    cutoffs = len(slab_molec) * [neigh_cutoff]
247
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True)
248
    nl.update(slab_molec)
249
    surf_h_vect = np.array([np.infty] * 3)
250
    for neigh_idx in nl.get_neighbors(h_idx)[0]:
251
        if neigh_idx < num_atoms_slab:
252
            dist = np.linalg.norm(slab_molec[neigh_idx].position -
253
                                  slab_molec[h_idx].position)
254
            if dist < np.linalg.norm(surf_h_vect):
255
                surf_h_vect = slab_molec[neigh_idx].position \
256
                              - slab_molec[h_idx].position
257
    if np.linalg.norm(surf_h_vect) != np.infty:
258
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
259
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
260
        return slab_molec
261

    
262

    
263
def dissociation(slab_molec, disso_atoms, num_atoms_slab):
264
    # TODO rethink
265
    # TODO multiple dissociation
266
    """Decides which H atoms to dissociate according to a list of atoms.
267

268
    Given a list of chemical symbols or atom indices it checks for every atom
269
    or any of its neighbor if it's a H and calls dissociate_h to try to carry
270
    out dissociation of that H. For atom indices, it checks both whether
271
    the atom index or its neighbors are H, for chemical symbols, it only checks
272
    if there is a neighbor H.
273
    @param slab_molec: The ase.Atoms object of the system adsorbate-slab.
274
    @param disso_atoms: The indices or chemical symbols of the atoms
275
    @param num_atoms_slab:
276
    @return:
277
    """
278
    from ase.neighborlist import natural_cutoffs, NeighborList
279
    molec = slab_molec[num_atoms_slab:]
280
    cutoffs = natural_cutoffs(molec)
281
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True)
282
    nl.update(molec)
283
    disso_structs = []
284
    for el in disso_atoms:
285
        if isinstance(el, int):
286
            if molec[el].symbol == 'H':
287
                disso_struct = dissociate_h(slab_molec, el + num_atoms_slab,
288
                                            num_atoms_slab)
289
                if disso_struct is not None:
290
                    disso_structs.append(disso_struct)
291
            else:
292
                for neigh_idx in nl.get_neighbors(el)[0]:
293
                    if molec[neigh_idx].symbol == 'H':
294
                        disso_struct = dissociate_h(slab_molec, neigh_idx +
295
                                                    num_atoms_slab,
296
                                                    num_atoms_slab)
297
                        if disso_struct is not None:
298
                            disso_structs.append(disso_struct)
299
        else:
300
            for atom in molec:
301
                if atom.symbol.lower() == el.lower():
302
                    for neigh_idx in nl.get_neighbors(atom.index)[0]:
303
                        if molec[neigh_idx].symbol == 'H':
304
                            disso_struct = dissociate_h(slab_molec, neigh_idx \
305
                                                        + num_atoms_slab,
306
                                                        num_atoms_slab)
307
                            if disso_struct is not None:
308
                                disso_structs.append(disso_struct)
309
    return disso_structs
310

    
311

    
312
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
313
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
314
                 coll_coeff):
315
    # TODO Rethink this function
316
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
317
    small rotations.
318

319
    @param molec: ase.Atoms object of the molecule to adsorb
320
    @param slab: ase.Atoms object of the surface on which to adsorb the
321
        molecule
322
    @param ctr_coord: The coordinates of the molecule to use as adsorption
323
        center.
324
    @param site_coord: The coordinates of the surface on which to adsorb the
325
        molecule
326
    @param num_pts: Number on which to sample Euler angles.
327
    @param min_coll_height: The lowermost height for which to detect a collision
328
    @param norm_vect: The vector perpendicular to the surface.
329
    @param slab_nghbs: Number of neigbors in the surface.
330
    @param molec_nghbs: Number of neighbors in the molecule.
331
    @param coll_coeff: The coefficient that multiplies the covalent radius of
332
        atoms resulting in a distance that two atoms being closer to that is
333
        considered as atomic collision.
334
    @return collision: bool, whether the structure generated has collisions
335
        between slab and adsorbate.
336
    """
337
    from copy import deepcopy
338
    slab_num_atoms = len(slab)
339
    collision = True
340
    max_corr = 6  # Should be an even number
341
    d_angle = 180 / ((max_corr / 2.0) * num_pts)
342
    num_corr = 0
343
    while collision and num_corr <= max_corr:
344
        k = num_corr * (-1) ** num_corr
345
        slab_molec = deepcopy(slab)
346
        molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
347
                           center=ctr_coord)
348
        add_adsorbate(slab_molec, molec, site_coord, ctr_coord, 2.5,
349
                      norm_vect=norm_vect)
350
        collision = check_collision(slab_molec, slab_num_atoms, min_coll_height,
351
                                    norm_vect, slab_nghbs, molec_nghbs,
352
                                    coll_coeff)
353
        num_corr += 1
354
    return slab_molec, collision
355

    
356

    
357
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
358
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
359
              disso_atoms):
360
    """Generates adsorbate-surface structures by sampling over Euler angles.
361

362
    This function generates a number of adsorbate-surface structures at
363
    different orientations of the adsorbate sampled at multiple Euler (zxz)
364
    angles.
365
    @param orig_molec: ase.Atoms object of the molecule to adsorb
366
    @param slab: ase.Atoms object of the surface on which to adsorb the molecule
367
    @param ctr_coord: The coordinates of the molecule to use as adsorption
368
        center.
369
    @param site_coord: The coordinates of the surface on which to adsorb the
370
        molecule
371
    @param num_pts: Number on which to sample Euler angles.
372
    @param min_coll_height: The lowermost height for which to detect a collision
373
    @param coll_coeff: The coefficient that multiplies the covalent radius of
374
        atoms resulting in a distance that two atoms being closer to that is
375
        considered as atomic collision.
376
    @param norm_vect: The vector perpendicular to the surface.
377
    @param slab_nghbs: Number of neigbors in the surface.
378
    @param molec_nghbs: Number of neighbors in the molecule.
379
    @param disso_atoms: List of atom types or atom numbers to try to dissociate.
380
    @return: list of ase.Atoms object conatining all the orientations of a given
381
        conformer
382
    """
383
    from copy import deepcopy
384
    slab_ads_list = []
385
    # rotation around z
386
    for alpha in np.arange(0, 360, 360 / num_pts):
387
        # rotation around x'
388
        for beta in np.arange(0, 180, 180 / num_pts):
389
            # rotation around z'
390
            for gamma in np.arange(0, 360, 360 / num_pts):
391
                molec = deepcopy(orig_molec)
392
                molec.euler_rotate(alpha, beta, gamma, center=ctr_coord)
393
                slab_molec, collision = correct_coll(molec, slab,
394
                                                     ctr_coord, site_coord,
395
                                                     num_pts, min_coll_height,
396
                                                     norm_vect,
397
                                                     slab_nghbs, molec_nghbs,
398
                                                     coll_coeff)
399
                if not collision:
400
                    slab_ads_list.append(slab_molec)
401
                    slab_ads_list.extend(dissociation(slab_molec, disso_atoms,
402
                                                      len(slab)))
403

    
404
    return slab_ads_list
405

    
406

    
407
def ads_chemcat(site, ctr, pts_angle):
408
    return "TO IMPLEMENT"
409

    
410

    
411
def adsorb_confs(conf_list, surf, molec_ctrs, sites, algo, num_pts, neigh_ctrs,
412
                 norm_vect, min_coll_height, coll_coeff, disso_atoms):
413
    """Generates a number of adsorbate-surface structure coordinates.
414

415
    Given a list of conformers, a surface, a list of atom indices (or list of
416
    list of indices) of both the surface and the adsorbate, it generates a
417
    number of adsorbate-surface structures for every possible combination of
418
    them at different orientations.
419
    @param conf_list: list of ase.Atoms of the different conformers
420
    @param surf: the ase.Atoms object of the surface
421
    @param molec_ctrs: the list atom indices of the adsorbate.
422
    @param sites: the list of atom indices of the surface.
423
    @param algo: the algorithm to use for the generation of adsorbates.
424
    @param num_pts: the number of points per angle orientation to sample
425
    @param neigh_ctrs: the indices of the neighboring atoms to the adsorption
426
        atoms.
427
    @param norm_vect: The vector perpendicular to the surface.
428
    @param min_coll_height: The lowermost height for which to detect a collision
429
    @param coll_coeff: The coefficient that multiplies the covalent radius of
430
        atoms resulting in a distance that two atoms being closer to that is
431
        considered as atomic collision.
432
    @param disso_atoms: List of atom types or atom numbers to try to dissociate.
433
    @return: list of ase.Atoms for the adsorbate-surface structures
434
    """
435
    from ase.neighborlist import natural_cutoffs, neighbor_list
436
    surf_ads_list = []
437
    sites_coords = get_atom_coords(surf, sites)
438
    if min_coll_height is False:
439
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
440
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
441
    else:
442
        surf_nghbs = 0
443
    for conf in conf_list:
444
        molec_ctr_coords = get_atom_coords(conf, molec_ctrs)
445
        molec_neigh_coords = get_atom_coords(conf, neigh_ctrs)
446
        if min_coll_height is False:
447
            conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
448
            molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
449
        else:
450
            molec_nghbs = 0
451
        for site in sites_coords:
452
            for ctr in molec_ctr_coords:
453
                if algo == 'euler':
454
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
455
                                                   num_pts, min_coll_height,
456
                                                   coll_coeff, norm_vect,
457
                                                   surf_nghbs, molec_nghbs,
458
                                                   disso_atoms))
459
                elif algo == 'chemcat':
460
                    surf_ads_list.extend(ads_chemcat(site, ctr, num_pts))
461
    return surf_ads_list
462

    
463

    
464
def run_screening(inp_vars):
465
    """Carry out the screening of adsorbate coordinates on a surface
466

467
    @param inp_vars: Calculation parameters from input file.
468
    """
469
    import os
470
    from modules.formats import read_coords, adapt_format
471
    from modules.calculation import run_calc
472

    
473
    if not os.path.isdir("isolated"):
474
        err = "'isolated' directory not found. It is needed in order to carry "
475
        "out the screening of structures to be adsorbed"
476
        logger.error(err)
477
        raise ValueError(err)
478

    
479
    conf_list = read_coords(inp_vars['code'], 'isolated', 'ase')
480
    logger.info(f"Found {len(conf_list)} structures of isolated conformers.")
481
    selected_confs = select_confs(conf_list, inp_vars['select_magns'],
482
                                  inp_vars['confs_per_magn'],
483
                                  inp_vars['code'])
484
    surf = adapt_format('ase', inp_vars['surf_file'])
485
    surf_ads_list = adsorb_confs(selected_confs, surf,
486
                                 inp_vars['molec_ads_ctrs'], inp_vars['sites'],
487
                                 inp_vars['ads_algo'],
488
                                 inp_vars['sample_points_per_angle'],
489
                                 inp_vars['molec_neigh_ctrs'],
490
                                 inp_vars['surf_norm_vect'],
491
                                 inp_vars['min_coll_height'],
492
                                 inp_vars['collision_threshold'],
493
                                 inp_vars['disso_atoms'])
494
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
495
                f'configurations, to carry out a calculation of.')
496
    run_calc('screening', inp_vars, surf_ads_list)