File VectorPartitionedTyings.h#

Methods to switch between storage types based on a mesh and DOFs that are partitioned in:

  • unknown DOFs

  • prescribed DOFs

  • dependent 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 VectorPartitionedTyings : public GooseFEM::Vector
#include <GooseFEM/VectorPartitionedTyings.h>

Class to switch between storage types.

In particular:

  • “nodevec”: nodal vectors [nnode, ndim].

  • ”elemvec”: nodal vectors stored per element [nelem, nne, ndim].

  • ”dofval”: DOF values [ndof].

  • ”dofval_u”: DOF values (Unknown), == dofval[iiu], [nnu].

  • ”dofval_p”: DOF values (Prescribed), == dofval[iiu], [nnp].

Public Functions

VectorPartitionedTyings() = default#
template<class E, class M>
inline VectorPartitionedTyings(const E &conn, const E &dofs, const M &Cdu, const M &Cdp, const M &Cdi)

Constructor.

Template Parameters
Parameters
inline size_t nnd() const
Returns

Number of dependent DOFs.

inline size_t nni() const
Returns

Number of independent DOFs.

inline size_t nnu() const
Returns

Number of independent unknown DOFs.

inline size_t nnp() const
Returns

Number of independent prescribed DOFs.

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

Dependent DOFs (list of DOF numbers).

Returns

Pointer.

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

Independent DOFs (list of DOF numbers).

Returns

Copy.

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

Independent unknown DOFs (list of DOF numbers).

Returns

Pointer.

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

Independent prescribed DOFs (list of DOF numbers).

Returns

Pointer.

template<class T>
inline void copy_p(const T &dofval_src, T &dofval_dest) const

Copy (part of) “dofval” to another “dofval”: dofval_dest[iip()] = dofval_src[iip()].

Parameters
  • dofval_src – DOF values, iip() updated, [ndof].

  • dofval_dest – DOF values, iip() updated, [ndof].

template<class T>
inline array_type::tensor<double, 1> AsDofs_i(const T &nodevec) const

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

Only the independent DOFs are retained.

Parameters

nodevec – nodal vectors [nnode, ndim].

Returns

dofval[iii()] [nni].

template<class T, class R>
inline void asDofs_i(const T &nodevec, R &dofval_i, bool apply_tyings = true) const

Same as InterpQuad_vector(), but writing to preallocated return.

Parameters
  • nodevec – [nnode, ndim].

  • dofval_i – [nni].

  • apply_tyings – If true the dependent DOFs are eliminated.

Private Functions

template<class T>
inline Eigen::VectorXd Eigen_asDofs_d(const T &nodevec) const#

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

Only the dependent DOFs are retained.

Parameters

nodevec – nodal vectors [nnode, ndim].

Returns

dofval[iid()] [nnd].

Private Members

array_type::tensor<size_t, 1> m_iiu#

See iiu().

array_type::tensor<size_t, 1> m_iip#

See iip().

array_type::tensor<size_t, 1> m_iii#

See iii().

array_type::tensor<size_t, 1> m_iid#

See iid().

size_t m_nnu#

See nnu().

size_t m_nnp#

See nnp().

size_t m_nni#

See nni().

size_t m_nnd#

See nnd().

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_Cdi#

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

Eigen::SparseMatrix<double> m_Cud#

Transpose of “m_Cdu”.

Eigen::SparseMatrix<double> m_Cpd#

Transpose of “m_Cdp”.

Eigen::SparseMatrix<double> m_Cid#

Transpose of “m_Cdi”.