PetscFEIntegrate#

Produce the integral for the given field for a chunk of elements by quadrature integration

Synopsis#

#include "petscfe.h" 
PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])

Not Collective

Input Parameters#

  • prob - The PetscDS specifying the discretizations and continuum functions

  • field - The field being integrated

  • Ne - The number of elements in the chunk

  • cgeom - The cell geometry for each cell in the chunk

  • coefficients - The array of FEM basis coefficients for the elements

  • probAux - The PetscDS specifying the auxiliary discretizations

  • coefficientsAux - The array of FEM auxiliary basis coefficients for the elements

Output Parameter#

  • integral - the integral for this field

Developer Notes#

The function name begins with PetscFE and yet the first argument is PetscDS and it has no PetscFE arguments.

See Also#

PetscFE, PetscDS, PetscFEIntegrateResidual(), PetscFEIntegrateBd()

Level#

intermediate

Location#

src/dm/dt/fe/interface/fe.c

Implementations#

PetscFEIntegrate_Basic() in src/dm/dt/fe/impls/basic/febasic.c


Index of all FE routines
Table of Contents for all manual pages
Index of all manual pages