1 |
/* |
2 |
* Copyright (c) 2008 The University of Notre Dame. All Rights Reserved. |
3 |
* |
4 |
* The University of Notre Dame grants you ("Licensee") a |
5 |
* non-exclusive, royalty free, license to use, modify and |
6 |
* redistribute this software in source and binary code form, provided |
7 |
* that the following conditions are met: |
8 |
* |
9 |
* 1. Acknowledgement of the program authors must be made in any |
10 |
* publication of scientific results based in part on use of the |
11 |
* program. An acceptable form of acknowledgement is citation of |
12 |
* the article in which the program was described (Matthew |
13 |
* A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
* J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
* Parallel Simulation Engine for Molecular Dynamics," |
16 |
* J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
* |
18 |
* 2. Redistributions of source code must retain the above copyright |
19 |
* notice, this list of conditions and the following disclaimer. |
20 |
* |
21 |
* 3. Redistributions in binary form must reproduce the above copyright |
22 |
* notice, this list of conditions and the following disclaimer in the |
23 |
* documentation and/or other materials provided with the |
24 |
* distribution. |
25 |
* |
26 |
* This software is provided "AS IS," without a warranty of any |
27 |
* kind. All express or implied conditions, representations and |
28 |
* warranties, including any implied warranty of merchantability, |
29 |
* fitness for a particular purpose or non-infringement, are hereby |
30 |
* excluded. The University of Notre Dame and its licensors shall not |
31 |
* be liable for any damages suffered by licensee as a result of |
32 |
* using, modifying or distributing the software or its |
33 |
* derivatives. In no event will the University of Notre Dame or its |
34 |
* licensors be liable for any lost revenue, profit or data, or for |
35 |
* direct, indirect, special, consequential, incidental or punitive |
36 |
* damages, however caused and regardless of the theory of liability, |
37 |
* arising out of the use of or inability to use software, even if the |
38 |
* University of Notre Dame has been advised of the possibility of |
39 |
* such damages. |
40 |
*/ |
41 |
#include <fstream> |
42 |
#include <iostream> |
43 |
#include "integrators/SMIPDForceManager.hpp" |
44 |
#include "math/CholeskyDecomposition.hpp" |
45 |
#include "utils/OOPSEConstant.hpp" |
46 |
#include "math/ConvexHull.hpp" |
47 |
#include "math/Triangle.hpp" |
48 |
|
49 |
namespace oopse { |
50 |
|
51 |
SMIPDForceManager::SMIPDForceManager(SimInfo* info) : ForceManager(info) { |
52 |
|
53 |
simParams = info->getSimParams(); |
54 |
veloMunge = new Velocitizer(info); |
55 |
|
56 |
// Create Hull, Convex Hull for now, other options later. |
57 |
|
58 |
surfaceMesh_ = new ConvexHull(); |
59 |
|
60 |
/* Check that the simulation has target pressure and target |
61 |
temperature set */ |
62 |
|
63 |
if (!simParams->haveTargetTemp()) { |
64 |
sprintf(painCave.errMsg, |
65 |
"SMIPDynamics error: You can't use the SMIPD integrator\n" |
66 |
" without a targetTemp!\n"); |
67 |
painCave.isFatal = 1; |
68 |
painCave.severity = OOPSE_ERROR; |
69 |
simError(); |
70 |
} else { |
71 |
targetTemp_ = simParams->getTargetTemp(); |
72 |
} |
73 |
|
74 |
if (!simParams->haveTargetPressure()) { |
75 |
sprintf(painCave.errMsg, |
76 |
"SMIPDynamics error: You can't use the SMIPD integrator\n" |
77 |
" without a targetPressure!\n"); |
78 |
painCave.isFatal = 1; |
79 |
simError(); |
80 |
} else { |
81 |
// Convert pressure from atm -> amu/(fs^2*Ang) |
82 |
targetPressure_ = simParams->getTargetPressure() / |
83 |
OOPSEConstant::pressureConvert; |
84 |
} |
85 |
|
86 |
if (simParams->getUsePeriodicBoundaryConditions()) { |
87 |
sprintf(painCave.errMsg, |
88 |
"SMIPDynamics error: You can't use the SMIPD integrator\n" |
89 |
" with periodic boundary conditions!\n"); |
90 |
painCave.isFatal = 1; |
91 |
simError(); |
92 |
} |
93 |
|
94 |
if (!simParams->haveViscosity()) { |
95 |
sprintf(painCave.errMsg, |
96 |
"SMIPDynamics error: You can't use the SMIPD integrator\n" |
97 |
" without a viscosity!\n"); |
98 |
painCave.isFatal = 1; |
99 |
painCave.severity = OOPSE_ERROR; |
100 |
simError(); |
101 |
}else{ |
102 |
viscosity_ = simParams->getViscosity(); |
103 |
} |
104 |
|
105 |
dt_ = simParams->getDt(); |
106 |
|
107 |
variance_ = 2.0 * OOPSEConstant::kb * targetTemp_ / dt_; |
108 |
|
109 |
// Build a vector of integrable objects to determine if the are |
110 |
// surface atoms |
111 |
Molecule* mol; |
112 |
StuntDouble* integrableObject; |
113 |
SimInfo::MoleculeIterator i; |
114 |
Molecule::IntegrableObjectIterator j; |
115 |
|
116 |
for (mol = info_->beginMolecule(i); mol != NULL; |
117 |
mol = info_->nextMolecule(i)) { |
118 |
for (integrableObject = mol->beginIntegrableObject(j); |
119 |
integrableObject != NULL; |
120 |
integrableObject = mol->nextIntegrableObject(j)) { |
121 |
localSites_.push_back(integrableObject); |
122 |
} |
123 |
} |
124 |
} |
125 |
|
126 |
void SMIPDForceManager::postCalculation(bool needStress){ |
127 |
SimInfo::MoleculeIterator i; |
128 |
Molecule::IntegrableObjectIterator j; |
129 |
Molecule* mol; |
130 |
StuntDouble* integrableObject; |
131 |
|
132 |
// Compute surface Mesh |
133 |
surfaceMesh_->computeHull(localSites_); |
134 |
|
135 |
// Get total area and number of surface stunt doubles |
136 |
RealType area = surfaceMesh_->getArea(); |
137 |
int nSurfaceSDs = surfaceMesh_->getNs(); |
138 |
std::vector<Triangle> sMesh = surfaceMesh_->getMesh(); |
139 |
int nTriangles = sMesh.size(); |
140 |
|
141 |
// Generate all of the necessary random forces |
142 |
std::vector<Vector3d> randNums = genTriangleForces(nTriangles, variance_); |
143 |
|
144 |
// Loop over the mesh faces and apply external pressure to each |
145 |
// of the faces |
146 |
std::vector<Triangle>::iterator face; |
147 |
std::vector<StuntDouble*>::iterator vertex; |
148 |
int thisFacet = 0; |
149 |
for (face = sMesh.begin(); face != sMesh.end(); ++face){ |
150 |
|
151 |
Triangle thisTriangle = *face; |
152 |
std::vector<StuntDouble*> vertexSDs = thisTriangle.getVertices(); |
153 |
RealType thisArea = thisTriangle.getArea(); |
154 |
Vector3d unitNormal = thisTriangle.getNormal(); |
155 |
unitNormal.normalize(); |
156 |
Vector3d centroid = thisTriangle.getCentroid(); |
157 |
Vector3d facetVel = thisTriangle.getFacetVelocity(); |
158 |
|
159 |
Mat3x3d hydroTensor = thisTriangle.computeHydrodynamicTensor(viscosity_); |
160 |
|
161 |
hydroTensor *= OOPSEConstant::viscoConvert; |
162 |
Mat3x3d S; |
163 |
CholeskyDecomposition(hydroTensor, S); |
164 |
|
165 |
Vector3d extPressure = -unitNormal*(targetPressure_ * thisArea)/OOPSEConstant::energyConvert; |
166 |
Vector3d randomForce = S * randNums[thisFacet++]; |
167 |
Vector3d dragForce = -hydroTensor * facetVel; |
168 |
|
169 |
Vector3d langevinForce = (extPressure + randomForce + dragForce); |
170 |
|
171 |
// Apply triangle force to stuntdouble vertices |
172 |
for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex){ |
173 |
if ((*vertex) != NULL){ |
174 |
Vector3d vertexForce = langevinForce / 3.0; |
175 |
(*vertex)->addFrc(vertexForce); |
176 |
if ((*vertex)->isDirectional()){ |
177 |
Vector3d vertexPos = (*vertex)->getPos(); |
178 |
Vector3d vertexCentroidVector = vertexPos - centroid; |
179 |
(*vertex)->addTrq(cross(vertexCentroidVector,vertexForce)); |
180 |
} |
181 |
} |
182 |
} |
183 |
} |
184 |
|
185 |
veloMunge->removeComDrift(); |
186 |
veloMunge->removeAngularDrift(); |
187 |
|
188 |
Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
189 |
currSnapshot->setVolume(surfaceMesh_->getVolume()); |
190 |
ForceManager::postCalculation(needStress); |
191 |
} |
192 |
|
193 |
|
194 |
std::vector<Vector3d> SMIPDForceManager::genTriangleForces(int nTriangles, |
195 |
RealType variance) |
196 |
{ |
197 |
|
198 |
// zero fill the random vector before starting: |
199 |
std::vector<Vector3d> gaussRand; |
200 |
gaussRand.resize(nTriangles); |
201 |
std::fill(gaussRand.begin(), gaussRand.end(), V3Zero); |
202 |
|
203 |
#ifdef IS_MPI |
204 |
if (worldRank == 0) { |
205 |
#endif |
206 |
RealType rx, ry, rz; |
207 |
for (int i = 0; i < nTriangles; i++) { |
208 |
rx = randNumGen_.randNorm(0.0, variance); |
209 |
ry = randNumGen_.randNorm(0.0, variance); |
210 |
rz = randNumGen_.randNorm(0.0, variance); |
211 |
gaussRand[i][0] = rx; |
212 |
gaussRand[i][1] = ry; |
213 |
gaussRand[i][2] = rz; |
214 |
} |
215 |
#ifdef IS_MPI |
216 |
} |
217 |
#endif |
218 |
|
219 |
// push these out to the other processors |
220 |
|
221 |
#ifdef IS_MPI |
222 |
if (worldRank == 0) { |
223 |
MPI_Bcast(&gaussRand[0], nTriangles*3 , MPI_REALTYPE, 0, MPI_COMM_WORLD); |
224 |
} else { |
225 |
MPI_Bcast(&gaussRand[0], nTriangles*3, MPI_REALTYPE, 0, MPI_COMM_WORLD); |
226 |
} |
227 |
#endif |
228 |
|
229 |
return gaussRand; |
230 |
} |
231 |
} |