root / ase / embed.py @ 20
Historique | Voir | Annoter | Télécharger (11,69 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, linkType = 'H'): |
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 |
|
43 |
self.embed(linkType)
|
44 |
return
|
45 |
|
46 |
#--- set the cluster ---
|
47 |
def set_cluster(self, atoms): |
48 |
import copy |
49 |
# set the min/max cluster dimensions
|
50 |
self.xyz_cl_min = atoms[0].get_position() |
51 |
self.xyz_cl_max = atoms[0].get_position() |
52 |
for atom in atoms: |
53 |
# assign the label "Cluster (10)" in atom.TAG
|
54 |
atom.set_tag(10)
|
55 |
xyz=atom.get_position() |
56 |
for i in xrange(3): |
57 |
# set the min/max cluster dimensions
|
58 |
if xyz[i] < self.xyz_cl_min[i]: |
59 |
self.xyz_cl_min[i] = xyz[i]
|
60 |
if xyz[i] > self.xyz_cl_max[i]: |
61 |
self.xyz_cl_max[i] = xyz[i]
|
62 |
|
63 |
# add self.d around min/max cluster dimensions
|
64 |
self.xyz_cl_min -= [self.d, self.d, self.d] |
65 |
self.xyz_cl_max += [self.d, self.d, self.d] |
66 |
# set the cluster for low and high level calculation
|
67 |
self.atoms_cluster = atoms
|
68 |
return
|
69 |
|
70 |
#--- set cluster by atom numbers ---
|
71 |
def set_cluster_by_numbers(self, numbers): |
72 |
cluster = Atoms() |
73 |
nat = len(self.atoms_system) |
74 |
for number in numbers: |
75 |
if nat > numbers:
|
76 |
raise RuntimeError("QMX: The number of the cluster atom ", number, "is bigger than the number of atoms") |
77 |
cluster.append(self.atoms_system[number-1]) |
78 |
self.set_cluster(cluster)
|
79 |
|
80 |
#--- set the system ---
|
81 |
def set_system(self, atoms): |
82 |
self.atoms_system = atoms
|
83 |
# assign the label "Cluster (10)" in atom.TAG
|
84 |
for atom in atoms: |
85 |
atom.set_tag(0)
|
86 |
# update search radius for link atoms
|
87 |
dx = 0
|
88 |
for atom in atoms: |
89 |
r = covalent_radii[atom.get_atomic_number()] |
90 |
if (r > dx):
|
91 |
dx = r |
92 |
self.d = dx * 2.1 |
93 |
return
|
94 |
|
95 |
#--- return cluster ---
|
96 |
def get_cluster(self): |
97 |
return self.atoms_cluster |
98 |
|
99 |
def get_system(self): |
100 |
return self.atoms_system |
101 |
|
102 |
#--- Embedding ---
|
103 |
def embed(self, linkType): |
104 |
# is the cluster and the host system definied ?
|
105 |
if self.atoms_cluster is None or self.atoms_system is None: |
106 |
return
|
107 |
self.find_cluster()
|
108 |
self.set_linkatoms(linkType)
|
109 |
print "link atoms found: ", len(self.linkatoms) |
110 |
if self.cell_cluster == "System": |
111 |
self.atoms_cluster.set_cell(self.atoms_system.get_cell()) |
112 |
elif self.cell_cluster == "Auto": |
113 |
positions = self.atoms_cluster.get_positions()
|
114 |
# build cell
|
115 |
cell = np.zeros((3, 3), float) |
116 |
#for all directions
|
117 |
for idir in xrange(3): |
118 |
#find the biggest dimensions of the cluster in x,y,z direction
|
119 |
l = positions[:, idir].max() - positions[:, idir].min() |
120 |
# calculate the box parameters (2*cluster + 60 pm)
|
121 |
l = (math.floor((2*l+6)/2.5)+1)*2.5 |
122 |
# apply cell parameters
|
123 |
cell[idir, idir] = l |
124 |
# set parameters to cluster
|
125 |
self.atoms_cluster.set_cell(cell)
|
126 |
# print information on the screen
|
127 |
print "dimension of box surrounding the cluster: ", |
128 |
for idir in xrange(3): |
129 |
print cell[idir, idir]*10, |
130 |
print " pm" |
131 |
else:
|
132 |
self.atoms_cluster.set_cell(self.cell_cluster) |
133 |
|
134 |
def find_cluster(self): |
135 |
# set tolerance
|
136 |
d = 0.001
|
137 |
#atoms
|
138 |
xyzs_cl=[] |
139 |
for atom_cl in self.atoms_cluster: |
140 |
xyzs_cl.append(atom_cl.get_position()) |
141 |
xyzs_sys=[] |
142 |
for atom_sys in self.atoms_system: |
143 |
xyzs_sys.append(atom_sys.get_position()) |
144 |
|
145 |
self.atom_map_sys_cl = np.zeros(len(self.atoms_system), int) |
146 |
self.atom_map_cl_sys = np.zeros(len(self.atoms_cluster), int) |
147 |
# loop over cluster atoms atom_sys
|
148 |
for iat_sys in xrange(len(self.atoms_system)): |
149 |
# get the coordinates of the system atom atom_sys
|
150 |
xyz_sys = xyzs_sys[iat_sys] |
151 |
# bSysOnly: no identical atom has been found
|
152 |
bSysOnly = True
|
153 |
# loop over system atoms atom_cl
|
154 |
for iat_cl in xrange(len(self.atoms_cluster)): |
155 |
# difference vector between both atoms
|
156 |
xyz_diff = np.abs(xyzs_sys[iat_sys]-xyzs_cl[iat_cl]) |
157 |
# identical atoms
|
158 |
if xyz_diff[0] < d and xyz_diff[1] < d and xyz_diff[2] < d: |
159 |
# set tag (CLUSTER+HOST: 10) to atom_sys
|
160 |
self.atoms_system[iat_sys].set_tag(10) |
161 |
# map the atom in the atom list
|
162 |
self.atom_map_sys_cl[iat_sys] = iat_cl
|
163 |
self.atom_map_cl_sys[iat_cl] = iat_sys
|
164 |
# atom has been identified
|
165 |
bSysOnly = False
|
166 |
break
|
167 |
if bSysOnly:
|
168 |
self.atom_map_sys_cl[iat_sys] = -1 |
169 |
|
170 |
def set_linkatoms(self, linkType, tol=15.): |
171 |
# local copies of xyz coordinates to avoid massive copying of xyz objects
|
172 |
xyzs_cl=[] |
173 |
for atom_cl in self.atoms_cluster: |
174 |
xyzs_cl.append(atom_cl.get_position()) |
175 |
xyzs_sys=[] |
176 |
for atom_sys in self.atoms_system: |
177 |
xyzs_sys.append(atom_sys.get_position()) |
178 |
# number of atoms in the cluster and the system
|
179 |
nat_cl=len(self.atoms_cluster) |
180 |
nat_sys=len(self.atoms_system) |
181 |
# system has pbc?
|
182 |
pbc = self.get_pbc()
|
183 |
# set the bond table
|
184 |
bonds = [] |
185 |
# set the 27 cell_vec, starting with the (0,0,0) vector for the unit cell
|
186 |
cells_L = [(0.0, 0.0, 0.0)] |
187 |
# get the cell vectors of the host system and build up a 3 by 3 supercell to search for neighbors in the surrounding
|
188 |
cell = self.atoms_system.get_cell()
|
189 |
if self.atoms_system.get_pbc().any(): |
190 |
for ix in xrange(-1, 2): |
191 |
for iy in xrange(-1, 2): |
192 |
for iz in xrange(-1, 2): |
193 |
if ix == 0 and iy == 0 and iz == 0: |
194 |
continue
|
195 |
cells_L.append(np.dot([ix, iy, iz], cell)) |
196 |
# save the radius of system atoms
|
197 |
rs_sys = [] |
198 |
for atom in self.atoms_system: |
199 |
rs_sys.append(covalent_radii[atom.get_atomic_number()]) |
200 |
# sum over cluster atoms (iat_cl)
|
201 |
for iat_cl in xrange(nat_cl): |
202 |
# get the cluster atom
|
203 |
atom_cl=self.atoms_cluster[iat_cl]
|
204 |
# ignore link atoms
|
205 |
if atom_cl.get_tag() == 50: |
206 |
continue
|
207 |
# xyz coordinates and covalent radius of the cluster atom iat_cl
|
208 |
xyz_cl = xyzs_cl[iat_cl] |
209 |
r_cl = covalent_radii[atom_cl.get_atomic_number()] |
210 |
|
211 |
# sum over system atoms (iat_sys)
|
212 |
for iat_sys in xrange(nat_sys): |
213 |
# avoid cluster atoms
|
214 |
if self.atoms_system[iat_sys].get_tag()==10: |
215 |
continue
|
216 |
# sum over all cell_L
|
217 |
for cell_L in cells_L: |
218 |
# xyz coordinates and covalent radius of the system atom iat_sys
|
219 |
xyz_sys = xyzs_sys[iat_sys]+cell_L |
220 |
# go only in distance self.d around the cluster
|
221 |
lcont = True
|
222 |
for i in xrange(3): |
223 |
if (xyz_sys[i] < self.xyz_cl_min[i] or |
224 |
xyz_sys[i] > self.xyz_cl_max[i]):
|
225 |
lcont = False
|
226 |
break
|
227 |
if not lcont: |
228 |
continue
|
229 |
# xyz coordinates and covalent radius of the system atom iat_sys
|
230 |
r_sys = rs_sys[iat_sys] |
231 |
# diff vector
|
232 |
xyz_diff = xyz_sys - xyz_cl |
233 |
# distance between the atoms
|
234 |
r = np.sqrt(np.dot(xyz_diff, xyz_diff)) |
235 |
# ratio of the distance to the sum of covalent radius
|
236 |
f = r / (r_cl + r_sys) |
237 |
if f <= (1+tol/100.) and f >= (1-2*tol/100.): |
238 |
s = cell_L, self.atom_map_cl_sys[iat_cl], iat_sys, r_cl
|
239 |
bonds.append(s) |
240 |
break
|
241 |
if f <= (1-2*tol/100.): |
242 |
raise RuntimeError("QMX: The cluster atom", iat_cl, " and the system atom", iat_sys, "came too close") |
243 |
|
244 |
r_h = covalent_radii[atomic_numbers[linkType]] |
245 |
for bond in bonds: |
246 |
cell_L, iat_cl_sys, iat_sys, r_cl = bond |
247 |
# assign the tags for the border atoms
|
248 |
atom_sys.set_tag(1)
|
249 |
atom_cl.set_tag(11)
|
250 |
#difference vector for the link atom, scaling
|
251 |
xyz_diff = xyzs_sys[iat_sys] + cell_L - xyzs_sys[iat_cl_sys] |
252 |
r = (r_cl + r_h) |
253 |
xyz_diff *= r / np.sqrt(np.dot(xyz_diff, xyz_diff)) |
254 |
# determine position of the link atom
|
255 |
xyz_diff += xyzs_sys[iat_cl_sys] |
256 |
# create link atom
|
257 |
atom = Atom(symbol=linkType, position=xyz_diff, tag=50)
|
258 |
# add atom to cluster
|
259 |
self.atoms_cluster.append(atom)
|
260 |
# add atom to the linkatoms
|
261 |
s = cell_L, iat_cl_sys, iat_sys, r, len(self.atoms_cluster)-1 |
262 |
self.linkatoms.append(s)
|
263 |
return
|
264 |
|
265 |
def set_positions(self, positions_new): |
266 |
# number of atoms
|
267 |
nat_sys=len(self.atoms_system) |
268 |
# go over all pairs of atoms
|
269 |
for iat_sys in xrange(nat_sys): |
270 |
xyz = positions_new[iat_sys] |
271 |
self.atoms_system[iat_sys].set_position(xyz)
|
272 |
iat_cl = self.atom_map_sys_cl[iat_sys]
|
273 |
if iat_cl > -1: |
274 |
self.atoms_cluster[iat_cl].set_position(xyz)
|
275 |
|
276 |
for cell_L, iat_cl_sys, iat_sys, r, iat in self.linkatoms: |
277 |
# determine position of the link atom
|
278 |
xyz_cl = positions_new[iat_cl_sys] |
279 |
xyz = positions_new[iat_sys] - xyz_cl + cell_L |
280 |
xyz *= r / np.sqrt(np.dot(xyz, xyz)) |
281 |
xyz += xyz_cl |
282 |
# update xyz coordinates of the cluster
|
283 |
self.atoms_cluster[iat].set_position(xyz)
|
284 |
|
285 |
def __getitem__(self, i): |
286 |
return self.atoms_system.__getitem__(i) |
287 |
|
288 |
def get_positions(self): |
289 |
return self.atoms_system.get_positions() |
290 |
|
291 |
def __add__(self, other): |
292 |
return self.atoms_system.__add__(other) |
293 |
|
294 |
def __delitem__(self, i): |
295 |
return self.atoms_system.__delitem__(i) |
296 |
|
297 |
def __len__(self): |
298 |
return self.atoms_system.__len__() |
299 |
|
300 |
def get_chemical_symbols(self, reduce=False): |
301 |
return self.atoms_system.get_chemical_symbols(reduce) |