Reference documentation for deal.II version 8.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
FEValues< dim, spacedim > Class Template Reference

#include <fe_values.h>

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

Public Member Functions

 FEValues (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
 
 FEValues (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
 
template<class DH , bool level_dof_access>
void reinit (const TriaIterator< DoFCellAccessor< DH, level_dof_access > > cell)
 
void reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell)
 
const Quadrature< dim > & get_quadrature () const
 
std::size_t memory_consumption () const
 
const FEValues< dim, spacedim > & get_present_fe_values () const
 
- Public Member Functions inherited from 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)
 
 ~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.")
 
 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 (ExcCannotInitializeField)
 
 DeclException0 (ExcInvalidUpdateFlag)
 
 DeclException0 (ExcFEDontMatch)
 
 DeclException0 (ExcFENotPrimitive)
 
const doubleshape_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
 
void get_function_values (const InputVector &fe_function, std::vector< number > &values) const
 
void get_function_values (const InputVector &fe_function, std::vector< Vector< number > > &values) const
 
void get_function_values (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< number > &values) const
 
void get_function_values (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< Vector< number > > &values) const
 
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
 
void get_function_gradients (const InputVector &fe_function, std::vector< Tensor< 1, spacedim > > &gradients) const
 
void get_function_gradients (const InputVector &fe_function, std::vector< std::vector< Tensor< 1, spacedim > > > &gradients) const
 
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
 
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
 
void get_function_grads (const InputVector &fe_function, std::vector< Tensor< 1, spacedim > > &gradients) const DEAL_II_DEPRECATED
 
void get_function_grads (const InputVector &fe_function, std::vector< std::vector< Tensor< 1, spacedim > > > &gradients) const DEAL_II_DEPRECATED
 
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
 
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
 
void get_function_hessians (const InputVector &fe_function, std::vector< Tensor< 2, spacedim > > &hessians) const
 
void get_function_hessians (const InputVector &fe_function, std::vector< std::vector< Tensor< 2, spacedim > > > &hessians, bool quadrature_points_fastest=false) const
 
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
 
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
 
void get_function_2nd_derivatives (const InputVector &, std::vector< Tensor< 2, spacedim > > &) const DEAL_II_DEPRECATED
 
void get_function_2nd_derivatives (const InputVector &, std::vector< std::vector< Tensor< 2, spacedim > > > &, bool=false) const DEAL_II_DEPRECATED
 
void get_function_laplacians (const InputVector &fe_function, std::vector< number > &laplacians) const
 
void get_function_laplacians (const InputVector &fe_function, std::vector< Vector< number > > &laplacians) const
 
void get_function_laplacians (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< number > &laplacians) const
 
void get_function_laplacians (const InputVector &fe_function, const VectorSlice< const std::vector< types::global_dof_index > > &indices, std::vector< Vector< number > > &laplacians) const
 
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
 
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
 
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
 
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
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (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)
 

Static Public Attributes

static const unsigned int integral_dimension = dim
 
- Static Public Attributes inherited from FEValuesBase< dim, spacedim >
static const unsigned int dimension
 
static const unsigned int space_dimension
 

Private Member Functions

void initialize (const UpdateFlags update_flags)
 
void do_reinit ()
 

Private Attributes

const Quadrature< dim > quadrature
 

Additional Inherited Members

- Public Attributes inherited from FEValuesBase< dim, spacedim >
const unsigned int n_quadrature_points
 
const unsigned int dofs_per_cell
 
- Protected Types inherited from FEValuesData< dim, spacedim >
typedef Table< 2, doubleShapeVector
 
typedef std::vector
< std::vector< Tensor
< 1, spacedim > > > 
GradientVector
 
typedef std::vector
< std::vector< Tensor
< 2, spacedim > > > 
HessianVector
 
- Protected Member Functions inherited from FEValuesBase< dim, spacedim >
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)
 
- Protected Member Functions inherited from FEValuesData< dim, spacedim >
void initialize (const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
 
- Protected Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (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 inherited from FEValuesBase< dim, spacedim >
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
 
- Protected Attributes inherited from FEValuesData< dim, spacedim >
ShapeVector shape_values
 
GradientVector shape_gradients
 
HessianVector shape_hessians
 
std::vector< doubleJxW_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 intshape_function_to_row_table
 
UpdateFlags update_flags
 

Detailed Description

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

Finite element evaluated in quadrature points of a cell.

This function implements the initialization routines for FEValuesBase, if values in quadrature points of a cell are needed. For further documentation see this class.

Author
Wolfgang Bangerth, 1998, Guido Kanschat, 2001

Definition at line 46 of file matrix_tools.h.

Constructor & Destructor Documentation

template<int dim, int spacedim>
FEValues< dim, spacedim >::FEValues ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Constructor. Gets cell independent data from mapping and finite element objects, matching the quadrature rule and update flags.

template<int dim, int spacedim>
FEValues< dim, spacedim >::FEValues ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Constructor. Uses MappingQ1 implicitly.

Member Function Documentation

template<int dim, int spacedim>
template<class DH , bool level_dof_access>
void FEValues< dim, spacedim >::reinit ( const TriaIterator< DoFCellAccessor< DH, level_dof_access > >  cell)

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a DoFHandler object", and the finite element associated with this object. It is assumed that the finite element used by the given cell is also the one used by this FEValues object.

template<int dim, int spacedim>
void FEValues< dim, spacedim >::reinit ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell)

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a Triangulation object", and the given finite element. Since iterators into triangulation alone only convey information about the geometry of a cell, but not about degrees of freedom possibly associated with this cell, you will not be able to call some functions of this class if they need information about degrees of freedom. These functions are, above all, the get_function_value/gradients/hessians/laplacians functions. If you want to call these functions, you have to call the reinit variants that take iterators into DoFHandler or other DoF handler type objects.

template<int dim, int spacedim>
const Quadrature<dim>& FEValues< dim, spacedim >::get_quadrature ( ) const

Return a reference to the copy of the quadrature formula stored by this object.

template<int dim, int spacedim>
std::size_t FEValues< dim, spacedim >::memory_consumption ( ) const

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

template<int dim, int spacedim>
const FEValues<dim,spacedim>& FEValues< dim, spacedim >::get_present_fe_values ( ) const

Return a reference to this very object.

Though it seems that it is not very useful, this function is there to provide capability to the hpFEValues class, in which case it provides the FEValues object for the present cell (remember that for hp finite elements, the actual FE object used may change from cell to cell, so we also need different FEValues objects for different cells; once you reinitialize the hpFEValues object for a specific cell, it retrieves the FEValues object for the FE on that cell and returns it through a function of the same name as this one; this function here therefore only provides the same interface so that one can templatize on FEValues/hpFEValues).

template<int dim, int spacedim>
void FEValues< dim, spacedim >::initialize ( const UpdateFlags  update_flags)
private

Do work common to the two constructors.

template<int dim, int spacedim>
void FEValues< dim, spacedim >::do_reinit ( )
private

The reinit() functions do only that part of the work that requires knowledge of the type of iterator. After setting present_cell(), they pass on to this function, which does the real work, and which is independent of the actual type of the cell iterator.

Member Data Documentation

template<int dim, int spacedim>
const unsigned int FEValues< dim, spacedim >::integral_dimension = dim
static

Dimension of the object over which we integrate. For the present class, this is equal to dim.

Definition at line 2435 of file fe_values.h.

template<int dim, int spacedim>
const Quadrature<dim> FEValues< dim, spacedim >::quadrature
private

Store a copy of the quadrature formula here.

Definition at line 2507 of file fe_values.h.


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