ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/restraints/RestraintForceManager.cpp
Revision: 1665
Committed: Tue Nov 22 20:38:56 2011 UTC (13 years, 5 months ago) by gezelter
File size: 17690 byte(s)
Log Message:
updated copyright notices

File Contents

# User Rev Content
1 cli2 1360 /*
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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 cli2 1360 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 cli2 1360 * 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 gezelter 1390 *
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 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 cli2 1360 */
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 gezelter 1390 #include "utils/PhysicalConstants.hpp"
50 cli2 1360 #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 gezelter 1390 namespace OpenMD {
59 cli2 1360
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 cli2 1364 "Restraint warning: If you use restraints without setting\n"
83     "\tstatusTime, no restraint data will be written to the rest\n"
84 cli2 1360 "\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 cli2 1407 std::vector<int> stuntDoubleIndex;
94    
95 cli2 1360 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 gezelter 1366 if (molIndex < 0) {
116 cli2 1360 sprintf(painCave.errMsg,
117 gezelter 1366 "Restraint Error: A molecular restraint was specified\n"
118     "\twith a molIndex that was less than 0\n");
119 cli2 1360 painCave.isFatal = 1;
120     simError();
121     }
122 gezelter 1366 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 cli2 1360
132 gezelter 1366 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 cli2 1407
141 gezelter 1366 // 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 gezelter 1398
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 cli2 1407
163 gezelter 1398 #endif
164    
165 cli2 1360 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 cli2 1364 if (stamp[i]->havePrint()) {
194     rest->setPrintRestraint(stamp[i]->getPrint());
195     }
196 cli2 1360
197     restraints_.push_back(rest);
198     mol->addProperty(new RestraintData("Restraint", rest));
199     restrainedMols_.push_back(mol);
200 gezelter 1398 #ifdef IS_MPI
201     }
202     #endif
203 cli2 1360 } 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 cli2 1364
221 cli2 1360 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 cli2 1407 stuntDoubleIndex.push_back(sd->getGlobalIntegrableObjectIndex());
239    
240 cli2 1360 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 cli2 1364 }
263     if (stamp[i]->havePrint()) {
264     rest->setPrintRestraint(stamp[i]->getPrint());
265     }
266    
267 cli2 1360 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 cli2 1407 if (simParam->getUseRestraints()) {
280 cli2 1360 std::string refFile = simParam->getRestraint_file();
281 cli2 1407 RestReader* rr = new RestReader(info, refFile, stuntDoubleIndex);
282 cli2 1360 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 gezelter 1464 void RestraintForceManager::calcForces(){
310 cli2 1360
311 gezelter 1464 ForceManager::calcForces();
312 cli2 1360 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     currSnapshot_->statData[Stats::LONG_RANGE_POTENTIAL] += restPot;
324     currSnapshot_->statData[Stats::VHARM] = restPot;
325    
326     //write out forces and current positions of restrained molecules
327     if (currSnapshot_->getTime() >= currRestTime_){
328     restOut->writeRest(restInfo_);
329     currRestTime_ += restTime_;
330     }
331     }
332    
333     RealType RestraintForceManager::doRestraints(RealType scalingFactor){
334     std::vector<Molecule*>::const_iterator rm;
335     GenericData* data;
336     Molecule::IntegrableObjectIterator ioi;
337     MolecularRestraint* mRest;
338     StuntDouble* sd;
339     RealType pTot;
340    
341     std::vector<StuntDouble*>::const_iterator ro;
342     ObjectRestraint* oRest;
343    
344     std::map<int, Restraint::RealPair> restInfo;
345    
346     unscaledPotential_ = 0.0;
347    
348     restInfo_.clear();
349    
350     for(rm=restrainedMols_.begin(); rm != restrainedMols_.end(); ++rm){
351    
352     // make sure this molecule (*rm) has a generic data for restraints:
353     data = (*rm)->getPropertyByName("Restraint");
354     if (data != NULL) {
355     // make sure we can reinterpret the generic data as restraint data:
356     RestraintData* restData= dynamic_cast<RestraintData*>(data);
357     if (restData != NULL) {
358     // make sure we can reinterpet the restraint data as a pointer to
359     // an MolecularRestraint:
360     mRest = dynamic_cast<MolecularRestraint*>(restData->getData());
361     if (mRest == NULL) {
362     sprintf( painCave.errMsg,
363     "Can not cast RestraintData to MolecularRestraint\n");
364 gezelter 1390 painCave.severity = OPENMD_ERROR;
365 cli2 1360 painCave.isFatal = 1;
366     simError();
367     }
368     } else {
369     sprintf( painCave.errMsg,
370     "Can not cast GenericData to RestraintData\n");
371 gezelter 1390 painCave.severity = OPENMD_ERROR;
372 cli2 1360 painCave.isFatal = 1;
373     simError();
374     }
375     } else {
376     sprintf( painCave.errMsg, "Can not find Restraint for RestrainedObject\n");
377 gezelter 1390 painCave.severity = OPENMD_ERROR;
378 cli2 1360 painCave.isFatal = 1;
379     simError();
380     }
381    
382     // phew. At this point, we should have the pointer to the
383     // correct MolecularRestraint in the variable mRest.
384    
385     Vector3d molCom = (*rm)->getCom();
386    
387     std::vector<Vector3d> struc;
388     std::vector<Vector3d> forces;
389    
390     for(sd = (*rm)->beginIntegrableObject(ioi); sd != NULL;
391 cli2 1407 sd = (*rm)->nextIntegrableObject(ioi)) {
392 cli2 1360 struc.push_back(sd->getPos());
393     }
394    
395     mRest->setScaleFactor(scalingFactor);
396     mRest->calcForce(struc, molCom);
397     forces = mRest->getRestraintForces();
398     int index = 0;
399    
400     for(sd = (*rm)->beginIntegrableObject(ioi); sd != NULL;
401     sd = (*rm)->nextIntegrableObject(ioi)) {
402     sd->addFrc(forces[index]);
403     struc.push_back(sd->getPos());
404     index++;
405     }
406    
407     unscaledPotential_ += mRest->getUnscaledPotential();
408    
409     restInfo = mRest->getRestraintInfo();
410 cli2 1364
411     // only collect data on restraints that we're going to print:
412     if (mRest->getPrintRestraint())
413     restInfo_.push_back(restInfo);
414 cli2 1360 }
415    
416     for(ro=restrainedObjs_.begin(); ro != restrainedObjs_.end(); ++ro){
417     // make sure this object (*ro) has a generic data for restraints:
418     data = (*ro)->getPropertyByName("Restraint");
419     if (data != NULL) {
420     // make sure we can reinterpret the generic data as restraint data:
421     RestraintData* restData= dynamic_cast<RestraintData*>(data);
422     if (restData != NULL) {
423     // make sure we can reinterpet the restraint data as a pointer to
424     // an ObjectRestraint:
425     oRest = dynamic_cast<ObjectRestraint*>(restData->getData());
426     if (oRest == NULL) {
427     sprintf( painCave.errMsg,
428     "Can not cast RestraintData to ObjectRestraint\n");
429 gezelter 1390 painCave.severity = OPENMD_ERROR;
430 cli2 1360 painCave.isFatal = 1;
431     simError();
432     }
433     } else {
434     sprintf( painCave.errMsg,
435     "Can not cast GenericData to RestraintData\n");
436 gezelter 1390 painCave.severity = OPENMD_ERROR;
437 cli2 1360 painCave.isFatal = 1;
438     simError();
439     }
440     } else {
441     sprintf( painCave.errMsg, "Can not find Restraint for RestrainedObject\n");
442 gezelter 1390 painCave.severity = OPENMD_ERROR;
443 cli2 1360 painCave.isFatal = 1;
444     simError();
445     }
446    
447     // phew. At this point, we should have the pointer to the
448     // correct Object restraint in the variable oRest.
449    
450     oRest->setScaleFactor(scalingFactor);
451    
452     Vector3d pos = (*ro)->getPos();
453    
454     if ( (*ro)->isDirectional() ) {
455    
456     // directional objects may have orientational restraints as well
457     // as positional, so get the rotation matrix first:
458    
459     RotMat3x3d A = (*ro)->getA();
460     oRest->calcForce(pos, A);
461     (*ro)->addFrc(oRest->getRestraintForce());
462     (*ro)->addTrq(oRest->getRestraintTorque());
463     } else {
464    
465     // plain vanilla positional restraints:
466    
467     oRest->calcForce(pos);
468     (*ro)->addFrc(oRest->getRestraintForce());
469     }
470    
471     unscaledPotential_ += oRest->getUnscaledPotential();
472    
473     restInfo = oRest->getRestraintInfo();
474 cli2 1364
475     // only collect data on restraints that we're going to print:
476     if (oRest->getPrintRestraint())
477     restInfo_.push_back(restInfo);
478 cli2 1360 }
479    
480     return unscaledPotential_ * scalingFactor;
481     }
482     }

Properties

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