OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
LangevinHullForceModifier.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 "integrators/LangevinHullForceModifier.hpp"
46
47#include <cmath>
48#include <fstream>
49#include <iostream>
50#include <map>
51#include <memory>
52#include <random>
53#include <string>
54
55#ifdef IS_MPI
56#include <mpi.h>
57#endif
58
59#include "brains/ForceModifier.hpp"
60#include "math/AlphaHull.hpp"
61#include "math/CholeskyDecomposition.hpp"
62#include "math/ConvexHull.hpp"
63#include "math/Triangle.hpp"
64#include "utils/Constants.hpp"
65
66using namespace std;
67
68namespace OpenMD {
69
70 LangevinHullForceModifier::LangevinHullForceModifier(SimInfo* info) :
71 ForceModifier {info} {
72 simParams_ = info->getSimParams();
73 veloMunge = std::make_unique<Velocitizer>(info);
74
75 // Create Hull, Convex Hull for now, other options later.
76
77 stringToEnumMap_["Convex"] = hullConvex;
78 stringToEnumMap_["AlphaShape"] = hullAlphaShape;
79 stringToEnumMap_["Unknown"] = hullUnknown;
80
81 const std::string ht = simParams_->getHULL_Method();
82
83 std::map<std::string, HullTypeEnum>::iterator iter;
84 iter = stringToEnumMap_.find(ht);
85 hullType_ = (iter == stringToEnumMap_.end()) ?
86 LangevinHullForceModifier::hullUnknown :
87 iter->second;
88
89 switch (hullType_) {
90 case hullConvex:
91 surfaceMesh_ = new ConvexHull();
92 break;
93 case hullAlphaShape:
94 surfaceMesh_ = new AlphaHull(simParams_->getAlpha());
95 break;
96 case hullUnknown:
97 default:
98 snprintf(
99 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
100 "LangevinHallForceManager: Unknown Hull_Method was requested!\n");
101 painCave.isFatal = 1;
102 simError();
103 break;
104 }
105
106 doThermalCoupling_ = true;
107 doPressureCoupling_ = true;
108
109 /* Check that the simulation has target pressure and target
110 temperature set */
111 if (!simParams_->haveTargetTemp()) {
112 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
113 "LangevinHullForceManager: no targetTemp (K) was set.\n"
114 "\tOpenMD is turning off the thermal coupling to the bath.\n");
115 painCave.isFatal = 0;
116 painCave.severity = OPENMD_INFO;
117 simError();
118 doThermalCoupling_ = false;
119 } else {
120 targetTemp_ = simParams_->getTargetTemp();
121
122 if (!simParams_->haveViscosity()) {
123 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
124 "LangevinHullForceManager: no viscosity was set.\n"
125 "\tOpenMD is turning off the thermal coupling to the bath.\n");
126 painCave.isFatal = 0;
127 painCave.severity = OPENMD_INFO;
128 simError();
129 doThermalCoupling_ = false;
130 } else {
131 viscosity_ = simParams_->getViscosity();
132 if (fabs(viscosity_) < 1e-6) {
133 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
134 "LangevinHullDynamics: The bath viscosity was set lower\n"
135 "\tthan 1e-6 poise. OpenMD is turning off the thermal\n"
136 "\tcoupling to the bath.\n");
137 painCave.isFatal = 0;
138 painCave.severity = OPENMD_INFO;
139 simError();
140 doThermalCoupling_ = false;
141 }
142 }
143 }
144 if (!simParams_->haveTargetPressure()) {
145 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
146 "LangevinHullForceManager: no targetPressure (atm) was set.\n"
147 "\tOpenMD is turning off the pressure coupling to the bath.\n");
148 painCave.isFatal = 0;
149 painCave.severity = OPENMD_INFO;
150 simError();
151 doPressureCoupling_ = false;
152 } else {
153 // Convert pressure from atm -> amu/(fs^2*Ang)
154 targetPressure_ =
155 simParams_->getTargetPressure() / Constants::pressureConvert;
156 }
157 if (simParams_->getUsePeriodicBoundaryConditions()) {
158 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
159 "LangevinHallForceManager: You can't use the Langevin Hull\n"
160 "\tintegrator with periodic boundary conditions turned on!\n");
161 painCave.isFatal = 1;
162 simError();
163 }
164
165 dt_ = simParams_->getDt();
166
167#ifdef IS_MPI
168 if (worldRank == 0) {
169#endif
170 if (doThermalCoupling_) {
171 randNumGen_ = info->getRandomNumberGenerator();
172
173 RealType stdDev = std::sqrt(2.0 * Constants::kb * targetTemp_ / dt_);
174
175 forceDistribution_ = std::normal_distribution<RealType>(0.0, stdDev);
176 }
177#ifdef IS_MPI
178 }
179#endif
180
181 // Build a vector of integrable objects to determine if the are
182 // surface atoms
183 Molecule* mol;
184 StuntDouble* sd;
185 SimInfo::MoleculeIterator i;
186 Molecule::IntegrableObjectIterator j;
187
188 for (mol = info_->beginMolecule(i); mol != NULL;
189 mol = info_->nextMolecule(i)) {
190 for (sd = mol->beginIntegrableObject(j); sd != NULL;
191 sd = mol->nextIntegrableObject(j)) {
192 localSites_.push_back(sd);
193 }
194 }
195
196 // We need to make an initial guess at the bounding box in order
197 // to compute long range forces in ForceMatrixDecomposition:
198
199 // Compute surface Mesh
200 surfaceMesh_->computeHull(localSites_);
201 }
202
203 LangevinHullForceModifier::~LangevinHullForceModifier() {
204 delete surfaceMesh_;
205 }
206
207 void LangevinHullForceModifier::modifyForces() {
208 int nTriangles, thisFacet;
209 RealType thisArea;
210 std::vector<Triangle> sMesh;
211 Triangle thisTriangle;
212 std::vector<Triangle>::iterator face;
213 std::vector<StuntDouble*> vertexSDs;
214 std::vector<StuntDouble*>::iterator vertex;
215
216 Vector3d unitNormal, centroid, facetVel;
217 Vector3d langevinForce, vertexForce;
218 Vector3d extPressure, randomForce, dragForce;
219
220 Mat3x3d hydroTensor, S;
221 std::vector<Vector3d> randNums;
222
223 // Compute surface Mesh
224 surfaceMesh_->computeHull(localSites_);
225 // Get number of surface stunt doubles
226 sMesh = surfaceMesh_->getMesh();
227 nTriangles = sMesh.size();
228
229 if (doThermalCoupling_) {
230 // Generate all of the necessary random forces
231 randNums = genTriangleForces(nTriangles);
232 }
233
234 // Loop over the mesh faces
235 thisFacet = 0;
236 for (face = sMesh.begin(); face != sMesh.end(); ++face) {
237 thisTriangle = *face;
238 vertexSDs = thisTriangle.getVertices();
239 thisArea = thisTriangle.getArea();
240 unitNormal = thisTriangle.getUnitNormal();
241 centroid = thisTriangle.getCentroid();
242 facetVel = thisTriangle.getFacetVelocity();
243
244 langevinForce = V3Zero;
245
246 if (doPressureCoupling_) {
247 extPressure = -unitNormal * (targetPressure_ * thisArea) /
248 Constants::energyConvert;
249 langevinForce += extPressure;
250 }
251
252 if (doThermalCoupling_) {
253 hydroTensor = thisTriangle.computeHydrodynamicTensor(viscosity_);
254 hydroTensor *= Constants::viscoConvert;
255 CholeskyDecomposition(hydroTensor, S);
256 randomForce = S * randNums[thisFacet++];
257 dragForce = -hydroTensor * facetVel;
258 langevinForce += randomForce + dragForce;
259 }
260
261 // Apply triangle force to stuntdouble vertices
262 for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex) {
263 if ((*vertex) != NULL) {
264 vertexForce = langevinForce / RealType(3.0);
265 (*vertex)->addFrc(vertexForce);
266 }
267 }
268 }
269
270 if (simParams_->getConserveLinearMomentum())
271 veloMunge->removeComDrift();
272 if (simParams_->getConserveAngularMomentum())
273 veloMunge->removeAngularDrift();
274
275 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
276 currSnapshot->setVolume(surfaceMesh_->getVolume());
277 currSnapshot->setHullVolume(surfaceMesh_->getVolume());
278 }
279
280 std::vector<Vector3d> LangevinHullForceModifier::genTriangleForces(
281 int nTriangles) {
282 // zero fill the random vector before starting:
283 std::vector<Vector3d> gaussRand(nTriangles);
284 std::fill(gaussRand.begin(), gaussRand.end(), V3Zero);
285
286#ifdef IS_MPI
287 if (worldRank == 0) {
288#endif
289 for (int i = 0; i < nTriangles; i++) {
290 gaussRand[i][0] = forceDistribution_(*randNumGen_);
291 gaussRand[i][1] = forceDistribution_(*randNumGen_);
292 gaussRand[i][2] = forceDistribution_(*randNumGen_);
293 }
294#ifdef IS_MPI
295 }
296 // push these out to the other processors
297 // Same command on all nodes:
298 MPI_Bcast(gaussRand[0].getArrayPointer(), nTriangles * 3, MPI_REALTYPE, 0,
299 MPI_COMM_WORLD);
300#endif
301
302 return gaussRand;
303 }
304} // namespace OpenMD
Compute alpha complex or alpha shape.
Definition AlphaHull.hpp:79
Abstract class for external ForceModifier classes.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
Definition SimInfo.hpp:248
The Snapshot class is a repository storing dynamic data during a Simulation.
Definition Snapshot.hpp:147
Snapshot * getCurrentSnapshot()
Returns the pointer of current snapshot.
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Triangle provides geometric data to OpenMD.
Definition Triangle.hpp:65
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.