12 #include "../config/config.hpp"
14 #ifdef MFEM_USE_CONDUIT
17 #include "../general/text.hpp"
18 #include <conduit_relay.hpp>
19 #include <conduit_blueprint.hpp>
24 using namespace conduit;
38 ConduitDataCollection::ConduitDataCollection(
const std::string& coll_name,
41 relay_protocol(
"hdf5")
50 const std::string& coll_name,
53 relay_protocol(
"hdf5")
56 MPI_Comm_rank(comm, &
myid);
76 MFEM_ABORT(
"Error creating directory: " << dir_name);
85 if (!blueprint::mesh::verify(n_mesh,verify_info))
87 MFEM_ABORT(
"Conduit Mesh Blueprint Verify Failed:\n"
88 << verify_info.to_json());
94 std::string
name = itr->first;
102 n_mesh[
"fields"][name]);
136 int num_domains = n_root[
"number_of_trees"].to_int();
141 MFEM_WARNING(
"num_procs must equal num_domains");
167 const std::string &main_toplogy_name,
180 std::string topo_name = main_toplogy_name;
184 topo_name = n_mesh[
"topologies"].schema().child_name(0);
187 MFEM_ASSERT(n_mesh.has_path(
"topologies/" + topo_name),
188 "Expected topology named \"" + topo_name +
"\" "
189 "(node is missing path \"topologies/" + topo_name +
"\")");
192 std::string coords_name =
193 n_mesh[
"topologies"][topo_name][
"coordset"].as_string();
196 MFEM_ASSERT(n_mesh.has_path(
"coordsets/" + coords_name),
197 "Expected topology named \"" + coords_name +
"\" "
198 "(node is missing path \"coordsets/" + coords_name +
"\")");
200 const Node &n_coordset = n_mesh[
"coordsets"][coords_name];
201 const Node &n_coordset_vals = n_coordset[
"values"];
204 int ndims = n_coordset_vals.number_of_children();
207 int num_verts = n_coordset_vals[0].dtype().number_of_elements();
209 const double *verts_ptr = NULL;
217 n_coordset_vals[0].dtype().is_double() &&
218 blueprint::mcarray::is_interleaved(n_coordset_vals) )
222 verts_ptr = n_coordset_vals[0].value();
229 NodeConstIterator itr = n_coordset_vals.children();
230 while (itr.has_next())
232 const Node &c_vals = itr.next();
233 std::string c_name = itr.name();
235 if ( c_vals.dtype().is_double() )
238 n_tmp[c_name].set_external(c_vals);
244 c_vals.to_double_array(n_tmp[c_name]);
253 n_tmp[
"z"].set(DataType::c_double(num_verts));
259 n_tmp[
"y"].set(DataType::c_double(num_verts));
262 Node &n_conv_coords_vals = n_conv[
"coordsets"][coords_name][
"values"];
263 blueprint::mcarray::to_interleaved(n_tmp,
265 verts_ptr = n_conv_coords_vals[0].value();
270 const Node &n_mesh_topo = n_mesh[
"topologies"][topo_name];
271 std::string mesh_ele_shape = n_mesh_topo[
"elements/shape"].as_string();
276 const Node &n_mesh_conn = n_mesh_topo[
"elements/connectivity"];
278 const int *elem_indices = NULL;
280 if (n_mesh_conn.dtype().is_int() &&
281 n_mesh_conn.is_compact() )
283 elem_indices = n_mesh_topo[
"elements/connectivity"].value();
287 Node &n_mesh_conn_conv=
288 n_conv[
"topologies"][topo_name][
"elements/connectivity"];
289 n_mesh_conn.to_int_array(n_mesh_conn_conv);
290 elem_indices = n_mesh_conn_conv.value();
294 n_mesh_topo[
"elements/connectivity"].dtype().number_of_elements();
295 num_mesh_ele = num_mesh_ele / num_idxs_per_ele;
298 const int *bndry_indices = NULL;
299 int num_bndry_ele = 0;
304 if ( n_mesh_topo.has_child(
"boundary_topology") )
306 std::string bndry_topo_name = n_mesh_topo[
"boundary_topology"].as_string();
315 if (n_mesh[
"topologies"].has_child(bndry_topo_name))
317 const Node &n_bndry_topo = n_mesh[
"topologies"][bndry_topo_name];
318 std::string bndry_ele_shape = n_bndry_topo[
"elements/shape"].as_string();
320 bndry_geo = ShapeNameToGeomType(bndry_ele_shape);
323 const Node &n_bndry_conn = n_bndry_topo[
"elements/connectivity"];
326 if ( n_bndry_conn.dtype().is_int() &&
327 n_bndry_conn.is_compact())
329 bndry_indices = n_bndry_conn.value();
333 Node &(n_bndry_conn_conv) =
334 n_conv[
"topologies"][bndry_topo_name][
"elements/connectivity"];
335 n_bndry_conn.to_int_array(n_bndry_conn_conv);
336 bndry_indices = (n_bndry_conn_conv).value();
341 n_bndry_topo[
"elements/connectivity"].dtype().number_of_elements();
342 num_bndry_ele = num_bndry_ele / num_idxs_per_bndry_ele;
350 const int *mesh_atts = NULL;
351 const int *bndry_atts = NULL;
360 std::string main_att_name =
"";
362 const Node &n_fields = n_mesh[
"fields"];
363 NodeConstIterator itr = n_fields.children();
365 while ( itr.has_next() && main_att_name ==
"" )
368 std::string fld_name = itr.name();
369 if ( fld_name.find(
"boundary") == std::string::npos &&
370 fld_name.find(
"_attribute") != std::string::npos )
372 main_att_name = fld_name;
376 if ( main_att_name !=
"" )
378 const Node &n_mesh_atts_vals = n_fields[main_att_name][
"values"];
381 if (n_mesh_atts_vals.dtype().is_int() &&
382 n_mesh_atts_vals.is_compact() )
384 mesh_atts = n_mesh_atts_vals.value();
388 Node &n_mesh_atts_vals_conv = n_conv[
"fields"][main_att_name][
"values"];
389 n_mesh_atts_vals.to_int_array(n_mesh_atts_vals_conv);
390 mesh_atts = n_mesh_atts_vals_conv.value();
402 std::string bnd_att_name =
"";
403 itr = n_fields.children();
405 while ( itr.has_next() && bnd_att_name ==
"" )
408 std::string fld_name = itr.name();
409 if ( fld_name.find(
"boundary") != std::string::npos &&
410 fld_name.find(
"_attribute") != std::string::npos )
412 bnd_att_name = fld_name;
416 if ( bnd_att_name !=
"" )
419 const Node &n_bndry_atts_vals =n_fields[bnd_att_name][
"values"];
422 if ( n_bndry_atts_vals.dtype().is_int() &&
423 n_bndry_atts_vals.is_compact())
425 bndry_atts = n_bndry_atts_vals.value();
429 Node &n_bndry_atts_vals_conv = n_conv[
"fields"][bnd_att_name][
"values"];
430 n_bndry_atts_vals.to_int_array(n_bndry_atts_vals_conv);
431 bndry_atts = n_bndry_atts_vals_conv.value();
455 const_cast<double*>(verts_ptr),
458 const_cast<int*>(elem_indices),
461 const_cast<int*>(mesh_atts),
464 const_cast<int*>(bndry_indices),
467 const_cast<int*>(bndry_atts),
473 if (n_mesh_topo.has_child(
"grid_function"))
475 std::string nodes_gf_name = n_mesh_topo[
"grid_function"].as_string();
478 const Node &n_mesh_gf = n_mesh[
"fields"][nodes_gf_name];
487 if (zero_copy && !n_conv.dtype().is_empty())
503 res =
new Mesh(*mesh,
true);
522 const double *vals_ptr = NULL;
528 if (n_field[
"values"].dtype().is_object())
530 vdim = n_field[
"values"].number_of_children();
535 if ( n_field[
"values"][0].dtype().is_double() )
538 if (n_field[
"values"].is_contiguous())
541 vals_ptr = n_field[
"values"].child(0).value();
544 else if (blueprint::mcarray::is_interleaved(n_field[
"values"]))
548 vals_ptr = n_field[
"values"].child(0).value();
554 blueprint::mcarray::to_contiguous(n_field[
"values"],
556 vals_ptr = n_conv[
"values"].child(0).value();
564 NodeConstIterator itr = n_field[
"values"].children();
565 while (itr.has_next())
567 const Node &c_vals = itr.next();
568 std::string c_name = itr.name();
570 if ( c_vals.dtype().is_double() )
573 n_tmp[c_name].set_external(c_vals);
579 c_vals.to_double_array(n_tmp[c_name]);
585 blueprint::mcarray::to_contiguous(n_tmp,
587 vals_ptr = n_conv[
"values"].child(0).value();
592 if (n_field[
"values"].dtype().is_double() &&
593 n_field[
"values"].is_compact())
595 vals_ptr = n_field[
"values"].value();
599 n_field[
"values"].to_double_array(n_conv[
"values"]);
600 vals_ptr = n_conv[
"values"].value();
604 if (zero_copy && !n_conv.dtype().is_empty())
611 std::string fec_name = n_field[
"basis"].as_string();
623 res =
new GridFunction(fes,const_cast<double*>(vals_ptr));
630 Vector vals_vec(const_cast<double*>(vals_ptr),fes->GetVSize());
646 const std::string &coordset_name,
647 const std::string &main_topology_name,
648 const std::string &boundary_topology_name)
652 MFEM_ASSERT(dim >= 1 && dim <= 3,
"invalid mesh dimension");
662 int num_vertices = mesh->
GetNV();
664 MFEM_ASSERT( ( stride == 3 *
sizeof(
double) ),
665 "Unexpected stride for Vertex");
667 Node &n_mesh_coords = n_mesh[
"coordsets"][coordset_name];
668 n_mesh_coords[
"type"] =
"explicit";
673 n_mesh_coords[
"values/x"].set_external(coords_ptr,
680 n_mesh_coords[
"values/y"].set_external(coords_ptr,
687 n_mesh_coords[
"values/z"].set_external(coords_ptr,
697 Node &n_topo = n_mesh[
"topologies"][main_topology_name];
699 n_topo[
"type"] =
"unstructured";
700 n_topo[
"coordset"] = coordset_name;
704 std::string ele_shape = ElementTypeToShapeName(ele_type);
706 n_topo[
"elements/shape"] = ele_shape;
710 if (gf_mesh_nodes != NULL)
712 n_topo[
"grid_function"] =
"mesh_nodes";
720 int num_ele = mesh->
GetNE();
723 int num_conn_idxs = num_ele * idxs_per_ele;
725 n_topo[
"elements/connectivity"].set(DataType::c_int(num_conn_idxs));
727 int *conn_ptr = n_topo[
"elements/connectivity"].value();
729 for (
int i=0; i < num_ele; i++)
734 memcpy(conn_ptr, ele_verts, idxs_per_ele *
sizeof(
int));
736 conn_ptr += idxs_per_ele;
739 if (gf_mesh_nodes != NULL)
742 n_mesh[
"fields/mesh_nodes"],
750 Node &n_mesh_att = n_mesh[
"fields/element_attribute"];
752 n_mesh_att[
"association"] =
"element";
753 n_mesh_att[
"topology"] = main_topology_name;
754 n_mesh_att[
"values"].set(DataType::c_int(num_ele));
756 int_array att_vals = n_mesh_att[
"values"].value();
757 for (
int i = 0; i < num_ele; i++)
769 n_topo[
"boundary_topology"] = boundary_topology_name;
771 Node &n_bndry_topo = n_mesh[
"topologies"][boundary_topology_name];
773 n_bndry_topo[
"type"] =
"unstructured";
774 n_bndry_topo[
"coordset"] = coordset_name;
778 std::string bndry_ele_shape = ElementTypeToShapeName(bndry_ele_type);
780 n_bndry_topo[
"elements/shape"] = bndry_ele_shape;
783 int num_bndry_ele = mesh->
GetNBE();
786 int num_bndry_conn_idxs = num_bndry_ele * bndry_idxs_per_ele;
788 n_bndry_topo[
"elements/connectivity"].set(DataType::c_int(num_bndry_conn_idxs));
790 int *bndry_conn_ptr = n_bndry_topo[
"elements/connectivity"].value();
792 for (
int i=0; i < num_bndry_ele; i++)
795 const int *bndry_ele_verts = bndry_ele->
GetVertices();
797 memcpy(bndry_conn_ptr, bndry_ele_verts, bndry_idxs_per_ele *
sizeof(
int));
799 bndry_conn_ptr += bndry_idxs_per_ele;
806 Node &n_bndry_mesh_att = n_mesh[
"fields/boundary_attribute"];
808 n_bndry_mesh_att[
"association"] =
"element";
809 n_bndry_mesh_att[
"topology"] = boundary_topology_name;
810 n_bndry_mesh_att[
"values"].set(DataType::c_int(num_bndry_ele));
812 int_array bndry_att_vals = n_bndry_mesh_att[
"values"].value();
813 for (
int i = 0; i < num_bndry_ele; i++)
824 const std::string &main_topology_name)
827 n_field[
"topology"] = main_topology_name;
834 n_field[
"values"].set_external(gf->
GetData(),
847 index_t stride =
sizeof(double) * entry_stride;
849 for (
int d = 0; d < vdim; d++)
851 std::ostringstream oss;
853 std::string comp_name = oss.str();
854 n_field[
"values"][comp_name].set_external(gf->
GetData(),
858 offset +=
sizeof(double) * vdim_stride;
889 const std::string &relay_protocol)
918 std::ostringstream oss;
936 const std::string &relay_protocol)
939 std::string root_proto =
"json";
941 if (relay_protocol ==
"hdf5")
948 Node &n_bp_idx = n_root[
"blueprint_index"];
950 blueprint::mesh::generate_index(n_mesh,
962 std::string gf_name = itr->first;
965 Node &idx_gf_ncomps = n_bp_idx[
"mesh/fields"][gf_name][
"number_of_components"];
968 if ( idx_gf_ncomps.to_int() != gf->
VectorDim() )
975 n_root[
"protocol/version"] =
"0.3.1";
979 n_root[
"number_of_files"] = num_domains;
980 n_root[
"number_of_trees"] = num_domains;
982 n_root[
"tree_pattern"] =
"";
991 const std::string &relay_protocol)
993 relay::io::save(n_mesh,
MeshFileName(domain_id, relay_protocol));
1003 std::string root_protocol =
"json";
1007 root_protocol =
"hdf5";
1011 relay::io::load(
RootFileName(), root_protocol, root_out);
1018 std::string root_json = root_out.to_json();
1020 int json_str_size = root_json.size() + 1;
1023 int mpi_status = MPI_Bcast((
void*)&json_str_size,
1029 if (mpi_status != MPI_SUCCESS)
1031 MFEM_ABORT(
"Broadcast of root file json string size failed");
1035 mpi_status = MPI_Bcast((
void*)root_json.c_str(),
1041 if (mpi_status != MPI_SUCCESS)
1043 MFEM_ABORT(
"Broadcast of root file json string failed");
1053 int json_str_size = -1;
1054 int mpi_status = MPI_Bcast(&json_str_size,
1060 if (mpi_status != MPI_SUCCESS)
1062 MFEM_ABORT(
"Broadcast of root file json string size failed");
1066 char *json_buff =
new char[json_str_size];
1067 mpi_status = MPI_Bcast(json_buff,
1073 if (mpi_status != MPI_SUCCESS)
1075 MFEM_ABORT(
"Broadcast of root file json string failed");
1079 Generator g(std::string(json_buff),
"json");
1082 delete [] json_buff;
1090 const std::string &relay_protocol)
1096 relay::io::load(
MeshFileName(domain_id, relay_protocol), n_mesh);
1100 if (!blueprint::mesh::verify(n_mesh,verify_info))
1102 MFEM_ABORT(
"Conduit Mesh Blueprint Verify Failed:\n"
1103 << verify_info.to_json());
1110 NodeConstIterator itr = n_mesh[
"fields"].children();
1112 std::string nodes_gf_name =
"";
1114 const Node &n_topo = n_mesh[
"topologies/main"];
1115 if (n_topo.has_child(
"grid_function"))
1117 nodes_gf_name = n_topo[
"grid_function"].as_string();
1120 while (itr.has_next())
1122 const Node &n_field = itr.next();
1123 std::string field_name = itr.name();
1127 if ( field_name != nodes_gf_name &&
1128 field_name.find(
"_attribute") == std::string::npos
1147 ConduitDataCollection::ElementTypeToShapeName(
Element::Type element_type)
1156 switch (element_type)
1173 ConduitDataCollection::ShapeNameToGeomType(
const std::string &shape_name)
1179 if (shape_name ==
"point")
1183 else if (shape_name ==
"line")
1187 else if (shape_name ==
"tri")
1191 else if (shape_name ==
"quad")
1195 else if (shape_name ==
"tet")
1199 else if (shape_name ==
"hex")
1205 MFEM_ABORT(
"Unsupported Element Shape: " << shape_name);
Ordering::Type GetOrdering() const
Return the ordering method.
int GetNDofs() const
Returns number of degrees of freedom.
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
const double * GetVertex(int i) const
Return pointer to vertex i's coordinates.
Class for grid function - Vector with associated FE space.
void SaveMeshAndFields(int domain_id, const conduit::Node &n_mesh, const std::string &file_protocol)
Saves all meshes and fields for the current cycle.
bool appendRankToFileName
Append rank to any output file names.
iterator begin()
Returns a begin iterator to the registered fields.
virtual void GetVertices(Array< int > &v) const =0
Returns element's vertices.
std::string RootFileName()
Returns blueprint root file name for the current cycle.
void DeleteAll()
Delete data owned by the DataCollection including field information.
int GetNBE() const
Returns number of boundary elements.
void MakeOwner(FiniteElementCollection *_fec)
Make the GridFunction the owner of fec and fes.
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.
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Element::Type GetElementType(int i) const
Returns the type of element i.
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
int GetNE() const
Returns number of elements.
bool own_data
Should the collection delete its mesh and fields.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Geometry::Type GetElementBaseGeometry(int i) const
std::string relay_protocol
virtual void Load(int cycle=0)
Load the collection based blueprint data.
ConduitDataCollection(const std::string &collection_name, Mesh *mesh=NULL)
Constructor. The collection name is used when saving the data.
MPI_Comm m_comm
Associated MPI communicator.
static Mesh * BlueprintMeshToMesh(const conduit::Node &n_mesh, const std::string &main_toplogy_name="", bool zero_copy=false)
Constructs and MFEM mesh from a Conduit Blueprint Description.
Mesh * mesh
The (common) mesh for the collected fields.
virtual ~ConduitDataCollection()
We will delete the mesh and fields if we own them.
int cycle
Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
std::string MeshFileName(int domain_id, const std::string &file_protocol="hdf5")
Returns the mesh file name for a given domain at the current cycle.
std::string to_padded_string(int i, int digits)
static const int NumVerts[NumGeom]
Type
Constants for the classes derived from Element.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
virtual void Save()
Save the collection and a Conduit blueprint root file.
void Register(const std::string &fname, T *field, bool own_data)
Register field field with name fname.
const Element * GetElement(int i) const
void LoadRootFile(conduit::Node &n_root_out)
Loads contents of the root field for the current cycle into n_root_out.
int GetVDim() const
Returns vector dimension.
FiniteElementSpace * FESpace()
int SpaceDimension() const
void SaveRootFile(int num_domains, const conduit::Node &n_mesh, const std::string &file_protocol)
Saves root file for the current cycle.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
void LoadMeshAndFields(int domain_id, const std::string &file_protocol)
Loads all meshes and fields of a given domain id for the current cycle.
int myid
MPI rank (in parallel)
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
virtual const char * Name() const
GFieldMap::const_iterator FieldMapConstIterator
void clear()
Clears the map of registered fields without reclaiming memory.
std::string MeshDirectoryName()
Returns the mesh output directory for the current cycle.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has '/' at the end...
Geometry::Type GetBdrElementBaseGeometry(int i) const
const FiniteElementCollection * FEColl() const
void GetNodes(Vector &node_coord) const
static void MeshToBlueprintMesh(Mesh *m, conduit::Node &out, const std::string &coordset_name="coords", const std::string &main_topology_name="main", const std::string &boundary_topology_name="boundary")
Describes a MFEM mesh using the mesh blueprint.
std::string name
Name of the collection, used as a directory name when saving.
std::string MeshFilePattern(const std::string &file_protocol="hdf5")
Returns the mesh file pattern for the current cycle.
static void GridFunctionToBlueprintField(GridFunction *gf, conduit::Node &out, const std::string &main_topology_name="main")
Describes a MFEM grid function using the mesh blueprint.
int num_procs
Number of MPI ranks (in parallel)
static GridFunction * BlueprintFieldToGridFunction(Mesh *mesh, const conduit::Node &n_field, bool zero_copy=false)
Constructs and MFEM Grid Function from a Conduit Blueprint Description.
Abstract data type element.
int GetAttribute(int i) const
Return the attribute of element i.
void SetProtocol(const std::string &protocol)
Set the Conduit relay i/o protocol to use.
const Element * GetBdrElement(int i) const