OGS
LocalToGlobalIndexMap.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#ifndef NDEBUG
7#include <iosfwd>
8#endif // NDEBUG
9
10#include <Eigen/Core>
11#include <vector>
12
14#include "MeshComponentMap.h"
15
16namespace MeshLib
17{
18class Node;
19class Element;
20} // namespace MeshLib
21
22namespace NumLib
23{
34{
35 // Enables using std::make_unique with private constructors from within
36 // member functions of LocalToGlobalIndexMap. Cf.
37 // https://web.archive.org/web/20210302195245/https://seanmiddleditch.com/enabling-make-unique-with-private-constructors/
39 {
40 explicit ConstructorTag() = default;
41 };
42
43public:
46
47public:
54 LocalToGlobalIndexMap(std::vector<MeshLib::MeshSubset>&& mesh_subsets,
55 NumLib::ComponentOrder const order);
56
65 LocalToGlobalIndexMap(std::vector<MeshLib::MeshSubset>&& mesh_subsets,
66 std::vector<int> const& vec_var_n_components,
67 NumLib::ComponentOrder const order);
68
79 std::vector<MeshLib::MeshSubset>&& mesh_subsets,
80 std::vector<int> const& vec_var_n_components,
81 std::vector<std::vector<MeshLib::Element*> const*> const&
82 vec_var_elements,
83 NumLib::ComponentOrder const order);
84
88 std::unique_ptr<LocalToGlobalIndexMap> deriveBoundaryConstrainedMap(
89 int const variable_id,
90 std::vector<int> const& component_ids,
91 MeshLib::MeshSubset&& new_mesh_subset) const;
92
97 std::unique_ptr<LocalToGlobalIndexMap> deriveBoundaryConstrainedMap(
98 MeshLib::MeshSubset&& new_mesh_subset) const;
99
102 std::size_t dofSizeWithGhosts() const;
103
107 std::size_t dofSizeWithoutGhosts() const;
108
109 std::size_t size() const;
110
111 int getNumberOfVariables() const;
112
113 int getNumberOfVariableComponents(int variable_id) const;
114
115 int getNumberOfGlobalComponents() const;
116
117 RowColumnIndices operator()(std::size_t const mesh_item_id,
118 const int global_component_id) const;
119
120 std::size_t getNumberOfElementDOF(std::size_t const mesh_item_id) const;
121
123 std::size_t const mesh_item_id) const;
124
125 std::vector<int> getElementVariableIDs(
126 std::size_t const mesh_item_id) const;
127
129 int const variable_id,
130 int const component_id) const;
131
133 int const global_component_id) const;
134
136 std::vector<GlobalIndexType> getGlobalIndices(
137 const MeshLib::Location& l) const;
138
140 std::vector<GlobalIndexType> const& getGhostIndices() const;
141
145 std::size_t const comp_id,
146 std::size_t const range_begin,
147 std::size_t const range_end) const;
148
149 MeshLib::MeshSubset const& getMeshSubset(int const variable_id,
150 int const component_id) const;
151
153 int const global_component_id) const;
154
157 int getGlobalComponent(int const variable_id, int const component_id) const;
158
165 explicit LocalToGlobalIndexMap(
166 std::vector<MeshLib::MeshSubset>&& mesh_subsets,
167 std::vector<int> const& global_component_ids,
168 std::vector<int> const& variable_component_offsets,
169 std::vector<MeshLib::Element*> const& elements,
170 NumLib::MeshComponentMap&& mesh_component_map,
171 ConstructorTag /*unused*/);
172
173private:
174 template <typename ElementIterator>
175 void findGlobalIndices(ElementIterator first, ElementIterator last,
176 std::vector<MeshLib::Node*> const& nodes,
177 std::size_t const mesh_id, const int comp_id,
178 const int comp_id_write);
179
180 template <typename ElementIterator>
182 ElementIterator first, ElementIterator last,
183 std::vector<MeshLib::Node*> const& nodes, std::size_t const mesh_id,
184 const int comp_id, const int comp_id_write);
185
187 std::vector<MeshLib::MeshSubset> _mesh_subsets;
189
190 using Table = Eigen::Matrix<LineIndex, Eigen::Dynamic, Eigen::Dynamic,
191 Eigen::RowMajor>;
192
197
200
201 std::vector<int> const _variable_component_offsets;
202#ifndef NDEBUG
204 friend std::ostream& operator<<(std::ostream& os,
205 LocalToGlobalIndexMap const& map);
206#endif // NDEBUG
207};
208
209} // namespace NumLib
210
211#ifndef NDEBUG
212template <>
213struct fmt::formatter<NumLib::LocalToGlobalIndexMap> : fmt::ostream_formatter
214{
215};
216#endif // NDEBUG
GlobalMatrix::IndexType GlobalIndexType
A subset of nodes on a single mesh.
Definition MeshSubset.h:17
RowColumnIndices operator()(std::size_t const mesh_item_id, const int global_component_id) const
std::size_t getNumberOfElementDOF(std::size_t const mesh_item_id) const
LocalToGlobalIndexMap(std::vector< MeshLib::MeshSubset > &&mesh_subsets, NumLib::ComponentOrder const order)
std::vector< MeshLib::MeshSubset > _mesh_subsets
A vector of mesh subsets for each process variables' components.
int getGlobalComponent(int const variable_id, int const component_id) const
NumLib::MeshComponentMap _mesh_component_map
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices
GlobalIndexType getGlobalIndex(MeshLib::Location const &l, int const variable_id, int const component_id) const
void findGlobalIndicesWithElementID(ElementIterator first, ElementIterator last, std::vector< MeshLib::Node * > const &nodes, std::size_t const mesh_id, const int comp_id, const int comp_id_write)
std::vector< int > const _variable_component_offsets
int getNumberOfVariableComponents(int variable_id) const
Eigen::Matrix< LineIndex, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > Table
std::unique_ptr< LocalToGlobalIndexMap > deriveBoundaryConstrainedMap(int const variable_id, std::vector< int > const &component_ids, MeshLib::MeshSubset &&new_mesh_subset) const
MeshLib::MeshSubset const & getMeshSubset(int const variable_id, int const component_id) const
std::size_t getNumberOfElementComponents(std::size_t const mesh_item_id) const
friend std::ostream & operator<<(std::ostream &os, LocalToGlobalIndexMap const &map)
Prints first rows of the table, every line, and the mesh component map.
void findGlobalIndices(ElementIterator first, ElementIterator last, std::vector< MeshLib::Node * > const &nodes, std::size_t const mesh_id, const int comp_id, const int comp_id_write)
GlobalIndexType getLocalIndex(MeshLib::Location const &l, std::size_t const comp_id, std::size_t const range_begin, std::size_t const range_end) const
std::vector< GlobalIndexType > const & getGhostIndices() const
Get ghost indices, forwarded from MeshComponentMap.
RowColumnIndices::LineIndex LineIndex
std::vector< GlobalIndexType > getGlobalIndices(const MeshLib::Location &l) const
Forwards the respective method from MeshComponentMap.
std::vector< int > getElementVariableIDs(std::size_t const mesh_item_id) const
Multidirectional mapping between mesh entities and degrees of freedom.
ComponentOrder
Ordering of components in global matrix/vector.
typename std::vector< GlobalIndexType > LineIndex