:orphan: # Compute_Lagrange_Basis_3D_Internal Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element. ## Synopsis ``` #include "petscdt.h" #include "petscdmmoab.h" PetscErrorCode Compute_Lagrange_Basis_3D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz, PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) ``` The routine is given the coordinates of the vertices of a hexahedra or tetrahedra. The routine evaluates the basis functions associated with each quadrature point provided, and their derivatives with respect to X, Y and Z. ## Notes Example Physical Element (HEX8) ```none 8------7 /| /| t s 5------6 | | / | | | | |/ | 4----|-3 0-------r |/ |/ 1------2 ``` ## Input Parameters - ***PetscInt nverts -*** the number of element vertices - ***PetscReal coords[3*nverts] -*** the physical coordinates of the vertices (in canonical numbering) - ***PetscInt npts -*** the number of evaluation points (quadrature points) - ***PetscReal quad[3*npts] -*** the evaluation points (quadrature points) in the reference space ## Output Parameters - ***PetscReal phypts[3*npts] -*** the evaluation points (quadrature points) transformed to the physical space - ***PetscReal jxw[npts] -*** the jacobian determinant * quadrature weight necessary for assembling discrete contributions - ***PetscReal phi[npts] -*** the bases evaluated at the specified quadrature points - ***PetscReal dphidx[npts] -*** the derivative of the bases wrt X-direction evaluated at the specified quadrature points - ***PetscReal dphidy[npts] -*** the derivative of the bases wrt Y-direction evaluated at the specified quadrature points - ***PetscReal dphidz[npts] -*** the derivative of the bases wrt Z-direction evaluated at the specified quadrature points - ***jacobian -*** jacobian - ***ijacobian -*** ijacobian - ***volume -*** volume ## Level advanced ## Location src/dm/impls/moab/dmmbfem.cxx --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/dm/impls/moab/dmmbfem.cxx) [Index of all DMMOAB routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)