PetscDSSetDynamicJacobian#
Set the pointwise dynamic Jacobian, \(dF/du_t\), function for given test and basis fields
Synopsis#
#include "petscds.h"
PetscErrorCode PetscDSSetDynamicJacobian(PetscDS ds, PetscInt f, PetscInt g, PetscPointJacFn *g0, PetscPointJacFn *g1, PetscPointJacFn *g2, PetscPointJacFn *g3)
Not Collective
Input Parameters#
ds - The
PetscDS
f - The test field number
g - The field number
g0 - integrand for the test and basis function term, see
PetscPointJacFn
g1 - integrand for the test function and basis function gradient term, see
PetscPointJacFn
g2 - integrand for the test function gradient and basis function term, see
PetscPointJacFn
g3 - integrand for the test function gradient and basis function gradient term, see
PetscPointJacFn
Note#
We are using a first order FEM model for the weak form:
\[
\int_\Omega \phi\, g_0(u, u_t, \nabla u, x, t) \psi + \phi\, {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi
+ \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
\]
See Also#
PetscDS
, PetscDSGetDynamicJacobian()
, PetscDSGetJacobian()
, PetscPointJacFn
Level#
intermediate
Location#
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages