File Matrix.h#

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 Matrix : public GooseFEM::MatrixBase<Matrix>
#include <GooseFEM/Matrix.h>

Sparse matrix.

See GooseFEM::Vector() for bookkeeping definitions.

Public Functions

Matrix() = default#
inline Matrix(const array_type::tensor<size_t, 2> &conn, const array_type::tensor<size_t, 2> &dofs)

Constructor.

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

Pointer to data.

inline void set(const array_type::tensor<size_t, 1> &rows, const array_type::tensor<size_t, 1> &cols, const array_type::tensor<double, 2> &matrix)

Overwrite matrix.

Parameters
  • rows – Row numbers [m].

  • cols – Column numbers [n].

  • matrix – Data entries matrix(i, j) for rows(i), cols(j) [m, n].

inline void add(const array_type::tensor<size_t, 1> &rows, const array_type::tensor<size_t, 1> &cols, const array_type::tensor<double, 2> &matrix)

Add matrix.

Parameters
  • rows – Row numbers [m].

  • cols – Column numbers [n].

  • matrix – Data entries matrix(i, j) for rows(i), cols(j) [m, n].

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 Eigen::VectorXd Eigen_AsDofs_nodevec(const T &nodevec) const#

Convert “nodevec” to “dofval” (overwrite entries that occur more than once).

Parameters

nodevec – input [nnode, ndim]

Returns

dofval output [ndof]

template<class T>
inline void Eigen_asNode_dofval_nodevec(const Eigen::VectorXd &dofval, T &nodevec) const#

Convert “dofval” to “nodevec” (overwrite entries that occur more than once).

Parameters

Private Members

friend MatrixBase< Matrix >
bool m_changed = true#

Signal changes to data.

Eigen::SparseMatrix<double> m_A#

The matrix.

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

Matrix entries.

Friends

friend class MatrixSolver
template<class D>
class MatrixBase
#include <GooseFEM/Matrix.h>

CRTP base class for a matrix.

Subclassed by GooseFEM::MatrixPartitionedBase< MatrixDiagonalPartitioned >, GooseFEM::MatrixPartitionedBase< MatrixPartitioned >, GooseFEM::MatrixPartitionedBase< MatrixPartitionedTyings >, GooseFEM::MatrixPartitionedBase< D >

Public Types

using derived_type = D

Underlying type.

Public Functions

inline size_t nelem() const

Number of elements.

Returns

Unsigned integer.

inline size_t nne() const

Number of nodes per element.

Returns

Unsigned integer.

inline size_t nnode() const

Number of nodes.

Returns

Unsigned integer.

inline size_t ndim() const

Number of dimensions.

Returns

Unsigned integer.

inline size_t ndof() const

Number of DOFs.

Returns

Unsigned integer.

inline const array_type::tensor<size_t, 2> &dofs() const

DOFs per node.

Returns

[nnode, ndim].

inline const array_type::tensor<size_t, 2> &conn() const

Connectivity.

Returns

[nelem, nne].

inline std::array<size_t, 1> shape_dofval() const

Shape of “dofval”.

Returns

[ndof].

inline std::array<size_t, 2> shape_nodevec() const

Shape of “nodevec”.

Returns

[nnode, ndim].

inline std::array<size_t, 3> shape_elemmat() const

Shape of “elemmat”.

Returns

[nelem, nne * ndim, nne * ndim].

template<class T>
inline void assemble(const T &elemmat)

Assemble from “elemmat”.

Parameters

elemmat – [nelem, nne * ndim, nne * ndim].

inline array_type::tensor<double, 2> Todense() const

Copy as dense matrix.

Returns

[ndof, ndof].

template<class T>
inline void todense(T &ret) const

Copy to dense matrix.

Parameters

ret – overwritten [ndof, ndof].

inline array_type::tensor<double, 2> Dot(const array_type::tensor<double, 2> &x) const

Dot-product \( b_i = A_{ij} x_j \).

Parameters

x – nodevec [nelem, ndim].

Returns

b nodevec [nelem, ndim].

inline array_type::tensor<double, 1> Dot(const array_type::tensor<double, 1> &x) const

Dot-product \( b_i = A_{ij} x_j \).

Parameters

x – dofval [ndof].

Returns

b dofval [ndof].

inline void dot(const array_type::tensor<double, 2> &x, array_type::tensor<double, 2> &b) const

Dot-product \( b_i = A_{ij} x_j \).

Parameters
inline void dot(const array_type::tensor<double, 1> &x, array_type::tensor<double, 1> &b) const

Dot-product \( b_i = A_{ij} x_j \).

Parameters
  • x – dofval [ndof].

  • b – (overwritten) dofval [ndof].

Protected Attributes

array_type::tensor<size_t, 2> m_conn#

Connectivity [nelem, nne].

array_type::tensor<size_t, 2> m_dofs#

DOF-numbers per node [nnode, ndim].

size_t m_nelem#

See nelem().

size_t m_nne#

See nne().

size_t m_nnode#

See nnode().

size_t m_ndim#

See ndim().

size_t m_ndof#

See ndof().

bool m_changed = true#

Signal changes to data.

Private Functions

inline auto derived_cast() -> derived_type&#
inline auto derived_cast() const -> const derived_type&#
template<class D>
class MatrixPartitionedBase : public GooseFEM::MatrixBase<D>
#include <GooseFEM/Matrix.h>

CRTP base class for a partitioned matrix.

Subclassed by GooseFEM::MatrixPartitionedTyingsBase< MatrixPartitionedTyings >, GooseFEM::MatrixPartitionedTyingsBase< D >

Public Types

using derived_type = D

Underlying type.

Public Functions

inline size_t nnu() const

Number of unknown DOFs.

Returns

Unsigned integer.

inline size_t nnp() const

Number of prescribed DOFs.

Returns

Unsigned integer.

inline const array_type::tensor<size_t, 1> &iiu() const

Unknown DOFs.

Returns

[nnu].

inline const array_type::tensor<size_t, 1> &iip() const

Prescribed DOFs.

Returns

[nnp].

inline array_type::tensor<double, 2> Reaction(const array_type::tensor<double, 2> &x, const array_type::tensor<double, 2> &b) const

Right-hand-size for corresponding to the prescribed DOFs:

\( b_p = A_{pu} * x_u + A_{pp} * x_p \)

and assemble them to the appropriate places in “nodevec”.

Parameters
Returns

Copy of b with \( b_p \) overwritten.

inline array_type::tensor<double, 1> Reaction(const array_type::tensor<double, 1> &x, const array_type::tensor<double, 1> &b) const

Same as Reaction(const array_type::tensor<double, 2>&, const array_type::tensor<double, 2>&), but of “dofval” input and output.

Parameters
  • x – “dofval” [ndof].

  • b – “dofval” [ndof].

Returns

Copy of b with \( b_p \) overwritten.

inline void reaction(const array_type::tensor<double, 2> &x, array_type::tensor<double, 2> &b) const

Same as Reaction(const array_type::tensor<double, 2>&, const array_type::tensor<double, 2>&), but inserting in-place.

Parameters
  • x – “nodevec” [nnode, ndim].

  • b – “nodevec” [nnode, ndim], \( b_p \) overwritten.

inline void reaction(const array_type::tensor<double, 1> &x, array_type::tensor<double, 1> &b) const

Same as Reaction(const array_type::tensor<double, 1>&, const array_type::tensor<double, 1>&), but inserting in-place.

Parameters
  • x – “dofval” [ndof].

  • b – “dofval” [ndof], \( b_p \) overwritten.

inline array_type::tensor<double, 1> Reaction_p(const array_type::tensor<double, 1> &x_u, const array_type::tensor<double, 1> &x_p) const

Same as Reaction(const array_type::tensor<double, 1>&, const array_type::tensor<double, 1>&), but with partitioned input and output.

Parameters
  • x_u – unknown “dofval” [nnu].

  • x_p – prescribed “dofval” [nnp].

Returns

b_p prescribed “dofval” [nnp].

inline void reaction_p(const array_type::tensor<double, 1> &x_u, const array_type::tensor<double, 1> &x_p, array_type::tensor<double, 1> &b_p) const

Same as Reaction_p(const array_type::tensor<double, 1>&, const array_type::tensor<double,

1>&), but writing to preallocated output.

Parameters
  • x_u – unknown “dofval” [nnu].

  • x_p – prescribed “dofval” [nnp].

  • b_p – (overwritten) prescribed “dofval” [nnp].

Protected Attributes

array_type::tensor<size_t, 1> m_iiu#

See iiu()

array_type::tensor<size_t, 1> m_iip#

See iip()

size_t m_nnu#

See nnu.

size_t m_nnp#

See nnp.

Private Functions

inline auto derived_cast() -> derived_type&#
inline auto derived_cast() const -> const derived_type&#
template<class D>
class MatrixPartitionedTyingsBase : public GooseFEM::MatrixPartitionedBase<D>
#include <GooseFEM/Matrix.h>

CRTP base class for a partitioned matrix with tying.

Public Types

using derived_type = D

Underlying type.

Public Functions

inline size_t nni() const

Number of independent DOFs.

Returns

Unsigned integer.

inline size_t nnd() const

Number of dependent DOFs.

Returns

Unsigned integer.

inline const array_type::tensor<size_t, 1> &iii() const

Independent DOFs.

Returns

[nnu].

inline const array_type::tensor<size_t, 1> &iid() const

Dependent DOFs.

Returns

[nnp].

Protected Attributes

array_type::tensor<size_t, 1> m_iii#

See iii()

array_type::tensor<size_t, 1> m_iid#

See iid()

size_t m_nni#

See nni.

size_t m_nnd#

See nnd.

Private Functions

inline auto derived_cast() -> derived_type&#
inline auto derived_cast() const -> const derived_type&#
template<class Solver = Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>>>
class MatrixSolver : public MatrixSolverBase<MatrixSolver<Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>>>>, public MatrixSolverSingleBase<MatrixSolver<Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>>>>
#include <GooseFEM/Matrix.h>

Solve \( x = A^{-1} b \), for A of the GooseFEM::Matrix() class.

You can solve for multiple right-hand-sides using one factorisation.

Public Functions

MatrixSolver() = default#

Private Functions

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

Compute inverse (evaluated by “solve”).

Private Members

friend MatrixSolverBase< MatrixSolver< Solver > >
friend MatrixSolverSingleBase< MatrixSolver< Solver > >
Solver m_solver#

Solver.

bool m_factor = true#

Signal to force factorization.

template<class D>
class MatrixSolverBase
#include <GooseFEM/Matrix.h>

CRTP base class for a solver class.

Public Types

using derived_type = D

Underlying type.

Public Functions

template<class M>
inline void solve(M &A, const array_type::tensor<double, 2> &b, array_type::tensor<double, 2> &x)

Solve \( x = A^{-1} b \).

Parameters
  • AGooseFEM (sparse) matrix, see e.g. GooseFEM::Matrix().

  • b – nodevec [nelem, ndim].

  • x – (overwritten) nodevec [nelem, ndim].

template<class M>
inline void solve(M &A, const array_type::tensor<double, 1> &b, array_type::tensor<double, 1> &x)

Solve \( x = A^{-1} b \).

Parameters

Private Functions

inline auto derived_cast() -> derived_type&#
inline auto derived_cast() const -> const derived_type&#
template<class D>
class MatrixSolverPartitionedBase
#include <GooseFEM/Matrix.h>

CRTP base class for a extra functions for a partitioned solver class.

Public Types

using derived_type = D

Underlying type.

Public Functions

template<class M>
inline array_type::tensor<double, 2> Solve(M &A, const array_type::tensor<double, 2> &b, const array_type::tensor<double, 2> &x)

Solve \( x = A^{-1} b \).

Parameters
Returns

x nodevec [nelem, ndim].

template<class M>
inline array_type::tensor<double, 1> Solve(M &A, const array_type::tensor<double, 1> &b, const array_type::tensor<double, 1> &x)

Solve \( x = A^{-1} b \).

Parameters
Returns

x dofval [ndof].

Private Functions

inline auto derived_cast() -> derived_type&#
inline auto derived_cast() const -> const derived_type&#
template<class D>
class MatrixSolverSingleBase
#include <GooseFEM/Matrix.h>

CRTP base class for a solver class.

Public Types

using derived_type = D

Underlying type.

Public Functions

template<class M>
inline array_type::tensor<double, 2> Solve(M &A, const array_type::tensor<double, 2> &b)

Solve \( x = A^{-1} b \).

Parameters
Returns

x nodevec [nelem, ndim].

template<class M>
inline array_type::tensor<double, 1> Solve(M &A, const array_type::tensor<double, 1> &b)

Solve \( x = A^{-1} b \).

Parameters
Returns

x dofval [ndof].

Private Functions

inline auto derived_cast() -> derived_type&#
inline auto derived_cast() const -> const derived_type&#