ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/restraints/RestraintForceManager.cpp
Revision: 1390
Committed: Wed Nov 25 20:02:06 2009 UTC (15 years, 5 months ago) by gezelter
File size: 17310 byte(s)
Log Message:
Almost all of the changes necessary to create OpenMD out of our old
project (OOPSE-4)

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

Properties

Name Value
svn:executable *