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)
 
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 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. More...
 
std::vector< GlobalIndexType > const & getGhostIndices () const
 Get ghost indices, forwarded from MeshComponentMap. More...
 
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 = Eigen::Matrix< LineIndex, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >
 

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. More...
 
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. More...
 

Member Typedef Documentation

◆ LineIndex

◆ RowColumnIndices

◆ Table

using NumLib::LocalToGlobalIndexMap::Table = Eigen::Matrix<LineIndex, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
private

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

121  : LocalToGlobalIndexMap(std::move(mesh_subsets),
122  std::vector<int>(mesh_subsets.size(), 1), order)
123 {
124 }
LocalToGlobalIndexMap(std::vector< MeshLib::MeshSubset > &&mesh_subsets, NumLib::ComponentOrder const order)

Referenced by deriveBoundaryConstrainedMap().

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

130  : _mesh_subsets(std::move(mesh_subsets)),
132  _variable_component_offsets(to_cumulative(vec_var_n_components))
133 {
134  // For each element of that MeshSubset save a line of global indices.
135  for (int variable_id = 0;
136  variable_id < static_cast<int>(vec_var_n_components.size());
137  ++variable_id)
138  {
139  for (int component_id = 0;
140  component_id < static_cast<int>(vec_var_n_components[variable_id]);
141  ++component_id)
142  {
143  auto const global_component_id =
144  getGlobalComponent(variable_id, component_id);
145 
146  auto const& ms = _mesh_subsets[global_component_id];
147  std::size_t const mesh_id = ms.getMeshID();
148 
149  findGlobalIndices(ms.elementsBegin(), ms.elementsEnd(),
150  ms.getNodes(), mesh_id, global_component_id,
151  global_component_id);
152  }
153  }
154 }
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< T > to_cumulative(std::vector< T > const &vec)

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

161  : _mesh_subsets(std::move(mesh_subsets)),
163  _variable_component_offsets(to_cumulative(vec_var_n_components))
164 {
165  assert(vec_var_n_components.size() == vec_var_elements.size());
166 
167  // For each element of that MeshSubset save a line of global indices.
168 
169  // _rows should be resized based on an element ID
170  std::size_t max_elem_id = 0;
171  for (std::vector<MeshLib::Element*> const* elements : vec_var_elements)
172  {
173  for (auto e : *elements)
174  {
175  max_elem_id = std::max(max_elem_id, e->getID());
176  }
177  }
178  _rows.resize(max_elem_id + 1, _mesh_subsets.size());
179 
180  for (int variable_id = 0;
181  variable_id < static_cast<int>(vec_var_n_components.size());
182  ++variable_id)
183  {
184  std::vector<MeshLib::Element*> const& var_elements =
185  *vec_var_elements[variable_id];
186  for (int component_id = 0;
187  component_id < static_cast<int>(vec_var_n_components[variable_id]);
188  ++component_id)
189  {
190  auto const global_component_id =
191  getGlobalComponent(variable_id, component_id);
192 
193  auto const& ms = _mesh_subsets[global_component_id];
194  std::size_t const mesh_id = ms.getMeshID();
195 
197  var_elements.cbegin(), var_elements.cend(), ms.getNodes(),
198  mesh_id, global_component_id, global_component_id);
199  }
200  }
201 }
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 203 of file LocalToGlobalIndexMap.cpp.

210  : _mesh_subsets(std::move(mesh_subsets)),
211  _mesh_component_map(std::move(mesh_component_map)),
212  _variable_component_offsets(variable_component_offsets)
213 {
214  // Each subset in the mesh_subsets represents a single component.
215  if (_mesh_subsets.size() != global_component_ids.size())
216  {
217  OGS_FATAL(
218  "Number of mesh subsets is not equal to number of components. "
219  "There are {:d} mesh subsets and {:d} components.",
220  _mesh_subsets.size(), global_component_ids.size());
221  }
222 
223  for (int i = 0; i < static_cast<int>(global_component_ids.size()); ++i)
224  {
225  auto const& ms = _mesh_subsets[i];
226 
227  // For all MeshSubset in mesh_subsets and each element of that
228  // MeshSubset save a line of global indices.
229  std::size_t const mesh_id = ms.getMeshID();
230 
231  findGlobalIndices(elements.cbegin(), elements.cend(), ms.getNodes(),
232  mesh_id, global_component_ids[i], i);
233  }
234 }
#define OGS_FATAL(...)
Definition: Error.h:26

References _mesh_subsets, findGlobalIndices(), and OGS_FATAL.

Member Function Documentation

◆ deriveBoundaryConstrainedMap() [1/2]

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

240 {
241  DBUG("Construct reduced local to global index map.");
242 
243  if (component_ids.empty())
244  {
245  OGS_FATAL("Expected non-empty vector of component ids.");
246  }
247 
248  // Elements of the new_mesh_subset's mesh.
249  std::vector<MeshLib::Element*> const& elements =
250  new_mesh_subset.getMesh().getElements();
251 
252  // Create a subset of the current mesh component map.
253  std::vector<int> global_component_ids;
254 
255  transform(cbegin(component_ids), cend(component_ids),
256  back_inserter(global_component_ids),
257  [&](auto const component_id)
258  { return getGlobalComponent(variable_id, component_id); });
259 
260  auto mesh_component_map = _mesh_component_map.getSubset(
261  _mesh_subsets, new_mesh_subset, global_component_ids);
262 
263  // Create copies of the new_mesh_subset for each of the global components.
264  // The last component is moved after the for-loop.
265  std::vector<MeshLib::MeshSubset> all_mesh_subsets;
266  for (int i = 0; i < static_cast<int>(global_component_ids.size()) - 1; ++i)
267  {
268  all_mesh_subsets.emplace_back(new_mesh_subset);
269  }
270  all_mesh_subsets.emplace_back(std::move(new_mesh_subset));
271 
272  return new LocalToGlobalIndexMap(
273  std::move(all_mesh_subsets), global_component_ids,
274  _variable_component_offsets, elements, std::move(mesh_component_map),
275  ConstructorTag{});
276 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Mesh const & getMesh() const
Definition: MeshSubset.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
MeshComponentMap getSubset(std::vector< MeshLib::MeshSubset > const &bulk_mesh_subsets, MeshLib::MeshSubset const &new_mesh_subset, std::vector< int > const &new_global_component_ids) const

References LocalToGlobalIndexMap(), _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::GenericNaturalBoundaryCondition< BoundaryConditionData, LocalAssemblerImplementation >::GenericNaturalBoundaryCondition(), ProcessLib::NormalTractionBoundaryCondition::NormalTractionBoundaryCondition< GlobalDim, LocalAssemblerImplementation >::NormalTractionBoundaryCondition(), ProcessLib::PrimaryVariableConstraintDirichletBoundaryCondition::PrimaryVariableConstraintDirichletBoundaryCondition(), ProcessLib::PythonBoundaryCondition::PythonBoundaryCondition(), ProcessLib::SolutionDependentDirichletBoundaryCondition::SolutionDependentDirichletBoundaryCondition(), ProcessLib::DeactivatedSubdomainDirichlet::config(), ProcessLib::DirichletBoundaryConditionWithinTimeInterval::config(), ProcessLib::createSourceTerm(), ProcessLib::createVariableDependentNeumannBoundaryCondition(), and ProcessLib::Output::doOutputAlways().

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

281 {
282  DBUG("Construct reduced local to global index map.");
283 
284  // Create a subset of the current mesh component map.
285  std::vector<int> global_component_ids;
286 
287  for (int i = 0; i < getNumberOfGlobalComponents(); ++i)
288  {
289  global_component_ids.push_back(i);
290  }
291 
292  // Elements of the new_mesh_subset's mesh.
293  std::vector<MeshLib::Element*> const& elements =
294  new_mesh_subset.getMesh().getElements();
295 
296  auto mesh_component_map = _mesh_component_map.getSubset(
297  _mesh_subsets, new_mesh_subset, global_component_ids);
298 
299  // Create copies of the new_mesh_subset for each of the global components.
300  // The last component is moved after the for-loop.
301  std::vector<MeshLib::MeshSubset> all_mesh_subsets;
302  for (int i = 0; i < static_cast<int>(global_component_ids.size()) - 1; ++i)
303  {
304  all_mesh_subsets.emplace_back(new_mesh_subset);
305  }
306  all_mesh_subsets.emplace_back(std::move(new_mesh_subset));
307 
308  return std::make_unique<LocalToGlobalIndexMap>(
309  std::move(all_mesh_subsets), global_component_ids,
310  _variable_component_offsets, elements, std::move(mesh_component_map),
311  ConstructorTag{});
312 }

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

315 {
317 }
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 319 of file LocalToGlobalIndexMap.cpp.

320 {
322 }
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 75 of file LocalToGlobalIndexMap.cpp.

79 {
80  _rows.resize(std::distance(first, last), _mesh_subsets.size());
81 
82  std::unordered_set<MeshLib::Node*> const set_nodes(nodes.begin(),
83  nodes.end());
84 
85  // For each element find the global indices for node/element
86  // components.
87  std::size_t elem_id = 0;
88  for (ElementIterator e = first; e != last; ++e, ++elem_id)
89  {
90  LineIndex indices;
91  indices.reserve((*e)->getNumberOfNodes());
92 
93  for (auto* n = (*e)->getNodes();
94  n < (*e)->getNodes() + (*e)->getNumberOfNodes();
95  ++n)
96  {
97  // Check if the element's node is in the given list of nodes.
98  if (set_nodes.find(*n) == set_nodes.end())
99  {
100  continue;
101  }
103  (*n)->getID());
104  auto const global_index =
106  if (global_index == std::numeric_limits<GlobalIndexType>::max())
107  {
108  continue;
109  }
110  indices.push_back(global_index);
111  }
112 
113  indices.shrink_to_fit();
114  _rows(elem_id, comp_id_write) = std::move(indices);
115  }
116 }
RowColumnIndices::LineIndex LineIndex
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().

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

44 {
45  std::unordered_set<MeshLib::Node*> const set_nodes(nodes.begin(),
46  nodes.end());
47 
48  // For each element find the global indices for node/element
49  // components.
50  for (ElementIterator e = first; e != last; ++e)
51  {
52  LineIndex indices;
53  indices.reserve((*e)->getNumberOfNodes());
54 
55  for (auto* n = (*e)->getNodes();
56  n < (*e)->getNodes() + (*e)->getNumberOfNodes();
57  ++n)
58  {
59  // Check if the element's node is in the given list of nodes.
60  if (set_nodes.find(*n) == set_nodes.end())
61  {
62  continue;
63  }
65  (*n)->getID());
66  indices.push_back(_mesh_component_map.getGlobalIndex(l, comp_id));
67  }
68 
69  indices.shrink_to_fit();
70  _rows((*e)->getID(), comp_id_write) = std::move(indices);
71  }
72 }

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

382 {
383  std::vector<int> vec;
384  for (int i = 0; i < getNumberOfVariables(); i++)
385  {
386  for (int j = 0; j < getNumberOfVariableComponents(i); j++)
387  {
388  auto comp_id = getGlobalComponent(i, j);
389  if (!_rows(mesh_item_id, comp_id).empty())
390  {
391  vec.push_back(i);
392  }
393  }
394  }
395  std::sort(vec.begin(), vec.end());
396  vec.erase(std::unique(vec.begin(), vec.end()), vec.end());
397 
398  return vec;
399 }
int getNumberOfVariableComponents(int variable_id) const

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

Referenced by ProcessLib::LIE::HydroMechanics::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerDataMatrix, LocalAssemblerDataMatrixNearFracture, LocalAssemblerDataFracture, GlobalDim, ConstructorArgs >::operator()(), and 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 424 of file LocalToGlobalIndexMap.cpp.

426 {
428 }
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

◆ getGlobalIndex() [1/2]

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

Definition at line 410 of file LocalToGlobalIndexMap.cpp.

412 {
413  return _mesh_component_map.getGlobalIndex(l, global_component_id);
414 }

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

◆ getGlobalIndex() [2/2]

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

◆ getGlobalIndices()

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

Forwards the respective method from MeshComponentMap.

Definition at line 417 of file LocalToGlobalIndexMap.cpp.

419 {
421 }
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 432 of file LocalToGlobalIndexMap.cpp.

435 {
436  return _mesh_component_map.getLocalIndex(l, comp_id, range_begin,
437  range_end);
438 }
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().

◆ getMeshSubset() [1/2]

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

Definition at line 446 of file LocalToGlobalIndexMap.cpp.

448 {
449  return _mesh_subsets[global_component_id];
450 }

References _mesh_subsets.

◆ getMeshSubset() [2/2]

◆ getNumberOfElementComponents()

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

Definition at line 366 of file LocalToGlobalIndexMap.cpp.

368 {
369  std::size_t n = 0;
370  for (Table::Index c = 0; c < _rows.cols(); ++c)
371  {
372  if (!_rows(mesh_item_id, c).empty())
373  {
374  n++;
375  }
376  }
377  return n;
378 }

References _rows, and MaterialPropertyLib::c.

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

Definition at line 353 of file LocalToGlobalIndexMap.cpp.

355 {
356  std::size_t ndof = 0;
357 
358  for (Table::Index c = 0; c < _rows.cols(); ++c)
359  {
360  ndof += _rows(mesh_item_id, c).size();
361  }
362 
363  return ndof;
364 }

References _rows, and MaterialPropertyLib::c.

Referenced by ProcessLib::BoundaryConditionAndSourceTerm::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerData, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::HydroMechanics::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerData, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::LIE::HydroMechanics::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerDataMatrix, LocalAssemblerDataMatrixNearFracture, LocalAssemblerDataFracture, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::LIE::SmallDeformation::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerDataMatrix, LocalAssemblerDataMatrixNearFracture, LocalAssemblerDataFracture, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::RichardsMechanics::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerData, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::LocalDataInitializer< LocalAssemblerInterface, SmallDeformationLocalAssembler, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::StokesFlow::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerData, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::TH2M::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerData, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::ThermoHydroMechanics::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerData, GlobalDim, ConstructorArgs >::operator()(), and ProcessLib::ThermoRichardsMechanics::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerData, GlobalDim, ConstructorArgs >::operator()().

◆ getNumberOfGlobalComponents()

◆ getNumberOfVariableComponents()

int NumLib::LocalToGlobalIndexMap::getNumberOfVariableComponents ( int  variable_id) const

Definition at line 329 of file LocalToGlobalIndexMap.cpp.

330 {
331  assert(variable_id < getNumberOfVariables());
332  return _variable_component_offsets[variable_id + 1] -
333  _variable_component_offsets[variable_id];
334 }

References _variable_component_offsets, and getNumberOfVariables().

Referenced by ProcessLib::ConstraintDirichletBoundaryCondition::ConstraintDirichletBoundaryCondition(), ProcessLib::GenericNaturalBoundaryCondition< BoundaryConditionData, LocalAssemblerImplementation >::GenericNaturalBoundaryCondition(), ProcessLib::NormalTractionBoundaryCondition::NormalTractionBoundaryCondition< GlobalDim, LocalAssemblerImplementation >::NormalTractionBoundaryCondition(), ProcessLib::PhaseFieldIrreversibleDamageOracleBoundaryCondition::PhaseFieldIrreversibleDamageOracleBoundaryCondition(), ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, IntegrationMethod, GlobalDim >::assemble(), ProcessLib::SourceTerms::Python::PythonSourceTermLocalAssembler< ShapeFunction, IntegrationMethod, GlobalDim >::assemble(), ProcessLib::checkParametersOfDirichletBoundaryCondition(), ProcessLib::ProcessVariable::createBoundaryConditionsForDeactivatedSubDomains(), ProcessLib::createPythonBoundaryCondition(), ProcessLib::createSourceTerm(), getElementVariableIDs(), ProcessLib::LIE::HydroMechanics::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerDataMatrix, LocalAssemblerDataMatrixNearFracture, LocalAssemblerDataFracture, GlobalDim, ConstructorArgs >::operator()(), ProcessLib::LIE::SmallDeformation::LocalDataInitializer< LocalAssemblerInterface, LocalAssemblerDataMatrix, LocalAssemblerDataMatrixNearFracture, LocalAssemblerDataFracture, GlobalDim, ConstructorArgs >::operator()(), and NumLib::transformVariableFromGlobalVector().

◆ getNumberOfVariables()

◆ operator()()

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

Definition at line 346 of file LocalToGlobalIndexMap.cpp.

348 {
349  return RowColumnIndices(_rows(mesh_item_id, component_id),
350  _columns(mesh_item_id, component_id));
351 }
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices

References _columns, and _rows.

◆ size()

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

Definition at line 341 of file LocalToGlobalIndexMap.cpp.

342 {
343  return _rows.rows();
344 }

References _rows.

Referenced by NumLib::LocalLinearLeastSquaresExtrapolator::calculateResiduals(), NumLib::getIndices(), and NumLib::getRowColumnIndices().

Friends And Related Function 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 453 of file LocalToGlobalIndexMap.cpp.

454 {
455  std::size_t const max_lines = 10;
456  std::size_t lines_printed = 0;
457 
458  os << "Rows of the local to global index map; " << map._rows.size()
459  << " rows\n";
460  for (std::size_t e = 0; e < map.size(); ++e)
461  {
462  os << "== e " << e << " ==\n";
463  for (int c = 0; c < map.getNumberOfGlobalComponents(); ++c)
464  {
465  auto const& line = map._rows(e, c);
466 
467  os << "c" << c << " { ";
468  std::copy(line.cbegin(), line.cend(),
469  std::ostream_iterator<std::size_t>(os, " "));
470  os << " }\n";
471  }
472 
473  if (lines_printed++ > max_lines)
474  {
475  os << "...\n";
476  break;
477  }
478  }
479 
480  os << "Mesh component map:\n" << map._mesh_component_map;
481  return os;
482 }
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37

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

A vector of mesh subsets for each process variables' components.

Definition at line 194 of file LocalToGlobalIndexMap.h.

Referenced by LocalToGlobalIndexMap(), deriveBoundaryConstrainedMap(), findGlobalIndices(), getMeshSubset(), and getNumberOfGlobalComponents().

◆ _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: