dockonsurf / modules / screening.py @ e1c5f171
Historique | Voir | Annoter | Télécharger (23,57 ko)
1 |
import logging |
---|---|
2 |
import numpy as np |
3 |
import ase |
4 |
|
5 |
logger = logging.getLogger('DockOnSurf')
|
6 |
|
7 |
|
8 |
def assign_prop(atoms: ase.Atoms, prop_name: str, prop_val): # TODO Needed? |
9 |
atoms.info[prop_name] = prop_val |
10 |
|
11 |
|
12 |
def select_confs(orig_conf_list: list, calc_dirs: list, magns: list, |
13 |
num_sel: int, code: str): |
14 |
"""Takes a list ase.Atoms and selects the most different magnitude-wise.
|
15 |
|
16 |
Given a list of ase.Atoms objects and a list of magnitudes, it selects a
|
17 |
number of the most different conformers according to every magnitude
|
18 |
specified.
|
19 |
|
20 |
@param orig_conf_list: list of ase.Atoms objects to select among.
|
21 |
@param calc_dirs: List of directories where to read the energies from.
|
22 |
@param magns: list of str with the names of the magnitudes to use for the
|
23 |
conformer selection. Supported magnitudes: 'energy', 'moi'.
|
24 |
@param num_sel: number of conformers to select for every of the magnitudes.
|
25 |
@param code: The code that generated the magnitude information.
|
26 |
Supported codes: See formats.py
|
27 |
@return: list of the selected ase.Atoms objects.
|
28 |
"""
|
29 |
from copy import deepcopy |
30 |
from modules.formats import collect_energies |
31 |
|
32 |
conf_list = deepcopy(orig_conf_list) |
33 |
conf_enrgs, mois, selected_ids = [], [], [] |
34 |
if num_sel >= len(conf_list): |
35 |
logger.warning('Number of conformers per magnitude is equal or larger '
|
36 |
'than the total number of conformers. Using all '
|
37 |
f'available conformers: {len(conf_list)}.')
|
38 |
return conf_list
|
39 |
|
40 |
# Read properties
|
41 |
if 'energy' in magns: |
42 |
conf_enrgs = collect_energies(calc_dirs, code, 'isolated')
|
43 |
if 'moi' in magns: |
44 |
mois = np.array([conf.get_moments_of_inertia() for conf in conf_list]) |
45 |
|
46 |
# Assign values
|
47 |
for i, conf in enumerate(conf_list): |
48 |
assign_prop(conf, 'idx', i)
|
49 |
if 'energy' in magns: |
50 |
assign_prop(conf, 'energy', conf_enrgs[i])
|
51 |
if 'moi' in magns: |
52 |
assign_prop(conf, 'moi', mois[i, 2]) |
53 |
|
54 |
# pick ids
|
55 |
for magn in magns: |
56 |
sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn])) |
57 |
if sorted_list[-1].info['idx'] not in selected_ids: |
58 |
selected_ids.append(sorted_list[-1].info['idx']) |
59 |
if num_sel > 1: |
60 |
for i in range(0, len(sorted_list) - 1, |
61 |
len(conf_list) // (num_sel - 1)): |
62 |
if sorted_list[i].info['idx'] not in selected_ids: |
63 |
selected_ids.append(sorted_list[i].info['idx'])
|
64 |
|
65 |
logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
|
66 |
return [conf_list[idx] for idx in selected_ids] |
67 |
|
68 |
|
69 |
def get_vect_angle(v1: list, v2: list, degrees=False): |
70 |
"""Computes the angle between two vectors.
|
71 |
|
72 |
@param v1: The first vector.
|
73 |
@param v2: The second vector.
|
74 |
@param degrees: Whether the result should be in radians (True) or in
|
75 |
degrees (False).
|
76 |
@return: The angle in radians if degrees = False, or in degrees if
|
77 |
degrees =True
|
78 |
"""
|
79 |
v1_u = v1 / np.linalg.norm(v1) |
80 |
v2_u = v2 / np.linalg.norm(v2) |
81 |
angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0)) |
82 |
return angle if not degrees else angle * 180 / np.pi |
83 |
|
84 |
|
85 |
def vect_avg(vects): |
86 |
"""Computes the element-wise mean of a set of vectors.
|
87 |
|
88 |
@param vects: list of lists-like: containing the vectors (num_vectors,
|
89 |
length_vector).
|
90 |
@return: vector average computed doing the element-wise mean.
|
91 |
"""
|
92 |
from utilities import try_command |
93 |
err = "vect_avg parameter vects must be a list-like, able to be converted" \
|
94 |
" np.array"
|
95 |
array = try_command(np.array, [(ValueError, err)], vects)
|
96 |
if len(array.shape) == 1: |
97 |
return array
|
98 |
else:
|
99 |
num_vects = array.shape[1]
|
100 |
return np.array([np.average(array[:, i]) for i in range(num_vects)]) |
101 |
|
102 |
|
103 |
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None): |
104 |
"""Gets the coordinates of the specified indices from a ase.Atoms object.
|
105 |
|
106 |
Given an ase.Atoms object and a list of atom indices specified in ctrs_list
|
107 |
it gets the coordinates of the specified atoms. If the element in the
|
108 |
ctrs_list is not an index but yet a list of indices, it computes the
|
109 |
element-wise mean of the coordinates of the atoms specified in the inner
|
110 |
list.
|
111 |
@param atoms: ase.Atoms object for which to obtain the coordinates of.
|
112 |
@param ctrs_list: list of (indices/list of indices) of the atoms for which
|
113 |
the coordinates should be extracted.
|
114 |
@return: np.ndarray of atomic coordinates.
|
115 |
"""
|
116 |
coords = [] |
117 |
err = "'ctrs_list' argument must be an integer, a list of integers or a " \
|
118 |
"list of lists of integers. Every integer must be in the range " \
|
119 |
"[0, num_atoms)"
|
120 |
if ctrs_list is None: |
121 |
ctrs_list = range(len(atoms)) |
122 |
elif isinstance(ctrs_list, int): |
123 |
if ctrs_list not in range(len(atoms)): |
124 |
logger.error(err) |
125 |
raise ValueError(err) |
126 |
return atoms[ctrs_list].position
|
127 |
for elem in ctrs_list: |
128 |
if isinstance(elem, list): |
129 |
coords.append(vect_avg([atoms[c].position for c in elem])) |
130 |
elif isinstance(elem, int): |
131 |
coords.append(atoms[elem].position) |
132 |
else:
|
133 |
|
134 |
logger.error(err) |
135 |
raise ValueError |
136 |
return np.array(coords)
|
137 |
|
138 |
|
139 |
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None, |
140 |
norm_vect=(0, 0, 1)): |
141 |
"""Add an adsorbate to a surface.
|
142 |
|
143 |
This function extends the functionality of ase.build.add_adsorbate
|
144 |
(https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
|
145 |
by enabling to change the z coordinate and the axis perpendicular to the
|
146 |
surface.
|
147 |
@param slab: ase.Atoms object containing the coordinates of the surface
|
148 |
@param adsorbate: ase.Atoms object containing the coordinates of the
|
149 |
adsorbate.
|
150 |
@param site_coord: The coordinates of the adsorption site on the surface.
|
151 |
@param ctr_coord: The coordinates of the adsorption center in the molecule.
|
152 |
@param height: The height above the surface where to adsorb.
|
153 |
@param offset: Offsets the adsorbate by a number of unit cells. Mostly
|
154 |
useful when adding more than one adsorbate.
|
155 |
@param norm_vect: The vector perpendicular to the surface.
|
156 |
"""
|
157 |
from copy import deepcopy |
158 |
info = slab.info.get('adsorbate_info', {})
|
159 |
pos = np.array([0.0, 0.0, 0.0]) # part of absolute coordinates |
160 |
spos = np.array([0.0, 0.0, 0.0]) # part relative to unit cell |
161 |
norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect) |
162 |
if offset is not None: |
163 |
spos += np.asarray(offset, float)
|
164 |
if isinstance(site_coord, str): |
165 |
# A site-name:
|
166 |
if 'sites' not in info: |
167 |
raise TypeError('If the atoms are not made by an ase.build ' |
168 |
'function, position cannot be a name.')
|
169 |
if site_coord not in info['sites']: |
170 |
raise TypeError('Adsorption site %s not supported.' % site_coord) |
171 |
spos += info['sites'][site_coord]
|
172 |
else:
|
173 |
pos += site_coord |
174 |
if 'cell' in info: |
175 |
cell = info['cell']
|
176 |
else:
|
177 |
cell = slab.get_cell() |
178 |
pos += np.dot(spos, cell) |
179 |
# Convert the adsorbate to an Atoms object
|
180 |
if isinstance(adsorbate, ase.Atoms): |
181 |
ads = deepcopy(adsorbate) |
182 |
elif isinstance(adsorbate, ase.Atom): |
183 |
ads = ase.Atoms([adsorbate]) |
184 |
else:
|
185 |
# Assume it is a string representing a single Atom
|
186 |
ads = ase.Atoms([ase.Atom(adsorbate)]) |
187 |
pos += height * norm_vect_u |
188 |
# Move adsorbate into position
|
189 |
ads.translate(pos - ctr_coord) |
190 |
# Attach the adsorbate
|
191 |
slab.extend(ads) |
192 |
|
193 |
|
194 |
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0, |
195 |
nn_molec=0, coll_coeff=0.9): |
196 |
"""Checks whether a slab and a molecule collide or not.
|
197 |
|
198 |
@param slab_molec: The system of adsorbate-slab for which to detect if there
|
199 |
are collisions.
|
200 |
@param nn_slab: Number of neigbors in the surface.
|
201 |
@param nn_molec: Number of neighbors in the molecule.
|
202 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
203 |
atoms resulting in a distance that two atoms being closer to that is
|
204 |
considered as atomic collision.
|
205 |
@param slab_num_atoms: Number of atoms of the bare slab.
|
206 |
@param min_height: The minimum height atoms can have to not be considered as
|
207 |
colliding.
|
208 |
@param vect: The vector perpendicular to the slab.
|
209 |
@return: bool, whether the surface and the molecule collide.
|
210 |
"""
|
211 |
from ase.neighborlist import natural_cutoffs, neighbor_list |
212 |
if min_height is not False: |
213 |
cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0], |
214 |
[-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]] |
215 |
if vect.tolist() in cart_axes: |
216 |
for atom in slab_molec[slab_num_atoms:]: |
217 |
for i, coord in enumerate(vect): |
218 |
if coord == 0: |
219 |
continue
|
220 |
if atom.position[i] * coord < min_height:
|
221 |
return True |
222 |
return False |
223 |
else:
|
224 |
slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff) |
225 |
slab_molec_nghbs = len(
|
226 |
neighbor_list("i", slab_molec, slab_molec_cutoffs))
|
227 |
if slab_molec_nghbs > nn_slab + nn_molec:
|
228 |
return True |
229 |
else:
|
230 |
return False |
231 |
|
232 |
|
233 |
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, neigh_cutoff=1): |
234 |
# TODO rethink
|
235 |
"""Tries to dissociate a H from the molecule and adsorbs it on the slab.
|
236 |
|
237 |
Tries to dissociate a H atom from the molecule and adsorb in on top of the
|
238 |
surface if the distance is shorter than two times the neigh_cutoff value.
|
239 |
@param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
|
240 |
@param h_idx: The index of the hydrogen atom to carry out adsorption of.
|
241 |
@param num_atoms_slab: The number of atoms of the slab without adsorbate.
|
242 |
@param neigh_cutoff: half the maximum distance between the surface and the
|
243 |
H for it to carry out dissociation.
|
244 |
@return: An ase.Atoms object of the system adsorbate-surface with H
|
245 |
"""
|
246 |
from copy import deepcopy |
247 |
from ase.neighborlist import NeighborList |
248 |
slab_molec = deepcopy(slab_molec_orig) |
249 |
cutoffs = len(slab_molec) * [neigh_cutoff]
|
250 |
nl = NeighborList(cutoffs, self_interaction=False, bothways=True) |
251 |
nl.update(slab_molec) |
252 |
surf_h_vect = np.array([np.infty] * 3)
|
253 |
for neigh_idx in nl.get_neighbors(h_idx)[0]: |
254 |
if neigh_idx < num_atoms_slab:
|
255 |
dist = np.linalg.norm(slab_molec[neigh_idx].position - |
256 |
slab_molec[h_idx].position) |
257 |
if dist < np.linalg.norm(surf_h_vect):
|
258 |
surf_h_vect = slab_molec[neigh_idx].position \ |
259 |
- slab_molec[h_idx].position |
260 |
if np.linalg.norm(surf_h_vect) != np.infty:
|
261 |
trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect) |
262 |
slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect |
263 |
return slab_molec
|
264 |
|
265 |
|
266 |
def dissociation(slab_molec, disso_atoms, num_atoms_slab): |
267 |
# TODO rethink
|
268 |
# TODO multiple dissociation
|
269 |
"""Decides which H atoms to dissociate according to a list of atoms.
|
270 |
|
271 |
Given a list of chemical symbols or atom indices it checks for every atom
|
272 |
or any of its neighbor if it's a H and calls dissociate_h to try to carry
|
273 |
out dissociation of that H. For atom indices, it checks both whether
|
274 |
the atom index or its neighbors are H, for chemical symbols, it only checks
|
275 |
if there is a neighbor H.
|
276 |
@param slab_molec: The ase.Atoms object of the system adsorbate-slab.
|
277 |
@param disso_atoms: The indices or chemical symbols of the atoms
|
278 |
@param num_atoms_slab:
|
279 |
@return:
|
280 |
"""
|
281 |
from ase.neighborlist import natural_cutoffs, NeighborList |
282 |
molec = slab_molec[num_atoms_slab:] |
283 |
cutoffs = natural_cutoffs(molec) |
284 |
nl = NeighborList(cutoffs, self_interaction=False, bothways=True) |
285 |
nl.update(molec) |
286 |
disso_structs = [] |
287 |
for el in disso_atoms: |
288 |
if isinstance(el, int): |
289 |
if molec[el].symbol == 'H': |
290 |
disso_struct = dissociate_h(slab_molec, el + num_atoms_slab, |
291 |
num_atoms_slab) |
292 |
if disso_struct is not None: |
293 |
disso_structs.append(disso_struct) |
294 |
else:
|
295 |
for neigh_idx in nl.get_neighbors(el)[0]: |
296 |
if molec[neigh_idx].symbol == 'H': |
297 |
disso_struct = dissociate_h(slab_molec, neigh_idx + |
298 |
num_atoms_slab, |
299 |
num_atoms_slab) |
300 |
if disso_struct is not None: |
301 |
disso_structs.append(disso_struct) |
302 |
else:
|
303 |
for atom in molec: |
304 |
if atom.symbol.lower() == el.lower():
|
305 |
for neigh_idx in nl.get_neighbors(atom.index)[0]: |
306 |
if molec[neigh_idx].symbol == 'H': |
307 |
disso_struct = dissociate_h(slab_molec, neigh_idx \ |
308 |
+ num_atoms_slab, |
309 |
num_atoms_slab) |
310 |
if disso_struct is not None: |
311 |
disso_structs.append(disso_struct) |
312 |
return disso_structs
|
313 |
|
314 |
|
315 |
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts, |
316 |
min_coll_height, norm_vect, slab_nghbs, molec_nghbs, |
317 |
coll_coeff, height=2.5):
|
318 |
# TODO Rethink this function
|
319 |
"""Tries to adsorb a molecule on a slab trying to avoid collisions by doing
|
320 |
small rotations.
|
321 |
|
322 |
@param molec: ase.Atoms object of the molecule to adsorb
|
323 |
@param slab: ase.Atoms object of the surface on which to adsorb the
|
324 |
molecule
|
325 |
@param ctr_coord: The coordinates of the molecule to use as adsorption
|
326 |
center.
|
327 |
@param site_coord: The coordinates of the surface on which to adsorb the
|
328 |
molecule
|
329 |
@param num_pts: Number on which to sample Euler angles.
|
330 |
@param min_coll_height: The lowermost height for which to detect a collision
|
331 |
@param norm_vect: The vector perpendicular to the surface.
|
332 |
@param slab_nghbs: Number of neigbors in the surface.
|
333 |
@param molec_nghbs: Number of neighbors in the molecule.
|
334 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
335 |
atoms resulting in a distance that two atoms being closer to that is
|
336 |
considered as atomic collision.
|
337 |
@param height: Height on which to try adsorption
|
338 |
@return collision: bool, whether the structure generated has collisions
|
339 |
between slab and adsorbate.
|
340 |
"""
|
341 |
from copy import deepcopy |
342 |
slab_num_atoms = len(slab)
|
343 |
collision = True
|
344 |
max_corr = 6 # Should be an even number |
345 |
d_angle = 180 / ((max_corr / 2.0) * num_pts) |
346 |
num_corr = 0
|
347 |
while collision and num_corr <= max_corr: |
348 |
k = num_corr * (-1) ** num_corr
|
349 |
slab_molec = deepcopy(slab) |
350 |
molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
|
351 |
center=ctr_coord) |
352 |
add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height, |
353 |
norm_vect=norm_vect) |
354 |
collision = check_collision(slab_molec, slab_num_atoms, min_coll_height, |
355 |
norm_vect, slab_nghbs, molec_nghbs, |
356 |
coll_coeff) |
357 |
num_corr += 1
|
358 |
return slab_molec, collision
|
359 |
|
360 |
|
361 |
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts, |
362 |
min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs, |
363 |
disso_atoms): |
364 |
"""Generates adsorbate-surface structures by sampling over Euler angles.
|
365 |
|
366 |
This function generates a number of adsorbate-surface structures at
|
367 |
different orientations of the adsorbate sampled at multiple Euler (zxz)
|
368 |
angles.
|
369 |
@param orig_molec: ase.Atoms object of the molecule to adsorb
|
370 |
@param slab: ase.Atoms object of the surface on which to adsorb the molecule
|
371 |
@param ctr_coord: The coordinates of the molecule to use as adsorption
|
372 |
center.
|
373 |
@param site_coord: The coordinates of the surface on which to adsorb the
|
374 |
molecule
|
375 |
@param num_pts: Number on which to sample Euler angles.
|
376 |
@param min_coll_height: The lowermost height for which to detect a collision
|
377 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
378 |
atoms resulting in a distance that two atoms being closer to that is
|
379 |
considered as atomic collision.
|
380 |
@param norm_vect: The vector perpendicular to the surface.
|
381 |
@param slab_nghbs: Number of neigbors in the surface.
|
382 |
@param molec_nghbs: Number of neighbors in the molecule.
|
383 |
@param disso_atoms: List of atom types or atom numbers to try to dissociate.
|
384 |
@return: list of ase.Atoms object conatining all the orientations of a given
|
385 |
conformer
|
386 |
"""
|
387 |
from copy import deepcopy |
388 |
slab_ads_list = [] |
389 |
# rotation around z
|
390 |
for alpha in np.arange(0, 360, 360 / num_pts): |
391 |
# rotation around x'
|
392 |
for beta in np.arange(0, 180, 180 / num_pts): |
393 |
# rotation around z'
|
394 |
for gamma in np.arange(0, 360, 360 / num_pts): |
395 |
molec = deepcopy(orig_molec) |
396 |
molec.euler_rotate(alpha, beta, gamma, center=ctr_coord) |
397 |
slab_molec, collision = correct_coll(molec, slab, |
398 |
ctr_coord, site_coord, |
399 |
num_pts, min_coll_height, |
400 |
norm_vect, |
401 |
slab_nghbs, molec_nghbs, |
402 |
coll_coeff) |
403 |
if not collision and \ |
404 |
not any([(slab_molec.positions == conf.positions).all() |
405 |
for conf in slab_ads_list]): |
406 |
slab_ads_list.append(slab_molec) |
407 |
slab_ads_list.extend(dissociation(slab_molec, disso_atoms, |
408 |
len(slab)))
|
409 |
|
410 |
return slab_ads_list
|
411 |
|
412 |
|
413 |
def ads_chemcat(site, ctr, pts_angle): |
414 |
return "TO IMPLEMENT" |
415 |
|
416 |
|
417 |
def adsorb_confs(conf_list, surf, molec_ctrs, sites, algo, num_pts, neigh_ctrs, |
418 |
norm_vect, min_coll_height, coll_coeff, disso_atoms): |
419 |
"""Generates a number of adsorbate-surface structure coordinates.
|
420 |
|
421 |
Given a list of conformers, a surface, a list of atom indices (or list of
|
422 |
list of indices) of both the surface and the adsorbate, it generates a
|
423 |
number of adsorbate-surface structures for every possible combination of
|
424 |
them at different orientations.
|
425 |
@param conf_list: list of ase.Atoms of the different conformers
|
426 |
@param surf: the ase.Atoms object of the surface
|
427 |
@param molec_ctrs: the list atom indices of the adsorbate.
|
428 |
@param sites: the list of atom indices of the surface.
|
429 |
@param algo: the algorithm to use for the generation of adsorbates.
|
430 |
@param num_pts: the number of points per angle orientation to sample
|
431 |
@param neigh_ctrs: the indices of the neighboring atoms to the adsorption
|
432 |
atoms.
|
433 |
@param norm_vect: The vector perpendicular to the surface.
|
434 |
@param min_coll_height: The lowermost height for which to detect a collision
|
435 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
436 |
atoms resulting in a distance that two atoms being closer to that is
|
437 |
considered as atomic collision.
|
438 |
@param disso_atoms: List of atom types or atom numbers to try to dissociate.
|
439 |
@return: list of ase.Atoms for the adsorbate-surface structures
|
440 |
"""
|
441 |
from ase.neighborlist import natural_cutoffs, neighbor_list |
442 |
surf_ads_list = [] |
443 |
sites_coords = get_atom_coords(surf, sites) |
444 |
if min_coll_height is False: |
445 |
surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff) |
446 |
surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs)) |
447 |
else:
|
448 |
surf_nghbs = 0
|
449 |
for conf in conf_list: |
450 |
molec_ctr_coords = get_atom_coords(conf, molec_ctrs) |
451 |
molec_neigh_coords = get_atom_coords(conf, neigh_ctrs) |
452 |
if min_coll_height is False: |
453 |
conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff) |
454 |
molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs)) |
455 |
else:
|
456 |
molec_nghbs = 0
|
457 |
for site in sites_coords: |
458 |
for ctr in molec_ctr_coords: |
459 |
if algo == 'euler': |
460 |
surf_ads_list.extend(ads_euler(conf, surf, ctr, site, |
461 |
num_pts, min_coll_height, |
462 |
coll_coeff, norm_vect, |
463 |
surf_nghbs, molec_nghbs, |
464 |
disso_atoms)) |
465 |
elif algo == 'chemcat': |
466 |
surf_ads_list.extend(ads_chemcat(site, ctr, num_pts)) |
467 |
return surf_ads_list
|
468 |
|
469 |
|
470 |
def run_screening(inp_vars): |
471 |
"""Carries out the screening of adsorbate structures on a surface.
|
472 |
|
473 |
@param inp_vars: Calculation parameters from input file.
|
474 |
"""
|
475 |
import os |
476 |
import random |
477 |
from modules.formats import collect_coords, adapt_format |
478 |
from modules.calculation import run_calc, check_finished_calcs |
479 |
|
480 |
logger.info('Carrying out procedures for the screening of adsorbate-surface'
|
481 |
' structures.')
|
482 |
if not os.path.isdir("isolated"): |
483 |
err = "'isolated' directory not found. It is needed in order to carry "
|
484 |
"out the screening of structures to be adsorbed"
|
485 |
logger.error(err) |
486 |
raise FileNotFoundError(err)
|
487 |
|
488 |
finished_calcs, unfinished_calcs = check_finished_calcs('isolated',
|
489 |
inp_vars['code'])
|
490 |
if len(unfinished_calcs) == 0: |
491 |
logger.info(f"Found {len(finished_calcs)} structures of isolated "
|
492 |
f"conformers.")
|
493 |
else:
|
494 |
logger.info(f"Found {len(finished_calcs)} structures of isolated "
|
495 |
f"conformers whose calculation finished normally.")
|
496 |
logger.warning(f"Found {len(unfinished_calcs)} calculations more that "
|
497 |
f"did not finish normally: {unfinished_calcs}. \n"
|
498 |
f"Using only the ones that finished normally: "
|
499 |
f"{finished_calcs}.")
|
500 |
|
501 |
conformer_atoms_list = collect_coords(finished_calcs, inp_vars['code'],
|
502 |
'isolated', inp_vars['special_atoms']) |
503 |
selected_confs = select_confs(conformer_atoms_list, finished_calcs, |
504 |
inp_vars['select_magns'],
|
505 |
inp_vars['confs_per_magn'],
|
506 |
inp_vars['code'])
|
507 |
surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms']) |
508 |
surf_ads_list = adsorb_confs(selected_confs, surf, |
509 |
inp_vars['molec_ads_ctrs'], inp_vars['sites'], |
510 |
inp_vars['ads_algo'],
|
511 |
inp_vars['sample_points_per_angle'],
|
512 |
inp_vars['molec_neigh_ctrs'],
|
513 |
inp_vars['surf_norm_vect'],
|
514 |
inp_vars['min_coll_height'],
|
515 |
inp_vars['collision_threshold'],
|
516 |
inp_vars['disso_atoms'])
|
517 |
if len(surf_ads_list) > inp_vars['max_structures']: |
518 |
surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
|
519 |
logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
|
520 |
f'configurations to carry out a calculation of.')
|
521 |
|
522 |
run_calc('screening', inp_vars, surf_ads_list)
|
523 |
logger.info('Finished the procedures for the screening of adsorbate-surface'
|
524 |
' structures section.')
|