28#define mkdir(dir, mode) _mkdir(dir)
36 const Mesh *mesh,
int myid)
39 const char path_delim =
'/';
40 std::string::size_type pos = 0;
48 pos = dir_name.find(path_delim, pos+1);
49 std::string subdir = dir_name.substr(0, pos);
52 err_flag = mkdir(subdir.c_str(), 0777);
53 err_flag = (err_flag && (errno != EEXIST)) ? 1 : 0;
55 if (
myid == 0 || pmesh == NULL)
57 err_flag = mkdir(subdir.c_str(), 0777);
58 err_flag = (err_flag && (errno != EEXIST)) ? 1 : 0;
62 while ( pos != std::string::npos );
67 MPI_Bcast(&err_flag, 1, MPI_INT, 0, pmesh->
GetComm());
78 std::string::size_type pos = collection_name.find_last_of(
'/');
79 if (pos == std::string::npos)
81 name = collection_name;
87 name = collection_name.substr(pos+1);
150 MPI_Comm_rank(comm, &
myid);
163 default: MFEM_ABORT(
"unknown format: " << fmt);
172 MFEM_VERIFY(!
compression,
"ZLib not enabled in MFEM build.");
194 MFEM_ABORT(
"this method is not implemented");
201 if (
error) {
return; }
227 MFEM_WARNING(
"Error creating directory: " << dir_name);
248 MFEM_WARNING(
"Error writing mesh to file: " << mesh_name);
270 std::string file_name = dir_name +
"/" + field_name;
283 (it->second)->
Save(field_file);
287 MFEM_WARNING(
"Error writing field to file: " << it->first);
296 (it->second)->
Save(q_field_file);
300 MFEM_WARNING(
"Error writing q-field to file: " << it->first);
381 const std::string& collection_name,
386 MPI_Comm_rank(comm, &
myid);
411 MPI_Comm_rank(comm, &
myid);
426 for (
int e=0; e<gf->
FESpace()->GetNE(); e++)
444 for (
int e=0; e<qf->
GetSpace()->GetNE(); e++)
450 LOD = std::max(LOD,locLOD);
458 std::ostringstream oss;
459 oss <<
"QF_" << qf_order <<
"_" << qf_vdim;
461 oss.str(), qf_order);
489 if (
myid != 0) {
return; }
494 std::ofstream root_file(root_name);
499 MFEM_WARNING(
"Error writing VisIt root file: " << root_name);
516 MFEM_WARNING(
"Cannot load parallel VisIt root file in serial.");
519 if (
m_comm == MPI_COMM_NULL)
521 MFEM_WARNING(
"Cannot load parallel VisIt root file without MPI"
530 MPI_Comm_size(
m_comm, &comm_size);
533 MFEM_WARNING(
"Processor number mismatch: VisIt root file: "
534 <<
num_procs <<
", MPI_comm: " << comm_size);
560 std::ifstream root_file(root_name);
561 std::stringstream buffer;
562 buffer << root_file.rdbuf();
566 MFEM_WARNING(
"Error reading the VisIt root file: " << root_name);
584 MFEM_WARNING(
"Unable to open mesh file: " << mesh_fname);
600 MFEM_WARNING(
"Reading parallel format in serial is not supported");
619 std::string fname = path_left + it->first + path_right;
625 MFEM_WARNING(
"Unable to open field file: " << fname);
631 if ((it->second).association ==
"nodes")
635 else if ((it->second).association ==
"elements" ||
636 (it->second).association ==
"quadrature")
644 if ((it->second).association ==
"nodes")
650 else if ((it->second).association ==
"elements" ||
651 (it->second).association ==
"quadrature")
657 MFEM_WARNING(
"Reading parallel format in serial is not supported");
667 std::string path_str =
671 picojson::object top, dsets,
main,
mesh, fields, field, mtags, ftags;
674 std::string file_ext_format =
".%0" + to_string(
pad_digits_rank) +
"d";
675 mtags[
"spatial_dim"] = picojson::value(to_string(
spatial_dim));
676 mtags[
"topo_dim"] = picojson::value(to_string(
topo_dim));
680 mesh[
"tags"] = picojson::value(mtags);
681 mesh[
"format"] = picojson::value(to_string(
format));
687 ftags[
"assoc"] = picojson::value((it->second).association);
688 ftags[
"comps"] = picojson::value(to_string((it->second).num_components));
689 ftags[
"lod"] = picojson::value(to_string((it->second).lod));
690 ftags[
"basis"] = picojson::value((it->second).basis);
691 ftags[
"order"] = picojson::value(to_string((it->second).order));
692 field[
"path"] = picojson::value(path_str + it->first + file_ext_format);
693 field[
"tags"] = picojson::value(ftags);
694 fields[it->first] = picojson::value(field);
697 main[
"cycle"] = picojson::value(
double(
cycle));
698 main[
"time"] = picojson::value(
time);
701 main[
"mesh"] = picojson::value(
mesh);
704 main[
"fields"] = picojson::value(fields);
707 dsets[
"main"] = picojson::value(
main);
708 top[
"dsets"] = picojson::value(dsets);
710 return picojson::value(top).serialize(
true);
715 picojson::value top, dsets,
main,
mesh, fields;
716 std::string parse_err = picojson::parse(top, json);
717 if (!parse_err.empty())
720 MFEM_WARNING(
"Unable to parse VisIt root data.");
725 dsets = top.get(
"dsets");
726 main = dsets.get(
"main");
727 cycle = int(
main.get(
"cycle").get<
double>());
728 time =
main.get(
"time").get<
double>();
729 if (
main.contains(
"time_step"))
735 fields =
main.get(
"fields");
740 std::string path =
mesh.get(
"path").get<std::string>();
741 size_t right_sep = path.rfind(
'_');
742 if (right_sep == std::string::npos)
745 MFEM_WARNING(
"Unable to parse VisIt root data.");
748 name = path.substr(0, right_sep);
750 if (
mesh.contains(
"format"))
757 to_int(
mesh.get(
"tags").get(
"max_lods").get<std::string>());
761 if (fields.is<picojson::object>())
763 picojson::object fields_obj = fields.get<picojson::object>();
764 for (picojson::object::iterator it = fields_obj.begin();
765 it != fields_obj.end(); ++it)
767 picojson::value tags = it->second.get(
"tags");
771 std::string basis =
"";
774 if (tags.contains(
"lod"))
776 lod =
to_int(tags.get(
"lod").get<std::string>());
779 if (tags.contains(
"basis"))
781 basis = tags.get(
"comps").get<std::string>();
784 if (tags.contains(
"order"))
786 order =
to_int(tags.get(
"comps").get<std::string>());
791 to_int(tags.get(
"comps").get<std::string>()),
825 MFEM_ASSERT(compression_level_ >= -1 && compression_level_ <= 9,
826 "Compression level must be between -1 and 9 (inclusive).");
852 const std::string& collection_name,
Mesh *mesh_)
876 const std::string &prefix)
878 return prefix +
".pvtu";
882 const std::string &prefix,
int rank)
899 MFEM_WARNING(
"Error creating directory: " << path);
909 if (
myid == 0 && !pvd_stream.is_open())
913 bool write_header =
true;
914 std::ifstream pvd_in;
915 if (
restart_mode && (pvd_in.open(pvdname,std::ios::binary),pvd_in.good()))
919 std::fstream::pos_type pos_begin = pvd_in.tellg();
920 std::fstream::pos_type pos_end = pos_begin;
922 std::regex regexp(
"timestep=\"([^[:space:]]+)\".*file=\"Cycle(\\d+)");
926 while (getline(pvd_in,line))
928 if (regex_search(line,match,regexp))
930 MFEM_ASSERT(match.size() == 3,
"Unable to parse DataSet");
931 double tvalue = std::stod(match[1]);
932 if (tvalue >=
GetTime()) {
break; }
933 int cvalue = std::stoi(match[2]);
935 " is too small for restart mode: trying to overwrite"
937 pos_end = pvd_in.tellg();
943 size_t count = pos_end - pos_begin;
946 write_header =
false;
947 std::vector<char> buf(count);
951 pvd_in.seekg(pos_begin);
952 pvd_in.read(buf.data(), count);
957 pvd_stream.open(pvdname,std::ios::out|std::ios::trunc|std::ios::binary);
958 pvd_stream.write(buf.data(), count);
961 pvd_stream.open(pvdname,std::ios::in|std::ios::out|std::ios::ate);
967 pvd_stream.open(pvdname,std::ios::out|std::ios::trunc);
968 pvd_stream <<
"<?xml version=\"1.0\"?>\n";
969 pvd_stream <<
"<VTKFile type=\"Collection\" version=\"2.2\"";
970 pvd_stream <<
" byte_order=\"" <<
VTKByteOrder() <<
"\">\n";
971 pvd_stream <<
"<Collection>" << std::endl;
989 "QuadratureFunction output is not supported for "
990 "ParaViewDataCollection on domain boundary!");
991 const std::string &field_name = qfield.first;
1007 pvtu_out <<
"<PPointData>\n";
1010 int vec_dim = field_it.second->VectorDim();
1012 <<
"\" Name=\"" << field_it.first
1013 <<
"\" NumberOfComponents=\"" << vec_dim <<
"\" "
1017 for (
auto &field_it : coeff_field_map)
1021 <<
"\" Name=\"" << field_it.first
1022 <<
"\" NumberOfComponents=\"" << vec_dim <<
"\" "
1025 for (
auto &field_it : vcoeff_field_map)
1027 int vec_dim = field_it.second->GetVDim();
1029 <<
"\" Name=\"" << field_it.first
1030 <<
"\" NumberOfComponents=\"" << vec_dim <<
"\" "
1033 pvtu_out <<
"</PPointData>\n";
1036 pvtu_out <<
"<PCellData>\n";
1037 pvtu_out <<
"\t<PDataArray type=\"Int32\" Name=\"" <<
"attribute"
1038 <<
"\" NumberOfComponents=\"1\""
1040 pvtu_out <<
"</PCellData>\n";
1046 pvd_stream <<
"<DataSet timestep=\"" <<
GetTime()
1047 <<
"\" group=\"\" part=\"" << 0 <<
"\" file=\""
1049 <<
"\" name=\"mesh\"/>\n";
1055 const std::string &q_field_name = q_field.first;
1059 std::ofstream pvtu_out(col_path +
"/" + q_fname);
1061 int vec_dim = q_field.second->GetVDim();
1062 pvtu_out <<
"<PPointData>\n";
1064 <<
"\" Name=\"" << q_field_name
1065 <<
"\" NumberOfComponents=\"" << vec_dim <<
"\" "
1068 pvtu_out <<
"</PPointData>\n";
1071 pvd_stream <<
"<DataSet timestep=\"" <<
GetTime()
1072 <<
"\" group=\"\" part=\"" << 0 <<
"\" file=\""
1073 << q_fname <<
"\" name=\"" << q_field_name <<
"\"/>\n";
1078 std::fstream::pos_type pos = pvd_stream.tellp();
1079 pvd_stream <<
"</Collection>\n";
1080 pvd_stream <<
"</VTKFile>" << std::endl;
1081 pvd_stream.seekp(pos);
1087 os <<
"<?xml version=\"1.0\"?>\n";
1088 os <<
"<VTKFile type=\"PUnstructuredGrid\"";
1089 os <<
" version =\"2.2\" byte_order=\"" <<
VTKByteOrder() <<
"\">\n";
1090 os <<
"<PUnstructuredGrid GhostLevel=\"0\">\n";
1092 os <<
"<PPoints>\n";
1094 os <<
" Name=\"Points\" NumberOfComponents=\"3\""
1096 os <<
"</PPoints>\n";
1099 os <<
"\t<PDataArray type=\"Int32\" ";
1100 os <<
" Name=\"connectivity\" NumberOfComponents=\"1\""
1102 os <<
"\t<PDataArray type=\"Int32\" ";
1103 os <<
" Name=\"offsets\" NumberOfComponents=\"1\""
1105 os <<
"\t<PDataArray type=\"UInt8\" ";
1106 os <<
" Name=\"types\" NumberOfComponents=\"1\""
1108 os <<
"</PCells>\n";
1112 const std::string &vtu_prefix)
1117 os <<
"<Piece Source=\"" << vtu_filename <<
"\"/>\n";
1119 os <<
"</PUnstructuredGrid>\n";
1120 os <<
"</VTKFile>\n";
1125 os <<
"<VTKFile type=\"UnstructuredGrid\"";
1128 os <<
" compressor=\"vtkZLibDataCompressor\"";
1130 os <<
" version=\"2.2\" byte_order=\"" <<
VTKByteOrder() <<
"\">\n";
1131 os <<
"<UnstructuredGrid>\n";
1136 os <<
"<PointData >\n";
1142 "GridFunction output is not supported for "
1143 "ParaViewDataCollection on domain boundary!");
1148 for (
const auto &kv : coeff_field_map)
1152 for (
const auto &kv : vcoeff_field_map)
1156 os <<
"</PointData>\n";
1159 os <<
"</UnstructuredGrid>\n";
1160 os <<
"</VTKFile>" << std::endl;
1169 std::vector<char> buf;
1170 int vec_dim = it->second->VectorDim();
1172 <<
"\" Name=\"" << it->first
1173 <<
"\" NumberOfComponents=\"" << vec_dim <<
"\" "
1182 it->second->GetValues(i, RefG->
RefPts, val, pmat);
1183 for (
int j = 0; j < val.
Size(); j++)
1196 it->second->GetVectorValues(i, RefG->
RefPts, vval, pmat);
1197 for (
int jj = 0; jj < vval.
Width(); jj++)
1199 for (
int ii = 0; ii < vval.
Height(); ii++)
1211 os <<
"</DataArray>" << std::endl;
1219 std::vector<char> buf;
1222 <<
"\" Name=\"" <<
name
1223 <<
"\" NumberOfComponents=\"" << vec_dim <<
"\""
1240 val = coeff.
Eval(*eltrans, ip);
1258 val = coeff.
Eval(*eltrans, ip);
1268 os <<
"</DataArray>" << std::endl;
1276 std::vector<char> buf;
1277 int vec_dim = coeff.
GetVDim();
1279 <<
"\" Name=\"" <<
name
1280 <<
"\" NumberOfComponents=\"" << vec_dim <<
"\""
1297 coeff.
Eval(val, *eltrans, ip);
1298 for (
int jj = 0; jj < val.
Size(); jj++)
1319 coeff.
Eval(val, *eltrans, ip);
1320 for (
int jj = 0; jj < val.
Size(); jj++)
1333 os <<
"</DataArray>" << std::endl;
1363 const std::string &collection_name,
Mesh *mesh)
1374void ParaViewHDFDataCollection::EnsureVTKHDF()
1381 MFEM_VERIFY(error_code == 0,
"Error creating directory " <<
prefix_path);
1385 bool use_mpi =
false;
1390#ifdef MFEM_PARALLEL_HDF5
1391 vtkhdf.reset(
new VTKHDF(fname, pmesh->GetComm(), {restart_mode, time}));
1393 MFEM_ABORT(
"Requires HDF5 library with parallel support enabled");
1404template <
typename FP_T>
1405void ParaViewHDFDataCollection::TSave()
1415 vtkhdf->DisableCompression();
1421 vtkhdf->SaveGridFunction<FP_T>(*field.second, field.first);
1423 vtkhdf->UpdateSteps(
time);
1433 default: MFEM_ABORT(
"Unsupported VTK format.");
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
int cycle
Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
real_t time
Physical time (for time-dependent simulations)
virtual void RegisterQField(const std::string &field_name, QuadratureFunction *qf)
Add a QuadratureFunction to the collection.
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.
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 SaveOneField(const FieldMapIterator &it)
Save one field to disk, assuming the collection directory exists.
void DeleteAll()
Delete data owned by the DataCollection including field information.
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.
virtual void SaveQField(const std::string &field_name)
Save one q-field, assuming the collection directory already exists.
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.
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.
std::string GetFieldFileName(const std::string &field_name) const
int myid
MPI rank (in parallel)
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)
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.
std::string GetMeshFileName() const
void DeleteData()
Delete data owned by the DataCollection keeping field information.
bool appendRankToFileName
Append rank to any output file names.
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.
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.
real_t GetTime() const
Get physical time (for time-dependent simulations)
Data type dense matrix using column-major storage.
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
virtual const char * Name() const
const NURBSExtension * GetNURBSext() const
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
const FiniteElementCollection * FEColl() const
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
RefinedGeometry * Refine(Geometry::Type Geom, int Times, int ETimes=1)
static int GetRefinementLevelFromElems(Geometry::Type geom, int Npts)
Get the Refinement level based on number of elements.
Class for grid function - Vector with associated FE space.
FiniteElementSpace * FESpace()
int VectorDim() const
Shortcut for calling FiniteElementSpace::GetVectorDim() on the underlying fes.
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Print the mesh to the given stream using the default MFEM mesh format.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
ElementTransformation * GetBdrElementTransformation(int i)
Returns a pointer to the transformation defining the i-th boundary element.
int SpaceDimension() const
Dimension of the physical space containing the mesh.
void PrintVTU(std::ostream &os, int ref=1, VTKFormat format=VTKFormat::ASCII, bool high_order_output=false, int compression_level=0, bool bdr_elements=false)
int GetNBE() const
Returns number of boundary elements.
Geometry::Type GetElementBaseGeometry(int i) const
Geometry::Type GetBdrElementBaseGeometry(int i) const
int GetOrder() const
If all KnotVector orders are identical, return that number. Otherwise, return NURBSFECollection::Vari...
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.
iterator begin()
Returns a begin iterator to the registered fields.
void DeleteData(bool own_data)
Clear all associations between names and fields.
void clear()
Clears the map of registered fields without reclaiming memory.
iterator find(const std::string &fname)
Returns an iterator to the field fname.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Class for parallel grid function.
Class for parallel meshes.
void ParPrint(std::ostream &out, const std::string &comments="") const
Abstract base class for ParaViewDataCollection and ParaViewHDFDataCollection.
bool IsBinaryFormat() const
Returns true if the output format is BINARY or BINARY32, false if ASCII.
ParaViewDataCollectionBase(const std::string &name, Mesh *mesh)
void SetLevelsOfDetail(int levels_of_detail_)
Set the refinement level.
void UseRestartMode(bool restart_mode_)
Enable or disable restart mode.
void SetBoundaryOutput(bool bdr_output_)
Configures collection to save only fields evaluated on boundaries of the mesh.
void SetHighOrderOutput(bool high_order_output_)
Sets whether or not to output the data as high-order elements (false by default).
void SetDataFormat(VTKFormat fmt)
Set the data format for the ParaView output files.
void SetCompressionLevel(int compression_level_)
Set the zlib compression level.
int GetCompressionLevel() const
If compression is enabled, return the compression level, else return 0.
void SaveVCoeffFieldVTU(std::ostream &out, int ref_, const std::string &name, VectorCoefficient &coeff)
std::string GenerateCollectionPath()
void WritePVTUHeader(std::ostream &out)
std::string GeneratePVDFileName()
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()
ParaViewDataCollection(const std::string &collection_name, Mesh *mesh_=nullptr)
Constructor. The collection name is used when saving the data.
void SaveCoeffFieldVTU(std::ostream &out, int ref_, const std::string &name, Coefficient &coeff)
std::string GenerateVTUFileName(const std::string &prefix, int rank)
const char * GetDataTypeString() const
std::string GeneratePVTUFileName(const std::string &prefix)
void Save() override
Save the collection.
ParaViewHDFDataCollection(const std::string &collection_name, Mesh *mesh_=nullptr)
Constructor. The collection name is used when saving the data.
~ParaViewHDFDataCollection()
Destructor.
void SetCompression(bool compression_) override
Enable or disable compression.
Represents values or vectors of values at quadrature points on a mesh.
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
int GetVDim() const
Get the vector dimension.
const IntegrationRule & GetIntRule(int idx) const
Get the IntegrationRule associated with entity (element or face) idx.
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Mesh * GetMesh() const
Returns the mesh.
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
virtual void Eval(Vector &V, ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the vector coefficient in the element described by T at the point ip, storing the result in ...
int Size() const
Returns the size of the vector.
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
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)
VisItDataCollection(const std::string &collection_name, Mesh *mesh_=NULL)
Constructor. The collection name is used when saving the data.
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.
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
void SetMesh(Mesh *new_mesh) override
Set/change the mesh associated with the collection.
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.
GeometryRefiner GlobGeometryRefiner
std::string to_padded_string(int i, int digits)
Convert an integer to a 0-padded string with the given number of digits.
void WriteBase64WithSizeAndClear(std::ostream &os, std::vector< char > &buf, int compression_level)
Encode in base 64 (and potentially compress) the given data, write it to the output stream (with a he...
VTKFormat
Data array format for VTK and VTU files.
@ ASCII
Data arrays will be written in ASCII format.
int to_int(const std::string &str)
Convert a string to an int.
void WriteBinaryOrASCII(std::ostream &os, std::vector< char > &buf, const T &val, const char *suffix, VTKFormat format)
Write either ASCII data to the stream or binary data to the buffer depending on the given format.
const char * VTKByteOrder()
Determine the byte order and return either "BigEndian" or "LittleEndian".
std::string VTKComponentLabels(int vdim)
Returns a string defining the component labels for vector-valued data arrays for use in XML VTU files...