Class GooseFEM::Element::QuadratureBaseCartesian#

template<class D>
class QuadratureBaseCartesian : public GooseFEM::Element::QuadratureBase<D>#

CRTP base class for interpolation and quadrature for a generic element in Cartesian coordinates.

Naming convention:

Public Types

using derived_type = D#

Underlying type.

Public Functions

template<class T>
inline void update_x(const T &x)#

Update the nodal positions.

This recomputes:

  • the shape functions,

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

  • the integration points volumes. Under the small deformations assumption this function should not be called.

Parameters

x – nodal coordinates (elemvec). Shape should match the earlier definition.

inline auto GradN() const -> const array_type::tensor<double, 4>&#

Shape function gradients (in global coordinates).

Returns

gradN stored per element, per integration point [nelem, nip, nne, ndim].

inline auto dV() const -> const array_type::tensor<double, 2>&#

Integration volume.

Returns

volume stored per element, per integration point [nelem, nip].

template<class T>
inline auto InterpQuad_vector(const T &elemvec) const -> array_type::tensor<double, 3>#

Interpolate element vector and evaluate at each quadrature point.

\( \vec{u}(\vec{x}_q) = N_i^e(\vec{x}) \vec{u}_i^e \)

Parameters

elemvec – nodal vector stored per element [nelem, nne, ndim].

Returns

qvector [nelem, nip, ndim].

template<class T, class R>
inline void interpQuad_vector(const T &elemvec, R &qvector) const#

Same as InterpQuad_vector(), but writing to preallocated return.

Parameters
template<class T>
inline auto GradN_vector(const T &elemvec) const -> array_type::tensor<double, 4>#

Element-by-element: dyadic product of the shape function gradients and a nodal vector.

Typical input: nodal displacements. Typical output: quadrature point strains. Within one element::

 for e in range(nelem):
     for q in range(nip):
         for m in range(nne):
             qtensor(e, q, i, j) += dNdx(e, q, m, i) * elemvec(e, m, j)
Note that the functions and their gradients are precomputed upon construction, or updated when calling update_x().

Parameters

elemvec – [nelem, nne, ndim]

Returns

qtensor [nelem, nip, tdim, tdim]

template<class T, class R>
inline void gradN_vector(const T &elemvec, R &qtensor) const#

Same as GradN_vector(), but writing to preallocated return.

Parameters
template<class T>
inline auto GradN_vector_T(const T &elemvec) const -> array_type::tensor<double, 4>#

The transposed output of GradN_vector().

Within one element::

 for e in range(nelem):
     for q in range(nip):
         for m in range(nne):
             qtensor(e, q, j, i) += dNdx(e, q, m, i) * elemvec(e, m, j)

Parameters

elemvec – [nelem, nne, ndim]

Returns

qtensor [nelem, nip, tdim, tdim]

template<class T, class R>
inline void gradN_vector_T(const T &elemvec, R &qtensor) const#

Same as GradN_vector_T(), but writing to preallocated return.

Parameters
template<class T>
inline auto SymGradN_vector(const T &elemvec) const -> array_type::tensor<double, 4>#

The symmetric output of GradN_vector().

Without one element::

 for e in range(nelem):
     for q in range(nip):
         for m in range(nne):
             qtensor(e, q, i, j) += 0.5 * dNdx(e, q, m, i) * elemvec(e, m, j)
             qtensor(e, q, j, i) += 0.5 * dNdx(e, q, m, i) * elemvec(e, m, j)

Parameters

elemvec – [nelem, nne, ndim]

Returns

qtensor [nelem, nip, tdim, tdim]

template<class T, class R>
inline void symGradN_vector(const T &elemvec, R &qtensor) const#

Same as SymGradN_vector(), but writing to preallocated return.

Parameters
template<class T>
inline auto Int_N_vector_dV(const T &qvector) const -> array_type::tensor<double, 3>#

Element-by-element: integral of a continuous vector-field.

\( \vec{f}_i^e = \int N_i^e(\vec{x}) \vec{f}(\vec{x}) d\Omega_e \)

which is integration numerically as follows

\( \vec{f}_i^e = \sum\limits_q N_i^e(\vec{x}_q) \vec{f}(\vec{x}_q) \)

Parameters

qvector – [nelem, nip. ndim]

Returns

elemvec [nelem, nne. ndim]

template<class T, class R>
inline void int_N_vector_dV(const T &qvector, R &elemvec) const#

Same as Int_N_vector_dV(), but writing to preallocated return.

Parameters
template<class T>
inline auto Int_N_scalar_NT_dV(const T &qscalar) const -> array_type::tensor<double, 3>#

Element-by-element: integral of the scalar product of the shape function with a scalar.

Within one one element::

 for e in range(nelem):
     for q in range(nip):
         for m in range(nne):
             for n in range(nne):
                 elemmat(e, m * ndim + i, n * ndim + i) +=
                     N(e, q, m) * qscalar(e, q) * N(e, q, n) * dV(e, q)
with i a tensor dimension. Note that the functions and their gradients are precomputed upon construction, or updated when calling update_x().

Parameters

qscalar – [nelem, nip]

Returns

elemmat [nelem, nne * ndim, nne * ndim]

template<class T, class R>
inline void int_N_scalar_NT_dV(const T &qscalar, R &elemmat) const#

Same as Int_N_scalar_NT_dV(), but writing to preallocated return.

Parameters
template<class T>
inline auto Int_gradN_dot_tensor2_dV(const T &qtensor) const -> array_type::tensor<double, 3>#

Element-by-element: integral of the dot product of the shape function gradients with a second order tensor.

Typical input: stress. Typical output: nodal force. Within one one element::

 for e in range(nelem):
     for q in range(nip):
         for m in range(nne):
             elemvec(e, m, j) += dNdx(e, q, m, i) * qtensor(e, q, i, j) * dV(e, q)
with i and j tensor dimensions. Note that the functions and their gradients are precomputed upon construction, or updated when calling update_x().

Parameters

qtensor – [nelem, nip, ndim, ndim]

Returns

elemvec [nelem, nne. ndim]

template<class T, class R>
inline void int_gradN_dot_tensor2_dV(const T &qtensor, R &elemvec) const#

Same as Int_gradN_dot_tensor2_dV(), but writing to preallocated return.

Parameters
template<class T>
inline auto Int_gradN_dot_tensor4_dot_gradNT_dV(const T &qtensor) const -> array_type::tensor<double, 3>#

Element-by-element: integral of the dot products of the shape function gradients with a fourth order tensor.

Typical input: stiffness tensor. Typical output: stiffness matrix. Within one one element::

 for e in range(nelem):
     for q in range(nip):
         for m in range(nne):
             for n in range(nne):
                 elemmat(e, m * ndim + j, n * ndim + k) +=
                     dNdx(e,q,m,i) * qtensor(e,q,i,j,k,l) * dNdx(e,q,n,l) * dV(e,q)
with i, j, k, and l tensor dimensions. Note that the functions and their gradients are precomputed upon construction, or updated when calling update_x().

Parameters

qtensor – [nelem, nip, ndim, ndim, ndim, ndim]

Returns

elemmat [nelem, nne * ndim, nne * ndim]

template<class T, class R>
inline void int_gradN_dot_tensor4_dot_gradNT_dV(const T &qtensor, R &elemmat) const#

Same as Int_gradN_dot_tensor4_dot_gradNT_dV(), but writing to preallocated return.

Parameters