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:
elemmat
: matrices stored per element, [nelem, nne * ndim, nne * ndim]elemvec
: nodal vectors stored per element, [nelem, nne, ndim]
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.
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< 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, 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#
-
double m_thick#
out-of-plane thickness
-
class QuadraturePlanar : public GooseFEM::Element::QuadratureBaseCartesian<QuadraturePlanar>
-
namespace Quad4
-
namespace Element