ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/restraints/RestraintForceManager.cpp
Revision: 1760
Committed: Thu Jun 21 19:26:46 2012 UTC (12 years, 10 months ago) by gezelter
File size: 17745 byte(s)
Log Message:
Some bugfixes (CholeskyDecomposition), more work on fluctuating charges,
migrating stats stuff into frameData

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

Properties

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