12 #ifndef MFEM_DATACOLLECTION 13 #define MFEM_DATACOLLECTION 15 #include "../config/config.hpp" 33 typedef std::map<std::string, T*>
MapType;
40 void Register(
const std::string& fname, T* field,
bool own_data)
80 bool Has(
const std::string& fname)
const 87 T*
Get(
const std::string& fname)
const 90 return it !=
field_map.end() ? it->second : NULL;
315 virtual void SetMesh(MPI_Comm comm,
Mesh *new_mesh);
370 virtual void SaveField(
const std::string &field_name);
372 virtual void SaveQField(
const std::string &q_field_name);
375 virtual void Load(
int cycle_ = 0);
461 virtual void SetMesh(MPI_Comm comm,
Mesh *new_mesh)
override;
492 virtual void Save()
override;
498 virtual void Load(
int cycle_ = 0)
override;
509 int levels_of_detail;
510 int compression_level;
511 std::fstream pvd_stream;
513 bool high_order_output;
549 virtual void Save()
override;
593 virtual void Load(
int cycle_ = 0)
override;
virtual void SaveMesh()
Save the mesh, creating the collection directory.
void SaveOneField(const FieldMapIterator &it)
Save one field to disk, assuming the collection directory exists.
void SetPrecision(int prec)
Set the precision (number of digits) used for the text output of doubles.
bool HasQField(const std::string &q_field_name) const
Check if a QuadratureFunction with the given name is in the collection.
virtual void SetPadDigitsRank(int digits)
Set the number of digits used for the MPI rank in filenames.
void WritePVTUHeader(std::ostream &out)
T * Get(const std::string &fname) const
Get a pointer to the field associated with name fname.
Class for grid function - Vector with associated FE space.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
int format
Output mesh format: see the Format enumeration.
void SetDataFormat(VTKFormat fmt)
int Error() const
Get the current error state.
bool appendRankToFileName
Append rank to any output file names.
iterator begin()
Returns a begin iterator to the registered fields.
int GetCompressionLevel() const
If compression is enabled, return the compression level, otherwise return 0.
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
void DeleteAll()
Delete data owned by the DataCollection including field information.
Helper class for ParaView visualization data.
iterator end()
Returns an end iterator to the registered fields.
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
GFieldMap::iterator FieldMapIterator
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
Format
Format constants to be used with SetFormat().
const std::string & GetPrefixPath() const
Get the path where the DataCollection will be saved.
Helper class for VisIt visualization data.
const char * GetDataTypeString() const
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.
QuadratureFunction * GetQField(const std::string &q_field_name)
Get a pointer to a QuadratureFunction in the collection.
QFieldMap::MapType QFieldMapType
const QFieldMapType & GetQFieldMap() const
Get a const reference to the internal q-field map.
bool own_data
Should the collection delete its mesh and fields.
VisItFieldInfo(std::string association_, int num_components_, int lod_=1)
virtual void SetCompression(bool comp)
Set the flag for use of gz compressed files.
int GetCycle() const
Get time cycle (for time-dependent simulations)
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
std::string GetMeshShortFileName() const
std::string GenerateCollectionPath()
virtual void SaveQField(const std::string &q_field_name)
Save one q-field, assuming the collection directory already exists.
MPI_Comm m_comm
Associated MPI communicator.
virtual void Load(int cycle_=0) override
Load the collection - not implemented in the ParaView writer.
double time
Physical time (for time-dependent simulations)
MapType::iterator iterator
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
iterator find(const std::string &fname)
Returns an iterator to the field fname.
void SetCompressionLevel(int compression_level_)
Set the zlib compression level.
MapType::const_iterator const_iterator
Mesh * mesh
The (common) mesh for the collected fields.
void SaveOneQField(const QFieldMapIterator &it)
Save one q-field to disk, assuming the collection directory exists.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
bool serial
Serial or parallel run? False iff mesh is a ParMesh.
int cycle
Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
VisItDataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
virtual void Save()
Save the collection to disk.
const_iterator end() const
Returns an end const iterator to the registered fields.
std::string GenerateVTUFileName(const std::string &prefix, int rank)
Data collection with VisIt I/O routines.
void SaveDataVTU(std::ostream &out, int ref)
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
int visit_max_levels_of_detail
std::string GetFieldFileName(const std::string &field_name) const
VTKFormat
Data array format for VTK and VTU files.
void SaveGFieldVTU(std::ostream &out, int ref_, const FieldMapIterator &it)
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
const MapType & GetMap() const
Returns a const reference to the underlying map.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
void WritePVTUFooter(std::ostream &out, const std::string &vtu_prefix)
void Register(const std::string &fname, T *field, bool own_data)
Register field field with name fname.
DataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Initialize the collection with its name and Mesh.
virtual void SetMesh(Mesh *new_mesh) override
Set/change the mesh associated with the collection.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
virtual void Load(int cycle_=0) override
Load the collection based on its VisIt data (described in its root file)
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
void SetHighOrderOutput(bool high_order_output_)
static const int pad_digits_default
Default value for pad_digits_*.
void SetCompression(bool compression_) override
void ResetError(int err_state=No_Error)
Reset the error state.
void SetTime(double t)
Set physical time (for time-dependent simulations)
QFieldMap::const_iterator QFieldMapConstIterator
const FieldMapType & GetFieldMap() const
Get a const reference to the internal field map.
std::string GeneratePVTUFileName(const std::string &prefix)
void Deregister(const std::string &fname, bool own_data)
Unregister association between field field and name fname.
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf) override
Add a quadrature function to the collection and update the root file.
void SetOwnData(bool o)
Set the ownership of collection data.
virtual void SetPadDigits(int digits)
Set the number of digits used for both the cycle and the MPI rank.
virtual ~VisItDataCollection()
We will delete the mesh and fields if we own them.
bool IsBinaryFormat() const
Returns true if the output format is BINARY or BINARY32, false if ASCII.
bool HasField(const std::string &field_name) const
Check if a grid function is part of the collection.
std::map< std::string, VisItFieldInfo > field_info_map
ParGridFunction * GetParField(const std::string &field_name)
Get a pointer to a parallel grid function in the collection.
void UseRestartMode(bool restart_mode_)
int NumFields() const
Returns the number of registered fields.
QFieldMap::iterator QFieldMapIterator
int precision
Precision (number of digits) used for the text output of doubles.
bool Has(const std::string &fname) const
Predicate to check if a field is associated with name fname.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
virtual void Save() override
Save the collection and a VisIt root file.
int myid
MPI rank (in parallel)
MPI_Comm GetComm() const
Return the associated MPI communicator or MPI_COMM_NULL.
virtual void SaveField(const std::string &field_name)
Save one field, assuming the collection directory already exists.
std::string GenerateVTUPath()
static const int precision_default
Default value for precision.
const_iterator find(const std::string &fname) const
Returns a const iterator to the field fname.
GFieldMap::MapType FieldMapType
void SetMaxLevelsOfDetail(int max_levels_of_detail)
Set VisIt parameter: maximum levels of detail for the MultiresControl.
virtual ~DataCollection()
Delete the mesh and fields if owned by the collection.
void SetTimeStep(double ts)
Set the simulation time step (for time-dependent simulations)
GFieldMap::const_iterator FieldMapConstIterator
std::string GeneratePVTUPath()
void clear()
Clears the map of registered fields without reclaiming memory.
void SaveRootFile()
Save a VisIt root file for the collection.
virtual void DeregisterQField(const std::string &field_name)
Remove a QuadratureFunction from the collection.
virtual void SetPadDigitsCycle(int digits)
Set the number of digits used for the cycle.
virtual void SetPadDigits(int digits) override
ParaViewDataCollection(const std::string &collection_name, mfem::Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has '/' at the end...
std::string GeneratePVDFileName()
Lightweight adaptor over an std::map from strings to pointer to T.
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
void SetLevelsOfDetail(int levels_of_detail_)
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
int visit_levels_of_detail
void LoadVisItRootFile(const std::string &root_name)
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
std::string GetMeshFileName() const
virtual void Save() override
std::string name
Name of the collection, used as a directory name when saving.
void DeleteAll()
Delete all data owned by VisItDataCollection including field data information.
void DeleteData(bool own_data)
Clear all associations between names and fields.
Class for parallel grid function.
double time_step
Time step i.e. delta_t (for time-dependent simulations)
int num_procs
Number of MPI ranks (in parallel)
Class for parallel meshes.
const std::string & GetCollectionName() const
Get the name of the collection.
const char * GetDataFormatString() const
Represents values or vectors of values at quadrature points on a mesh.
std::map< std::string, T * > MapType
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
const_iterator begin() const
Returns a begin const iterator to the registered fields.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
double GetTime() const
Get physical time (for time-dependent simulations)