ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/RestWriter.cpp
Revision: 2021
Committed: Tue Sep 23 15:25:08 2014 UTC (10 years, 7 months ago) by gezelter
File size: 9364 byte(s)
Log Message:
Took out a spurious debugging line

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date