ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpWriter.cpp
Revision: 1629
Committed: Wed Sep 14 21:15:17 2011 UTC (13 years, 7 months ago) by gezelter
File size: 16300 byte(s)
Log Message:
Merging changes from old branch into development branch

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 1390 * Copyright (c) 2009 The University of Notre Dame. All Rights Reserved.
3 gezelter 246 *
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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * 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 gezelter 1390 *
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, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40 gezelter 246 */
41    
42     #include "io/DumpWriter.hpp"
43     #include "primitives/Molecule.hpp"
44     #include "utils/simError.h"
45 tim 376 #include "io/basic_teebuf.hpp"
46 tim 615 #include "io/gzstream.hpp"
47     #include "io/Globals.hpp"
48    
49 gezelter 1437
50 gezelter 2 #ifdef IS_MPI
51     #include <mpi.h>
52 gezelter 246 #endif //is_mpi
53 gezelter 2
54 gezelter 1437 using namespace std;
55 gezelter 1390 namespace OpenMD {
56 gezelter 2
57 gezelter 507 DumpWriter::DumpWriter(SimInfo* info)
58     : info_(info), filename_(info->getDumpFileName()), eorFilename_(info->getFinalConfigFileName()){
59 tim 615
60     Globals* simParams = info->getSimParams();
61     needCompression_ = simParams->getCompressDumpFile();
62 chrisfen 726 needForceVector_ = simParams->getOutputForceVector();
63 gezelter 1629 needParticlePot_ = simParams->getOutputParticlePotential();
64 chuckv 791 createDumpFile_ = true;
65 tim 619 #ifdef HAVE_LIBZ
66 tim 615 if (needCompression_) {
67 tim 1024 filename_ += ".gz";
68     eorFilename_ += ".gz";
69 tim 615 }
70 tim 619 #endif
71 tim 615
72 tim 376 #ifdef IS_MPI
73    
74 tim 1024 if (worldRank == 0) {
75 tim 376 #endif // is_mpi
76 chuckv 791
77 tim 1024 dumpFile_ = createOStream(filename_);
78 tim 615
79 tim 1024 if (!dumpFile_) {
80     sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
81     filename_.c_str());
82     painCave.isFatal = 1;
83     simError();
84     }
85 tim 376
86     #ifdef IS_MPI
87    
88 tim 1024 }
89 tim 376
90     #endif // is_mpi
91    
92 tim 1024 }
93 tim 376
94    
95 gezelter 507 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename)
96     : info_(info), filename_(filename){
97 tim 615
98     Globals* simParams = info->getSimParams();
99     eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
100    
101     needCompression_ = simParams->getCompressDumpFile();
102 chrisfen 726 needForceVector_ = simParams->getOutputForceVector();
103 chuckv 791 createDumpFile_ = true;
104 tim 619 #ifdef HAVE_LIBZ
105 tim 615 if (needCompression_) {
106 tim 1024 filename_ += ".gz";
107     eorFilename_ += ".gz";
108 tim 615 }
109 tim 619 #endif
110 tim 615
111 gezelter 246 #ifdef IS_MPI
112 gezelter 2
113 tim 1024 if (worldRank == 0) {
114 gezelter 246 #endif // is_mpi
115 gezelter 2
116 chuckv 791
117 tim 1024 dumpFile_ = createOStream(filename_);
118 tim 615
119 tim 1024 if (!dumpFile_) {
120     sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
121     filename_.c_str());
122     painCave.isFatal = 1;
123     simError();
124     }
125 gezelter 2
126     #ifdef IS_MPI
127    
128 tim 1024 }
129 gezelter 2
130     #endif // is_mpi
131 gezelter 246
132 tim 1024 }
133 chuckv 791
134     DumpWriter::DumpWriter(SimInfo* info, const std::string& filename, bool writeDumpFile)
135 tim 1024 : info_(info), filename_(filename){
136 chuckv 791
137     Globals* simParams = info->getSimParams();
138     eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
139    
140     needCompression_ = simParams->getCompressDumpFile();
141     needForceVector_ = simParams->getOutputForceVector();
142 gezelter 1629 needParticlePot_ = simParams->getOutputParticlePotential();
143 chuckv 791
144     #ifdef HAVE_LIBZ
145     if (needCompression_) {
146     filename_ += ".gz";
147     eorFilename_ += ".gz";
148     }
149     #endif
150    
151     #ifdef IS_MPI
152    
153     if (worldRank == 0) {
154     #endif // is_mpi
155    
156     createDumpFile_ = writeDumpFile;
157     if (createDumpFile_) {
158     dumpFile_ = createOStream(filename_);
159    
160     if (!dumpFile_) {
161     sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
162     filename_.c_str());
163     painCave.isFatal = 1;
164     simError();
165     }
166     }
167     #ifdef IS_MPI
168    
169     }
170 tim 1024
171 chuckv 791
172     #endif // is_mpi
173    
174     }
175 gezelter 2
176 gezelter 507 DumpWriter::~DumpWriter() {
177 gezelter 2
178     #ifdef IS_MPI
179 gezelter 246
180     if (worldRank == 0) {
181 gezelter 2 #endif // is_mpi
182 chuckv 791 if (createDumpFile_){
183 tim 1024 writeClosing(*dumpFile_);
184 chuckv 791 delete dumpFile_;
185     }
186 gezelter 2 #ifdef IS_MPI
187 gezelter 246
188     }
189    
190 gezelter 2 #endif // is_mpi
191 gezelter 246
192 gezelter 507 }
193 gezelter 2
194 tim 1024 void DumpWriter::writeFrameProperties(std::ostream& os, Snapshot* s) {
195 gezelter 2
196 tim 1024 char buffer[1024];
197    
198     os << " <FrameData>\n";
199    
200     RealType currentTime = s->getTime();
201 gezelter 1437
202     if (isinf(currentTime) || isnan(currentTime)) {
203     sprintf( painCave.errMsg,
204     "DumpWriter detected a numerical error writing the time");
205     painCave.isFatal = 1;
206     simError();
207     }
208    
209 gezelter 1025 sprintf(buffer, " Time: %.10g\n", currentTime);
210 tim 1024 os << buffer;
211    
212 gezelter 246 Mat3x3d hmat;
213     hmat = s->getHmat();
214 gezelter 1437
215     for (unsigned int i = 0; i < 3; i++) {
216     for (unsigned int j = 0; j < 3; j++) {
217     if (isinf(hmat(i,j)) || isnan(hmat(i,j))) {
218     sprintf( painCave.errMsg,
219     "DumpWriter detected a numerical error writing the box");
220     painCave.isFatal = 1;
221     simError();
222     }
223     }
224     }
225    
226 tim 1024 sprintf(buffer, " Hmat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
227     hmat(0, 0), hmat(1, 0), hmat(2, 0),
228     hmat(0, 1), hmat(1, 1), hmat(2, 1),
229     hmat(0, 2), hmat(1, 2), hmat(2, 2));
230     os << buffer;
231 gezelter 2
232 tim 1024 RealType chi = s->getChi();
233     RealType integralOfChiDt = s->getIntegralOfChiDt();
234 gezelter 1437 if (isinf(chi) || isnan(chi) ||
235     isinf(integralOfChiDt) || isnan(integralOfChiDt)) {
236     sprintf( painCave.errMsg,
237     "DumpWriter detected a numerical error writing the thermostat");
238     painCave.isFatal = 1;
239     simError();
240     }
241 tim 1024 sprintf(buffer, " Thermostat: %.10g , %.10g\n", chi, integralOfChiDt);
242     os << buffer;
243 gezelter 246
244 tim 1024 Mat3x3d eta;
245     eta = s->getEta();
246 gezelter 1437
247     for (unsigned int i = 0; i < 3; i++) {
248     for (unsigned int j = 0; j < 3; j++) {
249     if (isinf(eta(i,j)) || isnan(eta(i,j))) {
250     sprintf( painCave.errMsg,
251     "DumpWriter detected a numerical error writing the barostat");
252     painCave.isFatal = 1;
253     simError();
254     }
255     }
256     }
257    
258 tim 1024 sprintf(buffer, " Barostat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
259     eta(0, 0), eta(1, 0), eta(2, 0),
260     eta(0, 1), eta(1, 1), eta(2, 1),
261     eta(0, 2), eta(1, 2), eta(2, 2));
262     os << buffer;
263 gezelter 246
264 tim 1024 os << " </FrameData>\n";
265 gezelter 507 }
266 gezelter 2
267 gezelter 507 void DumpWriter::writeFrame(std::ostream& os) {
268 gezelter 246
269 tim 1024 #ifdef IS_MPI
270     MPI_Status istatus;
271     #endif
272 gezelter 246
273     Molecule* mol;
274     StuntDouble* integrableObject;
275     SimInfo::MoleculeIterator mi;
276     Molecule::IntegrableObjectIterator ii;
277 gezelter 2
278 gezelter 246 #ifndef IS_MPI
279 tim 1024 os << " <Snapshot>\n";
280 gezelter 1025
281 tim 1024 writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot());
282 gezelter 2
283 tim 1024 os << " <StuntDoubles>\n";
284 gezelter 246 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
285 gezelter 2
286 xsun 1217
287     for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
288 tim 1024 integrableObject = mol->nextIntegrableObject(ii)) {
289 xsun 1217 os << prepareDumpLine(integrableObject);
290    
291 tim 1024 }
292     }
293     os << " </StuntDoubles>\n";
294 xsun 1217
295 tim 1024 os << " </Snapshot>\n";
296 gezelter 2
297 tim 1024 os.flush();
298     #else
299     //every node prepares the dump lines for integrable objects belong to itself
300     std::string buffer;
301     for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
302 xsun 1217
303    
304 tim 1024 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
305     integrableObject = mol->nextIntegrableObject(ii)) {
306 xsun 1217 buffer += prepareDumpLine(integrableObject);
307 gezelter 507 }
308 gezelter 2 }
309 tim 1024
310 gezelter 246 const int masterNode = 0;
311 gezelter 1629 int nProc;
312     MPI_Comm_size(MPI_COMM_WORLD, &nProc);
313 tim 1024 if (worldRank == masterNode) {
314     os << " <Snapshot>\n";
315     writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot());
316     os << " <StuntDoubles>\n";
317    
318 gezelter 1025 os << buffer;
319    
320 tim 1024 for (int i = 1; i < nProc; ++i) {
321 gezelter 2
322 tim 1024 // receive the length of the string buffer that was
323     // prepared by processor i
324 gezelter 2
325 gezelter 1629 MPI_Bcast(&i, 1, MPI_INT,masterNode,MPI_COMM_WORLD);
326 tim 1024 int recvLength;
327     MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &istatus);
328     char* recvBuffer = new char[recvLength];
329     if (recvBuffer == NULL) {
330     } else {
331 gezelter 1025 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD, &istatus);
332 tim 1024 os << recvBuffer;
333 gezelter 1313 delete [] recvBuffer;
334 tim 1024 }
335     }
336     os << " </StuntDoubles>\n";
337    
338     os << " </Snapshot>\n";
339 gezelter 507 os.flush();
340 gezelter 246 } else {
341 gezelter 1025 int sendBufferLength = buffer.size() + 1;
342 gezelter 1629 int myturn = 0;
343     for (int i = 1; i < nProc; ++i){
344     MPI_Bcast(&myturn,1, MPI_INT,masterNode,MPI_COMM_WORLD);
345     if (myturn == worldRank){
346     MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
347     MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
348     }
349     }
350 tim 1024 }
351 gezelter 2
352 tim 1024 #endif // is_mpi
353 gezelter 2
354 tim 1024 }
355 gezelter 2
356 tim 1024 std::string DumpWriter::prepareDumpLine(StuntDouble* integrableObject) {
357    
358     int index = integrableObject->getGlobalIntegrableObjectIndex();
359     std::string type("pv");
360     std::string line;
361     char tempBuffer[4096];
362 gezelter 2
363 tim 1024 Vector3d pos;
364     Vector3d vel;
365     pos = integrableObject->getPos();
366 gezelter 1437
367     if (isinf(pos[0]) || isnan(pos[0]) ||
368     isinf(pos[1]) || isnan(pos[1]) ||
369     isinf(pos[2]) || isnan(pos[2]) ) {
370     sprintf( painCave.errMsg,
371     "DumpWriter detected a numerical error writing the position"
372     " for object %d", index);
373     painCave.isFatal = 1;
374     simError();
375     }
376    
377 tim 1024 vel = integrableObject->getVel();
378 gezelter 1437
379     if (isinf(vel[0]) || isnan(vel[0]) ||
380     isinf(vel[1]) || isnan(vel[1]) ||
381     isinf(vel[2]) || isnan(vel[2]) ) {
382     sprintf( painCave.errMsg,
383     "DumpWriter detected a numerical error writing the velocity"
384     " for object %d", index);
385     painCave.isFatal = 1;
386     simError();
387     }
388    
389 gezelter 1025 sprintf(tempBuffer, "%18.10g %18.10g %18.10g %13e %13e %13e",
390 tim 1024 pos[0], pos[1], pos[2],
391     vel[0], vel[1], vel[2]);
392     line += tempBuffer;
393 gezelter 2
394 tim 1024 if (integrableObject->isDirectional()) {
395     type += "qj";
396     Quat4d q;
397     Vector3d ji;
398     q = integrableObject->getQ();
399 gezelter 1437
400     if (isinf(q[0]) || isnan(q[0]) ||
401     isinf(q[1]) || isnan(q[1]) ||
402     isinf(q[2]) || isnan(q[2]) ||
403     isinf(q[3]) || isnan(q[3]) ) {
404     sprintf( painCave.errMsg,
405     "DumpWriter detected a numerical error writing the quaternion"
406     " for object %d", index);
407     painCave.isFatal = 1;
408     simError();
409     }
410    
411 tim 1024 ji = integrableObject->getJ();
412 gezelter 1437
413     if (isinf(ji[0]) || isnan(ji[0]) ||
414     isinf(ji[1]) || isnan(ji[1]) ||
415     isinf(ji[2]) || isnan(ji[2]) ) {
416     sprintf( painCave.errMsg,
417     "DumpWriter detected a numerical error writing the angular"
418     " momentum for object %d", index);
419     painCave.isFatal = 1;
420     simError();
421     }
422    
423 gezelter 1025 sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e %13e",
424 tim 1024 q[0], q[1], q[2], q[3],
425     ji[0], ji[1], ji[2]);
426     line += tempBuffer;
427     }
428 gezelter 2
429 tim 1024 if (needForceVector_) {
430 gezelter 1437 type += "f";
431 tim 1024 Vector3d frc;
432 gezelter 1437
433 tim 1024 frc = integrableObject->getFrc();
434 gezelter 1437
435     if (isinf(frc[0]) || isnan(frc[0]) ||
436     isinf(frc[1]) || isnan(frc[1]) ||
437     isinf(frc[2]) || isnan(frc[2]) ) {
438     sprintf( painCave.errMsg,
439     "DumpWriter detected a numerical error writing the force"
440     " for object %d", index);
441     painCave.isFatal = 1;
442     simError();
443     }
444     sprintf(tempBuffer, " %13e %13e %13e",
445     frc[0], frc[1], frc[2]);
446 tim 1024 line += tempBuffer;
447 gezelter 1437
448     if (integrableObject->isDirectional()) {
449     type += "t";
450     Vector3d trq;
451    
452     trq = integrableObject->getTrq();
453    
454     if (isinf(trq[0]) || isnan(trq[0]) ||
455     isinf(trq[1]) || isnan(trq[1]) ||
456     isinf(trq[2]) || isnan(trq[2]) ) {
457     sprintf( painCave.errMsg,
458     "DumpWriter detected a numerical error writing the torque"
459     " for object %d", index);
460     painCave.isFatal = 1;
461     simError();
462     }
463    
464     sprintf(tempBuffer, " %13e %13e %13e",
465     trq[0], trq[1], trq[2]);
466     line += tempBuffer;
467     }
468 gezelter 246 }
469 gezelter 1629 if (needParticlePot_) {
470     type += "u";
471     RealType particlePot;
472    
473     particlePot = integrableObject->getParticlePot();
474    
475     if (isinf(particlePot) || isnan(particlePot)) {
476     sprintf( painCave.errMsg,
477     "DumpWriter detected a numerical error writing the particle "
478     " potential for object %d", index);
479     painCave.isFatal = 1;
480     simError();
481     }
482     sprintf(tempBuffer, " %13e", particlePot);
483     line += tempBuffer;
484     }
485 gezelter 1437
486 gezelter 1025 sprintf(tempBuffer, "%10d %7s %s\n", index, type.c_str(), line.c_str());
487 tim 1024 return std::string(tempBuffer);
488 gezelter 507 }
489 gezelter 2
490 gezelter 507 void DumpWriter::writeDump() {
491 tim 615 writeFrame(*dumpFile_);
492 gezelter 507 }
493 tim 376
494 gezelter 507 void DumpWriter::writeEor() {
495 tim 615 std::ostream* eorStream;
496 tim 376
497     #ifdef IS_MPI
498     if (worldRank == 0) {
499     #endif // is_mpi
500    
501 tim 615 eorStream = createOStream(eorFilename_);
502 tim 376
503     #ifdef IS_MPI
504     }
505     #endif // is_mpi
506    
507 tim 615 writeFrame(*eorStream);
508    
509     #ifdef IS_MPI
510     if (worldRank == 0) {
511     #endif // is_mpi
512 tim 1024 writeClosing(*eorStream);
513     delete eorStream;
514 tim 615 #ifdef IS_MPI
515     }
516     #endif // is_mpi
517    
518 gezelter 507 }
519 tim 376
520    
521 gezelter 507 void DumpWriter::writeDumpAndEor() {
522 tim 376 std::vector<std::streambuf*> buffers;
523 tim 615 std::ostream* eorStream;
524 tim 376 #ifdef IS_MPI
525     if (worldRank == 0) {
526     #endif // is_mpi
527    
528 tim 615 buffers.push_back(dumpFile_->rdbuf());
529 tim 376
530 tim 615 eorStream = createOStream(eorFilename_);
531 tim 376
532 tim 615 buffers.push_back(eorStream->rdbuf());
533 tim 376
534     #ifdef IS_MPI
535     }
536     #endif // is_mpi
537    
538     TeeBuf tbuf(buffers.begin(), buffers.end());
539     std::ostream os(&tbuf);
540    
541     writeFrame(os);
542 tim 615
543     #ifdef IS_MPI
544     if (worldRank == 0) {
545     #endif // is_mpi
546 tim 1024 writeClosing(*eorStream);
547     delete eorStream;
548 tim 615 #ifdef IS_MPI
549     }
550     #endif // is_mpi
551 tim 376
552 gezelter 507 }
553 tim 376
554 tim 1024 std::ostream* DumpWriter::createOStream(const std::string& filename) {
555 tim 619
556 tim 615 std::ostream* newOStream;
557 tim 619 #ifdef HAVE_LIBZ
558 tim 615 if (needCompression_) {
559 tim 1024 newOStream = new ogzstream(filename.c_str());
560 tim 615 } else {
561 tim 1024 newOStream = new std::ofstream(filename.c_str());
562 tim 615 }
563 tim 619 #else
564     newOStream = new std::ofstream(filename.c_str());
565     #endif
566 tim 1024 //write out MetaData first
567 gezelter 1390 (*newOStream) << "<OpenMD version=1>" << std::endl;
568 tim 1024 (*newOStream) << " <MetaData>" << std::endl;
569     (*newOStream) << info_->getRawMetaData();
570     (*newOStream) << " </MetaData>" << std::endl;
571 tim 615 return newOStream;
572 tim 1024 }
573 tim 376
574 tim 1024 void DumpWriter::writeClosing(std::ostream& os) {
575    
576 gezelter 1390 os << "</OpenMD>\n";
577 tim 1024 os.flush();
578     }
579    
580 gezelter 1390 }//end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date