Statistiques
| Révision :

root / ase / embed.py @ 15

Historique | Voir | Annoter | Télécharger (11,89 ko)

1
"""
2
This module is the EMBED module for ASE
3
implemented by T. Kerber
4

5
Torsten Kerber, ENS LYON: 2011, 07, 11
6

7
This work is supported by Award No. UK-C0017, made by King Abdullah
8
University of Science and Technology (KAUST)
9
"""
10

    
11
import math
12
from ase import Atom, Atoms
13
from ase.data import covalent_radii, atomic_numbers
14

    
15
import numpy as np
16

    
17
class Embed(Atoms):
18
    #--- constructor of the Embed class ---
19
    def __init__(self, system, cluster, cell_cluster = "Auto", cluster_pos = True):
20
        super(Embed, self).__init__()
21
        # define the atom map
22
        self.atom_map_sys_cl = []
23
        self.atom_map_cl_sys = []
24
        self.linkatoms = []
25
        # cluster dimensions
26
        self.xyz_cl_min = None
27
        self.xyz_cl_max = None
28
        # set the search radius for link atoms
29
        self.d = 10
30
        # define the systems for calculations
31
        if system is None or cluster is None:
32
            raise RuntimeError("Embed: system or cluster is not definied")
33
        self.set_system(system)
34
        if cluster_pos:
35
            self.set_cluster(cluster)
36
        else:
37
            self.set_cluster_by_numbers(cluster)
38
            
39
        # set the cell of the system
40
        self.set_cell(system.get_cell())
41
        self.cell_cluster = cell_cluster
42
        return
43

    
44
    #--- set the cluster ---
45
    def set_cluster(self, atoms):
46
        import copy
47
        # set the min/max cluster dimensions
48
        self.xyz_cl_min = atoms[0].get_position()
49
        self.xyz_cl_max = atoms[0].get_position()
50
        for atom in atoms:
51
            # assign the label "Cluster (10)" in atom.TAG
52
            atom.set_tag(10)
53
            xyz=atom.get_position()
54
            for i in xrange(3):
55
                # set the min/max cluster dimensions
56
                if xyz[i] < self.xyz_cl_min[i]:
57
                    self.xyz_cl_min[i] = xyz[i]
58
                if xyz[i] > self.xyz_cl_max[i]:
59
                    self.xyz_cl_max[i] = xyz[i]
60

    
61
        # add self.d around min/max cluster dimensions
62
        self.xyz_cl_min -= [self.d, self.d, self.d]
63
        self.xyz_cl_max += [self.d, self.d, self.d]
64
        # set the cluster for low and high level calculation
65
        self.atoms_cluster = atoms
66
        return
67
        
68
    #--- set cluster by atom numbers ---
69
    def set_cluster_by_numbers(self, numbers):
70
        cluster = Atoms()
71
        nat = len(self.atoms_system)
72
        for number in numbers:
73
            if nat > numbers:
74
                raise RuntimeError("QMX: The number of the cluster atom ", number, "is bigger than the number of atoms")
75
            cluster.append(self.atoms_system[number-1])
76
        self.set_cluster(cluster)
77

    
78
    #--- set the system ---
79
    def set_system(self, atoms):
80
        self.atoms_system = atoms
81
        # assign the label "Cluster (10)" in atom.TAG
82
        for atom in atoms:
83
            atom.set_tag(0)
84
        # update search radius for link atoms
85
        dx = 0
86
        for atom in atoms:
87
            r = covalent_radii[atom.get_atomic_number()]
88
            if (r > dx):
89
                dx = r
90
        self.d = dx * 2.1
91
        return
92

    
93
    #--- return cluster ---
94
    def get_cluster(self):
95
        return self.atoms_cluster
96

    
97
    def get_system(self):
98
        return self.atoms_system
99

    
100
    #--- Embedding ---
101
    def embed(self):
102
        # is the cluster and the host system definied ?
103
        if self.atoms_cluster is None or self.atoms_system is None:
104
            return
105
        self.find_cluster()
106
        self.set_linkatoms()
107
        print "link atoms found: ", len(self.linkatoms)
108
        if self.cell_cluster == "System":
109
            self.atoms_cluster.set_cell(self.atoms_system.get_cell())
110
        elif self.cell_cluster == "Auto":
111
            positions = self.atoms_cluster.get_positions()
112
            #find the biggest dimensions of the cluster in x,y,z direction
113
            l = 0
114
            for idir in xrange(3):
115
                l = max(l, positions[:, idir].max() - positions[:, idir].min())
116
            # calculate the box parameters (cluster + min 5 Ang)
117
            l = (math.floor(l/2.5)+1)*2.5 + 5.0
118
            # build cell
119
            cell = np.zeros((3, 3),  float)
120
            # apply cell parameters
121
            for idir in xrange(3):
122
                cell[idir, idir] = l
123
            # set parameters to cluster
124
            self.atoms_cluster.set_cell(cell)
125
            # print information on the screen
126
            print "length of box surrounding the cluster: ", 
127
            print l*10, 
128
            print "pm"
129
        else:
130
            self.atoms_cluster.set_cell(self.cell_cluster)
131

    
132
    def find_cluster(self):
133
        # set tolerance
134
        d = 0.001
135
        #atoms
136
        xyzs_cl=[]
137
        for atom_cl in self.atoms_cluster:
138
            xyzs_cl.append(atom_cl.get_position())
139
        xyzs_sys=[]
140
        for atom_sys in self.atoms_system:
141
            xyzs_sys.append(atom_sys.get_position())
142

    
143
        self.atom_map_sys_cl = np.zeros(len(self.atoms_system), int)
144
        self.atom_map_cl_sys = np.zeros(len(self.atoms_cluster), int)
145
        # loop over cluster atoms atom_sys
146
        for iat_sys in xrange(len(self.atoms_system)):
147
            # get the coordinates of the system atom atom_sys
148
            xyz_sys = xyzs_sys[iat_sys]
149
            # bSysOnly: no identical atom has been found
150
            bSysOnly = True
151
            # loop over system atoms atom_cl
152
            for iat_cl in xrange(len(self.atoms_cluster)):
153
                # difference vector between both atoms
154
                xyz_diff = np.abs(xyzs_sys[iat_sys]-xyzs_cl[iat_cl])
155
                # identical atoms
156
                if xyz_diff[0] < d and xyz_diff[1] < d and xyz_diff[2] < d:
157
                    # set tag (CLUSTER+HOST: 10) to atom_sys
158
                    self.atoms_system[iat_sys].set_tag(10)
159
                    # map the atom in the atom list
160
                    self.atom_map_sys_cl[iat_sys] = iat_cl
161
                    self.atom_map_cl_sys[iat_cl] = iat_sys
162
                    # atom has been identified
163
                    bSysOnly = False
164
                    break
165
            if bSysOnly:
166
                self.atom_map_sys_cl[iat_sys] = -1
167

    
168
    def set_linkatoms(self, tol=15., linkAtom=None, debug=False):
169
        # local copies of xyz coordinates to avoid massive copying of xyz objects
170
        xyzs_cl=[]
171
        for atom_cl in self.atoms_cluster:
172
            xyzs_cl.append(atom_cl.get_position())
173
        xyzs_sys=[]
174
        for atom_sys in self.atoms_system:
175
            xyzs_sys.append(atom_sys.get_position())
176
        # set the standard link atom
177
        if linkAtom is None:
178
            linkAtom ='H'
179
        # number of atoms in the cluster and the system
180
        nat_cl=len(self.atoms_cluster)
181
        nat_sys=len(self.atoms_system)
182
        # system has pbc?
183
        pbc = self.get_pbc()
184
        # set the bond table
185
        bonds = []
186
        # set the 27 cell_vec, starting with the (0,0,0) vector for the unit cell
187
        cells_L = [(0.0,  0.0,  0.0)]
188
        # get the cell vectors of the host system and build up a 3 by 3 supercell to search for neighbors in the surrounding
189
        cell = self.atoms_system.get_cell()
190
        if self.atoms_system.get_pbc().any():
191
            for ix in xrange(-1, 2):
192
                for iy in xrange(-1, 2):
193
                    for iz in xrange(-1, 2):
194
                        if ix == 0 and iy == 0 and iz == 0:
195
                            continue
196
                        cells_L.append(np.dot([ix, iy, iz], cell))
197
        # save the radius of system atoms
198
        rs_sys = []
199
        for atom in self.atoms_system:
200
            rs_sys.append(covalent_radii[atom.get_atomic_number()])
201
        # sum over cluster atoms (iat_cl)
202
        for iat_cl in xrange(nat_cl):
203
            # get the cluster atom
204
            atom_cl=self.atoms_cluster[iat_cl]
205
            # ignore link atoms
206
            if atom_cl.get_tag() == 50:
207
                continue
208
            # xyz coordinates and covalent radius of the cluster atom iat_cl
209
            xyz_cl = xyzs_cl[iat_cl]
210
            r_cl = covalent_radii[atom_cl.get_atomic_number()]
211

    
212
            # sum over system atoms (iat_sys)
213
            for iat_sys in xrange(nat_sys):
214
                # avoid cluster atoms
215
                if self.atoms_system[iat_sys].get_tag()==10:
216
                    continue
217
                # sum over all cell_L
218
                for cell_L in cells_L:
219
                    # xyz coordinates and covalent radius of the system atom iat_sys
220
                    xyz_sys = xyzs_sys[iat_sys]+cell_L
221
                    # go only in distance self.d around the cluster
222
                    lcont = True
223
                    for i in xrange(3):
224
                        if (xyz_sys[i] < self.xyz_cl_min[i] or
225
                            xyz_sys[i] > self.xyz_cl_max[i]):
226
                            lcont = False
227
                            break
228
                    if not lcont:
229
                        continue
230
                    # xyz coordinates and covalent radius of the system atom iat_sys
231
                    r_sys = rs_sys[iat_sys]
232
                    # diff vector
233
                    xyz_diff = xyz_sys - xyz_cl
234
                    # distance between the atoms
235
                    r = np.sqrt(np.dot(xyz_diff, xyz_diff))
236
                    # ratio of the distance to the sum of covalent radius
237
                    f = r / (r_cl + r_sys)
238
                    if debug:
239
                        print "Covalent radii = ",r_cl, r_sys
240
                        print "Distance ", f
241
                        print "tol = ",(1+tol/100.),(1-tol/100.),(1-2*tol/100.)
242
                    if f <= (1+tol/100.) and f >= (1-2*tol/100.):
243
                        s = cell_L, self.atom_map_cl_sys[iat_cl], iat_sys, r_cl
244
                        bonds.append(s)
245
                        break
246
                    if f <= (1-2*tol/100.):
247
                        raise RuntimeError("QMX: The cluster atom", iat_cl, " and the system atom", iat_sys, "came too close")
248

    
249
        r_h = covalent_radii[atomic_numbers[linkAtom]]
250
        for bond in bonds:
251
            cell_L, iat_cl_sys, iat_sys, r_cl = bond
252
            # assign the tags for the border atoms
253
            atom_sys.set_tag(1)
254
            atom_cl.set_tag(11)
255
            #difference vector for the link atom, scaling
256
            xyz_diff = xyzs_sys[iat_sys] + cell_L - xyzs_sys[iat_cl_sys]
257
            r = (r_cl + r_h)
258
            xyz_diff *=  r / np.sqrt(np.dot(xyz_diff, xyz_diff))
259
            # determine position of the link atom
260
            xyz_diff += xyzs_sys[iat_cl_sys]
261
            # create link atom
262
            atom = Atom(symbol=linkAtom, position=xyz_diff, tag=50)
263
            # add atom to cluster
264
            self.atoms_cluster.append(atom)
265
            # add atom to the linkatoms
266
            s = cell_L, iat_cl_sys, iat_sys, r, len(self.atoms_cluster)-1
267
            self.linkatoms.append(s)
268
        return
269

    
270
    def set_positions(self, positions_new):
271
        # number of atoms
272
        nat_sys=len(self.atoms_system)
273
        # go over all pairs of atoms
274
        for iat_sys in xrange(nat_sys):
275
            xyz = positions_new[iat_sys]
276
            self.atoms_system[iat_sys].set_position(xyz)
277
            iat_cl = self.atom_map_sys_cl[iat_sys]
278
            if iat_cl > -1:
279
                self.atoms_cluster[iat_cl].set_position(xyz)
280

    
281
        for cell_L, iat_cl_sys, iat_sys, r, iat in self.linkatoms:
282
            # determine position of the link atom
283
            xyz_cl = positions_new[iat_cl_sys]
284
            xyz = positions_new[iat_sys] - xyz_cl + cell_L
285
            xyz *= r / np.sqrt(np.dot(xyz, xyz))
286
            xyz += xyz_cl
287
            # update xyz coordinates of the cluster
288
            self.atoms_cluster[iat].set_position(xyz)
289

    
290
    def __getitem__(self, i):
291
        return self.atoms_system.__getitem__(i)
292

    
293
    def get_positions(self):
294
        return self.atoms_system.get_positions()
295

    
296
    def __add__(self, other):
297
        return self.atoms_system.__add__(other)
298

    
299
    def __delitem__(self, i):
300
        return self.atoms_system.__delitem__(i)
301

    
302
    def __len__(self):
303
        return self.atoms_system.__len__()
304
        
305
    def get_chemical_symbols(self, reduce=False):
306
        return self.atoms_system.get_chemical_symbols(reduce)