MFEM  v3.3.2
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 <spio/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  MFEM_VERIFY(num_elements > 0,
389  "TODO: processors with 0 " << mesh_name << " elements");
390 
391  const int element_size = !isBdry
392  ? mesh->GetElement(0)->GetNVertices()
394 
395  const int num_indices = num_elements * element_size;
396 
397  // Find the element shape
398  // Note: Assumes homogeneous elements, so only check the first element
399  const int geom =
400  isBdry ?
403  const std::string eltTypeStr =
404  !isBdry
405  ? getElementName( static_cast<Element::Type>(
406  mesh->GetElement(0)->GetType() ) )
407  : getElementName( static_cast<Element::Type>(
408  mesh->GetBdrElement(0)->GetType() ) );
409 
410  const std::string mesh_topo_str = "topologies/" + mesh_name;
411  const std::string mesh_attr_str = "fields/"+mesh_name+"_material_attribute";
412 
413  if ( !hasBP )
414  {
415  sidre::Group* topology_grp = bp_grp->createGroup(mesh_topo_str);
416 
417  // Add mesh topology
418  topology_grp->createViewString("type", "unstructured");
419  // Note: eltTypeStr comes form the mesh
420  topology_grp->createViewString("elements/shape", eltTypeStr);
421  topology_grp->createViewAndAllocate(
422  "elements/connectivity", sidre::INT_ID, num_indices);
423  topology_grp->createViewString("coordset", "coords");
424 
425  // If the mesh has nodes, set the name of the GridFunction holding the
426  // mesh nodes in the blueprint group.
427  if (!isBdry && mesh->GetNodes() != NULL)
428  {
429  topology_grp->createViewString("grid_function",m_meshNodesGFName);
430  }
431 
432  // Add material attribute field to blueprint
433  sidre::Group* attr_grp = bp_grp->createGroup(mesh_attr_str);
434  attr_grp->createViewString("association", "element");
435  attr_grp->createViewAndAllocate("values", sidre::INT_ID, num_elements);
436  attr_grp->createViewString("topology", mesh_name);
437  }
438 
439  // If rank 0, set up blueprint index for topologies group and material
440  // attribute field.
441  if (myid == 0)
442  {
443  const std::string bp_grp_path = bp_grp->getPathName();
444 
445  // Create blueprint index for topologies.
446  if (isBdry)
447  {
448  // "Shallow" copy the bp_grp view into the bp_index_grp sub-group.
449  // Note that the "topologies/mesh" sub-group has to exist, i.e. this
450  // method should be called first with mesh_name = "mesh".
451  bp_index_grp->getGroup("topologies/mesh")
452  ->copyView( bp_grp->getView("topologies/mesh/boundary_topology") );
453  }
454 
455  sidre::Group *bp_index_topo_grp =
456  bp_index_grp->createGroup(mesh_topo_str);
457  sidre::Group *topology_grp = bp_grp->getGroup(mesh_topo_str);
458 
459  bp_index_topo_grp->createViewString(
460  "path", bp_grp_path + "/" + mesh_topo_str);
461  bp_index_topo_grp->copyView( topology_grp->getView("type") );
462  bp_index_topo_grp->copyView( topology_grp->getView("coordset") );
463 
464  // If the mesh has nodes, set the name of the GridFunction holding the
465  // mesh nodes in the blueprint_index group.
466  if (!isBdry && mesh->GetNodes() != NULL)
467  {
468  bp_index_topo_grp->copyView(topology_grp->getView("grid_function"));
469  }
470 
471  // Create blueprint index for material attributes.
472  sidre::Group *bp_index_attr_grp =
473  bp_index_grp->createGroup(mesh_attr_str);
474  sidre::Group *attr_grp = bp_grp->getGroup(mesh_attr_str);
475 
476  bp_index_attr_grp->createViewString(
477  "path", bp_grp_path + "/" + mesh_attr_str );
478  bp_index_attr_grp->copyView( attr_grp->getView("association") );
479  bp_index_attr_grp->copyView( attr_grp->getView("topology") );
480 
481  int number_of_components = 1;
482  bp_index_attr_grp->createViewScalar("number_of_components",
483  number_of_components);
484  }
485 
486  // Finally, change ownership or copy the element arrays into Sidre
487  sidre::View* conn_view =
488  bp_grp->getGroup(mesh_topo_str)->getView("elements/connectivity");
489  sidre::View* attr_view =
490  bp_grp->getGroup(mesh_attr_str)->getView("values");
491  // The SidreDataCollection always owns these arrays:
492  Array<int> conn_array(conn_view->getData<int*>(), num_indices);
493  Array<int> attr_array(attr_view->getData<int*>(), num_elements);
494  if (!isBdry)
495  {
496  mesh->GetElementData(geom, conn_array, attr_array);
497  }
498  else
499  {
500  mesh->GetBdrElementData(geom, conn_array, attr_array);
501  }
502  MFEM_ASSERT(!conn_array.OwnsData(), "");
503  MFEM_ASSERT(!attr_array.OwnsData(), "");
504 }
505 
506 // private method
507 #ifdef MFEM_USE_MPI
508 void SidreDataCollection::createMeshBlueprintAdjacencies(bool hasBP)
509 {
510  ParMesh *pmesh = dynamic_cast<ParMesh*>(mesh);
511 
512  const int GRP_SZ = 25;
513  char group_str[GRP_SZ];
514 
515  // TODO(JRC): Separate this out into group hierarchy setup and data allocation
516  // stages like all of the other "createMeshBlueprint*" functions.
517  MFEM_VERIFY(hasBP == false, "The case hasBP == true is not supported yet!");
518 
519  sidre::Group* adjset_grp = NULL;
520  if (pmesh->GetNGroups() > 1)
521  {
522  adjset_grp = bp_grp->createGroup("adjsets/mesh");
523  adjset_grp->createViewString("association", "vertex");
524  adjset_grp->createViewString("topology", "mesh");
525 
526  if (myid == 0) { bp_index_grp->createGroup("adjsets"); }
527  }
528 
529  for (int gi = 1; gi < pmesh->GetNGroups(); ++gi)
530  {
531  int num_gneighbors = pmesh->gtopo.GetGroupSize(gi);
532  int num_gvertices = pmesh->GroupNVertices(gi);
533 
534  // Skip creation of empty groups
535  if (num_gneighbors > 1 && num_gvertices > 0)
536  {
537  std::snprintf(group_str, GRP_SZ, "groups/g%d_%d",
538  pmesh->gtopo.GetGroupMasterRank(gi),
539  pmesh->gtopo.GetGroupMasterGroup(gi));
540  sidre::Group* group_grp = adjset_grp->createGroup(group_str);
541 
542  sidre::View* gneighbors_view =
543  group_grp->createViewAndAllocate(
544  "neighbors", sidre::INT_ID, num_gneighbors - 1);
545  int* gneighbors_data = gneighbors_view->getData<int*>();
546 
547  // skip local domain when adding Blueprint neighbors
548  const int* gneighbors = pmesh->gtopo.GetGroup(gi);
549  for (int ni = 0, noff = 0; ni < num_gneighbors; ++ni)
550  {
551  if ( gneighbors[ni] == 0 )
552  {
553  noff++;
554  }
555  else
556  {
557  gneighbors_data[ni - noff] =
558  pmesh->gtopo.GetNeighborRank(gneighbors[ni]);
559  }
560  }
561 
562  sidre::View* gvertices_view =
563  group_grp->createViewAndAllocate(
564  "values", sidre::INT_ID, num_gvertices);
565  int* gvertices_data = gvertices_view->getData<int*>();
566 
567  for (int vi = 0; vi < num_gvertices; ++vi)
568  {
569  gvertices_data[vi] = pmesh->GroupVertex(gi, vi);
570  }
571  }
572  }
573 }
574 #endif
575 
576 // private method
577 void SidreDataCollection::verifyMeshBlueprint()
578 {
579  // Conduit will have a verify mesh blueprint capability in the future.
580  // Add call to that when it's available to check actual contents in sidre.
581 }
582 
584 {
585  DataCollection::SetMesh(new_mesh);
586 
587  // hasBP is used to indicate if the data currently in the blueprint should be
588  // used to replace the data in the mesh.
589  bool hasBP = bp_grp->getNumViews() > 0 || bp_grp->getNumGroups() > 0;
590  bool has_bnd_elts = (new_mesh->GetNBE() > 0);
591 
592  createMeshBlueprintStubs(hasBP);
593  createMeshBlueprintState(hasBP);
594  createMeshBlueprintCoordset(hasBP);
595 
596  GridFunction *nodes = new_mesh->GetNodes();
597 
598  // register the "mesh" topology in the blueprint.
599  createMeshBlueprintTopologies(hasBP, "mesh");
600 
601  if (has_bnd_elts)
602  {
603  // Set the "boundary_topology" of "mesh" to "boundary".
604  bp_grp->createViewString("topologies/mesh/boundary_topology", "boundary");
605 
606  // register the "boundary" topology in the blueprint.
607  createMeshBlueprintTopologies(hasBP, "boundary");
608  }
609 
610 #ifdef MFEM_USE_MPI
611  ParMesh *new_pmesh = dynamic_cast<ParMesh*>(new_mesh);
612  m_comm = new_pmesh ? new_pmesh->GetComm() : MPI_COMM_NULL;
613  if (new_pmesh)
614  {
615  createMeshBlueprintAdjacencies(hasBP);
616  }
617 #endif
618 
619  if (nodes)
620  {
621  // See the comment at the definition of 'hasBP' above.
622  if (hasBP)
623  {
624  // Get the bp mesh nodes name.
625  sidre::View *v_bp_nodes_name =
626  bp_grp->getView("topologies/mesh/grid_function");
627  std::string bp_nodes_name(v_bp_nodes_name->getString());
628 
629  // Check that the names match, e.g. when loading the collection.
630  MFEM_VERIFY(m_meshNodesGFName == bp_nodes_name,
631  "mismatch of requested and blueprint mesh nodes names");
632  // Support renaming bp_nodes_name --> m_meshNodesGFName ?
633  }
634 
635  if (m_owns_mesh_data)
636  {
637  // Make sure Sidre owns the data of the new_mesh's Nodes.
638  if (!GetNamedBuffer(m_meshNodesGFName))
639  {
640  int sz = new_mesh->GetNodalFESpace()->GetVSize();
641  double *gfData = AllocNamedBuffer(m_meshNodesGFName, sz)->getData();
642 
643  // See the comment at the definition of 'hasBP' above.
644  if (!hasBP)
645  {
646  MFEM_ASSERT(nodes->Size() == sz, "");
647  std::memcpy(gfData, nodes->GetData(), sizeof(double) * sz);
648  }
649  }
650  // Since the named buffer for m_meshNodesGFName exists, the call to
651  // RegisterField() below will replace the data of the nodes with the
652  // data from the named buffer.
653  }
654  else
655  {
656  // Make sure Sidre does not have a named buffer for m_meshNodesGFName.
657  MFEM_VERIFY(GetNamedBuffer(m_meshNodesGFName) == NULL, "");
658  }
659 
660  RegisterField(m_meshNodesGFName, nodes);
661 
662  if (own_data)
663  {
664  // Avoid double delete calls (for the nodes gf) when (own_data == true)
665  // and the new_mesh owns its Nodes --> take ownership from new_mesh.
666  // When new_mesh does not own its Nodes and (own_data == true), we can
667  // not take ownership --> verify that does not happen.
668  MFEM_VERIFY(new_mesh->OwnsNodes(), "mesh does not own its nodes, "
669  "can not take ownership");
670  new_mesh->SetNodesOwner(false);
671  }
672  }
673 }
674 
676 SetGroupPointers(axom::sidre::Group *global_grp,
677  axom::sidre::Group *domain_grp)
678 {
679  MFEM_VERIFY(domain_grp->hasGroup("blueprint"),
680  "Domain group does not contain a blueprint group.");
681  MFEM_VERIFY(global_grp->hasGroup("blueprint_index/" + name),
682  "Global group does not contain a blueprint index group.");
683 
684  bp_grp = domain_grp->getGroup("blueprint");
685  bp_index_grp = global_grp->getGroup("blueprint_index/" + name);
686  named_bufs_grp = domain_grp->getGroup("named_buffers");
687 }
688 
689 void SidreDataCollection::Load(const std::string& path,
690  const std::string& protocol)
691 {
693  // Reset DataStore?
694 
695 #ifdef MFEM_USE_MPI
696  if (m_comm != MPI_COMM_NULL)
697  {
698  axom::spio::IOManager reader(m_comm);
699  reader.read(bp_grp->getDataStore()->getRoot(), path);
700  }
701  else
702 #endif
703  {
704  bp_grp->load(path, protocol);
705  }
706 
707  // If the data collection created the datastore, it knows the layout of where
708  // the domain and global groups are, and can restore them after the Load().
709  //
710  // If the data collection did not create the datastore, the host code must
711  // reset these pointers after the load operation and also reset the state
712  // variables.
713  if (m_owns_datastore)
714  {
715  SetGroupPointers(m_datastore_ptr->getRoot()->getGroup(name + "_global"),
716  m_datastore_ptr->getRoot()->getGroup(name));
717 
719  }
720 }
721 
722 void SidreDataCollection::LoadExternalData(const std::string& path)
723 {
724 #ifdef MFEM_USE_MPI
725  if (m_comm != MPI_COMM_NULL)
726  {
727  axom::spio::IOManager reader(m_comm);
728  reader.loadExternalData(bp_grp->getDataStore()->getRoot(), path);
729  }
730  else
731 #endif
732  {
733  bp_grp->loadExternalData(path);
734  }
735 }
736 
738 {
739  SetTime( bp_grp->getView("state/time")->getData<double>() );
740  SetCycle( bp_grp->getView("state/cycle")->getData<int>() );
741  SetTimeStep( bp_grp->getView("state/time_step")->getData<double>() );
742 }
743 
745 {
746  bp_grp->getView("state/cycle")->setScalar(GetCycle());
747  bp_grp->getView("state/time")->setScalar(GetTime());
748  bp_grp->getView("state/time_step")->setScalar(GetTimeStep());
749 
750  if (myid == 0)
751  {
752  bp_index_grp->getView("state/cycle")->setScalar(GetCycle());
753  bp_index_grp->getView("state/time")->setScalar(time);
754  }
755 }
756 
758 {
759  verifyMeshBlueprint();
760  UpdateStateToDS();
761 }
762 
764 {
765  std::string filename = name;
766  std::string protocol = "sidre_hdf5";
767 
768  Save(filename, protocol);
769 }
770 
771 void SidreDataCollection::Save(const std::string& filename,
772  const std::string& protocol)
773 {
774  PrepareToSave();
775 
777 
778  std::string file_path = get_file_path(filename);
779 
780  sidre::Group * blueprint_indicies_grp = bp_index_grp->getParent();
781 #ifdef MFEM_USE_MPI
782  if (m_comm != MPI_COMM_NULL)
783  {
784  axom::spio::IOManager writer(m_comm);
785  sidre::DataStore *datastore = bp_grp->getDataStore();
786  writer.write(datastore->getRoot(), num_procs, file_path, protocol);
787  if (myid == 0)
788  {
789  if (protocol == "sidre_hdf5")
790  {
791  writer.writeGroupToRootFile(blueprint_indicies_grp,
792  file_path + ".root");
793  }
794  // Root file support only available in hdf5.
795  else
796  {
797  writer.write(blueprint_indicies_grp, 1,
798  file_path + ".root", protocol);
799  }
800  }
801  }
802  else
803 #endif
804  {
805  // If serial, use sidre group writer.
806  bp_grp->save(file_path, protocol);
807 
808  blueprint_indicies_grp
809  ->save(file_path + ".root", protocol);
810  }
811 }
812 
813 // private method
814 void SidreDataCollection::
815 addScalarBasedGridFunction(const std::string &field_name, GridFunction *gf,
816  const std::string &buffer_name,
817  axom::sidre::SidreLength offset)
818 {
819  sidre::Group* grp = bp_grp->getGroup("fields/" + field_name);
820  MFEM_ASSERT(grp != NULL, "field " << field_name << " does not exist");
821 
822  const int numDofs = gf->FESpace()->GetVSize();
823 
824  if (gf->GetData() == NULL)
825  {
826  AllocNamedBuffer(buffer_name, offset + numDofs);
827  // gf->data is set below.
828  }
829 
830  /*
831  * Mesh blueprint for a scalar-based grid function is of the form
832  * /fields/field_name/basis
833  * -- string value is GridFunction's FEC::Name
834  * /fields/field_name/values
835  * -- array of size numDofs
836  */
837 
838  // Make sure we have the View "values".
839  sidre::View *vv = alloc_view(grp, "values");
840 
841  // Describe and apply the "values" View.
842  // If the data store has buffer for field_name (e.g. AllocNamedBuffer was
843  // called, or it was loaded from file), use that buffer.
844  sidre::View *bv = named_buffers_grp()->getView(buffer_name);
845  if (bv)
846  {
847  MFEM_ASSERT(bv->hasBuffer() && bv->isDescribed(), "");
848 
849  // named buffers always have offset 0
850  MFEM_ASSERT(bv->getSchema().dtype().offset() == 0, "");
851  MFEM_ASSERT(bv->getNumElements() >= offset + numDofs, "");
852 
853  if (vv->isEmpty())
854  {
855  vv->attachBuffer(bv->getBuffer())
856  ->apply(sidre::DOUBLE_ID, numDofs, offset);
857  }
858 
859  gf->NewDataAndSize(vv->getData(), numDofs);
860  }
861  else
862  {
863  // If we are not managing the grid function's data,
864  // create a view with the external data
865  vv->setExternalDataPtr(sidre::DOUBLE_ID, numDofs, gf->GetData());
866  }
867  MFEM_ASSERT((numDofs > 0 && vv->isApplied()) ||
868  (numDofs == 0 && vv->isEmpty() && vv->isDescribed()),
869  "invalid View state");
870  MFEM_ASSERT(numDofs == 0 || vv->getData() == gf->GetData(),
871  "View data is different from GridFunction data");
872  MFEM_ASSERT(vv->getNumElements() == numDofs,
873  "View size is different from GridFunction size");
874 }
875 
876 // private method
877 void SidreDataCollection::
878 addVectorBasedGridFunction(const std::string& field_name, GridFunction *gf,
879  const std::string &buffer_name,
880  axom::sidre::SidreLength offset)
881 {
882  sidre::Group* grp = bp_grp->getGroup("fields/" + field_name);
883  MFEM_ASSERT(grp != NULL, "field " << field_name << " does not exist");
884 
885  const int FLD_SZ = 20;
886  char fidxName[FLD_SZ];
887 
888  int vdim = gf->FESpace()->GetVDim();
889  int ndof = gf->FESpace()->GetNDofs();
890  Ordering::Type ordering = gf->FESpace()->GetOrdering();
891 
892  if (gf->GetData() == NULL)
893  {
894  AllocNamedBuffer(buffer_name, offset + vdim*ndof);
895  // gf->data is set below.
896  }
897 
898  /*
899  * Mesh blueprint for a vector-based grid function is of the form
900  * /fields/field_name/basis
901  * -- string value is GridFunction's FEC::Name
902  * /fields/field_name/values/x0
903  * /fields/field_name/values/x1
904  * ...
905  * /fields/field_name/values/xn
906  * -- each coordinate is an array of size ndof
907  */
908 
909  // Get/create the Group "values".
910  sidre::Group *vg = alloc_group(grp, "values");
911 
912  // Create the Views "x0", "x1", etc inside the "values" Group, vg.
913  // If we have a named buffer for field_name, attach it to the Views;
914  // otherwise set the Views to use gf->GetData() as external data.
915  sidre::DataType dtype = sidre::DataType::c_double(ndof);
916  const int entry_stride = (ordering == Ordering::byNODES ? 1 : vdim);
917  const int vdim_stride = (ordering == Ordering::byNODES ? ndof : 1);
918  dtype.set_stride(dtype.stride()*entry_stride);
919 
920  sidre::View *bv = named_buffers_grp()->getView(buffer_name);
921  if (bv)
922  {
923  MFEM_ASSERT(bv->hasBuffer() && bv->isDescribed(), "");
924 
925  // named buffers always have offset 0
926  MFEM_ASSERT(bv->getSchema().dtype().offset() == 0, "");
927  dtype.set_offset(dtype.element_bytes()*offset);
928 
929  for (int d = 0; d < vdim; d++)
930  {
931  std::snprintf(fidxName, FLD_SZ, "x%d", d);
932  sidre::View *xv = alloc_view(vg, fidxName, dtype);
933  xv->attachBuffer(bv->getBuffer());
934  dtype.set_offset(dtype.offset() + dtype.element_bytes()*vdim_stride);
935  }
936 
937  gf->NewDataAndSize(bv->getData<double*>() + offset, vdim*ndof);
938  }
939  else
940  {
941  for (int d = 0; d < vdim; d++)
942  {
943  std::snprintf(fidxName, FLD_SZ, "x%d", d);
944  sidre::View *xv = alloc_view(vg, fidxName, dtype);
945  xv->setExternalDataPtr(gf->GetData());
946  dtype.set_offset(dtype.offset() + dtype.element_bytes()*vdim_stride);
947  }
948  }
949 
950 #ifdef MFEM_DEBUG
951  for (int d = 0; d < vdim; d++)
952  {
953  std::snprintf(fidxName, FLD_SZ, "x%d", d);
954  sidre::View *xv = vg->getView(fidxName);
955  MFEM_ASSERT((ndof > 0 && xv->isApplied()) ||
956  (ndof == 0 && xv->isEmpty() && xv->isDescribed()),
957  "invalid View state");
958  MFEM_ASSERT(ndof == 0 || xv->getData() == gf->GetData() + d*vdim_stride,
959  "View data is different from GridFunction data");
960  MFEM_ASSERT(xv->getNumElements() == ndof,
961  "View size is different from GridFunction size");
962  }
963 #endif
964 }
965 
966 // private method
967 // Should only be called on mpi rank 0 ( or if serial problem ).
968 void SidreDataCollection::
969 RegisterFieldInBPIndex(const std::string& field_name, GridFunction *gf)
970 {
971  sidre::Group *bp_field_grp = bp_grp->getGroup("fields/" + field_name);
972  sidre::Group *bp_index_field_grp =
973  bp_index_grp->createGroup("fields/" + field_name);
974 
975  bp_index_field_grp->createViewString( "path", bp_field_grp->getPathName() );
976  bp_index_field_grp->copyView( bp_field_grp->getView("topology") );
977  bp_index_field_grp->copyView( bp_field_grp->getView("basis") );
978 
979  // Note: The bp index requires GridFunction::VectorDim()
980  // since the GF might be scalar valued and have a vector basis
981  // (e.g. hdiv and hcurl spaces)
982  const int number_of_components = gf->VectorDim();
983  bp_index_field_grp->createViewScalar("number_of_components",
984  number_of_components);
985 }
986 
987 // private method
988 // Should only be called on mpi rank 0 ( or if serial problem ).
989 void SidreDataCollection::
990 DeregisterFieldInBPIndex(const std::string& field_name)
991 {
992  sidre::Group * fields_grp = bp_index_grp->getGroup("fields");
993  MFEM_VERIFY(fields_grp->hasGroup(field_name),
994  "No field exists in blueprint index with name " << name);
995 
996  // Note: This will destroy all orphaned views or buffer classes under this
997  // group also. If sidre owns this field data, the memory will be deleted
998  // unless it's referenced somewhere else in sidre.
999  fields_grp->destroyGroup(field_name);
1000 }
1001 
1002 void SidreDataCollection::RegisterField(const std::string &field_name,
1003  GridFunction *gf,
1004  const std::string &buffer_name,
1005  axom::sidre::SidreLength offset)
1006 {
1007  if ( field_name.empty() || buffer_name.empty() ||
1008  gf == NULL || gf->FESpace() == NULL )
1009  {
1010  return;
1011  }
1012 
1013  // Register field_name in the blueprint group.
1014  sidre::Group* f = bp_grp->getGroup("fields");
1015 
1016  if (f->hasGroup( field_name ))
1017  {
1018  // There are two possibilities:
1019  // 1. If HasField(field_name) is true - we are overwriting a field that
1020  // was previously registered.
1021  // 2. Otherwise, the field was loaded from a file, or defined outside of
1022  // the data collection.
1023  if (HasField(field_name))
1024  {
1025 #ifdef MFEM_DEBUG
1026  // Warn about overwriting field.
1027  // Skip warning when re-registering the nodal grid function
1028  if (field_name != m_meshNodesGFName)
1029  {
1030  MFEM_WARNING("field with the name '" << field_name<< "' is already "
1031  "registered, overwriting the old field");
1032  }
1033 #endif
1034  DeregisterField(field_name);
1035  }
1036  }
1037 
1038  sidre::Group* grp = f->createGroup( field_name );
1039 
1040  // Set the "basis" string using the gf's finite element space, overwrite if
1041  // necessary.
1042  sidre::View *v = alloc_view(grp, "basis");
1043  v->setString(gf->FESpace()->FEColl()->Name());
1044 
1045  // Set the topology of the GridFunction.
1046  // This is always 'mesh' except for a special case with the boundary material
1047  // attributes field.
1048  v = alloc_view(grp, "topology")->setString("mesh");
1049 
1050  MFEM_ASSERT(gf->Size() == gf->FESpace()->GetVSize(),
1051  "GridFunction size does not match FiniteElementSpace size");
1052 
1053  // Set the data views of the grid function
1054  // e.g. the number of coefficients per DoF -- either scalar-valued or
1055  // vector-valued
1056  bool const isScalarValued = (gf->FESpace()->GetVDim() == 1);
1057  if (isScalarValued)
1058  {
1059  // Set the View "<bp_grp>/fields/<field_name>/values"
1060  addScalarBasedGridFunction(field_name, gf, buffer_name, offset);
1061  }
1062  else // vector valued
1063  {
1064  // Set the Group "<bp_grp>/fields/<field_name>/values"
1065  addVectorBasedGridFunction(field_name, gf, buffer_name, offset);
1066  }
1067 
1068  // Register field_name in the blueprint_index group.
1069  if (myid == 0)
1070  {
1071  RegisterFieldInBPIndex(field_name, gf);
1072  }
1073 
1074  // Register field_name + gf in field_map.
1075  DataCollection::RegisterField(field_name, gf);
1076 }
1077 
1078 void SidreDataCollection::DeregisterField(const std::string& field_name)
1079 {
1080  // Deregister field_name from field_map.
1081  DataCollection::DeregisterField(field_name);
1082 
1083  sidre::Group * fields_grp = bp_grp->getGroup("fields");
1084  MFEM_VERIFY(fields_grp->hasGroup(field_name),
1085  "No field exists in blueprint with name " << field_name);
1086 
1087  // Delete field_name from the blueprint group.
1088 
1089  // Note: This will destroy all orphaned views or buffer classes under this
1090  // group also. If sidre owns this field data, the memory will be deleted
1091  // unless it's referenced somewhere else in sidre.
1092  fields_grp->destroyGroup(field_name);
1093 
1094  // Delete field_name from the blueprint_index group.
1095  if (myid == 0)
1096  {
1097  DeregisterFieldInBPIndex(field_name);
1098  }
1099 
1100  // Delete field_name from the named_buffers group, if allocated.
1101  FreeNamedBuffer(field_name);
1102 }
1103 
1104 // private method
1105 std::string SidreDataCollection::getElementName(Element::Type elementEnum)
1106 {
1107  // Note -- the mapping from Element::Type to string is based on
1108  // enum Element::Type { POINT, SEGMENT, TRIANGLE, QUADRILATERAL,
1109  // TETRAHEDRON, HEXAHEDRON };
1110  // Note: -- the string names are from conduit's blueprint
1111 
1112  switch (elementEnum)
1113  {
1114  case Element::POINT: return "point";
1115  case Element::SEGMENT: return "line";
1116  case Element::TRIANGLE: return "tri";
1117  case Element::QUADRILATERAL: return "quad";
1118  case Element::TETRAHEDRON: return "tet";
1119  case Element::HEXAHEDRON: return "hex";
1120  }
1121 
1122  return "unknown";
1123 }
1124 
1125 } // end namespace mfem
1126 
1127 #endif
axom::sidre::View * alloc_view(axom::sidre::Group *grp, const std::string &view_name)
int GetVSize() const
Definition: fespace.hpp:163
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:647
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:101
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:617
const Geometry::Type geom
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:113
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:614
bool HasField(const std::string &name) const
Check if a grid function is part of the collection.
virtual void DeregisterField(const std::string &field_name)
Remove a grid function from the collection.
SidreDataCollection(const std::string &collection_name, Mesh *the_mesh=NULL, bool owns_mesh_data=false)
Constructor that allocates and initializes a Sidre DataStore.
void LoadExternalData(const std::string &path)
Load external data after registering externally owned fields.
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
Definition: vector.hpp:121
double GetTimeStep() const
Get the simulation time step (for time-dependent simulations)
double time
Physical time (for time-dependent simulations)
Mesh * mesh
The (common) mesh for the collected fields.
int dim
Definition: ex3.cpp:47
bool serial
Serial or parallel run?
int GetCycle() const
Get time cycle (for time-dependent simulations)
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:912
const Element * GetElement(int i) const
Definition: mesh.hpp:672
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:153
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:309
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)
int SpaceDimension() const
Definition: mesh.hpp:642
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:117
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:684
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:611
void SetTimeStep(double ts)
Set the simulation time step (for time-dependent simulations)
virtual const char * Name() const
Definition: fe_coll.hpp:117
bool OwnsNodes() const
Return the mesh nodes ownership flag.
Definition: mesh.hpp:910
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:3310
void GetBdrElementData(int geom, Array< int > &bdr_elem_vtx, Array< int > &bdr_attr) const
Definition: mesh.hpp:657
virtual void Save()
Save the collection to file.
virtual int GetNVertices() const =0
void FreeNamedBuffer(const std::string &buffer_name)
Deallocate the named buffer buffer_name.
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:177
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5401
void SetComm(MPI_Comm comm)
Associate an MPI communicator with the collection.
virtual void SetMesh(Mesh *new_mesh)
Set/change the mesh associated with the collection.
void ChangeVertexDataOwnership(double *vertices, int len_vertices, bool zerocopy=false)
Set the internal Vertex array to point to the given vertices array without assuming ownership of the ...
Definition: mesh.cpp:2424
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
std::string name
Name of the collection, used as a directory name when saving.
int num_procs
Number of MPI ranks (in parallel)
Class for parallel meshes.
Definition: pmesh.hpp:29
int GetBdrElementBaseGeometry(int i=0) const
Definition: mesh.hpp:687
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:676
void GetElementData(const Array< Element * > &elem_array, int geom, Array< int > &elem_vtx, Array< int > &attr) const
Definition: mesh.cpp:6252
virtual int GetType() const =0
Returns element&#39;s type.
axom::sidre::Group * named_buffers_grp() const