ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/integrators/SMIPDForceManager.cpp
Revision: 3481
Committed: Wed Nov 26 14:26:17 2008 UTC (15 years, 7 months ago) by gezelter
File size: 8012 byte(s)
Log Message:
Cleaning up some cruft.

File Contents

# Content
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 "utils/OOPSEConstant.hpp"
45 #include "math/ConvexHull.hpp"
46 #include "math/Triangle.hpp"
47
48 namespace oopse {
49
50 SMIPDForceManager::SMIPDForceManager(SimInfo* info) : ForceManager(info) {
51
52 simParams = info->getSimParams();
53
54 // Create Hull, Convex Hull for now, other options later.
55
56 surfaceMesh_ = new ConvexHull();
57
58 /* Check that the simulation has target pressure and target
59 temperature set*/
60
61 if (!simParams->haveTargetTemp()) {
62 sprintf(painCave.errMsg,
63 "SMIPDynamics error: You can't use the SMIPD integrator\n"
64 " without a targetTemp!\n");
65 painCave.isFatal = 1;
66 painCave.severity = OOPSE_ERROR;
67 simError();
68 } else {
69 targetTemp_ = simParams->getTargetTemp();
70 }
71
72 if (!simParams->haveTargetPressure()) {
73 sprintf(painCave.errMsg,
74 "SMIPDynamics error: You can't use the SMIPD integrator\n"
75 " without a targetPressure!\n");
76 painCave.isFatal = 1;
77 simError();
78 } else {
79 // Convert pressure from atm -> amu/(fs^2*Ang)
80 targetPressure_ = simParams->getTargetPressure() /
81 OOPSEConstant::pressureConvert;
82 }
83
84 if (simParams->getUsePeriodicBoundaryConditions()) {
85 sprintf(painCave.errMsg,
86 "SMIPDynamics error: You can't use the SMIPD integrator\n"
87 " with periodic boundary conditions!\n");
88 painCave.isFatal = 1;
89 simError();
90 }
91
92 if (!simParams->haveViscosity()) {
93 sprintf(painCave.errMsg,
94 "SMIPDynamics error: You can't use the SMIPD integrator\n"
95 " without a viscosity!\n");
96 painCave.isFatal = 1;
97 painCave.severity = OOPSE_ERROR;
98 simError();
99 }else{
100 viscosity_ = simParams->getViscosity();
101 }
102
103 dt_ = simParams->getDt();
104
105 variance_ = 2.0 * OOPSEConstant::kb * targetTemp_ / dt_;
106
107 Molecule* mol;
108 StuntDouble* integrableObject;
109 SimInfo::MoleculeIterator i;
110 Molecule::IntegrableObjectIterator j;
111
112 // Build a vector of integrable objects to determine if the are
113 // surface atoms
114
115 for (mol = info_->beginMolecule(i); mol != NULL;
116 mol = info_->nextMolecule(i)) {
117 for (integrableObject = mol->beginIntegrableObject(j);
118 integrableObject != NULL;
119 integrableObject = mol->nextIntegrableObject(j)) {
120 localSites_.push_back(integrableObject);
121 }
122 }
123 }
124
125 void SMIPDForceManager::postCalculation(bool needStress){
126 SimInfo::MoleculeIterator i;
127 Molecule::IntegrableObjectIterator j;
128 Molecule* mol;
129 StuntDouble* integrableObject;
130
131 // Compute surface Mesh
132 surfaceMesh_->computeHull(localSites_);
133
134 // Get total area and number of surface stunt doubles
135 RealType area = surfaceMesh_->getArea();
136 int nSurfaceSDs = surfaceMesh_->getNs();
137 std::vector<Triangle> sMesh = surfaceMesh_->getMesh();
138 int nTriangles = sMesh.size();
139
140 // Generate all of the necessary random forces
141 std::vector<RealType> randNums = genTriangleForces(nTriangles, variance_);
142
143 // Loop over the mesh faces and apply random force to each of the faces
144 std::vector<Triangle>::iterator face;
145 std::vector<StuntDouble*>::iterator vertex;
146 int thisNumber = 0;
147 for (face = sMesh.begin(); face != sMesh.end(); ++face){
148
149 Triangle thisTriangle = *face;
150 std::vector<StuntDouble*> vertexSDs = thisTriangle.getVertices();
151 RealType thisArea = thisTriangle.getArea();
152 Vector3d unitNormal = thisTriangle.getNormal();
153 unitNormal.normalize();
154
155 Vector3d centroid = thisTriangle.getCentroid();
156 Vector3d facetVel = thisTriangle.getFacetVelocity();
157 RealType hydroLength = thisTriangle.getIncircleRadius() * 2.0 /
158 NumericConstant::PI;
159
160 // gamma is the drag coefficient normal to the face of the triangle
161 RealType gamma = viscosity_ * hydroLength *
162 OOPSEConstant::viscoConvert;
163
164 RealType extPressure = - (targetPressure_ * thisArea) /
165 OOPSEConstant::energyConvert;
166
167 RealType randomForce = randNums[thisNumber++] * sqrt(gamma);
168 RealType dragForce = -gamma * dot(facetVel, unitNormal);
169
170 Vector3d langevinForce = (extPressure + randomForce + dragForce) * unitNormal;
171
172 // Apply triangle force to stuntdouble vertices
173 for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex) {
174 if ((*vertex) != NULL) {
175 Vector3d vertexForce = langevinForce / 3.0;
176 (*vertex)->addFrc(vertexForce);
177 if ((*vertex)->isDirectional()) {
178 Vector3d vertexPos = (*vertex)->getPos();
179 Vector3d vertexCentroidVector = vertexPos - centroid;
180 (*vertex)->addTrq(cross(vertexCentroidVector,vertexForce));
181 }
182 }
183 }
184 }
185
186 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
187 currSnapshot->setVolume(surfaceMesh_->getVolume());
188 ForceManager::postCalculation(needStress);
189 }
190
191 std::vector<RealType> SMIPDForceManager::genTriangleForces(int nTriangles,
192 RealType variance) {
193
194 // zero fill the random vector before starting:
195 std::vector<RealType> gaussRand;
196 gaussRand.resize(nTriangles);
197 std::fill(gaussRand.begin(), gaussRand.end(), 0.0);
198
199 #ifdef IS_MPI
200 if (worldRank == 0) {
201 #endif
202 for (int i = 0; i < nTriangles; i++) {
203 gaussRand[i] = randNumGen_.randNorm(0.0, variance);
204 }
205 #ifdef IS_MPI
206 }
207 #endif
208
209 // push these out to the other processors
210
211 #ifdef IS_MPI
212 if (worldRank == 0) {
213 MPI_Bcast(&gaussRand[0], nTriangles, MPI_REALTYPE, 0, MPI_COMM_WORLD);
214 } else {
215 MPI_Bcast(&gaussRand[0], nTriangles, MPI_REALTYPE, 0, MPI_COMM_WORLD);
216 }
217 #endif
218
219 return gaussRand;
220 }
221 }

Properties

Name Value
svn:executable *