OGS
CollectAndInterpolateNodalDof.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
7
8namespace
9{
16 MeshLib::Element const& element, std::size_t const mesh_id,
17 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
18 int const variable, int const component, unsigned const num_nodes,
19 Eigen::Ref<Eigen::VectorXd> all_nodal_dof_for_this_component)
20{
21 bool dof_not_found = false;
22
23 for (unsigned element_node_id = 0; element_node_id < num_nodes;
24 ++element_node_id)
25 {
26 auto const& node = *element.getNode(element_node_id);
27 auto const node_id = node.getID();
29 node_id};
30 auto const dof_idx = dof_table.getGlobalIndex(loc, variable, component);
31
32 if (dof_idx == NumLib::MeshComponentMap::nop)
33 {
34 // We just skip this d.o.f. Actually we will also skip
35 // all other nodes of this mesh element for this (var,
36 // comp), because we assume that all linear nodes have a
37 // lower node id than any higher order node.
38 dof_not_found = true;
39 }
40 else
41 {
42 if (dof_not_found)
43 {
44 // We expect that for all mesh elements all linear
45 // nodes have a lower node id than any higher order
46 // node. I.e., there are no "holes" in the
47 // primary_variables_mat. We rely on there being no
48 // "holes" later on when interpolating the nodal
49 // d.o.f. to the integration points.
51 "This d.o.f. has been found in the d.o.f. "
52 "table, but before some d.o.f. has not been "
53 "found. Something has gone terribly wrong. "
54 "Some assumption in the implementation is "
55 "wrong.");
56 }
57 all_nodal_dof_for_this_component[element_node_id] = x[dof_idx];
58 }
59 }
60}
61
62} // namespace
63
65{
66
67Eigen::MatrixXd collectDofsToMatrix(
68 MeshLib::Element const& element,
69 std::size_t const mesh_id,
70 NumLib::LocalToGlobalIndexMap const& dof_table,
71 GlobalVector const& x)
72{
73 auto const num_var = dof_table.getNumberOfVariables();
74 auto const num_nodes = element.getNumberOfNodes();
75 auto const num_comp_total = dof_table.getNumberOfGlobalComponents();
76
77 Eigen::MatrixXd primary_variables_mat(num_nodes, num_comp_total);
78
79 for (int var = 0; var < num_var; ++var)
80 {
81 auto const num_comp = dof_table.getNumberOfVariableComponents(var);
82
83 for (int comp = 0; comp < num_comp; ++comp)
84 {
85 auto const global_component =
86 dof_table.getGlobalComponent(var, comp);
87 auto all_nodal_dof_for_this_component =
88 primary_variables_mat.col(global_component);
89
90 collectDofsToMatrixSingleComponentForSomeNodes(
91 element, mesh_id, dof_table, x, var, comp, num_nodes,
92 all_nodal_dof_for_this_component);
93 }
94 }
95
96 return primary_variables_mat;
97}
98
100 MeshLib::Element const& element, std::size_t const mesh_id,
101 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
102 int const variable, int const component)
103{
104 auto const num_nodes = element.getNumberOfBaseNodes();
105 Eigen::VectorXd primary_variables_vec(num_nodes);
106
107 collectDofsToMatrixSingleComponentForSomeNodes(
108 element, mesh_id, dof_table, x, variable, component, num_nodes,
109 primary_variables_vec);
110
111 return primary_variables_vec;
112}
113
114} // namespace ProcessLib::BoundaryConditionAndSourceTerm::Python
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenVector GlobalVector
std::size_t getID() const
virtual unsigned getNumberOfNodes() const =0
virtual unsigned getNumberOfBaseNodes() const =0
virtual const Node * getNode(unsigned idx) const =0
int getGlobalComponent(int const variable_id, int const component_id) const
GlobalIndexType getGlobalIndex(MeshLib::Location const &l, int const variable_id, int const component_id) const
int getNumberOfVariableComponents(int variable_id) const
static constexpr NUMLIB_EXPORT GlobalIndexType const nop
Eigen::MatrixXd collectDofsToMatrix(MeshLib::Element const &element, std::size_t const mesh_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x)
Eigen::VectorXd collectDofsToMatrixOnBaseNodesSingleComponent(MeshLib::Element const &element, std::size_t const mesh_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, int const variable, int const component)
void collectDofsToMatrixSingleComponentForSomeNodes(MeshLib::Element const &element, std::size_t const mesh_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, int const variable, int const component, unsigned const num_nodes, Eigen::Ref< Eigen::VectorXd > all_nodal_dof_for_this_component)