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