ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/integrators/SMIPDForceManager.cpp
(Generate patch)

Comparing trunk/OOPSE-4/src/integrators/SMIPDForceManager.cpp (file contents):
Revision 3481 by gezelter, Wed Nov 26 14:26:17 2008 UTC vs.
Revision 3550 by gezelter, Wed Oct 21 02:49:43 2009 UTC

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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines