![]() |
Reference documentation for deal.II version 8.1.0
|
#include <fe_values.h>
Classes | |
class | CellIterator |
Public Member Functions | |
FEValuesBase (const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe) | |
~FEValuesBase () | |
DeclException1 (ExcAccessToUninitializedField, char *,<< ("You are requesting information from an FEValues/FEFaceValues/FESubfaceValues ""object for which this kind of information has not been computed. What ""information these objects compute is determined by the update_* flags you ""pass to the constructor. Here, the operation you are attempting requires ""the <")<< arg1<< "> flag to be set, but it was apparently not specified upon construction.") | |
DeclException0 (ExcCannotInitializeField) | |
DeclException0 (ExcInvalidUpdateFlag) | |
DeclException0 (ExcFEDontMatch) | |
DeclException1 (ExcShapeFunctionNotPrimitive, int,<< "The shape function with index "<< arg1<< " is not primitive, i.e. it is vector-valued and "<< "has more than one non-zero vector component. This "<< "function cannot be called for these shape functions. "<< "Maybe you want to use the same function with the "<< "_component suffix?") | |
DeclException0 (ExcFENotPrimitive) | |
ShapeAccess Access to shape function values | |
const double & | shape_value (const unsigned int function_no, const unsigned int point_no) const |
double | shape_value_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const |
const Tensor< 1, spacedim > & | shape_grad (const unsigned int function_no, const unsigned int quadrature_point) const |
Tensor< 1, spacedim > | shape_grad_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const |
const Tensor< 2, spacedim > & | shape_hessian (const unsigned int function_no, const unsigned int point_no) const |
const Tensor< 2, spacedim > & | shape_2nd_derivative (const unsigned int function_no, const unsigned int point_no) const DEAL_II_DEPRECATED |
Tensor< 2, spacedim > | shape_hessian_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const |
Tensor< 2, spacedim > | shape_2nd_derivative_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const DEAL_II_DEPRECATED |
Access to values of global finite element fields | |
template<class InputVector , typename number > | |
void | get_function_values (const InputVector &fe_function, std::vector< number > &values) const |
template<class InputVector , typename number > | |
void | get_function_values (const InputVector &fe_function, std::vector< Vector< number > > &values) const |
template<class InputVector , typename number > | |
void | get_function_values (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< number > &values) const |
template<class InputVector , typename number > | |
void | get_function_values (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< Vector< number > > &values) const |
template<class InputVector > | |
void | get_function_values (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, VectorSlice< std::vector< std::vector< double > > > values, const bool quadrature_points_fastest) const |
Access to derivatives of global finite element fields | |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, std::vector< Tensor< 1, spacedim > > &gradients) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, std::vector< std::vector< Tensor< 1, spacedim > > > &gradients) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< Tensor< 1, spacedim > > &gradients) const |
template<class InputVector > | |
void | get_function_gradients (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, VectorSlice< std::vector< std::vector< Tensor< 1, spacedim > > > > gradients, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_grads (const InputVector &fe_function, std::vector< Tensor< 1, spacedim > > &gradients) const DEAL_II_DEPRECATED |
template<class InputVector > | |
void | get_function_grads (const InputVector &fe_function, std::vector< std::vector< Tensor< 1, spacedim > > > &gradients) const DEAL_II_DEPRECATED |
template<class InputVector > | |
void | get_function_grads (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< Tensor< 1, spacedim > > &gradients) const DEAL_II_DEPRECATED |
template<class InputVector > | |
void | get_function_grads (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< std::vector< Tensor< 1, spacedim > > > &gradients, bool quadrature_points_fastest=false) const DEAL_II_DEPRECATED |
Access to second derivatives (Hessian matrices and Laplacians) of global finite element fields | |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, std::vector< Tensor< 2, spacedim > > &hessians) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, std::vector< std::vector< Tensor< 2, spacedim > > > &hessians, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< Tensor< 2, spacedim > > &hessians) const |
template<class InputVector > | |
void | get_function_hessians (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, VectorSlice< std::vector< std::vector< Tensor< 2, spacedim > > > > hessians, bool quadrature_points_fastest=false) const |
template<class InputVector > | |
void | get_function_2nd_derivatives (const InputVector &, std::vector< Tensor< 2, spacedim > > &) const DEAL_II_DEPRECATED |
template<class InputVector > | |
void | get_function_2nd_derivatives (const InputVector &, std::vector< std::vector< Tensor< 2, spacedim > > > &, bool=false) const DEAL_II_DEPRECATED |
template<class InputVector , typename number > | |
void | get_function_laplacians (const InputVector &fe_function, std::vector< number > &laplacians) const |
template<class InputVector , typename number > | |
void | get_function_laplacians (const InputVector &fe_function, std::vector< Vector< number > > &laplacians) const |
template<class InputVector , typename number > | |
void | get_function_laplacians (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< number > &laplacians) const |
template<class InputVector , typename number > | |
void | get_function_laplacians (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< Vector< number > > &laplacians) const |
template<class InputVector , typename number > | |
void | get_function_laplacians (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< std::vector< number > > &laplacians, bool quadrature_points_fastest=false) const |
Geometry of the cell | |
const Point< spacedim > & | quadrature_point (const unsigned int i) const |
const std::vector< Point< spacedim > > & | get_quadrature_points () const |
double | JxW (const unsigned int quadrature_point) const |
const std::vector< double > & | get_JxW_values () const |
const DerivativeForm< 1, dim, spacedim > & | jacobian (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 1, dim, spacedim > > & | get_jacobians () const |
const DerivativeForm< 2, dim, spacedim > & | jacobian_grad (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 2, dim, spacedim > > & | get_jacobian_grads () const |
const DerivativeForm< 1, spacedim, dim > & | inverse_jacobian (const unsigned int quadrature_point) const |
const std::vector< DerivativeForm< 1, spacedim, dim > > & | get_inverse_jacobians () const |
const Point< spacedim > & | normal_vector (const unsigned int i) const |
const std::vector< Point< spacedim > > & | get_normal_vectors () const |
void | transform (std::vector< Tensor< 1, spacedim > > &transformed, const std::vector< Tensor< 1, dim > > &original, MappingType mapping) const |
const Point< spacedim > & | cell_normal_vector (const unsigned int i) const DEAL_II_DEPRECATED |
const std::vector< Point< spacedim > > & | get_cell_normal_vectors () const DEAL_II_DEPRECATED |
Extractors Methods to extract individual components | |
const FEValuesViews::Scalar< dim, spacedim > & | operator[] (const FEValuesExtractors::Scalar &scalar) const |
const FEValuesViews::Vector< dim, spacedim > & | operator[] (const FEValuesExtractors::Vector &vector) const |
const FEValuesViews::SymmetricTensor< 2, dim, spacedim > & | operator[] (const FEValuesExtractors::SymmetricTensor< 2 > &tensor) const |
const FEValuesViews::Tensor< 2, dim, spacedim > & | operator[] (const FEValuesExtractors::Tensor< 2 > &tensor) const |
Access to the raw data | |
const Mapping< dim, spacedim > & | get_mapping () const |
const FiniteElement< dim, spacedim > & | get_fe () const |
UpdateFlags | get_update_flags () const |
const Triangulation< dim, spacedim >::cell_iterator | get_cell () const |
CellSimilarity::Similarity | get_cell_similarity () const |
std::size_t | memory_consumption () const |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.") | |
DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1) | |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Public Attributes | |
const unsigned int | n_quadrature_points |
const unsigned int | dofs_per_cell |
Static Public Attributes | |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
Protected Member Functions | |
void | invalidate_present_cell () |
void | maybe_invalidate_previous_present_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell) |
UpdateFlags | compute_update_flags (const UpdateFlags update_flags) const |
void | check_cell_similarity (const typename Triangulation< dim, spacedim >::cell_iterator &cell) |
![]() | |
void | initialize (const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags) |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.") | |
DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1) | |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Protected Attributes | |
std::auto_ptr< const CellIteratorBase > | present_cell |
boost::signals2::connection | tria_listener |
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > | mapping |
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > | fe |
SmartPointer< typename Mapping< dim, spacedim >::InternalDataBase, FEValuesBase< dim, spacedim > > | mapping_data |
SmartPointer< typename Mapping< dim, spacedim >::InternalDataBase, FEValuesBase< dim, spacedim > > | fe_data |
CellSimilarity::Similarity | cell_similarity |
![]() | |
ShapeVector | shape_values |
GradientVector | shape_gradients |
HessianVector | shape_hessians |
std::vector< double > | JxW_values |
std::vector< DerivativeForm< 1, dim, spacedim > > | jacobians |
std::vector< DerivativeForm< 2, dim, spacedim > > | jacobian_grads |
std::vector< DerivativeForm< 1, spacedim, dim > > | inverse_jacobians |
std::vector< Point< spacedim > > | quadrature_points |
std::vector< Point< spacedim > > | normal_vectors |
std::vector< Tensor< 1, spacedim > > | boundary_forms |
std::vector< unsigned int > | shape_function_to_row_table |
UpdateFlags | update_flags |
Private Member Functions | |
FEValuesBase (const FEValuesBase &) | |
FEValuesBase & | operator= (const FEValuesBase &) |
Private Attributes | |
::internal::FEValuesViews::Cache< dim, spacedim > | fe_values_views_cache |
Friends | |
template<int , int > | |
class | FEValuesViews::Scalar |
template<int , int > | |
class | FEValuesViews::Vector |
template<int , int , int > | |
class | FEValuesViews::SymmetricTensor |
template<int , int , int > | |
class | FEValuesViews::Tensor |
Additional Inherited Members | |
![]() | |
typedef Table< 2, double > | ShapeVector |
typedef std::vector< std::vector< Tensor< 1, spacedim > > > | GradientVector |
typedef std::vector< std::vector< Tensor< 2, spacedim > > > | HessianVector |
FEValues, FEFaceValues and FESubfaceValues objects are interfaces to finite element and mapping classes on the one hand side, to cells and quadrature rules on the other side. They allow to evaluate values or derivatives of shape functions at the quadrature points of a quadrature formula when projected by a mapping from the unit cell onto a cell in real space. The reason for this abstraction is possible optimization: Depending on the type of finite element and mapping, some values can be computed once on the unit cell. Others must be computed on each cell, but maybe computation of several values at the same time offers ways for optimization. Since this interlay may be complex and depends on the actual finite element, it cannot be left to the applications programmer.
FEValues, FEFaceValues and FESubfaceValues provide only data handling: computations are left to objects of type Mapping and FiniteElement. These provide functions get_*_data
and fill_*_values
which are called by the constructor and reinit
functions of FEValues*
, respectively.
Usually, an object of FEValues*
is used in integration loops over all cells of a triangulation (or faces of cells). To take full advantage of the optimization features, it should be constructed before the loop so that information that does not depend on the location and shape of cells can be computed once and for all (this includes, for example, the values of shape functions at quadrature points for the most common elements: we can evaluate them on the unit cell and they will be the same when mapped to the real cell). Then, in the loop over all cells, it must be re-initialized for each grid cell to compute that part of the information that changes depending on the actual cell (for example, the gradient of shape functions equals the gradient on the unit cell – which can be computed once and for all – times the Jacobian matrix of the mapping between unit and real cell, which needs to be recomputed for each cell).
A typical piece of code, adding up local contributions to the Laplace matrix looks like this:
The individual functions used here are described below. Note that by design, the order of quadrature points used inside the FEValues object is the same as defined by the quadrature formula passed to the constructor of the FEValues object above.
The functions of this class fall into different cathegories:
shape_value(), shape_grad(), etc: return one of the values of this object at a time. These functions are inlined, so this is the suggested access to all finite element values. There should be no loss in performance with an optimizing compiler. If the finite element is vector valued, then these functions return the only non-zero component of the requested shape function. However, some finite elements have shape functions that have more than one non-zero component (we call them non-"primitive"), and in this case this set of functions will throw an exception since they cannot generate a useful result. Rather, use the next set of functions.
shape_value_component(), shape_grad_component(), etc: This is the same set of functions as above, except that for vector valued finite elements they return only one vector component. This is useful for elements of which shape functions have more than one non-zero component, since then the above functions cannot be used, and you have to walk over all (or only the non-zero) components of the shape function using this set of functions.
get_function_values(), get_function_gradients(), etc.: Compute a finite element function or its derivative in quadrature points.
The UpdateFlags object handed to the constructor is used to determine which of the data fields to compute. This way, it is possible to avoid expensive computations of useless derivatives. In the beginning, these flags are processed through the functions Mapping::update_once(), Mapping::update_each(), FiniteElement::update_once() FiniteElement::update_each(). All the results are bit-wise or'd and determine the fields actually computed. This enables Mapping and FiniteElement to schedule auxiliary data fields for updating. Still, it is recommended to give all needed update flags to FEValues.
The mechanisms by which this class works is also discussed on the page on The interplay of UpdateFlags, Mapping and FiniteElement in FEValues.
Definition at line 29 of file data_out.h.
FEValuesBase< dim, spacedim >::FEValuesBase | ( | const unsigned int | n_q_points, |
const unsigned int | dofs_per_cell, | ||
const UpdateFlags | update_flags, | ||
const Mapping< dim, spacedim > & | mapping, | ||
const FiniteElement< dim, spacedim > & | fe | ||
) |
Constructor. Set up the array sizes with n_q_points
quadrature points, dofs_per_cell
trial functions per cell and with the given pattern to update the fields when the reinit
function of the derived classes is called. The fields themselves are not set up, this must happen in the constructor of the derived class.
FEValuesBase< dim, spacedim >::~FEValuesBase | ( | ) |
Destructor.
|
private |
Copy constructor. Since objects of this class are not copyable, we make it private, and also do not implement it.
const double& FEValuesBase< dim, spacedim >::shape_value | ( | const unsigned int | function_no, |
const unsigned int | point_no | ||
) | const |
Value of a shape function at a quadrature point on the cell, face or subface selected the last time the reinit
function of the derived class was called.
If the shape function is vector-valued, then this returns the only non-zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_value_component() function.
function_no | 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. |
point_no | Number of the quadrature point at which function is to be evaluated |
double FEValuesBase< dim, spacedim >::shape_value_component | ( | const unsigned int | function_no, |
const unsigned int | point_no, | ||
const unsigned int | component | ||
) | const |
Compute one vector component of the value of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_value() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_value() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.
function_no | Number of the shape function to be evaluated |
point_no | Number of the quadrature point at which function is to be evaluated |
component | vector component to be evaluated |
const Tensor<1,spacedim>& FEValuesBase< dim, spacedim >::shape_grad | ( | const unsigned int | function_no, |
const unsigned int | quadrature_point | ||
) | const |
Compute the gradient of the function_no
th shape function at the quadrature_point
th quadrature point with respect to real cell coordinates. If you want to get the derivative in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component of the Tensor returned by this function. Since only a reference to the gradient's value is returned, there should be no major performance drawback.
If the shape function is vector-valued, then this returns the only non-zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then it will throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_grad_component() function.
The same holds for the arguments of this function as for the shape_value() function.
Tensor<1,spacedim> FEValuesBase< dim, spacedim >::shape_grad_component | ( | const unsigned int | function_no, |
const unsigned int | point_no, | ||
const unsigned int | component | ||
) | const |
Return one vector component of the gradient of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_grad() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_grad() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.
The same holds for the arguments of this function as for the shape_value_component() function.
const Tensor<2,spacedim>& FEValuesBase< dim, spacedim >::shape_hessian | ( | const unsigned int | function_no, |
const unsigned int | point_no | ||
) | const |
Second derivatives of the function_no
th shape function at the point_no
th quadrature point with respect to real cell coordinates. If you want to get the derivatives in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component. Since only a reference to the derivative values is returned, there should be no major performance drawback.
If the shape function is vector-valued, then this returns the only non-zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_grad_grad_component() function.
The same holds for the arguments of this function as for the shape_value() function.
const Tensor<2,spacedim>& FEValuesBase< dim, spacedim >::shape_2nd_derivative | ( | const unsigned int | function_no, |
const unsigned int | point_no | ||
) | const |
Tensor<2,spacedim> FEValuesBase< dim, spacedim >::shape_hessian_component | ( | const unsigned int | function_no, |
const unsigned int | point_no, | ||
const unsigned int | component | ||
) | const |
Return one vector component of the gradient of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_hessian() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_hessian() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.
The same holds for the arguments of this function as for the shape_value_component() function.
Tensor<2,spacedim> FEValuesBase< dim, spacedim >::shape_2nd_derivative_component | ( | const unsigned int | function_no, |
const unsigned int | point_no, | ||
const unsigned int | component | ||
) | const |
void FEValuesBase< dim, spacedim >::get_function_values | ( | const InputVector & | fe_function, |
std::vector< number > & | values | ||
) | const |
Returns the values of a finite element function restricted to the current cell, face or subface selected the last time the reinit
function of the derived class was called, at the quadrature points.
If the present cell is not active then values are interpolated to the current cell and point values are computed from that.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. To get values of multi-component elements, there is another get_function_values() below, returning a vector of vectors of results.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | values | The values of the function specified by fe_function at the quadrature points of the current cell. The object is assume to already have the correct size. |
values[q]
will contain the value of the field described by fe_function at the void FEValuesBase< dim, spacedim >::get_function_values | ( | const InputVector & | fe_function, |
std::vector< Vector< number > > & | values | ||
) | const |
This function does the same as the other get_function_values(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
values[q]
is a vector of values of the field described by fe_function at the values[q]
equals the number of components of the finite element, i.e. values[q](c)
returns the value of the void FEValuesBase< dim, spacedim >::get_function_values | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< number > & | values | ||
) | const |
Generate function values from an arbitrary vector.
This function offers the possibility to extract function values in quadrature points from vectors not corresponding to a whole discretization.
The vector indices
corresponds to the degrees of freedom on a single cell. Its length may even be a multiple of the number of dofs per cell. Then, the vectors in value
should allow for the same multiple of the components of the finite element.
You may want to use this function, if you want to access just a single block from a BlockVector, if you have a multi-level vector or if you already have a local representation of your finite element data.
void FEValuesBase< dim, spacedim >::get_function_values | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< Vector< number > > & | values | ||
) | const |
Generate vector function values from an arbitrary vector.
This function offers the possibility to extract function values in quadrature points from vectors not corresponding to a whole discretization.
The vector indices
corresponds to the degrees of freedom on a single cell. Its length may even be a multiple of the number of dofs per cell. Then, the vectors in value
should allow for the same multiple of the components of the finite element.
You may want to use this function, if you want to access just a single block from a BlockVector, if you have a multi-level vector or if you already have a local representation of your finite element data.
Since this function allows for fairly general combinations of argument sizes, be aware that the checks on the arguments may not detect errors.
void FEValuesBase< dim, spacedim >::get_function_values | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
VectorSlice< std::vector< std::vector< double > > > | values, | ||
const bool | quadrature_points_fastest | ||
) | const |
Generate vector function values from an arbitrary vector.
This function offers the possibility to extract function values in quadrature points from vectors not corresponding to a whole discretization.
The vector indices
corresponds to the degrees of freedom on a single cell. Its length may even be a multiple of the number of dofs per cell. Then, the vectors in value
should allow for the same multiple of the components of the finite element.
Depending on the value of the last argument, the outer vector of values
has either the length of the quadrature rule (quadrature_points_fastest == false
) or the length of components to be filled quadrature_points_fastest == true
. If p
is the current quadrature point number and i
is the vector component of the solution desired, the access to values
is values[p][i]
if quadrature_points_fastest == false
, and values[i][p]
otherwise.
You may want to use this function, if you want to access just a single block from a BlockVector, if you have a multi-level vector or if you already have a local representation of your finite element data.
Since this function allows for fairly general combinations of argument sizes, be aware that the checks on the arguments may not detect errors.
void FEValuesBase< dim, spacedim >::get_function_gradients | ( | const InputVector & | fe_function, |
std::vector< Tensor< 1, spacedim > > & | gradients | ||
) | const |
Compute the gradients of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's gradient instead of its value.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | gradients | The gradients of the function specified by fe_function at the quadrature points of the current cell. The gradients are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. |
gradients[q]
will contain the gradient of the field described by fe_function at the gradients[q][d]
represents the derivative in coordinate direction void FEValuesBase< dim, spacedim >::get_function_gradients | ( | const InputVector & | fe_function, |
std::vector< std::vector< Tensor< 1, spacedim > > > & | gradients | ||
) | const |
This function does the same as the other get_function_gradients(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
gradients[q]
is a vector of gradients of the field described by fe_function at the gradients[q]
equals the number of components of the finite element, i.e. gradients[q][c]
returns the gradient of the gradients[q][c][d]
is the derivative in coordinate direction void FEValuesBase< dim, spacedim >::get_function_gradients | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< Tensor< 1, spacedim > > & | gradients | ||
) | const |
Function gradient access with more flexibility. see get_function_values() with corresponding arguments.
void FEValuesBase< dim, spacedim >::get_function_gradients | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
VectorSlice< std::vector< std::vector< Tensor< 1, spacedim > > > > | gradients, | ||
bool | quadrature_points_fastest = false |
||
) | const |
Function gradient access with more flexibility. see get_function_values() with corresponding arguments.
void FEValuesBase< dim, spacedim >::get_function_grads | ( | const InputVector & | fe_function, |
std::vector< Tensor< 1, spacedim > > & | gradients | ||
) | const |
void FEValuesBase< dim, spacedim >::get_function_grads | ( | const InputVector & | fe_function, |
std::vector< std::vector< Tensor< 1, spacedim > > > & | gradients | ||
) | const |
void FEValuesBase< dim, spacedim >::get_function_grads | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< Tensor< 1, spacedim > > & | gradients | ||
) | const |
void FEValuesBase< dim, spacedim >::get_function_grads | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< std::vector< Tensor< 1, spacedim > > > & | gradients, | ||
bool | quadrature_points_fastest = false |
||
) | const |
void FEValuesBase< dim, spacedim >::get_function_hessians | ( | const InputVector & | fe_function, |
std::vector< Tensor< 2, spacedim > > & | hessians | ||
) | const |
Compute the tensor of second derivatives of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's second derivatives instead of its value.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | hessians | The Hessians of the function specified by fe_function at the quadrature points of the current cell. The Hessians are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. |
hessians[q]
will contain the Hessian of the field described by fe_function at the gradients[q][i][j]
represents the void FEValuesBase< dim, spacedim >::get_function_hessians | ( | const InputVector & | fe_function, |
std::vector< std::vector< Tensor< 2, spacedim > > > & | hessians, | ||
bool | quadrature_points_fastest = false |
||
) | const |
This function does the same as the other get_function_hessians(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
hessians[q]
is a vector of Hessians of the field described by fe_function at the hessians[q]
equals the number of components of the finite element, i.e. hessians[q][c]
returns the Hessian of the values[q][c][i][j]
is the void FEValuesBase< dim, spacedim >::get_function_hessians | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< Tensor< 2, spacedim > > & | hessians | ||
) | const |
Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.
void FEValuesBase< dim, spacedim >::get_function_hessians | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
VectorSlice< std::vector< std::vector< Tensor< 2, spacedim > > > > | hessians, | ||
bool | quadrature_points_fastest = false |
||
) | const |
Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.
void FEValuesBase< dim, spacedim >::get_function_2nd_derivatives | ( | const InputVector & | , |
std::vector< Tensor< 2, spacedim > > & | |||
) | const |
void FEValuesBase< dim, spacedim >::get_function_2nd_derivatives | ( | const InputVector & | , |
std::vector< std::vector< Tensor< 2, spacedim > > > & | , | ||
bool | = false |
||
) | const |
void FEValuesBase< dim, spacedim >::get_function_laplacians | ( | const InputVector & | fe_function, |
std::vector< number > & | laplacians | ||
) | const |
Compute the (scalar) Laplacian (i.e. the trace of the tensor of second derivatives) of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's second derivatives instead of its value.
This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.
[in] | fe_function | A vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell. |
[out] | laplacians | The Laplacians of the function specified by fe_function at the quadrature points of the current cell. The Laplacians are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. |
laplacians[q]
will contain the Laplacian of the field described by fe_function at the gradients[q][i][j]
represents the laplacians[q]=trace(hessians[q])
, where hessians
would be the output of the get_function_hessians() function.void FEValuesBase< dim, spacedim >::get_function_laplacians | ( | const InputVector & | fe_function, |
std::vector< Vector< number > > & | laplacians | ||
) | const |
This function does the same as the other get_function_laplacians(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.
laplacians[q]
is a vector of Laplacians of the field described by fe_function at the laplacians[q]
equals the number of components of the finite element, i.e. laplacians[q][c]
returns the Laplacian of the laplacians[q][c]=trace(hessians[q][c])
, where hessians
would be the output of the get_function_hessians() function. void FEValuesBase< dim, spacedim >::get_function_laplacians | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< number > & | laplacians | ||
) | const |
Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.
void FEValuesBase< dim, spacedim >::get_function_laplacians | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< Vector< number > > & | laplacians | ||
) | const |
Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.
void FEValuesBase< dim, spacedim >::get_function_laplacians | ( | const InputVector & | fe_function, |
const VectorSlice< const std::vector< types::global_dof_index > > & | indices, | ||
std::vector< std::vector< number > > & | laplacians, | ||
bool | quadrature_points_fastest = false |
||
) | const |
Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.
const Point<spacedim>& FEValuesBase< dim, spacedim >::quadrature_point | ( | const unsigned int | i | ) | const |
Position of the i
th quadrature point in real space.
const std::vector<Point<spacedim> >& FEValuesBase< dim, spacedim >::get_quadrature_points | ( | ) | const |
Return a pointer to the vector of quadrature points in real space.
double FEValuesBase< dim, spacedim >::JxW | ( | const unsigned int | quadrature_point | ) | const |
Mapped quadrature weight. If this object refers to a volume evaluation (i.e. the derived class is of type FEValues), then this is the Jacobi determinant times the weight of the *i
th unit quadrature point.
For surface evaluations (i.e. classes FEFaceValues or FESubfaceValues), it is the mapped surface element times the weight of the quadrature point.
You can think of the quantity returned by this function as the volume or surface element in the integral that we implement here by quadrature.
const std::vector<double>& FEValuesBase< dim, spacedim >::get_JxW_values | ( | ) | const |
Pointer to the array holding the values returned by JxW().
const DerivativeForm<1,dim,spacedim>& FEValuesBase< dim, spacedim >::jacobian | ( | const unsigned int | quadrature_point | ) | const |
Return the Jacobian of the transformation at the specified quadrature point, i.e.
const std::vector<DerivativeForm<1,dim,spacedim> >& FEValuesBase< dim, spacedim >::get_jacobians | ( | ) | const |
Pointer to the array holding the values returned by jacobian().
const DerivativeForm<2,dim,spacedim>& FEValuesBase< dim, spacedim >::jacobian_grad | ( | const unsigned int | quadrature_point | ) | const |
Return the second derivative of the transformation from unit to real cell, i.e. the first derivative of the Jacobian, at the specified quadrature point, i.e. .
const std::vector<DerivativeForm<2,dim,spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_grads | ( | ) | const |
Pointer to the array holding the values returned by jacobian_grads().
const DerivativeForm<1,spacedim,dim>& FEValuesBase< dim, spacedim >::inverse_jacobian | ( | const unsigned int | quadrature_point | ) | const |
Return the inverse Jacobian of the transformation at the specified quadrature point, i.e.
const std::vector<DerivativeForm<1,spacedim,dim> >& FEValuesBase< dim, spacedim >::get_inverse_jacobians | ( | ) | const |
Pointer to the array holding the values returned by inverse_jacobian().
const Point<spacedim>& FEValuesBase< dim, spacedim >::normal_vector | ( | const unsigned int | i | ) | const |
For a face, return the outward normal vector to the cell at the i
th quadrature point.
For a cell of codimension one, return the normal vector, as it is specified by the numbering of the vertices.
The length of the vector is normalized to one.
const std::vector<Point<spacedim> >& FEValuesBase< dim, spacedim >::get_normal_vectors | ( | ) | const |
Return the normal vectors at the quadrature points. For a face, these are the outward normal vectors to the cell. For a cell of codimension one, the orientation is given by the numbering of vertices.
void FEValuesBase< dim, spacedim >::transform | ( | std::vector< Tensor< 1, spacedim > > & | transformed, |
const std::vector< Tensor< 1, dim > > & | original, | ||
MappingType | mapping | ||
) | const |
Transform a set of vectors, one for each quadrature point. The mapping
can be any of the ones defined in MappingType.
const Point<spacedim>& FEValuesBase< dim, spacedim >::cell_normal_vector | ( | const unsigned int | i | ) | const |
Return the outward normal vector to the cell at the i
th quadrature point. The length of the vector is normalized to one.
const std::vector<Point<spacedim> >& FEValuesBase< dim, spacedim >::get_cell_normal_vectors | ( | ) | const |
Returns the vectors normal to the cell in each of the quadrature points.
const FEValuesViews::Scalar<dim,spacedim>& FEValuesBase< dim, spacedim >::operator[] | ( | const FEValuesExtractors::Scalar & | scalar | ) | const |
Create a view of the current FEValues object that represents a particular scalar component of the possibly vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
const FEValuesViews::Vector<dim,spacedim>& FEValuesBase< dim, spacedim >::operator[] | ( | const FEValuesExtractors::Vector & | vector | ) | const |
Create a view of the current FEValues object that represents a set of dim
scalar components (i.e. a vector) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
const FEValuesViews::SymmetricTensor<2,dim,spacedim>& FEValuesBase< dim, spacedim >::operator[] | ( | const FEValuesExtractors::SymmetricTensor< 2 > & | tensor | ) | const |
Create a view of the current FEValues object that represents a set of (dim*dim + dim)/2
scalar components (i.e. a symmetric 2nd order tensor) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
const FEValuesViews::Tensor<2,dim,spacedim>& FEValuesBase< dim, spacedim >::operator[] | ( | const FEValuesExtractors::Tensor< 2 > & | tensor | ) | const |
Create a view of the current FEValues object that represents a set of (dim*dim)
scalar components (i.e. a 2nd order tensor) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.
const Mapping<dim,spacedim>& FEValuesBase< dim, spacedim >::get_mapping | ( | ) | const |
Constant reference to the selected mapping object.
const FiniteElement<dim,spacedim>& FEValuesBase< dim, spacedim >::get_fe | ( | ) | const |
Constant reference to the selected finite element object.
UpdateFlags FEValuesBase< dim, spacedim >::get_update_flags | ( | ) | const |
Return the update flags set for this object.
const Triangulation<dim,spacedim>::cell_iterator FEValuesBase< dim, spacedim >::get_cell | ( | ) | const |
Return a triangulation iterator to the current cell.
CellSimilarity::Similarity FEValuesBase< dim, spacedim >::get_cell_similarity | ( | ) | const |
Return the relation of the current cell to the previous cell. This allows re-use of some cell data (like local matrices for equations with constant coefficients) if the result is CellSimilarity::translation
.
std::size_t FEValuesBase< dim, spacedim >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
|
protected |
A function that is connected to the triangulation in order to reset the stored 'present_cell' iterator to an invalid one whenever the triangulation is changed and the iterator consequently becomes invalid.
|
protected |
This function is called by the various reinit() functions in derived classes. Given the cell indicated by the argument, test whether we have to throw away the previously stored present_cell argument because it would require us to compare cells from different triangulations. In checking all this, also make sure that we have tria_listener connected to the triangulation to which we will set present_cell right after calling this function.
|
protected |
Initialize some update flags. Called from the initialize
functions of derived classes, which are in turn called from their constructors.
Basically, this function finds out using the finite element and mapping object already stored which flags need to be set to compute everything the user wants, as expressed through the flags passed as argument.
|
protected |
A function that checks whether the new cell is similar to the one previously used. Then, a significant amount of the data can be reused, e.g. the derivatives of the basis functions in real space, shape_grad.
|
private |
Copy operator. Since objects of this class are not copyable, we make it private, and also do not implement it.
Make the view classes friends of this class, since they access internal data.
Definition at line 2409 of file fe_values.h.
|
static |
Dimension in which this object operates.
Definition at line 1401 of file fe_values.h.
|
static |
Dimension of the space in which this object operates.
Definition at line 1406 of file fe_values.h.
const unsigned int FEValuesBase< dim, spacedim >::n_quadrature_points |
Number of quadrature points.
Definition at line 1411 of file fe_values.h.
const unsigned int FEValuesBase< dim, spacedim >::dofs_per_cell |
Number of shape functions per cell. If we use this base class to evaluate a finite element on faces of cells, this is still the number of degrees of freedom per cell, not per face.
Definition at line 1418 of file fe_values.h.
|
protected |
Store the cell selected last time the reinit() function was called. This is necessary for the get_function_*
functions as well as the functions of same name in the extractor classes.
Definition at line 2304 of file fe_values.h.
|
protected |
A signal connection we use to ensure we get informed whenever the triangulation changes. We need to know about that because it invalidates all cell iterators and, as part of that, the 'present_cell' iterator we keep around between subsequent calls to reinit() in order to compute the cell similarity.
Definition at line 2320 of file fe_values.h.
|
protected |
Storage for the mapping object.
Definition at line 2344 of file fe_values.h.
|
protected |
Store the finite element for later use.
Definition at line 2349 of file fe_values.h.
|
protected |
Internal data of mapping.
Definition at line 2355 of file fe_values.h.
|
protected |
Internal data of finite element.
Definition at line 2360 of file fe_values.h.
|
protected |
An enum variable that can store different states of the current cell in comparison to the previously visited cell. If wanted, additional states can be checked here and used in one of the methods used during reinit.
Definition at line 2377 of file fe_values.h.
|
private |
A cache for all possible FEValuesViews objects.
Definition at line 2403 of file fe_values.h.