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
|
def test_parameteric_constr():
import numpy as np
from ase.build import bulk
from ase.constraints import (
dict2constraint,
FixScaledParametricRelations,
FixCartesianParametricRelations,
)
from ase.calculators.emt import EMT
# Build the atoms object and attach a calculator
a = bulk("Ni", cubic=True)
a.calc = EMT()
# Get adjusted cell
cell = a.cell + 0.01
# Generate lattice constraint
param_lat = ["a"]
expr_lat = [
"a", "0", "0",
"0", "a", "0",
"0", "0", "a",
]
constr_lat = FixCartesianParametricRelations.from_expressions(
indices=[0, 1, 2],
params=param_lat,
expressions=expr_lat,
use_cell=True,
)
# Check expression generator
for const_expr, passed_expr in zip(constr_lat.expressions.flatten(), expr_lat):
assert const_expr == passed_expr
# Check adjust_cell
constr_lat.adjust_cell(a, cell)
# Check serialization and construction from dict
constr_lat_dict = constr_lat.todict()
dict2constraint(constr_lat_dict)
cell_diff = (cell - a.cell).flatten()
expected_cell_diff = np.array([0.01, 0.0, 0.0, 0.0, 0.01, 0.0, 0.0, 0.0, 0.01])
assert np.max(np.abs(cell_diff - expected_cell_diff)) < 1e-12
# Check adjust_stress
a.cell += 0.01
stress = a.get_stress().copy()
constr_lat.adjust_stress(a, stress)
stress_rat = stress / a.get_stress()
assert np.max(np.abs(stress_rat - np.array([1., 1., 1., 0., 0., 0.]))) < 1e-12
# Reset cell
a.cell -= 0.01
# Get adjusted cell/positions for the system
pos = a.get_positions().copy() + 0.01
# Generate proper atomic constraints
constr_atom = FixScaledParametricRelations(
[0, 1, 2, 3],
np.ndarray((12, 0)),
a.get_scaled_positions().flatten(),
)
# Check serialization and construction from dict
constr_atom_dict = constr_atom.todict()
dict2constraint(constr_atom_dict)
# Check adjust_positions
constr_atom.adjust_positions(a, pos)
assert np.max(np.abs(a.get_positions() - pos)) < 1e-12
# Check adjust_forces
assert np.max(np.abs(a.get_forces())) < 1e-12
# Check non-empty constraint
param_atom = ["dis"]
expr_atom = [
"dis", "dis", "dis",
"dis", "-0.5", "0.5",
"0.5", "dis", "0.5",
"0.5", "0.5", "dis",
]
constr_atom = FixScaledParametricRelations.from_expressions(
indices=[0, 1, 2, 3],
params=param_atom,
expressions=expr_atom,
)
# Restart position adjustment
pos += 0.01 * a.cell[0, 0]
# Check adjust_positions
constr_atom.adjust_positions(a, pos)
scaled_pos = a.cell.scaled_positions(pos)
pos_diff = (scaled_pos - a.get_scaled_positions()).flatten()
expected_pos_diff = np.array(
[0.01, 0.01, 0.01, 0.01, 0.0, 0.0, 0.0, 0.01, 0.0, 0.0, 0.0, 0.01]
)
assert np.max(np.abs(pos_diff - expected_pos_diff)) < 1e-12
# Check adjust_forces
a.set_positions(pos + 0.3)
forces = a.get_forces()
constr_atom.adjust_forces(a, forces)
forces_rat = forces / a.get_forces()
assert np.max(np.abs(forces_rat.flatten() / 100.0 - expected_pos_diff)) < 1e-12
# Check auto-remapping/expression generation, the -0.5 should now be 0.5
expr_atom[4] = "0.5"
current_expression = constr_atom.expressions.flatten()
for const_expr, passed_expr in zip(current_expression, expr_atom):
assert const_expr == passed_expr
# Check with Cartesian parametric constraints now
expr_atom = [
"dis", "dis", "dis",
"dis", "1.76", "1.76",
"1.76", "dis", "1.76",
"1.76", "1.76", "dis",
]
constr_atom = FixCartesianParametricRelations.from_expressions(
indices=[0, 1, 2, 3],
params=param_atom,
expressions=expr_atom,
)
# Restart position adjustment
a.set_positions(pos)
pos += 0.01
# Check adjust_positions
constr_atom.adjust_positions(a, pos)
pos_diff = (pos - a.get_positions()).flatten()
expected_pos_diff = np.array(
[0.01, 0.01, 0.01, 0.01, 0.0, 0.0, 0.0, 0.01, 0.0, 0.0, 0.0, 0.01]
)
assert np.max(np.abs(pos_diff - expected_pos_diff)) < 1e-12
# Check adjust_forces
a.set_positions(pos + 0.3)
forces = a.get_forces()
constr_atom.adjust_forces(a, forces)
forces_rat = forces / a.get_forces()
assert np.max(np.abs(forces_rat.flatten() / 100.0 - expected_pos_diff)) < 1e-12
|