13 #include "../general/text.hpp"
24 #define mkdir(dir, mode) _mkdir(dir)
32 const Mesh *mesh,
int myid)
35 const char path_delim =
'/';
36 std::string::size_type pos = 0;
44 pos = dir_name.find(path_delim, pos+1);
45 std::string subdir = dir_name.substr(0, pos);
48 err = mkdir(subdir.c_str(), 0777);
49 err = (err && (errno != EEXIST)) ? 1 : 0;
51 if (myid == 0 || pmesh == NULL)
53 err = mkdir(subdir.c_str(), 0777);
54 err = (err && (errno != EEXIST)) ? 1 : 0;
58 while ( pos != std::string::npos );
63 MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->
GetComm());
74 name = collection_name;
174 return (it !=
field_map.end()) ? it->second : NULL;
181 return (it !=
q_field_map.end()) ? it->second : NULL;
202 MFEM_ABORT(
"this method is not implemented");
209 if (
error) {
return; }
237 MFEM_WARNING(
"Error creating directory: " << dir_name);
241 std::string mesh_name = dir_name +
247 std::ofstream mesh_file(mesh_name.c_str());
251 if (pmesh &&
format == 1 )
263 MFEM_WARNING(
"Error writing mesh to file: " << mesh_name);
274 std::string file_name = dir_name +
"/" + field_name;
286 (it->second)->
Save(field_file);
290 MFEM_WARNING(
"Error writing field to file: " << it->first);
298 (it->second)->
Save(q_field_file);
302 MFEM_WARNING(
"Error writing q-field to file: " << it->first);
331 if (
own_data) {
delete it->second; }
337 if (
own_data) {
delete it->second; }
412 if (
myid != 0) {
return; }
416 std::ofstream root_file(root_name.c_str());
421 MFEM_WARNING(
"Error writting VisIt Root file: " << root_name);
452 std::ifstream root_file(root_name.c_str());
453 std::stringstream buffer;
454 buffer << root_file.rdbuf();
458 MFEM_WARNING(
"Error reading the VisIt Root file: " << root_name);
475 MFEM_WARNING(
"Unable to open mesh file: " << mesh_fname);
495 std::string fname = path_left + it->first + path_right;
496 std::ifstream file(fname.c_str());
500 MFEM_WARNING(
"Unable to open field file: " << fname);
515 picojson::object top, dsets,
main,
mesh, fields, field, mtags, ftags;
522 mesh[
"path"] = picojson::value(path_str + ((
format==0)?
"":
"p") +
"mesh" +
524 mesh[
"tags"] = picojson::value(mtags);
530 ftags[
"assoc"] = picojson::value((it->second).association);
531 ftags[
"comps"] = picojson::value(
to_string((it->second).num_components));
532 field[
"path"] = picojson::value(path_str + it->first + file_ext_format);
533 field[
"tags"] = picojson::value(ftags);
534 fields[it->first] = picojson::value(field);
537 main[
"cycle"] = picojson::value(
double(
cycle));
538 main[
"time"] = picojson::value(
time);
539 main[
"domains"] = picojson::value(
double(
num_procs));
540 main[
"mesh"] = picojson::value(mesh);
543 main[
"fields"] = picojson::value(fields);
546 dsets[
"main"] = picojson::value(main);
547 top[
"dsets"] = picojson::value(dsets);
549 return picojson::value(top).serialize(
true);
554 picojson::value top, dsets,
main,
mesh, fields;
555 std::string parse_err = picojson::parse(top, json);
556 if (!parse_err.empty())
559 MFEM_WARNING(
"Unable to parse visit root data.");
564 dsets = top.get(
"dsets");
565 main = dsets.get(
"main");
566 cycle = int(main.get(
"cycle").get<
double>());
567 time = main.get(
"time").get<
double>();
568 num_procs = int(main.get(
"domains").get<
double>());
569 mesh = main.get(
"mesh");
570 fields = main.get(
"fields");
575 std::string path = mesh.get(
"path").get<std::string>();
576 size_t right_sep = path.find(
'_');
577 if (right_sep == std::string::npos)
580 MFEM_WARNING(
"Unable to parse visit root data.");
583 name = path.substr(0, right_sep);
586 topo_dim =
to_int(mesh.get(
"tags").get(
"topo_dim").get<std::string>());
588 to_int(mesh.get(
"tags").get(
"max_lods").get<std::string>());
592 if (fields.is<picojson::object>())
594 picojson::object fields_obj = fields.get<picojson::object>();
595 for (picojson::object::iterator it = fields_obj.begin();
596 it != fields_obj.end(); ++it)
598 picojson::value tags = it->second.get(
"tags");
601 to_int(tags.get(
"comps").get<std::string>()));
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.
Class for grid function - Vector with associated FE space.
int format
Output mesh format: 0 - serial format (default), 1 - parallel format.
bool appendRankToFileName
Append rank to any output file names.
void DeleteAll()
Delete data owned by the DataCollection including field information.
std::string to_string(int i)
Helper class for VisIt visualization data.
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.
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.
int main(int argc, char *argv[])
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 time
Physical time (for time-dependent simulations)
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
FieldMapType::const_iterator FieldMapConstIterator
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.
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.
std::string to_padded_string(int i, int digits)
int to_int(const std::string &str)
std::string GetFieldFileName(const std::string &field_name)
int visit_max_levels_of_detail
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.
virtual void Load(int cycle_=0)
Load the collection based on its VisIt data (described in its root file)
int SpaceDimension() const
std::map< std::string, VisItFieldInfo > field_info_map
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.
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.
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.
double time_step
Time step i.e. delta_t (for time-dependent simulations)
int num_procs
Number of MPI ranks (in parallel)
void ParPrint(std::ostream &out) const
Save the mesh in a parallel mesh format.
Class for parallel meshes.
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.
virtual void Print(std::ostream &out=std::cout) const