DMProjectFunctionLabelLocal#
This projects the given function into the function space provided by the DM
, putting the coefficients in a local vector, setting values only for points in the given label.
Synopsis#
#include "petscdm.h"
#include "petscdmlabel.h"
#include "petscds.h"
PetscErrorCode DMProjectFunctionLabelLocal(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
Not Collective
Input Parameters#
Output Parameter#
localX - vector
Calling sequence of funcs
#
PetscErrorCode funcs(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
dim - The spatial dimension
x - The coordinates
Nc - The number of components
u - The output field values
ctx - optional user-defined function context
Developer Notes#
This API is specific to only particular usage of DM
The notes need to provide some information about what has to be provided to the DM
to be able to perform the computation.
See Also#
DM Basics, DM
, DMProjectFunction()
, DMProjectFunctionLocal()
, DMProjectFunctionLabel()
, DMComputeL2Diff()
Level#
developer
Location#
Implementations#
DMProjectFunctionLabelLocal_pforest(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], PetscErrorCode (**funcs) in src/dm/impls/forest/p4est/pforest.h
DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], 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