OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
ConstraintWriter.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include "io/ConstraintWriter.hpp"
46
47#include <algorithm>
48#include <iostream>
49#include <vector>
50
51#ifdef IS_MPI
52#include <mpi.h>
53#endif
54
55#include "utils/simError.h"
56
57namespace OpenMD {
58 ConstraintWriter::ConstraintWriter(SimInfo* info,
59 const std::string& filename) :
60 info_(info) {
61 // use a primary - secondary model, only the primary node writes
62 // to disk
63#ifdef IS_MPI
64 if (worldRank == 0) {
65#endif
66 output_.open(filename.c_str());
67
68 if (!output_) {
69 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
70 "Could not open %s for Constraint output\n", filename.c_str());
71 painCave.isFatal = 1;
72 simError();
73 }
74
75 output_ << "#time(fs)\t"
76 << "Index of atom 1\t"
77 << "Index of atom 2\tconstraint force" << std::endl;
78
79#ifdef IS_MPI
80 }
81#endif
82 }
83
84 ConstraintWriter::~ConstraintWriter() {
85#ifdef IS_MPI
86 if (worldRank == 0) {
87#endif
88 output_.close();
89#ifdef IS_MPI
90 }
91#endif
92 }
93
94 void ConstraintWriter::writeConstraintForces(
95 const std::list<ConstraintPair*>& constraints) {
96#ifndef IS_MPI
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;
104 }
105 }
106#else
107
108 const int primaryNode = 0;
109 int nproc;
110 int myNode;
111 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
112 MPI_Comm_rank(MPI_COMM_WORLD, &myNode);
113
114 std::vector<int> nConstraints(nproc, 0);
115 nConstraints[myNode] = constraints.size();
116
117 // do MPI_ALLREDUCE to exchange the total number of constraints:
118 MPI_Allreduce(MPI_IN_PLACE, &nConstraints[0], nproc, MPI_INT, MPI_SUM,
119 MPI_COMM_WORLD);
120
121 MPI_Status ierr;
122 int atom1, atom2, doPrint;
123 RealType force;
124
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);
137 }
138
139 } else {
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);
145
146 tmpData.atom1 = atom1;
147 tmpData.atom2 = atom2;
148 tmpData.constraintForce = force;
149 tmpData.printForce = (doPrint == 0) ? false : true;
150 constraintData.push_back(tmpData);
151 }
152 }
153 }
154
155 std::vector<ConstraintData>::iterator l;
156 for (l = constraintData.begin(); l != constraintData.end(); ++l) {
157 if (l->printForce) {
158 output_
159 << info_->getSnapshotManager()->getCurrentSnapshot()->getTime()
160 << "\t" << l->atom1 << "\t" << l->atom2 << "\t"
161 << l->constraintForce << std::endl;
162 }
163 }
164 } else {
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();
171
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,
175 MPI_COMM_WORLD);
176 MPI_Send(&printForce, 1, MPI_INT, primaryNode, 0, MPI_COMM_WORLD);
177 }
178 }
179#endif
180 }
181} // namespace OpenMD
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.