-
Notifications
You must be signed in to change notification settings - Fork 2
/
constrainmol.py
210 lines (178 loc) · 6.27 KB
/
constrainmol.py
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
import pyomo.environ as pyo
import parmed
import numpy as np
class ConstrainedMolecule(object):
def __init__(self, structure):
"""Initialize the ConstrainedMolecule from a parmed.Structure
Parameters
----------
structure: parmed.Structure
parmed structure with coordinates and bond lengths
Notes
-----
The pyomo model is created on initialization
"""
if not isinstance(structure, parmed.Structure):
raise TypeError(
f"structure: {structure} must be a parmed.Structure"
)
if len(structure.bonds) == 0:
raise ValueError(
f"structure: {structure} contains no bonds"
)
# Extract coordinates and bonds
xyz = structure.coordinates
constraints = {}
for bond in structure.bonds:
idx1 = bond.atom1.idx
idx2 = bond.atom2.idx
constraints[(idx1, idx2)] = bond.type.req
self.structure = parmed.structure.copy(structure)
self.model = self._create_model(xyz, constraints)
self.model_solved = False
def _create_model(self, xyz, constraints):
"""Create a pyomo model to make xyz satisfy bond constraints
Parameters
----------
xyz: np.ndarray, shape=(n_atoms, 3)
initial xyz coordinates
constraints: dict, keys=(idx1, idx2), values=bond length
dictionary with bond length constraints
Returns
-------
pyomo.ConcreteModel
the pyomo model
"""
# Create pyomo model
m = pyo.ConcreteModel()
# Create list of atom indicies
ids = range(xyz.shape[0])
m.atom_ids = pyo.Set(initialize=ids)
# Create target values
m.x_start = pyo.Param(
m.atom_ids,
initialize={idx: x for idx, x in zip(ids, xyz[:, 0])},
within=pyo.Reals,
mutable=True,
)
m.y_start = pyo.Param(
m.atom_ids,
initialize={idx: y for idx, y in zip(ids, xyz[:, 1])},
within=pyo.Reals,
mutable=True,
)
m.z_start = pyo.Param(
m.atom_ids,
initialize={idx: z for idx, z in zip(ids, xyz[:, 2])},
within=pyo.Reals,
mutable=True,
)
# Create position variables
m.x = pyo.Var(
m.atom_ids,
initialize={idx: x for idx, x in zip(ids, xyz[:, 0])},
within=pyo.Reals
)
m.y = pyo.Var(
m.atom_ids,
initialize={idx: y for idx, y in zip(ids, xyz[:, 1])},
within=pyo.Reals
)
m.z = pyo.Var(
m.atom_ids,
initialize={idx: z for idx, z in zip(ids, xyz[:, 2])},
within=pyo.Reals
)
# Create bond length parameter
m.bond_lengths = pyo.Param(
m.atom_ids, m.atom_ids, initialize=constraints
)
# Add bond length constraints
m.eq_calc_bound_length = pyo.Constraint(
m.atom_ids, m.atom_ids, rule=self._calc_bond_length
)
m.obj = pyo.Objective(
expr=sum(
(m.x[i] - m.x_start[i])**2 +
(m.y[i] - m.y_start[i])**2 +
(m.z[i] - m.z_start[i])**2
for i in m.atom_ids
)
)
return m
@staticmethod
def _calc_bond_length(m, i, j):
if (i, j) in m.bond_lengths.keys():
return (
(m.x[i] - m.x[j])**2 +
(m.y[i] - m.y[j])**2 +
(m.z[i] - m.z[j])**2 ==
m.bond_lengths[(i, j)]**2
)
else:
return pyo.Constraint.Skip
def solve(self, verbose=False):
"""Solve the pyomo model to find the constrained coordinates
Parameters
----------
verbose : boolean, optional, default=False
print the solver output to screen
Notes
-----
Updates ConstrainedMolecule.structure.coordinates if solve is successful
and sets ConstrainMolecule.model_solved to True
"""
result = pyo.SolverFactory('ipopt').solve(self.model, tee=verbose)
success = (
str(result["Solver"][0]["Termination condition"]) == 'optimal'
)
if not success:
raise ValueError(
"Optimal solution not found. You may want to run "
"'solve' with 'verbose=True' to see the solver output."
)
constrained_xyz = np.zeros((len(self.model.atom_ids), 3))
for idx in self.model.atom_ids:
constrained_xyz[idx, 0] = self.model.x[idx].value
constrained_xyz[idx, 1] = self.model.y[idx].value
constrained_xyz[idx, 2] = self.model.z[idx].value
self.structure.coordinates = constrained_xyz
self.model_solved = True
def update_xyz(self, xyz):
"""
Update the initial unconstrained coordinates
Parameters
----------
xyz: np.ndarray, shape=(n_atoms, 3)
the new coordinates
Returns
-------
"""
try:
xyz = np.array(xyz, dtype=np.float64)
except ValueError:
raise TypeError(f"New xyz: {xyz} must be array-like")
if xyz.shape != self.structure.coordinates.shape:
raise ValueError(
f"New xyz shape: {xyz.shape} does not match the coordinates "
f"used to instantiate the class: "
f"{self.structure.coordinates.shape}."
)
self.structure.coordinates = xyz
for idx in range(xyz.shape[0]):
self.model.x_start[idx] = xyz[idx, 0]
self.model.y_start[idx] = xyz[idx, 1]
self.model.z_start[idx] = xyz[idx, 2]
self.model.x[idx] = xyz[idx, 0]
self.model.y[idx] = xyz[idx, 1]
self.model.z[idx] = xyz[idx, 2]
self.model_solved = False
@property
def xyz(self):
return self.structure.coordinates
@xyz.setter
def xyz(self, xyz):
raise ValueError(
"Coordinates must be updated with ConstrainedMolecule.update_xyz. "
"See help(ConstrainedMolecule.update_xyz) for more information."
)