File ElementQuad4.h#
Quadrature for 4-noded quadrilateral element in 2d (GooseFEM::Mesh::ElementType::Quad4), 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 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:
elemmat
: matrices stored per element, [nelem, nne * ndim, nne * ndim]elemvec
: nodal vectors stored per element, [nelem, nne, ndim]
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.
Private Functions
-
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< Quadrature >
- friend QuadratureBaseCartesian< Quadrature >
-
size_t m_tdim = 2#
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, 1> m_w#
weight of each integration point [nip]
-
array_type::tensor<double, 2> m_xi#
-
array_type::tensor<double, 2> m_N#
-
array_type::tensor<double, 2> m_vol#
-
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
], withndim = 2
.
-
inline array_type::tensor<double, 1> w()#
Integration point weights.
- Returns:
Coordinates [nip].
-
inline size_t 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
], withndim = 2
.
-
inline array_type::tensor<double, 1> w()#
Integration point weights.
- Returns:
Coordinates [nip].
-
inline size_t 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
], withndim = 2
.
-
inline array_type::tensor<double, 1> w()#
Integration point weights.
- Returns:
Coordinates [nip].
-
inline size_t nip()#
-
class Quadrature : public GooseFEM::Element::QuadratureBaseCartesian<Quadrature>
-
namespace Quad4
-
namespace Element