ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/integrators/SMIPDForceManager.cpp
Revision: 3516
Committed: Wed Jul 22 15:00:21 2009 UTC (15 years, 9 months ago) by gezelter
File size: 8756 byte(s)
Log Message:
units for thermalConductivity, spelling error fixed

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 thermo = new Thermo(info);
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 (K)!\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 (atm)!\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->haveThermalConductivity()) {
95 sprintf(painCave.errMsg,
96 "SMIPDynamics error: You can't use the SMIPD integrator\n"
97 " without a thermalConductivity (W m^-1 K^-1)!\n");
98 painCave.isFatal = 1;
99 painCave.severity = OOPSE_ERROR;
100 simError();
101 }else{
102 thermalConductivity_ = simParams->getThermalConductivity() *
103 OOPSEConstant::thermalConductivityConvert;
104 }
105
106 if (!simParams->haveThermalLength()) {
107 sprintf(painCave.errMsg,
108 "SMIPDynamics error: You can't use the SMIPD integrator\n"
109 " without a thermalLength (Angstroms)!\n");
110 painCave.isFatal = 1;
111 painCave.severity = OOPSE_ERROR;
112 simError();
113 }else{
114 thermalLength_ = simParams->getThermalLength();
115 }
116
117 dt_ = simParams->getDt();
118
119 variance_ = 2.0 * OOPSEConstant::kb * targetTemp_ / dt_;
120
121 // Build a vector of integrable objects to determine if the are
122 // surface atoms
123 Molecule* mol;
124 StuntDouble* integrableObject;
125 SimInfo::MoleculeIterator i;
126 Molecule::IntegrableObjectIterator j;
127
128 for (mol = info_->beginMolecule(i); mol != NULL;
129 mol = info_->nextMolecule(i)) {
130 for (integrableObject = mol->beginIntegrableObject(j);
131 integrableObject != NULL;
132 integrableObject = mol->nextIntegrableObject(j)) {
133 localSites_.push_back(integrableObject);
134 }
135 }
136 }
137
138 void SMIPDForceManager::postCalculation(bool needStress){
139 SimInfo::MoleculeIterator i;
140 Molecule::IntegrableObjectIterator j;
141 Molecule* mol;
142 StuntDouble* integrableObject;
143
144 // Compute surface Mesh
145 surfaceMesh_->computeHull(localSites_);
146
147 // Get total area and number of surface stunt doubles
148 RealType area = surfaceMesh_->getArea();
149 int nSurfaceSDs = surfaceMesh_->getNs();
150 std::vector<Triangle> sMesh = surfaceMesh_->getMesh();
151 int nTriangles = sMesh.size();
152
153 // Generate all of the necessary random forces
154 std::vector<RealType> randNums = genTriangleForces(nTriangles, variance_);
155
156 RealType instaTemp = thermo->getTemperature();
157
158 // Loop over the mesh faces and apply external pressure to each
159 // of the faces
160 std::vector<Triangle>::iterator face;
161 std::vector<StuntDouble*>::iterator vertex;
162 int thisFacet = 0;
163 for (face = sMesh.begin(); face != sMesh.end(); ++face){
164
165 Triangle thisTriangle = *face;
166 std::vector<StuntDouble*> vertexSDs = thisTriangle.getVertices();
167 RealType thisArea = thisTriangle.getArea();
168 Vector3d unitNormal = thisTriangle.getNormal();
169 unitNormal.normalize();
170 Vector3d centroid = thisTriangle.getCentroid();
171 Vector3d facetVel = thisTriangle.getFacetVelocity();
172 RealType thisMass = thisTriangle.getFacetMass();
173
174 // gamma is the drag coefficient normal to the face of the triangle
175 RealType gamma = thermalConductivity_ * thisMass * thisArea
176 / (2.0 * thermalLength_ * OOPSEConstant::kB);
177
178 gamma *= fabs(1.0 - targetTemp_/instaTemp);
179
180 RealType extPressure = - (targetPressure_ * thisArea) /
181 OOPSEConstant::energyConvert;
182
183 RealType randomForce = randNums[thisFacet++] * sqrt(gamma);
184 RealType dragForce = -gamma * dot(facetVel, unitNormal);
185
186 Vector3d langevinForce = (extPressure + randomForce + dragForce) *
187 unitNormal;
188
189 // Apply triangle force to stuntdouble vertices
190 for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex){
191 if ((*vertex) != NULL){
192 Vector3d vertexForce = langevinForce / 3.0;
193 (*vertex)->addFrc(vertexForce);
194 if ((*vertex)->isDirectional()){
195 Vector3d vertexPos = (*vertex)->getPos();
196 Vector3d vertexCentroidVector = vertexPos - centroid;
197 (*vertex)->addTrq(cross(vertexCentroidVector,vertexForce));
198 }
199 }
200 }
201 }
202
203 veloMunge->removeComDrift();
204 veloMunge->removeAngularDrift();
205
206 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
207 currSnapshot->setVolume(surfaceMesh_->getVolume());
208 ForceManager::postCalculation(needStress);
209 }
210
211
212 std::vector<RealType> SMIPDForceManager::genTriangleForces(int nTriangles,
213 RealType variance)
214 {
215
216 // zero fill the random vector before starting:
217 std::vector<RealType> gaussRand;
218 gaussRand.resize(nTriangles);
219 std::fill(gaussRand.begin(), gaussRand.end(), 0.0);
220
221 #ifdef IS_MPI
222 if (worldRank == 0) {
223 #endif
224 for (int i = 0; i < nTriangles; i++) {
225 gaussRand[i] = randNumGen_.randNorm(0.0, variance);
226 }
227 #ifdef IS_MPI
228 }
229 #endif
230
231 // push these out to the other processors
232
233 #ifdef IS_MPI
234 if (worldRank == 0) {
235 MPI::COMM_WORLD.Bcast(&gaussRand[0], nTriangles, MPI::REALTYPE, 0);
236 } else {
237 MPI::COMM_WORLD.Bcast(&gaussRand[0], nTriangles, MPI::REALTYPE, 0);
238 }
239 #endif
240
241 return gaussRand;
242 }
243 }

Properties

Name Value
svn:executable *