OGS
NumLib::LocalToGlobalIndexMap Class Referencefinal

Detailed Description

Row and column indices in global linear algebra objects for each mesh item.

The row and column indices in global matrix and rhs vector for example, required for addition of local contributions from every mesh item (like node or cell) to global objects.

The number of rows should be equal to the number of mesh items and the number of columns should be equal to the number of the components on that mesh item.

Definition at line 33 of file LocalToGlobalIndexMap.h.

#include <LocalToGlobalIndexMap.h>

Collaboration diagram for NumLib::LocalToGlobalIndexMap:
[legend]

Classes

struct  ConstructorTag

Public Types

using RowColumnIndices = MathLib::RowColumnIndices<GlobalIndexType>
using LineIndex = RowColumnIndices::LineIndex

Public Member Functions

 LocalToGlobalIndexMap (std::vector< MeshLib::MeshSubset > &&mesh_subsets, NumLib::ComponentOrder const order)
 LocalToGlobalIndexMap (std::vector< MeshLib::MeshSubset > &&mesh_subsets, std::vector< int > const &vec_var_n_components, NumLib::ComponentOrder const order)
 LocalToGlobalIndexMap (std::vector< MeshLib::MeshSubset > &&mesh_subsets, std::vector< int > const &vec_var_n_components, std::vector< std::vector< MeshLib::Element * > const * > const &vec_var_elements, NumLib::ComponentOrder const order)
std::unique_ptr< LocalToGlobalIndexMapderiveBoundaryConstrainedMap (int const variable_id, std::vector< int > const &component_ids, MeshLib::MeshSubset &&new_mesh_subset) const
std::unique_ptr< LocalToGlobalIndexMapderiveBoundaryConstrainedMap (MeshLib::MeshSubset &&new_mesh_subset) const
std::size_t dofSizeWithGhosts () const
std::size_t dofSizeWithoutGhosts () const
std::size_t size () const
int getNumberOfVariables () const
int getNumberOfVariableComponents (int variable_id) const
int getNumberOfGlobalComponents () const
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
std::size_t getNumberOfElementComponents (std::size_t const mesh_item_id) const
std::vector< int > getElementVariableIDs (std::size_t const mesh_item_id) const
GlobalIndexType getGlobalIndex (MeshLib::Location const &l, int const variable_id, int const component_id) const
GlobalIndexType getGlobalIndex (MeshLib::Location const &l, int const global_component_id) const
std::vector< GlobalIndexTypegetGlobalIndices (const MeshLib::Location &l) const
 Forwards the respective method from MeshComponentMap.
std::vector< GlobalIndexType > const & getGhostIndices () const
 Get ghost indices, forwarded from MeshComponentMap.
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
MeshLib::MeshSubset const & getMeshSubset (int const variable_id, int const component_id) const
MeshLib::MeshSubset const & getMeshSubset (int const global_component_id) const
int getGlobalComponent (int const variable_id, int const component_id) const
 LocalToGlobalIndexMap (std::vector< MeshLib::MeshSubset > &&mesh_subsets, std::vector< int > const &global_component_ids, std::vector< int > const &variable_component_offsets, std::vector< MeshLib::Element * > const &elements, NumLib::MeshComponentMap &&mesh_component_map, ConstructorTag)

Private Types

using Table

Private Member Functions

template<typename ElementIterator>
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)
template<typename ElementIterator>
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)

Private Attributes

std::vector< MeshLib::MeshSubset_mesh_subsets
 A vector of mesh subsets for each process variables' components.
NumLib::MeshComponentMap _mesh_component_map
Table _rows
Table const & _columns = _rows
std::vector< int > const _variable_component_offsets

Friends

std::ostream & operator<< (std::ostream &os, LocalToGlobalIndexMap const &map)
 Prints first rows of the table, every line, and the mesh component map.

Member Typedef Documentation

◆ LineIndex

◆ RowColumnIndices

◆ Table

Initial value:
Eigen::Matrix<LineIndex, Eigen::Dynamic, Eigen::Dynamic,
Eigen::RowMajor>
RowColumnIndices::LineIndex LineIndex

Definition at line 190 of file LocalToGlobalIndexMap.h.

Constructor & Destructor Documentation

◆ LocalToGlobalIndexMap() [1/4]

NumLib::LocalToGlobalIndexMap::LocalToGlobalIndexMap ( std::vector< MeshLib::MeshSubset > && mesh_subsets,
NumLib::ComponentOrder const order )

Creates a MeshComponentMap internally and stores the global indices for each mesh element of the given mesh_subsets.

Attention
This constructor assumes the number of the given mesh subsets is equal to the number of variables, i.e. every variable has a single component.

Definition at line 98 of file LocalToGlobalIndexMap.cpp.

101 : LocalToGlobalIndexMap(std::move(mesh_subsets),
102 std::vector<int>(mesh_subsets.size(), 1), order)
103{
104}
LocalToGlobalIndexMap(std::vector< MeshLib::MeshSubset > &&mesh_subsets, NumLib::ComponentOrder const order)

References LocalToGlobalIndexMap(), and size().

Referenced by LocalToGlobalIndexMap(), and operator<<.

◆ LocalToGlobalIndexMap() [2/4]

NumLib::LocalToGlobalIndexMap::LocalToGlobalIndexMap ( std::vector< MeshLib::MeshSubset > && mesh_subsets,
std::vector< int > const & vec_var_n_components,
NumLib::ComponentOrder const order )

Creates a MeshComponentMap internally and stores the global indices for each mesh element of the given mesh_subsets.

Parameters
mesh_subsetsa vector of components
vec_var_n_componentsa vector of the number of variable components. The size of the vector should be equal to the number of variables. Sum of the entries should be equal to the size of the mesh_subsets.
ordertype of ordering values in a vector

Definition at line 106 of file LocalToGlobalIndexMap.cpp.

110 : _mesh_subsets(std::move(mesh_subsets)),
113{
114 // For each element of that MeshSubset save a line of global indices.
115 for (int variable_id = 0;
116 variable_id < static_cast<int>(vec_var_n_components.size());
117 ++variable_id)
118 {
119 for (int component_id = 0;
120 component_id < static_cast<int>(vec_var_n_components[variable_id]);
121 ++component_id)
122 {
123 auto const global_component_id =
124 getGlobalComponent(variable_id, component_id);
125
126 auto const& ms = _mesh_subsets[global_component_id];
127 std::size_t const mesh_id = ms.getMeshID();
128
129 findGlobalIndices(ms.elementsBegin(), ms.elementsEnd(),
130 ms.getNodes(), mesh_id, global_component_id,
131 global_component_id);
132 }
133 }
134}
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
std::vector< int > const _variable_component_offsets
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)
std::vector< ranges::range_value_t< R > > sizesToOffsets(R const &sizes)
Definition Algorithm.h:276

References _mesh_component_map, _mesh_subsets, _variable_component_offsets, findGlobalIndices(), and getGlobalComponent().

◆ LocalToGlobalIndexMap() [3/4]

NumLib::LocalToGlobalIndexMap::LocalToGlobalIndexMap ( std::vector< MeshLib::MeshSubset > && mesh_subsets,
std::vector< int > const & vec_var_n_components,
std::vector< std::vector< MeshLib::Element * > const * > const & vec_var_elements,
NumLib::ComponentOrder const order )

Creates a MeshComponentMap internally and stores the global indices for the given mesh elements

Parameters
mesh_subsetsa vector of components
vec_var_n_componentsa vector of the number of variable components. The size of the vector should be equal to the number of variables. Sum of the entries should be equal to the size of the mesh_subsets.
vec_var_elementsa vector of active mesh elements for each variable.
ordertype of ordering values in a vector

Definition at line 136 of file LocalToGlobalIndexMap.cpp.

141 : _mesh_subsets(std::move(mesh_subsets)),
144{
145 assert(vec_var_n_components.size() == vec_var_elements.size());
146
147 // For each element of that MeshSubset save a line of global indices.
148
149 // _rows should be resized based on an element ID
150 std::size_t max_elem_id = 0;
151 for (std::vector<MeshLib::Element*> const* elements : vec_var_elements)
152 {
153 for (auto e : *elements)
154 {
155 max_elem_id = std::max(max_elem_id, e->getID());
156 }
157 }
158 _rows.resize(max_elem_id + 1, _mesh_subsets.size());
159
160 for (int variable_id = 0;
161 variable_id < static_cast<int>(vec_var_n_components.size());
162 ++variable_id)
163 {
164 std::vector<MeshLib::Element*> const& var_elements =
165 *vec_var_elements[variable_id];
166 for (int component_id = 0;
167 component_id < static_cast<int>(vec_var_n_components[variable_id]);
168 ++component_id)
169 {
170 auto const global_component_id =
171 getGlobalComponent(variable_id, component_id);
172
173 auto const& ms = _mesh_subsets[global_component_id];
174 std::size_t const mesh_id = ms.getMeshID();
175
177 var_elements.cbegin(), var_elements.cend(), ms.getNodes(),
178 mesh_id, global_component_id, global_component_id);
179 }
180 }
181}
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)

References _mesh_component_map, _mesh_subsets, _rows, _variable_component_offsets, findGlobalIndicesWithElementID(), and getGlobalComponent().

◆ LocalToGlobalIndexMap() [4/4]

NumLib::LocalToGlobalIndexMap::LocalToGlobalIndexMap ( std::vector< MeshLib::MeshSubset > && mesh_subsets,
std::vector< int > const & global_component_ids,
std::vector< int > const & variable_component_offsets,
std::vector< MeshLib::Element * > const & elements,
NumLib::MeshComponentMap && mesh_component_map,
LocalToGlobalIndexMap::ConstructorTag  )
explicit

Private constructor (ensured by ConstructorTag) used by internally created local-to-global index maps. The mesh_component_map is passed as argument instead of being created by the constructor.

Attention
The passed mesh_component_map is in undefined state after this constructor.

Definition at line 183 of file LocalToGlobalIndexMap.cpp.

190 : _mesh_subsets(std::move(mesh_subsets)),
191 _mesh_component_map(std::move(mesh_component_map)),
192 _variable_component_offsets(variable_component_offsets)
193{
194 // Each subset in the mesh_subsets represents a single component.
195 if (_mesh_subsets.size() != global_component_ids.size())
196 {
197 OGS_FATAL(
198 "Number of mesh subsets is not equal to number of components. "
199 "There are {:d} mesh subsets and {:d} components.",
200 _mesh_subsets.size(), global_component_ids.size());
201 }
202
203 for (int i = 0; i < static_cast<int>(global_component_ids.size()); ++i)
204 {
205 auto const& ms = _mesh_subsets[i];
206
207 // For all MeshSubset in mesh_subsets and each element of that
208 // MeshSubset save a line of global indices.
209 std::size_t const mesh_id = ms.getMeshID();
210
211 findGlobalIndices(elements.cbegin(), elements.cend(), ms.getNodes(),
212 mesh_id, global_component_ids[i], i);
213 }
214}
#define OGS_FATAL(...)
Definition Error.h:19

References _mesh_component_map, _mesh_subsets, _variable_component_offsets, and findGlobalIndices().

Member Function Documentation

◆ deriveBoundaryConstrainedMap() [1/2]

std::unique_ptr< LocalToGlobalIndexMap > NumLib::LocalToGlobalIndexMap::deriveBoundaryConstrainedMap ( int const variable_id,
std::vector< int > const & component_ids,
MeshLib::MeshSubset && new_mesh_subset ) const

Derive a LocalToGlobalIndexMap constrained to the mesh subset and mesh subset's elements. A new mesh component map will be constructed using the passed mesh_subset for the given variable and component ids.

Definition at line 217 of file LocalToGlobalIndexMap.cpp.

221{
222 DBUG("Construct reduced local to global index map.");
223
224 if (component_ids.empty())
225 {
226 OGS_FATAL("Expected non-empty vector of component ids.");
227 }
228
229 // Elements of the new_mesh_subset's mesh.
230 std::vector<MeshLib::Element*> const& elements =
231 new_mesh_subset.getMesh().getElements();
232
233 // Create a subset of the current mesh component map.
234 std::vector<int> global_component_ids;
235
236 transform(cbegin(component_ids), cend(component_ids),
237 back_inserter(global_component_ids),
238 [&](auto const component_id)
239 { return getGlobalComponent(variable_id, component_id); });
240
241 auto mesh_component_map =
242 _mesh_component_map.getSubset(_mesh_subsets, new_mesh_subset);
243
244 // Create copies of the new_mesh_subset for each of the global components.
245 // The last component is moved after the for-loop.
246 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
247 for (int i = 0; i < static_cast<int>(global_component_ids.size()) - 1; ++i)
248 {
249 all_mesh_subsets.emplace_back(new_mesh_subset);
250 }
251 all_mesh_subsets.emplace_back(std::move(new_mesh_subset));
252
253 return std::make_unique<LocalToGlobalIndexMap>(
254 std::move(all_mesh_subsets), global_component_ids,
255 _variable_component_offsets, elements, std::move(mesh_component_map),
257}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
Mesh const & getMesh() const
Definition MeshSubset.h:83
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:100

References _mesh_component_map, _mesh_subsets, _variable_component_offsets, DBUG(), getGlobalComponent(), and OGS_FATAL.

Referenced by ProcessLib::ConstraintDirichletBoundaryCondition::ConstraintDirichletBoundaryCondition(), ProcessLib::DirichletBoundaryCondition::DirichletBoundaryCondition(), ProcessLib::PrimaryVariableConstraintDirichletBoundaryCondition::PrimaryVariableConstraintDirichletBoundaryCondition(), ProcessLib::PythonBoundaryCondition::PythonBoundaryCondition(), ProcessLib::SolutionDependentDirichletBoundaryCondition::SolutionDependentDirichletBoundaryCondition(), anonymous_namespace{ProcessOutputData.cpp}::computeDofTablesForSubmesh(), ProcessLib::DeactivatedSubdomainDirichlet::config(), ProcessLib::DirichletBoundaryConditionWithinTimeInterval::config(), ProcessLib::createBoundaryDOFTable(), ProcessLib::createReleaseNodalForce(), ProcessLib::createSourceTerm(), ProcessLib::createVariableDependentNeumannBoundaryCondition(), and ProcessLib::createWellboreCompensateNeumannBoundaryCondition().

◆ deriveBoundaryConstrainedMap() [2/2]

std::unique_ptr< LocalToGlobalIndexMap > NumLib::LocalToGlobalIndexMap::deriveBoundaryConstrainedMap ( MeshLib::MeshSubset && new_mesh_subset) const

Derive a LocalToGlobalIndexMap constrained to the mesh subset and mesh subset's elements. A new mesh component map will be constructed using the passed mesh_subset for all variables and components of the current LocalToGlobalIndexMap.

Definition at line 260 of file LocalToGlobalIndexMap.cpp.

262{
263 DBUG("Construct reduced local to global index map.");
264
265 // Create a subset of the current mesh component map.
266 std::vector<int> global_component_ids;
267
268 for (int i = 0; i < getNumberOfGlobalComponents(); ++i)
269 {
270 global_component_ids.push_back(i);
271 }
272
273 // Elements of the new_mesh_subset's mesh.
274 std::vector<MeshLib::Element*> const& elements =
275 new_mesh_subset.getMesh().getElements();
276
277 auto mesh_component_map =
278 _mesh_component_map.getSubset(_mesh_subsets, new_mesh_subset);
279
280 // Create copies of the new_mesh_subset for each of the global components.
281 // The last component is moved after the for-loop.
282 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
283 for (int i = 0; i < static_cast<int>(global_component_ids.size()) - 1; ++i)
284 {
285 all_mesh_subsets.emplace_back(new_mesh_subset);
286 }
287 all_mesh_subsets.emplace_back(std::move(new_mesh_subset));
288
289 return std::make_unique<LocalToGlobalIndexMap>(
290 std::move(all_mesh_subsets), global_component_ids,
291 _variable_component_offsets, elements, std::move(mesh_component_map),
293}

References _mesh_component_map, _mesh_subsets, _variable_component_offsets, DBUG(), and getNumberOfGlobalComponents().

◆ dofSizeWithGhosts()

std::size_t NumLib::LocalToGlobalIndexMap::dofSizeWithGhosts ( ) const

Returns total number of degrees of freedom including those located in the ghost nodes.

Definition at line 295 of file LocalToGlobalIndexMap.cpp.

296{
297 return _mesh_component_map.dofSizeWithGhosts();
298}

References _mesh_component_map.

◆ dofSizeWithoutGhosts()

std::size_t NumLib::LocalToGlobalIndexMap::dofSizeWithoutGhosts ( ) const

Returns total number of local degrees of freedom of the present rank, which does not count the unknowns associated with ghost nodes (for DDC with node-wise mesh partitioning).

Definition at line 300 of file LocalToGlobalIndexMap.cpp.

301{
302 return _mesh_component_map.dofSizeWithoutGhosts();
303}

References _mesh_component_map.

◆ findGlobalIndices()

template<typename ElementIterator>
void NumLib::LocalToGlobalIndexMap::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 )
private

Definition at line 55 of file LocalToGlobalIndexMap.cpp.

59{
60 _rows.resize(std::distance(first, last), _mesh_subsets.size());
61
62 std::unordered_set<MeshLib::Node*> const set_nodes(nodes.begin(),
63 nodes.end());
64
65 // For each element find the global indices for node/element
66 // components.
67 std::size_t elem_id = 0;
68 for (ElementIterator e = first; e != last; ++e, ++elem_id)
69 {
70 LineIndex indices;
71 indices.reserve((*e)->getNumberOfNodes());
72
73 for (auto* n = (*e)->getNodes();
74 n < (*e)->getNodes() + (*e)->getNumberOfNodes();
75 ++n)
76 {
77 // Check if the element's node is in the given list of nodes.
78 if (set_nodes.find(*n) == set_nodes.end())
79 {
80 continue;
81 }
82 MeshLib::Location l(mesh_id, MeshLib::MeshItemType::Node,
83 (*n)->getID());
84 auto const global_index =
85 _mesh_component_map.getGlobalIndex(l, comp_id);
86 if (global_index == std::numeric_limits<GlobalIndexType>::max())
87 {
88 continue;
89 }
90 indices.push_back(global_index);
91 }
92
93 indices.shrink_to_fit();
94 _rows(elem_id, comp_id_write) = std::move(indices);
95 }
96}

References _mesh_component_map, _mesh_subsets, _rows, and MeshLib::Node.

Referenced by LocalToGlobalIndexMap(), and LocalToGlobalIndexMap().

◆ findGlobalIndicesWithElementID()

template<typename ElementIterator>
void NumLib::LocalToGlobalIndexMap::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 )
private

Definition at line 20 of file LocalToGlobalIndexMap.cpp.

24{
25 std::unordered_set<MeshLib::Node*> const set_nodes(nodes.begin(),
26 nodes.end());
27
28 // For each element find the global indices for node/element
29 // components.
30 for (ElementIterator e = first; e != last; ++e)
31 {
32 LineIndex indices;
33 indices.reserve((*e)->getNumberOfNodes());
34
35 for (auto* n = (*e)->getNodes();
36 n < (*e)->getNodes() + (*e)->getNumberOfNodes();
37 ++n)
38 {
39 // Check if the element's node is in the given list of nodes.
40 if (set_nodes.find(*n) == set_nodes.end())
41 {
42 continue;
43 }
44 MeshLib::Location l(mesh_id, MeshLib::MeshItemType::Node,
45 (*n)->getID());
46 indices.push_back(_mesh_component_map.getGlobalIndex(l, comp_id));
47 }
48
49 indices.shrink_to_fit();
50 _rows((*e)->getID(), comp_id_write) = std::move(indices);
51 }
52}

References _mesh_component_map, _rows, and MeshLib::Node.

Referenced by LocalToGlobalIndexMap().

◆ getElementVariableIDs()

std::vector< int > NumLib::LocalToGlobalIndexMap::getElementVariableIDs ( std::size_t const mesh_item_id) const

Definition at line 361 of file LocalToGlobalIndexMap.cpp.

363{
364 std::vector<int> vec;
365 for (int i = 0; i < getNumberOfVariables(); i++)
366 {
367 for (int j = 0; j < getNumberOfVariableComponents(i); j++)
368 {
369 auto comp_id = getGlobalComponent(i, j);
370 if (!_rows(mesh_item_id, comp_id).empty())
371 {
372 vec.push_back(i);
373 }
374 }
375 }
376 std::sort(vec.begin(), vec.end());
377 vec.erase(std::unique(vec.begin(), vec.end()), vec.end());
378
379 return vec;
380}
int getNumberOfVariableComponents(int variable_id) const

References _rows, getGlobalComponent(), getNumberOfVariableComponents(), and getNumberOfVariables().

◆ getGhostIndices()

std::vector< GlobalIndexType > const & NumLib::LocalToGlobalIndexMap::getGhostIndices ( ) const

Get ghost indices, forwarded from MeshComponentMap.

Definition at line 405 of file LocalToGlobalIndexMap.cpp.

407{
408 return _mesh_component_map.getGhostIndices();
409}

References _mesh_component_map.

◆ getGlobalComponent()

int NumLib::LocalToGlobalIndexMap::getGlobalComponent ( int const variable_id,
int const component_id ) const

The global component id for the specific variable (like velocity) and a component (like x, or y, or z).

Definition at line 13 of file LocalToGlobalIndexMap.cpp.

15{
16 return _variable_component_offsets[variable_id] + component_id;
17}

References _variable_component_offsets.

Referenced by LocalToGlobalIndexMap(), LocalToGlobalIndexMap(), ProcessLib::BoundaryConditionAndSourceTerm::Python::collectDofsToMatrix(), ProcessLib::createPythonBoundaryCondition(), deriveBoundaryConstrainedMap(), getElementVariableIDs(), getGlobalIndex(), NumLib::getIndices(), and getMeshSubset().

◆ getGlobalIndex() [1/2]

GlobalIndexType NumLib::LocalToGlobalIndexMap::getGlobalIndex ( MeshLib::Location const & l,
int const global_component_id ) const

Definition at line 391 of file LocalToGlobalIndexMap.cpp.

393{
394 return _mesh_component_map.getGlobalIndex(l, global_component_id);
395}

References _mesh_component_map.

◆ getGlobalIndex() [2/2]

◆ getGlobalIndices()

std::vector< GlobalIndexType > NumLib::LocalToGlobalIndexMap::getGlobalIndices ( const MeshLib::Location & l) const

Forwards the respective method from MeshComponentMap.

Definition at line 398 of file LocalToGlobalIndexMap.cpp.

400{
401 return _mesh_component_map.getGlobalIndices(l);
402}

References _mesh_component_map.

◆ getLocalIndex()

GlobalIndexType NumLib::LocalToGlobalIndexMap::getLocalIndex ( MeshLib::Location const & l,
std::size_t const comp_id,
std::size_t const range_begin,
std::size_t const range_end ) const

Computes the index in a local (for DDC) vector for a given location and component; forwarded from MeshComponentMap.

Definition at line 413 of file LocalToGlobalIndexMap.cpp.

416{
417 return _mesh_component_map.getLocalIndex(l, comp_id, range_begin,
418 range_end);
419}

References _mesh_component_map.

Referenced by getIndexForComponentInSolutionVector().

◆ getMeshSubset() [1/2]

MeshLib::MeshSubset const & NumLib::LocalToGlobalIndexMap::getMeshSubset ( int const global_component_id) const

Definition at line 427 of file LocalToGlobalIndexMap.cpp.

429{
430 return _mesh_subsets[global_component_id];
431}

References _mesh_subsets.

◆ getMeshSubset() [2/2]

MeshLib::MeshSubset const & NumLib::LocalToGlobalIndexMap::getMeshSubset ( int const variable_id,
int const component_id ) const

◆ getNumberOfElementComponents()

std::size_t NumLib::LocalToGlobalIndexMap::getNumberOfElementComponents ( std::size_t const mesh_item_id) const

Definition at line 347 of file LocalToGlobalIndexMap.cpp.

349{
350 std::size_t n = 0;
351 for (Table::Index c = 0; c < _rows.cols(); ++c)
352 {
353 if (!_rows(mesh_item_id, c).empty())
354 {
355 n++;
356 }
357 }
358 return n;
359}

References _rows.

◆ getNumberOfElementDOF()

std::size_t NumLib::LocalToGlobalIndexMap::getNumberOfElementDOF ( std::size_t const mesh_item_id) const

Definition at line 334 of file LocalToGlobalIndexMap.cpp.

336{
337 std::size_t ndof = 0;
338
339 for (Table::Index c = 0; c < _rows.cols(); ++c)
340 {
341 ndof += _rows(mesh_item_id, c).size();
342 }
343
344 return ndof;
345}

References _rows.

◆ getNumberOfGlobalComponents()

◆ getNumberOfVariableComponents()

◆ getNumberOfVariables()

◆ operator()()

LocalToGlobalIndexMap::RowColumnIndices NumLib::LocalToGlobalIndexMap::operator() ( std::size_t const mesh_item_id,
const int global_component_id ) const

Definition at line 327 of file LocalToGlobalIndexMap.cpp.

329{
330 return RowColumnIndices(_rows(mesh_item_id, global_component_id),
331 _columns(mesh_item_id, global_component_id));
332}
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices

References _columns, and _rows.

◆ size()

◆ operator<<

std::ostream & operator<< ( std::ostream & os,
LocalToGlobalIndexMap const & map )
friend

Prints first rows of the table, every line, and the mesh component map.

Definition at line 434 of file LocalToGlobalIndexMap.cpp.

435{
436 std::size_t const max_lines = 10;
437 std::size_t lines_printed = 0;
438
439 os << "Rows of the local to global index map; " << map._rows.size()
440 << " rows\n";
441 for (std::size_t e = 0; e < map.size(); ++e)
442 {
443 os << "== e " << e << " ==\n";
444 for (int c = 0; c < map.getNumberOfGlobalComponents(); ++c)
445 {
446 auto const& line = map._rows(e, c);
447
448 os << "c" << c << " { ";
449 std::copy(line.cbegin(), line.cend(),
450 std::ostream_iterator<std::size_t>(os, " "));
451 os << " }\n";
452 }
453
454 if (lines_printed++ > max_lines)
455 {
456 os << "...\n";
457 break;
458 }
459 }
460
461 os << "Mesh component map:\n" << map._mesh_component_map;
462 return os;
463}

References LocalToGlobalIndexMap(), _mesh_component_map, _rows, getNumberOfGlobalComponents(), and size().

Member Data Documentation

◆ _columns

Table const& NumLib::LocalToGlobalIndexMap::_columns = _rows
private
See also
_rows

Definition at line 199 of file LocalToGlobalIndexMap.h.

Referenced by operator()().

◆ _mesh_component_map

◆ _mesh_subsets

std::vector<MeshLib::MeshSubset> NumLib::LocalToGlobalIndexMap::_mesh_subsets
private

◆ _rows

Table NumLib::LocalToGlobalIndexMap::_rows
private

Table contains for each element (first index) and each component (second index) a vector (LineIndex) of indices in the global stiffness matrix or vector

Definition at line 196 of file LocalToGlobalIndexMap.h.

Referenced by LocalToGlobalIndexMap(), findGlobalIndices(), findGlobalIndicesWithElementID(), getElementVariableIDs(), getNumberOfElementComponents(), getNumberOfElementDOF(), operator()(), operator<<, and size().

◆ _variable_component_offsets

std::vector<int> const NumLib::LocalToGlobalIndexMap::_variable_component_offsets
private

The documentation for this class was generated from the following files: