Enumerations | Functions

VectorTools Namespace Reference

Enumerations

enum  NormType {
  mean, L1_norm, L2_norm, Lp_norm,
  Linfty_norm, H1_seminorm, H1_norm, W1p_seminorm,
  W1p_norm, W1infty_seminorm, W1infty_norm
}

Functions

::ExceptionBaseExcInvalidBoundaryIndicator ()
::ExceptionBaseExcNonInterpolatingFE ()
::ExceptionBaseExcNoComponentSelected ()
Interpolation and projection
template<class VECTOR , class DH >
void interpolate (const Mapping< DH::dimension, DH::space_dimension > &mapping, const DH &dof, const Function< DH::space_dimension > &function, VECTOR &vec)
template<class VECTOR , class DH >
void interpolate (const DH &dof, const Function< DH::space_dimension > &function, VECTOR &vec)
template<int dim, class InVector , class OutVector , int spacedim>
void interpolate (const DoFHandler< dim, spacedim > &dof_1, const DoFHandler< dim, spacedim > &dof_2, const FullMatrix< double > &transfer, const InVector &data_1, OutVector &data_2)
template<int dim, int spacedim, template< int, int > class DH, class VECTOR >
void interpolate_to_different_mesh (const DH< dim, spacedim > &dof1, const VECTOR &u1, const DH< dim, spacedim > &dof2, VECTOR &u2)
template<int dim, int spacedim, template< int, int > class DH, class VECTOR >
void interpolate_to_different_mesh (const DH< dim, spacedim > &dof1, const VECTOR &u1, const DH< dim, spacedim > &dof2, const ConstraintMatrix &constraints, VECTOR &u2)
template<int dim, int spacedim, template< int, int > class DH, class VECTOR >
void interpolate_to_different_mesh (const InterGridMap< DH< dim, spacedim > > &intergridmap, const VECTOR &u1, const ConstraintMatrix &constraints, VECTOR &u2)
template<int dim, class VECTOR , int spacedim>
void project (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &constraints, const Quadrature< dim > &quadrature, const Function< spacedim > &function, VECTOR &vec, const bool enforce_zero_boundary=false, const Quadrature< dim-1 > &q_boundary=(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)), const bool project_to_boundary_first=false)
template<int dim, class VECTOR , int spacedim>
void project (const DoFHandler< dim, spacedim > &dof, const ConstraintMatrix &constraints, const Quadrature< dim > &quadrature, const Function< spacedim > &function, VECTOR &vec, const bool enforce_zero_boundary=false, const Quadrature< dim-1 > &q_boundary=(dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)), const bool project_to_boundary_first=false)
template<class DH >
void interpolate_boundary_values (const Mapping< DH::dimension, DH::space_dimension > &mapping, const DH &dof, const typename FunctionMap< DH::space_dimension >::type &function_map, std::map< unsigned int, double > &boundary_values, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
void interpolate_boundary_values (const Mapping< DH::dimension, DH::space_dimension > &mapping, const DH &dof, const types::boundary_id_t boundary_component, const Function< DH::space_dimension > &boundary_function, std::map< unsigned int, double > &boundary_values, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
void interpolate_boundary_values (const DH &dof, const types::boundary_id_t boundary_component, const Function< DH::space_dimension > &boundary_function, std::map< unsigned int, double > &boundary_values, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
void interpolate_boundary_values (const DH &dof, const typename FunctionMap< DH::space_dimension >::type &function_map, std::map< unsigned int, double > &boundary_values, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
void interpolate_boundary_values (const Mapping< DH::dimension, DH::space_dimension > &mapping, const DH &dof, const typename FunctionMap< DH::space_dimension >::type &function_map, ConstraintMatrix &constraints, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
void interpolate_boundary_values (const Mapping< DH::dimension, DH::space_dimension > &mapping, const DH &dof, const types::boundary_id_t boundary_component, const Function< DH::space_dimension > &boundary_function, ConstraintMatrix &constraints, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
void interpolate_boundary_values (const DH &dof, const types::boundary_id_t boundary_component, const Function< DH::space_dimension > &boundary_function, ConstraintMatrix &constraints, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
void interpolate_boundary_values (const DH &dof, const typename FunctionMap< DH::space_dimension >::type &function_map, ConstraintMatrix &constraints, const std::vector< bool > &component_mask=std::vector< bool >())
template<int dim, int spacedim>
void project_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_functions, const Quadrature< dim-1 > &q, std::map< unsigned int, double > &boundary_values, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
template<int dim, int spacedim>
void project_boundary_values (const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_function, const Quadrature< dim-1 > &q, std::map< unsigned int, double > &boundary_values, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
template<int dim, int spacedim>
void project_boundary_values (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_functions, const Quadrature< dim-1 > &q, ConstraintMatrix &constraints, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
template<int dim, int spacedim>
void project_boundary_values (const DoFHandler< dim, spacedim > &dof, const typename FunctionMap< spacedim >::type &boundary_function, const Quadrature< dim-1 > &q, ConstraintMatrix &constraints, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
template<int dim>
void project_boundary_values_curl_conforming (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim > &boundary_function, const types::boundary_id_t boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
template<int dim>
void project_boundary_values_curl_conforming (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim > &boundary_function, const types::boundary_id_t boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection)
template<int dim>
void project_boundary_values_div_conforming (const DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim > &boundary_function, const types::boundary_id_t boundary_component, ConstraintMatrix &constraints, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
template<int dim>
void project_boundary_values_div_conforming (const hp::DoFHandler< dim > &dof_handler, const unsigned int first_vector_component, const Function< dim > &boundary_function, const types::boundary_id_t boundary_component, ConstraintMatrix &constraints, const hp::MappingCollection< dim, dim > &mapping_collection=hp::StaticMappingQ1< dim >::mapping_collection)
template<int dim, template< int, int > class DH, int spacedim>
void compute_no_normal_flux_constraints (const DH< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id_t > &boundary_ids, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
Assembling of right hand sides
template<int dim, int spacedim>
void create_right_hand_side (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_right_hand_side (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_right_hand_side (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_right_hand_side (const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &orientation, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &orientation, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &orientation, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_point_source_vector (const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, const Point< dim > &orientation, Vector< double > &rhs_vector)
template<int dim, int spacedim>
void create_boundary_right_hand_side (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim-1 > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id_t > &boundary_indicators=std::set< types::boundary_id_t >())
template<int dim, int spacedim>
void create_boundary_right_hand_side (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim-1 > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id_t > &boundary_indicators=std::set< types::boundary_id_t >())
template<int dim, int spacedim>
void create_boundary_right_hand_side (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim-1 > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id_t > &boundary_indicators=std::set< types::boundary_id_t >())
template<int dim, int spacedim>
void create_boundary_right_hand_side (const hp::DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim-1 > &q, const Function< spacedim > &rhs, Vector< double > &rhs_vector, const std::set< types::boundary_id_t > &boundary_indicators=std::set< types::boundary_id_t >())
Evaluation of functions

and errors

template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim > *weight=0, const double exponent=2.)
template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim > &exact_solution, OutVector &difference, const Quadrature< dim > &q, const NormType &norm, const Function< spacedim > *weight=0, const double exponent=2.)
template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim > *weight=0, const double exponent=2.)
template<int dim, class InVector , class OutVector , int spacedim>
void integrate_difference (const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim > &exact_solution, OutVector &difference, const hp::QCollection< dim > &q, const NormType &norm, const Function< spacedim > *weight=0, const double exponent=2.)
template<int dim, class InVector , int spacedim>
void point_difference (const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim > &exact_solution, Vector< double > &difference, const Point< spacedim > &point)
template<int dim, class InVector , int spacedim>
void point_difference (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Function< spacedim > &exact_solution, Vector< double > &difference, const Point< spacedim > &point)
template<int dim, class InVector , int spacedim>
void point_value (const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point, Vector< double > &value)
template<int dim, class InVector , int spacedim>
void point_value (const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point, Vector< double > &value)
template<int dim, class InVector , int spacedim>
double point_value (const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point)
template<int dim, class InVector , int spacedim>
double point_value (const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point)
template<int dim, class InVector , int spacedim>
void point_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point, Vector< double > &value)
template<int dim, class InVector , int spacedim>
void point_value (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point, Vector< double > &value)
template<int dim, class InVector , int spacedim>
double point_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point)
template<int dim, class InVector , int spacedim>
double point_value (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point)
void subtract_mean_value (Vector< double > &v, const std::vector< bool > &p_select)
template<int dim, class InVector , int spacedim>
double compute_mean_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const InVector &v, const unsigned int component)
template<int dim, class InVector , int spacedim>
double compute_mean_value (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const InVector &v, const unsigned int component)

Detailed Description

Provide a namespace which offers some operations on vectors. Amoung these are assembling of standard vectors, integration of the difference of a finite element solution and a continuous function, interpolations and projections of continuous functions to the finite element space and other operations.

Note:
There exist two versions of almost each function. One with a Mapping argument and one without. If a code uses a mapping different from MappingQ1 the functions with mapping argument should be used. Code that uses only MappingQ1 may also use the functions without Mapping argument. Each of these latter functions create a MappingQ1 object and just call the respective functions with that object as mapping argument. The functions without Mapping argument still exist to ensure backward compatibility. Nevertheless it is advised to change the user's codes to store a specific Mapping object and to use the functions that take this Mapping object as argument. This gives the possibility to easily extend the user codes to work also on mappings of higher degree, this just by exchanging MappingQ1 by, for example, a MappingQ or another Mapping object of interest.

Description of operations

This collection of methods offers the following operations:

All functions use the finite element given to the DoFHandler object the last time that the degrees of freedom were distributed over the triangulation. Also, if access to an object describing the exact form of the boundary is needed, the pointer stored within the triangulation object is accessed.

Note:
Instantiations for this template are provided for some vector types, in particular Vector<float>, Vector<double>, BlockVector<float>, BlockVector<double>; others can be generated in application code (see the section on Template instantiations in the manual).
Author:
Wolfgang Bangerth, Ralf Hartmann, Guido Kanschat, 1998, 1999, 2000, 2001

Enumeration Type Documentation

Denote which norm/integral is to be computed by the integrate_difference() function of this class. The following possibilities are implemented:

Enumerator:
mean 

The function or difference of functions is integrated on each cell.

L1_norm 

The absolute value of the function is integrated.

L2_norm 

The square of the function is integrated and the the square root of the result is computed on each cell.

Lp_norm 

The absolute value to the pth power is integrated and the pth root is computed on each cell. The exponent p is the last parameter of the function.

Linfty_norm 

The maximum absolute value of the function.

H1_seminorm 

L2_norm of the gradient.

H1_norm 

The square of this norm is the square of the L2_norm plus the square of the H1_seminorm.

W1p_seminorm 

Lp_norm of the gradient.

W1p_norm 

same as H1_norm for Lp.

W1infty_seminorm 

Linfty_norm of the gradient.

W1infty_norm 

same as H1_norm for Linfty.

Definition at line 343 of file vectors.h.


Function Documentation

template<class VECTOR , class DH >
void VectorTools::interpolate ( const Mapping< DH::dimension, DH::space_dimension > &  mapping,
const DH &  dof,
const Function< DH::space_dimension > &  function,
VECTOR &  vec 
)

Compute the interpolation of function at the support points to the finite element space described by the Triangulation and FiniteElement object with which the given DoFHandler argument is initialized. It is assumed that the number of components of function matches that of the finite element used by dof.

Note that you may have to call hanging_nodes.distribute(vec) with the hanging nodes from space dof afterwards, to make the result continuous again.

The template argument DH may either be of type DoFHandler or hp::DoFHandler.

See the general documentation of this class for further information.

Todo:
The mapping argument should be replaced by a hp::MappingCollection in case of a hp::DoFHandler.
template<class VECTOR , class DH >
void VectorTools::interpolate ( const DH &  dof,
const Function< DH::space_dimension > &  function,
VECTOR &  vec 
)

Calls the interpolate() function above with mapping=MappingQ1<dim>@().

template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::interpolate ( const DoFHandler< dim, spacedim > &  dof_1,
const DoFHandler< dim, spacedim > &  dof_2,
const FullMatrix< double > &  transfer,
const InVector &  data_1,
OutVector &  data_2 
)

Interpolate different finite element spaces. The interpolation of vector data_1 is executed from the FE space represented by dof_1 to the vector data_2 on FE space dof_2. The interpolation on each cell is represented by the matrix transfer. Curved boundaries are neglected so far.

Note that you may have to call hanging_nodes.distribute(data_2) with the hanging nodes from space dof_2 afterwards, to make the result continuous again.

Note:
Instantiations for this template are provided for some vector types (see the general documentation of the class), but only the same vector for InVector and OutVector. Other combinations must be instantiated by hand.
template<int dim, int spacedim, template< int, int > class DH, class VECTOR >
void VectorTools::interpolate_to_different_mesh ( const DH< dim, spacedim > &  dof1,
const VECTOR &  u1,
const DH< dim, spacedim > &  dof2,
VECTOR &  u2 
)

Gives the interpolation of a dof1-function u1 to a dof2-function u2, where dof1 and dof2 represent different triangulations with a common coarse grid.

dof1 and dof2 need to have the same finite element discretization.

Note that for continuous elements on grids with hanging nodes (i.e. locally refined grids) this function does not give the expected output. Indeed, the resulting output vector does not necessarily respect continuity requirements at hanging nodes, due to local cellwise interpolation.

For this case (continuous elements on grids with hanging nodes), please use the interpolate_to_different_mesh function with an additional ConstraintMatrix argument, see below, or make the field conforming yourself by calling the ConstraintsMatrix::distribute function of your hanging node constraints object.

template<int dim, int spacedim, template< int, int > class DH, class VECTOR >
void VectorTools::interpolate_to_different_mesh ( const DH< dim, spacedim > &  dof1,
const VECTOR &  u1,
const DH< dim, spacedim > &  dof2,
const ConstraintMatrix constraints,
VECTOR &  u2 
)

Gives the interpolation of a dof1-function u1 to a dof2-function u2, where dof1 and dof2 represent different triangulations with a common coarse grid.

dof1 and dof2 need to have the same finite element discretization.

constraints is a hanging node constraints object corresponding to dof2. This object is particularly important when interpolating onto continuous elements on grids with hanging nodes (locally refined grids): Without it - due to cellwise interpolation - the resulting output vector does not necessarily respect continuity requirements at hanging nodes.

template<int dim, int spacedim, template< int, int > class DH, class VECTOR >
void VectorTools::interpolate_to_different_mesh ( const InterGridMap< DH< dim, spacedim > > &  intergridmap,
const VECTOR &  u1,
const ConstraintMatrix constraints,
VECTOR &  u2 
)

The same function as above, but takes an InterGridMap object directly as a parameter. Useful for interpolating several vectors at the same time.

intergridmap has to be initialized via InterGridMap::make_mapping pointing from a source DoFHandler to a destination DoFHandler.

template<int dim, class VECTOR , int spacedim>
void VectorTools::project ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const ConstraintMatrix constraints,
const Quadrature< dim > &  quadrature,
const Function< spacedim > &  function,
VECTOR &  vec,
const bool  enforce_zero_boundary = false,
const Quadrature< dim-1 > &  q_boundary = (dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)),
const bool  project_to_boundary_first = false 
)

Compute the projection of function to the finite element space.

By default, projection to the boundary and enforcement of zero boundary values are disabled. The ordering of arguments to this function is such that you need not give a second quadrature formula if you don't want to project to the boundary first, but that you must if you want to do so.

This function needs the mass matrix of the finite element space on the present grid. To this end, the mass matrix is assembled exactly using MatrixTools::create_mass_matrix. This function performs numerical quadrature using the given quadrature rule; you should therefore make sure that the given quadrature formula is also sufficient for the integration of the mass matrix.

See the general documentation of this class for further information.

In 1d, the default value of the boundary quadrature formula is an invalid object since integration on the boundary doesn't happen in 1d.

template<int dim, class VECTOR , int spacedim>
void VectorTools::project ( const DoFHandler< dim, spacedim > &  dof,
const ConstraintMatrix constraints,
const Quadrature< dim > &  quadrature,
const Function< spacedim > &  function,
VECTOR &  vec,
const bool  enforce_zero_boundary = false,
const Quadrature< dim-1 > &  q_boundary = (dim > 1?QGauss< dim-1 >(2):Quadrature< dim-1 >(0)),
const bool  project_to_boundary_first = false 
)

Calls the project() function, see above, with mapping=MappingQ1<dim>().

template<class DH >
void VectorTools::interpolate_boundary_values ( const Mapping< DH::dimension, DH::space_dimension > &  mapping,
const DH &  dof,
const typename FunctionMap< DH::space_dimension >::type &  function_map,
std::map< unsigned int, double > &  boundary_values,
const std::vector< bool > &  component_mask = std::vector< bool >() 
)

Prepare Dirichlet boundary conditions. Make up the list of degrees of freedom subject to Dirichlet boundary conditions and the values to be assigned to them, by interpolation around the boundary. If the boundary_values contained values before, the new ones are added, or the old ones overwritten if a node of the boundary part to be used was already in the map of boundary values.

The parameter boundary_component corresponds to the number boundary_indicator of the face. types::internal_face_boundary_id is an illegal value, since it is reserved for interior faces.

The flags in the last parameter, component_mask denote which components of the finite element space shall be interpolated. If it is left as specified by the default value (i.e. an empty array), all components are interpolated. If it is different from the default value, it is assumed that the number of entries equals the number of components in the boundary functions and the finite element, and those components in the given boundary function will be used for which the respective flag was set in the component mask.

It is assumed that the number of components of the function in boundary_function matches that of the finite element used by dof.

If the finite element used has shape functions that are non-zero in more than one component (in deal.II speak: they are non-primitive), then these components can presently not be used for interpolating boundary values. Thus, the elements in the component mask corresponding to the components of these non-primitive shape functions must be false.

See the general doc for more information.

template<class DH >
void VectorTools::interpolate_boundary_values ( const Mapping< DH::dimension, DH::space_dimension > &  mapping,
const DH &  dof,
const types::boundary_id_t  boundary_component,
const Function< DH::space_dimension > &  boundary_function,
std::map< unsigned int, double > &  boundary_values,
const std::vector< bool > &  component_mask = std::vector< bool >() 
)
Deprecated:
This function exists mainly for backward compatibility.

Same function as above, but taking only one pair of boundary indicator and corresponding boundary function. Calls the other function with remapped arguments.

template<class DH >
void VectorTools::interpolate_boundary_values ( const DH &  dof,
const types::boundary_id_t  boundary_component,
const Function< DH::space_dimension > &  boundary_function,
std::map< unsigned int, double > &  boundary_values,
const std::vector< bool > &  component_mask = std::vector< bool >() 
)

Calls the other interpolate_boundary_values() function, see above, with mapping=MappingQ1<dim>().

template<class DH >
void VectorTools::interpolate_boundary_values ( const DH &  dof,
const typename FunctionMap< DH::space_dimension >::type &  function_map,
std::map< unsigned int, double > &  boundary_values,
const std::vector< bool > &  component_mask = std::vector< bool >() 
)

Calls the other interpolate_boundary_values() function, see above, with mapping=MappingQ1<dim>().

template<int dim, int spacedim>
void VectorTools::project_boundary_values ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const typename FunctionMap< spacedim >::type &  boundary_functions,
const Quadrature< dim-1 > &  q,
std::map< unsigned int, double > &  boundary_values,
std::vector< unsigned int >  component_mapping = std::vector< unsigned int >() 
)

Project a function to the boundary of the domain, using the given quadrature formula for the faces. If the boundary_values contained values before, the new ones are added, or the old one overwritten if a node of the boundary part to be projected on already was in the variable.

If component_mapping is empty, it is assumed that the number of components of boundary_function matches that of the finite element used by dof.

In 1d, projection equals interpolation. Therefore, interpolate_boundary_values is called.

  • boundary_values: the result of this function, a map containing all indices of degrees of freedom at the boundary (as covered by the boundary parts in boundary_functions) and the computed dof value for this degree of freedom.
  • component_mapping: if the components in boundary_functions and dof do not coincide, this vector allows them to be remapped. If the vector is not empty, it has to have one entry for each component in dof. This entry is the component number in boundary_functions that should be used for this component in dof. By default, no remapping is applied.
template<int dim, int spacedim>
void VectorTools::project_boundary_values ( const DoFHandler< dim, spacedim > &  dof,
const typename FunctionMap< spacedim >::type &  boundary_function,
const Quadrature< dim-1 > &  q,
std::map< unsigned int, double > &  boundary_values,
std::vector< unsigned int >  component_mapping = std::vector< unsigned int >() 
)

Calls the project_boundary_values() function, see above, with mapping=MappingQ1<dim>().

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector 
)

Create a right hand side vector. Prior content of the given rhs_vector vector is deleted.

See the general documentation of this class for further information.

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector 
)

Calls the create_right_hand_side() function, see above, with mapping=MappingQ1<dim>().

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_right_hand_side ( const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Create a right hand side vector for a point source at point p. Prior content of the given rhs_vector vector is deleted.

See the general documentation of this class for further information.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Calls the create_point_source_vector() function, see above, with mapping=MappingQ1<dim>().

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects. The function uses the default Q1 mapping object. Note that if your hp::DoFHandler uses any active fe index other than zero, then you need to call the function above that provides a mapping object for each active fe index.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  orientation,
Vector< double > &  rhs_vector 
)

Create a right hand side vector for a point source at point p for vector-valued finite elements. Prior content of the given rhs_vector vector is deleted.

See the general documentation of this class for further information.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  orientation,
Vector< double > &  rhs_vector 
)

Calls the create_point_source_vector() function for vector-valued finite elements, see above, with mapping=MappingQ1<dim>().

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  orientation,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_point_source_vector ( const hp::DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
const Point< dim > &  orientation,
Vector< double > &  rhs_vector 
)

Like the previous set of functions, but for hp objects. The function uses the default Q1 mapping object. Note that if your hp::DoFHandler uses any active fe index other than zero, then you need to call the function above that provides a mapping object for each active fe index.

template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim-1 > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id_t > &  boundary_indicators = std::set< types::boundary_id_t >() 
)

Create a right hand side vector from boundary forces. Prior content of the given rhs_vector vector is deleted.

See the general documentation of this class for further information.

template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim-1 > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id_t > &  boundary_indicators = std::set< types::boundary_id_t >() 
)

Calls the create_boundary_right_hand_side() function, see above, with mapping=MappingQ1<dim>().

template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim-1 > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id_t > &  boundary_indicators = std::set< types::boundary_id_t >() 
)

Same as the set of functions above, but for hp objects.

template<int dim, int spacedim>
void VectorTools::create_boundary_right_hand_side ( const hp::DoFHandler< dim, spacedim > &  dof,
const hp::QCollection< dim-1 > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector,
const std::set< types::boundary_id_t > &  boundary_indicators = std::set< types::boundary_id_t >() 
)

Calls the create_boundary_right_hand_side() function, see above, with a single Q1 mapping as collection. This function therefore will only work if the only active fe index in use is zero.

template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::integrate_difference ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim > &  exact_solution,
OutVector &  difference,
const Quadrature< dim > &  q,
const NormType &  norm,
const Function< spacedim > *  weight = 0,
const double  exponent = 2. 
)

Compute the error of the finite element solution. Integrate the difference between a reference function which is given as a continuous function object, and a finite element function.

The value of exponent is used for computing $L^p$-norms and $W^{1,p}$-norms.

The additional argument weight allows to evaluate weighted norms. The weight function may be scalar, establishing a weight in the domain for all components equally. This may be used, for instance, to only integrates over parts of the domain.

The weight function may also be vector-valued, with as many components as the finite element function: Then, different components get different weights. A typical application is when the error with respect to only one or a subset of the solution variables is to be computed, in which the other components would have weight values equal to zero. The ComponentSelectFunction class is particularly useful for this purpose.

The weight function is expected to be positive, but negative values are not filtered. By default, no weighting function is given, i.e. weight=1 in the whole domain for all vector components uniformly. Note that one often wants to compute the error in only one component of a solution vector (e.g. for the pressure in the Stokes system, when the solution vector also contains the velocity components). In these cases, the weight should be a mask, i.e., be a vector function for which individual components are either zero or one. This can easily be achieved using the ComponentSelectFunction class.

It is assumed that the number of components of the function exact_solution matches that of the finite element used by dof.

See the general documentation of this class for more information.

Note:
Instantiations for this template are provided for some vector types (see the general documentation of the class), but only for InVectors as in the documentation of the class, OutVector only Vector<double> and Vector<float>.
template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::integrate_difference ( const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim > &  exact_solution,
OutVector &  difference,
const Quadrature< dim > &  q,
const NormType &  norm,
const Function< spacedim > *  weight = 0,
const double  exponent = 2. 
)

Calls the integrate_difference() function, see above, with mapping=MappingQ1<dim>().

template<int dim, class InVector , class OutVector , int spacedim>
void VectorTools::integrate_difference ( const hp::DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim > &  exact_solution,
OutVector &  difference,
const hp::QCollection< dim > &  q,
const NormType &  norm,
const Function< spacedim > *  weight = 0,
const double  exponent = 2. 
)

Calls the integrate_difference() function, see above, with mapping=MappingQ1<dim>().

template<int dim, class InVector , int spacedim>
void VectorTools::point_difference ( const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim > &  exact_solution,
Vector< double > &  difference,
const Point< spacedim > &  point 
)

Point error evaluation. Find the first cell containing the given point and compute the difference of a (possibly vector-valued) finite element function and a continuous function (with as many vector components as the finite element) at this point.

This is a wrapper function using a Q1-mapping for cell boundaries to call the other point_difference() function.

template<int dim, class InVector , int spacedim>
void VectorTools::point_difference ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Function< spacedim > &  exact_solution,
Vector< double > &  difference,
const Point< spacedim > &  point 
)

Point error evaluation. Find the first cell containing the given point and compute the difference of a (possibly vector-valued) finite element function and a continuous function (with as many vector components as the finite element) at this point.

Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.

template<int dim, class InVector , int spacedim>
void VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the (vector) value of this function through the last argument.

This is a wrapper function using a Q1-mapping for cell boundaries to call the other point_difference() function.

template<int dim, class InVector , int spacedim>
void VectorTools::point_value ( const hp::DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Same as above for hp.

template<int dim, class InVector , int spacedim>
double VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point 
)

Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the value of this function.

Compared with the other function of the same name, this is a wrapper function using a Q1-mapping for cells.

This function is used in the "Possibilities for extensions" part of the results section of step-3.

template<int dim, class InVector , int spacedim>
double VectorTools::point_value ( const hp::DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point 
)

Same as above for hp.

template<int dim, class InVector , int spacedim>
void VectorTools::point_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Evaluate a possibly vector-valued finite element function defined by the given DoFHandler and nodal vector at the given point, and return the (vector) value of this function through the last argument.

Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.

template<int dim, class InVector , int spacedim>
void VectorTools::point_value ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
)

Same as above for hp.

template<int dim, class InVector , int spacedim>
double VectorTools::point_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point 
)

Evaluate a scalar finite element function defined by the given DoFHandler and nodal vector at the given point, and return the value of this function.

Compared with the other function of the same name, this function uses an arbitrary mapping to evaluate the difference.

template<int dim, class InVector , int spacedim>
double VectorTools::point_value ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point 
)

Same as above for hp.

void VectorTools::subtract_mean_value ( Vector< double > &  v,
const std::vector< bool > &  p_select 
)

Mean value operations Subtract the (algebraic) mean value from a vector. This function is most frequently used as a mean-value filter for Stokes: The pressure in Stokes' equations with only Dirichlet boundaries for the velocities is only determined up to a constant. This function allows to subtract the mean value of the pressure. It is usually called in a preconditioner and generates updates with mean value zero. The mean value is computed as the mean value of the degrees of freedom values as given by the input vector; they are not weighted by the area of cells, i.e. the mean is computed as $\sum_i v_i$, rather than as $\int_\Omega v(x) = \int_\Omega \sum_i v_i \phi_i(x)$. The latter can be obtained from the VectorTools::compute_mean_function, however.

Apart from the vector v to operate on, this function takes a boolean mask that has a true entry for every component for which the mean value shall be computed and later subtracted. The argument is used to denote which components of the solution vector correspond to the pressure, and avoid touching all other components of the vector, such as the velocity components.

Note:
In the context of using this function to filter out the kernel of an operator (such as the null space of the Stokes operator that consists of the constant pressures), this function only makes sense for finite elements for which the null space indeed consists of the vector $(1,1,\ldots,1)^T$. This is the case for example for the usual Lagrange elements where the sum of all shape functions equals the function that is constant one. However, it is not true for some other functions: for example, for the FE_DGP element (another valid choice for the pressure in Stokes discretizations), the first shape function on each cell is constant while further elements are $L_2$ orthogonal to it (on the reference cell); consequently, the sum of all shape functions is not equal to one, and the vector that is associated with the constant mode is not equal to $(1,1,\ldots,1)^T$. For such elements, a different procedure has to be used when subtracting the mean value.
template<int dim, class InVector , int spacedim>
double VectorTools::compute_mean_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  quadrature,
const InVector &  v,
const unsigned int  component 
)

Compute the mean value of one component of the solution.

This function integrates the chosen component over the whole domain and returns the result, i.e. it computes $\int_\Omega [u_h(x)]_c \; dx$ where $c$ is the vector component and $u_h$ is the function representation of the nodal vector given as fourth argument. The integral is evaluated numerically using the quadrature formula given as third argument.

This function is used in the "Possibilities for extensions" part of the results section of step-3.

Note:
The function is most often used when solving a problem whose solution is only defined up to a constant, for example a pure Neumann problem or the pressure in a Stokes or Navier-Stokes problem. In both cases, subtracting the mean value as computed by the current function, from the nodal vector does not generally yield the desired result of a finite element function with mean value zero. In fact, it only works for Lagrangian elements. For all other elements, you will need to compute the mean value and subtract it right inside the evaluation routine.
template<int dim, class InVector , int spacedim>
double VectorTools::compute_mean_value ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  quadrature,
const InVector &  v,
const unsigned int  component 
)

Calls the other compute_mean_value() function, see above, with mapping=MappingQ1<dim>().

::ExceptionBase& VectorTools::ExcInvalidBoundaryIndicator ( ) [static]

Exception

::ExceptionBase& VectorTools::ExcNonInterpolatingFE ( ) [static]

Exception

::ExceptionBase& VectorTools::ExcNoComponentSelected ( ) [static]

Exception

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

deal.II documentation generated on Wed May 23 2012 06:08:32 by doxygen 1.7.3