Reference documentation for deal.II version GIT relicensing-426-g7976cfd195 2024-04-18 21:10:01+00:00
\(\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
polynomials_integrated_legendre_sz.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2017 - 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
17
19
21 : Polynomials::Polynomial<double>(get_coefficients(k))
22{}
23
24
25
26std::vector<double>
28{
29 std::vector<double> coefficients(k + 1);
30
31 // first two polynomials are hard-coded:
32 if (k == 0)
33 {
34 coefficients[0] = -1.;
35 return coefficients;
36 }
37 else if (k == 1)
38 {
39 coefficients[0] = 0.;
40 coefficients[1] = 1.;
41 return coefficients;
42 }
43
44 // General formula is:
45 // k*L_{k}(x) = (2*k-3)*x*L_{k-1} - (k-3)*L_{k-2}.
46 std::vector<double> coefficients_km2 = get_coefficients(k - 2);
47 std::vector<double> coefficients_km1 = get_coefficients(k - 1);
48
49 const double a = 1.0 / k;
50 const double b = 2.0 * k - 3.0;
51 const double c = k - 3.0;
52
53 // To maintain stability, delay the division (multiplication by a) until the
54 // end.
55 for (unsigned int i = 1; i <= k - 2; ++i)
56 {
57 coefficients[i] = b * coefficients_km1[i - 1] - c * coefficients_km2[i];
58 }
59
60 coefficients[0] = -c * coefficients_km2[0];
61 coefficients[k] = b * coefficients_km1[k - 1];
62 coefficients[k - 1] = b * coefficients_km1[k - 2];
63
64 for (double &coefficient : coefficients)
65 {
66 coefficient *= a;
67 }
68
69 return coefficients;
70}
71
72
73
74std::vector<Polynomials::Polynomial<double>>
76{
77 std::vector<Polynomials::Polynomial<double>> v;
78 v.reserve(degree + 1);
79 for (unsigned int i = 0; i <= degree; ++i)
80 {
81 v.push_back(IntegratedLegendreSZ(i));
82 }
83 return v;
84}
85
86
87
static std::vector< Polynomials::Polynomial< double > > generate_complete_basis(const unsigned int degree)
static std::vector< double > get_coefficients(const unsigned int k)
std::vector< double > coefficients
Definition polynomial.h:302
unsigned int degree() const
Definition polynomial.h:806
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503