petsc-3.12.5 2020-03-29
PetscDTGaussLobattoLegendreQuadrature
creates a set of the locations and weights of the Gauss-Lobatto-Legendre nodes of a given size on the domain [-1,1]
Synopsis
#include "petscdt.h"
PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w)
Not Collective
Input Parameter
Output Arguments
| x | - quadrature points
|
| w | - quadrature weights
|
Notes
For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
close enough to the desired solution
These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
See https://epubs.siam.org/doi/abs/10.1137/110855442 https://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes
See Also
PetscDTGaussQuadrature()
Level
intermediate
Location
src/dm/dt/interface/dt.c
Examples
src/ksp/ksp/examples/tutorials/ex68.c.html
src/ksp/ksp/examples/tutorials/ex69.c.html
src/ts/examples/tutorials/ex50.c.html
src/tao/unconstrained/examples/tutorials/spectraladjointassimilation.c.html
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages