Reference documentation for deal.II version 8.5.1
Classes | Public Types | Public Member Functions | Private Attributes | List of all members

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

Classes

struct  ShapeFunctionData
 

Public Types

typedef double value_type
 
typedef ::Tensor< 1, spacedim > gradient_type
 
typedef ::Tensor< 2, spacedim > hessian_type
 
typedef ::Tensor< 3, spacedim > third_derivative_type
 

Public Member Functions

 Scalar ()
 
 Scalar (const FEValuesBase< dim, spacedim > &fe_values_base, const unsigned int component)
 
Scalaroperator= (const Scalar< dim, spacedim > &)
 
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< typename ProductType< value_type, typename InputVector::value_type >::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_hessians (const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::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_third_derivatives (const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
 

Private Attributes

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

Detailed Description

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

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 137 of file fe_values.h.

Member Typedef Documentation

◆ value_type

template<int dim, int spacedim = dim>
typedef double FEValuesViews::Scalar< dim, spacedim >::value_type

A typedef 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 145 of file fe_values.h.

◆ gradient_type

template<int dim, int spacedim = dim>
typedef ::Tensor<1,spacedim> FEValuesViews::Scalar< dim, spacedim >::gradient_type

A typedef 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 152 of file fe_values.h.

◆ hessian_type

template<int dim, int spacedim = dim>
typedef ::Tensor<2,spacedim> FEValuesViews::Scalar< dim, spacedim >::hessian_type

A typedef 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 159 of file fe_values.h.

◆ third_derivative_type

template<int dim, int spacedim = dim>
typedef ::Tensor<3,spacedim> FEValuesViews::Scalar< dim, spacedim >::third_derivative_type

A typedef 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 166 of file fe_values.h.

Constructor & Destructor Documentation

◆ Scalar() [1/2]

template<int dim, int spacedim>
FEValuesViews::Scalar< dim, spacedim >::Scalar ( )

Default constructor. Creates an invalid object.

Definition at line 139 of file fe_values.cc.

◆ Scalar() [2/2]

template<int dim, int spacedim>
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 98 of file fe_values.cc.

Member Function Documentation

◆ operator=()

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

Copy operator. This is not a lightweight object so we don't allow copying and generate an exception if this function is called.

Definition at line 148 of file fe_values.cc.

◆ value()

template<int dim, int spacedim = dim>
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.

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::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.

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::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.

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::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.

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_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.

◆ get_function_values()

template<int dim, int spacedim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
std::vector< typename ProductType< value_type, typename InputVector::value_type >::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).

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

◆ get_function_gradients()

template<int dim, int spacedim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::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).

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

◆ get_function_hessians()

template<int dim, int spacedim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::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).

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

◆ get_function_laplacians()

template<int dim, int spacedim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
std::vector< typename ProductType< value_type, typename InputVector::value_type >::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).

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

◆ get_function_third_derivatives()

template<int dim, int spacedim>
template<class InputVector >
void FEValuesViews::Scalar< dim, spacedim >::get_function_third_derivatives ( const InputVector &  fe_function,
std::vector< typename ProductType< third_derivative_type, typename InputVector::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 1361 of file fe_values.cc.

Member Data Documentation

◆ fe_values

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

A pointer to the FEValuesBase object we operate on.

Definition at line 385 of file fe_values.h.

◆ component

template<int dim, int spacedim = dim>
const unsigned int FEValuesViews::Scalar< dim, spacedim >::component
private

The single scalar component this view represents of the FEValuesBase object.

Definition at line 391 of file fe_values.h.

◆ shape_function_data

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

Store the data about shape functions.

Definition at line 396 of file fe_values.h.


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