Classes | |
struct | ShapeFunctionData |
Public Types | |
typedef double | value_type |
typedef Tensor< 1, spacedim > | gradient_type |
typedef Tensor< 2, spacedim > | hessian_type |
Public Member Functions | |
Scalar () | |
Scalar (const FEValuesBase< dim, spacedim > &fe_values_base, const unsigned int component) | |
Scalar & | operator= (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 |
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_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 | 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.
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.
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>
.
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>
.
FEValuesViews::Scalar< dim, spacedim >::Scalar | ( | ) |
Default constructor. Creates an invalid object.
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).
Scalar& 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.
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.
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.
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.
void FEValuesViews::Scalar< dim, spacedim >::get_function_values | ( | const InputVector & | fe_function, | |
std::vector< value_type > & | values | |||
) | const [inline] |
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, at the quadrature points.
This function is the equivalent of the FEValuesBase::get_function_values function but it only works on the selected scalar component.
void FEValuesViews::Scalar< dim, spacedim >::get_function_gradients | ( | const InputVector & | fe_function, | |
std::vector< gradient_type > & | gradients | |||
) | const [inline] |
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, at the quadrature points.
This function is the equivalent of the FEValuesBase::get_function_gradients function but it only works on the selected scalar component.
void FEValuesViews::Scalar< dim, spacedim >::get_function_hessians | ( | const InputVector & | fe_function, | |
std::vector< hessian_type > & | hessians | |||
) | const [inline] |
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, at the quadrature points.
This function is the equivalent of the FEValuesBase::get_function_hessians function but it only works on the selected scalar component.
void FEValuesViews::Scalar< dim, spacedim >::get_function_laplacians | ( | const InputVector & | fe_function, | |
std::vector< value_type > & | laplacians | |||
) | const [inline] |
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, 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 scalar component.
const FEValuesBase<dim,spacedim>& FEValuesViews::Scalar< dim, spacedim >::fe_values [private] |
A reference to the FEValuesBase object we operator on.
const unsigned int FEValuesViews::Scalar< dim, spacedim >::component [private] |
The single scalar component this view represents of the FEValuesBase object.
std::vector<ShapeFunctionData> FEValuesViews::Scalar< dim, spacedim >::shape_function_data [private] |
Store the data about shape functions.