PetscDSSetBdJacobianPreconditioner#
Set the pointwise boundary Jacobian preconditioner function for given test and basis field
Synopsis#
#include "petscds.h"
PetscErrorCode PetscDSSetBdJacobianPreconditioner(PetscDS ds, PetscInt f, PetscInt g, void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
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
g1 - integrand for the test function and basis function gradient term
g2 - integrand for the test function gradient and basis function term
g3 - integrand for the test function gradient and basis function gradient term
Calling sequence of `g0’#
dim - the spatial dimension
Nf - the number of fields
NfAux - the number of auxiliary fields
uOff - the offset into u[] and u_t[] for each field
uOff_x - the offset into u_x[] for each field
u - each field evaluated at the current point
u_t - the time derivative of each field evaluated at the current point
u_x - the gradient of each field evaluated at the current point
aOff - the offset into a[] and a_t[] for each auxiliary field
aOff_x - the offset into a_x[] for each auxiliary field
a - each auxiliary field evaluated at the current point
a_t - the time derivative of each auxiliary field evaluated at the current point
a_x - the gradient of auxiliary each field evaluated at the current point
t - current time
u_tShift - the multiplier a for dF/dU_t
x - coordinates of the current point
n - normal at the current point
numConstants - number of constant parameters
constants - constant parameters
g0 - output values at the current point
Note#
g1
, g2
, and g3
have identical calling sequences to g0
and are omitted for brevity.
We are using a first order FEM model for the weak form: \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
See Also#
Level#
intermediate
Location#
Index of all DT routines
Table of Contents for all manual pages
Index of all manual pages