Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/function.h>
Public Member Functions | |
VectorFunctionFromTensorFunction (const TensorFunction< 1, dim, RangeNumberType > &tensor_function, const unsigned int selected_component=0, const unsigned int n_components=dim) | |
virtual | ~VectorFunctionFromTensorFunction ()=default |
virtual RangeNumberType | value (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_value (const Point< dim > &p, Vector< RangeNumberType > &values) const |
virtual void | vector_value_list (const std::vector< Point< dim > > &points, std::vector< Vector< RangeNumberType > > &value_list) const |
Public Member Functions inherited from Function< dim, RangeNumberType > | |
Function (const unsigned int n_components=1, const RangeNumberType initial_time=0.0) | |
virtual | ~Function ()=0 |
Function & | operator= (const Function &f) |
virtual void | value_list (const std::vector< Point< dim > > &points, std::vector< RangeNumberType > &values, const unsigned int component=0) const |
virtual void | vector_values (const std::vector< Point< dim > > &points, std::vector< std::vector< RangeNumberType > > &values) const |
virtual Tensor< 1, dim, RangeNumberType > | gradient (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_gradient (const Point< dim > &p, std::vector< Tensor< 1, dim, RangeNumberType > > &gradients) const |
virtual void | gradient_list (const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim, RangeNumberType > > &gradients, const unsigned int component=0) const |
virtual void | vector_gradients (const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim, RangeNumberType > > > &gradients) const |
virtual void | vector_gradient_list (const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim, RangeNumberType > > > &gradients) const |
virtual RangeNumberType | laplacian (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_laplacian (const Point< dim > &p, Vector< RangeNumberType > &values) const |
virtual void | laplacian_list (const std::vector< Point< dim > > &points, std::vector< RangeNumberType > &values, const unsigned int component=0) const |
virtual void | vector_laplacian_list (const std::vector< Point< dim > > &points, std::vector< Vector< RangeNumberType > > &values) const |
virtual SymmetricTensor< 2, dim, RangeNumberType > | hessian (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_hessian (const Point< dim > &p, std::vector< SymmetricTensor< 2, dim, RangeNumberType > > &values) const |
virtual void | hessian_list (const std::vector< Point< dim > > &points, std::vector< SymmetricTensor< 2, dim, RangeNumberType > > &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, RangeNumberType > > > &values) const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from FunctionTime< RangeNumberType > | |
FunctionTime (const RangeNumberType initial_time=RangeNumberType(0.0)) | |
virtual | ~FunctionTime ()=default |
RangeNumberType | get_time () const |
virtual void | set_time (const RangeNumberType new_time) |
virtual void | advance_time (const RangeNumberType 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 (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Attributes | |
const TensorFunction< 1, dim, RangeNumberType > & | tensor_function |
const unsigned int | selected_component |
Additional Inherited Members | |
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, RangeNumberType > | |
const unsigned int | n_components |
Static Public Attributes inherited from Function< dim, RangeNumberType > | |
static const unsigned int | dimension = dim |
This class is built as a means of translating the Tensor<1,dim, RangeNumberType>
values produced by objects of type TensorFunction and returning them as a multiple component version of the same thing as a Vector for use in, for example, the VectorTools::interpolate or the many other functions taking Function objects. It allows the user to place the desired components into an n_components
long vector starting at the selected_component
location in that vector and have all other components be 0.
For example: Say you created a class called
which extends the TensorFunction class and you have an object
of that class which you want to interpolate onto your mesh using the VectorTools::interpolate function, but the finite element you use for the DoFHandler object has 3 copies of a finite element with dim
components, for a total of 3*dim components. To interpolate onto that DoFHandler, you need an object of type Function that has 3*dim vector components. Creating such an object from the existing rhs
object is done using this piece of code:
where the dim
components of the tensor function are placed into the first dim
components of the function object.
Definition at line 794 of file function.h.
VectorFunctionFromTensorFunction< dim, RangeNumberType >::VectorFunctionFromTensorFunction | ( | const TensorFunction< 1, dim, RangeNumberType > & | tensor_function, |
const unsigned int | selected_component = 0 , |
||
const unsigned int | n_components = dim |
||
) |
Given a TensorFunction object that takes a Point
and returns a Tensor<1,dim, RangeNumberType>
value, convert this into an object that matches the Function<dim> interface.
By default, create a Vector object of the same size as tensor_function
returns, i.e., with dim
components.
tensor_function | The TensorFunction that will form one component of the resulting Vector Function object. |
n_components | The total number of vector components of the resulting TensorFunction object. |
selected_component | The first component that should be filled by the first argument. This should be such that the entire tensor_function fits inside the n_component length return vector. |
|
virtualdefault |
This destructor is defined as virtual so as to coincide with all other aspects of class.
|
virtual |
Return a single component of a vector-valued function at a given point.
Reimplemented from Function< dim, RangeNumberType >.
|
virtual |
Return all components of a vector-valued function at a given point.
values
shall have the right size beforehand, i.e. n_components.
Reimplemented from Function< dim, RangeNumberType >.
|
virtual |
Return all components of a vector-valued function at a list of points.
value_list
shall be the same size as points
and each element of the vector will be passed to vector_value() to evaluate the function
Reimplemented from Function< dim, RangeNumberType >.
|
private |
The TensorFunction object which we call when this class's vector_value() or vector_value_list() functions are called.
Definition at line 852 of file function.h.
|
private |
The first vector component whose value is to be filled by the given TensorFunction. The values will be placed in components selected_component to selected_component+dim-1 for a TensorFunction<1,dim, RangeNumberType>
object.
Definition at line 860 of file function.h.