Actual source code: petscgll.h

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: #ifndef __PETSCGLL_H
  3:  #include <petscsys.h>

  5: /*S
  6:     PetscGLL  - the locations, from [-1,1] and weights of the Gauss-Lobatto-Legendre nodes of a given size

  8:   Level: beginner

 10:   These values are usful in implementing spectral methods based on the Gauss-Lobatto-Legendre nodes

 12:   The array nodes[] contains the vertices of each node
 13:   The array weights[] are the integration weights

 15:   The mass matrix for the element corresponds to the diagona matrix whose entries are the weights[]

 17:   Developer Notes: This may eventually get merged into a more abstract or general object for managing 
 18:     integration schemes or discretization schemes.

 20:   References: XXXX

 22: .seealso: PetscGLLCreate(), PetscGLLDestroy(), PetscGLLView()
 23: S*/
 24: typedef struct {
 25:   PetscInt  n;
 26:   PetscReal *nodes;
 27:   PetscReal *weights;
 28: } PetscGLL;

 30: /*E
 31:   PetscGLLCreateType - algorithm used to compute the GLL nodes and weights

 33:   Level: beginner

 35: $  PETSCGLL_VIA_LINEARALGEBRA - compute the nodes via linear algebra
 36: $  PETSCGLL_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method

 38: .seealso: PetscGLL, PetscGLLCreate()
 39: E*/
 40: typedef enum {PETSCGLL_VIA_LINEARALGEBRA,PETSCGLL_VIA_NEWTON} PetscGLLCreateType;

 42: #endif
 43: PETSC_EXTERN PetscErrorCode PetscGLLCreate(PetscInt,PetscGLLCreateType,PetscGLL*);
 44: PETSC_EXTERN PetscErrorCode PetscGLLDestroy(PetscGLL*);
 45: PETSC_EXTERN PetscErrorCode PetscGLLView(PetscGLL*,PetscViewer);
 46: PETSC_EXTERN PetscErrorCode PetscGLLIntegrate(PetscGLL*,const PetscReal*,PetscReal*);
 47: PETSC_EXTERN PetscErrorCode PetscGLLElementLaplacianCreate(PetscGLL*,PetscReal***);
 48: PETSC_EXTERN PetscErrorCode PetscGLLElementLaplacianDestroy(PetscGLL*,PetscReal***);
 49: PETSC_EXTERN PetscErrorCode PetscGLLElementGradientCreate(PetscGLL*,PetscReal***,PetscReal***);
 50: PETSC_EXTERN PetscErrorCode PetscGLLElementGradientDestroy(PetscGLL*,PetscReal***,PetscReal***);
 51: PETSC_EXTERN PetscErrorCode PetscGLLElementAdvectionCreate(PetscGLL*,PetscReal***);
 52: PETSC_EXTERN PetscErrorCode PetscGLLElementAdvectionDestroy(PetscGLL*,PetscReal***);