12 #ifndef MFEM_DATACOLLECTION
13 #define MFEM_DATACOLLECTION
15 #include "../config/config.hpp"
211 virtual void SaveField(
const std::string &field_name);
213 virtual void SaveQField(
const std::string &q_field_name);
216 virtual void Load(
int cycle_ = 0);
289 virtual void Load(
int cycle_ = 0);
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.
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: 0 - serial format (default), 1 - parallel format.
double GetTime() const
Get physical time (for time-dependent simulations)
bool appendRankToFileName
Append rank to any output file names.
void DeleteAll()
Delete data owned by the DataCollection including field information.
bool HasQField(const std::string &q_field_name) const
Check if a QuadratureFunction with the given name is in the collection.
const std::string & GetCollectionName() const
Get the name of the collection.
std::map< std::string, GridFunction * > FieldMapType
Helper class for VisIt visualization data.
bool HasField(const std::string &name) const
Check if a grid function is part of the collection.
virtual void Save()
Save the collection and a VisIt root file.
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.
int Error() const
Get the current error state.
bool own_data
Should the collection delete its mesh and fields.
QFieldMapType q_field_map
QuadratureFunction * GetQField(const std::string &q_field_name)
Get a pointer to a QuadratureFunction in the collection.
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
QFieldMapType::iterator QFieldMapIterator
virtual void SaveQField(const std::string &q_field_name)
Save one q-field, assuming the collection directory already exists.
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
double time
Physical time (for time-dependent simulations)
void ResetError(int err=NO_ERROR)
Reset the error state.
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
FieldMapType::const_iterator FieldMapConstIterator
Mesh * mesh
The (common) mesh for the collected fields.
const QFieldMapType & GetQFieldMap() const
Get a const reference to the internal q-field map.
void SaveOneQField(const QFieldMapIterator &it)
Save one q-field to disk, assuming the collection directory exists.
void SetPadDigits(int digits)
Set the number of digits used for the cycle and MPI rank in filenames.
VisItFieldInfo(std::string _association, int _num_components)
bool serial
Serial or parallel run?
VisItDataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
const FieldMapType & GetFieldMap() const
Get a const reference to the internal field map.
int GetCycle() const
Get time cycle (for time-dependent simulations)
Data collection with VisIt I/O routines.
std::string GetFieldFileName(const std::string &field_name)
int visit_max_levels_of_detail
void SetFormat(int fmt)
Set the desired output mesh format: 0 - serial format (default), 1 - parallel format.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
DataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Initialize the collection with its name and Mesh.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
static const int pad_digits_default
Default value for pad_digits.
void SetTime(double t)
Set physical time (for time-dependent simulations)
const std::string & GetPrefixPath() const
Get the path where the DataCollection will be saved.
virtual void Load(int cycle_=0)
Load the collection based on its VisIt data (described in its root file)
void SetOwnData(bool o)
Set the ownership of collection data.
virtual ~VisItDataCollection()
We will delete the mesh and fields if we own them.
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.
FieldMapType::iterator FieldMapIterator
The fields and their names (used when saving)
int precision
Precision (number of digits) used for the text output of doubles.
int myid
MPI rank (in parallel)
virtual void SaveField(const std::string &field_name)
Save one field, assuming the collection directory already exists.
static const int precision_default
Default value for precision.
QFieldMapType::const_iterator QFieldMapConstIterator
void SetMaxLevelsOfDetail(int max_levels_of_detail)
Set VisIt parameter: maximum levels of detail for the MultiresControl.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
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)
virtual void RegisterQField(const std::string &q_field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
virtual void DeregisterQField(const std::string &field_name)
Remove a QuadratureFunction from the collection.
void SaveRootFile()
Save a VisIt root file for the collection.
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
int pad_digits
Number of digits used for the cycle and MPI rank in filenames.
void LoadVisItRootFile(const std::string &root_name)
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
std::string name
Name of the collection, used as a directory name when saving.
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
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 representing a function through its values (scalar or vector) at quadrature points...
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
std::map< std::string, QuadratureFunction * > QFieldMapType