OGS
ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim > Struct Template Reference

Detailed Description

template<typename BcOrStData, typename ShapeFunction, typename LowerOrderShapeFunction, int GlobalDim>
struct ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >

Common parts of the implementation of Python boundary conditions and source terms.

The implementation of Python boundary conditions and source terms are almost the same, however, their interfaces differ slightly. This struct contains all common parts of the local assemblers.

Definition at line 34 of file BcAndStLocalAssemblerImpl.h.

#include <BcAndStLocalAssemblerImpl.h>

Collaboration diagram for ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >:
[legend]

Public Types

using Traits
 

Public Member Functions

 BcAndStLocalAssemblerImpl (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method_, bool is_axially_symmetric, BcOrStData const &data)
 
void assemble (std::size_t const boundary_element_id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const t, GlobalVector const &x, GlobalVector &b, GlobalMatrix *const Jac) const
 

Public Attributes

BcOrStData const & bc_or_st_data
 
MeshLib::Element const & element
 
NumLib::GenericIntegrationMethod const & integration_method
 
std::vector< typename Traits::NsAndWeight > const nss_and_weights
 

Private Types

using ShapeMatrixPolicy = typename Traits::ShapeMatrixPolicy
 

Private Member Functions

std::array< double, 3 > interpolateCoords (typename Traits::NsAndWeight const &ns_and_weight, NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatrixPolicy > const &fe) const
 
void assembleLocalRhs (double const flux, typename Traits::NsAndWeight const &ns_and_weight, Eigen::VectorXd &local_rhs) const
 
void assembleLocalJacobian (std::vector< double > const &dFlux, typename Traits::NsAndWeight const &ns_and_weight, Eigen::MatrixXd &local_Jac) const
 

Member Typedef Documentation

◆ ShapeMatrixPolicy

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
using ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::ShapeMatrixPolicy = typename Traits::ShapeMatrixPolicy
private

Definition at line 41 of file BcAndStLocalAssemblerImpl.h.

◆ Traits

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
using ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::Traits
Initial value:
NsAndWeightsTraits<ShapeFunction, LowerOrderShapeFunction, GlobalDim>

Definition at line 37 of file BcAndStLocalAssemblerImpl.h.

Constructor & Destructor Documentation

◆ BcAndStLocalAssemblerImpl()

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::BcAndStLocalAssemblerImpl ( MeshLib::Element const & e,
NumLib::GenericIntegrationMethod const & integration_method_,
bool is_axially_symmetric,
BcOrStData const & data )
inline

Definition at line 44 of file BcAndStLocalAssemblerImpl.h.

49 : bc_or_st_data(data),
50 element(e),
51 integration_method(integration_method_),
53 computeNsAndWeights<ShapeFunction, LowerOrderShapeFunction,
54 GlobalDim>(e, is_axially_symmetric,
56 {
57 }
auto computeNsAndWeights(MeshLib::Element const &element, bool const is_axially_symmetric, IntegrationMethod const &integration_method)

Member Function Documentation

◆ assemble()

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
void ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assemble ( std::size_t const boundary_element_id,
NumLib::LocalToGlobalIndexMap const & dof_table_boundary,
double const t,
GlobalVector const & x,
GlobalVector & b,
GlobalMatrix *const Jac ) const
inline

Assembles the local rhs and (possibly) Jacobian for this BC or ST and adds them to the global rhs vector and Jacobian matrix, respectively.

Definition at line 61 of file BcAndStLocalAssemblerImpl.h.

65 {
66 // Variables for input data
67 Eigen::MatrixXd const primary_variables_mat = ProcessLib::
70 dof_table_boundary, x);
71
72 // Variables for output data
73 auto const& indices_this_component =
74 dof_table_boundary(boundary_element_id,
76 .rows;
77 std::vector<GlobalIndexType> indices_all_components_for_Jac;
78
79 auto const num_dof_this_component = indices_this_component.size();
80
81 Eigen::VectorXd local_rhs =
82 Eigen::VectorXd::Zero(num_dof_this_component);
83 Eigen::MatrixXd local_Jac;
84
85 if (Jac)
86 {
87 indices_all_components_for_Jac =
88 NumLib::getIndices(boundary_element_id, dof_table_boundary);
89
90 auto const num_dof_all_components =
91 indices_all_components_for_Jac.size();
92
93 local_Jac = Eigen::MatrixXd::Zero(num_dof_this_component,
94 num_dof_all_components);
95 }
96
97 // Variables for temporary data
98 auto const num_comp_total =
99 dof_table_boundary.getNumberOfGlobalComponents();
100 std::vector<double> prim_vars_data(num_comp_total);
101 auto prim_vars = MathLib::toVector(prim_vars_data);
102
103 // Variables for the integration loop
104 unsigned const num_integration_points =
107 ShapeFunction, ShapeMatrixPolicy>(element);
108
109 for (unsigned ip = 0; ip < num_integration_points; ip++)
110 {
111 auto const& ns_and_weight = nss_and_weights[ip];
112 auto const coords = interpolateCoords(ns_and_weight, fe);
113
115 primary_variables_mat,
117 nss_and_weights[ip], prim_vars);
118
119 auto const [flag, flux, dFlux] =
120 bc_or_st_data.getFlagAndFluxAndDFlux(t, coords, prim_vars_data);
121
122 if (!flag)
123 {
124 // No flux value for this integration point. Skip assembly of
125 // the entire element.
126 return;
127 }
128
129 assembleLocalRhs(flux, ns_and_weight, local_rhs);
130
131 if (static_cast<int>(dFlux.size()) != num_comp_total)
132 {
133 // This strict check is technically mandatory only if a Jacobian
134 // is assembled. However, it is done as a consistency check also
135 // for cases without Jacobian assembly.
136 OGS_FATAL(
137 "The Python BC or ST must return the derivative of the "
138 "flux w.r.t. each component of each primary variable. "
139 "{:d} components expected. {:d} components returned from "
140 "Python."
141 "The expected number of components depends on multiple "
142 "factors. In the case of vectorial primary variables (such "
143 "as displacement), it will depend on the space dimension. "
144 "In the case of staggered coupling, the number of "
145 "components can be different for each staggered process "
146 "and different from monolithic coupling.",
147 num_comp_total, dFlux.size());
148 }
149
150 if (Jac)
151 {
152 assembleLocalJacobian(dFlux, ns_and_weight, local_Jac);
153 }
154 }
155
156 b.add(indices_this_component, local_rhs);
157
158 if (Jac)
159 {
161 indices_this_component, indices_all_components_for_Jac};
162 Jac->add(rci, local_Jac);
163 }
164 }
#define OGS_FATAL(...)
Definition Error.h:26
int add(IndexType row, IndexType col, double val)
Definition EigenMatrix.h:86
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:76
std::size_t getID() const
Get id of the mesh.
Definition Mesh.h:121
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
constexpr ranges::views::view_closure coords
Definition Mesh.h:232
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatricesType > createIsoparametricFiniteElement(MeshLib::Element const &e)
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)
void assembleLocalRhs(double const flux, typename Traits::NsAndWeight const &ns_and_weight, Eigen::VectorXd &local_rhs) const
std::array< double, 3 > interpolateCoords(typename Traits::NsAndWeight const &ns_and_weight, NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatrixPolicy > const &fe) const
void assembleLocalJacobian(std::vector< double > const &dFlux, typename Traits::NsAndWeight const &ns_and_weight, Eigen::MatrixXd &local_Jac) const
MeshLib::Mesh const & bc_or_st_mesh
The domain, where this BC or ST will be applied.
Definition BcOrStData.h:46
std::vector< std::reference_wrapper< ProcessVariable > > const & all_process_variables_for_this_process
Definition BcOrStData.h:51

References MathLib::EigenMatrix::add(), MathLib::EigenVector::add(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::all_process_variables_for_this_process, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assembleLocalJacobian(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assembleLocalRhs(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::bc_or_st_data, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::bc_or_st_mesh, ProcessLib::BoundaryConditionAndSourceTerm::Python::collectDofsToMatrix(), NumLib::createIsoparametricFiniteElement(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::element, MeshLib::Mesh::getID(), NumLib::getIndices(), NumLib::LocalToGlobalIndexMap::getNumberOfGlobalComponents(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::global_component_id, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::integration_method, ProcessLib::BoundaryConditionAndSourceTerm::Python::interpolate(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::interpolateCoords(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::nss_and_weights, OGS_FATAL, NumLib::LocalToGlobalIndexMap::size(), and MathLib::toVector().

◆ assembleLocalJacobian()

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
void ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assembleLocalJacobian ( std::vector< double > const & dFlux,
typename Traits::NsAndWeight const & ns_and_weight,
Eigen::MatrixXd & local_Jac ) const
inlineprivate

Assembles the term \(N_\mathrm{this}^T \cdot \mathrm{dFlux}_c \cdot \mathrm{weight} \cdot N_c \) for each component \(c\) of each primary variable and adds it to the passed local Jacobian.

The shape functions \(N_\alpha\) of different primary variables \(\alpha\) might differ, e.g., for Taylor-Hood elements.

The index this refers to the primary variable to which this boundary condition or source term belongs. The index \(c\) (called other in the source code in the method body) runs over all components of all primary variables of the ProcessLib::Process to which this boundary condition or source term belongs.

Definition at line 205 of file BcAndStLocalAssemblerImpl.h.

209 {
210 auto const& pv_refs =
212 auto const weight = ns_and_weight.weight();
213
214 auto const this_shp_fct_order = bc_or_st_data.shape_function_order;
215 auto const N_this = ns_and_weight.N(this_shp_fct_order);
216
217 Eigen::Index left = 0;
218
219 for (auto pv_ref : pv_refs)
220 {
221 auto const& pv = pv_ref.get();
222 auto const other_num_comp = pv.getNumberOfGlobalComponents();
223 auto const other_shp_fct_order = pv.getShapeFunctionOrder();
224 auto const N_other = ns_and_weight.N(other_shp_fct_order);
225 auto const width = N_other.size();
226
227 for (decltype(+other_num_comp) other_comp = 0;
228 other_comp < other_num_comp;
229 ++other_comp)
230 {
231 // The assignment -= takes into account the sign convention
232 // of 1st-order in time in ODE systems in OpenGeoSys.
233 local_Jac.middleCols(left, width).noalias() -=
234 N_this.transpose() * (dFlux[other_comp] * weight) * N_other;
235
236 left += width;
237 }
238 }
239
240 assert(left == local_Jac.cols());
241 }

References ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::all_process_variables_for_this_process, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::bc_or_st_data, ProcessLib::BoundaryConditionAndSourceTerm::Python::NsAndWeight< ShapeFunction, LowerOrderShapeFunction, ShapeMatrix, LowerOrderShapeMatrix >::N(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::shape_function_order, and ProcessLib::BoundaryConditionAndSourceTerm::Python::NsAndWeight< ShapeFunction, LowerOrderShapeFunction, ShapeMatrix, LowerOrderShapeMatrix >::weight().

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

◆ assembleLocalRhs()

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
void ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assembleLocalRhs ( double const flux,
typename Traits::NsAndWeight const & ns_and_weight,
Eigen::VectorXd & local_rhs ) const
inlineprivate

◆ interpolateCoords()

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
std::array< double, 3 > ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::interpolateCoords ( typename Traits::NsAndWeight const & ns_and_weight,
NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatrixPolicy > const & fe ) const
inlineprivate

Determines the coordinates that the point associated with the passed shape matrix by interpolating the element's node coordinates

Note
We use the shape function of higher order to interpolate coordinates.

Definition at line 172 of file BcAndStLocalAssemblerImpl.h.

176 {
177 auto const& N_higher = ns_and_weight.NHigherOrder();
178 return fe.interpolateCoordinates(N_higher);
179 }
std::array< double, 3 > interpolateCoordinates(typename ShapeMatrices::ShapeType const &N) const
Interpolates the coordinates of the element with the given shape matrix.

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::interpolateCoordinates(), and ProcessLib::BoundaryConditionAndSourceTerm::Python::NsAndWeight< ShapeFunction, LowerOrderShapeFunction, ShapeMatrix, LowerOrderShapeMatrix >::NHigherOrder().

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

Member Data Documentation

◆ bc_or_st_data

◆ element

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
MeshLib::Element const& ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::element

◆ integration_method

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
NumLib::GenericIntegrationMethod const& ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::integration_method

◆ nss_and_weights

template<typename BcOrStData , typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
std::vector<typename Traits::NsAndWeight> const ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::nss_and_weights

The documentation for this struct was generated from the following file: