13 #include "../mesh/nurbs.hpp"
14 #include "../general/text.hpp"
25 #define mkdir(dir, mode) _mkdir(dir)
33 const Mesh *mesh,
int myid)
36 const char path_delim =
'/';
37 std::string::size_type pos = 0;
45 pos = dir_name.find(path_delim, pos+1);
46 std::string subdir = dir_name.substr(0, pos);
49 err = mkdir(subdir.c_str(), 0777);
50 err = (err && (errno != EEXIST)) ? 1 : 0;
52 if (myid == 0 || pmesh == NULL)
54 err = mkdir(subdir.c_str(), 0777);
55 err = (err && (errno != EEXIST)) ? 1 : 0;
59 while ( pos != std::string::npos );
64 MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->
GetComm());
75 std::string::size_type pos = collection_name.find_last_of(
'/');
76 if (pos == std::string::npos)
78 name = collection_name;
84 name = collection_name.substr(pos+1);
147 MPI_Comm_rank(comm, &
myid);
160 default: MFEM_ABORT(
"unknown format: " << fmt);
168 #ifdef MFEM_USE_GZSTREAM
169 MFEM_ASSERT(!
compression,
"GZStream not enabled in MFEM build.");
191 MFEM_ABORT(
"this method is not implemented");
198 if (
error) {
return; }
226 MFEM_WARNING(
"Error creating directory: " << dir_name);
232 ofgzstream mesh_file(mesh_name.c_str(), mode);
248 MFEM_WARNING(
"Error writing mesh to file: " << mesh_name);
270 std::string file_name = dir_name +
"/" + field_name;
284 (it->second)->
Save(field_file);
288 MFEM_WARNING(
"Error writing field to file: " << it->first);
297 (it->second)->
Save(q_field_file);
301 MFEM_WARNING(
"Error writing q-field to file: " << it->first);
382 const std::string& collection_name,
387 MPI_Comm_rank(comm, &
myid);
412 MPI_Comm_rank(comm, &
myid);
463 if (
myid != 0) {
return; }
468 std::ofstream root_file(root_name.c_str());
473 MFEM_WARNING(
"Error writing VisIt root file: " << root_name);
490 MFEM_WARNING(
"Cannot load parallel VisIt root file in serial.");
493 if (
m_comm == MPI_COMM_NULL)
495 MFEM_WARNING(
"Cannot load parallel VisIt root file without MPI"
504 MPI_Comm_size(
m_comm, &comm_size);
507 MFEM_WARNING(
"Processor number mismatch: VisIt root file: "
508 <<
num_procs <<
", MPI_comm: " << comm_size);
534 std::ifstream root_file(root_name.c_str());
535 std::stringstream buffer;
536 buffer << root_file.rdbuf();
540 MFEM_WARNING(
"Error reading the VisIt root file: " << root_name);
556 MFEM_WARNING(
"Unable to open mesh file: " << mesh_fname);
572 MFEM_WARNING(
"Reading parallel format in serial is not supported");
591 std::string fname = path_left + it->first + path_right;
597 MFEM_WARNING(
"Unable to open field file: " << fname);
613 MFEM_WARNING(
"Reading parallel format in serial is not supported");
623 std::string path_str =
627 picojson::object top, dsets,
main,
mesh, fields, field, mtags, ftags;
636 mesh[
"tags"] = picojson::value(mtags);
643 ftags[
"assoc"] = picojson::value((it->second).association);
644 ftags[
"comps"] = picojson::value(
to_string((it->second).num_components));
646 field[
"path"] = picojson::value(path_str + it->first + file_ext_format);
647 field[
"tags"] = picojson::value(ftags);
648 fields[it->first] = picojson::value(field);
651 main[
"cycle"] = picojson::value(
double(
cycle));
652 main[
"time"] = picojson::value(
time);
653 main[
"time_step"] = picojson::value(
time_step);
654 main[
"domains"] = picojson::value(
double(
num_procs));
655 main[
"mesh"] = picojson::value(mesh);
658 main[
"fields"] = picojson::value(fields);
661 dsets[
"main"] = picojson::value(main);
662 top[
"dsets"] = picojson::value(dsets);
664 return picojson::value(top).serialize(
true);
669 picojson::value top, dsets,
main,
mesh, fields;
670 std::string parse_err = picojson::parse(top, json);
671 if (!parse_err.empty())
674 MFEM_WARNING(
"Unable to parse VisIt root data.");
679 dsets = top.get(
"dsets");
680 main = dsets.get(
"main");
681 cycle = int(main.get(
"cycle").get<
double>());
682 time = main.get(
"time").get<
double>();
683 if (main.contains(
"time_step"))
685 time_step = main.get(
"time_step").get<
double>();
687 num_procs = int(main.get(
"domains").get<
double>());
688 mesh = main.get(
"mesh");
689 fields = main.get(
"fields");
694 std::string path = mesh.get(
"path").get<std::string>();
695 size_t right_sep = path.find(
'_');
696 if (right_sep == std::string::npos)
699 MFEM_WARNING(
"Unable to parse VisIt root data.");
702 name = path.substr(0, right_sep);
704 if (mesh.contains(
"format"))
709 topo_dim =
to_int(mesh.get(
"tags").get(
"topo_dim").get<std::string>());
711 to_int(mesh.get(
"tags").get(
"max_lods").get<std::string>());
715 if (fields.is<picojson::object>())
717 picojson::object fields_obj = fields.get<picojson::object>();
718 for (picojson::object::iterator it = fields_obj.begin();
719 it != fields_obj.end(); ++it)
721 picojson::value tags = it->second.get(
"tags");
724 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.
virtual void Print(std::ostream &out=mfem::out) const
Class for grid function - Vector with associated FE space.
int format
Output mesh format: see the Format enumeration.
bool appendRankToFileName
Append rank to any output file names.
iterator begin()
Returns a begin iterator to the registered fields.
void DeleteAll()
Delete data owned by the DataCollection including field information.
std::string to_string(int i)
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
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Helper class for VisIt visualization data.
virtual void Save()
Save the collection and a VisIt root file.
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.
void SetCompression(bool comp)
Set the flag for use of gz compressed files.
int main(int argc, char *argv[])
void ParseVisItRootString(const std::string &json)
Read in a VisIt root file in JSON format.
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.
double time
Physical time (for time-dependent simulations)
std::map< std::string, VisItFieldInfo >::iterator FieldInfoMapIterator
iterator find(const std::string &fname)
Returns an iterator to the field fname.
std::string GetMeshShortFileName() const
Mesh * mesh
The (common) mesh for the collected fields.
int GetNE() const
Returns number of elements in the mesh.
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.
std::string to_padded_string(int i, int digits)
int to_int(const std::string &str)
int visit_max_levels_of_detail
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
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.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
void SetLevelsOfDetail(int levels_of_detail)
Set VisIt parameter: default levels of detail for the MultiresControl.
static const int pad_digits_default
Default value for pad_digits_*.
FiniteElementSpace * FESpace()
virtual void Load(int cycle_=0)
Load the collection based on its VisIt data (described in its root file)
int GetOrder() const
If all orders are identical, return that number. Otherwise, return NURBSFECollection::VariableOrder.
int SpaceDimension() const
std::map< std::string, VisItFieldInfo > field_info_map
const NURBSExtension * GetNURBSext() const
QFieldMap::iterator QFieldMapIterator
int precision
Precision (number of digits) used for the text output of doubles.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
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.
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.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual ~DataCollection()
Delete the mesh and fields if owned by the collection.
std::string GetFieldFileName(const std::string &field_name) const
void clear()
Clears the map of registered fields without reclaiming memory.
void SaveRootFile()
Save a VisIt root file for the collection.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has '/' at the end...
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i'th element.
virtual void Load(int cycle_=0)
Load the collection. Not implemented in the base class DataCollection.
std::string GetMeshFileName() const
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 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.
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)
void ParPrint(std::ostream &out) const
Save the mesh in a parallel mesh format.
Class for parallel meshes.
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.