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