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 "math/CholeskyDecomposition.hpp" |
45 |
#include "utils/OOPSEConstant.hpp" |
46 |
#include "hydrodynamics/Sphere.hpp" |
47 |
#include "hydrodynamics/Ellipsoid.hpp" |
48 |
#include "utils/ElementsTable.hpp" |
49 |
#include "math/ConvexHull.hpp" |
50 |
#include "math/Triangle.hpp" |
51 |
|
52 |
|
53 |
namespace oopse { |
54 |
|
55 |
SMIPDForceManager::SMIPDForceManager(SimInfo* info) : ForceManager(info), forceTolerance_(1e-6), maxIterNum_(4) { |
56 |
simParams = info->getSimParams(); |
57 |
veloMunge = new Velocitizer(info); |
58 |
|
59 |
// Create Hull, Convex Hull for now, other options later. |
60 |
surfaceMesh_ = new ConvexHull(); |
61 |
|
62 |
|
63 |
|
64 |
/* Check that the simulation has target pressure and target |
65 |
temperature set*/ |
66 |
|
67 |
if (!simParams->haveTargetTemp()) { |
68 |
sprintf(painCave.errMsg, "You can't use the SMIPDynamics integrator without a targetTemp!\n"); |
69 |
painCave.isFatal = 1; |
70 |
painCave.severity = OOPSE_ERROR; |
71 |
simError(); |
72 |
} else { |
73 |
targetTemp_ = simParams->getTargetTemp(); |
74 |
} |
75 |
|
76 |
if (!simParams->haveTargetPressure()) { |
77 |
sprintf(painCave.errMsg, "SMIPDynamics error: You can't use the SMIPD integrator\n" |
78 |
" without a targetPressure!\n"); |
79 |
|
80 |
painCave.isFatal = 1; |
81 |
simError(); |
82 |
} else { |
83 |
targetPressure_ = simParams->getTargetPressure()/OOPSEConstant::pressureConvert; |
84 |
} |
85 |
|
86 |
|
87 |
if (simParams->getUsePeriodicBoundaryConditions()) { |
88 |
sprintf(painCave.errMsg, "SMIPDynamics error: You can't use the SMIPD integrator\n" |
89 |
" with periodic boundary conditions !\n"); |
90 |
|
91 |
painCave.isFatal = 1; |
92 |
simError(); |
93 |
} |
94 |
|
95 |
if (!simParams->haveViscosity()) { |
96 |
sprintf(painCave.errMsg, "You can't use SMIPDynamics without a viscosity!\n"); |
97 |
painCave.isFatal = 1; |
98 |
painCave.severity = OOPSE_ERROR; |
99 |
simError(); |
100 |
} |
101 |
|
102 |
|
103 |
|
104 |
|
105 |
//Compute initial hull |
106 |
/* |
107 |
surfaceMesh_->computeHull(localSites_); |
108 |
Area0_ = surfaceMesh_->getArea(); |
109 |
int nTriangles = surfaceMesh_->getNMeshElements(); |
110 |
// variance_ = 2.0 * OOPSEConstant::kb*simParams->getTargetTemp()/simParams->getDt(); |
111 |
gamma_0_ = (NumericConstant::PI * targetPressure_* targetPressure_ * Area0_ * Area0_ * simParams->getDt()) / |
112 |
(4.0 * nTriangles * nTriangles* OOPSEConstant::kb*simParams->getTargetTemp()); |
113 |
//RealType eta0 = gamma_0/ |
114 |
*/ |
115 |
|
116 |
// Build the hydroProp map: |
117 |
std::map<std::string, HydroProp*> hydroPropMap; |
118 |
|
119 |
Molecule* mol; |
120 |
StuntDouble* integrableObject; |
121 |
SimInfo::MoleculeIterator i; |
122 |
Molecule::IntegrableObjectIterator j; |
123 |
bool needHydroPropFile = false; |
124 |
|
125 |
for (mol = info->beginMolecule(i); mol != NULL; |
126 |
mol = info->nextMolecule(i)) { |
127 |
for (integrableObject = mol->beginIntegrableObject(j); |
128 |
integrableObject != NULL; |
129 |
integrableObject = mol->nextIntegrableObject(j)) { |
130 |
|
131 |
if (integrableObject->isRigidBody()) { |
132 |
RigidBody* rb = static_cast<RigidBody*>(integrableObject); |
133 |
if (rb->getNumAtoms() > 1) needHydroPropFile = true; |
134 |
} |
135 |
|
136 |
} |
137 |
} |
138 |
|
139 |
|
140 |
if (needHydroPropFile) { |
141 |
if (simParams->haveHydroPropFile()) { |
142 |
hydroPropMap = parseFrictionFile(simParams->getHydroPropFile()); |
143 |
} else { |
144 |
sprintf( painCave.errMsg, |
145 |
"HydroPropFile must be set to a file name if SMIPDynamics\n" |
146 |
"\tis specified for rigidBodies which contain more than one atom\n" |
147 |
"\tTo create a HydroPropFile, run the \"Hydro\" program.\n\n" |
148 |
"\tFor use with SMIPD, the default viscosity in Hydro should be\n" |
149 |
"\tset to 1.0 because the friction and random forces will be\n" |
150 |
"\tdynamically re-set assuming this is true.\n"); |
151 |
painCave.severity = OOPSE_ERROR; |
152 |
painCave.isFatal = 1; |
153 |
simError(); |
154 |
} |
155 |
|
156 |
for (mol = info->beginMolecule(i); mol != NULL; |
157 |
mol = info->nextMolecule(i)) { |
158 |
for (integrableObject = mol->beginIntegrableObject(j); |
159 |
integrableObject != NULL; |
160 |
integrableObject = mol->nextIntegrableObject(j)) { |
161 |
|
162 |
std::map<std::string, HydroProp*>::iterator iter = hydroPropMap.find(integrableObject->getType()); |
163 |
if (iter != hydroPropMap.end()) { |
164 |
hydroProps_.push_back(iter->second); |
165 |
} else { |
166 |
sprintf( painCave.errMsg, |
167 |
"Can not find resistance tensor for atom [%s]\n", integrableObject->getType().c_str()); |
168 |
painCave.severity = OOPSE_ERROR; |
169 |
painCave.isFatal = 1; |
170 |
simError(); |
171 |
} |
172 |
} |
173 |
} |
174 |
} else { |
175 |
|
176 |
std::map<std::string, HydroProp*> hydroPropMap; |
177 |
for (mol = info->beginMolecule(i); mol != NULL; |
178 |
mol = info->nextMolecule(i)) { |
179 |
for (integrableObject = mol->beginIntegrableObject(j); |
180 |
integrableObject != NULL; |
181 |
integrableObject = mol->nextIntegrableObject(j)) { |
182 |
Shape* currShape = NULL; |
183 |
|
184 |
if (integrableObject->isAtom()){ |
185 |
Atom* atom = static_cast<Atom*>(integrableObject); |
186 |
AtomType* atomType = atom->getAtomType(); |
187 |
if (atomType->isGayBerne()) { |
188 |
DirectionalAtomType* dAtomType = dynamic_cast<DirectionalAtomType*>(atomType); |
189 |
GenericData* data = dAtomType->getPropertyByName("GayBerne"); |
190 |
if (data != NULL) { |
191 |
GayBerneParamGenericData* gayBerneData = dynamic_cast<GayBerneParamGenericData*>(data); |
192 |
|
193 |
if (gayBerneData != NULL) { |
194 |
GayBerneParam gayBerneParam = gayBerneData->getData(); |
195 |
currShape = new Ellipsoid(V3Zero, |
196 |
gayBerneParam.GB_l / 2.0, |
197 |
gayBerneParam.GB_d / 2.0, |
198 |
Mat3x3d::identity()); |
199 |
} else { |
200 |
sprintf( painCave.errMsg, |
201 |
"Can not cast GenericData to GayBerneParam\n"); |
202 |
painCave.severity = OOPSE_ERROR; |
203 |
painCave.isFatal = 1; |
204 |
simError(); |
205 |
} |
206 |
} else { |
207 |
sprintf( painCave.errMsg, "Can not find Parameters for GayBerne\n"); |
208 |
painCave.severity = OOPSE_ERROR; |
209 |
painCave.isFatal = 1; |
210 |
simError(); |
211 |
} |
212 |
} else { |
213 |
if (atomType->isLennardJones()){ |
214 |
GenericData* data = atomType->getPropertyByName("LennardJones"); |
215 |
if (data != NULL) { |
216 |
LJParamGenericData* ljData = dynamic_cast<LJParamGenericData*>(data); |
217 |
if (ljData != NULL) { |
218 |
LJParam ljParam = ljData->getData(); |
219 |
currShape = new Sphere(atom->getPos(), 2.0); |
220 |
} else { |
221 |
sprintf( painCave.errMsg, |
222 |
"Can not cast GenericData to LJParam\n"); |
223 |
painCave.severity = OOPSE_ERROR; |
224 |
painCave.isFatal = 1; |
225 |
simError(); |
226 |
} |
227 |
} |
228 |
} else { |
229 |
int aNum = etab.GetAtomicNum((atom->getType()).c_str()); |
230 |
if (aNum != 0) { |
231 |
currShape = new Sphere(atom->getPos(), 2.0); |
232 |
} else { |
233 |
sprintf( painCave.errMsg, |
234 |
"Could not find atom type in default element.txt\n"); |
235 |
painCave.severity = OOPSE_ERROR; |
236 |
painCave.isFatal = 1; |
237 |
simError(); |
238 |
} |
239 |
} |
240 |
} |
241 |
} |
242 |
HydroProp* currHydroProp = currShape->getHydroProp(1.0,simParams->getTargetTemp()); |
243 |
std::map<std::string, HydroProp*>::iterator iter = hydroPropMap.find(integrableObject->getType()); |
244 |
if (iter != hydroPropMap.end()) |
245 |
hydroProps_.push_back(iter->second); |
246 |
else { |
247 |
currHydroProp->complete(); |
248 |
hydroPropMap.insert(std::map<std::string, HydroProp*>::value_type(integrableObject->getType(), currHydroProp)); |
249 |
hydroProps_.push_back(currHydroProp); |
250 |
} |
251 |
} |
252 |
} |
253 |
} |
254 |
|
255 |
/* Compute hull first time through to get the area of t=0*/ |
256 |
|
257 |
//Build a vector of integrable objects to determine if the are surface atoms |
258 |
for (mol = info_->beginMolecule(i); mol != NULL; mol = info_->nextMolecule(i)) { |
259 |
for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; |
260 |
integrableObject = mol->nextIntegrableObject(j)) { |
261 |
localSites_.push_back(integrableObject); |
262 |
} |
263 |
} |
264 |
|
265 |
|
266 |
} |
267 |
|
268 |
std::map<std::string, HydroProp*> SMIPDForceManager::parseFrictionFile(const std::string& filename) { |
269 |
std::map<std::string, HydroProp*> props; |
270 |
std::ifstream ifs(filename.c_str()); |
271 |
if (ifs.is_open()) { |
272 |
|
273 |
} |
274 |
|
275 |
const unsigned int BufferSize = 65535; |
276 |
char buffer[BufferSize]; |
277 |
while (ifs.getline(buffer, BufferSize)) { |
278 |
HydroProp* currProp = new HydroProp(buffer); |
279 |
props.insert(std::map<std::string, HydroProp*>::value_type(currProp->getName(), currProp)); |
280 |
} |
281 |
|
282 |
return props; |
283 |
} |
284 |
|
285 |
void SMIPDForceManager::postCalculation(bool needStress){ |
286 |
SimInfo::MoleculeIterator i; |
287 |
Molecule::IntegrableObjectIterator j; |
288 |
Molecule* mol; |
289 |
StuntDouble* integrableObject; |
290 |
RealType mass; |
291 |
Vector3d pos; |
292 |
Vector3d frc; |
293 |
Mat3x3d A; |
294 |
Mat3x3d Atrans; |
295 |
Vector3d Tb; |
296 |
Vector3d ji; |
297 |
unsigned int index = 0; |
298 |
int fdf; |
299 |
|
300 |
fdf = 0; |
301 |
|
302 |
/*Compute surface Mesh*/ |
303 |
surfaceMesh_->computeHull(localSites_); |
304 |
|
305 |
/* Get area and number of surface stunt doubles and compute new variance */ |
306 |
RealType area = surfaceMesh_->getArea(); |
307 |
int nSurfaceSDs = surfaceMesh_->getNs(); |
308 |
|
309 |
|
310 |
std::vector<Triangle> sMesh = surfaceMesh_->getMesh(); |
311 |
int nTriangles = sMesh.size(); |
312 |
|
313 |
|
314 |
|
315 |
/* Compute variance for random forces */ |
316 |
|
317 |
|
318 |
|
319 |
|
320 |
std::vector<RealType> randNums = genTriangleForces(nTriangles, 1.0); |
321 |
|
322 |
/* Loop over the mesh faces and apply random force to each of the faces*/ |
323 |
|
324 |
|
325 |
std::vector<Triangle>::iterator face; |
326 |
std::vector<StuntDouble*>::iterator vertex; |
327 |
int thisNumber = 0; |
328 |
for (face = sMesh.begin(); face != sMesh.end(); ++face){ |
329 |
|
330 |
Triangle thisTriangle = *face; |
331 |
std::vector<StuntDouble*> vertexSDs = thisTriangle.getVertices(); |
332 |
RealType thisArea = thisTriangle.getArea(); |
333 |
// RealType sigma_t = sqrt(NumericConstant::PI/2.0)*((targetPressure_)*thisArea) /OOPSEConstant::energyConvert; |
334 |
// gamma_t_ = (NumericConstant::PI * targetPressure_* targetPressure_ * thisArea * thisArea * simParams->getDt()) /(4.0 * OOPSEConstant::kB*simParams->getTargetTemp()); |
335 |
|
336 |
/* Get Random Force */ |
337 |
Vector3d unitNormal = thisTriangle.getNormal(); |
338 |
unitNormal.normalize(); |
339 |
//Vector3d randomForce = -randNums[thisNumber] * sigma_t * unitNormal; |
340 |
Vector3d centroid = thisTriangle.getCentroid(); |
341 |
Vector3d facetVel = thisTriangle.getFacetVelocity(); |
342 |
RealType hydroLength = thisTriangle.getIncircleRadius()*2.0/3.14; |
343 |
|
344 |
RealType f_normal = simParams->getViscosity()*hydroLength*1.439326479e4; |
345 |
RealType extPressure = -(targetPressure_ * thisArea)/OOPSEConstant::energyConvert; |
346 |
RealType randomForce = randNums[thisNumber++] * sqrt(2.0 * f_normal * OOPSEConstant::kb*targetTemp_/simParams->getDt()); |
347 |
|
348 |
RealType dragForce = -f_normal * dot(facetVel, unitNormal); |
349 |
|
350 |
|
351 |
Vector3d langevinForce = (extPressure + randomForce + dragForce) * unitNormal; |
352 |
|
353 |
// Vector3d dragForce = - gamma_t_ * dot(facetVel, unitNormal) * unitNormal / OOPSEConstant::energyConvert; |
354 |
|
355 |
// std::cout << " " << randomForce << " " << f_normal << std::endl; |
356 |
|
357 |
|
358 |
for (vertex = vertexSDs.begin(); vertex != vertexSDs.end(); ++vertex){ |
359 |
if ((*vertex) != NULL){ |
360 |
// mass = integrableObject->getMass(); |
361 |
Vector3d vertexForce = langevinForce/3.0; |
362 |
(*vertex)->addFrc(vertexForce); |
363 |
|
364 |
if ((*vertex)->isDirectional()){ |
365 |
|
366 |
Vector3d vertexPos = (*vertex)->getPos(); |
367 |
Vector3d vertexCentroidVector = vertexPos - centroid; |
368 |
(*vertex)->addTrq(cross(vertexCentroidVector,vertexForce)); |
369 |
} |
370 |
} |
371 |
} |
372 |
} |
373 |
|
374 |
/* Now loop over all surface particles and apply the drag*/ |
375 |
/* |
376 |
std::vector<StuntDouble*> surfaceSDs = surfaceMesh_->getSurfaceAtoms(); |
377 |
for (vertex = surfaceSDs.begin(); vertex != surfaceSDs.end(); ++vertex){ |
378 |
integrableObject = *vertex; |
379 |
mass = integrableObject->getMass(); |
380 |
if (integrableObject->isDirectional()){ |
381 |
|
382 |
// preliminaries for directional objects: |
383 |
|
384 |
A = integrableObject->getA(); |
385 |
Atrans = A.transpose(); |
386 |
Vector3d rcrLab = Atrans * hydroProps_[index]->getCOR(); |
387 |
//apply random force and torque at center of resistance |
388 |
Mat3x3d I = integrableObject->getI(); |
389 |
Vector3d omegaBody; |
390 |
|
391 |
// What remains contains velocity explicitly, but the velocity required |
392 |
// is at the full step: v(t + h), while we have initially the velocity |
393 |
// at the half step: v(t + h/2). We need to iterate to converge the |
394 |
// friction force and friction torque vectors. |
395 |
|
396 |
// this is the velocity at the half-step: |
397 |
|
398 |
Vector3d vel =integrableObject->getVel(); |
399 |
Vector3d angMom = integrableObject->getJ(); |
400 |
|
401 |
//estimate velocity at full-step using everything but friction forces: |
402 |
|
403 |
frc = integrableObject->getFrc(); |
404 |
Vector3d velStep = vel + (dt2_ /mass * OOPSEConstant::energyConvert) * frc; |
405 |
|
406 |
Tb = integrableObject->lab2Body(integrableObject->getTrq()); |
407 |
Vector3d angMomStep = angMom + (dt2_ * OOPSEConstant::energyConvert) * Tb; |
408 |
|
409 |
Vector3d omegaLab; |
410 |
Vector3d vcdLab; |
411 |
Vector3d vcdBody; |
412 |
Vector3d frictionForceBody; |
413 |
Vector3d frictionForceLab(0.0); |
414 |
Vector3d oldFFL; // used to test for convergence |
415 |
Vector3d frictionTorqueBody(0.0); |
416 |
Vector3d oldFTB; // used to test for convergence |
417 |
Vector3d frictionTorqueLab; |
418 |
RealType fdot; |
419 |
RealType tdot; |
420 |
|
421 |
//iteration starts here: |
422 |
|
423 |
for (int k = 0; k < maxIterNum_; k++) { |
424 |
|
425 |
if (integrableObject->isLinear()) { |
426 |
int linearAxis = integrableObject->linearAxis(); |
427 |
int l = (linearAxis +1 )%3; |
428 |
int m = (linearAxis +2 )%3; |
429 |
omegaBody[l] = angMomStep[l] /I(l, l); |
430 |
omegaBody[m] = angMomStep[m] /I(m, m); |
431 |
|
432 |
} else { |
433 |
omegaBody[0] = angMomStep[0] /I(0, 0); |
434 |
omegaBody[1] = angMomStep[1] /I(1, 1); |
435 |
omegaBody[2] = angMomStep[2] /I(2, 2); |
436 |
} |
437 |
|
438 |
omegaLab = Atrans * omegaBody; |
439 |
|
440 |
// apply friction force and torque at center of resistance |
441 |
|
442 |
vcdLab = velStep + cross(omegaLab, rcrLab); |
443 |
vcdBody = A * vcdLab; |
444 |
frictionForceBody = -(hydroProps_[index]->getXitt() * vcdBody + hydroProps_[index]->getXirt() * omegaBody); |
445 |
oldFFL = frictionForceLab; |
446 |
frictionForceLab = Atrans * frictionForceBody; |
447 |
oldFTB = frictionTorqueBody; |
448 |
frictionTorqueBody = -(hydroProps_[index]->getXitr() * vcdBody + hydroProps_[index]->getXirr() * omegaBody); |
449 |
frictionTorqueLab = Atrans * frictionTorqueBody; |
450 |
|
451 |
// re-estimate velocities at full-step using friction forces: |
452 |
|
453 |
velStep = vel + (dt2_ / mass * OOPSEConstant::energyConvert) * (frc + frictionForceLab); |
454 |
angMomStep = angMom + (dt2_ * OOPSEConstant::energyConvert) * (Tb + frictionTorqueBody); |
455 |
|
456 |
// check for convergence (if the vectors have converged, fdot and tdot will both be 1.0): |
457 |
|
458 |
fdot = dot(frictionForceLab, oldFFL) / frictionForceLab.lengthSquare(); |
459 |
tdot = dot(frictionTorqueBody, oldFTB) / frictionTorqueBody.lengthSquare(); |
460 |
|
461 |
if (fabs(1.0 - fdot) <= forceTolerance_ && fabs(1.0 - tdot) <= forceTolerance_) |
462 |
break; // iteration ends here |
463 |
} |
464 |
|
465 |
integrableObject->addFrc(frictionForceLab); |
466 |
integrableObject->addTrq(frictionTorqueLab + cross(rcrLab, frictionForceLab)); |
467 |
|
468 |
|
469 |
} else { |
470 |
//spherical atom |
471 |
|
472 |
// What remains contains velocity explicitly, but the velocity required |
473 |
// is at the full step: v(t + h), while we have initially the velocity |
474 |
// at the half step: v(t + h/2). We need to iterate to converge the |
475 |
// friction force vector. |
476 |
|
477 |
// this is the velocity at the half-step: |
478 |
|
479 |
Vector3d vel =integrableObject->getVel(); |
480 |
|
481 |
//estimate velocity at full-step using everything but friction forces: |
482 |
|
483 |
frc = integrableObject->getFrc(); |
484 |
Vector3d velStep = vel + (dt2_ / mass * OOPSEConstant::energyConvert) * frc; |
485 |
|
486 |
Vector3d frictionForce(0.0); |
487 |
Vector3d oldFF; // used to test for convergence |
488 |
RealType fdot; |
489 |
|
490 |
//iteration starts here: |
491 |
|
492 |
for (int k = 0; k < maxIterNum_; k++) { |
493 |
|
494 |
oldFF = frictionForce; |
495 |
frictionForce = -hydroProps_[index]->getXitt() * velStep; |
496 |
//frictionForce = -gamma_t*velStep; |
497 |
// re-estimate velocities at full-step using friction forces: |
498 |
|
499 |
velStep = vel + (dt2_ / mass * OOPSEConstant::energyConvert) * (frc + frictionForce); |
500 |
|
501 |
// check for convergence (if the vector has converged, fdot will be 1.0): |
502 |
|
503 |
fdot = dot(frictionForce, oldFF) / frictionForce.lengthSquare(); |
504 |
|
505 |
if (fabs(1.0 - fdot) <= forceTolerance_) |
506 |
break; // iteration ends here |
507 |
} |
508 |
|
509 |
integrableObject->addFrc(frictionForce); |
510 |
|
511 |
|
512 |
} |
513 |
|
514 |
|
515 |
} |
516 |
*/ |
517 |
Snapshot* currSnapshot = info_->getSnapshotManager()->getCurrentSnapshot(); |
518 |
currSnapshot->setVolume(surfaceMesh_->getVolume()); |
519 |
ForceManager::postCalculation(needStress); |
520 |
} |
521 |
|
522 |
void SMIPDForceManager::genRandomForceAndTorque(Vector3d& force, Vector3d& torque, unsigned int index, RealType variance) { |
523 |
|
524 |
|
525 |
Vector<RealType, 6> Z; |
526 |
Vector<RealType, 6> generalForce; |
527 |
|
528 |
|
529 |
Z[0] = randNumGen_.randNorm(0, variance); |
530 |
Z[1] = randNumGen_.randNorm(0, variance); |
531 |
Z[2] = randNumGen_.randNorm(0, variance); |
532 |
Z[3] = randNumGen_.randNorm(0, variance); |
533 |
Z[4] = randNumGen_.randNorm(0, variance); |
534 |
Z[5] = randNumGen_.randNorm(0, variance); |
535 |
|
536 |
generalForce = hydroProps_[index]->getS()*Z; |
537 |
|
538 |
force[0] = generalForce[0]; |
539 |
force[1] = generalForce[1]; |
540 |
force[2] = generalForce[2]; |
541 |
torque[0] = generalForce[3]; |
542 |
torque[1] = generalForce[4]; |
543 |
torque[2] = generalForce[5]; |
544 |
|
545 |
} |
546 |
std::vector<RealType> SMIPDForceManager::genTriangleForces(int nTriangles, RealType variance) { |
547 |
|
548 |
// zero fill the random vector before starting: |
549 |
std::vector<RealType> gaussRand; |
550 |
gaussRand.resize(nTriangles); |
551 |
std::fill(gaussRand.begin(), gaussRand.end(), 0.0); |
552 |
|
553 |
|
554 |
|
555 |
#ifdef IS_MPI |
556 |
if (worldRank == 0) { |
557 |
#endif |
558 |
for (int i = 0; i < nTriangles; i++) { |
559 |
//gaussRand[i] = fabs(randNumGen_.randNorm(0.0, 1.0)); |
560 |
gaussRand[i] = randNumGen_.randNorm(0.0, 1.0); |
561 |
} |
562 |
#ifdef IS_MPI |
563 |
} |
564 |
#endif |
565 |
|
566 |
// push these out to the other processors |
567 |
|
568 |
#ifdef IS_MPI |
569 |
if (worldRank == 0) { |
570 |
MPI_Bcast(&gaussRand[0], nTriangles, MPI_REAL, 0, MPI_COMM_WORLD); |
571 |
} else { |
572 |
MPI_Bcast(&gaussRand[0], nTriangles, MPI_REAL, 0, MPI_COMM_WORLD); |
573 |
} |
574 |
#endif |
575 |
|
576 |
for (int i = 0; i < nTriangles; i++) { |
577 |
gaussRand[i] = gaussRand[i] * variance; |
578 |
} |
579 |
|
580 |
return gaussRand; |
581 |
} |
582 |
|
583 |
|
584 |
|
585 |
|
586 |
|
587 |
} |