OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
RestWriter.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/RestWriter.hpp"
46
47#include <iostream>
48#include <sstream>
49#include <string>
50
51#ifdef IS_MPI
52#include <mpi.h>
53#endif
54
56#include "utils/simError.h"
57
58namespace OpenMD {
59 RestWriter::RestWriter(SimInfo* info, const std::string& filename,
60 std::vector<Restraint*> restraints) :
61 info_(info) {
62 std::vector<Restraint*>::const_iterator resti;
63
64 createRestFile_ = false;
65
66#ifdef IS_MPI
67 MPI_Status* istatus = NULL;
68#endif
69
70 int printAny = 0;
71 for (resti = restraints.begin(); resti != restraints.end(); ++resti) {
72 if ((*resti)->getPrintRestraint()) { printAny = 1; }
73 }
74
75#ifdef IS_MPI
76 MPI_Allreduce(MPI_IN_PLACE, &printAny, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
77#endif
78
79 if (printAny) createRestFile_ = true;
80
81#ifdef IS_MPI
82 if (worldRank == 0) {
83#endif
84
85 if (createRestFile_) {
86 output_ = new std::ofstream(filename.c_str());
87
88 if (!output_) {
89 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
90 "Could not open %s for restraint output.\n",
91 filename.c_str());
92 painCave.isFatal = 1;
93 simError();
94 }
95 }
96
97#ifdef IS_MPI
98 }
99#endif // is_mpi
100
101#ifndef IS_MPI
102
103 if (createRestFile_) (*output_) << "#time\t";
104
105 for (resti = restraints.begin(); resti != restraints.end(); ++resti) {
106 if ((*resti)->getPrintRestraint()) {
107 std::string myName = (*resti)->getRestraintName();
108 int myType = (*resti)->getRestraintType();
109
110 (*output_) << myName << ":";
111
112 if (myType & Restraint::rtDisplacement)
113 (*output_) << "\tPosition(angstroms)\tEnergy(kcal/mol)";
114
115 if (myType & Restraint::rtAbsoluteZ)
116 (*output_) << "\tPosition(angstroms)\tEnergy(kcal/mol)";
117
118 if (myType & Restraint::rtTwist)
119 (*output_) << "\tTwistAngle(radians)\tEnergy(kcal/mol)";
120
121 if (myType & Restraint::rtSwingX)
122 (*output_) << "\tSwingXAngle(radians)\tEnergy(kcal/mol)";
123
124 if (myType & Restraint::rtSwingY)
125 (*output_) << "\tSwingYAngle(radians)\tEnergy(kcal/mol)";
126 }
127 }
128
129 if (createRestFile_) (*output_) << "\n";
130 if (createRestFile_) (*output_).flush();
131
132#else
133
134 std::string buffer;
135
136 for (resti = restraints.begin(); resti != restraints.end(); ++resti) {
137 if ((*resti)->getPrintRestraint()) {
138 std::string myName = (*resti)->getRestraintName();
139 int myType = (*resti)->getRestraintType();
140
141 buffer += (myName + ":");
142
143 if (myType & Restraint::rtDisplacement)
144 buffer += "\tPosition(angstroms)\tEnergy(kcal/mol)";
145
146 if (myType & Restraint::rtAbsoluteZ)
147 buffer += "\tPosition(angstroms)\tEnergy(kcal/mol)";
148
149 if (myType & Restraint::rtTwist)
150 buffer += "\tTwistAngle(radians)\tEnergy(kcal/mol)";
151
152 if (myType & Restraint::rtSwingX)
153 buffer += "\tSwingXAngle(radians)\tEnergy(kcal/mol)";
154
155 if (myType & Restraint::rtSwingY)
156 buffer += "\tSwingYAngle(radians)\tEnergy(kcal/mol)";
157
158 buffer += "\n";
159 }
160 }
161
162 const int primaryNode = 0;
163
164 if (worldRank == primaryNode) {
165 if (createRestFile_) (*output_) << "#time\t";
166 if (createRestFile_) (*output_) << buffer;
167
168 int nProc;
169 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
170
171 for (int i = 1; i < nProc; ++i) {
172 // receive the length of the string buffer that was
173 // prepared by processor i
174
175 int recvLength;
176 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
177 char* recvBuffer = new char[recvLength];
178 if (recvBuffer == NULL) {
179 } else {
180 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
181 istatus);
182 if (createRestFile_) (*output_) << recvBuffer;
183 delete[] recvBuffer;
184 }
185 }
186 if (createRestFile_) (*output_).flush();
187 } else {
188 int sendBufferLength = buffer.size() + 1;
189 MPI_Send(&sendBufferLength, 1, MPI_INT, primaryNode, 0, MPI_COMM_WORLD);
190 MPI_Send((void*)buffer.c_str(), sendBufferLength, MPI_CHAR, primaryNode,
191 0, MPI_COMM_WORLD);
192 }
193
194#endif // is_mpi
195 }
196
197 void RestWriter::writeRest(
198 std::vector<std::map<int, Restraint::RealPair>> restInfo) {
199#ifdef IS_MPI
200 MPI_Status* istatus = NULL;
201#endif
202
203#ifndef IS_MPI
204 if (createRestFile_)
205 (*output_)
206 << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
207
208 // output some information about the molecules
209 std::vector<std::map<int, Restraint::RealPair>>::const_iterator i;
210 std::map<int, Restraint::RealPair>::const_iterator j;
211
212 if (createRestFile_) {
213 for (i = restInfo.begin(); i != restInfo.end(); ++i) {
214 for (j = (*i).begin(); j != (*i).end(); ++j) {
215 (*output_) << "\t" << (j->second).first << "\t" << (j->second).second;
216 }
217 (*output_) << std::endl;
218 }
219 (*output_).flush();
220 }
221#else
222 std::string buffer, first, second;
223 std::stringstream ss;
224
225 std::vector<std::map<int, Restraint::RealPair>>::const_iterator i;
226 std::map<int, Restraint::RealPair>::const_iterator j;
227
228 if (createRestFile_) {
229 for (i = restInfo.begin(); i != restInfo.end(); ++i) {
230 for (j = (*i).begin(); j != (*i).end(); ++j) {
231 ss.clear();
232 ss << (j->second).first;
233 ss >> first;
234 ss.clear();
235 ss << (j->second).second;
236 ss >> second;
237 buffer += ("\t" + first + "\t" + second);
238 }
239 buffer += "\n";
240 }
241 }
242
243 const int primaryNode = 0;
244
245 if (createRestFile_) {
246 if (worldRank == primaryNode) {
247 (*output_)
248 << info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
249 (*output_) << buffer;
250
251 int nProc;
252 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
253 for (int i = 1; i < nProc; ++i) {
254 // receive the length of the string buffer that was
255 // prepared by processor i
256
257 int recvLength;
258 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, istatus);
259 char* recvBuffer = new char[recvLength];
260 if (recvBuffer == NULL) {
261 } else {
262 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD,
263 istatus);
264 if (createRestFile_) (*output_) << recvBuffer;
265
266 delete[] recvBuffer;
267 }
268 }
269 (*output_).flush();
270 } else {
271 int sendBufferLength = buffer.size() + 1;
272 MPI_Send(&sendBufferLength, 1, MPI_INT, primaryNode, 0, MPI_COMM_WORLD);
273 MPI_Send((void*)buffer.c_str(), sendBufferLength, MPI_CHAR, primaryNode,
274 0, MPI_COMM_WORLD);
275 }
276 }
277#endif // is_mpi
278 }
279
280 RestWriter::~RestWriter() {
281#ifdef IS_MPI
282
283 if (worldRank == 0) {
284#endif // is_mpi
285 if (createRestFile_) {
286 writeClosing(*output_);
287 delete output_;
288 }
289#ifdef IS_MPI
290 }
291#endif // is_mpi
292 }
293
294 void RestWriter::writeClosing(std::ostream& os) { os.flush(); }
295
296} // namespace OpenMD
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.