Reference documentation for deal.II version 9.3.3
|
#include <deal.II/fe/fe_values.h>
Classes | |
struct | OutputType |
struct | ShapeFunctionData |
Public Types | |
using | value_type = double |
using | gradient_type = ::Tensor< 1, spacedim > |
using | hessian_type = ::Tensor< 2, spacedim > |
using | third_derivative_type = ::Tensor< 3, spacedim > |
template<typename Number > | |
using | solution_value_type = typename ProductType< Number, value_type >::type |
template<typename Number > | |
using | solution_gradient_type = typename ProductType< Number, gradient_type >::type |
template<typename Number > | |
using | solution_laplacian_type = typename ProductType< Number, value_type >::type |
template<typename Number > | |
using | solution_hessian_type = typename ProductType< Number, hessian_type >::type |
template<typename Number > | |
using | solution_third_derivative_type = typename ProductType< Number, third_derivative_type >::type |
Public Member Functions | |
Scalar () | |
Scalar (const FEValuesBase< dim, spacedim > &fe_values_base, const unsigned int component) | |
Scalar (const Scalar< dim, spacedim > &)=delete | |
Scalar (Scalar< dim, spacedim > &&)=default | |
~Scalar ()=default | |
Scalar & | operator= (const Scalar< dim, spacedim > &)=delete |
Scalar & | operator= (Scalar< dim, spacedim > &&) noexcept=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 |
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< solution_value_type< typename InputVector::value_type > > &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<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, std::vector< solution_gradient_type< typename InputVector::value_type > > &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<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, std::vector< solution_hessian_type< typename InputVector::value_type > > &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<class InputVector > | |
void | get_function_laplacians (const InputVector &fe_function, std::vector< solution_laplacian_type< typename InputVector::value_type > > &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<class InputVector > | |
void | get_function_third_derivatives (const InputVector &fe_function, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &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 | |
const SmartPointer< const FEValuesBase< dim, spacedim > > | fe_values |
const unsigned int | component |
std::vector< ShapeFunctionData > | shape_function_data |
A class representing a view to a single scalar component of a possibly vector-valued finite element. Views are discussed in the Handling vector valued problems module.
You get an object of this type if you apply a FEValuesExtractors::Scalar to an FEValues, FEFaceValues or FESubfaceValues object.
Definition at line 146 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::value_type = double |
An alias for the data type of values of the view this class represents. Since we deal with a single components, the value type is a scalar double.
Definition at line 154 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::gradient_type = ::Tensor<1, spacedim> |
An alias for the type of gradients of the view this class represents. Here, for a scalar component of the finite element, the gradient is a Tensor<1,dim>
.
Definition at line 161 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::hessian_type = ::Tensor<2, spacedim> |
An alias for the type of second derivatives of the view this class represents. Here, for a scalar component of the finite element, the Hessian is a Tensor<2,dim>
.
Definition at line 168 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::third_derivative_type = ::Tensor<3, spacedim> |
An alias for the type of third derivatives of the view this class represents. Here, for a scalar component of the finite element, the Third derivative is a Tensor<3,dim>
.
Definition at line 175 of file fe_values.h.
using FEValuesViews::Scalar< 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 scalar components of a finite element field whose degrees of freedom are described by a vector with elements of type Number
.
Definition at line 184 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::solution_gradient_type = typename ProductType<Number, gradient_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 scalar components of a finite element field whose degrees of freedom are described by a vector with elements of type Number
.
Definition at line 193 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::solution_laplacian_type = typename ProductType<Number, value_type>::type |
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 scalar components of a finite element field whose degrees of freedom are described by a vector with elements of type Number
.
Definition at line 203 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::solution_hessian_type = typename ProductType<Number, hessian_type>::type |
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 scalar components of a finite element field whose degrees of freedom are described by a vector with elements of type Number
.
Definition at line 213 of file fe_values.h.
using FEValuesViews::Scalar< dim, spacedim >::solution_third_derivative_type = typename ProductType<Number, third_derivative_type>::type |
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 scalar components of a finite element field whose degrees of freedom are described by a vector with elements of type Number
.
Definition at line 223 of file fe_values.h.
FEValuesViews::Scalar< dim, spacedim >::Scalar |
Default constructor. Creates an invalid object.
Definition at line 182 of file fe_values.cc.
FEValuesViews::Scalar< dim, spacedim >::Scalar | ( | const FEValuesBase< dim, spacedim > & | fe_values_base, |
const unsigned int | component | ||
) |
Constructor for an object that represents a single scalar component of a FEValuesBase object (or of one of the classes derived from FEValuesBase).
Definition at line 145 of file fe_values.cc.
|
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.
|
default |
Move constructor.
|
default |
Destructor.
|
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.
|
defaultnoexcept |
Move assignment operator.
value_type FEValuesViews::Scalar< dim, spacedim >::value | ( | const unsigned int | shape_function, |
const unsigned int | q_point | ||
) | const |
Return the value of the vector component selected by this view, for the shape function and quadrature point selected by the arguments.
shape_function | Number 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_point | Number of the quadrature point at which function is to be evaluated. |
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_type FEValuesViews::Scalar< dim, spacedim >::gradient | ( | const unsigned int | shape_function, |
const unsigned int | q_point | ||
) | const |
Return the gradient (a tensor of rank 1) of the vector component selected by this view, for the shape function and quadrature point selected by the arguments.
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_type FEValuesViews::Scalar< 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 component selected by this view, for the shape function and quadrature point selected by the arguments.
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_type FEValuesViews::Scalar< 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 component selected by this view, for the shape function and quadrature point selected by the arguments.
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. void FEValuesViews::Scalar< dim, spacedim >::get_function_values | ( | const InputVector & | fe_function, |
std::vector< solution_value_type< typename InputVector::value_type > > & | values | ||
) | const |
Return the values 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_values 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 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).
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 1547 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_values_from_local_dof_values | ( | const InputVector & | dof_values, |
std::vector< solution_value_type< typename InputVector::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
(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.
[in] | dof_values | A 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] | values | A vector of values of the given finite element field, at the quadrature points on the current object. |
InputVector | The InputVector type must allow creation of an ArrayView object from it; this is satisfied by the std::vector class, among others. |
Definition at line 1578 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_gradients | ( | const InputVector & | fe_function, |
std::vector< solution_gradient_type< typename InputVector::value_type > > & | gradients | ||
) | const |
Return the gradients 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_gradients 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 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).
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 1602 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_gradients_from_local_dof_values | ( | const InputVector & | dof_values, |
std::vector< solution_gradient_type< typename InputVector::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 1632 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_hessians | ( | const InputVector & | fe_function, |
std::vector< solution_hessian_type< typename InputVector::value_type > > & | hessians | ||
) | const |
Return the Hessians 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_hessians 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 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).
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 1656 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_hessians_from_local_dof_values | ( | const InputVector & | dof_values, |
std::vector< solution_hessian_type< typename InputVector::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 1686 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_laplacians | ( | const InputVector & | fe_function, |
std::vector< solution_laplacian_type< typename InputVector::value_type > > & | laplacians | ||
) | const |
Return the Laplacians 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. 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 scalar component.
The data type stored by the output vector must be what you get when you multiply the Laplacians 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).
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 1710 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_laplacians_from_local_dof_values | ( | const InputVector & | dof_values, |
std::vector< solution_laplacian_type< typename InputVector::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 1740 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::get_function_third_derivatives | ( | const InputVector & | fe_function, |
std::vector< solution_third_derivative_type< typename InputVector::value_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).
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 1764 of file fe_values.cc.
void FEValuesViews::Scalar< dim, spacedim >::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 |
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 1795 of file fe_values.cc.
|
private |
A pointer to the FEValuesBase object we operate on.
Definition at line 628 of file fe_values.h.
|
private |
The single scalar component this view represents of the FEValuesBase object.
Definition at line 634 of file fe_values.h.
|
private |
Store the data about shape functions.
Definition at line 639 of file fe_values.h.