xref: /petsc/src/dm/impls/moab/dmmbfem.cxx (revision 362febeeeb69b91ebadcb4b2dc0a22cb6dfc4097)
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   PetscReal det = DMatrix_Determinant_2x2_Internal(inmat);
1563d025dbSVijay Mahadevan   if (outmat) {
1663d025dbSVijay Mahadevan     outmat[0] = inmat[3] / det;
1763d025dbSVijay Mahadevan     outmat[1] = -inmat[1] / det;
1863d025dbSVijay Mahadevan     outmat[2] = -inmat[2] / det;
1963d025dbSVijay Mahadevan     outmat[3] = inmat[0] / det;
2063d025dbSVijay Mahadevan   }
2163d025dbSVijay Mahadevan   if (determinant) *determinant = det;
22*362febeeSStefano Zampini   return 0;
2363d025dbSVijay Mahadevan }
2463d025dbSVijay Mahadevan 
25181f196bSVijay Mahadevan static inline PetscReal DMatrix_Determinant_3x3_Internal(const PetscReal inmat[3 * 3])
2663d025dbSVijay Mahadevan {
2763d025dbSVijay Mahadevan   return   inmat[0] * (inmat[8] * inmat[4] - inmat[7] * inmat[5])
2863d025dbSVijay Mahadevan            - inmat[3] * (inmat[8] * inmat[1] - inmat[7] * inmat[2])
2963d025dbSVijay Mahadevan            + inmat[6] * (inmat[5] * inmat[1] - inmat[4] * inmat[2]);
3063d025dbSVijay Mahadevan }
3163d025dbSVijay Mahadevan 
3263d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_3x3_Internal (const PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
3363d025dbSVijay Mahadevan {
34181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_3x3_Internal(inmat);
3563d025dbSVijay Mahadevan   if (outmat) {
3663d025dbSVijay Mahadevan     outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det;
3763d025dbSVijay Mahadevan     outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det;
3863d025dbSVijay Mahadevan     outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det;
3963d025dbSVijay Mahadevan     outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det;
4063d025dbSVijay Mahadevan     outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det;
4163d025dbSVijay Mahadevan     outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det;
4263d025dbSVijay Mahadevan     outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det;
4363d025dbSVijay Mahadevan     outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det;
4463d025dbSVijay Mahadevan     outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det;
4563d025dbSVijay Mahadevan   }
4663d025dbSVijay Mahadevan   if (determinant) *determinant = det;
47*362febeeSStefano Zampini   return 0;
4863d025dbSVijay Mahadevan }
4963d025dbSVijay Mahadevan 
50181f196bSVijay Mahadevan inline PetscReal DMatrix_Determinant_4x4_Internal(PetscReal inmat[4 * 4])
51181f196bSVijay Mahadevan {
52181f196bSVijay Mahadevan   return
53181f196bSVijay Mahadevan     inmat[0 + 0 * 4] * (
54181f196bSVijay Mahadevan       inmat[1 + 1 * 4] * (inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4])
55181f196bSVijay Mahadevan       - inmat[1 + 2 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4])
56181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4]))
57181f196bSVijay Mahadevan     - inmat[0 + 1 * 4] * (
58181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * (inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4])
59181f196bSVijay Mahadevan       - inmat[1 + 2 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4])
60181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4]))
61181f196bSVijay Mahadevan     + inmat[0 + 2 * 4] * (
62181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4])
63181f196bSVijay Mahadevan       - inmat[1 + 1 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4])
64181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4]))
65181f196bSVijay Mahadevan     - inmat[0 + 3 * 4] * (
66181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4])
67181f196bSVijay Mahadevan       - inmat[1 + 1 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4])
68181f196bSVijay Mahadevan       + inmat[1 + 2 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4]));
69181f196bSVijay Mahadevan }
70181f196bSVijay Mahadevan 
71181f196bSVijay Mahadevan inline PetscErrorCode DMatrix_Invert_4x4_Internal (PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
72181f196bSVijay Mahadevan {
73181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_4x4_Internal(inmat);
74181f196bSVijay Mahadevan   if (outmat) {
75181f196bSVijay 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;
76181f196bSVijay 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;
77181f196bSVijay 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;
78181f196bSVijay 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;
79181f196bSVijay 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;
80181f196bSVijay 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;
81181f196bSVijay 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;
82181f196bSVijay 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;
83181f196bSVijay 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;
84181f196bSVijay 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;
85181f196bSVijay 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;
86181f196bSVijay 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;
87181f196bSVijay 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;
88181f196bSVijay 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;
89181f196bSVijay 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;
90181f196bSVijay 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;
91181f196bSVijay Mahadevan   }
92181f196bSVijay Mahadevan   if (determinant) *determinant = det;
93*362febeeSStefano Zampini   return 0;
94181f196bSVijay Mahadevan }
95181f196bSVijay Mahadevan 
96cab5ea25SPierre Jolivet /*@C
9797b73a88SSatish Balay   Compute_Lagrange_Basis_1D_Internal - Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element.
9863d025dbSVijay Mahadevan 
9963d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a linear or quadratic edge element.
10063d025dbSVijay Mahadevan 
10163d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
10263d025dbSVijay Mahadevan   and their derivatives with respect to X.
10363d025dbSVijay Mahadevan 
10463d025dbSVijay Mahadevan   Notes:
10563d025dbSVijay Mahadevan 
10663d025dbSVijay Mahadevan   Example Physical Element
107a2b725a8SWilliam Gropp .vb
10863d025dbSVijay Mahadevan     1-------2        1----3----2
10963d025dbSVijay Mahadevan       EDGE2             EDGE3
110a2b725a8SWilliam Gropp .ve
11163d025dbSVijay Mahadevan 
11263d025dbSVijay Mahadevan   Input Parameter:
113a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
114a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
115a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
116a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
11763d025dbSVijay Mahadevan 
11863d025dbSVijay Mahadevan   Output Parameter:
119a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
120a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
121a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
122a2b725a8SWilliam Gropp -  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
12363d025dbSVijay Mahadevan 
124edc382c3SSatish Balay   Level: advanced
125edc382c3SSatish Balay 
12663d025dbSVijay Mahadevan @*/
12763d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_1D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
128181f196bSVijay Mahadevan                                                   PetscReal *jxw, PetscReal *phi, PetscReal *dphidx,
129181f196bSVijay Mahadevan                                                   PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
13063d025dbSVijay Mahadevan {
13163d025dbSVijay Mahadevan   int             i, j;
13263d025dbSVijay Mahadevan   PetscErrorCode  ierr;
13363d025dbSVijay Mahadevan 
134181f196bSVijay Mahadevan   PetscFunctionBegin;
135181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 9);
136181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 10);
137181f196bSVijay Mahadevan   PetscValidPointer(volume, 11);
138181f196bSVijay Mahadevan   if (phypts) {
139580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr);
140181f196bSVijay Mahadevan   }
14163d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
142580bdb30SBarry Smith     ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr);
14363d025dbSVijay Mahadevan   }
14463d025dbSVijay Mahadevan   if (nverts == 2) { /* Linear Edge */
14563d025dbSVijay Mahadevan 
1462da392ccSBarry Smith     for (j = 0; j < npts; j++) {
147a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
148181f196bSVijay Mahadevan       const PetscReal r = quad[j];
14963d025dbSVijay Mahadevan 
150181f196bSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r);
151181f196bSVijay Mahadevan       phi[1 + offset] = (       r);
15263d025dbSVijay Mahadevan 
153181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 1.0 };
15463d025dbSVijay Mahadevan 
155181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
15663d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
157181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
158181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
159181f196bSVijay Mahadevan         if (phypts) {
160181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
161181f196bSVijay Mahadevan         }
16263d025dbSVijay Mahadevan       }
16363d025dbSVijay Mahadevan 
16463d025dbSVijay Mahadevan       /* invert the jacobian */
165181f196bSVijay Mahadevan       *volume = jacobian[0];
166181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
167181f196bSVijay Mahadevan       jxw[j] *= *volume;
16863d025dbSVijay Mahadevan 
16963d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
17063d025dbSVijay Mahadevan       for (i = 0; i < nverts; i++) {
171181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
17263d025dbSVijay Mahadevan       }
17363d025dbSVijay Mahadevan     }
1742da392ccSBarry Smith   } else if (nverts == 3) { /* Quadratic Edge */
17563d025dbSVijay Mahadevan 
1762da392ccSBarry Smith     for (j = 0; j < npts; j++) {
177a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
178181f196bSVijay Mahadevan       const PetscReal r = quad[j];
17963d025dbSVijay Mahadevan 
180181f196bSVijay Mahadevan       phi[0 + offset] = 1.0 + r * ( 2.0 * r - 3.0);
181181f196bSVijay Mahadevan       phi[1 + offset] = 4.0 * r * ( 1.0 - r);
182181f196bSVijay Mahadevan       phi[2 + offset] = r * ( 2.0 * r - 1.0);
18363d025dbSVijay Mahadevan 
184181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r), 4.0 * r - 1.0};
18563d025dbSVijay Mahadevan 
186181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
18763d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
188181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
189181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
190181f196bSVijay Mahadevan         if (phypts) {
191181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
192181f196bSVijay Mahadevan         }
19363d025dbSVijay Mahadevan       }
19463d025dbSVijay Mahadevan 
19563d025dbSVijay Mahadevan       /* invert the jacobian */
196181f196bSVijay Mahadevan       *volume = jacobian[0];
197181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
198181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
19963d025dbSVijay Mahadevan 
20063d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
20163d025dbSVijay Mahadevan       for (i = 0; i < nverts; i++) {
202181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
20363d025dbSVijay Mahadevan       }
20463d025dbSVijay Mahadevan     }
2052da392ccSBarry Smith   } else 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);
20663d025dbSVijay Mahadevan   PetscFunctionReturn(0);
20763d025dbSVijay Mahadevan }
20863d025dbSVijay Mahadevan 
209cab5ea25SPierre Jolivet /*@C
21097b73a88SSatish Balay   Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.
21163d025dbSVijay Mahadevan 
21263d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a quadrangle or triangle.
21363d025dbSVijay Mahadevan 
21463d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
21563d025dbSVijay Mahadevan   and their derivatives with respect to X and Y.
21663d025dbSVijay Mahadevan 
21763d025dbSVijay Mahadevan   Notes:
21863d025dbSVijay Mahadevan 
21963d025dbSVijay Mahadevan   Example Physical Element (QUAD4)
220a2b725a8SWilliam Gropp .vb
22163d025dbSVijay Mahadevan     4------3        s
22263d025dbSVijay Mahadevan     |      |        |
22363d025dbSVijay Mahadevan     |      |        |
22463d025dbSVijay Mahadevan     |      |        |
22563d025dbSVijay Mahadevan     1------2        0-------r
226a2b725a8SWilliam Gropp .ve
22763d025dbSVijay Mahadevan 
22863d025dbSVijay Mahadevan   Input Parameter:
229a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
230a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
231a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
232a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
23363d025dbSVijay Mahadevan 
23463d025dbSVijay Mahadevan   Output Parameter:
235a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
236a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
237a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
238a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
239a2b725a8SWilliam Gropp -  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
24063d025dbSVijay Mahadevan 
241edc382c3SSatish Balay   Level: advanced
242edc382c3SSatish Balay 
24363d025dbSVijay Mahadevan @*/
24463d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_2D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
245181f196bSVijay Mahadevan                                                   PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy,
246181f196bSVijay Mahadevan                                                   PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
24763d025dbSVijay Mahadevan {
248a86ed394SVijay Mahadevan   PetscInt       i, j, k;
24963d025dbSVijay Mahadevan   PetscErrorCode ierr;
25063d025dbSVijay Mahadevan 
25163d025dbSVijay Mahadevan   PetscFunctionBegin;
252181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 10);
253181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 11);
254181f196bSVijay Mahadevan   PetscValidPointer(volume, 12);
255580bdb30SBarry Smith   ierr = PetscArrayzero(phi, npts);CHKERRQ(ierr);
256181f196bSVijay Mahadevan   if (phypts) {
257580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr);
258181f196bSVijay Mahadevan   }
25963d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
260580bdb30SBarry Smith     ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr);
261580bdb30SBarry Smith     ierr = PetscArrayzero(dphidy, npts * nverts);CHKERRQ(ierr);
26263d025dbSVijay Mahadevan   }
26363d025dbSVijay Mahadevan   if (nverts == 4) { /* Linear Quadrangle */
26463d025dbSVijay Mahadevan 
26563d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
26663d025dbSVijay Mahadevan     {
267a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
268181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
269181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
27063d025dbSVijay Mahadevan 
27163d025dbSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r) * ( 1.0 - s);
27263d025dbSVijay Mahadevan       phi[1 + offset] =         r   * ( 1.0 - s);
27363d025dbSVijay Mahadevan       phi[2 + offset] =         r   *         s;
27463d025dbSVijay Mahadevan       phi[3 + offset] = ( 1.0 - r) *         s;
27563d025dbSVijay Mahadevan 
276181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
277181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
27863d025dbSVijay Mahadevan 
279580bdb30SBarry Smith       ierr = PetscArrayzero(jacobian, 4);CHKERRQ(ierr);
280580bdb30SBarry Smith       ierr = PetscArrayzero(ijacobian, 4);CHKERRQ(ierr);
28163d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
282181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
28363d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
28463d025dbSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
28563d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
28663d025dbSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
287181f196bSVijay Mahadevan         if (phypts) {
288181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
289181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
290181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
291181f196bSVijay Mahadevan         }
29263d025dbSVijay Mahadevan       }
29363d025dbSVijay Mahadevan 
29463d025dbSVijay Mahadevan       /* invert the jacobian */
295181f196bSVijay Mahadevan       ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
296181f196bSVijay 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);
29763d025dbSVijay Mahadevan 
298181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
29963d025dbSVijay Mahadevan 
300181f196bSVijay Mahadevan       /*  Let us compute the bases derivatives by scaling with inverse jacobian. */
301181f196bSVijay Mahadevan       if (dphidx) {
30263d025dbSVijay Mahadevan         for (i = 0; i < nverts; i++) {
303a86ed394SVijay Mahadevan           for (k = 0; k < 2; ++k) {
30463d025dbSVijay Mahadevan             if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0];
30563d025dbSVijay Mahadevan             if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1];
306181f196bSVijay Mahadevan           } /* for k=[0..2) */
307181f196bSVijay Mahadevan         } /* for i=[0..nverts) */
308181f196bSVijay Mahadevan       } /* if (dphidx) */
30963d025dbSVijay Mahadevan     }
3102da392ccSBarry Smith   } else if (nverts == 3) { /* Linear triangle */
3112da392ccSBarry Smith     const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1];
31263d025dbSVijay Mahadevan 
313580bdb30SBarry Smith     ierr = PetscArrayzero(jacobian, 4);CHKERRQ(ierr);
314580bdb30SBarry Smith     ierr = PetscArrayzero(ijacobian, 4);CHKERRQ(ierr);
31563d025dbSVijay Mahadevan 
31663d025dbSVijay Mahadevan     /* Jacobian is constant */
317181f196bSVijay Mahadevan     jacobian[0] = (coords[0 * 3 + 0] - x2); jacobian[1] = (coords[1 * 3 + 0] - x2);
318181f196bSVijay Mahadevan     jacobian[2] = (coords[0 * 3 + 1] - y2); jacobian[3] = (coords[1 * 3 + 1] - y2);
31963d025dbSVijay Mahadevan 
32063d025dbSVijay Mahadevan     /* invert the jacobian */
321181f196bSVijay Mahadevan     ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
3222da392ccSBarry Smith     if (*volume < 1e-12) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangular element has zero volume: %g. Degenerate element or invalid connectivity\n", (double)*volume);
323181f196bSVijay Mahadevan 
324181f196bSVijay Mahadevan     const PetscReal Dx[3] = { ijacobian[0], ijacobian[2], - ijacobian[0] - ijacobian[2] };
325181f196bSVijay Mahadevan     const PetscReal Dy[3] = { ijacobian[1], ijacobian[3], - ijacobian[1] - ijacobian[3] };
32663d025dbSVijay Mahadevan 
3272da392ccSBarry Smith     for (j = 0; j < npts; j++) {
328a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
329181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
330181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
33163d025dbSVijay Mahadevan 
332181f196bSVijay Mahadevan       if (jxw) jxw[j] *= 0.5;
33363d025dbSVijay Mahadevan 
334181f196bSVijay Mahadevan       /* const PetscReal Ni[3]  = { r, s, 1.0 - r - s }; */
335181f196bSVijay Mahadevan       const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s;
336181f196bSVijay Mahadevan       const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s;
337181f196bSVijay Mahadevan       if (phypts) {
338181f196bSVijay Mahadevan         phypts[offset + 0] = phipts_x;
339181f196bSVijay Mahadevan         phypts[offset + 1] = phipts_y;
340181f196bSVijay Mahadevan       }
34163d025dbSVijay Mahadevan 
342110fc3b0SBarry Smith       /* \phi_0 = (b.y - c.y) x + (b.x - c.x) y + c.x b.y - b.x c.y */
343181f196bSVijay Mahadevan       phi[0 + offset] = (  ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2));
344110fc3b0SBarry Smith       /* \phi_1 = (c.y - a.y) x + (a.x - c.x) y + c.x a.y - a.x c.y */
345181f196bSVijay Mahadevan       phi[1 + offset] = (  ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2));
34663d025dbSVijay Mahadevan       phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset];
34763d025dbSVijay Mahadevan 
34863d025dbSVijay Mahadevan       if (dphidx) {
349181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
350181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
351181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
35263d025dbSVijay Mahadevan       }
35363d025dbSVijay Mahadevan 
35463d025dbSVijay Mahadevan       if (dphidy) {
355181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
356181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
357181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
35863d025dbSVijay Mahadevan       }
35963d025dbSVijay Mahadevan 
36063d025dbSVijay Mahadevan     }
3612da392ccSBarry Smith   } else 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);
36263d025dbSVijay Mahadevan   PetscFunctionReturn(0);
36363d025dbSVijay Mahadevan }
36463d025dbSVijay Mahadevan 
365cab5ea25SPierre Jolivet /*@C
36697b73a88SSatish Balay   Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element.
36763d025dbSVijay Mahadevan 
36863d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a hexahedra or tetrahedra.
36963d025dbSVijay Mahadevan 
37063d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
37163d025dbSVijay Mahadevan   and their derivatives with respect to X, Y and Z.
37263d025dbSVijay Mahadevan 
37363d025dbSVijay Mahadevan   Notes:
37463d025dbSVijay Mahadevan 
37563d025dbSVijay Mahadevan   Example Physical Element (HEX8)
376a2b725a8SWilliam Gropp .vb
37763d025dbSVijay Mahadevan       8------7
37863d025dbSVijay Mahadevan      /|     /|        t  s
37963d025dbSVijay Mahadevan     5------6 |        | /
38063d025dbSVijay Mahadevan     | |    | |        |/
38163d025dbSVijay Mahadevan     | 4----|-3        0-------r
38263d025dbSVijay Mahadevan     |/     |/
38363d025dbSVijay Mahadevan     1------2
384a2b725a8SWilliam Gropp .ve
38563d025dbSVijay Mahadevan 
38663d025dbSVijay Mahadevan   Input Parameter:
387a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
388a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
389a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
390a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
39163d025dbSVijay Mahadevan 
39263d025dbSVijay Mahadevan   Output Parameter:
393a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
394a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
395a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
396a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
397a2b725a8SWilliam Gropp .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
398a2b725a8SWilliam Gropp -  PetscReal dphidz[npts] -    the derivative of the bases wrt Z-direction evaluated at the specified quadrature points
39963d025dbSVijay Mahadevan 
400edc382c3SSatish Balay   Level: advanced
401edc382c3SSatish Balay 
40263d025dbSVijay Mahadevan @*/
40363d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_3D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
404181f196bSVijay Mahadevan                                                   PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz,
405181f196bSVijay Mahadevan                                                   PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
40663d025dbSVijay Mahadevan {
407a86ed394SVijay Mahadevan   PetscInt       i, j, k;
40863d025dbSVijay Mahadevan   PetscErrorCode ierr;
40963d025dbSVijay Mahadevan 
41063d025dbSVijay Mahadevan   PetscFunctionBegin;
411181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 11);
412181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 12);
413181f196bSVijay Mahadevan   PetscValidPointer(volume, 13);
4142da392ccSBarry Smith 
415580bdb30SBarry Smith   ierr = PetscArrayzero(phi, npts);CHKERRQ(ierr);
416181f196bSVijay Mahadevan   if (phypts) {
417580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr);
418181f196bSVijay Mahadevan   }
41963d025dbSVijay Mahadevan   if (dphidx) {
420580bdb30SBarry Smith     ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr);
421580bdb30SBarry Smith     ierr = PetscArrayzero(dphidy, npts * nverts);CHKERRQ(ierr);
422580bdb30SBarry Smith     ierr = PetscArrayzero(dphidz, npts * nverts);CHKERRQ(ierr);
42363d025dbSVijay Mahadevan   }
42463d025dbSVijay Mahadevan 
42563d025dbSVijay Mahadevan   if (nverts == 8) { /* Linear Hexahedra */
42663d025dbSVijay Mahadevan 
4272da392ccSBarry Smith     for (j = 0; j < npts; j++) {
428a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
429181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
430181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
431181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
43263d025dbSVijay Mahadevan 
433a86ed394SVijay Mahadevan       phi[offset + 0] = ( 1.0 - r) * ( 1.0 - s) * ( 1.0 - t); /* 0,0,0 */
434a86ed394SVijay Mahadevan       phi[offset + 1] = (       r) * ( 1.0 - s) * ( 1.0 - t); /* 1,0,0 */
435a86ed394SVijay Mahadevan       phi[offset + 2] = (       r) * (       s) * ( 1.0 - t); /* 1,1,0 */
436a86ed394SVijay Mahadevan       phi[offset + 3] = ( 1.0 - r) * (       s) * ( 1.0 - t); /* 0,1,0 */
437a86ed394SVijay Mahadevan       phi[offset + 4] = ( 1.0 - r) * ( 1.0 - s) * (       t); /* 0,0,1 */
438a86ed394SVijay Mahadevan       phi[offset + 5] = (       r) * ( 1.0 - s) * (       t); /* 1,0,1 */
439a86ed394SVijay Mahadevan       phi[offset + 6] = (       r) * (       s) * (       t); /* 1,1,1 */
440a86ed394SVijay Mahadevan       phi[offset + 7] = ( 1.0 - r) * (       s) * (       t); /* 0,1,1 */
44163d025dbSVijay Mahadevan 
442181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s) * ( 1.0 - t),
44363d025dbSVijay Mahadevan                                        ( 1.0 - s) * ( 1.0 - t),
444a86ed394SVijay Mahadevan                                        (       s) * ( 1.0 - t),
445a86ed394SVijay Mahadevan                                      - (       s) * ( 1.0 - t),
446a86ed394SVijay Mahadevan                                      - ( 1.0 - s) * (       t),
447a86ed394SVijay Mahadevan                                        ( 1.0 - s) * (       t),
448a86ed394SVijay Mahadevan                                        (       s) * (       t),
449a86ed394SVijay Mahadevan                                      - (       s) * (       t)
45063d025dbSVijay Mahadevan                                     };
45163d025dbSVijay Mahadevan 
452181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r) * ( 1.0 - t),
453a86ed394SVijay Mahadevan                                        - (       r) * ( 1.0 - t),
454a86ed394SVijay Mahadevan                                          (       r) * ( 1.0 - t),
45563d025dbSVijay Mahadevan                                          ( 1.0 - r) * ( 1.0 - t),
456a86ed394SVijay Mahadevan                                        - ( 1.0 - r) * (       t),
457a86ed394SVijay Mahadevan                                        - (       r) * (       t),
458a86ed394SVijay Mahadevan                                          (       r) * (       t),
459a86ed394SVijay Mahadevan                                          ( 1.0 - r) * (       t)
46063d025dbSVijay Mahadevan                                       };
46163d025dbSVijay Mahadevan 
462181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r) * ( 1.0 - s),
463a86ed394SVijay Mahadevan                                         - (       r) * ( 1.0 - s),
464a86ed394SVijay Mahadevan                                         - (       r) * (       s),
465a86ed394SVijay Mahadevan                                         - ( 1.0 - r) * (       s),
46663d025dbSVijay Mahadevan                                           ( 1.0 - r) * ( 1.0 - s),
467a86ed394SVijay Mahadevan                                           (       r) * ( 1.0 - s),
468a86ed394SVijay Mahadevan                                           (       r) * (       s),
469a86ed394SVijay Mahadevan                                           ( 1.0 - r) * (       s)
47063d025dbSVijay Mahadevan                                      };
47163d025dbSVijay Mahadevan 
472580bdb30SBarry Smith       ierr = PetscArrayzero(jacobian, 9);CHKERRQ(ierr);
473580bdb30SBarry Smith       ierr = PetscArrayzero(ijacobian, 9);CHKERRQ(ierr);
47463d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
475181f196bSVijay Mahadevan         const PetscReal* vertex = coords + i * 3;
47663d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
47763d025dbSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
47863d025dbSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
47963d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
48063d025dbSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
48163d025dbSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
48263d025dbSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
48363d025dbSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
48463d025dbSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
485181f196bSVijay Mahadevan         if (phypts) {
486181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertex[0];
487181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertex[1];
488181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertex[2];
489181f196bSVijay Mahadevan         }
49063d025dbSVijay Mahadevan       }
49163d025dbSVijay Mahadevan 
49263d025dbSVijay Mahadevan       /* invert the jacobian */
493181f196bSVijay Mahadevan       ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
494a86ed394SVijay 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);
49563d025dbSVijay Mahadevan 
496a86ed394SVijay Mahadevan       if (jxw) jxw[j] *= (*volume);
49763d025dbSVijay Mahadevan 
49863d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
49963d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
500a86ed394SVijay Mahadevan         for (k = 0; k < 3; ++k) {
50163d025dbSVijay Mahadevan           if (dphidx) dphidx[i + offset] += dNi_dxi[i]   * ijacobian[0 * 3 + k];
50263d025dbSVijay Mahadevan           if (dphidy) dphidy[i + offset] += dNi_deta[i]  * ijacobian[1 * 3 + k];
50363d025dbSVijay Mahadevan           if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k];
50463d025dbSVijay Mahadevan         }
50563d025dbSVijay Mahadevan       }
50663d025dbSVijay Mahadevan     }
5072da392ccSBarry Smith   } else if (nverts == 4) { /* Linear Tetrahedra */
5082da392ccSBarry Smith     PetscReal       Dx[4]={0,0,0,0}, Dy[4]={0,0,0,0}, Dz[4]={0,0,0,0};
5092da392ccSBarry Smith     const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2];
51063d025dbSVijay Mahadevan 
511580bdb30SBarry Smith     ierr = PetscArrayzero(jacobian, 9);CHKERRQ(ierr);
512580bdb30SBarry Smith     ierr = PetscArrayzero(ijacobian, 9);CHKERRQ(ierr);
51363d025dbSVijay Mahadevan 
514181f196bSVijay Mahadevan     /* compute the jacobian */
515181f196bSVijay Mahadevan     jacobian[0] = coords[1 * 3 + 0] - x0;  jacobian[1] = coords[2 * 3 + 0] - x0; jacobian[2] = coords[3 * 3 + 0] - x0;
516181f196bSVijay Mahadevan     jacobian[3] = coords[1 * 3 + 1] - y0;  jacobian[4] = coords[2 * 3 + 1] - y0; jacobian[5] = coords[3 * 3 + 1] - y0;
517181f196bSVijay Mahadevan     jacobian[6] = coords[1 * 3 + 2] - z0;  jacobian[7] = coords[2 * 3 + 2] - z0; jacobian[8] = coords[3 * 3 + 2] - z0;
51863d025dbSVijay Mahadevan 
51963d025dbSVijay Mahadevan     /* invert the jacobian */
520181f196bSVijay Mahadevan     ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
5212da392ccSBarry Smith     if (*volume < 1e-8) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral element has zero volume: %g. Degenerate element or invalid connectivity\n", (double)*volume);
52263d025dbSVijay Mahadevan 
523181f196bSVijay Mahadevan     if (dphidx) {
524181f196bSVijay Mahadevan       Dx[0] =   ( coords[1 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
525181f196bSVijay Mahadevan                  - coords[1 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
5262da392ccSBarry Smith                  - coords[1 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume;
527181f196bSVijay Mahadevan       Dx[1] = - ( coords[1 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
528181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
5292da392ccSBarry Smith                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume;
530181f196bSVijay Mahadevan       Dx[2] =   ( coords[1 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
531181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
5322da392ccSBarry Smith                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume;
533181f196bSVijay Mahadevan       Dx[3] =  - ( Dx[0] + Dx[1] + Dx[2]);
534181f196bSVijay Mahadevan       Dy[0] =   ( coords[0 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
535181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
5362da392ccSBarry Smith                  + coords[0 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume;
537181f196bSVijay Mahadevan       Dy[1] = - ( coords[0 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
538181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
5392da392ccSBarry Smith                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume;
540181f196bSVijay Mahadevan       Dy[2] =   ( coords[0 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
541181f196bSVijay Mahadevan                  - coords[0 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
5422da392ccSBarry Smith                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume;
543181f196bSVijay Mahadevan       Dy[3] =  - ( Dy[0] + Dy[1] + Dy[2]);
544181f196bSVijay Mahadevan       Dz[0] =   ( coords[0 + 1 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3])
545181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3])
5462da392ccSBarry Smith                  + coords[0 + 3 * 3] * (coords[1 + 2 * 3] - coords[1 + 1 * 3])) / *volume;
547181f196bSVijay Mahadevan       Dz[1] = - ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3])
548181f196bSVijay Mahadevan                   + coords[0 + 2 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3])
5492da392ccSBarry Smith                   - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 2 * 3])) / *volume;
550181f196bSVijay Mahadevan       Dz[2] =   ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3])
551181f196bSVijay Mahadevan                  + coords[0 + 1 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3])
5522da392ccSBarry Smith                  - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 1 * 3])) / *volume;
553181f196bSVijay Mahadevan       Dz[3] =  - ( Dz[0] + Dz[1] + Dz[2]);
554181f196bSVijay Mahadevan     }
55563d025dbSVijay Mahadevan 
5562da392ccSBarry Smith     for (j = 0; j < npts; j++) {
557a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
558181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
559181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
560181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
56163d025dbSVijay Mahadevan 
562181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
56363d025dbSVijay Mahadevan 
56463d025dbSVijay Mahadevan       phi[offset + 0] = 1.0 - r - s - t;
56563d025dbSVijay Mahadevan       phi[offset + 1] = r;
56663d025dbSVijay Mahadevan       phi[offset + 2] = s;
56763d025dbSVijay Mahadevan       phi[offset + 3] = t;
56863d025dbSVijay Mahadevan 
569181f196bSVijay Mahadevan       if (phypts) {
570181f196bSVijay Mahadevan         for (i = 0; i < nverts; ++i) {
571181f196bSVijay Mahadevan           const PetscScalar* vertices = coords + i * 3;
572181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
573181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
574181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
575181f196bSVijay Mahadevan         }
576181f196bSVijay Mahadevan       }
577181f196bSVijay Mahadevan 
578181f196bSVijay Mahadevan       /* Now set the derivatives */
57963d025dbSVijay Mahadevan       if (dphidx) {
580181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
581181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
582181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
583181f196bSVijay Mahadevan         dphidx[3 + offset] = Dx[3];
58463d025dbSVijay Mahadevan 
585181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
586181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
587181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
588181f196bSVijay Mahadevan         dphidy[3 + offset] = Dy[3];
58963d025dbSVijay Mahadevan 
590181f196bSVijay Mahadevan         dphidz[0 + offset] = Dz[0];
591181f196bSVijay Mahadevan         dphidz[1 + offset] = Dz[1];
592181f196bSVijay Mahadevan         dphidz[2 + offset] = Dz[2];
593181f196bSVijay Mahadevan         dphidz[3 + offset] = Dz[3];
59463d025dbSVijay Mahadevan       }
59563d025dbSVijay Mahadevan 
59663d025dbSVijay Mahadevan     } /* Tetrahedra -- ends */
5972da392ccSBarry Smith   } else 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);
59863d025dbSVijay Mahadevan   PetscFunctionReturn(0);
59963d025dbSVijay Mahadevan }
60063d025dbSVijay Mahadevan 
601cab5ea25SPierre Jolivet /*@C
60297b73a88SSatish Balay   DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
60363d025dbSVijay Mahadevan 
60463d025dbSVijay Mahadevan   The routine takes the coordinates of the vertices of an element and computes basis functions associated with
60563d025dbSVijay Mahadevan   each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
60663d025dbSVijay Mahadevan 
60763d025dbSVijay Mahadevan   Input Parameter:
608a2b725a8SWilliam Gropp +  PetscInt  nverts,           the number of element vertices
60963d025dbSVijay Mahadevan .  PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering)
61063d025dbSVijay Mahadevan .  PetscInt  npts,             the number of evaluation points (quadrature points)
611a2b725a8SWilliam Gropp -  PetscReal quad[3*npts],     the evaluation points (quadrature points) in the reference space
61263d025dbSVijay Mahadevan 
61363d025dbSVijay Mahadevan   Output Parameter:
614a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts],   the evaluation points (quadrature points) transformed to the physical space
61563d025dbSVijay Mahadevan .  PetscReal jxw[npts],        the jacobian determinant * quadrature weight necessary for assembling discrete contributions
61663d025dbSVijay Mahadevan .  PetscReal fe_basis[npts],        the bases values evaluated at the specified quadrature points
617a2b725a8SWilliam 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
61863d025dbSVijay Mahadevan 
619edc382c3SSatish Balay   Level: advanced
620edc382c3SSatish Balay 
62163d025dbSVijay Mahadevan @*/
622181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMComputeBasis(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature,
623181f196bSVijay Mahadevan                                      PetscReal *phypts, PetscReal *jacobian_quadrature_weight_product,
624181f196bSVijay Mahadevan                                      PetscReal *fe_basis, PetscReal **fe_basis_derivatives)
62563d025dbSVijay Mahadevan {
62663d025dbSVijay Mahadevan   PetscErrorCode  ierr;
627b21a1e88SVijay Mahadevan   PetscInt        npoints,idim;
62863d025dbSVijay Mahadevan   bool            compute_der;
62963d025dbSVijay Mahadevan   const PetscReal *quadpts, *quadwts;
630181f196bSVijay Mahadevan   PetscReal       jacobian[9], ijacobian[9], volume;
63163d025dbSVijay Mahadevan 
63263d025dbSVijay Mahadevan   PetscFunctionBegin;
63363d025dbSVijay Mahadevan   PetscValidPointer(coordinates, 3);
63478dc7ee3SMatthew G. Knepley   PetscValidHeaderSpecific(quadrature, PETSCQUADRATURE_CLASSID, 4);
63563d025dbSVijay Mahadevan   PetscValidPointer(fe_basis, 7);
63663d025dbSVijay Mahadevan   compute_der = (fe_basis_derivatives != NULL);
63763d025dbSVijay Mahadevan 
63863d025dbSVijay Mahadevan   /* Get the quadrature points and weights for the given quadrature rule */
639b21a1e88SVijay Mahadevan   ierr = PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts);CHKERRQ(ierr);
640b21a1e88SVijay Mahadevan   if (idim != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%D) vs quadrature (%D)\n",idim,dim);
641181f196bSVijay Mahadevan   if (jacobian_quadrature_weight_product) {
642580bdb30SBarry Smith     ierr = PetscArraycpy(jacobian_quadrature_weight_product, quadwts, npoints);CHKERRQ(ierr);
643181f196bSVijay Mahadevan   }
64463d025dbSVijay Mahadevan 
64563d025dbSVijay Mahadevan   switch (dim) {
64663d025dbSVijay Mahadevan   case 1:
64763d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, npoints, quadpts, phypts,
64863d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
649181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
650181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
65163d025dbSVijay Mahadevan     break;
65263d025dbSVijay Mahadevan   case 2:
65363d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, npoints, quadpts, phypts,
65463d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
65563d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
656181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[1] : NULL),
657181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
65863d025dbSVijay Mahadevan     break;
65963d025dbSVijay Mahadevan   case 3:
66063d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, npoints, quadpts, phypts,
66163d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
66263d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
66363d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[1] : NULL),
664181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[2] : NULL),
665181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
66663d025dbSVijay Mahadevan     break;
66763d025dbSVijay Mahadevan   default:
66863d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
66963d025dbSVijay Mahadevan   }
67063d025dbSVijay Mahadevan   PetscFunctionReturn(0);
67163d025dbSVijay Mahadevan }
67263d025dbSVijay Mahadevan 
673cab5ea25SPierre Jolivet /*@C
67497b73a88SSatish Balay   DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given
67563d025dbSVijay Mahadevan   dimension and polynomial order (deciphered from number of element vertices).
67663d025dbSVijay Mahadevan 
67763d025dbSVijay Mahadevan   Input Parameter:
67863d025dbSVijay Mahadevan 
679a2b725a8SWilliam Gropp +  PetscInt  dim   -   the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
680a2b725a8SWilliam Gropp -  PetscInt nverts -   the number of vertices in the physical element
68163d025dbSVijay Mahadevan 
68263d025dbSVijay Mahadevan   Output Parameter:
683a2b725a8SWilliam Gropp .  PetscQuadrature quadrature -  the quadrature object with default settings to integrate polynomials defined over the element
68463d025dbSVijay Mahadevan 
685edc382c3SSatish Balay   Level: advanced
686edc382c3SSatish Balay 
68763d025dbSVijay Mahadevan @*/
688181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature)
68963d025dbSVijay Mahadevan {
69063d025dbSVijay Mahadevan   PetscReal       *w, *x;
691b21a1e88SVijay Mahadevan   PetscInt        nc=1;
69263d025dbSVijay Mahadevan   PetscErrorCode  ierr;
69363d025dbSVijay Mahadevan 
69463d025dbSVijay Mahadevan   PetscFunctionBegin;
69563d025dbSVijay Mahadevan   /* Create an appropriate quadrature rule to sample basis */
69663d025dbSVijay Mahadevan   switch (dim)
69763d025dbSVijay Mahadevan   {
69863d025dbSVijay Mahadevan   case 1:
69963d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
700e6a796c3SToby Isaac     ierr = PetscDTStroudConicalQuadrature(1, nc, nverts, 0, 1.0, quadrature);CHKERRQ(ierr);
70163d025dbSVijay Mahadevan     break;
70263d025dbSVijay Mahadevan   case 2:
70363d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
70463d025dbSVijay Mahadevan     if (nverts == 3) {
705a86ed394SVijay Mahadevan       const PetscInt order = 2;
706a86ed394SVijay Mahadevan       const PetscInt npoints = (order == 2 ? 3 : 6);
70763d025dbSVijay Mahadevan       ierr = PetscMalloc2(npoints * 2, &x, npoints, &w);CHKERRQ(ierr);
708181f196bSVijay Mahadevan       if (npoints == 3) {
70963d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0;
71063d025dbSVijay Mahadevan         x[3] = x[4] = 2.0 / 3.0;
71163d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 1.0 / 3.0;
7122da392ccSBarry Smith       } else if (npoints == 6) {
71363d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = 0.44594849091597;
71463d025dbSVijay Mahadevan         x[3] = x[4] = 0.10810301816807;
71563d025dbSVijay Mahadevan         x[5] = 0.44594849091597;
71663d025dbSVijay Mahadevan         x[6] = x[7] = x[8] = 0.09157621350977;
71763d025dbSVijay Mahadevan         x[9] = x[10] = 0.81684757298046;
71863d025dbSVijay Mahadevan         x[11] = 0.09157621350977;
71963d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 0.22338158967801;
720181f196bSVijay Mahadevan         w[3] = w[4] = w[5] = 0.10995174365532;
7212da392ccSBarry Smith       } else SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %D", npoints);
72263d025dbSVijay Mahadevan       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr);
72363d025dbSVijay Mahadevan       ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr);
724b21a1e88SVijay Mahadevan       ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr);
725e6a796c3SToby Isaac       /* ierr = PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */
7262da392ccSBarry Smith     } else {
727b21a1e88SVijay Mahadevan       ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr);
72863d025dbSVijay Mahadevan     }
72963d025dbSVijay Mahadevan     break;
73063d025dbSVijay Mahadevan   case 3:
73163d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
73263d025dbSVijay Mahadevan     if (nverts == 4) {
733a86ed394SVijay Mahadevan       PetscInt order;
734a86ed394SVijay Mahadevan       const PetscInt npoints = 4; // Choose between 4 and 10
735181f196bSVijay Mahadevan       ierr = PetscMalloc2(npoints * 3, &x, npoints, &w);CHKERRQ(ierr);
736181f196bSVijay Mahadevan       if (npoints == 4) { /*  KEAST1, K1,  N=4, O=4 */
737181f196bSVijay Mahadevan         const PetscReal x_4[12] = { 0.5854101966249685, 0.1381966011250105, 0.1381966011250105,
738181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.1381966011250105,
739181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.5854101966249685,
740181f196bSVijay Mahadevan                                     0.1381966011250105, 0.5854101966249685, 0.1381966011250105
741181f196bSVijay Mahadevan                                   };
742580bdb30SBarry Smith         ierr = PetscArraycpy(x, x_4, 12);CHKERRQ(ierr);
743181f196bSVijay Mahadevan 
744181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0;
745181f196bSVijay Mahadevan         order = 4;
7462da392ccSBarry Smith       } else if (npoints == 10) { /*  KEAST3, K3  N=10, O=10 */
747181f196bSVijay Mahadevan         const PetscReal x_10[30] = { 0.5684305841968444, 0.1438564719343852, 0.1438564719343852,
748181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.1438564719343852,
749181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.5684305841968444,
750181f196bSVijay Mahadevan                                      0.1438564719343852, 0.5684305841968444, 0.1438564719343852,
751181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.5000000000000000,
752181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.5000000000000000,
753181f196bSVijay Mahadevan                                      0.5000000000000000, 0.5000000000000000, 0.0000000000000000,
754181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.0000000000000000,
755181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.0000000000000000,
756181f196bSVijay Mahadevan                                      0.0000000000000000, 0.0000000000000000, 0.5000000000000000
757181f196bSVijay Mahadevan                                    };
758580bdb30SBarry Smith         ierr = PetscArraycpy(x, x_10, 30);CHKERRQ(ierr);
759181f196bSVijay Mahadevan 
760181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 0.2177650698804054;
761181f196bSVijay Mahadevan         w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631;
762181f196bSVijay Mahadevan         order = 10;
7632da392ccSBarry Smith       } else SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %D", npoints);
764181f196bSVijay Mahadevan       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr);
765181f196bSVijay Mahadevan       ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr);
766b21a1e88SVijay Mahadevan       ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr);
767e6a796c3SToby Isaac       /* ierr = PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */
7682da392ccSBarry Smith     } else {
769b21a1e88SVijay Mahadevan       ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr);
77063d025dbSVijay Mahadevan     }
77163d025dbSVijay Mahadevan     break;
77263d025dbSVijay Mahadevan   default:
77363d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
77463d025dbSVijay Mahadevan   }
77563d025dbSVijay Mahadevan   PetscFunctionReturn(0);
77663d025dbSVijay Mahadevan }
77763d025dbSVijay Mahadevan 
778181f196bSVijay Mahadevan /* Compute Jacobians */
779181f196bSVijay Mahadevan PetscErrorCode ComputeJacobian_Internal (const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts,
780a86ed394SVijay Mahadevan   PetscReal *jacobian, PetscReal *ijacobian, PetscReal* dvolume)
781181f196bSVijay Mahadevan {
782a86ed394SVijay Mahadevan   PetscInt       i;
7832417220eSVijay Mahadevan   PetscReal      volume=1.0;
784181f196bSVijay Mahadevan   PetscErrorCode ierr;
785181f196bSVijay Mahadevan 
786181f196bSVijay Mahadevan   PetscFunctionBegin;
787181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
788181f196bSVijay Mahadevan   PetscValidPointer(quad, 4);
789181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 5);
790580bdb30SBarry Smith   ierr = PetscArrayzero(jacobian, dim * dim);CHKERRQ(ierr);
791181f196bSVijay Mahadevan   if (ijacobian) {
792580bdb30SBarry Smith     ierr = PetscArrayzero(ijacobian, dim * dim);CHKERRQ(ierr);
793181f196bSVijay Mahadevan   }
794181f196bSVijay Mahadevan   if (phypts) {
795580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, /*npts=1 * */ 3);CHKERRQ(ierr);
796181f196bSVijay Mahadevan   }
797181f196bSVijay Mahadevan 
798181f196bSVijay Mahadevan   if (dim == 1) {
799181f196bSVijay Mahadevan     const PetscReal& r = quad[0];
800181f196bSVijay Mahadevan     if (nverts == 2) { /* Linear Edge */
801181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 1.0 };
802181f196bSVijay Mahadevan 
803181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
804181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
805181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
806181f196bSVijay Mahadevan       }
8072da392ccSBarry Smith     } else if (nverts == 3) { /* Quadratic Edge */
808181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r), 4.0 * r - 1.0};
809181f196bSVijay Mahadevan 
810181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
811181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
812181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
813181f196bSVijay Mahadevan       }
8142da392ccSBarry Smith     } else {
815181f196bSVijay 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);
816181f196bSVijay Mahadevan     }
817181f196bSVijay Mahadevan 
818181f196bSVijay Mahadevan     if (ijacobian) {
819181f196bSVijay Mahadevan       /* invert the jacobian */
820181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
821181f196bSVijay Mahadevan     }
8222da392ccSBarry Smith   } else if (dim == 2) {
823181f196bSVijay Mahadevan 
824181f196bSVijay Mahadevan     if (nverts == 4) { /* Linear Quadrangle */
825181f196bSVijay Mahadevan       const PetscReal& r = quad[0];
826181f196bSVijay Mahadevan       const PetscReal& s = quad[1];
827181f196bSVijay Mahadevan 
828181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
829181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
830181f196bSVijay Mahadevan 
831181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
832181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
833181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]  * vertices[0];
834181f196bSVijay Mahadevan         jacobian[2] += dNi_dxi[i]  * vertices[1];
835181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
836181f196bSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
837181f196bSVijay Mahadevan       }
8382da392ccSBarry Smith     } else if (nverts == 3) { /* Linear triangle */
839181f196bSVijay Mahadevan       const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1];
840181f196bSVijay Mahadevan 
841181f196bSVijay Mahadevan       /* Jacobian is constant */
842181f196bSVijay Mahadevan       jacobian[0] = (coordinates[0 * 3 + 0] - x2); jacobian[1] = (coordinates[1 * 3 + 0] - x2);
843181f196bSVijay Mahadevan       jacobian[2] = (coordinates[0 * 3 + 1] - y2); jacobian[3] = (coordinates[1 * 3 + 1] - y2);
8442da392ccSBarry Smith     } else 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);
845181f196bSVijay Mahadevan 
846181f196bSVijay Mahadevan     /* invert the jacobian */
847181f196bSVijay Mahadevan     if (ijacobian) {
848a86ed394SVijay Mahadevan       ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr);
849181f196bSVijay Mahadevan     }
8502da392ccSBarry Smith   } else {
851181f196bSVijay Mahadevan 
852181f196bSVijay Mahadevan     if (nverts == 8) { /* Linear Hexahedra */
853181f196bSVijay Mahadevan       const PetscReal &r = quad[0];
854181f196bSVijay Mahadevan       const PetscReal &s = quad[1];
855181f196bSVijay Mahadevan       const PetscReal &t = quad[2];
856181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s) * ( 1.0 - t),
857181f196bSVijay Mahadevan                                        ( 1.0 - s) * ( 1.0 - t),
858a86ed394SVijay Mahadevan                                        (       s) * ( 1.0 - t),
859a86ed394SVijay Mahadevan                                      - (       s) * ( 1.0 - t),
860a86ed394SVijay Mahadevan                                      - ( 1.0 - s) * (       t),
861a86ed394SVijay Mahadevan                                        ( 1.0 - s) * (       t),
862a86ed394SVijay Mahadevan                                        (       s) * (       t),
863a86ed394SVijay Mahadevan                                      - (       s) * (       t)
864181f196bSVijay Mahadevan                                     };
865181f196bSVijay Mahadevan 
866181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r) * ( 1.0 - t),
867a86ed394SVijay Mahadevan                                        - (       r) * ( 1.0 - t),
868a86ed394SVijay Mahadevan                                          (       r) * ( 1.0 - t),
869181f196bSVijay Mahadevan                                          ( 1.0 - r) * ( 1.0 - t),
870a86ed394SVijay Mahadevan                                        - ( 1.0 - r) * (       t),
871a86ed394SVijay Mahadevan                                        - (       r) * (       t),
872a86ed394SVijay Mahadevan                                          (       r) * (       t),
873a86ed394SVijay Mahadevan                                          ( 1.0 - r) * (       t)
874181f196bSVijay Mahadevan                                       };
875181f196bSVijay Mahadevan 
876181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r) * ( 1.0 - s),
877a86ed394SVijay Mahadevan                                         - (       r) * ( 1.0 - s),
878a86ed394SVijay Mahadevan                                         - (       r) * (       s),
879a86ed394SVijay Mahadevan                                         - ( 1.0 - r) * (       s),
880181f196bSVijay Mahadevan                                           ( 1.0 - r) * ( 1.0 - s),
881a86ed394SVijay Mahadevan                                           (       r) * ( 1.0 - s),
882a86ed394SVijay Mahadevan                                           (       r) * (       s),
883a86ed394SVijay Mahadevan                                           ( 1.0 - r) * (       s)
884181f196bSVijay Mahadevan                                      };
885a86ed394SVijay Mahadevan 
886181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
887181f196bSVijay Mahadevan         const PetscReal* vertex = coordinates + i * 3;
888181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
889181f196bSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
890181f196bSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
891181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
892181f196bSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
893181f196bSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
894181f196bSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
895181f196bSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
896181f196bSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
897181f196bSVijay Mahadevan       }
8982da392ccSBarry Smith     } else if (nverts == 4) { /* Linear Tetrahedra */
899181f196bSVijay Mahadevan       const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2];
900181f196bSVijay Mahadevan 
901181f196bSVijay Mahadevan       /* compute the jacobian */
902181f196bSVijay Mahadevan       jacobian[0] = coordinates[1 * 3 + 0] - x0;  jacobian[1] = coordinates[2 * 3 + 0] - x0; jacobian[2] = coordinates[3 * 3 + 0] - x0;
903181f196bSVijay Mahadevan       jacobian[3] = coordinates[1 * 3 + 1] - y0;  jacobian[4] = coordinates[2 * 3 + 1] - y0; jacobian[5] = coordinates[3 * 3 + 1] - y0;
904181f196bSVijay Mahadevan       jacobian[6] = coordinates[1 * 3 + 2] - z0;  jacobian[7] = coordinates[2 * 3 + 2] - z0; jacobian[8] = coordinates[3 * 3 + 2] - z0;
9052da392ccSBarry Smith     } else 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);
906181f196bSVijay Mahadevan 
907181f196bSVijay Mahadevan     if (ijacobian) {
908181f196bSVijay Mahadevan       /* invert the jacobian */
909a86ed394SVijay Mahadevan       ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr);
910181f196bSVijay Mahadevan     }
911181f196bSVijay Mahadevan 
912181f196bSVijay Mahadevan   }
913a86ed394SVijay 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);
914a86ed394SVijay Mahadevan   if (dvolume) *dvolume = volume;
915181f196bSVijay Mahadevan   PetscFunctionReturn(0);
916181f196bSVijay Mahadevan }
917181f196bSVijay Mahadevan 
918181f196bSVijay Mahadevan PetscErrorCode FEMComputeBasis_JandF(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quadrature, PetscReal *phypts,
919181f196bSVijay Mahadevan                                      PetscReal *phibasis, PetscReal *jacobian, PetscReal *ijacobian, PetscReal* volume)
920181f196bSVijay Mahadevan {
921181f196bSVijay Mahadevan   PetscErrorCode ierr;
922181f196bSVijay Mahadevan 
923181f196bSVijay Mahadevan   PetscFunctionBegin;
924181f196bSVijay Mahadevan   switch (dim) {
925181f196bSVijay Mahadevan     case 1:
926181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts,
927181f196bSVijay Mahadevan             NULL, phibasis, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
928181f196bSVijay Mahadevan       break;
929181f196bSVijay Mahadevan     case 2:
930181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts,
931181f196bSVijay Mahadevan             NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
932181f196bSVijay Mahadevan       break;
933181f196bSVijay Mahadevan     case 3:
934181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts,
935181f196bSVijay Mahadevan             NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
936181f196bSVijay Mahadevan       break;
937181f196bSVijay Mahadevan     default:
938181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
939181f196bSVijay Mahadevan   }
940181f196bSVijay Mahadevan   PetscFunctionReturn(0);
941181f196bSVijay Mahadevan }
942181f196bSVijay Mahadevan 
943cab5ea25SPierre Jolivet /*@C
94497b73a88SSatish Balay   DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the
945a86ed394SVijay Mahadevan   canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration,
946a86ed394SVijay Mahadevan   the basis function at the parametric point is also evaluated optionally.
947a86ed394SVijay Mahadevan 
948a86ed394SVijay Mahadevan   Input Parameter:
949a2b725a8SWilliam Gropp +  PetscInt  dim -         the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
950a2b725a8SWilliam Gropp .  PetscInt nverts -       the number of vertices in the physical element
951a2b725a8SWilliam Gropp .  PetscReal coordinates - the coordinates of vertices in the physical element
952a2b725a8SWilliam Gropp -  PetscReal[3] xphy -     the coordinates of physical point for which natural coordinates (in reference frame) are sought
953a86ed394SVijay Mahadevan 
954a86ed394SVijay Mahadevan   Output Parameter:
955a2b725a8SWilliam Gropp +  PetscReal[3] natparam - the natural coordinates (in reference frame) corresponding to xphy
956a2b725a8SWilliam Gropp -  PetscReal[nverts] phi - the basis functions evaluated at the natural coordinates (natparam)
957a86ed394SVijay Mahadevan 
958edc382c3SSatish Balay   Level: advanced
959edc382c3SSatish Balay 
960a86ed394SVijay Mahadevan @*/
961181f196bSVijay Mahadevan PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal* xphy, PetscReal* natparam, PetscReal* phi)
962181f196bSVijay Mahadevan {
963a86ed394SVijay Mahadevan   /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */
964181f196bSVijay Mahadevan   const PetscReal tol = 1.0e-10;
965181f196bSVijay Mahadevan   const PetscInt  max_iterations = 10;
966181f196bSVijay Mahadevan   const PetscReal error_tol_sqr = tol*tol;
967181f196bSVijay Mahadevan   PetscReal       phibasis[8], jacobian[9], ijacobian[9], volume;
968181f196bSVijay Mahadevan   PetscReal       phypts[3] = {0.0, 0.0, 0.0};
969181f196bSVijay Mahadevan   PetscReal       delta[3] = {0.0, 0.0, 0.0};
970181f196bSVijay Mahadevan   PetscErrorCode  ierr;
971181f196bSVijay Mahadevan   PetscInt        iters=0;
972181f196bSVijay Mahadevan   PetscReal       error=1.0;
973181f196bSVijay Mahadevan 
974181f196bSVijay Mahadevan   PetscFunctionBegin;
975181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
976181f196bSVijay Mahadevan   PetscValidPointer(xphy, 4);
977181f196bSVijay Mahadevan   PetscValidPointer(natparam, 5);
978181f196bSVijay Mahadevan 
979580bdb30SBarry Smith   ierr = PetscArrayzero(jacobian, dim * dim);CHKERRQ(ierr);
980580bdb30SBarry Smith   ierr = PetscArrayzero(ijacobian, dim * dim);CHKERRQ(ierr);
981580bdb30SBarry Smith   ierr = PetscArrayzero(phibasis, nverts);CHKERRQ(ierr);
982181f196bSVijay Mahadevan 
983a86ed394SVijay Mahadevan   /* zero initial guess */
984a86ed394SVijay Mahadevan   natparam[0] = natparam[1] = natparam[2] = 0.0;
985181f196bSVijay Mahadevan 
986a86ed394SVijay Mahadevan   /* Compute delta = evaluate( xi) - x */
987a86ed394SVijay Mahadevan   ierr = FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume);CHKERRQ(ierr);
988a86ed394SVijay Mahadevan 
989a86ed394SVijay Mahadevan   error = 0.0;
990a86ed394SVijay Mahadevan   switch(dim) {
991a86ed394SVijay Mahadevan     case 3:
992181f196bSVijay Mahadevan       delta[2] = phypts[2] - xphy[2];
993a86ed394SVijay Mahadevan       error += (delta[2]*delta[2]);
994a86ed394SVijay Mahadevan     case 2:
995a86ed394SVijay Mahadevan       delta[1] = phypts[1] - xphy[1];
996a86ed394SVijay Mahadevan       error += (delta[1]*delta[1]);
997a86ed394SVijay Mahadevan     case 1:
998a86ed394SVijay Mahadevan       delta[0] = phypts[0] - xphy[0];
999a86ed394SVijay Mahadevan       error += (delta[0]*delta[0]);
1000a86ed394SVijay Mahadevan       break;
1001a86ed394SVijay Mahadevan   }
1002a86ed394SVijay Mahadevan 
1003181f196bSVijay Mahadevan   while (error > error_tol_sqr) {
1004181f196bSVijay Mahadevan 
1005181f196bSVijay Mahadevan     if (++iters > max_iterations)
1006181f196bSVijay 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]);
1007181f196bSVijay Mahadevan 
1008181f196bSVijay Mahadevan     /* Compute natparam -= J.inverse() * delta */
1009181f196bSVijay Mahadevan     switch(dim) {
1010181f196bSVijay Mahadevan       case 1:
1011181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0];
1012181f196bSVijay Mahadevan         break;
1013181f196bSVijay Mahadevan       case 2:
1014181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1];
1015181f196bSVijay Mahadevan         natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1];
1016181f196bSVijay Mahadevan         break;
1017181f196bSVijay Mahadevan       case 3:
1018181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2];
1019181f196bSVijay Mahadevan         natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2];
1020181f196bSVijay Mahadevan         natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2];
1021181f196bSVijay Mahadevan         break;
1022181f196bSVijay Mahadevan     }
1023181f196bSVijay Mahadevan 
1024181f196bSVijay Mahadevan     /* Compute delta = evaluate( xi) - x */
1025a86ed394SVijay Mahadevan     ierr = FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume);CHKERRQ(ierr);
1026181f196bSVijay Mahadevan 
1027a86ed394SVijay Mahadevan     error = 0.0;
1028a86ed394SVijay Mahadevan     switch(dim) {
1029a86ed394SVijay Mahadevan       case 3:
1030181f196bSVijay Mahadevan         delta[2] = phypts[2] - xphy[2];
1031a86ed394SVijay Mahadevan         error += (delta[2]*delta[2]);
1032a86ed394SVijay Mahadevan       case 2:
1033a86ed394SVijay Mahadevan         delta[1] = phypts[1] - xphy[1];
1034a86ed394SVijay Mahadevan         error += (delta[1]*delta[1]);
1035a86ed394SVijay Mahadevan       case 1:
1036a86ed394SVijay Mahadevan         delta[0] = phypts[0] - xphy[0];
1037a86ed394SVijay Mahadevan         error += (delta[0]*delta[0]);
1038a86ed394SVijay Mahadevan         break;
1039a86ed394SVijay Mahadevan     }
1040181f196bSVijay Mahadevan   }
1041181f196bSVijay Mahadevan   if (phi) {
1042580bdb30SBarry Smith     ierr = PetscArraycpy(phi, phibasis, nverts);CHKERRQ(ierr);
1043181f196bSVijay Mahadevan   }
1044181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1045181f196bSVijay Mahadevan }
1046