PetscSimplePointFn#

A prototype of a simple pointwise function that can be passed to, for example, DMPlexTransformExtrudeSetNormalFunction()

Synopsis#

PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscSimplePointFn(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx);

Calling Sequence#

  • dim - The coordinate dimension of the original mesh (usually a surface)

  • time - The current time, or 0.

  • x - The location of the current normal, in the coordinate space of the original mesh

  • r - The layer number of this point

  • u - The user provides the computed normal on output

  • ctx - An optional user context, this context may be obtained by the calling code with DMGetApplicationContext()

Developer Note#

The handling of ctx in the use of such functions may not be ideal since the context is not provided when the function pointer is provided with, for example, DMSwarmSetCoordinateFunction()

See Also#

PetscPointFn, DMPlexTransformExtrudeSetNormalFunction(), DMSwarmSetCoordinateFunction()

Level#

beginner

Location#

include/petscdstypes.h

Examples#

src/snes/tutorials/ex36.c
src/snes/tutorials/ex71.c
src/ts/tutorials/ex76.c
src/ts/tutorials/ex47.c
src/ts/tutorials/ex46.c
src/snes/tutorials/ex27.c


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