ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpWriter.cpp
(Generate patch)

Comparing:
trunk/src/io/DumpWriter.cpp (file contents), Revision 1025 by gezelter, Wed Aug 30 20:33:44 2006 UTC vs.
branches/development/src/io/DumpWriter.cpp (file contents), Revision 1465 by chuckv, Fri Jul 9 23:08:25 2010 UTC

# Line 1 | Line 1
1   /*
2 < * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
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. Acknowledgement of the program authors must be made in any
10 < *    publication of scientific results based in part on use of the
11 < *    program.  An acceptable form of acknowledgement is citation of
12 < *    the article in which the program was described (Matthew
13 < *    A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 < *    J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 < *    Parallel Simulation Engine for Molecular Dynamics,"
16 < *    J. Comput. Chem. 26, pp. 252-271 (2005))
17 < *
18 < * 2. Redistributions of source code must retain the above copyright
9 > * 1. Redistributions of source code must retain the above copyright
10   *    notice, this list of conditions and the following disclaimer.
11   *
12 < * 3. Redistributions in binary form must reproduce the above copyright
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.
# Line 37 | Line 28
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]  Vardeman & Gezelter, in progress (2009).                        
40   */
41  
42   #include "io/DumpWriter.hpp"
# Line 46 | Line 46
46   #include "io/gzstream.hpp"
47   #include "io/Globals.hpp"
48  
49 +
50   #ifdef IS_MPI
51   #include <mpi.h>
52   #endif //is_mpi
53  
54 < namespace oopse {
54 > using namespace std;
55 > namespace OpenMD {
56  
57    DumpWriter::DumpWriter(SimInfo* info)
58      : info_(info), filename_(info->getDumpFileName()), eorFilename_(info->getFinalConfigFileName()){
# Line 194 | Line 196 | namespace oopse {
196      os << "    <FrameData>\n";
197  
198      RealType currentTime = s->getTime();
199 +
200 +    if (isinf(currentTime) || isnan(currentTime)) {      
201 +      sprintf( painCave.errMsg,
202 +               "DumpWriter detected a numerical error writing the time");      
203 +      painCave.isFatal = 1;
204 +      simError();
205 +    }
206 +    
207      sprintf(buffer, "        Time: %.10g\n", currentTime);
208      os << buffer;
209  
210      Mat3x3d hmat;
211      hmat = s->getHmat();
212 +
213 +    for (unsigned int i = 0; i < 3; i++) {
214 +      for (unsigned int j = 0; j < 3; j++) {
215 +        if (isinf(hmat(i,j)) || isnan(hmat(i,j))) {      
216 +          sprintf( painCave.errMsg,
217 +                   "DumpWriter detected a numerical error writing the box");
218 +          painCave.isFatal = 1;
219 +          simError();
220 +        }        
221 +      }
222 +    }
223 +    
224      sprintf(buffer, "        Hmat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
225              hmat(0, 0), hmat(1, 0), hmat(2, 0),
226              hmat(0, 1), hmat(1, 1), hmat(2, 1),
# Line 207 | Line 229 | namespace oopse {
229  
230      RealType chi = s->getChi();
231      RealType integralOfChiDt = s->getIntegralOfChiDt();
232 +    if (isinf(chi) || isnan(chi) ||
233 +        isinf(integralOfChiDt) || isnan(integralOfChiDt)) {      
234 +      sprintf( painCave.errMsg,
235 +               "DumpWriter detected a numerical error writing the thermostat");
236 +      painCave.isFatal = 1;
237 +      simError();
238 +    }
239      sprintf(buffer, "  Thermostat: %.10g , %.10g\n", chi, integralOfChiDt);
240      os << buffer;
241  
242      Mat3x3d eta;
243      eta = s->getEta();
244 +
245 +    for (unsigned int i = 0; i < 3; i++) {
246 +      for (unsigned int j = 0; j < 3; j++) {
247 +        if (isinf(eta(i,j)) || isnan(eta(i,j))) {      
248 +          sprintf( painCave.errMsg,
249 +                   "DumpWriter detected a numerical error writing the barostat");
250 +          painCave.isFatal = 1;
251 +          simError();
252 +        }        
253 +      }
254 +    }
255 +
256      sprintf(buffer, "    Barostat: {{ %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }, { %.10g, %.10g, %.10g }}\n",
257              eta(0, 0), eta(1, 0), eta(2, 0),
258              eta(0, 1), eta(1, 1), eta(2, 1),
# Line 240 | Line 281 | namespace oopse {
281      os << "    <StuntDoubles>\n";
282      for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
283  
284 <      for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
284 >      
285 >      for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;  
286             integrableObject = mol->nextIntegrableObject(ii)) {  
287 <        os << prepareDumpLine(integrableObject);
288 <
287 >          os << prepareDumpLine(integrableObject);
288 >          
289        }
290      }    
291      os << "    </StuntDoubles>\n";
292 <
292 >    
293      os << "  </Snapshot>\n";
294  
295      os.flush();
# Line 255 | Line 297 | namespace oopse {
297      //every node prepares the dump lines for integrable objects belong to itself
298      std::string buffer;
299      for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
300 <      
300 >
301 >
302        for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
303             integrableObject = mol->nextIntegrableObject(ii)) {  
304 <        buffer += prepareDumpLine(integrableObject);
304 >          buffer += prepareDumpLine(integrableObject);
305        }
306      }
307      
# Line 285 | Line 328 | namespace oopse {
328          } else {
329            MPI_Recv(recvBuffer, recvLength, MPI_CHAR, i, 0, MPI_COMM_WORLD, &istatus);
330            os << recvBuffer;
331 <          delete recvBuffer;
331 >          delete [] recvBuffer;
332          }
333        }
334        os << "    </StuntDoubles>\n";
# Line 312 | Line 355 | namespace oopse {
355      Vector3d pos;
356      Vector3d vel;
357      pos = integrableObject->getPos();
358 +
359 +    if (isinf(pos[0]) || isnan(pos[0]) ||
360 +        isinf(pos[1]) || isnan(pos[1]) ||
361 +        isinf(pos[2]) || isnan(pos[2]) ) {      
362 +      sprintf( painCave.errMsg,
363 +               "DumpWriter detected a numerical error writing the position"
364 +               " for object %d", index);      
365 +      painCave.isFatal = 1;
366 +      simError();
367 +    }
368 +
369      vel = integrableObject->getVel();          
370 +
371 +    if (isinf(vel[0]) || isnan(vel[0]) ||
372 +        isinf(vel[1]) || isnan(vel[1]) ||
373 +        isinf(vel[2]) || isnan(vel[2]) ) {      
374 +      sprintf( painCave.errMsg,
375 +               "DumpWriter detected a numerical error writing the velocity"
376 +               " for object %d", index);      
377 +      painCave.isFatal = 1;
378 +      simError();
379 +    }
380 +
381      sprintf(tempBuffer, "%18.10g %18.10g %18.10g %13e %13e %13e",
382              pos[0], pos[1], pos[2],
383              vel[0], vel[1], vel[2]);                    
# Line 323 | Line 388 | namespace oopse {
388        Quat4d q;
389        Vector3d ji;
390        q = integrableObject->getQ();
391 +
392 +      if (isinf(q[0]) || isnan(q[0]) ||
393 +          isinf(q[1]) || isnan(q[1]) ||
394 +          isinf(q[2]) || isnan(q[2]) ||
395 +          isinf(q[3]) || isnan(q[3]) ) {      
396 +        sprintf( painCave.errMsg,
397 +                 "DumpWriter detected a numerical error writing the quaternion"
398 +                 " for object %d", index);      
399 +        painCave.isFatal = 1;
400 +        simError();
401 +      }
402 +
403        ji = integrableObject->getJ();
404 +
405 +      if (isinf(ji[0]) || isnan(ji[0]) ||
406 +          isinf(ji[1]) || isnan(ji[1]) ||
407 +          isinf(ji[2]) || isnan(ji[2]) ) {      
408 +        sprintf( painCave.errMsg,
409 +                 "DumpWriter detected a numerical error writing the angular"
410 +                 " momentum for object %d", index);      
411 +        painCave.isFatal = 1;
412 +        simError();
413 +      }
414 +
415        sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e %13e",
416                q[0], q[1], q[2], q[3],
417                ji[0], ji[1], ji[2]);
# Line 331 | Line 419 | namespace oopse {
419      }
420  
421      if (needForceVector_) {
422 <      type += "ft";
422 >      type += "f";
423        Vector3d frc;
424 <      Vector3d trq;
424 >
425        frc = integrableObject->getFrc();
426 <      trq = integrableObject->getTrq();
427 <              
428 <      sprintf(tempBuffer, " %13e %13e %13e %13e %13e %13e",
429 <              frc[0], frc[1], frc[2],
430 <              trq[0], trq[1], trq[2]);
426 >
427 >      if (isinf(frc[0]) || isnan(frc[0]) ||
428 >          isinf(frc[1]) || isnan(frc[1]) ||
429 >          isinf(frc[2]) || isnan(frc[2]) ) {      
430 >        sprintf( painCave.errMsg,
431 >                 "DumpWriter detected a numerical error writing the force"
432 >                 " for object %d", index);      
433 >        painCave.isFatal = 1;
434 >        simError();
435 >      }
436 >      sprintf(tempBuffer, " %13e %13e %13e",
437 >              frc[0], frc[1], frc[2]);
438        line += tempBuffer;
439 +      
440 +      if (integrableObject->isDirectional()) {
441 +        type += "t";
442 +        Vector3d trq;
443 +        
444 +        trq = integrableObject->getTrq();
445 +        
446 +        if (isinf(trq[0]) || isnan(trq[0]) ||
447 +            isinf(trq[1]) || isnan(trq[1]) ||
448 +            isinf(trq[2]) || isnan(trq[2]) ) {      
449 +          sprintf( painCave.errMsg,
450 +                   "DumpWriter detected a numerical error writing the torque"
451 +                   " for object %d", index);      
452 +          painCave.isFatal = 1;
453 +          simError();
454 +        }
455 +        
456 +        sprintf(tempBuffer, " %13e %13e %13e",
457 +                trq[0], trq[1], trq[2]);
458 +        line += tempBuffer;
459 +      }      
460      }
461 <        
461 >    
462      sprintf(tempBuffer, "%10d %7s %s\n", index, type.c_str(), line.c_str());
463      return std::string(tempBuffer);
464    }
# Line 424 | Line 540 | namespace oopse {
540      newOStream = new std::ofstream(filename.c_str());
541   #endif
542      //write out MetaData first
543 <    (*newOStream) << "<OOPSE version=4>" << std::endl;
543 >    (*newOStream) << "<OpenMD version=1>" << std::endl;
544      (*newOStream) << "  <MetaData>" << std::endl;
545      (*newOStream) << info_->getRawMetaData();
546      (*newOStream) << "  </MetaData>" << std::endl;
# Line 433 | Line 549 | namespace oopse {
549  
550    void DumpWriter::writeClosing(std::ostream& os) {
551  
552 <    os << "</OOPSE>\n";
552 >    os << "</OpenMD>\n";
553      os.flush();
554    }
555  
556 < }//end namespace oopse
556 > }//end namespace OpenMD

Comparing:
trunk/src/io/DumpWriter.cpp (property svn:keywords), Revision 1025 by gezelter, Wed Aug 30 20:33:44 2006 UTC vs.
branches/development/src/io/DumpWriter.cpp (property svn:keywords), Revision 1465 by chuckv, Fri Jul 9 23:08:25 2010 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines