PetscDSSetBdJacobianPreconditioner#
Set the pointwise boundary Jacobian preconditioner function for given test and basis field that constructs the matrix used to construct the preconditioner
Synopsis#
#include "petscds.h"
PetscErrorCode PetscDSSetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, PetscBdPointJacFn *g0, PetscBdPointJacFn *g1, PetscBdPointJacFn *g2, PetscBdPointJacFn *g3)
Not Collective; No Fortran Support
Input Parameters#
ds - The
PetscDS
f - The test field number
g - The field number
g0 - integrand for the test and basis function term, see
PetscBdPointJacFn
g1 - integrand for the test function and basis function gradient term, see
PetscBdPointJacFn
g2 - integrand for the test function gradient and basis function term, see
PetscBdPointJacFn
g3 - integrand for the test function gradient and basis function gradient term, see
PetscBdPointJacFn
Note#
We are using a first order FEM model for the weak form:
Developer Note#
The name is confusing since the function computes a matrix used to construct the preconditioner, not a preconditioner.
See Also#
PetscDS
, PetscBdPointJacFn
, PetscDSGetBdJacobianPreconditioner()
Level#
intermediate
Location#
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages