Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/function.h>
Public Member Functions | |
ComponentSelectFunction (const unsigned int selected, const RangeNumberType value, const unsigned int n_components) | |
ComponentSelectFunction (const unsigned int selected, const unsigned int n_components) | |
ComponentSelectFunction (const std::pair< unsigned int, unsigned int > &selected, const unsigned int n_components) | |
virtual void | substitute_function_value_with (const Functions::ConstantFunction< dim, RangeNumberType > &f) |
virtual void | vector_value (const Point< dim > &p, Vector< RangeNumberType > &return_value) const override |
virtual void | vector_value_list (const std::vector< Point< dim >> &points, std::vector< Vector< RangeNumberType >> &values) const override |
std::size_t | memory_consumption () const |
Public Member Functions inherited from Functions::ConstantFunction< dim, RangeNumberType > | |
ConstantFunction (const RangeNumberType value, const unsigned int n_components=1) | |
ConstantFunction (const std::vector< RangeNumberType > &values) | |
ConstantFunction (const Vector< RangeNumberType > &values) | |
ConstantFunction (const RangeNumberType *begin_ptr, const unsigned int n_components) | |
virtual RangeNumberType | value (const Point< dim > &p, const unsigned int component=0) const override |
virtual void | value_list (const std::vector< Point< dim >> &points, std::vector< RangeNumberType > &return_values, const unsigned int component=0) const override |
virtual Tensor< 1, dim, RangeNumberType > | gradient (const Point< dim > &p, const unsigned int component=0) const override |
virtual void | vector_gradient (const Point< dim > &p, std::vector< Tensor< 1, dim, RangeNumberType >> &gradients) const override |
virtual void | gradient_list (const std::vector< Point< dim >> &points, std::vector< Tensor< 1, dim, RangeNumberType >> &gradients, const unsigned int component=0) const override |
virtual void | vector_gradient_list (const std::vector< Point< dim >> &points, std::vector< std::vector< Tensor< 1, dim, RangeNumberType >>> &gradients) const override |
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 |
Function & | operator= (const Function &f) |
virtual void | vector_values (const std::vector< Point< dim >> &points, std::vector< std::vector< RangeNumberType >> &values) const |
virtual void | vector_gradients (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 () |
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 Attributes | |
const std::pair< unsigned int, unsigned int > | selected_components |
Protected Attributes inherited from Functions::ConstantFunction< dim, RangeNumberType > | |
std::vector< RangeNumberType > | function_value_vector |
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 ::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 is a constant vector-valued function, in which one or more components of the vector have a constant value and all other components are zero. It is especially useful as a weight function for VectorTools::integrate_difference, where it allows to integrate only one or a few vector components, rather than the entire vector-valued solution. In other words, it acts as a component mask with a single component selected (see the the glossary entry on component masks). See the step-20 tutorial program for a detailed explanation and a use case.
Definition at line 544 of file function.h.
ComponentSelectFunction< dim, RangeNumberType >::ComponentSelectFunction | ( | const unsigned int | selected, |
const RangeNumberType | value, | ||
const unsigned int | n_components | ||
) |
Constructor if only a single component shall be non-zero. Arguments denote the component selected, the value for that component and the total number of vector components.
ComponentSelectFunction< dim, RangeNumberType >::ComponentSelectFunction | ( | const unsigned int | selected, |
const unsigned int | n_components | ||
) |
Constructor. As before, but the value for the selected component is assumed to be one. In essence, this function then works as a mask.
ComponentSelectFunction< dim, RangeNumberType >::ComponentSelectFunction | ( | const std::pair< unsigned int, unsigned int > & | selected, |
const unsigned int | n_components | ||
) |
Constructor if multiple components shall have non-zero, unit values (i.e. this should be a mask for multiple components). The first argument denotes a half-open interval of components (for example std::pair(0,dim) for the first dim components), and the second argument is the total number of vector components.
|
virtual |
Substitute function value with value of a ConstantFunction<dim, RangeNumberType>
object and keep the current selection pattern.
This is useful if you want to have different values in different components since the provided constructors of ComponentSelectFunction<dim, RangeNumberType>
class can only have same value for all components.
f
from its beginning. So the number of components of f
cannot be less than the calling object.
|
overridevirtual |
Return the value of the function at the given point for all components.
Reimplemented from Functions::ConstantFunction< dim, RangeNumberType >.
|
overridevirtual |
Set values
to the point values of the function at the points
, for all components. It is assumed that values
already has the right size, i.e. the same size as the points
array.
Reimplemented from Functions::ConstantFunction< dim, RangeNumberType >.
std::size_t ComponentSelectFunction< dim, RangeNumberType >::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.
|
protected |
Half-open interval of the indices of selected components.
Definition at line 621 of file function.h.