petsc-3.8.4 2018-03-24
Report Typos and Errors

PetscDSGetExactSolution

Get the pointwise exact solution function for a given test field

Synopsis

#include "petscds.h" 
PetscErrorCode PetscDSGetExactSolution(PetscDS prob, PetscInt f, PetscErrorCode (**sol)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx))
Not collective

Input Parameters

prob - The PetscDS
f - The test field number

Output Parameter

exactSol -exact solution for the test field

Note: The calling sequence for the solution functions is given by

sol(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx)

dim - the spatial dimension
t - current time
x - coordinates of the current point
Nc - the number of field components
u - the solution field evaluated at the current point
ctx - a user context

See Also

PetscDSSetExactSolution()

Level:intermediate
Location:
src/dm/dt/interface/dtds.c
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages