Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/function_lib.h>
Public Member Functions | |
InterpolatedTensorProductGridData (const std::array< std::vector< double >, dim > &coordinate_values, const Table< dim, double > &data_values) | |
virtual double | value (const Point< dim > &p, const unsigned int component=0) const override |
virtual Tensor< 1, dim > | gradient (const Point< dim > &p, const unsigned int component=0) const override |
Public Member Functions inherited from Function< dim > | |
Function (const unsigned int n_components=1, const time_type initial_time=0.0) | |
virtual | ~Function () override=0 |
Function & | operator= (const Function &f) |
virtual void | vector_value (const Point< dim > &p, Vector< double > &values) const |
virtual void | value_list (const std::vector< Point< dim >> &points, std::vector< double > &values, const unsigned int component=0) const |
virtual void | vector_value_list (const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) const |
virtual void | vector_values (const std::vector< Point< dim >> &points, std::vector< std::vector< double >> &values) const |
virtual void | vector_gradient (const Point< dim > &p, std::vector< Tensor< 1, dim, double >> &gradients) const |
virtual void | gradient_list (const std::vector< Point< dim >> &points, std::vector< Tensor< 1, dim, double >> &gradients, const unsigned int component=0) const |
virtual void | vector_gradients (const std::vector< Point< dim >> &points, std::vector< std::vector< Tensor< 1, dim, double >>> &gradients) const |
virtual void | vector_gradient_list (const std::vector< Point< dim >> &points, std::vector< std::vector< Tensor< 1, dim, double >>> &gradients) const |
virtual double | laplacian (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_laplacian (const Point< dim > &p, Vector< double > &values) const |
virtual void | laplacian_list (const std::vector< Point< dim >> &points, std::vector< double > &values, const unsigned int component=0) const |
virtual void | vector_laplacian_list (const std::vector< Point< dim >> &points, std::vector< Vector< double >> &values) const |
virtual SymmetricTensor< 2, dim, double > | hessian (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_hessian (const Point< dim > &p, std::vector< SymmetricTensor< 2, dim, double >> &values) const |
virtual void | hessian_list (const std::vector< Point< dim >> &points, std::vector< SymmetricTensor< 2, dim, double >> &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, double >>> &values) const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from FunctionTime< Number > | |
FunctionTime (const Number initial_time=Number(0.0)) | |
virtual | ~FunctionTime ()=default |
Number | get_time () const |
virtual void | set_time (const Number new_time) |
virtual void | advance_time (const Number delta_t) |
Public Member Functions inherited from Subscriptor | |
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) |
Protected Member Functions | |
TableIndices< dim > | table_index_of_point (const Point< dim > &p) const |
Protected Attributes | |
const std::array< std::vector< double >, dim > | coordinate_values |
const Table< dim, double > | data_values |
Additional Inherited Members | |
Public Types inherited from Function< dim > | |
using | time_type = typename FunctionTime< typename numbers::NumberTraits< double >::real_type >::time_type |
Public Types inherited from FunctionTime< Number > | |
using | time_type = Number |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Public Attributes inherited from Function< dim > | |
const unsigned int | n_components |
Static Public Attributes inherited from Function< dim > | |
static const unsigned int | dimension |
A scalar function that computes its values by (bi-, tri-)linear interpolation from a set of point data that are arranged on a possibly non-uniform tensor product mesh. In other words, considering the three- dimensional case, let there be points \(x_0,\ldots, x_{K-1}\), \(y_0,\ldots,y_{L-1}\), \(z_1,\ldots,z_{M-1}\), and data \(d_{klm}\) defined at point \((x_k,y_l,z_m)^T\), then evaluating the function at a point \(\mathbf x=(x,y,z)\) will find the box so that \(x_k\le x\le x_{k+1}, y_l\le y\le y_{l+1}, z_m\le z\le z_{m+1}\), and do a trilinear interpolation of the data on this cell. Similar operations are done in lower dimensions.
This class is most often used for either evaluating coefficients or right hand sides that are provided experimentally at a number of points inside the domain, or for comparing outputs of a solution on a finite element mesh against previously obtained data defined on a grid.
If a point is requested outside the box defined by the end points of the coordinate arrays, then the function is assumed to simply extend by constant values beyond the last data point in each coordinate direction. (The class does not throw an error if a point lies outside the box since it frequently happens that a point lies just outside the box by an amount on the order of numerical roundoff.)
Definition at line 1426 of file function_lib.h.
Functions::InterpolatedTensorProductGridData< dim >::InterpolatedTensorProductGridData | ( | const std::array< std::vector< double >, dim > & | coordinate_values, |
const Table< dim, double > & | data_values | ||
) |
Constructor.
coordinate_values | An array of dim arrays. Each of the inner arrays contains the coordinate values \(x_0,\ldots, x_{K-1}\) and similarly for the other coordinate directions. These arrays need not have the same size. Obviously, we need dim such arrays for a dim- dimensional function object. The coordinate values within this array are assumed to be strictly ascending to allow for efficient lookup. |
data_values | A dim-dimensional table of data at each of the mesh points defined by the coordinate arrays above. Note that the Table class has a number of conversion constructors that allow converting other data types into a table where you specify this argument. |
Definition at line 2479 of file function_lib.cc.
|
overridevirtual |
Compute the value of the function set by bilinear interpolation of the given data set.
p | The point at which the function is to be evaluated. |
component | The vector component. Since this function is scalar, only zero is a valid argument here. |
Reimplemented from Function< dim >.
Definition at line 2544 of file function_lib.cc.
|
overridevirtual |
Compute the gradient of the function defined by bilinear interpolation of the given data set.
p | The point at which the function gradient is to be evaluated. |
component | The vector component. Since this function is scalar, only zero is a valid argument here. |
Reimplemented from Function< dim >.
Definition at line 2575 of file function_lib.cc.
|
protected |
Find the index in the table of the rectangle containing an input point
Definition at line 2507 of file function_lib.cc.
|
protected |
The set of coordinate values in each of the coordinate directions.
Definition at line 1484 of file function_lib.h.
|
protected |
The data that is to be interpolated.
Definition at line 1489 of file function_lib.h.