ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/integrators/LangevinHullForceManager.cpp
Revision: 1855
Committed: Tue Apr 2 18:31:51 2013 UTC (12 years, 1 month ago) by gezelter
File size: 9353 byte(s)
Log Message:
Fixed a bunch of bugs in CubicSpline debug sections, ForceMatrix Decomp 
computing costhetas for non-existent rotation matrices, Hidden accumulator counts, and SMIPD non-coupling

File Contents

# User Rev Content
1 gezelter 1629 /*
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 gezelter 1850 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
39 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 1629 */
42     #include <fstream>
43     #include <iostream>
44     #include "integrators/LangevinHullForceManager.hpp"
45     #include "utils/PhysicalConstants.hpp"
46     #include "math/ConvexHull.hpp"
47     #include "math/AlphaHull.hpp"
48     #include "math/Triangle.hpp"
49     #include "math/CholeskyDecomposition.hpp"
50     #ifdef IS_MPI
51     #include <mpi.h>
52     #endif
53    
54     namespace OpenMD {
55    
56     LangevinHullForceManager::LangevinHullForceManager(SimInfo* info) : ForceManager(info) {
57    
58     simParams = info->getSimParams();
59     veloMunge = new Velocitizer(info);
60    
61     // Create Hull, Convex Hull for now, other options later.
62    
63     stringToEnumMap_["Convex"] = hullConvex;
64     stringToEnumMap_["AlphaShape"] = hullAlphaShape;
65     stringToEnumMap_["Unknown"] = hullUnknown;
66    
67     const std::string ht = simParams->getHULL_Method();
68    
69    
70    
71     std::map<std::string, HullTypeEnum>::iterator iter;
72     iter = stringToEnumMap_.find(ht);
73     hullType_ = (iter == stringToEnumMap_.end()) ? LangevinHullForceManager::hullUnknown : iter->second;
74     if (hullType_ == hullUnknown) {
75     std::cerr << "WARNING! Hull Type Unknown!\n";
76     }
77    
78     switch(hullType_) {
79     case hullConvex :
80     surfaceMesh_ = new ConvexHull();
81     break;
82     case hullAlphaShape :
83     surfaceMesh_ = new AlphaHull(simParams->getAlpha());
84     break;
85     case hullUnknown :
86     default :
87     break;
88     }
89     /* Check that the simulation has target pressure and target
90     temperature set */
91    
92     if (!simParams->haveTargetTemp()) {
93     sprintf(painCave.errMsg,
94     "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
95     "\twithout a targetTemp (K)!\n");
96     painCave.isFatal = 1;
97     painCave.severity = OPENMD_ERROR;
98     simError();
99     } else {
100     targetTemp_ = simParams->getTargetTemp();
101     }
102    
103     if (!simParams->haveTargetPressure()) {
104     sprintf(painCave.errMsg,
105     "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
106     "\twithout a targetPressure (atm)!\n");
107     painCave.isFatal = 1;
108     simError();
109     } else {
110     // Convert pressure from atm -> amu/(fs^2*Ang)
111     targetPressure_ = simParams->getTargetPressure() /
112     PhysicalConstants::pressureConvert;
113     }
114    
115     if (simParams->getUsePeriodicBoundaryConditions()) {
116     sprintf(painCave.errMsg,
117     "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
118     "\twith periodic boundary conditions!\n");
119     painCave.isFatal = 1;
120     simError();
121     }
122    
123     if (!simParams->haveViscosity()) {
124     sprintf(painCave.errMsg,
125     "LangevinHullDynamics error: You can't use the Langevin Hull integrator\n"
126     "\twithout a viscosity!\n");
127     painCave.isFatal = 1;
128     painCave.severity = OPENMD_ERROR;
129     simError();
130     }else{
131     viscosity_ = simParams->getViscosity();
132     }
133    
134 gezelter 1855 doThermalCoupling_ = true;
135     if ( fabs(viscosity_) < 1e-6 ) {
136     sprintf(painCave.errMsg,
137     "LangevinHullDynamics: The bath viscosity was set lower than\n"
138     "\t1e-6 poise. OpenMD is turning off the thermal coupling to\n"
139     "\t the bath.\n");
140     painCave.isFatal = 0;
141     painCave.severity = OPENMD_INFO;
142     simError();
143     doThermalCoupling_ = false;
144     }
145    
146 gezelter 1629 dt_ = simParams->getDt();
147    
148     variance_ = 2.0 * PhysicalConstants::kb * targetTemp_ / dt_;
149    
150     // Build a vector of integrable objects to determine if the are
151     // surface atoms
152     Molecule* mol;
153 gezelter 1767 StuntDouble* sd;
154 gezelter 1629 SimInfo::MoleculeIterator i;
155     Molecule::IntegrableObjectIterator j;
156    
157     for (mol = info_->beginMolecule(i); mol != NULL;
158     mol = info_->nextMolecule(i)) {
159 gezelter 1767 for (sd = mol->beginIntegrableObject(j);
160     sd != NULL;
161     sd = mol->nextIntegrableObject(j)) {
162     localSites_.push_back(sd);
163 gezelter 1629 }
164     }
165     }
166    
167     void LangevinHullForceManager::postCalculation(){
168    
169     // Compute surface Mesh
170     surfaceMesh_->computeHull(localSites_);
171    
172     // Get total area and number of surface stunt doubles
173     RealType area = surfaceMesh_->getArea();
174     std::vector<Triangle> sMesh = surfaceMesh_->getMesh();
175     int nTriangles = sMesh.size();
176    
177     // Generate all of the necessary random forces
178     std::vector<Vector3d> randNums = genTriangleForces(nTriangles, variance_);
179    
180     // Loop over the mesh faces and apply external pressure to each
181     // of the faces
182     std::vector<Triangle>::iterator face;
183     std::vector<StuntDouble*>::iterator vertex;
184     int thisFacet = 0;
185     for (face = sMesh.begin(); face != sMesh.end(); ++face){
186     Triangle thisTriangle = *face;
187     std::vector<StuntDouble*> vertexSDs = thisTriangle.getVertices();
188     RealType thisArea = thisTriangle.getArea();
189     Vector3d unitNormal = thisTriangle.getUnitNormal();
190     //unitNormal.normalize();
191     Vector3d centroid = thisTriangle.getCentroid();
192     Vector3d facetVel = thisTriangle.getFacetVelocity();
193     RealType thisMass = thisTriangle.getFacetMass();
194     Mat3x3d hydroTensor = thisTriangle.computeHydrodynamicTensor(viscosity_);
195    
196     hydroTensor *= PhysicalConstants::viscoConvert;
197     Mat3x3d S;
198     CholeskyDecomposition(hydroTensor, S);
199 gezelter 1760
200 gezelter 1629 Vector3d extPressure = -unitNormal * (targetPressure_ * thisArea) /
201     PhysicalConstants::energyConvert;
202    
203 gezelter 1855 Vector3d randomForce(V3Zero);
204     Vector3d dragForce(V3Zero);
205     if (doThermalCoupling_) {
206     randomForce = S * randNums[thisFacet++];
207     dragForce = -hydroTensor * facetVel;
208     }
209 gezelter 1629
210     Vector3d langevinForce = (extPressure + randomForce + dragForce);
211    
212     // Apply triangle force to stuntdouble vertices
213     for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex){
214     if ((*vertex) != NULL){
215 gezelter 1668 Vector3d vertexForce = langevinForce / RealType(3.0);
216 gezelter 1629 (*vertex)->addFrc(vertexForce);
217     }
218     }
219     }
220    
221     veloMunge->removeComDrift();
222     veloMunge->removeAngularDrift();
223    
224     Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
225     currSnapshot->setVolume(surfaceMesh_->getVolume());
226     ForceManager::postCalculation();
227     }
228    
229    
230     std::vector<Vector3d> LangevinHullForceManager::genTriangleForces(int nTriangles,
231     RealType variance)
232     {
233    
234     // zero fill the random vector before starting:
235     std::vector<Vector3d> gaussRand;
236     gaussRand.resize(nTriangles);
237     std::fill(gaussRand.begin(), gaussRand.end(), V3Zero);
238    
239     #ifdef IS_MPI
240     if (worldRank == 0) {
241     #endif
242     for (int i = 0; i < nTriangles; i++) {
243     gaussRand[i][0] = randNumGen_.randNorm(0.0, variance);
244     gaussRand[i][1] = randNumGen_.randNorm(0.0, variance);
245     gaussRand[i][2] = randNumGen_.randNorm(0.0, variance);
246     }
247     #ifdef IS_MPI
248     }
249     #endif
250    
251     // push these out to the other processors
252    
253     #ifdef IS_MPI
254     if (worldRank == 0) {
255     MPI::COMM_WORLD.Bcast(&gaussRand[0], nTriangles*3, MPI::REALTYPE, 0);
256     } else {
257     MPI::COMM_WORLD.Bcast(&gaussRand[0], nTriangles*3, MPI::REALTYPE, 0);
258     }
259     #endif
260    
261     return gaussRand;
262     }
263     }

Properties

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