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 40 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 197 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 105 of file LocalToGlobalIndexMap.cpp.

108 : LocalToGlobalIndexMap(std::move(mesh_subsets),
109 std::vector<int>(mesh_subsets.size(), 1), order)
110{
111}
LocalToGlobalIndexMap(std::vector< MeshLib::MeshSubset > &&mesh_subsets, NumLib::ComponentOrder const order)

◆ 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 113 of file LocalToGlobalIndexMap.cpp.

117 : _mesh_subsets(std::move(mesh_subsets)),
120{
121 // For each element of that MeshSubset save a line of global indices.
122 for (int variable_id = 0;
123 variable_id < static_cast<int>(vec_var_n_components.size());
124 ++variable_id)
125 {
126 for (int component_id = 0;
127 component_id < static_cast<int>(vec_var_n_components[variable_id]);
128 ++component_id)
129 {
130 auto const global_component_id =
131 getGlobalComponent(variable_id, component_id);
132
133 auto const& ms = _mesh_subsets[global_component_id];
134 std::size_t const mesh_id = ms.getMeshID();
135
136 findGlobalIndices(ms.elementsBegin(), ms.elementsEnd(),
137 ms.getNodes(), mesh_id, global_component_id,
138 global_component_id);
139 }
140 }
141}
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:283

References _mesh_subsets, 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 143 of file LocalToGlobalIndexMap.cpp.

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

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

References _mesh_subsets, findGlobalIndices(), and OGS_FATAL.

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 224 of file LocalToGlobalIndexMap.cpp.

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

References _mesh_component_map, _mesh_subsets, _variable_component_offsets, DBUG(), getGlobalComponent(), NumLib::MeshComponentMap::getSubset(), 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::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 267 of file LocalToGlobalIndexMap.cpp.

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

References _mesh_component_map, _mesh_subsets, _variable_component_offsets, DBUG(), getNumberOfGlobalComponents(), and NumLib::MeshComponentMap::getSubset().

◆ 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 302 of file LocalToGlobalIndexMap.cpp.

303{
305}
std::size_t dofSizeWithGhosts() const
The number of dofs including the those located in the ghost nodes.

References _mesh_component_map, and NumLib::MeshComponentMap::dofSizeWithGhosts().

◆ 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 307 of file LocalToGlobalIndexMap.cpp.

308{
310}
std::size_t dofSizeWithoutGhosts() const

References _mesh_component_map, and NumLib::MeshComponentMap::dofSizeWithoutGhosts().

Referenced by NumLib::LocalLinearLeastSquaresExtrapolator::extrapolate().

◆ 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 62 of file LocalToGlobalIndexMap.cpp.

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

References _mesh_component_map, _mesh_subsets, _rows, NumLib::MeshComponentMap::getGlobalIndex(), 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 27 of file LocalToGlobalIndexMap.cpp.

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

References _mesh_component_map, _rows, NumLib::MeshComponentMap::getGlobalIndex(), and MeshLib::Node.

Referenced by LocalToGlobalIndexMap().

◆ getElementVariableIDs()

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

Definition at line 368 of file LocalToGlobalIndexMap.cpp.

370{
371 std::vector<int> vec;
372 for (int i = 0; i < getNumberOfVariables(); i++)
373 {
374 for (int j = 0; j < getNumberOfVariableComponents(i); j++)
375 {
376 auto comp_id = getGlobalComponent(i, j);
377 if (!_rows(mesh_item_id, comp_id).empty())
378 {
379 vec.push_back(i);
380 }
381 }
382 }
383 std::sort(vec.begin(), vec.end());
384 vec.erase(std::unique(vec.begin(), vec.end()), vec.end());
385
386 return vec;
387}
int getNumberOfVariableComponents(int variable_id) const

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

Referenced by ProcessLib::LIE::SmallDeformation::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerDataMatrix, LocalAssemblerDataMatrixNearFracture, LocalAssemblerDataFracture, GlobalDim, ConstructorArgs >::operator()().

◆ getGhostIndices()

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

Get ghost indices, forwarded from MeshComponentMap.

Definition at line 412 of file LocalToGlobalIndexMap.cpp.

414{
416}
std::vector< GlobalIndexType > const & getGhostIndices() const
Get ghost indices (for DDC).

References _mesh_component_map, and NumLib::MeshComponentMap::getGhostIndices().

Referenced by NumLib::LocalLinearLeastSquaresExtrapolator::extrapolate().

◆ 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 20 of file LocalToGlobalIndexMap.cpp.

22{
23 return _variable_component_offsets[variable_id] + component_id;
24}

References _variable_component_offsets.

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

◆ getGlobalIndex() [1/2]

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

Definition at line 398 of file LocalToGlobalIndexMap.cpp.

400{
401 return _mesh_component_map.getGlobalIndex(l, global_component_id);
402}

References _mesh_component_map, and NumLib::MeshComponentMap::getGlobalIndex().

◆ getGlobalIndex() [2/2]

◆ getGlobalIndices()

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

Forwards the respective method from MeshComponentMap.

Definition at line 405 of file LocalToGlobalIndexMap.cpp.

407{
409}
std::vector< GlobalIndexType > getGlobalIndices(const Location &l) const

References _mesh_component_map, and NumLib::MeshComponentMap::getGlobalIndices().

◆ 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 420 of file LocalToGlobalIndexMap.cpp.

423{
424 return _mesh_component_map.getLocalIndex(l, comp_id, range_begin,
425 range_end);
426}
GlobalIndexType getLocalIndex(Location const &l, int const comp_id, std::size_t const range_begin, std::size_t const range_end) const

References _mesh_component_map, and NumLib::MeshComponentMap::getLocalIndex().

Referenced by getIndexForComponentInSolutionVector().

◆ getMeshSubset() [1/2]

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

Definition at line 434 of file LocalToGlobalIndexMap.cpp.

436{
437 return _mesh_subsets[global_component_id];
438}

References _mesh_subsets.

◆ getMeshSubset() [2/2]

◆ getNumberOfElementComponents()

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

Definition at line 354 of file LocalToGlobalIndexMap.cpp.

356{
357 std::size_t n = 0;
358 for (Table::Index c = 0; c < _rows.cols(); ++c)
359 {
360 if (!_rows(mesh_item_id, c).empty())
361 {
362 n++;
363 }
364 }
365 return n;
366}

References _rows.

Referenced by ProcessLib::LIE::SmallDeformation::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerDataMatrix, LocalAssemblerDataMatrixNearFracture, LocalAssemblerDataFracture, GlobalDim, ConstructorArgs >::operator()().

◆ getNumberOfElementDOF()

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

◆ 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 334 of file LocalToGlobalIndexMap.cpp.

336{
337 return RowColumnIndices(_rows(mesh_item_id, global_component_id),
338 _columns(mesh_item_id, global_component_id));
339}
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices

References _columns, and _rows.

◆ size()

Friends And Related Symbol Documentation

◆ 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 441 of file LocalToGlobalIndexMap.cpp.

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

Member Data Documentation

◆ _columns

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

Definition at line 206 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 203 of file LocalToGlobalIndexMap.h.

Referenced by LocalToGlobalIndexMap(), findGlobalIndices(), findGlobalIndicesWithElementID(), getElementVariableIDs(), getNumberOfElementComponents(), getNumberOfElementDOF(), 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: