MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
sidredatacollection.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 
13 #include "../config/config.hpp"
14 
15 #ifdef MFEM_USE_SIDRE
16 
17 #include "fem.hpp"
18 
19 #ifdef MFEM_USE_MPI
20 #include <sidre/IOManager.hpp>
21 #endif
22 
23 #include <string>
24 #include <iomanip> // for setw, setfill
25 #include <cstdio> // for snprintf()
26 
27 namespace sidre = axom::sidre;
28 
29 namespace mfem
30 {
31 
32 // Constructor that will automatically create the sidre data store and necessary
33 // data groups for domain and global data.
34 SidreDataCollection::SidreDataCollection(const std::string& collection_name,
35  Mesh * the_mesh, bool own_mesh_data)
36  : mfem::DataCollection(collection_name, the_mesh),
37  m_owns_datastore(true),
38  m_owns_mesh_data(own_mesh_data),
39  m_meshNodesGFName("mesh_nodes")
40 {
41  m_datastore_ptr = new sidre::DataStore();
42 
43  sidre::Group * global_grp =
44  m_datastore_ptr->getRoot()->createGroup(collection_name + "_global");
45  sidre::Group * domain_grp =
46  m_datastore_ptr->getRoot()->createGroup(collection_name);
47 
48  m_bp_grp = domain_grp->createGroup("blueprint");
49  // Currently only rank 0 adds anything to bp_index.
50  m_bp_index_grp = global_grp->createGroup("blueprint_index/" + name);
51 
52  m_named_bufs_grp = domain_grp->createGroup("named_buffers");
53 
54  if (the_mesh)
55  {
56  SetMesh(the_mesh);
57  }
58 #ifdef MFEM_USE_MPI
59  else
60  {
61  m_comm = MPI_COMM_NULL;
62  }
63 #endif
64 }
65 
66 // Second constructor that allows external code to specify data groups to place
67 // domain and global data in.
68 
69 // TODO - Conduit will have the capability to generate a blueprint index group
70 // in the future. When this is available, all the blueprint index code can be
71 // removed from the data collection class.
72 SidreDataCollection::SidreDataCollection(const std::string& collection_name,
73  axom::sidre::Group* bp_index_grp,
74  axom::sidre::Group* domain_grp,
75  bool own_mesh_data)
76  : mfem::DataCollection(collection_name),
77  m_owns_datastore(false),
78  m_owns_mesh_data(own_mesh_data),
79  m_meshNodesGFName("mesh_nodes"),
80  m_datastore_ptr(NULL),
81  m_bp_index_grp(bp_index_grp)
82 {
83  m_bp_grp = domain_grp->createGroup("blueprint");
84 
85  m_named_bufs_grp = domain_grp->createGroup("named_buffers");
86 
87 #ifdef MFEM_USE_MPI
88  m_comm = MPI_COMM_NULL;
89 #endif
90 }
91 
93 {
94  attr_map.DeleteData(true);
95  if (m_owns_datastore)
96  {
97  delete m_datastore_ptr;
98  }
99 }
100 
101 #ifdef MFEM_USE_MPI
102 void SidreDataCollection::SetComm(MPI_Comm comm)
103 {
104  m_comm = comm;
105  serial = false;
106  appendRankToFileName = true;
107  MPI_Comm_rank(m_comm, &myid);
108  MPI_Comm_size(m_comm, &num_procs);
109 }
110 #endif
111 
112 // protected method
114 {
115  MFEM_ASSERT(m_named_bufs_grp != NULL,
116  "No group 'named_buffers' in data collection. Verify that"
117  " SetMesh was called to set the mesh in the data collection.");
118  return m_named_bufs_grp;
119 }
120 
121 // protected method
122 axom::sidre::View *
123 SidreDataCollection::alloc_view(axom::sidre::Group *grp,
124  const std::string &view_name)
125 {
126  MFEM_ASSERT(grp, "Group pointer is NULL");
127  sidre::View *v = NULL;
128 
129  if (! grp->hasView(view_name) )
130  {
131  v = grp->createView(view_name);
132  MFEM_ASSERT(v, "error allocating View " << view_name
133  << " in group " << grp->getPathName());
134  }
135  else
136  {
137  v = grp->getView(view_name);
138  }
139 
140  return v;
141 }
142 
143 // protected method
144 axom::sidre::View *
145 SidreDataCollection::alloc_view(axom::sidre::Group *grp,
146  const std::string &view_name,
147  const axom::sidre::DataType &dtype)
148 {
149  MFEM_ASSERT(grp, "Group pointer is NULL");
150  sidre::View *v = NULL;
151 
152  if (! grp->hasView(view_name))
153  {
154  v = grp->createView(view_name, dtype);
155  MFEM_ASSERT(v, "error allocating View " << view_name
156  << " in group " << grp->getPathName());
157  }
158  else
159  {
160  v = grp->getView(view_name);
161  MFEM_ASSERT(v->getSchema().dtype().equals(dtype), "");
162  }
163  return v;
164 }
165 
166 // protected method
167 axom::sidre::Group *
168 SidreDataCollection::alloc_group(axom::sidre::Group *grp,
169  const std::string &group_name)
170 {
171  MFEM_ASSERT(grp, "Group pointer is NULL");
172  sidre::Group *g = NULL;
173 
174  if (! grp->hasGroup(group_name) )
175  {
176  g = grp->createGroup(group_name);
177  MFEM_ASSERT(g, "error allocating Group " << group_name
178  << " in group " << grp->getPathName());
179  }
180  else
181  {
182  g = grp->getGroup(group_name);
183  }
184  return g;
185 }
186 
187 // protected method
188 std::string
189 SidreDataCollection::get_file_path(const std::string &filename) const
190 {
191  std::stringstream fNameSstr;
192 
193  // Note: If non-empty, prefix_path has a separator ('/') at the end
194  fNameSstr << prefix_path << filename;
195 
196  if (GetCycle() >= 0)
197  {
198  fNameSstr << "_" << std::setfill('0') << std::setw(pad_digits_cycle)
199  << GetCycle();
200  }
201 
202  return fNameSstr.str();
203 }
204 
205 axom::sidre::View *
206 SidreDataCollection::AllocNamedBuffer(const std::string& buffer_name,
207  axom::sidre::SidreLength sz,
208  axom::sidre::TypeID type)
209 {
210  sz = std::max(sz, sidre::SidreLength(0));
211  sidre::Group *f = named_buffers_grp();
212  sidre::View *v = NULL;
213 
214  if (! f->hasView(buffer_name) )
215  {
216  // create a buffer view
217  v = f->createViewAndAllocate(buffer_name, type, sz);
218  }
219  else
220  {
221  v = f->getView(buffer_name);
222  MFEM_ASSERT(v->getTypeID() == type, "type does not match existing type");
223 
224  // Here v is the view holding the buffer in the named_buffers group, so
225  // its size is the full size of the buffer.
226 
227  // check if we need to resize.
228  if (!v->isApplied() || v->getNumElements() < sz)
229  {
230  // resize, even if the buffer has more than 1 View.
231  // v->reallocate(sz); // this will not work for more than 1 view.
232  sidre::DataType dtype(v->getSchema().dtype());
233  dtype.set_number_of_elements(sz);
234  f->destroyViewAndData(buffer_name);
235  v = f->createViewAndAllocate(buffer_name, dtype);
236  }
237  }
238  MFEM_ASSERT(v && v->isApplied(), "allocation failed");
239  return v;
240 }
241 
242 // private method
243 void SidreDataCollection::createMeshBlueprintStubs(bool hasBP)
244 {
245  if (!hasBP)
246  {
247  m_bp_grp->createGroup("state");
248  m_bp_grp->createGroup("coordsets");
249  m_bp_grp->createGroup("topologies");
250  m_bp_grp->createGroup("fields");
251  }
252 
253  // If rank is 0, set up blueprint index state group.
254  if (myid == 0)
255  {
256  m_bp_index_grp->createGroup("state");
257  m_bp_index_grp->createGroup("coordsets");
258  m_bp_index_grp->createGroup("topologies");
259  m_bp_index_grp->createGroup("fields");
260  }
261 }
262 
263 // private method
264 void SidreDataCollection::createMeshBlueprintState(bool hasBP)
265 {
266  if (!hasBP)
267  {
268  // Set up blueprint state group.
269  m_bp_grp->createViewScalar("state/cycle", 0);
270  m_bp_grp->createViewScalar("state/time", 0.);
271  m_bp_grp->createViewScalar("state/domain", myid);
272  m_bp_grp->createViewScalar("state/time_step", 0.);
273  }
274 
275  // If rank is 0, set up blueprint index state group.
276  if (myid == 0)
277  {
278  m_bp_index_grp->createViewScalar("state/cycle", 0);
279  m_bp_index_grp->createViewScalar("state/time", 0.);
280  m_bp_index_grp->createViewScalar("state/number_of_domains", num_procs);
281  }
282 }
283 
284 // private method
285 void SidreDataCollection::createMeshBlueprintCoordset(bool hasBP)
286 {
287  int dim = mesh->SpaceDimension();
288  MFEM_ASSERT(dim >= 1 && dim <= 3, "invalid mesh dimension");
289 
290  // Assuming mfem::Vertex has the layout of a double array.
291  const int NUM_COORDS = sizeof(mfem::Vertex) / sizeof(double);
292 
293  const int num_vertices = mesh->GetNV();
294  const int coordset_len = NUM_COORDS * num_vertices;
295 
296  // Add blueprint if not present
297  if ( !hasBP )
298  {
299  m_bp_grp->createViewString("coordsets/coords/type", "explicit");
300 
301  sidre::DataType dtype =
302  sidre::DataType::c_double(num_vertices);
303  const size_t stride = dtype.stride();
304  dtype.set_stride(stride*NUM_COORDS);
305 
306  // Set up views for x, y, z values
307  sidre::View *vx, *vy = NULL, *vz = NULL;
308  vx = m_bp_grp->createView("coordsets/coords/values/x", dtype);
309 
310  if (dim >= 2)
311  {
312  dtype.set_offset(dtype.offset() + stride);
313  vy = m_bp_grp->createView("coordsets/coords/values/y", dtype);
314  }
315  if (dim >= 3)
316  {
317  dtype.set_offset(dtype.offset() + stride);
318  vz = m_bp_grp->createView("coordsets/coords/values/z", dtype);
319  }
320 
321  if (m_owns_mesh_data)
322  {
323  // Allocate buffer for coord values.
324  sidre::Buffer* coordbuf =
325  AllocNamedBuffer("vertex_coords", coordset_len)->getBuffer();
326 
327  vx->attachBuffer(coordbuf);
328  if (dim >= 2) { vy->attachBuffer(coordbuf); }
329  if (dim >= 3) { vz->attachBuffer(coordbuf); }
330  }
331  else
332  {
333  double *coordbuf = mesh->GetVertex(0);
334 
335  vx->setExternalDataPtr(coordbuf);
336  if (dim >= 2) { vy->setExternalDataPtr(coordbuf); }
337  if (dim >= 3) { vz->setExternalDataPtr(coordbuf); }
338  }
339 
340  }
341 
342  // If rank 0, set up blueprint index for coordinate set.
343  if (myid == 0)
344  {
345  m_bp_index_grp->createViewString(
346  "coordsets/coords/path", m_bp_grp->getPathName() + "/coordsets/coords");
347 
348  m_bp_index_grp->getGroup("coordsets/coords")->copyView(
349  m_bp_grp->getView("coordsets/coords/type") );
350 
351  m_bp_index_grp->createViewString(
352  "coordsets/coords/coord_system/type", "cartesian");
353 
354  // These are empty views, their existence in the group tree is used to
355  // define the number of dims
356  m_bp_index_grp->createView("coordsets/coords/coord_system/axes/x");
357 
358  if (dim >= 2)
359  {
360  m_bp_index_grp->createView("coordsets/coords/coord_system/axes/y");
361  }
362 
363  if (dim == 3)
364  {
365  m_bp_index_grp->createView("coordsets/coords/coord_system/axes/z");
366  }
367  }
368 
369  if (m_owns_mesh_data)
370  {
371  double *coord_values = GetNamedBuffer("vertex_coords")->getData();
372  // Change ownership of the mesh vertex data to sidre
373  mesh->ChangeVertexDataOwnership(coord_values, coordset_len, hasBP);
374  }
375 }
376 
377 // private method
378 void SidreDataCollection::
379 createMeshBlueprintTopologies(bool hasBP, const std::string& mesh_name)
380 {
381  const bool isBdry = (mesh_name == "boundary");
382 
383  const int num_elements = !isBdry
384  ? mesh->GetNE()
385  : mesh->GetNBE();
386 
387  const std::string mesh_topo_str = "topologies/" + mesh_name;
388  const std::string mesh_attr_str = mesh_name + "_material_attribute";
389 
390  int element_size = 0;
391  int num_indices = 0;
392  int geom = 0;
393  std::string eltTypeStr = "point";
394 
395  if (num_elements > 0)
396  {
397  mfem::Element *elem = !isBdry
398  ? mesh->GetElement(0)
399  : mesh->GetBdrElement(0);
400  element_size = elem->GetNVertices();
401 
402  num_indices = num_elements * element_size;
403 
404  // Find the element shape
405  // Note: Assumes homogeneous elements, so only check the first element
406  geom = elem->GetGeometryType();
407  eltTypeStr = getElementName(elem->GetType());
408  }
409 
410  // Create the blueprint "topology" group, if not present
411  if ( !hasBP )
412  {
413  sidre::Group* topology_grp = m_bp_grp->createGroup(mesh_topo_str);
414 
415  topology_grp->createViewString("type", "unstructured");
416  topology_grp->createViewString("elements/shape", eltTypeStr);
417  topology_grp->createViewAndAllocate(
418  "elements/connectivity", sidre::INT_ID, num_indices);
419  topology_grp->createViewString("coordset", "coords");
420 
421  // If the mesh has nodes, set the name of the GridFunction holding the
422  // mesh nodes in the blueprint group.
423  if (!isBdry && mesh->GetNodes() != NULL)
424  {
425  topology_grp->createViewString("grid_function",m_meshNodesGFName);
426  }
427  }
428 
429  // Add the mesh's attributes as an attribute field
430  RegisterAttributeField(mesh_attr_str, isBdry);
431 
432  // Change ownership or copy the element arrays into Sidre
433  if (num_elements > 0)
434  {
435  sidre::View* conn_view =
436  m_bp_grp->getGroup(mesh_topo_str)->getView("elements/connectivity");
437 
438  // The SidreDataCollection always owns these arrays:
439  Array<int> conn_array(conn_view->getData<int*>(), num_indices);
440  Array<int>* attr_array = attr_map.Get(mesh_attr_str);
441  if (!isBdry)
442  {
443  mesh->GetElementData(geom, conn_array, *attr_array);
444  }
445  else
446  {
447  mesh->GetBdrElementData(geom, conn_array, *attr_array);
448  }
449  MFEM_ASSERT(!conn_array.OwnsData(), "");
450  MFEM_ASSERT(!attr_array->OwnsData(), "");
451  }
452 
453  // If rank 0, set up blueprint index for topologies group
454  if (myid == 0)
455  {
456  const std::string m_bp_grp_path = m_bp_grp->getPathName();
457 
458  if (isBdry)
459  {
460  // "Shallow" copy the m_bp_grp view into the m_bp_index_grp sub-group.
461  // Note that the "topologies/mesh" sub-group has to exist, i.e. this
462  // method should be called first with mesh_name = "mesh".
463  m_bp_index_grp->getGroup("topologies/mesh")
464  ->copyView( m_bp_grp->getView("topologies/mesh/boundary_topology") );
465  }
466 
467  sidre::Group *bp_index_topo_grp =
468  m_bp_index_grp->createGroup(mesh_topo_str);
469  sidre::Group *topology_grp = m_bp_grp->getGroup(mesh_topo_str);
470 
471  bp_index_topo_grp->createViewString(
472  "path", m_bp_grp_path + "/" + mesh_topo_str);
473  bp_index_topo_grp->copyView( topology_grp->getView("type") );
474  bp_index_topo_grp->copyView( topology_grp->getView("coordset") );
475 
476  // If the mesh has nodes, set the name of the GridFunction holding the
477  // mesh nodes in the blueprint_index group.
478  if (!isBdry && mesh->GetNodes() != NULL)
479  {
480  bp_index_topo_grp->copyView(topology_grp->getView("grid_function"));
481  }
482  }
483 }
484 
485 // private method
486 #ifdef MFEM_USE_MPI
487 void SidreDataCollection::createMeshBlueprintAdjacencies(bool hasBP)
488 {
489  ParMesh *pmesh = dynamic_cast<ParMesh*>(mesh);
490 
491  const int GRP_SZ = 25;
492  char group_str[GRP_SZ];
493 
494  // TODO(JRC): Separate this out into group hierarchy setup and data allocation
495  // stages like all of the other "createMeshBlueprint*" functions.
496  MFEM_VERIFY(hasBP == false, "The case hasBP == true is not supported yet!");
497 
498  sidre::Group* adjset_grp = NULL;
499  if (pmesh->GetNGroups() > 1)
500  {
501  adjset_grp = m_bp_grp->createGroup("adjsets/mesh");
502  adjset_grp->createViewString("association", "vertex");
503  adjset_grp->createViewString("topology", "mesh");
504 
505  if (myid == 0) { m_bp_index_grp->createGroup("adjsets"); }
506  }
507 
508  for (int gi = 1; gi < pmesh->GetNGroups(); ++gi)
509  {
510  int num_gneighbors = pmesh->gtopo.GetGroupSize(gi);
511  int num_gvertices = pmesh->GroupNVertices(gi);
512 
513  // Skip creation of empty groups
514  if (num_gneighbors > 1 && num_gvertices > 0)
515  {
516  std::snprintf(group_str, GRP_SZ, "groups/g%d_%d",
517  pmesh->gtopo.GetGroupMasterRank(gi),
518  pmesh->gtopo.GetGroupMasterGroup(gi));
519  sidre::Group* group_grp = adjset_grp->createGroup(group_str);
520 
521  sidre::View* gneighbors_view =
522  group_grp->createViewAndAllocate(
523  "neighbors", sidre::INT_ID, num_gneighbors - 1);
524  int* gneighbors_data = gneighbors_view->getData<int*>();
525 
526  // skip local domain when adding Blueprint neighbors
527  const int* gneighbors = pmesh->gtopo.GetGroup(gi);
528  for (int ni = 0, noff = 0; ni < num_gneighbors; ++ni)
529  {
530  if ( gneighbors[ni] == 0 )
531  {
532  noff++;
533  }
534  else
535  {
536  gneighbors_data[ni - noff] =
537  pmesh->gtopo.GetNeighborRank(gneighbors[ni]);
538  }
539  }
540 
541  sidre::View* gvertices_view =
542  group_grp->createViewAndAllocate(
543  "values", sidre::INT_ID, num_gvertices);
544  int* gvertices_data = gvertices_view->getData<int*>();
545 
546  for (int vi = 0; vi < num_gvertices; ++vi)
547  {
548  gvertices_data[vi] = pmesh->GroupVertex(gi, vi);
549  }
550  }
551  }
552 }
553 #endif
554 
555 // private method
556 void SidreDataCollection::verifyMeshBlueprint()
557 {
558  // Conduit will have a verify mesh blueprint capability in the future.
559  // Add call to that when it's available to check actual contents in sidre.
560 }
561 
562 
564 {
565  // check if this rank has any boundary elements
566  int hasBndElts = mesh->GetNBE() > 0 ? 1 : 0;
567 
568 #ifdef MFEM_USE_MPI
569  // check if any rank has boundary elements
570  ParMesh *pmesh = dynamic_cast<ParMesh*>(mesh);
571  if (pmesh)
572  {
573  int hasBndElts_g;
574  MPI_Allreduce(&hasBndElts, &hasBndElts_g, 1,
575  MPI_INT, MPI_MAX,pmesh->GetComm());
576 
577  hasBndElts = hasBndElts_g;
578  }
579 #endif
580 
581  return hasBndElts > 0? true : false;
582 }
583 
585 {
586  DataCollection::SetMesh(new_mesh);
587 
588  // hasBP is used to indicate if the data currently in the blueprint should be
589  // used to replace the data in the mesh.
590  bool hasBP = m_bp_grp->getNumViews() > 0 || m_bp_grp->getNumGroups() > 0;
591  bool has_bnd_elts = HasBoundaryMesh();
592 
593  createMeshBlueprintStubs(hasBP);
594  createMeshBlueprintState(hasBP);
595  createMeshBlueprintCoordset(hasBP);
596 
597  GridFunction *nodes = new_mesh->GetNodes();
598 
599  // register the "mesh" topology in the blueprint.
600  createMeshBlueprintTopologies(hasBP, "mesh");
601 
602  if (has_bnd_elts)
603  {
604  // Set the "boundary_topology" of "mesh" to "boundary".
605  m_bp_grp->createViewString("topologies/mesh/boundary_topology", "boundary");
606 
607  // register the "boundary" topology in the blueprint.
608  createMeshBlueprintTopologies(hasBP, "boundary");
609  }
610 
611 #ifdef MFEM_USE_MPI
612  ParMesh *new_pmesh = dynamic_cast<ParMesh*>(new_mesh);
613  m_comm = new_pmesh ? new_pmesh->GetComm() : MPI_COMM_NULL;
614  if (new_pmesh)
615  {
616  createMeshBlueprintAdjacencies(hasBP);
617  }
618 #endif
619 
620  if (nodes)
621  {
622  // See the comment at the definition of 'hasBP' above.
623  if (hasBP)
624  {
625  // Get the bp mesh nodes name.
626  sidre::View *v_bp_nodes_name =
627  m_bp_grp->getView("topologies/mesh/grid_function");
628  std::string bp_nodes_name(v_bp_nodes_name->getString());
629 
630  // Check that the names match, e.g. when loading the collection.
631  MFEM_VERIFY(m_meshNodesGFName == bp_nodes_name,
632  "mismatch of requested and blueprint mesh nodes names");
633  // Support renaming bp_nodes_name --> m_meshNodesGFName ?
634  }
635 
636  if (m_owns_mesh_data)
637  {
638  // Make sure Sidre owns the data of the new_mesh's Nodes.
639  if (!GetNamedBuffer(m_meshNodesGFName))
640  {
641  int sz = new_mesh->GetNodalFESpace()->GetVSize();
642  double *gfData = AllocNamedBuffer(m_meshNodesGFName, sz)->getData();
643 
644  // See the comment at the definition of 'hasBP' above.
645  if (!hasBP)
646  {
647  MFEM_ASSERT(nodes->Size() == sz, "");
648  std::memcpy(gfData, nodes->GetData(), sizeof(double) * sz);
649  }
650  }
651  // Since the named buffer for m_meshNodesGFName exists, the call to
652  // RegisterField() below will replace the data of the nodes with the
653  // data from the named buffer.
654  }
655  else
656  {
657  // Make sure Sidre does not have a named buffer for m_meshNodesGFName.
658  MFEM_VERIFY(GetNamedBuffer(m_meshNodesGFName) == NULL, "");
659  }
660 
661  RegisterField(m_meshNodesGFName, nodes);
662 
663  if (own_data)
664  {
665  // Avoid double delete calls (for the nodes gf) when (own_data == true)
666  // and the new_mesh owns its Nodes --> take ownership from new_mesh.
667  // When new_mesh does not own its Nodes and (own_data == true), we can
668  // not take ownership --> verify that does not happen.
669  MFEM_VERIFY(new_mesh->OwnsNodes(), "mesh does not own its nodes, "
670  "can not take ownership");
671  new_mesh->SetNodesOwner(false);
672  }
673  }
674 }
675 
676 #ifdef MFEM_USE_MPI
677 void SidreDataCollection::SetMesh(MPI_Comm comm, Mesh *new_mesh)
678 {
679  // use SidreDataCollection's custom SetMesh, then set MPI info
680  SetMesh(new_mesh);
681 
682  m_comm = comm;
683  MPI_Comm_rank(comm, &myid);
684  MPI_Comm_size(comm, &num_procs);
685 }
686 #endif
687 
689 SetGroupPointers(axom::sidre::Group *bp_index_grp,
690  axom::sidre::Group *domain_grp)
691 {
692  MFEM_VERIFY(domain_grp->hasGroup("blueprint"),
693  "Domain group does not contain a blueprint group.");
694 
695  m_bp_grp = domain_grp->getGroup("blueprint");
696  m_bp_index_grp = bp_index_grp;
697  m_named_bufs_grp = domain_grp->getGroup("named_buffers");
698 }
699 
700 void SidreDataCollection::Load(const std::string& path,
701  const std::string& protocol)
702 {
704  // Reset DataStore?
705 
706 #ifdef MFEM_USE_MPI
707  if (m_comm != MPI_COMM_NULL)
708  {
709  axom::sidre::IOManager reader(m_comm);
710  reader.read(m_bp_grp->getDataStore()->getRoot(), path);
711  }
712  else
713 #endif
714  {
715  m_bp_grp->load(path, protocol);
716  }
717 
718  // If the data collection created the datastore, it knows the layout of where
719  // the domain and global groups are, and can restore them after the Load().
720  //
721  // If the data collection did not create the datastore, the host code must
722  // reset these pointers after the load operation and also reset the state
723  // variables.
724  if (m_owns_datastore)
725  {
726  SetGroupPointers(m_datastore_ptr->getRoot()->getGroup(name + "_global"),
727  m_datastore_ptr->getRoot()->getGroup(name));
728 
730  }
731 }
732 
733 void SidreDataCollection::LoadExternalData(const std::string& path)
734 {
735 #ifdef MFEM_USE_MPI
736  if (m_comm != MPI_COMM_NULL)
737  {
738  axom::sidre::IOManager reader(m_comm);
739  reader.loadExternalData(m_bp_grp->getDataStore()->getRoot(), path);
740  }
741  else
742 #endif
743  {
744  m_bp_grp->loadExternalData(path);
745  }
746 }
747 
749 {
750  SetTime( m_bp_grp->getView("state/time")->getData<double>() );
751  SetCycle( m_bp_grp->getView("state/cycle")->getData<int>() );
752  SetTimeStep( m_bp_grp->getView("state/time_step")->getData<double>() );
753 }
754 
756 {
757  m_bp_grp->getView("state/cycle")->setScalar(GetCycle());
758  m_bp_grp->getView("state/time")->setScalar(GetTime());
759  m_bp_grp->getView("state/time_step")->setScalar(GetTimeStep());
760 
761  if (myid == 0)
762  {
763  m_bp_index_grp->getView("state/cycle")->setScalar(GetCycle());
764  m_bp_index_grp->getView("state/time")->setScalar(time);
765  }
766 }
767 
769 {
770  verifyMeshBlueprint();
771  UpdateStateToDS();
772 }
773 
775 {
776  std::string filename = name;
777  std::string protocol = "sidre_hdf5";
778 
779  Save(filename, protocol);
780 }
781 
782 void SidreDataCollection::Save(const std::string& filename,
783  const std::string& protocol)
784 {
785  PrepareToSave();
786 
788 
789  std::string file_path = get_file_path(filename);
790 
791  sidre::Group * blueprint_indicies_grp = m_bp_index_grp->getParent();
792 #ifdef MFEM_USE_MPI
793  if (m_comm != MPI_COMM_NULL)
794  {
795  axom::sidre::IOManager writer(m_comm);
796  sidre::DataStore *datastore = m_bp_grp->getDataStore();
797  writer.write(datastore->getRoot(), num_procs, file_path, protocol);
798  if (myid == 0)
799  {
800  if (protocol == "sidre_hdf5")
801  {
802  writer.writeGroupToRootFile(blueprint_indicies_grp,
803  file_path + ".root");
804  }
805  // Root file support only available in hdf5.
806  else
807  {
808  writer.write(blueprint_indicies_grp, 1,
809  file_path + ".root", protocol);
810  }
811  }
812  }
813  else
814 #endif
815  {
816  // If serial, use sidre group writer.
817  m_bp_grp->save(file_path, protocol);
818 
819  blueprint_indicies_grp
820  ->save(file_path + ".root", protocol);
821  }
822 }
823 
824 // private method
825 void SidreDataCollection::
826 addScalarBasedGridFunction(const std::string &field_name, GridFunction *gf,
827  const std::string &buffer_name,
828  axom::sidre::SidreLength offset)
829 {
830  sidre::Group* grp = m_bp_grp->getGroup("fields/" + field_name);
831  MFEM_ASSERT(grp != NULL, "field " << field_name << " does not exist");
832 
833  const int numDofs = gf->FESpace()->GetVSize();
834 
835  if (gf->GetData() == NULL)
836  {
837  AllocNamedBuffer(buffer_name, offset + numDofs);
838  // gf->data is set below.
839  }
840 
841  /*
842  * Mesh blueprint for a scalar-based grid function is of the form
843  * /fields/field_name/basis
844  * -- string value is GridFunction's FEC::Name
845  * /fields/field_name/values
846  * -- array of size numDofs
847  */
848 
849  // Make sure we have the View "values".
850  sidre::View *vv = alloc_view(grp, "values");
851 
852  // Describe and apply the "values" View.
853  // If the data store has buffer for field_name (e.g. AllocNamedBuffer was
854  // called, or it was loaded from file), use that buffer.
855  sidre::View *bv = named_buffers_grp()->getView(buffer_name);
856  if (bv)
857  {
858  MFEM_ASSERT(bv->hasBuffer() && bv->isDescribed(), "");
859 
860  // named buffers always have offset 0
861  MFEM_ASSERT(bv->getSchema().dtype().offset() == 0, "");
862  MFEM_ASSERT(bv->getNumElements() >= offset + numDofs, "");
863 
864  if (vv->isEmpty())
865  {
866  vv->attachBuffer(bv->getBuffer())
867  ->apply(sidre::DOUBLE_ID, numDofs, offset);
868  }
869 
870  gf->NewDataAndSize(vv->getData(), numDofs);
871  }
872  else
873  {
874  // If we are not managing the grid function's data,
875  // create a view with the external data
876  vv->setExternalDataPtr(sidre::DOUBLE_ID, numDofs, gf->GetData());
877  }
878  MFEM_ASSERT((numDofs > 0 && vv->isApplied()) ||
879  (numDofs == 0 && vv->isEmpty() && vv->isDescribed()),
880  "invalid View state");
881  MFEM_ASSERT(numDofs == 0 || vv->getData() == gf->GetData(),
882  "View data is different from GridFunction data");
883  MFEM_ASSERT(vv->getNumElements() == numDofs,
884  "View size is different from GridFunction size");
885 }
886 
887 // private method
888 void SidreDataCollection::
889 addVectorBasedGridFunction(const std::string& field_name, GridFunction *gf,
890  const std::string &buffer_name,
891  axom::sidre::SidreLength offset)
892 {
893  sidre::Group* grp = m_bp_grp->getGroup("fields/" + field_name);
894  MFEM_ASSERT(grp != NULL, "field " << field_name << " does not exist");
895 
896  const int FLD_SZ = 20;
897  char fidxName[FLD_SZ];
898 
899  int vdim = gf->FESpace()->GetVDim();
900  int ndof = gf->FESpace()->GetNDofs();
901  Ordering::Type ordering = gf->FESpace()->GetOrdering();
902 
903  if (gf->GetData() == NULL)
904  {
905  AllocNamedBuffer(buffer_name, offset + vdim*ndof);
906  // gf->data is set below.
907  }
908 
909  /*
910  * Mesh blueprint for a vector-based grid function is of the form
911  * /fields/field_name/basis
912  * -- string value is GridFunction's FEC::Name
913  * /fields/field_name/values/x0
914  * /fields/field_name/values/x1
915  * ...
916  * /fields/field_name/values/xn
917  * -- each coordinate is an array of size ndof
918  */
919 
920  // Get/create the Group "values".
921  sidre::Group *vg = alloc_group(grp, "values");
922 
923  // Create the Views "x0", "x1", etc inside the "values" Group, vg.
924  // If we have a named buffer for field_name, attach it to the Views;
925  // otherwise set the Views to use gf->GetData() as external data.
926  sidre::DataType dtype = sidre::DataType::c_double(ndof);
927  const int entry_stride = (ordering == Ordering::byNODES ? 1 : vdim);
928  const int vdim_stride = (ordering == Ordering::byNODES ? ndof : 1);
929  dtype.set_stride(dtype.stride()*entry_stride);
930 
931  sidre::View *bv = named_buffers_grp()->getView(buffer_name);
932  if (bv)
933  {
934  MFEM_ASSERT(bv->hasBuffer() && bv->isDescribed(), "");
935 
936  // named buffers always have offset 0
937  MFEM_ASSERT(bv->getSchema().dtype().offset() == 0, "");
938  dtype.set_offset(dtype.element_bytes()*offset);
939 
940  for (int d = 0; d < vdim; d++)
941  {
942  std::snprintf(fidxName, FLD_SZ, "x%d", d);
943  sidre::View *xv = alloc_view(vg, fidxName, dtype);
944  xv->attachBuffer(bv->getBuffer());
945  dtype.set_offset(dtype.offset() + dtype.element_bytes()*vdim_stride);
946  }
947 
948  gf->NewDataAndSize(bv->getData<double*>() + offset, vdim*ndof);
949  }
950  else
951  {
952  for (int d = 0; d < vdim; d++)
953  {
954  std::snprintf(fidxName, FLD_SZ, "x%d", d);
955  sidre::View *xv = alloc_view(vg, fidxName, dtype);
956  xv->setExternalDataPtr(gf->GetData());
957  dtype.set_offset(dtype.offset() + dtype.element_bytes()*vdim_stride);
958  }
959  }
960 
961 #ifdef MFEM_DEBUG
962  for (int d = 0; d < vdim; d++)
963  {
964  std::snprintf(fidxName, FLD_SZ, "x%d", d);
965  sidre::View *xv = vg->getView(fidxName);
966  MFEM_ASSERT((ndof > 0 && xv->isApplied()) ||
967  (ndof == 0 && xv->isEmpty() && xv->isDescribed()),
968  "invalid View state");
969  MFEM_ASSERT(ndof == 0 || xv->getData() == gf->GetData() + d*vdim_stride,
970  "View data is different from GridFunction data");
971  MFEM_ASSERT(xv->getNumElements() == ndof,
972  "View size is different from GridFunction size");
973  }
974 #endif
975 }
976 
977 // private method
978 // Should only be called on mpi rank 0 ( or if serial problem ).
979 void SidreDataCollection::
980 RegisterFieldInBPIndex(const std::string& field_name, GridFunction *gf)
981 {
982  sidre::Group *bp_field_grp = m_bp_grp->getGroup("fields/" + field_name);
983  sidre::Group *bp_index_field_grp =
984  m_bp_index_grp->createGroup("fields/" + field_name);
985 
986  bp_index_field_grp->createViewString( "path", bp_field_grp->getPathName() );
987  bp_index_field_grp->copyView( bp_field_grp->getView("topology") );
988  bp_index_field_grp->copyView( bp_field_grp->getView("basis") );
989 
990  // Note: The bp index requires GridFunction::VectorDim()
991  // since the GF might be scalar valued and have a vector basis
992  // (e.g. hdiv and hcurl spaces)
993  const int number_of_components = gf->VectorDim();
994  bp_index_field_grp->createViewScalar("number_of_components",
995  number_of_components);
996 }
997 
998 // private method
999 // Should only be called on mpi rank 0 ( or if serial problem ).
1000 void SidreDataCollection::
1001 DeregisterFieldInBPIndex(const std::string& field_name)
1002 {
1003  sidre::Group * fields_grp = m_bp_index_grp->getGroup("fields");
1004  MFEM_VERIFY(fields_grp->hasGroup(field_name),
1005  "No field exists in blueprint index with name " << name);
1006 
1007  // Note: This will destroy all orphaned views or buffer classes under this
1008  // group also. If sidre owns this field data, the memory will be deleted
1009  // unless it's referenced somewhere else in sidre.
1010  fields_grp->destroyGroup(field_name);
1011 }
1012 
1013 void SidreDataCollection::RegisterField(const std::string &field_name,
1014  GridFunction *gf,
1015  const std::string &buffer_name,
1016  axom::sidre::SidreLength offset)
1017 {
1018  if ( field_name.empty() || buffer_name.empty() ||
1019  gf == NULL || gf->FESpace() == NULL )
1020  {
1021  return;
1022  }
1023 
1024  // Register field_name in the blueprint group.
1025  sidre::Group* f = m_bp_grp->getGroup("fields");
1026 
1027  if (f->hasGroup( field_name ))
1028  {
1029  // There are two possibilities:
1030  // 1. If HasField(field_name) is true - we are overwriting a field that
1031  // was previously registered.
1032  // 2. Otherwise, the field was loaded from a file, or defined outside of
1033  // the data collection.
1034  if (HasField(field_name))
1035  {
1036 #ifdef MFEM_DEBUG
1037  // Warn about overwriting field.
1038  // Skip warning when re-registering the nodal grid function
1039  if (field_name != m_meshNodesGFName)
1040  {
1041  MFEM_WARNING("field with the name '" << field_name<< "' is already "
1042  "registered, overwriting the old field");
1043  }
1044 #endif
1045  DeregisterField(field_name);
1046  }
1047  }
1048 
1049  sidre::Group* grp = f->createGroup( field_name );
1050 
1051  // Set the "basis" string using the gf's finite element space, overwrite if
1052  // necessary.
1053  sidre::View *v = alloc_view(grp, "basis");
1054  v->setString(gf->FESpace()->FEColl()->Name());
1055 
1056  // Set the topology of the GridFunction.
1057  // This is always 'mesh' except for a special case with the boundary material
1058  // attributes field.
1059  v = alloc_view(grp, "topology")->setString("mesh");
1060 
1061  MFEM_ASSERT(gf->Size() == gf->FESpace()->GetVSize(),
1062  "GridFunction size does not match FiniteElementSpace size");
1063 
1064  // Set the data views of the grid function
1065  // e.g. the number of coefficients per DoF -- either scalar-valued or
1066  // vector-valued
1067  bool const isScalarValued = (gf->FESpace()->GetVDim() == 1);
1068  if (isScalarValued)
1069  {
1070  // Set the View "<m_bp_grp>/fields/<field_name>/values"
1071  addScalarBasedGridFunction(field_name, gf, buffer_name, offset);
1072  }
1073  else // vector valued
1074  {
1075  // Set the Group "<m_bp_grp>/fields/<field_name>/values"
1076  addVectorBasedGridFunction(field_name, gf, buffer_name, offset);
1077  }
1078 
1079  // Register field_name in the blueprint_index group.
1080  if (myid == 0)
1081  {
1082  RegisterFieldInBPIndex(field_name, gf);
1083  }
1084 
1085  // Register field_name + gf in field_map.
1086  DataCollection::RegisterField(field_name, gf);
1087 }
1088 
1089 void SidreDataCollection::RegisterAttributeField(const std::string& attr_name,
1090  bool is_bdry)
1091 {
1092  MFEM_ASSERT(
1093  mesh != NULL,
1094  "Need to set mesh before registering attributes in SidreDataCollection.");
1095 
1096  // Register attr_name in the blueprint group.
1097  sidre::Group* f = m_bp_grp->getGroup("fields");
1098  if (f->hasGroup( attr_name ))
1099  {
1100  bool isAttr = attr_map.Has(attr_name);
1101  bool isFld = field_map.Has(attr_name);
1102 
1103  if (isAttr)
1104  {
1105  MFEM_WARNING("field with the name '" << attr_name << "' is already "
1106  " registered as an attribute, overwriting old values.");
1107  DeregisterAttributeField(attr_name);
1108  }
1109  else if (isFld)
1110  {
1111  MFEM_WARNING("field with the name '" << attr_name << "' is already "
1112  " registered as a field, skipping register attribute.");
1113  return;
1114  }
1115  }
1116 
1117  // Generate sidre views and groups for this mesh attribute and allocate space
1118  addIntegerAttributeField(attr_name, is_bdry);
1119 
1120  if (myid == 0)
1121  {
1122  RegisterAttributeFieldInBPIndex(attr_name);
1123  }
1124 
1125  // Register new attribute array with attr_map
1126  sidre::View* a =
1127  m_bp_grp->getGroup("fields")->getGroup(attr_name)->getView("values");
1128  Array<int>* attr = new Array<int>(a->getData<int*>(), a->getNumElements());
1129 
1130  attr_map.Register(attr_name, attr, true);
1131 }
1132 
1133 void SidreDataCollection::RegisterAttributeFieldInBPIndex(
1134  const std::string& attr_name)
1135 {
1136  const std::string m_bp_grp_path = m_bp_grp->getPathName();
1137 
1138  MFEM_ASSERT(m_bp_grp->getGroup("fields") != NULL,
1139  "Mesh blueprint does not have 'fields' group");
1140  MFEM_ASSERT(m_bp_index_grp->getGroup("fields") != NULL,
1141  "Mesh blueprint index does not have 'fields' group");
1142 
1143  // get the BP attr group
1144  sidre::Group* attr_grp =
1145  m_bp_grp->getGroup("fields")->getGroup(attr_name);
1146 
1147  // create blueprint index for this attribute
1148  sidre::Group *bp_index_attr_grp =
1149  m_bp_index_grp->getGroup("fields")->createGroup(attr_name);
1150 
1151  bp_index_attr_grp->createViewString("path", attr_grp->getPathName() );
1152  bp_index_attr_grp->copyView( attr_grp->getView("association") );
1153  bp_index_attr_grp->copyView( attr_grp->getView("topology") );
1154  bp_index_attr_grp->createViewScalar("number_of_components", 1);
1155 }
1156 
1157 void SidreDataCollection::DeregisterAttributeField(const std::string& attr_name)
1158 {
1159  attr_map.Deregister(name, true);
1160 
1161  sidre::Group * attr_grp = m_bp_grp->getGroup("fields");
1162  MFEM_VERIFY(attr_grp->hasGroup(attr_name),
1163  "No field exists in blueprint with name " << attr_name);
1164 
1165  // Delete attr_name from the blueprint group.
1166 
1167  // Note: This will destroy all orphaned views or buffer classes under this
1168  // group also. If sidre owns this field data, the memory will be deleted
1169  // unless it's referenced somewhere else in sidre.
1170  attr_grp->destroyGroup(attr_name);
1171 
1172  // Delete field_name from the blueprint_index group.
1173  if (myid == 0)
1174  {
1175  DeregisterAttributeFieldInBPIndex(attr_name);
1176  }
1177 
1178  // Delete field_name from the named_buffers group, if allocated.
1179  FreeNamedBuffer(attr_name);
1180 }
1181 
1182 void SidreDataCollection::DeregisterAttributeFieldInBPIndex(
1183  const std::string& attr_name)
1184 {
1185  sidre::Group * fields_grp = m_bp_index_grp->getGroup("fields");
1186  MFEM_VERIFY(fields_grp->hasGroup(attr_name),
1187  "No attribute exists in blueprint index with name " << attr_name);
1188 
1189  // Note: This will destroy all orphaned views or buffer classes under this
1190  // group also. If sidre owns this field data, the memory will be deleted
1191  // unless it's referenced somewhere else in sidre.
1192  fields_grp->destroyGroup(attr_name);
1193 }
1194 
1195 void SidreDataCollection::
1196 addIntegerAttributeField(const std::string& attr_name, bool is_bdry)
1197 {
1198  sidre::Group* fld_grp = m_bp_grp->getGroup("fields");
1199  MFEM_ASSERT(fld_grp != NULL, "'fields' group does not exist");
1200 
1201  const int num_elem = is_bdry? mesh->GetNBE() : mesh->GetNE();
1202  std::string topo_name = is_bdry ? "boundary" : "mesh";
1203 
1204  sidre::Group* attr_grp = fld_grp->createGroup(attr_name);
1205  attr_grp->createViewString("association", "element");
1206  attr_grp->createViewAndAllocate("values", sidre::INT_ID, num_elem);
1207  attr_grp->createViewString("topology", topo_name);
1208 }
1209 
1210 void SidreDataCollection::DeregisterField(const std::string& field_name)
1211 {
1212  // Deregister field_name from field_map.
1213  DataCollection::DeregisterField(field_name);
1214 
1215  sidre::Group * fields_grp = m_bp_grp->getGroup("fields");
1216  MFEM_VERIFY(fields_grp->hasGroup(field_name),
1217  "No field exists in blueprint with name " << field_name);
1218 
1219  // Delete field_name from the blueprint group.
1220 
1221  // Note: This will destroy all orphaned views or buffer classes under this
1222  // group also. If sidre owns this field data, the memory will be deleted
1223  // unless it's referenced somewhere else in sidre.
1224  fields_grp->destroyGroup(field_name);
1225 
1226  // Delete field_name from the blueprint_index group.
1227  if (myid == 0)
1228  {
1229  DeregisterFieldInBPIndex(field_name);
1230  }
1231 
1232  // Delete field_name from the named_buffers group, if allocated.
1233  FreeNamedBuffer(field_name);
1234 }
1235 
1236 // private method
1237 std::string SidreDataCollection::getElementName(Element::Type elementEnum)
1238 {
1239  // Note -- the mapping from Element::Type to string is based on
1240  // enum Element::Type { POINT, SEGMENT, TRIANGLE, QUADRILATERAL,
1241  // TETRAHEDRON, HEXAHEDRON };
1242  // Note: -- the string names are from conduit's blueprint
1243 
1244  switch (elementEnum)
1245  {
1246  case Element::POINT: return "point";
1247  case Element::SEGMENT: return "line";
1248  case Element::TRIANGLE: return "tri";
1249  case Element::QUADRILATERAL: return "quad";
1250  case Element::TETRAHEDRON: return "tet";
1251  case Element::HEXAHEDRON: return "hex";
1252  case Element::WEDGE:
1253  default: ;
1254  }
1255 
1256  return "unknown";
1257 }
1258 
1259 } // end namespace mfem
1260 
1261 #endif
axom::sidre::View * alloc_view(axom::sidre::Group *grp, const std::string &view_name)
Geometry::Type GetGeometryType() const
Definition: element.hpp:52
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:347
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:719
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:131
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
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.
Definition: mesh.hpp:679
const Geometry::Type geom
Definition: ex1.cpp:40
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:676
bool HasField(const std::string &name) const
Check if a grid function is part of 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.
bool Has(const std::string &fname) const
Predicate to check if a field is associated with name fname.
axom::sidre::View * AllocNamedBuffer(const std::string &buffer_name, axom::sidre::SidreLength sz, axom::sidre::TypeID type=axom::sidre::DOUBLE_ID)
Return newly allocated or existing named buffer for buffer_name.
Data type for vertex.
Definition: vertex.hpp:21
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.
axom::sidre::View * GetNamedBuffer(const std::string &buffer_name) const
Get a pointer to the sidre::View holding the named buffer for buffer_name.
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:159
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
MPI_Comm m_comm
Associated MPI communicator.
double time
Physical time (for time-dependent simulations)
Mesh * mesh
The (common) mesh for the collected fields.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
int dim
Definition: ex3.cpp:48
bool serial
Serial or parallel run? False iff mesh is a ParMesh.
Type
Ordering methods:
Definition: fespace.hpp:30
int GetCycle() const
Get time cycle (for time-dependent simulations)
void RegisterAttributeField(const std::string &name, bool is_bdry)
Registers an attribute field in the Sidre DataStore.
Type
Constants for the classes derived from Element.
Definition: element.hpp:41
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.
Definition: mesh.hpp:1033
void Register(const std::string &fname, T *field, bool own_data)
Register field field with name fname.
const Element * GetElement(int i) const
Definition: mesh.hpp:744
void SetGroupPointers(axom::sidre::Group *global_grp, axom::sidre::Group *domain_grp)
Reset the domain and global datastore group pointers.
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:336
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:430
void UpdateStateToDS()
Updates the data store&#39;s cycle, time, and time-step variables with the values from the SidreDataColle...
void SetTime(double t)
Set physical time (for time-dependent simulations)
void Deregister(const std::string &fname, bool own_data)
Unregister association between field field and name fname.
int SpaceDimension() const
Definition: mesh.hpp:714
void UpdateStateFromDS()
Updates the DataCollection&#39;s cycle, time, and time-step variables with the values from the data store...
MPI_Comm GetComm() const
Definition: pmesh.hpp:223
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.
axom::sidre::Group * alloc_group(axom::sidre::Group *grp, const std::string &group_name)
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:673
void SetTimeStep(double ts)
Set the simulation time step (for time-dependent simulations)
virtual const char * Name() const
Definition: fe_coll.hpp:53
bool OwnsNodes() const
Return the mesh nodes ownership flag.
Definition: mesh.hpp:1031
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:3639
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
Definition: mesh.hpp:729
virtual void Save()
Save the collection to file.
std::string prefix_path
A path where the directory with results is saved. If not empty, it has &#39;/&#39; at the end...
virtual int GetNVertices() const =0
void FreeNamedBuffer(const std::string &buffer_name)
Deallocate the named buffer buffer_name.
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:361
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5837
void SetComm(MPI_Comm comm)
Associate an MPI communicator with the collection.
T * Get(const std::string &fname) const
Get a pointer to the field associated with name fname.
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 ...
Definition: mesh.cpp:2778
std::string name
Name of the collection, used as a directory name when saving.
void DeleteData(bool own_data)
Clear all associations between names and fields.
int num_procs
Number of MPI ranks (in parallel)
Class for parallel meshes.
Definition: pmesh.hpp:32
Abstract data type element.
Definition: element.hpp:28
virtual Type GetType() const =0
Returns element&#39;s type.
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:748
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.cpp:7117
void DeregisterAttributeField(const std::string &name)
axom::sidre::Group * named_buffers_grp() const