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