File ElementHex8.h#

Quadrature for 8-noded hexahedral element in 3d (GooseFEM::Mesh::ElementType::Hex8), in a Cartesian coordinate system.

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

Element quadrature and interpolation.

namespace Hex8#

8-noded hexahedral element in 3d (GooseFEM::Mesh::ElementType::Hex8).

class Quadrature : public GooseFEM::Element::QuadratureBaseCartesian<Quadrature>
#include <GooseFEM/ElementHex8.h>

Interpolation and quadrature.

Fixed dimensions:

  • ndim = 3: number of dimensions.

  • nne = 8: number of nodes per element.

Naming convention:

Public Functions

Quadrature() = default#
template<class T>
inline Quadrature(const T &x)

Constructor: use default Gauss integration.

The following is pre-computed during construction:

  • the shape functions,

  • the shape function gradients (in local and global) coordinates,

  • the integration points volumes. They can be reused without any cost. They only have to be recomputed when the nodal position changes (note that they are assumed to be constant under a small-strain assumption). In that case use update_x() to update the nodal positions and to recompute the above listed quantities.

Parameters

x – nodal coordinates (elemvec).

template<class T, class X, class W>
inline Quadrature(const T &x, const X &xi, const W &w)

Constructor with custom integration.

The following is pre-computed during construction:

  • the shape functions,

  • the shape function gradients (in local and global) coordinates,

  • the integration points volumes. They can be reused without any cost. They only have to be recomputed when the nodal position changes (note that they are assumed to be constant under a small-strain assumption). In that case use update_x() to update the nodal positions and to recompute the above listed quantities.

Parameters
  • x – nodal coordinates (elemvec).

  • xi – Integration point coordinates (local coordinates) [nip].

  • w – Integration point weights [nip].

Private Functions

template<class T, class R>
inline void int_N_scalar_NT_dV_impl(const T &qscalar, R &elemmat) const#
template<class T, class R>
inline void int_gradN_dot_tensor2_dV_impl(const T &qtensor, R &elemvec) const#

Private Members

friend QuadratureBase< Quadrature >
friend QuadratureBaseCartesian< Quadrature >
size_t m_tdim = 3#

Dynamic alias of s_tdim (remove in C++17)

size_t m_nelem#

Number of elements.

size_t m_nip#

Number of integration points per element.

array_type::tensor<double, 3> m_x#

nodal positions stored per element [nelem, nne, ndim]

array_type::tensor<double, 1> m_w#

weight of each integration point [nip]

array_type::tensor<double, 2> m_xi#

local coordinate per integration point [nip, ndim]

array_type::tensor<double, 2> m_N#

shape functions [nip, nne]

array_type::tensor<double, 3> m_dNxi#

local shape func grad [nip, nne, ndim]

array_type::tensor<double, 4> m_dNx#

global shape func grad [nelem, nip, nne, ndim]

array_type::tensor<double, 2> m_vol#

integration point volume [nelem, nip]

Private Static Attributes

static constexpr size_t s_nne = 8#

Number of nodes per element.

static constexpr size_t s_ndim = 3#

Number of dimensions for nodal vectors.

static constexpr size_t s_tdim = 3#

Number of dimensions for tensors.

namespace Gauss#

gauss quadrature: quadrature points such that integration is exact for these bi-linear elements::

Functions

inline size_t nip()#

Number of integration points:

 nip = nne = 8

Returns

unsigned int

inline array_type::tensor<double, 2> xi()#

Integration point coordinates (local coordinates).

Returns

Coordinates [nip, ndim], with ndim = 3.

inline array_type::tensor<double, 1> w()#

Integration point weights.

Returns

Coordinates [nip].

namespace Nodal#

nodal quadrature: quadrature points coincide with the nodes.

The order is the same as in the connectivity.

Functions

inline size_t nip()#

Number of integration points:

 nip = nne = 8

Returns

unsigned int

inline array_type::tensor<double, 2> xi()#

Integration point coordinates (local coordinates).

Returns

Coordinates [nip, ndim], with ndim = 3.

inline array_type::tensor<double, 1> w()#

Integration point weights.

Returns

Coordinates [nip].