PetscFVIntegrateRHSFunction#
Produce the cell residual vector for a chunk of elements by quadrature integration
Synopsis#
#include "petscfv.h" 
PetscErrorCode PetscFVIntegrateRHSFunction(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[])
Not Collective
Input Parameters#
- fvm - The - PetscFVobject for the field being integrated
- prob - The - PetscDSspecifying the discretizations and continuum functions
- field - The field being integrated 
- Nf - The number of faces in the chunk 
- fgeom - The face geometry for each face in the chunk 
- neighborVol - The volume for each pair of cells in the chunk 
- uL - The state from the cell on the left 
- uR - The state from the cell on the right 
Output Parameters#
- fluxL - the left fluxes for each face 
- fluxR - the right fluxes for each face 
See Also#
Level#
developer
Location#
Index of all FV routines
Table of Contents for all manual pages
Index of all manual pages