template<typename ShapeFunction, typename LowerOrderShapeFunction, typename IntegrationMethod, int GlobalDim>
class ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, IntegrationMethod, GlobalDim >
Definition at line 21 of file PythonBoundaryConditionLocalAssembler.h.
|
| 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 |
|
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 |
|
template<typename ShapeFunction , typename LowerOrderShapeFunction , typename IntegrationMethod , int GlobalDim>
Definition at line 77 of file PythonBoundaryConditionLocalAssembler.h.
79 {
80 using HigherOrderMeshElement = typename ShapeFunction::MeshElement;
81
83
85 HigherOrderMeshElement>::coordinates[local_node_id]}};
86
87 bool const is_axially_symmetric = false;
89 LowerOrderShapeFunction,
90 typename Traits::LowerOrderShapeMatrixPolicy, GlobalDim,
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)
MeshLib::Element const & element
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().
template<typename ShapeFunction , typename LowerOrderShapeFunction , typename IntegrationMethod , int GlobalDim>
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 {
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.
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)
BcOrStData const & bc_or_st_data
MeshLib::Mesh const & bc_or_st_mesh
The domain, where this BC or ST will be applied.
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_.