Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Classes | Public Types | Public Member Functions | Private Attributes | List of all members

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

Classes

struct  OutputType
 
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 >
 

Public Member Functions

 Vector ()
 
 Vector (const FEValuesBase< dim, spacedim > &fe_values_base, const unsigned int first_vector_component)
 
Vectoroperator= (const Vector< dim, spacedim > &)=delete
 
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<class InputVector >
void get_function_values (const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
 
template<class InputVector >
void get_function_values_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::value_type > &values) const
 
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
 
template<class InputVector >
void get_function_gradients_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::gradient_type > &gradients) const
 
template<class InputVector >
void get_function_symmetric_gradients (const InputVector &fe_function, std::vector< typename ProductType< symmetric_gradient_type, typename InputVector::value_type >::type > &symmetric_gradients) const
 
template<class InputVector >
void get_function_symmetric_gradients_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::symmetric_gradient_type > &symmetric_gradients) const
 
template<class InputVector >
void get_function_divergences (const InputVector &fe_function, std::vector< typename ProductType< divergence_type, typename InputVector::value_type >::type > &divergences) const
 
template<class InputVector >
void get_function_divergences_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::divergence_type > &divergences) const
 
template<class InputVector >
void get_function_curls (const InputVector &fe_function, std::vector< typename ProductType< curl_type, typename InputVector::value_type >::type > &curls) const
 
template<class InputVector >
void get_function_curls_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::curl_type > &curls) const
 
template<class InputVector >
void get_function_hessians (const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
 
template<class InputVector >
void get_function_hessians_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::hessian_type > &hessians) const
 
template<class InputVector >
void get_function_laplacians (const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
 
template<class InputVector >
void get_function_laplacians_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::laplacian_type > &laplacians) const
 
template<class InputVector >
void get_function_third_derivatives (const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
 
template<class InputVector >
void get_function_third_derivatives_from_local_dof_values (const InputVector &dof_values, std::vector< typename OutputType< typename InputVector::value_type >::third_derivative_type > &third_derivatives) const
 

Private Attributes

const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
 
const 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 module.

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 579 of file fe_values.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 587 of file fe_values.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 597 of file fe_values.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 609 of file fe_values.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 616 of file fe_values.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 624 of file fe_values.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 631 of file fe_values.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 638 of file fe_values.h.

Constructor & Destructor Documentation

◆ Vector() [1/2]

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

Default constructor. Creates an invalid object.

Definition at line 261 of file fe_values.cc.

◆ Vector() [2/2]

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 188 of file fe_values.cc.

Member Function Documentation

◆ operator=()

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.

◆ 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<class InputVector >
void Vector< dim, spacedim >::get_function_values ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< value_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 1855 of file fe_values.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< typename OutputType< typename InputVector< dim, spacedim >::value_type >::value_type > &  values 
) const

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 1886 of file fe_values.cc.

◆ get_function_gradients()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_gradients ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< gradient_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 1911 of file fe_values.cc.

◆ get_function_gradients_from_local_dof_values()

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

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 1942 of file fe_values.cc.

◆ get_function_symmetric_gradients()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_symmetric_gradients ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< symmetric_gradient_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 1967 of file fe_values.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< typename OutputType< typename InputVector< dim, spacedim >::value_type >::symmetric_gradient_type > &  symmetric_gradients 
) const

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 1998 of file fe_values.cc.

◆ get_function_divergences()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_divergences ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< divergence_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 2022 of file fe_values.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< typename OutputType< typename InputVector< dim, spacedim >::value_type >::divergence_type > &  divergences 
) const

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 2054 of file fe_values.cc.

◆ get_function_curls()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_curls ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< curl_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 2079 of file fe_values.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< typename OutputType< typename InputVector< dim, spacedim >::value_type >::curl_type > &  curls 
) const

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 2110 of file fe_values.cc.

◆ get_function_hessians()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_hessians ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< hessian_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 2135 of file fe_values.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< typename OutputType< typename InputVector< dim, spacedim >::value_type >::hessian_type > &  hessians 
) const

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 2166 of file fe_values.cc.

◆ get_function_laplacians()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_laplacians ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< value_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 2191 of file fe_values.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< typename OutputType< typename InputVector< dim, spacedim >::value_type >::laplacian_type > &  laplacians 
) const

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 2227 of file fe_values.cc.

◆ get_function_third_derivatives()

template<int dim, int spacedim>
template<class InputVector >
void Vector< dim, spacedim >::get_function_third_derivatives ( const InputVector< dim, spacedim > &  fe_function,
std::vector< typename ProductType< third_derivative_type, typename InputVector< dim, spacedim >::value_type >::type > &  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 2255 of file fe_values.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< typename OutputType< typename InputVector< dim, spacedim >::value_type >::third_derivative_type > &  third_derivatives 
) const

Same as above, but using a vector of local degree-of-freedom values.

The dof_values vector must have a length equal to number of DoFs on a cell, and each entry dof_values[i] is the value of the local DoF i. The fundamental prerequisite for the InputVector is that it must be possible to create an ArrayView from it; this is satisfied by the std::vector class.

The DoF values typically would be obtained in the following way:

Vector<double> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution, local_dof_values);

or, for a generic Number type,

std::vector<Number> local_dof_values(cell->get_fe().dofs_per_cell);
cell->get_dof_values(solution,
local_dof_values.begin(),
local_dof_values.end());

Definition at line 2286 of file fe_values.cc.

Member Data Documentation

◆ fe_values

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

A pointer to the FEValuesBase object we operate on.

Definition at line 1213 of file fe_values.h.

◆ first_vector_component

template<int dim, int spacedim = dim>
const 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 1219 of file fe_values.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 1224 of file fe_values.h.


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