ReferenceCCMAAlgorithm.cpp 20.7 KB
Newer Older
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

/* Portions copyright (c) 2006-2009 Stanford University and Simbios.
 * 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>

#include "../SimTKUtilities/SimTKOpenMMCommon.h"
#include "../SimTKUtilities/SimTKOpenMMUtilities.h"
#include "../SimTKUtilities/SimTKOpenMMLog.h"
31
#include "ReferenceCCMAAlgorithm.h"
32
#include "ReferenceDynamics.h"
33
#include "quern.h"
34
35
36
37
#include "openmm/Vec3.h"
#include <map>

using std::map;
38
using std::pair;
39
40
41
42
43
44
using std::vector;
using std::set;
using OpenMM::Vec3;

/**---------------------------------------------------------------------------------------

45
   ReferenceCCMAAlgorithm constructor
46

47
48
49
50
51
52
53
         @param numberOfAtoms    number of atoms
         @param numberOfConstraints      number of constraints
         @param atomIndices              atom indices for contraints
         @param distance                 distances for constraints
         @param masses                   atom masses
         @param angles                   angle force field terms
         @param tolerance                constraint tolerance
54
55
56

   --------------------------------------------------------------------------------------- */

57
ReferenceCCMAAlgorithm::ReferenceCCMAAlgorithm( int numberOfAtoms,
58
59
60
                                                  int numberOfConstraints,
                                                  int** atomIndices,
                                                  RealOpenMM* distance,
61
                                                  RealOpenMM* masses,
62
                                                  vector<AngleInfo>& angles,
63
64
65
66
                                                  RealOpenMM tolerance){

   // ---------------------------------------------------------------------------------------

67
   // static const char* methodName = "\nReferenceCCMAAlgorithm::ReferenceCCMAAlgorithm";
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

   static const RealOpenMM zero        =  0.0;
   static const RealOpenMM one         =  1.0;
   static const RealOpenMM two         =  2.0;
   static const RealOpenMM three       =  3.0;
   static const RealOpenMM oneM        = -1.0;

   static const int threeI             =  3;

   // ---------------------------------------------------------------------------------------

   _numberOfConstraints        = numberOfConstraints;
   _atomIndices                = atomIndices;
   _distance                   = distance;

   _maximumNumberOfIterations  = 150;
   _tolerance                  = tolerance;
   _hasInitializedMasses       = false;

   // work arrays

   if (_numberOfConstraints > 0) {
       _r_ij                       = SimTKOpenMMUtilities::allocateTwoDRealOpenMMArray( numberOfConstraints, threeI, NULL,
                                                                                        1, zero, "r_ij" );

       _d_ij2                      = SimTKOpenMMUtilities::allocateOneDRealOpenMMArray( numberOfConstraints, NULL, 1, zero, "dij_2" );
       _distanceTolerance          = SimTKOpenMMUtilities::allocateOneDRealOpenMMArray( numberOfConstraints, NULL, 1, zero, "distanceTolerance" );
       _reducedMasses              = SimTKOpenMMUtilities::allocateOneDRealOpenMMArray( numberOfConstraints, NULL, 1, zero, "reducedMasses" );
   }
97
98
99
100
101
102
103
104
105
106
107
108
109
   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;
110
               }
111
               double scale;
112
113
114
115
               int atomj0 = _atomIndices[j][0];
               int atomj1 = _atomIndices[j][1];
               int atomk0 = _atomIndices[k][0];
               int atomk1 = _atomIndices[k][1];
116
117
               RealOpenMM invMass0 = one/masses[atomj0];
               RealOpenMM invMass1 = one/masses[atomj1];
118
               int atoma, atomb, atomc;
119
120
               if (atomj0 == atomk0) {
                   atoma = atomj1;
121
122
                   atomb = atomj0;
                   atomc = atomk1;
123
124
125
126
                   scale = invMass0/(invMass0+invMass1);
               }
               else if (atomj1 == atomk1) {
                   atoma = atomj0;
127
128
                   atomb = atomj1;
                   atomc = atomk0;
129
130
131
132
                   scale = invMass1/(invMass0+invMass1);
               }
               else if (atomj0 == atomk1) {
                   atoma = atomj1;
133
134
                   atomb = atomj0;
                   atomc = atomk0;
135
136
137
138
                   scale = invMass0/(invMass0+invMass1);
               }
               else if (atomj1 == atomk0) {
                   atoma = atomj0;
139
140
                   atomb = atomj1;
                   atomc = atomk1;
141
142
                   scale = invMass1/(invMass0+invMass1);
               }
143
               else
144
145
                   continue; // These constraints are not connected.

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

148
149
150
151
152
               bool foundConstraint = false;
               for (int other = 0; other < numberOfConstraints; other++) {
                   if ((_atomIndices[other][0] == atoma && _atomIndices[other][1] == atomc) || (_atomIndices[other][0] == atomc && _atomIndices[other][1] == atoma)) {
                       double d1 = _distance[j];
                       double d2 = _distance[k];
153
                       double d3 = _distance[other];
154
155
                       matrix[j].push_back(pair<int, double>(k, scale*(d1*d1+d2*d2-d3*d3)/(2.0*d1*d2)));
                       foundConstraint = true;
156
157
158
                       break;
                   }
               }
159
160
161
162
163
164
165
166
167
168
169
170
               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;
                       }
                   }
               }
171
172
173
           }
       }

174
175
176
177
178
179
180
181
182
183
184
       // 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);
185
186
           }
       }
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
       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);
       for (int i = 0; i < numberOfConstraints; i++) {
           // Extract column i of the inverse matrix.

           for (int j = 0; j < numberOfConstraints; j++)
               rhs[j] = (i == j ? 1.0 : 0.0);
           QUERN_multiply_with_q_transpose(numberOfConstraints, qRowStart, qColIndex, qValue, &rhs[0]);
           QUERN_solve_with_r(numberOfConstraints, rRowStart, rColIndex, rValue, &rhs[0], &rhs[0]);
           for (int j = 0; j < numberOfConstraints; j++) {
               double value = rhs[j]*_distance[i]/_distance[j];
203
               if (FABS((RealOpenMM)value) > 0.02)
204
205
206
207
208
209
                   _matrix[j].push_back(pair<int, RealOpenMM>(i, (RealOpenMM) value));
           }
       }
       QUERN_free_result(qRowStart, qColIndex, qValue);
       QUERN_free_result(rRowStart, rColIndex, rValue);
   }
210
211
212
213
}

/**---------------------------------------------------------------------------------------

214
   ReferenceCCMAAlgorithm destructor
215
216
217

   --------------------------------------------------------------------------------------- */

218
ReferenceCCMAAlgorithm::~ReferenceCCMAAlgorithm( ){
219
220
221

   // ---------------------------------------------------------------------------------------

222
   // static const char* methodName = "\nReferenceCCMAAlgorithm::~ReferenceCCMAAlgorithm";
223
224
225
226
227
228
229
230
231
232

   // ---------------------------------------------------------------------------------------

    if (_numberOfConstraints > 0) {
       SimTKOpenMMUtilities::freeTwoDRealOpenMMArray( _r_ij,  "r_ij" );

       SimTKOpenMMUtilities::freeOneDRealOpenMMArray( _d_ij2, "d_ij2" );
       SimTKOpenMMUtilities::freeOneDRealOpenMMArray( _distanceTolerance, "distanceTolerance" );
       SimTKOpenMMUtilities::freeOneDRealOpenMMArray( _reducedMasses, "reducedMasses" );
    }
233
234
//    for (unsigned int i = 0; i < _matrices.size(); i++)
//        SimTKOpenMMUtilities::freeTwoDRealOpenMMArray(_matrices[i], "");
235
236
237
238
239
240
241
242
243
244
}

/**---------------------------------------------------------------------------------------

   Get number of constraints

   @return number of constraints

   --------------------------------------------------------------------------------------- */

245
int ReferenceCCMAAlgorithm::getNumberOfConstraints( void ) const {
246
247
248

   // ---------------------------------------------------------------------------------------

249
   // static const char* methodName = "\nReferenceCCMAAlgorithm::getNumberOfConstraints";
250
251
252
253
254
255
256
257
258
259
260
261
262
263

   // ---------------------------------------------------------------------------------------

   return _numberOfConstraints;
}

/**---------------------------------------------------------------------------------------

   Get maximum number of iterations

   @return maximum number of iterations

   --------------------------------------------------------------------------------------- */

264
int ReferenceCCMAAlgorithm::getMaximumNumberOfIterations( void ) const {
265
266
267

   // ---------------------------------------------------------------------------------------

268
   // static const char* methodName = "\nReferenceCCMAAlgorithm::getMaximumNumberOfIterations";
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284

   // ---------------------------------------------------------------------------------------

   return _maximumNumberOfIterations;
}

/**---------------------------------------------------------------------------------------

   Set maximum number of iterations

   @param maximumNumberOfIterations   new maximum number of iterations

   @return ReferenceDynamics::DefaultReturn

   --------------------------------------------------------------------------------------- */

285
int ReferenceCCMAAlgorithm::setMaximumNumberOfIterations( int maximumNumberOfIterations ){
286
287
288

   // ---------------------------------------------------------------------------------------

289
   // static const char* methodName = "\nReferenceCCMAAlgorithm::setMaximumNumberOfIterations";
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305

   // ---------------------------------------------------------------------------------------

   _maximumNumberOfIterations = maximumNumberOfIterations;

   return ReferenceDynamics::DefaultReturn;
}

/**---------------------------------------------------------------------------------------

   Get tolerance

   @return tolerance

   --------------------------------------------------------------------------------------- */

306
RealOpenMM ReferenceCCMAAlgorithm::getTolerance( void ) const {
307
308
309

   // ---------------------------------------------------------------------------------------

310
   // static const char* methodName = "\nReferenceCCMAAlgorithm::getTolerance";
311
312
313
314
315
316
317
318
319
320
321
322
323
324

   // ---------------------------------------------------------------------------------------

   return _tolerance;
}

/**---------------------------------------------------------------------------------------

   Set tolerance

   @param tolerance new tolerance

   --------------------------------------------------------------------------------------- */

325
void ReferenceCCMAAlgorithm::setTolerance( RealOpenMM tolerance ){
326
327
328
329


   // ---------------------------------------------------------------------------------------

330
   // static const char* methodName = "\nReferenceCCMAAlgorithm::setTolerance";
331
332
333
334
335
336
337
338

   // ---------------------------------------------------------------------------------------

   _tolerance = tolerance;;
}

/**---------------------------------------------------------------------------------------

339
   Apply CCMA algorithm
340
341
342
343
344
345
346
347
348
349
350

   @param numberOfAtoms    number of atoms
   @param atomCoordinates  atom coordinates
   @param atomCoordinatesP atom coordinates prime
   @param inverseMasses    1/mass

   @return ReferenceDynamics::DefaultReturn if converge; else
    return ReferenceDynamics::ErrorReturn

   --------------------------------------------------------------------------------------- */

351
int ReferenceCCMAAlgorithm::apply( int numberOfAtoms, RealOpenMM** atomCoordinates,
352
353
354
355
356
                                         RealOpenMM** atomCoordinatesP,
                                         RealOpenMM* inverseMasses ){

   // ---------------------------------------------------------------------------------------

357
   static const char* methodName = "\nReferenceCCMAAlgorithm::apply";
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409

   static const RealOpenMM zero        =  0.0;
   static const RealOpenMM one         =  1.0;
   static const RealOpenMM two         =  2.0;
   static const RealOpenMM half        =  0.5;

   static const RealOpenMM epsilon6    = (RealOpenMM) 1.0e-06;

   // ---------------------------------------------------------------------------------------

   // temp arrays

   RealOpenMM** r_ij                   = _r_ij;
   RealOpenMM* d_ij2                   = _d_ij2;
   RealOpenMM* distanceTolerance       = _distanceTolerance;
   RealOpenMM* reducedMasses           = _reducedMasses;

   // calculate reduced masses on 1st pass

   if( !_hasInitializedMasses ){
      _hasInitializedMasses = true;
      for( int ii = 0; ii < _numberOfConstraints; ii++ ){
         int atomI          = _atomIndices[ii][0];
         int atomJ          = _atomIndices[ii][1];
         reducedMasses[ii]  = half/( inverseMasses[atomI] + inverseMasses[atomJ] );
      }
   }

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

   RealOpenMM tolerance     = getTolerance();
              tolerance    *= two;
   for( int ii = 0; ii < _numberOfConstraints; ii++ ){

      int atomI   = _atomIndices[ii][0];
      int atomJ   = _atomIndices[ii][1];
      for( int jj = 0; jj < 3; jj++ ){
         r_ij[ii][jj] = atomCoordinates[atomI][jj] - atomCoordinates[atomJ][jj];
      }
      d_ij2[ii]              = DOT3( r_ij[ii], r_ij[ii] );
      distanceTolerance[ii]  = d_ij2[ii]*tolerance;
      if( distanceTolerance[ii] > zero ){
         distanceTolerance[ii] = one/distanceTolerance[ii];
      }
   }
   RealOpenMM lowerTol = one-two*getTolerance()+getTolerance()*getTolerance();
   RealOpenMM upperTol = one+two*getTolerance()+getTolerance()*getTolerance();

   // main loop

   int iterations           = 0;
   int numberConverged      = 0;
410
411
   vector<RealOpenMM> constraintDelta(_numberOfConstraints);
   vector<RealOpenMM> tempDelta(_numberOfConstraints);
412
   while (iterations < getMaximumNumberOfIterations()) {
413
414
415
416
417
418
419
420
421
422
423
424
425
      numberConverged  = 0;
      for( int ii = 0; ii < _numberOfConstraints; ii++ ){

         int atomI   = _atomIndices[ii][0];
         int atomJ   = _atomIndices[ii][1];

         RealOpenMM rp_ij[3];
         for( int jj = 0; jj < 3; jj++ ){
            rp_ij[jj] = atomCoordinatesP[atomI][jj] - atomCoordinatesP[atomJ][jj];
         }
         RealOpenMM rp2  = DOT3( rp_ij, rp_ij );
         RealOpenMM dist2= _distance[ii]*_distance[ii];
         RealOpenMM diff = dist2 - rp2;
426
         constraintDelta[ii] = zero;
427
428
429
430
431
432
         RealOpenMM rrpr  = DOT3(  rp_ij, r_ij[ii] );
         if( rrpr <  d_ij2[ii]*epsilon6 ){
             std::stringstream message;
             message << iterations <<" "<<atomI<<" "<<atomJ<< " Error: sign of rrpr < 0?\n";
             SimTKOpenMMLog::printMessage( message );
         } else {
433
             constraintDelta[ii] = reducedMasses[ii]*diff/rrpr;
434
435
436
437
438
         }
         if (rp2 >= lowerTol*dist2 && rp2 <= upperTol*dist2) {
            numberConverged++;
         }
      }
439
440
441
      if( numberConverged == _numberOfConstraints )
         break;
      iterations++;
442

443
444
445
446
447
      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];
448
                  sum += element.second*constraintDelta[element.first];
449
              }
450
              tempDelta[i] = sum;
451
          }
452
          constraintDelta = tempDelta;
453
454
455
456
457
458
      }
      for( int ii = 0; ii < _numberOfConstraints; ii++ ){

         int atomI   = _atomIndices[ii][0];
         int atomJ   = _atomIndices[ii][1];
         for( int jj = 0; jj < 3; jj++ ){
459
            RealOpenMM dr                = constraintDelta[ii]*r_ij[ii][jj];
460
461
462
463
464
465
            atomCoordinatesP[atomI][jj] += inverseMasses[atomI]*dr;
            atomCoordinatesP[atomJ][jj] -= inverseMasses[atomJ]*dr;
         }
      }
   }

466
   return (numberConverged == _numberOfConstraints ? ReferenceDynamics::DefaultReturn : ReferenceDynamics::ErrorReturn);
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481

}

/**---------------------------------------------------------------------------------------

   Report any violated constriants

   @param numberOfAtoms    number of atoms
   @param atomCoordinates  atom coordinates
   @param message          report

   @return number of violated constraints

   --------------------------------------------------------------------------------------- */

482
int ReferenceCCMAAlgorithm::reportCCMA( int numberOfAtoms, RealOpenMM** atomCoordinates,
483
484
485
486
                                          std::stringstream& message ){

   // ---------------------------------------------------------------------------------------

487
   static const char* methodName = "\nReferenceCCMAAlgorithm::reportCCMA";
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529

   static const RealOpenMM zero        =  0.0;
   static const RealOpenMM one         =  1.0;
   static const RealOpenMM two         =  2.0;
   static const RealOpenMM three       =  3.0;
   static const RealOpenMM oneM        = -1.0;
   static const RealOpenMM half        =  0.5;

   // ---------------------------------------------------------------------------------------

   int numberOfConstraints = getNumberOfConstraints();

   // loop over constraints calculating distance and comparing to
   // expected distance -- report any contraints that are violated

   int numberConverged  = 0;
   RealOpenMM tolerance = getTolerance();
   for( int ii = 0; ii < _numberOfConstraints; ii++ ){

      int atomI   = _atomIndices[ii][0];
      int atomJ   = _atomIndices[ii][1];

      RealOpenMM rp2  = zero;
      for( int jj = 0; jj < 3; jj++ ){
         rp2 += (atomCoordinates[atomI][jj] - atomCoordinates[atomJ][jj])*(atomCoordinates[atomI][jj] - atomCoordinates[atomJ][jj]);
      }
      RealOpenMM diff = FABS( rp2 - (_distance[ii]*_distance[ii]) );
      if( diff > tolerance ){
         message << ii << " constraint violated: " << atomI << " " << atomJ << "] d=" << SQRT( rp2 ) << " " << rp2 << " d0=" << _distance[ii];
         message << " diff=" << diff;
         message << " [" << atomCoordinates[atomI][0] << " " << atomCoordinates[atomI][1] << " " << atomCoordinates[atomI][2] << "] ";
         message << " [" << atomCoordinates[atomJ][0] << " " << atomCoordinates[atomJ][1] << " " << atomCoordinates[atomJ][2] << "] ";
         message << "\n";
      } else {
         numberConverged++;
      }
   }

   return (numberOfConstraints-numberConverged);

}