ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/integrators/LangevinHullForceManager.cpp
Revision: 1629
Committed: Wed Sep 14 21:15:17 2011 UTC (13 years, 7 months ago) by gezelter
File size: 8970 byte(s)
Log Message:
Merging changes from old branch into development branch

File Contents

# Content
1 /*
2 * Copyright (c) 2008, 2009, 2010 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. Redistributions of source code must retain the above copyright
10 * notice, this list of conditions and the following disclaimer.
11 *
12 * 2. Redistributions in binary form must reproduce the above copyright
13 * notice, this list of conditions and the following disclaimer in the
14 * documentation and/or other materials provided with the
15 * distribution.
16 *
17 * This software is provided "AS IS," without a warranty of any
18 * kind. All express or implied conditions, representations and
19 * warranties, including any implied warranty of merchantability,
20 * fitness for a particular purpose or non-infringement, are hereby
21 * excluded. The University of Notre Dame and its licensors shall not
22 * be liable for any damages suffered by licensee as a result of
23 * using, modifying or distributing the software or its
24 * derivatives. In no event will the University of Notre Dame or its
25 * licensors be liable for any lost revenue, profit or data, or for
26 * direct, indirect, special, consequential, incidental or punitive
27 * damages, however caused and regardless of the theory of liability,
28 * arising out of the use of or inability to use software, even if the
29 * University of Notre Dame has been advised of the possibility of
30 * such damages.
31 *
32 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 * research, please cite the appropriate papers when you publish your
34 * work. Good starting points are:
35 *
36 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 * [4] Vardeman & Gezelter, in progress (2009).
40 */
41 #include <fstream>
42 #include <iostream>
43 #include "integrators/LangevinHullForceManager.hpp"
44 #include "utils/PhysicalConstants.hpp"
45 #include "math/ConvexHull.hpp"
46 #include "math/AlphaHull.hpp"
47 #include "math/Triangle.hpp"
48 #include "math/CholeskyDecomposition.hpp"
49 #ifdef IS_MPI
50 #include <mpi.h>
51 #endif
52
53 namespace OpenMD {
54
55 LangevinHullForceManager::LangevinHullForceManager(SimInfo* info) : ForceManager(info) {
56
57 simParams = info->getSimParams();
58 veloMunge = new Velocitizer(info);
59
60 // Create Hull, Convex Hull for now, other options later.
61
62 stringToEnumMap_["Convex"] = hullConvex;
63 stringToEnumMap_["AlphaShape"] = hullAlphaShape;
64 stringToEnumMap_["Unknown"] = hullUnknown;
65
66 const std::string ht = simParams->getHULL_Method();
67
68
69
70 std::map<std::string, HullTypeEnum>::iterator iter;
71 iter = stringToEnumMap_.find(ht);
72 hullType_ = (iter == stringToEnumMap_.end()) ? LangevinHullForceManager::hullUnknown : iter->second;
73 if (hullType_ == hullUnknown) {
74 std::cerr << "WARNING! Hull Type Unknown!\n";
75 }
76
77 switch(hullType_) {
78 case hullConvex :
79 surfaceMesh_ = new ConvexHull();
80 break;
81 case hullAlphaShape :
82 surfaceMesh_ = new AlphaHull(simParams->getAlpha());
83 break;
84 case hullUnknown :
85 default :
86 break;
87 }
88 /* Check that the simulation has target pressure and target
89 temperature set */
90
91 if (!simParams->haveTargetTemp()) {
92 sprintf(painCave.errMsg,
93 "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
94 "\twithout a targetTemp (K)!\n");
95 painCave.isFatal = 1;
96 painCave.severity = OPENMD_ERROR;
97 simError();
98 } else {
99 targetTemp_ = simParams->getTargetTemp();
100 }
101
102 if (!simParams->haveTargetPressure()) {
103 sprintf(painCave.errMsg,
104 "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
105 "\twithout a targetPressure (atm)!\n");
106 painCave.isFatal = 1;
107 simError();
108 } else {
109 // Convert pressure from atm -> amu/(fs^2*Ang)
110 targetPressure_ = simParams->getTargetPressure() /
111 PhysicalConstants::pressureConvert;
112 }
113
114 if (simParams->getUsePeriodicBoundaryConditions()) {
115 sprintf(painCave.errMsg,
116 "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
117 "\twith periodic boundary conditions!\n");
118 painCave.isFatal = 1;
119 simError();
120 }
121
122 if (!simParams->haveViscosity()) {
123 sprintf(painCave.errMsg,
124 "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
125 "\twithout a viscosity!\n");
126 painCave.isFatal = 1;
127 painCave.severity = OPENMD_ERROR;
128 simError();
129 }else{
130 viscosity_ = simParams->getViscosity();
131 }
132
133 dt_ = simParams->getDt();
134
135 variance_ = 2.0 * PhysicalConstants::kb * targetTemp_ / dt_;
136
137 // Build a vector of integrable objects to determine if the are
138 // surface atoms
139 Molecule* mol;
140 StuntDouble* integrableObject;
141 SimInfo::MoleculeIterator i;
142 Molecule::IntegrableObjectIterator j;
143
144 for (mol = info_->beginMolecule(i); mol != NULL;
145 mol = info_->nextMolecule(i)) {
146 for (integrableObject = mol->beginIntegrableObject(j);
147 integrableObject != NULL;
148 integrableObject = mol->nextIntegrableObject(j)) {
149 localSites_.push_back(integrableObject);
150 }
151 }
152 }
153
154 void LangevinHullForceManager::postCalculation(){
155 SimInfo::MoleculeIterator i;
156 Molecule::IntegrableObjectIterator j;
157 Molecule* mol;
158 StuntDouble* integrableObject;
159
160 // Compute surface Mesh
161 surfaceMesh_->computeHull(localSites_);
162
163 // Get total area and number of surface stunt doubles
164 RealType area = surfaceMesh_->getArea();
165 std::vector<Triangle> sMesh = surfaceMesh_->getMesh();
166 int nTriangles = sMesh.size();
167
168 // Generate all of the necessary random forces
169 std::vector<Vector3d> randNums = genTriangleForces(nTriangles, variance_);
170
171 // Loop over the mesh faces and apply external pressure to each
172 // of the faces
173 std::vector<Triangle>::iterator face;
174 std::vector<StuntDouble*>::iterator vertex;
175 int thisFacet = 0;
176 for (face = sMesh.begin(); face != sMesh.end(); ++face){
177 Triangle thisTriangle = *face;
178 std::vector<StuntDouble*> vertexSDs = thisTriangle.getVertices();
179 RealType thisArea = thisTriangle.getArea();
180 Vector3d unitNormal = thisTriangle.getUnitNormal();
181 //unitNormal.normalize();
182 Vector3d centroid = thisTriangle.getCentroid();
183 Vector3d facetVel = thisTriangle.getFacetVelocity();
184 RealType thisMass = thisTriangle.getFacetMass();
185 Mat3x3d hydroTensor = thisTriangle.computeHydrodynamicTensor(viscosity_);
186
187 hydroTensor *= PhysicalConstants::viscoConvert;
188 Mat3x3d S;
189 CholeskyDecomposition(hydroTensor, S);
190
191 Vector3d extPressure = -unitNormal * (targetPressure_ * thisArea) /
192 PhysicalConstants::energyConvert;
193
194 Vector3d randomForce = S * randNums[thisFacet++];
195 Vector3d dragForce = -hydroTensor * facetVel;
196
197 Vector3d langevinForce = (extPressure + randomForce + dragForce);
198
199 // Apply triangle force to stuntdouble vertices
200 for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex){
201 if ((*vertex) != NULL){
202 Vector3d vertexForce = langevinForce / 3.0;
203 (*vertex)->addFrc(vertexForce);
204 }
205 }
206 }
207
208 veloMunge->removeComDrift();
209 veloMunge->removeAngularDrift();
210
211 Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
212 currSnapshot->setVolume(surfaceMesh_->getVolume());
213 ForceManager::postCalculation();
214 }
215
216
217 std::vector<Vector3d> LangevinHullForceManager::genTriangleForces(int nTriangles,
218 RealType variance)
219 {
220
221 // zero fill the random vector before starting:
222 std::vector<Vector3d> gaussRand;
223 gaussRand.resize(nTriangles);
224 std::fill(gaussRand.begin(), gaussRand.end(), V3Zero);
225
226 #ifdef IS_MPI
227 if (worldRank == 0) {
228 #endif
229 for (int i = 0; i < nTriangles; i++) {
230 gaussRand[i][0] = randNumGen_.randNorm(0.0, variance);
231 gaussRand[i][1] = randNumGen_.randNorm(0.0, variance);
232 gaussRand[i][2] = randNumGen_.randNorm(0.0, variance);
233 }
234 #ifdef IS_MPI
235 }
236 #endif
237
238 // push these out to the other processors
239
240 #ifdef IS_MPI
241 if (worldRank == 0) {
242 MPI::COMM_WORLD.Bcast(&gaussRand[0], nTriangles*3, MPI::REALTYPE, 0);
243 } else {
244 MPI::COMM_WORLD.Bcast(&gaussRand[0], nTriangles*3, MPI::REALTYPE, 0);
245 }
246 #endif
247
248 return gaussRand;
249 }
250 }

Properties

Name Value
svn:eol-style native
svn:executable *