OGS
ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim > Class Template Referencefinal

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >

Definition at line 28 of file RobinBoundaryConditionLocalAssembler.h.

#include <RobinBoundaryConditionLocalAssembler.h>

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

Public Member Functions

 RobinBoundaryConditionLocalAssembler (MeshLib::Element const &e, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, RobinBoundaryConditionData const &data)
 
void assemble (std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const t, std::vector< GlobalVector * > const &xs, int const process_id, GlobalMatrix *K, GlobalVector &b, GlobalMatrix *Jac) override
 
- Public Member Functions inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >
 GenericNaturalBoundaryConditionLocalAssembler (MeshLib::Element const &e, bool is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method)
 
- Public Member Functions inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface
virtual ~GenericNaturalBoundaryConditionLocalAssemblerInterface ()=default
 

Public Attributes

 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
 

Private Types

using Base
 
using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
 

Private Attributes

RobinBoundaryConditionData const & _data
 
Base::NodalMatrixType _local_K
 
Base::NodalVectorType _local_rhs
 

Additional Inherited Members

- Protected Types inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >
using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
 
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
- Protected Attributes inherited from ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >
NumLib::GenericIntegrationMethod const & _integration_method
 
std::vector< NAndWeight, Eigen::aligned_allocator< NAndWeight > > const _ns_and_weights
 
MeshLib::Element const & _element
 

Member Typedef Documentation

◆ Base

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::Base
private
Initial value:
GenericNaturalBoundaryConditionLocalAssembler(MeshLib::Element const &e, bool is_axially_symmetric, NumLib::GenericIntegrationMethod const &integration_method)

Definition at line 32 of file RobinBoundaryConditionLocalAssembler.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 34 of file RobinBoundaryConditionLocalAssembler.h.

Constructor & Destructor Documentation

◆ RobinBoundaryConditionLocalAssembler()

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

Definition at line 37 of file RobinBoundaryConditionLocalAssembler.h.

43 : Base(e, is_axially_symmetric, integration_method),
44 _data(data),
45 _local_K(local_matrix_size, local_matrix_size),
46 _local_rhs(local_matrix_size)
47 {
48 }
GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim > Base

Member Function Documentation

◆ assemble()

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

Implements ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface.

Definition at line 50 of file RobinBoundaryConditionLocalAssembler.h.

55 {
56 _local_K.setZero();
57 _local_rhs.setZero();
58
59 auto& x = *xs[process_id];
60
61 auto const indices = NumLib::getIndices(id, dof_table_boundary);
62 auto const local_x = x.get(indices);
63 auto const u =
65 local_x, ShapeFunction::NPOINTS);
66
67 unsigned const n_integration_points =
69
70 typename Base::NodalVectorType const alpha =
72 .template topRows<ShapeFunction::MeshElement::n_all_nodes>();
73 typename Base::NodalVectorType const u_0 =
75 .template topRows<ShapeFunction::MeshElement::n_all_nodes>();
76
77 for (unsigned ip = 0; ip < n_integration_points; ++ip)
78 {
79 auto const& ip_data = Base::_ns_and_weights[ip];
80 auto const& N = ip_data.N;
81 auto const& w = ip_data.weight;
82
83 ParameterLib::SpatialPosition const position{
84 std::nullopt, Base::_element.getID(), ip,
88 Base::_element, N))};
89
90 double integral_measure = 1.0;
92 {
93 integral_measure = (*_data.integral_measure)(t, position)[0];
94 }
95
96 double const a = alpha.dot(N) * w * integral_measure;
97
98 // The local K matrix is used for both, the Newton and Picard
99 // methods.
100 _local_K.noalias() += N.transpose() * N * a;
101
102 if (Jac != nullptr)
103 {
104 _local_rhs.noalias() -= N.transpose() * (u - u_0).dot(N) * a;
105 }
106 else
107 {
108 _local_rhs.noalias() += N.transpose() * (u_0.dot(N) * a);
109 }
110 }
111
112 b.add(indices, _local_rhs);
113 if (Jac != nullptr)
114 {
116 indices),
117 _local_K);
118 }
119 else if (K != nullptr)
120 {
122 indices),
123 _local_K);
124 }
125 else
126 {
127 OGS_FATAL(
128 "In Robin boundary condition assembler, both, the Jacobian and "
129 "the K-matrices are null, but one matrix must be provided.");
130 }
131 }
#define OGS_FATAL(...)
Definition Error.h:26
int add(IndexType row, IndexType col, double val)
Definition EigenMatrix.h:87
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:76
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
std::vector< NAndWeight, Eigen::aligned_allocator< NAndWeight > > const _ns_and_weights
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
virtual Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > getNodalValuesOnElement(MeshLib::Element const &element, double const t) const
Returns a matrix of values for all nodes of the given element.
Definition Parameter.h:164
ParameterLib::Parameter< double > const *const integral_measure
ParameterLib::Parameter< double > const & alpha
ParameterLib::Parameter< double > const & u_0

References ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_data, ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_element, ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_local_K, ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_local_rhs, ProcessLib::GenericNaturalBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_ns_and_weights, MathLib::EigenMatrix::add(), MathLib::EigenVector::add(), ProcessLib::RobinBoundaryConditionData::alpha, MeshLib::Element::getID(), NumLib::getIndices(), ParameterLib::Parameter< T >::getNodalValuesOnElement(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::RobinBoundaryConditionData::integral_measure, NumLib::interpolateCoordinates(), OGS_FATAL, MathLib::toVector(), and ProcessLib::RobinBoundaryConditionData::u_0.

Member Data Documentation

◆ _data

template<typename ShapeFunction , int GlobalDim>
RobinBoundaryConditionData const& ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_data
private

◆ _local_K

template<typename ShapeFunction , int GlobalDim>
Base::NodalMatrixType ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_local_K
private

◆ _local_rhs

template<typename ShapeFunction , int GlobalDim>
Base::NodalVectorType ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::_local_rhs
private

◆ EIGEN_MAKE_ALIGNED_OPERATOR_NEW

template<typename ShapeFunction , int GlobalDim>
ProcessLib::RobinBoundaryConditionLocalAssembler< ShapeFunction, GlobalDim >::EIGEN_MAKE_ALIGNED_OPERATOR_NEW

Definition at line 140 of file RobinBoundaryConditionLocalAssembler.h.


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