#include "petscfe.h" PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)Not collective
fem = The PetscFE object for the field being integrated | - . prob - The PetscDS specifing the discretizations and continuum functions | |
fieldI | - The test field being integrated | |
fieldJ | - The basis field being integrated | |
Ne | - The number of elements in the chunk | |
geom | - 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 specifing the auxiliary discretizations | |
coefficientsAux | - The array of FEM auxiliary basis coefficients for the elements |
Output Parameter
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)*/ PetscErrorCode PetscFEIntegrateJacobian(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFECellGeom *geom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar elemMat[]) { PetscErrorCode ierr;
PetscFunctionBegin; PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); if (fem->ops->integratejacobian) {ierr = (*fem->ops->integratejacobian)(fem, prob, fieldI, fieldJ, Ne, geom, coefficients, coefficients_t, probAux, coefficientsAux, elemMat);CHKERRQ(ierr);} PetscFunctionReturn(0); }
#undef __FUNCT__ #define __FUNCT__ "PetscFEIntegrateBdJacobian" /*C PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
Not collective
fem = The PetscFE object for the field being integrated | - . prob - The PetscDS specifing the discretizations and continuum functions | |
fieldI | - The test field being integrated | |
fieldJ | - The basis field being integrated | |
Ne | - The number of elements in the chunk | |
geom | - 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 specifing the auxiliary discretizations | |
coefficientsAux | - The array of FEM auxiliary basis coefficients for the elements |
Output Parameter
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)*/ PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFECellGeom *geom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar elemMat[]) { PetscErrorCode ierr;
PetscFunctionBegin; PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); if (fem->ops->integratebdjacobian) {ierr = (*fem->ops->integratebdjacobian)(fem, prob, fieldI, fieldJ, Ne, geom, coefficients, coefficients_t, probAux, coefficientsAux, elemMat);CHKERRQ(ierr);} PetscFunctionReturn(0); }
#undef __FUNCT__ #define __FUNCT__ "PetscFERefine" /*@ PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more sparsity). It is also used to create an interpolation between regularly refined meshes.
Level:developer
Location:src/dm/dt/interface/dtfe.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages