1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMMOAB/DMMoabFEMComputeBasis.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>DMMoabFEMComputeBasis</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.14.5 2021-03-03</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.14.5 v3.14.5 docs/manualpages/DMMOAB/DMMoabFEMComputeBasis.html "><small>Report Typos and Errors</small></a></div>
<A NAME="DMMoabFEMComputeBasis"><H1>DMMoabFEMComputeBasis</H1></A>
Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
<H3><FONT COLOR="#CC3333">Synopsis</FONT></H3>
<PRE>
#include "petscdt.h"
#include "petscdmmoab.h"
<A HREF="../Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../DMMOAB/DMMoabFEMComputeBasis.html#DMMoabFEMComputeBasis">DMMoabFEMComputeBasis</A>(const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> dim, const <A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> nverts, const <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> *coordinates, const <A HREF="../FE/PetscQuadrature.html#PetscQuadrature">PetscQuadrature</A> quadrature,
<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> *phypts, <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> *jacobian_quadrature_weight_product,
<A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> *fe_basis, <A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> **fe_basis_derivatives)
</PRE>
The routine takes the coordinates of the vertices of an element and computes basis functions associated with
each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
<P>
<H3><FONT COLOR="#CC3333">Input Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> nverts, the number of element vertices</B></TD>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> coords[3*nverts], the physical coordinates of the vertices (in canonical numbering)</B></TD>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscInt.html#PetscInt">PetscInt</A> npts, the number of evaluation points (quadrature points)</B></TD>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> quad[3*npts], the evaluation points (quadrature points) in the reference space</B></TD></TABLE>
<P>
<H3><FONT COLOR="#CC3333">Output Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> phypts[3*npts], the evaluation points (quadrature points) transformed to the physical space</B></TD>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> jxw[npts], the jacobian determinant * quadrature weight necessary for assembling discrete contributions</B></TD>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> fe_basis[npts], the bases values evaluated at the specified quadrature points</B></TD>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B><A HREF="../Sys/PetscReal.html#PetscReal">PetscReal</A> fe_basis_derivatives[dim][npts], the derivative of the bases wrt (X,Y,Z)</B></TD><TD>- directions (depending on the dimension) evaluated at the specified quadrature points
</TD></TR></TABLE>
<P>
<P>
<P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>advanced<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/dm/impls/moab/dmmbfem.cxx.html#DMMoabFEMComputeBasis">src/dm/impls/moab/dmmbfem.cxx</A>
<BR><A HREF="./index.html">Index of all DMMOAB routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|