OGS
ProcessLib::BoundaryConditionAndSourceTerm::Python Namespace Reference

Classes

struct  BcAndStLocalAssemblerImpl
struct  BcOrStData
struct  FlagAndFluxAndDFlux
struct  NsAndWeight
struct  NsAndWeight< ShapeFunction, ShapeFunction, ShapeMatrix, ShapeMatrix >
struct  NsAndWeightsTraits
 Collects common type aliases needed when working with NsAndWeight. More...

Functions

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)
template<typename NsAndWeight>
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)
template<typename ShapeFunction, typename LowerOrderShapeFunction, int GlobalDim, typename IntegrationMethod>
auto computeNsAndWeights (MeshLib::Element const &element, bool const is_axially_symmetric, IntegrationMethod const &integration_method)

Function Documentation

◆ collectDofsToMatrix()

Eigen::MatrixXd ProcessLib::BoundaryConditionAndSourceTerm::Python::collectDofsToMatrix ( MeshLib::Element const & element,
std::size_t const mesh_id,
NumLib::LocalToGlobalIndexMap const & dof_table,
GlobalVector const & x )

Collects the degrees of freedom of the passed element from the passed global vector into a matrix.

The dimensions of the returned matrix are #nodes x #total_components, i.e., each row of the matrix will contain all degrees of freedom at a specific node of the passed element.

The order of nodes is determined by the passed element. The order of components is determined by the passed dof_table.

Note
OGS currently has implemented Lagrange finite elements only, i.e., all degrees of freedom are nodal degrees of freedom.
In the case of Taylor-Hood elements, some components are defined only on base nodes. In that case, the returned matrix might contain some uninitialized entries. It is the responsibility of the caller to handle the returned matrix correctly.

Definition at line 67 of file CollectAndInterpolateNodalDof.cpp.

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}

References NumLib::LocalToGlobalIndexMap::getGlobalComponent(), NumLib::LocalToGlobalIndexMap::getNumberOfGlobalComponents(), MeshLib::Element::getNumberOfNodes(), NumLib::LocalToGlobalIndexMap::getNumberOfVariableComponents(), and NumLib::LocalToGlobalIndexMap::getNumberOfVariables().

Referenced by ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assemble().

◆ collectDofsToMatrixOnBaseNodesSingleComponent()

Eigen::VectorXd ProcessLib::BoundaryConditionAndSourceTerm::Python::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 )

Does the same as collectDofsToMatrix(), just for a single component and only on the base nodes of the passed mesh element.

Definition at line 99 of file CollectAndInterpolateNodalDof.cpp.

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}

References MeshLib::Element::getNumberOfBaseNodes().

Referenced by ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::interpolate().

◆ computeNsAndWeights()

template<typename ShapeFunction, typename LowerOrderShapeFunction, int GlobalDim, typename IntegrationMethod>
auto ProcessLib::BoundaryConditionAndSourceTerm::Python::computeNsAndWeights ( MeshLib::Element const & element,
bool const is_axially_symmetric,
IntegrationMethod const & integration_method )

Computes shape matrices and integration weights for all integration points in the given mesh element.

Note
This is an extension of GenericNaturalBoundaryConditionLocalAssembler::initNsAndWeights().

Definition at line 142 of file NsAndWeight.h.

145{
146 using Traits =
148 using VecOfNsAndWeight = std::vector<typename Traits::NsAndWeight>;
149
150 VecOfNsAndWeight nss_and_weights;
151 nss_and_weights.reserve(integration_method.getNumberOfPoints());
152
153 auto sms =
155 typename Traits::ShapeMatrixPolicy, GlobalDim,
157 element, is_axially_symmetric, integration_method);
158
159 if constexpr (std::is_same_v<ShapeFunction, LowerOrderShapeFunction>)
160 {
161 static_assert(ShapeFunction::ORDER < 2,
162 "We do not expect higher order shape functions here. "
163 "Something must have gone terribly wrong.");
164
165 for (unsigned ip = 0; ip < sms.size(); ++ip)
166 {
167 auto& sm = sms[ip];
168 double const w =
169 sm.detJ * sm.integralMeasure *
170 integration_method.getWeightedPoint(ip).getWeight();
171
172 nss_and_weights.emplace_back(std::move(sm.N), w);
173 }
174 }
175 else
176 {
177 auto sms_lower = NumLib::initShapeMatrices<
178 LowerOrderShapeFunction,
179 typename Traits::LowerOrderShapeMatrixPolicy, GlobalDim,
180 NumLib::ShapeMatrixType::N>(element, is_axially_symmetric,
181 integration_method);
182
183 for (unsigned ip = 0; ip < sms.size(); ++ip)
184 {
185 auto& sm = sms[ip];
186
187 // Note: we use det(J) of the higher order shape function. For
188 // linear geometries (i.e., higher order elements but with flat
189 // surfaces) there should be no difference.
190 double const w =
191 sm.detJ * sm.integralMeasure *
192 integration_method.getWeightedPoint(ip).getWeight();
193
194 nss_and_weights.emplace_back(std::move(sm.N),
195 std::move(sms_lower[ip].N), w);
196 }
197 }
198
199 return nss_and_weights;
200}
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
@ N_J
calculates N, dNdr, J, and detJ
Collects common type aliases needed when working with NsAndWeight.

References NumLib::initShapeMatrices(), NumLib::N, and NumLib::N_J.

Referenced by ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::BcAndStLocalAssemblerImpl().

◆ interpolate()

template<typename NsAndWeight>
void ProcessLib::BoundaryConditionAndSourceTerm::Python::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 )

Interpolates the passed matrix of degrees of freedom to the interpolated_primary_variables via the passed shape matrices in ns_and_weight.

The interpolation function order for each variable is determined by the passed ProcessVariable vector.

For the layout of primary_variables_mat see collectDofsToMatrix().

Definition at line 56 of file CollectAndInterpolateNodalDof.h.

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}

References ProcessLib::BoundaryConditionAndSourceTerm::Python::NsAndWeight< ShapeFunction, LowerOrderShapeFunction, ShapeMatrix, LowerOrderShapeMatrix >::N().

Referenced by ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assemble().