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

Detailed Description

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

Definition at line 21 of file PythonBoundaryConditionLocalAssembler.h.

#include <PythonBoundaryConditionLocalAssembler.h>

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

Public Member Functions

 PythonBoundaryConditionLocalAssembler (MeshLib::Element const &e, std::size_t const, bool is_axially_symmetric, unsigned const integration_order, 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
 
virtual double interpolate (unsigned const local_node_id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, GlobalVector const &x, int const var, int const comp) const =0
 
- Public Member Functions inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface
virtual ~GenericNaturalBoundaryConditionLocalAssemblerInterface ()=default
 
virtual void assemble (std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const t, std::vector< GlobalVector * > const &x, int const process_id, GlobalMatrix &K, GlobalVector &b, GlobalMatrix *Jac)=0
 

Private Types

using LocAsmImpl = ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< PythonBcData, ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >
 
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 , typename IntegrationMethod , int GlobalDim>
using ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::LocAsmImpl = ProcessLib::BoundaryConditionAndSourceTerm::Python:: BcAndStLocalAssemblerImpl<PythonBcData, ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim>
private

Definition at line 24 of file PythonBoundaryConditionLocalAssembler.h.

◆ Traits

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

Definition at line 28 of file PythonBoundaryConditionLocalAssembler.h.

Constructor & Destructor Documentation

◆ PythonBoundaryConditionLocalAssembler()

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

Definition at line 31 of file PythonBoundaryConditionLocalAssembler.h.

37 : impl_{e, is_axially_symmetric, integration_order, data}
38 {
39 }

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , typename LowerOrderShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, 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 41 of file PythonBoundaryConditionLocalAssembler.h.

46 {
47 impl_.assemble(boundary_element_id, dof_table_boundary, t,
48 *xs[process_id], b, Jac);
49 }
static const double t
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

References ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::assemble(), ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::impl_, and MathLib::t.

◆ computeLowerOrderShapeMatrix()

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

Definition at line 77 of file PythonBoundaryConditionLocalAssembler.h.

79 {
80 using HigherOrderMeshElement = typename ShapeFunction::MeshElement;
81
82 assert(local_node_id < impl_.element.getNumberOfNodes());
83
84 std::array natural_coordss{MathLib::Point3d{NumLib::NaturalCoordinates<
85 HigherOrderMeshElement>::coordinates[local_node_id]}};
86
87 bool const is_axially_symmetric = false; // does not matter for N
88 auto const shape_matrices = NumLib::computeShapeMatrices<
89 LowerOrderShapeFunction,
90 typename Traits::LowerOrderShapeMatrixPolicy, GlobalDim,
91 NumLib::ShapeMatrixType::N>(impl_.element, is_axially_symmetric,
92 natural_coordss);
93 return shape_matrices.front().N;
94 }
virtual unsigned getNumberOfNodes() const =0
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::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::element, MeshLib::Element::getNumberOfNodes(), ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::impl_, and NumLib::N.

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

◆ interpolate()

template<typename ShapeFunction , typename LowerOrderShapeFunction , typename IntegrationMethod , int GlobalDim>
double ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, 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 51 of file PythonBoundaryConditionLocalAssembler.h.

55 {
56 if constexpr (ShapeFunction::ORDER < 2 ||
57 LowerOrderShapeFunction::ORDER > 1)
58 {
59 return std::numeric_limits<double>::quiet_NaN();
60 }
61 else
62 {
63 auto const N = computeLowerOrderShapeMatrix(local_node_id);
64
65 auto const nodal_values_base_node =
70 dof_table_boundary, x, var, comp);
71
72 return N * nodal_values_base_node;
73 }
74 }
std::size_t getID() const
Get id of the mesh.
Definition: Mesh.h:115
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)
MeshLib::Mesh const & bc_or_st_mesh
The domain, where this BC or ST will be applied.
Definition: BcOrStData.h:46

References ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::bc_or_st_data, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::bc_or_st_mesh, ProcessLib::BoundaryConditionAndSourceTerm::Python::collectDofsToMatrixOnBaseNodesSingleComponent(), ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::computeLowerOrderShapeMatrix(), ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::element, MeshLib::Mesh::getID(), and ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >::impl_.

Member Data Documentation

◆ impl_


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