Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Private Attributes | List of all members
FunctionFromFunctionObjects< dim, RangeNumberType > Class Template Reference

#include <deal.II/base/function.h>

Inheritance diagram for FunctionFromFunctionObjects< dim, RangeNumberType >:
[legend]

Public Member Functions

 FunctionFromFunctionObjects (const unsigned int n_components=1, const double initial_time=0)
 
 FunctionFromFunctionObjects (const std::vector< std::function< RangeNumberType(const Point< dim > &)>> &values, const double initial_time=0.0)
 
 FunctionFromFunctionObjects (const std::vector< std::function< RangeNumberType(const Point< dim > &)>> &values, const std::vector< std::function< Tensor< 1, dim, RangeNumberType >(const Point< dim > &)>> &gradients, const double initial_time=0.0)
 
virtual RangeNumberType value (const Point< dim > &p, const unsigned int component=0) const override
 
virtual Tensor< 1, dim, RangeNumberType > gradient (const Point< dim > &p, const unsigned int component=0) const override
 
void set_function_values (const std::vector< std::function< RangeNumberType(const Point< dim > &)>> &values)
 
void set_function_gradients (const std::vector< std::function< Tensor< 1, dim, RangeNumberType >(const Point< dim > &)>> &gradients)
 
- Public Member Functions inherited from Function< dim, RangeNumberType >
 Function (const unsigned int n_components=1, const time_type initial_time=0.0)
 
virtual ~Function () override=0
 
Functionoperator= (const Function &f)
 
virtual void vector_value (const Point< dim > &p, Vector< RangeNumberType > &values) const
 
virtual void value_list (const std::vector< Point< dim >> &points, std::vector< RangeNumberType > &values, const unsigned int component=0) const
 
virtual void vector_value_list (const std::vector< Point< dim >> &points, std::vector< Vector< RangeNumberType >> &values) const
 
virtual void vector_values (const std::vector< Point< dim >> &points, std::vector< std::vector< RangeNumberType >> &values) 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< numbers::NumberTraits< RangeNumberType >::real_type >
 FunctionTime (const numbers::NumberTraits< RangeNumberType >::real_type initial_time=numbers::NumberTraits< RangeNumberType >::real_type(0.0))
 
virtual ~FunctionTime ()=default
 
numbers::NumberTraits< RangeNumberType >::real_type get_time () const
 
virtual void set_time (const numbers::NumberTraits< RangeNumberType >::real_type new_time)
 
virtual void advance_time (const numbers::NumberTraits< RangeNumberType >::real_type delta_t)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (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)
 

Private Attributes

std::vector< std::function< RangeNumberType(const Point< dim > &)> > function_values
 
std::vector< std::function< Tensor< 1, dim, RangeNumberType >const Point< dim > &)> > function_gradients
 

Additional Inherited Members

- Public Types inherited from Function< dim, RangeNumberType >
using time_type = typename FunctionTime< typename numbers::NumberTraits< RangeNumberType >::real_type >::time_type
 
- Public Types inherited from FunctionTime< numbers::NumberTraits< RangeNumberType >::real_type >
using time_type = numbers::NumberTraits< RangeNumberType >::real_type
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (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
 

Detailed Description

template<int dim, typename RangeNumberType = double>
class FunctionFromFunctionObjects< dim, RangeNumberType >

This class is similar to the ScalarFunctionFromFunctionObject and VectorFunctionFromFunctionObject classes in that it allows for the easy conversion of a vector of function objects to something that satisfies the interface of the Function base class.

The difference is that here the Function object generated may be vector valued, and you can specify the gradients of the function. The number of vector components is deduced from the size of the vector in the constructor.

To be more concrete, let us consider the following example:

RangeNumberType
first_component(const Point<2> &p)
{
return 1.0;
}
RangeNumberType
second_component(const Point<2> &p)
{
return 2.0;
}
zero_gradient(const Point<2> &) {
}
custom_function({&first_component, &second_component},
{&zero_gradient, &zero_gradient});
Author
Luca Heltai, 2019

Definition at line 869 of file function.h.

Constructor & Destructor Documentation

◆ FunctionFromFunctionObjects() [1/3]

template<int dim, typename RangeNumberType = double>
FunctionFromFunctionObjects< dim, RangeNumberType >::FunctionFromFunctionObjects ( const unsigned int  n_components = 1,
const double  initial_time = 0 
)
explicit

Default constructor.

This constructor does not initialize the internal methods. To have a usable function, you need to call at least the set_function_values() method. If you need also the gradients of the solution, then you must also call the set_function_gradients() method.

◆ FunctionFromFunctionObjects() [2/3]

template<int dim, typename RangeNumberType = double>
FunctionFromFunctionObjects< dim, RangeNumberType >::FunctionFromFunctionObjects ( const std::vector< std::function< RangeNumberType(const Point< dim > &)>> &  values,
const double  initial_time = 0.0 
)

Constructor for functions of which you only know the values.

The resulting function will have a number of components equal to the size of the vector values. A call to the FunctionFromFunctionObject::gradient() method will trigger an exception, unless you first call the set_function_gradients() method.

◆ FunctionFromFunctionObjects() [3/3]

template<int dim, typename RangeNumberType = double>
FunctionFromFunctionObjects< dim, RangeNumberType >::FunctionFromFunctionObjects ( const std::vector< std::function< RangeNumberType(const Point< dim > &)>> &  values,
const std::vector< std::function< Tensor< 1, dim, RangeNumberType >(const Point< dim > &)>> &  gradients,
const double  initial_time = 0.0 
)

Constructor for functions of which you know both the values and the gradients.

The resulting function will have a number of components equal to the size of the vector values. If the size of values and gradients does not match, an exception is triggered.

Member Function Documentation

◆ value()

template<int dim, typename RangeNumberType = double>
virtual RangeNumberType FunctionFromFunctionObjects< dim, RangeNumberType >::value ( const Point< dim > &  p,
const unsigned int  component = 0 
) const
overridevirtual

Return the value of the function at the given point. Unless there is only one component (i.e. the function is scalar), you should state the component you want to have evaluated; it defaults to zero, i.e. the first component.

Reimplemented from Function< dim, RangeNumberType >.

◆ gradient()

template<int dim, typename RangeNumberType = double>
virtual Tensor<1, dim, RangeNumberType> FunctionFromFunctionObjects< dim, RangeNumberType >::gradient ( const Point< dim > &  p,
const unsigned int  component = 0 
) const
overridevirtual

Return the gradient of the function at the given point. Unless there is only one component (i.e. the function is scalar), you should state the component you want to have evaluated; it defaults to zero, i.e. the first component.

Reimplemented from Function< dim, RangeNumberType >.

◆ set_function_values()

template<int dim, typename RangeNumberType = double>
void FunctionFromFunctionObjects< dim, RangeNumberType >::set_function_values ( const std::vector< std::function< RangeNumberType(const Point< dim > &)>> &  values)

Reset the function values of this object. An assertion is thrown if the size of the values parameter does not match the number of components of this object.

◆ set_function_gradients()

template<int dim, typename RangeNumberType = double>
void FunctionFromFunctionObjects< dim, RangeNumberType >::set_function_gradients ( const std::vector< std::function< Tensor< 1, dim, RangeNumberType >(const Point< dim > &)>> &  gradients)

Reset the function gradients of this object. An assertion is thrown if the size of the gradients parameter does not match the number of components of this object.

Member Data Documentation

◆ function_values

template<int dim, typename RangeNumberType = double>
std::vector<std::function<RangeNumberType(const Point<dim> &)> > FunctionFromFunctionObjects< dim, RangeNumberType >::function_values
private

The actual function values.

Definition at line 958 of file function.h.

◆ function_gradients

template<int dim, typename RangeNumberType = double>
std::vector< std::function<Tensor<1, dim, RangeNumberType>const Point<dim> &)> > FunctionFromFunctionObjects< dim, RangeNumberType >::function_gradients
private

The actual function gradients.

Definition at line 965 of file function.h.


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