chebyshev
moment_kinetics.chebyshev
— Modulemoment_kinetics.chebyshev.chebyshev_base_info
— TypeChebyshev pseudospectral discretization
moment_kinetics.calculus.elementwise_derivative!
— Methodelementwise_derivative!(coord, ff, adv_fac, spectral::chebyshev_info)
Chebyshev transform f to get Chebyshev spectral coefficients and use them to calculate f'.
Note: Chebyshev derivative does not make use of upwinding information within each element.
moment_kinetics.calculus.elementwise_derivative!
— Methodelementwise_derivative!(coord, ff, chebyshev::chebyshev_info)
Chebyshev transform f to get Chebyshev spectral coefficients and use them to calculate f'.
moment_kinetics.chebyshev.cheb_derivative_matrix_elementwise!
— Methodderivative matrix for Gauss-Lobatto points using the analytical specification from Chapter 8.2 from Trefethen 1994 https://people.maths.ox.ac.uk/trefethen/8all.pdf full list of Chapters may be obtained here https://people.maths.ox.ac.uk/trefethen/pdetext.html
moment_kinetics.chebyshev.cheb_derivative_matrix_elementwise_radau_by_FFT!
— MethodDerivative matrix for Chebyshev-Radau grid using the FFT. Note that a similar function could be constructed for the Chebyshev-Lobatto grid, if desired.
moment_kinetics.chebyshev.chebyshev_backward_transform!
— Methodmoment_kinetics.chebyshev.chebyshev_derivative_single_element!
— Methodmoment_kinetics.chebyshev.chebyshev_forward_transform!
— Methodtakes the real function ff on a Chebyshev grid in z (domain [-1, 1]), which corresponds to the domain [π, 2π] in variable theta = ArcCos(z). interested in functions of form f(z) = sumn cn Tn(z) using Tn(cos(theta)) = cos(ntheta) and z = cos(theta) gives f(z) = sumn cn cos(ntheta) thus a Chebyshev transform is equivalent to a discrete cosine transform doing this directly turns out to be slower than extending the domain from [0, 2pi] and using the fact that f(z) must be even (as cosines are all even) on this extended domain, can do a standard complex-to-complex fft fext is an array used to store f(theta) on the extended grid theta ∈ [0,2π) ff is f(theta) on the grid [π,2π] the Chebyshev coefficients of ff are calculated and stored in chebyf n is the number of grid points on the Chebyshev-Gauss-Lobatto grid transform is the plan for the complex-to-complex, in-place fft
moment_kinetics.chebyshev.chebyshev_radau_backward_transform!
— Methodmoment_kinetics.chebyshev.chebyshev_spectral_derivative!
— Methoduse Chebyshev basis to compute the first derivative of f
moment_kinetics.chebyshev.chebyshevmoments
— Methodcompute and return modified Chebyshev moments of the first kind: ∫dx Tᵢ(x) over range [-1,1]
moment_kinetics.chebyshev.chebyshevpoints
— Methodreturns the Chebyshev-Gauss-Lobatto grid points on an n point grid
moment_kinetics.chebyshev.clenshaw_curtis_weights
— Methodreturns wgts array containing the integration weights associated with all grid points for Clenshaw-Curtis quadrature
moment_kinetics.chebyshev.scaled_chebyshev_grid
— Methodinitialize chebyshev grid scaled to interval [-boxlength/2, boxlength/2] we no longer pass the boxlength to this function, but instead pass precomputed arrays elementscale and element_shift that are needed to compute the grid.
ngrid – number of points per element (including boundary points) nelementlocal – number of elements in the local (distributed memory MPI) grid n – total number of points in the local grid (excluding duplicate points) elementscale – the scale factor in the transform from the coordinates where the element limits are -1, 1 to the coordinate where the limits are Aj = coord.grid[imin[j]-1] and Bj = coord.grid[imax[j]] elementscale = 0.5*(Bj - Aj) elementshift – the centre of the element in the extended grid coordinate element_shift = 0.5*(Aj + Bj) imin – the array of minimum indices of each element on the extended grid. By convention, the duplicated points are not included, so for element index j > 1 the lower boundary point is actually imin[j] - 1 imax – the array of maximum indices of each element on the extended grid.
moment_kinetics.chebyshev.setup_chebyshev_pseudospectral
— Methodcreate arrays needed for explicit Chebyshev pseudospectral treatment and create the plans for the forward and backward fast Fourier transforms
moment_kinetics.chebyshev.update_df_chebyshev!
— Methodcompute the Chebyshev spectral coefficients of the spatial derivative of f
moment_kinetics.chebyshev.update_fcheby!
— MethodChebyshev transform f to get Chebyshev spectral coefficients