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 74 of file CollectAndInterpolateNodalDof.cpp.

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
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}
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)

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 106 of file CollectAndInterpolateNodalDof.cpp.

110{
111 auto const num_nodes = element.getNumberOfBaseNodes();
112 Eigen::VectorXd primary_variables_vec(num_nodes);
113
115 element, mesh_id, dof_table, x, variable, component, num_nodes,
116 primary_variables_vec);
117
118 return primary_variables_vec;
119}

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 149 of file NsAndWeight.h.

152{
153 using Traits =
155 using VecOfNsAndWeight = std::vector<typename Traits::NsAndWeight>;
156
157 VecOfNsAndWeight nss_and_weights;
158 nss_and_weights.reserve(integration_method.getNumberOfPoints());
159
160 auto sms =
161 NumLib::initShapeMatrices<ShapeFunction,
162 typename Traits::ShapeMatrixPolicy, GlobalDim,
164 element, is_axially_symmetric, integration_method);
165
166 if constexpr (std::is_same_v<ShapeFunction, LowerOrderShapeFunction>)
167 {
168 static_assert(ShapeFunction::ORDER < 2,
169 "We do not expect higher order shape functions here. "
170 "Something must have gone terribly wrong.");
171
172 for (unsigned ip = 0; ip < sms.size(); ++ip)
173 {
174 auto& sm = sms[ip];
175 double const w =
176 sm.detJ * sm.integralMeasure *
177 integration_method.getWeightedPoint(ip).getWeight();
178
179 nss_and_weights.emplace_back(std::move(sm.N), w);
180 }
181 }
182 else
183 {
184 auto sms_lower = NumLib::initShapeMatrices<
185 LowerOrderShapeFunction,
186 typename Traits::LowerOrderShapeMatrixPolicy, GlobalDim,
187 NumLib::ShapeMatrixType::N>(element, is_axially_symmetric,
188 integration_method);
189
190 for (unsigned ip = 0; ip < sms.size(); ++ip)
191 {
192 auto& sm = sms[ip];
193
194 // Note: we use det(J) of the higher order shape function. For
195 // linear geometries (i.e., higher order elements but with flat
196 // surfaces) there should be no difference.
197 double const w =
198 sm.detJ * sm.integralMeasure *
199 integration_method.getWeightedPoint(ip).getWeight();
200
201 nss_and_weights.emplace_back(std::move(sm.N),
202 std::move(sms_lower[ip].N), w);
203 }
204 }
205
206 return nss_and_weights;
207}
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.

◆ 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 63 of file CollectAndInterpolateNodalDof.h.

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}

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

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