ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/DumpWriter.cpp
Revision: 1993
Committed: Tue Apr 29 17:32:31 2014 UTC (11 years ago) by gezelter
File size: 25278 byte(s)
Log Message:
Added sitePotentials for the Choi et al. potential-frequency maps for nitriles

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

Properties

Name Value
svn:keywords Author Id Revision Date