3 #ifndef LOCAL_STATIC_CONDENSATION_HPP
4 #define LOCAL_STATIC_CONDENSATION_HPP
28 const Eigen::MatrixXd & Perm,
43 std::tuple<Eigen::MatrixXd, Eigen::VectorXd, Eigen::MatrixXd, Eigen::VectorXd>
44 compute(
const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT)
47 Eigen::MatrixXd AT =
m_Perm * lsT.first *
m_Perm.transpose();
48 Eigen::VectorXd bT =
m_Perm * lsT.second;
50 Eigen::MatrixXd AT_sys, AT_sc;
51 Eigen::VectorXd bT_sys, bT_sc;
60 Eigen::VectorXd b2 = bT.tail(
m_dim_sc);
62 Eigen::PartialPivLU<Eigen::MatrixXd> A22pivlu(A22);
63 Eigen::MatrixXd A22inv_A21 = A22pivlu.solve(A21);
64 Eigen::VectorXd A22inv_b2 = A22pivlu.solve(b2);
65 AT_sys = A11 - A12 * A22inv_A21;
66 bT_sys = b1 - A12 * A22inv_b2;
73 AT_sc = Eigen::MatrixXd::Zero(0,0);
74 bT_sc = Eigen::VectorXd::Zero(0);
77 return std::make_tuple(AT_sys, bT_sys, AT_sc, bT_sc);
std::vector< size_t > globalDOFs_sc()
Returns global DOFs that are statically condensend.
Definition: local_static_condensation.hpp:86
std::vector< size_t > m_globalDOFs_sc
Definition: local_static_condensation.hpp:93
size_t dim_sc()
Returns the number of DOFs that are statically condensed.
Definition: local_static_condensation.hpp:88
size_t m_dim_sc
Definition: local_static_condensation.hpp:95
std::tuple< Eigen::MatrixXd, Eigen::VectorXd, Eigen::MatrixXd, Eigen::VectorXd > compute(const std::pair< Eigen::MatrixXd, Eigen::VectorXd > &lsT)
Compute the local static condensation.
Definition: local_static_condensation.hpp:44
LocalStaticCondensation(const Eigen::MatrixXd &Perm, const std::vector< size_t > &globalDOFs_sys, const std::vector< size_t > &globalDOFs_sc)
Constructor.
Definition: local_static_condensation.hpp:27
size_t dim_sys()
Returns the number of DOFs that are not statically condensed.
Definition: local_static_condensation.hpp:83
std::vector< size_t > m_globalDOFs_sys
Definition: local_static_condensation.hpp:92
std::vector< size_t > globalDOFs_sys()
Returns global DOFs that are not statically condensend.
Definition: local_static_condensation.hpp:81
size_t m_dim_sys
Definition: local_static_condensation.hpp:94
Eigen::MatrixXd m_Perm
Definition: local_static_condensation.hpp:88
Definition: ddr-magnetostatics.hpp:40
Structure to store information for, and perform, local static condensation.
Definition: local_static_condensation.hpp:25