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 <array>
11#include <cstdint>
12#include <dolfinx/mesh/MeshTags.h>
13#include <hdf5.h>
14#include <mpi.h>
15#include <pugixml.hpp>
16#include <span>
17#include <string>
18#include <utility>
19#include <variant>
20#include <vector>
21
22namespace pugi
23{
24class xml_node;
25} // namespace pugi
26
27namespace dolfinx
28{
29
30namespace mesh
31{
32template <std::floating_point T>
33class Geometry;
34template <std::floating_point T>
35class Mesh;
36class Topology;
37} // namespace mesh
38
41{
42
47template <std::floating_point U>
48void add_mesh(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
49 const mesh::Mesh<U>& mesh, const std::string& path_prefix);
50
61template <std::floating_point U>
62void add_topology_data(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
63 const std::string& path_prefix,
64 const mesh::Topology& topology,
65 const mesh::Geometry<U>& geometry, int cell_dim,
66 std::span<const std::int32_t> entities);
67
69template <std::floating_point U>
70void add_geometry_data(MPI_Comm comm, pugi::xml_node& xml_node, hid_t h5_id,
71 const std::string& path_prefix,
73
80std::pair<std::variant<std::vector<float>, std::vector<double>>,
81 std::array<std::size_t, 2>>
82read_geometry_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node& node);
83
90std::pair<std::vector<std::int64_t>, std::array<std::size_t, 2>>
91read_topology_data(MPI_Comm comm, hid_t h5_id, const pugi::xml_node& node);
92
94template <typename T, std::floating_point U>
95void add_meshtags(MPI_Comm comm, const mesh::MeshTags<T>& meshtags,
96 const mesh::Geometry<U>& geometry, pugi::xml_node& xml_node,
97 hid_t h5_id, const std::string& name)
98{
99 spdlog::info("XDMF: add meshtags ({})", name.c_str());
100 // Get mesh
101 const int dim = meshtags.dim();
102 std::shared_ptr<const common::IndexMap> entity_map
103 = meshtags.topology()->index_map(dim);
104 if (!entity_map)
105 {
106 throw std::runtime_error("Missing entities. Did you forget to call "
107 "dolfinx::mesh::Topology::create_entities?");
108 }
109 const std::int32_t num_local_entities = entity_map->size_local();
110
111 // Find number of tagged entities in local range
112 auto it = std::ranges::lower_bound(meshtags.indices(), num_local_entities);
113 const int num_active_entities = std::distance(meshtags.indices().begin(), it);
114
115 const std::string path_prefix = "/MeshTags/" + name;
117 comm, xml_node, h5_id, path_prefix, *meshtags.topology(), geometry, dim,
118 std::span<const std::int32_t>(meshtags.indices().data(),
119 num_active_entities));
120
121 // Add attribute node with values
122 pugi::xml_node attribute_node = xml_node.append_child("Attribute");
123 assert(attribute_node);
124 attribute_node.append_attribute("Name") = name.c_str();
125 attribute_node.append_attribute("AttributeType") = "Scalar";
126 attribute_node.append_attribute("Center") = "Cell";
127
128 std::int64_t global_num_values = 0;
129 const std::int64_t local_num_values = num_active_entities;
130 MPI_Allreduce(&local_num_values, &global_num_values, 1, MPI_INT64_T, MPI_SUM,
131 comm);
132 const std::int64_t num_local = num_active_entities;
133 std::int64_t offset = 0;
134 MPI_Exscan(&num_local, &offset, 1, MPI_INT64_T, MPI_SUM, comm);
135 const bool use_mpi_io = (dolfinx::MPI::size(comm) > 1);
136 xdmf_utils::add_data_item(
137 attribute_node, h5_id, path_prefix + std::string("/Values"),
138 std::span<const T>(meshtags.values().data(), num_active_entities), offset,
139 {global_num_values, 1}, "", use_mpi_io);
140}
141} // namespace io::xdmf_mesh
142} // namespace dolfinx
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:34
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:41
int size(MPI_Comm comm)
Definition MPI.cpp:72
Geometry data structures and algorithms.
Definition BoundingBoxTree.h:22
Low-level methods for reading XDMF files.
Definition xdmf_mesh.h:41
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:95
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, const 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
void add_geometry_data(MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, const std::string &path_prefix, const mesh::Geometry< U > &geometry)
Add Geometry xml node.
Definition xdmf_mesh.cpp:160
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
Top-level namespace.
Definition defines.h:12