| 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, 234107 (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 |  | 
| 314 | #ifndef IS_MPI | 
| 315 | os << "  <Snapshot>\n"; | 
| 316 |  | 
| 317 | writeFrameProperties(os, info_->getSnapshotManager()->getCurrentSnapshot()); | 
| 318 |  | 
| 319 | os << "    <StuntDoubles>\n"; | 
| 320 | for (mol = info_->beginMolecule(mi); mol != NULL; | 
| 321 | mol = info_->nextMolecule(mi)) { | 
| 322 |  | 
| 323 | for (sd = mol->beginIntegrableObject(ii); sd != NULL; | 
| 324 | sd = mol->nextIntegrableObject(ii)) { | 
| 325 | os << prepareDumpLine(sd); | 
| 326 |  | 
| 327 | } | 
| 328 | } | 
| 329 | os << "    </StuntDoubles>\n"; | 
| 330 |  | 
| 331 | if (doSiteData_) { | 
| 332 | os << "    <SiteData>\n"; | 
| 333 | for (mol = info_->beginMolecule(mi); mol != NULL; | 
| 334 | 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* 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 |  | 
| 362 | const int masterNode = 0; | 
| 363 | int worldRank = MPI::COMM_WORLD.Get_rank(); | 
| 364 | int nProc = MPI::COMM_WORLD.Get_size(); | 
| 365 |  | 
| 366 | if (worldRank == masterNode) { | 
| 367 | os << "  <Snapshot>\n"; | 
| 368 | writeFrameProperties(os, | 
| 369 | info_->getSnapshotManager()->getCurrentSnapshot()); | 
| 370 | os << "    <StuntDoubles>\n"; | 
| 371 | } | 
| 372 |  | 
| 373 | //every node prepares the dump lines for integrable objects belong to itself | 
| 374 | std::string buffer; | 
| 375 | for (mol = info_->beginMolecule(mi); mol != NULL; | 
| 376 | mol = info_->nextMolecule(mi)) { | 
| 377 | for (sd = mol->beginIntegrableObject(ii); sd != NULL; | 
| 378 | sd = mol->nextIntegrableObject(ii)) { | 
| 379 | buffer += prepareDumpLine(sd); | 
| 380 | } | 
| 381 | } | 
| 382 |  | 
| 383 | if (worldRank == masterNode) { | 
| 384 | os << buffer; | 
| 385 |  | 
| 386 | for (int i = 1; i < nProc; ++i) { | 
| 387 | // tell processor i to start sending us data: | 
| 388 | MPI::COMM_WORLD.Bcast(&i, 1, MPI::INT, masterNode); | 
| 389 |  | 
| 390 | // receive the length of the string buffer that was | 
| 391 | // prepared by processor i: | 
| 392 | int recvLength; | 
| 393 | MPI::COMM_WORLD.Recv(&recvLength, 1, MPI::INT, i, MPI::ANY_TAG, | 
| 394 | istatus); | 
| 395 |  | 
| 396 | // create a buffer to receive the data | 
| 397 | char* recvBuffer = new char[recvLength]; | 
| 398 | if (recvBuffer == NULL) { | 
| 399 | } else { | 
| 400 | // receive the data: | 
| 401 | MPI::COMM_WORLD.Recv(recvBuffer, recvLength, MPI::CHAR, i, | 
| 402 | MPI::ANY_TAG, istatus); | 
| 403 | // send it to the file: | 
| 404 | os << recvBuffer; | 
| 405 | // get rid of the receive buffer: | 
| 406 | delete [] recvBuffer; | 
| 407 | } | 
| 408 | } | 
| 409 | } else { | 
| 410 | int sendBufferLength = buffer.size() + 1; | 
| 411 | int myturn = 0; | 
| 412 | for (int i = 1; i < nProc; ++i){ | 
| 413 | // wait for the master node to call our number: | 
| 414 | MPI::COMM_WORLD.Bcast(&myturn, 1, MPI::INT, masterNode); | 
| 415 | if (myturn == worldRank){ | 
| 416 | // send the length of our buffer: | 
| 417 | MPI::COMM_WORLD.Send(&sendBufferLength, 1, MPI::INT, masterNode, 0); | 
| 418 |  | 
| 419 | // send our buffer: | 
| 420 | MPI::COMM_WORLD.Send((void *)buffer.c_str(), sendBufferLength, | 
| 421 | MPI::CHAR, masterNode, 0); | 
| 422 |  | 
| 423 | } | 
| 424 | } | 
| 425 | } | 
| 426 |  | 
| 427 | if (worldRank == masterNode) { | 
| 428 | os << "    </StuntDoubles>\n"; | 
| 429 | } | 
| 430 |  | 
| 431 | if (doSiteData_) { | 
| 432 | if (worldRank == masterNode) { | 
| 433 | os << "    <SiteData>\n"; | 
| 434 | } | 
| 435 | buffer.clear(); | 
| 436 | for (mol = info_->beginMolecule(mi); mol != NULL; | 
| 437 | mol = info_->nextMolecule(mi)) { | 
| 438 |  | 
| 439 | for (sd = mol->beginIntegrableObject(ii); sd != NULL; | 
| 440 | sd = mol->nextIntegrableObject(ii)) { | 
| 441 |  | 
| 442 | int ioIndex = sd->getGlobalIntegrableObjectIndex(); | 
| 443 | // do one for the IO itself | 
| 444 | buffer += prepareSiteLine(sd, ioIndex, 0); | 
| 445 |  | 
| 446 | if (sd->isRigidBody()) { | 
| 447 |  | 
| 448 | RigidBody* rb = static_cast<RigidBody*>(sd); | 
| 449 | int siteIndex = 0; | 
| 450 | for (Atom* atom = rb->beginAtom(ai); atom != NULL; | 
| 451 | atom = rb->nextAtom(ai)) { | 
| 452 | buffer += prepareSiteLine(atom, ioIndex, siteIndex); | 
| 453 | siteIndex++; | 
| 454 | } | 
| 455 | } | 
| 456 | } | 
| 457 | } | 
| 458 |  | 
| 459 | if (worldRank == masterNode) { | 
| 460 | os << buffer; | 
| 461 |  | 
| 462 | for (int i = 1; i < nProc; ++i) { | 
| 463 |  | 
| 464 | // tell processor i to start sending us data: | 
| 465 | MPI::COMM_WORLD.Bcast(&i, 1, MPI::INT, masterNode); | 
| 466 |  | 
| 467 | // receive the length of the string buffer that was | 
| 468 | // prepared by processor i: | 
| 469 | int recvLength; | 
| 470 | MPI::COMM_WORLD.Recv(&recvLength, 1, MPI::INT, i, MPI::ANY_TAG, | 
| 471 | istatus); | 
| 472 |  | 
| 473 | // create a buffer to receive the data | 
| 474 | char* recvBuffer = new char[recvLength]; | 
| 475 | if (recvBuffer == NULL) { | 
| 476 | } else { | 
| 477 | // receive the data: | 
| 478 | MPI::COMM_WORLD.Recv(recvBuffer, recvLength, MPI::CHAR, i, | 
| 479 | MPI::ANY_TAG, istatus); | 
| 480 | // send it to the file: | 
| 481 | os << recvBuffer; | 
| 482 | // get rid of the receive buffer: | 
| 483 | delete [] recvBuffer; | 
| 484 | } | 
| 485 | } | 
| 486 | } else { | 
| 487 | int sendBufferLength = buffer.size() + 1; | 
| 488 | int myturn = 0; | 
| 489 | for (int i = 1; i < nProc; ++i){ | 
| 490 | // wait for the master node to call our number: | 
| 491 | MPI::COMM_WORLD.Bcast(&myturn, 1, MPI::INT, masterNode); | 
| 492 | if (myturn == worldRank){ | 
| 493 | // send the length of our buffer: | 
| 494 | MPI::COMM_WORLD.Send(&sendBufferLength, 1, MPI::INT, masterNode, 0); | 
| 495 | // send our buffer: | 
| 496 | MPI::COMM_WORLD.Send((void *)buffer.c_str(), sendBufferLength, | 
| 497 | MPI::CHAR, masterNode, 0); | 
| 498 | } | 
| 499 | } | 
| 500 | } | 
| 501 |  | 
| 502 | if (worldRank == masterNode) { | 
| 503 | os << "    </SiteData>\n"; | 
| 504 | } | 
| 505 | } | 
| 506 |  | 
| 507 | if (worldRank == masterNode) { | 
| 508 | os << "  </Snapshot>\n"; | 
| 509 | os.flush(); | 
| 510 | } | 
| 511 |  | 
| 512 | #endif // is_mpi | 
| 513 |  | 
| 514 | } | 
| 515 |  | 
| 516 | std::string DumpWriter::prepareDumpLine(StuntDouble* sd) { | 
| 517 |  | 
| 518 | int index = sd->getGlobalIntegrableObjectIndex(); | 
| 519 | std::string type("pv"); | 
| 520 | std::string line; | 
| 521 | char tempBuffer[4096]; | 
| 522 |  | 
| 523 | Vector3d pos; | 
| 524 | Vector3d vel; | 
| 525 | pos = sd->getPos(); | 
| 526 |  | 
| 527 | if (isinf(pos[0]) || isnan(pos[0]) || | 
| 528 | isinf(pos[1]) || isnan(pos[1]) || | 
| 529 | isinf(pos[2]) || isnan(pos[2]) ) { | 
| 530 | sprintf( painCave.errMsg, | 
| 531 | "DumpWriter detected a numerical error writing the position" | 
| 532 | " for object %d", index); | 
| 533 | painCave.isFatal = 1; | 
| 534 | simError(); | 
| 535 | } | 
| 536 |  | 
| 537 | vel = sd->getVel(); | 
| 538 |  | 
| 539 | if (isinf(vel[0]) || isnan(vel[0]) || | 
| 540 | isinf(vel[1]) || isnan(vel[1]) || | 
| 541 | isinf(vel[2]) || isnan(vel[2]) ) { | 
| 542 | sprintf( painCave.errMsg, | 
| 543 | "DumpWriter detected a numerical error writing the velocity" | 
| 544 | " for object %d", index); | 
| 545 | painCave.isFatal = 1; | 
| 546 | simError(); | 
| 547 | } | 
| 548 |  | 
| 549 | sprintf(tempBuffer, "%18.10g %18.10g %18.10g %13e %13e %13e", | 
| 550 | pos[0], pos[1], pos[2], | 
| 551 | vel[0], vel[1], vel[2]); | 
| 552 | line += tempBuffer; | 
| 553 |  | 
| 554 | if (sd->isDirectional()) { | 
| 555 | type += "qj"; | 
| 556 | Quat4d q; | 
| 557 | Vector3d ji; | 
| 558 | q = sd->getQ(); | 
| 559 |  | 
| 560 | if (isinf(q[0]) || isnan(q[0]) || | 
| 561 | isinf(q[1]) || isnan(q[1]) || | 
| 562 | isinf(q[2]) || isnan(q[2]) || | 
| 563 | isinf(q[3]) || isnan(q[3]) ) { | 
| 564 | sprintf( painCave.errMsg, | 
| 565 | "DumpWriter detected a numerical error writing the quaternion" | 
| 566 | " for object %d", index); | 
| 567 | painCave.isFatal = 1; | 
| 568 | simError(); | 
| 569 | } | 
| 570 |  | 
| 571 | ji = sd->getJ(); | 
| 572 |  | 
| 573 | if (isinf(ji[0]) || isnan(ji[0]) || | 
| 574 | isinf(ji[1]) || isnan(ji[1]) || | 
| 575 | isinf(ji[2]) || isnan(ji[2]) ) { | 
| 576 | sprintf( painCave.errMsg, | 
| 577 | "DumpWriter detected a numerical error writing the angular" | 
| 578 | " momentum for object %d", index); | 
| 579 | painCave.isFatal = 1; | 
| 580 | simError(); | 
| 581 | } | 
| 582 |  | 
| 583 | sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e %13e", | 
| 584 | q[0], q[1], q[2], q[3], | 
| 585 | ji[0], ji[1], ji[2]); | 
| 586 | line += tempBuffer; | 
| 587 | } | 
| 588 |  | 
| 589 | if (needForceVector_) { | 
| 590 | type += "f"; | 
| 591 | Vector3d frc = sd->getFrc(); | 
| 592 | if (isinf(frc[0]) || isnan(frc[0]) || | 
| 593 | isinf(frc[1]) || isnan(frc[1]) || | 
| 594 | isinf(frc[2]) || isnan(frc[2]) ) { | 
| 595 | sprintf( painCave.errMsg, | 
| 596 | "DumpWriter detected a numerical error writing the force" | 
| 597 | " for object %d", index); | 
| 598 | painCave.isFatal = 1; | 
| 599 | simError(); | 
| 600 | } | 
| 601 | sprintf(tempBuffer, " %13e %13e %13e", | 
| 602 | frc[0], frc[1], frc[2]); | 
| 603 | line += tempBuffer; | 
| 604 |  | 
| 605 | if (sd->isDirectional()) { | 
| 606 | type += "t"; | 
| 607 | Vector3d trq = sd->getTrq(); | 
| 608 | if (isinf(trq[0]) || isnan(trq[0]) || | 
| 609 | isinf(trq[1]) || isnan(trq[1]) || | 
| 610 | isinf(trq[2]) || isnan(trq[2]) ) { | 
| 611 | sprintf( painCave.errMsg, | 
| 612 | "DumpWriter detected a numerical error writing the torque" | 
| 613 | " for object %d", index); | 
| 614 | painCave.isFatal = 1; | 
| 615 | simError(); | 
| 616 | } | 
| 617 | sprintf(tempBuffer, " %13e %13e %13e", | 
| 618 | trq[0], trq[1], trq[2]); | 
| 619 | line += tempBuffer; | 
| 620 | } | 
| 621 | } | 
| 622 |  | 
| 623 | sprintf(tempBuffer, "%10d %7s %s\n", index, type.c_str(), line.c_str()); | 
| 624 | return std::string(tempBuffer); | 
| 625 | } | 
| 626 |  | 
| 627 | std::string DumpWriter::prepareSiteLine(StuntDouble* sd, int ioIndex, int siteIndex) { | 
| 628 | int storageLayout = info_->getSnapshotManager()->getStorageLayout(); | 
| 629 |  | 
| 630 | std::string id; | 
| 631 | std::string type; | 
| 632 | std::string line; | 
| 633 | char tempBuffer[4096]; | 
| 634 |  | 
| 635 | if (sd->isRigidBody()) { | 
| 636 | sprintf(tempBuffer, "%10d           ", ioIndex); | 
| 637 | id = std::string(tempBuffer); | 
| 638 | } else { | 
| 639 | sprintf(tempBuffer, "%10d %10d", ioIndex, siteIndex); | 
| 640 | id = std::string(tempBuffer); | 
| 641 | } | 
| 642 |  | 
| 643 | if (needFlucQ_) { | 
| 644 | if (storageLayout & DataStorage::dslFlucQPosition) { | 
| 645 | type += "c"; | 
| 646 | RealType fqPos = sd->getFlucQPos(); | 
| 647 | if (isinf(fqPos) || isnan(fqPos) ) { | 
| 648 | sprintf( painCave.errMsg, | 
| 649 | "DumpWriter detected a numerical error writing the" | 
| 650 | " fluctuating charge for object %s", id.c_str()); | 
| 651 | painCave.isFatal = 1; | 
| 652 | simError(); | 
| 653 | } | 
| 654 | sprintf(tempBuffer, " %13e ", fqPos); | 
| 655 | line += tempBuffer; | 
| 656 | } | 
| 657 |  | 
| 658 | if (storageLayout & DataStorage::dslFlucQVelocity) { | 
| 659 | type += "w"; | 
| 660 | RealType fqVel = sd->getFlucQVel(); | 
| 661 | if (isinf(fqVel) || isnan(fqVel) ) { | 
| 662 | sprintf( painCave.errMsg, | 
| 663 | "DumpWriter detected a numerical error writing the" | 
| 664 | " fluctuating charge velocity for object %s", id.c_str()); | 
| 665 | painCave.isFatal = 1; | 
| 666 | simError(); | 
| 667 | } | 
| 668 | sprintf(tempBuffer, " %13e ", fqVel); | 
| 669 | line += tempBuffer; | 
| 670 | } | 
| 671 |  | 
| 672 | if (needForceVector_) { | 
| 673 | if (storageLayout & DataStorage::dslFlucQForce) { | 
| 674 | type += "g"; | 
| 675 | RealType fqFrc = sd->getFlucQFrc(); | 
| 676 | if (isinf(fqFrc) || isnan(fqFrc) ) { | 
| 677 | sprintf( painCave.errMsg, | 
| 678 | "DumpWriter detected a numerical error writing the" | 
| 679 | " fluctuating charge force for object %s", id.c_str()); | 
| 680 | painCave.isFatal = 1; | 
| 681 | simError(); | 
| 682 | } | 
| 683 | sprintf(tempBuffer, " %13e ", fqFrc); | 
| 684 | line += tempBuffer; | 
| 685 | } | 
| 686 | } | 
| 687 | } | 
| 688 |  | 
| 689 | if (needElectricField_) { | 
| 690 | if (storageLayout & DataStorage::dslElectricField) { | 
| 691 | type += "e"; | 
| 692 | Vector3d eField= sd->getElectricField(); | 
| 693 | if (isinf(eField[0]) || isnan(eField[0]) || | 
| 694 | isinf(eField[1]) || isnan(eField[1]) || | 
| 695 | isinf(eField[2]) || isnan(eField[2]) ) { | 
| 696 | sprintf( painCave.errMsg, | 
| 697 | "DumpWriter detected a numerical error writing the electric" | 
| 698 | " field for object %s", id.c_str()); | 
| 699 | painCave.isFatal = 1; | 
| 700 | simError(); | 
| 701 | } | 
| 702 | sprintf(tempBuffer, " %13e %13e %13e", | 
| 703 | eField[0], eField[1], eField[2]); | 
| 704 | line += tempBuffer; | 
| 705 | } | 
| 706 | } | 
| 707 |  | 
| 708 |  | 
| 709 | if (needParticlePot_) { | 
| 710 | if (storageLayout & DataStorage::dslParticlePot) { | 
| 711 | type += "u"; | 
| 712 | RealType particlePot = sd->getParticlePot(); | 
| 713 | if (isinf(particlePot) || isnan(particlePot)) { | 
| 714 | sprintf( painCave.errMsg, | 
| 715 | "DumpWriter detected a numerical error writing the particle " | 
| 716 | " potential for object %s", id.c_str()); | 
| 717 | painCave.isFatal = 1; | 
| 718 | simError(); | 
| 719 | } | 
| 720 | sprintf(tempBuffer, " %13e", particlePot); | 
| 721 | line += tempBuffer; | 
| 722 | } | 
| 723 | } | 
| 724 |  | 
| 725 | sprintf(tempBuffer, "%s %7s %s\n", id.c_str(), type.c_str(), line.c_str()); | 
| 726 | return std::string(tempBuffer); | 
| 727 | } | 
| 728 |  | 
| 729 | void DumpWriter::writeDump() { | 
| 730 | writeFrame(*dumpFile_); | 
| 731 | } | 
| 732 |  | 
| 733 | void DumpWriter::writeEor() { | 
| 734 | std::ostream* eorStream; | 
| 735 |  | 
| 736 | #ifdef IS_MPI | 
| 737 | if (worldRank == 0) { | 
| 738 | #endif // is_mpi | 
| 739 |  | 
| 740 | eorStream = createOStream(eorFilename_); | 
| 741 |  | 
| 742 | #ifdef IS_MPI | 
| 743 | } | 
| 744 | #endif | 
| 745 |  | 
| 746 | writeFrame(*eorStream); | 
| 747 |  | 
| 748 | #ifdef IS_MPI | 
| 749 | if (worldRank == 0) { | 
| 750 | #endif | 
| 751 |  | 
| 752 | writeClosing(*eorStream); | 
| 753 | delete eorStream; | 
| 754 |  | 
| 755 | #ifdef IS_MPI | 
| 756 | } | 
| 757 | #endif // is_mpi | 
| 758 |  | 
| 759 | } | 
| 760 |  | 
| 761 |  | 
| 762 | void DumpWriter::writeDumpAndEor() { | 
| 763 | std::vector<std::streambuf*> buffers; | 
| 764 | std::ostream* eorStream; | 
| 765 | #ifdef IS_MPI | 
| 766 | if (worldRank == 0) { | 
| 767 | #endif // is_mpi | 
| 768 |  | 
| 769 | buffers.push_back(dumpFile_->rdbuf()); | 
| 770 |  | 
| 771 | eorStream = createOStream(eorFilename_); | 
| 772 |  | 
| 773 | buffers.push_back(eorStream->rdbuf()); | 
| 774 |  | 
| 775 | #ifdef IS_MPI | 
| 776 | } | 
| 777 | #endif // is_mpi | 
| 778 |  | 
| 779 | TeeBuf tbuf(buffers.begin(), buffers.end()); | 
| 780 | std::ostream os(&tbuf); | 
| 781 |  | 
| 782 | writeFrame(os); | 
| 783 |  | 
| 784 | #ifdef IS_MPI | 
| 785 | if (worldRank == 0) { | 
| 786 | #endif // is_mpi | 
| 787 | writeClosing(*eorStream); | 
| 788 | delete eorStream; | 
| 789 | #ifdef IS_MPI | 
| 790 | } | 
| 791 | #endif // is_mpi | 
| 792 |  | 
| 793 | } | 
| 794 |  | 
| 795 | std::ostream* DumpWriter::createOStream(const std::string& filename) { | 
| 796 |  | 
| 797 | std::ostream* newOStream; | 
| 798 | #ifdef HAVE_ZLIB | 
| 799 | if (needCompression_) { | 
| 800 | newOStream = new ogzstream(filename.c_str()); | 
| 801 | } else { | 
| 802 | newOStream = new std::ofstream(filename.c_str()); | 
| 803 | } | 
| 804 | #else | 
| 805 | newOStream = new std::ofstream(filename.c_str()); | 
| 806 | #endif | 
| 807 | //write out MetaData first | 
| 808 | (*newOStream) << "<OpenMD version=2>" << std::endl; | 
| 809 | (*newOStream) << "  <MetaData>" << std::endl; | 
| 810 | (*newOStream) << info_->getRawMetaData(); | 
| 811 | (*newOStream) << "  </MetaData>" << std::endl; | 
| 812 | return newOStream; | 
| 813 | } | 
| 814 |  | 
| 815 | void DumpWriter::writeClosing(std::ostream& os) { | 
| 816 |  | 
| 817 | os << "</OpenMD>\n"; | 
| 818 | os.flush(); | 
| 819 | } | 
| 820 |  | 
| 821 | }//end namespace OpenMD |