#include "petscfe.h" PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])Not collective
ds | - The PetscDS specifying the discretizations and continuum functions | |
jtype | - The type of matrix pointwise functions that should be used | |
key | - The (label+value, fieldI*Nf + fieldJ) being integrated | |
s | - The side of the cell being integrated, 0 for negative and 1 for positive | |
Ne | - The number of elements in the chunk | |
fgeom | - The face geometry for each cell in the chunk | |
coefficients | - The array of FEM basis coefficients for the elements for the Jacobian evaluation point | |
coefficients_t | - The array of FEM basis time derivative coefficients for the elements | |
probAux | - The PetscDS specifying the auxiliary discretizations | |
coefficientsAux | - The array of FEM auxiliary basis coefficients for the elements | |
t | - The time | |
u_tShift | - A multiplier for the dF/du_t term (as opposed to the dF/du term) |
Output Parameter
elemMat | - the element matrices for the Jacobian from each element |
Loop over batch of elements (e):
Loop over element matrix entries (f,fc,g,gc --> i,j):
Loop over quadrature points (q):
Make u_q and gradU_q (loops over fields,Nb,Ncomp)
elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
+ \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
+ \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
+ \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)