File: test_nose_hoover_chain.py

package info (click to toggle)
python-ase 3.26.0-3
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 15,484 kB
  • sloc: python: 148,112; xml: 2,728; makefile: 110; javascript: 47
file content (309 lines) | stat: -rw-r--r-- 9,506 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
# fmt: off
from __future__ import annotations

from copy import deepcopy

import numpy as np
import pytest

import ase.build
import ase.units
from ase import Atoms
from ase.md.nose_hoover_chain import (
    MTKNPT,
    IsotropicMTKBarostat,
    IsotropicMTKNPT,
    MTKBarostat,
    NoseHooverChainNVT,
    NoseHooverChainThermostat,
)
from ase.md.velocitydistribution import MaxwellBoltzmannDistribution, Stationary


@pytest.fixture
def hcp_Cu() -> Atoms:
    atoms = ase.build.bulk(
        "Cu", crystalstructure='hcp', a=2.53, c=4.11
    ).repeat(2)
    return atoms


@pytest.mark.parametrize("tchain", [1, 3])
@pytest.mark.parametrize("tloop", [1, 3])
def test_thermostat_round_trip(hcp_Cu: Atoms, tchain: int, tloop: int):
    atoms = hcp_Cu.copy()

    timestep = 1.0 * ase.units.fs
    thermostat = NoseHooverChainThermostat(
        num_atoms_global=len(atoms),
        masses=atoms.get_masses()[:, None],
        temperature_K=1000,
        tdamp=100 * timestep,
        tchain=tchain,
        tloop=tloop,
    )

    rng = np.random.default_rng(0)
    p = rng.standard_normal(size=(len(atoms), 3))

    # Forward `n` steps and backward `n` steps with`, which should go back to
    # the initial state.
    n = 1000
    p_start = p.copy()
    eta_start = thermostat._eta.copy()
    p_eta_start = thermostat._p_eta.copy()

    for _ in range(n):
        p = thermostat.integrate_nhc(p, timestep)
    assert not np.allclose(p, p_start, atol=1e-6)
    assert not np.allclose(thermostat._eta, eta_start, atol=1e-6)
    assert not np.allclose(thermostat._p_eta, p_eta_start, atol=1e-6)

    for _ in range(n):
        p = thermostat.integrate_nhc(p, -timestep)
    assert np.allclose(p, p_start, atol=1e-6)

    # These values are apparently very machine-dependent:
    assert np.allclose(thermostat._eta, eta_start, atol=1e-5)
    assert np.allclose(thermostat._p_eta, p_eta_start, atol=1e-4)


@pytest.mark.parametrize("tchain", [1, 3])
@pytest.mark.parametrize("tloop", [1, 3])
def test_thermostat_truncation_error(hcp_Cu: Atoms, tchain: int, tloop: int):
    """Compare thermostat integration with delta by n steps and delta/2 by 2n
    steps. The difference between the two results should decrease with delta
    until reaching rounding error.
    """
    atoms = hcp_Cu.copy()

    delta0 = 1e-1
    n = 100
    m = 10

    list_p_diff = []
    list_eta_diff = []
    list_p_eta_diff = []
    for i in range(m):
        delta = delta0 * (2 ** -i)
        thermostat = NoseHooverChainThermostat(
            num_atoms_global=len(atoms),
            masses=atoms.get_masses()[:, None],
            temperature_K=1000,
            tdamp=100 * delta0,
            tchain=tchain,
            tloop=tloop,
        )

        rng = np.random.default_rng(0)
        p = rng.standard_normal(size=(len(atoms), 3))
        thermostat._eta = rng.standard_normal(size=(tchain, ))
        thermostat._p_eta = rng.standard_normal(size=(tchain, ))

        thermostat1 = deepcopy(thermostat)
        p1 = p.copy()
        for _ in range(n):
            p1 = thermostat1.integrate_nhc(p1, delta)

        thermostat2 = deepcopy(thermostat)
        p2 = p.copy()
        for _ in range(2 * n):
            p2 = thermostat2.integrate_nhc(p2, delta / 2)

        # O(delta^3) truncation error
        list_p_diff.append(np.linalg.norm(p1 - p2))
        list_eta_diff.append(
            np.linalg.norm(thermostat1._eta - thermostat2._eta)
        )
        list_p_eta_diff.append(
            np.linalg.norm(thermostat1._p_eta - thermostat2._p_eta)
        )

    print(np.array(list_p_diff))
    print(np.array(list_eta_diff))
    print(np.array(list_p_eta_diff))

    # Check that the differences decrease with delta until reaching rounding
    # error.
    eps = 1e-12
    for i in range(1, m):
        assert (
            (list_p_diff[i] < eps)
            or (list_p_diff[i] < list_p_diff[i - 1])
        )
        assert (
            (list_eta_diff[i] < eps)
            or (list_eta_diff[i] < list_eta_diff[i - 1])
        )
        assert (
            (list_p_eta_diff[i] < eps)
            or (list_p_eta_diff[i] < list_p_eta_diff[i - 1])
        )


@pytest.mark.parametrize("pchain", [1, 3])
@pytest.mark.parametrize("ploop", [1, 3])
def test_isotropic_barostat(asap3, hcp_Cu: Atoms, pchain: int, ploop: int):
    atoms = hcp_Cu.copy()
    atoms.calc = asap3.EMT()

    timestep = 1.0 * ase.units.fs
    barostat = IsotropicMTKBarostat(
        num_atoms_global=len(atoms),
        temperature_K=1000,
        pdamp=1000 * timestep,
        pchain=pchain,
        ploop=ploop,
    )

    rng = np.random.default_rng(0)
    p_eps = float(rng.standard_normal())

    # Forward `n` steps and backward `n` steps with`, which should go back to
    # the initial state.
    n = 1000
    p_eps_start = p_eps
    xi_start = barostat._xi.copy()
    p_xi_start = barostat._p_xi.copy()
    for _ in range(n):
        p_eps = barostat.integrate_nhc_baro(p_eps, timestep)
    assert not np.allclose(p_eps, p_eps_start, atol=1e-6)
    assert not np.allclose(barostat._xi, xi_start, atol=1e-6)
    assert not np.allclose(barostat._p_xi, p_xi_start, atol=1e-6)

    for _ in range(n):
        p_eps = barostat.integrate_nhc_baro(p_eps, -timestep)
    assert np.allclose(p_eps, p_eps_start, atol=1e-6)
    assert np.allclose(barostat._xi, xi_start, atol=1e-6)
    assert np.allclose(barostat._p_xi, p_xi_start, atol=1e-6)


@pytest.mark.parametrize("pchain", [1, 3])
@pytest.mark.parametrize("ploop", [1, 3])
def test_anisotropic_barostat(asap3, hcp_Cu: Atoms, pchain: int, ploop: int):
    atoms = hcp_Cu.copy()
    atoms.calc = asap3.EMT()

    timestep = 1.0 * ase.units.fs
    barostat = MTKBarostat(
        num_atoms_global=len(atoms),
        temperature_K=1000,
        pdamp=1000 * timestep,
        pchain=pchain,
    )

    rng = np.random.default_rng(0)
    p_g = rng.standard_normal((3, 3))
    p_g = 0.5 * (p_g + p_g.T)

    n = 1000
    p_g_start = p_g.copy()
    xi_start = barostat._xi.copy()
    p_xi_start = barostat._p_xi.copy()

    for _ in range(n):
        p_g = barostat.integrate_nhc_baro(p_g, timestep)
    # extended variables should be updated by n * timestep
    assert not np.allclose(p_g, p_g_start, atol=1e-6)
    assert not np.allclose(barostat._xi, xi_start, atol=1e-6)
    assert not np.allclose(barostat._p_xi, p_xi_start, atol=1e-6)

    for _ in range(n):
        p_g = barostat.integrate_nhc_baro(p_g, -timestep)
    # Now, the extended variables should be back to the initial state
    assert np.allclose(p_g, p_g_start, atol=1e-6)
    assert np.allclose(barostat._xi, xi_start, atol=1e-6)
    assert np.allclose(barostat._p_xi, p_xi_start, atol=1e-6)


@pytest.mark.parametrize("tchain", [1, 3])
def test_nose_hoover_chain_nvt(asap3, tchain: int):
    atoms = ase.build.bulk("Cu").repeat((2, 2, 2))
    atoms.calc = asap3.EMT()

    temperature_K = 300
    rng = np.random.default_rng(0)
    MaxwellBoltzmannDistribution(
        atoms,
        temperature_K=temperature_K, force_temp=True, rng=rng
    )
    Stationary(atoms)

    timestep = 1.0 * ase.units.fs
    md = NoseHooverChainNVT(
        atoms,
        timestep=timestep,
        temperature_K=temperature_K,
        tdamp=100 * timestep,
        tchain=tchain,
    )
    conserved_energy1 = md.get_conserved_energy()
    md.run(100)
    conserved_energy2 = md.get_conserved_energy()
    assert np.allclose(np.sum(atoms.get_momenta(), axis=0), 0.0)
    assert np.isclose(conserved_energy1, conserved_energy2, atol=1e-3)


@pytest.mark.parametrize("tchain", [1, 3])
@pytest.mark.parametrize("pchain", [1, 3])
def test_isotropic_mtk_npt(asap3, hcp_Cu: Atoms, tchain: int, pchain: int):
    atoms = hcp_Cu.copy()
    atoms.calc = asap3.EMT()

    temperature_K = 300
    rng = np.random.default_rng(0)
    MaxwellBoltzmannDistribution(
        atoms,
        temperature_K=temperature_K, force_temp=True, rng=rng
    )
    Stationary(atoms)

    timestep = 1.0 * ase.units.fs
    md = IsotropicMTKNPT(
        atoms,
        timestep=timestep,
        temperature_K=temperature_K,
        pressure_au=10.0 * ase.units.GPa,
        tdamp=100 * timestep,
        pdamp=1000 * timestep,
        tchain=tchain,
        pchain=pchain,
    )

    conserved_energy1 = md.get_conserved_energy()
    md.run(100)
    conserved_energy2 = md.get_conserved_energy()
    assert np.allclose(np.sum(atoms.get_momenta(), axis=0), 0.0)
    assert np.isclose(conserved_energy1, conserved_energy2, atol=1e-3)


@pytest.mark.parametrize("tchain", [1, 3])
@pytest.mark.parametrize("pchain", [1, 3])
def test_anisotropic_npt(asap3, hcp_Cu: Atoms, tchain: int, pchain: int):
    atoms = hcp_Cu.copy()
    atoms.calc = asap3.EMT()

    temperature_K = 300
    rng = np.random.default_rng(0)
    MaxwellBoltzmannDistribution(
        atoms,
        temperature_K=temperature_K, force_temp=True, rng=rng
    )
    Stationary(atoms)

    timestep = 1.0 * ase.units.fs
    md = MTKNPT(
        atoms,
        timestep=timestep,
        temperature_K=temperature_K,
        pressure_au=10.0 * ase.units.GPa,
        tdamp=100 * timestep,
        pdamp=1000 * timestep,
    )
    conserved_energy1 = md.get_conserved_energy()
    positions1 = atoms.get_positions().copy()
    md.run(100)
    conserved_energy2 = md.get_conserved_energy()
    assert np.allclose(np.sum(atoms.get_momenta(), axis=0), 0.0)
    assert np.isclose(conserved_energy1, conserved_energy2, atol=1e-3)
    assert not np.allclose(atoms.get_positions(), positions1, atol=1e-6)