OGS
|
Classes | |
struct | BcAndStLocalAssemblerImpl |
struct | BcOrStData |
struct | FlagAndFluxAndDFlux |
struct | NsAndWeight |
struct | NsAndWeight< ShapeFunction, ShapeFunction, ShapeMatrix, ShapeMatrix > |
struct | NsAndWeightsTraits |
Collects common type aliases needed when working with NsAndWeight. More... | |
Functions | |
Eigen::MatrixXd | collectDofsToMatrix (MeshLib::Element const &element, std::size_t const mesh_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x) |
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) |
template<typename NsAndWeight > | |
void | interpolate (Eigen::MatrixXd const &primary_variables_mat, std::vector< std::reference_wrapper< ProcessVariable > > const &pv_refs, NsAndWeight const &ns_and_weight, Eigen::Ref< Eigen::VectorXd > interpolated_primary_variables) |
template<typename ShapeFunction , typename LowerOrderShapeFunction , int GlobalDim, typename IntegrationMethod > | |
auto | computeNsAndWeights (MeshLib::Element const &element, bool const is_axially_symmetric, IntegrationMethod const &integration_method) |
Eigen::MatrixXd ProcessLib::BoundaryConditionAndSourceTerm::Python::collectDofsToMatrix | ( | MeshLib::Element const & | element, |
std::size_t const | mesh_id, | ||
NumLib::LocalToGlobalIndexMap const & | dof_table, | ||
GlobalVector const & | x ) |
Collects the degrees of freedom of the passed element from the passed global vector into a matrix.
The dimensions of the returned matrix are #nodes x #total_components, i.e., each row of the matrix will contain all degrees of freedom at a specific node of the passed element.
The order of nodes is determined by the passed element. The order of components is determined by the passed dof_table.
Definition at line 74 of file CollectAndInterpolateNodalDof.cpp.
References NumLib::LocalToGlobalIndexMap::getGlobalComponent(), NumLib::LocalToGlobalIndexMap::getNumberOfGlobalComponents(), MeshLib::Element::getNumberOfNodes(), NumLib::LocalToGlobalIndexMap::getNumberOfVariableComponents(), and NumLib::LocalToGlobalIndexMap::getNumberOfVariables().
Eigen::VectorXd ProcessLib::BoundaryConditionAndSourceTerm::Python::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 ) |
Does the same as collectDofsToMatrix(), just for a single component and only on the base nodes of the passed mesh element.
Definition at line 106 of file CollectAndInterpolateNodalDof.cpp.
References MeshLib::Element::getNumberOfBaseNodes().
Referenced by ProcessLib::PythonBoundaryConditionLocalAssembler< ShapeFunction, LowerOrderShapeFunction, GlobalDim >::interpolate().
auto ProcessLib::BoundaryConditionAndSourceTerm::Python::computeNsAndWeights | ( | MeshLib::Element const & | element, |
bool const | is_axially_symmetric, | ||
IntegrationMethod const & | integration_method ) |
Computes shape matrices and integration weights for all integration points in the given mesh element.
Definition at line 149 of file NsAndWeight.h.
References NumLib::initShapeMatrices(), NumLib::N, and NumLib::N_J.
void ProcessLib::BoundaryConditionAndSourceTerm::Python::interpolate | ( | Eigen::MatrixXd const & | primary_variables_mat, |
std::vector< std::reference_wrapper< ProcessVariable > > const & | pv_refs, | ||
NsAndWeight const & | ns_and_weight, | ||
Eigen::Ref< Eigen::VectorXd > | interpolated_primary_variables ) |
Interpolates the passed matrix of degrees of freedom to the interpolated_primary_variables
via the passed shape matrices in ns_and_weight
.
The interpolation function order for each variable is determined by the passed ProcessVariable vector.
For the layout of primary_variables_mat
see collectDofsToMatrix().
Definition at line 63 of file CollectAndInterpolateNodalDof.h.