xref: /petsc/src/dm/impls/moab/dmmbfem.cxx (revision a2b725a8db0d6bf6cc2a1c6df7dd8029aadfff6e)
163d025dbSVijay Mahadevan 
263d025dbSVijay Mahadevan #include <petscconf.h>
363d025dbSVijay Mahadevan #include <petscdt.h>            /*I "petscdt.h" I*/
463d025dbSVijay Mahadevan #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
563d025dbSVijay Mahadevan 
663d025dbSVijay Mahadevan /* Utility functions */
763d025dbSVijay Mahadevan static inline PetscReal DMatrix_Determinant_2x2_Internal ( const PetscReal inmat[2 * 2] )
863d025dbSVijay Mahadevan {
963d025dbSVijay Mahadevan   return  inmat[0] * inmat[3] - inmat[1] * inmat[2];
1063d025dbSVijay Mahadevan }
1163d025dbSVijay Mahadevan 
1263d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_2x2_Internal (const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant)
1363d025dbSVijay Mahadevan {
1463d025dbSVijay Mahadevan   if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 2x2 inversion.");
1563d025dbSVijay Mahadevan   PetscReal det = DMatrix_Determinant_2x2_Internal(inmat);
1663d025dbSVijay Mahadevan   if (outmat) {
1763d025dbSVijay Mahadevan     outmat[0] = inmat[3] / det;
1863d025dbSVijay Mahadevan     outmat[1] = -inmat[1] / det;
1963d025dbSVijay Mahadevan     outmat[2] = -inmat[2] / det;
2063d025dbSVijay Mahadevan     outmat[3] = inmat[0] / det;
2163d025dbSVijay Mahadevan   }
2263d025dbSVijay Mahadevan   if (determinant) *determinant = det;
2363d025dbSVijay Mahadevan   PetscFunctionReturn(0);
2463d025dbSVijay Mahadevan }
2563d025dbSVijay Mahadevan 
26181f196bSVijay Mahadevan static inline PetscReal DMatrix_Determinant_3x3_Internal ( const PetscReal inmat[3 * 3] )
2763d025dbSVijay Mahadevan {
2863d025dbSVijay Mahadevan   return   inmat[0] * (inmat[8] * inmat[4] - inmat[7] * inmat[5])
2963d025dbSVijay Mahadevan            - inmat[3] * (inmat[8] * inmat[1] - inmat[7] * inmat[2])
3063d025dbSVijay Mahadevan            + inmat[6] * (inmat[5] * inmat[1] - inmat[4] * inmat[2]);
3163d025dbSVijay Mahadevan }
3263d025dbSVijay Mahadevan 
3363d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_3x3_Internal (const PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
3463d025dbSVijay Mahadevan {
3563d025dbSVijay Mahadevan   if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 3x3 inversion.");
36181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_3x3_Internal(inmat);
3763d025dbSVijay Mahadevan   if (outmat) {
3863d025dbSVijay Mahadevan     outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det;
3963d025dbSVijay Mahadevan     outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det;
4063d025dbSVijay Mahadevan     outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det;
4163d025dbSVijay Mahadevan     outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det;
4263d025dbSVijay Mahadevan     outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det;
4363d025dbSVijay Mahadevan     outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det;
4463d025dbSVijay Mahadevan     outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det;
4563d025dbSVijay Mahadevan     outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det;
4663d025dbSVijay Mahadevan     outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det;
4763d025dbSVijay Mahadevan   }
4863d025dbSVijay Mahadevan   if (determinant) *determinant = det;
4963d025dbSVijay Mahadevan   PetscFunctionReturn(0);
5063d025dbSVijay Mahadevan }
5163d025dbSVijay Mahadevan 
52181f196bSVijay Mahadevan inline PetscReal DMatrix_Determinant_4x4_Internal ( PetscReal inmat[4 * 4] )
53181f196bSVijay Mahadevan {
54181f196bSVijay Mahadevan   return
55181f196bSVijay Mahadevan     inmat[0 + 0 * 4] * (
56181f196bSVijay Mahadevan       inmat[1 + 1 * 4] * ( inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4] )
57181f196bSVijay Mahadevan       - inmat[1 + 2 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4] )
58181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4] ) )
59181f196bSVijay Mahadevan     - inmat[0 + 1 * 4] * (
60181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * ( inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4] )
61181f196bSVijay Mahadevan       - inmat[1 + 2 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4] )
62181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4] ) )
63181f196bSVijay Mahadevan     + inmat[0 + 2 * 4] * (
64181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4] )
65181f196bSVijay Mahadevan       - inmat[1 + 1 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4] )
66181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4] ) )
67181f196bSVijay Mahadevan     - inmat[0 + 3 * 4] * (
68181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4] )
69181f196bSVijay Mahadevan       - inmat[1 + 1 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4] )
70181f196bSVijay Mahadevan       + inmat[1 + 2 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4] ) );
71181f196bSVijay Mahadevan }
72181f196bSVijay Mahadevan 
73181f196bSVijay Mahadevan inline PetscErrorCode DMatrix_Invert_4x4_Internal (PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
74181f196bSVijay Mahadevan {
75181f196bSVijay Mahadevan   if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 4x4 inversion.");
76181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_4x4_Internal(inmat);
77181f196bSVijay Mahadevan   if (outmat) {
78181f196bSVijay Mahadevan     outmat[0] =  (inmat[5] * inmat[10] * inmat[15] + inmat[6] * inmat[11] * inmat[13] + inmat[7] * inmat[9] * inmat[14] - inmat[5] * inmat[11] * inmat[14] - inmat[6] * inmat[9] * inmat[15] - inmat[7] * inmat[10] * inmat[13]) / det;
79181f196bSVijay Mahadevan     outmat[1] =  (inmat[1] * inmat[11] * inmat[14] + inmat[2] * inmat[9] * inmat[15] + inmat[3] * inmat[10] * inmat[13] - inmat[1] * inmat[10] * inmat[15] - inmat[2] * inmat[11] * inmat[13] - inmat[3] * inmat[9] * inmat[14]) / det;
80181f196bSVijay Mahadevan     outmat[2] =  (inmat[1] * inmat[6] * inmat[15] + inmat[2] * inmat[7] * inmat[13] + inmat[3] * inmat[5] * inmat[14] - inmat[1] * inmat[7] * inmat[14] - inmat[2] * inmat[5] * inmat[15] - inmat[3] * inmat[6] * inmat[13]) / det;
81181f196bSVijay Mahadevan     outmat[3] =  (inmat[1] * inmat[7] * inmat[10] + inmat[2] * inmat[5] * inmat[11] + inmat[3] * inmat[6] * inmat[9] - inmat[1] * inmat[6] * inmat[11] - inmat[2] * inmat[7] * inmat[9] - inmat[3] * inmat[5] * inmat[10]) / det;
82181f196bSVijay Mahadevan     outmat[4] =  (inmat[4] * inmat[11] * inmat[14] + inmat[6] * inmat[8] * inmat[15] + inmat[7] * inmat[10] * inmat[12] - inmat[4] * inmat[10] * inmat[15] - inmat[6] * inmat[11] * inmat[12] - inmat[7] * inmat[8] * inmat[14]) / det;
83181f196bSVijay Mahadevan     outmat[5] =  (inmat[0] * inmat[10] * inmat[15] + inmat[2] * inmat[11] * inmat[12] + inmat[3] * inmat[8] * inmat[14] - inmat[0] * inmat[11] * inmat[14] - inmat[2] * inmat[8] * inmat[15] - inmat[3] * inmat[10] * inmat[12]) / det;
84181f196bSVijay Mahadevan     outmat[6] =  (inmat[0] * inmat[7] * inmat[14] + inmat[2] * inmat[4] * inmat[15] + inmat[3] * inmat[6] * inmat[12] - inmat[0] * inmat[6] * inmat[15] - inmat[2] * inmat[7] * inmat[12] - inmat[3] * inmat[4] * inmat[14]) / det;
85181f196bSVijay Mahadevan     outmat[7] =  (inmat[0] * inmat[6] * inmat[11] + inmat[2] * inmat[7] * inmat[8] + inmat[3] * inmat[4] * inmat[10] - inmat[0] * inmat[7] * inmat[10] - inmat[2] * inmat[4] * inmat[11] - inmat[3] * inmat[6] * inmat[8]) / det;
86181f196bSVijay Mahadevan     outmat[8] =  (inmat[4] * inmat[9] * inmat[15] + inmat[5] * inmat[11] * inmat[12] + inmat[7] * inmat[8] * inmat[13] - inmat[4] * inmat[11] * inmat[13] - inmat[5] * inmat[8] * inmat[15] - inmat[7] * inmat[9] * inmat[12]) / det;
87181f196bSVijay Mahadevan     outmat[9] =  (inmat[0] * inmat[11] * inmat[13] + inmat[1] * inmat[8] * inmat[15] + inmat[3] * inmat[9] * inmat[12] - inmat[0] * inmat[9] * inmat[15] - inmat[1] * inmat[11] * inmat[12] - inmat[3] * inmat[8] * inmat[13]) / det;
88181f196bSVijay Mahadevan     outmat[10] = (inmat[0] * inmat[5] * inmat[15] + inmat[1] * inmat[7] * inmat[12] + inmat[3] * inmat[4] * inmat[13] - inmat[0] * inmat[7] * inmat[13] - inmat[1] * inmat[4] * inmat[15] - inmat[3] * inmat[5] * inmat[12]) / det;
89181f196bSVijay Mahadevan     outmat[11] = (inmat[0] * inmat[7] * inmat[9] + inmat[1] * inmat[4] * inmat[11] + inmat[3] * inmat[5] * inmat[8] - inmat[0] * inmat[5] * inmat[11] - inmat[1] * inmat[7] * inmat[8] - inmat[3] * inmat[4] * inmat[9]) / det;
90181f196bSVijay Mahadevan     outmat[12] = (inmat[4] * inmat[10] * inmat[13] + inmat[5] * inmat[8] * inmat[14] + inmat[6] * inmat[9] * inmat[12] - inmat[4] * inmat[9] * inmat[14] - inmat[5] * inmat[10] * inmat[12] - inmat[6] * inmat[8] * inmat[13]) / det;
91181f196bSVijay Mahadevan     outmat[13] = (inmat[0] * inmat[9] * inmat[14] + inmat[1] * inmat[10] * inmat[12] + inmat[2] * inmat[8] * inmat[13] - inmat[0] * inmat[10] * inmat[13] - inmat[1] * inmat[8] * inmat[14] - inmat[2] * inmat[9] * inmat[12]) / det;
92181f196bSVijay Mahadevan     outmat[14] = (inmat[0] * inmat[6] * inmat[13] + inmat[1] * inmat[4] * inmat[14] + inmat[2] * inmat[5] * inmat[12] - inmat[0] * inmat[5] * inmat[14] - inmat[1] * inmat[6] * inmat[12] - inmat[2] * inmat[4] * inmat[13]) / det;
93181f196bSVijay Mahadevan     outmat[15] = (inmat[0] * inmat[5] * inmat[10] + inmat[1] * inmat[6] * inmat[8] + inmat[2] * inmat[4] * inmat[9] - inmat[0] * inmat[6] * inmat[9] - inmat[1] * inmat[4] * inmat[10] - inmat[2] * inmat[5] * inmat[8]) / det;
94181f196bSVijay Mahadevan   }
95181f196bSVijay Mahadevan   if (determinant) *determinant = det;
96181f196bSVijay Mahadevan   PetscFunctionReturn(0);
97181f196bSVijay Mahadevan }
98181f196bSVijay Mahadevan 
9963d025dbSVijay Mahadevan 
10063d025dbSVijay Mahadevan /*@
10197b73a88SSatish Balay   Compute_Lagrange_Basis_1D_Internal - Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element.
10263d025dbSVijay Mahadevan 
10363d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a linear or quadratic edge element.
10463d025dbSVijay Mahadevan 
10563d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
10663d025dbSVijay Mahadevan   and their derivatives with respect to X.
10763d025dbSVijay Mahadevan 
10863d025dbSVijay Mahadevan   Notes:
10963d025dbSVijay Mahadevan 
11063d025dbSVijay Mahadevan   Example Physical Element
111*a2b725a8SWilliam Gropp .vb
11263d025dbSVijay Mahadevan     1-------2        1----3----2
11363d025dbSVijay Mahadevan       EDGE2             EDGE3
114*a2b725a8SWilliam Gropp .ve
11563d025dbSVijay Mahadevan 
11663d025dbSVijay Mahadevan   Input Parameter:
117*a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
118*a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
119*a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
120*a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
12163d025dbSVijay Mahadevan 
12263d025dbSVijay Mahadevan   Output Parameter:
123*a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
124*a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
125*a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
126*a2b725a8SWilliam Gropp -  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
12763d025dbSVijay Mahadevan 
128edc382c3SSatish Balay   Level: advanced
129edc382c3SSatish Balay 
13063d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 1-D
13163d025dbSVijay Mahadevan @*/
13263d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_1D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
133181f196bSVijay Mahadevan     PetscReal *jxw, PetscReal *phi, PetscReal *dphidx,
134181f196bSVijay Mahadevan     PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
13563d025dbSVijay Mahadevan {
13663d025dbSVijay Mahadevan   int i, j;
13763d025dbSVijay Mahadevan   PetscErrorCode  ierr;
13863d025dbSVijay Mahadevan 
139181f196bSVijay Mahadevan   PetscFunctionBegin;
140181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 9);
141181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 10);
142181f196bSVijay Mahadevan   PetscValidPointer(volume, 11);
143181f196bSVijay Mahadevan   if (phypts) {
14463d025dbSVijay Mahadevan     ierr = PetscMemzero(phypts, npts * 3 * sizeof(PetscReal));CHKERRQ(ierr);
145181f196bSVijay Mahadevan   }
14663d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
14763d025dbSVijay Mahadevan     ierr = PetscMemzero(dphidx, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr);
14863d025dbSVijay Mahadevan   }
14963d025dbSVijay Mahadevan   if (nverts == 2) { /* Linear Edge */
15063d025dbSVijay Mahadevan 
15163d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
15263d025dbSVijay Mahadevan     {
153a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
154181f196bSVijay Mahadevan       const PetscReal r = quad[j];
15563d025dbSVijay Mahadevan 
156181f196bSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r );
157181f196bSVijay Mahadevan       phi[1 + offset] = (       r );
15863d025dbSVijay Mahadevan 
159181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 1.0 };
16063d025dbSVijay Mahadevan 
161181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
16263d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
163181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
164181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
165181f196bSVijay Mahadevan         if (phypts) {
166181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
167181f196bSVijay Mahadevan         }
16863d025dbSVijay Mahadevan       }
16963d025dbSVijay Mahadevan 
17063d025dbSVijay Mahadevan       /* invert the jacobian */
171181f196bSVijay Mahadevan       *volume = jacobian[0];
172181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
173181f196bSVijay Mahadevan       jxw[j] *= *volume;
17463d025dbSVijay Mahadevan 
17563d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
17663d025dbSVijay Mahadevan       for ( i = 0; i < nverts; i++ ) {
177181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
17863d025dbSVijay Mahadevan       }
17963d025dbSVijay Mahadevan 
18063d025dbSVijay Mahadevan     }
18163d025dbSVijay Mahadevan   }
18263d025dbSVijay Mahadevan   else if (nverts == 3) { /* Quadratic Edge */
18363d025dbSVijay Mahadevan 
18463d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
18563d025dbSVijay Mahadevan     {
186a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
187181f196bSVijay Mahadevan       const PetscReal r = quad[j];
18863d025dbSVijay Mahadevan 
189181f196bSVijay Mahadevan       phi[0 + offset] = 1.0 + r * ( 2.0 * r - 3.0 );
190181f196bSVijay Mahadevan       phi[1 + offset] = 4.0 * r * ( 1.0 - r );
191181f196bSVijay Mahadevan       phi[2 + offset] = r * ( 2.0 * r - 1.0 );
19263d025dbSVijay Mahadevan 
193181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r ), 4.0 * r - 1.0};
19463d025dbSVijay Mahadevan 
195181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
19663d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
197181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
198181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
199181f196bSVijay Mahadevan         if (phypts) {
200181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
201181f196bSVijay Mahadevan         }
20263d025dbSVijay Mahadevan       }
20363d025dbSVijay Mahadevan 
20463d025dbSVijay Mahadevan       /* invert the jacobian */
205181f196bSVijay Mahadevan       *volume = jacobian[0];
206181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
207181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
20863d025dbSVijay Mahadevan 
20963d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
21063d025dbSVijay Mahadevan       for ( i = 0; i < nverts; i++ ) {
211181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
21263d025dbSVijay Mahadevan       }
21363d025dbSVijay Mahadevan     }
21463d025dbSVijay Mahadevan   }
21563d025dbSVijay Mahadevan   else {
21663d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support EDGE2 and EDGE3 basis evaluations in 1-D : %D", nverts);
21763d025dbSVijay Mahadevan   }
21863d025dbSVijay Mahadevan #if 0
21963d025dbSVijay Mahadevan   /* verify if the computed basis functions are consistent */
22063d025dbSVijay Mahadevan   for ( j = 0; j < npts; j++ ) {
22163d025dbSVijay Mahadevan     PetscScalar phisum = 0, dphixsum = 0;
22263d025dbSVijay Mahadevan     for ( i = 0; i < nverts; i++ ) {
22363d025dbSVijay Mahadevan       phisum += phi[i + j * nverts];
22463d025dbSVijay Mahadevan       if (dphidx) dphixsum += dphidx[i + j * nverts];
22563d025dbSVijay Mahadevan     }
22663d025dbSVijay Mahadevan     PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D = %g, %g, %g\n", j, phisum, dphixsum);
22763d025dbSVijay Mahadevan   }
22863d025dbSVijay Mahadevan #endif
22963d025dbSVijay Mahadevan   PetscFunctionReturn(0);
23063d025dbSVijay Mahadevan }
23163d025dbSVijay Mahadevan 
23263d025dbSVijay Mahadevan 
23363d025dbSVijay Mahadevan /*@
23497b73a88SSatish Balay   Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.
23563d025dbSVijay Mahadevan 
23663d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a quadrangle or triangle.
23763d025dbSVijay Mahadevan 
23863d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
23963d025dbSVijay Mahadevan   and their derivatives with respect to X and Y.
24063d025dbSVijay Mahadevan 
24163d025dbSVijay Mahadevan   Notes:
24263d025dbSVijay Mahadevan 
24363d025dbSVijay Mahadevan   Example Physical Element (QUAD4)
244*a2b725a8SWilliam Gropp .vb
24563d025dbSVijay Mahadevan     4------3        s
24663d025dbSVijay Mahadevan     |      |        |
24763d025dbSVijay Mahadevan     |      |        |
24863d025dbSVijay Mahadevan     |      |        |
24963d025dbSVijay Mahadevan     1------2        0-------r
250*a2b725a8SWilliam Gropp .ve
25163d025dbSVijay Mahadevan 
25263d025dbSVijay Mahadevan   Input Parameter:
253*a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
254*a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
255*a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
256*a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
25763d025dbSVijay Mahadevan 
25863d025dbSVijay Mahadevan   Output Parameter:
259*a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
260*a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
261*a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
262*a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
263*a2b725a8SWilliam Gropp -  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
26463d025dbSVijay Mahadevan 
265edc382c3SSatish Balay   Level: advanced
266edc382c3SSatish Balay 
26763d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 2-D
26863d025dbSVijay Mahadevan @*/
26963d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_2D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
270181f196bSVijay Mahadevan     PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy,
271181f196bSVijay Mahadevan     PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
27263d025dbSVijay Mahadevan {
273a86ed394SVijay Mahadevan   PetscInt i, j, k;
27463d025dbSVijay Mahadevan   PetscErrorCode   ierr;
27563d025dbSVijay Mahadevan 
27663d025dbSVijay Mahadevan   PetscFunctionBegin;
277181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 10);
278181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 11);
279181f196bSVijay Mahadevan   PetscValidPointer(volume, 12);
28063d025dbSVijay Mahadevan   ierr = PetscMemzero(phi, npts * sizeof(PetscReal));CHKERRQ(ierr);
281181f196bSVijay Mahadevan   if (phypts) {
28263d025dbSVijay Mahadevan     ierr = PetscMemzero(phypts, npts * 3 * sizeof(PetscReal));CHKERRQ(ierr);
283181f196bSVijay Mahadevan   }
28463d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
28563d025dbSVijay Mahadevan     ierr = PetscMemzero(dphidx, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr);
28663d025dbSVijay Mahadevan     ierr = PetscMemzero(dphidy, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr);
28763d025dbSVijay Mahadevan   }
28863d025dbSVijay Mahadevan   if (nverts == 4) { /* Linear Quadrangle */
28963d025dbSVijay Mahadevan 
29063d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
29163d025dbSVijay Mahadevan     {
292a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
293181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
294181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
29563d025dbSVijay Mahadevan 
29663d025dbSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r ) * ( 1.0 - s );
29763d025dbSVijay Mahadevan       phi[1 + offset] =         r   * ( 1.0 - s );
29863d025dbSVijay Mahadevan       phi[2 + offset] =         r   *         s;
29963d025dbSVijay Mahadevan       phi[3 + offset] = ( 1.0 - r ) *         s;
30063d025dbSVijay Mahadevan 
301181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
302181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
30363d025dbSVijay Mahadevan 
30463d025dbSVijay Mahadevan       ierr = PetscMemzero(jacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr);
30563d025dbSVijay Mahadevan       ierr = PetscMemzero(ijacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr);
30663d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
307181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
30863d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
30963d025dbSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
31063d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
31163d025dbSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
312181f196bSVijay Mahadevan         if (phypts) {
313181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
314181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
315181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
316181f196bSVijay Mahadevan         }
31763d025dbSVijay Mahadevan       }
31863d025dbSVijay Mahadevan 
31963d025dbSVijay Mahadevan       /* invert the jacobian */
320181f196bSVijay Mahadevan       ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
321181f196bSVijay Mahadevan       if ( *volume < 1e-12 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Quadrangular element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
32263d025dbSVijay Mahadevan 
323181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
32463d025dbSVijay Mahadevan 
325181f196bSVijay Mahadevan       /*  Let us compute the bases derivatives by scaling with inverse jacobian. */
326181f196bSVijay Mahadevan       if (dphidx) {
32763d025dbSVijay Mahadevan         for ( i = 0; i < nverts; i++ ) {
328a86ed394SVijay Mahadevan           for ( k = 0; k < 2; ++k) {
32963d025dbSVijay Mahadevan             if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0];
33063d025dbSVijay Mahadevan             if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1];
331181f196bSVijay Mahadevan           } /* for k=[0..2) */
332181f196bSVijay Mahadevan         } /* for i=[0..nverts) */
333181f196bSVijay Mahadevan       } /* if (dphidx) */
33463d025dbSVijay Mahadevan     }
33563d025dbSVijay Mahadevan   }
33663d025dbSVijay Mahadevan   else if (nverts == 3) { /* Linear triangle */
33763d025dbSVijay Mahadevan 
33863d025dbSVijay Mahadevan     ierr = PetscMemzero(jacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr);
33963d025dbSVijay Mahadevan     ierr = PetscMemzero(ijacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr);
34063d025dbSVijay Mahadevan 
341181f196bSVijay Mahadevan     const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1];
342181f196bSVijay Mahadevan 
34363d025dbSVijay Mahadevan     /* Jacobian is constant */
344181f196bSVijay Mahadevan     jacobian[0] = (coords[0 * 3 + 0] - x2); jacobian[1] = (coords[1 * 3 + 0] - x2);
345181f196bSVijay Mahadevan     jacobian[2] = (coords[0 * 3 + 1] - y2); jacobian[3] = (coords[1 * 3 + 1] - y2);
34663d025dbSVijay Mahadevan 
34763d025dbSVijay Mahadevan     /* invert the jacobian */
348181f196bSVijay Mahadevan     ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
349181f196bSVijay Mahadevan     if ( *volume < 1e-12 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangular element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
350181f196bSVijay Mahadevan 
351181f196bSVijay Mahadevan     const PetscReal Dx[3] = { ijacobian[0], ijacobian[2], - ijacobian[0] - ijacobian[2] };
352181f196bSVijay Mahadevan     const PetscReal Dy[3] = { ijacobian[1], ijacobian[3], - ijacobian[1] - ijacobian[3] };
35363d025dbSVijay Mahadevan 
35463d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
35563d025dbSVijay Mahadevan     {
356a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
357181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
358181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
35963d025dbSVijay Mahadevan 
360181f196bSVijay Mahadevan       if (jxw) jxw[j] *= 0.5;
36163d025dbSVijay Mahadevan 
362181f196bSVijay Mahadevan       /* const PetscReal Ni[3]  = { r, s, 1.0 - r - s }; */
363181f196bSVijay Mahadevan       const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s;
364181f196bSVijay Mahadevan       const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s;
365181f196bSVijay Mahadevan       if (phypts) {
366181f196bSVijay Mahadevan         phypts[offset + 0] = phipts_x;
367181f196bSVijay Mahadevan         phypts[offset + 1] = phipts_y;
368181f196bSVijay Mahadevan       }
36963d025dbSVijay Mahadevan 
370181f196bSVijay Mahadevan       /* \phi_0 = (b.y − c.y) x + (b.x − c.x) y + c.x b.y − b.x c.y */
371181f196bSVijay Mahadevan       phi[0 + offset] = (  ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2) );
372181f196bSVijay Mahadevan       /* \phi_1 = (c.y − a.y) x + (a.x − c.x) y + c.x a.y − a.x c.y */
373181f196bSVijay Mahadevan       phi[1 + offset] = (  ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2) );
37463d025dbSVijay Mahadevan       phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset];
37563d025dbSVijay Mahadevan 
37663d025dbSVijay Mahadevan       if (dphidx) {
377181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
378181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
379181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
38063d025dbSVijay Mahadevan       }
38163d025dbSVijay Mahadevan 
38263d025dbSVijay Mahadevan       if (dphidy) {
383181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
384181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
385181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
38663d025dbSVijay Mahadevan       }
38763d025dbSVijay Mahadevan 
38863d025dbSVijay Mahadevan     }
38963d025dbSVijay Mahadevan   }
39063d025dbSVijay Mahadevan   else {
39163d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support QUAD4 and TRI3 basis evaluations in 2-D : %D", nverts);
39263d025dbSVijay Mahadevan   }
39363d025dbSVijay Mahadevan #if 0
39463d025dbSVijay Mahadevan   /* verify if the computed basis functions are consistent */
39563d025dbSVijay Mahadevan   for ( j = 0; j < npts; j++ ) {
39663d025dbSVijay Mahadevan     PetscScalar phisum = 0, dphixsum = 0, dphiysum = 0;
39763d025dbSVijay Mahadevan     for ( i = 0; i < nverts; i++ ) {
39863d025dbSVijay Mahadevan       phisum += phi[i + j * nverts];
39963d025dbSVijay Mahadevan       if (dphidx) dphixsum += dphidx[i + j * nverts];
40063d025dbSVijay Mahadevan       if (dphidy) dphiysum += dphidy[i + j * nverts];
40163d025dbSVijay Mahadevan     }
40263d025dbSVijay Mahadevan     PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D = %g, %g, %g\n", j, phisum, dphixsum, dphiysum);
40363d025dbSVijay Mahadevan   }
40463d025dbSVijay Mahadevan #endif
40563d025dbSVijay Mahadevan   PetscFunctionReturn(0);
40663d025dbSVijay Mahadevan }
40763d025dbSVijay Mahadevan 
40863d025dbSVijay Mahadevan 
40963d025dbSVijay Mahadevan /*@
41097b73a88SSatish Balay   Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element.
41163d025dbSVijay Mahadevan 
41263d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a hexahedra or tetrahedra.
41363d025dbSVijay Mahadevan 
41463d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
41563d025dbSVijay Mahadevan   and their derivatives with respect to X, Y and Z.
41663d025dbSVijay Mahadevan 
41763d025dbSVijay Mahadevan   Notes:
41863d025dbSVijay Mahadevan 
41963d025dbSVijay Mahadevan   Example Physical Element (HEX8)
420*a2b725a8SWilliam Gropp .vb
42163d025dbSVijay Mahadevan       8------7
42263d025dbSVijay Mahadevan      /|     /|        t  s
42363d025dbSVijay Mahadevan     5------6 |        | /
42463d025dbSVijay Mahadevan     | |    | |        |/
42563d025dbSVijay Mahadevan     | 4----|-3        0-------r
42663d025dbSVijay Mahadevan     |/     |/
42763d025dbSVijay Mahadevan     1------2
428*a2b725a8SWilliam Gropp .ve
42963d025dbSVijay Mahadevan 
43063d025dbSVijay Mahadevan   Input Parameter:
431*a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
432*a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
433*a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
434*a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
43563d025dbSVijay Mahadevan 
43663d025dbSVijay Mahadevan   Output Parameter:
437*a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
438*a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
439*a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
440*a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
441*a2b725a8SWilliam Gropp .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
442*a2b725a8SWilliam Gropp -  PetscReal dphidz[npts] -    the derivative of the bases wrt Z-direction evaluated at the specified quadrature points
44363d025dbSVijay Mahadevan 
444edc382c3SSatish Balay   Level: advanced
445edc382c3SSatish Balay 
44663d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 3-D
44763d025dbSVijay Mahadevan @*/
44863d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_3D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
449181f196bSVijay Mahadevan     PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz,
450181f196bSVijay Mahadevan     PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
45163d025dbSVijay Mahadevan {
452a86ed394SVijay Mahadevan   PetscInt i, j, k;
45363d025dbSVijay Mahadevan   PetscErrorCode ierr;
45463d025dbSVijay Mahadevan 
45563d025dbSVijay Mahadevan   PetscFunctionBegin;
456181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 11);
457181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 12);
458181f196bSVijay Mahadevan   PetscValidPointer(volume, 13);
45963d025dbSVijay Mahadevan   /* Reset arrays. */
46063d025dbSVijay Mahadevan   ierr = PetscMemzero(phi, npts * sizeof(PetscReal));CHKERRQ(ierr);
461181f196bSVijay Mahadevan   if (phypts) {
46263d025dbSVijay Mahadevan     ierr = PetscMemzero(phypts, npts * 3 * sizeof(PetscReal));CHKERRQ(ierr);
463181f196bSVijay Mahadevan   }
46463d025dbSVijay Mahadevan   if (dphidx) {
46563d025dbSVijay Mahadevan     ierr = PetscMemzero(dphidx, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr);
46663d025dbSVijay Mahadevan     ierr = PetscMemzero(dphidy, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr);
46763d025dbSVijay Mahadevan     ierr = PetscMemzero(dphidz, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr);
46863d025dbSVijay Mahadevan   }
46963d025dbSVijay Mahadevan 
47063d025dbSVijay Mahadevan   if (nverts == 8) { /* Linear Hexahedra */
47163d025dbSVijay Mahadevan 
47263d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
47363d025dbSVijay Mahadevan     {
474a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
475181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
476181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
477181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
47863d025dbSVijay Mahadevan 
479a86ed394SVijay Mahadevan       phi[offset + 0] = ( 1.0 - r ) * ( 1.0 - s ) * ( 1.0 - t ); /* 0,0,0 */
480a86ed394SVijay Mahadevan       phi[offset + 1] = (       r ) * ( 1.0 - s ) * ( 1.0 - t ); /* 1,0,0 */
481a86ed394SVijay Mahadevan       phi[offset + 2] = (       r ) * (       s ) * ( 1.0 - t ); /* 1,1,0 */
482a86ed394SVijay Mahadevan       phi[offset + 3] = ( 1.0 - r ) * (       s ) * ( 1.0 - t ); /* 0,1,0 */
483a86ed394SVijay Mahadevan       phi[offset + 4] = ( 1.0 - r ) * ( 1.0 - s ) * (       t ); /* 0,0,1 */
484a86ed394SVijay Mahadevan       phi[offset + 5] = (       r ) * ( 1.0 - s ) * (       t ); /* 1,0,1 */
485a86ed394SVijay Mahadevan       phi[offset + 6] = (       r ) * (       s ) * (       t ); /* 1,1,1 */
486a86ed394SVijay Mahadevan       phi[offset + 7] = ( 1.0 - r ) * (       s ) * (       t ); /* 0,1,1 */
48763d025dbSVijay Mahadevan 
488181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s ) * ( 1.0 - t ),
48963d025dbSVijay Mahadevan                                        ( 1.0 - s ) * ( 1.0 - t ),
490a86ed394SVijay Mahadevan                                        (       s ) * ( 1.0 - t ),
491a86ed394SVijay Mahadevan                                      - (       s ) * ( 1.0 - t ),
492a86ed394SVijay Mahadevan                                      - ( 1.0 - s ) * (       t ),
493a86ed394SVijay Mahadevan                                        ( 1.0 - s ) * (       t ),
494a86ed394SVijay Mahadevan                                        (       s ) * (       t ),
495a86ed394SVijay Mahadevan                                      - (       s ) * (       t )
49663d025dbSVijay Mahadevan                                     };
49763d025dbSVijay Mahadevan 
498181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r ) * ( 1.0 - t ),
499a86ed394SVijay Mahadevan                                        - (       r ) * ( 1.0 - t ),
500a86ed394SVijay Mahadevan                                          (       r ) * ( 1.0 - t ),
50163d025dbSVijay Mahadevan                                          ( 1.0 - r ) * ( 1.0 - t ),
502a86ed394SVijay Mahadevan                                        - ( 1.0 - r ) * (       t ),
503a86ed394SVijay Mahadevan                                        - (       r ) * (       t ),
504a86ed394SVijay Mahadevan                                          (       r ) * (       t ),
505a86ed394SVijay Mahadevan                                          ( 1.0 - r ) * (       t )
50663d025dbSVijay Mahadevan                                       };
50763d025dbSVijay Mahadevan 
508181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r ) * ( 1.0 - s ),
509a86ed394SVijay Mahadevan                                         - (       r ) * ( 1.0 - s ),
510a86ed394SVijay Mahadevan                                         - (       r ) * (       s ),
511a86ed394SVijay Mahadevan                                         - ( 1.0 - r ) * (       s ),
51263d025dbSVijay Mahadevan                                           ( 1.0 - r ) * ( 1.0 - s ),
513a86ed394SVijay Mahadevan                                           (       r ) * ( 1.0 - s ),
514a86ed394SVijay Mahadevan                                           (       r ) * (       s ),
515a86ed394SVijay Mahadevan                                           ( 1.0 - r ) * (       s )
51663d025dbSVijay Mahadevan                                      };
51763d025dbSVijay Mahadevan 
51863d025dbSVijay Mahadevan       ierr = PetscMemzero(jacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr);
51963d025dbSVijay Mahadevan       ierr = PetscMemzero(ijacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr);
52063d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
521181f196bSVijay Mahadevan         const PetscReal* vertex = coords + i * 3;
52263d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
52363d025dbSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
52463d025dbSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
52563d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
52663d025dbSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
52763d025dbSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
52863d025dbSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
52963d025dbSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
53063d025dbSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
531181f196bSVijay Mahadevan         if (phypts) {
532181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertex[0];
533181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertex[1];
534181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertex[2];
535181f196bSVijay Mahadevan         }
53663d025dbSVijay Mahadevan       }
53763d025dbSVijay Mahadevan 
53863d025dbSVijay Mahadevan       /* invert the jacobian */
539181f196bSVijay Mahadevan       ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
540a86ed394SVijay Mahadevan       if ( *volume < 1e-8 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Hexahedral element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
54163d025dbSVijay Mahadevan 
542a86ed394SVijay Mahadevan       if (jxw) jxw[j] *= (*volume);
54363d025dbSVijay Mahadevan 
54463d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
54563d025dbSVijay Mahadevan       for ( i = 0; i < nverts; ++i ) {
546a86ed394SVijay Mahadevan         for ( k = 0; k < 3; ++k) {
54763d025dbSVijay Mahadevan           if (dphidx) dphidx[i + offset] += dNi_dxi[i]   * ijacobian[0 * 3 + k];
54863d025dbSVijay Mahadevan           if (dphidy) dphidy[i + offset] += dNi_deta[i]  * ijacobian[1 * 3 + k];
54963d025dbSVijay Mahadevan           if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k];
55063d025dbSVijay Mahadevan         }
55163d025dbSVijay Mahadevan       }
55263d025dbSVijay Mahadevan     }
55363d025dbSVijay Mahadevan   }
55463d025dbSVijay Mahadevan   else if (nverts == 4) { /* Linear Tetrahedra */
55563d025dbSVijay Mahadevan 
55663d025dbSVijay Mahadevan     ierr = PetscMemzero(jacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr);
55763d025dbSVijay Mahadevan     ierr = PetscMemzero(ijacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr);
55863d025dbSVijay Mahadevan 
559181f196bSVijay Mahadevan     const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2];
560181f196bSVijay Mahadevan 
561181f196bSVijay Mahadevan     /* compute the jacobian */
562181f196bSVijay Mahadevan     jacobian[0] = coords[1 * 3 + 0] - x0;  jacobian[1] = coords[2 * 3 + 0] - x0; jacobian[2] = coords[3 * 3 + 0] - x0;
563181f196bSVijay Mahadevan     jacobian[3] = coords[1 * 3 + 1] - y0;  jacobian[4] = coords[2 * 3 + 1] - y0; jacobian[5] = coords[3 * 3 + 1] - y0;
564181f196bSVijay Mahadevan     jacobian[6] = coords[1 * 3 + 2] - z0;  jacobian[7] = coords[2 * 3 + 2] - z0; jacobian[8] = coords[3 * 3 + 2] - z0;
56563d025dbSVijay Mahadevan 
56663d025dbSVijay Mahadevan     /* invert the jacobian */
567181f196bSVijay Mahadevan     ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
568a86ed394SVijay Mahadevan     if ( *volume < 1e-8 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
56963d025dbSVijay Mahadevan 
5706ea892caSVijay Mahadevan     PetscReal Dx[4]={0,0,0,0}, Dy[4]={0,0,0,0}, Dz[4]={0,0,0,0};
571181f196bSVijay Mahadevan     if (dphidx) {
572181f196bSVijay Mahadevan       Dx[0] =   ( coords[1 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
573181f196bSVijay Mahadevan                  - coords[1 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
574181f196bSVijay Mahadevan                  - coords[1 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3] )
575181f196bSVijay Mahadevan                 ) / *volume;
576181f196bSVijay Mahadevan       Dx[1] = - ( coords[1 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
577181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
578181f196bSVijay Mahadevan                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3] )
579181f196bSVijay Mahadevan                 ) / *volume;
580181f196bSVijay Mahadevan       Dx[2] =   ( coords[1 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
581181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
582181f196bSVijay Mahadevan                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3] )
583181f196bSVijay Mahadevan                 ) / *volume;
584181f196bSVijay Mahadevan       Dx[3] =  - ( Dx[0] + Dx[1] + Dx[2] );
585181f196bSVijay Mahadevan       Dy[0] =   ( coords[0 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
586181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
587181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3] )
588181f196bSVijay Mahadevan                 ) / *volume;
589181f196bSVijay Mahadevan       Dy[1] = - ( coords[0 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
590181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
591181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3] )
592181f196bSVijay Mahadevan                 ) / *volume;
593181f196bSVijay Mahadevan       Dy[2] =   ( coords[0 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
594181f196bSVijay Mahadevan                  - coords[0 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
595181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3] )
596181f196bSVijay Mahadevan                 ) / *volume;
597181f196bSVijay Mahadevan       Dy[3] =  - ( Dy[0] + Dy[1] + Dy[2] );
598181f196bSVijay Mahadevan       Dz[0] =   ( coords[0 + 1 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3] )
599181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3] )
600181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * (coords[1 + 2 * 3] - coords[1 + 1 * 3] )
601181f196bSVijay Mahadevan                 ) / *volume;
602181f196bSVijay Mahadevan       Dz[1] = - ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3] )
603181f196bSVijay Mahadevan                   + coords[0 + 2 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3] )
604181f196bSVijay Mahadevan                   - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 2 * 3] )
605181f196bSVijay Mahadevan                 ) / *volume;
606181f196bSVijay Mahadevan       Dz[2] =   ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3] )
607181f196bSVijay Mahadevan                  + coords[0 + 1 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3] )
608181f196bSVijay Mahadevan                  - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 1 * 3] )
609181f196bSVijay Mahadevan                 ) / *volume;
610181f196bSVijay Mahadevan       Dz[3] =  - ( Dz[0] + Dz[1] + Dz[2] );
611181f196bSVijay Mahadevan     }
61263d025dbSVijay Mahadevan 
61363d025dbSVijay Mahadevan     for ( j = 0; j < npts; j++ )
61463d025dbSVijay Mahadevan     {
615a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
616181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
617181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
618181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
61963d025dbSVijay Mahadevan 
620181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
62163d025dbSVijay Mahadevan 
62263d025dbSVijay Mahadevan       phi[offset + 0] = 1.0 - r - s - t;
62363d025dbSVijay Mahadevan       phi[offset + 1] = r;
62463d025dbSVijay Mahadevan       phi[offset + 2] = s;
62563d025dbSVijay Mahadevan       phi[offset + 3] = t;
62663d025dbSVijay Mahadevan 
627181f196bSVijay Mahadevan       if (phypts) {
628181f196bSVijay Mahadevan         for (i = 0; i < nverts; ++i) {
629181f196bSVijay Mahadevan           const PetscScalar* vertices = coords + i * 3;
630181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
631181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
632181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
633181f196bSVijay Mahadevan         }
634181f196bSVijay Mahadevan       }
635181f196bSVijay Mahadevan 
636181f196bSVijay Mahadevan       /* Now set the derivatives */
63763d025dbSVijay Mahadevan       if (dphidx) {
638181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
639181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
640181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
641181f196bSVijay Mahadevan         dphidx[3 + offset] = Dx[3];
64263d025dbSVijay Mahadevan 
643181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
644181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
645181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
646181f196bSVijay Mahadevan         dphidy[3 + offset] = Dy[3];
64763d025dbSVijay Mahadevan 
648181f196bSVijay Mahadevan         dphidz[0 + offset] = Dz[0];
649181f196bSVijay Mahadevan         dphidz[1 + offset] = Dz[1];
650181f196bSVijay Mahadevan         dphidz[2 + offset] = Dz[2];
651181f196bSVijay Mahadevan         dphidz[3 + offset] = Dz[3];
65263d025dbSVijay Mahadevan       }
65363d025dbSVijay Mahadevan 
65463d025dbSVijay Mahadevan     } /* Tetrahedra -- ends */
65563d025dbSVijay Mahadevan   }
65663d025dbSVijay Mahadevan   else
65763d025dbSVijay Mahadevan   {
65863d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support HEX8 and TET4 basis evaluations in 3-D : %D", nverts);
65963d025dbSVijay Mahadevan   }
66063d025dbSVijay Mahadevan #if 0
66163d025dbSVijay Mahadevan   /* verify if the computed basis functions are consistent */
66263d025dbSVijay Mahadevan   for ( j = 0; j < npts; j++ ) {
66363d025dbSVijay Mahadevan     PetscScalar phisum = 0, dphixsum = 0, dphiysum = 0, dphizsum = 0;
66463d025dbSVijay Mahadevan     const int offset = j * nverts;
66563d025dbSVijay Mahadevan     for ( i = 0; i < nverts; i++ ) {
66663d025dbSVijay Mahadevan       phisum += phi[i + offset];
66763d025dbSVijay Mahadevan       if (dphidx) dphixsum += dphidx[i + offset];
66863d025dbSVijay Mahadevan       if (dphidy) dphiysum += dphidy[i + offset];
66963d025dbSVijay Mahadevan       if (dphidz) dphizsum += dphidz[i + offset];
67063d025dbSVijay Mahadevan       if (dphidx) PetscPrintf(PETSC_COMM_WORLD, "\t Values [%d]: [JxW] [phi, dphidx, dphidy, dphidz] = %g, %g, %g, %g, %g\n", j, jxw[j], phi[i + offset], dphidx[i + offset], dphidy[i + offset], dphidz[i + offset]);
67163d025dbSVijay Mahadevan     }
67263d025dbSVijay Mahadevan     if (dphidx) PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D (%g, %g, %g) = %g, %g, %g, %g\n", j, quad[3 * j + 0], quad[3 * j + 1], quad[3 * j + 2], phisum, dphixsum, dphiysum, dphizsum);
67363d025dbSVijay Mahadevan   }
67463d025dbSVijay Mahadevan #endif
67563d025dbSVijay Mahadevan   PetscFunctionReturn(0);
67663d025dbSVijay Mahadevan }
67763d025dbSVijay Mahadevan 
67863d025dbSVijay Mahadevan 
67963d025dbSVijay Mahadevan 
68063d025dbSVijay Mahadevan /*@
68197b73a88SSatish Balay   DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
68263d025dbSVijay Mahadevan 
68363d025dbSVijay Mahadevan   The routine takes the coordinates of the vertices of an element and computes basis functions associated with
68463d025dbSVijay Mahadevan   each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
68563d025dbSVijay Mahadevan 
68663d025dbSVijay Mahadevan   Input Parameter:
687*a2b725a8SWilliam Gropp +  PetscInt  nverts,           the number of element vertices
68863d025dbSVijay Mahadevan .  PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering)
68963d025dbSVijay Mahadevan .  PetscInt  npts,             the number of evaluation points (quadrature points)
690*a2b725a8SWilliam Gropp -  PetscReal quad[3*npts],     the evaluation points (quadrature points) in the reference space
69163d025dbSVijay Mahadevan 
69263d025dbSVijay Mahadevan   Output Parameter:
693*a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts],   the evaluation points (quadrature points) transformed to the physical space
69463d025dbSVijay Mahadevan .  PetscReal jxw[npts],        the jacobian determinant * quadrature weight necessary for assembling discrete contributions
69563d025dbSVijay Mahadevan .  PetscReal fe_basis[npts],        the bases values evaluated at the specified quadrature points
696*a2b725a8SWilliam Gropp -  PetscReal fe_basis_derivatives[dim][npts],  the derivative of the bases wrt (X,Y,Z)-directions (depending on the dimension) evaluated at the specified quadrature points
69763d025dbSVijay Mahadevan 
698edc382c3SSatish Balay   Level: advanced
699edc382c3SSatish Balay 
70063d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 3-D
70163d025dbSVijay Mahadevan @*/
702181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMComputeBasis ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature,
703181f196bSVijay Mahadevan                                        PetscReal *phypts, PetscReal *jacobian_quadrature_weight_product,
704181f196bSVijay Mahadevan                                        PetscReal *fe_basis, PetscReal **fe_basis_derivatives)
70563d025dbSVijay Mahadevan {
70663d025dbSVijay Mahadevan   PetscErrorCode  ierr;
707b21a1e88SVijay Mahadevan   PetscInt        npoints,idim;
70863d025dbSVijay Mahadevan   bool            compute_der;
70963d025dbSVijay Mahadevan   const PetscReal *quadpts, *quadwts;
710181f196bSVijay Mahadevan   PetscReal       jacobian[9], ijacobian[9], volume;
71163d025dbSVijay Mahadevan 
71263d025dbSVijay Mahadevan   PetscFunctionBegin;
71363d025dbSVijay Mahadevan   PetscValidPointer(coordinates, 3);
71463d025dbSVijay Mahadevan   PetscValidHeaderSpecific(quadrature, PETSC_OBJECT_CLASSID, 4);
71563d025dbSVijay Mahadevan   PetscValidPointer(fe_basis, 7);
71663d025dbSVijay Mahadevan   compute_der = (fe_basis_derivatives != NULL);
71763d025dbSVijay Mahadevan 
71863d025dbSVijay Mahadevan   /* Get the quadrature points and weights for the given quadrature rule */
719b21a1e88SVijay Mahadevan   ierr = PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts);CHKERRQ(ierr);
720b21a1e88SVijay Mahadevan   if (idim != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%D) vs quadrature (%D)\n",idim,dim);
721181f196bSVijay Mahadevan   if (jacobian_quadrature_weight_product) {
72263d025dbSVijay Mahadevan     ierr = PetscMemcpy(jacobian_quadrature_weight_product, quadwts, npoints * sizeof(PetscReal));CHKERRQ(ierr);
723181f196bSVijay Mahadevan   }
72463d025dbSVijay Mahadevan 
72563d025dbSVijay Mahadevan   switch (dim) {
72663d025dbSVijay Mahadevan   case 1:
72763d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, npoints, quadpts, phypts,
72863d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
729181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
730181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
73163d025dbSVijay Mahadevan     break;
73263d025dbSVijay Mahadevan   case 2:
73363d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, npoints, quadpts, phypts,
73463d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
73563d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
736181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[1] : NULL),
737181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
73863d025dbSVijay Mahadevan     break;
73963d025dbSVijay Mahadevan   case 3:
74063d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, npoints, quadpts, phypts,
74163d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
74263d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
74363d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[1] : NULL),
744181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[2] : NULL),
745181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
74663d025dbSVijay Mahadevan     break;
74763d025dbSVijay Mahadevan   default:
74863d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
74963d025dbSVijay Mahadevan   }
75063d025dbSVijay Mahadevan   PetscFunctionReturn(0);
75163d025dbSVijay Mahadevan }
75263d025dbSVijay Mahadevan 
75363d025dbSVijay Mahadevan 
75463d025dbSVijay Mahadevan 
75563d025dbSVijay Mahadevan /*@
75697b73a88SSatish Balay   DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given
75763d025dbSVijay Mahadevan   dimension and polynomial order (deciphered from number of element vertices).
75863d025dbSVijay Mahadevan 
75963d025dbSVijay Mahadevan   Input Parameter:
76063d025dbSVijay Mahadevan 
761*a2b725a8SWilliam Gropp +  PetscInt  dim   -   the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
762*a2b725a8SWilliam Gropp -  PetscInt nverts -   the number of vertices in the physical element
76363d025dbSVijay Mahadevan 
76463d025dbSVijay Mahadevan   Output Parameter:
765*a2b725a8SWilliam Gropp .  PetscQuadrature quadrature -  the quadrature object with default settings to integrate polynomials defined over the element
76663d025dbSVijay Mahadevan 
767edc382c3SSatish Balay   Level: advanced
768edc382c3SSatish Balay 
76963d025dbSVijay Mahadevan .keywords: DMMoab, Quadrature, PetscDT
77063d025dbSVijay Mahadevan @*/
771181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMCreateQuadratureDefault ( const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature )
77263d025dbSVijay Mahadevan {
77363d025dbSVijay Mahadevan   PetscReal *w, *x;
774b21a1e88SVijay Mahadevan   PetscInt nc=1;
77563d025dbSVijay Mahadevan   PetscErrorCode  ierr;
77663d025dbSVijay Mahadevan 
77763d025dbSVijay Mahadevan   PetscFunctionBegin;
77863d025dbSVijay Mahadevan   /* Create an appropriate quadrature rule to sample basis */
77963d025dbSVijay Mahadevan   switch (dim)
78063d025dbSVijay Mahadevan   {
78163d025dbSVijay Mahadevan   case 1:
78263d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
783b21a1e88SVijay Mahadevan     ierr = PetscDTGaussJacobiQuadrature(1, nc, nverts, 0, 1.0, quadrature);CHKERRQ(ierr);
78463d025dbSVijay Mahadevan     break;
78563d025dbSVijay Mahadevan   case 2:
78663d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
78763d025dbSVijay Mahadevan     if (nverts == 3) {
788a86ed394SVijay Mahadevan       const PetscInt order = 2;
789a86ed394SVijay Mahadevan       const PetscInt npoints = (order == 2 ? 3 : 6);
79063d025dbSVijay Mahadevan       ierr = PetscMalloc2(npoints * 2, &x, npoints, &w);CHKERRQ(ierr);
791181f196bSVijay Mahadevan       if (npoints == 3) {
79263d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0;
79363d025dbSVijay Mahadevan         x[3] = x[4] = 2.0 / 3.0;
79463d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 1.0 / 3.0;
795181f196bSVijay Mahadevan       }
796181f196bSVijay Mahadevan       else if (npoints == 6) {
79763d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = 0.44594849091597;
79863d025dbSVijay Mahadevan         x[3] = x[4] = 0.10810301816807;
79963d025dbSVijay Mahadevan         x[5] = 0.44594849091597;
80063d025dbSVijay Mahadevan         x[6] = x[7] = x[8] = 0.09157621350977;
80163d025dbSVijay Mahadevan         x[9] = x[10] = 0.81684757298046;
80263d025dbSVijay Mahadevan         x[11] = 0.09157621350977;
80363d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 0.22338158967801;
804181f196bSVijay Mahadevan         w[3] = w[4] = w[5] = 0.10995174365532;
805181f196bSVijay Mahadevan       }
806181f196bSVijay Mahadevan       else {
807181f196bSVijay Mahadevan         SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %D", npoints);
80863d025dbSVijay Mahadevan       }
80963d025dbSVijay Mahadevan       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr);
81063d025dbSVijay Mahadevan       ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr);
811b21a1e88SVijay Mahadevan       ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr);
812181f196bSVijay Mahadevan       /* ierr = PetscDTGaussJacobiQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */
81363d025dbSVijay Mahadevan     }
81463d025dbSVijay Mahadevan     else {
815b21a1e88SVijay Mahadevan       ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr);
81663d025dbSVijay Mahadevan     }
81763d025dbSVijay Mahadevan     break;
81863d025dbSVijay Mahadevan   case 3:
81963d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
82063d025dbSVijay Mahadevan     if (nverts == 4) {
821a86ed394SVijay Mahadevan       PetscInt order;
822a86ed394SVijay Mahadevan       const PetscInt npoints = 4; // Choose between 4 and 10
823181f196bSVijay Mahadevan       ierr = PetscMalloc2(npoints * 3, &x, npoints, &w);CHKERRQ(ierr);
824181f196bSVijay Mahadevan       if (npoints == 4) { /*  KEAST1, K1,  N=4, O=4 */
825181f196bSVijay Mahadevan         const PetscReal x_4[12] = { 0.5854101966249685, 0.1381966011250105, 0.1381966011250105,
826181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.1381966011250105,
827181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.5854101966249685,
828181f196bSVijay Mahadevan                                     0.1381966011250105, 0.5854101966249685, 0.1381966011250105
829181f196bSVijay Mahadevan                                   };
830181f196bSVijay Mahadevan         ierr = PetscMemcpy(x, x_4, 12 * sizeof(PetscReal));CHKERRQ(ierr);
831181f196bSVijay Mahadevan 
832181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0;
833181f196bSVijay Mahadevan         order = 4;
834181f196bSVijay Mahadevan       }
835181f196bSVijay Mahadevan       else if (npoints == 10) { /*  KEAST3, K3  N=10, O=10 */
836181f196bSVijay Mahadevan         const PetscReal x_10[30] = { 0.5684305841968444, 0.1438564719343852, 0.1438564719343852,
837181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.1438564719343852,
838181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.5684305841968444,
839181f196bSVijay Mahadevan                                      0.1438564719343852, 0.5684305841968444, 0.1438564719343852,
840181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.5000000000000000,
841181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.5000000000000000,
842181f196bSVijay Mahadevan                                      0.5000000000000000, 0.5000000000000000, 0.0000000000000000,
843181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.0000000000000000,
844181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.0000000000000000,
845181f196bSVijay Mahadevan                                      0.0000000000000000, 0.0000000000000000, 0.5000000000000000
846181f196bSVijay Mahadevan                                    };
847181f196bSVijay Mahadevan         ierr = PetscMemcpy(x, x_10, 30 * sizeof(PetscReal));CHKERRQ(ierr);
848181f196bSVijay Mahadevan 
849181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 0.2177650698804054;
850181f196bSVijay Mahadevan         w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631;
851181f196bSVijay Mahadevan         order = 10;
852181f196bSVijay Mahadevan       }
853181f196bSVijay Mahadevan       else {
854181f196bSVijay Mahadevan         SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %D", npoints);
855181f196bSVijay Mahadevan       }
856181f196bSVijay Mahadevan       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr);
857181f196bSVijay Mahadevan       ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr);
858b21a1e88SVijay Mahadevan       ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr);
859181f196bSVijay Mahadevan       /* ierr = PetscDTGaussJacobiQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */
86063d025dbSVijay Mahadevan     }
86163d025dbSVijay Mahadevan     else {
862b21a1e88SVijay Mahadevan       ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr);
86363d025dbSVijay Mahadevan     }
86463d025dbSVijay Mahadevan     break;
86563d025dbSVijay Mahadevan   default:
86663d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
86763d025dbSVijay Mahadevan   }
86863d025dbSVijay Mahadevan   PetscFunctionReturn(0);
86963d025dbSVijay Mahadevan }
87063d025dbSVijay Mahadevan 
871181f196bSVijay Mahadevan /* Compute Jacobians */
872181f196bSVijay Mahadevan PetscErrorCode ComputeJacobian_Internal ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts,
873a86ed394SVijay Mahadevan   PetscReal *jacobian, PetscReal *ijacobian, PetscReal* dvolume)
874181f196bSVijay Mahadevan {
875a86ed394SVijay Mahadevan   PetscInt i;
8762417220eSVijay Mahadevan   PetscReal volume=1.0;
877181f196bSVijay Mahadevan   PetscErrorCode ierr;
878181f196bSVijay Mahadevan 
879181f196bSVijay Mahadevan   PetscFunctionBegin;
880181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
881181f196bSVijay Mahadevan   PetscValidPointer(quad, 4);
882181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 5);
883181f196bSVijay Mahadevan   ierr = PetscMemzero(jacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr);
884181f196bSVijay Mahadevan   if (ijacobian) {
885181f196bSVijay Mahadevan     ierr = PetscMemzero(ijacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr);
886181f196bSVijay Mahadevan   }
887181f196bSVijay Mahadevan   if (phypts) {
888181f196bSVijay Mahadevan     ierr = PetscMemzero(phypts, /*npts=1 * */ 3 * sizeof(PetscReal));CHKERRQ(ierr);
889181f196bSVijay Mahadevan   }
890181f196bSVijay Mahadevan 
891181f196bSVijay Mahadevan   if (dim == 1) {
892181f196bSVijay Mahadevan 
893181f196bSVijay Mahadevan     const PetscReal& r = quad[0];
894181f196bSVijay Mahadevan     if (nverts == 2) { /* Linear Edge */
895181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 1.0 };
896181f196bSVijay Mahadevan 
897181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
898181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
899181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
900181f196bSVijay Mahadevan       }
901181f196bSVijay Mahadevan     }
902181f196bSVijay Mahadevan     else if (nverts == 3) { /* Quadratic Edge */
903181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r ), 4.0 * r - 1.0};
904181f196bSVijay Mahadevan 
905181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
906181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
907181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
908181f196bSVijay Mahadevan       }
909181f196bSVijay Mahadevan     }
910181f196bSVijay Mahadevan     else {
911181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 1-D entity vertices are invalid. Currently only support EDGE2 and EDGE3 basis evaluations in 1-D : %D", nverts);
912181f196bSVijay Mahadevan     }
913181f196bSVijay Mahadevan 
914181f196bSVijay Mahadevan     if (ijacobian) {
915181f196bSVijay Mahadevan       /* invert the jacobian */
916181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
917181f196bSVijay Mahadevan     }
918181f196bSVijay Mahadevan 
919181f196bSVijay Mahadevan   }
920181f196bSVijay Mahadevan   else if (dim == 2) {
921181f196bSVijay Mahadevan 
922181f196bSVijay Mahadevan     if (nverts == 4) { /* Linear Quadrangle */
923181f196bSVijay Mahadevan       const PetscReal& r = quad[0];
924181f196bSVijay Mahadevan       const PetscReal& s = quad[1];
925181f196bSVijay Mahadevan 
926181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
927181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
928181f196bSVijay Mahadevan 
929181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
930181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
931181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]  * vertices[0];
932181f196bSVijay Mahadevan         jacobian[2] += dNi_dxi[i]  * vertices[1];
933181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
934181f196bSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
935181f196bSVijay Mahadevan       }
936181f196bSVijay Mahadevan     }
937181f196bSVijay Mahadevan     else if (nverts == 3) { /* Linear triangle */
938181f196bSVijay Mahadevan       const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1];
939181f196bSVijay Mahadevan 
940181f196bSVijay Mahadevan       /* Jacobian is constant */
941181f196bSVijay Mahadevan       jacobian[0] = (coordinates[0 * 3 + 0] - x2); jacobian[1] = (coordinates[1 * 3 + 0] - x2);
942181f196bSVijay Mahadevan       jacobian[2] = (coordinates[0 * 3 + 1] - y2); jacobian[3] = (coordinates[1 * 3 + 1] - y2);
943181f196bSVijay Mahadevan     }
944181f196bSVijay Mahadevan     else {
945181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 2-D entity vertices are invalid. Currently only support QUAD4 and TRI3 basis evaluations in 2-D : %D", nverts);
946181f196bSVijay Mahadevan     }
947181f196bSVijay Mahadevan 
948181f196bSVijay Mahadevan     /* invert the jacobian */
949181f196bSVijay Mahadevan     if (ijacobian) {
950a86ed394SVijay Mahadevan       ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr);
951181f196bSVijay Mahadevan     }
952181f196bSVijay Mahadevan 
953181f196bSVijay Mahadevan   }
954181f196bSVijay Mahadevan   else {
955181f196bSVijay Mahadevan 
956181f196bSVijay Mahadevan     if (nverts == 8) { /* Linear Hexahedra */
957181f196bSVijay Mahadevan       const PetscReal& r = quad[0];
958181f196bSVijay Mahadevan       const PetscReal& s = quad[1];
959181f196bSVijay Mahadevan       const PetscReal& t = quad[2];
960181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s ) * ( 1.0 - t ),
961181f196bSVijay Mahadevan                                        ( 1.0 - s ) * ( 1.0 - t ),
962a86ed394SVijay Mahadevan                                        (       s ) * ( 1.0 - t ),
963a86ed394SVijay Mahadevan                                      - (       s ) * ( 1.0 - t ),
964a86ed394SVijay Mahadevan                                      - ( 1.0 - s ) * (       t ),
965a86ed394SVijay Mahadevan                                        ( 1.0 - s ) * (       t ),
966a86ed394SVijay Mahadevan                                        (       s ) * (       t ),
967a86ed394SVijay Mahadevan                                      - (       s ) * (       t )
968181f196bSVijay Mahadevan                                     };
969181f196bSVijay Mahadevan 
970181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r ) * ( 1.0 - t ),
971a86ed394SVijay Mahadevan                                        - (       r ) * ( 1.0 - t ),
972a86ed394SVijay Mahadevan                                          (       r ) * ( 1.0 - t ),
973181f196bSVijay Mahadevan                                          ( 1.0 - r ) * ( 1.0 - t ),
974a86ed394SVijay Mahadevan                                        - ( 1.0 - r ) * (       t ),
975a86ed394SVijay Mahadevan                                        - (       r ) * (       t ),
976a86ed394SVijay Mahadevan                                          (       r ) * (       t ),
977a86ed394SVijay Mahadevan                                          ( 1.0 - r ) * (       t )
978181f196bSVijay Mahadevan                                       };
979181f196bSVijay Mahadevan 
980181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r ) * ( 1.0 - s ),
981a86ed394SVijay Mahadevan                                         - (       r ) * ( 1.0 - s ),
982a86ed394SVijay Mahadevan                                         - (       r ) * (       s ),
983a86ed394SVijay Mahadevan                                         - ( 1.0 - r ) * (       s ),
984181f196bSVijay Mahadevan                                           ( 1.0 - r ) * ( 1.0 - s ),
985a86ed394SVijay Mahadevan                                           (       r ) * ( 1.0 - s ),
986a86ed394SVijay Mahadevan                                           (       r ) * (       s ),
987a86ed394SVijay Mahadevan                                           ( 1.0 - r ) * (       s )
988181f196bSVijay Mahadevan                                      };
989a86ed394SVijay Mahadevan 
990181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
991181f196bSVijay Mahadevan         const PetscReal* vertex = coordinates + i * 3;
992181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
993181f196bSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
994181f196bSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
995181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
996181f196bSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
997181f196bSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
998181f196bSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
999181f196bSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
1000181f196bSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
1001181f196bSVijay Mahadevan       }
1002181f196bSVijay Mahadevan 
1003181f196bSVijay Mahadevan     }
1004181f196bSVijay Mahadevan     else if (nverts == 4) { /* Linear Tetrahedra */
1005181f196bSVijay Mahadevan       const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2];
1006181f196bSVijay Mahadevan 
1007181f196bSVijay Mahadevan       /* compute the jacobian */
1008181f196bSVijay Mahadevan       jacobian[0] = coordinates[1 * 3 + 0] - x0;  jacobian[1] = coordinates[2 * 3 + 0] - x0; jacobian[2] = coordinates[3 * 3 + 0] - x0;
1009181f196bSVijay Mahadevan       jacobian[3] = coordinates[1 * 3 + 1] - y0;  jacobian[4] = coordinates[2 * 3 + 1] - y0; jacobian[5] = coordinates[3 * 3 + 1] - y0;
1010181f196bSVijay Mahadevan       jacobian[6] = coordinates[1 * 3 + 2] - z0;  jacobian[7] = coordinates[2 * 3 + 2] - z0; jacobian[8] = coordinates[3 * 3 + 2] - z0;
1011181f196bSVijay Mahadevan     } /* Tetrahedra -- ends */
1012181f196bSVijay Mahadevan     else
1013181f196bSVijay Mahadevan     {
1014181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 3-D entity vertices are invalid. Currently only support HEX8 and TET4 basis evaluations in 3-D : %D", nverts);
1015181f196bSVijay Mahadevan     }
1016181f196bSVijay Mahadevan 
1017181f196bSVijay Mahadevan     if (ijacobian) {
1018181f196bSVijay Mahadevan       /* invert the jacobian */
1019a86ed394SVijay Mahadevan       ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr);
1020181f196bSVijay Mahadevan     }
1021181f196bSVijay Mahadevan 
1022181f196bSVijay Mahadevan   }
1023a86ed394SVijay Mahadevan   if ( volume < 1e-12 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Element has zero volume: %g. Degenerate element or invalid connectivity\n", volume);
1024a86ed394SVijay Mahadevan   if (dvolume) *dvolume = volume;
1025181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1026181f196bSVijay Mahadevan }
1027181f196bSVijay Mahadevan 
1028181f196bSVijay Mahadevan 
1029181f196bSVijay Mahadevan PetscErrorCode FEMComputeBasis_JandF ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quadrature, PetscReal *phypts,
1030181f196bSVijay Mahadevan                                        PetscReal *phibasis, PetscReal *jacobian, PetscReal *ijacobian, PetscReal* volume  )
1031181f196bSVijay Mahadevan {
1032181f196bSVijay Mahadevan   PetscErrorCode  ierr;
1033181f196bSVijay Mahadevan 
1034181f196bSVijay Mahadevan   PetscFunctionBegin;
1035181f196bSVijay Mahadevan   switch (dim) {
1036181f196bSVijay Mahadevan     case 1:
1037181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts,
1038181f196bSVijay Mahadevan             NULL, phibasis, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
1039181f196bSVijay Mahadevan       break;
1040181f196bSVijay Mahadevan     case 2:
1041181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts,
1042181f196bSVijay Mahadevan             NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
1043181f196bSVijay Mahadevan       break;
1044181f196bSVijay Mahadevan     case 3:
1045181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts,
1046181f196bSVijay Mahadevan             NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
1047181f196bSVijay Mahadevan       break;
1048181f196bSVijay Mahadevan     default:
1049181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
1050181f196bSVijay Mahadevan   }
1051181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1052181f196bSVijay Mahadevan }
1053181f196bSVijay Mahadevan 
1054181f196bSVijay Mahadevan 
1055181f196bSVijay Mahadevan 
1056a86ed394SVijay Mahadevan /*@
105797b73a88SSatish Balay   DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the
1058a86ed394SVijay Mahadevan   canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration,
1059a86ed394SVijay Mahadevan   the basis function at the parametric point is also evaluated optionally.
1060a86ed394SVijay Mahadevan 
1061a86ed394SVijay Mahadevan   Input Parameter:
1062*a2b725a8SWilliam Gropp +  PetscInt  dim -         the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
1063*a2b725a8SWilliam Gropp .  PetscInt nverts -       the number of vertices in the physical element
1064*a2b725a8SWilliam Gropp .  PetscReal coordinates - the coordinates of vertices in the physical element
1065*a2b725a8SWilliam Gropp -  PetscReal[3] xphy -     the coordinates of physical point for which natural coordinates (in reference frame) are sought
1066a86ed394SVijay Mahadevan 
1067a86ed394SVijay Mahadevan   Output Parameter:
1068*a2b725a8SWilliam Gropp +  PetscReal[3] natparam - the natural coordinates (in reference frame) corresponding to xphy
1069*a2b725a8SWilliam Gropp -  PetscReal[nverts] phi - the basis functions evaluated at the natural coordinates (natparam)
1070a86ed394SVijay Mahadevan 
1071edc382c3SSatish Balay   Level: advanced
1072edc382c3SSatish Balay 
1073a86ed394SVijay Mahadevan .keywords: DMMoab, Mapping, FEM
1074a86ed394SVijay Mahadevan @*/
1075181f196bSVijay Mahadevan PetscErrorCode DMMoabPToRMapping( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal* xphy, PetscReal* natparam, PetscReal* phi)
1076181f196bSVijay Mahadevan {
1077a86ed394SVijay Mahadevan   /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */
1078181f196bSVijay Mahadevan   const PetscReal tol = 1.0e-10;
1079181f196bSVijay Mahadevan   const PetscInt max_iterations = 10;
1080181f196bSVijay Mahadevan   const PetscReal error_tol_sqr = tol*tol;
1081181f196bSVijay Mahadevan   PetscReal phibasis[8], jacobian[9], ijacobian[9], volume;
1082181f196bSVijay Mahadevan   PetscReal phypts[3] = {0.0, 0.0, 0.0};
1083181f196bSVijay Mahadevan   PetscReal delta[3] = {0.0, 0.0, 0.0};
1084181f196bSVijay Mahadevan   PetscErrorCode  ierr;
1085181f196bSVijay Mahadevan   PetscInt iters=0;
1086181f196bSVijay Mahadevan   PetscReal error=1.0;
1087181f196bSVijay Mahadevan 
1088181f196bSVijay Mahadevan   PetscFunctionBegin;
1089181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
1090181f196bSVijay Mahadevan   PetscValidPointer(xphy, 4);
1091181f196bSVijay Mahadevan   PetscValidPointer(natparam, 5);
1092181f196bSVijay Mahadevan 
1093181f196bSVijay Mahadevan   ierr = PetscMemzero(jacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr);
1094181f196bSVijay Mahadevan   ierr = PetscMemzero(ijacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr);
1095181f196bSVijay Mahadevan   ierr = PetscMemzero(phibasis, nverts * sizeof(PetscReal));CHKERRQ(ierr);
1096181f196bSVijay Mahadevan 
1097a86ed394SVijay Mahadevan   /* zero initial guess */
1098a86ed394SVijay Mahadevan   natparam[0] = natparam[1] = natparam[2] = 0.0;
1099181f196bSVijay Mahadevan 
1100a86ed394SVijay Mahadevan   /* Compute delta = evaluate( xi ) - x */
1101a86ed394SVijay Mahadevan   ierr = FEMComputeBasis_JandF ( dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume );CHKERRQ(ierr);
1102a86ed394SVijay Mahadevan 
1103a86ed394SVijay Mahadevan   error = 0.0;
1104a86ed394SVijay Mahadevan   switch(dim) {
1105a86ed394SVijay Mahadevan     case 3:
1106181f196bSVijay Mahadevan       delta[2] = phypts[2] - xphy[2];
1107a86ed394SVijay Mahadevan       error += (delta[2]*delta[2]);
1108a86ed394SVijay Mahadevan     case 2:
1109a86ed394SVijay Mahadevan       delta[1] = phypts[1] - xphy[1];
1110a86ed394SVijay Mahadevan       error += (delta[1]*delta[1]);
1111a86ed394SVijay Mahadevan     case 1:
1112a86ed394SVijay Mahadevan       delta[0] = phypts[0] - xphy[0];
1113a86ed394SVijay Mahadevan       error += (delta[0]*delta[0]);
1114a86ed394SVijay Mahadevan       break;
1115a86ed394SVijay Mahadevan   }
1116a86ed394SVijay Mahadevan 
1117a86ed394SVijay Mahadevan #if 0
1118a86ed394SVijay Mahadevan   PetscInfo4(NULL, "Current point in physical space : (%g, %g, %g) and error : %g\n", xphy[0], xphy[1], xphy[2], error);
1119a86ed394SVijay Mahadevan #endif
1120181f196bSVijay Mahadevan 
1121181f196bSVijay Mahadevan   while (error > error_tol_sqr) {
1122181f196bSVijay Mahadevan 
1123181f196bSVijay Mahadevan     if(++iters > max_iterations)
1124181f196bSVijay Mahadevan       SETERRQ3(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Maximum Newton iterations (10) reached. Current point in reference space : (%g, %g, %g)", natparam[0], natparam[1], natparam[2]);
1125181f196bSVijay Mahadevan 
1126181f196bSVijay Mahadevan     /* Compute natparam -= J.inverse() * delta */
1127181f196bSVijay Mahadevan     switch(dim) {
1128181f196bSVijay Mahadevan       case 1:
1129181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0];
1130181f196bSVijay Mahadevan         break;
1131181f196bSVijay Mahadevan       case 2:
1132181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1];
1133181f196bSVijay Mahadevan         natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1];
1134181f196bSVijay Mahadevan         break;
1135181f196bSVijay Mahadevan       case 3:
1136181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2];
1137181f196bSVijay Mahadevan         natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2];
1138181f196bSVijay Mahadevan         natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2];
1139181f196bSVijay Mahadevan         break;
1140181f196bSVijay Mahadevan     }
1141181f196bSVijay Mahadevan 
1142181f196bSVijay Mahadevan     /* Compute delta = evaluate( xi ) - x */
1143a86ed394SVijay Mahadevan     ierr = FEMComputeBasis_JandF ( dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume );CHKERRQ(ierr);
1144181f196bSVijay Mahadevan 
1145a86ed394SVijay Mahadevan     error = 0.0;
1146a86ed394SVijay Mahadevan     switch(dim) {
1147a86ed394SVijay Mahadevan       case 3:
1148181f196bSVijay Mahadevan         delta[2] = phypts[2] - xphy[2];
1149a86ed394SVijay Mahadevan         error += (delta[2]*delta[2]);
1150a86ed394SVijay Mahadevan       case 2:
1151a86ed394SVijay Mahadevan         delta[1] = phypts[1] - xphy[1];
1152a86ed394SVijay Mahadevan         error += (delta[1]*delta[1]);
1153a86ed394SVijay Mahadevan       case 1:
1154a86ed394SVijay Mahadevan         delta[0] = phypts[0] - xphy[0];
1155a86ed394SVijay Mahadevan         error += (delta[0]*delta[0]);
1156a86ed394SVijay Mahadevan         break;
1157a86ed394SVijay Mahadevan     }
1158a86ed394SVijay Mahadevan #if 0
1159a86ed394SVijay Mahadevan     PetscInfo5(NULL, "Newton [%D] : Current point in reference space : (%g, %g, %g) and error : %g\n", iters, natparam[0], natparam[1], natparam[2], error);
1160a86ed394SVijay Mahadevan #endif
1161181f196bSVijay Mahadevan   }
1162181f196bSVijay Mahadevan   if (phi) {
1163181f196bSVijay Mahadevan     ierr = PetscMemcpy(phi, phibasis, nverts * sizeof(PetscReal));CHKERRQ(ierr);
1164181f196bSVijay Mahadevan   }
1165181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1166181f196bSVijay Mahadevan }
1167181f196bSVijay Mahadevan 
1168