petsc-3.14.6 2021-03-30
PetscGaussLobattoLegendreElementAdvectionCreate
computes the advection operator for a single 1d GLL element
Synopsis
#include "petscdt.h"
PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
Not Collective
Input Parameter
| n | - the number of GLL nodes
|
| nodes | - the GLL nodes
|
| weights | - the GLL weightss
|
Output Parameter
| AA | - the stiffness element
|
Notes
Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy()
This is the same as the Gradient operator multiplied by the diagonal mass matrix
You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
See Also
PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy()
Level
beginner
Location
src/dm/dt/interface/dt.c
Examples
src/ts/tutorials/ex50.c.html
src/tao/unconstrained/tutorials/spectraladjointassimilation.c.html
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages