HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
local_static_condensation.hpp
Go to the documentation of this file.
1 // Structure to store information for and perform local static condensation
2 
3 #ifndef LOCAL_STATIC_CONDENSATION_HPP
4 #define LOCAL_STATIC_CONDENSATION_HPP
5 
6 #include <Eigen/Dense>
7 
8 namespace HArDCore3D
9 {
10 
17 
25  {
28  const Eigen::MatrixXd & Perm,
29  const std::vector<size_t> & globalDOFs_sys,
30  const std::vector<size_t> & globalDOFs_sc
31  )
32  : m_Perm(Perm),
35  m_dim_sys(globalDOFs_sys.size()),
36  m_dim_sc(globalDOFs_sc.size())
37  {
38  // do nothing
39  };
40 
42 
43  std::tuple<Eigen::MatrixXd, Eigen::VectorXd, Eigen::MatrixXd, Eigen::VectorXd>
44  compute(const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT)
45  {
46  // We permute the DOFs that are statically condensed
47  Eigen::MatrixXd AT = m_Perm * lsT.first * m_Perm.transpose();
48  Eigen::VectorXd bT = m_Perm * lsT.second;
49 
50  Eigen::MatrixXd AT_sys, AT_sc;
51  Eigen::VectorXd bT_sys, bT_sc;
52 
53  if (m_dim_sc > 0){
54  // Extract 4 blocks of AT, bT for static condensation
55  Eigen::MatrixXd A11 = AT.topLeftCorner(m_dim_sys, m_dim_sys);
56  Eigen::MatrixXd A12 = AT.topRightCorner(m_dim_sys, m_dim_sc);
57  Eigen::MatrixXd A21 = AT.bottomLeftCorner(m_dim_sc, m_dim_sys);
58  Eigen::MatrixXd A22 = AT.bottomRightCorner(m_dim_sc, m_dim_sc);
59  Eigen::VectorXd b1 = bT.head(m_dim_sys);
60  Eigen::VectorXd b2 = bT.tail(m_dim_sc);
61  // Create condensed system (AT_reduced, bT_reduced) and SC recovery operator (RT, cT)
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;
67  AT_sc = -A22inv_A21;
68  bT_sc = A22inv_b2;
69  }else{
70  // No static condensation
71  AT_sys = AT;
72  bT_sys = bT;
73  AT_sc = Eigen::MatrixXd::Zero(0,0);
74  bT_sc = Eigen::VectorXd::Zero(0);
75  }
76 
77  return std::make_tuple(AT_sys, bT_sys, AT_sc, bT_sc);
78  };
79 
81  inline std::vector<size_t> globalDOFs_sys() { return m_globalDOFs_sys;};
83  inline size_t dim_sys() { return m_dim_sys;};
84 
86  inline std::vector<size_t> globalDOFs_sc() { return m_globalDOFs_sc;};
88  inline size_t dim_sc() { return m_dim_sc;};
89 
90  // Members
91  Eigen::MatrixXd m_Perm;
92  std::vector<size_t> m_globalDOFs_sys;
93  std::vector<size_t> m_globalDOFs_sc;
94  size_t m_dim_sys;
95  size_t m_dim_sc;
96  };
97 
98 
100 
101 } // end of namespace HArDCore3D
102 
103 #endif
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