PetscDualSpaceApplyDefault#
Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
Synopsis#
#include "petscfe.h"
PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
Input Parameters#
sp - The PetscDualSpace object
f - The basis functional index
time - The time
cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
Nc - The number of components for the function
func - The input function
ctx - A context for the function
Output Parameter#
value - The output value
Note: The calling sequence for the callback func is given by#
PetscInt numComponents, PetscScalar values[], void *ctx)
and the idea is to evaluate the functional as an integral
n(f) = int dx n(x) . f(x)
where both n and f have Nc components.
See Also#
Level#
beginner
Location#
src/dm/dt/dualspace/interface/dualspace.c
Index of all DUALSPACE routines
Table of Contents for all manual pages
Index of all manual pages