Reference documentation for deal.II version 9.6.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_bessel.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2010 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_function_bessel_h
16#define dealii_function_bessel_h
17
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/point.h>
23
25
26namespace Functions
27{
33 template <int dim>
34 class Bessel1 : public Function<dim>
35 {
36 public:
40 Bessel1(const unsigned int order,
41 const double wave_number,
42 const Point<dim> center = Point<dim>());
43
44 virtual double
45 value(const Point<dim> &points,
46 const unsigned int component = 0) const override;
47
48 virtual void
49 value_list(const std::vector<Point<dim>> &points,
50 std::vector<double> &values,
51 const unsigned int component = 0) const override;
52
53 virtual Tensor<1, dim>
54 gradient(const Point<dim> &p,
55 const unsigned int component = 0) const override;
56
57 virtual void
58 gradient_list(const std::vector<Point<dim>> &points,
59 std::vector<Tensor<1, dim>> &gradients,
60 const unsigned int component = 0) const override;
61
62 private:
63 unsigned int order;
66 };
67} // namespace Functions
68
70
71#endif
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< double > &values, const unsigned int component=0) const override
virtual double value(const Point< dim > &points, const unsigned int component=0) const override
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const override
Bessel1(const unsigned int order, const double wave_number, const Point< dim > center=Point< dim >())
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504