ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/DumpWriter.cpp
Revision: 1793
Committed: Fri Aug 31 21:16:10 2012 UTC (12 years, 8 months ago) by gezelter
File size: 20520 byte(s)
Log Message:
Cleaning up some warning messages on linux

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, 24107 (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 #include "io/DumpWriter.hpp"
44 #include "primitives/Molecule.hpp"
45 #include "utils/simError.h"
46 #include "io/basic_teebuf.hpp"
47 #ifdef HAVE_ZLIB
48 #include "io/gzstream.hpp"
49 #endif
50 #include "io/Globals.hpp"
51
52 #ifdef _MSC_VER
53 #define isnan(x) _isnan((x))
54 #define isinf(x) (!_finite(x) && !_isnan(x))
55 #endif
56
57 #ifdef IS_MPI
58 #include <mpi.h>
59 #endif
60
61 using namespace std;
62 namespace OpenMD {
63
64 DumpWriter::DumpWriter(SimInfo* info)
65 : info_(info), filename_(info->getDumpFileName()), eorFilename_(info->getFinalConfigFileName()){
66
67 Globals* simParams = info->getSimParams();
68 needCompression_ = simParams->getCompressDumpFile();
69 needForceVector_ = simParams->getOutputForceVector();
70 needParticlePot_ = simParams->getOutputParticlePotential();
71 needFlucQ_ = simParams->getOutputFluctuatingCharges();
72 needElectricField_ = simParams->getOutputElectricField();
73
74 if (needParticlePot_ || needFlucQ_ || needElectricField_) {
75 doSiteData_ = true;
76 } else {
77 doSiteData_ = false;
78 }
79
80 createDumpFile_ = true;
81 #ifdef HAVE_LIBZ
82 if (needCompression_) {
83 filename_ += ".gz";
84 eorFilename_ += ".gz";
85 }
86 #endif
87
88 #ifdef IS_MPI
89
90 if (worldRank == 0) {
91 #endif // is_mpi
92
93 dumpFile_ = createOStream(filename_);
94
95 if (!dumpFile_) {
96 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
97 filename_.c_str());
98 painCave.isFatal = 1;
99 simError();
100 }
101
102 #ifdef IS_MPI
103
104 }
105
106 #endif // is_mpi
107
108 }
109
110
111 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename)
112 : info_(info), filename_(filename){
113
114 Globals* simParams = info->getSimParams();
115 eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
116
117 needCompression_ = simParams->getCompressDumpFile();
118 needForceVector_ = simParams->getOutputForceVector();
119 needParticlePot_ = simParams->getOutputParticlePotential();
120 needFlucQ_ = simParams->getOutputFluctuatingCharges();
121 needElectricField_ = simParams->getOutputElectricField();
122
123 if (needParticlePot_ || needFlucQ_ || needElectricField_) {
124 doSiteData_ = true;
125 } else {
126 doSiteData_ = false;
127 }
128
129 createDumpFile_ = true;
130 #ifdef HAVE_LIBZ
131 if (needCompression_) {
132 filename_ += ".gz";
133 eorFilename_ += ".gz";
134 }
135 #endif
136
137 #ifdef IS_MPI
138
139 if (worldRank == 0) {
140 #endif // is_mpi
141
142
143 dumpFile_ = createOStream(filename_);
144
145 if (!dumpFile_) {
146 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
147 filename_.c_str());
148 painCave.isFatal = 1;
149 simError();
150 }
151
152 #ifdef IS_MPI
153
154 }
155
156 #endif // is_mpi
157
158 }
159
160 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename, bool writeDumpFile)
161 : info_(info), filename_(filename){
162
163 Globals* simParams = info->getSimParams();
164 eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
165
166 needCompression_ = simParams->getCompressDumpFile();
167 needForceVector_ = simParams->getOutputForceVector();
168 needParticlePot_ = simParams->getOutputParticlePotential();
169 needFlucQ_ = simParams->getOutputFluctuatingCharges();
170 needElectricField_ = simParams->getOutputElectricField();
171
172 if (needParticlePot_ || needFlucQ_ || needElectricField_) {
173 doSiteData_ = true;
174 } else {
175 doSiteData_ = false;
176 }
177
178 #ifdef HAVE_LIBZ
179 if (needCompression_) {
180 filename_ += ".gz";
181 eorFilename_ += ".gz";
182 }
183 #endif
184
185 #ifdef IS_MPI
186
187 if (worldRank == 0) {
188 #endif // is_mpi
189
190 createDumpFile_ = writeDumpFile;
191 if (createDumpFile_) {
192 dumpFile_ = createOStream(filename_);
193
194 if (!dumpFile_) {
195 sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n",
196 filename_.c_str());
197 painCave.isFatal = 1;
198 simError();
199 }
200 }
201 #ifdef IS_MPI
202
203 }
204
205
206 #endif // is_mpi
207
208 }
209
210 DumpWriter::~DumpWriter() {
211
212 #ifdef IS_MPI
213
214 if (worldRank == 0) {
215 #endif // is_mpi
216 if (createDumpFile_){
217 writeClosing(*dumpFile_);
218 delete dumpFile_;
219 }
220 #ifdef IS_MPI
221
222 }
223
224 #endif // is_mpi
225
226 }
227
228 void DumpWriter::writeFrameProperties(std::ostream& os, Snapshot* s) {
229
230 char buffer[1024];
231
232 os << " <FrameData>\n";
233
234 RealType currentTime = s->getTime();
235
236 if (isinf(currentTime) || isnan(currentTime)) {
237 sprintf( painCave.errMsg,
238 "DumpWriter detected a numerical error writing the time");
239 painCave.isFatal = 1;
240 simError();
241 }
242
243 sprintf(buffer, " Time: %.10g\n", currentTime);
244 os << buffer;
245
246 Mat3x3d hmat;
247 hmat = s->getHmat();
248
249 for (unsigned int i = 0; i < 3; i++) {
250 for (unsigned int j = 0; j < 3; j++) {
251 if (isinf(hmat(i,j)) || isnan(hmat(i,j))) {
252 sprintf( painCave.errMsg,
253 "DumpWriter detected a numerical error writing the box");
254 painCave.isFatal = 1;
255 simError();
256 }
257 }
258 }
259
260 sprintf(buffer, " Hmat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
261 hmat(0, 0), hmat(1, 0), hmat(2, 0),
262 hmat(0, 1), hmat(1, 1), hmat(2, 1),
263 hmat(0, 2), hmat(1, 2), hmat(2, 2));
264 os << buffer;
265
266 pair<RealType, RealType> thermostat = s->getThermostat();
267
268 if (isinf(thermostat.first) || isnan(thermostat.first) ||
269 isinf(thermostat.second) || isnan(thermostat.second)) {
270 sprintf( painCave.errMsg,
271 "DumpWriter detected a numerical error writing the thermostat");
272 painCave.isFatal = 1;
273 simError();
274 }
275 sprintf(buffer, " Thermostat: %.10g , %.10g\n", thermostat.first,
276 thermostat.second);
277 os << buffer;
278
279 Mat3x3d eta;
280 eta = s->getBarostat();
281
282 for (unsigned int i = 0; i < 3; i++) {
283 for (unsigned int j = 0; j < 3; j++) {
284 if (isinf(eta(i,j)) || isnan(eta(i,j))) {
285 sprintf( painCave.errMsg,
286 "DumpWriter detected a numerical error writing the barostat");
287 painCave.isFatal = 1;
288 simError();
289 }
290 }
291 }
292
293 sprintf(buffer, " Barostat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
294 eta(0, 0), eta(1, 0), eta(2, 0),
295 eta(0, 1), eta(1, 1), eta(2, 1),
296 eta(0, 2), eta(1, 2), eta(2, 2));
297 os << buffer;
298
299 os << " </FrameData>\n";
300 }
301
302 void DumpWriter::writeFrame(std::ostream& os) {
303
304 #ifdef IS_MPI
305 MPI_Status istatus;
306 #endif
307
308 Molecule* mol;
309 StuntDouble* sd;
310 SimInfo::MoleculeIterator mi;
311 Molecule::IntegrableObjectIterator ii;
312 RigidBody::AtomIterator ai;
313 Atom* atom;
314
315 #ifndef IS_MPI
316 os << " <Snapshot>\n";
317
318 writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot());
319
320 os << " <StuntDoubles>\n";
321 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
322
323
324 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
325 sd = mol->nextIntegrableObject(ii)) {
326 os << prepareDumpLine(sd);
327
328 }
329 }
330 os << " </StuntDoubles>\n";
331
332 if (doSiteData_) {
333 os << " <SiteData>\n";
334 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
335
336 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
337 sd = mol->nextIntegrableObject(ii)) {
338
339 int ioIndex = sd->getGlobalIntegrableObjectIndex();
340 // do one for the IO itself
341 os << prepareSiteLine(sd, ioIndex, 0);
342
343 if (sd->isRigidBody()) {
344
345 RigidBody* rb = static_cast<RigidBody*>(sd);
346 int siteIndex = 0;
347 for (atom = rb->beginAtom(ai); atom != NULL;
348 atom = rb->nextAtom(ai)) {
349 os << prepareSiteLine(atom, ioIndex, siteIndex);
350 siteIndex++;
351 }
352 }
353 }
354 }
355 os << " </SiteData>\n";
356 }
357 os << " </Snapshot>\n";
358
359 os.flush();
360 #else
361 //every node prepares the dump lines for integrable objects belong to itself
362 std::string buffer;
363 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
364
365
366 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
367 sd = mol->nextIntegrableObject(ii)) {
368 buffer += prepareDumpLine(sd);
369 }
370 }
371
372 const int masterNode = 0;
373 int nProc;
374 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
375 if (worldRank == masterNode) {
376 os << " <Snapshot>\n";
377 writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot());
378 os << " <StuntDoubles>\n";
379
380 os << buffer;
381
382 for (int i = 1; i < nProc; ++i) {
383
384 // receive the length of the string buffer that was
385 // prepared by processor i
386
387 MPI_Bcast(&i, 1, MPI_INT,masterNode,MPI_COMM_WORLD);
388 int recvLength;
389 MPI_Recv(&recvLength, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &istatus);
390 char* recvBuffer = new char[recvLength];
391 if (recvBuffer == NULL) {
392 } else {
393 MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD, &istatus);
394 os << recvBuffer;
395 delete [] recvBuffer;
396 }
397 }
398 os << " </StuntDoubles>\n";
399
400 os << " </Snapshot>\n";
401 os.flush();
402 } else {
403 int sendBufferLength = buffer.size() + 1;
404 int myturn = 0;
405 for (int i = 1; i < nProc; ++i){
406 MPI_Bcast(&myturn,1, MPI_INT,masterNode,MPI_COMM_WORLD);
407 if (myturn == worldRank){
408 MPI_Send(&sendBufferLength, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
409 MPI_Send((void *)buffer.c_str(), sendBufferLength, MPI_CHAR, masterNode, 0, MPI_COMM_WORLD);
410 }
411 }
412 }
413
414 #endif // is_mpi
415
416 }
417
418 std::string DumpWriter::prepareDumpLine(StuntDouble* sd) {
419
420 int index = sd->getGlobalIntegrableObjectIndex();
421 std::string type("pv");
422 std::string line;
423 char tempBuffer[4096];
424
425 Vector3d pos;
426 Vector3d vel;
427 pos = sd->getPos();
428
429 if (isinf(pos[0]) || isnan(pos[0]) ||
430 isinf(pos[1]) || isnan(pos[1]) ||
431 isinf(pos[2]) || isnan(pos[2]) ) {
432 sprintf( painCave.errMsg,
433 "DumpWriter detected a numerical error writing the position"
434 " for object %d", index);
435 painCave.isFatal = 1;
436 simError();
437 }
438
439 vel = sd->getVel();
440
441 if (isinf(vel[0]) || isnan(vel[0]) ||
442 isinf(vel[1]) || isnan(vel[1]) ||
443 isinf(vel[2]) || isnan(vel[2]) ) {
444 sprintf( painCave.errMsg,
445 "DumpWriter detected a numerical error writing the velocity"
446 " for object %d", index);
447 painCave.isFatal = 1;
448 simError();
449 }
450
451 sprintf(tempBuffer, "%18.10g %18.10g %18.10g %13e %13e %13e",
452 pos[0], pos[1], pos[2],
453 vel[0], vel[1], vel[2]);
454 line += tempBuffer;
455
456 if (sd->isDirectional()) {
457 type += "qj";
458 Quat4d q;
459 Vector3d ji;
460 q = sd->getQ();
461
462 if (isinf(q[0]) || isnan(q[0]) ||
463 isinf(q[1]) || isnan(q[1]) ||
464 isinf(q[2]) || isnan(q[2]) ||
465 isinf(q[3]) || isnan(q[3]) ) {
466 sprintf( painCave.errMsg,
467 "DumpWriter detected a numerical error writing the quaternion"
468 " for object %d", index);
469 painCave.isFatal = 1;
470 simError();
471 }
472
473 ji = sd->getJ();
474
475 if (isinf(ji[0]) || isnan(ji[0]) ||
476 isinf(ji[1]) || isnan(ji[1]) ||
477 isinf(ji[2]) || isnan(ji[2]) ) {
478 sprintf( painCave.errMsg,
479 "DumpWriter detected a numerical error writing the angular"
480 " momentum for object %d", index);
481 painCave.isFatal = 1;
482 simError();
483 }
484
485 sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e %13e",
486 q[0], q[1], q[2], q[3],
487 ji[0], ji[1], ji[2]);
488 line += tempBuffer;
489 }
490
491 if (needForceVector_) {
492 type += "f";
493 Vector3d frc = sd->getFrc();
494 if (isinf(frc[0]) || isnan(frc[0]) ||
495 isinf(frc[1]) || isnan(frc[1]) ||
496 isinf(frc[2]) || isnan(frc[2]) ) {
497 sprintf( painCave.errMsg,
498 "DumpWriter detected a numerical error writing the force"
499 " for object %d", index);
500 painCave.isFatal = 1;
501 simError();
502 }
503 sprintf(tempBuffer, " %13e %13e %13e",
504 frc[0], frc[1], frc[2]);
505 line += tempBuffer;
506
507 if (sd->isDirectional()) {
508 type += "t";
509 Vector3d trq = sd->getTrq();
510 if (isinf(trq[0]) || isnan(trq[0]) ||
511 isinf(trq[1]) || isnan(trq[1]) ||
512 isinf(trq[2]) || isnan(trq[2]) ) {
513 sprintf( painCave.errMsg,
514 "DumpWriter detected a numerical error writing the torque"
515 " for object %d", index);
516 painCave.isFatal = 1;
517 simError();
518 }
519 sprintf(tempBuffer, " %13e %13e %13e",
520 trq[0], trq[1], trq[2]);
521 line += tempBuffer;
522 }
523 }
524
525 sprintf(tempBuffer, "%10d %7s %s\n", index, type.c_str(), line.c_str());
526 return std::string(tempBuffer);
527 }
528
529 std::string DumpWriter::prepareSiteLine(StuntDouble* sd, int ioIndex, int siteIndex) {
530
531
532 std::string id;
533 std::string type;
534 std::string line;
535 char tempBuffer[4096];
536
537 if (sd->isRigidBody()) {
538 sprintf(tempBuffer, "%10d ", ioIndex);
539 id = std::string(tempBuffer);
540 } else {
541 sprintf(tempBuffer, "%10d %10d", ioIndex, siteIndex);
542 id = std::string(tempBuffer);
543 }
544
545 if (needFlucQ_) {
546 type += "cw";
547 RealType fqPos = sd->getFlucQPos();
548 if (isinf(fqPos) || isnan(fqPos) ) {
549 sprintf( painCave.errMsg,
550 "DumpWriter detected a numerical error writing the"
551 " fluctuating charge for object %s", id.c_str());
552 painCave.isFatal = 1;
553 simError();
554 }
555 sprintf(tempBuffer, " %13e ", fqPos);
556 line += tempBuffer;
557
558 RealType fqVel = sd->getFlucQVel();
559 if (isinf(fqVel) || isnan(fqVel) ) {
560 sprintf( painCave.errMsg,
561 "DumpWriter detected a numerical error writing the"
562 " fluctuating charge velocity for object %s", id.c_str());
563 painCave.isFatal = 1;
564 simError();
565 }
566 sprintf(tempBuffer, " %13e ", fqVel);
567 line += tempBuffer;
568
569 if (needForceVector_) {
570 type += "g";
571 RealType fqFrc = sd->getFlucQFrc();
572 if (isinf(fqFrc) || isnan(fqFrc) ) {
573 sprintf( painCave.errMsg,
574 "DumpWriter detected a numerical error writing the"
575 " fluctuating charge force for object %s", id.c_str());
576 painCave.isFatal = 1;
577 simError();
578 }
579 sprintf(tempBuffer, " %13e ", fqFrc);
580 line += tempBuffer;
581 }
582 }
583
584 if (needElectricField_) {
585 type += "e";
586 Vector3d eField= sd->getElectricField();
587 if (isinf(eField[0]) || isnan(eField[0]) ||
588 isinf(eField[1]) || isnan(eField[1]) ||
589 isinf(eField[2]) || isnan(eField[2]) ) {
590 sprintf( painCave.errMsg,
591 "DumpWriter detected a numerical error writing the electric"
592 " field for object %s", id.c_str());
593 painCave.isFatal = 1;
594 simError();
595 }
596 sprintf(tempBuffer, " %13e %13e %13e",
597 eField[0], eField[1], eField[2]);
598 line += tempBuffer;
599 }
600
601
602 if (needParticlePot_) {
603 type += "u";
604 RealType particlePot = sd->getParticlePot();
605 if (isinf(particlePot) || isnan(particlePot)) {
606 sprintf( painCave.errMsg,
607 "DumpWriter detected a numerical error writing the particle "
608 " potential for object %s", id.c_str());
609 painCave.isFatal = 1;
610 simError();
611 }
612 sprintf(tempBuffer, " %13e", particlePot);
613 line += tempBuffer;
614 }
615
616
617 sprintf(tempBuffer, "%s %7s %s\n", id.c_str(), type.c_str(), line.c_str());
618 return std::string(tempBuffer);
619 }
620
621 void DumpWriter::writeDump() {
622 writeFrame(*dumpFile_);
623 }
624
625 void DumpWriter::writeEor() {
626 std::ostream* eorStream;
627
628 #ifdef IS_MPI
629 if (worldRank == 0) {
630 #endif // is_mpi
631
632 eorStream = createOStream(eorFilename_);
633
634 #ifdef IS_MPI
635 }
636 #endif // is_mpi
637
638 writeFrame(*eorStream);
639
640 #ifdef IS_MPI
641 if (worldRank == 0) {
642 #endif // is_mpi
643 writeClosing(*eorStream);
644 delete eorStream;
645 #ifdef IS_MPI
646 }
647 #endif // is_mpi
648
649 }
650
651
652 void DumpWriter::writeDumpAndEor() {
653 std::vector<std::streambuf*> buffers;
654 std::ostream* eorStream;
655 #ifdef IS_MPI
656 if (worldRank == 0) {
657 #endif // is_mpi
658
659 buffers.push_back(dumpFile_->rdbuf());
660
661 eorStream = createOStream(eorFilename_);
662
663 buffers.push_back(eorStream->rdbuf());
664
665 #ifdef IS_MPI
666 }
667 #endif // is_mpi
668
669 TeeBuf tbuf(buffers.begin(), buffers.end());
670 std::ostream os(&tbuf);
671
672 writeFrame(os);
673
674 #ifdef IS_MPI
675 if (worldRank == 0) {
676 #endif // is_mpi
677 writeClosing(*eorStream);
678 delete eorStream;
679 #ifdef IS_MPI
680 }
681 #endif // is_mpi
682
683 }
684
685 std::ostream* DumpWriter::createOStream(const std::string& filename) {
686
687 std::ostream* newOStream;
688 #ifdef HAVE_ZLIB
689 if (needCompression_) {
690 newOStream = new ogzstream(filename.c_str());
691 } else {
692 newOStream = new std::ofstream(filename.c_str());
693 }
694 #else
695 newOStream = new std::ofstream(filename.c_str());
696 #endif
697 //write out MetaData first
698 (*newOStream) << "<OpenMD version=2>" << std::endl;
699 (*newOStream) << " <MetaData>" << std::endl;
700 (*newOStream) << info_->getRawMetaData();
701 (*newOStream) << " </MetaData>" << std::endl;
702 return newOStream;
703 }
704
705 void DumpWriter::writeClosing(std::ostream& os) {
706
707 os << "</OpenMD>\n";
708 os.flush();
709 }
710
711 }//end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date