OGS
CollectAndInterpolateNodalDof.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
8#include <Eigen/Core>
9
10namespace MeshLib
11{
12class Element;
13} // namespace MeshLib
14
16{
34Eigen::MatrixXd collectDofsToMatrix(
35 MeshLib::Element const& element,
36 std::size_t const mesh_id,
37 NumLib::LocalToGlobalIndexMap const& dof_table,
38 GlobalVector const& x);
39
43 MeshLib::Element const& element, std::size_t const mesh_id,
44 NumLib::LocalToGlobalIndexMap const& dof_table, GlobalVector const& x,
45 int const variable, int const component);
46
55template <typename NsAndWeight>
57 Eigen::MatrixXd const& primary_variables_mat,
58 std::vector<std::reference_wrapper<ProcessVariable>> const& pv_refs,
59 NsAndWeight const& ns_and_weight,
60 Eigen::Ref<Eigen::VectorXd>
61 interpolated_primary_variables)
62{
63 Eigen::Index component_flattened = 0;
64
65 // We assume that all_process_variables_for_this_process have the same
66 // order as the d.o.f. table. Therefore we can iterate over
67 // all_process_variables_for_this_process.
68 for (auto pv_ref : pv_refs)
69 {
70 auto const& pv = pv_ref.get();
71 auto const num_comp = pv.getNumberOfGlobalComponents();
72 auto const shp_fct_order = pv.getShapeFunctionOrder();
73 auto const N = ns_and_weight.N(shp_fct_order);
74
75 for (auto comp = decltype(num_comp){0}; comp < num_comp; ++comp)
76 {
77 // This computation assumes that there are no "holes" in the
78 // primary_variables_mat. I.e., all nodal d.o.f. for a certain
79 // (var, comp) must be stored contiguously in the respective column
80 // of primary_variables_mat.
81 interpolated_primary_variables[component_flattened] =
82 N *
83 primary_variables_mat.col(component_flattened).head(N.size());
84 component_flattened++;
85 }
86 }
87}
88} // namespace ProcessLib::BoundaryConditionAndSourceTerm::Python
MathLib::EigenVector GlobalVector
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:44