45#include "io/ConstraintWriter.hpp"
55#include "utils/simError.h"
58 ConstraintWriter::ConstraintWriter(SimInfo* info,
59 const std::string& filename) :
66 output_.open(filename.c_str());
69 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
70 "Could not open %s for Constraint output\n", filename.c_str());
75 output_ <<
"#time(fs)\t"
76 <<
"Index of atom 1\t"
77 <<
"Index of atom 2\tconstraint force" << std::endl;
84 ConstraintWriter::~ConstraintWriter() {
94 void ConstraintWriter::writeConstraintForces(
95 const std::list<ConstraintPair*>& constraints) {
97 std::list<ConstraintPair*>::const_iterator i;
98 for (i = constraints.begin(); i != constraints.end(); ++i) {
99 if ((*i)->getPrintForce()) {
100 output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime()
101 <<
"\t" << (*i)->getConsElem1()->getGlobalIndex() <<
"\t"
102 << (*i)->getConsElem2()->getGlobalIndex() <<
"\t"
103 << (*i)->getConstraintForce() << std::endl;
108 const int primaryNode = 0;
111 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
112 MPI_Comm_rank(MPI_COMM_WORLD, &myNode);
114 std::vector<int> nConstraints(nproc, 0);
115 nConstraints[myNode] = constraints.size();
118 MPI_Allreduce(MPI_IN_PLACE, &nConstraints[0], nproc, MPI_INT, MPI_SUM,
122 int atom1, atom2, doPrint;
125 if (myNode == primaryNode) {
126 std::vector<ConstraintData> constraintData;
127 ConstraintData tmpData;
128 for (
int i = 0; i < nproc; ++i) {
129 if (i == primaryNode) {
130 std::list<ConstraintPair*>::const_iterator j;
131 for (j = constraints.begin(); j != constraints.end(); ++j) {
132 tmpData.atom1 = (*j)->getConsElem1()->getGlobalIndex();
133 tmpData.atom2 = (*j)->getConsElem2()->getGlobalIndex();
134 tmpData.constraintForce = (*j)->getConstraintForce();
135 tmpData.printForce = (*j)->getPrintForce();
136 constraintData.push_back(tmpData);
140 for (
int k = 0; k < nConstraints[i]; ++k) {
141 MPI_Recv(&atom1, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &ierr);
142 MPI_Recv(&atom2, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &ierr);
143 MPI_Recv(&force, 1, MPI_REALTYPE, i, 0, MPI_COMM_WORLD, &ierr);
144 MPI_Recv(&doPrint, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &ierr);
146 tmpData.atom1 = atom1;
147 tmpData.atom2 = atom2;
148 tmpData.constraintForce = force;
149 tmpData.printForce = (doPrint == 0) ?
false : true;
150 constraintData.push_back(tmpData);
155 std::vector<ConstraintData>::iterator l;
156 for (l = constraintData.begin(); l != constraintData.end(); ++l) {
159 << info_->getSnapshotManager()->getCurrentSnapshot()->getTime()
160 <<
"\t" << l->atom1 <<
"\t" << l->atom2 <<
"\t"
161 << l->constraintForce << std::endl;
165 std::list<ConstraintPair*>::const_iterator j;
166 for (j = constraints.begin(); j != constraints.end(); ++j) {
167 int atom1 = (*j)->getConsElem1()->getGlobalIndex();
168 int atom2 = (*j)->getConsElem2()->getGlobalIndex();
169 RealType constraintForce = (*j)->getConstraintForce();
170 int printForce = (int)(*j)->getPrintForce();
172 MPI_Send(&atom1, 1, MPI_INT, primaryNode, 0, MPI_COMM_WORLD);
173 MPI_Send(&atom2, 1, MPI_INT, primaryNode, 0, MPI_COMM_WORLD);
174 MPI_Send(&constraintForce, 1, MPI_REALTYPE, primaryNode, 0,
176 MPI_Send(&printForce, 1, MPI_INT, primaryNode, 0, MPI_COMM_WORLD);
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.