Statistiques
| Révision :

root / ase / embed.py @ 12

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

1
import math
2
from ase import Atom, Atoms
3
from ase.data import covalent_radii, atomic_numbers
4

    
5
import numpy as np
6

    
7
class Embed(Atoms):
8
    #--- constructor of the Embed class ---
9
    def __init__(self, system=None, cluster=None, cell_cluster = "Auto"):
10
        super(Embed, self).__init__()
11
        # define the atom map
12
        self.atom_map_sys_cl = []
13
        self.atom_map_cl_sys = []
14
        self.linkatoms = []
15
        # cluster dimensions
16
        self.xyz_cl_min = None
17
        self.xyz_cl_max = None
18
        # set the search radius for link atoms
19
        self.d = 10
20
        # define the systems for calculations
21
        self.set_system(system)
22
        self.set_cluster(cluster)
23
        # set the cell of the system
24
        self.set_cell(system.get_cell())
25
        self.cell_cluster = cell_cluster
26
        return
27

    
28
    #--- set the cluster ---
29
    def set_cluster(self, atoms):
30
        import copy
31
        # set the min/max cluster dimensions
32
        self.xyz_cl_min = atoms[0].get_position()
33
        self.xyz_cl_max = atoms[0].get_position()
34
        for atom in atoms:
35
            # assign the label "Cluster (10)" in atom.TAG
36
            atom.set_tag(10)
37
            xyz=atom.get_position()
38
            for i in xrange(3):
39
                # set the min/max cluster dimensions
40
                if xyz[i] < self.xyz_cl_min[i]:
41
                    self.xyz_cl_min[i] = xyz[i]
42
                if xyz[i] > self.xyz_cl_max[i]:
43
                    self.xyz_cl_max[i] = xyz[i]
44

    
45
        # add self.d around min/max cluster dimensions
46
        self.xyz_cl_min -= [self.d, self.d, self.d]
47
        self.xyz_cl_max += [self.d, self.d, self.d]
48
        # set the cluster for low and high level calculation
49
        self.atoms_cluster = atoms
50
        return
51

    
52
    #--- set the system ---
53
    def set_system(self, atoms):
54
        self.atoms_system = atoms
55
        # assign the label "Cluster (10)" in atom.TAG
56
        for atom in atoms:
57
            atom.set_tag(0)
58
        # update search radius for link atoms
59
        dx = 0
60
        for atom in atoms:
61
            r = covalent_radii[atom.get_atomic_number()]
62
            if (r > dx):
63
                dx = r
64
        self.d = dx * 2.1
65
        return
66

    
67
    #--- return cluster ---
68
    def get_cluster(self):
69
        return self.atoms_cluster
70

    
71
    def get_system(self):
72
        return self.atoms_system
73

    
74
    #--- Embedding ---
75
    def embed(self):
76
        # is the cluster and the host system definied ?
77
        if self.atoms_cluster is None or self.atoms_system is None:
78
            return
79
        self.find_cluster()
80
        self.set_linkatoms()
81
        print "link atoms found: ", len(self.linkatoms)
82
        if self.cell_cluster == "System":
83
            self.atoms_cluster.set_cell(self.atoms_system.get_cell())
84
        elif self.cell_cluster == "Auto":
85
            positions = self.atoms_cluster.get_positions()
86
            #find the biggest dimensions of the cluster in x,y,z direction
87
            l = 0
88
            for idir in xrange(3):
89
                l = max(l, positions[:, idir].max() - positions[:, idir].min())
90
            # calculate the box parameters (cluster + min 5 Ang)
91
            l = (math.floor(l/2.5)+1)*2.5 + 5.0
92
            # build cell
93
            cell = np.zeros((3, 3),  float)
94
            # apply cell parameters
95
            for idir in xrange(3):
96
                cell[idir, idir] = l
97
            # set parameters to cluster
98
            self.atoms_cluster.set_cell(cell)
99
            # print information on the screen
100
            print "length of box surrounding the cluster: ", 
101
            print l*10, 
102
            print "pm"
103
        else:
104
            self.atoms_cluster.set_cell(self.cell_cluster)
105

    
106
    def find_cluster(self):
107
        # set tolerance
108
        d = 0.001
109
        #atoms
110
        xyzs_cl=[]
111
        for atom_cl in self.atoms_cluster:
112
            xyzs_cl.append(atom_cl.get_position())
113
        xyzs_sys=[]
114
        for atom_sys in self.atoms_system:
115
            xyzs_sys.append(atom_sys.get_position())
116

    
117
        self.atom_map_sys_cl = np.zeros(len(self.atoms_system), int)
118
        self.atom_map_cl_sys = np.zeros(len(self.atoms_cluster), int)
119
        # loop over cluster atoms atom_sys
120
        for iat_sys in xrange(len(self.atoms_system)):
121
            # get the coordinates of the system atom atom_sys
122
            xyz_sys = xyzs_sys[iat_sys]
123
            # bSysOnly: no identical atom has been found
124
            bSysOnly = True
125
            # loop over system atoms atom_cl
126
            for iat_cl in xrange(len(self.atoms_cluster)):
127
                # difference vector between both atoms
128
                xyz_diff = np.abs(xyzs_sys[iat_sys]-xyzs_cl[iat_cl])
129
                # identical atoms
130
                if xyz_diff[0] < d and xyz_diff[1] < d and xyz_diff[2] < d:
131
                    # set tag (CLUSTER+HOST: 10) to atom_sys
132
                    self.atoms_system[iat_sys].set_tag(10)
133
                    # map the atom in the atom list
134
                    self.atom_map_sys_cl[iat_sys] = iat_cl
135
                    self.atom_map_cl_sys[iat_cl] = iat_sys
136
                    # atom has been identified
137
                    bSysOnly = False
138
                    break
139
            if bSysOnly:
140
                self.atom_map_sys_cl[iat_sys] = -1
141

    
142
    def set_linkatoms(self, tol=15., linkAtom=None, debug=False):
143
        # local copies of xyz coordinates to avoid massive copying of xyz objects
144
        xyzs_cl=[]
145
        for atom_cl in self.atoms_cluster:
146
            xyzs_cl.append(atom_cl.get_position())
147
        xyzs_sys=[]
148
        for atom_sys in self.atoms_system:
149
            xyzs_sys.append(atom_sys.get_position())
150
        # set the standard link atom
151
        if linkAtom is None:
152
            linkAtom ='H'
153
        # number of atoms in the cluster and the system
154
        nat_cl=len(self.atoms_cluster)
155
        nat_sys=len(self.atoms_system)
156
        # system has pbc?
157
        pbc = self.get_pbc()
158
        # set the bond table
159
        bonds = []
160
        # set the 27 cell_vec, starting with the (0,0,0) vector for the unit cell
161
        cells_L = [(0.0,  0.0,  0.0)]
162
        # get the cell vectors of the host system and build up a 3 by 3 supercell to search for neighbors in the surrounding
163
        cell = self.atoms_system.get_cell()
164
        if self.atoms_system.get_pbc().any():
165
            for ix in xrange(-1, 2):
166
                for iy in xrange(-1, 2):
167
                    for iz in xrange(-1, 2):
168
                        if ix == 0 and iy == 0 and iz == 0:
169
                            continue
170
                        cells_L.append(np.dot([ix, iy, iz], cell))
171
        # save the radius of system atoms
172
        rs_sys = []
173
        for atom in self.atoms_system:
174
            rs_sys.append(covalent_radii[atom.get_atomic_number()])
175
        # sum over cluster atoms (iat_cl)
176
        for iat_cl in xrange(nat_cl):
177
            # get the cluster atom
178
            atom_cl=self.atoms_cluster[iat_cl]
179
            # ignore link atoms
180
            if atom_cl.get_tag() == 50:
181
                continue
182
            # xyz coordinates and covalent radius of the cluster atom iat_cl
183
            xyz_cl = xyzs_cl[iat_cl]
184
            r_cl = covalent_radii[atom_cl.get_atomic_number()]
185

    
186
            # sum over system atoms (iat_sys)
187
            for iat_sys in xrange(nat_sys):
188
                # avoid cluster atoms
189
                if self.atoms_system[iat_sys].get_tag()==10:
190
                    continue
191
                # sum over all cell_L
192
                for cell_L in cells_L:
193
                    # xyz coordinates and covalent radius of the system atom iat_sys
194
                    xyz_sys = xyzs_sys[iat_sys]+cell_L
195
                    # go only in distance self.d around the cluster
196
                    lcont = True
197
                    for i in xrange(3):
198
                        if (xyz_sys[i] < self.xyz_cl_min[i] or
199
                            xyz_sys[i] > self.xyz_cl_max[i]):
200
                            lcont = False
201
                            break
202
                    if not lcont:
203
                        continue
204
                    # xyz coordinates and covalent radius of the system atom iat_sys
205
                    r_sys = rs_sys[iat_sys]
206
                    # diff vector
207
                    xyz_diff = xyz_sys - xyz_cl
208
                    # distance between the atoms
209
                    r = np.sqrt(np.dot(xyz_diff, xyz_diff))
210
                    # ratio of the distance to the sum of covalent radius
211
                    f = r / (r_cl + r_sys)
212
                    if debug:
213
                        print "Covalent radii = ",r_cl, r_sys
214
                        print "Distance ", f
215
                        print "tol = ",(1+tol/100.),(1-tol/100.),(1-2*tol/100.)
216
                    if f <= (1+tol/100.) and f >= (1-2*tol/100.):
217
                        s = cell_L, self.atom_map_cl_sys[iat_cl], iat_sys, r_cl
218
                        bonds.append(s)
219
                        break
220
                    if f <= (1-2*tol/100.):
221
                        raise RuntimeError("QMX: The cluster atom", iat_cl, " and the system atom", iat_sys, "came too close")
222

    
223
        r_h = covalent_radii[atomic_numbers[linkAtom]]
224
        for bond in bonds:
225
            cell_L, iat_cl_sys, iat_sys, r_cl = bond
226
            # assign the tags for the border atoms
227
            atom_sys.set_tag(1)
228
            atom_cl.set_tag(11)
229
            #difference vector for the link atom, scaling
230
            xyz_diff = xyzs_sys[iat_sys] + cell_L - xyzs_sys[iat_cl_sys]
231
            r = (r_cl + r_h)
232
            xyz_diff *=  r / np.sqrt(np.dot(xyz_diff, xyz_diff))
233
            # determine position of the link atom
234
            xyz_diff += xyzs_sys[iat_cl_sys]
235
            # create link atom
236
            atom = Atom(symbol=linkAtom, position=xyz_diff, tag=50)
237
            # add atom to cluster
238
            self.atoms_cluster.append(atom)
239
            # add atom to the linkatoms
240
            s = cell_L, iat_cl_sys, iat_sys, r, len(self.atoms_cluster)-1
241
            self.linkatoms.append(s)
242
        return
243

    
244
    def set_positions(self, positions_new):
245
        # number of atoms
246
        nat_sys=len(self.atoms_system)
247
        # go over all pairs of atoms
248
        for iat_sys in xrange(nat_sys):
249
            xyz = positions_new[iat_sys]
250
            self.atoms_system[iat_sys].set_position(xyz)
251
            iat_cl = self.atom_map_sys_cl[iat_sys]
252
            if iat_cl > -1:
253
                self.atoms_cluster[iat_cl].set_position(xyz)
254

    
255
        for cell_L, iat_cl_sys, iat_sys, r, iat in self.linkatoms:
256
            # determine position of the link atom
257
            xyz_cl = positions_new[iat_cl_sys]
258
            xyz = positions_new[iat_sys] - xyz_cl + cell_L
259
            xyz *= r / np.sqrt(np.dot(xyz, xyz))
260
            xyz += xyz_cl
261
            # update xyz coordinates of the cluster
262
            self.atoms_cluster[iat].set_position(xyz)
263

    
264
    def __getitem__(self, i):
265
        return self.atoms_system.__getitem__(i)
266

    
267
    def get_positions(self):
268
        return self.atoms_system.get_positions()
269

    
270
    def __add__(self, other):
271
        return self.atoms_system.__add__(other)
272

    
273
    def __delitem__(self, i):
274
        return self.atoms_system.__delitem__(i)
275

    
276
    def __len__(self):
277
        return self.atoms_system.__len__()
278
        
279
    def get_chemical_symbols(self, reduce=False):
280
        return self.atoms_system.get_chemical_symbols(reduce)