ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/DumpWriter.cpp
Revision: 1983
Committed: Tue Apr 15 20:36:19 2014 UTC (11 years ago) by gezelter
File size: 24422 byte(s)
Log Message:
Fixes for ConstraintWriter in parallel, some cppcheck cleanup
starting preparation for 2.2 release

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

Properties

Name Value
svn:keywords Author Id Revision Date