petsc-3.10.5 2019-03-28
PetscGLLCreate
creates a set of the locations and weights of the Gauss-Lobatto-Legendre (GLL) nodes of a given size on the domain [-1,1]
Synopsis
PetscErrorCode PetscGLLCreate(PetscInt n,PetscGLLCreateType type,PetscGLL *gll)
Not Collective
Input Parameter
Output Parameter
gll -the nodes
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 http://epubs.siam.org/doi/abs/10.1137/110855442 http://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes
See Also
PetscGLL, PetscGLLDestroy(), PetscGLLView(), PetscGLLIntegrate(), PetscGLLElementLaplacianCreate(), PetscGLLElementLaplacianDestroy(),
PetscGLLElementGradientCreate(), PetscGLLElementGradientDestroy(), PetscGLLElementAdvectionCreate(), PetscGLLElementAdvectionDestroy()
Level
beginner
Location
src/sys/classes/gll/petscgll.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 Sys routines
Table of Contents for all manual pages
Index of all manual pages