VectorTools Class Reference
[Numerical algorithms]

List of all members.

Classes

class  ExcInvalidBoundaryIndicator
class  ExcNoComponentSelected
class  ExcNonInterpolatingFE

Public Types

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

Static Public Member Functions

Interpolation and projection



template<class VECTOR , class DH >
static 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 >
static void interpolate (const DH &dof, const Function< DH::space_dimension > &function, VECTOR &vec)
template<int dim, class InVector , class OutVector , int spacedim>
static 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, class VECTOR , int spacedim>
static 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>
static 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 >
static 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 >
static void interpolate_boundary_values (const Mapping< DH::dimension, DH::space_dimension > &mapping, const DH &dof, const unsigned char 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 >
static void interpolate_boundary_values (const DH &dof, const unsigned char 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 >
static 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 >
static 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 >
static void interpolate_boundary_values (const Mapping< DH::dimension, DH::space_dimension > &mapping, const DH &dof, const unsigned char boundary_component, const Function< DH::space_dimension > &boundary_function, ConstraintMatrix &constraints, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
static void interpolate_boundary_values (const DH &dof, const unsigned char boundary_component, const Function< DH::space_dimension > &boundary_function, ConstraintMatrix &constraints, const std::vector< bool > &component_mask=std::vector< bool >())
template<class DH >
static 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>
static 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>
static 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>
static 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>
static 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, template< int, int > class DH, int spacedim>
static void compute_no_normal_flux_constraints (const DH< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< unsigned char > &boundary_ids, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
Assembling of right hand sides



template<int dim, int spacedim>
static 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>
static 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>
static 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>
static 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>
static 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>
static void create_point_source_vector (const DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
template<int dim, int spacedim>
static 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>
static void create_point_source_vector (const hp::DoFHandler< dim, spacedim > &dof, const Point< spacedim > &p, Vector< double > &rhs_vector)
template<int dim, int spacedim>
static 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< unsigned char > &boundary_indicators=std::set< unsigned char >())
template<int dim, int spacedim>
static 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< unsigned char > &boundary_indicators=std::set< unsigned char >())
template<int dim, int spacedim>
static 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< unsigned char > &boundary_indicators=std::set< unsigned char >())
template<int dim, int spacedim>
static 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< unsigned char > &boundary_indicators=std::set< unsigned char >())
Evaluation of functions

and errors



template<int dim, class InVector , class OutVector , int spacedim>
static 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>
static 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>
static 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>
static 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>
static 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>
static 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>
static 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>
static double point_value (const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point)
template<int dim, class InVector , int spacedim>
static 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>
static double point_value (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InVector &fe_function, const Point< spacedim > &point)



static void subtract_mean_value (Vector< double > &v, const std::vector< bool > &p_select)
template<int dim, class InVector , int spacedim>
static 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>
static double compute_mean_value (const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &quadrature, const InVector &v, const unsigned int component)

Detailed Description

Provide a class 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

Member Enumeration 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.


Member Function Documentation

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

Compute the interpolation of function at the support points to the finite element space. 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 >
static void VectorTools::interpolate ( const DH &  dof,
const Function< DH::space_dimension > &  function,
VECTOR &  vec 
) [inline, static]

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

template<int dim, class InVector , class OutVector , int spacedim>
static 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 
) [inline, static]

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, class VECTOR , int spacedim>
static 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 
) [inline, static]

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>
static 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 
) [inline, static]

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

template<class DH >
static 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 >() 
) [inline, static]

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. 255 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 >
static void VectorTools::interpolate_boundary_values ( const Mapping< DH::dimension, DH::space_dimension > &  mapping,
const DH &  dof,
const unsigned char  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 >() 
) [inline, static]
Deprecated:
This function is there 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 >
static void VectorTools::interpolate_boundary_values ( const DH &  dof,
const unsigned char  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 >() 
) [inline, static]

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

template<class DH >
static 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 >() 
) [inline, static]

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

template<class DH >
static 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,
ConstraintMatrix constraints,
const std::vector< bool > &  component_mask = std::vector< bool >() 
) [inline, static]

Insert the (algebraic) constraints due to Dirichlet boundary conditions to the ConstraintMatrix. This function makes 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 ConstraintMatrix constraints contained values or other constraints 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 list of constraints. This is handled by using inhomogeneous constraints. Please note that when combining adaptive meshes and this kind of constraints, the Dirichlet conditions should be set first, and then completed by hanging node constraints, in order to make sure that the discretization remains consistent.

The parameter boundary_component corresponds to the number boundary_indicator of the face. 255 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 >
static void VectorTools::interpolate_boundary_values ( const Mapping< DH::dimension, DH::space_dimension > &  mapping,
const DH &  dof,
const unsigned char  boundary_component,
const Function< DH::space_dimension > &  boundary_function,
ConstraintMatrix constraints,
const std::vector< bool > &  component_mask = std::vector< bool >() 
) [inline, static]
Deprecated:
This function is there 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 >
static void VectorTools::interpolate_boundary_values ( const DH &  dof,
const unsigned char  boundary_component,
const Function< DH::space_dimension > &  boundary_function,
ConstraintMatrix constraints,
const std::vector< bool > &  component_mask = std::vector< bool >() 
) [inline, static]

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

template<class DH >
static void VectorTools::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 >() 
) [inline, static]

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

template<int dim, int spacedim>
static 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 >() 
) [inline, static]

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>
static 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 >() 
) [inline, static]

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

template<int dim, int spacedim>
static 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,
ConstraintMatrix constraints,
std::vector< unsigned int component_mapping = std::vector< unsigned int >() 
) [inline, static]

Project a function to the boundary of the domain, using the given quadrature formula for the faces. If the ConstraintMatrix constraints contained values or other constraints 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 list of constraints. This is handled by using inhomogeneous constraints. Please note that when combining adaptive meshes and this kind of constraints, the Dirichlet conditions should be set first, and then completed by hanging node constraints, in order to make sure that the discretization remains consistent.

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.

  • 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>
static void VectorTools::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 >() 
) [inline, static]

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

template<int dim, template< int, int > class DH, int spacedim>
static void VectorTools::compute_no_normal_flux_constraints ( const DH< dim, spacedim > &  dof_handler,
const unsigned int  first_vector_component,
const std::set< unsigned char > &  boundary_ids,
ConstraintMatrix constraints,
const Mapping< dim, spacedim > &  mapping = StaticMappingQ1< dim >::mapping 
) [inline, static]

Compute the constraints that correspond to boundary conditions of the form $\vec n \cdot \vec u=0$, i.e. no normal flux if $\vec u$ is a vector-valued quantity. These conditions have exactly the form handled by the ConstraintMatrix class, so instead of creating a map between boundary degrees of freedom and corresponding value, we here create a list of constraints that are written into a ConstraintMatrix. This object may already have some content, for example from hanging node constraints, that remains untouched. These constraints have to be applied to the linear system like any other such constraints, i.e. you have to condense the linear system with the constraints before solving, and you have to distribute the solution vector afterwards.

The use of this function is explained in more detail in step-31. It doesn't make much sense in 1d, so the function throws an exception in that case.

The second argument of this function denotes the first vector component in the finite element that corresponds to the vector function that you want to constrain. For example, if we were solving a Stokes equation in 2d and the finite element had components $(u,v,p)$, then first_vector_component would be zero. On the other hand, if we solved the Maxwell equations in 3d and the finite element has components $(E_x,E_y,E_z,B_x,B_y,B_z)$ and we want the boundary condition $\vec n\cdot \vec B=0$, then first_vector_component would be 3. Vectors are implicitly assumed to have exactly dim components that are ordered in the same way as we usually order the coordinate directions, i.e. $x$-, $y$-, and finally $z$-component. The function assumes, but can't check, that the vector components in the range [first_vector_component,first_vector_component+dim) come from the same base finite element. For example, in the Stokes example above, it would not make sense to use a FESystem<dim>(FE_Q<dim>(2), 1, FE_Q<dim>(1), dim) (note that the first velocity vector component is a $Q_2$ element, whereas all the other ones are $Q_1$ elements) as there would be points on the boundary where the $x$-velocity is defined but no corresponding $y$- or $z$-velocities.

The third argument denotes the set of boundary indicators on which the boundary condition is to be enforced. Note that, as explained below, this is one of the few functions where it makes a difference where we call the function multiple times with only one boundary indicator, or whether we call the function onces with the whole set of boundary indicators at once.

The last argument is denoted to compute the normal vector $\vec n$ at the boundary points.

Computing constraints in 2d

Computing these constraints requires some smarts. The main question revolves around the question what the normal vector is. Consider the following situation:

no_normal_flux_1.png

Here, we have two cells that use a bilinear mapping (i.e. MappingQ1). Consequently, for each of the cells, the normal vector is perpendicular to the straight edge. If the two edges at the top and right are meant to approximate a curved boundary (as indicated by the dashed line), then neither of the two computed normal vectors are equal to the exact normal vector (though they approximate it as the mesh is refined further). What is worse, if we constrain $\vec n \cdot \vec u=0$ at the common vertex with the normal vector from both cells, then we constrain the vector $\vec u$ with respect to two linearly independent vectors; consequently, the constraint would be $\vec u=0$ at this point (i.e. all components of the vector), which is not what we wanted.

To deal with this situation, the algorithm works in the following way: at each point where we want to constrain $\vec u$, we first collect all normal vectors that adjacent cells might compute at this point. We then do not constrain $\vec n \cdot \vec u=0$ for each of these normal vectors but only for the average of the normal vectors. In the example above, we therefore record only a single constraint $\vec n \cdot \vec {\bar u}=0$, where $\vec {\bar u}$ is the average of the two indicated normal vectors.

Unfortunately, this is not quite enough. Consider the situation here:

no_normal_flux_2.png

If again the top and right edges approximate a curved boundary, and the left boundary a separate boundary (for example straight) so that the exact boundary has indeed a corner at the top left vertex, then the above construction would not work: here, we indeed want the constraint that $\vec u$ at this point (because the normal velocities with respect to both the left normal as well as the top normal vector should be zero), not that the velocity in the direction of the average normal vector is zero.

Consequently, we use the following heuristic to determine whether all normal vectors computed at one point are to be averaged: if two normal vectors for the same point are computed on different cells, then they are to be averaged. This covers the first example above. If they are computed from the same cell, then the fact that they are different is considered indication that they come from different parts of the boundary that might be joined by a real corner, and must not be averaged.

There is one problem with this scheme. If, for example, the same domain we have considered above, is discretized with the following mesh, then we get into trouble:

no_normal_flux_2.png

Here, the algorithm assumes that the boundary does not have a corner at the point where faces $F1$ and $F2$ join because at that point there are two different normal vectors computed from different cells. If you intend for there to be a corner of the exact boundary at this point, the only way to deal with this is to assign the two parts of the boundary different boundary indicators and call this function twice, once for each boundary indicators; doing so will yield only one normal vector at this point per invocation (because we consider only one boundary part at a time), with the result that the normal vectors will not be averaged.

Computing constraints in 3d

The situation is more complicated in 3d. Consider the following case where we want to compute the constraints at the marked vertex:

no_normal_flux_4.png

Here, we get four different normal vectors, one from each of the four faces that meet at the vertex. Even though they may form a complete set of vectors, it is not our intent to constrain all components of the vector field at this point. Rather, we would like to still allow tangential flow, where the term "tangential" has to be suitably defined.

In a case like this, the algorithm proceeds as follows: for each cell that has computed two tangential vectors at this point, we compute the unconstrained direction as the outer product of the two tangential vectors (if necessary multiplied by minus one). We then average these tangential vectors. Finally, we compute constraints for the two directions perpendicular to this averaged tangential direction.

There are cases where one cell contributes two tangential directions and another one only one; for example, this would happen if both top and front faces of the left cell belong to the boundary selected whereas only the top face of the right cell belongs to it. This case is not currently implemented.

Results

Because it makes for good pictures, here are two images of vector fields on a circle and on a sphere to which the constraints computed by this function have been applied:

no_normal_flux_5.png
no_normal_flux_6.png

The vectors fields are not physically reasonable but the tangentiality constraint is clearly enforced. The fact that the vector fields are zero at some points on the boundary is an artifact of the way it is created, it is not constrained to be zero at these points.

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

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>
static void VectorTools::create_right_hand_side ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  q,
const Function< spacedim > &  rhs,
Vector< double > &  rhs_vector 
) [inline, static]

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

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

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

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

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

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

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>
static void VectorTools::create_point_source_vector ( const DoFHandler< dim, spacedim > &  dof,
const Point< spacedim > &  p,
Vector< double > &  rhs_vector 
) [inline, static]

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

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

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

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

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>
static 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< unsigned char > &  boundary_indicators = std::set< unsigned char >() 
) [inline, static]

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>
static 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< unsigned char > &  boundary_indicators = std::set< unsigned char >() 
) [inline, static]

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

template<int dim, int spacedim>
static 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< unsigned char > &  boundary_indicators = std::set< unsigned char >() 
) [inline, static]

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

template<int dim, int spacedim>
static 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< unsigned char > &  boundary_indicators = std::set< unsigned char >() 
) [inline, static]

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>
static 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. 
) [inline, static]

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.

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>
static 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. 
) [inline, static]

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

template<int dim, class InVector , class OutVector , int spacedim>
static void VectorTools::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. 
) [inline, static]
template<int dim, class InVector , class OutVector , int spacedim>
static 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. 
) [inline, static]

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

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

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>
static 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 
) [inline, static]

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>
static void VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point,
Vector< double > &  value 
) [inline, static]

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>
static double VectorTools::point_value ( const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point 
) [inline, static]

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.

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

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>
static double VectorTools::point_value ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  dof,
const InVector &  fe_function,
const Point< spacedim > &  point 
) [inline, static]

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.

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

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 degree 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).

Apart from the vector v to operate on, this function takes a bit vector. This 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.

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

Compute the mean value of one component of the solution.

This function integrates the chosen component over the whole domain and returns the result.

Subtracting this mean value from the node 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. Therefore, it is necessary to compute the mean value and subtract it in the evaluation routine.

So far, this is needed only in the error evaluation for Stokes with complete Dirichlet boundaries for the velocities.

template<int dim, class InVector , int spacedim>
static double VectorTools::compute_mean_value ( const DoFHandler< dim, spacedim > &  dof,
const Quadrature< dim > &  quadrature,
const InVector &  v,
const unsigned int  component 
) [inline, static]

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


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

deal.II documentation generated on Mon Nov 23 22:58:17 2009 by doxygen 1.6.1