ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/integrators/LangevinHullForceManager.cpp
Revision: 1874
Committed: Wed May 15 15:09:35 2013 UTC (11 years, 11 months ago) by gezelter
File size: 10164 byte(s)
Log Message:
Fixed a bunch of cppcheck warnings.

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 gezelter 1858 using namespace std;
55 gezelter 1629 namespace OpenMD {
56 gezelter 1858
57     LangevinHullForceManager::LangevinHullForceManager(SimInfo* info) :
58     ForceManager(info) {
59 gezelter 1866
60 gezelter 1629 simParams = info->getSimParams();
61     veloMunge = new Velocitizer(info);
62    
63     // Create Hull, Convex Hull for now, other options later.
64    
65     stringToEnumMap_["Convex"] = hullConvex;
66     stringToEnumMap_["AlphaShape"] = hullAlphaShape;
67     stringToEnumMap_["Unknown"] = hullUnknown;
68    
69     const std::string ht = simParams->getHULL_Method();
70    
71     std::map<std::string, HullTypeEnum>::iterator iter;
72     iter = stringToEnumMap_.find(ht);
73 gezelter 1858 hullType_ = (iter == stringToEnumMap_.end()) ?
74     LangevinHullForceManager::hullUnknown : iter->second;
75    
76 gezelter 1629 switch(hullType_) {
77     case hullConvex :
78     surfaceMesh_ = new ConvexHull();
79     break;
80     case hullAlphaShape :
81     surfaceMesh_ = new AlphaHull(simParams->getAlpha());
82     break;
83     case hullUnknown :
84     default :
85 gezelter 1858 sprintf(painCave.errMsg,
86     "LangevinHallForceManager: Unknown Hull_Method was requested!\n");
87     painCave.isFatal = 1;
88     simError();
89 gezelter 1629 break;
90     }
91 gezelter 1858
92     doThermalCoupling_ = true;
93     doPressureCoupling_ = true;
94    
95 gezelter 1629 /* Check that the simulation has target pressure and target
96 gezelter 1858 temperature set */
97 gezelter 1629 if (!simParams->haveTargetTemp()) {
98     sprintf(painCave.errMsg,
99 gezelter 1858 "LangevinHullForceManager: no targetTemp (K) was set.\n"
100     "\tOpenMD is turning off the thermal coupling to the bath.\n");
101     painCave.isFatal = 0;
102     painCave.severity = OPENMD_INFO;
103 gezelter 1629 simError();
104 gezelter 1858 doThermalCoupling_ = false;
105 gezelter 1629 } else {
106     targetTemp_ = simParams->getTargetTemp();
107 gezelter 1866
108     if (!simParams->haveViscosity()) {
109     sprintf(painCave.errMsg,
110     "LangevinHullForceManager: no viscosity was set.\n"
111     "\tOpenMD is turning off the thermal coupling to the bath.\n");
112     painCave.isFatal = 0;
113     painCave.severity = OPENMD_INFO;
114     simError();
115     doThermalCoupling_ = false;
116     }else{
117     viscosity_ = simParams->getViscosity();
118     if ( fabs(viscosity_) < 1e-6 ) {
119     sprintf(painCave.errMsg,
120     "LangevinHullDynamics: The bath viscosity was set lower\n"
121     "\tthan 1e-6 poise. OpenMD is turning off the thermal\n"
122     "\tcoupling to the bath.\n");
123     painCave.isFatal = 0;
124     painCave.severity = OPENMD_INFO;
125     simError();
126     doThermalCoupling_ = false;
127     }
128     }
129 gezelter 1629 }
130     if (!simParams->haveTargetPressure()) {
131     sprintf(painCave.errMsg,
132 gezelter 1858 "LangevinHullForceManager: no targetPressure (atm) was set.\n"
133     "\tOpenMD is turning off the pressure coupling to the bath.\n");
134     painCave.isFatal = 0;
135     painCave.severity = OPENMD_INFO;
136 gezelter 1629 simError();
137 gezelter 1858 doPressureCoupling_ = false;
138 gezelter 1629 } else {
139     // Convert pressure from atm -> amu/(fs^2*Ang)
140     targetPressure_ = simParams->getTargetPressure() /
141     PhysicalConstants::pressureConvert;
142     }
143     if (simParams->getUsePeriodicBoundaryConditions()) {
144     sprintf(painCave.errMsg,
145 gezelter 1858 "LangevinHallForceManager: You can't use the Langevin Hull\n"
146     "\tintegrator with periodic boundary conditions turned on!\n");
147 gezelter 1629 painCave.isFatal = 1;
148     simError();
149     }
150 gezelter 1855
151 gezelter 1629 dt_ = simParams->getDt();
152 gezelter 1866
153 gezelter 1858 if (doThermalCoupling_)
154     variance_ = 2.0 * PhysicalConstants::kb * targetTemp_ / dt_;
155 gezelter 1629
156     // Build a vector of integrable objects to determine if the are
157     // surface atoms
158     Molecule* mol;
159 gezelter 1767 StuntDouble* sd;
160 gezelter 1629 SimInfo::MoleculeIterator i;
161     Molecule::IntegrableObjectIterator j;
162 gezelter 1858
163 gezelter 1629 for (mol = info_->beginMolecule(i); mol != NULL;
164     mol = info_->nextMolecule(i)) {
165 gezelter 1767 for (sd = mol->beginIntegrableObject(j);
166     sd != NULL;
167     sd = mol->nextIntegrableObject(j)) {
168     localSites_.push_back(sd);
169 gezelter 1629 }
170 gezelter 1858 }
171    
172     // We need to make an initial guess at the bounding box in order
173     // to compute long range forces in ForceMatrixDecomposition:
174    
175     // Compute surface Mesh
176     surfaceMesh_->computeHull(localSites_);
177 gezelter 1629 }
178 gezelter 1867
179     LangevinHullForceManager::~LangevinHullForceManager() {
180     delete surfaceMesh_;
181 gezelter 1868 delete veloMunge;
182 gezelter 1867 }
183 gezelter 1858
184 gezelter 1629 void LangevinHullForceManager::postCalculation(){
185    
186 gezelter 1858 int nTriangles, thisFacet;
187     RealType area, thisArea, thisMass;
188     vector<Triangle> sMesh;
189     Triangle thisTriangle;
190     vector<Triangle>::iterator face;
191     vector<StuntDouble*> vertexSDs;
192     vector<StuntDouble*>::iterator vertex;
193    
194     Vector3d unitNormal, centroid, facetVel;
195     Vector3d langevinForce, vertexForce;
196     Vector3d extPressure, randomForce, dragForce;
197    
198     Mat3x3d hydroTensor, S;
199     vector<Vector3d> randNums;
200    
201 gezelter 1629 // Compute surface Mesh
202     surfaceMesh_->computeHull(localSites_);
203     // Get total area and number of surface stunt doubles
204 gezelter 1858 area = surfaceMesh_->getArea();
205     sMesh = surfaceMesh_->getMesh();
206     nTriangles = sMesh.size();
207 gezelter 1629
208 gezelter 1858 if (doThermalCoupling_) {
209     // Generate all of the necessary random forces
210     randNums = genTriangleForces(nTriangles, variance_);
211     }
212    
213     // Loop over the mesh faces
214     thisFacet = 0;
215 gezelter 1629 for (face = sMesh.begin(); face != sMesh.end(); ++face){
216 gezelter 1858 thisTriangle = *face;
217     vertexSDs = thisTriangle.getVertices();
218     thisArea = thisTriangle.getArea();
219     unitNormal = thisTriangle.getUnitNormal();
220     centroid = thisTriangle.getCentroid();
221     facetVel = thisTriangle.getFacetVelocity();
222     thisMass = thisTriangle.getFacetMass();
223 gezelter 1760
224 gezelter 1858 langevinForce = V3Zero;
225 gezelter 1629
226 gezelter 1858 if (doPressureCoupling_) {
227     extPressure = -unitNormal * (targetPressure_ * thisArea) /
228     PhysicalConstants::energyConvert;
229     langevinForce += extPressure;
230     }
231    
232 gezelter 1855 if (doThermalCoupling_) {
233 gezelter 1858 hydroTensor = thisTriangle.computeHydrodynamicTensor(viscosity_);
234     hydroTensor *= PhysicalConstants::viscoConvert;
235     CholeskyDecomposition(hydroTensor, S);
236 gezelter 1855 randomForce = S * randNums[thisFacet++];
237     dragForce = -hydroTensor * facetVel;
238 gezelter 1858 langevinForce += randomForce + dragForce;
239 gezelter 1855 }
240 gezelter 1629
241     // Apply triangle force to stuntdouble vertices
242     for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex){
243     if ((*vertex) != NULL){
244 gezelter 1858 vertexForce = langevinForce / RealType(3.0);
245 gezelter 1629 (*vertex)->addFrc(vertexForce);
246     }
247     }
248     }
249    
250     veloMunge->removeComDrift();
251     veloMunge->removeAngularDrift();
252    
253     Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
254 gezelter 1858 currSnapshot->setVolume(surfaceMesh_->getVolume());
255 gezelter 1856 currSnapshot->setHullVolume(surfaceMesh_->getVolume());
256 gezelter 1629 ForceManager::postCalculation();
257     }
258    
259 gezelter 1858 vector<Vector3d> LangevinHullForceManager::genTriangleForces(int nTriangles,
260     RealType var) {
261 gezelter 1629 // zero fill the random vector before starting:
262 gezelter 1858 vector<Vector3d> gaussRand;
263 gezelter 1629 gaussRand.resize(nTriangles);
264     std::fill(gaussRand.begin(), gaussRand.end(), V3Zero);
265    
266     #ifdef IS_MPI
267     if (worldRank == 0) {
268     #endif
269     for (int i = 0; i < nTriangles; i++) {
270 gezelter 1858 gaussRand[i][0] = randNumGen_.randNorm(0.0, var);
271     gaussRand[i][1] = randNumGen_.randNorm(0.0, var);
272     gaussRand[i][2] = randNumGen_.randNorm(0.0, var);
273 gezelter 1629 }
274     #ifdef IS_MPI
275     }
276     #endif
277    
278     // push these out to the other processors
279    
280     #ifdef IS_MPI
281 gezelter 1874 // Same command on all nodes:
282     MPI::COMM_WORLD.Bcast(&gaussRand[0], nTriangles*3, MPI::REALTYPE, 0);
283 gezelter 1629 #endif
284    
285     return gaussRand;
286     }
287     }

Properties

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