| fem | - The PetscFE object for the field being integrated
|
| prob | - The PetscDS specifying the discretizations and continuum functions
|
| jtype | - The type of matrix pointwise functions that should be used
|
| fieldI | - The test field being integrated
|
| fieldJ | - The basis field being integrated
|
| Ne | - The number of elements in the chunk
|
| cgeom | - The cell 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)
|
elemMat -the element matrices for the Jacobian from each element
Note
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)
See Also
PetscFEIntegrateResidual()
Level
intermediate
Location
src/dm/dt/fe/interface/fe.c
Implementations
src/dm/dt/fe/impls/basic/febasic.c:574:PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
Index of all FE routines
Table of Contents for all manual pages
Index of all manual pages