13 #include "../config/config.hpp"
20 #include <spio/IOManager.hpp>
27 namespace sidre = asctoolkit::sidre;
35 Mesh * the_mesh,
bool own_mesh_data)
37 m_owns_datastore(true),
38 m_owns_mesh_data(own_mesh_data),
39 m_meshNodesGFName(
"mesh_nodes")
41 m_datastore_ptr =
new sidre::DataStore();
43 sidre::DataGroup * global_grp =
44 m_datastore_ptr->getRoot()->createGroup(collection_name +
"_global");
45 sidre::DataGroup * domain_grp =
46 m_datastore_ptr->getRoot()->createGroup(collection_name);
48 bp_grp = domain_grp->createGroup(
"blueprint");
50 bp_index_grp = global_grp->createGroup(
"blueprint_index/" +
name);
52 named_bufs_grp = domain_grp->createGroup(
"named_buffers");
61 m_comm = MPI_COMM_NULL;
73 asctoolkit::sidre::DataGroup* global_grp,
74 asctoolkit::sidre::DataGroup* domain_grp,
77 m_owns_datastore(false),
78 m_owns_mesh_data(own_mesh_data),
79 m_meshNodesGFName(
"mesh_nodes"),
82 bp_grp = domain_grp->createGroup(
"blueprint");
85 bp_index_grp = global_grp->createGroup(
"blueprint_index/" +
name);
87 named_bufs_grp = domain_grp->createGroup(
"named_buffers");
90 m_comm = MPI_COMM_NULL;
98 delete m_datastore_ptr;
108 MPI_Comm_rank(m_comm, &
myid);
116 MFEM_ASSERT(named_bufs_grp != NULL,
117 "No group 'named_buffers' in data collection. Verify that"
118 " SetMesh was called to set the mesh in the data collection.");
119 return named_bufs_grp;
123 asctoolkit::sidre::DataView *
125 const std::string &view_name)
127 MFEM_ASSERT(grp,
"DataGroup pointer is NULL");
128 sidre::DataView *v = grp->getView(view_name);
131 v = grp->createView(view_name);
132 MFEM_ASSERT(v,
"error allocating DataView " << view_name
133 <<
" in group " << grp->getPathName());
139 asctoolkit::sidre::DataView *
141 const std::string &view_name,
142 const asctoolkit::sidre::DataType &dtype)
144 MFEM_ASSERT(grp,
"DataGroup pointer is NULL");
145 sidre::DataView *v = grp->getView(view_name);
148 v = grp->createView(view_name, dtype);
149 MFEM_ASSERT(v,
"error allocating DataView " << view_name
150 <<
" in group " << grp->getPathName());
154 MFEM_ASSERT(v->getSchema().dtype().equals(dtype),
"");
160 asctoolkit::sidre::DataGroup *
162 const std::string &group_name)
164 MFEM_ASSERT(grp,
"DataGroup pointer is NULL");
165 sidre::DataGroup *g = grp->getGroup(group_name);
168 g = grp->createGroup(group_name);
169 MFEM_ASSERT(g,
"error allocating DataGroup " << group_name
170 <<
" in group " << grp->getPathName());
179 std::stringstream fNameSstr;
186 fNameSstr <<
"_" << std::setfill(
'0') << std::setw(
pad_digits)
190 return fNameSstr.str();
193 asctoolkit::sidre::DataView *
195 asctoolkit::sidre::SidreLength sz,
196 asctoolkit::sidre::TypeID type)
198 sz = std::max(sz, sidre::SidreLength(0));
200 sidre::DataView *v = f->getView(buffer_name);
204 v = f->createViewAndAllocate(buffer_name, type, sz);
208 MFEM_ASSERT(v->getTypeID() == type,
"type does not match existing type");
214 if (!v->isApplied() || v->getNumElements() < sz)
218 sidre::DataType dtype(v->getSchema().dtype());
219 dtype.set_number_of_elements(sz);
220 f->destroyViewAndData(buffer_name);
221 v = f->createViewAndAllocate(buffer_name, dtype);
224 MFEM_ASSERT(v && v->isApplied(),
"allocation failed");
229 void SidreDataCollection::createMeshBlueprintStubs(
bool hasBP)
233 bp_grp->createGroup(
"state");
234 bp_grp->createGroup(
"coordsets");
235 bp_grp->createGroup(
"topologies");
236 bp_grp->createGroup(
"fields");
242 bp_index_grp->createGroup(
"state");
243 bp_index_grp->createGroup(
"coordsets");
244 bp_index_grp->createGroup(
"topologies");
245 bp_index_grp->createGroup(
"fields");
250 void SidreDataCollection::createMeshBlueprintState(
bool hasBP)
255 bp_grp->createViewScalar(
"state/cycle", 0);
256 bp_grp->createViewScalar(
"state/time", 0.);
257 bp_grp->createViewScalar(
"state/domain",
myid);
258 bp_grp->createViewScalar(
"state/time_step", 0.);
264 bp_index_grp->createViewScalar(
"state/cycle", 0);
265 bp_index_grp->createViewScalar(
"state/time", 0.);
266 bp_index_grp->createViewScalar(
"state/number_of_domains",
num_procs);
271 void SidreDataCollection::createMeshBlueprintCoordset(
bool hasBP)
274 MFEM_ASSERT(dim >= 1 && dim <= 3,
"invalid mesh dimension");
277 const int NUM_COORDS =
sizeof(
mfem::Vertex) /
sizeof(
double);
280 const int coordset_len = NUM_COORDS * num_vertices;
285 bp_grp->createViewString(
"coordsets/coords/type",
"explicit");
287 sidre::DataType dtype =
288 sidre::DataType::c_double(num_vertices);
289 const size_t stride = dtype.stride();
290 dtype.set_stride(stride*NUM_COORDS);
293 sidre::DataView *vx, *vy = NULL, *vz = NULL;
294 vx = bp_grp->createView(
"coordsets/coords/values/x", dtype);
298 dtype.set_offset(dtype.offset() + stride);
299 vy = bp_grp->createView(
"coordsets/coords/values/y", dtype);
303 dtype.set_offset(dtype.offset() + stride);
304 vz = bp_grp->createView(
"coordsets/coords/values/z", dtype);
307 if (m_owns_mesh_data)
310 sidre::DataBuffer* coordbuf =
313 vx->attachBuffer(coordbuf);
314 if (dim >= 2) { vy->attachBuffer(coordbuf); }
315 if (dim >= 3) { vz->attachBuffer(coordbuf); }
321 vx->setExternalDataPtr(coordbuf);
322 if (dim >= 2) { vy->setExternalDataPtr(coordbuf); }
323 if (dim >= 3) { vz->setExternalDataPtr(coordbuf); }
331 bp_index_grp->createViewString(
332 "coordsets/coords/path", bp_grp->getPathName() +
"/coordsets/coords");
334 bp_index_grp->getGroup(
"coordsets/coords")->copyView(
335 bp_grp->getView(
"coordsets/coords/type") );
337 bp_index_grp->createViewString(
338 "coordsets/coords/coord_system/type",
"cartesian");
342 bp_index_grp->createView(
"coordsets/coords/coord_system/axes/x");
346 bp_index_grp->createView(
"coordsets/coords/coord_system/axes/y");
351 bp_index_grp->createView(
"coordsets/coords/coord_system/axes/z");
355 if (m_owns_mesh_data)
357 double *coord_values =
GetNamedBuffer(
"vertex_coords")->getData();
364 void SidreDataCollection::
365 createMeshBlueprintTopologies(
bool hasBP,
const std::string& mesh_name)
367 const bool isBdry = (mesh_name ==
"boundary");
369 const int num_elements = !isBdry
373 MFEM_VERIFY(num_elements > 0,
374 "TODO: processors with 0 " << mesh_name <<
" elements");
376 const int element_size = !isBdry
380 const int num_indices = num_elements * element_size;
388 const std::string eltTypeStr =
390 ? getElementName( static_cast<Element::Type>(
392 : getElementName( static_cast<Element::Type>(
393 mesh->GetBdrElement(0)->GetType() ) );
395 const std::string mesh_topo_str =
"topologies/" + mesh_name;
396 const std::string mesh_attr_str =
"fields/"+mesh_name+
"_material_attribute";
400 sidre::DataGroup* topology_grp = bp_grp->createGroup(mesh_topo_str);
403 topology_grp->createViewString(
"type",
"unstructured");
405 topology_grp->createViewString(
"elements/shape", eltTypeStr);
406 topology_grp->createViewAndAllocate(
407 "elements/connectivity", sidre::INT_ID, num_indices);
408 topology_grp->createViewString(
"coordset",
"coords");
414 topology_grp->createViewString(
"grid_function",m_meshNodesGFName);
418 sidre::DataGroup* attr_grp = bp_grp->createGroup(mesh_attr_str);
419 attr_grp->createViewString(
"association",
"element");
420 attr_grp->createViewAndAllocate(
"values", sidre::INT_ID, num_elements);
421 attr_grp->createViewString(
"topology", mesh_name);
428 const std::string bp_grp_path = bp_grp->getPathName();
436 bp_index_grp->getGroup(
"topologies/mesh")
437 ->copyView( bp_grp->getView(
"topologies/mesh/boundary_topology") );
440 sidre::DataGroup *bp_index_topo_grp =
441 bp_index_grp->createGroup(mesh_topo_str);
442 sidre::DataGroup *topology_grp = bp_grp->getGroup(mesh_topo_str);
444 bp_index_topo_grp->createViewString(
445 "path", bp_grp_path +
"/" + mesh_topo_str);
446 bp_index_topo_grp->copyView( topology_grp->getView(
"type") );
447 bp_index_topo_grp->copyView( topology_grp->getView(
"coordset") );
453 bp_index_topo_grp->copyView(topology_grp->getView(
"grid_function"));
457 sidre::DataGroup *bp_index_attr_grp =
458 bp_index_grp->createGroup(mesh_attr_str);
459 sidre::DataGroup *attr_grp = bp_grp->getGroup(mesh_attr_str);
461 bp_index_attr_grp->createViewString(
462 "path", bp_grp_path +
"/" + mesh_attr_str );
463 bp_index_attr_grp->copyView( attr_grp->getView(
"association") );
464 bp_index_attr_grp->copyView( attr_grp->getView(
"topology") );
466 int number_of_components = 1;
467 bp_index_attr_grp->createViewScalar(
"number_of_components",
468 number_of_components);
472 sidre::DataView* conn_view =
473 bp_grp->getGroup(mesh_topo_str)->getView(
"elements/connectivity");
474 sidre::DataView* attr_view =
475 bp_grp->getGroup(mesh_attr_str)->getView(
"values");
477 Array<int> conn_array(conn_view->getData<
int*>(), num_indices);
478 Array<int> attr_array(attr_view->getData<
int*>(), num_elements);
487 MFEM_ASSERT(!conn_array.OwnsData(),
"");
488 MFEM_ASSERT(!attr_array.OwnsData(),
"");
492 void SidreDataCollection::verifyMeshBlueprint()
504 m_comm = pmesh ? pmesh->
GetComm() : MPI_COMM_NULL;
509 bool hasBP = bp_grp->getNumViews() > 0 || bp_grp->getNumGroups() > 0;
510 bool has_bnd_elts = (new_mesh->
GetNBE() > 0);
512 createMeshBlueprintStubs(hasBP);
513 createMeshBlueprintState(hasBP);
514 createMeshBlueprintCoordset(hasBP);
519 createMeshBlueprintTopologies(hasBP,
"mesh");
524 bp_grp->createViewString(
"topologies/mesh/boundary_topology",
"boundary");
527 createMeshBlueprintTopologies(hasBP,
"boundary");
536 sidre::DataView *v_bp_nodes_name =
537 bp_grp->getView(
"topologies/mesh/grid_function");
538 std::string bp_nodes_name(v_bp_nodes_name->getString());
541 MFEM_VERIFY(m_meshNodesGFName == bp_nodes_name,
542 "mismatch of requested and blueprint mesh nodes names");
546 if (m_owns_mesh_data)
557 MFEM_ASSERT(nodes->
Size() == sz,
"");
558 std::memcpy(gfData, nodes->
GetData(),
sizeof(double) * sz);
579 MFEM_VERIFY(new_mesh->
OwnsNodes(),
"mesh does not own its nodes, "
580 "can not take ownership");
588 asctoolkit::sidre::DataGroup *domain_grp)
590 MFEM_VERIFY(domain_grp->hasGroup(
"blueprint"),
591 "Domain group does not contain a blueprint group.");
592 MFEM_VERIFY(global_grp->hasGroup(
"blueprint_index/" +
name),
593 "Global group does not contain a blueprint index group.");
595 bp_grp = domain_grp->getGroup(
"blueprint");
596 bp_index_grp = global_grp->getGroup(
"blueprint_index/" +
name);
597 named_bufs_grp = domain_grp->getGroup(
"named_buffers");
601 const std::string& protocol)
607 if (m_comm != MPI_COMM_NULL)
609 asctoolkit::spio::IOManager reader(m_comm);
610 reader.read(bp_grp->getDataStore()->getRoot(), path);
615 bp_grp->load(path, protocol);
624 if (m_owns_datastore)
627 m_datastore_ptr->getRoot()->getGroup(
name));
636 if (m_comm != MPI_COMM_NULL)
638 asctoolkit::spio::IOManager reader(m_comm);
639 reader.loadExternalData(bp_grp->getDataStore()->getRoot(), path);
644 bp_grp->loadExternalData(path);
650 SetTime( bp_grp->getView(
"state/time")->getData<
double>() );
651 SetCycle( bp_grp->getView(
"state/cycle")->getData<
int>() );
652 SetTimeStep( bp_grp->getView(
"state/time_step")->getData<
double>() );
657 bp_grp->getView(
"state/cycle")->setScalar(
GetCycle());
658 bp_grp->getView(
"state/time")->setScalar(
GetTime());
659 bp_grp->getView(
"state/time_step")->setScalar(
GetTimeStep());
663 bp_index_grp->getView(
"state/cycle")->setScalar(
GetCycle());
664 bp_index_grp->getView(
"state/time")->setScalar(
time);
670 verifyMeshBlueprint();
676 std::string filename =
name;
677 std::string protocol =
"sidre_hdf5";
679 Save(filename, protocol);
683 const std::string& protocol)
691 sidre::DataGroup * blueprint_indicies_grp = bp_index_grp->getParent();
693 if (m_comm != MPI_COMM_NULL)
695 asctoolkit::spio::IOManager writer(m_comm);
696 sidre::DataStore *datastore = bp_grp->getDataStore();
697 writer.write(datastore->getRoot(),
num_procs, file_path, protocol);
700 if (protocol ==
"sidre_hdf5")
702 writer.writeGroupToRootFile(blueprint_indicies_grp,
703 file_path +
".root");
708 writer.write(blueprint_indicies_grp, 1,
709 file_path +
".root", protocol);
717 bp_grp->save(file_path, protocol);
719 blueprint_indicies_grp
720 ->save(file_path +
".root", protocol);
725 void SidreDataCollection::
726 addScalarBasedGridFunction(
const std::string &field_name,
GridFunction *gf,
727 const std::string &buffer_name,
728 asctoolkit::sidre::SidreLength offset)
730 sidre::DataGroup* grp = bp_grp->getGroup(
"fields/" + field_name);
731 MFEM_ASSERT(grp != NULL,
"field " << field_name <<
" does not exist");
750 sidre::DataView *vv =
alloc_view(grp,
"values");
758 MFEM_ASSERT(bv->hasBuffer() && bv->isDescribed(),
"");
761 MFEM_ASSERT(bv->getSchema().dtype().offset() == 0,
"");
762 MFEM_ASSERT(bv->getNumElements() >= offset + numDofs,
"");
766 vv->attachBuffer(bv->getBuffer())
767 ->apply(sidre::DOUBLE_ID, numDofs, offset);
776 vv->setExternalDataPtr(sidre::DOUBLE_ID, numDofs, gf->
GetData());
778 MFEM_ASSERT((numDofs > 0 && vv->isApplied()) ||
779 (numDofs == 0 && vv->isEmpty() && vv->isDescribed()),
780 "invlid DataView state");
781 MFEM_ASSERT(numDofs == 0 || vv->getData() == gf->
GetData(),
782 "DataView data is different from GridFunction data");
783 MFEM_ASSERT(vv->getNumElements() == numDofs,
784 "DataView size is different from GridFunction size");
788 void SidreDataCollection::
789 addVectorBasedGridFunction(
const std::string& field_name, GridFunction *gf,
790 const std::string &buffer_name,
791 asctoolkit::sidre::SidreLength offset)
793 sidre::DataGroup* grp = bp_grp->getGroup(
"fields/" + field_name);
794 MFEM_ASSERT(grp != NULL,
"field " << field_name <<
" does not exist");
796 const int FLD_SZ = 20;
797 char fidxName[FLD_SZ];
799 int vdim = gf->FESpace()->GetVDim();
800 int ndof = gf->FESpace()->GetNDofs();
803 if (gf->GetData() == NULL)
826 sidre::DataType dtype = sidre::DataType::c_double(ndof);
829 dtype.set_stride(dtype.stride()*entry_stride);
834 MFEM_ASSERT(bv->hasBuffer() && bv->isDescribed(),
"");
837 MFEM_ASSERT(bv->getSchema().dtype().offset() == 0,
"");
838 dtype.set_offset(dtype.element_bytes()*offset);
840 for (
int d = 0; d < vdim; d++)
842 std::snprintf(fidxName, FLD_SZ,
"x%d", d);
843 sidre::DataView *xv =
alloc_view(vg, fidxName, dtype);
844 xv->attachBuffer(bv->getBuffer());
845 dtype.set_offset(dtype.offset() + dtype.element_bytes()*vdim_stride);
848 gf->NewDataAndSize(bv->getData<
double*>() + offset, vdim*ndof);
852 for (
int d = 0; d < vdim; d++)
854 std::snprintf(fidxName, FLD_SZ,
"x%d", d);
855 sidre::DataView *xv =
alloc_view(vg, fidxName, dtype);
856 xv->setExternalDataPtr(gf->GetData());
857 dtype.set_offset(dtype.offset() + dtype.element_bytes()*vdim_stride);
862 for (
int d = 0; d < vdim; d++)
864 std::snprintf(fidxName, FLD_SZ,
"x%d", d);
865 sidre::DataView *xv = vg->getView(fidxName);
866 MFEM_ASSERT((ndof > 0 && xv->isApplied()) ||
867 (ndof == 0 && xv->isEmpty() && xv->isDescribed()),
868 "invlid DataView state");
869 MFEM_ASSERT(ndof == 0 || xv->getData() == gf->GetData() + d*vdim_stride,
870 "DataView data is different from GridFunction data");
871 MFEM_ASSERT(xv->getNumElements() == ndof,
872 "DataView size is different from GridFunction size");
879 void SidreDataCollection::
880 RegisterFieldInBPIndex(
const std::string& field_name, GridFunction *gf)
882 sidre::DataGroup *bp_field_grp = bp_grp->getGroup(
"fields/" + field_name);
883 sidre::DataGroup *bp_index_field_grp =
884 bp_index_grp->createGroup(
"fields/" + field_name);
886 bp_index_field_grp->createViewString(
"path", bp_field_grp->getPathName() );
887 bp_index_field_grp->copyView( bp_field_grp->getView(
"topology") );
888 bp_index_field_grp->copyView( bp_field_grp->getView(
"basis") );
893 const int number_of_components = gf->VectorDim();
894 bp_index_field_grp->createViewScalar(
"number_of_components",
895 number_of_components);
900 void SidreDataCollection::
901 DeregisterFieldInBPIndex(
const std::string& field_name)
903 sidre::DataGroup * fields_grp = bp_index_grp->getGroup(
"fields");
904 MFEM_VERIFY(fields_grp->hasGroup(field_name),
905 "No field exists in blueprint index with name " <<
name);
910 fields_grp->destroyGroup(field_name);
915 const std::string &buffer_name,
916 asctoolkit::sidre::SidreLength offset)
918 if ( field_name.empty() || buffer_name.empty() ||
919 gf == NULL || gf->
FESpace() == NULL )
925 sidre::DataGroup* f = bp_grp->getGroup(
"fields");
927 if (f->hasGroup( field_name ))
937 MFEM_WARNING(
"field with the name '" << field_name<<
"' is already "
938 "registered, overwriting the old field"));
943 sidre::DataGroup* grp = f->createGroup( field_name );
947 sidre::DataView *v =
alloc_view(grp,
"basis");
953 v =
alloc_view(grp,
"topology")->setString(
"mesh");
956 "GridFunction size does not match FiniteElementSpace size");
965 addScalarBasedGridFunction(field_name, gf, buffer_name, offset);
970 addVectorBasedGridFunction(field_name, gf, buffer_name, offset);
976 RegisterFieldInBPIndex(field_name, gf);
988 sidre::DataGroup * fields_grp = bp_grp->getGroup(
"fields");
989 MFEM_VERIFY(fields_grp->hasGroup(field_name),
990 "No field exists in blueprint with name " << field_name);
997 fields_grp->destroyGroup(field_name);
1002 DeregisterFieldInBPIndex(field_name);
1010 std::string SidreDataCollection::getElementName(
Element::Type elementEnum)
1017 switch (elementEnum)
asctoolkit::sidre::DataView * AllocNamedBuffer(const std::string &buffer_name, asctoolkit::sidre::SidreLength sz, asctoolkit::sidre::TypeID type=asctoolkit::sidre::DOUBLE_ID)
Return newly allocated or existing named buffer for buffer_name.
asctoolkit::sidre::DataGroup * alloc_group(asctoolkit::sidre::DataGroup *grp, const std::string &group_name)
const double * GetVertex(int i) const
Return pointer to vertex i's coordinates.
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Class for grid function - Vector with associated FE space.
void Load(const std::string &path, const std::string &protocol)
Load the Sidre DataStore from file.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
double GetTime() const
Get physical time (for time-dependent simulations)
bool appendRankToFileName
Append rank to any output file names.
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
void DeleteAll()
Delete data owned by the DataCollection including field information.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Register a GridFunction in the Sidre DataStore.
int GetNBE() const
Returns number of boundary elements.
asctoolkit::sidre::DataView * GetNamedBuffer(const std::string &buffer_name) const
Get a pointer to the sidre::DataView holding the named buffer for buffer_name.
asctoolkit::sidre::DataView * alloc_view(asctoolkit::sidre::DataGroup *grp, const std::string &view_name)
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
bool HasField(const std::string &name) const
Check if a grid function is part of the collection.
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
SidreDataCollection(const std::string &collection_name, Mesh *the_mesh=NULL, bool owns_mesh_data=false)
Constructor that allocates and initializes a Sidre DataStore.
void LoadExternalData(const std::string &path)
Load external data after registering externally owned fields.
virtual void DeregisterField(const std::string &field_name)
De-register field_name from the SidreDataCollection.
bool own_data
Should the collection delete its mesh and fields.
asctoolkit::sidre::DataGroup * named_buffers_grp() const
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
double time
Physical time (for time-dependent simulations)
Mesh * mesh
The (common) mesh for the collected fields.
bool serial
Serial or parallel run?
int GetCycle() const
Get time cycle (for time-dependent simulations)
Type
Constants for the classes derived from Element.
static int create_directory(const std::string &dir_name, const Mesh *mesh, int myid)
void SetNodesOwner(bool nodes_owner)
Set the mesh nodes ownership flag.
const Element * GetElement(int i) const
void SetGroupPointers(asctoolkit::sidre::DataGroup *global_grp, asctoolkit::sidre::DataGroup *domain_grp)
Reset the domain and global datastore group pointers.
int GetVDim() const
Returns vector dimension.
FiniteElementSpace * FESpace()
void UpdateStateToDS()
Updates the data store's cycle, time, and time-step variables with the values from the SidreDataColle...
void SetTime(double t)
Set physical time (for time-dependent simulations)
int SpaceDimension() const
void UpdateStateFromDS()
Updates the DataCollection's cycle, time, and time-step variables with the values from the data store...
virtual ~SidreDataCollection()
Delete all owned data.
int myid
MPI rank (in parallel)
std::string get_file_path(const std::string &filename) const
virtual void PrepareToSave()
Prepare the DataStore for writing.
int GetElementBaseGeometry(int i=0) const
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
void SetTimeStep(double ts)
Set the simulation time step (for time-dependent simulations)
virtual const char * Name() const
bool OwnsNodes() const
Return the mesh nodes ownership flag.
const FiniteElementSpace * GetNodalFESpace() const
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
virtual void Save()
Save the collection to file.
virtual int GetNVertices() const =0
void FreeNamedBuffer(const std::string &buffer_name)
Deallocate the named buffer buffer_name.
int pad_digits
Number of digits used for the cycle and MPI rank in filenames.
const FiniteElementCollection * FEColl() const
void GetNodes(Vector &node_coord) const
void SetComm(MPI_Comm comm)
Associate an MPI communicator with the collection.
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
void ChangeVertexDataOwnership(double *vertices, int len_vertices, bool zerocopy=false)
Set the internal Vertex array to point to the given vertices array without assuming ownership of the ...
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.
int num_procs
Number of MPI ranks (in parallel)
Class for parallel meshes.
int GetBdrElementBaseGeometry(int i=0) const
const Element * GetBdrElement(int i) const
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
virtual int GetType() const =0
Returns element's type.