root / ase / constraints.py
Historique | Voir | Annoter | Télécharger (21,74 ko)
1 | 1 | tkerber | from math import sqrt |
---|---|---|---|
2 | 1 | tkerber | |
3 | 1 | tkerber | import numpy as np |
4 | 1 | tkerber | |
5 | 1 | tkerber | __all__ = ['FixCartesian', 'FixBondLength', 'FixedMode', 'FixConstraintSingle', |
6 | 1 | tkerber | 'FixAtoms', 'UnitCellFilter', 'FixScaled', 'StrainFilter', |
7 | 1 | tkerber | 'FixedPlane', 'Filter', 'FixConstraint', 'FixedLine', |
8 | 1 | tkerber | 'FixBondLengths']
|
9 | 1 | tkerber | |
10 | 1 | tkerber | |
11 | 1 | tkerber | def slice2enlist(s): |
12 | 1 | tkerber | """Convert a slice object into a list of (new, old) tuples."""
|
13 | 1 | tkerber | if isinstance(s, (list, tuple)): |
14 | 1 | tkerber | return enumerate(s) |
15 | 1 | tkerber | if s.step == None: |
16 | 1 | tkerber | step = 1
|
17 | 1 | tkerber | else:
|
18 | 1 | tkerber | step = s.step |
19 | 1 | tkerber | if s.start == None: |
20 | 1 | tkerber | start = 0
|
21 | 1 | tkerber | else:
|
22 | 1 | tkerber | start = s.start |
23 | 1 | tkerber | return enumerate(range(start, s.stop, step)) |
24 | 1 | tkerber | |
25 | 1 | tkerber | class FixConstraint: |
26 | 1 | tkerber | """Base class for classes that fix one or more atoms in some way."""
|
27 | 1 | tkerber | |
28 | 1 | tkerber | def index_shuffle(self, ind): |
29 | 1 | tkerber | """Change the indices.
|
30 | 1 | tkerber |
|
31 | 1 | tkerber | When the ordering of the atoms in the Atoms object changes,
|
32 | 1 | tkerber | this method can be called to shuffle the indices of the
|
33 | 1 | tkerber | constraints.
|
34 | 1 | tkerber |
|
35 | 1 | tkerber | ind -- List or tuple of indices.
|
36 | 1 | tkerber |
|
37 | 1 | tkerber | """
|
38 | 1 | tkerber | raise NotImplementedError |
39 | 1 | tkerber | |
40 | 1 | tkerber | def repeat(self, m, n): |
41 | 1 | tkerber | """ basic method to multiply by m, needs to know the length
|
42 | 1 | tkerber | of the underlying atoms object for the assignment of
|
43 | 1 | tkerber | multiplied constraints to work.
|
44 | 1 | tkerber | """
|
45 | 1 | tkerber | raise NotImplementedError |
46 | 1 | tkerber | |
47 | 1 | tkerber | class FixConstraintSingle(FixConstraint): |
48 | 1 | tkerber | "Base class for classes that fix a single atom."
|
49 | 1 | tkerber | |
50 | 1 | tkerber | def index_shuffle(self, ind): |
51 | 1 | tkerber | "The atom index must be stored as self.a."
|
52 | 1 | tkerber | newa = -1 # Signal error |
53 | 1 | tkerber | for new, old in slice2enlist(ind): |
54 | 1 | tkerber | if old == self.a: |
55 | 1 | tkerber | newa = new |
56 | 1 | tkerber | break
|
57 | 1 | tkerber | if newa == -1: |
58 | 1 | tkerber | raise IndexError('Constraint not part of slice') |
59 | 1 | tkerber | self.a = newa
|
60 | 1 | tkerber | |
61 | 1 | tkerber | class FixAtoms(FixConstraint): |
62 | 1 | tkerber | """Constraint object for fixing some chosen atoms."""
|
63 | 1 | tkerber | def __init__(self, indices=None, mask=None): |
64 | 1 | tkerber | """Constrain chosen atoms.
|
65 | 1 | tkerber |
|
66 | 1 | tkerber | Parameters
|
67 | 1 | tkerber | ----------
|
68 | 1 | tkerber | indices : list of int
|
69 | 1 | tkerber | Indices for those atoms that should be constrained.
|
70 | 1 | tkerber | mask : list of bool
|
71 | 1 | tkerber | One boolean per atom indicating if the atom should be
|
72 | 1 | tkerber | constrained or not.
|
73 | 1 | tkerber |
|
74 | 1 | tkerber | Examples
|
75 | 1 | tkerber | --------
|
76 | 1 | tkerber | Fix all Copper atoms:
|
77 | 1 | tkerber |
|
78 | 1 | tkerber | >>> c = FixAtoms(mask=[s == 'Cu' for s in atoms.get_chemical_symbols()])
|
79 | 1 | tkerber | >>> atoms.set_constraint(c)
|
80 | 1 | tkerber |
|
81 | 1 | tkerber | Fix all atoms with z-coordinate less than 1.0 Angstrom:
|
82 | 1 | tkerber |
|
83 | 1 | tkerber | >>> c = FixAtoms(mask=atoms.positions[:, 2] < 1.0)
|
84 | 1 | tkerber | >>> atoms.set_constraint(c)
|
85 | 1 | tkerber | """
|
86 | 1 | tkerber | |
87 | 1 | tkerber | if indices is None and mask is None: |
88 | 1 | tkerber | raise ValueError('Use "indices" or "mask".') |
89 | 1 | tkerber | if indices is not None and mask is not None: |
90 | 1 | tkerber | raise ValueError('Use only one of "indices" and "mask".') |
91 | 1 | tkerber | |
92 | 1 | tkerber | if mask is not None: |
93 | 1 | tkerber | self.index = np.asarray(mask, bool) |
94 | 1 | tkerber | else:
|
95 | 1 | tkerber | # Check for duplicates
|
96 | 1 | tkerber | srt = np.sort(indices) |
97 | 1 | tkerber | for i in range(len(indices)-1): |
98 | 1 | tkerber | if srt[i] == srt[i+1]: |
99 | 1 | tkerber | raise ValueError("FixAtoms: The indices array contained duplicates. Perhaps you wanted to specify a mask instead, but forgot the mask= keyword.") |
100 | 1 | tkerber | self.index = np.asarray(indices, int) |
101 | 1 | tkerber | |
102 | 1 | tkerber | if self.index.ndim != 1: |
103 | 1 | tkerber | raise ValueError('Wrong argument to FixAtoms class!') |
104 | 1 | tkerber | |
105 | 1 | tkerber | def adjust_positions(self, old, new): |
106 | 1 | tkerber | new[self.index] = old[self.index] |
107 | 1 | tkerber | |
108 | 1 | tkerber | def adjust_forces(self, positions, forces): |
109 | 1 | tkerber | forces[self.index] = 0.0 |
110 | 1 | tkerber | |
111 | 1 | tkerber | def index_shuffle(self, ind): |
112 | 1 | tkerber | # See docstring of superclass
|
113 | 1 | tkerber | if self.index.dtype == bool: |
114 | 1 | tkerber | self.index = self.index[ind] |
115 | 1 | tkerber | else:
|
116 | 1 | tkerber | index = [] |
117 | 1 | tkerber | for new, old in slice2enlist(ind): |
118 | 1 | tkerber | if old in self.index: |
119 | 1 | tkerber | index.append(new) |
120 | 1 | tkerber | if len(index) == 0: |
121 | 1 | tkerber | raise IndexError('All indices in FixAtoms not part of slice') |
122 | 1 | tkerber | self.index = np.asarray(index, int) |
123 | 1 | tkerber | |
124 | 1 | tkerber | def copy(self): |
125 | 1 | tkerber | if self.index.dtype == bool: |
126 | 1 | tkerber | return FixAtoms(mask=self.index.copy()) |
127 | 1 | tkerber | else:
|
128 | 1 | tkerber | return FixAtoms(indices=self.index.copy()) |
129 | 1 | tkerber | |
130 | 1 | tkerber | def __repr__(self): |
131 | 1 | tkerber | if self.index.dtype == bool: |
132 | 1 | tkerber | return 'FixAtoms(mask=%s)' % ints2string(self.index.astype(int)) |
133 | 1 | tkerber | return 'FixAtoms(indices=%s)' % ints2string(self.index) |
134 | 1 | tkerber | |
135 | 1 | tkerber | def repeat(self, m, n): |
136 | 1 | tkerber | i0 = 0
|
137 | 1 | tkerber | l = len(self.index) |
138 | 1 | tkerber | natoms = 0
|
139 | 1 | tkerber | if isinstance(m, int): |
140 | 1 | tkerber | m = (m, m, m) |
141 | 1 | tkerber | index_new = [] |
142 | 1 | tkerber | for m2 in range(m[2]): |
143 | 1 | tkerber | for m1 in range(m[1]): |
144 | 1 | tkerber | for m0 in range(m[0]): |
145 | 1 | tkerber | i1 = i0 + n |
146 | 1 | tkerber | if self.index.dtype == bool: |
147 | 1 | tkerber | index_new += self.index
|
148 | 1 | tkerber | else:
|
149 | 1 | tkerber | index_new += [i+natoms for i in self.index] |
150 | 1 | tkerber | i0 = i1 |
151 | 1 | tkerber | natoms += n |
152 | 1 | tkerber | if self.index.dtype == bool: |
153 | 1 | tkerber | self.index = np.asarray(index_new, bool) |
154 | 1 | tkerber | else:
|
155 | 1 | tkerber | self.index = np.asarray(index_new, int) |
156 | 1 | tkerber | return self |
157 | 1 | tkerber | |
158 | 1 | tkerber | def ints2string(x, threshold=10): |
159 | 1 | tkerber | """Convert ndarray of ints to string."""
|
160 | 1 | tkerber | if len(x) <= threshold: |
161 | 1 | tkerber | return str(x.tolist()) |
162 | 1 | tkerber | return str(x[:threshold].tolist())[:-1] + ', ...]' |
163 | 1 | tkerber | |
164 | 1 | tkerber | class FixBondLengths(FixConstraint): |
165 | 1 | tkerber | def __init__(self, pairs, iterations=10): |
166 | 1 | tkerber | self.constraints = [FixBondLength(a1, a2)
|
167 | 1 | tkerber | for a1, a2 in pairs] |
168 | 1 | tkerber | self.iterations = iterations
|
169 | 1 | tkerber | |
170 | 1 | tkerber | def adjust_positions(self, old, new): |
171 | 1 | tkerber | for i in range(self.iterations): |
172 | 1 | tkerber | for constraint in self.constraints: |
173 | 1 | tkerber | constraint.adjust_positions(old, new) |
174 | 1 | tkerber | |
175 | 1 | tkerber | def adjust_forces(self, positions, forces): |
176 | 1 | tkerber | for i in range(self.iterations): |
177 | 1 | tkerber | for constraint in self.constraints: |
178 | 1 | tkerber | constraint.adjust_forces(positions, forces) |
179 | 1 | tkerber | |
180 | 1 | tkerber | def copy(self): |
181 | 1 | tkerber | return FixBondLengths([constraint.indices for constraint in self.constraints]) |
182 | 1 | tkerber | |
183 | 1 | tkerber | class FixBondLength(FixConstraint): |
184 | 1 | tkerber | """Constraint object for fixing a bond length."""
|
185 | 1 | tkerber | def __init__(self, a1, a2): |
186 | 1 | tkerber | """Fix distance between atoms with indices a1 and a2."""
|
187 | 1 | tkerber | self.indices = [a1, a2]
|
188 | 1 | tkerber | |
189 | 1 | tkerber | def adjust_positions(self, old, new): |
190 | 1 | tkerber | p1, p2 = old[self.indices]
|
191 | 1 | tkerber | d = p2 - p1 |
192 | 1 | tkerber | p = sqrt(np.dot(d, d)) |
193 | 1 | tkerber | q1, q2 = new[self.indices]
|
194 | 1 | tkerber | d = q2 - q1 |
195 | 1 | tkerber | q = sqrt(np.dot(d, d)) |
196 | 1 | tkerber | d *= 0.5 * (p - q) / q
|
197 | 1 | tkerber | new[self.indices] = (q1 - d, q2 + d)
|
198 | 1 | tkerber | |
199 | 1 | tkerber | def adjust_forces(self, positions, forces): |
200 | 1 | tkerber | d = np.subtract.reduce(positions[self.indices])
|
201 | 1 | tkerber | d2 = np.dot(d, d) |
202 | 1 | tkerber | d *= 0.5 * np.dot(np.subtract.reduce(forces[self.indices]), d) / d2 |
203 | 1 | tkerber | forces[self.indices] += (-d, d)
|
204 | 1 | tkerber | |
205 | 1 | tkerber | def index_shuffle(self, ind): |
206 | 1 | tkerber | 'Shuffle the indices of the two atoms in this constraint'
|
207 | 1 | tkerber | newa = [-1, -1] # Signal error |
208 | 1 | tkerber | for new, old in slice2enlist(ind): |
209 | 1 | tkerber | for i, a in enumerate(self.indices): |
210 | 1 | tkerber | if old == a:
|
211 | 1 | tkerber | newa[i] = new |
212 | 1 | tkerber | if newa[0] == -1 or newa[1] == -1: |
213 | 1 | tkerber | raise IndexError('Constraint not part of slice') |
214 | 1 | tkerber | self.indices = newa
|
215 | 1 | tkerber | |
216 | 1 | tkerber | def copy(self): |
217 | 1 | tkerber | return FixBondLength(*self.indices) |
218 | 1 | tkerber | |
219 | 1 | tkerber | def __repr__(self): |
220 | 1 | tkerber | return 'FixBondLength(%d, %d)' % tuple(self.indices) |
221 | 1 | tkerber | |
222 | 1 | tkerber | class FixedMode(FixConstraint): |
223 | 1 | tkerber | """Constrain atoms to move along directions orthogonal to
|
224 | 1 | tkerber | a given mode only."""
|
225 | 1 | tkerber | |
226 | 1 | tkerber | def __init__(self,indices,mode): |
227 | 1 | tkerber | if indices is None: |
228 | 1 | tkerber | raise ValueError('Use "indices".') |
229 | 1 | tkerber | if indices is not None: |
230 | 1 | tkerber | self.index = np.asarray(indices, int) |
231 | 1 | tkerber | self.mode = (np.asarray(mode) / np.sqrt((mode **2).sum())).reshape(-1) |
232 | 1 | tkerber | |
233 | 1 | tkerber | def adjust_positions(self, oldpositions, newpositions): |
234 | 1 | tkerber | newpositions = newpositions.ravel() |
235 | 1 | tkerber | oldpositions = oldpositions.ravel() |
236 | 1 | tkerber | step = newpositions - oldpositions |
237 | 1 | tkerber | newpositions -= self.mode * np.dot(step, self.mode) |
238 | 1 | tkerber | newpositions = newpositions.reshape(-1,3) |
239 | 1 | tkerber | oldpositions = oldpositions.reshape(-1,3) |
240 | 1 | tkerber | |
241 | 1 | tkerber | def adjust_forces(self, positions, forces): |
242 | 1 | tkerber | forces = forces.ravel() |
243 | 1 | tkerber | forces -= self.mode * np.dot(forces, self.mode) |
244 | 1 | tkerber | forces = forces.reshape(-1 ,3) |
245 | 1 | tkerber | |
246 | 1 | tkerber | def copy(self): |
247 | 1 | tkerber | return FixedMode(self.index.copy(), self.mode) |
248 | 1 | tkerber | |
249 | 1 | tkerber | def __repr__(self): |
250 | 1 | tkerber | return 'FixedMode(%d, %s)' % (ints2string(self.index), self.mode.tolist()) |
251 | 1 | tkerber | |
252 | 1 | tkerber | class FixedPlane(FixConstraintSingle): |
253 | 1 | tkerber | """Constrain an atom *a* to move in a given plane only.
|
254 | 1 | tkerber |
|
255 | 1 | tkerber | The plane is defined by its normal: *direction*."""
|
256 | 1 | tkerber | |
257 | 1 | tkerber | def __init__(self, a, direction): |
258 | 1 | tkerber | self.a = a
|
259 | 1 | tkerber | self.dir = np.asarray(direction) / sqrt(np.dot(direction, direction))
|
260 | 1 | tkerber | |
261 | 1 | tkerber | def adjust_positions(self, oldpositions, newpositions): |
262 | 1 | tkerber | step = newpositions[self.a] - oldpositions[self.a] |
263 | 1 | tkerber | newpositions[self.a] -= self.dir * np.dot(step, self.dir) |
264 | 1 | tkerber | |
265 | 1 | tkerber | def adjust_forces(self, positions, forces): |
266 | 1 | tkerber | forces[self.a] -= self.dir * np.dot(forces[self.a], self.dir) |
267 | 1 | tkerber | |
268 | 1 | tkerber | def copy(self): |
269 | 1 | tkerber | return FixedPlane(self.a, self.dir) |
270 | 1 | tkerber | |
271 | 1 | tkerber | def __repr__(self): |
272 | 1 | tkerber | return 'FixedPlane(%d, %s)' % (self.a, self.dir.tolist()) |
273 | 1 | tkerber | |
274 | 1 | tkerber | |
275 | 1 | tkerber | class FixedLine(FixConstraintSingle): |
276 | 1 | tkerber | """Constrain an atom *a* to move on a given line only.
|
277 | 1 | tkerber |
|
278 | 1 | tkerber | The line is defined by its *direction*."""
|
279 | 1 | tkerber | |
280 | 1 | tkerber | def __init__(self, a, direction): |
281 | 1 | tkerber | self.a = a
|
282 | 1 | tkerber | self.dir = np.asarray(direction) / sqrt(np.dot(direction, direction))
|
283 | 1 | tkerber | |
284 | 1 | tkerber | def adjust_positions(self, oldpositions, newpositions): |
285 | 1 | tkerber | step = newpositions[self.a] - oldpositions[self.a] |
286 | 1 | tkerber | x = np.dot(step, self.dir)
|
287 | 1 | tkerber | newpositions[self.a] = oldpositions[self.a] + x * self.dir |
288 | 1 | tkerber | |
289 | 1 | tkerber | def adjust_forces(self, positions, forces): |
290 | 1 | tkerber | forces[self.a] = self.dir * np.dot(forces[self.a], self.dir) |
291 | 1 | tkerber | |
292 | 1 | tkerber | def copy(self): |
293 | 1 | tkerber | return FixedLine(self.a, self.dir) |
294 | 1 | tkerber | |
295 | 1 | tkerber | def __repr__(self): |
296 | 1 | tkerber | return 'FixedLine(%d, %s)' % (self.a, self.dir.tolist()) |
297 | 1 | tkerber | |
298 | 1 | tkerber | class FixCartesian(FixConstraintSingle): |
299 | 1 | tkerber | "Fix an atom in the directions of the cartesian coordinates."
|
300 | 1 | tkerber | def __init__(self, a, mask=[1,1,1]): |
301 | 1 | tkerber | self.a=a
|
302 | 1 | tkerber | self.mask = -(np.array(mask)-1) |
303 | 1 | tkerber | |
304 | 1 | tkerber | def adjust_positions(self, old, new): |
305 | 1 | tkerber | step = new - old |
306 | 1 | tkerber | step[self.a] *= self.mask |
307 | 1 | tkerber | new = old + step |
308 | 1 | tkerber | |
309 | 1 | tkerber | def adjust_forces(self, positions, forces): |
310 | 1 | tkerber | forces[self.a] *= self.mask |
311 | 1 | tkerber | |
312 | 1 | tkerber | def copy(self): |
313 | 1 | tkerber | return FixCartesian(self.a, self.mask) |
314 | 1 | tkerber | |
315 | 1 | tkerber | def __repr__(self): |
316 | 1 | tkerber | return 'FixCartesian(indice=%s mask=%s)' % (self.a, self.mask) |
317 | 1 | tkerber | |
318 | 1 | tkerber | class fix_cartesian(FixCartesian): |
319 | 1 | tkerber | "Backwards compatibility for FixCartesian."
|
320 | 1 | tkerber | def __init__(self, a, mask=[1,1,1]): |
321 | 1 | tkerber | import warnings |
322 | 1 | tkerber | super(fix_cartesian, self).__init__(a, mask) |
323 | 1 | tkerber | warnings.warn('fix_cartesian is deprecated. Please use FixCartesian' \
|
324 | 1 | tkerber | ' instead.', DeprecationWarning, stacklevel=2) |
325 | 1 | tkerber | |
326 | 1 | tkerber | class FixScaled(FixConstraintSingle): |
327 | 1 | tkerber | "Fix an atom in the directions of the unit vectors."
|
328 | 1 | tkerber | def __init__(self, cell, a, mask=[1,1,1]): |
329 | 1 | tkerber | self.cell = cell
|
330 | 1 | tkerber | self.a = a
|
331 | 1 | tkerber | self.mask = np.array(mask)
|
332 | 1 | tkerber | |
333 | 1 | tkerber | def adjust_positions(self, old, new): |
334 | 1 | tkerber | scaled_old = np.linalg.solve(self.cell.T, old.T).T
|
335 | 1 | tkerber | scaled_new = np.linalg.solve(self.cell.T, new.T).T
|
336 | 1 | tkerber | for n in range(3): |
337 | 1 | tkerber | if self.mask[n]: |
338 | 1 | tkerber | scaled_new[self.a, n] = scaled_old[self.a, n] |
339 | 1 | tkerber | new[self.a] = np.dot(scaled_new, self.cell)[self.a] |
340 | 1 | tkerber | |
341 | 1 | tkerber | def adjust_forces(self, positions, forces): |
342 | 1 | tkerber | scaled_forces = np.linalg.solve(self.cell.T, forces.T).T
|
343 | 1 | tkerber | scaled_forces[self.a] *= -(self.mask-1) |
344 | 1 | tkerber | forces[self.a] = np.dot(scaled_forces, self.cell)[self.a] |
345 | 1 | tkerber | |
346 | 1 | tkerber | def copy(self): |
347 | 1 | tkerber | return FixScaled(self.cell ,self.a, self.mask) |
348 | 1 | tkerber | |
349 | 1 | tkerber | def __repr__(self): |
350 | 1 | tkerber | return 'FixScaled(indice=%s mask=%s)' % (self.a, self.mask) |
351 | 1 | tkerber | |
352 | 1 | tkerber | class fix_scaled(FixScaled): |
353 | 1 | tkerber | "Backwards compatibility for FixScaled."
|
354 | 1 | tkerber | def __init__(self, cell, a, mask=[1,1,1]): |
355 | 1 | tkerber | import warnings |
356 | 1 | tkerber | super(fix_scaled, self).__init__(cell, a, mask) |
357 | 1 | tkerber | warnings.warn('fix_scaled is deprecated. Please use FixScaled ' \
|
358 | 1 | tkerber | 'instead.', DeprecationWarning, stacklevel=2) |
359 | 1 | tkerber | |
360 | 1 | tkerber | class Filter: |
361 | 1 | tkerber | """Subset filter class."""
|
362 | 1 | tkerber | def __init__(self, atoms, indices=None, mask=None): |
363 | 1 | tkerber | """Filter atoms.
|
364 | 1 | tkerber |
|
365 | 1 | tkerber | This filter can be used to hide degrees of freedom in an Atoms
|
366 | 1 | tkerber | object.
|
367 | 1 | tkerber |
|
368 | 1 | tkerber | Parameters
|
369 | 1 | tkerber | ----------
|
370 | 1 | tkerber | indices : list of int
|
371 | 1 | tkerber | Indices for those atoms that should remain visible.
|
372 | 1 | tkerber | mask : list of bool
|
373 | 1 | tkerber | One boolean per atom indicating if the atom should remain
|
374 | 1 | tkerber | visible or not.
|
375 | 1 | tkerber | """
|
376 | 1 | tkerber | |
377 | 1 | tkerber | self.atoms = atoms
|
378 | 1 | tkerber | self.constraints = []
|
379 | 1 | tkerber | |
380 | 1 | tkerber | if indices is None and mask is None: |
381 | 1 | tkerber | raise ValuError('Use "indices" or "mask".') |
382 | 1 | tkerber | if indices is not None and mask is not None: |
383 | 1 | tkerber | raise ValuError('Use only one of "indices" and "mask".') |
384 | 1 | tkerber | |
385 | 1 | tkerber | if mask is not None: |
386 | 1 | tkerber | self.index = np.asarray(mask, bool) |
387 | 1 | tkerber | self.n = self.index.sum() |
388 | 1 | tkerber | else:
|
389 | 1 | tkerber | self.index = np.asarray(indices, int) |
390 | 1 | tkerber | self.n = len(self.index) |
391 | 1 | tkerber | |
392 | 1 | tkerber | def get_cell(self): |
393 | 1 | tkerber | """Returns the computational cell.
|
394 | 1 | tkerber |
|
395 | 1 | tkerber | The computational cell is the same as for the original system.
|
396 | 1 | tkerber | """
|
397 | 1 | tkerber | return self.atoms.get_cell() |
398 | 1 | tkerber | |
399 | 1 | tkerber | def get_pbc(self): |
400 | 1 | tkerber | """Returns the periodic boundary conditions.
|
401 | 1 | tkerber |
|
402 | 1 | tkerber | The boundary conditions are the same as for the original system.
|
403 | 1 | tkerber | """
|
404 | 1 | tkerber | return self.atoms.get_pbc() |
405 | 1 | tkerber | |
406 | 1 | tkerber | def get_positions(self): |
407 | 1 | tkerber | "Return the positions of the visible atoms."
|
408 | 1 | tkerber | return self.atoms.get_positions()[self.index] |
409 | 1 | tkerber | |
410 | 1 | tkerber | def set_positions(self, positions): |
411 | 1 | tkerber | "Set the positions of the visible atoms."
|
412 | 1 | tkerber | pos = self.atoms.get_positions()
|
413 | 1 | tkerber | pos[self.index] = positions
|
414 | 1 | tkerber | self.atoms.set_positions(pos)
|
415 | 1 | tkerber | |
416 | 1 | tkerber | def get_momenta(self): |
417 | 1 | tkerber | "Return the momenta of the visible atoms."
|
418 | 1 | tkerber | return self.atoms.get_momenta()[self.index] |
419 | 1 | tkerber | |
420 | 1 | tkerber | def set_momenta(self, momenta): |
421 | 1 | tkerber | "Set the momenta of the visible atoms."
|
422 | 1 | tkerber | mom = self.atoms.get_momenta()
|
423 | 1 | tkerber | mom[self.index] = momenta
|
424 | 1 | tkerber | self.atoms.set_momenta(mom)
|
425 | 1 | tkerber | |
426 | 1 | tkerber | def get_atomic_numbers(self): |
427 | 1 | tkerber | "Return the atomic numbers of the visible atoms."
|
428 | 1 | tkerber | return self.atoms.get_atomic_numbers()[self.index] |
429 | 1 | tkerber | |
430 | 1 | tkerber | def set_atomic_numbers(self, atomic_numbers): |
431 | 1 | tkerber | "Set the atomic numbers of the visible atoms."
|
432 | 1 | tkerber | z = self.atoms.get_atomic_numbers()
|
433 | 1 | tkerber | z[self.index] = atomic_numbers
|
434 | 1 | tkerber | self.atoms.set_atomic_numbers(z)
|
435 | 1 | tkerber | |
436 | 1 | tkerber | def get_tags(self): |
437 | 1 | tkerber | "Return the tags of the visible atoms."
|
438 | 1 | tkerber | return self.atoms.get_tags()[self.index] |
439 | 1 | tkerber | |
440 | 1 | tkerber | def set_tags(self, tags): |
441 | 1 | tkerber | "Set the tags of the visible atoms."
|
442 | 1 | tkerber | tg = self.atoms.get_tags()
|
443 | 1 | tkerber | tg[self.index] = tags
|
444 | 1 | tkerber | self.atoms.set_tags(tg)
|
445 | 1 | tkerber | |
446 | 1 | tkerber | def get_forces(self, *args, **kwargs): |
447 | 1 | tkerber | return self.atoms.get_forces(*args, **kwargs)[self.index] |
448 | 1 | tkerber | |
449 | 1 | tkerber | def get_stress(self): |
450 | 1 | tkerber | return self.atoms.get_stress() |
451 | 1 | tkerber | |
452 | 1 | tkerber | def get_stresses(self): |
453 | 1 | tkerber | return self.atoms.get_stresses()[self.index] |
454 | 1 | tkerber | |
455 | 1 | tkerber | def get_masses(self): |
456 | 1 | tkerber | return self.atoms.get_masses()[self.index] |
457 | 1 | tkerber | |
458 | 1 | tkerber | def get_potential_energy(self): |
459 | 1 | tkerber | """Calculate potential energy.
|
460 | 1 | tkerber |
|
461 | 1 | tkerber | Returns the potential energy of the full system.
|
462 | 1 | tkerber | """
|
463 | 1 | tkerber | return self.atoms.get_potential_energy() |
464 | 1 | tkerber | |
465 | 1 | tkerber | def get_calculator(self): |
466 | 1 | tkerber | """Returns the calculator.
|
467 | 1 | tkerber |
|
468 | 1 | tkerber | WARNING: The calculator is unaware of this filter, and sees a
|
469 | 1 | tkerber | different number of atoms.
|
470 | 1 | tkerber | """
|
471 | 1 | tkerber | return self.atoms.get_calculator() |
472 | 1 | tkerber | |
473 | 1 | tkerber | def has(self, name): |
474 | 1 | tkerber | """Check for existance of array."""
|
475 | 1 | tkerber | return self.atoms.has(name) |
476 | 1 | tkerber | |
477 | 1 | tkerber | def __len__(self): |
478 | 1 | tkerber | "Return the number of movable atoms."
|
479 | 1 | tkerber | return self.n |
480 | 1 | tkerber | |
481 | 1 | tkerber | def __getitem__(self, i): |
482 | 1 | tkerber | "Return an atom."
|
483 | 1 | tkerber | return self.atoms[self.index[i]] |
484 | 1 | tkerber | |
485 | 1 | tkerber | |
486 | 1 | tkerber | class StrainFilter: |
487 | 1 | tkerber | """Modify the supercell while keeping the scaled positions fixed.
|
488 | 1 | tkerber |
|
489 | 1 | tkerber | Presents the strain of the supercell as the generalized positions,
|
490 | 1 | tkerber | and the global stress tensor (times the volume) as the generalized
|
491 | 1 | tkerber | force.
|
492 | 1 | tkerber |
|
493 | 1 | tkerber | This filter can be used to relax the unit cell until the stress is
|
494 | 1 | tkerber | zero. If MDMin is used for this, the timestep (dt) to be used
|
495 | 1 | tkerber | depends on the system size. 0.01/x where x is a typical dimension
|
496 | 1 | tkerber | seems like a good choice.
|
497 | 1 | tkerber |
|
498 | 1 | tkerber | The stress and strain are presented as 6-vectors, the order of the
|
499 | 1 | tkerber | components follow the standard engingeering practice: xx, yy, zz,
|
500 | 1 | tkerber | yz, xz, xy.
|
501 | 1 | tkerber |
|
502 | 1 | tkerber | """
|
503 | 1 | tkerber | def __init__(self, atoms, mask=None): |
504 | 1 | tkerber | """Create a filter applying a homogeneous strain to a list of atoms.
|
505 | 1 | tkerber |
|
506 | 1 | tkerber | The first argument, atoms, is the atoms object.
|
507 | 1 | tkerber |
|
508 | 1 | tkerber | The optional second argument, mask, is a list of six booleans,
|
509 | 1 | tkerber | indicating which of the six independent components of the
|
510 | 1 | tkerber | strain that are allowed to become non-zero. It defaults to
|
511 | 1 | tkerber | [1,1,1,1,1,1].
|
512 | 1 | tkerber |
|
513 | 1 | tkerber | """
|
514 | 1 | tkerber | |
515 | 1 | tkerber | self.atoms = atoms
|
516 | 1 | tkerber | self.strain = np.zeros(6) |
517 | 1 | tkerber | |
518 | 1 | tkerber | if mask is None: |
519 | 1 | tkerber | self.mask = np.ones(6) |
520 | 1 | tkerber | else:
|
521 | 1 | tkerber | self.mask = np.array(mask)
|
522 | 1 | tkerber | |
523 | 1 | tkerber | self.origcell = atoms.get_cell()
|
524 | 1 | tkerber | |
525 | 1 | tkerber | def get_positions(self): |
526 | 1 | tkerber | return self.strain.reshape((2, 3)) |
527 | 1 | tkerber | |
528 | 1 | tkerber | def set_positions(self, new): |
529 | 1 | tkerber | new = new.ravel() * self.mask
|
530 | 1 | tkerber | eps = np.array([[1.0 + new[0], 0.5 * new[5], 0.5 * new[4]], |
531 | 1 | tkerber | [0.5 * new[5], 1.0 + new[1], 0.5 * new[3]], |
532 | 1 | tkerber | [0.5 * new[4], 0.5 * new[3], 1.0 + new[2]]]) |
533 | 1 | tkerber | |
534 | 1 | tkerber | self.atoms.set_cell(np.dot(self.origcell, eps), scale_atoms=True) |
535 | 1 | tkerber | self.strain[:] = new
|
536 | 1 | tkerber | |
537 | 1 | tkerber | def get_forces(self): |
538 | 1 | tkerber | stress = self.atoms.get_stress()
|
539 | 1 | tkerber | return -self.atoms.get_volume() * (stress * self.mask).reshape((2, 3)) |
540 | 1 | tkerber | |
541 | 1 | tkerber | def get_potential_energy(self): |
542 | 1 | tkerber | return self.atoms.get_potential_energy() |
543 | 1 | tkerber | |
544 | 1 | tkerber | def has(self, x): |
545 | 1 | tkerber | return self.atoms.has(x) |
546 | 1 | tkerber | |
547 | 1 | tkerber | def __len__(self): |
548 | 1 | tkerber | return 2 |
549 | 1 | tkerber | |
550 | 1 | tkerber | class UnitCellFilter: |
551 | 1 | tkerber | """Modify the supercell and the atom positions. """
|
552 | 1 | tkerber | def __init__(self, atoms, mask=None): |
553 | 1 | tkerber | """Create a filter that returns the atomic forces and unit
|
554 | 1 | tkerber | cell stresses together, so they can simultaneously be
|
555 | 1 | tkerber | minimized.
|
556 | 1 | tkerber |
|
557 | 1 | tkerber | The first argument, atoms, is the atoms object.
|
558 | 1 | tkerber |
|
559 | 1 | tkerber | The optional second argument, mask, is a list of booleans,
|
560 | 1 | tkerber | indicating which of the six independent
|
561 | 1 | tkerber | components of the strain are relaxed.
|
562 | 1 | tkerber | 1, True = relax to zero
|
563 | 1 | tkerber | 0, False = fixed, ignore this component
|
564 | 1 | tkerber |
|
565 | 1 | tkerber | use atom Constraints, e.g. FixAtoms, to control relaxation of
|
566 | 1 | tkerber | the atoms.
|
567 | 1 | tkerber |
|
568 | 1 | tkerber | #this should be equivalent to the StrainFilter
|
569 | 1 | tkerber | >>> atoms = Atoms(...)
|
570 | 1 | tkerber | >>> atoms.set_constraint(FixAtoms(mask=[True for atom in atoms]))
|
571 | 1 | tkerber | >>> ucf = UCFilter(atoms)
|
572 | 1 | tkerber |
|
573 | 1 | tkerber | You should not attach this UCFilter object to a
|
574 | 1 | tkerber | trajectory. Instead, create a trajectory for the atoms, and
|
575 | 1 | tkerber | attach it to an optimizer like this:
|
576 | 1 | tkerber |
|
577 | 1 | tkerber | >>> atoms = Atoms(...)
|
578 | 1 | tkerber | >>> ucf = UCFilter(atoms)
|
579 | 1 | tkerber | >>> qn = QuasiNewton(ucf)
|
580 | 1 | tkerber | >>> traj = PickleTrajectory('TiO2.traj','w',atoms)
|
581 | 1 | tkerber | >>> qn.attach(traj)
|
582 | 1 | tkerber | >>> qn.run(fmax=0.05)
|
583 | 1 | tkerber |
|
584 | 1 | tkerber | Helpful conversion table
|
585 | 1 | tkerber | ========================
|
586 | 1 | tkerber | 0.05 eV/A^3 = 8 GPA
|
587 | 1 | tkerber | 0.003 eV/A^3 = 0.48 GPa
|
588 | 1 | tkerber | 0.0006 eV/A^3 = 0.096 GPa
|
589 | 1 | tkerber | 0.0003 eV/A^3 = 0.048 GPa
|
590 | 1 | tkerber | 0.0001 eV/A^3 = 0.02 GPa
|
591 | 1 | tkerber | """
|
592 | 1 | tkerber | |
593 | 1 | tkerber | self.atoms = atoms
|
594 | 1 | tkerber | self.strain = np.zeros(6) |
595 | 1 | tkerber | |
596 | 1 | tkerber | if mask is None: |
597 | 1 | tkerber | self.mask = np.ones(6) |
598 | 1 | tkerber | else:
|
599 | 1 | tkerber | self.mask = np.array(mask)
|
600 | 1 | tkerber | |
601 | 1 | tkerber | self.origcell = atoms.get_cell()
|
602 | 1 | tkerber | |
603 | 1 | tkerber | def get_positions(self): |
604 | 1 | tkerber | '''
|
605 | 1 | tkerber | this returns an array with shape (natoms+2,3).
|
606 | 1 | tkerber |
|
607 | 1 | tkerber | the first natoms rows are the positions of the atoms, the last
|
608 | 1 | tkerber | two rows are the strains associated with the unit cell
|
609 | 1 | tkerber | '''
|
610 | 1 | tkerber | |
611 | 1 | tkerber | atom_positions = self.atoms.get_positions()
|
612 | 1 | tkerber | strains = self.strain.reshape((2, 3)) |
613 | 1 | tkerber | |
614 | 1 | tkerber | natoms = len(self.atoms) |
615 | 1 | tkerber | all_pos = np.zeros((natoms+2,3),np.float) |
616 | 1 | tkerber | all_pos[0:natoms,:] = atom_positions
|
617 | 1 | tkerber | all_pos[natoms:,:] = strains |
618 | 1 | tkerber | |
619 | 1 | tkerber | return all_pos
|
620 | 1 | tkerber | |
621 | 1 | tkerber | def set_positions(self, new): |
622 | 1 | tkerber | '''
|
623 | 1 | tkerber | new is an array with shape (natoms+2,3).
|
624 | 1 | tkerber |
|
625 | 1 | tkerber | the first natoms rows are the positions of the atoms, the last
|
626 | 1 | tkerber | two rows are the strains used to change the cell shape.
|
627 | 1 | tkerber |
|
628 | 1 | tkerber | The atom positions are set first, then the unit cell is
|
629 | 1 | tkerber | changed keeping the atoms in their scaled positions.
|
630 | 1 | tkerber | '''
|
631 | 1 | tkerber | |
632 | 1 | tkerber | natoms = len(self.atoms) |
633 | 1 | tkerber | |
634 | 1 | tkerber | atom_positions = new[0:natoms,:]
|
635 | 1 | tkerber | self.atoms.set_positions(atom_positions)
|
636 | 1 | tkerber | |
637 | 1 | tkerber | new = new[natoms:,:] #this is only the strains
|
638 | 1 | tkerber | new = new.ravel() * self.mask
|
639 | 1 | tkerber | eps = np.array([[1.0 + new[0], 0.5 * new[5], 0.5 * new[4]], |
640 | 1 | tkerber | [0.5 * new[5], 1.0 + new[1], 0.5 * new[3]], |
641 | 1 | tkerber | [0.5 * new[4], 0.5 * new[3], 1.0 + new[2]]]) |
642 | 1 | tkerber | |
643 | 1 | tkerber | self.atoms.set_cell(np.dot(self.origcell, eps), scale_atoms=True) |
644 | 1 | tkerber | self.strain[:] = new
|
645 | 1 | tkerber | |
646 | 1 | tkerber | def get_forces(self,apply_constraint=False): |
647 | 1 | tkerber | '''
|
648 | 1 | tkerber | returns an array with shape (natoms+2,3) of the atomic forces
|
649 | 1 | tkerber | and unit cell stresses.
|
650 | 1 | tkerber |
|
651 | 1 | tkerber | the first natoms rows are the forces on the atoms, the last
|
652 | 1 | tkerber | two rows are the stresses on the unit cell, which have been
|
653 | 1 | tkerber | reshaped to look like "atomic forces". i.e.,
|
654 | 1 | tkerber |
|
655 | 1 | tkerber | f[-2] = -vol*[sxx,syy,szz]*mask[0:3]
|
656 | 1 | tkerber | f[-1] = -vol*[syz, sxz, sxy]*mask[3:]
|
657 | 1 | tkerber |
|
658 | 1 | tkerber | apply_constraint is an argument expected by ase
|
659 | 1 | tkerber | '''
|
660 | 1 | tkerber | |
661 | 1 | tkerber | stress = self.atoms.get_stress()
|
662 | 1 | tkerber | atom_forces = self.atoms.get_forces()
|
663 | 1 | tkerber | |
664 | 1 | tkerber | natoms = len(self.atoms) |
665 | 1 | tkerber | all_forces = np.zeros((natoms+2,3),np.float) |
666 | 1 | tkerber | all_forces[0:natoms,:] = atom_forces
|
667 | 1 | tkerber | |
668 | 1 | tkerber | vol = self.atoms.get_volume()
|
669 | 1 | tkerber | stress_forces = -vol * (stress * self.mask).reshape((2, 3)) |
670 | 1 | tkerber | all_forces[natoms:,:] = stress_forces |
671 | 1 | tkerber | return all_forces
|
672 | 1 | tkerber | |
673 | 1 | tkerber | def get_potential_energy(self): |
674 | 1 | tkerber | return self.atoms.get_potential_energy() |
675 | 1 | tkerber | |
676 | 1 | tkerber | def has(self, x): |
677 | 1 | tkerber | return self.atoms.has(x) |
678 | 1 | tkerber | |
679 | 1 | tkerber | def __len__(self): |
680 | 1 | tkerber | return (2 + len(self.atoms)) |