petsc-3.5.4 2015-05-23
Report Typos and Errors

DMPlexTSSetRHSFunctionLocal

set a local residual evaluation function

Synopsis

#include "petscdmplex.h" 
#include "petscts.h" 
PetscErrorCode DMPlexTSSetRHSFunctionLocal(DM dm, void (*riemann)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx), void *ctx)
Logically Collective

Input Arguments

dm - DM to associate callback with
riemann - Riemann solver
ctx - optional context for Riemann solve

Calling sequence for riemann

riemann(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)

x - The coordinates at a point on the interface
n - The normal vector to the interface
uL - The state vector to the left of the interface
uR - The state vector to the right of the interface
flux - output array of flux through the interface
ctx - optional user context

See Also

DMTSSetRHSFunctionLocal()

Level:beginner
Location:
src/ts/utils/dmplexts.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/ts/examples/tutorials/ex11.c.html