Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/function_lib.h>
Public Member Functions | |
CutOffFunctionBase (const double radius=1., const Point< dim > center=Point< dim >(), const unsigned int n_components=1, const unsigned int select=CutOffFunctionBase< dim >::no_component, const bool integrate_to_one=false, const double unitary_integral_value=1.0) | |
virtual | ~CutOffFunctionBase ()=default |
void | new_center (const Point< dim > &p) |
void | new_radius (const double r) |
virtual void | set_center (const Point< dim > &p) |
virtual void | set_radius (const double r) |
const Point< dim > & | get_center () const |
double | get_radius () const |
bool | integrates_to_one () const |
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 double | value (const Point< dim > &p, const unsigned int component=0) const |
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 Tensor< 1, dim, double > | gradient (const Point< dim > &p, const unsigned int component=0) 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) |
Static Public Attributes | |
static const unsigned int | no_component = numbers::invalid_unsigned_int |
Static Public Attributes inherited from Function< dim > | |
static const unsigned int | dimension |
Protected Attributes | |
Point< dim > | center |
double | radius |
const unsigned int | selected |
bool | integrate_to_one |
const double | unitary_integral_value |
double | rescaling |
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 |
Base function for cut-off function. This class stores the center and the radius of the supporting ball of a cut-off function. It also stores the number of the non-zero component, if the function is vector-valued.
This class can also be used for approximated Dirac delta functions. These are special cut-off functions whose integral is always equal to one, independently of the radius of the supporting ball.
Definition at line 924 of file function_lib.h.
Functions::CutOffFunctionBase< dim >::CutOffFunctionBase | ( | const double | radius = 1. , |
const Point< dim > | center = Point<dim>() , |
||
const unsigned int | n_components = 1 , |
||
const unsigned int | select = CutOffFunctionBase<dim>::no_component , |
||
const bool | integrate_to_one = false , |
||
const double | unitary_integral_value = 1.0 |
||
) |
Constructor.
[in] | radius | Radius of the ball |
[in] | center | Center of the ball |
[in] | n_components | Number of components of this function object |
[in] | select | If this is different from CutOffFunctionBase<dim>::no_component, then the function will be non-zero for this component only |
[in] | integrate_to_one | Rescale the value of the function whenever a new radius is set, to guarantee that the integral is equal to one |
[in] | unitary_integral_value | Value of the integral when the radius is equal to 1.0. Derived classes will need to supply this value, to guarantee that the rescaling is performed correctly. |
Definition at line 30 of file function_lib_cutoff.cc.
|
virtualdefault |
Virtual destructor.
void Functions::CutOffFunctionBase< dim >::new_center | ( | const Point< dim > & | p | ) |
Move the center of the ball to new point p
.
Definition at line 53 of file function_lib_cutoff.cc.
void Functions::CutOffFunctionBase< dim >::new_radius | ( | const double | r | ) |
Set the radius of the ball to r
.
Definition at line 80 of file function_lib_cutoff.cc.
|
virtual |
Set the center of the ball to the point p
.
Reimplemented in Functions::CutOffFunctionTensorProduct< dim >.
Definition at line 62 of file function_lib_cutoff.cc.
|
virtual |
Set the radius of the ball to r
Reimplemented in Functions::CutOffFunctionTensorProduct< dim >.
Definition at line 89 of file function_lib_cutoff.cc.
const Point< dim > & Functions::CutOffFunctionBase< dim >::get_center | ( | ) | const |
Return the center stored in this object.
Definition at line 71 of file function_lib_cutoff.cc.
double Functions::CutOffFunctionBase< dim >::get_radius | ( | ) | const |
Return the radius stored in this object.
Definition at line 104 of file function_lib_cutoff.cc.
bool Functions::CutOffFunctionBase< dim >::integrates_to_one | ( | ) | const |
Return a boolean indicating whether this function integrates to one.
Definition at line 113 of file function_lib_cutoff.cc.
|
static |
Value used in the constructor of this and derived classes to denote that no component is selected.
Definition at line 931 of file function_lib.h.
|
protected |
Center of the integration ball.
Definition at line 1011 of file function_lib.h.
|
protected |
Radius of the ball.
Definition at line 1016 of file function_lib.h.
|
protected |
Component selected. If no_component
, the function is the same in all components.
Definition at line 1022 of file function_lib.h.
|
protected |
Flag that controls whether we rescale the value when the radius changes.
Definition at line 1027 of file function_lib.h.
|
protected |
The reference integral value. Derived classes should specify what their integral is when radius
= 1.0.
Definition at line 1033 of file function_lib.h.
|
protected |
Current rescaling to apply the cut-off function.
Definition at line 1038 of file function_lib.h.