PetscDSGetObjective#

Get the pointwise objective function for a given test field that was provided with PetscDSSetObjective()

Synopsis#

#include "petscds.h" 
PetscErrorCode PetscDSGetObjective(PetscDS ds, PetscInt f, PetscPointFn **obj)

Not Collective

Input Parameters#

  • ds - The PetscDS

  • f - The test field number

Output Parameter#

  • obj - integrand for the test function term, see PetscPointFn

Note#

We are using a first order FEM model for the weak form: \( \int_\Omega \phi\,\mathrm{obj}(u, u_t, \nabla u, x, t)\)

See Also#

PetscPointFn, PetscDS, PetscDSSetObjective(), PetscDSGetResidual()

Level#

intermediate

Location#

src/dm/dt/interface/dtds.c


Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages