21 #define mkdir(dir, mode) _mkdir(dir)
37 string out_str = ss.str();
38 out_str = out_str.substr(out_str.find_first_not_of(
" \t"));
45 oss << setw(digits) << setfill(
'0') << i;
52 stringstream(str) >> i;
60 name = collection_name;
68 pad_digits = pad_digits_default;
74 name = collection_name;
91 pad_digits = pad_digits_default;
102 if (HasField(field_name))
103 return field_map[field_name];
117 err = mkdir(dir_name.c_str(), 0777);
118 err = (err && (errno != EEXIST)) ? 1 : 0;
121 if (myid == 0 || pmesh == NULL)
123 err = mkdir(dir_name.c_str(), 0777);
124 err = (err && (errno != EEXIST)) ? 1 : 0;
126 MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->
GetComm());
131 MPI_Bcast(&err, 1, MPI_INT, 0, pmesh->
GetComm());
137 MFEM_WARNING(
"Error creating directory: " << dir_name);
143 mesh_name = dir_name +
"/mesh";
146 ofstream mesh_file(mesh_name.c_str());
147 mesh->Print(mesh_file);
151 MFEM_WARNING(
"Error writing mesh to file: " << mesh_name);
156 for (map<string,GridFunction*>::iterator it = field_map.begin();
157 it != field_map.end(); ++it)
160 field_name = dir_name +
"/" + it->first;
162 field_name = dir_name +
"/" + it->first +
"." +
164 ofstream field_file(field_name.c_str());
165 (it->second)->Save(field_file);
169 MFEM_WARNING(
"Error writting field to file: " << field_name);
179 for (map<string,GridFunction*>::iterator it = field_map.begin();
180 it != field_map.end(); ++it)
200 for (map<string,GridFunction*>::iterator it = field_map.begin();
201 it != field_map.end(); ++it)
255 ofstream root_file(root_name.c_str());
260 MFEM_WARNING(
"Error writting VisIt Root file: " << root_name);
285 ifstream root_file(root_name.c_str());
287 buffer << root_file.rdbuf();
291 MFEM_WARNING(
"Error reading the VisIt Root file: " << root_name);
303 ifstream file(mesh_fname.c_str());
307 MFEM_WARNING(
"Unable to open mesh file: " << mesh_fname);
323 for (map<string,VisItFieldInfo>::iterator it =
field_info_map.begin();
326 string fname = path_left + it->first + path_right;
327 ifstream file(fname.c_str());
331 MFEM_WARNING(
"Unable to open field file: " << fname);
346 picojson::object top, dsets,
main,
mesh, fields, field, mtags, ftags;
353 mesh[
"path"] = picojson::value(path_str +
"mesh" + file_ext_format);
354 mesh[
"tags"] = picojson::value(mtags);
357 for (map<string,VisItFieldInfo>::iterator it =
field_info_map.begin();
360 ftags[
"assoc"] = picojson::value((it->second).association);
361 ftags[
"comps"] = picojson::value(
to_string((it->second).num_components));
362 field[
"path"] = picojson::value(path_str + it->first + file_ext_format);
363 field[
"tags"] = picojson::value(ftags);
364 fields[it->first] = picojson::value(field);
367 main[
"cycle"] = picojson::value(
double(
cycle));
368 main[
"time"] = picojson::value(
time);
369 main[
"domains"] = picojson::value(
double(
num_procs));
370 main[
"mesh"] = picojson::value(mesh);
372 main[
"fields"] = picojson::value(fields);
374 dsets[
"main"] = picojson::value(main);
375 top[
"dsets"] = picojson::value(dsets);
377 return picojson::value(top).serialize(
true);
382 picojson::value top, dsets,
main,
mesh, fields;
383 string parse_err = picojson::parse(top, json);
384 if (!parse_err.empty())
387 MFEM_WARNING(
"Unable to parse visit root data.");
392 dsets = top.get(
"dsets");
393 main = dsets.get(
"main");
394 cycle = int(main.get(
"cycle").get<
double>());
395 time = main.get(
"time").get<
double>();
396 num_procs = int(main.get(
"domains").get<
double>());
397 mesh = main.get(
"mesh");
398 fields = main.get(
"fields");
403 string path = mesh.get(
"path").get<
string>();
404 size_t right_sep = path.find(
'_');
405 if (right_sep == string::npos)
408 MFEM_WARNING(
"Unable to parse visit root data.");
411 name = path.substr(0, right_sep);
414 topo_dim =
to_int(mesh.get(
"tags").get(
"topo_dim").get<
string>());
416 to_int(mesh.get(
"tags").get(
"max_lods").get<
string>());
420 if (fields.is<picojson::object>())
422 picojson::object fields_obj = fields.get<picojson::object>();
423 for (picojson::object::iterator it = fields_obj.begin();
424 it != fields_obj.end(); ++it)
426 picojson::value tags = it->second.get(
"tags");
429 to_int(tags.get(
"comps").get<
string>()));
std::map< std::string, GridFunction * > field_map
The fields and their names (used when saving)
Class for grid function - Vector with associated FE space.
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
string to_padded_string(int i, int digits)
void DeleteAll()
Delete data owned by the DataCollection including field information.
Helper class for VisIt visualization data.
virtual void Save()
Save the collection and a VisIt root file.
bool own_data
Should the collection delete its mesh and fields.
double time
Physical time (for time-dependent simulations)
Mesh * mesh
The (common) mesh for the collected fields.
bool serial
Serial or parallel run? If false, append rank (myid) to file names.
void ParseVisItRootString(std::string json)
Read in a VisIt root file in JSON format.
int visit_max_levels_of_detail
GridFunction * GetField(const char *field_name)
Get a pointer to a grid function in the collection.
void DeleteData()
Delete data owned by the DataCollection keeping field information.
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
void LoadVisItRootFile(std::string root_name)
int main(int argc, char *argv[])
int myid
MPI rank (in parallel)
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.
VisItDataCollection(const char *collection_name)
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection.
DataCollection(const char *collection_name)
Create an empty collection with the given name.
int pad_digits
Number of digits used for the cycle and MPI rank in filenames.
std::string name
Name of the collection, used as a directory name when saving.
int num_procs
Number of MPI ranks (in parallel)
Class for parallel meshes.
std::string GetVisItRootString()
Prepare the VisIt root file in JSON format for the current collection.