"platforms/reference/tests/TestReferenceCustomCVForce.cpp" did not exist on "bc025a7be08f8fdb445fca7cf9debacc5aaa1692"
ReferenceCCMAAlgorithm.cpp 14 KB
Newer Older
1

2
/* Portions copyright (c) 2006-2015 Stanford University and Simbios.
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
 * Contributors: Peter Eastman, Pande Group
 *
 * Permission is hereby granted, free of charge, to any person obtaining
 * a copy of this software and associated documentation files (the
 * "Software"), to deal in the Software without restriction, including
 * without limitation the rights to use, copy, modify, merge, publish,
 * distribute, sublicense, and/or sell copies of the Software, and to
 * permit persons to whom the Software is furnished to do so, subject
 * to the following conditions:
 *
 * The above copyright notice and this permission notice shall be included
 * in all copies or substantial portions of the Software.
 *
 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
 * OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
 * MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
 * IN NO EVENT SHALL THE AUTHORS, CONTRIBUTORS OR COPYRIGHT HOLDERS BE
 * LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION
 * OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION
 * WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
 */

#include <string.h>
#include <sstream>

28
#include "SimTKOpenMMUtilities.h"
29
#include "ReferenceCCMAAlgorithm.h"
30
#include "ReferenceDynamics.h"
31
#include "quern.h"
32
#include "openmm/OpenMMException.h"
33
#include "openmm/Vec3.h"
34
#include "openmm/internal/ThreadPool.h"
35
36
37
#include <map>

using std::map;
38
using std::pair;
39
40
using std::vector;
using std::set;
41
using namespace OpenMM;
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
// This class extracts columns from the inverse matrix one at a time.  It is done in parallel,
// since this can be very slow.

class ExtractMatrixTask : public ThreadPool::Task {
public:
    ExtractMatrixTask(int numConstraints, vector<vector<pair<int, RealOpenMM> > >& matrix, const vector<RealOpenMM>& distance, RealOpenMM elementCutoff,
                      const int* qRowStart, const int* qColIndex, const int* rRowStart, const int* rColIndex, const double* qValue, const double* rValue) :
                numConstraints(numConstraints), matrix(matrix), distance(distance), elementCutoff(elementCutoff), qRowStart(qRowStart), qColIndex(qColIndex),
                rRowStart(rRowStart), rColIndex(rColIndex), qValue(qValue), rValue(rValue) {
    }
    
    void execute(ThreadPool& pool, int threadIndex) {
        vector<double> rhs(numConstraints);
        for (int i = threadIndex; i < numConstraints; i += pool.getNumThreads()) {
            // Extract column i of the inverse matrix.

            for (int j = 0; j < numConstraints; j++)
                rhs[j] = (i == j ? 1.0 : 0.0);
            QUERN_multiply_with_q_transpose(numConstraints, qRowStart, qColIndex, qValue, &rhs[0]);
            QUERN_solve_with_r(numConstraints, rRowStart, rColIndex, rValue, &rhs[0], &rhs[0]);
            for (int j = 0; j < numConstraints; j++) {
                double value = rhs[j]*distance[j]/distance[i];
                if (FABS((RealOpenMM) value) > elementCutoff)
                    matrix[i].push_back(pair<int, RealOpenMM>(j, (RealOpenMM) value));
            }
        }
    }
private:
    int numConstraints;
    vector<vector<pair<int, RealOpenMM> > >& matrix;
    const vector<RealOpenMM>& distance;
    RealOpenMM elementCutoff;
    const int *qRowStart, *qColIndex, *rRowStart, *rColIndex;
    const double *qValue, *rValue;
};

79
ReferenceCCMAAlgorithm::ReferenceCCMAAlgorithm(int numberOfAtoms,
80
81
82
83
84
85
                                               int numberOfConstraints,
                                               const vector<pair<int, int> >& atomIndices,
                                               const vector<RealOpenMM>& distance,
                                               vector<RealOpenMM>& masses,
                                               vector<AngleInfo>& angles,
                                               RealOpenMM elementCutoff) {
86
    _numberOfConstraints = numberOfConstraints;
87
    _elementCutoff = elementCutoff;
88
89
    _atomIndices = atomIndices;
    _distance = distance;
90

91
92
    _maximumNumberOfIterations = 150;
    _hasInitializedMasses = false;
93

94
    // work arrays
95
96

    if (_numberOfConstraints > 0) {
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
        _r_ij.resize(numberOfConstraints);
        _d_ij2 = SimTKOpenMMUtilities::allocateOneDRealOpenMMArray(numberOfConstraints, NULL, 1, 0.0, "dij_2");
        _distanceTolerance = SimTKOpenMMUtilities::allocateOneDRealOpenMMArray(numberOfConstraints, NULL, 1, 0.0, "distanceTolerance");
        _reducedMasses = SimTKOpenMMUtilities::allocateOneDRealOpenMMArray(numberOfConstraints, NULL, 1, 0.0, "reducedMasses");
    }
    if (numberOfConstraints > 0)
    {
        // Compute the constraint coupling matrix

        vector<vector<int> > atomAngles(numberOfAtoms);
        for (int i = 0; i < (int) angles.size(); i++)
            atomAngles[angles[i].atom2].push_back(i);
        vector<vector<pair<int, double> > > matrix(numberOfConstraints);
        for (int j = 0; j < numberOfConstraints; j++) {
            for (int k = 0; k < numberOfConstraints; k++) {
                if (j == k) {
                    matrix[j].push_back(pair<int, double>(j, 1.0));
                    continue;
                }
                double scale;
                int atomj0 = _atomIndices[j].first;
                int atomj1 = _atomIndices[j].second;
                int atomk0 = _atomIndices[k].first;
                int atomk1 = _atomIndices[k].second;
                RealOpenMM invMass0 = 1/masses[atomj0];
                RealOpenMM invMass1 = 1/masses[atomj1];
                int atoma, atomb, atomc;
                if (atomj0 == atomk0) {
                    atoma = atomj1;
                    atomb = atomj0;
                    atomc = atomk1;
                    scale = invMass0/(invMass0+invMass1);
                }
                else if (atomj1 == atomk1) {
                    atoma = atomj0;
                    atomb = atomj1;
                    atomc = atomk0;
                    scale = invMass1/(invMass0+invMass1);
                }
                else if (atomj0 == atomk1) {
                    atoma = atomj1;
                    atomb = atomj0;
                    atomc = atomk0;
                    scale = invMass0/(invMass0+invMass1);
                }
                else if (atomj1 == atomk0) {
                    atoma = atomj0;
                    atomb = atomj1;
                    atomc = atomk1;
                    scale = invMass1/(invMass0+invMass1);
                }
                else
                    continue; // These constraints are not connected.

                // Look for a third constraint forming a triangle with these two.

                bool foundConstraint = false;
                for (int other = 0; other < numberOfConstraints; other++) {
                    if ((_atomIndices[other].first == atoma && _atomIndices[other].second == atomc) || (_atomIndices[other].first == atomc && _atomIndices[other].second == atoma)) {
                        double d1 = _distance[j];
                        double d2 = _distance[k];
                        double d3 = _distance[other];
                        matrix[j].push_back(pair<int, double>(k, scale*(d1*d1+d2*d2-d3*d3)/(2.0*d1*d2)));
                        foundConstraint = true;
                        break;
                    }
                }
                if (!foundConstraint) {
                    // We didn't find one, so look for an angle force field term.

                    const vector<int>& angleCandidates = atomAngles[atomb];
                    for (vector<int>::const_iterator iter = angleCandidates.begin(); iter != angleCandidates.end(); iter++) {
                        const AngleInfo& angle = angles[*iter];
                        if ((angle.atom1 == atoma && angle.atom3 == atomc) || (angle.atom3 == atoma && angle.atom1 == atomc)) {
                            matrix[j].push_back(pair<int, double>(k, scale*cos(angle.angle)));
                            break;
                        }
                    }
                }
            }
        }

        // Invert it using QR.

        vector<int> matrixRowStart;
        vector<int> matrixColIndex;
        vector<double> matrixValue;
        for (int i = 0; i < numberOfConstraints; i++) {
            matrixRowStart.push_back(matrixValue.size());
            for (int j = 0; j < (int) matrix[i].size(); j++) {
                pair<int, double> element = matrix[i][j];
                matrixColIndex.push_back(element.first);
                matrixValue.push_back(element.second);
            }
        }
        matrixRowStart.push_back(matrixValue.size());
        int *qRowStart, *qColIndex, *rRowStart, *rColIndex;
        double *qValue, *rValue;
        QUERN_compute_qr(numberOfConstraints, numberOfConstraints, &matrixRowStart[0], &matrixColIndex[0], &matrixValue[0], NULL,
                &qRowStart, &qColIndex, &qValue, &rRowStart, &rColIndex, &rValue);
        vector<double> rhs(numberOfConstraints);
        _matrix.resize(numberOfConstraints);
199
200
201
202
        ThreadPool threads;
        ExtractMatrixTask task(numberOfConstraints, _matrix, _distance, _elementCutoff, qRowStart, qColIndex, rRowStart, rColIndex, qValue, rValue);
        threads.execute(task);
        threads.waitForThreads();
203
204
        QUERN_free_result(qRowStart, qColIndex, qValue);
        QUERN_free_result(rRowStart, rColIndex, rValue);
205
206
207
    }
}

208
209
210
211
212
213
ReferenceCCMAAlgorithm::~ReferenceCCMAAlgorithm() {
    if (_numberOfConstraints > 0) {
        SimTKOpenMMUtilities::freeOneDRealOpenMMArray(_d_ij2, "d_ij2");
        SimTKOpenMMUtilities::freeOneDRealOpenMMArray(_distanceTolerance, "distanceTolerance");
        SimTKOpenMMUtilities::freeOneDRealOpenMMArray(_reducedMasses, "reducedMasses");
    }
214
215
}

216
217
int ReferenceCCMAAlgorithm::getNumberOfConstraints() const {
    return _numberOfConstraints;
218
219
}

220
221
int ReferenceCCMAAlgorithm::getMaximumNumberOfIterations() const {
    return _maximumNumberOfIterations;
222
223
}

224
225
void ReferenceCCMAAlgorithm::setMaximumNumberOfIterations(int maximumNumberOfIterations) {
    _maximumNumberOfIterations = maximumNumberOfIterations;
226
227
}

228
void ReferenceCCMAAlgorithm::apply(vector<RealVec>& atomCoordinates,
229
                                         vector<RealVec>& atomCoordinatesP,
230
231
                                         vector<RealOpenMM>& inverseMasses, RealOpenMM tolerance) {
    applyConstraints(atomCoordinates, atomCoordinatesP, inverseMasses, false, tolerance);
232
233
}

234
void ReferenceCCMAAlgorithm::applyToVelocities(std::vector<OpenMM::RealVec>& atomCoordinates,
235
236
               std::vector<OpenMM::RealVec>& velocities, std::vector<RealOpenMM>& inverseMasses, RealOpenMM tolerance) {
    applyConstraints(atomCoordinates, velocities, inverseMasses, true, tolerance);
237
238
}

239
void ReferenceCCMAAlgorithm::applyConstraints(vector<RealVec>& atomCoordinates,
240
                                         vector<RealVec>& atomCoordinatesP,
241
                                         vector<RealOpenMM>& inverseMasses, bool constrainingVelocities, RealOpenMM tolerance) {
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
    // temp arrays

    vector<RealVec>& r_ij = _r_ij;
    RealOpenMM* d_ij2 = _d_ij2;
    RealOpenMM* reducedMasses = _reducedMasses;

    // calculate reduced masses on 1st pass

    if (!_hasInitializedMasses) {
        _hasInitializedMasses = true;
        for (int ii = 0; ii < _numberOfConstraints; ii++) {
           int atomI = _atomIndices[ii].first;
           int atomJ = _atomIndices[ii].second;
           reducedMasses[ii] = 0.5/(inverseMasses[atomI] + inverseMasses[atomJ]);
        }
    }
258

259
    // setup: r_ij for each (i,j) constraint
260

261
262
263
264
265
266
    for (int ii = 0; ii < _numberOfConstraints; ii++) {
        int atomI = _atomIndices[ii].first;
        int atomJ = _atomIndices[ii].second;
        r_ij[ii] = atomCoordinates[atomI] - atomCoordinates[atomJ];
        d_ij2[ii] = r_ij[ii].dot(r_ij[ii]);
    }
267
268
    RealOpenMM lowerTol = 1-2*tolerance+tolerance*tolerance;
    RealOpenMM upperTol = 1+2*tolerance+tolerance*tolerance;
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284

    // main loop

    int iterations = 0;
    int numberConverged = 0;
    vector<RealOpenMM> constraintDelta(_numberOfConstraints);
    vector<RealOpenMM> tempDelta(_numberOfConstraints);
    while (iterations < getMaximumNumberOfIterations()) {
        numberConverged  = 0;
        for (int ii = 0; ii < _numberOfConstraints; ii++) {
            int atomI = _atomIndices[ii].first;
            int atomJ = _atomIndices[ii].second;
            RealVec rp_ij = atomCoordinatesP[atomI] - atomCoordinatesP[atomJ];
            if (constrainingVelocities) {
                RealOpenMM rrpr = rp_ij.dot(r_ij[ii]);
                constraintDelta[ii] = -2*reducedMasses[ii]*rrpr/d_ij2[ii];
285
                if (fabs(constraintDelta[ii]) <= tolerance)
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
                    numberConverged++;
            }
            else {
                RealOpenMM rp2  = rp_ij.dot(rp_ij);
                RealOpenMM dist2 = _distance[ii]*_distance[ii];
                RealOpenMM diff = dist2 - rp2;
                constraintDelta[ii] = 0;
                RealOpenMM rrpr = DOT3(rp_ij, r_ij[ii]);
                constraintDelta[ii] = reducedMasses[ii]*diff/rrpr;
                if (rp2 >= lowerTol*dist2 && rp2 <= upperTol*dist2)
                    numberConverged++;
            }
        }
        if (numberConverged == _numberOfConstraints)
            break;
        iterations++;

        if (_matrix.size() > 0) {
            for (int i = 0; i < _numberOfConstraints; i++) {
                RealOpenMM sum = 0.0;
                for (int j = 0; j < (int) _matrix[i].size(); j++) {
                    pair<int, RealOpenMM> element = _matrix[i][j];
                    sum += element.second*constraintDelta[element.first];
                }
                tempDelta[i] = sum;
            }
            constraintDelta = tempDelta;
        }
        for (int ii = 0; ii < _numberOfConstraints; ii++) {
            int atomI = _atomIndices[ii].first;
            int atomJ = _atomIndices[ii].second;
            RealVec dr = r_ij[ii]*constraintDelta[ii];
            atomCoordinatesP[atomI] += dr*inverseMasses[atomI];
            atomCoordinatesP[atomJ] -= dr*inverseMasses[atomJ];
        }
    }
322
}
323
324
325
326

const vector<vector<pair<int, RealOpenMM> > >& ReferenceCCMAAlgorithm::getMatrix() const {
    return _matrix;
}