petsc-3.7.5 2017-01-01
Report Typos and Errors

PetscDSGetJacobian

Get the pointwise Jacobian function for given test and basis field

Synopsis

#include "petscds.h" 
PetscErrorCode PetscDSGetJacobian(PetscDS prob, 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[], PetscScalar g0[]),
                                  void (**g1)(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[], PetscScalar g1[]),
                                  void (**g2)(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[], PetscScalar g2[]),
                                  void (**g3)(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[], PetscScalar g3[]))
Not collective

Input Parameters

prob - The PetscDS
f - The test field number
g - The field number

Output Parameters

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

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

The calling sequence for the callbacks g0, g1, g2 and g3 is given by

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, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])

dim - the spatial dimension
Nf - the number of 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
g0 - output values at the current point

See Also

PetscDSSetJacobian()

Level:intermediate
Location:
src/dm/dt/interface/dtds.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages