NoseHooverIntegrator.cpp 18.8 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
31
/* -------------------------------------------------------------------------- *
 *                                   OpenMM                                   *
 * -------------------------------------------------------------------------- *
 * This is part of the OpenMM molecular simulation toolkit originating from   *
 * Simbios, the NIH National Center for Physics-Based Simulation of           *
 * Biological Structures at Stanford, funded under the NIH Roadmap for        *
 * Medical Research, grant U54 GM072970. See https://simtk.org.               *
 *                                                                            *
 * Portions copyright (c) 2019 Stanford University and the Authors.           *
 * Authors: Andreas Krämer and Andrew C. Simmonett                            *
 * Contributors:                                                              *
 *                                                                            *
 * 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.                                     *
 * -------------------------------------------------------------------------- */

32
#include "openmm/NoseHooverIntegrator.h"
33
34
35
#include "openmm/Context.h"
#include "openmm/Force.h"
#include "openmm/System.h"
36
#include "openmm/NoseHooverChain.h"
37
38
39
#include "openmm/OpenMMException.h"
#include "openmm/CMMotionRemover.h"
#include "openmm/internal/ContextImpl.h"
40
#include "openmm/internal/AssertionUtilities.h"
41
#include "openmm/kernels.h"
42
#include <iostream>
43
#include <string>
44
#include <algorithm>
45
46
47
48
49

using namespace OpenMM;
using std::string;
using std::vector;

50
NoseHooverIntegrator::NoseHooverIntegrator(double stepSize):
51
52
    forcesAreValid(false),
    hasSubsystemThermostats_(true)
53
{
54
55
56
    setStepSize(stepSize);
    setConstraintTolerance(1e-5);
}
57
NoseHooverIntegrator::NoseHooverIntegrator(double temperature, double collisionFrequency, double stepSize,
58
59
60
61
                                           int chainLength, int numMTS, int numYoshidaSuzuki) : 
    forcesAreValid(false),
    hasSubsystemThermostats_(false) {

62
63
    setStepSize(stepSize);
    setConstraintTolerance(1e-5);
64
    addThermostat(temperature, collisionFrequency, chainLength, numMTS, numYoshidaSuzuki);
65
}
66

67
NoseHooverIntegrator::~NoseHooverIntegrator() {}
68

69
std::pair<double, double> NoseHooverIntegrator::propagateChain(std::pair<double, double> kineticEnergy, int chainID) {
70
    return nhcKernel.getAs<NoseHooverChainKernel>().propagateChain(*context, noseHooverChains.at(chainID), kineticEnergy, getStepSize());
71
72
}

73
74

int NoseHooverIntegrator::addThermostat(double temperature, double collisionFrequency,
75
                                        int chainLength, int numMTS, int numYoshidaSuzuki) {
76
    hasSubsystemThermostats_ = false;
77
78
    return addSubsystemThermostat(std::vector<int>(), std::vector<std::pair<int, int>>(), temperature,
                                  collisionFrequency, temperature, collisionFrequency, chainLength, numMTS, numYoshidaSuzuki);
79
80
}

81
int NoseHooverIntegrator::addSubsystemThermostat(const std::vector<int>& thermostatedParticles,
82
                                                 const std::vector< std::pair< int, int> > &thermostatedPairs,
83
84
                                                 double temperature, double collisionFrequency,
                                                 double relativeTemperature, double relativeCollisionFrequency,
85
                                                 int chainLength, int numMTS, int numYoshidaSuzuki) {
86
    int chainID = noseHooverChains.size();
87
88
89
90
91
92
93
94
    // check if one thermostat already applies to all atoms or pairs
    if ( (chainID > 0) && (noseHooverChains[0].getThermostatedAtoms().size()*noseHooverChains[0].getThermostatedPairs().size() == 0) ) {
        throw OpenMMException(
            "Cannot add thermostat, since one of the thermostats already in the integrator applies to all particles. "
            "To manually add thermostats, use the constructor that takes only the "
            "step size as an input argument."
        );
    }
95
96
97
98
99
100
101
    int nDOF = 0; // is set in initializeThermostats()
    noseHooverChains.emplace_back(temperature, relativeTemperature,
                                 collisionFrequency, relativeCollisionFrequency,
                                 nDOF, chainLength, numMTS,
                                 numYoshidaSuzuki, chainID,
                                 thermostatedParticles, thermostatedPairs);
    return chainID;
102
103
}

104
const NoseHooverChain& NoseHooverIntegrator::getThermostat(int chainID) const {
105
106
107
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
    return noseHooverChains[chainID];
}
108

109
void NoseHooverIntegrator::initializeThermostats(const System &system) {
110

111
112
113
114
115
    allAtoms.clear();
    allPairs.clear();
    std::set<int> allAtomsSet;
    for (int atom = 0; atom < system.getNumParticles(); ++atom) allAtomsSet.insert(atom);

116
117
118
    for (auto &thermostat : noseHooverChains) {
        const auto &thermostatedParticles = thermostat.getThermostatedAtoms();
        const auto &thermostatedPairs = thermostat.getThermostatedPairs();
119

120
121
122
        for( auto& pair : thermostatedPairs ) {
            allAtomsSet.erase(pair.first);
            allAtomsSet.erase(pair.second);
123
            allPairs.emplace_back(pair.first, pair.second, thermostat.getRelativeTemperature());
124
125
        }

126
        // figure out the number of DOFs
127
        int nDOF = 3*(thermostatedParticles.size() + thermostatedPairs.size());
128
129
130
131
        for (int constraintNum = 0; constraintNum < system.getNumConstraints(); constraintNum++) {
            int particle1, particle2;
            double distance;
            system.getConstraintParameters(constraintNum, particle1, particle2, distance);
132
133
134
135
136
            bool particle1_in_thermostatedParticles = ((std::find(thermostatedParticles.begin(),
                                                                  thermostatedParticles.end(), particle1)
                                                                    != thermostatedParticles.end())) ||
                                                      (std::find_if(thermostatedPairs.begin(),
                                                                    thermostatedPairs.end(),
137
138
                                                                    [&particle1](const std::pair<int, int>& pair){
                                                                           return pair.first == particle1 || pair.second == particle1;})
139
140
141
142
143
144
                                                                      != thermostatedPairs.end());
            bool particle2_in_thermostatedParticles = ((std::find(thermostatedParticles.begin(),
                                                                  thermostatedParticles.end(), particle2)
                                                                    != thermostatedParticles.end())) ||
                                                      (std::find_if(thermostatedPairs.begin(),
                                                                    thermostatedPairs.end(),
145
146
                                                                    [&particle2](const std::pair<int, int>& pair){
                                                                           return pair.first == particle2 || pair.second == particle2;})
147
                                                                      != thermostatedPairs.end());
148
149
150
151
152
153
154
155
156
            if ((system.getParticleMass(particle1) > 0) && (system.getParticleMass(particle2) > 0)){
                if ((particle1_in_thermostatedParticles && !particle2_in_thermostatedParticles) ||
                     (!particle1_in_thermostatedParticles && particle2_in_thermostatedParticles)){
                    throw OpenMMException("Cannot add only one of particles " + std::to_string(particle1) + " and " + std::to_string(particle2)
                                            + " to NoseHooverChain, because they are connected by a constraint.");
                }
                if (particle1_in_thermostatedParticles && particle2_in_thermostatedParticles){
                    nDOF -= 1;
                }
157
            }
158
159
160
161
        }

        // remove 3 degrees of freedom from thermostats that act on absolute motions
        int numForces = system.getNumForces();
162
        if (thermostatedPairs.size() == 0){
163
164
            for (int forceNum = 0; forceNum < numForces; ++forceNum) {
                if (dynamic_cast<const CMMotionRemover*>(&system.getForce(forceNum))) nDOF -= 3;
165
166
            }
        }
167

168
        // set number of DoFs for chain 
169
        thermostat.setNumDegreesOfFreedom(nDOF);
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


    for (int chain1 = 0; chain1 < noseHooverChains.size(); ++chain1){
        const auto& nhc = noseHooverChains[chain1];

       // make sure that thermostats do not overlap
        for (int chain2 = 0; chain2 < chain1; ++chain2){
            const auto& other_nhc = noseHooverChains[chain2];
            for (auto &particle: nhc.getThermostatedAtoms()){
                bool isParticleInOtherChain = (std::find(other_nhc.getThermostatedAtoms().begin(),
                                                         other_nhc.getThermostatedAtoms().end(),
                                                         particle) != other_nhc.getThermostatedAtoms().end()) ||
                                              (std::find_if(other_nhc.getThermostatedPairs().begin(),
                                                            other_nhc.getThermostatedPairs().end(),
                                               [&particle](const std::pair<int, int>& pair){ return pair.first == particle || pair.second == particle;})
                                                  != other_nhc.getThermostatedPairs().end());
                if (isParticleInOtherChain){
                    throw OpenMMException("Found particle " + std::to_string(particle) + "in a different NoseHooverChain, "
                                          "but particles can only be thermostated by one thermostat.");
                }
            }
            for (auto &pair: nhc.getThermostatedPairs()){
                bool isParticleInOtherChain = (std::find(other_nhc.getThermostatedAtoms().begin(),
                                                         other_nhc.getThermostatedAtoms().end(),
                                                         pair.first) != other_nhc.getThermostatedAtoms().end()) ||
                                              (std::find(other_nhc.getThermostatedAtoms().begin(),
                                                         other_nhc.getThermostatedAtoms().end(),
                                                         pair.second) != other_nhc.getThermostatedAtoms().end()) ||
                                              (std::find_if(other_nhc.getThermostatedPairs().begin(),
                                                            other_nhc.getThermostatedPairs().end(),
                                               [&pair](const std::pair<int, int>& other_pair){
                                                        return pair.first == other_pair.first || pair.first == other_pair.second ||
                                                               pair.second == other_pair.first || pair.second == other_pair.second;})
                                                  != other_nhc.getThermostatedPairs().end());
                if (isParticleInOtherChain){
                    throw OpenMMException("Found pair " + std::to_string(pair.first) + "," +
                                          std::to_string(pair.second) + " in a different NoseHooverChain, "
                                          "but particles can only be thermostated by one thermostat.");
                }
            }
211
        }
212

213
214
215
        // make sure that massless particles are not thermostated
        for(auto particle: nhc.getThermostatedAtoms()){
            double mass = system.getParticleMass(particle);
216
            if (mass == 0.0) {
217
                throw OpenMMException("Found a particle with no mass (" + std::to_string(particle) + ") in a thermostat. Massless particles cannot be thermostated.");
218
219
220
            }
        }
    }
221
222

    for(const auto& atom : allAtomsSet) allAtoms.push_back(atom);
223
224
}

225
double NoseHooverIntegrator::getTemperature(int chainID) const {
226
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
227
    return noseHooverChains[chainID].getTemperature();
228
229
}

230
void NoseHooverIntegrator::setTemperature(double temperature, int chainID){
231
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
232
    noseHooverChains[chainID].setTemperature(temperature);
233
234
235

}

236
double NoseHooverIntegrator::getRelativeTemperature(int chainID) const {
237
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
238
    return noseHooverChains[chainID].getRelativeTemperature();
239
240
241
}

void NoseHooverIntegrator::setRelativeTemperature(double temperature, int chainID){
242
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
243
    noseHooverChains[chainID].setRelativeTemperature(temperature);
244
245
246
247

}

double NoseHooverIntegrator::getCollisionFrequency(int chainID) const {
248
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
249
    return noseHooverChains[chainID].getCollisionFrequency();
250
251
252
}

void NoseHooverIntegrator::setCollisionFrequency(double frequency, int chainID){
253
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
254
    noseHooverChains[chainID].setCollisionFrequency(frequency);
255
256
}

257
double NoseHooverIntegrator::getRelativeCollisionFrequency(int chainID) const {
258
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
259
    return noseHooverChains[chainID].getRelativeCollisionFrequency();
260
261
262
}

void NoseHooverIntegrator::setRelativeCollisionFrequency(double frequency, int chainID){
263
    ASSERT_VALID_INDEX(chainID, noseHooverChains);
264
    noseHooverChains[chainID].setRelativeCollisionFrequency(frequency);
265
266
267
}

double NoseHooverIntegrator::computeKineticEnergy() {
268
    double kE = 0.0;
269
    if(noseHooverChains.size() > 0) {
270
        for (const auto &nhc: noseHooverChains){
271
            kE += nhcKernel.getAs<NoseHooverChainKernel>().computeMaskedKineticEnergy(*context, nhc, true).first;
272
273
274
275
276
        }
    } else {
        kE = vvKernel.getAs<IntegrateVelocityVerletStepKernel>().computeKineticEnergy(*context, *this);
    }
    return kE;
277
278
}

279
double NoseHooverIntegrator::computeHeatBathEnergy() {
280
281
    double energy = 0;
    for(auto &nhc : noseHooverChains) {
282
        if (context && (nhc.getNumDegreesOfFreedom() > 0)) {
Andy Simmonett's avatar
Andy Simmonett committed
283
284
            energy += nhcKernel.getAs<NoseHooverChainKernel>().computeHeatBathEnergy(*context, nhc);
        }
285
286
287
288
    }
    return energy;
}

289
void NoseHooverIntegrator::initialize(ContextImpl& contextRef) {
290
291
    if (owner != NULL && &contextRef.getOwner() != owner)
        throw OpenMMException("This Integrator is already bound to a context");
292

293
    context = &contextRef;
294
    const System& system = context->getSystem();
295
296
297
    owner = &contextRef.getOwner();
    vvKernel = context->getPlatform().createKernel(IntegrateVelocityVerletStepKernel::Name(), contextRef);
    vvKernel.getAs<IntegrateVelocityVerletStepKernel>().initialize(contextRef.getSystem(), *this);
298
299
    nhcKernel = context->getPlatform().createKernel(NoseHooverChainKernel::Name(), contextRef);
    nhcKernel.getAs<NoseHooverChainKernel>().initialize();
300
    forcesAreValid = false;
301
302

    // check for drude particles and build the Nose-Hoover Chains
303
    for (auto& thermostat: noseHooverChains){
304
        // if there are no thermostated particles or pairs in the lists this is a regular thermostat for the whole (non-Drude) system
305
306
        if ( (thermostat.getThermostatedAtoms().size() == 0) && (thermostat.getThermostatedPairs().size() == 0) ){
            std::vector<int> thermostatedParticles;
307
308
309
            for(int particle = 0; particle < system.getNumParticles(); ++particle) {
                double mass = system.getParticleMass(particle);
                if ( (mass > 0) && (mass < 0.8) ){
310
                    std::cout << "Warning: Found particles with mass between 0.0 and 0.8 dalton. Did you mean to make a DrudeNoseHooverIntegrator instead? "
311
312
313
                                 "The thermostat you are about to use will not treat these particles as Drude particles!" << std::endl;
                }
                if(system.getParticleMass(particle) > 0) {
314
                   thermostatedParticles.push_back(particle);
315
316
                }
            }
317
            thermostat.setThermostatedAtoms(thermostatedParticles);
318
319
320
        }
    }

321
    initializeThermostats(system);
322
323
}

324
void NoseHooverIntegrator::cleanup() {
325
    vvKernel = Kernel();
326
    nhcKernel = Kernel();
327
328
}

329
vector<string> NoseHooverIntegrator::getKernelNames() {
330
    std::vector<std::string> names;
331
    names.push_back(NoseHooverChainKernel::Name());
332
333
334
335
    names.push_back(IntegrateVelocityVerletStepKernel::Name());
    return names;
}

336
void NoseHooverIntegrator::step(int steps) {
337
338
    if (context == NULL)
        throw OpenMMException("This Integrator is not bound to a context!");
339
    std::pair<double, double> scale, kineticEnergy;
340
341
    for (int i = 0; i < steps; ++i) {
        context->updateContextState();
342
        for(auto &nhc : noseHooverChains) {
343
            kineticEnergy = nhcKernel.getAs<NoseHooverChainKernel>().computeMaskedKineticEnergy(*context, nhc, false);
344
345
            scale = nhcKernel.getAs<NoseHooverChainKernel>().propagateChain(*context, nhc, kineticEnergy, getStepSize());
            nhcKernel.getAs<NoseHooverChainKernel>().scaleVelocities(*context, nhc, scale);
346
        }
347
        vvKernel.getAs<IntegrateVelocityVerletStepKernel>().execute(*context, *this, forcesAreValid);
348
        for(auto &nhc : noseHooverChains) {
349
            kineticEnergy = nhcKernel.getAs<NoseHooverChainKernel>().computeMaskedKineticEnergy(*context, nhc, false);
350
351
            scale = nhcKernel.getAs<NoseHooverChainKernel>().propagateChain(*context, nhc, kineticEnergy, getStepSize());
            nhcKernel.getAs<NoseHooverChainKernel>().scaleVelocities(*context, nhc, scale);
352
        }
353
354
    }
}