Namespace GooseFEM::Element::Quad4#

namespace Quad4

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

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

Interpolation and quadrature.

Fixed dimensions:

  • ndim = 2: number of dimensions.

  • nne = 4: number of nodes per element.

Naming convention:

Public Functions

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

class QuadratureAxisymmetric : public GooseFEM::Element::QuadratureBaseCartesian<QuadratureAxisymmetric>
#include <GooseFEM/ElementQuad4Axisymmetric.h>

Interpolation and quadrature.

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

template<class T>
inline QuadratureAxisymmetric(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 QuadratureAxisymmetric(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].

inline const array_type::tensor<double, 6> &B() const

Get the B-matrix (shape function gradients) (in global coordinates).

Note that the functions and their gradients are precomputed upon construction, or updated when calling update_x().

Returns

B matrix stored per element, per integration point [nelem, nne, tdim, tdim, tdim]

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

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

namespace Gauss

Gauss quadrature: quadrature points such that integration is exact for this bi-linear element:

 + ----------- +
 |             |
 |   3     2   |
 |             |
 |   0     1   |
 |             |
 + ----------- +

Functions

inline size_t nip()

Number of integration points:

 nip = nne = 4

Returns

unsigned int

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

Integration point coordinates (local coordinates).

Returns

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

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

Integration point weights.

Returns

Coordinates [nip].

namespace MidPoint

midpoint quadrature: quadrature points in the middle of the element::

 + ------- +
 |         |
 |    0    |
 |         |
 + ------- +

Functions

inline size_t nip()

Number of integration points::

 nip = 1

Returns

unsigned int

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

Integration point coordinates (local coordinates).

Returns

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

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:

 3 -- 2
 |    |
 0 -- 1

Functions

inline size_t nip()

Number of integration points::

 nip = nne = 4

Returns

unsigned int

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

Integration point coordinates (local coordinates).

Returns

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

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

Integration point weights.

Returns

Coordinates [nip].