OGS
CollectAndInterpolateNodalDof.h
Go to the documentation of this file.
1
11#pragma once
12
15#include <Eigen/Core>
16
17namespace MeshLib
18{
19class Element;
20} // namespace MeshLib
21
23{
41Eigen::MatrixXd collectDofsToMatrix(
42 MeshLib::Element const& element,
43 std::size_t const mesh_id,
44 NumLib::LocalToGlobalIndexMap const& dof_table,
45 GlobalVector const& x);
46
50 MeshLib::Element const& element, std::size_t const mesh_id,
51 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
52 int const variable, int const component);
53
62template <typename NsAndWeight>
64 Eigen::MatrixXd const& primary_variables_mat,
65 std::vector<std::reference_wrapper<ProcessVariable>> const& pv_refs,
66 NsAndWeight const& ns_and_weight,
67 Eigen::Ref<Eigen::VectorXd>
68 interpolated_primary_variables)
69{
70 Eigen::Index component_flattened = 0;
71
72 // We assume that all_process_variables_for_this_process have the same
73 // order as the d.o.f. table. Therefore we can iterate over
74 // all_process_variables_for_this_process.
75 for (auto pv_ref : pv_refs)
76 {
77 auto const& pv = pv_ref.get();
78 auto const num_comp = pv.getNumberOfGlobalComponents();
79 auto const shp_fct_order = pv.getShapeFunctionOrder();
80 auto const N = ns_and_weight.N(shp_fct_order);
81
82 for (auto comp = decltype(num_comp){0}; comp < num_comp; ++comp)
83 {
84 // This computation assumes that there are no "holes" in the
85 // primary_variables_mat. I.e., all nodal d.o.f. for a certain
86 // (var, comp) must be stored contiguously in the respective column
87 // of primary_variables_mat.
88 interpolated_primary_variables[component_flattened] =
89 N *
90 primary_variables_mat.col(component_flattened).head(N.size());
91 component_flattened++;
92 }
93 }
94}
95} // namespace ProcessLib::BoundaryConditionAndSourceTerm::Python
Global vector based on Eigen vector.
Definition EigenVector.h:25
Eigen::MatrixXd collectDofsToMatrix(MeshLib::Element const &element, std::size_t const mesh_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x)
void interpolate(Eigen::MatrixXd const &primary_variables_mat, std::vector< std::reference_wrapper< ProcessVariable > > const &pv_refs, NsAndWeight const &ns_and_weight, Eigen::Ref< Eigen::VectorXd > interpolated_primary_variables)
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)
Eigen::Ref< const Eigen::RowVectorXd > N(unsigned order) const
Definition NsAndWeight.h:51