OGS
CollectAndInterpolateNodalDof.cpp
Go to the documentation of this file.
1
12
14
15namespace
16{
23 MeshLib::Element const& element, std::size_t const mesh_id,
24 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
25 int const variable, int const component, unsigned const num_nodes,
26 Eigen::Ref<Eigen::VectorXd> all_nodal_dof_for_this_component)
27{
28 bool dof_not_found = false;
29
30 for (unsigned element_node_id = 0; element_node_id < num_nodes;
31 ++element_node_id)
32 {
33 auto const& node = *element.getNode(element_node_id);
34 auto const node_id = node.getID();
36 node_id};
37 auto const dof_idx = dof_table.getGlobalIndex(loc, variable, component);
38
39 if (dof_idx == NumLib::MeshComponentMap::nop)
40 {
41 // We just skip this d.o.f. Actually we will also skip
42 // all other nodes of this mesh element for this (var,
43 // comp), because we assume that all linear nodes have a
44 // lower node id than any higher order node.
45 dof_not_found = true;
46 }
47 else
48 {
49 if (dof_not_found)
50 {
51 // We expect that for all mesh elements all linear
52 // nodes have a lower node id than any higher order
53 // node. I.e., there are no "holes" in the
54 // primary_variables_mat. We rely on there being no
55 // "holes" later on when interpolating the nodal
56 // d.o.f. to the integration points.
58 "This d.o.f. has been found in the d.o.f. "
59 "table, but before some d.o.f. has not been "
60 "found. Something has gone terribly wrong. "
61 "Some assumption in the implementation is "
62 "wrong.");
63 }
64 all_nodal_dof_for_this_component[element_node_id] = x[dof_idx];
65 }
66 }
67}
68
69} // namespace
70
72{
73
74Eigen::MatrixXd collectDofsToMatrix(
75 MeshLib::Element const& element,
76 std::size_t const mesh_id,
77 NumLib::LocalToGlobalIndexMap const& dof_table,
78 GlobalVector const& x)
79{
80 auto const num_var = dof_table.getNumberOfVariables();
81 auto const num_nodes = element.getNumberOfNodes();
82 auto const num_comp_total = dof_table.getNumberOfGlobalComponents();
83
84 Eigen::MatrixXd primary_variables_mat(num_nodes, num_comp_total);
85
86 for (int var = 0; var < num_var; ++var)
87 {
88 auto const num_comp = dof_table.getNumberOfVariableComponents(var);
89
90 for (int comp = 0; comp < num_comp; ++comp)
91 {
92 auto const global_component =
93 dof_table.getGlobalComponent(var, comp);
94 auto all_nodal_dof_for_this_component =
95 primary_variables_mat.col(global_component);
96
97 collectDofsToMatrixSingleComponentForSomeNodes(
98 element, mesh_id, dof_table, x, var, comp, num_nodes,
99 all_nodal_dof_for_this_component);
100 }
101 }
102
103 return primary_variables_mat;
104}
105
107 MeshLib::Element const& element, std::size_t const mesh_id,
108 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
109 int const variable, int const component)
110{
111 auto const num_nodes = element.getNumberOfBaseNodes();
112 Eigen::VectorXd primary_variables_vec(num_nodes);
113
114 collectDofsToMatrixSingleComponentForSomeNodes(
115 element, mesh_id, dof_table, x, variable, component, num_nodes,
116 primary_variables_vec);
117
118 return primary_variables_vec;
119}
120
121} // namespace ProcessLib::BoundaryConditionAndSourceTerm::Python
#define OGS_FATAL(...)
Definition Error.h:26
Global vector based on Eigen vector.
Definition EigenVector.h:25
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)