DMPlexComputeBdJacobianSingleByLabel#

Compute the local boundary Jacobian for terms matching the input label

Synopsis#

#include "petscdmplex.h"   
PetscErrorCode DMPlexComputeBdJacobianSingleByLabel(DM dm, PetscWeakForm wf, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, IS facetIS, Vec locX, Vec locX_t, PetscReal t, DMField coordField, PetscReal X_tShift, Mat Jac, Mat JacP)

Not collective

Input Parameters#

  • dm - The output DM

  • wf - The PetscWeakForm holding forms on this boundary

  • label - The DMLabel indicating what faces should be integrated over

  • numValues - The number of label values

  • values - The array of label values

  • fieldI - The test field for these integrals

  • facetIS - The IS giving the set of possible faces to integrate over (intersected with the label)

  • locX - The local solution

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

  • t - The time

  • coordField - The DMField object with coordinates for these faces

  • X_tShift - The multiplier for dF/dxdot

Output Parameters#

  • Jac - The local Jacobian

  • JacP - The local Jacobian preconditioner

See Also#

DMPlexComputeBdJacobianSingle(), DMPlexComputeJacobianByKey(), 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