MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
conduitdatacollection.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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#include "../config/config.hpp"
13
14#ifdef MFEM_USE_CONDUIT
15
16#include "fem.hpp"
17#include "../general/text.hpp"
18#include <conduit_relay.hpp>
19#include <conduit_blueprint.hpp>
20
21#include <string>
22#include <sstream>
23
24using namespace conduit;
25
26namespace mfem
27{
28
29//---------------------------------------------------------------------------//
30// class ConduitDataCollection implementation
31//---------------------------------------------------------------------------//
32
33//------------------------------
34// begin public methods
35//------------------------------
36
37//---------------------------------------------------------------------------//
39 Mesh *mesh)
40 : DataCollection(coll_name, mesh),
41 relay_protocol("hdf5")
42{
43 appendRankToFileName = true; // always include rank in file names
44 cycle = 0; // always include cycle in directory names
45}
46
47#ifdef MFEM_USE_MPI
48//---------------------------------------------------------------------------//
50 const std::string& coll_name,
51 Mesh *mesh)
52 : DataCollection(coll_name, mesh),
53 relay_protocol("hdf5")
54{
55 m_comm = comm;
56 MPI_Comm_rank(comm, &myid);
57 MPI_Comm_size(comm, &num_procs);
58 appendRankToFileName = true; // always include rank in file names
59 cycle = 0; // always include cycle in directory names
60}
61#endif
62
63//---------------------------------------------------------------------------//
68
69//---------------------------------------------------------------------------//
71{
72 std::string dir_name = MeshDirectoryName();
73 int err_ = create_directory(dir_name, mesh, myid);
74 if (err_)
75 {
76 MFEM_ABORT("Error creating directory: " << dir_name);
77 }
78
79 Node n_mesh;
80 // future? If moved into Mesh class
81 // mesh->toConduitBlueprint(n_mesh);
83
84 Node verify_info;
85 if (!blueprint::mesh::verify(n_mesh,verify_info))
86 {
87 MFEM_ABORT("Conduit Mesh Blueprint Verify Failed:\n"
88 << verify_info.to_json());
89 }
90
91 // wrap all grid functions
93 for ( itr = field_map.begin(); itr != field_map.end(); itr++)
94 {
95 std::string name = itr->first;
96 GridFunction *gf = itr->second;
97 // don't save mesh nodes twice ...
98 if ( gf != mesh->GetNodes())
99 {
100 // future? If moved into GridFunction class
101 //gf->toConduitBlueprint(n_mesh["fields"][it->first]);
103 n_mesh["fields"][name]);
104 }
105 }
106
107 // wrap all quadrature functions
109 for ( qf_itr = q_field_map.begin(); qf_itr != q_field_map.end(); qf_itr++)
110 {
111 std::string name = qf_itr->first;
112 QuadratureFunction *qf = qf_itr->second;
114 n_mesh["fields"][name]);
115 }
116
117 // save mesh data
119 n_mesh,
121
122 if (myid == 0)
123 {
124 // save root file
126 n_mesh,
128 }
129}
130
131//---------------------------------------------------------------------------//
133{
134 DeleteAll();
135 this->cycle = cycle;
136
137 // Note: We aren't currently using much info from the root file ...
138 // with cycle, we can use implicit mfem conduit file layout
139
140 Node n_root;
141 LoadRootFile(n_root);
142 relay_protocol = n_root["protocol/name"].as_string();
143
144 // for MPI case, we assume that we have # of mpi tasks
145 // == number of domains
146
147 int num_domains = n_root["number_of_trees"].to_int();
148
149 if (num_procs != num_domains)
150 {
152 MFEM_WARNING("num_procs must equal num_domains");
153 return;
154 }
155
156 // load the mesh and fields
158
159 // TODO: am I properly wielding this?
160 own_data = true;
161
162}
163
164//---------------------------------------------------------------------------//
165void
166ConduitDataCollection::SetProtocol(const std::string &protocol)
167{
168 relay_protocol = protocol;
169}
170
171// Conduit data type id for the MFEM precision
172constexpr conduit::index_t mfem_precision_conduit_id =
173#if defined(MFEM_USE_DOUBLE)
174 CONDUIT_NATIVE_DOUBLE_ID;
175#elif defined(MFEM_USE_SINGLE)
177#else
178#error Unknown MFEM precision
179#endif
180
181//------------------------------
182// begin static public methods
183//------------------------------
184
185//---------------------------------------------------------------------------//
188 const std::string &main_toplogy_name,
189 bool zero_copy)
190{
191 // n_conv holds converted data (when necessary for mfem api)
192 // if n_conv is used ( !n_conv.dtype().empty() ) we
193 // now that some data allocation was necessary, so we
194 // can't return a mesh that zero copies the conduit data
195 Node n_conv;
196
197 //
198 // we need to find the topology and its coordset.
199 //
200
201 std::string topo_name = main_toplogy_name;
202 // if topo name is not set, look for first topology
203 if (topo_name == "")
204 {
205 topo_name = n_mesh["topologies"].schema().child_name(0);
206 }
207
208 MFEM_ASSERT(n_mesh.has_path("topologies/" + topo_name),
209 "Expected topology named \"" + topo_name + "\" "
210 "(node is missing path \"topologies/" + topo_name + "\")");
211
212 // find the coord set
213 std::string coords_name =
214 n_mesh["topologies"][topo_name]["coordset"].as_string();
215
216
217 MFEM_ASSERT(n_mesh.has_path("coordsets/" + coords_name),
218 "Expected topology named \"" + coords_name + "\" "
219 "(node is missing path \"coordsets/" + coords_name + "\")");
220
221 const Node &n_coordset = n_mesh["coordsets"][coords_name];
222 const Node &n_coordset_vals = n_coordset["values"];
223
224 // get the number of dims of the coordset
225 int ndims = n_coordset_vals.number_of_children();
226
227 // get the number of points
228 int num_verts = n_coordset_vals[0].dtype().number_of_elements();
229 // get vals for points
230 const real_t *verts_ptr = NULL;
231
232 // the mfem mesh constructor needs coords with interleaved (aos) type
233 // ordering, even for 1d + 2d we always need 3 real_t (double/float) b/c it
234 // uses Array<Vertex> and Vertex is a pod of 3 real_t. we check for this
235 // case, if we don't have it we convert the data
236
237 if (ndims == 3 &&
238 n_coordset_vals[0].dtype().id() == mfem_precision_conduit_id &&
239 blueprint::mcarray::is_interleaved(n_coordset_vals) )
240 {
241 // already interleaved mcarray of 3 real_t (double/float),
242 // return ptr to beginning
243 verts_ptr = n_coordset_vals[0].value();
244 }
245 else
246 {
247 Node n_tmp;
248 // check all vals, if we don't have real_t (double/float) convert
249 // to real_t
250 NodeConstIterator itr = n_coordset_vals.children();
251 while (itr.has_next())
252 {
253 const Node &c_vals = itr.next();
254 std::string c_name = itr.name();
255
256 if ( c_vals.dtype().id() == mfem_precision_conduit_id )
257 {
258 // zero copy current coords
259 n_tmp[c_name].set_external(c_vals);
260 }
261 else
262 {
263 // convert
264 c_vals.to_data_type(mfem_precision_conduit_id, n_tmp[c_name]);
265 }
266 }
267
268 // check if we need to add extra dims to get
269 // proper interleaved array
270 if (ndims < 3)
271 {
272 // add dummy z
273 n_tmp["z"].set(DataType(mfem_precision_conduit_id, num_verts));
274 }
275
276 if (ndims < 2)
277 {
278 // add dummy y
279 n_tmp["y"].set(DataType(mfem_precision_conduit_id, num_verts));
280 }
281
282 Node &n_conv_coords_vals = n_conv["coordsets"][coords_name]["values"];
283 blueprint::mcarray::to_interleaved(n_tmp,
284 n_conv_coords_vals);
285 verts_ptr = n_conv_coords_vals[0].value();
286 }
287
288
289
290 const Node &n_mesh_topo = n_mesh["topologies"][topo_name];
291 std::string mesh_ele_shape = n_mesh_topo["elements/shape"].as_string();
292
293 mfem::Geometry::Type mesh_geo = ShapeNameToGeomType(mesh_ele_shape);
294 int num_idxs_per_ele = Geometry::NumVerts[mesh_geo];
295
296 const Node &n_mesh_conn = n_mesh_topo["elements/connectivity"];
297
298 const int *elem_indices = NULL;
299 // mfem requires ints, we could have int64s, etc convert if necessary
300 if (n_mesh_conn.dtype().is_int() &&
301 n_mesh_conn.is_compact() )
302 {
303 elem_indices = n_mesh_topo["elements/connectivity"].value();
304 }
305 else
306 {
307 Node &n_mesh_conn_conv=
308 n_conv["topologies"][topo_name]["elements/connectivity"];
309 n_mesh_conn.to_int_array(n_mesh_conn_conv);
310 elem_indices = n_mesh_conn_conv.value();
311 }
312
313 int num_mesh_ele =
314 n_mesh_topo["elements/connectivity"].dtype().number_of_elements();
315 num_mesh_ele = num_mesh_ele / num_idxs_per_ele;
316
317
318 const int *bndry_indices = NULL;
319 int num_bndry_ele = 0;
320 // init to something b/c the mesh constructor will use this for a
321 // table lookup, even if we don't have boundary info.
323
324 if ( n_mesh_topo.has_child("boundary_topology") )
325 {
326 std::string bndry_topo_name = n_mesh_topo["boundary_topology"].as_string();
327
328 // In VisIt, we encountered a case were a mesh specified a boundary
329 // topology, but the boundary topology was omitted from the blueprint
330 // index, so it's data could not be obtained.
331 //
332 // This guard prevents an error in that case, allowing the mesh to be
333 // created without boundary info
334
335 if (n_mesh["topologies"].has_child(bndry_topo_name))
336 {
337 const Node &n_bndry_topo = n_mesh["topologies"][bndry_topo_name];
338 std::string bndry_ele_shape = n_bndry_topo["elements/shape"].as_string();
339
340 bndry_geo = ShapeNameToGeomType(bndry_ele_shape);
341 int num_idxs_per_bndry_ele = Geometry::NumVerts[mesh_geo];
342
343 const Node &n_bndry_conn = n_bndry_topo["elements/connectivity"];
344
345 // mfem requires ints, we could have int64s, etc convert if necessary
346 if ( n_bndry_conn.dtype().is_int() &&
347 n_bndry_conn.is_compact())
348 {
349 bndry_indices = n_bndry_conn.value();
350 }
351 else
352 {
353 Node &n_bndry_conn_conv =
354 n_conv["topologies"][bndry_topo_name]["elements/connectivity"];
355 n_bndry_conn.to_int_array(n_bndry_conn_conv);
356 bndry_indices = (n_bndry_conn_conv).value();
357
358 }
359
360 num_bndry_ele =
361 n_bndry_topo["elements/connectivity"].dtype().number_of_elements();
362 num_bndry_ele = num_bndry_ele / num_idxs_per_bndry_ele;
363 }
364 }
365 else
366 {
367 // Skipping Boundary Element Data
368 }
369
370 const int *mesh_atts = NULL;
371 const int *bndry_atts = NULL;
372
373 // These variables are used in debug code below.
374 // int num_mesh_atts_entires = 0;
375 // int num_bndry_atts_entires = 0;
376
377 // the attribute fields could have several names
378 // for the element attributes check for first occurrence of field with
379 // name containing "_attribute", that doesn't contain "boundary"
380 std::string main_att_name = "";
381
382 const Node &n_fields = n_mesh["fields"];
383 NodeConstIterator itr = n_fields.children();
384
385 while ( itr.has_next() && main_att_name == "" )
386 {
387 itr.next();
388 std::string fld_name = itr.name();
389 if ( fld_name.find("boundary") == std::string::npos &&
390 fld_name.find("_attribute") != std::string::npos )
391 {
392 main_att_name = fld_name;
393 }
394 }
395
396 if ( main_att_name != "" )
397 {
398 const Node &n_mesh_atts_vals = n_fields[main_att_name]["values"];
399
400 // mfem requires ints, we could have int64s, etc convert if necessary
401 if (n_mesh_atts_vals.dtype().is_int() &&
402 n_mesh_atts_vals.is_compact() )
403 {
404 mesh_atts = n_mesh_atts_vals.value();
405 }
406 else
407 {
408 Node &n_mesh_atts_vals_conv = n_conv["fields"][main_att_name]["values"];
409 n_mesh_atts_vals.to_int_array(n_mesh_atts_vals_conv);
410 mesh_atts = n_mesh_atts_vals_conv.value();
411 }
412
413 // num_mesh_atts_entires = n_mesh_atts_vals.dtype().number_of_elements();
414 }
415 else
416 {
417 // Skipping Mesh Attribute Data
418 }
419
420 // for the boundary attributes check for first occurrence of field with
421 // name containing "_attribute", that also contains "boundary"
422 std::string bnd_att_name = "";
423 itr = n_fields.children();
424
425 while ( itr.has_next() && bnd_att_name == "" )
426 {
427 itr.next();
428 std::string fld_name = itr.name();
429 if ( fld_name.find("boundary") != std::string::npos &&
430 fld_name.find("_attribute") != std::string::npos )
431 {
432 bnd_att_name = fld_name;
433 }
434 }
435
436 if ( bnd_att_name != "" )
437 {
438 // Info: "Getting Boundary Attribute Data"
439 const Node &n_bndry_atts_vals =n_fields[bnd_att_name]["values"];
440
441 // mfem requires ints, we could have int64s, etc convert if necessary
442 if ( n_bndry_atts_vals.dtype().is_int() &&
443 n_bndry_atts_vals.is_compact())
444 {
445 bndry_atts = n_bndry_atts_vals.value();
446 }
447 else
448 {
449 Node &n_bndry_atts_vals_conv = n_conv["fields"][bnd_att_name]["values"];
450 n_bndry_atts_vals.to_int_array(n_bndry_atts_vals_conv);
451 bndry_atts = n_bndry_atts_vals_conv.value();
452 }
453
454 // num_bndry_atts_entires = n_bndry_atts_vals.dtype().number_of_elements();
455
456 }
457 else
458 {
459 // Skipping Boundary Attribute Data
460 }
461
462 // Info: "Number of Vertices: " << num_verts << endl
463 // << "Number of Mesh Elements: " << num_mesh_ele << endl
464 // << "Number of Boundary Elements: " << num_bndry_ele << endl
465 // << "Number of Mesh Attribute Entries: "
466 // << num_mesh_atts_entires << endl
467 // << "Number of Boundary Attribute Entries: "
468 // << num_bndry_atts_entires << endl);
469
470 // Construct MFEM Mesh Object with externally owned data
471 // Note: if we don't have a gf, we need to provide the proper space dim
472 // if nodes gf is attached later, it resets the space dim based
473 // on the gf's fes.
474 Mesh *mesh = new Mesh(// from coordset
475 const_cast<real_t*>(verts_ptr),
476 num_verts,
477 // from topology
478 const_cast<int*>(elem_indices),
479 mesh_geo,
480 // from mesh_attribute field
481 const_cast<int*>(mesh_atts),
482 num_mesh_ele,
483 // from boundary topology
484 const_cast<int*>(bndry_indices),
485 bndry_geo,
486 // from boundary_attribute field
487 const_cast<int*>(bndry_atts),
488 num_bndry_ele,
489 ndims, // dim
490 ndims); // space dim
491
492 // Attach Nodes Grid Function, if it exists
493 if (n_mesh_topo.has_child("grid_function"))
494 {
495 std::string nodes_gf_name = n_mesh_topo["grid_function"].as_string();
496
497 // fetch blueprint field for the nodes gf
498 const Node &n_mesh_gf = n_mesh["fields"][nodes_gf_name];
499 // create gf
501 n_mesh_gf);
502 // attach to mesh
503 mesh->NewNodes(*nodes,true);
504 }
505
506
507 if (zero_copy && !n_conv.dtype().is_empty())
508 {
509 //Info: "Cannot zero-copy since data conversions were necessary"
510 zero_copy = false;
511 }
512
513 Mesh *res = NULL;
514
515 if (zero_copy)
516 {
517 res = mesh;
518 }
519 else
520 {
521 // the mesh above contains references to external data, to get a
522 // copy independent of the conduit data, we use:
523 res = new Mesh(*mesh,true);
524 delete mesh;
525 }
526
527 return res;
528}
529
530//---------------------------------------------------------------------------//
533 const Node &n_field,
534 bool zero_copy)
535{
536 // n_conv holds converted data (when necessary for mfem api)
537 // if n_conv is used ( !n_conv.dtype().empty() ) we
538 // know that some data allocation was necessary, so we
539 // can't return a gf that zero copies the conduit data
540 Node n_conv;
541
542 const real_t *vals_ptr = NULL;
543
544 int vdim = 1;
545
547
548 if (n_field["values"].dtype().is_object())
549 {
550 vdim = n_field["values"].number_of_children();
551
552 // need to check that we have real_t (double/float) and
553 // cover supported layouts
554
555 if ( n_field["values"][0].dtype().id() == mfem_precision_conduit_id )
556 {
557 // check for contig
558 if (n_field["values"].is_contiguous())
559 {
560 // conduit mcarray contig == mfem byNODES
561 vals_ptr = n_field["values"].child(0).value();
562 }
563 // check for interleaved
564 else if (blueprint::mcarray::is_interleaved(n_field["values"]))
565 {
566 // conduit mcarray interleaved == mfem byVDIM
567 ordering = Ordering::byVDIM;
568 vals_ptr = n_field["values"].child(0).value();
569 }
570 else
571 {
572 // for mcarray generic case -- default to byNODES
573 // and provide values w/ contiguous (soa) ordering
574 blueprint::mcarray::to_contiguous(n_field["values"],
575 n_conv["values"]);
576 vals_ptr = n_conv["values"].child(0).value();
577 }
578 }
579 else // convert to real_t (double/float) and use contig
580 {
581 Node n_tmp;
582 // check all vals, if we don't have real_t (double/float) convert
583 // to real_t
584 NodeConstIterator itr = n_field["values"].children();
585 while (itr.has_next())
586 {
587 const Node &c_vals = itr.next();
588 std::string c_name = itr.name();
589
590 if ( c_vals.dtype().id() == mfem_precision_conduit_id )
591 {
592 // zero copy current coords
593 n_tmp[c_name].set_external(c_vals);
594 }
595 else
596 {
597 // convert
598 c_vals.to_data_type(mfem_precision_conduit_id, n_tmp[c_name]);
599 }
600 }
601
602 // for mcarray generic case -- default to byNODES
603 // and provide values w/ contiguous (soa) ordering
604 blueprint::mcarray::to_contiguous(n_tmp,
605 n_conv["values"]);
606 vals_ptr = n_conv["values"].child(0).value();
607 }
608 }
609 else
610 {
611 if (n_field["values"].dtype().id() == mfem_precision_conduit_id &&
612 n_field["values"].is_compact())
613 {
614 vals_ptr = n_field["values"].value();
615 }
616 else
617 {
618 n_field["values"].to_data_type(mfem_precision_conduit_id,
619 n_conv["values"]);
620 vals_ptr = n_conv["values"].value();
621 }
622 }
623
624 if (zero_copy && !n_conv.dtype().is_empty())
625 {
626 //Info: "Cannot zero-copy since data conversions were necessary"
627 zero_copy = false;
628 }
629
630 // we need basis name to create the proper mfem fec
631 std::string fec_name = n_field["basis"].as_string();
632
633 GridFunction *res = NULL;
635 fec_name.c_str());
637 fec,
638 vdim,
639 ordering);
640
641 if (zero_copy)
642 {
643 res = new GridFunction(fes,const_cast<real_t*>(vals_ptr));
644 }
645 else
646 {
647 // copy case, this constructor will alloc the space for the GF data
648 res = new GridFunction(fes);
649 // create an mfem vector that wraps the conduit data
650 Vector vals_vec(const_cast<real_t*>(vals_ptr),fes->GetVSize());
651 // copy values into the result
652 (*res) = vals_vec;
653 }
654
655 // TODO: I believe the GF already has ownership of fes, so this should be all
656 // we need to do to avoid leaking objs created here?
657 res->MakeOwner(fec);
658
659 return res;
660}
661
662//---------------------------------------------------------------------------//
665 const Node &n_field,
666 bool zero_copy)
667{
668 // n_conv holds converted data (when necessary for mfem api)
669 // if n_conv is used ( !n_conv.dtype().empty() ) we
670 // know that some data allocation was necessary, so we
671 // can't return a qf that zero copies the conduit data
672 Node n_conv;
673
674 const real_t *vals_ptr = NULL;
675 int vdim = 1;
676
677 if (n_field["values"].dtype().is_object())
678 {
679 vdim = n_field["values"].number_of_children();
680
681 // need to check that we have real_t (double/float) and
682 // cover supported layouts
683 if ( n_field["values"][0].dtype().id() == mfem_precision_conduit_id )
684 {
685 // quad funcs use what mfem calls byVDIM
686 // and what conduit calls interleaved
687 // check for interleaved
688 if (blueprint::mcarray::is_interleaved(n_field["values"]))
689 {
690 // conduit mcarray interleaved == mfem byVDIM
691 vals_ptr = n_field["values"].child(0).value();
692 }
693 else
694 {
695 // for mcarray generic case -- default to byVDIM
696 // aka interleaved
697 blueprint::mcarray::to_interleaved(n_field["values"],
698 n_conv["values"]);
699 vals_ptr = n_conv["values"].child(0).value();
700 }
701 }
702 else // convert to real_t (double/float) and use interleaved
703 {
704 Node n_tmp;
705 // check all vals, if we don't have real_t (double/float) convert
706 // to real_t
707 NodeConstIterator itr = n_field["values"].children();
708 while (itr.has_next())
709 {
710 const Node &c_vals = itr.next();
711 std::string c_name = itr.name();
712
713 if ( c_vals.dtype().id() == mfem_precision_conduit_id )
714 {
715 // zero copy current coords
716 n_tmp[c_name].set_external(c_vals);
717 }
718 else
719 {
720 // convert
721 c_vals.to_data_type(mfem_precision_conduit_id, n_tmp[c_name]);
722 }
723 }
724
725 // for mcarray generic case -- default to byVDIM
726 // aka interleaved
727 blueprint::mcarray::to_interleaved(n_tmp,
728 n_conv["values"]);
729 vals_ptr = n_conv["values"].child(0).value();
730 }
731 }
732 else // scalar case
733 {
734 if (n_field["values"].dtype().id() == mfem_precision_conduit_id &&
735 n_field["values"].is_compact())
736 {
737 vals_ptr = n_field["values"].value();
738 }
739 else
740 {
741 n_field["values"].to_data_type(mfem_precision_conduit_id,
742 n_conv["values"]);
743 vals_ptr = n_conv["values"].value();
744 }
745 }
746
747 if (zero_copy && !n_conv.dtype().is_empty())
748 {
749 //Info: "Cannot zero-copy since data conversions were necessary"
750 zero_copy = false;
751 }
752
753 // we need basis name to create the proper mfem quad space and quad func
754 // the pattern used to encode the quad space params is:
755 // QF_{ORDER}_{VDIM}
756 // ORDER is the degree of the polynomials for the quad rule
757 // VDIM is the number of components at each quad point (scalar, vector, etc)
758
759 int qf_order = 0;
760 int qf_vdim = 0;
761 std::string qf_name = n_field["basis"].as_string();
762 const char *qf_name_cstr = qf_name.c_str();
763 if (!strncmp(qf_name_cstr, "QF_", 3))
764 {
765 // parse {ORDER}
766 qf_order = atoi(qf_name_cstr + 3);
767 // find second `_`
768 const char *qf_vdim_cstr = strstr(qf_name_cstr+3,"_");
769 if (qf_vdim_cstr == NULL)
770 {
771 MFEM_ABORT("Error parsing quadrature function description string: "
772 << qf_name << std::endl
773 << "Expected: QF_{ORDER}_{VDIM}");
774 }
775 // parse {VDIM}
776 qf_vdim = atoi(qf_vdim_cstr+1);
777 }
778 else
779 {
780 MFEM_ABORT("Error parsing quadrature function description string: "
781 << qf_name << std::endl
782 << "Expected: QF_{ORDER}_{VDIM}");
783 }
784 MFEM_VERIFY(qf_vdim == vdim, "vector dimension mismatch: vdim = " << vdim
785 << ", qf_vdim = " << qf_vdim);
786
787 mfem::QuadratureSpace *quad_space = new mfem::QuadratureSpace(mesh, qf_order);
789
790 if (zero_copy)
791 {
792 res->SetSpace(quad_space, const_cast<real_t*>(vals_ptr), vdim);
793 res->SetOwnsSpace(true);
794 }
795 else
796 {
797 res->SetSpace(quad_space, vdim);
798 res->SetOwnsSpace(true);
799 // copy case, this constructor will alloc the space for the quad data
800 // create an mfem vector that wraps the conduit data
801 Vector vals_vec(const_cast<real_t*>(vals_ptr),res->Size());
802 // copy values into the result
803 (*res) = vals_vec;
804 }
805
806 return res;
807}
808
809
810
811//---------------------------------------------------------------------------//
812void
814 Node &n_mesh,
815 const std::string &coordset_name,
816 const std::string &main_topology_name,
817 const std::string &boundary_topology_name,
818 const std::string &main_adjset_name)
819{
820 int dim = mesh->SpaceDimension();
821
822 MFEM_ASSERT(dim >= 1 && dim <= 3, "invalid mesh dimension");
823
824 ////////////////////////////////////////////
825 // Setup main coordset
826 ////////////////////////////////////////////
827
828 // Assumes mfem::Vertex has the layout of a real_t (double/float) array.
829
830 // this logic assumes an mfem vertex is always 3 real_t (double/float) wide
831 int stride = sizeof(mfem::Vertex);
832 int num_vertices = mesh->GetNV();
833
834 MFEM_ASSERT( ( stride == 3 * sizeof(real_t) ),
835 "Unexpected stride for Vertex");
836
837 Node &n_mesh_coords = n_mesh["coordsets"][coordset_name];
838 n_mesh_coords["type"] = "explicit";
839
840
841 real_t *coords_ptr = mesh->GetVertex(0);
842
843 n_mesh_coords["values/x"].set_external(coords_ptr,
844 num_vertices,
845 0,
846 stride);
847
848 if (dim >= 2)
849 {
850 n_mesh_coords["values/y"].set_external(coords_ptr,
851 num_vertices,
852 sizeof(real_t),
853 stride);
854 }
855 if (dim >= 3)
856 {
857 n_mesh_coords["values/z"].set_external(coords_ptr,
858 num_vertices,
859 sizeof(real_t) * 2,
860 stride);
861 }
862
863 ////////////////////////////////////////////
864 // Setup main topo
865 ////////////////////////////////////////////
866
867 Node &n_topo = n_mesh["topologies"][main_topology_name];
868
869 n_topo["type"] = "unstructured";
870 n_topo["coordset"] = coordset_name;
871
872 Element::Type ele_type = mesh->GetElementType(0);
873
874 std::string ele_shape = ElementTypeToShapeName(ele_type);
875
876 n_topo["elements/shape"] = ele_shape;
877
878 GridFunction *gf_mesh_nodes = mesh->GetNodes();
879
880 if (gf_mesh_nodes != NULL)
881 {
882 n_topo["grid_function"] = "mesh_nodes";
883 }
884
885 // connectivity
886 // TODO: generic case, i don't think we can zero-copy (mfem allocs
887 // an array per element) so we alloc our own contig array and
888 // copy out. Some other cases (sidre) may actually have contig
889 // allocation but I am not sure how to detect this case from mfem
890 int num_ele = mesh->GetNE();
891 int geom = mesh->GetTypicalElementGeometry();
892 int idxs_per_ele = Geometry::NumVerts[geom];
893 int num_conn_idxs = num_ele * idxs_per_ele;
894
895 n_topo["elements/connectivity"].set(DataType::c_int(num_conn_idxs));
896
897 int *conn_ptr = n_topo["elements/connectivity"].value();
898
899 for (int i=0; i < num_ele; i++)
900 {
901 const Element *ele = mesh->GetElement(i);
902 const int *ele_verts = ele->GetVertices();
903
904 memcpy(conn_ptr, ele_verts, idxs_per_ele * sizeof(int));
905
906 conn_ptr += idxs_per_ele;
907 }
908
909 if (gf_mesh_nodes != NULL)
910 {
911 GridFunctionToBlueprintField(gf_mesh_nodes,
912 n_mesh["fields/mesh_nodes"],
913 main_topology_name);
914 }
915
916 ////////////////////////////////////////////
917 // Setup mesh attribute
918 ////////////////////////////////////////////
919
920 Node &n_mesh_att = n_mesh["fields/element_attribute"];
921
922 n_mesh_att["association"] = "element";
923 n_mesh_att["topology"] = main_topology_name;
924 n_mesh_att["values"].set(DataType::c_int(num_ele));
925
926 int_array att_vals = n_mesh_att["values"].value();
927 for (int i = 0; i < num_ele; i++)
928 {
929 att_vals[i] = mesh->GetAttribute(i);
930 }
931
932 ////////////////////////////////////////////
933 // Setup bndry topo "boundary"
934 ////////////////////////////////////////////
935
936 // guard vs if we have boundary elements
938 {
939 n_topo["boundary_topology"] = boundary_topology_name;
940
941 Node &n_bndry_topo = n_mesh["topologies"][boundary_topology_name];
942
943 n_bndry_topo["type"] = "unstructured";
944 n_bndry_topo["coordset"] = coordset_name;
945
946 int num_bndry_ele = mesh->GetNBE();
947
948 Element *BE0 = NULL; // representative boundary element
949 if (num_bndry_ele > 0) { BE0 = mesh->GetBdrElement(0); }
950
951 // must initialize this to something, pick POINT if no boundary elements
952 Element::Type bndry_ele_type = (BE0) ? BE0->GetType() : Element::POINT;
953 std::string bndry_ele_shape = ElementTypeToShapeName(bndry_ele_type);
954 n_bndry_topo["elements/shape"] = bndry_ele_shape;
955
956 // must initialize this to something, pick POINT if no boundary elements
957 int bndry_geom = (BE0) ? BE0->GetGeometryType() : Geometry::POINT;
958 int bndry_idxs_per_ele = Geometry::NumVerts[bndry_geom];
959 int num_bndry_conn_idxs = num_bndry_ele * bndry_idxs_per_ele;
960
961 n_bndry_topo["elements/connectivity"].set(DataType::c_int(num_bndry_conn_idxs));
962
963 int *bndry_conn_ptr = n_bndry_topo["elements/connectivity"].value();
964
965 for (int i=0; i < num_bndry_ele; i++)
966 {
967 const Element *bndry_ele = mesh->GetBdrElement(i);
968 const int *bndry_ele_verts = bndry_ele->GetVertices();
969
970 memcpy(bndry_conn_ptr, bndry_ele_verts, bndry_idxs_per_ele * sizeof(int));
971
972 bndry_conn_ptr += bndry_idxs_per_ele;
973 }
974
975 ////////////////////////////////////////////
976 // Setup bndry mesh attribute
977 ////////////////////////////////////////////
978
979 Node &n_bndry_mesh_att = n_mesh["fields/boundary_attribute"];
980
981 n_bndry_mesh_att["association"] = "element";
982 n_bndry_mesh_att["topology"] = boundary_topology_name;
983 n_bndry_mesh_att["values"].set(DataType::c_int(num_bndry_ele));
984
985 int_array bndry_att_vals = n_bndry_mesh_att["values"].value();
986 for (int i = 0; i < num_bndry_ele; i++)
987 {
988 bndry_att_vals[i] = mesh->GetBdrAttribute(i);
989 }
990 }
991
992 ////////////////////////////////////////////
993 // Setup adjsets
994 ////////////////////////////////////////////
995
996#ifdef MFEM_USE_MPI
997 ParMesh *pmesh = dynamic_cast<ParMesh*>(mesh);
998 if (pmesh)
999 {
1000 ////////////////////////////////////////////
1001 // Setup main adjset
1002 ////////////////////////////////////////////
1003
1004 Node &n_adjset = n_mesh["adjsets"][main_adjset_name];
1005
1006 n_adjset["association"] = "vertex";
1007 n_adjset["topology"] = main_topology_name;
1008 n_adjset["groups"].set(DataType::object());
1009
1010 const GroupTopology &pmesh_gtopo = pmesh->gtopo;
1011 const int local_rank = pmesh->GetMyRank();
1012 const int num_groups = pmesh_gtopo.NGroups();
1013 // NOTE: skip the first group since its the local-only group
1014 for (int i = 1; i < num_groups; i++)
1015 {
1016 const int num_group_nbrs = pmesh_gtopo.GetGroupSize(i);
1017 const int *group_nbrs = pmesh_gtopo.GetGroup(i);
1018 const int num_group_verts = pmesh->GroupNVertices(i);
1019
1020 // NOTE: 'neighbor' values are local to this processor, but Blueprint
1021 // expects global domain identifiers, so we collapse this layer of
1022 // indirection
1023 Array<int> group_ranks(num_group_nbrs);
1024 std::string group_name = "group";
1025 {
1026 for (int j = 0; j < num_group_nbrs; j++)
1027 {
1028 group_ranks[j] = pmesh_gtopo.GetNeighborRank(group_nbrs[j]);
1029 }
1030 group_ranks.Sort();
1031 for (int j = 0; j < num_group_nbrs; j++)
1032 {
1033 group_name += "_" + std::to_string(group_ranks[j]);
1034 }
1035
1036 // NOTE: Blueprint only wants remote ranks in its neighbor list,
1037 // so we remove the local rank after the canonicalized Blueprint
1038 // group name is formed
1039 group_ranks.DeleteFirst(local_rank);
1040 }
1041 Node &n_group = n_adjset["groups"][group_name];
1042
1043 n_group["neighbors"].set(group_ranks.GetData(), group_ranks.Size());
1044 n_group["values"].set(DataType::c_int(num_group_verts));
1045
1046 int_array group_vals = n_group["values"].value();
1047 for (int j = 0; j < num_group_verts; j++)
1048 {
1049 group_vals[j] = pmesh->GroupVertex(i, j);
1050 }
1051 }
1052
1053 // NOTE: We don't create an adjset for face neighbor data because
1054 // these faces aren't listed in the 'boundary_topology_name' topology
1055 // (this topology only covers the faces between 'main_topology_name'
1056 // elements and void). To include a face neighbor data adjset, this
1057 // function would need to export a topology with either (1) all faces
1058 // in the mesh topology or (2) all boundary faces, including neighbors.
1059
1060 ////////////////////////////////////////////
1061 // Setup distributed state
1062 ////////////////////////////////////////////
1063
1064 Node &n_domid = n_mesh["state/domain_id"];
1065 n_domid.set(local_rank);
1066 }
1067#endif
1068}
1069
1070//---------------------------------------------------------------------------//
1071void
1073 Node &n_field,
1074 const std::string &main_topology_name)
1075{
1076 n_field["basis"] = gf->FESpace()->FEColl()->Name();
1077 n_field["topology"] = main_topology_name;
1078
1079 int vdim = gf->FESpace()->GetVDim();
1080 int ndofs = gf->FESpace()->GetNDofs();
1081
1082 if (vdim == 1) // scalar case
1083 {
1084 n_field["values"].set_external(const_cast<real_t *>(gf->HostRead()),
1085 ndofs);
1086 }
1087 else // vector case
1088 {
1089 // deal with striding of all components
1090
1091 Ordering::Type ordering = gf->FESpace()->GetOrdering();
1092
1093 int entry_stride = (ordering == Ordering::byNODES ? 1 : vdim);
1094 int vdim_stride = (ordering == Ordering::byNODES ? ndofs : 1);
1095
1096 index_t offset = 0;
1097 index_t stride = sizeof(real_t) * entry_stride;
1098
1099 for (int d = 0; d < vdim; d++)
1100 {
1101 std::ostringstream oss;
1102 oss << "v" << d;
1103 std::string comp_name = oss.str();
1104 n_field["values"][comp_name].set_external(const_cast<real_t *>(gf->HostRead()),
1105 ndofs,
1106 offset,
1107 stride);
1108 offset += sizeof(real_t) * vdim_stride;
1109 }
1110 }
1111
1112}
1113
1114//---------------------------------------------------------------------------//
1115void
1118 Node &n_field,
1119 const std::string &main_topology_name)
1120{
1121 // For quadrature functions, use basis pattern:
1122 // QF_{ORDER}_{VDIM}
1123
1124 int qf_vdim = qf->GetVDim();
1125 int qf_order = qf->GetSpace()->GetOrder();
1126 int qf_size = qf->GetSpace()->GetSize();
1127
1128 {
1129 std::ostringstream oss;
1130 oss << "QF_" << qf_order << "_" << qf_vdim;
1131
1132 n_field["basis"] = oss.str();
1133 n_field["topology"] = main_topology_name;
1134 }
1135
1136 if (qf_vdim == 1) // scalar case
1137 {
1138 n_field["values"].set_external(const_cast<real_t *>(qf->HostRead()),
1139 qf_size);
1140 }
1141 else // vector case
1142 {
1143 // deal with striding of all components
1144 // quadrature functions are always byVDIM
1145 // or what conduit calls interleaved
1146
1147 index_t offset = 0;
1148 index_t stride = sizeof(real_t) * qf_vdim;
1149
1150 for (int d = 0; d < qf_vdim; d++)
1151 {
1152 std::ostringstream oss;
1153 oss << "v" << d;
1154 std::string comp_name = oss.str();
1155 n_field["values"][comp_name].set_external(const_cast<real_t *>(qf->HostRead()),
1156 qf_size,
1157 offset,
1158 stride);
1159 offset += sizeof(real_t);
1160 }
1161 }
1162
1163}
1164
1165
1166
1167//------------------------------
1168// end static public methods
1169//------------------------------
1170
1171//------------------------------
1172// end public methods
1173//------------------------------
1174
1175//------------------------------
1176// begin protected methods
1177//------------------------------
1178
1179//---------------------------------------------------------------------------//
1180std::string
1182{
1183 std::string res = prefix_path + name + "_" +
1185 ".root";
1186 return res;
1187}
1188
1189//---------------------------------------------------------------------------//
1190std::string
1192 const std::string &relay_protocol_)
1193{
1194 std::string res = prefix_path +
1195 name +
1196 "_" +
1198 "/domain_" +
1199 to_padded_string(domain_id, pad_digits_rank) +
1200 "." +
1201 relay_protocol_;
1202
1203 return res;
1204}
1205
1206//---------------------------------------------------------------------------//
1207std::string
1209{
1210 std::string res = prefix_path +
1211 name +
1212 "_" +
1214 return res;
1215}
1216
1217//---------------------------------------------------------------------------//
1218std::string
1219ConduitDataCollection::MeshFilePattern(const std::string &relay_protocol_)
1220{
1221 std::ostringstream oss;
1222 oss << name
1223 << "_"
1225 << "/domain_%0"
1227 << "d."
1228 << relay_protocol_;
1229
1230 return oss.str();
1231}
1232
1233
1234//---------------------------------------------------------------------------//
1235void
1237 const Node &n_mesh,
1238 const std::string &relay_protocol_)
1239{
1240 // default to json root file, except for hdf5 case
1241 std::string root_proto = "json";
1242
1243 if (relay_protocol_ == "hdf5")
1244 {
1245 root_proto = relay_protocol_;
1246 }
1247
1248 Node n_root;
1249 // create blueprint index
1250 Node &n_bp_idx = n_root["blueprint_index"];
1251
1252 blueprint::mesh::generate_index(n_mesh,
1253 "",
1254 num_domains,
1255 n_bp_idx["mesh"]);
1256
1257 // there are cases where the data backing the gf fields doesn't
1258 // accurately represent the number of components in physical space,
1259 // so we loop over all gfs and fix those that are incorrect
1260
1262 for ( itr = field_map.begin(); itr != field_map.end(); itr++)
1263 {
1264 std::string gf_name = itr->first;
1265 GridFunction *gf = itr->second;
1266
1267 Node &idx_gf_ncomps = n_bp_idx["mesh/fields"][gf_name]["number_of_components"];
1268 // check that the number_of_components in the index matches what we expect
1269 // correct if necessary
1270 if ( idx_gf_ncomps.to_int() != gf->VectorDim() )
1271 {
1272 idx_gf_ncomps = gf->VectorDim();
1273 }
1274 }
1275 // add extra header info
1276 n_root["protocol/name"] = relay_protocol_;
1277 n_root["protocol/version"] = "0.3.1";
1278
1279
1280 // we will save one file per domain, so trees == files
1281 n_root["number_of_files"] = num_domains;
1282 n_root["number_of_trees"] = num_domains;
1283 n_root["file_pattern"] = MeshFilePattern(relay_protocol_);
1284 n_root["tree_pattern"] = "";
1285
1286 // Add the time, time step, and cycle
1287 n_root["blueprint_index/mesh/state/time"] = time;
1288 n_root["blueprint_index/mesh/state/time_step"] = time_step;
1289 n_root["blueprint_index/mesh/state/cycle"] = cycle;
1290
1291 relay::io::save(n_root, RootFileName(), root_proto);
1292}
1293
1294//---------------------------------------------------------------------------//
1295void
1297 const Node &n_mesh,
1298 const std::string &relay_protocol_)
1299{
1300 relay::io::save(n_mesh, MeshFileName(domain_id, relay_protocol_));
1301}
1302
1303//---------------------------------------------------------------------------//
1304void
1306{
1307 if (myid == 0)
1308 {
1309 // assume root file is json, unless hdf5 is specified
1310 std::string root_protocol = "json";
1311
1312 if ( relay_protocol.find("hdf5") != std::string::npos )
1313 {
1314 root_protocol = "hdf5";
1315 }
1316
1317
1318 relay::io::load(RootFileName(), root_protocol, root_out);
1319#ifdef MFEM_USE_MPI
1320 // broadcast contents of root file other ranks
1321 // (conduit relay mpi would simplify, but we would need to link another
1322 // lib for mpi case)
1323
1324 // create json string
1325 std::string root_json = root_out.to_json();
1326 // string size +1 for null term
1327 int json_str_size = root_json.size() + 1;
1328
1329 // broadcast json string buffer size
1330 int mpi_status = MPI_Bcast((void*)&json_str_size, // ptr
1331 1, // size
1332 MPI_INT, // type
1333 0, // root
1334 m_comm); // comm
1335
1336 if (mpi_status != MPI_SUCCESS)
1337 {
1338 MFEM_ABORT("Broadcast of root file json string size failed");
1339 }
1340
1341 // broadcast json string
1342 mpi_status = MPI_Bcast((void*)root_json.c_str(), // ptr
1343 json_str_size, // size
1344 MPI_CHAR, // type
1345 0, // root
1346 m_comm); // comm
1347
1348 if (mpi_status != MPI_SUCCESS)
1349 {
1350 MFEM_ABORT("Broadcast of root file json string failed");
1351 }
1352
1353#endif
1354 }
1355
1356#ifdef MFEM_USE_MPI
1357 else
1358 {
1359 // recv json string buffer size via broadcast
1360 int json_str_size = -1;
1361 int mpi_status = MPI_Bcast(&json_str_size, // ptr
1362 1, // size
1363 MPI_INT, // type
1364 0, // root
1365 m_comm); // comm
1366
1367 if (mpi_status != MPI_SUCCESS)
1368 {
1369 MFEM_ABORT("Broadcast of root file json string size failed");
1370 }
1371
1372 // recv json string buffer via broadcast
1373 char *json_buff = new char[json_str_size];
1374 mpi_status = MPI_Bcast(json_buff, // ptr
1375 json_str_size, // size
1376 MPI_CHAR, // type
1377 0, // root
1378 m_comm); // comm
1379
1380 if (mpi_status != MPI_SUCCESS)
1381 {
1382 MFEM_ABORT("Broadcast of root file json string failed");
1383 }
1384
1385 // reconstruct root file contents
1386 Generator g(std::string(json_buff),"json");
1387 g.walk(root_out);
1388 // cleanup temp buffer
1389 delete [] json_buff;
1390 }
1391#endif
1392}
1393
1394//---------------------------------------------------------------------------//
1395void
1397 const std::string &relay_protocol_)
1398{
1399 // Note: This path doesn't use any info from the root file
1400 // it uses the implicit mfem ConduitDataCollection layout
1401
1402 Node n_mesh;
1403 relay::io::load( MeshFileName(domain_id, relay_protocol_), n_mesh);
1404
1405
1406 Node verify_info;
1407 if (!blueprint::mesh::verify(n_mesh,verify_info))
1408 {
1409 MFEM_ABORT("Conduit Mesh Blueprint Verify Failed:\n"
1410 << verify_info.to_json());
1411 }
1412
1413 mesh = BlueprintMeshToMesh(n_mesh);
1414
1415 field_map.clear();
1416
1417 NodeConstIterator itr = n_mesh["fields"].children();
1418
1419 std::string nodes_gf_name = "";
1420
1421 const Node &n_topo = n_mesh["topologies/main"];
1422 if (n_topo.has_child("grid_function"))
1423 {
1424 nodes_gf_name = n_topo["grid_function"].as_string();
1425 }
1426
1427 while (itr.has_next())
1428 {
1429 const Node &n_field = itr.next();
1430 std::string field_name = itr.name();
1431
1432 // skip mesh nodes gf since they are already processed
1433 // skip attribute fields, they aren't grid functions
1434 if ( field_name != nodes_gf_name &&
1435 field_name.find("_attribute") == std::string::npos
1436 )
1437 {
1439 field_map.Register(field_name, gf, true);
1440 }
1441 }
1442}
1443
1444//------------------------------
1445// end protected methods
1446//------------------------------
1447
1448//------------------------------
1449// begin static private methods
1450//------------------------------
1451
1452//---------------------------------------------------------------------------//
1453std::string
1454ConduitDataCollection::ElementTypeToShapeName(Element::Type element_type)
1455{
1456 // Adapted from SidreDataCollection
1457
1458 // Note -- the mapping from Element::Type to string is based on
1459 // enum Element::Type { POINT, SEGMENT, TRIANGLE, QUADRILATERAL,
1460 // TETRAHEDRON, HEXAHEDRON };
1461 // Note: -- the string names are from conduit's blueprint
1462
1463 switch (element_type)
1464 {
1465 case Element::POINT: return "point";
1466 case Element::SEGMENT: return "line";
1467 case Element::TRIANGLE: return "tri";
1468 case Element::QUADRILATERAL: return "quad";
1469 case Element::TETRAHEDRON: return "tet";
1470 case Element::HEXAHEDRON: return "hex";
1471 case Element::WEDGE:
1472 default: ;
1473 }
1474
1475 return "unknown";
1476}
1477
1478//---------------------------------------------------------------------------//
1480ConduitDataCollection::ShapeNameToGeomType(const std::string &shape_name)
1481{
1482 // Note: must init to something to avoid invalid memory access
1483 // in the mfem mesh constructor
1485
1486 if (shape_name == "point")
1487 {
1489 }
1490 else if (shape_name == "line")
1491 {
1493 }
1494 else if (shape_name == "tri")
1495 {
1497 }
1498 else if (shape_name == "quad")
1499 {
1501 }
1502 else if (shape_name == "tet")
1503 {
1505 }
1506 else if (shape_name == "hex")
1507 {
1509 }
1510 else
1511 {
1512 MFEM_ABORT("Unsupported Element Shape: " << shape_name);
1513 }
1514
1515 return res;
1516}
1517
1518//------------------------------
1519// end static private methods
1520//------------------------------
1521
1522} // end namespace mfem
1523
1524#endif
void DeleteFirst(const T &el)
Delete the first entry with value == 'el'.
Definition array.hpp:990
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
Definition array.hpp:312
int Size() const
Return the logical size of the array.
Definition array.hpp:166
T * GetData()
Returns the data.
Definition array.hpp:140
static Mesh * BlueprintMeshToMesh(const conduit::Node &n_mesh, const std::string &main_toplogy_name="", bool zero_copy=false)
Constructs and MFEM mesh from a Conduit Blueprint Description.
std::string RootFileName()
Returns blueprint root file name for the current cycle.
void SaveRootFile(int num_domains, const conduit::Node &n_mesh, const std::string &file_protocol)
Saves root file for the current cycle.
void LoadMeshAndFields(int domain_id, const std::string &file_protocol)
Loads all meshes and fields of a given domain id for the current cycle.
static GridFunction * BlueprintFieldToGridFunction(Mesh *mesh, const conduit::Node &n_field, bool zero_copy=false)
Constructs and MFEM Grid Function from a Conduit Blueprint Description.
virtual void Load(int cycle=0)
Load the collection based blueprint data.
static void MeshToBlueprintMesh(Mesh *m, conduit::Node &out, const std::string &coordset_name="coords", const std::string &main_topology_name="main", const std::string &boundary_topology_name="boundary", const std::string &main_adjset_name="main_adjset")
Describes a MFEM mesh using the mesh blueprint.
std::string MeshDirectoryName()
Returns the mesh output directory for the current cycle.
void SetProtocol(const std::string &protocol)
Set the Conduit relay i/o protocol to use.
std::string MeshFileName(int domain_id, const std::string &file_protocol="hdf5")
Returns the mesh file name for a given domain at the current cycle.
static void GridFunctionToBlueprintField(GridFunction *gf, conduit::Node &out, const std::string &main_topology_name="main")
Describes a MFEM grid function using the mesh blueprint.
std::string MeshFilePattern(const std::string &file_protocol="hdf5")
Returns the mesh file pattern for the current cycle.
static QuadratureFunction * BlueprintFieldToQuadratureFunction(Mesh *mesh, const conduit::Node &n_field, bool zero_copy=false)
Constructs and MFEM Quadrature Function from a Conduit Blueprint Description.
static void QuadratureFunctionToBlueprintField(QuadratureFunction *qf, conduit::Node &out, const std::string &main_topology_name="main")
Describes a MFEM quadrature function using the mesh blueprint.
virtual void Save()
Save the collection and a Conduit blueprint root file.
void SaveMeshAndFields(int domain_id, const conduit::Node &n_mesh, const std::string &file_protocol)
Saves all meshes and fields for the current cycle.
virtual ~ConduitDataCollection()
We will delete the mesh and fields if we own them.
void LoadRootFile(conduit::Node &n_root_out)
Loads contents of the root field for the current cycle into n_root_out.
ConduitDataCollection(const std::string &collection_name, Mesh *mesh=NULL)
Constructor. The collection name is used when saving the data.
int cycle
Time cycle; for time-dependent simulations cycle >= 0, otherwise = -1.
real_t time
Physical time (for time-dependent simulations)
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)
void DeleteAll()
Delete data owned by the DataCollection including field information.
GFieldMap::const_iterator FieldMapConstIterator
QFieldMap::const_iterator QFieldMapConstIterator
int pad_digits_cycle
Number of digits used for the cycle and MPI rank in filenames.
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)
real_t time_step
Time step i.e. delta_t (for time-dependent simulations)
int num_procs
Number of MPI ranks (in parallel)
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.
Abstract data type element.
Definition element.hpp:29
Geometry::Type GetGeometryType() const
Definition element.hpp:55
virtual void GetVertices(Array< int > &v) const =0
Get the indices defining the vertices.
virtual Type GetType() const =0
Returns element's type.
Type
Constants for the classes derived from Element.
Definition element.hpp:41
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name.
Definition fe_coll.cpp:124
virtual const char * Name() const
Definition fe_coll.hpp:79
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition fespace.hpp:823
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:854
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:856
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:826
int GetVDim() const
Returns the vector dimension of the finite element space.
Definition fespace.hpp:819
static const int NumVerts[NumGeom]
Definition geom.hpp:53
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec_owned and fes.
Definition gridfunc.hpp:123
FiniteElementSpace * FESpace()
int VectorDim() const
Shortcut for calling FiniteElementSpace::GetVectorDim() on the underlying fes.
Definition gridfunc.cpp:347
int GetNeighborRank(int i) const
Return the MPI rank of neighbor 'i'.
const int * GetGroup(int g) const
Return a pointer to a list of neighbors for a given group. Neighbor 0 is the local processor.
int GetGroupSize(int g) const
Get the number of processors in a group.
int NGroups() const
Return the number of groups.
Mesh data type.
Definition mesh.hpp:65
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition mesh.cpp:7962
int GetAttribute(int i) const
Return the attribute of element i.
Definition mesh.hpp:1484
Geometry::Type GetTypicalElementGeometry() const
If the local mesh is not empty, return GetElementGeometry(0); otherwise, return a typical Geometry pr...
Definition mesh.cpp:1628
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition mesh.hpp:1490
const Element * GetElement(int i) const
Return pointer to the i'th element object.
Definition mesh.hpp:1434
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1377
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
Definition mesh.hpp:1449
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1309
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:9627
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition mesh.hpp:1374
int GetNBE() const
Returns number of boundary elements.
Definition mesh.hpp:1380
void NewNodes(GridFunction &nodes, bool make_owner=false)
Replace the internal node GridFunction with the given GridFunction.
Definition mesh.cpp:9654
virtual bool HasBoundaryElements() const
Checks if the mesh has boundary elements.
Definition mesh.hpp:1336
const real_t * GetVertex(int i) const
Return pointer to vertex i's coordinates.
Definition mesh.hpp:1416
void Register(const std::string &fname, T *field, bool own_data)
Register field field with name fname.
iterator end()
Returns an end iterator to the registered fields.
iterator begin()
Returns a begin iterator to the registered fields.
void clear()
Clears the map of registered fields without reclaiming memory.
Type
Ordering methods:
Definition ordering.hpp:17
Class for parallel meshes.
Definition pmesh.hpp:34
int GetMyRank() const
Definition pmesh.hpp:407
int GroupVertex(int group, int i) const
Accessors for entities within a shared group structure.
Definition pmesh.hpp:493
GroupTopology gtopo
Definition pmesh.hpp:459
int GroupNVertices(int group) const
Definition pmesh.hpp:477
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
Definition qfunction.hpp:94
void SetOwnsSpace(bool own)
Set the QuadratureSpaceBase ownership flag.
int GetVDim() const
Get the vector dimension.
Definition qfunction.hpp:87
void SetSpace(QuadratureSpaceBase *qspace_, int vdim_=-1)
Change the QuadratureSpaceBase and optionally the vector dimension.
int GetSize() const
Return the total number of quadrature points.
Definition qspace.hpp:107
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Definition qspace.hpp:110
Class representing the storage layout of a QuadratureFunction.
Definition qspace.hpp:164
Vector data type.
Definition vector.hpp:82
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:524
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
Data type for vertex.
Definition vertex.hpp:23
int dim
Definition ex24.cpp:53
constexpr conduit::index_t mfem_precision_conduit_id
std::string to_padded_string(int i, int digits)
Convert an integer to a 0-padded string with the given number of digits.
Definition text.hpp:54
float real_t
Definition config.hpp:46
std::array< int, NCMesh::MaxFaceNodes > nodes