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