File TyingsPeriodic.h#

Tools to store and apply nodal/DOF tyings.

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.

namespace Tyings#

Tools to store and apply nodal/DOF tyings.

class Control
#include <GooseFEM/TyingsPeriodic.h>

Add control nodes to an existing system.

Public Functions

Control() = default#
template<class C, class N>
inline Control(const C &coor, const N &dofs)

Constructor.

Template Parameters:
Parameters:
  • coor – Nodal coordinates [nnode, ndim].

  • dofs – DOF-numbers per node [nnode, ndim].

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

Nodal coordinates, for the system with control nodes added to it.

Returns:

[nnode + ndim, ndim], with nnode the number of nodes of the original system.

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

DOF-numbers per node, for the system with control nodes added to it.

Returns:

[nnode + ndim, ndim], with nnode the number of nodes of the original system.

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

DOF-numbers of each control node.

Returns:

[ndim, ndim].

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

Node-numbers of the control nodes.

Returns:

[ndim].

Private Members

array_type::tensor<double, 2> m_coor#

See coor().

array_type::tensor<size_t, 2> m_dofs#

See dofs().

array_type::tensor<size_t, 2> m_control_dofs#

See controlDofs().

array_type::tensor<size_t, 1> m_control_nodes#

See controlNodes().

class Periodic
#include <GooseFEM/TyingsPeriodic.h>

Nodal tyings per periodic boundary conditions.

The idea is that the displacement of all DOFs of a node are tied to another node and to the average displacement gradient. The latter is applied/measured using ‘virtual’ control nodes.

Consider the DOF list \( u \) renumbered such that it is split up in independent and dependent DOFs as follows

\( u = \begin{bmatrix} u_i \\ u_d \end{bmatrix}\)

whereby the independent DOFs are furthermore split up in unknown and prescribed nodes as follows

\( u_i = \begin{bmatrix} u_u \\ u_p \end{bmatrix}\)

such that

\( u = \begin{bmatrix} u_u \\ u_p \\ u_d \end{bmatrix}\)

Todo:

Document how the DOFs are tied to the control nodes, and what the has to do with the mean.

Public Functions

Periodic() = default#
template<class C, class D, class S, class T>
inline Periodic(const C &coor, const D &dofs, const S &control_dofs, const T &nodal_tyings)

Constructor.

Template Parameters:
Parameters:
  • coor – Nodal coordinates [nnode, ndim].

  • dofs – DOF-numbers per node [nnode, ndim].

  • control_dofs – DOF-numbers per control node [ndim, ndim].

  • nodal_tyings – List of nodal tyings, see nodal_tyings(). [ntyings, 2].

template<class C, class D, class S, class T, class U>
inline Periodic(const C &coor, const D &dofs, const S &control_dofs, const T &nodal_tyings, const U &iip)

Constructor.

Template Parameters:
Parameters:
  • coor – Nodal coordinates [nnode, ndim].

  • dofs – DOF-numbers per node [nnode, ndim].

  • control_dofs – DOF-numbers per control node [ndim, ndim].

  • nodal_tyings – List of nodal tyings, see nodal_tyings(). [ntyings, 2].

  • iip – List of prescribed DOF-numbers.

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, 2> &dofs() const
Returns:

DOF-numbers per node, as used internally (after renumbering), [nnode, ndim].

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

DOF-numbers for each control node, as used internally (after renumbering), [ndim, ndim].

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

Return the applied nodal tyings.

Per tying (row) two node numbers are specified, according to the convention (independent, dependent).

Returns:

[ntyings, 2].

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

Dependent DOFs.

Returns:

List of DOF numbers.

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

Independent DOFs.

Returns:

List of DOF numbers.

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

Independent unknown DOFs.

Returns:

List of DOF numbers.

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

Independent prescribed DOFs.

Returns:

List of DOF numbers.

inline Eigen::SparseMatrix<double> Cdi() const

Return tying matrix such as to get the dependent DOFs \( u_d \) from the independent DOFs \( u_i \) as follows.

\( u_d = C_{di} u_i \)

Note that this can be further partitioned in

\( u_d = C_{du} u_u + C_{dp} u_p \)

See Cdu() and Cdp().

Returns:

Sparse matrix.

inline Eigen::SparseMatrix<double> Cdu() const

Unknown part of the partitioned tying matrix, see Cdi().

Returns:

Sparse matrix.

inline Eigen::SparseMatrix<double> Cdp() const

Prescribed part of the partitioned tying matrix, see Cdi().

Returns:

Sparse matrix.

Private Members

size_t m_nnu#

See nnu().

size_t m_nnp#

See nnp().

size_t m_nni#

See nni().

size_t m_nnd#

See nnd().

size_t m_ndim#

Number of dimensions.

size_t m_nties#

Number of nodal ties.

array_type::tensor<size_t, 2> m_dofs#

See dofs().

array_type::tensor<size_t, 2> m_control#

See control().

array_type::tensor<size_t, 2> m_tyings#

See nodal_tyings().

array_type::tensor<double, 2> m_coor#

Nodal coordinates [nnode, ndim].