ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/restraints/RestraintForceManager.cpp
Revision: 1398
Committed: Tue Dec 8 22:17:02 2009 UTC (15 years, 4 months ago) by gezelter
File size: 17562 byte(s)
Log Message:
Parallel bugfix in RestraintForceManager
Reverted back to hydrodynamics on triangular plates for SMIPDForceManager
Removed thermalLength and thermalConductivity keywords from Globals.
Bug tracking in openmdformat (not yet resolved).

File Contents

# Content
1 /*
2 * Copyright (c) 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
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
42 #include <cmath>
43 #include "restraints/RestraintForceManager.hpp"
44 #include "restraints/MolecularRestraint.hpp"
45 #include "restraints/ObjectRestraint.hpp"
46 #include "io/RestReader.hpp"
47 #include "utils/simError.h"
48 #include "utils/PhysicalConstants.hpp"
49 #include "utils/StringUtils.hpp"
50 #include "selection/SelectionEvaluator.hpp"
51 #include "selection/SelectionManager.hpp"
52 #ifdef IS_MPI
53 #include <mpi.h>
54 #endif
55
56
57 namespace OpenMD {
58
59 RestraintForceManager::RestraintForceManager(SimInfo* info): ForceManager(info) {
60
61 // order of affairs:
62 //
63 // 1) create restraints from the restraintStamps found in the MD
64 // file.
65 //
66 // 2) Create RestraintReader to parse the input files for the ideal
67 // structures. This reader will set reference structures, and will
68 // calculate molecular centers of mass, etc.
69 //
70 // 3) sit around and wait for calcForces to be called. When it comes,
71 // call the normal force manager calcForces, then loop through the
72 // restrained objects and do their restraint forces.
73
74 currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
75 Globals* simParam = info_->getSimParams();
76
77 if (simParam->haveStatusTime()){
78 restTime_ = simParam->getStatusTime();
79 } else {
80 sprintf(painCave.errMsg,
81 "Restraint warning: If you use restraints without setting\n"
82 "\tstatusTime, no restraint data will be written to the rest\n"
83 "\tfile.\n");
84 painCave.isFatal = 0;
85 simError();
86 restTime_ = simParam->getRunTime();
87 }
88
89 int nRestraintStamps = simParam->getNRestraintStamps();
90 std::vector<RestraintStamp*> stamp = simParam->getRestraintStamps();
91
92 for (int i = 0; i < nRestraintStamps; i++){
93
94 std::string myType = toUpperCopy(stamp[i]->getType());
95
96 if (myType.compare("MOLECULAR")==0){
97
98 int molIndex;
99 std::vector<Vector3d> ref;
100 Vector3d refCom;
101
102 if (!stamp[i]->haveMolIndex()) {
103 sprintf(painCave.errMsg,
104 "Restraint Error: A molecular restraint was specified\n"
105 "\twithout providing a value for molIndex.\n");
106 painCave.isFatal = 1;
107 simError();
108 } else {
109 molIndex = stamp[i]->getMolIndex();
110 }
111
112 if (molIndex < 0) {
113 sprintf(painCave.errMsg,
114 "Restraint Error: A molecular restraint was specified\n"
115 "\twith a molIndex that was less than 0\n");
116 painCave.isFatal = 1;
117 simError();
118 }
119 if (molIndex >= info_->getNGlobalMolecules()) {
120 sprintf(painCave.errMsg,
121 "Restraint Error: A molecular restraint was specified with\n"
122 "\ta molIndex that was greater than the total number of molecules\n");
123 painCave.isFatal = 1;
124 simError();
125 }
126
127 Molecule* mol = info_->getMoleculeByGlobalIndex(molIndex);
128
129 if (mol == NULL) {
130 #ifdef IS_MPI
131 // getMoleculeByGlobalIndex returns a NULL in parallel if
132 // this proc doesn't have the molecule. Do a quick check to
133 // make sure another processor is supposed to have it.
134
135 int myrank = MPI::COMM_WORLD.Get_rank();
136 if (info_->getMolToProc(molIndex) == myrank) {
137 // If we were supposed to have it but got a null, then freak out.
138 #endif
139
140 sprintf(painCave.errMsg,
141 "Restraint Error: A molecular restraint was specified, but\n"
142 "\tno molecule was found with global index %d.\n",
143 molIndex);
144 painCave.isFatal = 1;
145 simError();
146
147 #ifdef IS_MPI
148 }
149 #endif
150 }
151
152
153 #ifdef IS_MPI
154 // only handle this molecular restraint if this processor owns the
155 // molecule
156 int myrank = MPI::COMM_WORLD.Get_rank();
157 if (info_->getMolToProc(molIndex) == myrank) {
158 #endif
159
160 MolecularRestraint* rest = new MolecularRestraint();
161
162 std::string restPre("mol_");
163 std::stringstream restName;
164 restName << restPre << molIndex;
165 rest->setRestraintName(restName.str());
166
167 if (stamp[i]->haveDisplacementSpringConstant()) {
168 rest->setDisplacementForceConstant(stamp[i]->getDisplacementSpringConstant());
169 }
170 if (stamp[i]->haveTwistSpringConstant()) {
171 rest->setTwistForceConstant(stamp[i]->getTwistSpringConstant());
172 }
173 if (stamp[i]->haveSwingXSpringConstant()) {
174 rest->setSwingXForceConstant(stamp[i]->getSwingXSpringConstant());
175 }
176 if (stamp[i]->haveSwingYSpringConstant()) {
177 rest->setSwingYForceConstant(stamp[i]->getSwingYSpringConstant());
178 }
179 if (stamp[i]->haveRestrainedTwistAngle()) {
180 rest->setRestrainedTwistAngle(stamp[i]->getRestrainedTwistAngle() * M_PI/180.0);
181 }
182 if (stamp[i]->haveRestrainedSwingYAngle()) {
183 rest->setRestrainedSwingYAngle(stamp[i]->getRestrainedSwingYAngle() * M_PI/180.0);
184 }
185 if (stamp[i]->haveRestrainedSwingXAngle()) {
186 rest->setRestrainedSwingXAngle(stamp[i]->getRestrainedSwingXAngle() * M_PI/180.0);
187 }
188 if (stamp[i]->havePrint()) {
189 rest->setPrintRestraint(stamp[i]->getPrint());
190 }
191
192 restraints_.push_back(rest);
193 mol->addProperty(new RestraintData("Restraint", rest));
194 restrainedMols_.push_back(mol);
195 #ifdef IS_MPI
196 }
197 #endif
198 } else if (myType.compare("OBJECT") == 0) {
199
200 std::string objectSelection;
201
202 if (!stamp[i]->haveObjectSelection()) {
203 sprintf(painCave.errMsg,
204 "Restraint Error: An object restraint was specified\n"
205 "\twithout providing a selection script in the\n"
206 "\tobjectSelection variable.\n");
207 painCave.isFatal = 1;
208 simError();
209 } else {
210 objectSelection = stamp[i]->getObjectSelection();
211 }
212
213 SelectionEvaluator evaluator(info);
214 SelectionManager seleMan(info);
215
216 evaluator.loadScriptString(objectSelection);
217 seleMan.setSelectionSet(evaluator.evaluate());
218 int selectionCount = seleMan.getSelectionCount();
219
220 sprintf(painCave.errMsg,
221 "Restraint Info: The specified restraint objectSelection,\n"
222 "\t\t%s\n"
223 "\twill result in %d integrable objects being\n"
224 "\trestrained.\n", objectSelection.c_str(), selectionCount);
225 painCave.isFatal = 0;
226 simError();
227
228 int selei;
229 StuntDouble* sd;
230
231 for (sd = seleMan.beginSelected(selei); sd != NULL;
232 sd = seleMan.nextSelected(selei)) {
233
234 ObjectRestraint* rest = new ObjectRestraint();
235
236 if (stamp[i]->haveDisplacementSpringConstant()) {
237 rest->setDisplacementForceConstant(stamp[i]->getDisplacementSpringConstant());
238 }
239 if (stamp[i]->haveTwistSpringConstant()) {
240 rest->setTwistForceConstant(stamp[i]->getTwistSpringConstant());
241 }
242 if (stamp[i]->haveSwingXSpringConstant()) {
243 rest->setSwingXForceConstant(stamp[i]->getSwingXSpringConstant());
244 }
245 if (stamp[i]->haveSwingYSpringConstant()) {
246 rest->setSwingYForceConstant(stamp[i]->getSwingYSpringConstant());
247 }
248 if (stamp[i]->haveRestrainedTwistAngle()) {
249 rest->setRestrainedTwistAngle(stamp[i]->getRestrainedTwistAngle());
250 }
251 if (stamp[i]->haveRestrainedSwingXAngle()) {
252 rest->setRestrainedSwingXAngle(stamp[i]->getRestrainedSwingXAngle());
253 }
254 if (stamp[i]->haveRestrainedSwingYAngle()) {
255 rest->setRestrainedSwingYAngle(stamp[i]->getRestrainedSwingYAngle());
256 }
257 if (stamp[i]->havePrint()) {
258 rest->setPrintRestraint(stamp[i]->getPrint());
259 }
260
261 restraints_.push_back(rest);
262 sd->addProperty(new RestraintData("Restraint", rest));
263 restrainedObjs_.push_back(sd);
264 }
265
266 }
267 }
268
269 // ThermodynamicIntegration subclasses RestraintForceManager, and there
270 // are times when it won't use restraints at all, so only open the
271 // restraint file if we are actually using restraints:
272
273 if (simParam->getUseRestraints()) {
274 std::string refFile = simParam->getRestraint_file();
275 RestReader* rr = new RestReader(info, refFile);
276
277 rr->readReferenceStructure();
278 }
279
280 restOutput_ = getPrefix(info_->getFinalConfigFileName()) + ".rest";
281 restOut = new RestWriter(info_, restOutput_.c_str(), restraints_);
282
283 if(!restOut){
284 sprintf(painCave.errMsg, "Restraint error: Failed to create RestWriter\n");
285 painCave.isFatal = 1;
286 simError();
287 }
288
289 // todo: figure out the scale factor. Right now, just scale them all to 1
290 std::vector<Restraint*>::const_iterator resti;
291 for(resti=restraints_.begin(); resti != restraints_.end(); ++resti){
292 (*resti)->setScaleFactor(1.0);
293 }
294 }
295
296 RestraintForceManager::~RestraintForceManager(){
297 if (restOut)
298 delete restOut;
299 }
300
301 void RestraintForceManager::init() {
302 currRestTime_ = currSnapshot_->getTime();
303 }
304
305 void RestraintForceManager::calcForces(bool needPotential, bool needStress){
306
307 ForceManager::calcForces(needPotential, needStress);
308 RealType restPot_local, restPot;
309
310 restPot_local = doRestraints(1.0);
311
312 #ifdef IS_MPI
313 MPI::COMM_WORLD.Allreduce(&restPot_local, &restPot, 1,
314 MPI::REALTYPE, MPI::SUM);
315 #else
316 restPot = restPot_local;
317 #endif
318 currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
319 currSnapshot_->statData[Stats::LONG_RANGE_POTENTIAL] += restPot;
320 currSnapshot_->statData[Stats::VHARM] = restPot;
321
322 //write out forces and current positions of restrained molecules
323 if (currSnapshot_->getTime() >= currRestTime_){
324 restOut->writeRest(restInfo_);
325 currRestTime_ += restTime_;
326 }
327 }
328
329 RealType RestraintForceManager::doRestraints(RealType scalingFactor){
330 std::vector<Molecule*>::const_iterator rm;
331 GenericData* data;
332 Molecule::IntegrableObjectIterator ioi;
333 MolecularRestraint* mRest;
334 StuntDouble* sd;
335 RealType pTot;
336
337 std::vector<StuntDouble*>::const_iterator ro;
338 ObjectRestraint* oRest;
339
340 std::map<int, Restraint::RealPair> restInfo;
341
342 unscaledPotential_ = 0.0;
343
344 restInfo_.clear();
345
346 for(rm=restrainedMols_.begin(); rm != restrainedMols_.end(); ++rm){
347
348 // make sure this molecule (*rm) has a generic data for restraints:
349 data = (*rm)->getPropertyByName("Restraint");
350 if (data != NULL) {
351 // make sure we can reinterpret the generic data as restraint data:
352 RestraintData* restData= dynamic_cast<RestraintData*>(data);
353 if (restData != NULL) {
354 // make sure we can reinterpet the restraint data as a pointer to
355 // an MolecularRestraint:
356 mRest = dynamic_cast<MolecularRestraint*>(restData->getData());
357 if (mRest == NULL) {
358 sprintf( painCave.errMsg,
359 "Can not cast RestraintData to MolecularRestraint\n");
360 painCave.severity = OPENMD_ERROR;
361 painCave.isFatal = 1;
362 simError();
363 }
364 } else {
365 sprintf( painCave.errMsg,
366 "Can not cast GenericData to RestraintData\n");
367 painCave.severity = OPENMD_ERROR;
368 painCave.isFatal = 1;
369 simError();
370 }
371 } else {
372 sprintf( painCave.errMsg, "Can not find Restraint for RestrainedObject\n");
373 painCave.severity = OPENMD_ERROR;
374 painCave.isFatal = 1;
375 simError();
376 }
377
378 // phew. At this point, we should have the pointer to the
379 // correct MolecularRestraint in the variable mRest.
380
381 Vector3d molCom = (*rm)->getCom();
382
383 std::vector<Vector3d> struc;
384 std::vector<Vector3d> forces;
385
386 for(sd = (*rm)->beginIntegrableObject(ioi); sd != NULL;
387 sd = (*rm)->nextIntegrableObject(ioi)) {
388 struc.push_back(sd->getPos());
389 }
390
391 mRest->setScaleFactor(scalingFactor);
392 mRest->calcForce(struc, molCom);
393 forces = mRest->getRestraintForces();
394 int index = 0;
395
396 for(sd = (*rm)->beginIntegrableObject(ioi); sd != NULL;
397 sd = (*rm)->nextIntegrableObject(ioi)) {
398 sd->addFrc(forces[index]);
399 struc.push_back(sd->getPos());
400 index++;
401 }
402
403 unscaledPotential_ += mRest->getUnscaledPotential();
404
405 restInfo = mRest->getRestraintInfo();
406
407 // only collect data on restraints that we're going to print:
408 if (mRest->getPrintRestraint())
409 restInfo_.push_back(restInfo);
410 }
411
412 for(ro=restrainedObjs_.begin(); ro != restrainedObjs_.end(); ++ro){
413 // make sure this object (*ro) has a generic data for restraints:
414 data = (*ro)->getPropertyByName("Restraint");
415 if (data != NULL) {
416 // make sure we can reinterpret the generic data as restraint data:
417 RestraintData* restData= dynamic_cast<RestraintData*>(data);
418 if (restData != NULL) {
419 // make sure we can reinterpet the restraint data as a pointer to
420 // an ObjectRestraint:
421 oRest = dynamic_cast<ObjectRestraint*>(restData->getData());
422 if (oRest == NULL) {
423 sprintf( painCave.errMsg,
424 "Can not cast RestraintData to ObjectRestraint\n");
425 painCave.severity = OPENMD_ERROR;
426 painCave.isFatal = 1;
427 simError();
428 }
429 } else {
430 sprintf( painCave.errMsg,
431 "Can not cast GenericData to RestraintData\n");
432 painCave.severity = OPENMD_ERROR;
433 painCave.isFatal = 1;
434 simError();
435 }
436 } else {
437 sprintf( painCave.errMsg, "Can not find Restraint for RestrainedObject\n");
438 painCave.severity = OPENMD_ERROR;
439 painCave.isFatal = 1;
440 simError();
441 }
442
443 // phew. At this point, we should have the pointer to the
444 // correct Object restraint in the variable oRest.
445
446 oRest->setScaleFactor(scalingFactor);
447
448 Vector3d pos = (*ro)->getPos();
449
450 if ( (*ro)->isDirectional() ) {
451
452 // directional objects may have orientational restraints as well
453 // as positional, so get the rotation matrix first:
454
455 RotMat3x3d A = (*ro)->getA();
456 oRest->calcForce(pos, A);
457 (*ro)->addFrc(oRest->getRestraintForce());
458 (*ro)->addTrq(oRest->getRestraintTorque());
459 } else {
460
461 // plain vanilla positional restraints:
462
463 oRest->calcForce(pos);
464 (*ro)->addFrc(oRest->getRestraintForce());
465 }
466
467 unscaledPotential_ += oRest->getUnscaledPotential();
468
469 restInfo = oRest->getRestraintInfo();
470
471 // only collect data on restraints that we're going to print:
472 if (oRest->getPrintRestraint())
473 restInfo_.push_back(restInfo);
474 }
475
476 return unscaledPotential_ * scalingFactor;
477 }
478 }

Properties

Name Value
svn:executable *