PetscDSSetJacobianPreconditioner#

Set the pointwise Jacobian function for given test and basis fields that constructs the matrix used to compute the preconditioner. If this is missing, the system matrix is used to build the preconditioner.

Synopsis#

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 \]

Developer Note#

The name is confusing since the function computes a matrix used to construct the preconditioner, not a preconditioner.

See Also#

PetscDS, PetscDSGetJacobianPreconditioner(), PetscDSSetJacobian(), PetscPointJacFn

Level#

intermediate

Location#

src/dm/dt/interface/dtds.c

Examples#

src/snes/tutorials/ex64.c
src/snes/tutorials/ex69.c
src/snes/tutorials/ex62.c


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