Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/de/ded/xdmf__mesh_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
xdmf_mesh.h
1// Copyright (C) 2012-2018 Chris N. Richardson and Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "xdmf_utils.h"
10#include <array>
11#include <concepts>
12#include <cstdint>
13#include <dolfinx/common/MPI.h>
14#include <dolfinx/mesh/MeshTags.h>
15#include <hdf5.h>
16#include <mpi.h>
17#include <pugixml.hpp>
18#include <span>
19#include <string>
20#include <utility>
21#include <variant>
22#include <vector>
23
24namespace pugi
25{
26class xml_node;
27} // namespace pugi
28
29namespace dolfinx
30{
31
32namespace mesh
33{
34template <std::floating_point T>
35class Geometry;
36template <std::floating_point T>
37class Mesh;
38class Topology;
39} // namespace mesh
40
42namespace io::xdmf_mesh
43{
44
49template <std::floating_point U>
50void add_mesh(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
51 const mesh::Mesh<U>& mesh, const std::string& path_prefix);
52
63template <std::floating_point U>
64void add_topology_data(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
65 std::string path_prefix, const mesh::Topology& topology,
66 const mesh::Geometry<U>& geometry, int cell_dim,
67 std::span<const std::int32_t> entities);
68
70template <std::floating_point U>
71void add_geometry_data(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
72 std::string path_prefix,
73 const mesh::Geometry<U>& geometry);
74
81std::pair<std::variant<std::vector<float>, std::vector<double>>,
82 std::array<std::size_t, 2>>
83read_geometry_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node& node);
84
91std::pair<std::vector<std::int64_t>, std::array<std::size_t, 2>>
92read_topology_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node& node);
93
95template <typename T, std::floating_point U>
96void add_meshtags(MPI_Comm comm, const mesh::MeshTags<T>& meshtags,
97 const mesh::Geometry<U>& geometry, pugi::xml_node& xml_node,
98 hid_t h5_id, const std::string& name)
99{
100 LOG(INFO) << "XDMF: add meshtags (" << name << ")";
101 // Get mesh
102 const int dim = meshtags.dim();
103 std::shared_ptr<const common::IndexMap> entity_map
104 = meshtags.topology()->index_map(dim);
105 if (!entity_map)
106 {
107 throw std::runtime_error("Missing entities. Did you forget to call "
108 "dolfinx::mesh::Topology::create_entities?");
109 }
110 const std::int32_t num_local_entities = entity_map->size_local();
111
112 // Find number of tagged entities in local range
113 auto it = std::lower_bound(meshtags.indices().begin(),
114 meshtags.indices().end(), num_local_entities);
115 const int num_active_entities = std::distance(meshtags.indices().begin(), it);
116
117 const std::string path_prefix = "/MeshTags/" + name;
119 comm, xml_node, h5_id, path_prefix, *meshtags.topology(), geometry, dim,
120 std::span<const std::int32_t>(meshtags.indices().data(),
121 num_active_entities));
122
123 // Add attribute node with values
124 pugi::xml_node attribute_node = xml_node.append_child("Attribute");
125 assert(attribute_node);
126 attribute_node.append_attribute("Name") = name.c_str();
127 attribute_node.append_attribute("AttributeType") = "Scalar";
128 attribute_node.append_attribute("Center") = "Cell";
129
130 std::int64_t global_num_values = 0;
131 const std::int64_t local_num_values = num_active_entities;
132 MPI_Allreduce(&local_num_values, &global_num_values, 1, MPI_INT64_T, MPI_SUM,
133 comm);
134 const std::int64_t num_local = num_active_entities;
135 std::int64_t offset = 0;
136 MPI_Exscan(&num_local, &offset, 1, MPI_INT64_T, MPI_SUM, comm);
137 const bool use_mpi_io = (dolfinx::MPI::size(comm) > 1);
138 xdmf_utils::add_data_item(
139 attribute_node, h5_id, path_prefix + std::string("/Values"),
140 std::span<const T>(meshtags.values().data(), num_active_entities), offset,
141 {global_num_values, 1}, "", use_mpi_io);
142}
143} // namespace io::xdmf_mesh
144} // namespace dolfinx
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:33
MeshTags associate values with mesh topology entities.
Definition MeshTags.h:33
std::span< const std::int32_t > indices() const
Definition MeshTags.h:101
std::span< const T > values() const
Values attached to topology entities.
Definition MeshTags.h:104
int dim() const
Return topological dimension of tagged entities.
Definition MeshTags.h:107
std::shared_ptr< const Topology > topology() const
Return topology.
Definition MeshTags.h:110
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relat...
Definition Topology.h:44
int size(MPI_Comm comm)
Definition MPI.cpp:72
std::pair< std::vector< std::int64_t >, std::array< std::size_t, 2 > > read_topology_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node &node)
Read topology (cell connectivity) data.
Definition xdmf_mesh.cpp:294
std::pair< std::variant< std::vector< float >, std::vector< double > >, std::array< std::size_t, 2 > > read_geometry_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node &node)
Read geometry (coordinate) data.
Definition xdmf_mesh.cpp:254
void add_meshtags(MPI_Comm comm, const mesh::MeshTags< T > &meshtags, const mesh::Geometry< U > &geometry, pugi::xml_node &xml_node, hid_t h5_id, const std::string &name)
Add mesh tags to XDMF file.
Definition xdmf_mesh.h:96
void add_geometry_data(MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, std::string path_prefix, const mesh::Geometry< U > &geometry)
Add Geometry xml node.
Definition xdmf_mesh.cpp:160
void add_mesh(MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, const mesh::Mesh< U > &mesh, const std::string &path_prefix)
Definition xdmf_mesh.cpp:214
void add_topology_data(MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, std::string path_prefix, const mesh::Topology &topology, const mesh::Geometry< U > &geometry, int cell_dim, std::span< const std::int32_t > entities)
Definition xdmf_mesh.cpp:23
Top-level namespace.
Definition defines.h:12