FEValuesBase< dim, spacedim > Class Template Reference
[Finite element access/FEValues classes]

Inheritance diagram for FEValuesBase< dim, spacedim >:
Inheritance graph
[legend]

List of all members.

Classes

class  CellIteratorBase
class  ExcAccessToUninitializedField
class  ExcCannotInitializeField
class  ExcFEDontMatch
class  ExcFENotPrimitive
class  ExcInvalidUpdateFlag
class  ExcShapeFunctionNotPrimitive

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 ()
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 Tensor< 2, spacedim > & jacobian (const unsigned int quadrature_point) const
const std::vector< Tensor
< 2, spacedim > > & 
get_jacobians () const
const Tensor< 3, spacedim > & jacobian_grad (const unsigned int quadrature_point) const
const std::vector< Tensor
< 3, spacedim > > & 
get_jacobian_grads () const
const Tensor< 2, spacedim > & inverse_jacobian (const unsigned int quadrature_point) const
const std::vector< Tensor
< 2, spacedim > > & 
get_inverse_jacobians () const
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
const std::vector< Point
< spacedim > > & 
get_cell_normal_vectors () const
const Point< spacedim > & cell_normal_vector (const unsigned int i) const
CellSimilarity::Similarity get_cell_similarity () const
unsigned int memory_consumption () const
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
ShapeAccess Access to shape function values



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, 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
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
FunctionAccess Access to values of global finite element functions



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< unsigned int > > &indices, std::vector< number > &values) const
template<class InputVector , typename number >
void get_function_values (const InputVector &fe_function, const VectorSlice< const std::vector< unsigned int > > &indices, std::vector< Vector< number > > &values) const
template<class InputVector , typename number >
void get_function_values (const InputVector &fe_function, const VectorSlice< const std::vector< unsigned int > > &indices, std::vector< std::vector< number > > &values, const bool quadrature_points_fastest) const
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, std::vector< Tensor< 1, spacedim > > &gradients) const
template<class InputVector >
void get_function_grads (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_grads (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< unsigned int > > &indices, std::vector< Tensor< 1, spacedim > > &gradients) const
template<class InputVector >
void get_function_grads (const InputVector &fe_function, const VectorSlice< const std::vector< unsigned int > > &indices, std::vector< Tensor< 1, spacedim > > &gradients) const
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, const VectorSlice< const std::vector< unsigned int > > &indices, 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, const VectorSlice< const std::vector< unsigned int > > &indices, std::vector< std::vector< Tensor< 1, spacedim > > > &gradients, bool quadrature_points_fastest=false) const
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< unsigned int > > &indices, std::vector< Tensor< 2, spacedim > > &hessians) const
template<class InputVector >
void get_function_hessians (const InputVector &fe_function, const VectorSlice< const std::vector< unsigned int > > &indices, 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
template<class InputVector >
void get_function_2nd_derivatives (const InputVector &, std::vector< std::vector< Tensor< 2, spacedim > > > &, bool=false) const
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< unsigned int > > &indices, std::vector< number > &laplacians) const
template<class InputVector , typename number >
void get_function_laplacians (const InputVector &fe_function, const VectorSlice< const std::vector< unsigned int > > &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< unsigned int > > &indices, std::vector< std::vector< number > > &laplacians, bool quadrature_points_fastest=false) const

Public Attributes

const unsigned int n_quadrature_points
const unsigned int dofs_per_cell

Protected Member Functions

UpdateFlags compute_update_flags (const UpdateFlags update_flags) const
void check_cell_similarity (const typename Triangulation< dim, spacedim >::cell_iterator &cell)

Protected Attributes

std::auto_ptr< const
CellIteratorBase
present_cell
const SmartPointer< const
Mapping< dim, spacedim > > 
mapping
const SmartPointer< const
FiniteElement< dim, spacedim > > 
fe
SmartPointer< typename Mapping
< dim, spacedim >
::InternalDataBase > 
mapping_data
SmartPointer< typename Mapping
< dim, spacedim >
::InternalDataBase > 
fe_data
CellSimilarity::Similarity cell_similarity

Private Member Functions

 FEValuesBase (const FEValuesBase &)
FEValuesBaseoperator= (const FEValuesBase &)

Private Attributes

internal::FEValuesViews::Cache
< dim, spacedim > 
fe_values_views_cache

Friends

class FEValuesViews::Scalar
class FEValuesViews::Vector

Detailed Description

template<int dim, int spacedim>
class FEValuesBase< dim, spacedim >

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.

General usage

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:

 FEValues values (mapping, finite_element, quadrature, flags);
 for (cell = dof_handler.begin_active();
      cell != dof_handler.end();
      ++cell)
   {
     values.reinit(cell);
     for (unsigned int q=0; q<quadrature.size(); ++q)
       for (unsigned int i=0; i<finite_element.dofs_per_cell; ++i)
         for (unsigned int j=0; j<finite_element.dofs_per_cell; ++j)
         A(i,j) += fe_values.shape_value(i,q) *
                   fe_values.shape_value(j,q) *
                   fe_values.JxW(q);
     ...
   }

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.

Member functions

The functions of this class fall into different cathegories:

UpdateFlags

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.

Author:
Wolfgang Bangerth, 1998, 2003, Guido Kanschat, 2001

Constructor & Destructor Documentation

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

template<int dim, int spacedim>
FEValuesBase< dim, spacedim >::~FEValuesBase (  ) 

Destructor.

template<int dim, int spacedim>
FEValuesBase< dim, spacedim >::FEValuesBase ( const FEValuesBase< dim, spacedim > &   )  [private]

Copy constructor. Since objects of this class are not copyable, we make it private, and also do not implement it.


Member Function Documentation

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

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

template<int dim, int spacedim>
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 computed
  • point_no Number of the quadrature point at which function is to be computed
template<int dim, int spacedim>
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 computed
  • point_no Number of the quadrature point at which function is to be computed
  • component vector component to be computed
template<int dim, int spacedim>
const Tensor<1,spacedim>& FEValuesBase< dim, spacedim >::shape_grad ( const unsigned int  function,
const unsigned int  quadrature_point 
) const

Compute the gradient of the ith shape function at the jth 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. 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 throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_grad_component() function.

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

template<int dim, int spacedim>
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_noth shape function at the point_noth 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.

template<int dim, int spacedim>
const Tensor<2,spacedim>& FEValuesBase< dim, spacedim >::shape_2nd_derivative ( const unsigned int  function_no,
const unsigned int  point_no 
) const
template<int dim, int spacedim>
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.

template<int dim, int spacedim>
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
template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
std::vector< number > &  values 
) const [inline]

Returns the values of the finite element function characterized by fe_function restricted to the 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 an active one the interpolated function values are returned.

To get values of multi-component elements, there is another get_function_values() below, returning a vector of vectors of results.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. If it is a vector-valued one, then use the other get_function_values() function.

The function assumes that the values object already has the correct size.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
std::vector< Vector< number > > &  values 
) const [inline]

Access to vector valued finite element functions.

This function does the same as the other get_function_values(), but applied to multi-component elements.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< number > &  values 
) const [inline]

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.

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< Vector< number > > &  values 
) const [inline]

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.

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< std::vector< number > > &  values,
const bool  quadrature_points_fastest 
) const [inline]

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.

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
std::vector< Tensor< 1, spacedim > > &  gradients 
) const [inline]

Compute the gradients of the finite element function characterized by fe_function restricted to cell at the quadrature points.

If the present cell is not an active one the interpolated function values are returned.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. If it is a vector-valued one, then use the other get_function_gradients() function.

The function assumes that the gradients object already has the right size.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

The output are the gradients of the function represented by these DoF values, as computed in real space (as opposed to on the unit cell).

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_grads ( const InputVector &  fe_function,
std::vector< Tensor< 1, spacedim > > &  gradients 
) const [inline]
template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
std::vector< std::vector< Tensor< 1, spacedim > > > &  gradients 
) const [inline]

Compute the gradients of the finite element function characterized by fe_function restricted to cell at the quadrature points.

If the present cell is not an active one the interpolated function values are returned.

The function assumes that the gradients object already has the right size.

This function does the same as the other get_function_values(), but applied to multi-component elements.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

The output are the gradients of the function represented by these DoF values, as computed in real space (as opposed to on the unit cell).

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_grads ( const InputVector &  fe_function,
std::vector< std::vector< Tensor< 1, spacedim > > > &  gradients 
) const [inline]
template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< Tensor< 1, spacedim > > &  gradients 
) const [inline]

Function gradient access with more flexibility. see get_function_values() with corresponding arguments.

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_grads ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< Tensor< 1, spacedim > > &  gradients 
) const [inline]
template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< std::vector< Tensor< 1, spacedim > > > &  gradients,
bool  quadrature_points_fastest = false 
) const [inline]

Function gradient access with more flexibility. see get_function_values() with corresponding arguments.

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_grads ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< std::vector< Tensor< 1, spacedim > > > &  gradients,
bool  quadrature_points_fastest = false 
) const [inline]
template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
std::vector< Tensor< 2, spacedim > > &  hessians 
) const [inline]

Compute the tensor of second derivatives of the finite element function characterized by fe_function restricted to cell at the quadrature points.

The function assumes that the hessians object already has the correct size.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. If it is a vector-valued one, then use the other get_function_hessians() function.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

The output are the second derivatives of the function represented by these DoF values, as computed in real space (as opposed to on the unit cell).

template<int dim, int spacedim>
template<class InputVector >
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 [inline]

Compute the tensor of second derivatives of the finite element function characterized by fe_function restricted to cell at the quadrature points.

The function assumes that the hessians object already has the right size.

This function does the same as the other one with the same name, but applies to vector-valued finite elements.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

The output are the second derivatives of the function represented by these DoF values, as computed in real space (as opposed to on the unit cell).

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< Tensor< 2, spacedim > > &  hessians 
) const [inline]

Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< std::vector< Tensor< 2, spacedim > > > &  hessians,
bool  quadrature_points_fastest = false 
) const [inline]

Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_2nd_derivatives ( const InputVector &  ,
std::vector< Tensor< 2, spacedim > > &   
) const [inline]
template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_2nd_derivatives ( const InputVector &  ,
std::vector< std::vector< Tensor< 2, spacedim > > > &  ,
bool  = false 
) const [inline]
template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
std::vector< number > &  laplacians 
) const [inline]

Compute the (scalar) Laplacian of the finite element function characterized by fe_function restricted to cell at the quadrature points. The Laplacian output vector is equivalent to getting trace(hessians), where hessian would be the output of the get_function_hessians() function.

The function assumes that the laplacians object already has the correct size.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. If it is a vector-valued one, then use the other get_function_laplacians() function.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

The output are the traces of the second derivatives (i.e. Laplacians) of the function represented by these DoF values, as computed in real space (as opposed to on the unit cell).

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
std::vector< Vector< number > > &  laplacians 
) const [inline]

Compute the (scalar) Laplacian of the finite element function characterized by fe_function restricted to cell at the quadrature points. The Laplacian output vector is equivalent to getting trace(hessians), with hessian corresponding to the output of the get_function_hessians() function.

The function assumes that the laplacians object already has the correct size.

This function does the same as the other one with the same name, but applies to vector-valued finite elements.

The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the sequential PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DofHandler object with which this FEValues object was last initialized.

The output are the traces of the second derivatives (i.e. Laplacians) of the function represented by these DoF values, as computed in real space (as opposed to on the unit cell).

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< number > &  laplacians 
) const [inline]

Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< Vector< number > > &  laplacians 
) const [inline]

Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.

template<int dim, int spacedim>
template<class InputVector , typename number >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
const VectorSlice< const std::vector< unsigned int > > &  indices,
std::vector< std::vector< number > > &  laplacians,
bool  quadrature_points_fastest = false 
) const [inline]

Access to the second derivatives of a function with more flexibility. see get_function_values() with corresponding arguments.

template<int dim, int spacedim>
const Point<spacedim>& FEValuesBase< dim, spacedim >::quadrature_point ( const unsigned int  i  )  const

Position of the ith quadrature point in real space.

template<int dim, int spacedim>
const std::vector<Point<spacedim> >& FEValuesBase< dim, spacedim >::get_quadrature_points (  )  const

Return a pointer to the vector of quadrature points.

template<int dim, int spacedim>
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 *ith 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 $dx, ds$ in the integral that we implement here by quadrature.

template<int dim, int spacedim>
const std::vector<double>& FEValuesBase< dim, spacedim >::get_JxW_values (  )  const

Pointer to the array holding the values returned by JxW().

template<int dim, int spacedim>
const Tensor<2,spacedim>& FEValuesBase< dim, spacedim >::jacobian ( const unsigned int  quadrature_point  )  const

Return the Jacobian of the transformation at the specified quadrature point, i.e. $J_{ij}=dx_i/d\hat x_j$

template<int dim, int spacedim>
const std::vector<Tensor<2,spacedim> >& FEValuesBase< dim, spacedim >::get_jacobians (  )  const

Pointer to the array holding the values returned by jacobian().

template<int dim, int spacedim>
const Tensor<3,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. $G_{ijk}=dJ_{jk}/d\hat x_i$.

template<int dim, int spacedim>
const std::vector<Tensor<3,spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_grads (  )  const

Pointer to the array holding the values returned by jacobian_grads().

template<int dim, int spacedim>
const Tensor<2,spacedim>& 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. $J_{ij}=d\hat x_i/dx_j$

template<int dim, int spacedim>
const std::vector<Tensor<2,spacedim> >& FEValuesBase< dim, spacedim >::get_inverse_jacobians (  )  const

Pointer to the array holding the values returned by inverse_jacobian().

template<int dim, int spacedim>
const Mapping<dim,spacedim>& FEValuesBase< dim, spacedim >::get_mapping (  )  const

Constant reference to the selected mapping object.

template<int dim, int spacedim>
const FiniteElement<dim,spacedim>& FEValuesBase< dim, spacedim >::get_fe (  )  const

Constant reference to the selected finite element object.

template<int dim, int spacedim>
UpdateFlags FEValuesBase< dim, spacedim >::get_update_flags (  )  const

Return the update flags set for this object.

template<int dim, int spacedim>
const Triangulation<dim,spacedim>::cell_iterator FEValuesBase< dim, spacedim >::get_cell (  )  const

Return a triangulation iterator to the current cell.

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

template<int dim, int spacedim>
const Point<spacedim>& FEValuesBase< dim, spacedim >::cell_normal_vector ( const unsigned int  i  )  const

Return the outward normal vector to the cell at the ith quadrature point. The length of the vector is normalized to one.

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

template<int dim, int spacedim>
unsigned int FEValuesBase< dim, spacedim >::memory_consumption (  )  const

Determine an estimate for the memory consumption (in bytes) of this object.

Reimplemented in FEValues< dim, spacedim >, FEFaceValuesBase< dim, spacedim >, and FEValues< dim >.

template<int dim, int spacedim>
UpdateFlags FEValuesBase< dim, spacedim >::compute_update_flags ( const UpdateFlags  update_flags  )  const [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.

template<int dim, int spacedim>
void FEValuesBase< dim, spacedim >::check_cell_similarity ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell  )  [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.

template<int dim, int spacedim>
FEValuesBase& FEValuesBase< dim, spacedim >::operator= ( const FEValuesBase< dim, spacedim > &   )  [private]

Copy operator. Since objects of this class are not copyable, we make it private, and also do not implement it.

Reimplemented from Subscriptor.


Friends And Related Function Documentation

template<int dim, int spacedim>
friend class FEValuesViews::Scalar [friend]

Make the view classes friends of this class, since they access internal data.

template<int dim, int spacedim>
friend class FEValuesViews::Vector [friend]

Member Data Documentation

template<int dim, int spacedim>
const unsigned int FEValuesBase< dim, spacedim >::n_quadrature_points

Number of quadrature points.

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

template<int dim, int spacedim>
std::auto_ptr<const CellIteratorBase> FEValuesBase< dim, spacedim >::present_cell [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.

template<int dim, int spacedim>
const SmartPointer<const Mapping<dim,spacedim> > FEValuesBase< dim, spacedim >::mapping [protected]

Storage for the mapping object.

template<int dim, int spacedim>
const SmartPointer<const FiniteElement<dim,spacedim> > FEValuesBase< dim, spacedim >::fe [protected]

Store the finite element for later use.

template<int dim, int spacedim>
SmartPointer<typename Mapping<dim,spacedim>::InternalDataBase> FEValuesBase< dim, spacedim >::mapping_data [protected]

Internal data of mapping.

template<int dim, int spacedim>
SmartPointer<typename Mapping<dim,spacedim>::InternalDataBase> FEValuesBase< dim, spacedim >::fe_data [protected]

Internal data of finite element.

template<int dim, int spacedim>
CellSimilarity::Similarity FEValuesBase< dim, spacedim >::cell_similarity [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.

template<int dim, int spacedim>
internal::FEValuesViews::Cache<dim,spacedim> FEValuesBase< dim, spacedim >::fe_values_views_cache [private]

A cache for all possible FEValuesViews objects.


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

deal.II documentation generated on Mon Nov 23 22:57:44 2009 by doxygen 1.6.1