File MatrixPartitionedTyings.h#

Sparse matrix that is partitioned in:

  • unknown DOFs

  • prescribed DOFs

  • tied DOFs

Copyright

Copyright 2017. Tom de Geus. All rights reserved.

License: This project is released under the GNU Public License (GPLv3).

namespace GooseFEM

Toolbox to perform finite element computations.

class MatrixPartitionedTyings : public GooseFEM::MatrixPartitionedTyingsBase<MatrixPartitionedTyings>
#include <GooseFEM/MatrixPartitionedTyings.h>

Sparse matrix from with dependent DOFs are eliminated, and the remaining (small) independent system is partitioned in an unknown and a prescribed part.

In particular:

\( A_{ii} = \begin{bmatrix} A_{uu} & A_{up} \\ A_{pu} & A_{pp} \end{bmatrix} \)

See VectorPartitionedTyings() for bookkeeping definitions.

Public Functions

MatrixPartitionedTyings() = default#
inline MatrixPartitionedTyings(const array_type::tensor<size_t, 2> &conn, const array_type::tensor<size_t, 2> &dofs, const Eigen::SparseMatrix<double> &Cdu, const Eigen::SparseMatrix<double> &Cdp)

Constructor.

Parameters
inline const Eigen::SparseMatrix<double> &data_uu() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_up() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_pu() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_pp() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_ud() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_pd() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_du() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_dp() const

Pointer to data.

inline const Eigen::SparseMatrix<double> &data_dd() const

Pointer to data.

Protected Attributes

Eigen::SparseMatrix<double> m_Auu#

The matrix.

Eigen::SparseMatrix<double> m_Aup#

The matrix.

Eigen::SparseMatrix<double> m_Apu#

The matrix.

Eigen::SparseMatrix<double> m_App#

The matrix.

Eigen::SparseMatrix<double> m_Aud#

The matrix.

Eigen::SparseMatrix<double> m_Apd#

The matrix.

Eigen::SparseMatrix<double> m_Adu#

The matrix.

Eigen::SparseMatrix<double> m_Adp#

The matrix.

Eigen::SparseMatrix<double> m_Add#

The matrix.

Eigen::SparseMatrix<double> m_ACuu#

// The matrix for which the tyings have been applied.

Eigen::SparseMatrix<double> m_ACup#

// The matrix for which the tyings have been applied.

Eigen::SparseMatrix<double> m_ACpu#

// The matrix for which the tyings have been applied.

Eigen::SparseMatrix<double> m_ACpp#

// The matrix for which the tyings have been applied.

std::vector<Eigen::Triplet<double>> m_Tuu#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tup#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tpu#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tpp#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tud#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tpd#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tdu#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tdp#

Matrix entries.

std::vector<Eigen::Triplet<double>> m_Tdd#

Matrix entries.

Eigen::SparseMatrix<double> m_Cdu#

Tying matrix, see Tyings::Periodic::Cdu().

Eigen::SparseMatrix<double> m_Cdp#

Tying matrix, see Tyings::Periodic::Cdp().

Eigen::SparseMatrix<double> m_Cud#

Transpose of “m_Cdu”.

Eigen::SparseMatrix<double> m_Cpd#

Transpose of “m_Cdp”.

Private Functions

template<class T>
inline void assemble_impl(const T &elemmat)#
template<class T>
inline void todense_impl(T &ret) const#
template<class T>
inline void dot_nodevec_impl(const T &x, T &b) const#
template<class T>
inline void dot_dofval_impl(const T &x, T &b) const#
template<class T>
inline void reaction_nodevec_impl(const T &x, T &b) const#
template<class T>
inline void reaction_dofval_impl(const T &x, T &b) const#
inline void reaction_p_impl(const array_type::tensor<double, 1> &x_u, const array_type::tensor<double, 1> &x_p, array_type::tensor<double, 1> &b_p) const#
inline Eigen::VectorXd AsDofs_u(const array_type::tensor<double, 1> &dofval) const#
inline Eigen::VectorXd AsDofs_u(const array_type::tensor<double, 2> &nodevec) const#
inline Eigen::VectorXd AsDofs_p(const array_type::tensor<double, 1> &dofval) const#
inline Eigen::VectorXd AsDofs_p(const array_type::tensor<double, 2> &nodevec) const#
inline Eigen::VectorXd AsDofs_d(const array_type::tensor<double, 1> &dofval) const#
inline Eigen::VectorXd AsDofs_d(const array_type::tensor<double, 2> &nodevec) const#

Private Members

friend MatrixBase< MatrixPartitionedTyings >
friend MatrixPartitionedBase< MatrixPartitionedTyings >
friend MatrixPartitionedTyingsBase< MatrixPartitionedTyings >

Friends

friend class MatrixPartitionedTyingsSolver
template<class Solver = Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>>>
class MatrixPartitionedTyingsSolver : public MatrixSolverBase<MatrixPartitionedTyingsSolver<Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>>>>, public MatrixSolverPartitionedBase<MatrixPartitionedTyingsSolver<Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>>>>
#include <GooseFEM/MatrixPartitionedTyings.h>

Solver for MatrixPartitionedTyings().

This solver class can be used to solve for multiple right-hand-sides using one factorisation.

Solving proceeds as follows:

\( A' = A_{ii} + A_{id} * C_{di} + C_{di}^T * A_{di} + C_{di}^T * A_{dd} * C_{di} \)

\( b' = b_i + C_{di}^T * b_d \)

\( x_u = A'_{uu} \ ( b'_u - A'_{up} * x_p ) \)

\( x_i = \begin{bmatrix} x_u \\ x_p \end{bmatrix} \)

\( x_d = C_{di} * x_i \)

Public Functions

MatrixPartitionedTyingsSolver() = default#
inline array_type::tensor<double, 1> Solve_u(MatrixPartitionedTyings &A, const array_type::tensor<double, 1> &b_u, const array_type::tensor<double, 1> &b_d, const array_type::tensor<double, 1> &x_p)

Same as Solve(MatrixPartitionedTyings&, const array_type::tensor<double, 2>&, const

array_type::tensor<double, 2>&), but with partitioned input and output.

Parameters
  • A – sparse matrix, see MatrixPartitionedTyings().

  • b_u – unknown dofval [nnu].

  • b_d – dependent dofval [nnd].

  • x_p – prescribed dofval [nnp]

Returns

x_u unknown dofval [nnu].

inline void solve_u(MatrixPartitionedTyings &A, const array_type::tensor<double, 1> &b_u, const array_type::tensor<double, 1> &b_d, const array_type::tensor<double, 1> &x_p, array_type::tensor<double, 1> &x_u)

Same as Solve_u(MatrixPartitionedTyings&, const array_type::tensor<double, 1>&, constarray_type::tensor<double, 1>&, const array_type::tensor<double, 1>&), but writing to pre-allocated output.

Parameters
  • A – sparse matrix, see MatrixPartitionedTyings().

  • b_u – unknown dofval [nnu].

  • b_d – dependent dofval [nnd].

  • x_p – prescribed dofval [nnp]

  • x_u – (overwritten) unknown dofval [nnu].

Private Functions

template<class T>
inline void solve_nodevec_impl(MatrixPartitionedTyings &A, const T &b, T &x)#
template<class T>
inline void solve_dofval_impl(MatrixPartitionedTyings &A, const T &b, T &x)#
inline void factorize(MatrixPartitionedTyings &A)#

compute inverse (evaluated by “solve”)

Private Members

friend MatrixSolverBase< MatrixPartitionedTyingsSolver< Solver > >
friend MatrixSolverPartitionedBase< MatrixPartitionedTyingsSolver< Solver > >
Solver m_solver#

solver

bool m_factor = true#

signal to force factorization