root / ase / vibrations.py @ 17
Historique | Voir | Annoter | Télécharger (15,57 ko)
1 | 1 | tkerber | # -*- coding: utf-8 -*-
|
---|---|---|---|
2 | 1 | tkerber | |
3 | 1 | tkerber | """Vibrational modes."""
|
4 | 1 | tkerber | |
5 | 1 | tkerber | import pickle |
6 | 1 | tkerber | from math import sin, pi, sqrt |
7 | 1 | tkerber | from os import remove |
8 | 1 | tkerber | from os.path import isfile |
9 | 1 | tkerber | import sys |
10 | 1 | tkerber | |
11 | 1 | tkerber | import numpy as np |
12 | 1 | tkerber | |
13 | 1 | tkerber | import ase.units as units |
14 | 1 | tkerber | from ase.io.trajectory import PickleTrajectory |
15 | 1 | tkerber | from ase.parallel import rank, barrier, paropen |
16 | 1 | tkerber | |
17 | 1 | tkerber | |
18 | 1 | tkerber | class Vibrations: |
19 | 1 | tkerber | """Class for calculating vibrational modes using finite difference.
|
20 | 1 | tkerber |
|
21 | 1 | tkerber | The vibrational modes are calculated from a finite difference
|
22 | 1 | tkerber | approximation of the Hessian matrix.
|
23 | 1 | tkerber |
|
24 | 1 | tkerber | The *summary()*, *get_energies()* and *get_frequencies()*
|
25 | 1 | tkerber | methods all take an optional *method* keyword. Use
|
26 | 1 | tkerber | method='Frederiksen' to use the method described in:
|
27 | 1 | tkerber |
|
28 | 1 | tkerber | T. Frederiksen, M. Paulsson, M. Brandbyge, A. P. Jauho:
|
29 | 1 | tkerber | "Inelastic transport theory from first-principles: methodology
|
30 | 1 | tkerber | and applications for nanoscale devices",
|
31 | 1 | tkerber | Phys. Rev. B 75, 205413 (2007)
|
32 | 1 | tkerber |
|
33 | 1 | tkerber | atoms: Atoms object
|
34 | 1 | tkerber | The atoms to work on.
|
35 | 1 | tkerber | indices: list of int
|
36 | 1 | tkerber | List of indices of atoms to vibrate. Default behavior is
|
37 | 1 | tkerber | to vibrate all atoms.
|
38 | 1 | tkerber | name: str
|
39 | 1 | tkerber | Name to use for files.
|
40 | 1 | tkerber | delta: float
|
41 | 1 | tkerber | Magnitude of displacements.
|
42 | 1 | tkerber | nfree: int
|
43 | 1 | tkerber | Number of displacements per atom and cartesian coordinate,
|
44 | 1 | tkerber | 2 and 4 are supported. Default is 2 which will displace
|
45 | 1 | tkerber | each atom +delta and -delta for each cartesian coordinate.
|
46 | 1 | tkerber |
|
47 | 1 | tkerber | Example:
|
48 | 1 | tkerber |
|
49 | 1 | tkerber | >>> from ase import Atoms
|
50 | 1 | tkerber | >>> from ase.calculators.emt import EMT
|
51 | 1 | tkerber | >>> from ase.optimizers import BFGS
|
52 | 1 | tkerber | >>> from ase.vibrations import Vibrations
|
53 | 1 | tkerber | >>> n2 = Atoms('N2', [(0, 0, 0), (0, 0, 1.1)],
|
54 | 1 | tkerber | ... calculator=EMT())
|
55 | 1 | tkerber | >>> BFGS(n2).run(fmax=0.01)
|
56 | 1 | tkerber | BFGS: 0 19:16:06 0.042171 2.9357
|
57 | 1 | tkerber | BFGS: 1 19:16:07 0.104197 3.9270
|
58 | 1 | tkerber | BFGS: 2 19:16:07 0.000963 0.4142
|
59 | 1 | tkerber | BFGS: 3 19:16:07 0.000027 0.0698
|
60 | 1 | tkerber | BFGS: 4 19:16:07 0.000000 0.0010
|
61 | 1 | tkerber | >>> vib = Vibrations(n2)
|
62 | 1 | tkerber | >>> vib.run()
|
63 | 1 | tkerber | >>> vib.summary()
|
64 | 1 | tkerber | ---------------------
|
65 | 1 | tkerber | # meV cm^-1
|
66 | 1 | tkerber | ---------------------
|
67 | 1 | tkerber | 0 0.0i 0.0i
|
68 | 1 | tkerber | 1 0.0i 0.0i
|
69 | 1 | tkerber | 2 0.0i 0.0i
|
70 | 1 | tkerber | 3 1.6 13.1
|
71 | 1 | tkerber | 4 1.6 13.1
|
72 | 1 | tkerber | 5 232.7 1877.2
|
73 | 1 | tkerber | ---------------------
|
74 | 1 | tkerber | Zero-point energy: 0.118 eV
|
75 | 1 | tkerber | Thermodynamic properties at 298.00 K
|
76 | 1 | tkerber | Enthalpy: 0.050 eV
|
77 | 1 | tkerber | Entropy : 0.648 meV/K
|
78 | 1 | tkerber | T*S : 0.193 eV
|
79 | 1 | tkerber | E->G : -0.025 eV
|
80 | 1 | tkerber |
|
81 | 1 | tkerber | >>> vib.write_mode(-1) # write last mode to trajectory file
|
82 | 1 | tkerber |
|
83 | 1 | tkerber | """
|
84 | 1 | tkerber | def __init__(self, atoms, indices=None, name='vib', delta=0.01, nfree=2): |
85 | 1 | tkerber | assert nfree in [2, 4] |
86 | 1 | tkerber | self.atoms = atoms
|
87 | 1 | tkerber | if indices is None: |
88 | 1 | tkerber | indices = range(len(atoms)) |
89 | 1 | tkerber | self.indices = np.asarray(indices)
|
90 | 1 | tkerber | self.name = name
|
91 | 1 | tkerber | self.delta = delta
|
92 | 1 | tkerber | self.nfree = nfree
|
93 | 1 | tkerber | self.H = None |
94 | 1 | tkerber | self.ir = None |
95 | 1 | tkerber | |
96 | 1 | tkerber | def run(self): |
97 | 1 | tkerber | """Run the vibration calculations.
|
98 | 1 | tkerber |
|
99 | 1 | tkerber | This will calculate the forces for 6 displacements per atom
|
100 | 1 | tkerber | ±x, ±y, ±z. Only those calculations that are not already done
|
101 | 1 | tkerber | will be started. Be aware that an interrupted calculation may
|
102 | 1 | tkerber | produce an empty file (ending with .pckl), which must be deleted
|
103 | 1 | tkerber | before restarting the job. Otherwise the forces will not be
|
104 | 1 | tkerber | calculated for that displacement."""
|
105 | 1 | tkerber | filename = self.name + '.eq.pckl' |
106 | 1 | tkerber | if not isfile(filename): |
107 | 1 | tkerber | barrier() |
108 | 1 | tkerber | forces = self.atoms.get_forces()
|
109 | 1 | tkerber | if self.ir: |
110 | 1 | tkerber | dipole = self.calc.get_dipole_moment(self.atoms) |
111 | 1 | tkerber | if rank == 0: |
112 | 1 | tkerber | fd = open(filename, 'w') |
113 | 1 | tkerber | if self.ir: |
114 | 1 | tkerber | pickle.dump([forces, dipole], fd) |
115 | 1 | tkerber | sys.stdout.write( |
116 | 1 | tkerber | 'Writing %s, dipole moment = (%.6f %.6f %.6f)\n' %
|
117 | 1 | tkerber | (filename, dipole[0], dipole[1], dipole[2])) |
118 | 1 | tkerber | else:
|
119 | 1 | tkerber | pickle.dump(forces, fd) |
120 | 1 | tkerber | sys.stdout.write('Writing %s\n' % filename)
|
121 | 1 | tkerber | fd.close() |
122 | 1 | tkerber | sys.stdout.flush() |
123 | 1 | tkerber | |
124 | 1 | tkerber | p = self.atoms.positions.copy()
|
125 | 1 | tkerber | for a in self.indices: |
126 | 1 | tkerber | for i in range(3): |
127 | 1 | tkerber | for sign in [-1, 1]: |
128 | 1 | tkerber | for ndis in range(1, self.nfree//2+1): |
129 | 1 | tkerber | filename = ('%s.%d%s%s.pckl' %
|
130 | 1 | tkerber | (self.name, a, 'xyz'[i], ndis*' +-'[sign])) |
131 | 1 | tkerber | if isfile(filename):
|
132 | 1 | tkerber | continue
|
133 | 1 | tkerber | barrier() |
134 | 1 | tkerber | self.atoms.positions[a, i] = (p[a, i] +
|
135 | 1 | tkerber | ndis * sign * self.delta)
|
136 | 1 | tkerber | forces = self.atoms.get_forces()
|
137 | 1 | tkerber | if self.ir: |
138 | 1 | tkerber | dipole = self.calc.get_dipole_moment(self.atoms) |
139 | 1 | tkerber | if rank == 0: |
140 | 1 | tkerber | fd = open(filename, 'w') |
141 | 1 | tkerber | if self.ir: |
142 | 1 | tkerber | pickle.dump([forces, dipole], fd) |
143 | 1 | tkerber | sys.stdout.write( |
144 | 1 | tkerber | 'Writing %s, ' % filename +
|
145 | 1 | tkerber | 'dipole moment = (%.6f %.6f %.6f)\n' %
|
146 | 1 | tkerber | (dipole[0], dipole[1], dipole[2])) |
147 | 1 | tkerber | else:
|
148 | 1 | tkerber | pickle.dump(forces, fd) |
149 | 1 | tkerber | sys.stdout.write('Writing %s\n' % filename)
|
150 | 1 | tkerber | fd.close() |
151 | 1 | tkerber | sys.stdout.flush() |
152 | 1 | tkerber | self.atoms.positions[a, i] = p[a, i]
|
153 | 1 | tkerber | self.atoms.set_positions(p)
|
154 | 1 | tkerber | |
155 | 1 | tkerber | def clean(self): |
156 | 1 | tkerber | if isfile(self.name + '.eq.pckl'): |
157 | 1 | tkerber | remove(self.name + '.eq.pckl') |
158 | 1 | tkerber | |
159 | 1 | tkerber | for a in self.indices: |
160 | 1 | tkerber | for i in 'xyz': |
161 | 1 | tkerber | for sign in '-+': |
162 | 1 | tkerber | for ndis in range(1, self.nfree/2+1): |
163 | 1 | tkerber | name = '%s.%d%s%s.pckl' % (self.name, a, i, ndis*sign) |
164 | 1 | tkerber | if isfile(name):
|
165 | 1 | tkerber | remove(name) |
166 | 1 | tkerber | |
167 | 1 | tkerber | def read(self, method='standard', direction='central'): |
168 | 1 | tkerber | self.method = method.lower()
|
169 | 1 | tkerber | self.direction = direction.lower()
|
170 | 1 | tkerber | assert self.method in ['standard', 'frederiksen'] |
171 | 1 | tkerber | assert self.direction in ['central', 'forward', 'backward'] |
172 | 1 | tkerber | |
173 | 1 | tkerber | n = 3 * len(self.indices) |
174 | 1 | tkerber | H = np.empty((n, n)) |
175 | 1 | tkerber | r = 0
|
176 | 1 | tkerber | if direction != 'central': |
177 | 1 | tkerber | feq = pickle.load(open(self.name + '.eq.pckl')) |
178 | 1 | tkerber | for a in self.indices: |
179 | 1 | tkerber | for i in 'xyz': |
180 | 1 | tkerber | name = '%s.%d%s' % (self.name, a, i) |
181 | 1 | tkerber | fminus = pickle.load(open(name + '-.pckl')) |
182 | 1 | tkerber | fplus = pickle.load(open(name + '+.pckl')) |
183 | 1 | tkerber | if self.method == 'frederiksen': |
184 | 1 | tkerber | fminus[a] -= fminus.sum(0)
|
185 | 1 | tkerber | fplus[a] -= fplus.sum(0)
|
186 | 1 | tkerber | if self.nfree == 4: |
187 | 1 | tkerber | fminusminus = pickle.load(open(name + '--.pckl')) |
188 | 1 | tkerber | fplusplus = pickle.load(open(name + '++.pckl')) |
189 | 1 | tkerber | if self.method == 'frederiksen': |
190 | 1 | tkerber | fminusminus[a] -= fminusminus.sum(0)
|
191 | 1 | tkerber | fplusplus[a] -= fplusplus.sum(0)
|
192 | 1 | tkerber | if self.direction == 'central': |
193 | 1 | tkerber | if self.nfree == 2: |
194 | 1 | tkerber | H[r] = .5 * (fminus - fplus)[self.indices].ravel() |
195 | 1 | tkerber | else:
|
196 | 1 | tkerber | H[r] = H[r] = (-fminusminus + |
197 | 1 | tkerber | 8 * fminus -
|
198 | 1 | tkerber | 8 * fplus +
|
199 | 1 | tkerber | fplusplus)[self.indices].ravel() / 12.0 |
200 | 1 | tkerber | elif self.direction == 'forward': |
201 | 1 | tkerber | H[r] = (feq - fplus)[self.indices].ravel()
|
202 | 1 | tkerber | else: # self.direction == 'backward': |
203 | 1 | tkerber | H[r] = (fminus - feq)[self.indices].ravel()
|
204 | 1 | tkerber | H[r] /= 2 * self.delta |
205 | 1 | tkerber | r += 1
|
206 | 1 | tkerber | H += H.copy().T |
207 | 1 | tkerber | self.H = H
|
208 | 1 | tkerber | m = self.atoms.get_masses()
|
209 | 1 | tkerber | self.im = np.repeat(m[self.indices]**-0.5, 3) |
210 | 1 | tkerber | omega2, modes = np.linalg.eigh(self.im[:, None] * H * self.im) |
211 | 1 | tkerber | self.modes = modes.T.copy()
|
212 | 1 | tkerber | |
213 | 1 | tkerber | # Conversion factor:
|
214 | 1 | tkerber | s = units._hbar * 1e10 / sqrt(units._e * units._amu)
|
215 | 1 | tkerber | self.hnu = s * omega2.astype(complex)**0.5 |
216 | 1 | tkerber | |
217 | 1 | tkerber | def get_energies(self, method='standard', direction='central'): |
218 | 1 | tkerber | """Get vibration energies in eV."""
|
219 | 1 | tkerber | if (self.H is None or method.lower() != self.method or |
220 | 1 | tkerber | direction.lower() != self.direction):
|
221 | 1 | tkerber | self.read(method, direction)
|
222 | 1 | tkerber | return self.hnu |
223 | 1 | tkerber | |
224 | 1 | tkerber | def get_frequencies(self, method='standard', direction='central'): |
225 | 1 | tkerber | """Get vibration frequencies in cm^-1."""
|
226 | 1 | tkerber | s = 0.01 * units._e / units._c / units._hplanck
|
227 | 1 | tkerber | return s * self.get_energies(method, direction) |
228 | 1 | tkerber | |
229 | 1 | tkerber | def summary(self, method='standard', direction='central', T=298., |
230 | 1 | tkerber | threshold=10, freq=None, log=sys.stdout): |
231 | 1 | tkerber | """Print a summary of the frequencies and derived thermodynamic
|
232 | 1 | tkerber | properties. The equations for the calculation of the enthalpy and
|
233 | 1 | tkerber | entropy diverge for zero frequencies and a threshold value is used
|
234 | 1 | tkerber | to ignore extremely low frequencies (default = 10 cm^-1).
|
235 | 1 | tkerber |
|
236 | 1 | tkerber | Parameters:
|
237 | 1 | tkerber |
|
238 | 1 | tkerber | method : string
|
239 | 1 | tkerber | Can be 'standard'(default) or 'Frederiksen'.
|
240 | 1 | tkerber | direction: string
|
241 | 1 | tkerber | Direction for finite differences. Can be one of 'central'
|
242 | 1 | tkerber | (default), 'forward', 'backward'.
|
243 | 1 | tkerber | T : float
|
244 | 1 | tkerber | Temperature in K at which thermodynamic properties are calculated.
|
245 | 1 | tkerber | threshold : float
|
246 | 1 | tkerber | Threshold value for low frequencies (default 10 cm^-1).
|
247 | 1 | tkerber | freq : numpy array
|
248 | 1 | tkerber | Optional. Can be used to create a summary on a set of known
|
249 | 1 | tkerber | frequencies.
|
250 | 1 | tkerber | destination : string or object with a .write() method
|
251 | 1 | tkerber | Where to print the summary info. Default is to print to
|
252 | 1 | tkerber | standard output. If another string is given, it creates that
|
253 | 1 | tkerber | text file and writes the output to it. If a file object (or any
|
254 | 1 | tkerber | object with a .write(string) function) is given, it writes to that
|
255 | 1 | tkerber | file.
|
256 | 1 | tkerber |
|
257 | 1 | tkerber | Notes:
|
258 | 1 | tkerber |
|
259 | 1 | tkerber | The enthalpy and entropy calculations are very sensitive to low
|
260 | 1 | tkerber | frequencies and the threshold value must be used with caution."""
|
261 | 1 | tkerber | |
262 | 1 | tkerber | if isinstance(log, str): |
263 | 1 | tkerber | log = paropen(log, 'a')
|
264 | 1 | tkerber | write = log.write |
265 | 1 | tkerber | |
266 | 1 | tkerber | s = 0.01 * units._e / units._c / units._hplanck
|
267 | 1 | tkerber | if freq != None: |
268 | 1 | tkerber | hnu = freq / s |
269 | 1 | tkerber | else:
|
270 | 1 | tkerber | hnu = self.get_energies(method, direction)
|
271 | 1 | tkerber | write('---------------------\n')
|
272 | 1 | tkerber | write(' # meV cm^-1\n')
|
273 | 1 | tkerber | write('---------------------\n')
|
274 | 1 | tkerber | for n, e in enumerate(hnu): |
275 | 1 | tkerber | if e.imag != 0: |
276 | 1 | tkerber | c = 'i'
|
277 | 1 | tkerber | e = e.imag |
278 | 1 | tkerber | else:
|
279 | 1 | tkerber | c = ' '
|
280 | 1 | tkerber | write('%3d %6.1f%s %7.1f%s\n' % (n, 1000 * e, c, s * e, c)) |
281 | 1 | tkerber | write('---------------------\n')
|
282 | 1 | tkerber | write('Zero-point energy: %.3f eV\n' %
|
283 | 1 | tkerber | self.get_zero_point_energy(freq=freq))
|
284 | 1 | tkerber | write('Thermodynamic properties at %.2f K\n' % T)
|
285 | 1 | tkerber | write('Enthalpy: %.3f eV\n' % self.get_enthalpy(method=method, |
286 | 1 | tkerber | direction=direction, |
287 | 1 | tkerber | T=T, |
288 | 1 | tkerber | threshold=threshold, |
289 | 1 | tkerber | freq=freq)) |
290 | 1 | tkerber | write('Entropy : %.3f meV/K\n' %
|
291 | 1 | tkerber | (1E3 * self.get_entropy(method=method, |
292 | 1 | tkerber | direction=direction, |
293 | 1 | tkerber | T=T, |
294 | 1 | tkerber | threshold=threshold, |
295 | 1 | tkerber | freq=freq))) |
296 | 1 | tkerber | write('T*S : %.3f eV\n' %
|
297 | 1 | tkerber | (T * self.get_entropy(method=method,
|
298 | 1 | tkerber | direction=direction, |
299 | 1 | tkerber | T=T, |
300 | 1 | tkerber | threshold=threshold, |
301 | 1 | tkerber | freq=freq))) |
302 | 1 | tkerber | write('E->G : %.3f eV\n' %
|
303 | 1 | tkerber | (self.get_zero_point_energy(freq=freq) +
|
304 | 1 | tkerber | self.get_enthalpy(method=method,
|
305 | 1 | tkerber | direction=direction, |
306 | 1 | tkerber | T=T, |
307 | 1 | tkerber | threshold=threshold, |
308 | 1 | tkerber | freq=freq) - |
309 | 1 | tkerber | T * self.get_entropy(method=method,
|
310 | 1 | tkerber | direction=direction, |
311 | 1 | tkerber | T=T, |
312 | 1 | tkerber | threshold=threshold, |
313 | 1 | tkerber | freq=freq))) |
314 | 1 | tkerber | |
315 | 1 | tkerber | def get_zero_point_energy(self, freq=None): |
316 | 1 | tkerber | if freq is None: |
317 | 1 | tkerber | return 0.5 * self.hnu.real.sum() |
318 | 1 | tkerber | else:
|
319 | 1 | tkerber | s = 0.01 * units._e / units._c / units._hplanck
|
320 | 1 | tkerber | return 0.5 * freq.real.sum() / s |
321 | 1 | tkerber | |
322 | 1 | tkerber | def get_enthalpy(self, method='standard', direction='central', T=298.0, |
323 | 1 | tkerber | threshold=10, freq=None): |
324 | 1 | tkerber | H = 0.0
|
325 | 1 | tkerber | if freq is None: |
326 | 1 | tkerber | freq = self.get_frequencies(method=method, direction=direction)
|
327 | 1 | tkerber | for f in freq: |
328 | 1 | tkerber | if f.imag == 0 and f.real >= threshold: |
329 | 1 | tkerber | # The formula diverges for f->0
|
330 | 1 | tkerber | x = (f.real * 100 * units._hplanck * units._c) / units._k
|
331 | 1 | tkerber | H += units._k / units._e * x / (np.exp(x / T) - 1)
|
332 | 1 | tkerber | return H
|
333 | 1 | tkerber | |
334 | 1 | tkerber | def get_entropy(self, method='standard', direction='central', T=298.0, |
335 | 1 | tkerber | threshold=10, freq=None): |
336 | 1 | tkerber | S = 0.0
|
337 | 1 | tkerber | if freq is None: |
338 | 1 | tkerber | freq=self.get_frequencies(method=method, direction=direction)
|
339 | 1 | tkerber | for f in freq: |
340 | 1 | tkerber | if f.imag == 0 and f.real >= threshold: |
341 | 1 | tkerber | # The formula diverges for f->0
|
342 | 1 | tkerber | x = (f.real * 100 * units._hplanck * units._c) / units._k
|
343 | 1 | tkerber | S += (units._k / units._e * (((x / T) / (np.exp(x / T) - 1)) -
|
344 | 1 | tkerber | np.log(1 - np.exp(-x / T))))
|
345 | 1 | tkerber | return S
|
346 | 1 | tkerber | |
347 | 1 | tkerber | def get_mode(self, n): |
348 | 1 | tkerber | mode = np.zeros((len(self.atoms), 3)) |
349 | 1 | tkerber | mode[self.indices] = (self.modes[n] * self.im).reshape((-1, 3)) |
350 | 1 | tkerber | return mode
|
351 | 1 | tkerber | |
352 | 1 | tkerber | def write_mode(self, n, kT=units.kB * 300, nimages=30): |
353 | 1 | tkerber | """Write mode to trajectory file."""
|
354 | 1 | tkerber | mode = self.get_mode(n) * sqrt(kT / abs(self.hnu[n])) |
355 | 1 | tkerber | p = self.atoms.positions.copy()
|
356 | 1 | tkerber | n %= 3 * len(self.indices) |
357 | 1 | tkerber | traj = PickleTrajectory('%s.%d.traj' % (self.name, n), 'w') |
358 | 1 | tkerber | calc = self.atoms.get_calculator()
|
359 | 1 | tkerber | self.atoms.set_calculator()
|
360 | 1 | tkerber | for x in np.linspace(0, 2 * pi, nimages, endpoint=False): |
361 | 1 | tkerber | self.atoms.set_positions(p + sin(x) * mode)
|
362 | 1 | tkerber | traj.write(self.atoms)
|
363 | 1 | tkerber | self.atoms.set_positions(p)
|
364 | 1 | tkerber | self.atoms.set_calculator(calc)
|
365 | 1 | tkerber | traj.close() |
366 | 1 | tkerber | |
367 | 1 | tkerber | def write_jmol(self): |
368 | 1 | tkerber | """Writes file for viewing of the modes with jmol."""
|
369 | 1 | tkerber | |
370 | 1 | tkerber | fd = open(self.name + '.xyz', 'w') |
371 | 1 | tkerber | symbols = self.atoms.get_chemical_symbols()
|
372 | 1 | tkerber | f = self.get_frequencies()
|
373 | 1 | tkerber | for n in range(3 * len(self.atoms)): |
374 | 1 | tkerber | fd.write('%6d\n' % len(self.atoms)) |
375 | 1 | tkerber | if f[n].imag != 0: |
376 | 1 | tkerber | c = 'i'
|
377 | 1 | tkerber | f[n] = f[n].imag |
378 | 1 | tkerber | else:
|
379 | 1 | tkerber | c = ' '
|
380 | 1 | tkerber | fd.write('Mode #%d, f = %.1f%s cm^-1' % (n, f[n], c))
|
381 | 1 | tkerber | if self.ir: |
382 | 1 | tkerber | fd.write(', I = %.4f (D/Å)^2 amu^-1.\n' % self.intensities[n]) |
383 | 1 | tkerber | else:
|
384 | 1 | tkerber | fd.write('.\n')
|
385 | 1 | tkerber | mode = self.get_mode(n)
|
386 | 1 | tkerber | for i, pos in enumerate(self.atoms.positions): |
387 | 1 | tkerber | fd.write('%2s %12.5f %12.5f %12.5f %12.5f %12.5f %12.5f \n' %
|
388 | 1 | tkerber | (symbols[i], pos[0], pos[1], pos[2], |
389 | 1 | tkerber | mode[i,0], mode[i,1], mode[i,2])) |
390 | 1 | tkerber | fd.close() |