dockonsurf / modules / screening.py @ 587dca22
Historique | Voir | Annoter | Télécharger (43,79 ko)
1 |
import logging |
---|---|
2 |
import numpy as np |
3 |
import ase |
4 |
|
5 |
logger = logging.getLogger('DockOnSurf')
|
6 |
|
7 |
|
8 |
def assign_prop(atoms: ase.Atoms, prop_name: str, prop_val): # TODO Needed? |
9 |
atoms.info[prop_name] = prop_val |
10 |
|
11 |
|
12 |
def select_confs(orig_conf_list: list, calc_dirs: list, magns: list, |
13 |
num_sel: int, code: str): |
14 |
"""Takes a list ase.Atoms and selects the most different magnitude-wise.
|
15 |
|
16 |
Given a list of ase.Atoms objects and a list of magnitudes, it selects a
|
17 |
number of the most different conformers according to every magnitude
|
18 |
specified.
|
19 |
|
20 |
@param orig_conf_list: list of ase.Atoms objects to select among.
|
21 |
@param calc_dirs: List of directories where to read the energies from.
|
22 |
@param magns: list of str with the names of the magnitudes to use for the
|
23 |
conformer selection. Supported magnitudes: 'energy', 'moi'.
|
24 |
@param num_sel: number of conformers to select for every of the magnitudes.
|
25 |
@param code: The code that generated the magnitude information.
|
26 |
Supported codes: See formats.py
|
27 |
@return: list of the selected ase.Atoms objects.
|
28 |
"""
|
29 |
from copy import deepcopy |
30 |
from modules.formats import collect_energies |
31 |
|
32 |
conf_list = deepcopy(orig_conf_list) |
33 |
conf_enrgs, mois, selected_ids = [], [], [] |
34 |
if num_sel >= len(conf_list): |
35 |
logger.warning('Number of conformers per magnitude is equal or larger '
|
36 |
'than the total number of conformers. Using all '
|
37 |
f'available conformers: {len(conf_list)}.')
|
38 |
return conf_list
|
39 |
|
40 |
# Read properties
|
41 |
if 'energy' in magns: |
42 |
conf_enrgs = collect_energies(calc_dirs, code, 'isolated')
|
43 |
if 'moi' in magns: |
44 |
mois = np.array([conf.get_moments_of_inertia() for conf in conf_list]) |
45 |
|
46 |
# Assign values
|
47 |
for i, conf in enumerate(conf_list): |
48 |
assign_prop(conf, 'idx', i)
|
49 |
if 'energy' in magns: |
50 |
assign_prop(conf, 'energy', conf_enrgs[i])
|
51 |
if 'moi' in magns: |
52 |
assign_prop(conf, 'moi', mois[i, 2]) |
53 |
|
54 |
# pick ids
|
55 |
for magn in magns: |
56 |
sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn])) |
57 |
if sorted_list[-1].info['idx'] not in selected_ids: |
58 |
selected_ids.append(sorted_list[-1].info['idx']) |
59 |
if num_sel > 1: |
60 |
for i in range(0, len(sorted_list) - 1, |
61 |
len(conf_list) // (num_sel - 1)): |
62 |
if sorted_list[i].info['idx'] not in selected_ids: |
63 |
selected_ids.append(sorted_list[i].info['idx'])
|
64 |
|
65 |
logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
|
66 |
return [conf_list[idx] for idx in selected_ids] |
67 |
|
68 |
|
69 |
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True): |
70 |
"""Computes the angle between two vectors.
|
71 |
|
72 |
@param v1: The first vector.
|
73 |
@param v2: The second vector.
|
74 |
@param ref: Orthogonal vector to both v1 and v2,
|
75 |
along which the sign of the rotation is defined (i.e. positive if
|
76 |
counterclockwise angle when facing ref)
|
77 |
@param degrees: Whether the result should be in radians (True) or in
|
78 |
degrees (False).
|
79 |
@return: The angle in radians if degrees = False, or in degrees if
|
80 |
degrees =True
|
81 |
"""
|
82 |
v1_u = v1 / np.linalg.norm(v1) |
83 |
v2_u = v2 / np.linalg.norm(v2) |
84 |
angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0)) |
85 |
if ref is not None: |
86 |
# Give sign according to ref direction
|
87 |
angle *= (1 if np.dot(np.cross(v1, v2), ref) >= 0 else -1) |
88 |
|
89 |
return angle if not degrees else angle * 180 / np.pi |
90 |
|
91 |
|
92 |
def vect_avg(vects): |
93 |
"""Computes the element-wise mean of a set of vectors.
|
94 |
|
95 |
@param vects: list of lists-like: containing the vectors (num_vectors,
|
96 |
length_vector).
|
97 |
@return: vector average computed doing the element-wise mean.
|
98 |
"""
|
99 |
from utilities import try_command |
100 |
err = "vect_avg parameter vects must be a list-like, able to be converted" \
|
101 |
" np.array"
|
102 |
array = try_command(np.array, [(ValueError, err)], vects)
|
103 |
if len(array.shape) == 1: |
104 |
return array
|
105 |
else:
|
106 |
num_vects = array.shape[1]
|
107 |
return np.array([np.average(array[:, i]) for i in range(num_vects)]) |
108 |
|
109 |
|
110 |
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None): |
111 |
"""Gets the coordinates of the specified indices from a ase.Atoms object.
|
112 |
|
113 |
Given an ase.Atoms object and a list of atom indices specified in ctrs_list
|
114 |
it gets the coordinates of the specified atoms. If the element in the
|
115 |
ctrs_list is not an index but yet a list of indices, it computes the
|
116 |
element-wise mean of the coordinates of the atoms specified in the inner
|
117 |
list.
|
118 |
@param atoms: ase.Atoms object for which to obtain the coordinates of.
|
119 |
@param ctrs_list: list of (indices/list of indices) of the atoms for which
|
120 |
the coordinates should be extracted.
|
121 |
@return: np.ndarray of atomic coordinates.
|
122 |
"""
|
123 |
coords = [] |
124 |
err = "'ctrs_list' argument must be an integer, a list of integers or a " \
|
125 |
"list of lists of integers. Every integer must be in the range " \
|
126 |
"[0, num_atoms)"
|
127 |
if ctrs_list is None: |
128 |
ctrs_list = range(len(atoms)) |
129 |
elif isinstance(ctrs_list, int): |
130 |
if ctrs_list not in range(len(atoms)): |
131 |
logger.error(err) |
132 |
raise ValueError(err) |
133 |
return atoms[ctrs_list].position
|
134 |
for elem in ctrs_list: |
135 |
if isinstance(elem, list): |
136 |
coords.append(vect_avg([atoms[c].position for c in elem])) |
137 |
elif isinstance(elem, int): |
138 |
coords.append(atoms[elem].position) |
139 |
else:
|
140 |
logger.error(err) |
141 |
raise ValueError |
142 |
return np.array(coords)
|
143 |
|
144 |
|
145 |
def compute_norm_vect(atoms, idxs, cell): |
146 |
"""Computes the local normal vector of a surface at a given site.
|
147 |
|
148 |
Given an ase.Atoms object and a site defined as a linear combination of
|
149 |
atoms it computes the vector perpendicular to the surface, considering the
|
150 |
local environment of the site.
|
151 |
@param atoms: ase.Atoms object of the surface.
|
152 |
@param idxs: list or int: Index or list of indices of the atom/s that define
|
153 |
the site
|
154 |
@param cell: Unit cell. A 3x3 matrix (the three unit cell vectors)
|
155 |
@return: numpy.ndarray of the coordinates of the vector locally
|
156 |
perpendicular to the surface.
|
157 |
"""
|
158 |
from ASANN import coordination_numbers as coord_nums |
159 |
if isinstance(idxs, list): |
160 |
atm_vect = [-np.round(coord_nums(atoms.get_scaled_positions(), |
161 |
pbc=np.any(cell), |
162 |
cell_vectors=cell)[3][i], 2) |
163 |
for i in idxs] |
164 |
norm_vec = vect_avg([vect / np.linalg.norm(vect) for vect in atm_vect]) |
165 |
elif isinstance(idxs, int): |
166 |
norm_vec = -coord_nums(atoms.get_scaled_positions(), |
167 |
pbc=np.any(cell), |
168 |
cell_vectors=cell)[3][idxs]
|
169 |
else:
|
170 |
err = "'idxs' must be either an int or a list"
|
171 |
logger.error(err) |
172 |
raise ValueError(err) |
173 |
norm_vec = np.round(norm_vec, 2) / np.linalg.norm(np.round(norm_vec, 2)) |
174 |
logger.info(f"The perpendicular vector to the surface at site '{idxs}' is "
|
175 |
f"{norm_vec}")
|
176 |
return norm_vec
|
177 |
|
178 |
|
179 |
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None, |
180 |
norm_vect=(0, 0, 1)): |
181 |
"""Add an adsorbate to a surface.
|
182 |
|
183 |
This function extends the functionality of ase.build.add_adsorbate
|
184 |
(https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
|
185 |
by enabling to change the z coordinate and the axis perpendicular to the
|
186 |
surface.
|
187 |
@param slab: ase.Atoms object containing the coordinates of the surface
|
188 |
@param adsorbate: ase.Atoms object containing the coordinates of the
|
189 |
adsorbate.
|
190 |
@param site_coord: The coordinates of the adsorption site on the surface.
|
191 |
@param ctr_coord: The coordinates of the adsorption center in the molecule.
|
192 |
@param height: The height above the surface where to adsorb.
|
193 |
@param offset: Offsets the adsorbate by a number of unit cells. Mostly
|
194 |
useful when adding more than one adsorbate.
|
195 |
@param norm_vect: The vector perpendicular to the surface.
|
196 |
"""
|
197 |
from copy import deepcopy |
198 |
info = slab.info.get('adsorbate_info', {})
|
199 |
pos = np.array([0.0, 0.0, 0.0]) # part of absolute coordinates |
200 |
spos = np.array([0.0, 0.0, 0.0]) # part relative to unit cell |
201 |
norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect) |
202 |
if offset is not None: |
203 |
spos += np.asarray(offset, float)
|
204 |
if isinstance(site_coord, str): |
205 |
# A site-name:
|
206 |
if 'sites' not in info: |
207 |
raise TypeError('If the atoms are not made by an ase.build ' |
208 |
'function, position cannot be a name.')
|
209 |
if site_coord not in info['sites']: |
210 |
raise TypeError('Adsorption site %s not supported.' % site_coord) |
211 |
spos += info['sites'][site_coord]
|
212 |
else:
|
213 |
pos += site_coord |
214 |
if 'cell' in info: |
215 |
cell = info['cell']
|
216 |
else:
|
217 |
cell = slab.get_cell() |
218 |
pos += np.dot(spos, cell) |
219 |
# Convert the adsorbate to an Atoms object
|
220 |
if isinstance(adsorbate, ase.Atoms): |
221 |
ads = deepcopy(adsorbate) |
222 |
elif isinstance(adsorbate, ase.Atom): |
223 |
ads = ase.Atoms([adsorbate]) |
224 |
else:
|
225 |
# Assume it is a string representing a single Atom
|
226 |
ads = ase.Atoms([ase.Atom(adsorbate)]) |
227 |
pos += height * norm_vect_u |
228 |
# Move adsorbate into position
|
229 |
ads.translate(pos - ctr_coord) |
230 |
# Attach the adsorbate
|
231 |
slab.extend(ads) |
232 |
|
233 |
|
234 |
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0, |
235 |
nn_molec=0, coll_coeff=1.2): |
236 |
"""Checks whether a slab and a molecule collide or not.
|
237 |
|
238 |
@param slab_molec: The system of adsorbate-slab for which to detect if there
|
239 |
are collisions.
|
240 |
@param nn_slab: Number of neigbors in the surface.
|
241 |
@param nn_molec: Number of neighbors in the molecule.
|
242 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
243 |
atoms resulting in a distance that two atoms being closer to that is
|
244 |
considered as atomic collision.
|
245 |
@param slab_num_atoms: Number of atoms of the bare slab.
|
246 |
@param min_height: The minimum height atoms can have to not be considered as
|
247 |
colliding.
|
248 |
@param vect: The vector perpendicular to the slab.
|
249 |
@return: bool, whether the surface and the molecule collide.
|
250 |
"""
|
251 |
from ase.neighborlist import natural_cutoffs, neighbor_list |
252 |
|
253 |
# Check structure overlap by height
|
254 |
if min_height is not False: |
255 |
cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0], |
256 |
[-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]] |
257 |
if vect.tolist() not in cart_axes: |
258 |
err_msg = "'min_coll_height' option is only implemented for " \
|
259 |
"'surf_norm_vect' to be one of the x, y or z axes. "
|
260 |
logger.error(err_msg) |
261 |
raise ValueError(err_msg) |
262 |
for atom in slab_molec[slab_num_atoms:]: |
263 |
for i, coord in enumerate(vect): |
264 |
if coord == 0: |
265 |
continue
|
266 |
if atom.position[i] * coord < min_height * coord:
|
267 |
return True |
268 |
|
269 |
# Check structure overlap by sphere collision
|
270 |
if coll_coeff is not False: |
271 |
slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff) |
272 |
slab_molec_nghbs = len(
|
273 |
neighbor_list("i", slab_molec, slab_molec_cutoffs))
|
274 |
if slab_molec_nghbs > nn_slab + nn_molec:
|
275 |
return True |
276 |
|
277 |
return False |
278 |
|
279 |
|
280 |
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts, |
281 |
min_coll_height, norm_vect, slab_nghbs, molec_nghbs, |
282 |
coll_coeff, height=2.5):
|
283 |
# TODO Rethink this function
|
284 |
"""Tries to adsorb a molecule on a slab trying to avoid collisions by doing
|
285 |
small rotations.
|
286 |
|
287 |
@param molec: ase.Atoms object of the molecule to adsorb
|
288 |
@param slab: ase.Atoms object of the surface on which to adsorb the
|
289 |
molecule
|
290 |
@param ctr_coord: The coordinates of the molecule to use as adsorption
|
291 |
center.
|
292 |
@param site_coord: The coordinates of the surface on which to adsorb the
|
293 |
molecule
|
294 |
@param num_pts: Number on which to sample Euler angles.
|
295 |
@param min_coll_height: The lowermost height for which to detect a collision
|
296 |
@param norm_vect: The vector perpendicular to the surface.
|
297 |
@param slab_nghbs: Number of neigbors in the surface.
|
298 |
@param molec_nghbs: Number of neighbors in the molecule.
|
299 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
300 |
atoms resulting in a distance that two atoms being closer to that is
|
301 |
considered as atomic collision.
|
302 |
@param height: Height on which to try adsorption
|
303 |
@return collision: bool, whether the structure generated has collisions
|
304 |
between slab and adsorbate.
|
305 |
"""
|
306 |
from copy import deepcopy |
307 |
slab_num_atoms = len(slab)
|
308 |
slab_molec = [] |
309 |
collision = True
|
310 |
max_corr = 6 # Should be an even number |
311 |
d_angle = 180 / ((max_corr / 2.0) * num_pts) |
312 |
num_corr = 0
|
313 |
while collision and num_corr <= max_corr: |
314 |
k = num_corr * (-1) ** num_corr
|
315 |
slab_molec = deepcopy(slab) |
316 |
molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
|
317 |
center=ctr_coord) |
318 |
add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height, |
319 |
norm_vect=norm_vect) |
320 |
collision = check_collision(slab_molec, slab_num_atoms, min_coll_height, |
321 |
norm_vect, slab_nghbs, molec_nghbs, |
322 |
coll_coeff) |
323 |
num_corr += 1
|
324 |
return slab_molec, collision
|
325 |
|
326 |
|
327 |
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor, |
328 |
neigh_cutoff=1):
|
329 |
# TODO rethink
|
330 |
"""Tries to dissociate a H from the molecule and adsorbs it on the slab.
|
331 |
|
332 |
Tries to dissociate a H atom from the molecule and adsorb in on top of the
|
333 |
surface if the distance is shorter than two times the neigh_cutoff value.
|
334 |
@param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
|
335 |
@param h_idx: The index of the hydrogen atom to carry out adsorption of.
|
336 |
@param num_atoms_slab: The number of atoms of the slab without adsorbate.
|
337 |
@param h_acceptor: List of atom types or atom numbers that are H-acceptors.
|
338 |
@param neigh_cutoff: half the maximum distance between the surface and the
|
339 |
H for it to carry out dissociation.
|
340 |
@return: An ase.Atoms object of the system adsorbate-surface with H
|
341 |
"""
|
342 |
from copy import deepcopy |
343 |
from ase.neighborlist import NeighborList |
344 |
slab_molec = deepcopy(slab_molec_orig) |
345 |
cutoffs = len(slab_molec) * [neigh_cutoff]
|
346 |
nl = NeighborList(cutoffs, self_interaction=False, bothways=True, skin=0.0) |
347 |
nl.update(slab_molec) |
348 |
surf_h_vect = np.array([np.infty] * 3)
|
349 |
if h_acceptor == 'all': |
350 |
h_acceptor = list(range(num_atoms_slab)) |
351 |
for neigh_idx in nl.get_neighbors(h_idx)[0]: |
352 |
for elm in h_acceptor: |
353 |
if isinstance(elm, int): |
354 |
if neigh_idx == elm and neigh_idx < num_atoms_slab: |
355 |
dist = np.linalg.norm(slab_molec[neigh_idx].position - |
356 |
slab_molec[h_idx].position) |
357 |
if dist < np.linalg.norm(surf_h_vect):
|
358 |
surf_h_vect = slab_molec[neigh_idx].position \ |
359 |
- slab_molec[h_idx].position |
360 |
else:
|
361 |
if slab_molec[neigh_idx].symbol == elm \
|
362 |
and neigh_idx < num_atoms_slab:
|
363 |
dist = np.linalg.norm(slab_molec[neigh_idx].position - |
364 |
slab_molec[h_idx].position) |
365 |
if dist < np.linalg.norm(surf_h_vect):
|
366 |
surf_h_vect = slab_molec[neigh_idx].position \ |
367 |
- slab_molec[h_idx].position |
368 |
|
369 |
if np.linalg.norm(surf_h_vect) != np.infty:
|
370 |
trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect) |
371 |
slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect |
372 |
return slab_molec
|
373 |
|
374 |
|
375 |
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab): |
376 |
# TODO multiple dissociation
|
377 |
"""Decides which H atoms to dissociate according to a list of atoms.
|
378 |
|
379 |
Given a list of chemical symbols or atom indices it checks for every atom
|
380 |
or any of its neighbor if it's a H and calls dissociate_h to try to carry
|
381 |
out dissociation of that H. For atom indices, it checks both whether
|
382 |
the atom index or its neighbors are H, for chemical symbols, it only checks
|
383 |
if there is a neighbor H.
|
384 |
@param slab_molec: The ase.Atoms object of the system adsorbate-slab.
|
385 |
@param h_donor: List of atom types or atom numbers that are H-donors.
|
386 |
@param h_acceptor: List of atom types or atom numbers that are H-acceptors.
|
387 |
@param num_atoms_slab: Number of atoms of the bare slab.
|
388 |
@return:
|
389 |
"""
|
390 |
from ase.neighborlist import natural_cutoffs, NeighborList |
391 |
molec = slab_molec[num_atoms_slab:] |
392 |
cutoffs = natural_cutoffs(molec) |
393 |
nl = NeighborList(cutoffs, self_interaction=False, bothways=True) |
394 |
nl.update(molec) |
395 |
disso_structs = [] |
396 |
for el in h_donor: |
397 |
if isinstance(el, int): |
398 |
if molec[el].symbol == 'H': |
399 |
disso_struct = dissociate_h(slab_molec, el + num_atoms_slab, |
400 |
num_atoms_slab, h_acceptor) |
401 |
if disso_struct is not None: |
402 |
disso_structs.append(disso_struct) |
403 |
else:
|
404 |
for neigh_idx in nl.get_neighbors(el)[0]: |
405 |
if molec[neigh_idx].symbol == 'H': |
406 |
disso_struct = dissociate_h(slab_molec, neigh_idx + |
407 |
num_atoms_slab, |
408 |
num_atoms_slab, h_acceptor) |
409 |
if disso_struct is not None: |
410 |
disso_structs.append(disso_struct) |
411 |
else:
|
412 |
for atom in molec: |
413 |
if atom.symbol.lower() == el.lower():
|
414 |
for neigh_idx in nl.get_neighbors(atom.index)[0]: |
415 |
if molec[neigh_idx].symbol == 'H': |
416 |
disso_struct = dissociate_h(slab_molec, neigh_idx |
417 |
+ num_atoms_slab, |
418 |
num_atoms_slab, |
419 |
h_acceptor) |
420 |
if disso_struct is not None: |
421 |
disso_structs.append(disso_struct) |
422 |
return disso_structs
|
423 |
|
424 |
|
425 |
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts, |
426 |
min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs, |
427 |
h_donor, h_acceptor): |
428 |
"""Generates adsorbate-surface structures by sampling over Euler angles.
|
429 |
|
430 |
This function generates a number of adsorbate-surface structures at
|
431 |
different orientations of the adsorbate sampled at multiple Euler (zxz)
|
432 |
angles.
|
433 |
@param orig_molec: ase.Atoms object of the molecule to adsorb.
|
434 |
@param slab: ase.Atoms object of the surface on which to adsorb the
|
435 |
molecule.
|
436 |
@param ctr_coord: The coordinates of the molecule to use as adsorption
|
437 |
center.
|
438 |
@param site_coord: The coordinates of the surface on which to adsorb the
|
439 |
molecule.
|
440 |
@param num_pts: Number on which to sample Euler angles.
|
441 |
@param min_coll_height: The lowest height for which to detect a collision.
|
442 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
443 |
atoms resulting in a distance that two atoms being closer to that is
|
444 |
considered as atomic collision.
|
445 |
@param norm_vect: The vector perpendicular to the surface.
|
446 |
@param slab_nghbs: Number of neigbors in the surface.
|
447 |
@param molec_nghbs: Number of neighbors in the molecule.
|
448 |
@param h_donor: List of atom types or atom numbers that are H-donors.
|
449 |
@param h_acceptor: List of atom types or atom numbers that are H-acceptors.
|
450 |
@return: list of ase.Atoms object conatining all the orientations of a given
|
451 |
conformer.
|
452 |
"""
|
453 |
from copy import deepcopy |
454 |
slab_ads_list = [] |
455 |
# rotation around z
|
456 |
for alpha in np.arange(0, 360, 360 / num_pts): |
457 |
# rotation around x'
|
458 |
for beta in np.arange(0, 180, 180 / num_pts): |
459 |
# rotation around z'
|
460 |
for gamma in np.arange(0, 360, 360 / num_pts): |
461 |
molec = deepcopy(orig_molec) |
462 |
molec.euler_rotate(alpha, beta, gamma, center=ctr_coord) |
463 |
slab_molec, collision = correct_coll(molec, slab, |
464 |
ctr_coord, site_coord, |
465 |
num_pts, min_coll_height, |
466 |
norm_vect, |
467 |
slab_nghbs, molec_nghbs, |
468 |
coll_coeff) |
469 |
if not collision and not any([np.allclose(slab_molec.positions, |
470 |
conf.positions) |
471 |
for conf in slab_ads_list]): |
472 |
slab_ads_list.append(slab_molec) |
473 |
if h_donor is not False: |
474 |
slab_ads_list.extend(dissociation(slab_molec, h_donor, |
475 |
h_acceptor, |
476 |
len(slab)))
|
477 |
|
478 |
return slab_ads_list
|
479 |
|
480 |
|
481 |
def chemcat_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf, |
482 |
ctr2_surf, bond_vector, bond_angle_target, |
483 |
dihedral_angle_target=None, mol_dihedral_angle_target=None): |
484 |
"""Performs translation and rotation of an adsorbate defined by an
|
485 |
adsorption bond length, direction, angle and dihedral angle
|
486 |
|
487 |
Carles modification of chemcat's transform_adsorbate to work with
|
488 |
coordinates instead of ase.Atom
|
489 |
Parameters:
|
490 |
molecule (ase.Atoms): The molecule to adsorb.
|
491 |
|
492 |
surf (ase.Atoms): The surface ontop of which to adsorb.
|
493 |
|
494 |
ctr1_mol (int/list): The position of the adsorption center in the
|
495 |
molecule that will be bound to the surface.
|
496 |
|
497 |
ctr2_mol (int/list): The position of a second center of the
|
498 |
adsorbate used to define the adsorption bond angle, and the dihedral
|
499 |
adsorption angle.
|
500 |
|
501 |
ctr3_mol (int/list): The position of a third center in the
|
502 |
adsorbate used to define the adsorbate dihedral angle.
|
503 |
|
504 |
ctr1_surf (int/list): The position of the site on the surface that
|
505 |
will be bound to the molecule.
|
506 |
|
507 |
ctr2_surf (int/list): The position of a second center of the
|
508 |
surface used to define the dihedral adsorption angle.
|
509 |
|
510 |
bond_vector (numpy.ndarray): The adsorption bond desired.
|
511 |
Details: offset = vect(atom1_surf, atom1_mol)
|
512 |
|
513 |
bond_angle_target (float or int): The adsorption bond angle desired (in
|
514 |
degrees).
|
515 |
Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
|
516 |
|
517 |
dihedral_angle_target (float or int): The dihedral adsorption angle
|
518 |
desired (in degrees).
|
519 |
Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
|
520 |
atom1_mol, atom2_mol)
|
521 |
The rotation vector is facing the adsorbate from the surface
|
522 |
(i.e. counterclockwise rotation when facing the surface (i.e.
|
523 |
view from top))
|
524 |
|
525 |
mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
|
526 |
desired (in degrees).
|
527 |
Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
|
528 |
atom2_mol, atom3_mol)
|
529 |
The rotation vector is facing atom2_mol from atom1_mol
|
530 |
|
531 |
Returns:
|
532 |
None (the `molecule` object is modified in-place)
|
533 |
"""
|
534 |
vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf, |
535 |
ctr1_surf) |
536 |
vect_inter = get_atom_coords(molecule, ctr1_mol) \ |
537 |
- get_atom_coords(surf, ctr1_surf) |
538 |
vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule, |
539 |
ctr1_mol) |
540 |
vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule, |
541 |
ctr2_mol) |
542 |
|
543 |
# Check if dihedral angles can be defined
|
544 |
do_dihedral = dihedral_angle_target is not None |
545 |
do_mol_dihedral = mol_dihedral_angle_target is not None |
546 |
dihedral_use_mol2 = False
|
547 |
# Check if vect_surf and vect_inter are not aligned
|
548 |
if np.allclose(np.cross(vect_surf, vect_inter), 0): |
549 |
logger.warning( |
550 |
"Surface atoms are incompatible with adsorption "
|
551 |
"direction/bond. An adsorption dihedral angle cannot be defined.")
|
552 |
do_dihedral = False
|
553 |
# Check if requested bond angle is not flat
|
554 |
if np.isclose((bond_angle_target + 90) % 180 - 90, 0): |
555 |
logger.warning("Requested bond angle is flat. Only a single dihedral "
|
556 |
"angle can be defined (ctr2_surf, ctr1_surf, ctr2_mol, "
|
557 |
"ctr3_mol).")
|
558 |
do_mol_dihedral = False
|
559 |
dihedral_use_mol2 = True
|
560 |
logger.warning("Dihedral adsorption angle rotation will be perfomed "
|
561 |
"with (ctr2_surf, ctr1_surf, ctr2_mol, ctr3_mol).")
|
562 |
# Check if vect_mol and vect2_mol are not aligned
|
563 |
if np.allclose(np.cross(vect_mol, vect2_mol), 0): |
564 |
logger.warning("Adsorbates atoms are aligned. An adsorbate dihedral "
|
565 |
"angle cannot be defined.")
|
566 |
do_mol_dihedral = False
|
567 |
if not do_dihedral: |
568 |
logger.warning("Dihedral adsorption angle rotation will not be "
|
569 |
"performed.")
|
570 |
if not do_mol_dihedral: |
571 |
logger.warning("Adsorbate dihedral angle rotation will not be "
|
572 |
"performed.")
|
573 |
|
574 |
###########################
|
575 |
# Translation #
|
576 |
###########################
|
577 |
|
578 |
# Compute and apply translation of adsorbate
|
579 |
translation = bond_vector - vect_inter |
580 |
molecule.translate(translation) |
581 |
|
582 |
# Update adsorption bond
|
583 |
vect_inter = get_atom_coords(molecule, ctr1_mol) - \ |
584 |
get_atom_coords(surf, ctr1_surf) |
585 |
|
586 |
# Check if translation was successful
|
587 |
if np.allclose(vect_inter, bond_vector):
|
588 |
pass # print("Translation successfully applied (error: ~ {:.5g} unit " |
589 |
# "length)".format(np.linalg.norm(vect_inter - bond_vector)))
|
590 |
else:
|
591 |
err = 'An unknown error occured during the translation'
|
592 |
logger.error(err) |
593 |
raise AssertionError(err) |
594 |
|
595 |
###########################
|
596 |
# Bond angle rotation #
|
597 |
###########################
|
598 |
|
599 |
# Compute rotation vector
|
600 |
rotation_vector = np.cross(-vect_inter, vect_mol) |
601 |
if np.allclose(rotation_vector, 0, atol=1e-5): |
602 |
# If molecular bonds are aligned, any vector orthogonal to vect_inter
|
603 |
# can be used Such vector can be found as the orthogonal rejection of
|
604 |
# either X-axis, Y-axis or Z-axis onto vect_inter (since they cannot
|
605 |
# be all aligned)
|
606 |
non_aligned_vector = np.zeros(3)
|
607 |
# Select the most orthogonal axis (lowest dot product):
|
608 |
non_aligned_vector[np.argmin(np.abs(vect_inter))] = 1
|
609 |
rotation_vector = non_aligned_vector - np.dot(non_aligned_vector, |
610 |
vect_inter) / np.dot( |
611 |
vect_inter, vect_inter) * vect_inter |
612 |
|
613 |
# Retrieve current bond angle
|
614 |
bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector) |
615 |
|
616 |
# Apply rotation to reach desired bond_angle
|
617 |
molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector, |
618 |
center=get_atom_coords(molecule, ctr1_mol)) |
619 |
|
620 |
# Update molecular bonds
|
621 |
vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule, |
622 |
ctr1_mol) |
623 |
vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule, |
624 |
ctr2_mol) |
625 |
|
626 |
# Check if rotation was successful
|
627 |
bond_angle = get_vect_angle(-vect_inter, vect_mol) |
628 |
if np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90, 0, |
629 |
atol=1e-3) and np.allclose(get_atom_coords(molecule, ctr1_mol) |
630 |
- get_atom_coords(surf, |
631 |
ctr1_surf), |
632 |
vect_inter): |
633 |
pass # print("Rotation successfully applied (error: {:.5f}°)".format( |
634 |
# (bond_angle - bond_angle_target + 90) % 180 - 90))
|
635 |
else:
|
636 |
err = 'An unknown error occured during the rotation'
|
637 |
logger.error(err) |
638 |
raise AssertionError(err) |
639 |
|
640 |
###########################
|
641 |
# Dihedral angle rotation #
|
642 |
###########################
|
643 |
|
644 |
# Perform dihedral rotation if possible
|
645 |
if do_dihedral:
|
646 |
# Retrieve current dihedral angle (by computing the angle between the
|
647 |
# vector rejection of vect_surf and vect_mol onto vect_inter)
|
648 |
vect_inter_inner = np.dot(vect_inter, vect_inter) |
649 |
vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \ |
650 |
vect_inter_inner * vect_inter |
651 |
if dihedral_use_mol2:
|
652 |
vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \ |
653 |
vect_inter_inner * vect_inter |
654 |
else:
|
655 |
vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \ |
656 |
vect_inter_inner * vect_inter |
657 |
dihedral_angle_ini = get_vect_angle(vect_surf_reject, vect_mol_reject, |
658 |
vect_inter) |
659 |
|
660 |
# Apply dihedral rotation along vect_inter
|
661 |
molecule.rotate(dihedral_angle_target - dihedral_angle_ini, |
662 |
v=vect_inter, center=get_atom_coords(molecule, |
663 |
ctr1_mol)) |
664 |
|
665 |
# Update molecular bonds
|
666 |
vect_mol = get_atom_coords(molecule, ctr2_mol) - \ |
667 |
get_atom_coords(molecule, ctr1_mol) |
668 |
vect2_mol = get_atom_coords(molecule, ctr3_mol) - \ |
669 |
get_atom_coords(molecule, ctr2_mol) |
670 |
|
671 |
# Check if rotation was successful
|
672 |
# Check dihedral rotation
|
673 |
if dihedral_use_mol2:
|
674 |
vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \ |
675 |
vect_inter_inner * vect_inter |
676 |
else:
|
677 |
vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \ |
678 |
vect_inter_inner * vect_inter |
679 |
dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject, |
680 |
vect_inter) |
681 |
# Check bond rotation is unmodified
|
682 |
bond_angle = get_vect_angle(-vect_inter, vect_mol) |
683 |
if np.isclose((dihedral_angle - dihedral_angle_target + 90) % 180 - 90, |
684 |
0, atol=1e-3) \ |
685 |
and np.isclose((bond_angle - bond_angle_target + 90) % |
686 |
180 - 90, 0, atol=1e-5) \ |
687 |
and np.allclose(get_atom_coords(molecule, ctr1_mol)
|
688 |
- get_atom_coords(surf, ctr1_surf), |
689 |
vect_inter): |
690 |
pass # print( "Dihedral rotation successfully applied (error: { |
691 |
# :.5f}°)".format((dihedral_angle - dihedral_angle_target + 90) %
|
692 |
# 180 - 90))
|
693 |
else:
|
694 |
err = 'An unknown error occured during the dihedral rotation'
|
695 |
logger.error(err) |
696 |
raise AssertionError(err) |
697 |
|
698 |
#####################################
|
699 |
# Adsorbate dihedral angle rotation #
|
700 |
#####################################
|
701 |
|
702 |
# Perform adsorbate dihedral rotation if possible
|
703 |
if do_mol_dihedral:
|
704 |
# Retrieve current adsorbate dihedral angle (by computing the angle
|
705 |
# between the orthogonal rejection of vect_inter and vect2_mol onto
|
706 |
# vect_mol)
|
707 |
vect_mol_inner = np.dot(vect_mol, vect_mol) |
708 |
bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \ |
709 |
vect_mol_inner * vect_mol |
710 |
bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \ |
711 |
vect_mol_inner * vect_mol |
712 |
dihedral_angle_ini = get_vect_angle(bond_inter_reject, |
713 |
bond2_mol_reject, vect_mol) |
714 |
|
715 |
# Apply dihedral rotation along vect_mol
|
716 |
molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini, |
717 |
v=vect_mol, center=get_atom_coords(molecule, ctr1_mol)) |
718 |
|
719 |
# Update molecular bonds
|
720 |
vect_mol = get_atom_coords(molecule, ctr2_mol) \ |
721 |
- get_atom_coords(molecule, ctr1_mol) |
722 |
vect2_mol = get_atom_coords(molecule, ctr3_mol) \ |
723 |
- get_atom_coords(molecule, ctr2_mol) |
724 |
|
725 |
# Check if rotation was successful
|
726 |
# Check adsorbate dihedral rotation
|
727 |
vect_mol_inner = np.dot(vect_mol, vect_mol) |
728 |
bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \ |
729 |
vect_mol_inner * vect_mol |
730 |
mol_dihedral_angle = get_vect_angle(bond_inter_reject, |
731 |
bond2_mol_reject, vect_mol) |
732 |
# Check dihedral rotation
|
733 |
vect_inter_inner = np.dot(vect_inter, vect_inter) |
734 |
vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \ |
735 |
vect_inter_inner * vect_inter |
736 |
vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \ |
737 |
vect_inter_inner * vect_inter |
738 |
dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject, |
739 |
vect_inter) |
740 |
# Check bond rotation is unmodified
|
741 |
bond_angle = get_vect_angle(-vect_inter, vect_mol) |
742 |
if np.isclose((mol_dihedral_angle - mol_dihedral_angle_target + 90) % |
743 |
180 - 90, 0, atol=1e-3) \ |
744 |
and np.isclose((dihedral_angle -
|
745 |
dihedral_angle_target + 90) % 180 - 90, 0, |
746 |
atol=1e-5) \
|
747 |
and np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90, |
748 |
0, atol=1e-5) \ |
749 |
and np.allclose(get_atom_coords(molecule, ctr1_mol) -
|
750 |
get_atom_coords(surf, ctr1_surf), |
751 |
vect_inter): |
752 |
pass # print( |
753 |
# "Adsorbate dihedral rotation successfully applied (error:
|
754 |
# {:.5f}°)".format((mol_dihedral_angle - mol_dihedral_angle_target
|
755 |
# + 90) % 180 - 90))
|
756 |
else:
|
757 |
err = 'An unknown error occured during the adsorbate dihedral ' \
|
758 |
'rotation'
|
759 |
logger.error(err) |
760 |
raise AssertionError(err) |
761 |
|
762 |
|
763 |
def ads_chemcat(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf, |
764 |
ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect, |
765 |
slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel): |
766 |
"""Generates adsorbate-surface structures by sampling over chemcat angles.
|
767 |
|
768 |
@param orig_molec: ase.Atoms object of the molecule to adsorb (adsorbate).
|
769 |
@param slab: ase.Atoms object of the surface on which to adsorb the molecule
|
770 |
@param ctr1_mol: The index/es of the center in the adsorbate to use as
|
771 |
adsorption center.
|
772 |
@param ctr2_mol: The index/es of the center in the adsorbate to use for the
|
773 |
definition of the surf-adsorbate angle, surf-adsorbate dihedral angle
|
774 |
and adsorbate dihedral angle.
|
775 |
@param ctr3_mol: The index/es of the center in the adsorbate to use for the
|
776 |
definition of the adsorbate dihedral angle.
|
777 |
@param ctr1_surf: The index/es of the center in the surface to use as
|
778 |
adsorption center.
|
779 |
@param ctr2_surf: The index/es of the center in the surface to use for the
|
780 |
definition of the surf-adsorbate dihedral angle.
|
781 |
@param num_pts: Number on which to sample Euler angles.
|
782 |
@param min_coll_height: The lowest height for which to detect a collision
|
783 |
@param coll_coeff: The coefficient that multiplies the covalent radius of
|
784 |
atoms resulting in a distance that two atoms being closer to that is
|
785 |
considered as atomic collision.
|
786 |
@param norm_vect: The vector perpendicular to the surface.
|
787 |
@param slab_nghbs: Number of neigbors in the surface.
|
788 |
@param molec_nghbs: Number of neighbors in the molecule.
|
789 |
@param h_donor: List of atom types or atom numbers that are H-donors.
|
790 |
@param h_acceptor: List of atom types or atom numbers that are H-acceptors.
|
791 |
@param max_hel: Maximum value for sampling the helicopter
|
792 |
(surf-adsorbate dihedral) angle.
|
793 |
@return: list of ase.Atoms object conatining all the orientations of a given
|
794 |
conformer.
|
795 |
"""
|
796 |
from copy import deepcopy |
797 |
slab_ads_list = [] |
798 |
# Rotation over bond angle # TODO Check sampling
|
799 |
for alpha in np.arange(90, 180, 90 / min(num_pts, 2)): |
800 |
# Rotation over surf-adsorbate dihedral angle
|
801 |
for beta in np.arange(0, max_hel, max_hel / num_pts): |
802 |
# Rotation over adsorbate bond dihedral angle
|
803 |
for gamma in np.arange(0, 360, 360 / num_pts): |
804 |
new_molec = deepcopy(orig_molec) |
805 |
chemcat_rotate(new_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, |
806 |
ctr1_surf, ctr2_surf, norm_vect, alpha, |
807 |
beta, gamma) |
808 |
site_coords = get_atom_coords(slab, ctr1_surf) |
809 |
ctr_coords = get_atom_coords(new_molec, ctr1_mol) |
810 |
slab_molec, collision = correct_coll(new_molec, slab, |
811 |
ctr_coords, site_coords, |
812 |
num_pts, min_coll_height, |
813 |
norm_vect, slab_nghbs, |
814 |
molec_nghbs, coll_coeff) |
815 |
if not collision and \ |
816 |
not any([np.allclose(slab_molec.positions, |
817 |
conf.positions) |
818 |
for conf in slab_ads_list]): |
819 |
slab_ads_list.append(slab_molec) |
820 |
if h_donor is not False: |
821 |
slab_ads_list.extend(dissociation(slab_molec, h_donor, |
822 |
h_acceptor, |
823 |
len(slab)))
|
824 |
|
825 |
return slab_ads_list
|
826 |
|
827 |
|
828 |
def adsorb_confs(conf_list, surf, inp_vars): |
829 |
"""Generates a number of adsorbate-surface structure coordinates.
|
830 |
|
831 |
Given a list of conformers, a surface, a list of atom indices (or list of
|
832 |
list of indices) of both the surface and the adsorbate, it generates a
|
833 |
number of adsorbate-surface structures for every possible combination of
|
834 |
them at different orientations.
|
835 |
@param conf_list: list of ase.Atoms of the different conformers
|
836 |
@param surf: the ase.Atoms object of the surface
|
837 |
@param inp_vars: Calculation parameters from input file.
|
838 |
@return: list of ase.Atoms for the adsorbate-surface structures
|
839 |
"""
|
840 |
from ase.neighborlist import natural_cutoffs, neighbor_list |
841 |
molec_ctrs = inp_vars['molec_ctrs']
|
842 |
sites = inp_vars['sites']
|
843 |
angles = inp_vars['set_angles']
|
844 |
num_pts = inp_vars['sample_points_per_angle']
|
845 |
norm_vect = inp_vars['surf_norm_vect']
|
846 |
min_coll_height = inp_vars['min_coll_height']
|
847 |
coll_coeff = inp_vars['collision_threshold']
|
848 |
h_donor = inp_vars['h_donor']
|
849 |
h_acceptor = inp_vars['h_acceptor']
|
850 |
|
851 |
if inp_vars['pbc_cell'] is not False: |
852 |
surf.set_pbc(True)
|
853 |
surf.set_cell(inp_vars['pbc_cell'])
|
854 |
|
855 |
surf_ads_list = [] |
856 |
sites_coords = get_atom_coords(surf, sites) |
857 |
if coll_coeff is not False: |
858 |
surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff) |
859 |
surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs)) |
860 |
else:
|
861 |
surf_nghbs = 0
|
862 |
for i, conf in enumerate(conf_list): |
863 |
molec_ctr_coords = get_atom_coords(conf, molec_ctrs) |
864 |
if inp_vars['pbc_cell'] is not False: |
865 |
conf.set_pbc(True)
|
866 |
conf.set_cell(inp_vars['pbc_cell'])
|
867 |
if coll_coeff is not False: |
868 |
conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff) |
869 |
molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs)) |
870 |
else:
|
871 |
molec_nghbs = 0
|
872 |
for s, site in enumerate(sites_coords): |
873 |
if norm_vect == 'auto': |
874 |
norm_vect = compute_norm_vect(surf, sites[s], |
875 |
inp_vars['pbc_cell'])
|
876 |
for c, ctr in enumerate(molec_ctr_coords): |
877 |
if angles == 'euler': |
878 |
surf_ads_list.extend(ads_euler(conf, surf, ctr, site, |
879 |
num_pts, min_coll_height, |
880 |
coll_coeff, norm_vect, |
881 |
surf_nghbs, molec_nghbs, |
882 |
h_donor, h_acceptor)) |
883 |
elif angles == 'chemcat': |
884 |
mol_ctr1 = molec_ctrs[c] |
885 |
mol_ctr2 = inp_vars["molec_ctrs2"][c]
|
886 |
mol_ctr3 = inp_vars["molec_ctrs3"][c]
|
887 |
surf_ctr1 = sites[s] |
888 |
surf_ctr2 = inp_vars["surf_ctrs2"][s]
|
889 |
max_h = inp_vars["max_helic_angle"]
|
890 |
surf_ads_list.extend(ads_chemcat(conf, surf, mol_ctr1, |
891 |
mol_ctr2, mol_ctr3, |
892 |
surf_ctr1, surf_ctr2, |
893 |
num_pts, min_coll_height, |
894 |
coll_coeff, norm_vect, |
895 |
surf_nghbs, molec_nghbs, |
896 |
h_donor, h_acceptor, |
897 |
max_h)) |
898 |
return surf_ads_list
|
899 |
|
900 |
|
901 |
def run_screening(inp_vars): |
902 |
"""Carries out the screening of adsorbate structures on a surface.
|
903 |
|
904 |
@param inp_vars: Calculation parameters from input file.
|
905 |
"""
|
906 |
import os |
907 |
import random |
908 |
from modules.formats import collect_coords, adapt_format |
909 |
from modules.calculation import run_calc, check_finished_calcs |
910 |
|
911 |
logger.info('Carrying out procedures for the screening of adsorbate-surface'
|
912 |
' structures.')
|
913 |
if not os.path.isdir("isolated"): |
914 |
err = "'isolated' directory not found. It is needed in order to carry "
|
915 |
"out the screening of structures to be adsorbed"
|
916 |
logger.error(err) |
917 |
raise FileNotFoundError(err)
|
918 |
|
919 |
finished_calcs, unfinished_calcs = check_finished_calcs('isolated',
|
920 |
inp_vars['code'])
|
921 |
logger.info(f"Found {len(finished_calcs)} structures of isolated "
|
922 |
f"conformers whose calculation finished normally.")
|
923 |
if len(unfinished_calcs) != 0: |
924 |
logger.warning(f"Found {len(unfinished_calcs)} calculations more that "
|
925 |
f"did not finish normally: {unfinished_calcs}. \n"
|
926 |
f"Using only the ones that finished normally: "
|
927 |
f"{finished_calcs}.")
|
928 |
|
929 |
conformer_atoms_list = collect_coords(finished_calcs, inp_vars['code'],
|
930 |
'isolated', inp_vars['special_atoms']) |
931 |
selected_confs = select_confs(conformer_atoms_list, finished_calcs, |
932 |
inp_vars['select_magns'],
|
933 |
inp_vars['confs_per_magn'],
|
934 |
inp_vars['code'])
|
935 |
surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms']) |
936 |
surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars) |
937 |
if len(surf_ads_list) > inp_vars['max_structures']: |
938 |
surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
|
939 |
logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
|
940 |
f'configurations to carry out a calculation of.')
|
941 |
|
942 |
run_calc('screening', inp_vars, surf_ads_list)
|
943 |
logger.info('Finished the procedures for the screening of adsorbate-surface'
|
944 |
' structures section.')
|