Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/function_derivative.h>
Public Member Functions | |
FunctionDerivative (const Function< dim > &f, const Point< dim > &direction, const double h=1.e-6) | |
FunctionDerivative (const Function< dim > &f, const std::vector< Point< dim > > &direction, const double h=1.e-6) | |
void | set_formula (typename AutoDerivativeFunction< dim >::DifferenceFormula formula=AutoDerivativeFunction< dim >::Euler) |
void | set_h (const double h) |
virtual double | value (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_value (const Point< dim > &p, Vector< double > &value) const |
virtual void | value_list (const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from AutoDerivativeFunction< dim > | |
AutoDerivativeFunction (const double h, const unsigned int n_components=1, const double initial_time=0.0) | |
virtual | ~AutoDerivativeFunction ()=default |
void | set_formula (const DifferenceFormula formula=Euler) |
void | set_h (const double h) |
virtual Tensor< 1, dim > | gradient (const Point< dim > &p, const unsigned int component=0) const |
virtual void | vector_gradient (const Point< dim > &p, std::vector< Tensor< 1, dim > > &gradients) const |
virtual void | gradient_list (const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const |
virtual void | vector_gradient_list (const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const |
Public Member Functions inherited from Function< dim > | |
Function (const unsigned int n_components=1, const double initial_time=0.0) | |
virtual | ~Function ()=0 |
Function & | operator= (const Function &f) |
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 (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 Function< dim > & | f |
double | h |
AutoDerivativeFunction< dim >::DifferenceFormula | formula |
std::vector< Tensor< 1, dim > > | incr |
Additional Inherited Members | |
Public Types inherited from AutoDerivativeFunction< dim > | |
enum | DifferenceFormula { Euler, UpwindEuler, FourthOrder } |
Static Public Member Functions inherited from AutoDerivativeFunction< dim > | |
static DifferenceFormula | get_formula_of_order (const unsigned int ord) |
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 |
Derivative of a function object. The value access functions of this class return the directional derivative of a function with respect to a direction provided on construction. If b
is the vector, the derivative b . grad f
is computed. This derivative is evaluated directly, not by computing the gradient of f
and its scalar product with b
.
The derivative is computed numerically, using one of the provided difference formulas (see set_formula
for available schemes). Experimenting with h
and the difference scheme may be necessary to obtain sufficient results.
Definition at line 44 of file function_derivative.h.
FunctionDerivative< dim >::FunctionDerivative | ( | const Function< dim > & | f, |
const Point< dim > & | direction, | ||
const double | h = 1.e-6 |
||
) |
Constructor. Provided are the functions to compute derivatives of, the direction vector of the differentiation and the step size h
of the difference formula.
Definition at line 25 of file function_derivative.cc.
FunctionDerivative< dim >::FunctionDerivative | ( | const Function< dim > & | f, |
const std::vector< Point< dim > > & | direction, | ||
const double | h = 1.e-6 |
||
) |
Constructor. Provided are the functions to compute derivatives of and the direction vector of the differentiation in each quadrature point and the difference step size.
This is the constructor for a variable velocity field. Most probably, a new object of FunctionDerivative
has to be constructed for each set of quadrature points.
The number of quadrature point must still be the same, when values are accessed.
Definition at line 40 of file function_derivative.cc.
void FunctionDerivative< dim >::set_formula | ( | typename AutoDerivativeFunction< dim >::DifferenceFormula | formula = AutoDerivativeFunction<dim>::Euler | ) |
Choose the difference formula. This is set to the default in the constructor.
Formulas implemented right now are first order backward Euler (UpwindEuler
), second order symmetric Euler (Euler
) and a symmetric fourth order formula (FourthOrder
).
Definition at line 58 of file function_derivative.cc.
void FunctionDerivative< dim >::set_h | ( | const double | h | ) |
Change the base step size of the difference formula
Definition at line 80 of file function_derivative.cc.
|
virtual |
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 >.
Definition at line 91 of file function_derivative.cc.
|
virtual |
Return all components of a vector-valued function at a given point.
values
shall have the right size beforehand, i.e. n_components.
The default implementation will call value() for each component.
Reimplemented from Function< dim >.
Definition at line 116 of file function_derivative.cc.
|
virtual |
Set values
to the point values of the specified component of the function at the points
. It is assumed that values
already has the right size, i.e. the same size as the points
array.
By default, this function repeatedly calls value() for each point separately, to fill the output array.
Reimplemented from Function< dim >.
Definition at line 158 of file function_derivative.cc.
std::size_t FunctionDerivative< dim >::memory_consumption | ( | ) | const |
Return an estimate for the memory consumption, in bytes, of this object. This is not exact (but will usually be close) because calculating the memory usage of trees (e.g., std::map
) is difficult.
Definition at line 224 of file function_derivative.cc.
|
private |
Function for differentiation.
Definition at line 108 of file function_derivative.h.
|
private |
Step size of the difference formula.
Definition at line 113 of file function_derivative.h.
|
private |
Difference formula.
Definition at line 118 of file function_derivative.h.
|
private |
Helper object. Contains the increment vector for the formula.
Definition at line 123 of file function_derivative.h.