Classes | |
struct | ShapeFunctionData |
Public Types | |
typedef Tensor< 1, spacedim > | value_type |
typedef Tensor< 2, spacedim > | gradient_type |
typedef SymmetricTensor < 2, spacedim > | symmetric_gradient_type |
typedef double | divergence_type |
typedef Tensor< 3, spacedim > | hessian_type |
Public Member Functions | |
Vector () | |
Vector (const FEValuesBase< dim, spacedim > &fe_values_base, const unsigned int first_vector_component) | |
Vector & | operator= (const Vector< 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 |
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 |
hessian_type | hessian (const unsigned int shape_function, const unsigned int q_point) const |
template<class InputVector > | |
void | get_function_values (const InputVector &fe_function, std::vector< value_type > &values) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, std::vector< gradient_type > &gradients) const |
template<class InputVector > | |
void | get_function_symmetric_gradients (const InputVector &fe_function, std::vector< symmetric_gradient_type > &symmetric_gradients) const |
template<class InputVector > | |
void | get_function_divergences (const InputVector &fe_function, std::vector< divergence_type > &divergences) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, std::vector< hessian_type > &hessians) const |
template<class InputVector > | |
void | get_function_laplacians (const InputVector &fe_function, std::vector< value_type > &laplacians) const |
Private Attributes | |
const FEValuesBase< dim, spacedim > & | fe_values |
const unsigned int | first_vector_component |
std::vector< ShapeFunctionData > | shape_function_data |
dim
components forming a vector part of a vector-valued finite element. Views are discussed in the Handling vector valued problems module. typedef Tensor<1,spacedim> FEValuesViews::Vector< dim, spacedim >::value_type |
A typedef 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>.
typedef Tensor<2,spacedim> FEValuesViews::Vector< dim, spacedim >::gradient_type |
A typedef 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>
.
typedef SymmetricTensor<2,spacedim> FEValuesViews::Vector< dim, spacedim >::symmetric_gradient_type |
A typedef 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>
.
typedef double FEValuesViews::Vector< dim, spacedim >::divergence_type |
A typedef 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.
typedef Tensor<3,spacedim> FEValuesViews::Vector< dim, spacedim >::hessian_type |
A typedef 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>
.
FEValuesViews::Vector< dim, spacedim >::Vector | ( | ) |
Default constructor. Creates an invalid object.
FEValuesViews::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.
Vector& FEValuesViews::Vector< dim, spacedim >::operator= | ( | const Vector< 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.
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.
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.
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 , where
represents the
dim
components selected from the FEValuesBase object, and is the location of the
-th quadrature point.
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.
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.
void FEValuesViews::Vector< dim, spacedim >::get_function_values | ( | const InputVector< dim, spacedim > & | fe_function, | |
std::vector< value_type > & | values | |||
) | const [inline] |
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, at the quadrature points.
This function is the equivalent of the FEValuesBase::get_function_values function but it only works on the selected vector components.
void FEValuesViews::Vector< dim, spacedim >::get_function_gradients | ( | const InputVector< dim, spacedim > & | fe_function, | |
std::vector< gradient_type > & | gradients | |||
) | const [inline] |
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, at the quadrature points.
This function is the equivalent of the FEValuesBase::get_function_gradients function but it only works on the selected vector components.
void FEValuesViews::Vector< dim, spacedim >::get_function_symmetric_gradients | ( | const InputVector< dim, spacedim > & | fe_function, | |
std::vector< symmetric_gradient_type > & | symmetric_gradients | |||
) | const [inline] |
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, at the quadrature points.
There is no equivalent function such as FEValuesBase::get_function_gradients in the FEValues classes but the information can be obtained from FEValuesBase::get_function_gradients, of course.
void FEValuesViews::Vector< dim, spacedim >::get_function_divergences | ( | const InputVector< dim, spacedim > & | fe_function, | |
std::vector< divergence_type > & | divergences | |||
) | const [inline] |
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, at the quadrature points.
There is no equivalent function such as FEValuesBase::get_function_gradients in the FEValues classes but the information can be obtained from FEValuesBase::get_function_gradients, of course.
void FEValuesViews::Vector< dim, spacedim >::get_function_hessians | ( | const InputVector< dim, spacedim > & | fe_function, | |
std::vector< hessian_type > & | hessians | |||
) | const [inline] |
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, at the quadrature points.
This function is the equivalent of the FEValuesBase::get_function_hessians function but it only works on the selected vector components.
void FEValuesViews::Vector< dim, spacedim >::get_function_laplacians | ( | const InputVector< dim, spacedim > & | fe_function, | |
std::vector< value_type > & | laplacians | |||
) | const [inline] |
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, at the quadrature points. 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.
const FEValuesBase<dim,spacedim>& FEValuesViews::Vector< dim, spacedim >::fe_values [private] |
A reference to the FEValuesBase object we operator on.
const unsigned int FEValuesViews::Vector< dim, spacedim >::first_vector_component [private] |
The first component of the vector this view represents of the FEValuesBase object.
std::vector<ShapeFunctionData> FEValuesViews::Vector< dim, spacedim >::shape_function_data [private] |
Store the data about shape functions.