#include "petscds.h" PetscErrorCode PetscDSGetDynamicJacobian(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 Many br
prob | - The PetscDS Many br | |
f | - The test field number Many br | |
g | - The field number Many br |
g0 | - integrand for the test and basis function term Many br | |
g1 | - integrand for the test function and basis function gradient term Many br | |
g2 | - integrand for the test function gradient and basis function term Many br | |
g3 | - integrand for the test function gradient and basis function gradient term Many br |
\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 Many br
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 Many br | |
Nf | - the number of fields Many br | |
uOff | - the offset into u[] and u_t[] for each field Many br | |
uOff_x | - the offset into u_x[] for each field Many br | |
u | - each field evaluated at the current point Many br | |
u_t | - the time derivative of each field evaluated at the current point Many br | |
u_x | - the gradient of each field evaluated at the current point Many br | |
aOff | - the offset into a[] and a_t[] for each auxiliary field Many br | |
aOff_x | - the offset into a_x[] for each auxiliary field Many br | |
a | - each auxiliary field evaluated at the current point Many br | |
a_t | - the time derivative of each auxiliary field evaluated at the current point Many br | |
a_x | - the gradient of auxiliary each field evaluated at the current point Many br | |
t | - current time Many br | |
u_tShift | - the multiplier a for dF/dU_t Many br | |
x | - coordinates of the current point Many br | |
g0 | - output values at the current point Many br |
Many br
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