Line data Source code
1 1 : // Distributed under the MIT License.
2 : // See LICENSE.txt for details.
3 :
4 : /*!
5 : * \file
6 : * Declares functionality to retrieve spectral quantities associated with
7 : * a particular choice of basis functions and quadrature.
8 : */
9 :
10 : #pragma once
11 :
12 : #include <cstddef>
13 : #include <iosfwd>
14 : #include <limits>
15 : #include <utility>
16 :
17 : #include "NumericalAlgorithms/Spectral/Basis.hpp"
18 : #include "NumericalAlgorithms/Spectral/Quadrature.hpp"
19 :
20 : /// \cond
21 : class Matrix;
22 : class DataVector;
23 : template <size_t>
24 : class Mesh;
25 : namespace Options {
26 : class Option;
27 : template <typename T>
28 : struct create_from_yaml;
29 : } // namespace Options
30 : /// \endcond
31 :
32 : /*!
33 : * \ingroup SpectralGroup
34 : * \brief Functionality associated with a particular choice of basis functions
35 : * and quadrature for spectral operations.
36 : *
37 : * \details The functions in this namespace provide low-level access to
38 : * collocation points, quadrature weights and associated matrices, such as for
39 : * differentiation and interpolation. They are available in two versions: either
40 : * templated directly on the enum cases of the Spectral::Basis and
41 : * Spectral::Quadrature types, or taking a one-dimensional Mesh as their
42 : * argument.
43 : *
44 : * \note Generally you should prefer working with a Mesh. Use its
45 : * Mesh::slice_through() method to retrieve the mesh in a particular dimension:
46 : * \snippet Test_Spectral.cpp get_points_for_mesh
47 : *
48 : *
49 : * Most algorithms in this namespace are adapted from \cite Kopriva.
50 : */
51 : namespace Spectral {
52 : namespace detail {
53 : constexpr size_t minimum_number_of_points(const Basis /*basis*/,
54 : const Quadrature quadrature) {
55 : // NOLINTNEXTLINE(bugprone-branch-clone)
56 : if (quadrature == Quadrature::Gauss) {
57 : return 1;
58 : // NOLINTNEXTLINE(bugprone-branch-clone)
59 : } else if (quadrature == Quadrature::GaussLobatto) {
60 : return 2;
61 : // NOLINTNEXTLINE(bugprone-branch-clone)
62 : } else if (quadrature == Quadrature::CellCentered) {
63 : return 1;
64 : // NOLINTNEXTLINE(bugprone-branch-clone)
65 : } else if (quadrature == Quadrature::FaceCentered) {
66 : return 2;
67 : } else if (quadrature == Quadrature::Equiangular) {
68 : return 1;
69 : }
70 : return std::numeric_limits<size_t>::max();
71 : }
72 : } // namespace detail
73 :
74 : /*!
75 : * \brief Minimum number of possible collocation points for a quadrature type.
76 : *
77 : * \details Since Gauss-Lobatto quadrature has points on the domain boundaries
78 : * it must have at least two collocation points. Gauss quadrature can have only
79 : * one collocation point.
80 : *
81 : * \details For `CellCentered` the minimum number of points is 1, while for
82 : * `FaceCentered` it is 2.
83 : */
84 : template <Basis basis, Quadrature quadrature>
85 1 : constexpr size_t minimum_number_of_points =
86 : detail::minimum_number_of_points(basis, quadrature);
87 :
88 : /*!
89 : * \brief Maximum number of allowed collocation points.
90 : *
91 : * \details We choose a limit of 24 FD grid points because for DG-subcell the
92 : * number of points in an element is `2 * number_dg_points - 1` for cell
93 : * centered, and `2 * number_dg_points` for face-centered. Because there is no
94 : * way of generically retrieving the maximum number of grid points for a non-FD
95 : * basis, we need to hard-code both values here. If the number of grid points is
96 : * increased for the non-FD bases, it should also be increased for the FD basis.
97 : * Note that for good task-based parallelization 24 grid points is already a
98 : * fairly large number.
99 : */
100 : template <Basis basis>
101 1 : constexpr size_t maximum_number_of_points =
102 : basis == Basis::FiniteDifference ? 40 : 20;
103 :
104 : /*!
105 : * \brief Compute the function values of the basis function \f$\Phi_k(x)\f$
106 : * (zero-indexed).
107 : */
108 : template <Basis BasisType, typename T>
109 1 : T compute_basis_function_value(size_t k, const T& x);
110 :
111 : /*!
112 : * \brief Compute the inverse of the weight function \f$w(x)\f$ w.r.t. which
113 : * the basis functions are orthogonal. See the description of
114 : * `quadrature_weights(size_t)` for details.
115 : * This is arbitrarily set to 1 for FiniteDifference basis, to integrate
116 : * using the midpoint method (see `quadrature_weights (size_t)` for details).
117 : */
118 : template <Basis>
119 1 : DataVector compute_inverse_weight_function_values(const DataVector&);
120 :
121 : /*!
122 : * \brief Compute the normalization square of the basis function \f$\Phi_k\f$
123 : * (zero-indexed), i.e. the weighted definite integral over its square.
124 : */
125 : template <Basis BasisType>
126 1 : double compute_basis_function_normalization_square(size_t k);
127 :
128 : /*!
129 : * \brief Compute the collocation points and weights associated to the
130 : * basis and quadrature.
131 : *
132 : * \details This function is expected to return the tuple
133 : * \f$(\xi_k,w_k)\f$ where the \f$\xi_k\f$ are the collocation
134 : * points and the \f$w_k\f$ are defined in the description of
135 : * `quadrature_weights(size_t)`.
136 : *
137 : * \warning for a `FiniteDifference` basis or `CellCentered` and `FaceCentered`
138 : * quadratures, the weights are defined to integrate with the midpoint method
139 : */
140 : template <Basis BasisType, Quadrature QuadratureType>
141 1 : std::pair<DataVector, DataVector> compute_collocation_points_and_weights(
142 : size_t num_points);
143 :
144 : /*!
145 : * \brief Collocation points
146 : * \param num_points The number of collocation points
147 : */
148 : template <Basis BasisType, Quadrature QuadratureType>
149 1 : const DataVector& collocation_points(size_t num_points);
150 :
151 : /*!
152 : * \brief Collocation points for a one-dimensional mesh.
153 : *
154 : * \see collocation_points(size_t)
155 : */
156 1 : const DataVector& collocation_points(const Mesh<1>& mesh);
157 :
158 : /*!
159 : * \brief Weights to compute definite integrals.
160 : *
161 : * \details These are the coefficients to contract with the nodal
162 : * function values \f$f_k\f$ to approximate the definite integral \f$I[f]=\int
163 : * f(x)\mathrm{d}x\f$.
164 : *
165 : * Note that the term _quadrature_ also often refers to the quantity
166 : * \f$Q[f]=\int f(x)w(x)\mathrm{d}x\approx \sum_k f_k w_k\f$. Here, \f$w(x)\f$
167 : * denotes the basis-specific weight function w.r.t. to which the basis
168 : * functions \f$\Phi_k\f$ are orthogonal, i.e \f$\int\Phi_i(x)\Phi_j(x)w(x)=0\f$
169 : * for \f$i\neq j\f$. The weights \f$w_k\f$ approximate this inner product. To
170 : * approximate the definite integral \f$I[f]\f$ we must employ the
171 : * coefficients \f$\frac{w_k}{w(\xi_k)}\f$ instead, where the \f$\xi_k\f$ are
172 : * the collocation points. These are the coefficients this function returns.
173 : * Only for a unit weight function \f$w(x)=1\f$, i.e. a Legendre basis, is
174 : * \f$I[f]=Q[f]\f$ so this function returns the \f$w_k\f$ identically.
175 : *
176 : * For a `FiniteDifference` basis or `CellCentered` and `FaceCentered`
177 : * quadratures, the interpretation of the quadrature weights in term
178 : * of an approximation to \f$I(q)\f$ remains correct, but its explanation
179 : * in terms of orthonormal basis is not, i.e. we set \f$w_k\f$ to the grid
180 : * spacing at each point, and the inverse weight \f$\frac{1}{w(\xi_k)}=1\f$ to
181 : * recover the midpoint method for definite integrals.
182 : *
183 : * \param num_points The number of collocation points
184 : */
185 : template <Basis BasisType, Quadrature QuadratureType>
186 1 : const DataVector& quadrature_weights(size_t num_points);
187 :
188 : /*!
189 : * \brief Quadrature weights for a one-dimensional mesh.
190 : *
191 : * \see quadrature_weights(size_t)
192 : */
193 1 : const DataVector& quadrature_weights(const Mesh<1>& mesh);
194 :
195 : /// @{
196 : /*!
197 : * \brief %Matrix used to compute the derivative of a function.
198 : *
199 : * \details For a function represented by the nodal coefficients \f$u_j\f$ a
200 : * matrix multiplication with the differentiation matrix \f$D_{ij}\f$ gives the
201 : * coefficients of the function's derivative. Since \f$u(x)\f$ is expanded in
202 : * Lagrange polynomials \f$u(x)=\sum_j u_j l_j(x)\f$ the differentiation matrix
203 : * is computed as \f$D_{ij}=l_j^\prime(\xi_i)\f$ where the \f$\xi_i\f$ are the
204 : * collocation points.
205 : *
206 : * The finite difference matrix uses summation by parts operators,
207 : * \f$D_{2-1}, D_{4-2}, D_{4-3}\f$, and \f$D_{6-5}\f$ from \cite Diener2005tn.
208 : *
209 : * \param num_points The number of collocation points
210 : */
211 : template <Basis BasisType, Quadrature QuadratureType>
212 1 : const Matrix& differentiation_matrix(size_t num_points);
213 : template <Basis BasisType, Quadrature QuadratureType>
214 1 : const Matrix& differentiation_matrix_transpose(size_t num_points);
215 : /// @}
216 :
217 : /// @{
218 : /*!
219 : * \brief Differentiation matrix for a one-dimensional mesh.
220 : *
221 : * \see differentiation_matrix(size_t)
222 : */
223 1 : const Matrix& differentiation_matrix(const Mesh<1>& mesh);
224 1 : const Matrix& differentiation_matrix_transpose(const Mesh<1>& mesh);
225 : /// @}
226 :
227 : /*!
228 : * \brief %Matrix used to compute the divergence of the flux in weak form.
229 : *
230 : * This is the transpose of the differentiation matrix multiplied by quadrature
231 : * weights that appear in DG integrals:
232 : *
233 : * \begin{equation}
234 : * \frac{D^T_{ij}} \frac{w_j}{w_i}
235 : * \end{equation}
236 : *
237 : * \param num_points The number of collocation points
238 : */
239 : template <Basis BasisType, Quadrature QuadratureType>
240 1 : const Matrix& weak_flux_differentiation_matrix(size_t num_points);
241 :
242 : /*!
243 : * \brief %Matrix used to compute the divergence of the flux in weak form.
244 : *
245 : * \see weak_flux_differentiation_matrix(size_t)
246 : */
247 1 : const Matrix& weak_flux_differentiation_matrix(const Mesh<1>& mesh);
248 :
249 : /*!
250 : * \brief %Matrix used to perform an indefinite integral of a function over the
251 : * logical grid. The left boundary condition is such that the integral is 0 at
252 : * \f$\xi=-1\f$
253 : *
254 : * Currently only Legendre and Chebyshev polynomials are implemented, but we
255 : * provide a derivation for how to compute the indefinite integration matrix for
256 : * general Jacobi polynomials.
257 : *
258 : * #### Legendre Polynomials
259 : * The Legendre polynomials have the identity:
260 : *
261 : * \f{align*}{
262 : * P_n(x) = \frac{1}{2n+1}\frac{d}{dx}\left(P_{n+1}(x) - P_{n-1}(x)\right)
263 : * \f}
264 : *
265 : * The goal is to evaluate the integral of a function \f$u\f$ expanded in terms
266 : * of Legendre polynomials as
267 : *
268 : * \f{align*}{
269 : * u(x) = \sum_{i=0}^{N} c_i P_i(x)
270 : * \f}
271 : *
272 : * We similarly expand the indefinite integral of \f$u\f$ as
273 : *
274 : * \f{align*}{
275 : * \left.\int u(y) dy\right\rvert_{x}=&\sum_{i=0}^N \tilde{c}_i P_i(x) \\
276 : * =&\left.\int\sum_{i=1}^{N}\frac{c_i}{2i+1}
277 : * \left(P_{i+1}(y)-P_{i-1}(y)\right)dy\right\rvert_{x}
278 : * + \tilde{c}_0 P_0(x) \\
279 : * =&\sum_{i=1}^{N}\left(\frac{c_{i-1}}{2i-1} - \frac{c_{i+1}}{2i+3}\right)
280 : * P_i(x) + \tilde{c}_0 P_0(x)
281 : * \f}
282 : *
283 : * Thus we get that for \f$i>0\f$
284 : *
285 : * \f{align*}{
286 : * \tilde{c}_i=\frac{c_{i-1}}{2i-1}-\frac{c_{i+1}}{2i+3}
287 : * \f}
288 : *
289 : * and \f$\tilde{c}_0\f$ is a constant of integration, which we choose such that
290 : * the integral is 0 at the left boundary of the domain (\f$x=-1\f$). The
291 : * condition for this is:
292 : *
293 : * \f{align*}{
294 : * \tilde{c}_0=\sum_{i=1}^{N}(-1)^{i+1}\tilde{c}_i
295 : * \f}
296 : *
297 : * The matrix returned by this function is the product of the tridiagonal matrix
298 : * for the \f$\tilde{c}_i\f$ and the matrix for the boundary condition.
299 : *
300 : * #### Chebyshev Polynomials
301 : *
302 : * A similar derivation leads to the relations:
303 : *
304 : * \f{align*}{
305 : * \tilde{c}_i=&\frac{c_{i-1}-c_{i+1}}{2i},&\mathrm{if}\;i>1 \\
306 : * \tilde{c}_1=&c_0 - \frac{c_2}{2},&\mathrm{if}\;i=1 \\
307 : * \f}
308 : *
309 : * We again have:
310 : *
311 : * \f{align*}{
312 : * \tilde{c}_0=\sum_{i=1}^N(-1)^{i+1}\tilde{c}_i
313 : * \f}
314 : *
315 : * These are then used to define the indefinite integration matrix.
316 : *
317 : * #### Jacobi Polynomials
318 : *
319 : * For general Jacobi polynomials \f$P^{(\alpha,\beta)}_n(x)\f$ given by
320 : *
321 : * \f{align*}{
322 : * (1-x)^\alpha(1+x)^\beta P^{(\alpha,\beta)}_n(x)=\frac{(-1)^n}{2^n n!}
323 : * \frac{d^n}{dx^n}\left[(1-x)^{\alpha+n}(1+x)^{\beta+n}\right]
324 : * \f}
325 : *
326 : * we have that
327 : *
328 : * \f{align*}{
329 : * P^{(\alpha,\beta)}_n(x)=\frac{d}{dx}\left(
330 : * b^{(\alpha,\beta)}_{n-1,n}P^{(\alpha,\beta)}_{n-1}(x) +
331 : * b^{(\alpha,\beta)}_{n,n}P^{(\alpha,\beta)}_n(x) +
332 : * b^{(\alpha,\beta)}_{n+1,n}P^{(\alpha,\beta)}_{n+1}(x)
333 : * \right)
334 : * \f}
335 : *
336 : * where
337 : *
338 : * \f{align*}{
339 : * b^{(\alpha,\beta)}_{n-1,n}=&-\frac{1}{n+\alpha+\beta}
340 : * a^{(\alpha,\beta)}_{n-1,n} \\
341 : * b^{(\alpha,\beta)}_{n,n}=&-\frac{2}{\alpha+\beta}
342 : * a^{(\alpha,\beta)}_{n,n} \\
343 : * b^{(\alpha,\beta)}_{n+1,n}=&\frac{1}{n+1}
344 : * a^{(\alpha,\beta)}_{n+1,n} \\
345 : * a^{(\alpha,\beta)}_{n-1,n}=&\frac{2(n+\alpha)(n+\beta)}
346 : * {(2n+\alpha+\beta+1)(2n+\alpha+\beta)} \\
347 : * a^{(\alpha,\beta)}_{n,n}=&-\frac{\alpha^2-\beta^2}
348 : * {(2n+\alpha+\beta+2)(2n+\alpha+\beta)} \\
349 : * a^{(\alpha,\beta)}_{n-1,n}=&\frac{2(n+1)(n+\alpha+\beta+1)}
350 : * {(2n+\alpha+\beta+2)(2n+\alpha+\beta+1)}
351 : * \f}
352 : *
353 : * Following the same derivation we get that
354 : *
355 : * \f{align*}{
356 : * \tilde{c}_i=c_{i+1}b^{(\alpha,\beta)}_{i,i+1}
357 : * +c_i b^{(\alpha,\beta)}_{i,i}
358 : * +c_{i-1}b^{(\alpha,\beta)}_{i,i-1}
359 : * \f}
360 : *
361 : * and the boundary condition is
362 : *
363 : * \f{align*}{
364 : * \tilde{c}_0=\sum_{i=1}^N(-1)^{i+1}
365 : * \frac{\Gamma(i+\alpha+1)}{i!\Gamma(\alpha+1)} \tilde{c}_i
366 : * \f}
367 : *
368 : * where \f$\Gamma(x)\f$ is the Gamma function.
369 : */
370 : template <Basis BasisType, Quadrature QuadratureType>
371 1 : const Matrix& integration_matrix(size_t num_points);
372 :
373 : /*!
374 : * \brief Indefinite integration matrix for a one-dimensional mesh.
375 : *
376 : * \see integration_matrix(size_t)
377 : */
378 1 : const Matrix& integration_matrix(const Mesh<1>& mesh);
379 :
380 : /*!
381 : * \brief %Matrix used to interpolate to the \p target_points.
382 : *
383 : * \warning For each target point located outside of the logical coordinate
384 : * interval covered by `BasisType` (often \f$[-1,1]\f$), the resulting matrix
385 : * performs polynomial extrapolation instead of interpolation. The extrapolation
386 : * will be correct but may suffer from reduced accuracy, especially for
387 : * higher-order polynomials (i.e., larger values of `num_points`).
388 : *
389 : * \param num_points The number of collocation points
390 : * \param target_points The points to interpolate to
391 : */
392 : template <Basis BasisType, Quadrature QuadratureType, typename T>
393 1 : Matrix interpolation_matrix(size_t num_points, const T& target_points);
394 :
395 : /*!
396 : * \brief Interpolation matrix to the \p target_points for a one-dimensional
397 : * mesh.
398 : *
399 : * \see interpolation_matrix(size_t, const T&)
400 : */
401 : template <typename T>
402 1 : Matrix interpolation_matrix(const Mesh<1>& mesh, const T& target_points);
403 :
404 : /// @{
405 : /*!
406 : * \brief Matrices that interpolate to the lower and upper boundaries of the
407 : * element.
408 : *
409 : * Assumes that the logical coordinates are \f$[-1, 1]\f$. The first element of
410 : * the pair interpolates to \f$\xi=-1\f$ and the second to \f$\xi=1\f$. These
411 : * are just the Lagrange interpolating polynomials evaluated at \f$\xi=\pm1\f$.
412 : * For Gauss-Lobatto points the only non-zero element is at the boundaries
413 : * and is one and so is not implemented.
414 : *
415 : * \warning This can only be called with Gauss points.
416 : */
417 1 : const std::pair<Matrix, Matrix>& boundary_interpolation_matrices(
418 : const Mesh<1>& mesh);
419 :
420 : template <Basis BasisType, Quadrature QuadratureType>
421 1 : const std::pair<Matrix, Matrix>& boundary_interpolation_matrices(
422 : size_t num_points);
423 : /// @}
424 :
425 : /// @{
426 : /*!
427 : * \brief Interpolates values from the boundary into the volume, which is needed
428 : * when applying time derivative or Bjorhus-type boundary conditions in a
429 : * discontinuous Galerkin scheme using Gauss points.
430 : *
431 : * Assumes that the logical coordinates are \f$[-1, 1]\f$.
432 : * The interpolation is done by assuming the time derivative correction is zero
433 : * on interior nodes. With a nodal Lagrange polynomial basis this means that
434 : * only the \f$\ell^{\mathrm{Gauss-Lobatto}}_{0}\f$ and
435 : * \f$\ell^{\mathrm{Gauss-Lobatto}}_{N}\f$ polynomials/basis functions
436 : * contribute to the correction. In order to interpolate the correction from the
437 : * nodes at the boundary, the Gauss-Lobatto Lagrange polynomials must be
438 : * evaluated at the Gauss grid points. The returned pair of `DataVector`s stores
439 : *
440 : * \f{align*}{
441 : * &\ell^{\mathrm{Gauss-Lobatto}}_{0}(\xi_j^{\mathrm{Gauss}}), \\
442 : * &\ell^{\mathrm{Gauss-Lobatto}}_{N}(\xi_j^{\mathrm{Gauss}}).
443 : * \f}
444 : *
445 : * This is a different correction from lifting. Lifting is done using the mass
446 : * matrix, which is an integral over the basis functions, while here we use
447 : * interpolation.
448 : *
449 : * \warning This can only be called with Gauss points.
450 : */
451 1 : const std::pair<DataVector, DataVector>& boundary_interpolation_term(
452 : const Mesh<1>& mesh);
453 :
454 : template <Basis BasisType, Quadrature QuadratureType>
455 1 : const std::pair<DataVector, DataVector>& boundary_interpolation_term(
456 : size_t num_points);
457 : /// @}
458 :
459 : /// @{
460 : /*!
461 : * \brief Terms used during the lifting portion of a discontinuous Galerkin
462 : * scheme when using Gauss points.
463 : *
464 : * Assumes that the logical coordinates are \f$[-1, 1]\f$. The first element of
465 : * the pair is the Lagrange polyonmials evaluated at \f$\xi=-1\f$ divided by the
466 : * weights and the second at \f$\xi=1\f$. Specifically,
467 : *
468 : * \f{align*}{
469 : * \frac{\ell_j(\xi=\pm1)}{w_j}
470 : * \f}
471 : *
472 : * \warning This can only be called with Gauss points.
473 : */
474 1 : const std::pair<DataVector, DataVector>& boundary_lifting_term(
475 : const Mesh<1>& mesh);
476 :
477 : template <Basis BasisType, Quadrature QuadratureType>
478 1 : const std::pair<DataVector, DataVector>& boundary_lifting_term(
479 : size_t num_points);
480 : /// @}
481 :
482 : /*!
483 : * \brief %Matrix used to transform from the spectral coefficients (modes) of a
484 : * function to its nodal coefficients. Also referred to as the _Vandermonde
485 : * matrix_.
486 : *
487 : * \details The Vandermonde matrix is computed as
488 : * \f$\mathcal{V}_{ij}=\Phi_j(\xi_i)\f$ where the \f$\Phi_j(x)\f$ are the
489 : * spectral basis functions used in the modal expansion
490 : * \f$u(x)=\sum_j \widetilde{u}_j\Phi_j(x)\f$, e.g. normalized Legendre
491 : * polynomials, and the \f$\xi_i\f$ are the collocation points used to construct
492 : * the interpolating Lagrange polynomials in the nodal expansion
493 : * \f$u(x)=\sum_j u_j l_j(x)\f$. Then the Vandermonde matrix arises as
494 : * \f$u(\xi_i)=u_i=\sum_j \widetilde{u}_j\Phi_j(\xi_i)=\sum_j
495 : * \mathcal{V}_{ij}\widetilde{u}_j\f$.
496 : *
497 : * \param num_points The number of collocation points
498 :
499 : * \see nodal_to_modal_matrix(size_t)
500 : */
501 : template <Basis BasisType, Quadrature QuadratureType>
502 1 : const Matrix& modal_to_nodal_matrix(size_t num_points);
503 :
504 : /*!
505 : * \brief Transformation matrix from modal to nodal coefficients for a
506 : * one-dimensional mesh.
507 : *
508 : * \see modal_to_nodal_matrix(size_t)
509 : */
510 1 : const Matrix& modal_to_nodal_matrix(const Mesh<1>& mesh);
511 :
512 : /*!
513 : * \brief %Matrix used to transform from the nodal coefficients of a function to
514 : * its spectral coefficients (modes). Also referred to as the inverse
515 : * _Vandermonde matrix_.
516 : *
517 : * \details This is the inverse to the Vandermonde matrix \f$\mathcal{V}\f$
518 : * computed in modal_to_nodal_matrix(size_t). It can be computed
519 : * analytically for Gauss quadrature by evaluating
520 : * \f$\sum_j\mathcal{V}^{-1}_{ij}u_j=\widetilde{u}_i=
521 : * \frac{(u,\Phi_i)}{\gamma_i}\f$
522 : * for a Lagrange basis function \f$u(x)=l_k(x)\f$ to find
523 : * \f$\mathcal{V}^{-1}_{ij}=\mathcal{V}_{ji}\frac{w_j}{\gamma_i}\f$ where the
524 : * \f$w_j\f$ are the Gauss quadrature weights and \f$\gamma_i\f$ is the norm
525 : * square of the spectral basis function \f$\Phi_i\f$.
526 : *
527 : * \param num_points The number of collocation points
528 : *
529 : * \see modal_to_nodal_matrix(size_t)
530 : */
531 : template <Basis BasisType, Quadrature QuadratureType>
532 1 : const Matrix& nodal_to_modal_matrix(size_t num_points);
533 :
534 : /*!
535 : * \brief Transformation matrix from nodal to modal coefficients for a
536 : * one-dimensional mesh.
537 : *
538 : * \see nodal_to_modal_matrix(size_t)
539 : */
540 1 : const Matrix& nodal_to_modal_matrix(const Mesh<1>& mesh);
541 :
542 : /*!
543 : * \brief %Matrix used to linearize a function.
544 : *
545 : * \details Filters out all except the lowest two modes by applying
546 : * \f$\mathcal{V}^{-1}\cdot\mathrm{diag}(1,1,0,0,...)\cdot\mathcal{V}\f$ to the
547 : * nodal coefficients, where \f$\mathcal{V}\f$ is the Vandermonde matrix
548 : * computed in `modal_to_nodal_matrix(size_t)`.
549 : *
550 : * \param num_points The number of collocation points
551 : *
552 : * \see modal_to_nodal_matrix(size_t)
553 : * \see nodal_to_modal_matrix(size_t)
554 : */
555 : template <Basis BasisType, Quadrature QuadratureType>
556 1 : const Matrix& linear_filter_matrix(size_t num_points);
557 :
558 : /*!
559 : * \brief Linear filter matrix for a one-dimensional mesh.
560 : *
561 : * \see linear_filter_matrix(size_t)
562 : */
563 1 : const Matrix& linear_filter_matrix(const Mesh<1>& mesh);
564 :
565 : } // namespace Spectral
|