OGS
LocalToGlobalIndexMap.h
Go to the documentation of this file.
1
11#pragma once
12
13#ifndef NDEBUG
14#include <iosfwd>
15#endif // NDEBUG
16
17#include <Eigen/Core>
18#include <vector>
19
21#include "MeshComponentMap.h"
22
23namespace MeshLib
24{
25class Node;
26class Element;
27} // namespace MeshLib
28
29namespace NumLib
30{
41{
42 // Enables using std::make_unique with private constructors from within
43 // member functions of LocalToGlobalIndexMap. Cf.
44 // https://web.archive.org/web/20210302195245/https://seanmiddleditch.com/enabling-make-unique-with-private-constructors/
46 {
47 explicit ConstructorTag() = default;
48 };
49
50public:
53
54public:
61 LocalToGlobalIndexMap(std::vector<MeshLib::MeshSubset>&& mesh_subsets,
62 NumLib::ComponentOrder const order);
63
72 LocalToGlobalIndexMap(std::vector<MeshLib::MeshSubset>&& mesh_subsets,
73 std::vector<int> const& vec_var_n_components,
74 NumLib::ComponentOrder const order);
75
86 std::vector<MeshLib::MeshSubset>&& mesh_subsets,
87 std::vector<int> const& vec_var_n_components,
88 std::vector<std::vector<MeshLib::Element*> const*> const&
89 vec_var_elements,
90 NumLib::ComponentOrder const order);
91
95 std::unique_ptr<LocalToGlobalIndexMap> deriveBoundaryConstrainedMap(
96 int const variable_id,
97 std::vector<int> const& component_ids,
98 MeshLib::MeshSubset&& new_mesh_subset) const;
99
104 std::unique_ptr<LocalToGlobalIndexMap> deriveBoundaryConstrainedMap(
105 MeshLib::MeshSubset&& new_mesh_subset) const;
106
109 std::size_t dofSizeWithGhosts() const;
110
114 std::size_t dofSizeWithoutGhosts() const;
115
116 std::size_t size() const;
117
118 int getNumberOfVariables() const;
119
120 int getNumberOfVariableComponents(int variable_id) const;
121
122 int getNumberOfGlobalComponents() const;
123
124 RowColumnIndices operator()(std::size_t const mesh_item_id,
125 const int global_component_id) const;
126
127 std::size_t getNumberOfElementDOF(std::size_t const mesh_item_id) const;
128
130 std::size_t const mesh_item_id) const;
131
132 std::vector<int> getElementVariableIDs(
133 std::size_t const mesh_item_id) const;
134
136 int const variable_id,
137 int const component_id) const;
138
140 int const global_component_id) const;
141
143 std::vector<GlobalIndexType> getGlobalIndices(
144 const MeshLib::Location& l) const;
145
147 std::vector<GlobalIndexType> const& getGhostIndices() const;
148
152 std::size_t const comp_id,
153 std::size_t const range_begin,
154 std::size_t const range_end) const;
155
156 MeshLib::MeshSubset const& getMeshSubset(int const variable_id,
157 int const component_id) const;
158
160 int const global_component_id) const;
161
164 int getGlobalComponent(int const variable_id, int const component_id) const;
165
172 explicit LocalToGlobalIndexMap(
173 std::vector<MeshLib::MeshSubset>&& mesh_subsets,
174 std::vector<int> const& global_component_ids,
175 std::vector<int> const& variable_component_offsets,
176 std::vector<MeshLib::Element*> const& elements,
177 NumLib::MeshComponentMap&& mesh_component_map,
178 ConstructorTag /*unused*/);
179
180private:
181 template <typename ElementIterator>
182 void findGlobalIndices(ElementIterator first, ElementIterator last,
183 std::vector<MeshLib::Node*> const& nodes,
184 std::size_t const mesh_id, const int comp_id,
185 const int comp_id_write);
186
187 template <typename ElementIterator>
189 ElementIterator first, ElementIterator last,
190 std::vector<MeshLib::Node*> const& nodes, std::size_t const mesh_id,
191 const int comp_id, const int comp_id_write);
192
194 std::vector<MeshLib::MeshSubset> _mesh_subsets;
196
197 using Table = Eigen::Matrix<LineIndex, Eigen::Dynamic, Eigen::Dynamic,
198 Eigen::RowMajor>;
199
204
207
208 std::vector<int> const _variable_component_offsets;
209#ifndef NDEBUG
211 friend std::ostream& operator<<(std::ostream& os,
212 LocalToGlobalIndexMap const& map);
213#endif // NDEBUG
214};
215
216} // namespace NumLib
GlobalMatrix::IndexType GlobalIndexType
A subset of nodes on a single mesh.
Definition MeshSubset.h:26
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
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< IDX_TYPE > LineIndex