petsc-3.10.5 2019-03-28
Report Typos and Errors

PetscDualSpaceApplyFVM

Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.

Synopsis

#include "petscfe.h" 
PetscErrorCode PetscDualSpaceApplyFVM(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFVCellGeom *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 currently just use the centroid
Nc - The number of components for the function
func - The input function
ctx - A context for the function

Output Parameter

value -The output value (scalar)

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)

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

PetscDualSpaceCreate()

Level

developer

Location

src/dm/dt/dualspace/interface/dualspace.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages