45#include "integrators/LangevinHullForceModifier.hpp"
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"
70 LangevinHullForceModifier::LangevinHullForceModifier(
SimInfo* info) :
72 simParams_ = info->getSimParams();
73 veloMunge = std::make_unique<Velocitizer>(info);
77 stringToEnumMap_[
"Convex"] = hullConvex;
78 stringToEnumMap_[
"AlphaShape"] = hullAlphaShape;
79 stringToEnumMap_[
"Unknown"] = hullUnknown;
81 const std::string ht = simParams_->getHULL_Method();
83 std::map<std::string, HullTypeEnum>::iterator iter;
84 iter = stringToEnumMap_.find(ht);
85 hullType_ = (iter == stringToEnumMap_.end()) ?
86 LangevinHullForceModifier::hullUnknown :
94 surfaceMesh_ =
new AlphaHull(simParams_->getAlpha());
99 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
100 "LangevinHallForceManager: Unknown Hull_Method was requested!\n");
101 painCave.isFatal = 1;
106 doThermalCoupling_ =
true;
107 doPressureCoupling_ =
true;
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;
118 doThermalCoupling_ =
false;
120 targetTemp_ = simParams_->getTargetTemp();
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;
129 doThermalCoupling_ =
false;
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;
140 doThermalCoupling_ =
false;
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;
151 doPressureCoupling_ =
false;
155 simParams_->getTargetPressure() / Constants::pressureConvert;
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;
165 dt_ = simParams_->getDt();
168 if (worldRank == 0) {
170 if (doThermalCoupling_) {
171 randNumGen_ = info->getRandomNumberGenerator();
173 RealType stdDev = std::sqrt(2.0 * Constants::kb * targetTemp_ / dt_);
175 forceDistribution_ = std::normal_distribution<RealType>(0.0, stdDev);
185 SimInfo::MoleculeIterator i;
186 Molecule::IntegrableObjectIterator j;
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);
200 surfaceMesh_->computeHull(localSites_);
203 LangevinHullForceModifier::~LangevinHullForceModifier() {
207 void LangevinHullForceModifier::modifyForces() {
208 int nTriangles, thisFacet;
210 std::vector<Triangle> sMesh;
212 std::vector<Triangle>::iterator face;
213 std::vector<StuntDouble*> vertexSDs;
214 std::vector<StuntDouble*>::iterator vertex;
216 Vector3d unitNormal, centroid, facetVel;
217 Vector3d langevinForce, vertexForce;
218 Vector3d extPressure, randomForce, dragForce;
221 std::vector<Vector3d> randNums;
224 surfaceMesh_->computeHull(localSites_);
226 sMesh = surfaceMesh_->getMesh();
227 nTriangles = sMesh.size();
229 if (doThermalCoupling_) {
231 randNums = genTriangleForces(nTriangles);
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();
244 langevinForce = V3Zero;
246 if (doPressureCoupling_) {
247 extPressure = -unitNormal * (targetPressure_ * thisArea) /
248 Constants::energyConvert;
249 langevinForce += extPressure;
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;
262 for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex) {
263 if ((*vertex) != NULL) {
264 vertexForce = langevinForce / RealType(3.0);
265 (*vertex)->addFrc(vertexForce);
270 if (simParams_->getConserveLinearMomentum())
271 veloMunge->removeComDrift();
272 if (simParams_->getConserveAngularMomentum())
273 veloMunge->removeAngularDrift();
276 currSnapshot->setVolume(surfaceMesh_->getVolume());
277 currSnapshot->setHullVolume(surfaceMesh_->getVolume());
280 std::vector<Vector3d> LangevinHullForceModifier::genTriangleForces(
283 std::vector<Vector3d> gaussRand(nTriangles);
284 std::fill(gaussRand.begin(), gaussRand.end(), V3Zero);
287 if (worldRank == 0) {
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_);
298 MPI_Bcast(gaussRand[0].getArrayPointer(), nTriangles * 3, MPI_REALTYPE, 0,
Compute alpha complex or alpha shape.
Abstract class for external ForceModifier classes.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
The Snapshot class is a repository storing dynamic data during a Simulation.
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.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.