Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches

#include <deal.II/fe/fe_values_views.h>

Classes

struct  ShapeFunctionData
 

Public Types

using value_type = ::Tensor<1, spacedim>
 
using gradient_type = ::Tensor<2, spacedim>
 
using symmetric_gradient_type = ::SymmetricTensor<2, spacedim>
 
using divergence_type = double
 
using curl_type = typename ::internal::CurlType<spacedim>::type
 
using hessian_type = ::Tensor<3, spacedim>
 
using third_derivative_type = ::Tensor<4, spacedim>
 
template<typename Number >
using solution_value_type = typename ProductType<Number, value_type>::type
 
template<typename Number >
using solution_gradient_type
 
template<typename Number >
using solution_symmetric_gradient_type
 
template<typename Number >
using solution_divergence_type
 
template<typename Number >
using solution_laplacian_type
 
template<typename Number >
using solution_curl_type = typename ProductType<Number, curl_type>::type
 
template<typename Number >
using solution_hessian_type
 
template<typename Number >
using solution_third_derivative_type
 

Public Member Functions

 Vector ()
 
 Vector (const FEValuesBase< dim, spacedim > &fe_values_base, const unsigned int first_vector_component)
 
 Vector (const Vector< dim, spacedim > &)=delete
 
 Vector (Vector< dim, spacedim > &&)=default
 
 ~Vector ()=default
 
Vectoroperator= (const Vector< dim, spacedim > &)=delete
 
Vectoroperator= (Vector< dim, spacedim > &&)=default
 
value_type value (const unsigned int shape_function, const unsigned int q_point) const
 
gradient_type gradient (const unsigned int shape_function, const unsigned int q_point) const
 
symmetric_gradient_type symmetric_gradient (const unsigned int shape_function, const unsigned int q_point) const
 
divergence_type divergence (const unsigned int shape_function, const unsigned int q_point) const
 
curl_type curl (const unsigned int shape_function, const unsigned int q_point) const
 
hessian_type hessian (const unsigned int shape_function, const unsigned int q_point) const
 
third_derivative_type third_derivative (const unsigned int shape_function, const unsigned int q_point) const
 
template<typename Number >
void get_function_values (const ReadVector< Number > &fe_function, std::vector< solution_value_type< Number > > &values) const
 
template<class InputVector >
void get_function_values_from_local_dof_values (const InputVector &dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
 
template<typename Number >
void get_function_gradients (const ReadVector< Number > &fe_function, std::vector< solution_gradient_type< Number > > &gradients) const
 
template<class InputVector >
void get_function_gradients_from_local_dof_values (const InputVector &dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
 
template<typename Number >
void get_function_symmetric_gradients (const ReadVector< Number > &fe_function, std::vector< solution_symmetric_gradient_type< Number > > &symmetric_gradients) const
 
template<class InputVector >
void get_function_symmetric_gradients_from_local_dof_values (const InputVector &dof_values, std::vector< solution_symmetric_gradient_type< typename InputVector::value_type > > &symmetric_gradients) const
 
template<typename Number >
void get_function_divergences (const ReadVector< Number > &fe_function, std::vector< solution_divergence_type< Number > > &divergences) const
 
template<class InputVector >
void get_function_divergences_from_local_dof_values (const InputVector &dof_values, std::vector< solution_divergence_type< typename InputVector::value_type > > &divergences) const
 
template<typename Number >
void get_function_curls (const ReadVector< Number > &fe_function, std::vector< solution_curl_type< Number > > &curls) const
 
template<class InputVector >
void get_function_curls_from_local_dof_values (const InputVector &dof_values, std::vector< solution_curl_type< typename InputVector::value_type > > &curls) const
 
template<typename Number >
void get_function_hessians (const ReadVector< Number > &fe_function, std::vector< solution_hessian_type< Number > > &hessians) const
 
template<class InputVector >
void get_function_hessians_from_local_dof_values (const InputVector &dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
 
template<typename Number >
void get_function_laplacians (const ReadVector< Number > &fe_function, std::vector< solution_laplacian_type< Number > > &laplacians) const
 
template<class InputVector >
void get_function_laplacians_from_local_dof_values (const InputVector &dof_values, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
 
template<typename Number >
void get_function_third_derivatives (const ReadVector< Number > &fe_function, std::vector< solution_third_derivative_type< Number > > &third_derivatives) const
 
template<class InputVector >
void get_function_third_derivatives_from_local_dof_values (const InputVector &dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
 

Private Attributes

SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
 
unsigned int first_vector_component
 
std::vector< ShapeFunctionDatashape_function_data
 

Detailed Description

template<int dim, int spacedim = dim>
class FEValuesViews::Vector< dim, spacedim >

A class representing a view to a set of spacedim components forming a vector part of a vector-valued finite element. Views are discussed in the Handling vector valued problems topic.

Note that in the current context, a vector is meant in the sense physics uses it: it has spacedim components that behave in specific ways under coordinate system transformations. Examples include velocity or displacement fields. This is opposed to how mathematics uses the word "vector" (and how we use this word in other contexts in the library, for example in the Vector class), where it really stands for a collection of numbers. An example of this latter use of the word could be the set of concentrations of chemical species in a flame; however, these are really just a collection of scalar variables, since they do not change if the coordinate system is rotated, unlike the components of a velocity vector, and consequently, this class should not be used for this context.

This class allows to query the value, gradient and divergence of (components of) shape functions and solutions representing vectors. The gradient of a vector \(d_{k}, 0\le k<\text{dim}\) is defined as \(S_{ij} = \frac{\partial d_{i}}{\partial x_j}, 0\le i,j<\text{dim}\).

You get an object of this type if you apply a FEValuesExtractors::Vector to an FEValues, FEFaceValues or FESubfaceValues object.

Definition at line 597 of file fe_values_views.h.

Member Typedef Documentation

◆ value_type

template<int dim, int spacedim = dim>
using FEValuesViews::Vector< dim, spacedim >::value_type = ::Tensor<1, spacedim>

An alias for the data type of values of the view this class represents. Since we deal with a set of dim components, the value type is a Tensor<1,spacedim>.

Definition at line 605 of file fe_values_views.h.

◆ gradient_type

template<int dim, int spacedim = dim>
using FEValuesViews::Vector< dim, spacedim >::gradient_type = ::Tensor<2, spacedim>

An alias for the type of gradients of the view this class represents. Here, for a set of dim components of the finite element, the gradient is a Tensor<2,spacedim>.

See the general documentation of this class for how exactly the gradient of a vector is defined.

Definition at line 615 of file fe_values_views.h.

◆ symmetric_gradient_type

template<int dim, int spacedim = dim>
using FEValuesViews::Vector< dim, spacedim >::symmetric_gradient_type = ::SymmetricTensor<2, spacedim>

An alias for the type of symmetrized gradients of the view this class represents. Here, for a set of dim components of the finite element, the symmetrized gradient is a SymmetricTensor<2,spacedim>.

The symmetric gradient of a vector field \(\mathbf v\) is defined as \(\varepsilon(\mathbf v)=\frac 12 (\nabla \mathbf v + \nabla \mathbf v^T)\).

Definition at line 627 of file fe_values_views.h.

◆ divergence_type

template<int dim, int spacedim = dim>
using FEValuesViews::Vector< dim, spacedim >::divergence_type = double

An alias for the type of the divergence of the view this class represents. Here, for a set of dim components of the finite element, the divergence of course is a scalar.

Definition at line 634 of file fe_values_views.h.

◆ curl_type

template<int dim, int spacedim = dim>
using FEValuesViews::Vector< dim, spacedim >::curl_type = typename ::internal::CurlType<spacedim>::type

An alias for the type of the curl of the view this class represents. Here, for a set of spacedim=2 components of the finite element, the curl is a Tensor<1, 1>. For spacedim=3 it is a Tensor<1, dim>.

Definition at line 642 of file fe_values_views.h.

◆ hessian_type

template<int dim, int spacedim = dim>
using FEValuesViews::Vector< dim, spacedim >::hessian_type = ::Tensor<3, spacedim>

An alias for the type of second derivatives of the view this class represents. Here, for a set of dim components of the finite element, the Hessian is a Tensor<3,dim>.

Definition at line 649 of file fe_values_views.h.

◆ third_derivative_type

template<int dim, int spacedim = dim>
using FEValuesViews::Vector< dim, spacedim >::third_derivative_type = ::Tensor<4, spacedim>

An alias for the type of third derivatives of the view this class represents. Here, for a set of dim components of the finite element, the third derivative is a Tensor<4,dim>.

Definition at line 656 of file fe_values_views.h.

◆ solution_value_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_value_type = typename ProductType<Number, value_type>::type

An alias for the data type of the product of a Number and the values of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 665 of file fe_values_views.h.

◆ solution_gradient_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_gradient_type
Initial value:
typename internal::ProductTypeImpl< std::decay_t< T >, std::decay_t< U > >::type type

An alias for the data type of the product of a Number and the gradients of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 674 of file fe_values_views.h.

◆ solution_symmetric_gradient_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_symmetric_gradient_type
Initial value:

An alias for the data type of the product of a Number and the symmetric gradients of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 684 of file fe_values_views.h.

◆ solution_divergence_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_divergence_type
Initial value:

An alias for the data type of the product of a Number and the divergences of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 694 of file fe_values_views.h.

◆ solution_laplacian_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_laplacian_type
Initial value:

An alias for the data type of the product of a Number and the laplacians of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 704 of file fe_values_views.h.

◆ solution_curl_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_curl_type = typename ProductType<Number, curl_type>::type

An alias for the data type of the product of a Number and the curls of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 714 of file fe_values_views.h.

◆ solution_hessian_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_hessian_type
Initial value:

An alias for the data type of the product of a Number and the hessians of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 723 of file fe_values_views.h.

◆ solution_third_derivative_type

template<int dim, int spacedim = dim>
template<typename Number >
using FEValuesViews::Vector< dim, spacedim >::solution_third_derivative_type
Initial value:

An alias for the data type of the product of a Number and the third derivatives of the view this class provides. This is the data type of vector components of a finite element field whose degrees of freedom are described by a vector with elements of type Number.

Definition at line 733 of file fe_values_views.h.

Constructor & Destructor Documentation

◆ Vector() [1/4]

template<int dim, int spacedim>
Vector< dim, spacedim >::Vector ( )

Default constructor. Creates an invalid object.

Definition at line 184 of file fe_values_views.cc.

◆ Vector() [2/4]

template<int dim, int spacedim>
Vector< dim, spacedim >::Vector ( const FEValuesBase< dim, spacedim > & fe_values_base,
const unsigned int first_vector_component )

Constructor for an object that represents dim components of a FEValuesBase object (or of one of the classes derived from FEValuesBase), representing a vector-valued variable.

The second argument denotes the index of the first component of the selected vector.

Definition at line 114 of file fe_values_views.cc.

◆ Vector() [3/4]

template<int dim, int spacedim = dim>
FEValuesViews::Vector< dim, spacedim >::Vector ( const Vector< dim, spacedim > & )
delete

Copy constructor. This is not a lightweight object so we don't allow copying and generate a compile-time error if this function is called.

◆ Vector() [4/4]

template<int dim, int spacedim = dim>
FEValuesViews::Vector< dim, spacedim >::Vector ( Vector< dim, spacedim > && )
default

Move constructor.

◆ ~Vector()

template<int dim, int spacedim = dim>
FEValuesViews::Vector< dim, spacedim >::~Vector ( )
default

Destructor.

Member Function Documentation

◆ operator=() [1/2]

template<int dim, int spacedim = dim>
Vector & FEValuesViews::Vector< dim, spacedim >::operator= ( const Vector< dim, spacedim > & )
delete

Copy operator. This is not a lightweight object so we don't allow copying and generate a compile-time error if this function is called.

◆ operator=() [2/2]

template<int dim, int spacedim = dim>
Vector & FEValuesViews::Vector< dim, spacedim >::operator= ( Vector< dim, spacedim > && )
default

Move assignment operator.

◆ value()

template<int dim, int spacedim = dim>
value_type FEValuesViews::Vector< dim, spacedim >::value ( const unsigned int shape_function,
const unsigned int q_point ) const

Return the value of the vector components selected by this view, for the shape function and quadrature point selected by the arguments. Here, since the view represents a vector-valued part of the FEValues object with dim components, the return type is a tensor of rank 1 with dim components.

Parameters
shape_functionNumber of the shape function to be evaluated. Note that this number runs from zero to dofs_per_cell, even in the case of an FEFaceValues or FESubfaceValues object.
q_pointNumber of the quadrature point at which function is to be evaluated.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ gradient()

template<int dim, int spacedim = dim>
gradient_type FEValuesViews::Vector< dim, spacedim >::gradient ( const unsigned int shape_function,
const unsigned int q_point ) const

Return the gradient (a tensor of rank 2) of the vector component selected by this view, for the shape function and quadrature point selected by the arguments.

See the general documentation of this class for how exactly the gradient of a vector is defined.

Note
The meaning of the arguments is as documented for the value() function.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ symmetric_gradient()

template<int dim, int spacedim = dim>
symmetric_gradient_type FEValuesViews::Vector< dim, spacedim >::symmetric_gradient ( const unsigned int shape_function,
const unsigned int q_point ) const

Return the symmetric gradient (a symmetric tensor of rank 2) of the vector component selected by this view, for the shape function and quadrature point selected by the arguments.

The symmetric gradient is defined as \(\frac 12 [(\nabla \phi_i(x_q)) + (\nabla \phi_i(x_q))^T]\), where \(\phi_i\) represents the dim components selected from the FEValuesBase object, and \(x_q\) is the location of the \(q\)-th quadrature point.

Note
The meaning of the arguments is as documented for the value() function.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ divergence()

template<int dim, int spacedim = dim>
divergence_type FEValuesViews::Vector< dim, spacedim >::divergence ( const unsigned int shape_function,
const unsigned int q_point ) const

Return the scalar divergence of the vector components selected by this view, for the shape function and quadrature point selected by the arguments.

Note
The meaning of the arguments is as documented for the value() function.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ curl()

template<int dim, int spacedim = dim>
curl_type FEValuesViews::Vector< dim, spacedim >::curl ( const unsigned int shape_function,
const unsigned int q_point ) const

Return the vector curl of the vector components selected by this view, for the shape function and quadrature point selected by the arguments. For 1d this function does not make any sense. Thus it is not implemented for spacedim=1. In 2d the curl is defined as

\begin{equation*} \operatorname{curl}(u) \dealcoloneq \frac{du_2}{dx} -\frac{du_1}{dy}, \end{equation*}

whereas in 3d it is given by

\begin{equation*} \operatorname{curl}(u) \dealcoloneq \left( \begin{array}{c} \frac{du_3}{dy}-\frac{du_2}{dz}\\ \frac{du_1}{dz}-\frac{du_3}{dx}\\ \frac{du_2}{dx}-\frac{du_1}{dy} \end{array} \right). \end{equation*}

Note
The meaning of the arguments is as documented for the value() function.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ hessian()

template<int dim, int spacedim = dim>
hessian_type FEValuesViews::Vector< dim, spacedim >::hessian ( const unsigned int shape_function,
const unsigned int q_point ) const

Return the Hessian (the tensor of rank 2 of all second derivatives) of the vector components selected by this view, for the shape function and quadrature point selected by the arguments.

Note
The meaning of the arguments is as documented for the value() function.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ third_derivative()

template<int dim, int spacedim = dim>
third_derivative_type FEValuesViews::Vector< dim, spacedim >::third_derivative ( const unsigned int shape_function,
const unsigned int q_point ) const

Return the tensor of rank 3 of all third derivatives of the vector components selected by this view, for the shape function and quadrature point selected by the arguments.

Note
The meaning of the arguments is as documented for the value() function.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_function_values()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_values ( const ReadVector< Number > & fe_function,
std::vector< solution_value_type< Number > > & values ) const

Return the values of the selected vector components of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

This function is the equivalent of the FEValuesBase::get_function_values function but it only works on the selected vector components.

The data type stored by the output vector must be what you get when you multiply the values of shape functions (i.e., value_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 622 of file fe_values_views.cc.

◆ get_function_values_from_local_dof_values()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_values_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_value_type< typename InputVector< dim, spacedim >::value_type > > & values ) const

Same as above, but using a vector of local degree-of-freedom values. In other words, instead of extracting the nodal values of the degrees of freedom located on the current cell from a global vector associated with a DoFHandler object (as the function above does), this function instead takes these local nodal values through its first argument. A typical way to obtain such a vector is by calling code such as

cell->get_dof_values (dof_values, local_dof_values);

(See DoFCellAccessor::get_dof_values() for more information on this function.) The point of the current function is then that one could modify these local values first, for example by applying a limiter or by ensuring that all nodal values are positive, before evaluating the finite element field that corresponds to these local values on the current cell. Another application is where one wants to postprocess the solution on a cell into a different finite element space on every cell, without actually creating a corresponding DoFHandler – in that case, all one would compute is a local representation of that postprocessed function, characterized by its nodal values; this function then allows the evaluation of that representation at quadrature points.

Parameters
[in]dof_valuesA vector of local nodal values. This vector must have a length equal to number of DoFs on the current cell, and must be ordered in the same order as degrees of freedom are numbered on the reference cell.
[out]valuesA vector of values of the given finite element field, at the quadrature points on the current object.
Template Parameters
InputVectorThe InputVector type must allow creation of an ArrayView object from it; this is satisfied by the std::vector class, among others.

Definition at line 648 of file fe_values_views.cc.

◆ get_function_gradients()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_gradients ( const ReadVector< Number > & fe_function,
std::vector< solution_gradient_type< Number > > & gradients ) const

Return the gradients of the selected vector components of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

This function is the equivalent of the FEValuesBase::get_function_gradients function but it only works on the selected vector components.

The data type stored by the output vector must be what you get when you multiply the gradients of shape functions (i.e., gradient_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 672 of file fe_values_views.cc.

◆ get_function_gradients_from_local_dof_values()

template<int dim, int spacedim>
template<typename InputVector >
void Vector< dim, spacedim >::get_function_gradients_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_gradient_type< typename InputVector< dim, spacedim >::value_type > > & gradients ) const

This function relates to get_function_gradients() in the same way as get_function_values_from_local_dof_values() relates to get_function_values(). See the documentation of get_function_values_from_local_dof_values() for more information.

Definition at line 700 of file fe_values_views.cc.

◆ get_function_symmetric_gradients()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_symmetric_gradients ( const ReadVector< Number > & fe_function,
std::vector< solution_symmetric_gradient_type< Number > > & symmetric_gradients ) const

Return the symmetrized gradients of the selected vector components of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

The symmetric gradient of a vector field \(\mathbf v\) is defined as \(\varepsilon(\mathbf v)=\frac 12 (\nabla \mathbf v + \nabla \mathbf v^T)\).

Note
There is no equivalent function such as FEValuesBase::get_function_symmetric_gradients in the FEValues classes but the information can be obtained from FEValuesBase::get_function_gradients, of course.

The data type stored by the output vector must be what you get when you multiply the symmetric gradients of shape functions (i.e., symmetric_gradient_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 724 of file fe_values_views.cc.

◆ get_function_symmetric_gradients_from_local_dof_values()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_symmetric_gradients_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_symmetric_gradient_type< typename InputVector< dim, spacedim >::value_type > > & symmetric_gradients ) const

This function relates to get_function_symmetric_gradients() in the same way as get_function_values_from_local_dof_values() relates to get_function_values(). See the documentation of get_function_values_from_local_dof_values() for more information.

Definition at line 753 of file fe_values_views.cc.

◆ get_function_divergences()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_divergences ( const ReadVector< Number > & fe_function,
std::vector< solution_divergence_type< Number > > & divergences ) const

Return the divergence of the selected vector components of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

There is no equivalent function such as FEValuesBase::get_function_divergences in the FEValues classes but the information can be obtained from FEValuesBase::get_function_gradients, of course.

The data type stored by the output vector must be what you get when you multiply the divergences of shape functions (i.e., divergence_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 778 of file fe_values_views.cc.

◆ get_function_divergences_from_local_dof_values()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_divergences_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_divergence_type< typename InputVector< dim, spacedim >::value_type > > & divergences ) const

This function relates to get_function_divergences() in the same way as get_function_values_from_local_dof_values() relates to get_function_values(). See the documentation of get_function_values_from_local_dof_values() for more information.

Definition at line 807 of file fe_values_views.cc.

◆ get_function_curls()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_curls ( const ReadVector< Number > & fe_function,
std::vector< solution_curl_type< Number > > & curls ) const

Return the curl of the selected vector components of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

There is no equivalent function such as FEValuesBase::get_function_curls in the FEValues classes but the information can be obtained from FEValuesBase::get_function_gradients, of course.

The data type stored by the output vector must be what you get when you multiply the curls of shape functions (i.e., curl_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 831 of file fe_values_views.cc.

◆ get_function_curls_from_local_dof_values()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_curls_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_curl_type< typename InputVector< dim, spacedim >::value_type > > & curls ) const

This function relates to get_function_curls() in the same way as get_function_values_from_local_dof_values() relates to get_function_values(). See the documentation of get_function_values_from_local_dof_values() for more information.

Definition at line 859 of file fe_values_views.cc.

◆ get_function_hessians()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_hessians ( const ReadVector< Number > & fe_function,
std::vector< solution_hessian_type< Number > > & hessians ) const

Return the Hessians of the selected vector components of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

This function is the equivalent of the FEValuesBase::get_function_hessians function but it only works on the selected vector components.

The data type stored by the output vector must be what you get when you multiply the Hessians of shape functions (i.e., hessian_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 883 of file fe_values_views.cc.

◆ get_function_hessians_from_local_dof_values()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_hessians_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_hessian_type< typename InputVector< dim, spacedim >::value_type > > & hessians ) const

This function relates to get_function_hessians() in the same way as get_function_values_from_local_dof_values() relates to get_function_values(). See the documentation of get_function_values_from_local_dof_values() for more information.

Definition at line 911 of file fe_values_views.cc.

◆ get_function_laplacians()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_laplacians ( const ReadVector< Number > & fe_function,
std::vector< solution_laplacian_type< Number > > & laplacians ) const

Return the Laplacians of the selected vector components of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called. The Laplacians are the trace of the Hessians.

This function is the equivalent of the FEValuesBase::get_function_laplacians function but it only works on the selected vector components.

The data type stored by the output vector must be what you get when you multiply the Laplacians of shape functions (i.e., laplacian_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 935 of file fe_values_views.cc.

◆ get_function_laplacians_from_local_dof_values()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_laplacians_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_laplacian_type< typename InputVector< dim, spacedim >::value_type > > & laplacians ) const

This function relates to get_function_laplacians() in the same way as get_function_values_from_local_dof_values() relates to get_function_values(). See the documentation of get_function_values_from_local_dof_values() for more information.

Definition at line 968 of file fe_values_views.cc.

◆ get_function_third_derivatives()

template<int dim, int spacedim>
template<typename Number >
void Vector< dim, spacedim >::get_function_third_derivatives ( const ReadVector< Number > & fe_function,
std::vector< solution_third_derivative_type< Number > > & third_derivatives ) const

Return the third derivatives of the selected scalar component of the finite element function characterized by fe_function at the quadrature points of the cell, face or subface selected the last time the reinit function of the FEValues object was called.

This function is the equivalent of the FEValuesBase::get_function_third_derivatives function but it only works on the selected scalar component.

The data type stored by the output vector must be what you get when you multiply the third derivatives of shape functions (i.e., third_derivative_type) times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_third_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 995 of file fe_values_views.cc.

◆ get_function_third_derivatives_from_local_dof_values()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_third_derivatives_from_local_dof_values ( const InputVector< dim, spacedim > & dof_values,
std::vector< solution_third_derivative_type< typename InputVector< dim, spacedim >::value_type > > & third_derivatives ) const

This function relates to get_function_third_derivatives() in the same way as get_function_values_from_local_dof_values() relates to get_function_values(). See the documentation of get_function_values_from_local_dof_values() for more information.

Definition at line 1024 of file fe_values_views.cc.

Member Data Documentation

◆ fe_values

template<int dim, int spacedim = dim>
SmartPointer<const FEValuesBase<dim, spacedim> > FEValuesViews::Vector< dim, spacedim >::fe_values
private

A pointer to the FEValuesBase object we operate on.

Definition at line 1274 of file fe_values_views.h.

◆ first_vector_component

template<int dim, int spacedim = dim>
unsigned int FEValuesViews::Vector< dim, spacedim >::first_vector_component
private

The first component of the vector this view represents of the FEValuesBase object.

Definition at line 1280 of file fe_values_views.h.

◆ shape_function_data

template<int dim, int spacedim = dim>
std::vector<ShapeFunctionData> FEValuesViews::Vector< dim, spacedim >::shape_function_data
private

Store the data about shape functions.

Definition at line 1285 of file fe_values_views.h.


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