Quadratures
ClimaCore.Quadratures.QuadratureStyle
— TypeClimaCore.Quadratures.GLL
— TypeGLL{Nq}()
Gauss-Legendre-Lobatto quadrature using Nq
quadrature points.
ClimaCore.Quadratures.GL
— TypeGL{Nq}()
Gauss-Legendre quadrature using Nq
quadrature points.
ClimaCore.Quadratures.Uniform
— TypeUniform{Nq}()
Uniformly-spaced quadrature.
ClimaCore.Quadratures.degrees_of_freedom
— Functiondegrees_of_freedom(QuadratureStyle) -> Int
Returns the degreesoffreedom of the QuadratureStyle
concrete type
ClimaCore.Quadratures.polynomial_degree
— Functionpolynomial_degree(QuadratureStyle) -> Int
Returns the polynomial degree of the QuadratureStyle
concrete type
ClimaCore.Quadratures.quadrature_points
— Functionpoints, weights = quadrature_points(::Type{FT}, quadrature_style)
The points and weights of the quadrature rule in floating point type FT
.
ClimaCore.Quadratures.barycentric_weights
— Functionbarycentric_weights(x::SVector{Nq}) where {Nq}
The barycentric weights associated with the array of point locations x
:
\[w_j = \frac{1}{\prod_{k \ne j} (x_i - x_j)}\]
See [14], equation 3.2.
ClimaCore.Quadratures.interpolation_matrix
— Functioninterpolation_matrix(x::SVector, r::SVector{Nq})
The matrix which interpolates the Lagrange polynomial of degree Nq-1
through the points r
, to points x
. The matrix coefficients are computed using the Barycentric formula of [14], section 4:
\[I_{ij} = \begin{cases} 1 & \text{if } x_i = r_j, \\ 0 & \text{if } x_i = r_k \text{ for } k \ne j, \\ \frac{\displaystyle \frac{w_j}{x_i - r_j}}{\displaystyle \sum_k \frac{w_k}{x_i - r_k}} & \text{otherwise,} \end{cases}\]
where $w_j$ are the barycentric weights, see barycentric_weights
.
ClimaCore.Quadratures.differentiation_matrix
— Functiondifferentiation_matrix(r::SVector{Nq, T}) where {Nq, T}
The spectral differentiation matrix for the Lagrange polynomial of degree Nq-1
interpolating at points r
.
The matrix coefficients are computed using the [14], section 9.3:
\[D_{ij} = \begin{cases} \displaystyle \frac{w_j}{w_i (x_i - x_j)} &\text{ if } i \ne j \\ -\sum_{k \ne j} D_{kj} &\text{ if } i = j \end{cases}\]
where $w_j$ are the barycentric weights, see barycentric_weights
.
differentiation_matrix(FT, quadstyle::QuadratureStyle)
The spectral differentiation matrix at the quadrature points of quadstyle
, using floating point types FT
.
ClimaCore.Quadratures.orthonormal_poly
— FunctionV = orthonormal_poly(points, quad)
V_{ij}
contains the j-1
th Legendre polynomial evaluated at points[i]
. i.e. it is the mapping from the modal to the nodal representation.