petsc-3.14.6 2021-03-30
DMProjectFunctionLocal
This projects the given function into the function space provided, putting the coefficients in a local vector.
Synopsis
#include "petscdm.h"
#include "petscdmlabel.h"
#include "petscds.h"
PetscErrorCode DMProjectFunctionLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
Not collective
Input Parameters
| dm | - The DM
|
| time | - The time
|
| funcs | - The coordinate functions to evaluate, one per field
|
| ctxs | - Optional array of contexts to pass to each coordinate function. ctxs itself may be null.
|
| mode | - The insertion mode for values
|
Output Parameter
Calling sequence of func
func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
| dim | - The spatial dimension
|
| x | - The coordinates
|
| Nf | - The number of fields
|
| u | - The output field values
|
| ctx | - optional user-defined function context
|
See Also
DMProjectFunction(), DMProjectFunctionLabel(), DMComputeL2Diff()
Level
developer
Location
src/dm/interface/dm.c
Examples
src/dm/field/tutorials/ex1.c.html
src/snes/tutorials/ex12.c.html
src/snes/tutorials/ex77.c.html
src/ts/tutorials/ex48.c.html
Implementations
DMProjectFunctionLocal_pforest(DM dm, PetscReal time, PetscErrorCode (**funcs) in src/dm/impls/forest/p4est/pforest.c
DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs) in src/dm/impls/plex/plexproject.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages