|
| FEFieldFunction (const DoFHandlerType &dh, const VectorType &data_vector, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping) |
|
void | set_active_cell (const typename DoFHandlerType::active_cell_iterator &newcell) |
|
virtual void | vector_value (const Point< dim > &p, Vector< typename VectorType::value_type > &values) const override |
|
virtual VectorType::value_type | value (const Point< dim > &p, const unsigned int component=0) const override |
|
virtual void | value_list (const std::vector< Point< dim >> &points, std::vector< typename VectorType::value_type > &values, const unsigned int component=0) const override |
|
virtual void | vector_value_list (const std::vector< Point< dim >> &points, std::vector< Vector< typename VectorType::value_type >> &values) const override |
|
virtual void | vector_gradient (const Point< dim > &p, std::vector< Tensor< 1, dim, typename VectorType::value_type >> &gradients) const override |
|
virtual Tensor< 1, dim, typename VectorType::value_type > | gradient (const Point< dim > &p, const unsigned int component=0) const override |
|
virtual void | vector_gradient_list (const std::vector< Point< dim >> &p, std::vector< std::vector< Tensor< 1, dim, typename VectorType::value_type >>> &gradients) const override |
|
virtual void | gradient_list (const std::vector< Point< dim >> &p, std::vector< Tensor< 1, dim, typename VectorType::value_type >> &gradients, const unsigned int component=0) const override |
|
virtual VectorType::value_type | laplacian (const Point< dim > &p, const unsigned int component=0) const override |
|
virtual void | vector_laplacian (const Point< dim > &p, Vector< typename VectorType::value_type > &values) const override |
|
virtual void | laplacian_list (const std::vector< Point< dim >> &points, std::vector< typename VectorType::value_type > &values, const unsigned int component=0) const override |
|
virtual void | vector_laplacian_list (const std::vector< Point< dim >> &points, std::vector< Vector< typename VectorType::value_type >> &values) const override |
|
unsigned int | compute_point_locations (const std::vector< Point< dim >> &points, std::vector< typename DoFHandlerType::active_cell_iterator > &cells, std::vector< std::vector< Point< dim >>> &qpoints, std::vector< std::vector< unsigned int >> &maps) const |
|
| Function (const unsigned int n_components=1, const time_type initial_time=0.0) |
|
| Function (const Function &f)=default |
|
virtual | ~Function () override=0 |
|
Function & | operator= (const Function &f) |
|
virtual void | vector_value (const Point< dim > &p, Vector< Vector< double > ::value_type > &values) const |
|
virtual void | value_list (const std::vector< Point< dim >> &points, std::vector< Vector< double > ::value_type > &values, const unsigned int component=0) const |
|
virtual void | vector_value_list (const std::vector< Point< dim >> &points, std::vector< Vector< Vector< double > ::value_type >> &values) const |
|
virtual void | vector_values (const std::vector< Point< dim >> &points, std::vector< std::vector< Vector< double > ::value_type >> &values) const |
|
virtual void | vector_gradient (const Point< dim > &p, std::vector< Tensor< 1, dim, Vector< double > ::value_type >> &gradients) const |
|
virtual void | gradient_list (const std::vector< Point< dim >> &points, std::vector< Tensor< 1, dim, Vector< double > ::value_type >> &gradients, const unsigned int component=0) const |
|
virtual void | vector_gradients (const std::vector< Point< dim >> &points, std::vector< std::vector< Tensor< 1, dim, Vector< double > ::value_type >>> &gradients) const |
|
virtual void | vector_gradient_list (const std::vector< Point< dim >> &points, std::vector< std::vector< Tensor< 1, dim, Vector< double > ::value_type >>> &gradients) const |
|
virtual void | vector_laplacian (const Point< dim > &p, Vector< Vector< double > ::value_type > &values) const |
|
virtual void | laplacian_list (const std::vector< Point< dim >> &points, std::vector< Vector< double > ::value_type > &values, const unsigned int component=0) const |
|
virtual void | vector_laplacian_list (const std::vector< Point< dim >> &points, std::vector< Vector< Vector< double > ::value_type >> &values) const |
|
virtual SymmetricTensor< 2, dim, Vector< double > ::value_type > | hessian (const Point< dim > &p, const unsigned int component=0) const |
|
virtual void | vector_hessian (const Point< dim > &p, std::vector< SymmetricTensor< 2, dim, Vector< double > ::value_type >> &values) const |
|
virtual void | hessian_list (const std::vector< Point< dim >> &points, std::vector< SymmetricTensor< 2, dim, Vector< double > ::value_type >> &values, const unsigned int component=0) const |
|
virtual void | vector_hessian_list (const std::vector< Point< dim >> &points, std::vector< std::vector< SymmetricTensor< 2, dim, Vector< double > ::value_type >>> &values) const |
|
std::size_t | memory_consumption () const |
|
| FunctionTime (const numbers::NumberTraits< Vector< double > ::value_type >::real_type initial_time=numbers::NumberTraits< Vector< double > ::value_type >::real_type(0.0)) |
|
virtual | ~FunctionTime ()=default |
|
numbers::NumberTraits< Vector< double > ::value_type >::real_type | get_time () const |
|
virtual void | set_time (const numbers::NumberTraits< Vector< double > ::value_type >::real_type new_time) |
|
virtual void | advance_time (const numbers::NumberTraits< Vector< double > ::value_type >::real_type delta_t) |
|
| Subscriptor () |
|
| Subscriptor (const Subscriptor &) |
|
| Subscriptor (Subscriptor &&) noexcept |
|
virtual | ~Subscriptor () |
|
Subscriptor & | operator= (const Subscriptor &) |
|
Subscriptor & | operator= (Subscriptor &&) noexcept |
|
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
|
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
|
unsigned int | n_subscriptions () const |
|
template<typename StreamType > |
void | list_subscribers (StreamType &stream) const |
|
void | list_subscribers () const |
|
template<class Archive > |
void | serialize (Archive &ar, const unsigned int version) |
|
template<int dim, typename DoFHandlerType = DoFHandler<dim>, typename VectorType = Vector<double>>
class Functions::FEFieldFunction< dim, DoFHandlerType, VectorType >
This is an interpolation function for the given dof handler and the given solution vector. The points at which this function can be evaluated MUST be inside the domain of the dof handler, but except from this, no other requirement is given. This function is rather slow, as it needs to construct a quadrature object for the point (or set of points) where you want to evaluate your finite element function. In order to do so, it needs to find out where the points lie.
If you know in advance in which cell your points lie, you can accelerate things a bit, by calling set_active_cell before asking for values or gradients of the function. If you don't do this, and your points don't lie in the cell that is currently stored, the function GridTools::find_cell_around_point is called to find out where the point is. You can specify an optional mapping to use when looking for points in the grid. If you don't do so, this function uses a Q1 mapping.
Once the FEFieldFunction knows where the points lie, it creates a quadrature formula for those points, and calls FEValues::get_function_values or FEValues::get_function_gradients with the given quadrature points.
If you only need the quadrature points but not the values of the finite element function (you might want this for the adjoint interpolation), you can also use the function compute_point_locations
alone.
An example of how to use this function is the following:
FEFieldFunction<dim> fe_function_1 (dh_1, solution_1);
fe_function_1, solution_2);
The snippet of code above will work assuming that the second triangulation is entirely included in the first one.
FEFieldFunction is designed to be an easy way to get the results of your computations across different, possibly non matching, grids. No knowledge of the location of the points is assumed in this class, which makes it rely entirely on the GridTools::find_active_cell_around_point utility for its job. However the class can be fed an "educated guess" of where the points that will be computed actually are by using the FEFieldFunction::set_active_cell method, so if you have a smart way to tell where your points are, you will save a lot of computational time by letting this class know.
When using this class with a parallel distributed triangulation object and evaluating the solution at a particular point, not every processor will own the cell at which the solution is evaluated. Rather, it may be that the cell in which this point is found is in fact a ghost or artificial cell (see GlossArtificialCell and GlossGhostCell). The solution can be evaluated on ghost cells, but for artificial cells we have no access to the solution there and functions that evaluate the solution at such a point will trigger an exception of type VectorTools::ExcPointNotAvailableHere.
To deal with this situation, you will want to use code as follows when, for example, evaluating the solution at the origin (here using a parallel TrilinosWrappers vector to hold the solution):
solution_function (dof_handler, solution);
double solution_at_origin;
bool point_found = true;
try
{
solution_at_origin = solution_function.value (origin);
}
{
point_found = false;
}
if (point_found == true)
...do something...;
- Author
- Luca Heltai, 2006, Markus Buerg, 2012, Wolfgang Bangerth, 2013
Definition at line 166 of file fe_field_function.h.