DMPlexComputeJacobianByKey#

Compute the local Jacobian for terms matching the input key

Synopsis#

#include "petscdmplex.h"   
PetscErrorCode DMPlexComputeJacobianByKey(DM dm, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Mat Jac, Mat JacP, void *user)

Collective

Input Parameters#

  • dm - The output DM

  • key - The PetscFormKey indicating what should be integrated

  • cellIS - The IS give a set of cells to integrate over

  • t - The time

  • X_tShift - The multiplier for the Jacobian with respect to \(X_t\)

  • locX - The local solution

  • locX_t - The time derivative of the local solution, or NULL for time-independent problems

  • user - An optional user context, passed to the pointwise functions

Output Parameters#

  • Jac - The local Jacobian

  • JacP - The local Jacobian preconditioner

See Also#

DMPlexComputeResidualByKey(), DMPlexComputeResidualHybridByKey(), DMPlexComputeJacobianHybridByKey(), PetscFormKey

Level#

developer

Location#

src/dm/impls/plex/plexfem.c


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