ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/io/DumpWriter.cpp
Revision: 2983
Committed: Wed Aug 30 20:33:44 2006 UTC (17 years, 10 months ago) by gezelter
File size: 11789 byte(s)
Log Message:
fixed some bugs in DumpWriter, eliminated some extra printing of
debugging information

File Contents

# Content
1 /*
2 * Copyright (c) 2005 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. Acknowledgement of the program authors must be made in any
10 * publication of scientific results based in part on use of the
11 * program. An acceptable form of acknowledgement is citation of
12 * the article in which the program was described (Matthew
13 * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 * Parallel Simulation Engine for Molecular Dynamics,"
16 * J. Comput. Chem. 26, pp. 252-271 (2005))
17 *
18 * 2. Redistributions of source code must retain the above copyright
19 * notice, this list of conditions and the following disclaimer.
20 *
21 * 3. Redistributions in binary form must reproduce the above copyright
22 * notice, this list of conditions and the following disclaimer in the
23 * documentation and/or other materials provided with the
24 * distribution.
25 *
26 * This software is provided "AS IS," without a warranty of any
27 * kind. All express or implied conditions, representations and
28 * warranties, including any implied warranty of merchantability,
29 * fitness for a particular purpose or non-infringement, are hereby
30 * excluded. The University of Notre Dame and its licensors shall not
31 * be liable for any damages suffered by licensee as a result of
32 * using, modifying or distributing the software or its
33 * derivatives. In no event will the University of Notre Dame or its
34 * licensors be liable for any lost revenue, profit or data, or for
35 * direct, indirect, special, consequential, incidental or punitive
36 * damages, however caused and regardless of the theory of liability,
37 * arising out of the use of or inability to use software, even if the
38 * University of Notre Dame has been advised of the possibility of
39 * such damages.
40 */
41
42 #include "io/DumpWriter.hpp"
43 #include "primitives/Molecule.hpp"
44 #include "utils/simError.h"
45 #include "io/basic_teebuf.hpp"
46 #include "io/gzstream.hpp"
47 #include "io/Globals.hpp"
48
49 #ifdef IS_MPI
50 #include <mpi.h>
51 #endif //is_mpi
52
53 namespace oopse {
54
55 DumpWriter::DumpWriter(SimInfo* info)
56 : info_(info), filename_(info->getDumpFileName()), eorFilename_(info->getFinalConfigFileName()){
57
58 Globals* simParams = info->getSimParams();
59 needCompression_ = simParams->getCompressDumpFile();
60 needForceVector_ = simParams->getOutputForceVector();
61 createDumpFile_ = true;
62 #ifdef HAVE_LIBZ
63 if (needCompression_) {
64 filename_ += ".gz";
65 eorFilename_ += ".gz";
66 }
67 #endif
68
69 #ifdef IS_MPI
70
71 if (worldRank == 0) {
72 #endif // is_mpi
73
74 dumpFile_ = createOStream(filename_);
75
76 if (!dumpFile_) {
77 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
78 filename_.c_str());
79 painCave.isFatal = 1;
80 simError();
81 }
82
83 #ifdef IS_MPI
84
85 }
86
87 #endif // is_mpi
88
89 }
90
91
92 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename)
93 : info_(info), filename_(filename){
94
95 Globals* simParams = info->getSimParams();
96 eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
97
98 needCompression_ = simParams->getCompressDumpFile();
99 needForceVector_ = simParams->getOutputForceVector();
100 createDumpFile_ = true;
101 #ifdef HAVE_LIBZ
102 if (needCompression_) {
103 filename_ += ".gz";
104 eorFilename_ += ".gz";
105 }
106 #endif
107
108 #ifdef IS_MPI
109
110 if (worldRank == 0) {
111 #endif // is_mpi
112
113
114 dumpFile_ = createOStream(filename_);
115
116 if (!dumpFile_) {
117 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
118 filename_.c_str());
119 painCave.isFatal = 1;
120 simError();
121 }
122
123 #ifdef IS_MPI
124
125 }
126
127 #endif // is_mpi
128
129 }
130
131 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename, bool writeDumpFile)
132 : info_(info), filename_(filename){
133
134 Globals* simParams = info->getSimParams();
135 eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
136
137 needCompression_ = simParams->getCompressDumpFile();
138 needForceVector_ = simParams->getOutputForceVector();
139
140 #ifdef HAVE_LIBZ
141 if (needCompression_) {
142 filename_ += ".gz";
143 eorFilename_ += ".gz";
144 }
145 #endif
146
147 #ifdef IS_MPI
148
149 if (worldRank == 0) {
150 #endif // is_mpi
151
152 createDumpFile_ = writeDumpFile;
153 if (createDumpFile_) {
154 dumpFile_ = createOStream(filename_);
155
156 if (!dumpFile_) {
157 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
158 filename_.c_str());
159 painCave.isFatal = 1;
160 simError();
161 }
162 }
163 #ifdef IS_MPI
164
165 }
166
167
168 #endif // is_mpi
169
170 }
171
172 DumpWriter::~DumpWriter() {
173
174 #ifdef IS_MPI
175
176 if (worldRank == 0) {
177 #endif // is_mpi
178 if (createDumpFile_){
179 writeClosing(*dumpFile_);
180 delete dumpFile_;
181 }
182 #ifdef IS_MPI
183
184 }
185
186 #endif // is_mpi
187
188 }
189
190 void DumpWriter::writeFrameProperties(std::ostream& os, Snapshot* s) {
191
192 char buffer[1024];
193
194 os << " <FrameData>\n";
195
196 RealType currentTime = s->getTime();
197 sprintf(buffer, " Time: %.10g\n", currentTime);
198 os << buffer;
199
200 Mat3x3d hmat;
201 hmat = s->getHmat();
202 sprintf(buffer, " Hmat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
203 hmat(0, 0), hmat(1, 0), hmat(2, 0),
204 hmat(0, 1), hmat(1, 1), hmat(2, 1),
205 hmat(0, 2), hmat(1, 2), hmat(2, 2));
206 os << buffer;
207
208 RealType chi = s->getChi();
209 RealType integralOfChiDt = s->getIntegralOfChiDt();
210 sprintf(buffer, " Thermostat: %.10g , %.10g\n", chi, integralOfChiDt);
211 os << buffer;
212
213 Mat3x3d eta;
214 eta = s->getEta();
215 sprintf(buffer, " Barostat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
216 eta(0, 0), eta(1, 0), eta(2, 0),
217 eta(0, 1), eta(1, 1), eta(2, 1),
218 eta(0, 2), eta(1, 2), eta(2, 2));
219 os << buffer;
220
221 os << " </FrameData>\n";
222 }
223
224 void DumpWriter::writeFrame(std::ostream& os) {
225
226 #ifdef IS_MPI
227 MPI_Status istatus;
228 #endif
229
230 Molecule* mol;
231 StuntDouble* integrableObject;
232 SimInfo::MoleculeIterator mi;
233 Molecule::IntegrableObjectIterator ii;
234
235 #ifndef IS_MPI
236 os << " <Snapshot>\n";
237
238 writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot());
239
240 os << " <StuntDoubles>\n";
241 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
242
243 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
244 integrableObject = mol->nextIntegrableObject(ii)) {
245 os << prepareDumpLine(integrableObject);
246
247 }
248 }
249 os << " </StuntDoubles>\n";
250
251 os << " </Snapshot>\n";
252
253 os.flush();
254 #else
255 //every node prepares the dump lines for integrable objects belong to itself
256 std::string buffer;
257 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
258
259 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
260 integrableObject = mol->nextIntegrableObject(ii)) {
261 buffer += prepareDumpLine(integrableObject);
262 }
263 }
264
265 const int masterNode = 0;
266
267 if (worldRank == masterNode) {
268 os << " <Snapshot>\n";
269 writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot());
270 os << " <StuntDoubles>\n";
271
272 os << buffer;
273
274 int nProc;
275 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
276 for (int i = 1; i < nProc; ++i) {
277
278 // receive the length of the string buffer that was
279 // prepared by processor i
280
281 int recvLength;
282 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &istatus);
283 char* recvBuffer = new char[recvLength];
284 if (recvBuffer == NULL) {
285 } else {
286 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD, &istatus);
287 os << recvBuffer;
288 delete recvBuffer;
289 }
290 }
291 os << " </StuntDoubles>\n";
292
293 os << " </Snapshot>\n";
294 os.flush();
295 } else {
296 int sendBufferLength = buffer.size() + 1;
297 MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
298 MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
299 }
300
301 #endif // is_mpi
302
303 }
304
305 std::string DumpWriter::prepareDumpLine(StuntDouble* integrableObject) {
306
307 int index = integrableObject->getGlobalIntegrableObjectIndex();
308 std::string type("pv");
309 std::string line;
310 char tempBuffer[4096];
311
312 Vector3d pos;
313 Vector3d vel;
314 pos = integrableObject->getPos();
315 vel = integrableObject->getVel();
316 sprintf(tempBuffer, "%18.10g %18.10g %18.10g %13e %13e %13e",
317 pos[0], pos[1], pos[2],
318 vel[0], vel[1], vel[2]);
319 line += tempBuffer;
320
321 if (integrableObject->isDirectional()) {
322 type += "qj";
323 Quat4d q;
324 Vector3d ji;
325 q = integrableObject->getQ();
326 ji = integrableObject->getJ();
327 sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e %13e",
328 q[0], q[1], q[2], q[3],
329 ji[0], ji[1], ji[2]);
330 line += tempBuffer;
331 }
332
333 if (needForceVector_) {
334 type += "ft";
335 Vector3d frc;
336 Vector3d trq;
337 frc = integrableObject->getFrc();
338 trq = integrableObject->getTrq();
339
340 sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e",
341 frc[0], frc[1], frc[2],
342 trq[0], trq[1], trq[2]);
343 line += tempBuffer;
344 }
345
346 sprintf(tempBuffer, "%10d %7s %s\n", index, type.c_str(), line.c_str());
347 return std::string(tempBuffer);
348 }
349
350 void DumpWriter::writeDump() {
351 writeFrame(*dumpFile_);
352 }
353
354 void DumpWriter::writeEor() {
355 std::ostream* eorStream;
356
357 #ifdef IS_MPI
358 if (worldRank == 0) {
359 #endif // is_mpi
360
361 eorStream = createOStream(eorFilename_);
362
363 #ifdef IS_MPI
364 }
365 #endif // is_mpi
366
367 writeFrame(*eorStream);
368
369 #ifdef IS_MPI
370 if (worldRank == 0) {
371 #endif // is_mpi
372 writeClosing(*eorStream);
373 delete eorStream;
374 #ifdef IS_MPI
375 }
376 #endif // is_mpi
377
378 }
379
380
381 void DumpWriter::writeDumpAndEor() {
382 std::vector<std::streambuf*> buffers;
383 std::ostream* eorStream;
384 #ifdef IS_MPI
385 if (worldRank == 0) {
386 #endif // is_mpi
387
388 buffers.push_back(dumpFile_->rdbuf());
389
390 eorStream = createOStream(eorFilename_);
391
392 buffers.push_back(eorStream->rdbuf());
393
394 #ifdef IS_MPI
395 }
396 #endif // is_mpi
397
398 TeeBuf tbuf(buffers.begin(), buffers.end());
399 std::ostream os(&tbuf);
400
401 writeFrame(os);
402
403 #ifdef IS_MPI
404 if (worldRank == 0) {
405 #endif // is_mpi
406 writeClosing(*eorStream);
407 delete eorStream;
408 #ifdef IS_MPI
409 }
410 #endif // is_mpi
411
412 }
413
414 std::ostream* DumpWriter::createOStream(const std::string& filename) {
415
416 std::ostream* newOStream;
417 #ifdef HAVE_LIBZ
418 if (needCompression_) {
419 newOStream = new ogzstream(filename.c_str());
420 } else {
421 newOStream = new std::ofstream(filename.c_str());
422 }
423 #else
424 newOStream = new std::ofstream(filename.c_str());
425 #endif
426 //write out MetaData first
427 (*newOStream) << "<OOPSE version=4>" << std::endl;
428 (*newOStream) << " <MetaData>" << std::endl;
429 (*newOStream) << info_->getRawMetaData();
430 (*newOStream) << " </MetaData>" << std::endl;
431 return newOStream;
432 }
433
434 void DumpWriter::writeClosing(std::ostream& os) {
435
436 os << "</OOPSE>\n";
437 os.flush();
438 }
439
440 }//end namespace oopse