DOLFINx 0.10.0.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 <algorithm>
11#include <array>
12#include <concepts>
13#include <cstdint>
14#include <dolfinx/common/MPI.h>
15#include <dolfinx/mesh/MeshTags.h>
16#include <hdf5.h>
17#include <mpi.h>
18#include <pugixml.hpp>
19#include <span>
20#include <string>
21#include <utility>
22#include <variant>
23#include <vector>
24
25namespace pugi
26{
27class xml_node;
28} // namespace pugi
29
30namespace dolfinx
31{
32
33namespace mesh
34{
35template <std::floating_point T>
36class Geometry;
37template <std::floating_point T>
38class Mesh;
39class Topology;
40} // namespace mesh
41
43namespace io::xdmf_mesh
44{
45
50template <std::floating_point U>
51void add_mesh(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
52 const mesh::Mesh<U>& mesh, const std::string& path_prefix);
53
64template <std::floating_point U>
65void add_topology_data(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
66 std::string path_prefix, const mesh::Topology& topology,
67 const mesh::Geometry<U>& geometry, int cell_dim,
68 std::span<const std::int32_t> entities);
69
71template <std::floating_point U>
72void add_geometry_data(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
73 std::string path_prefix,
74 const mesh::Geometry<U>& geometry);
75
82std::pair<std::variant<std::vector<float>, std::vector<double>>,
83 std::array<std::size_t, 2>>
84read_geometry_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node& node);
85
92std::pair<std::vector<std::int64_t>, std::array<std::size_t, 2>>
93read_topology_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node& node);
94
96template <typename T, std::floating_point U>
97void add_meshtags(MPI_Comm comm, const mesh::MeshTags<T>& meshtags,
98 const mesh::Geometry<U>& geometry, pugi::xml_node& xml_node,
99 hid_t h5_id, const std::string& name)
100{
101 spdlog::info("XDMF: add meshtags ({})", name.c_str());
102 // Get mesh
103 const int dim = meshtags.dim();
104 std::shared_ptr<const common::IndexMap> entity_map
105 = meshtags.topology()->index_map(dim);
106 if (!entity_map)
107 {
108 throw std::runtime_error("Missing entities. Did you forget to call "
109 "dolfinx::mesh::Topology::create_entities?");
110 }
111 const std::int32_t num_local_entities = entity_map->size_local();
112
113 // Find number of tagged entities in local range
114 auto it = std::ranges::lower_bound(meshtags.indices(), 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:34
MeshTags associate values with mesh topology entities.
Definition utils.h:26
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:97
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