ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/DumpWriter.cpp
Revision: 1217
Committed: Wed Jan 23 21:23:02 2008 UTC (17 years, 3 months ago) by xsun
File size: 11810 byte(s)
Log Message:
Nothing important

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
244 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
245 integrableObject = mol->nextIntegrableObject(ii)) {
246 os << prepareDumpLine(integrableObject);
247
248 }
249 }
250 os << " </StuntDoubles>\n";
251
252 os << " </Snapshot>\n";
253
254 os.flush();
255 #else
256 //every node prepares the dump lines for integrable objects belong to itself
257 std::string buffer;
258 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
259
260
261 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
262 integrableObject = mol->nextIntegrableObject(ii)) {
263 buffer += prepareDumpLine(integrableObject);
264 }
265 }
266
267 const int masterNode = 0;
268
269 if (worldRank == masterNode) {
270 os << " <Snapshot>\n";
271 writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot());
272 os << " <StuntDoubles>\n";
273
274 os << buffer;
275
276 int nProc;
277 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
278 for (int i = 1; i < nProc; ++i) {
279
280 // receive the length of the string buffer that was
281 // prepared by processor i
282
283 int recvLength;
284 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &istatus);
285 char* recvBuffer = new char[recvLength];
286 if (recvBuffer == NULL) {
287 } else {
288 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD, &istatus);
289 os << recvBuffer;
290 delete recvBuffer;
291 }
292 }
293 os << " </StuntDoubles>\n";
294
295 os << " </Snapshot>\n";
296 os.flush();
297 } else {
298 int sendBufferLength = buffer.size() + 1;
299 MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
300 MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
301 }
302
303 #endif // is_mpi
304
305 }
306
307 std::string DumpWriter::prepareDumpLine(StuntDouble* integrableObject) {
308
309 int index = integrableObject->getGlobalIntegrableObjectIndex();
310 std::string type("pv");
311 std::string line;
312 char tempBuffer[4096];
313
314 Vector3d pos;
315 Vector3d vel;
316 pos = integrableObject->getPos();
317 vel = integrableObject->getVel();
318 sprintf(tempBuffer, "%18.10g %18.10g %18.10g %13e %13e %13e",
319 pos[0], pos[1], pos[2],
320 vel[0], vel[1], vel[2]);
321 line += tempBuffer;
322
323 if (integrableObject->isDirectional()) {
324 type += "qj";
325 Quat4d q;
326 Vector3d ji;
327 q = integrableObject->getQ();
328 ji = integrableObject->getJ();
329 sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e %13e",
330 q[0], q[1], q[2], q[3],
331 ji[0], ji[1], ji[2]);
332 line += tempBuffer;
333 }
334
335 if (needForceVector_) {
336 type += "ft";
337 Vector3d frc;
338 Vector3d trq;
339 frc = integrableObject->getFrc();
340 trq = integrableObject->getTrq();
341
342 sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e",
343 frc[0], frc[1], frc[2],
344 trq[0], trq[1], trq[2]);
345 line += tempBuffer;
346 }
347
348 sprintf(tempBuffer, "%10d %7s %s\n", index, type.c_str(), line.c_str());
349 return std::string(tempBuffer);
350 }
351
352 void DumpWriter::writeDump() {
353 writeFrame(*dumpFile_);
354 }
355
356 void DumpWriter::writeEor() {
357 std::ostream* eorStream;
358
359 #ifdef IS_MPI
360 if (worldRank == 0) {
361 #endif // is_mpi
362
363 eorStream = createOStream(eorFilename_);
364
365 #ifdef IS_MPI
366 }
367 #endif // is_mpi
368
369 writeFrame(*eorStream);
370
371 #ifdef IS_MPI
372 if (worldRank == 0) {
373 #endif // is_mpi
374 writeClosing(*eorStream);
375 delete eorStream;
376 #ifdef IS_MPI
377 }
378 #endif // is_mpi
379
380 }
381
382
383 void DumpWriter::writeDumpAndEor() {
384 std::vector<std::streambuf*> buffers;
385 std::ostream* eorStream;
386 #ifdef IS_MPI
387 if (worldRank == 0) {
388 #endif // is_mpi
389
390 buffers.push_back(dumpFile_->rdbuf());
391
392 eorStream = createOStream(eorFilename_);
393
394 buffers.push_back(eorStream->rdbuf());
395
396 #ifdef IS_MPI
397 }
398 #endif // is_mpi
399
400 TeeBuf tbuf(buffers.begin(), buffers.end());
401 std::ostream os(&tbuf);
402
403 writeFrame(os);
404
405 #ifdef IS_MPI
406 if (worldRank == 0) {
407 #endif // is_mpi
408 writeClosing(*eorStream);
409 delete eorStream;
410 #ifdef IS_MPI
411 }
412 #endif // is_mpi
413
414 }
415
416 std::ostream* DumpWriter::createOStream(const std::string& filename) {
417
418 std::ostream* newOStream;
419 #ifdef HAVE_LIBZ
420 if (needCompression_) {
421 newOStream = new ogzstream(filename.c_str());
422 } else {
423 newOStream = new std::ofstream(filename.c_str());
424 }
425 #else
426 newOStream = new std::ofstream(filename.c_str());
427 #endif
428 //write out MetaData first
429 (*newOStream) << "<OOPSE version=4>" << std::endl;
430 (*newOStream) << " <MetaData>" << std::endl;
431 (*newOStream) << info_->getRawMetaData();
432 (*newOStream) << " </MetaData>" << std::endl;
433 return newOStream;
434 }
435
436 void DumpWriter::writeClosing(std::ostream& os) {
437
438 os << "</OOPSE>\n";
439 os.flush();
440 }
441
442 }//end namespace oopse