dockonsurf / modules / screening.py @ 5f3f4b69
Historique | Voir | Annoter | Télécharger (12,16 ko)
1 |
import logging |
---|---|
2 |
import numpy as np |
3 |
import ase |
4 |
|
5 |
logger = logging.getLogger('DockOnSurf')
|
6 |
|
7 |
|
8 |
def get_vect_angle(v1, v2, degrees=False): |
9 |
"""Computes the angle between two vectors.
|
10 |
|
11 |
@param v1: The first vector.
|
12 |
@param v2: The second vector.
|
13 |
@param degrees: Whether the result should be in radians (True) or in
|
14 |
degrees (False).
|
15 |
@return: The angle in radians if degrees = False, or in degrees if
|
16 |
degrees =True
|
17 |
"""
|
18 |
v1_u = v1 / np.linalg.norm(v1) |
19 |
v2_u = v2 / np.linalg.norm(v2) |
20 |
angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0)) |
21 |
return angle if not degrees else angle * 180 / np.pi |
22 |
|
23 |
|
24 |
def vect_avg(vects): |
25 |
"""Computes the element-wise mean of a set of vectors.
|
26 |
|
27 |
@param vects: list of lists-like: containing the vectors (num_vectors,
|
28 |
length_vector).
|
29 |
@return: vector average computed doing the element-wise mean.
|
30 |
"""
|
31 |
from utilities import try_command |
32 |
err = "vect_avg parameter vects must be a list-like, able to be converted" \
|
33 |
" np.array"
|
34 |
array = try_command(np.array, [(ValueError, err)], vects)
|
35 |
if len(array.shape) == 1: |
36 |
return array
|
37 |
else:
|
38 |
num_vects = array.shape[1]
|
39 |
return np.array([np.average(array[:, i]) for i in range(num_vects)]) |
40 |
|
41 |
|
42 |
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None): |
43 |
"""Gets the coordinates of the specified indices from a ase.Atoms object.
|
44 |
|
45 |
Given an ase.Atoms object and a list of atom indices specified in ctrs_list
|
46 |
it gets the coordinates of the specified atoms. If the element in the
|
47 |
ctrs_list is not an index but yet a list of indices, it computes the
|
48 |
element-wise mean of the coordinates of the atoms specified in the inner
|
49 |
list.
|
50 |
@param atoms: ase.Atoms object for which to obtain the coordinates of.
|
51 |
@param ctrs_list: list of (indices/list of indices) of the atoms for which
|
52 |
the coordinates should be extracted.
|
53 |
@return: np.ndarray of atomic coordinates.
|
54 |
"""
|
55 |
coords = [] |
56 |
err = "'ctrs_list' argument must be an integer, a list of integers or a " \
|
57 |
"list of lists of integers. Every integer must be in the range " \
|
58 |
"[0, num_atoms)"
|
59 |
if ctrs_list is None: |
60 |
ctrs_list = range(len(atoms)) |
61 |
elif isinstance(ctrs_list, int): |
62 |
if ctrs_list not in range(len(atoms)): |
63 |
logger.error(err) |
64 |
raise ValueError(err) |
65 |
return atoms[ctrs_list].position
|
66 |
for elem in ctrs_list: |
67 |
if isinstance(elem, list): |
68 |
coords.append(vect_avg([atoms[c].position for c in elem])) |
69 |
elif isinstance(elem, int): |
70 |
coords.append(atoms[elem].position) |
71 |
else:
|
72 |
|
73 |
logger.error(err) |
74 |
raise ValueError |
75 |
return np.array(coords)
|
76 |
|
77 |
|
78 |
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None, |
79 |
norm_vect=(0, 0, 1)): |
80 |
"""Add an adsorbate to a surface.
|
81 |
|
82 |
This function extends the functionality of ase.build.add_adsorbate
|
83 |
(https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
|
84 |
by enabling to change the z coordinate and the axis perpendicular to the
|
85 |
surface.
|
86 |
@param slab: ase.Atoms object containing the coordinates of the surface
|
87 |
@param adsorbate: ase.Atoms object containing the coordinates of the
|
88 |
adsorbate.
|
89 |
@param site_coord: The coordinates of the adsorption site on the surface.
|
90 |
@param ctr_coord: The coordinates of the adsorption center in the molecule.
|
91 |
@param height: The height above the surface where to adsorb.
|
92 |
@param offset: Offsets the adsorbate by a number of unit cells. Mostly
|
93 |
useful when adding more than one adsorbate.
|
94 |
@param norm_vect: The vector perpendicular to the surface.
|
95 |
"""
|
96 |
from copy import deepcopy |
97 |
info = slab.info.get('adsorbate_info', {})
|
98 |
pos = np.array([0.0, 0.0, 0.0]) # part of absolute coordinates |
99 |
spos = np.array([0.0, 0.0, 0.0]) # part relative to unit cell |
100 |
norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect) |
101 |
if offset is not None: |
102 |
spos += np.asarray(offset, float)
|
103 |
if isinstance(site_coord, str): |
104 |
# A site-name:
|
105 |
if 'sites' not in info: |
106 |
raise TypeError('If the atoms are not made by an ase.build ' |
107 |
'function, position cannot be a name.')
|
108 |
if site_coord not in info['sites']: |
109 |
raise TypeError('Adsorption site %s not supported.' % site_coord) |
110 |
spos += info['sites'][site_coord]
|
111 |
else:
|
112 |
pos += site_coord |
113 |
if 'cell' in info: |
114 |
cell = info['cell']
|
115 |
else:
|
116 |
cell = slab.get_cell() |
117 |
pos += np.dot(spos, cell) |
118 |
# Convert the adsorbate to an Atoms object
|
119 |
if isinstance(adsorbate, ase.Atoms): |
120 |
ads = deepcopy(adsorbate) |
121 |
elif isinstance(adsorbate, ase.Atom): |
122 |
ads = ase.Atoms([adsorbate]) |
123 |
else:
|
124 |
# Assume it is a string representing a single Atom
|
125 |
ads = ase.Atoms([ase.Atom(adsorbate)]) |
126 |
pos += height * norm_vect_u |
127 |
# Move adsorbate into position
|
128 |
ads.translate(pos - ctr_coord) |
129 |
# Attach the adsorbate
|
130 |
slab.extend(ads) |
131 |
|
132 |
|
133 |
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab, |
134 |
nn_molec, coll_coeff): |
135 |
"""Checks whether a slab and a molecule collide or not.
|
136 |
|
137 |
@param slab_molec: The system of adsorbate-slab for which to detect if there
|
138 |
are collisions.
|
139 |
@param nn_slab: Number of neigbors in the surface.
|
140 |
@param nn_molec: Number of neighbors in the molecule.
|
141 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
142 |
atoms resulting in a distance that two atoms being closer to that is
|
143 |
considered as atomic collision.
|
144 |
@param slab_num_atoms: Number of atoms of the bare slab.
|
145 |
@param min_height: The minimum height atoms can have to not be considered as
|
146 |
colliding.
|
147 |
@param vect: The vector perpendicular to the slab.
|
148 |
@return: bool, whether the surface and the molecule collide.
|
149 |
"""
|
150 |
from ase.neighborlist import natural_cutoffs, neighbor_list |
151 |
if (vect == np.array([1.0, 0.0, 0.0])).all \ |
152 |
or (vect == np.array([0.0, 1.0, 0.0])).all \ |
153 |
or (vect == np.array([0.0, 0.0, 1.0])).all: |
154 |
for atom in slab_molec[slab_num_atoms:]: |
155 |
if np.linalg.norm(atom.position * vect) < min_height:
|
156 |
return True |
157 |
slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff) |
158 |
slab_molec_nghbs = len(neighbor_list("i", slab_molec, slab_molec_cutoffs)) |
159 |
if slab_molec_nghbs > nn_slab + nn_molec:
|
160 |
return True |
161 |
else:
|
162 |
return False |
163 |
|
164 |
|
165 |
|
166 |
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts, |
167 |
coll_bttm, coll_coeff, norm_vect, slab_nghbs, molec_nghbs): |
168 |
"""Generates adsorbate-surface structures by sampling over Euler angles.
|
169 |
|
170 |
This function generates a number of adsorbate-surface structures at
|
171 |
different orientations of the adsorbate sampled at multiple Euler (zxz)
|
172 |
angles.
|
173 |
@param orig_molec: ase.Atoms object of the molecule to adsorb
|
174 |
@param slab: ase.Atoms object of the surface on which to adsorb the molecule
|
175 |
@param ctr_coord: The coordinates of the molecule to use as adsorption
|
176 |
center.
|
177 |
@param site_coord: The coordinates of the surface on which to adsorb the
|
178 |
molecule
|
179 |
@param num_pts: Number on which to sample Euler angles.
|
180 |
@param coll_bttm: The lowermost height for which to detect a collision
|
181 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
182 |
atoms resulting in a distance that two atoms being closer to that is
|
183 |
considered as atomic collision.
|
184 |
@param norm_vect: The vector perpendicular to the surface.
|
185 |
@param slab_nghbs: Number of neigbors in the surface.
|
186 |
@param molec_nghbs: Number of neighbors in the molecule.
|
187 |
@return: list of ase.Atoms object conatining all the orientations of a given
|
188 |
conformer
|
189 |
"""
|
190 |
from copy import deepcopy |
191 |
slab_ads_list = [] |
192 |
# rotation around z
|
193 |
for alpha in np.arange(0, 360, 360 / num_pts): |
194 |
# rotation around x'
|
195 |
for beta in np.arange(0, 180, 180 / num_pts): |
196 |
# rotation around z'
|
197 |
for gamma in np.arange(0, 360, 360 / num_pts): |
198 |
molec = deepcopy(orig_molec) |
199 |
molec.euler_rotate(alpha, beta, gamma, center=ctr_coord) |
200 |
slab_molec, collision = correct_coll(molec, slab, |
201 |
ctr_coord, site_coord, num_pts, coll_bttm, norm_vect, |
202 |
slab_nghbs, molec_nghbs, coll_coeff) |
203 |
if not collision: |
204 |
slab_ads_list.append(slab_molec) |
205 |
|
206 |
return slab_ads_list
|
207 |
|
208 |
|
209 |
def ads_chemcat(site, ctr, pts_angle): |
210 |
return "TO IMPLEMENT" |
211 |
|
212 |
|
213 |
def adsorb_confs(conf_list, surf, ads_ctrs, sites, algo, num_pts, neigh_ctrs, |
214 |
norm_vect, coll_bttm): |
215 |
"""Generates a number of adsorbate-surface structure coordinates.
|
216 |
|
217 |
Given a list of conformers, a surface, a list of atom indices (or list of
|
218 |
list of indices) of both the surface and the adsorbate, it generates a
|
219 |
number of adsorbate-surface structures for every possible combination of
|
220 |
them at different orientations.
|
221 |
@param conf_list: list of ase.Atoms of the different conformers
|
222 |
@param surf: the ase.Atoms object of the surface
|
223 |
@param ads_ctrs: the list atom indices of the adsorbate.
|
224 |
@param sites: the list of atom indices of the surface.
|
225 |
@param algo: the algorithm to use for the generation of adsorbates.
|
226 |
@param num_pts: the number of points per angle orientation to sample
|
227 |
@param neigh_ctrs: the indices of the neighboring atoms to the adsorption
|
228 |
atoms.
|
229 |
@param norm_vect: The vector perpendicular to the surface.
|
230 |
@param coll_bttm: The lowermost height for which to detect a collision
|
231 |
@return: list of ase.Atoms for the adsorbate-surface structures
|
232 |
"""
|
233 |
surf_ads_list = [] |
234 |
sites_coords = get_atom_coords(surf, sites) |
235 |
for conf in conf_list: |
236 |
molec_ctr_coords = get_atom_coords(conf, ads_ctrs) |
237 |
molec_neigh_coords = get_atom_coords(conf, neigh_ctrs) |
238 |
for site in sites_coords: |
239 |
for i, molec_ctr in enumerate(molec_ctr_coords): |
240 |
if algo == 'euler': |
241 |
surf_ads_list.extend(ads_euler(conf, surf, molec_ctr, site, |
242 |
num_pts, norm_vect, |
243 |
coll_bttm, |
244 |
molec_neigh_coords[i])) |
245 |
elif algo == 'chemcat': |
246 |
surf_ads_list.extend(ads_chemcat(site, molec_ctr, num_pts)) |
247 |
return surf_ads_list
|
248 |
|
249 |
|
250 |
def run_screening(inp_vars): |
251 |
"""Carry out the screening of adsorbate coordinates on a surface
|
252 |
|
253 |
@param inp_vars: Calculation parameters from input file.
|
254 |
"""
|
255 |
import os |
256 |
from modules.formats import read_coords, read_energies, \ |
257 |
rdkit_mol_to_ase_atoms, adapt_format |
258 |
from modules.clustering import get_rmsd, clustering |
259 |
from modules.isolated import get_moments_of_inertia |
260 |
from modules.calculation import run_calc |
261 |
|
262 |
if not os.path.isdir("isolated"): |
263 |
err = "'isolated' directory not found. It is needed in order to carry "
|
264 |
"out the screening of structures to be adsorbed"
|
265 |
logger.error(err) |
266 |
raise ValueError(err) |
267 |
|
268 |
conf_list = read_coords(inp_vars['code'], 'isolated', 'rdkit') |
269 |
# TODO Implement neighbors algorithm / align molecule to surface
|
270 |
# neigh_list = get_neighbors(conf_list[0], inp_vars['molec_ads_ctrs'])
|
271 |
conf_enrgs = read_energies(inp_vars['code'], 'isolated') |
272 |
mois = np.array([get_moments_of_inertia(conf)[0] for conf in conf_list]) |
273 |
rmsd_mtx = get_rmsd(conf_list) |
274 |
exemplars = clustering(rmsd_mtx) # TODO Kind of Clustering
|
275 |
clust_conf_list = [conf_list[idx] for idx in exemplars] |
276 |
conf_list = [rdkit_mol_to_ase_atoms(conf) for conf in clust_conf_list] |
277 |
surf = adapt_format('ase', inp_vars['surf_file']) |
278 |
surf_ads_list = adsorb_confs(conf_list, surf, inp_vars['molec_ads_ctrs'],
|
279 |
inp_vars['sites'], inp_vars['ads_algo'], |
280 |
inp_vars['sample_points_per_angle'],
|
281 |
inp_vars['molec_neigh_ctrs'],
|
282 |
inp_vars['surf_norm_vect'],
|
283 |
inp_vars['collision_bottom'])
|
284 |
run_calc('screening', inp_vars, surf_ads_list)
|