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

dockonsurf / modules / screening.py @ e8bebcca

Historique | Voir | Annoter | Télécharger (42,23 ko)

1 e07c09eb Carles
import logging
2 f3d1e601 Carles Marti
import numpy as np
3 f3d1e601 Carles Marti
import ase
4 e07c09eb Carles
5 e07c09eb Carles
logger = logging.getLogger('DockOnSurf')
6 e07c09eb Carles
7 e07c09eb Carles
8 bfe93f0d Carles Marti
def assign_prop(atoms: ase.Atoms, prop_name: str, prop_val):  # TODO Needed?
9 bfe93f0d Carles Marti
    atoms.info[prop_name] = prop_val
10 bfe93f0d Carles Marti
11 bfe93f0d Carles Marti
12 e1c5f171 Carles Marti
def select_confs(orig_conf_list: list, calc_dirs: list, magns: list,
13 e1c5f171 Carles Marti
                 num_sel: int, code: str):
14 bfe93f0d Carles Marti
    """Takes a list ase.Atoms and selects the most different magnitude-wise.
15 bfe93f0d Carles Marti

16 bfe93f0d Carles Marti
    Given a list of ase.Atoms objects and a list of magnitudes, it selects a
17 bfe93f0d Carles Marti
    number of the most different conformers according to every magnitude
18 bfe93f0d Carles Marti
    specified.
19 1e9e784d Carles Marti
     
20 bfe93f0d Carles Marti
    @param orig_conf_list: list of ase.Atoms objects to select among.
21 1e9e784d Carles Marti
    @param calc_dirs: List of directories where to read the energies from.
22 bfe93f0d Carles Marti
    @param magns: list of str with the names of the magnitudes to use for the
23 1e9e784d Carles Marti
        conformer selection. Supported magnitudes: 'energy', 'moi'.
24 bfe93f0d Carles Marti
    @param num_sel: number of conformers to select for every of the magnitudes.
25 bfe93f0d Carles Marti
    @param code: The code that generated the magnitude information.
26 bfe93f0d Carles Marti
         Supported codes: See formats.py
27 bfe93f0d Carles Marti
    @return: list of the selected ase.Atoms objects.
28 bfe93f0d Carles Marti
    """
29 bfe93f0d Carles Marti
    from copy import deepcopy
30 fd2384fc Carles Marti
    from modules.formats import collect_energies
31 bfe93f0d Carles Marti
32 bfe93f0d Carles Marti
    conf_list = deepcopy(orig_conf_list)
33 e1c5f171 Carles Marti
    conf_enrgs, mois, selected_ids = [], [], []
34 bfe93f0d Carles Marti
    if num_sel >= len(conf_list):
35 bfe93f0d Carles Marti
        logger.warning('Number of conformers per magnitude is equal or larger '
36 bfe93f0d Carles Marti
                       'than the total number of conformers. Using all '
37 695dcff8 Carles Marti
                       f'available conformers: {len(conf_list)}.')
38 bfe93f0d Carles Marti
        return conf_list
39 bfe93f0d Carles Marti
40 bfe93f0d Carles Marti
    # Read properties
41 bfe93f0d Carles Marti
    if 'energy' in magns:
42 1e9e784d Carles Marti
        conf_enrgs = collect_energies(calc_dirs, code, 'isolated')
43 bfe93f0d Carles Marti
    if 'moi' in magns:
44 bfe93f0d Carles Marti
        mois = np.array([conf.get_moments_of_inertia() for conf in conf_list])
45 bfe93f0d Carles Marti
46 bfe93f0d Carles Marti
    # Assign values
47 bfe93f0d Carles Marti
    for i, conf in enumerate(conf_list):
48 bfe93f0d Carles Marti
        assign_prop(conf, 'idx', i)
49 bfe93f0d Carles Marti
        if 'energy' in magns:
50 bfe93f0d Carles Marti
            assign_prop(conf, 'energy', conf_enrgs[i])
51 bfe93f0d Carles Marti
        if 'moi' in magns:
52 bfe93f0d Carles Marti
            assign_prop(conf, 'moi', mois[i, 2])
53 bfe93f0d Carles Marti
54 bfe93f0d Carles Marti
    # pick ids
55 bfe93f0d Carles Marti
    for magn in magns:
56 bfe93f0d Carles Marti
        sorted_list = sorted(conf_list, key=lambda conf: abs(conf.info[magn]))
57 bfe93f0d Carles Marti
        if sorted_list[-1].info['idx'] not in selected_ids:
58 bfe93f0d Carles Marti
            selected_ids.append(sorted_list[-1].info['idx'])
59 bfe93f0d Carles Marti
        if num_sel > 1:
60 bfe93f0d Carles Marti
            for i in range(0, len(sorted_list) - 1,
61 bfe93f0d Carles Marti
                           len(conf_list) // (num_sel - 1)):
62 bfe93f0d Carles Marti
                if sorted_list[i].info['idx'] not in selected_ids:
63 bfe93f0d Carles Marti
                    selected_ids.append(sorted_list[i].info['idx'])
64 bfe93f0d Carles Marti
65 695dcff8 Carles Marti
    logger.info(f'Selected {len(selected_ids)} conformers for adsorption.')
66 bfe93f0d Carles Marti
    return [conf_list[idx] for idx in selected_ids]
67 bfe93f0d Carles Marti
68 bfe93f0d Carles Marti
69 7dd94df7 Carles Marti
def get_vect_angle(v1: list, v2: list, ref=None, degrees=True):
70 b4aef3d7 Carles Marti
    """Computes the angle between two vectors.
71 b4aef3d7 Carles Marti

72 b4aef3d7 Carles Marti
    @param v1: The first vector.
73 b4aef3d7 Carles Marti
    @param v2: The second vector.
74 b516a1d6 Carles Marti
    @param ref: Orthogonal vector to both v1 and v2,
75 b516a1d6 Carles Marti
        along which the sign of the rotation is defined (i.e. positive if
76 b516a1d6 Carles Marti
        counterclockwise angle when facing ref)
77 b4aef3d7 Carles Marti
    @param degrees: Whether the result should be in radians (True) or in
78 b4aef3d7 Carles Marti
        degrees (False).
79 b4aef3d7 Carles Marti
    @return: The angle in radians if degrees = False, or in degrees if
80 b4aef3d7 Carles Marti
        degrees =True
81 b4aef3d7 Carles Marti
    """
82 b4aef3d7 Carles Marti
    v1_u = v1 / np.linalg.norm(v1)
83 b4aef3d7 Carles Marti
    v2_u = v2 / np.linalg.norm(v2)
84 b4aef3d7 Carles Marti
    angle = np.arccos(np.clip(np.dot(v1_u, v2_u), -1.0, 1.0))
85 b516a1d6 Carles Marti
    if ref is not None:
86 b516a1d6 Carles Marti
        # Give sign according to ref direction
87 b516a1d6 Carles Marti
        angle *= (1 if np.dot(np.cross(v1, v2), ref) >= 0 else -1)
88 b516a1d6 Carles Marti
89 b4aef3d7 Carles Marti
    return angle if not degrees else angle * 180 / np.pi
90 b4aef3d7 Carles Marti
91 b4aef3d7 Carles Marti
92 12a12bdd Carles Marti
def vect_avg(vects):
93 a5cc42ff Carles Marti
    """Computes the element-wise mean of a set of vectors.
94 12a12bdd Carles Marti

95 dadc6016 Carles Marti
    @param vects: list of lists-like: containing the vectors (num_vectors,
96 dadc6016 Carles Marti
        length_vector).
97 a5cc42ff Carles Marti
    @return: vector average computed doing the element-wise mean.
98 12a12bdd Carles Marti
    """
99 12a12bdd Carles Marti
    from utilities import try_command
100 12a12bdd Carles Marti
    err = "vect_avg parameter vects must be a list-like, able to be converted" \
101 12a12bdd Carles Marti
          " np.array"
102 dadc6016 Carles Marti
    array = try_command(np.array, [(ValueError, err)], vects)
103 dadc6016 Carles Marti
    if len(array.shape) == 1:
104 dadc6016 Carles Marti
        return array
105 12a12bdd Carles Marti
    else:
106 dadc6016 Carles Marti
        num_vects = array.shape[1]
107 dadc6016 Carles Marti
        return np.array([np.average(array[:, i]) for i in range(num_vects)])
108 12a12bdd Carles Marti
109 12a12bdd Carles Marti
110 8fdff302 Carles Marti
def get_atom_coords(atoms: ase.Atoms, ctrs_list=None):
111 8fdff302 Carles Marti
    """Gets the coordinates of the specified indices from a ase.Atoms object.
112 a5cc42ff Carles Marti

113 a5cc42ff Carles Marti
    Given an ase.Atoms object and a list of atom indices specified in ctrs_list
114 a5cc42ff Carles Marti
    it gets the coordinates of the specified atoms. If the element in the
115 a5cc42ff Carles Marti
    ctrs_list is not an index but yet a list of indices, it computes the
116 a5cc42ff Carles Marti
    element-wise mean of the coordinates of the atoms specified in the inner
117 a5cc42ff Carles Marti
    list.
118 a5cc42ff Carles Marti
    @param atoms: ase.Atoms object for which to obtain the coordinates of.
119 a5cc42ff Carles Marti
    @param ctrs_list: list of (indices/list of indices) of the atoms for which
120 a5cc42ff Carles Marti
                      the coordinates should be extracted.
121 a5cc42ff Carles Marti
    @return: np.ndarray of atomic coordinates.
122 a5cc42ff Carles Marti
    """
123 12a12bdd Carles Marti
    coords = []
124 8fdff302 Carles Marti
    err = "'ctrs_list' argument must be an integer, a list of integers or a " \
125 8fdff302 Carles Marti
          "list of lists of integers. Every integer must be in the range " \
126 8fdff302 Carles Marti
          "[0, num_atoms)"
127 8fdff302 Carles Marti
    if ctrs_list is None:
128 8fdff302 Carles Marti
        ctrs_list = range(len(atoms))
129 8fdff302 Carles Marti
    elif isinstance(ctrs_list, int):
130 8fdff302 Carles Marti
        if ctrs_list not in range(len(atoms)):
131 8fdff302 Carles Marti
            logger.error(err)
132 8fdff302 Carles Marti
            raise ValueError(err)
133 8fdff302 Carles Marti
        return atoms[ctrs_list].position
134 8fdff302 Carles Marti
    for elem in ctrs_list:
135 12a12bdd Carles Marti
        if isinstance(elem, list):
136 8fdff302 Carles Marti
            coords.append(vect_avg([atoms[c].position for c in elem]))
137 12a12bdd Carles Marti
        elif isinstance(elem, int):
138 12a12bdd Carles Marti
            coords.append(atoms[elem].position)
139 12a12bdd Carles Marti
        else:
140 12a12bdd Carles Marti
            logger.error(err)
141 12a12bdd Carles Marti
            raise ValueError
142 12a12bdd Carles Marti
    return np.array(coords)
143 f3d1e601 Carles Marti
144 f3d1e601 Carles Marti
145 dadc6016 Carles Marti
def add_adsorbate(slab, adsorbate, site_coord, ctr_coord, height, offset=None,
146 1d22a086 Carles Marti
                  norm_vect=(0, 0, 1)):
147 1d22a086 Carles Marti
    """Add an adsorbate to a surface.
148 1d22a086 Carles Marti

149 1d22a086 Carles Marti
    This function extends the functionality of ase.build.add_adsorbate
150 1d22a086 Carles Marti
    (https://wiki.fysik.dtu.dk/ase/ase/build/surface.html#ase.build.add_adsorbate)
151 1d22a086 Carles Marti
    by enabling to change the z coordinate and the axis perpendicular to the
152 1d22a086 Carles Marti
    surface.
153 1d22a086 Carles Marti
    @param slab: ase.Atoms object containing the coordinates of the surface
154 1d22a086 Carles Marti
    @param adsorbate: ase.Atoms object containing the coordinates of the
155 1d22a086 Carles Marti
        adsorbate.
156 dadc6016 Carles Marti
    @param site_coord: The coordinates of the adsorption site on the surface.
157 dadc6016 Carles Marti
    @param ctr_coord: The coordinates of the adsorption center in the molecule.
158 dadc6016 Carles Marti
    @param height: The height above the surface where to adsorb.
159 1d22a086 Carles Marti
    @param offset: Offsets the adsorbate by a number of unit cells. Mostly
160 1d22a086 Carles Marti
        useful when adding more than one adsorbate.
161 1d22a086 Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
162 1d22a086 Carles Marti
    """
163 36d92f4f Carles Marti
    from copy import deepcopy
164 1d22a086 Carles Marti
    info = slab.info.get('adsorbate_info', {})
165 dadc6016 Carles Marti
    pos = np.array([0.0, 0.0, 0.0])  # part of absolute coordinates
166 1d22a086 Carles Marti
    spos = np.array([0.0, 0.0, 0.0])  # part relative to unit cell
167 dadc6016 Carles Marti
    norm_vect_u = np.array(norm_vect) / np.linalg.norm(norm_vect)
168 1d22a086 Carles Marti
    if offset is not None:
169 1d22a086 Carles Marti
        spos += np.asarray(offset, float)
170 dadc6016 Carles Marti
    if isinstance(site_coord, str):
171 1d22a086 Carles Marti
        # A site-name:
172 1d22a086 Carles Marti
        if 'sites' not in info:
173 1d22a086 Carles Marti
            raise TypeError('If the atoms are not made by an ase.build '
174 1d22a086 Carles Marti
                            'function, position cannot be a name.')
175 dadc6016 Carles Marti
        if site_coord not in info['sites']:
176 dadc6016 Carles Marti
            raise TypeError('Adsorption site %s not supported.' % site_coord)
177 dadc6016 Carles Marti
        spos += info['sites'][site_coord]
178 1d22a086 Carles Marti
    else:
179 dadc6016 Carles Marti
        pos += site_coord
180 1d22a086 Carles Marti
    if 'cell' in info:
181 1d22a086 Carles Marti
        cell = info['cell']
182 1d22a086 Carles Marti
    else:
183 1d22a086 Carles Marti
        cell = slab.get_cell()
184 1d22a086 Carles Marti
    pos += np.dot(spos, cell)
185 1d22a086 Carles Marti
    # Convert the adsorbate to an Atoms object
186 1d22a086 Carles Marti
    if isinstance(adsorbate, ase.Atoms):
187 36d92f4f Carles Marti
        ads = deepcopy(adsorbate)
188 1d22a086 Carles Marti
    elif isinstance(adsorbate, ase.Atom):
189 1d22a086 Carles Marti
        ads = ase.Atoms([adsorbate])
190 1d22a086 Carles Marti
    else:
191 1d22a086 Carles Marti
        # Assume it is a string representing a single Atom
192 1d22a086 Carles Marti
        ads = ase.Atoms([ase.Atom(adsorbate)])
193 dadc6016 Carles Marti
    pos += height * norm_vect_u
194 1d22a086 Carles Marti
    # Move adsorbate into position
195 dadc6016 Carles Marti
    ads.translate(pos - ctr_coord)
196 1d22a086 Carles Marti
    # Attach the adsorbate
197 1d22a086 Carles Marti
    slab.extend(ads)
198 1d22a086 Carles Marti
199 1d22a086 Carles Marti
200 a4b57124 Carles Marti
def check_collision(slab_molec, slab_num_atoms, min_height, vect, nn_slab=0,
201 e8bebcca Carles Marti
                    nn_molec=0, coll_coeff=1.2):
202 5f3f4b69 Carles Marti
    """Checks whether a slab and a molecule collide or not.
203 5f3f4b69 Carles Marti

204 5f3f4b69 Carles Marti
    @param slab_molec: The system of adsorbate-slab for which to detect if there
205 5f3f4b69 Carles Marti
        are collisions.
206 5f3f4b69 Carles Marti
    @param nn_slab: Number of neigbors in the surface.
207 5f3f4b69 Carles Marti
    @param nn_molec: Number of neighbors in the molecule.
208 5f3f4b69 Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
209 5f3f4b69 Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
210 5f3f4b69 Carles Marti
        considered as atomic collision.
211 5f3f4b69 Carles Marti
    @param slab_num_atoms: Number of atoms of the bare slab.
212 5f3f4b69 Carles Marti
    @param min_height: The minimum height atoms can have to not be considered as
213 5f3f4b69 Carles Marti
        colliding.
214 5f3f4b69 Carles Marti
    @param vect: The vector perpendicular to the slab.
215 5f3f4b69 Carles Marti
    @return: bool, whether the surface and the molecule collide.
216 e8bebcca Carles Marti
    """
217 5f3f4b69 Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
218 e8bebcca Carles Marti
219 e8bebcca Carles Marti
    # Check structure overlap by height
220 5fb01677 Carles Marti
    if min_height is not False:
221 a4b57124 Carles Marti
        cart_axes = [[1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0],
222 a4b57124 Carles Marti
                     [-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, -1.0]]
223 e8bebcca Carles Marti
        if vect.tolist() not in cart_axes:
224 e8bebcca Carles Marti
            err_msg = "'min_coll_height' option is only implemented for " \
225 e8bebcca Carles Marti
                      "'surf_norm_vect' to be one of the x, y or z axes. "
226 e8bebcca Carles Marti
            logger.error(err_msg)
227 e8bebcca Carles Marti
            raise ValueError(err_msg)
228 e8bebcca Carles Marti
        for atom in slab_molec[slab_num_atoms:]:
229 e8bebcca Carles Marti
            for i, coord in enumerate(vect):
230 e8bebcca Carles Marti
                if coord == 0:
231 e8bebcca Carles Marti
                    continue
232 e8bebcca Carles Marti
                if atom.position[i] * coord < min_height * coord:
233 e8bebcca Carles Marti
                    return True
234 e8bebcca Carles Marti
235 e8bebcca Carles Marti
    # Check structure overlap by sphere collision
236 e8bebcca Carles Marti
    if coll_coeff is not False:
237 5fb01677 Carles Marti
        slab_molec_cutoffs = natural_cutoffs(slab_molec, mult=coll_coeff)
238 b4b2f307 Carles Marti
        slab_molec_nghbs = len(
239 b4b2f307 Carles Marti
            neighbor_list("i", slab_molec, slab_molec_cutoffs))
240 5fb01677 Carles Marti
        if slab_molec_nghbs > nn_slab + nn_molec:
241 5fb01677 Carles Marti
            return True
242 e8bebcca Carles Marti
243 e8bebcca Carles Marti
    return False
244 e8bebcca Carles Marti
245 e8bebcca Carles Marti
246 e8bebcca Carles Marti
def correct_coll(molec, slab, ctr_coord, site_coord, num_pts,
247 e8bebcca Carles Marti
                 min_coll_height, norm_vect, slab_nghbs, molec_nghbs,
248 e8bebcca Carles Marti
                 coll_coeff, height=2.5):
249 e8bebcca Carles Marti
    # TODO Rethink this function
250 e8bebcca Carles Marti
    """Tries to adsorb a molecule on a slab trying to avoid collisions by doing
251 e8bebcca Carles Marti
    small rotations.
252 e8bebcca Carles Marti

253 e8bebcca Carles Marti
    @param molec: ase.Atoms object of the molecule to adsorb
254 e8bebcca Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the
255 e8bebcca Carles Marti
        molecule
256 e8bebcca Carles Marti
    @param ctr_coord: The coordinates of the molecule to use as adsorption
257 e8bebcca Carles Marti
        center.
258 e8bebcca Carles Marti
    @param site_coord: The coordinates of the surface on which to adsorb the
259 e8bebcca Carles Marti
        molecule
260 e8bebcca Carles Marti
    @param num_pts: Number on which to sample Euler angles.
261 e8bebcca Carles Marti
    @param min_coll_height: The lowermost height for which to detect a collision
262 e8bebcca Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
263 e8bebcca Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
264 e8bebcca Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
265 e8bebcca Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
266 e8bebcca Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
267 e8bebcca Carles Marti
        considered as atomic collision.
268 e8bebcca Carles Marti
    @param height: Height on which to try adsorption
269 e8bebcca Carles Marti
    @return collision: bool, whether the structure generated has collisions
270 e8bebcca Carles Marti
        between slab and adsorbate.
271 e8bebcca Carles Marti
    """
272 e8bebcca Carles Marti
    from copy import deepcopy
273 e8bebcca Carles Marti
    slab_num_atoms = len(slab)
274 e8bebcca Carles Marti
    slab_molec = []
275 e8bebcca Carles Marti
    collision = True
276 e8bebcca Carles Marti
    max_corr = 6  # Should be an even number
277 e8bebcca Carles Marti
    d_angle = 180 / ((max_corr / 2.0) * num_pts)
278 e8bebcca Carles Marti
    num_corr = 0
279 e8bebcca Carles Marti
    while collision and num_corr <= max_corr:
280 e8bebcca Carles Marti
        k = num_corr * (-1) ** num_corr
281 e8bebcca Carles Marti
        slab_molec = deepcopy(slab)
282 e8bebcca Carles Marti
        molec.euler_rotate(k * d_angle, k * d_angle / 2, k * d_angle,
283 e8bebcca Carles Marti
                           center=ctr_coord)
284 e8bebcca Carles Marti
        add_adsorbate(slab_molec, molec, site_coord, ctr_coord, height,
285 e8bebcca Carles Marti
                      norm_vect=norm_vect)
286 e8bebcca Carles Marti
        collision = check_collision(slab_molec, slab_num_atoms, min_coll_height,
287 e8bebcca Carles Marti
                                    norm_vect, slab_nghbs, molec_nghbs,
288 e8bebcca Carles Marti
                                    coll_coeff)
289 e8bebcca Carles Marti
        num_corr += 1
290 e8bebcca Carles Marti
    return slab_molec, collision
291 5f3f4b69 Carles Marti
292 5f3f4b69 Carles Marti
293 c25aa299 Carles Marti
def dissociate_h(slab_molec_orig, h_idx, num_atoms_slab, h_acceptor,
294 c25aa299 Carles Marti
                 neigh_cutoff=1):
295 b4b2f307 Carles Marti
    # TODO rethink
296 91ae8d86 Carles Marti
    """Tries to dissociate a H from the molecule and adsorbs it on the slab.
297 b4b2f307 Carles Marti

298 91ae8d86 Carles Marti
    Tries to dissociate a H atom from the molecule and adsorb in on top of the
299 91ae8d86 Carles Marti
    surface if the distance is shorter than two times the neigh_cutoff value.
300 b4b2f307 Carles Marti
    @param slab_molec_orig: The ase.Atoms object of the system adsorbate-slab.
301 b4b2f307 Carles Marti
    @param h_idx: The index of the hydrogen atom to carry out adsorption of.
302 b4b2f307 Carles Marti
    @param num_atoms_slab: The number of atoms of the slab without adsorbate.
303 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
304 b4b2f307 Carles Marti
    @param neigh_cutoff: half the maximum distance between the surface and the
305 b4b2f307 Carles Marti
        H for it to carry out dissociation.
306 b4b2f307 Carles Marti
    @return: An ase.Atoms object of the system adsorbate-surface with H
307 b4b2f307 Carles Marti
    """
308 b4b2f307 Carles Marti
    from copy import deepcopy
309 b4b2f307 Carles Marti
    from ase.neighborlist import NeighborList
310 b4b2f307 Carles Marti
    slab_molec = deepcopy(slab_molec_orig)
311 b4b2f307 Carles Marti
    cutoffs = len(slab_molec) * [neigh_cutoff]
312 c25aa299 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True, skin=0.0)
313 b4b2f307 Carles Marti
    nl.update(slab_molec)
314 b4b2f307 Carles Marti
    surf_h_vect = np.array([np.infty] * 3)
315 c25aa299 Carles Marti
    if h_acceptor == 'all':
316 c25aa299 Carles Marti
        h_acceptor = list(range(num_atoms_slab))
317 b4b2f307 Carles Marti
    for neigh_idx in nl.get_neighbors(h_idx)[0]:
318 c25aa299 Carles Marti
        for elm in h_acceptor:
319 c25aa299 Carles Marti
            if isinstance(elm, int):
320 c25aa299 Carles Marti
                if neigh_idx == elm and neigh_idx < num_atoms_slab:
321 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
322 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
323 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
324 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
325 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
326 c25aa299 Carles Marti
            else:
327 c25aa299 Carles Marti
                if slab_molec[neigh_idx].symbol == elm \
328 c25aa299 Carles Marti
                        and neigh_idx < num_atoms_slab:
329 c25aa299 Carles Marti
                    dist = np.linalg.norm(slab_molec[neigh_idx].position -
330 c25aa299 Carles Marti
                                          slab_molec[h_idx].position)
331 c25aa299 Carles Marti
                    if dist < np.linalg.norm(surf_h_vect):
332 c25aa299 Carles Marti
                        surf_h_vect = slab_molec[neigh_idx].position \
333 c25aa299 Carles Marti
                                      - slab_molec[h_idx].position
334 c25aa299 Carles Marti
335 b4b2f307 Carles Marti
    if np.linalg.norm(surf_h_vect) != np.infty:
336 b4b2f307 Carles Marti
        trans_vect = surf_h_vect - surf_h_vect / np.linalg.norm(surf_h_vect)
337 b4b2f307 Carles Marti
        slab_molec[h_idx].position = slab_molec[h_idx].position + trans_vect
338 b4b2f307 Carles Marti
        return slab_molec
339 b4b2f307 Carles Marti
340 b4b2f307 Carles Marti
341 c25aa299 Carles Marti
def dissociation(slab_molec, h_donor, h_acceptor, num_atoms_slab):
342 b4b2f307 Carles Marti
    # TODO multiple dissociation
343 b4b2f307 Carles Marti
    """Decides which H atoms to dissociate according to a list of atoms.
344 b4b2f307 Carles Marti

345 b4b2f307 Carles Marti
    Given a list of chemical symbols or atom indices it checks for every atom
346 b4b2f307 Carles Marti
    or any of its neighbor if it's a H and calls dissociate_h to try to carry
347 b4b2f307 Carles Marti
    out dissociation of that H. For atom indices, it checks both whether
348 b4b2f307 Carles Marti
    the atom index or its neighbors are H, for chemical symbols, it only checks
349 b4b2f307 Carles Marti
    if there is a neighbor H.
350 b4b2f307 Carles Marti
    @param slab_molec: The ase.Atoms object of the system adsorbate-slab.
351 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
352 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
353 c25aa299 Carles Marti
    @param num_atoms_slab: Number of atoms of the bare slab.
354 b4b2f307 Carles Marti
    @return:
355 b4b2f307 Carles Marti
    """
356 b4b2f307 Carles Marti
    from ase.neighborlist import natural_cutoffs, NeighborList
357 b4b2f307 Carles Marti
    molec = slab_molec[num_atoms_slab:]
358 b4b2f307 Carles Marti
    cutoffs = natural_cutoffs(molec)
359 b4b2f307 Carles Marti
    nl = NeighborList(cutoffs, self_interaction=False, bothways=True)
360 b4b2f307 Carles Marti
    nl.update(molec)
361 b4b2f307 Carles Marti
    disso_structs = []
362 c25aa299 Carles Marti
    for el in h_donor:
363 b4b2f307 Carles Marti
        if isinstance(el, int):
364 b4b2f307 Carles Marti
            if molec[el].symbol == 'H':
365 b4b2f307 Carles Marti
                disso_struct = dissociate_h(slab_molec, el + num_atoms_slab,
366 c25aa299 Carles Marti
                                            num_atoms_slab, h_acceptor)
367 b4b2f307 Carles Marti
                if disso_struct is not None:
368 b4b2f307 Carles Marti
                    disso_structs.append(disso_struct)
369 b4b2f307 Carles Marti
            else:
370 b4b2f307 Carles Marti
                for neigh_idx in nl.get_neighbors(el)[0]:
371 b4b2f307 Carles Marti
                    if molec[neigh_idx].symbol == 'H':
372 b4b2f307 Carles Marti
                        disso_struct = dissociate_h(slab_molec, neigh_idx +
373 b4b2f307 Carles Marti
                                                    num_atoms_slab,
374 c25aa299 Carles Marti
                                                    num_atoms_slab, h_acceptor)
375 b4b2f307 Carles Marti
                        if disso_struct is not None:
376 b4b2f307 Carles Marti
                            disso_structs.append(disso_struct)
377 b4b2f307 Carles Marti
        else:
378 b4b2f307 Carles Marti
            for atom in molec:
379 b4b2f307 Carles Marti
                if atom.symbol.lower() == el.lower():
380 b4b2f307 Carles Marti
                    for neigh_idx in nl.get_neighbors(atom.index)[0]:
381 b4b2f307 Carles Marti
                        if molec[neigh_idx].symbol == 'H':
382 5261a07f Carles Marti
                            disso_struct = dissociate_h(slab_molec, neigh_idx
383 b4b2f307 Carles Marti
                                                        + num_atoms_slab,
384 c25aa299 Carles Marti
                                                        num_atoms_slab,
385 c25aa299 Carles Marti
                                                        h_acceptor)
386 b4b2f307 Carles Marti
                            if disso_struct is not None:
387 b4b2f307 Carles Marti
                                disso_structs.append(disso_struct)
388 b4b2f307 Carles Marti
    return disso_structs
389 b4b2f307 Carles Marti
390 b4b2f307 Carles Marti
391 3ab0865c Carles Marti
def ads_euler(orig_molec, slab, ctr_coord, site_coord, num_pts,
392 b4b2f307 Carles Marti
              min_coll_height, coll_coeff, norm_vect, slab_nghbs, molec_nghbs,
393 c25aa299 Carles Marti
              h_donor, h_acceptor):
394 3ab0865c Carles Marti
    """Generates adsorbate-surface structures by sampling over Euler angles.
395 3ab0865c Carles Marti

396 3ab0865c Carles Marti
    This function generates a number of adsorbate-surface structures at
397 3ab0865c Carles Marti
    different orientations of the adsorbate sampled at multiple Euler (zxz)
398 3ab0865c Carles Marti
    angles.
399 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb.
400 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the
401 5261a07f Carles Marti
        molecule.
402 3ab0865c Carles Marti
    @param ctr_coord: The coordinates of the molecule to use as adsorption
403 3ab0865c Carles Marti
        center.
404 3ab0865c Carles Marti
    @param site_coord: The coordinates of the surface on which to adsorb the
405 5261a07f Carles Marti
        molecule.
406 3ab0865c Carles Marti
    @param num_pts: Number on which to sample Euler angles.
407 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision.
408 3ab0865c Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
409 3ab0865c Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
410 3ab0865c Carles Marti
        considered as atomic collision.
411 3ab0865c Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
412 3ab0865c Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
413 3ab0865c Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
414 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
415 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
416 3ab0865c Carles Marti
    @return: list of ase.Atoms object conatining all the orientations of a given
417 5261a07f Carles Marti
        conformer.
418 3ab0865c Carles Marti
    """
419 3ab0865c Carles Marti
    from copy import deepcopy
420 3ab0865c Carles Marti
    slab_ads_list = []
421 3ab0865c Carles Marti
    # rotation around z
422 3ab0865c Carles Marti
    for alpha in np.arange(0, 360, 360 / num_pts):
423 3ab0865c Carles Marti
        # rotation around x'
424 3ab0865c Carles Marti
        for beta in np.arange(0, 180, 180 / num_pts):
425 3ab0865c Carles Marti
            # rotation around z'
426 3ab0865c Carles Marti
            for gamma in np.arange(0, 360, 360 / num_pts):
427 3ab0865c Carles Marti
                molec = deepcopy(orig_molec)
428 3ab0865c Carles Marti
                molec.euler_rotate(alpha, beta, gamma, center=ctr_coord)
429 3ab0865c Carles Marti
                slab_molec, collision = correct_coll(molec, slab,
430 5fb01677 Carles Marti
                                                     ctr_coord, site_coord,
431 5fb01677 Carles Marti
                                                     num_pts, min_coll_height,
432 5fb01677 Carles Marti
                                                     norm_vect,
433 5fb01677 Carles Marti
                                                     slab_nghbs, molec_nghbs,
434 5fb01677 Carles Marti
                                                     coll_coeff)
435 5261a07f Carles Marti
                if not collision and not any([np.allclose(slab_molec.positions,
436 5261a07f Carles Marti
                                                          conf.positions)
437 5261a07f Carles Marti
                                              for conf in slab_ads_list]):
438 3ab0865c Carles Marti
                    slab_ads_list.append(slab_molec)
439 c25aa299 Carles Marti
                    if h_donor is not False:
440 c25aa299 Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
441 c25aa299 Carles Marti
                                                          h_acceptor,
442 c25aa299 Carles Marti
                                                          len(slab)))
443 3ab0865c Carles Marti
444 3ab0865c Carles Marti
    return slab_ads_list
445 f3d1e601 Carles Marti
446 f3d1e601 Carles Marti
447 7dd94df7 Carles Marti
def chemcat_rotate(molecule, surf, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
448 7dd94df7 Carles Marti
                   ctr2_surf, bond_vector, bond_angle_target,
449 7dd94df7 Carles Marti
                   dihedral_angle_target=None, mol_dihedral_angle_target=None):
450 7dd94df7 Carles Marti
    """Performs translation and rotation of an adsorbate defined by an
451 7dd94df7 Carles Marti
    adsorption bond length, direction, angle and dihedral angle
452 7dd94df7 Carles Marti

453 7dd94df7 Carles Marti
    Carles modification of chemcat's transform_adsorbate to work with
454 7dd94df7 Carles Marti
    coordinates instead of ase.Atom
455 7dd94df7 Carles Marti
    Parameters:
456 7dd94df7 Carles Marti
        molecule (ase.Atoms): The molecule to adsorb.
457 7dd94df7 Carles Marti

458 7dd94df7 Carles Marti
        surf (ase.Atoms): The surface ontop of which to adsorb.
459 7dd94df7 Carles Marti

460 7dd94df7 Carles Marti
        ctr1_mol (int/list): The position of the adsorption center in the
461 7dd94df7 Carles Marti
        molecule that will be bound to the surface.
462 7dd94df7 Carles Marti

463 7dd94df7 Carles Marti
        ctr2_mol (int/list): The position of a second center of the
464 7dd94df7 Carles Marti
        adsorbate used to define the adsorption bond angle, and the dihedral
465 7dd94df7 Carles Marti
        adsorption angle.
466 7dd94df7 Carles Marti

467 7dd94df7 Carles Marti
        ctr3_mol (int/list): The position of a third center in the
468 7dd94df7 Carles Marti
        adsorbate used to define the adsorbate dihedral angle.
469 7dd94df7 Carles Marti

470 7dd94df7 Carles Marti
        ctr1_surf (int/list): The position of the site on the surface that
471 7dd94df7 Carles Marti
        will be bound to the molecule.
472 7dd94df7 Carles Marti

473 7dd94df7 Carles Marti
        ctr2_surf (int/list): The position of a second center of the
474 7dd94df7 Carles Marti
        surface used to define the dihedral adsorption angle.
475 7dd94df7 Carles Marti

476 7dd94df7 Carles Marti
        bond_vector (numpy.ndarray): The adsorption bond desired.
477 7dd94df7 Carles Marti
            Details: offset = vect(atom1_surf, atom1_mol)
478 7dd94df7 Carles Marti

479 7dd94df7 Carles Marti
        bond_angle_target (float or int): The adsorption bond angle desired (in
480 7dd94df7 Carles Marti
            degrees).
481 7dd94df7 Carles Marti
            Details: bond_angle_target = angle(atom1_surf, atom1_mol, atom2_mol)
482 7dd94df7 Carles Marti

483 7dd94df7 Carles Marti
        dihedral_angle_target (float or int): The dihedral adsorption angle
484 7dd94df7 Carles Marti
            desired (in degrees).
485 7dd94df7 Carles Marti
            Details: dihedral_angle_target = dihedral(atom2_surf, atom1_surf,
486 7dd94df7 Carles Marti
            atom1_mol, atom2_mol)
487 7dd94df7 Carles Marti
                The rotation vector is facing the adsorbate from the surface
488 7dd94df7 Carles Marti
                (i.e. counterclockwise rotation when facing the surface (i.e.
489 7dd94df7 Carles Marti
                view from top))
490 7dd94df7 Carles Marti

491 7dd94df7 Carles Marti
        mol_dihedral_angle_target (float or int): The adsorbate dihedral angle
492 7dd94df7 Carles Marti
            desired (in degrees).
493 7dd94df7 Carles Marti
            Details: mol_dihedral_angle_target = dihedral(atom1_surf, atom1_mol,
494 7dd94df7 Carles Marti
            atom2_mol, atom3_mol)
495 7dd94df7 Carles Marti
                The rotation vector is facing atom2_mol from atom1_mol
496 7dd94df7 Carles Marti

497 7dd94df7 Carles Marti
    Returns:
498 7dd94df7 Carles Marti
        None (the `molecule` object is modified in-place)
499 7dd94df7 Carles Marti
    """
500 7dd94df7 Carles Marti
    vect_surf = get_atom_coords(surf, ctr2_surf) - get_atom_coords(surf,
501 7dd94df7 Carles Marti
                                                                   ctr1_surf)
502 7dd94df7 Carles Marti
    vect_inter = get_atom_coords(molecule, ctr1_mol) - get_atom_coords(surf,
503 7dd94df7 Carles Marti
                                                                       ctr1_surf)
504 7dd94df7 Carles Marti
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
505 7dd94df7 Carles Marti
                                                                     ctr1_mol)
506 7dd94df7 Carles Marti
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
507 7dd94df7 Carles Marti
                                                                      ctr2_mol)
508 7dd94df7 Carles Marti
509 7dd94df7 Carles Marti
    # Check if dihedral angles can be defined
510 7dd94df7 Carles Marti
    do_dihedral = dihedral_angle_target is not None
511 7dd94df7 Carles Marti
    do_mol_dihedral = mol_dihedral_angle_target is not None
512 7dd94df7 Carles Marti
    dihedral_use_mol2 = False
513 7dd94df7 Carles Marti
    # Check if vect_surf and vect_inter are not aligned
514 7dd94df7 Carles Marti
    if np.allclose(np.cross(vect_surf, vect_inter), 0):
515 7dd94df7 Carles Marti
        logger.warning(
516 7dd94df7 Carles Marti
            "Surface atoms are incompatible with adsorption "
517 7dd94df7 Carles Marti
            "direction/bond. An adsorption dihedral angle cannot be defined.")
518 7dd94df7 Carles Marti
        do_dihedral = False
519 7dd94df7 Carles Marti
    # Check if requested bond angle is not flat
520 7dd94df7 Carles Marti
    if np.isclose((bond_angle_target + 90) % 180 - 90, 0):
521 7dd94df7 Carles Marti
        logger.warning("Requested bond angle is flat. Only a single dihedral "
522 7dd94df7 Carles Marti
                       "angle can be defined (ctr2_surf, ctr1_surf, ctr2_mol, "
523 7dd94df7 Carles Marti
                       "ctr3_mol).")
524 7dd94df7 Carles Marti
        do_mol_dihedral = False
525 7dd94df7 Carles Marti
        dihedral_use_mol2 = True
526 7dd94df7 Carles Marti
        logger.warning("Dihedral adsorption angle rotation will be perfomed "
527 7dd94df7 Carles Marti
                       "with (ctr2_surf, ctr1_surf, ctr2_mol, ctr3_mol).")
528 7dd94df7 Carles Marti
    # Check if vect_mol and vect2_mol are not aligned
529 7dd94df7 Carles Marti
    if np.allclose(np.cross(vect_mol, vect2_mol), 0):
530 7dd94df7 Carles Marti
        logger.warning("Adsorbates atoms are aligned. An adsorbate dihedral "
531 7dd94df7 Carles Marti
                       "angle cannot be defined.")
532 7dd94df7 Carles Marti
        do_mol_dihedral = False
533 7dd94df7 Carles Marti
    if not do_dihedral:
534 7dd94df7 Carles Marti
        logger.warning("Dihedral adsorption angle rotation will not be "
535 7dd94df7 Carles Marti
                       "performed.")
536 7dd94df7 Carles Marti
    if not do_mol_dihedral:
537 7dd94df7 Carles Marti
        logger.warning("Adsorbate dihedral angle rotation will not be "
538 7dd94df7 Carles Marti
                       "performed.")
539 7dd94df7 Carles Marti
540 7dd94df7 Carles Marti
    ###########################
541 7dd94df7 Carles Marti
    #       Translation       #
542 7dd94df7 Carles Marti
    ###########################
543 7dd94df7 Carles Marti
544 7dd94df7 Carles Marti
    # Compute and apply translation of adsorbate
545 7dd94df7 Carles Marti
    translation = bond_vector - vect_inter
546 7dd94df7 Carles Marti
    molecule.translate(translation)
547 7dd94df7 Carles Marti
548 7dd94df7 Carles Marti
    # Update adsorption bond
549 7dd94df7 Carles Marti
    vect_inter = get_atom_coords(molecule, ctr1_mol) - get_atom_coords(surf,
550 7dd94df7 Carles Marti
                                                                       ctr1_surf)
551 7dd94df7 Carles Marti
552 7dd94df7 Carles Marti
    # Check if translation was successful
553 7dd94df7 Carles Marti
    if np.allclose(vect_inter, bond_vector):
554 7dd94df7 Carles Marti
        pass  # print("Translation successfully applied (error: ~ {:.5g} unit "
555 7dd94df7 Carles Marti
        # "length)".format(np.linalg.norm(vect_inter - bond_vector)))
556 7dd94df7 Carles Marti
    else:
557 7dd94df7 Carles Marti
        err = 'An unknown error occured during the translation'
558 7dd94df7 Carles Marti
        logger.error(err)
559 7dd94df7 Carles Marti
        raise AssertionError(err)
560 7dd94df7 Carles Marti
561 7dd94df7 Carles Marti
    ###########################
562 7dd94df7 Carles Marti
    #   Bond angle rotation   #
563 7dd94df7 Carles Marti
    ###########################
564 7dd94df7 Carles Marti
565 7dd94df7 Carles Marti
    # Compute rotation vector
566 7dd94df7 Carles Marti
    rotation_vector = np.cross(-vect_inter, vect_mol)
567 7dd94df7 Carles Marti
    if np.allclose(rotation_vector, 0, atol=1e-5):
568 7dd94df7 Carles Marti
        # If molecular bonds are aligned, any vector orthogonal to vect_inter
569 7dd94df7 Carles Marti
        # can be used Such vector can be found as the orthogonal rejection of
570 7dd94df7 Carles Marti
        # either X-axis, Y-axis or Z-axis onto vect_inter (since they cannot
571 7dd94df7 Carles Marti
        # be all aligned)
572 7dd94df7 Carles Marti
        non_aligned_vector = np.zeros(3)
573 7dd94df7 Carles Marti
        # Select the most orthogonal axis (lowest dot product):
574 7dd94df7 Carles Marti
        non_aligned_vector[np.argmin(np.abs(vect_inter))] = 1
575 7dd94df7 Carles Marti
        rotation_vector = non_aligned_vector - np.dot(non_aligned_vector,
576 7dd94df7 Carles Marti
                                                      vect_inter) / np.dot(
577 7dd94df7 Carles Marti
            vect_inter, vect_inter) * vect_inter
578 7dd94df7 Carles Marti
579 7dd94df7 Carles Marti
    # Retrieve current bond angle
580 7dd94df7 Carles Marti
    bond_angle_ini = get_vect_angle(-vect_inter, vect_mol, rotation_vector)
581 7dd94df7 Carles Marti
582 7dd94df7 Carles Marti
    # Apply rotation to reach desired bond_angle
583 7dd94df7 Carles Marti
    molecule.rotate(bond_angle_target - bond_angle_ini, v=rotation_vector,
584 7dd94df7 Carles Marti
                    center=get_atom_coords(molecule, ctr1_mol))
585 7dd94df7 Carles Marti
586 7dd94df7 Carles Marti
    # Update molecular bonds
587 7dd94df7 Carles Marti
    vect_mol = get_atom_coords(molecule, ctr2_mol) - get_atom_coords(molecule,
588 7dd94df7 Carles Marti
                                                                     ctr1_mol)
589 7dd94df7 Carles Marti
    vect2_mol = get_atom_coords(molecule, ctr3_mol) - get_atom_coords(molecule,
590 7dd94df7 Carles Marti
                                                                      ctr2_mol)
591 7dd94df7 Carles Marti
592 7dd94df7 Carles Marti
    # Check if rotation was successful
593 7dd94df7 Carles Marti
    bond_angle = get_vect_angle(-vect_inter, vect_mol)
594 7dd94df7 Carles Marti
    if np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90, 0,
595 7dd94df7 Carles Marti
                  atol=1e-3) and np.allclose(get_atom_coords(molecule, ctr1_mol)
596 7dd94df7 Carles Marti
                                             - get_atom_coords(surf,
597 7dd94df7 Carles Marti
                                                               ctr1_surf),
598 7dd94df7 Carles Marti
                                             vect_inter):
599 7dd94df7 Carles Marti
        pass  # print("Rotation successfully applied (error: {:.5f}°)".format(
600 7dd94df7 Carles Marti
        # (bond_angle - bond_angle_target + 90) % 180 - 90))
601 7dd94df7 Carles Marti
    else:
602 7dd94df7 Carles Marti
        err = 'An unknown error occured during the rotation'
603 7dd94df7 Carles Marti
        logger.error(err)
604 7dd94df7 Carles Marti
        raise AssertionError(err)
605 7dd94df7 Carles Marti
606 7dd94df7 Carles Marti
    ###########################
607 7dd94df7 Carles Marti
    # Dihedral angle rotation #
608 7dd94df7 Carles Marti
    ###########################
609 7dd94df7 Carles Marti
610 7dd94df7 Carles Marti
    # Perform dihedral rotation if possible
611 7dd94df7 Carles Marti
    if do_dihedral:
612 7dd94df7 Carles Marti
        # Retrieve current dihedral angle (by computing the angle between the
613 7dd94df7 Carles Marti
        # vector rejection of vect_surf and vect_mol onto vect_inter)
614 7dd94df7 Carles Marti
        vect_inter_inner = np.dot(vect_inter, vect_inter)
615 7dd94df7 Carles Marti
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
616 7dd94df7 Carles Marti
                           vect_inter_inner * vect_inter
617 7dd94df7 Carles Marti
        if dihedral_use_mol2:
618 7dd94df7 Carles Marti
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
619 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
620 7dd94df7 Carles Marti
        else:
621 7dd94df7 Carles Marti
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
622 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
623 7dd94df7 Carles Marti
        dihedral_angle_ini = get_vect_angle(vect_surf_reject, vect_mol_reject,
624 7dd94df7 Carles Marti
                                            vect_inter)
625 7dd94df7 Carles Marti
626 7dd94df7 Carles Marti
        # Apply dihedral rotation along vect_inter
627 7dd94df7 Carles Marti
        molecule.rotate(dihedral_angle_target - dihedral_angle_ini,
628 7dd94df7 Carles Marti
                        v=vect_inter, center=get_atom_coords(molecule,
629 7dd94df7 Carles Marti
                                                             ctr1_mol))
630 7dd94df7 Carles Marti
631 7dd94df7 Carles Marti
        # Update molecular bonds
632 7dd94df7 Carles Marti
        vect_mol = get_atom_coords(molecule, ctr2_mol) - \
633 7dd94df7 Carles Marti
                   get_atom_coords(molecule, ctr1_mol)
634 7dd94df7 Carles Marti
        vect2_mol = get_atom_coords(molecule, ctr3_mol) - \
635 7dd94df7 Carles Marti
                    get_atom_coords(molecule, ctr2_mol)
636 7dd94df7 Carles Marti
637 7dd94df7 Carles Marti
        # Check if rotation was successful
638 7dd94df7 Carles Marti
        # Check dihedral rotation
639 7dd94df7 Carles Marti
        if dihedral_use_mol2:
640 7dd94df7 Carles Marti
            vect_mol_reject = vect2_mol - np.dot(vect2_mol, vect_inter) / \
641 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
642 7dd94df7 Carles Marti
        else:
643 7dd94df7 Carles Marti
            vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
644 7dd94df7 Carles Marti
                              vect_inter_inner * vect_inter
645 7dd94df7 Carles Marti
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
646 7dd94df7 Carles Marti
                                        vect_inter)
647 7dd94df7 Carles Marti
        # Check bond rotation is unmodified
648 7dd94df7 Carles Marti
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
649 7dd94df7 Carles Marti
        if np.isclose((dihedral_angle - dihedral_angle_target + 90) % 180 - 90,
650 5261a07f Carles Marti
                      0, atol=1e-3) \
651 5261a07f Carles Marti
                and np.isclose((bond_angle - bond_angle_target + 90) %
652 5261a07f Carles Marti
                               180 - 90, 0, atol=1e-5) \
653 c25aa299 Carles Marti
                and np.allclose(get_atom_coords(molecule, ctr1_mol)
654 c25aa299 Carles Marti
                                - get_atom_coords(surf, ctr1_surf),
655 c25aa299 Carles Marti
                                vect_inter):
656 7dd94df7 Carles Marti
            pass  # print( "Dihedral rotation successfully applied (error: {
657 7dd94df7 Carles Marti
            # :.5f}°)".format((dihedral_angle - dihedral_angle_target + 90) %
658 7dd94df7 Carles Marti
            # 180 - 90))
659 7dd94df7 Carles Marti
        else:
660 7dd94df7 Carles Marti
            err = 'An unknown error occured during the dihedral rotation'
661 7dd94df7 Carles Marti
            logger.error(err)
662 7dd94df7 Carles Marti
            raise AssertionError(err)
663 7dd94df7 Carles Marti
664 7dd94df7 Carles Marti
    #####################################
665 7dd94df7 Carles Marti
    # Adsorbate dihedral angle rotation #
666 7dd94df7 Carles Marti
    #####################################
667 7dd94df7 Carles Marti
668 7dd94df7 Carles Marti
    # Perform adsorbate dihedral rotation if possible
669 7dd94df7 Carles Marti
    if do_mol_dihedral:
670 7dd94df7 Carles Marti
        # Retrieve current adsorbate dihedral angle (by computing the angle
671 7dd94df7 Carles Marti
        # between the orthogonal rejection of vect_inter and vect2_mol onto
672 7dd94df7 Carles Marti
        # vect_mol)
673 7dd94df7 Carles Marti
        vect_mol_inner = np.dot(vect_mol, vect_mol)
674 7dd94df7 Carles Marti
        bond_inter_reject = -vect_inter - np.dot(-vect_inter, vect_mol) / \
675 5261a07f Carles Marti
            vect_mol_inner * vect_mol
676 7dd94df7 Carles Marti
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
677 5261a07f Carles Marti
            vect_mol_inner * vect_mol
678 7dd94df7 Carles Marti
        dihedral_angle_ini = get_vect_angle(bond_inter_reject,
679 7dd94df7 Carles Marti
                                            bond2_mol_reject, vect_mol)
680 7dd94df7 Carles Marti
681 7dd94df7 Carles Marti
        # Apply dihedral rotation along vect_mol
682 7dd94df7 Carles Marti
        molecule.rotate(mol_dihedral_angle_target - dihedral_angle_ini,
683 7dd94df7 Carles Marti
                        v=vect_mol, center=get_atom_coords(molecule, ctr1_mol))
684 7dd94df7 Carles Marti
685 7dd94df7 Carles Marti
        # Update molecular bonds
686 7dd94df7 Carles Marti
        vect_mol = get_atom_coords(molecule, ctr2_mol) \
687 5261a07f Carles Marti
            - get_atom_coords(molecule, ctr1_mol)
688 7dd94df7 Carles Marti
        vect2_mol = get_atom_coords(molecule, ctr3_mol) \
689 5261a07f Carles Marti
            - get_atom_coords(molecule, ctr2_mol)
690 7dd94df7 Carles Marti
691 7dd94df7 Carles Marti
        # Check if rotation was successful
692 7dd94df7 Carles Marti
        # Check adsorbate dihedral rotation
693 7dd94df7 Carles Marti
        vect_mol_inner = np.dot(vect_mol, vect_mol)
694 7dd94df7 Carles Marti
        bond2_mol_reject = vect2_mol - np.dot(vect2_mol, vect_mol) / \
695 5261a07f Carles Marti
            vect_mol_inner * vect_mol
696 7dd94df7 Carles Marti
        mol_dihedral_angle = get_vect_angle(bond_inter_reject,
697 7dd94df7 Carles Marti
                                            bond2_mol_reject, vect_mol)
698 7dd94df7 Carles Marti
        # Check dihedral rotation
699 7dd94df7 Carles Marti
        vect_inter_inner = np.dot(vect_inter, vect_inter)
700 7dd94df7 Carles Marti
        vect_surf_reject = vect_surf - np.dot(vect_surf, vect_inter) / \
701 5261a07f Carles Marti
            vect_inter_inner * vect_inter
702 7dd94df7 Carles Marti
        vect_mol_reject = vect_mol - np.dot(vect_mol, vect_inter) / \
703 5261a07f Carles Marti
            vect_inter_inner * vect_inter
704 7dd94df7 Carles Marti
        dihedral_angle = get_vect_angle(vect_surf_reject, vect_mol_reject,
705 7dd94df7 Carles Marti
                                        vect_inter)
706 7dd94df7 Carles Marti
        # Check bond rotation is unmodified
707 7dd94df7 Carles Marti
        bond_angle = get_vect_angle(-vect_inter, vect_mol)
708 7dd94df7 Carles Marti
        if np.isclose((mol_dihedral_angle - mol_dihedral_angle_target + 90) %
709 7dd94df7 Carles Marti
                      180 - 90, 0, atol=1e-3) \
710 7dd94df7 Carles Marti
                and np.isclose((dihedral_angle -
711 7dd94df7 Carles Marti
                                dihedral_angle_target + 90) % 180 - 90, 0,
712 7dd94df7 Carles Marti
                               atol=1e-5) \
713 7dd94df7 Carles Marti
                and np.isclose((bond_angle - bond_angle_target + 90) % 180 - 90,
714 7dd94df7 Carles Marti
                               0, atol=1e-5) \
715 7dd94df7 Carles Marti
                and np.allclose(get_atom_coords(molecule, ctr1_mol) -
716 7dd94df7 Carles Marti
                                get_atom_coords(surf, ctr1_surf),
717 7dd94df7 Carles Marti
                                vect_inter):
718 7dd94df7 Carles Marti
            pass  # print(
719 7dd94df7 Carles Marti
            # "Adsorbate dihedral rotation successfully applied (error:
720 7dd94df7 Carles Marti
            # {:.5f}°)".format((mol_dihedral_angle - mol_dihedral_angle_target
721 7dd94df7 Carles Marti
            # + 90) % 180 - 90))
722 7dd94df7 Carles Marti
        else:
723 7dd94df7 Carles Marti
            err = 'An unknown error occured during the adsorbate dihedral ' \
724 7dd94df7 Carles Marti
                  'rotation'
725 7dd94df7 Carles Marti
            logger.error(err)
726 7dd94df7 Carles Marti
            raise AssertionError(err)
727 7dd94df7 Carles Marti
728 7dd94df7 Carles Marti
729 7dd94df7 Carles Marti
def ads_chemcat(orig_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol, ctr1_surf,
730 7dd94df7 Carles Marti
                ctr2_surf, num_pts, min_coll_height, coll_coeff, norm_vect,
731 c25aa299 Carles Marti
                slab_nghbs, molec_nghbs, h_donor, h_acceptor, max_hel):
732 5261a07f Carles Marti
    """Generates adsorbate-surface structures by sampling over chemcat angles.
733 5261a07f Carles Marti

734 5261a07f Carles Marti
    @param orig_molec: ase.Atoms object of the molecule to adsorb (adsorbate).
735 5261a07f Carles Marti
    @param slab: ase.Atoms object of the surface on which to adsorb the molecule
736 5261a07f Carles Marti
    @param ctr1_mol: The index/es of the center in the adsorbate to use as
737 5261a07f Carles Marti
        adsorption center.
738 5261a07f Carles Marti
    @param ctr2_mol: The index/es of the center in the adsorbate to use for the
739 5261a07f Carles Marti
        definition of the surf-adsorbate angle, surf-adsorbate dihedral angle
740 5261a07f Carles Marti
        and adsorbate dihedral angle.
741 5261a07f Carles Marti
    @param ctr3_mol: The index/es of the center in the adsorbate to use for the
742 5261a07f Carles Marti
        definition of the adsorbate dihedral angle.
743 5261a07f Carles Marti
    @param ctr1_surf: The index/es of the center in the surface to use as
744 5261a07f Carles Marti
        adsorption center.
745 5261a07f Carles Marti
    @param ctr2_surf: The index/es of the center in the surface to use for the
746 5261a07f Carles Marti
        definition of the surf-adsorbate dihedral angle.
747 5261a07f Carles Marti
    @param num_pts: Number on which to sample Euler angles.
748 5261a07f Carles Marti
    @param min_coll_height: The lowest height for which to detect a collision
749 5261a07f Carles Marti
    @param coll_coeff: The coefficient that multiplies the covalent radius of
750 5261a07f Carles Marti
        atoms resulting in a distance that two atoms being closer to that is
751 5261a07f Carles Marti
        considered as atomic collision.
752 5261a07f Carles Marti
    @param norm_vect: The vector perpendicular to the surface.
753 5261a07f Carles Marti
    @param slab_nghbs: Number of neigbors in the surface.
754 5261a07f Carles Marti
    @param molec_nghbs: Number of neighbors in the molecule.
755 c25aa299 Carles Marti
    @param h_donor: List of atom types or atom numbers that are H-donors.
756 c25aa299 Carles Marti
    @param h_acceptor: List of atom types or atom numbers that are H-acceptors.
757 c25aa299 Carles Marti
    @param max_hel: Maximum value for sampling the helicopter
758 5261a07f Carles Marti
        (surf-adsorbate dihedral) angle.
759 5261a07f Carles Marti
    @return: list of ase.Atoms object conatining all the orientations of a given
760 5261a07f Carles Marti
        conformer.
761 5261a07f Carles Marti
    """
762 7dd94df7 Carles Marti
    from copy import deepcopy
763 7dd94df7 Carles Marti
    slab_ads_list = []
764 7dd94df7 Carles Marti
    # Rotation over bond angle
765 7dd94df7 Carles Marti
    for alpha in np.arange(90, 180, 90 / min(num_pts, 2)):
766 7dd94df7 Carles Marti
        # Rotation over surf-adsorbate dihedral angle
767 c25aa299 Carles Marti
        for beta in np.arange(0, max_hel, max_hel / num_pts):
768 7dd94df7 Carles Marti
            # Rotation over adsorbate bond dihedral angle
769 7dd94df7 Carles Marti
            for gamma in np.arange(0, 360, 360 / num_pts):
770 7dd94df7 Carles Marti
                new_molec = deepcopy(orig_molec)
771 7dd94df7 Carles Marti
                chemcat_rotate(new_molec, slab, ctr1_mol, ctr2_mol, ctr3_mol,
772 7dd94df7 Carles Marti
                               ctr1_surf, ctr2_surf, norm_vect, alpha,
773 7dd94df7 Carles Marti
                               beta, gamma)
774 7dd94df7 Carles Marti
                site_coords = get_atom_coords(slab, ctr1_surf)
775 7dd94df7 Carles Marti
                ctr_coords = get_atom_coords(new_molec, ctr1_mol)
776 7dd94df7 Carles Marti
                slab_molec, collision = correct_coll(new_molec, slab,
777 7dd94df7 Carles Marti
                                                     ctr_coords, site_coords,
778 7dd94df7 Carles Marti
                                                     num_pts, min_coll_height,
779 7dd94df7 Carles Marti
                                                     norm_vect, slab_nghbs,
780 e8bebcca Carles Marti
                                                     molec_nghbs, coll_coeff)
781 7dd94df7 Carles Marti
                if not collision and \
782 5261a07f Carles Marti
                        not any([np.allclose(slab_molec.positions,
783 5261a07f Carles Marti
                                             conf.positions)
784 7dd94df7 Carles Marti
                                 for conf in slab_ads_list]):
785 7dd94df7 Carles Marti
                    slab_ads_list.append(slab_molec)
786 c25aa299 Carles Marti
                    if h_donor is not False:
787 c25aa299 Carles Marti
                        slab_ads_list.extend(dissociation(slab_molec, h_donor,
788 c25aa299 Carles Marti
                                                          h_acceptor,
789 c25aa299 Carles Marti
                                                          len(slab)))
790 7dd94df7 Carles Marti
791 7dd94df7 Carles Marti
    return slab_ads_list
792 f3d1e601 Carles Marti
793 f3d1e601 Carles Marti
794 7dd94df7 Carles Marti
def adsorb_confs(conf_list, surf, inp_vars):
795 a5cc42ff Carles Marti
    """Generates a number of adsorbate-surface structure coordinates.
796 a5cc42ff Carles Marti

797 a5cc42ff Carles Marti
    Given a list of conformers, a surface, a list of atom indices (or list of
798 a5cc42ff Carles Marti
    list of indices) of both the surface and the adsorbate, it generates a
799 a5cc42ff Carles Marti
    number of adsorbate-surface structures for every possible combination of
800 a5cc42ff Carles Marti
    them at different orientations.
801 a5cc42ff Carles Marti
    @param conf_list: list of ase.Atoms of the different conformers
802 a5cc42ff Carles Marti
    @param surf: the ase.Atoms object of the surface
803 7dd94df7 Carles Marti
    @param inp_vars: Calculation parameters from input file.
804 a5cc42ff Carles Marti
    @return: list of ase.Atoms for the adsorbate-surface structures
805 a5cc42ff Carles Marti
    """
806 bb55f47c Carles Marti
    from ase.neighborlist import natural_cutoffs, neighbor_list
807 7dd94df7 Carles Marti
    molec_ctrs = inp_vars['molec_ctrs']
808 7dd94df7 Carles Marti
    sites = inp_vars['sites']
809 7dd94df7 Carles Marti
    angles = inp_vars['set_angles']
810 7dd94df7 Carles Marti
    num_pts = inp_vars['sample_points_per_angle']
811 7dd94df7 Carles Marti
    norm_vect = inp_vars['surf_norm_vect']
812 7dd94df7 Carles Marti
    min_coll_height = inp_vars['min_coll_height']
813 7dd94df7 Carles Marti
    coll_coeff = inp_vars['collision_threshold']
814 c25aa299 Carles Marti
    h_donor = inp_vars['h_donor']
815 c25aa299 Carles Marti
    h_acceptor = inp_vars['h_acceptor']
816 7dd94df7 Carles Marti
817 bc703cab Carles Marti
    if inp_vars['pbc_cell'] is not False:
818 bc703cab Carles Marti
        surf.set_pbc(True)
819 bc703cab Carles Marti
        surf.set_cell(inp_vars['pbc_cell'])
820 bc703cab Carles Marti
821 f3d1e601 Carles Marti
    surf_ads_list = []
822 f3d1e601 Carles Marti
    sites_coords = get_atom_coords(surf, sites)
823 e8bebcca Carles Marti
    if coll_coeff is not False:
824 5fb01677 Carles Marti
        surf_cutoffs = natural_cutoffs(surf, mult=coll_coeff)
825 5fb01677 Carles Marti
        surf_nghbs = len(neighbor_list("i", surf, surf_cutoffs))
826 5fb01677 Carles Marti
    else:
827 5fb01677 Carles Marti
        surf_nghbs = 0
828 7dd94df7 Carles Marti
    for i, conf in enumerate(conf_list):
829 bb55f47c Carles Marti
        molec_ctr_coords = get_atom_coords(conf, molec_ctrs)
830 bc703cab Carles Marti
        if inp_vars['pbc_cell'] is not False:
831 bc703cab Carles Marti
            conf.set_pbc(True)
832 bc703cab Carles Marti
            conf.set_cell(inp_vars['pbc_cell'])
833 e8bebcca Carles Marti
        if coll_coeff is not False:
834 5fb01677 Carles Marti
            conf_cutoffs = natural_cutoffs(conf, mult=coll_coeff)
835 5fb01677 Carles Marti
            molec_nghbs = len(neighbor_list("i", conf, conf_cutoffs))
836 5fb01677 Carles Marti
        else:
837 5fb01677 Carles Marti
            molec_nghbs = 0
838 7dd94df7 Carles Marti
        for s, site in enumerate(sites_coords):
839 7dd94df7 Carles Marti
            for c, ctr in enumerate(molec_ctr_coords):
840 7dd94df7 Carles Marti
                if angles == 'euler':
841 bb55f47c Carles Marti
                    surf_ads_list.extend(ads_euler(conf, surf, ctr, site,
842 5fb01677 Carles Marti
                                                   num_pts, min_coll_height,
843 bb55f47c Carles Marti
                                                   coll_coeff, norm_vect,
844 b4b2f307 Carles Marti
                                                   surf_nghbs, molec_nghbs,
845 c25aa299 Carles Marti
                                                   h_donor, h_acceptor))
846 7dd94df7 Carles Marti
                elif angles == 'chemcat':
847 7dd94df7 Carles Marti
                    mol_ctr1 = molec_ctrs[c]
848 7dd94df7 Carles Marti
                    mol_ctr2 = inp_vars["molec_ctrs2"][c]
849 7dd94df7 Carles Marti
                    mol_ctr3 = inp_vars["molec_ctrs3"][c]
850 7dd94df7 Carles Marti
                    surf_ctr1 = sites[s]
851 7dd94df7 Carles Marti
                    surf_ctr2 = inp_vars["surf_ctrs2"][s]
852 a98d4172 Carles Marti
                    max_h = inp_vars["max_helic_angle"]
853 7dd94df7 Carles Marti
                    surf_ads_list.extend(ads_chemcat(conf, surf, mol_ctr1,
854 7dd94df7 Carles Marti
                                                     mol_ctr2, mol_ctr3,
855 7dd94df7 Carles Marti
                                                     surf_ctr1, surf_ctr2,
856 7dd94df7 Carles Marti
                                                     num_pts, min_coll_height,
857 7dd94df7 Carles Marti
                                                     coll_coeff, norm_vect,
858 7dd94df7 Carles Marti
                                                     surf_nghbs, molec_nghbs,
859 c25aa299 Carles Marti
                                                     h_donor, h_acceptor,
860 c25aa299 Carles Marti
                                                     max_h))
861 f3d1e601 Carles Marti
    return surf_ads_list
862 f3d1e601 Carles Marti
863 f3d1e601 Carles Marti
864 4614bb6a Carles
def run_screening(inp_vars):
865 91ae8d86 Carles Marti
    """Carries out the screening of adsorbate structures on a surface.
866 e07c09eb Carles

867 e07c09eb Carles
    @param inp_vars: Calculation parameters from input file.
868 e07c09eb Carles
    """
869 e07c09eb Carles
    import os
870 57e3a8c7 Carles Marti
    import random
871 fd2384fc Carles Marti
    from modules.formats import collect_coords, adapt_format
872 cf8fe0e3 Carles Marti
    from modules.calculation import run_calc, check_finished_calcs
873 e07c09eb Carles
874 76f4ac19 Carles Marti
    logger.info('Carrying out procedures for the screening of adsorbate-surface'
875 76f4ac19 Carles Marti
                ' structures.')
876 e07c09eb Carles
    if not os.path.isdir("isolated"):
877 e07c09eb Carles
        err = "'isolated' directory not found. It is needed in order to carry "
878 e07c09eb Carles
        "out the screening of structures to be adsorbed"
879 e07c09eb Carles
        logger.error(err)
880 cf8fe0e3 Carles Marti
        raise FileNotFoundError(err)
881 e07c09eb Carles
882 cf8fe0e3 Carles Marti
    finished_calcs, unfinished_calcs = check_finished_calcs('isolated',
883 cf8fe0e3 Carles Marti
                                                            inp_vars['code'])
884 1a1164e0 Carles Marti
    logger.info(f"Found {len(finished_calcs)} structures of isolated "
885 1a1164e0 Carles Marti
                f"conformers whose calculation finished normally.")
886 1a1164e0 Carles Marti
    if len(unfinished_calcs) != 0:
887 cf8fe0e3 Carles Marti
        logger.warning(f"Found {len(unfinished_calcs)} calculations more that "
888 76f4ac19 Carles Marti
                       f"did not finish normally: {unfinished_calcs}. \n"
889 76f4ac19 Carles Marti
                       f"Using only the ones that finished normally: "
890 76f4ac19 Carles Marti
                       f"{finished_calcs}.")
891 cf8fe0e3 Carles Marti
892 fd2384fc Carles Marti
    conformer_atoms_list = collect_coords(finished_calcs, inp_vars['code'],
893 fd2384fc Carles Marti
                                          'isolated', inp_vars['special_atoms'])
894 1e9e784d Carles Marti
    selected_confs = select_confs(conformer_atoms_list, finished_calcs,
895 fd2384fc Carles Marti
                                  inp_vars['select_magns'],
896 bfe93f0d Carles Marti
                                  inp_vars['confs_per_magn'],
897 bfe93f0d Carles Marti
                                  inp_vars['code'])
898 90819cc3 Carles Marti
    surf = adapt_format('ase', inp_vars['surf_file'], inp_vars['special_atoms'])
899 7dd94df7 Carles Marti
    surf_ads_list = adsorb_confs(selected_confs, surf, inp_vars)
900 7d97341d Carles Marti
    if len(surf_ads_list) > inp_vars['max_structures']:
901 57e3a8c7 Carles Marti
        surf_ads_list = random.sample(surf_ads_list, inp_vars['max_structures'])
902 bfe93f0d Carles Marti
    logger.info(f'Generated {len(surf_ads_list)} adsorbate-surface atomic '
903 d9167fea Carles Marti
                f'configurations to carry out a calculation of.')
904 d9167fea Carles Marti
905 f3d1e601 Carles Marti
    run_calc('screening', inp_vars, surf_ads_list)
906 14f39d2a Carles Marti
    logger.info('Finished the procedures for the screening of adsorbate-surface'
907 14f39d2a Carles Marti
                ' structures section.')