OGS
ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim > Class Template Referencefinal

Detailed Description

template<typename ShapeFunction, typename LowerOrderShapeFunction, int GlobalDim>
class ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >

Definition at line 21 of file PythonBoundaryConditionLocalAssembler.h.

#include <PythonBoundaryConditionLocalAssembler.h>

Inheritance diagram for ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >:
[legend]

Public Member Functions

 PythonBoundaryConditionLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, PythonBcData const &data)
 
void assemble (std::size_t const boundary_element_id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const t, std::vector< GlobalVector * > const &xs, int const process_id, GlobalMatrix &, GlobalVector &b, GlobalMatrix *const Jac) override
 
double interpolate (unsigned const local_node_id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, GlobalVector const &x, int const var, int const comp) const override
 
- Public Member Functions inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface
virtual ~GenericNaturalBoundaryConditionLocalAssemblerInterface ()=default
 

Private Types

using LocAsmImpl
 
using Traits = typename LocAsmImpl::Traits
 

Private Member Functions

Traits::LowerOrderShapeMatrix computeLowerOrderShapeMatrix (unsigned const local_node_id) const
 

Private Attributes

LocAsmImpl const impl_
 

Member Typedef Documentation

◆ LocAsmImpl

template<typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
using ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::LocAsmImpl
private
Initial value:
ProcessLib::BoundaryConditionAndSourceTerm::Python::
BcAndStLocalAssemblerImpl<PythonBcData, ShapeFunction,
LowerOrderShapeFunction, GlobalDim>

Definition at line 24 of file PythonBoundaryConditionLocalAssembler.h.

◆ Traits

template<typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
using ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::Traits = typename LocAsmImpl::Traits
private

Definition at line 27 of file PythonBoundaryConditionLocalAssembler.h.

Constructor & Destructor Documentation

◆ PythonBoundaryConditionLocalAssembler()

template<typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::PythonBoundaryConditionLocalAssembler ( MeshLib::Element const & e,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool is_axially_symmetric,
PythonBcData const & data )
inline

Definition at line 30 of file PythonBoundaryConditionLocalAssembler.h.

36 : impl_{e, integration_method, is_axially_symmetric, data}
37 {
38 }

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
void ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::assemble ( std::size_t const boundary_element_id,
NumLib::LocalToGlobalIndexMap const & dof_table_boundary,
double const t,
std::vector< GlobalVector * > const & xs,
int const process_id,
GlobalMatrix & ,
GlobalVector & b,
GlobalMatrix *const Jac )
inlineoverridevirtual

Implements ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface.

Definition at line 40 of file PythonBoundaryConditionLocalAssembler.h.

45 {
46 impl_.assemble(boundary_element_id, dof_table_boundary, t,
47 *xs[process_id], b, Jac);
48 }

References ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::impl_.

◆ computeLowerOrderShapeMatrix()

template<typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
Traits::LowerOrderShapeMatrix ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::computeLowerOrderShapeMatrix ( unsigned const local_node_id) const
inlineprivate

Definition at line 76 of file PythonBoundaryConditionLocalAssembler.h.

78 {
79 using HigherOrderMeshElement = typename ShapeFunction::MeshElement;
80
81 assert(local_node_id < impl_.element.getNumberOfNodes());
82
83 std::array natural_coordss{MathLib::Point3d{NumLib::NaturalCoordinates<
84 HigherOrderMeshElement>::coordinates[local_node_id]}};
85
86 bool const is_axially_symmetric = false; // does not matter for N
87 auto const shape_matrices = NumLib::computeShapeMatrices<
88 LowerOrderShapeFunction,
89 typename Traits::LowerOrderShapeMatrixPolicy, GlobalDim,
90 NumLib::ShapeMatrixType::N>(impl_.element, is_axially_symmetric,
91 natural_coordss);
92 return shape_matrices.front().N;
93 }
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)

References NumLib::computeShapeMatrices(), ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::impl_, and NumLib::N.

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

◆ interpolate()

template<typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim>
double ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::interpolate ( unsigned const local_node_id,
NumLib::LocalToGlobalIndexMap const & dof_table_boundary,
GlobalVector const & x,
int const var,
int const comp ) const
inlineoverridevirtual

Interpolates the given component of the given variable to the given local_node_id.

The local_node_id is the number of the node within the current boundary element.

Implements ProcessLib::PythonBoundaryConditionLocalAssemblerInterface.

Definition at line 50 of file PythonBoundaryConditionLocalAssembler.h.

54 {
55 if constexpr (ShapeFunction::ORDER < 2 ||
56 LowerOrderShapeFunction::ORDER > 1)
57 {
58 return std::numeric_limits<double>::quiet_NaN();
59 }
60 else
61 {
62 auto const N = computeLowerOrderShapeMatrix(local_node_id);
63
64 auto const nodal_values_base_node =
67 impl_.element,
68 impl_.bc_or_st_data.bc_or_st_mesh.getID(),
69 dof_table_boundary, x, var, comp);
70
71 return N * nodal_values_base_node;
72 }
73 }
Traits::LowerOrderShapeMatrix computeLowerOrderShapeMatrix(unsigned const local_node_id) const
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)

References ProcessLib::BoundaryConditionAndSourceTerm::Python::collectDofsToMatrixOnBaseNodesSingleComponent(), ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::computeLowerOrderShapeMatrix(), and ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::impl_.

Member Data Documentation

◆ impl_


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