OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
RestraintForceModifier.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include "restraints/RestraintForceModifier.hpp"
46
47#include <config.h>
48
49#include <cmath>
50#include <memory>
51
52#ifdef IS_MPI
53#include <mpi.h>
54#endif
55
56#include "brains/ForceModifier.hpp"
57#include "io/RestReader.hpp"
58#include "restraints/MolecularRestraint.hpp"
59#include "restraints/ObjectRestraint.hpp"
60#include "selection/SelectionEvaluator.hpp"
61#include "selection/SelectionManager.hpp"
62#include "utils/Constants.hpp"
63#include "utils/MemoryUtils.hpp"
64#include "utils/StringUtils.hpp"
65#include "utils/simError.h"
66
67namespace OpenMD {
68
69 RestraintForceModifier::RestraintForceModifier(SimInfo* info) :
70 ForceModifier {info} {
71 // order of affairs:
72 //
73 // 1) create restraints from the restraintStamps found in the MD
74 // file.
75 //
76 // 2) Create RestraintReader to parse the input files for the ideal
77 // structures. This reader will set reference structures, and will
78 // calculate molecular centers of mass, etc.
79 //
80 // 3) sit around and wait for calcForces to be called. When it comes,
81 // call the normal force manager calcForces, then loop through the
82 // restrained objects and do their restraint forces.
83
84 Globals* simParam = info_->getSimParams();
85 currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
86
87 currRestTime_ = currSnapshot_->getTime();
88
89 if (simParam->haveStatusTime()) {
90 restTime_ = simParam->getStatusTime();
91 } else {
92 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
93 "Restraint warning: If you use restraints without setting\n"
94 "\tstatusTime, no restraint data will be written to the rest\n"
95 "\tfile.\n");
96 painCave.isFatal = 0;
97 simError();
98 restTime_ = simParam->getRunTime();
99 }
100
101 int nRestraintStamps = simParam->getNRestraintStamps();
102 std::vector<RestraintStamp*> stamp = simParam->getRestraintStamps();
103
104 std::vector<int> stuntDoubleIndex;
105
106 for (int i = 0; i < nRestraintStamps; i++) {
107 std::string myType = toUpperCopy(stamp[i]->getType());
108
109 if (myType.compare("MOLECULAR") == 0) {
110 int molIndex(-1);
111 Vector3d refCom;
112
113 if (!stamp[i]->haveMolIndex()) {
114 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
115 "Restraint Error: A molecular restraint was specified\n"
116 "\twithout providing a value for molIndex.\n");
117 painCave.isFatal = 1;
118 simError();
119 } else {
120 molIndex = stamp[i]->getMolIndex();
121 }
122
123 if (molIndex < 0) {
124 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
125 "Restraint Error: A molecular restraint was specified\n"
126 "\twith a molIndex that was less than 0\n");
127 painCave.isFatal = 1;
128 simError();
129 }
130 if (molIndex >= info_->getNGlobalMolecules()) {
131 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
132 "Restraint Error: A molecular restraint was specified with\n"
133 "\ta molIndex that was greater than the total number of "
134 "molecules\n");
135 painCave.isFatal = 1;
136 simError();
137 }
138
139 Molecule* mol = info_->getMoleculeByGlobalIndex(molIndex);
140
141 if (mol == NULL) {
142#ifdef IS_MPI
143 // getMoleculeByGlobalIndex returns a NULL in parallel if
144 // this proc doesn't have the molecule. Do a quick check to
145 // make sure another processor is supposed to have it.
146
147 int myrank;
148 MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
149
150 if (info_->getMolToProc(molIndex) == myrank) {
151 // If we were supposed to have it but got a null, then freak out.
152#endif
153
154 snprintf(
155 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
156 "Restraint Error: A molecular restraint was specified, but\n"
157 "\tno molecule was found with global index %d.\n",
158 molIndex);
159 painCave.isFatal = 1;
160 simError();
161
162#ifdef IS_MPI
163 }
164#endif
165 }
166
167#ifdef IS_MPI
168 // only handle this molecular restraint if this processor owns the
169 // molecule
170 int myrank;
171 MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
172 if (info_->getMolToProc(molIndex) == myrank) {
173#endif
174
175 MolecularRestraint* rest = new MolecularRestraint();
176
177 std::string restPre("mol_");
178 std::stringstream restName;
179 restName << restPre << molIndex;
180 rest->setRestraintName(restName.str());
181
182 if (stamp[i]->haveDisplacementSpringConstant()) {
183 rest->setDisplacementForceConstant(
184 stamp[i]->getDisplacementSpringConstant());
185 }
186 if (stamp[i]->haveAbsoluteSpringConstant()) {
187 rest->setAbsoluteForceConstant(
188 stamp[i]->getAbsoluteSpringConstant());
189 }
190 if (stamp[i]->haveTwistSpringConstant()) {
191 rest->setTwistForceConstant(stamp[i]->getTwistSpringConstant());
192 }
193 if (stamp[i]->haveSwingXSpringConstant()) {
194 rest->setSwingXForceConstant(stamp[i]->getSwingXSpringConstant());
195 }
196 if (stamp[i]->haveSwingYSpringConstant()) {
197 rest->setSwingYForceConstant(stamp[i]->getSwingYSpringConstant());
198 }
199 if (stamp[i]->haveAbsolutePositionZ()) {
200 rest->setAbsolutePositionZ(stamp[i]->getAbsolutePositionZ());
201 }
202 if (stamp[i]->haveRestrainedTwistAngle()) {
203 rest->setRestrainedTwistAngle(stamp[i]->getRestrainedTwistAngle() *
204 Constants::PI / 180.0);
205 }
206 if (stamp[i]->haveRestrainedSwingYAngle()) {
207 rest->setRestrainedSwingYAngle(
208 stamp[i]->getRestrainedSwingYAngle() * Constants::PI / 180.0);
209 }
210 if (stamp[i]->haveRestrainedSwingXAngle()) {
211 rest->setRestrainedSwingXAngle(
212 stamp[i]->getRestrainedSwingXAngle() * Constants::PI / 180.0);
213 }
214 if (stamp[i]->havePrint()) {
215 rest->setPrintRestraint(stamp[i]->getPrint());
216 }
217
218 restraints_.push_back(rest);
219 mol->addProperty(std::make_shared<RestraintData>("Restraint", rest));
220 restrainedMols_.push_back(mol);
221#ifdef IS_MPI
222 }
223#endif
224 } else if (myType.compare("OBJECT") == 0) {
225 std::string objectSelection;
226
227 if (!stamp[i]->haveObjectSelection()) {
228 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
229 "Restraint Error: An object restraint was specified\n"
230 "\twithout providing a selection script in the\n"
231 "\tobjectSelection variable.\n");
232 painCave.isFatal = 1;
233 simError();
234 } else {
235 objectSelection = stamp[i]->getObjectSelection();
236 }
237
238 SelectionEvaluator evaluator(info);
239 SelectionManager seleMan(info);
240
241 evaluator.loadScriptString(objectSelection);
242 seleMan.setSelectionSet(evaluator.evaluate());
243 int selectionCount = seleMan.getSelectionCount();
244
245#ifdef IS_MPI
246 MPI_Allreduce(MPI_IN_PLACE, &selectionCount, 1, MPI_INT, MPI_SUM,
247 MPI_COMM_WORLD);
248#endif
249
250 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
251 "Restraint Info: The specified restraint objectSelection,\n"
252 "\t\t%s\n"
253 "\twill result in %d integrable objects being\n"
254 "\trestrained.\n",
255 objectSelection.c_str(), selectionCount);
256 painCave.severity = OPENMD_INFO;
257 painCave.isFatal = 0;
258 simError();
259
260 int selei;
261 StuntDouble* sd;
262
263 for (sd = seleMan.beginSelected(selei); sd != NULL;
264 sd = seleMan.nextSelected(selei)) {
265 stuntDoubleIndex.push_back(sd->getGlobalIntegrableObjectIndex());
266
267 ObjectRestraint* rest = new ObjectRestraint();
268
269 if (stamp[i]->haveDisplacementSpringConstant()) {
270 rest->setDisplacementForceConstant(
271 stamp[i]->getDisplacementSpringConstant());
272 }
273 if (stamp[i]->haveAbsoluteSpringConstant()) {
274 rest->setAbsoluteForceConstant(
275 stamp[i]->getAbsoluteSpringConstant());
276 }
277 if (stamp[i]->haveTwistSpringConstant()) {
278 rest->setTwistForceConstant(stamp[i]->getTwistSpringConstant());
279 }
280 if (stamp[i]->haveSwingXSpringConstant()) {
281 rest->setSwingXForceConstant(stamp[i]->getSwingXSpringConstant());
282 }
283 if (stamp[i]->haveSwingYSpringConstant()) {
284 rest->setSwingYForceConstant(stamp[i]->getSwingYSpringConstant());
285 }
286 if (stamp[i]->haveAbsolutePositionZ()) {
287 rest->setAbsolutePositionZ(stamp[i]->getAbsolutePositionZ());
288 }
289 if (stamp[i]->haveRestrainedTwistAngle()) {
290 rest->setRestrainedTwistAngle(stamp[i]->getRestrainedTwistAngle());
291 }
292 if (stamp[i]->haveRestrainedSwingXAngle()) {
293 rest->setRestrainedSwingXAngle(
294 stamp[i]->getRestrainedSwingXAngle());
295 }
296 if (stamp[i]->haveRestrainedSwingYAngle()) {
297 rest->setRestrainedSwingYAngle(
298 stamp[i]->getRestrainedSwingYAngle());
299 }
300 if (stamp[i]->havePrint()) {
301 rest->setPrintRestraint(stamp[i]->getPrint());
302 }
303
304 restraints_.push_back(rest);
305 sd->addProperty(std::make_shared<RestraintData>("Restraint", rest));
306 restrainedObjs_.push_back(sd);
307 }
308 }
309 }
310
311 // ThermodynamicIntegration subclasses RestraintForceManager, and there
312 // are times when it won't use restraints at all, so only open the
313 // restraint file if we are actually using restraints:
314 if (simParam->getUseRestraints()) {
315 std::string refFile = simParam->getRestraint_file();
316 RestReader* rr = new RestReader(info, refFile, stuntDoubleIndex);
317 rr->readReferenceStructure();
318 delete rr;
319 }
320
321 restOutput_ = getPrefix(info_->getFinalConfigFileName()) + ".rest";
322 restOut = new RestWriter(info_, restOutput_.c_str(), restraints_);
323 if (!restOut) {
324 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
325 "Restraint error: Failed to create RestWriter\n");
326 painCave.isFatal = 1;
327 simError();
328 }
329
330 // todo: figure out the scale factor. Right now, just scale them all to 1
331 std::vector<Restraint*>::const_iterator resti;
332 for (resti = restraints_.begin(); resti != restraints_.end(); ++resti) {
333 (*resti)->setScaleFactor(1.0);
334 }
335 }
336
337 RestraintForceModifier::~RestraintForceModifier() {
338 Utils::deletePointers(restraints_);
339
340 delete restOut;
341 }
342
343 void RestraintForceModifier::modifyForces() {
344 RealType restPot(0.0);
345
346 restPot = doRestraints(1.0);
347
348#ifdef IS_MPI
349 MPI_Allreduce(MPI_IN_PLACE, &restPot, 1, MPI_REALTYPE, MPI_SUM,
350 MPI_COMM_WORLD);
351#endif
352
353 currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
354 RealType rp = currSnapshot_->getRestraintPotential();
355 currSnapshot_->setRestraintPotential(rp + restPot);
356
357 RealType pe = currSnapshot_->getPotentialEnergy();
358 currSnapshot_->setRawPotential(pe);
359 currSnapshot_->setPotentialEnergy(pe + restPot);
360
361 // write out forces and current positions of restrained molecules
362 if (currSnapshot_->getTime() >= currRestTime_) {
363 restOut->writeRest(restInfo_);
364 currRestTime_ += restTime_;
365 }
366 }
367
368 RealType RestraintForceModifier::doRestraints(RealType scalingFactor) {
369 std::vector<Molecule*>::const_iterator rm;
370 std::shared_ptr<GenericData> data;
371 Molecule::IntegrableObjectIterator ioi;
372 MolecularRestraint* mRest = NULL;
373 ObjectRestraint* oRest = NULL;
374 StuntDouble* sd;
375
376 std::vector<StuntDouble*>::const_iterator ro;
377
378 std::map<int, Restraint::RealPair> restInfo;
379
380 unscaledPotential_ = 0.0;
381
382 restInfo_.clear();
383
384 for (rm = restrainedMols_.begin(); rm != restrainedMols_.end(); ++rm) {
385 // make sure this molecule (*rm) has a generic data for restraints:
386 data = (*rm)->getPropertyByName("Restraint");
387 if (data != nullptr) {
388 // make sure we can reinterpret the generic data as restraint data:
389 std::shared_ptr<RestraintData> restData =
390 std::dynamic_pointer_cast<RestraintData>(data);
391 if (restData != nullptr) {
392 // make sure we can reinterpet the restraint data as a
393 // MolecularRestraint
394 mRest = dynamic_cast<MolecularRestraint*>(restData->getData());
395 if (mRest == NULL) {
396 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
397 "Can not cast RestraintData to MolecularRestraint\n");
398 painCave.severity = OPENMD_ERROR;
399 painCave.isFatal = 1;
400 simError();
401 }
402 } else {
403 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
404 "Can not cast GenericData to RestraintData\n");
405 painCave.severity = OPENMD_ERROR;
406 painCave.isFatal = 1;
407 simError();
408 }
409 } else {
410 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
411 "Can not find Restraint for RestrainedObject\n");
412 painCave.severity = OPENMD_ERROR;
413 painCave.isFatal = 1;
414 simError();
415 }
416
417 // phew. At this point, we should have the pointer to the
418 // correct MolecularRestraint in the variable mRest.
419
420 Vector3d molCom = (*rm)->getCom();
421
422 std::vector<Vector3d> struc;
423 std::vector<Vector3d> forces;
424
425 for (sd = (*rm)->beginIntegrableObject(ioi); sd != NULL;
426 sd = (*rm)->nextIntegrableObject(ioi)) {
427 struc.push_back(sd->getPos());
428 }
429
430 mRest->setScaleFactor(scalingFactor);
431 mRest->calcForce(struc, molCom);
432 forces = mRest->getRestraintForces();
433 int index = 0;
434
435 for (sd = (*rm)->beginIntegrableObject(ioi); sd != NULL;
436 sd = (*rm)->nextIntegrableObject(ioi)) {
437 sd->addFrc(forces[index]);
438 struc.push_back(sd->getPos());
439 index++;
440 }
441
442 unscaledPotential_ += mRest->getUnscaledPotential();
443
444 // only collect data on restraints that we're going to print:
445 if (mRest->getPrintRestraint()) {
446 restInfo = mRest->getRestraintInfo();
447 restInfo_.push_back(restInfo);
448 }
449 }
450
451 for (ro = restrainedObjs_.begin(); ro != restrainedObjs_.end(); ++ro) {
452 // make sure this object (*ro) has a generic data for restraints:
453 data = (*ro)->getPropertyByName("Restraint");
454 if (data != NULL) {
455 // make sure we can reinterpret the generic data as restraint data:
456 std::shared_ptr<RestraintData> restData =
457 std::dynamic_pointer_cast<RestraintData>(data);
458 if (restData != nullptr) {
459 // make sure we can reinterpet the restraint data as a pointer to
460 // an ObjectRestraint:
461 oRest = dynamic_cast<ObjectRestraint*>(restData->getData());
462 if (oRest == NULL) {
463 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
464 "Can not cast RestraintData to ObjectRestraint\n");
465 painCave.severity = OPENMD_ERROR;
466 painCave.isFatal = 1;
467 simError();
468 }
469 } else {
470 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
471 "Can not cast GenericData to RestraintData\n");
472 painCave.severity = OPENMD_ERROR;
473 painCave.isFatal = 1;
474 simError();
475 }
476 } else {
477 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
478 "Can not find Restraint for RestrainedObject\n");
479 painCave.severity = OPENMD_ERROR;
480 painCave.isFatal = 1;
481 simError();
482 }
483
484 // phew. At this point, we should have the pointer to the
485 // correct Object restraint in the variable oRest.
486 oRest->setScaleFactor(scalingFactor);
487
488 Vector3d pos = (*ro)->getPos();
489
490 if ((*ro)->isDirectional()) {
491 // directional objects may have orientational restraints as well
492 // as positional, so get the rotation matrix first:
493
494 RotMat3x3d A = (*ro)->getA();
495 oRest->calcForce(pos, A);
496 (*ro)->addFrc(oRest->getRestraintForce());
497 (*ro)->addTrq(oRest->getRestraintTorque());
498
499 } else {
500 // plain vanilla positional restraints:
501
502 oRest->calcForce(pos);
503 (*ro)->addFrc(oRest->getRestraintForce());
504 }
505
506 unscaledPotential_ += oRest->getUnscaledPotential();
507
508 // only collect data on restraints that we're going to print:
509 if (oRest->getPrintRestraint()) {
510 restInfo = oRest->getRestraintInfo();
511 restInfo_.push_back(restInfo);
512 }
513 }
514
515 return unscaledPotential_ * scalingFactor;
516 }
517} // namespace OpenMD
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
std::string getPrefix(const std::string &str)