OGS
|
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>
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 |
|
private |
Definition at line 41 of file BcAndStLocalAssemblerImpl.h.
using ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::Traits |
Definition at line 37 of file BcAndStLocalAssemblerImpl.h.
|
inline |
Definition at line 44 of file BcAndStLocalAssemblerImpl.h.
|
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.
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().
|
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.
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().
|
inlineprivate |
Assembles the term \(N^T \cdot \mathrm{flux} \cdot \mathrm{weight} \) and adds it to the passed local rhs vector.
Definition at line 183 of file BcAndStLocalAssemblerImpl.h.
References 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().
|
inlineprivate |
Determines the coordinates that the point associated with the passed shape matrix by interpolating the element's node coordinates
Definition at line 172 of file BcAndStLocalAssemblerImpl.h.
References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::interpolateCoordinates(), and ProcessLib::BoundaryConditionAndSourceTerm::Python::NsAndWeight< ShapeFunction, LowerOrderShapeFunction, ShapeMatrix, LowerOrderShapeMatrix >::NHigherOrder().
BcOrStData const& ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::bc_or_st_data |
Definition at line 244 of file BcAndStLocalAssemblerImpl.h.
Referenced by ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assemble(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assembleLocalJacobian(), and ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assembleLocalRhs().
MeshLib::Element const& ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::element |
NumLib::GenericIntegrationMethod const& ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::integration_method |
std::vector<typename Traits::NsAndWeight> const ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::nss_and_weights |