Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Member Functions | Private Member Functions | Private Attributes | List of all members
Functions::Spherical< dim > Class Template Reference

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

Inheritance diagram for Functions::Spherical< dim >:
[legend]

Public Member Functions

 Spherical (const Point< dim > &center=Point< dim >(), const unsigned int n_components=1)
 
virtual double value (const Point< dim > &point, const unsigned int component=0) const override
 
virtual Tensor< 1, dim > gradient (const Point< dim > &p, const unsigned int component=0) const override
 
virtual SymmetricTensor< 2, dim > hessian (const Point< dim > &p, const unsigned int component=0) const override
 
std::size_t memory_consumption () const
 
Tensor< 1, 3 > gradient (const Point< 3 > &p_, const unsigned int component) const
 
SymmetricTensor< 2, 3 > hessian (const Point< 3 > &p_, const unsigned int component) const
 
- Public Member Functions inherited from Function< dim >
 Function (const unsigned int n_components=1, const time_type initial_time=0.0)
 
 Function (const Function &f)=default
 
virtual ~Function () override=0
 
Functionoperator= (const Function &f)
 
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 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 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< numbers::NumberTraits< double >::real_type >
 FunctionTime (const numbers::NumberTraits< double >::real_type initial_time=numbers::NumberTraits< double >::real_type(0.0))
 
virtual ~FunctionTime ()=default
 
numbers::NumberTraits< double >::real_type get_time () const
 
virtual void set_time (const numbers::NumberTraits< double >::real_type new_time)
 
virtual void advance_time (const numbers::NumberTraits< double >::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 Member Functions

virtual double svalue (const std::array< double, dim > &sp, const unsigned int component) const
 
virtual std::array< double, dim > sgradient (const std::array< double, dim > &sp, const unsigned int component) const
 
virtual std::array< double, 6 > shessian (const std::array< double, dim > &sp, const unsigned int component) const
 

Private Attributes

const Tensor< 1, dim > coordinate_system_offset
 

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< numbers::NumberTraits< double >::real_type >
using time_type = numbers::NumberTraits< double >::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 >
const unsigned int n_components
 
- Static Public Attributes inherited from Function< dim >
static const unsigned int dimension
 

Detailed Description

template<int dim>
class Functions::Spherical< dim >

An abstract base class for a scalar-valued function \(f=f(r,\theta,\phi)\) defined in spherical coordinates. This class wraps transformation of values, gradients and hessians from spherical coordinates to the Cartesian coordinate system used by the Function base class. Therefore derived classes only need to implement those functions in spherical coordinates (specifically svalue(), sgradient() and shessian() ). The convention for angles is the same as in GeometricUtilities::Coordinates.

Note
This function is currently only implemented for dim==3 .
Author
Denis Davydov, 2017

Definition at line 45 of file function_spherical.h.

Constructor & Destructor Documentation

◆ Spherical()

template<int dim>
Functions::Spherical< dim >::Spherical ( const Point< dim > &  center = Point<dim>(),
const unsigned int  n_components = 1 
)

Constructor which should be provided with center defining the origin of the coordinate system.

Note that components of this function are treated as entirely separate quantities – not as the components of a vector that will be re-interpreted in a different coordinate system.

Definition at line 156 of file function_spherical.cc.

Member Function Documentation

◆ value()

template<int dim>
double Functions::Spherical< dim >::value ( const Point< dim > &  point,
const unsigned int  component = 0 
) const
overridevirtual

Return the value of the function at the given point.

This function converts the given point to spherical coordinates, calls svalue() with it, and returns the result.

Reimplemented from Function< dim >.

Definition at line 168 of file function_spherical.cc.

◆ gradient() [1/2]

template<int dim>
Tensor< 1, dim > Functions::Spherical< dim >::gradient ( const Point< dim > &  p,
const unsigned int  component = 0 
) const
overridevirtual

Return the gradient with respect to the Cartesian coordinates at point p.

This function converts the given point to spherical coordinates, calls sgradient() with it, and converts the result into Cartesian coordinates.

Reimplemented from Function< dim >.

Definition at line 181 of file function_spherical.cc.

◆ hessian() [1/2]

template<int dim>
SymmetricTensor< 2, dim > Functions::Spherical< dim >::hessian ( const Point< dim > &  p,
const unsigned int  component = 0 
) const
overridevirtual

Return the Hessian with respect to the Cartesian coordinates at point p.

This function converts the given point to spherical coordinates, calls sgradient and shessian() with it, and converts the result into Cartesian coordinates.

Reimplemented from Function< dim >.

Definition at line 237 of file function_spherical.cc.

◆ memory_consumption()

template<int dim>
std::size_t Functions::Spherical< dim >::memory_consumption

Definition at line 308 of file function_spherical.cc.

◆ svalue()

template<int dim>
double Functions::Spherical< dim >::svalue ( const std::array< double, dim > &  sp,
const unsigned int  component 
) const
privatevirtual

Return the value at point sp. Here, sp is provided in spherical coordinates.

Definition at line 317 of file function_spherical.cc.

◆ sgradient()

template<int dim>
std::array< double, dim > Functions::Spherical< dim >::sgradient ( const std::array< double, dim > &  sp,
const unsigned int  component 
) const
privatevirtual

Return the gradient in spherical coordinates.

The returned object should contain derivatives in the following order: \(\{ f_{,r},\, f_{,\theta},\, f_{,\phi}\}\).

Definition at line 328 of file function_spherical.cc.

◆ shessian()

template<int dim>
std::array< double, 6 > Functions::Spherical< dim >::shessian ( const std::array< double, dim > &  sp,
const unsigned int  component 
) const
privatevirtual

Return the Hessian in spherical coordinates.

The returned object should contain derivatives in the following order: \(\{ f_{,rr},\, f_{,\theta\theta},\, f_{,\phi\phi},\, f_{,r\theta},\, f_{,r\phi},\, f_{,\theta\phi}\}\).

Definition at line 339 of file function_spherical.cc.

◆ gradient() [2/2]

Tensor< 1, 3 > Functions::Spherical< 3 >::gradient ( const Point< 3 > &  p_,
const unsigned int  component 
) const

Definition at line 193 of file function_spherical.cc.

◆ hessian() [2/2]

SymmetricTensor< 2, 3 > Functions::Spherical< 3 >::hessian ( const Point< 3 > &  p_,
const unsigned int  component 
) const

Definition at line 248 of file function_spherical.cc.

Member Data Documentation

◆ coordinate_system_offset

template<int dim>
const Tensor<1, dim> Functions::Spherical< dim >::coordinate_system_offset
private

A vector from the origin to the center of spherical coordinate system.

Definition at line 127 of file function_spherical.h.


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