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 | |
| ::ExceptionBase & | ExcInvalidBoundaryIndicator () |
| ::ExceptionBase & | ExcNonInterpolatingFE () |
| ::ExceptionBase & | ExcNoComponentSelected () |
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) |
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.
This collection of methods offers the following operations:
Interpolation: assign each degree of freedom in the vector to be the value of the function given as argument. This is identical to saying that the resulting finite element function (which is isomorphic to the output vector) has exact function values in all support points of trial functions. The support point of a trial function is the point where its value equals one, e.g. for linear trial functions the support points are four corners of an element. This function therefore relies on the assumption that a finite element is used for which the degrees of freedom are function values (Lagrange elements) rather than gradients, normal derivatives, second derivatives, etc (Hermite elements, quintic Argyris element, etc.).
It seems inevitable that some values of the vector to be created are set twice or even more than that. The reason is that we have to loop over all cells and get the function values for each of the trial functions located thereon. This applies also to the functions located on faces and corners which we thus visit more than once. While setting the value in the vector is not an expensive operation, the evaluation of the given function may be, taking into account that a virtual function has to be called.
Projection: compute the L2-projection of the given function onto the finite element space, i.e. if f is the function to be projected, compute fh in Vh such that (fh,vh)=(f,vh) for all discrete test functions vh. This is done through the solution of the linear system of equations M v = f where M is the mass matrix
and
. The solution vector
then is the nodal representation of the projection fh. The project() functions are used in the step-21 and step-23 tutorial programs.
In order to get proper results, it be may necessary to treat boundary conditions right. Below are listed some cases where this may be needed. If needed, this is done by L2-projection of the trace of the given function onto the finite element space restricted to the boundary of the domain, then taking this information and using it to eliminate the boundary nodes from the mass matrix of the whole domain, using the MatrixTools::apply_boundary_values() function. The projection of the trace of the function to the boundary is done with the VectorTools::project_boundary_values() (see below) function, which is called with a map of boundary functions FunctioMap in which all boundary indicators from zero to types::internal_face_boundary_id-1 (types::internal_face_boundary_id is used for other purposes, see the Triangulation class documentation) point to the function to be projected. The projection to the boundary takes place using a second quadrature formula on the boundary given to the project() function. The first quadrature formula is used to compute the right hand side and for numerical quadrature of the mass matrix.
The projection of the boundary values first, then eliminating them from the global system of equations is not needed usually. It may be necessary if you want to enforce special restrictions on the boundary values of the projected function, for example in time dependent problems: you may want to project the initial values but need consistency with the boundary values for later times. Since the latter are projected onto the boundary in each time step, it is necessary that we also project the boundary values of the initial values, before projecting them to the whole domain.
Obviously, the results of the two schemes for projection are different. Usually, when projecting to the boundary first, the L2-norm of the difference between original function and projection over the whole domain will be larger (factors of five have been observed) while the L2-norm of the error integrated over the boundary should of course be less. The reverse should also hold if no projection to the boundary is performed.
The selection whether the projection to the boundary first is needed is done with the project_to_boundary_first flag passed to the function. If false is given, the additional quadrature formula for faces is ignored.
You should be aware of the fact that if no projection to the boundary is requested, a function with zero boundary values may not have zero boundary values after projection. There is a flag for this especially important case, which tells the function to enforce zero boundary values on the respective boundary parts. Since enforced zero boundary values could also have been reached through projection, but are more economically obtain using other methods, the project_to_boundary_first flag is ignored if the enforce_zero_boundary flag is set.
The solution of the linear system is presently done using a simple CG method without preconditioning and without multigrid. This is clearly not too efficient, but sufficient in many cases and simple to implement. This detail may change in the future.
Creation of right hand side vectors: The create_right_hand_side() function computes the vector
. This is the same as what the MatrixCreator::create_* functions which take a right hand side do, but without assembling a matrix.
Creation of right hand side vectors for point sources: The create_point_source_vector() function computes the vector
.
Creation of boundary right hand side vectors: The create_boundary_right_hand_side() function computes the vector
. This is the right hand side contribution of boundary forces when having inhomogeneous Neumann boundary values in Laplace's equation or other second order operators. This function also takes an optional argument denoting over which parts of the boundary the integration shall extend. If the default argument is used, it is applied to all boundaries.
Interpolation of boundary values: The MatrixTools::apply_boundary_values() function takes a list of boundary nodes and their values. You can get such a list by interpolation of a boundary function using the interpolate_boundary_values() function. To use it, you have to specify a list of pairs of boundary indicators (of type types::boundary_id_t; see the section in the documentation of the Triangulation class for more details) and the according functions denoting the dirichlet boundary values of the nodes on boundary faces with this boundary indicator.
Usually, all other boundary conditions, such as inhomogeneous Neumann values or mixed boundary conditions are handled in the weak formulation. No attempt is made to include these into the process of matrix and vector assembly therefore.
Within this function, boundary values are interpolated, i.e. a node is given the point value of the boundary function. In some cases, it may be necessary to use the L2-projection of the boundary function or any other method. For this purpose we refer to the project_boundary_values() function below.
You should be aware that the boundary function may be evaluated at nodes on the interior of faces. These, however, need not be on the true boundary, but rather are on the approximation of the boundary represented by the mapping of the unit cell to the real cell. Since this mapping will in most cases not be the exact one at the face, the boundary function is evaluated at points which are not on the boundary and you should make sure that the returned values are reasonable in some sense anyway.
In 1d the situation is a bit different since there faces (i.e. vertices) have no boundary indicator. It is assumed that if the boundary indicator zero is given in the list of boundary functions, the left boundary point is to be interpolated while the right boundary point is associated with the boundary index 1 in the map. The respective boundary functions are then evaluated at the place of the respective boundary point.
Projection of boundary values: The project_boundary_values() function acts similar to the interpolate_boundary_values() function, apart from the fact that it does not get the nodal values of boundary nodes by interpolation but rather through the L2-projection of the trace of the function to the boundary.
The projection takes place on all boundary parts with boundary indicators listed in the map (FunctioMap::FunctionMap) of boundary functions. These boundary parts may or may not be continuous. For these boundary parts, the mass matrix is assembled using the MatrixTools::create_boundary_mass_matrix() function, as well as the appropriate right hand side. Then the resulting system of equations is solved using a simple CG method (without preconditioning), which is in most cases sufficient for the present purpose.
Computing errors: The function integrate_difference() performs the calculation of the error between a given (continuous) reference function and the finite element solution in different norms. The integration is performed using a given quadrature formula and assumes that the given finite element objects equals that used for the computation of the solution.
The result is stored in a vector (named difference), where each entry equals the given norm of the difference on a cell. The order of entries is the same as a cell_iterator takes when started with begin_active and promoted with the ++ operator.
This data, one number per active cell, can be used to generate graphical output by directly passing it to the DataOut class through the DataOut::add_data_vector function. Alternatively, it can be interpolated to the nodal points of a finite element field using the DoFTools::distribute_cell_to_dof_vector function.
Presently, there is the possibility to compute the following values from the difference, on each cell: mean, L1_norm, L2_norm, Linfty_norm, H1_seminorm and H1_norm, see VectorTools::NormType. For the mean difference value, the reference function minus the numerical solution is computed, not the other way round.
The infinity norm of the difference on a given cell returns the maximum absolute value of the difference at the quadrature points given by the quadrature formula parameter. This will in some cases not be too good an approximation, since for example the Gauss quadrature formulae do not evaluate the difference at the end or corner points of the cells. You may want to choose a quadrature formula with more quadrature points or one with another distribution of the quadrature points in this case. You should also take into account the superconvergence properties of finite elements in some points: for example in 1D, the standard finite element method is a collocation method and should return the exact value at nodal points. Therefore, the trapezoidal rule should always return a vanishing L-infinity error. Conversely, in 2D the maximum L-infinity error should be located at the vertices or at the center of the cell, which would make it plausible to use the Simpson quadrature rule. On the other hand, there may be superconvergence at Gauss integration points. These examples are not intended as a rule of thumb, rather they are thought to illustrate that the use of the wrong quadrature formula may show a significantly wrong result and care should be taken to chose the right formula.
The H1 seminorm is the L2 norm of the gradient of the difference. The square of the full H1 norm is the sum of the square of seminorm and the square of the L2 norm.
To get the global L1 error, you have to sum up the entries in difference, e.g. using Vector::l1_norm() function. For the global L2 difference, you have to sum up the squares of the entries and take the root of the sum, e.g. using Vector::l2_norm(). These two operations represent the l1 and l2 norms of the vectors, but you need not take the absolute value of each entry, since the cellwise norms are already positive.
To get the global mean difference, simply sum up the elements as above. To get the
norm, take the maximum of the vector elements, e.g. using the Vector::linfty_norm() function.
For the global H1 norm and seminorm, the same rule applies as for the L2 norm: compute the l2 norm of the cell error vector.
Note that, in the codimension one case, if you ask for a norm that requires the computation of a gradient, then the provided function is automatically projected along the curve, and the difference is only computed on the tangential part of the gradient, since no information is available on the normal component of the gradient anyway.
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.
Vector<float>, Vector<double>, BlockVector<float>, BlockVector<double>; others can be generated in application code (see the section on Template instantiations in the manual).Denote which norm/integral is to be computed by the integrate_difference() function of this class. The following possibilities are implemented:
| 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. |
| 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.
mapping argument should be replaced by a hp::MappingCollection in case of a hp::DoFHandler. | void VectorTools::interpolate | ( | const DH & | dof, |
| const Function< DH::space_dimension > & | function, | ||
| VECTOR & | vec | ||
| ) |
Calls the interpolate() function above with mapping=MappingQ1<dim>@().
| 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.
| 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.
| 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.
| 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.
| 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.
| 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>().
| 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.
| 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 >() |
||
| ) |
Same function as above, but taking only one pair of boundary indicator and corresponding boundary function. Calls the other function with remapped arguments.
| 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>().
| 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>().
| 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. | 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>().
| 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.
| 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>().
| 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.
| 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.
| 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.
| 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>().
| 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.
| 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.
| 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.
| 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>().
| 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.
| 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.
| 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.
| 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>().
| 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.
| 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.
| 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
-norms and
-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.
| 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>().
| 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>().
| 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.
| 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.
| 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.
| 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.
| 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.
| double VectorTools::point_value | ( | const hp::DoFHandler< dim, spacedim > & | dof, |
| const InVector & | fe_function, | ||
| const Point< spacedim > & | point | ||
| ) |
Same as above for hp.
| 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.
| 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.
| 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.
| 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
, rather than as
. 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.
. 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
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
. For such elements, a different procedure has to be used when subtracting the mean value. | 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
where
is the vector component and
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.
| 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
documentation generated on Wed May 23 2012 06:08:32 by
doxygen
1.7.3