GooseFEM 1.4.1.dev2+g78f16df
Loading...
Searching...
No Matches
GooseFEM::Element::QuadratureBaseCartesian< D > Class Template Reference

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

#include <GooseFEM/Element.h>

Inheritance diagram for GooseFEM::Element::QuadratureBaseCartesian< D >:
GooseFEM::Element::QuadratureBase< D >

Public Types

using derived_type = D
 Underlying type.
 
- Public Types inherited from GooseFEM::Element::QuadratureBase< D >
using derived_type = D
 Underlying type.
 

Public Member Functions

template<class T >
void update_x (const T &x)
 Update the nodal positions.
 
auto GradN () const -> const array_type::tensor< double, 4 > &
 Shape function gradients (in global coordinates).
 
auto dV () const -> const array_type::tensor< double, 2 > &
 Integration volume.
 
template<class T >
auto InterpQuad_vector (const T &elemvec) const -> array_type::tensor< double, 3 >
 Interpolate element vector and evaluate at each quadrature point.
 
template<class T , class R >
void interpQuad_vector (const T &elemvec, R &qvector) const
 Same as InterpQuad_vector(), but writing to preallocated return.
 
template<class T >
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.
 
template<class T , class R >
void gradN_vector (const T &elemvec, R &qtensor) const
 Same as GradN_vector(), but writing to preallocated return.
 
template<class T >
auto GradN_vector_T (const T &elemvec) const -> array_type::tensor< double, 4 >
 The transposed output of GradN_vector().
 
template<class T , class R >
void gradN_vector_T (const T &elemvec, R &qtensor) const
 Same as GradN_vector_T(), but writing to preallocated return.
 
template<class T >
auto SymGradN_vector (const T &elemvec) const -> array_type::tensor< double, 4 >
 The symmetric output of GradN_vector().
 
template<class T , class R >
void symGradN_vector (const T &elemvec, R &qtensor) const
 Same as SymGradN_vector(), but writing to preallocated return.
 
template<class T >
auto Int_N_vector_dV (const T &qvector) const -> array_type::tensor< double, 3 >
 Element-by-element: integral of a continuous vector-field.
 
template<class T , class R >
void int_N_vector_dV (const T &qvector, R &elemvec) const
 Same as Int_N_vector_dV(), but writing to preallocated return.
 
template<class T >
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.
 
template<class T , class R >
void int_N_scalar_NT_dV (const T &qscalar, R &elemmat) const
 Same as Int_N_scalar_NT_dV(), but writing to preallocated return.
 
template<class T >
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.
 
template<class T , class R >
void int_gradN_dot_tensor2_dV (const T &qtensor, R &elemvec) const
 Same as Int_gradN_dot_tensor2_dV(), but writing to preallocated return.
 
template<class T >
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.
 
template<class T , class R >
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.
 
- Public Member Functions inherited from GooseFEM::Element::QuadratureBase< D >
auto nelem () const
 Number of elements.
 
auto nne () const
 Number of nodes per element.
 
auto ndim () const
 Number of dimensions for node vectors.
 
auto tdim () const
 Number of dimensions for integration point tensors.
 
auto nip () const
 Number of integration points.
 
template<class T , class R >
void asTensor (const T &arg, R &ret) const
 Convert "qscalar" to "qtensor" of certain rank.
 
template<size_t rank, class T >
auto AsTensor (const T &arg) const
 Convert "qscalar" to "qtensor" of certain rank.
 
template<class T >
auto AsTensor (size_t rank, const T &arg) const
 Convert "qscalar" to "qtensor" of certain rank.
 
auto shape_elemvec () const -> std::array< size_t, 3 >
 Get the shape of an "elemvec".
 
auto shape_elemvec (size_t arg) const -> std::array< size_t, 3 >
 Get the shape of an "elemvec".
 
auto shape_elemmat () const -> std::array< size_t, 3 >
 Get the shape of an "elemmat".
 
template<size_t rank = 0>
auto shape_qtensor () const -> std::array< size_t, 2+rank >
 Get the shape of a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
auto shape_qtensor (size_t rank) const -> std::vector< size_t >
 Get the shape of a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
template<size_t trank>
auto shape_qtensor (size_t rank, size_t arg) const -> std::array< size_t, 2+trank >
 Get the shape of a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
auto shape_qtensor (size_t rank, size_t arg) const -> std::vector< size_t >
 Get the shape of a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
auto shape_qscalar () const -> std::array< size_t, 2 >
 Get the shape of a "qscalar" (a "qtensor" of rank 0)
 
auto shape_qvector () const -> std::array< size_t, 3 >
 Get the shape of a "qvector" (a "qtensor" of rank 1)
 
auto shape_qvector (size_t arg) const -> std::array< size_t, 3 >
 Get the shape of a "qvector" (a "qtensor" of rank 1)
 
template<class R >
auto allocate_elemvec () const
 Get an allocated array_type::tensor to store a "elemvec".
 
template<class R >
auto allocate_elemvec (R val) const
 Get an allocated and initialised xt::xarray to store a "elemvec".
 
template<class R >
auto allocate_elemmat () const
 Get an allocated array_type::tensor to store a "elemmat".
 
template<class R >
auto allocate_elemmat (R val) const
 Get an allocated and initialised xt::xarray to store a "elemmat".
 
template<size_t rank = 0, class R >
auto allocate_qtensor () const
 Get an allocated array_type::tensor to store a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
template<size_t rank = 0, class R >
auto allocate_qtensor (R val) const
 Get an allocated and initialised array_type::tensor to store a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
template<class R >
auto allocate_qtensor (size_t rank) const
 Get an allocated xt::xarray to store a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
template<class R >
auto allocate_qtensor (size_t rank, R val) const
 Get an allocated and initialised xt::xarray to store a "qtensor" of a certain rank (0 = scalar, 1, vector, 2 = 2nd-order tensor, etc.).
 
template<class R >
auto allocate_qscalar () const
 Get an allocated array_type::tensor to store a "qscalar" (a "qtensor" of rank 0).
 
template<class R >
auto allocate_qscalar (R val) const
 Get an allocated and initialised xt::xarray to store a "qscalar" (a "qtensor" of rank 0).
 

Protected Member Functions

void compute_dN ()
 Update the shape function gradients (called when the nodal positions are updated).
 

Friends

class QuadratureBase< D >
 

Detailed Description

template<class D>
class GooseFEM::Element::QuadratureBaseCartesian< D >

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

Naming convention:

Definition at line 546 of file Element.h.

Member Typedef Documentation

◆ derived_type

template<class D >
using GooseFEM::Element::QuadratureBaseCartesian< D >::derived_type = D

Underlying type.

Definition at line 551 of file Element.h.

Member Function Documentation

◆ compute_dN()

template<class D >
void GooseFEM::Element::QuadratureBaseCartesian< D >::compute_dN ( )
inlineprotected

Update the shape function gradients (called when the nodal positions are updated).

Definition at line 874 of file Element.h.

◆ dV()

Integration volume.

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

Definition at line 584 of file Element.h.

◆ GradN()

template<class D >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::GradN ( ) const -> const array_type::tensor<double, 4>&
inline

Shape function gradients (in global coordinates).

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

Definition at line 575 of file Element.h.

◆ GradN_vector()

template<class D >
template<class T >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::GradN_vector ( const T & elemvec) const -> array_type::tensor<double, 4>
inline

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]

Definition at line 635 of file Element.h.

◆ gradN_vector()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::gradN_vector ( const T & elemvec,
R & qtensor ) const
inline

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

Parameters
elemvec[nelem, nne, ndim]
qtensoroverwritten [nelem, nip, tdim, tdim]

Definition at line 649 of file Element.h.

◆ GradN_vector_T()

template<class D >
template<class T >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::GradN_vector_T ( const T & elemvec) const -> array_type::tensor<double, 4>
inline

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]

Definition at line 667 of file Element.h.

◆ gradN_vector_T()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::gradN_vector_T ( const T & elemvec,
R & qtensor ) const
inline

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

Parameters
elemvec[nelem, nne, ndim]
qtensoroverwritten [nelem, nip, tdim, tdim]

Definition at line 681 of file Element.h.

◆ Int_gradN_dot_tensor2_dV()

template<class D >
template<class T >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::Int_gradN_dot_tensor2_dV ( const T & qtensor) const -> array_type::tensor<double, 3>
inline

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]

Definition at line 808 of file Element.h.

◆ int_gradN_dot_tensor2_dV()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::int_gradN_dot_tensor2_dV ( const T & qtensor,
R & elemvec ) const
inline

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

Parameters
qtensor[nelem, nip, ndim, ndim]
elemvecoverwritten [nelem, nne. ndim]

Definition at line 822 of file Element.h.

◆ Int_gradN_dot_tensor4_dot_gradNT_dV()

template<class D >
template<class T >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::Int_gradN_dot_tensor4_dot_gradNT_dV ( const T & qtensor) const -> array_type::tensor<double, 3>
inline

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]

Definition at line 850 of file Element.h.

◆ int_gradN_dot_tensor4_dot_gradNT_dV()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::int_gradN_dot_tensor4_dot_gradNT_dV ( const T & qtensor,
R & elemmat ) const
inline

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

Parameters
qtensor[nelem, nip, ndim, ndim, ndim, ndim]
elemmatoverwritten [nelem, nne * ndim, nne * ndim]

Definition at line 865 of file Element.h.

◆ Int_N_scalar_NT_dV()

template<class D >
template<class T >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::Int_N_scalar_NT_dV ( const T & qscalar) const -> array_type::tensor<double, 3>
inline

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]

Definition at line 771 of file Element.h.

◆ int_N_scalar_NT_dV()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::int_N_scalar_NT_dV ( const T & qscalar,
R & elemmat ) const
inline

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

Parameters
qscalar[nelem, nip]
elemmatoverwritten [nelem, nne * ndim, nne * ndim]

Definition at line 785 of file Element.h.

◆ Int_N_vector_dV()

template<class D >
template<class T >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::Int_N_vector_dV ( const T & qvector) const -> array_type::tensor<double, 3>
inline

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]

Definition at line 732 of file Element.h.

◆ int_N_vector_dV()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::int_N_vector_dV ( const T & qvector,
R & elemvec ) const
inline

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

Parameters
qvector[nelem, nip. ndim]
elemvecoverwritten [nelem, nne. ndim]

Definition at line 747 of file Element.h.

◆ InterpQuad_vector()

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

Interpolate element vector and evaluate at each quadrature point.

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

Parameters
elemvecnodal vector stored per element [nelem, nne, ndim].
Returns
qvector [nelem, nip, ndim].

Definition at line 598 of file Element.h.

◆ interpQuad_vector()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::interpQuad_vector ( const T & elemvec,
R & qvector ) const
inline

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

Parameters
elemvecnodal vector stored per element [nelem, nne, ndim].
qvector[nelem, nip, ndim].

Definition at line 613 of file Element.h.

◆ SymGradN_vector()

template<class D >
template<class T >
auto GooseFEM::Element::QuadratureBaseCartesian< D >::SymGradN_vector ( const T & elemvec) const -> array_type::tensor<double, 4>
inline

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]

Definition at line 700 of file Element.h.

◆ symGradN_vector()

template<class D >
template<class T , class R >
void GooseFEM::Element::QuadratureBaseCartesian< D >::symGradN_vector ( const T & elemvec,
R & qtensor ) const
inline

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

Parameters
elemvec[nelem, nne, ndim]
qtensoroverwritten [nelem, nip, tdim, tdim]

Definition at line 714 of file Element.h.

◆ update_x()

template<class D >
template<class T >
void GooseFEM::Element::QuadratureBaseCartesian< D >::update_x ( const T & x)
inline

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
xnodal coordinates (elemvec). Shape should match the earlier definition.

Definition at line 564 of file Element.h.

Friends And Related Symbol Documentation

◆ QuadratureBase< D >

template<class D >
friend class QuadratureBase< D >
friend

Definition at line 885 of file Element.h.


The documentation for this class was generated from the following file: