![]() |
MFEM v4.9.0
Finite element discretization library
|
Low-level class for writing VTKHDF data (for use in ParaView). More...
#include <vtkhdf.hpp>
Classes | |
| struct | Restart |
| Helper used in VTKHDF::VTKHDF for enabling/disabling restart mode. More... | |
Public Member Functions | |
| VTKHDF (const std::string &filename, Restart restart=Restart::Disabled()) | |
| Create a new VTKHDF file for serial I/O. | |
| VTKHDF (const std::string &filename, MPI_Comm comm_, Restart restart=Restart::Disabled()) | |
| Create a new VTKHDF file for parallel I/O. | |
| void | UpdateSteps (real_t t) |
| Update the time step data after saving the mesh and grid functions. | |
| void | DisableCompression () |
| Disable zlib compression. | |
| void | EnableCompression (int level=6) |
| Enable zlib compression at the specified level. | |
| template<typename FP_T = real_t> | |
| void | SaveMesh (const Mesh &mesh, bool high_order=true, int ref=-1) |
| Save the mesh, appending as a new time step. | |
| template<typename FP_T = real_t> | |
| void | SaveGridFunction (const GridFunction &gf, const std::string &name) |
| Save the grid function with the given name, appending as a new time step. | |
| void | Flush () |
| Flush the file. | |
| ~VTKHDF () | |
Not copyable or movable. | |
| VTKHDF (const VTKHDF &)=delete | |
| VTKHDF (VTKHDF &&)=delete | |
| VTKHDF & | operator= (const VTKHDF &)=delete |
| VTKHDF & | operator= (VTKHDF &&)=delete |
Low-level class for writing VTKHDF data (for use in ParaView).
Users should typically use ParaViewHDFDataCollection, Mesh::SaveVTKHDF, or GridFunction::SaveVTKHDF instead.
Definition at line 36 of file vtkhdf.hpp.
| mfem::VTKHDF::VTKHDF | ( | const std::string & | filename, |
| Restart | restart = Restart::Disabled() ) |
Create a new VTKHDF file for serial I/O.
If restart is enabled, then the file (if it exists) will be opened, and time steps before the given time will be preserved. If restart is not enabled, the file (if it exists) will be deleted, and a new file will be created.
Definition at line 460 of file vtkhdf.cpp.
| mfem::VTKHDF::VTKHDF | ( | const std::string & | filename, |
| MPI_Comm | comm_, | ||
| Restart | restart = Restart::Disabled() ) |
Create a new VTKHDF file for parallel I/O.
If restart is enabled, then the file (if it exists) will be opened, and time steps before the given time will be preserved. If restart is not enabled, the file (if it exists) will be deleted, and a new file will be created.
If using restart mode, the file must have previously been saved with the same number of MPI ranks.
Definition at line 482 of file vtkhdf.cpp.
|
delete |
|
delete |
| mfem::VTKHDF::~VTKHDF | ( | ) |
Definition at line 794 of file vtkhdf.cpp.
|
inline |
Disable zlib compression.
Definition at line 312 of file vtkhdf.hpp.
|
inline |
Enable zlib compression at the specified level.
level must be between 0 and 9, in increasing order of compression.
Definition at line 317 of file vtkhdf.hpp.
| void mfem::VTKHDF::Flush | ( | ) |
| template void mfem::VTKHDF::SaveGridFunction< double > | ( | const GridFunction & | gf, |
| const std::string & | name ) |
Save the grid function with the given name, appending as a new time step.
Definition at line 753 of file vtkhdf.cpp.
| template void mfem::VTKHDF::SaveMesh< double > | ( | const Mesh & | mesh, |
| bool | high_order = true, | ||
| int | ref = -1 ) |
Save the mesh, appending as a new time step.
If high_order is true, ref determines the polynomial degree of the mesh elements (-1 indicates the same as the mesh nodes). If high_order is false, the elements are uniformly subdivided according to ref.
Definition at line 534 of file vtkhdf.cpp.
| void mfem::VTKHDF::UpdateSteps | ( | real_t | t | ) |
Update the time step data after saving the mesh and grid functions.
Definition at line 507 of file vtkhdf.cpp.