Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
function_spherical.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2016 - 2020 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_function_spherical_h
17#define dealii_function_spherical_h
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/point.h>
23
24#include <array>
25
27
28namespace Functions
29{
43 template <int dim>
44 class Spherical : public Function<dim>
45 {
46 public:
56 const unsigned int n_components = 1);
57
64 virtual double
65 value(const Point<dim> & point,
66 const unsigned int component = 0) const override;
67
75 virtual Tensor<1, dim>
76 gradient(const Point<dim> & p,
77 const unsigned int component = 0) const override;
78
87 hessian(const Point<dim> & p,
88 const unsigned int component = 0) const override;
89
93 virtual std::size_t
94 memory_consumption() const override;
95
96 private:
101 virtual double
102 svalue(const std::array<double, dim> &sp,
103 const unsigned int component) const;
104
111 virtual std::array<double, dim>
112 sgradient(const std::array<double, dim> &sp,
113 const unsigned int component) const;
114
122 virtual std::array<double, 6>
123 shessian(const std::array<double, dim> &sp,
124 const unsigned int component) const;
125
130 };
131} // namespace Functions
132
134
135#endif
const unsigned int n_components
Definition function.h:164
virtual SymmetricTensor< 2, dim > hessian(const Point< dim > &p, const unsigned int component=0) const override
const Tensor< 1, dim > coordinate_system_offset
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual std::array< double, dim > sgradient(const std::array< double, dim > &sp, const unsigned int component) const
virtual double value(const Point< dim > &point, const unsigned int component=0) const override
virtual std::size_t memory_consumption() const override
virtual std::array< double, 6 > shessian(const std::array< double, dim > &sp, const unsigned int component) const
virtual double svalue(const std::array< double, dim > &sp, const unsigned int component) const
Definition point.h:112
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
Point< 3 > center