ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/restraints/RestraintForceManager.cpp
Revision: 1767
Committed: Fri Jul 6 22:01:58 2012 UTC (12 years, 9 months ago) by gezelter
File size: 17747 byte(s)
Log Message:
Various fixes required to compile OpenMD with the MS Visual C++ compiler

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

Properties

Name Value
svn:executable *
svn:keywords Author Id Revision Date