petsc-3.10.5 2019-03-28
Report Typos and Errors

Compute_Lagrange_Basis_2D_Internal

Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.

Synopsis

#include "petscdt.h" 
#include "petscdmmoab.h"   
PetscErrorCode Compute_Lagrange_Basis_2D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
    PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy,
    PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
The routine is given the coordinates of the vertices of a quadrangle or triangle.

The routine evaluates the basis functions associated with each quadrature point provided, and their derivatives with respect to X and Y.

Notes

Example Physical Element (QUAD4)

4------3 s | | | | | | | | | 1------2 0-------r

Input Parameter

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 Parameter

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

Keywords

DMMoab, FEM, 2-D

Level

advanced

Location

src/dm/impls/moab/dmmbfem.cxx
Index of all DMMOAB routines
Table of Contents for all manual pages
Index of all manual pages