File ElementQuad4Planar.h#

Quadrature for 4-noded quadrilateral element in 2d (GooseFEM::Mesh::ElementType::Quad4), in a Cartesian coordinate system.

The different with ElementQuad4.h is that here the tensors live in 3d and are assumed plane strain.

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 Quad4

4-noded quadrilateral element in 2d (GooseFEM::Mesh::ElementType::Quad4).

class QuadraturePlanar : public GooseFEM::Element::QuadratureBaseCartesian<QuadraturePlanar>
#include <GooseFEM/ElementQuad4Planar.h>

Interpolation and quadrature.

Similar to Element::Quad4::Quadrature with the only different that quadrature point tensors are 3d (“plane strain”) while the mesh is 2d.

Fixed dimensions:

  • ndim = 2: number of dimensions.

  • tdim = 3: number of dimensions or tensor.

  • nne = 4: number of nodes per element.

Naming convention:

Public Functions

QuadraturePlanar() = default#
template<class T>
inline QuadraturePlanar(const T &x, double thick = 1.0)

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).

  • thick – out-of-plane thickness (incorporated in the element volume).

template<class T, class X, class W>
inline QuadraturePlanar(const T &x, const X &xi, const W &w, double thick = 1.0)

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].

  • thick – out-of-plane thickness (incorporated in the element volume).

Private Functions

template<class T, class R>
inline void gradN_vector_impl(const T &elemvec, R &qtensor) const#
template<class T, class R>
inline void gradN_vector_T_impl(const T &elemvec, R &qtensor) const#
template<class T, class R>
inline void symGradN_vector_impl(const T &elemvec, R &qtensor) const#
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#
inline void compute_dN_impl()#

Private Members

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

Number of dimensions for tensors.

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]

double m_thick#

out-of-plane thickness

Private Static Attributes

static constexpr size_t s_nne = 4#

Number of nodes per element.

static constexpr size_t s_ndim = 2#

Number of dimensions for nodal vectors.

static constexpr size_t s_tdim = 3#

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