Statistiques
| Révision :

root / ase / calculators / tip3p.py @ 1

Historique | Voir | Annoter | Télécharger (5,83 ko)

1 1 tkerber
"""TIP3P potential, constraints and dynamics."""
2 1 tkerber
from math import pi, sin, cos
3 1 tkerber
4 1 tkerber
import numpy as np
5 1 tkerber
6 1 tkerber
import ase.units as units
7 1 tkerber
from ase.parallel import world
8 1 tkerber
from ase.md.md import MolecularDynamics
9 1 tkerber
10 1 tkerber
qH = 0.417
11 1 tkerber
sigma0 = 3.15061
12 1 tkerber
epsilon0 = 0.1521 * units.kcal / units.mol
13 1 tkerber
rOH = 0.9572
14 1 tkerber
thetaHOH = 104.52 / 180 * pi
15 1 tkerber
16 1 tkerber
17 1 tkerber
class TIP3P:
18 1 tkerber
    def __init__(self, rc=9.0, width=1.0):
19 1 tkerber
        self.energy = None
20 1 tkerber
        self.forces = None
21 1 tkerber
22 1 tkerber
        self.rc1 = rc - width
23 1 tkerber
        self.rc2 = rc
24 1 tkerber
25 1 tkerber
    def get_spin_polarized(self):
26 1 tkerber
        return False
27 1 tkerber
28 1 tkerber
    def update(self, atoms):
29 1 tkerber
        if (self.energy is None or
30 1 tkerber
            len(self.numbers) != len(atoms) or
31 1 tkerber
            (self.numbers != atoms.get_atomic_numbers()).any()):
32 1 tkerber
            self.calculate(atoms)
33 1 tkerber
        elif ((self.positions != atoms.get_positions()).any() or
34 1 tkerber
              (self.pbc != atoms.get_pbc()).any() or
35 1 tkerber
              (self.cell != atoms.get_cell()).any()):
36 1 tkerber
            self.calculate(atoms)
37 1 tkerber
38 1 tkerber
    def calculation_required(self, atoms, quantities):
39 1 tkerber
        if len(quantities) == 0:
40 1 tkerber
            return False
41 1 tkerber
42 1 tkerber
        return (self.energy is None or
43 1 tkerber
                len(self.numbers) != len(atoms) or
44 1 tkerber
                (self.numbers != atoms.get_atomic_numbers()).any() or
45 1 tkerber
                (self.positions != atoms.get_positions()).any() or
46 1 tkerber
                (self.pbc != atoms.get_pbc()).any() or
47 1 tkerber
                (self.cell != atoms.get_cell()).any())
48 1 tkerber
49 1 tkerber
    def get_potential_energy(self, atoms):
50 1 tkerber
        self.update(atoms)
51 1 tkerber
        return self.energy
52 1 tkerber
53 1 tkerber
    def get_forces(self, atoms):
54 1 tkerber
        self.update(atoms)
55 1 tkerber
        return self.forces.copy()
56 1 tkerber
57 1 tkerber
    def get_stress(self, atoms):
58 1 tkerber
        raise NotImplementedError
59 1 tkerber
60 1 tkerber
    def calculate(self, atoms):
61 1 tkerber
        self.positions = atoms.get_positions().copy()
62 1 tkerber
        self.cell = atoms.get_cell().copy()
63 1 tkerber
        self.pbc = atoms.get_pbc().copy()
64 1 tkerber
        natoms = len(atoms)
65 1 tkerber
        nH2O = natoms // 3
66 1 tkerber
67 1 tkerber
        assert self.pbc.all()
68 1 tkerber
        C = self.cell.diagonal()
69 1 tkerber
        assert not (self.cell - np.diag(C)).any()
70 1 tkerber
        assert (C >= 2 * self.rc2).all()
71 1 tkerber
        self.numbers = atoms.get_atomic_numbers()
72 1 tkerber
        Z = self.numbers.reshape((-1, 3))
73 1 tkerber
        assert (Z[:, 1:] == 1).all() and (Z[:, 0] == 8).all()
74 1 tkerber
75 1 tkerber
        R = self.positions.reshape((nH2O, 3, 3))
76 1 tkerber
        RO = R[:, 0]
77 1 tkerber
78 1 tkerber
        self.energy = 0.0
79 1 tkerber
        self.forces = np.zeros((natoms, 3))
80 1 tkerber
81 1 tkerber
        if world is None:
82 1 tkerber
            mya = range(nH2O - 1)
83 1 tkerber
        else:
84 1 tkerber
            rank = world.rank
85 1 tkerber
            size = world.size
86 1 tkerber
            assert nH2O // (2 * size) == 0
87 1 tkerber
            mynH2O = nH2O // 2 // size
88 1 tkerber
            mya = (range(rank * n, (rank + 1) * n) +
89 1 tkerber
                   range((size - rank - 1) * n, (size - rank) * n))
90 1 tkerber
91 1 tkerber
        q = np.empty(3)
92 1 tkerber
        q[:] = qH * (units.Hartree * units.Bohr)**0.5
93 1 tkerber
        q[0] *= -2
94 1 tkerber
95 1 tkerber
        for a in mya:
96 1 tkerber
            DOO = (RO[a + 1:] - RO[a] + 0.5 * C) % C - 0.5 * C
97 1 tkerber
            dOO = (DOO**2).sum(axis=1)**0.5
98 1 tkerber
            x1 = dOO > self.rc1
99 1 tkerber
            x2 = dOO < self.rc2
100 1 tkerber
            f = np.zeros(nH2O - a - 1)
101 1 tkerber
            f[x2] = 1.0
102 1 tkerber
            dfdd = np.zeros(nH2O - a - 1)
103 1 tkerber
            x12 = np.logical_and(x1, x2)
104 1 tkerber
            d = (dOO[x12] - self.rc1) / (self.rc2 - self.rc1)
105 1 tkerber
            f[x12] -= d**2 * (3.0 - 2.0 * d)
106 1 tkerber
            dfdd[x12] -= 6.0 / (self.rc2 - self.rc1) * d * (1.0 - d)
107 1 tkerber
108 1 tkerber
            y = (sigma0 / dOO)**6
109 1 tkerber
            y2 = y**2
110 1 tkerber
            e = 4 * epsilon0 * (y2 - y)
111 1 tkerber
            self.energy += np.dot(e, f)
112 1 tkerber
            dedd = 24 * epsilon0 * (2 * y2 - y) / dOO * f - e * dfdd
113 1 tkerber
            F = (dedd / dOO)[:, np.newaxis] * DOO
114 1 tkerber
            self.forces[(a + 1) * 3::3] += F
115 1 tkerber
            self.forces[a * 3] -= F.sum(axis=0)
116 1 tkerber
117 1 tkerber
            for i in range(3):
118 1 tkerber
                D = (R[a + 1:] - R[a, i] + 0.5 * C) % C - 0.5 * C
119 1 tkerber
                d = (D**2).sum(axis=2)**0.5
120 1 tkerber
                e = q[i] * q / d
121 1 tkerber
                self.energy += np.dot(f, e).sum()
122 1 tkerber
                F = (e / d**2 * f[:, np.newaxis])[:, :, np.newaxis] * D
123 1 tkerber
                F[:, 0] -= (e.sum(axis=1) * dfdd / dOO)[:, np.newaxis] * DOO
124 1 tkerber
                self.forces[(a + 1) * 3:] += F.reshape((-1, 3))
125 1 tkerber
                self.forces[a * 3 + i] -= F.sum(axis=0).sum(axis=0)
126 1 tkerber
127 1 tkerber
        if world is not None:
128 1 tkerber
            self.energy = world.sum(self.energy)
129 1 tkerber
            world.sum(self.forces)
130 1 tkerber
131 1 tkerber
132 1 tkerber
class H2OConstraint:
133 1 tkerber
    """Constraint object for a rigid H2O molecule."""
134 1 tkerber
    def __init__(self, r=rOH, theta=thetaHOH, iterations=23, masses=None):
135 1 tkerber
        self.r = r
136 1 tkerber
        self.theta = theta
137 1 tkerber
        self.iterations = iterations
138 1 tkerber
        self.m = masses
139 1 tkerber
140 1 tkerber
    def set_masses(self, masses):
141 1 tkerber
        self.m = masses
142 1 tkerber
143 1 tkerber
    def adjust_positions(self, old, new):
144 1 tkerber
        bonds = [(0, 1, self.r), (0, 2, self.r)]
145 1 tkerber
        if self.theta:
146 1 tkerber
            bonds.append((1, 2, sin(self.theta / 2) * self.r * 2))
147 1 tkerber
        for iter in range(self.iterations):
148 1 tkerber
            for i, j, r in bonds:
149 1 tkerber
                D = old[i::3] - old[j::3]
150 1 tkerber
                m1 = self.m[i]
151 1 tkerber
                m2 = self.m[j]
152 1 tkerber
                a = new[i::3]
153 1 tkerber
                b = new[j::3]
154 1 tkerber
                B = a - b
155 1 tkerber
                x = (D**2).sum(axis=1)
156 1 tkerber
                y = (D * B).sum(axis=1)
157 1 tkerber
                z = (B**2).sum(axis=1) - r**2
158 1 tkerber
                k = m1 * m2 / (m1 + m2) * ((y**2 - x * z)**0.5 - y) / x
159 1 tkerber
                k.shape = (-1, 1)
160 1 tkerber
                a += k / m1 * D
161 1 tkerber
                b -= k / m2 * D
162 1 tkerber
163 1 tkerber
    def adjust_forces(self, positions, forces):
164 1 tkerber
        pass
165 1 tkerber
166 1 tkerber
    def copy(self):
167 1 tkerber
        return H2OConstraint(self.r, self.theta, self.iterations, self.m)
168 1 tkerber
169 1 tkerber
170 1 tkerber
class Verlet(MolecularDynamics):
171 1 tkerber
    def step(self, f):
172 1 tkerber
        atoms = self.atoms
173 1 tkerber
        m = atoms.get_masses()[:, np.newaxis]
174 1 tkerber
        v = self.atoms.get_velocities()
175 1 tkerber
        r0 = atoms.get_positions()
176 1 tkerber
        r = r0 + self.dt * v + self.dt**2 * f / m
177 1 tkerber
        atoms.set_positions(r)
178 1 tkerber
        r = atoms.get_positions()
179 1 tkerber
        v = (r - r0) / self.dt
180 1 tkerber
        self.atoms.set_velocities(v)
181 1 tkerber
        return atoms.get_forces()