petsc-3.7.3 2016-08-01
Report Typos and Errors

PetscDualSpaceApply

Apply a functional from the dual space basis to an input function

Synopsis

#include "petscfe.h" 
PetscErrorCode PetscDualSpaceApply(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFECellGeom *geom, PetscInt numComp, 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
geom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
numComp - The number of components for the function
func - The input function
ctx - A context for the function

Output Parameter

value -numComp output values

Note: The calling sequence for the callback func is given by

func(PetscInt dim, PetscReal time, const PetscReal x[],
     PetscInt numComponents, PetscScalar values[], void *ctx)

See Also

PetscDualSpaceCreate()

Level:developer
Location:
src/dm/dt/interface/dtfe.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages