xref: /petsc/src/dm/impls/moab/dmmbfem.cxx (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 */
79371c9d4SSatish Balay static inline PetscReal DMatrix_Determinant_2x2_Internal(const PetscReal inmat[2 * 2]) {
863d025dbSVijay Mahadevan   return inmat[0] * inmat[3] - inmat[1] * inmat[2];
963d025dbSVijay Mahadevan }
1063d025dbSVijay Mahadevan 
119371c9d4SSatish Balay static inline PetscErrorCode DMatrix_Invert_2x2_Internal(const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant) {
1263d025dbSVijay Mahadevan   PetscReal det = DMatrix_Determinant_2x2_Internal(inmat);
1363d025dbSVijay Mahadevan   if (outmat) {
1463d025dbSVijay Mahadevan     outmat[0] = inmat[3] / det;
1563d025dbSVijay Mahadevan     outmat[1] = -inmat[1] / det;
1663d025dbSVijay Mahadevan     outmat[2] = -inmat[2] / det;
1763d025dbSVijay Mahadevan     outmat[3] = inmat[0] / det;
1863d025dbSVijay Mahadevan   }
1963d025dbSVijay Mahadevan   if (determinant) *determinant = det;
20362febeeSStefano Zampini   return 0;
2163d025dbSVijay Mahadevan }
2263d025dbSVijay Mahadevan 
239371c9d4SSatish Balay static inline PetscReal DMatrix_Determinant_3x3_Internal(const PetscReal inmat[3 * 3]) {
249371c9d4SSatish Balay   return inmat[0] * (inmat[8] * inmat[4] - inmat[7] * inmat[5]) - inmat[3] * (inmat[8] * inmat[1] - inmat[7] * inmat[2]) + inmat[6] * (inmat[5] * inmat[1] - inmat[4] * inmat[2]);
2563d025dbSVijay Mahadevan }
2663d025dbSVijay Mahadevan 
279371c9d4SSatish Balay static inline PetscErrorCode DMatrix_Invert_3x3_Internal(const PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant) {
28181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_3x3_Internal(inmat);
2963d025dbSVijay Mahadevan   if (outmat) {
3063d025dbSVijay Mahadevan     outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det;
3163d025dbSVijay Mahadevan     outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det;
3263d025dbSVijay Mahadevan     outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det;
3363d025dbSVijay Mahadevan     outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det;
3463d025dbSVijay Mahadevan     outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det;
3563d025dbSVijay Mahadevan     outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det;
3663d025dbSVijay Mahadevan     outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det;
3763d025dbSVijay Mahadevan     outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det;
3863d025dbSVijay Mahadevan     outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det;
3963d025dbSVijay Mahadevan   }
4063d025dbSVijay Mahadevan   if (determinant) *determinant = det;
41362febeeSStefano Zampini   return 0;
4263d025dbSVijay Mahadevan }
4363d025dbSVijay Mahadevan 
449371c9d4SSatish Balay inline PetscReal DMatrix_Determinant_4x4_Internal(PetscReal inmat[4 * 4]) {
459371c9d4SSatish Balay   return inmat[0 + 0 * 4] * (inmat[1 + 1 * 4] * (inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4]) - inmat[1 + 2 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4]) + inmat[1 + 3 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4])) - inmat[0 + 1 * 4] * (inmat[1 + 0 * 4] * (inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4]) - inmat[1 + 2 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4]) + inmat[1 + 3 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4])) + inmat[0 + 2 * 4] * (inmat[1 + 0 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4]) - inmat[1 + 1 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4]) + inmat[1 + 3 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4])) - inmat[0 + 3 * 4] * (inmat[1 + 0 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4]) - inmat[1 + 1 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4]) + inmat[1 + 2 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4]));
46181f196bSVijay Mahadevan }
47181f196bSVijay Mahadevan 
489371c9d4SSatish Balay inline PetscErrorCode DMatrix_Invert_4x4_Internal(PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant) {
49181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_4x4_Internal(inmat);
50181f196bSVijay Mahadevan   if (outmat) {
51181f196bSVijay 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;
52181f196bSVijay 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;
53181f196bSVijay 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;
54181f196bSVijay 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;
55181f196bSVijay 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;
56181f196bSVijay 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;
57181f196bSVijay 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;
58181f196bSVijay 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;
59181f196bSVijay 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;
60181f196bSVijay 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;
61181f196bSVijay 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;
62181f196bSVijay 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;
63181f196bSVijay 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;
64181f196bSVijay 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;
65181f196bSVijay 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;
66181f196bSVijay 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;
67181f196bSVijay Mahadevan   }
68181f196bSVijay Mahadevan   if (determinant) *determinant = det;
69362febeeSStefano Zampini   return 0;
70181f196bSVijay Mahadevan }
71181f196bSVijay Mahadevan 
72cab5ea25SPierre Jolivet /*@C
7397b73a88SSatish Balay   Compute_Lagrange_Basis_1D_Internal - Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element.
7463d025dbSVijay Mahadevan 
7563d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a linear or quadratic edge element.
7663d025dbSVijay Mahadevan 
7763d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
7863d025dbSVijay Mahadevan   and their derivatives with respect to X.
7963d025dbSVijay Mahadevan 
8063d025dbSVijay Mahadevan   Notes:
8163d025dbSVijay Mahadevan 
8263d025dbSVijay Mahadevan   Example Physical Element
83a2b725a8SWilliam Gropp .vb
8463d025dbSVijay Mahadevan     1-------2        1----3----2
8563d025dbSVijay Mahadevan       EDGE2             EDGE3
86a2b725a8SWilliam Gropp .ve
8763d025dbSVijay Mahadevan 
88d8d19677SJose E. Roman   Input Parameters:
89a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
90a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
91a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
92a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
9363d025dbSVijay Mahadevan 
94d8d19677SJose E. Roman   Output Parameters:
95a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
96a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
97a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
986b867d5aSJose E. Roman .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
9967b8a455SSatish Balay .  jacobian -                  jacobian
10067b8a455SSatish Balay .  ijacobian -                 ijacobian
10167b8a455SSatish Balay -  volume -                    volume
10263d025dbSVijay Mahadevan 
103edc382c3SSatish Balay   Level: advanced
104edc382c3SSatish Balay 
10563d025dbSVijay Mahadevan @*/
1069371c9d4SSatish Balay PetscErrorCode Compute_Lagrange_Basis_1D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) {
10763d025dbSVijay Mahadevan   int i, j;
10863d025dbSVijay Mahadevan 
109181f196bSVijay Mahadevan   PetscFunctionBegin;
110181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 9);
111181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 10);
112181f196bSVijay Mahadevan   PetscValidPointer(volume, 11);
113*48a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
11463d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
1159566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
11663d025dbSVijay Mahadevan   }
11763d025dbSVijay Mahadevan   if (nverts == 2) { /* Linear Edge */
11863d025dbSVijay Mahadevan 
1192da392ccSBarry Smith     for (j = 0; j < npts; j++) {
120a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
121181f196bSVijay Mahadevan       const PetscReal r      = quad[j];
12263d025dbSVijay Mahadevan 
123181f196bSVijay Mahadevan       phi[0 + offset] = (1.0 - r);
124181f196bSVijay Mahadevan       phi[1 + offset] = (r);
12563d025dbSVijay Mahadevan 
126181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2] = {-1.0, 1.0};
12763d025dbSVijay Mahadevan 
128181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
12963d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
130181f196bSVijay Mahadevan         const PetscReal *vertices = coords + i * 3;
131181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
1329371c9d4SSatish Balay         if (phypts) { phypts[3 * j + 0] += phi[i + offset] * vertices[0]; }
13363d025dbSVijay Mahadevan       }
13463d025dbSVijay Mahadevan 
13563d025dbSVijay Mahadevan       /* invert the jacobian */
136181f196bSVijay Mahadevan       *volume      = jacobian[0];
137181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
138181f196bSVijay Mahadevan       jxw[j] *= *volume;
13963d025dbSVijay Mahadevan 
14063d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
14163d025dbSVijay Mahadevan       for (i = 0; i < nverts; i++) {
142181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
14363d025dbSVijay Mahadevan       }
14463d025dbSVijay Mahadevan     }
1452da392ccSBarry Smith   } else if (nverts == 3) { /* Quadratic Edge */
14663d025dbSVijay Mahadevan 
1472da392ccSBarry Smith     for (j = 0; j < npts; j++) {
148a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
149181f196bSVijay Mahadevan       const PetscReal r      = quad[j];
15063d025dbSVijay Mahadevan 
151181f196bSVijay Mahadevan       phi[0 + offset] = 1.0 + r * (2.0 * r - 3.0);
152181f196bSVijay Mahadevan       phi[1 + offset] = 4.0 * r * (1.0 - r);
153181f196bSVijay Mahadevan       phi[2 + offset] = r * (2.0 * r - 1.0);
15463d025dbSVijay Mahadevan 
155181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3] = {4 * r - 3.0, 4 * (1.0 - 2.0 * r), 4.0 * r - 1.0};
15663d025dbSVijay Mahadevan 
157181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
15863d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
159181f196bSVijay Mahadevan         const PetscReal *vertices = coords + i * 3;
160181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
1619371c9d4SSatish Balay         if (phypts) { phypts[3 * j + 0] += phi[i + offset] * vertices[0]; }
16263d025dbSVijay Mahadevan       }
16363d025dbSVijay Mahadevan 
16463d025dbSVijay Mahadevan       /* invert the jacobian */
165181f196bSVijay Mahadevan       *volume      = jacobian[0];
166181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
167181f196bSVijay Mahadevan       if (jxw) 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     }
17463a3b9bcSJacob Faibussowitsch   } else SETERRQ(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 : %" PetscInt_FMT, nverts);
17563d025dbSVijay Mahadevan   PetscFunctionReturn(0);
17663d025dbSVijay Mahadevan }
17763d025dbSVijay Mahadevan 
178cab5ea25SPierre Jolivet /*@C
17997b73a88SSatish Balay   Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.
18063d025dbSVijay Mahadevan 
18163d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a quadrangle or triangle.
18263d025dbSVijay Mahadevan 
18363d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
18463d025dbSVijay Mahadevan   and their derivatives with respect to X and Y.
18563d025dbSVijay Mahadevan 
18663d025dbSVijay Mahadevan   Notes:
18763d025dbSVijay Mahadevan 
18863d025dbSVijay Mahadevan   Example Physical Element (QUAD4)
189a2b725a8SWilliam Gropp .vb
19063d025dbSVijay Mahadevan     4------3        s
19163d025dbSVijay Mahadevan     |      |        |
19263d025dbSVijay Mahadevan     |      |        |
19363d025dbSVijay Mahadevan     |      |        |
19463d025dbSVijay Mahadevan     1------2        0-------r
195a2b725a8SWilliam Gropp .ve
19663d025dbSVijay Mahadevan 
197d8d19677SJose E. Roman   Input Parameters:
198a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
199a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
200a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
201a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
20263d025dbSVijay Mahadevan 
203d8d19677SJose E. Roman   Output Parameters:
204a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
205a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
206a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
207a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
2086b867d5aSJose E. Roman .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
20967b8a455SSatish Balay .  jacobian -                  jacobian
21067b8a455SSatish Balay .  ijacobian -                 ijacobian
21167b8a455SSatish Balay -  volume -                    volume
21263d025dbSVijay Mahadevan 
213edc382c3SSatish Balay   Level: advanced
214edc382c3SSatish Balay 
21563d025dbSVijay Mahadevan @*/
2169371c9d4SSatish Balay PetscErrorCode Compute_Lagrange_Basis_2D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) {
217a86ed394SVijay Mahadevan   PetscInt i, j, k;
21863d025dbSVijay Mahadevan 
21963d025dbSVijay Mahadevan   PetscFunctionBegin;
220181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 10);
221181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 11);
222181f196bSVijay Mahadevan   PetscValidPointer(volume, 12);
2239566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phi, npts));
224*48a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
22563d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
2269566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
2279566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidy, npts * nverts));
22863d025dbSVijay Mahadevan   }
22963d025dbSVijay Mahadevan   if (nverts == 4) { /* Linear Quadrangle */
23063d025dbSVijay Mahadevan 
2319371c9d4SSatish Balay     for (j = 0; j < npts; j++) {
232a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
233181f196bSVijay Mahadevan       const PetscReal r      = quad[0 + j * 2];
234181f196bSVijay Mahadevan       const PetscReal s      = quad[1 + j * 2];
23563d025dbSVijay Mahadevan 
23663d025dbSVijay Mahadevan       phi[0 + offset] = (1.0 - r) * (1.0 - s);
23763d025dbSVijay Mahadevan       phi[1 + offset] = r * (1.0 - s);
23863d025dbSVijay Mahadevan       phi[2 + offset] = r * s;
23963d025dbSVijay Mahadevan       phi[3 + offset] = (1.0 - r) * s;
24063d025dbSVijay Mahadevan 
241181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = {-1.0 + s, 1.0 - s, s, -s};
242181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = {-1.0 + r, -r, r, 1.0 - r};
24363d025dbSVijay Mahadevan 
2449566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(jacobian, 4));
2459566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ijacobian, 4));
24663d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
247181f196bSVijay Mahadevan         const PetscReal *vertices = coords + i * 3;
24863d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
24963d025dbSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
25063d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
25163d025dbSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
252181f196bSVijay Mahadevan         if (phypts) {
253181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
254181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
255181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
256181f196bSVijay Mahadevan         }
25763d025dbSVijay Mahadevan       }
25863d025dbSVijay Mahadevan 
25963d025dbSVijay Mahadevan       /* invert the jacobian */
2609566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
26108401ef6SPierre Jolivet       PetscCheck(*volume >= 1e-12, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Quadrangular element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
26263d025dbSVijay Mahadevan 
263181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
26463d025dbSVijay Mahadevan 
265181f196bSVijay Mahadevan       /*  Let us compute the bases derivatives by scaling with inverse jacobian. */
266181f196bSVijay Mahadevan       if (dphidx) {
26763d025dbSVijay Mahadevan         for (i = 0; i < nverts; i++) {
268a86ed394SVijay Mahadevan           for (k = 0; k < 2; ++k) {
26963d025dbSVijay Mahadevan             if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0];
27063d025dbSVijay Mahadevan             if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1];
271181f196bSVijay Mahadevan           } /* for k=[0..2) */
272181f196bSVijay Mahadevan         }   /* for i=[0..nverts) */
273181f196bSVijay Mahadevan       }     /* if (dphidx) */
27463d025dbSVijay Mahadevan     }
2752da392ccSBarry Smith   } else if (nverts == 3) { /* Linear triangle */
2762da392ccSBarry Smith     const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1];
27763d025dbSVijay Mahadevan 
2789566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(jacobian, 4));
2799566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ijacobian, 4));
28063d025dbSVijay Mahadevan 
28163d025dbSVijay Mahadevan     /* Jacobian is constant */
2829371c9d4SSatish Balay     jacobian[0] = (coords[0 * 3 + 0] - x2);
2839371c9d4SSatish Balay     jacobian[1] = (coords[1 * 3 + 0] - x2);
2849371c9d4SSatish Balay     jacobian[2] = (coords[0 * 3 + 1] - y2);
2859371c9d4SSatish Balay     jacobian[3] = (coords[1 * 3 + 1] - y2);
28663d025dbSVijay Mahadevan 
28763d025dbSVijay Mahadevan     /* invert the jacobian */
2889566063dSJacob Faibussowitsch     PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
28908401ef6SPierre Jolivet     PetscCheck(*volume >= 1e-12, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangular element has zero volume: %g. Degenerate element or invalid connectivity", (double)*volume);
290181f196bSVijay Mahadevan 
291181f196bSVijay Mahadevan     const PetscReal Dx[3] = {ijacobian[0], ijacobian[2], -ijacobian[0] - ijacobian[2]};
292181f196bSVijay Mahadevan     const PetscReal Dy[3] = {ijacobian[1], ijacobian[3], -ijacobian[1] - ijacobian[3]};
29363d025dbSVijay Mahadevan 
2942da392ccSBarry Smith     for (j = 0; j < npts; j++) {
295a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
296181f196bSVijay Mahadevan       const PetscReal r      = quad[0 + j * 2];
297181f196bSVijay Mahadevan       const PetscReal s      = quad[1 + j * 2];
29863d025dbSVijay Mahadevan 
299181f196bSVijay Mahadevan       if (jxw) jxw[j] *= 0.5;
30063d025dbSVijay Mahadevan 
301181f196bSVijay Mahadevan       /* const PetscReal Ni[3]  = { r, s, 1.0 - r - s }; */
302181f196bSVijay Mahadevan       const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s;
303181f196bSVijay Mahadevan       const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s;
304181f196bSVijay Mahadevan       if (phypts) {
305181f196bSVijay Mahadevan         phypts[offset + 0] = phipts_x;
306181f196bSVijay Mahadevan         phypts[offset + 1] = phipts_y;
307181f196bSVijay Mahadevan       }
30863d025dbSVijay Mahadevan 
309110fc3b0SBarry Smith       /* \phi_0 = (b.y - c.y) x + (b.x - c.x) y + c.x b.y - b.x c.y */
310181f196bSVijay Mahadevan       phi[0 + offset] = (ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2));
311110fc3b0SBarry Smith       /* \phi_1 = (c.y - a.y) x + (a.x - c.x) y + c.x a.y - a.x c.y */
312181f196bSVijay Mahadevan       phi[1 + offset] = (ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2));
31363d025dbSVijay Mahadevan       phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset];
31463d025dbSVijay Mahadevan 
31563d025dbSVijay Mahadevan       if (dphidx) {
316181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
317181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
318181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
31963d025dbSVijay Mahadevan       }
32063d025dbSVijay Mahadevan 
32163d025dbSVijay Mahadevan       if (dphidy) {
322181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
323181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
324181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
32563d025dbSVijay Mahadevan       }
32663d025dbSVijay Mahadevan     }
32763a3b9bcSJacob Faibussowitsch   } else SETERRQ(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 : %" PetscInt_FMT, nverts);
32863d025dbSVijay Mahadevan   PetscFunctionReturn(0);
32963d025dbSVijay Mahadevan }
33063d025dbSVijay Mahadevan 
331cab5ea25SPierre Jolivet /*@C
33297b73a88SSatish Balay   Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element.
33363d025dbSVijay Mahadevan 
33463d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a hexahedra or tetrahedra.
33563d025dbSVijay Mahadevan 
33663d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
33763d025dbSVijay Mahadevan   and their derivatives with respect to X, Y and Z.
33863d025dbSVijay Mahadevan 
33963d025dbSVijay Mahadevan   Notes:
34063d025dbSVijay Mahadevan 
34163d025dbSVijay Mahadevan   Example Physical Element (HEX8)
342a2b725a8SWilliam Gropp .vb
34363d025dbSVijay Mahadevan       8------7
34463d025dbSVijay Mahadevan      /|     /|        t  s
34563d025dbSVijay Mahadevan     5------6 |        | /
34663d025dbSVijay Mahadevan     | |    | |        |/
34763d025dbSVijay Mahadevan     | 4----|-3        0-------r
34863d025dbSVijay Mahadevan     |/     |/
34963d025dbSVijay Mahadevan     1------2
350a2b725a8SWilliam Gropp .ve
35163d025dbSVijay Mahadevan 
352d8d19677SJose E. Roman   Input Parameters:
353a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
354a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
355a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
356a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
35763d025dbSVijay Mahadevan 
358d8d19677SJose E. Roman   Output Parameters:
359a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
360a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
361a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
362a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
363a2b725a8SWilliam Gropp .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
3646b867d5aSJose E. Roman .  PetscReal dphidz[npts] -    the derivative of the bases wrt Z-direction evaluated at the specified quadrature points
36567b8a455SSatish Balay .  jacobian -                  jacobian
36667b8a455SSatish Balay .  ijacobian -                 ijacobian
36767b8a455SSatish Balay -  volume -                    volume
36863d025dbSVijay Mahadevan 
369edc382c3SSatish Balay   Level: advanced
370edc382c3SSatish Balay 
37163d025dbSVijay Mahadevan @*/
3729371c9d4SSatish Balay PetscErrorCode Compute_Lagrange_Basis_3D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz, PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) {
373a86ed394SVijay Mahadevan   PetscInt i, j, k;
37463d025dbSVijay Mahadevan 
37563d025dbSVijay Mahadevan   PetscFunctionBegin;
376181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 11);
377181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 12);
378181f196bSVijay Mahadevan   PetscValidPointer(volume, 13);
3792da392ccSBarry Smith 
3809566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phi, npts));
381*48a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
38263d025dbSVijay Mahadevan   if (dphidx) {
3839566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
3849566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidy, npts * nverts));
3859566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidz, npts * nverts));
38663d025dbSVijay Mahadevan   }
38763d025dbSVijay Mahadevan 
38863d025dbSVijay Mahadevan   if (nverts == 8) { /* Linear Hexahedra */
38963d025dbSVijay Mahadevan 
3902da392ccSBarry Smith     for (j = 0; j < npts; j++) {
391a86ed394SVijay Mahadevan       const PetscInt   offset = j * nverts;
392181f196bSVijay Mahadevan       const PetscReal &r      = quad[j * 3 + 0];
393181f196bSVijay Mahadevan       const PetscReal &s      = quad[j * 3 + 1];
394181f196bSVijay Mahadevan       const PetscReal &t      = quad[j * 3 + 2];
39563d025dbSVijay Mahadevan 
396a86ed394SVijay Mahadevan       phi[offset + 0] = (1.0 - r) * (1.0 - s) * (1.0 - t); /* 0,0,0 */
397a86ed394SVijay Mahadevan       phi[offset + 1] = (r) * (1.0 - s) * (1.0 - t);       /* 1,0,0 */
398a86ed394SVijay Mahadevan       phi[offset + 2] = (r) * (s) * (1.0 - t);             /* 1,1,0 */
399a86ed394SVijay Mahadevan       phi[offset + 3] = (1.0 - r) * (s) * (1.0 - t);       /* 0,1,0 */
400a86ed394SVijay Mahadevan       phi[offset + 4] = (1.0 - r) * (1.0 - s) * (t);       /* 0,0,1 */
401a86ed394SVijay Mahadevan       phi[offset + 5] = (r) * (1.0 - s) * (t);             /* 1,0,1 */
402a86ed394SVijay Mahadevan       phi[offset + 6] = (r) * (s) * (t);                   /* 1,1,1 */
403a86ed394SVijay Mahadevan       phi[offset + 7] = (1.0 - r) * (s) * (t);             /* 0,1,1 */
40463d025dbSVijay Mahadevan 
4059371c9d4SSatish Balay       const PetscReal dNi_dxi[8] = {-(1.0 - s) * (1.0 - t), (1.0 - s) * (1.0 - t), (s) * (1.0 - t), -(s) * (1.0 - t), -(1.0 - s) * (t), (1.0 - s) * (t), (s) * (t), -(s) * (t)};
40663d025dbSVijay Mahadevan 
4079371c9d4SSatish Balay       const PetscReal dNi_deta[8] = {-(1.0 - r) * (1.0 - t), -(r) * (1.0 - t), (r) * (1.0 - t), (1.0 - r) * (1.0 - t), -(1.0 - r) * (t), -(r) * (t), (r) * (t), (1.0 - r) * (t)};
40863d025dbSVijay Mahadevan 
4099371c9d4SSatish Balay       const PetscReal dNi_dzeta[8] = {-(1.0 - r) * (1.0 - s), -(r) * (1.0 - s), -(r) * (s), -(1.0 - r) * (s), (1.0 - r) * (1.0 - s), (r) * (1.0 - s), (r) * (s), (1.0 - r) * (s)};
41063d025dbSVijay Mahadevan 
4119566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(jacobian, 9));
4129566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ijacobian, 9));
41363d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
414181f196bSVijay Mahadevan         const PetscReal *vertex = coords + i * 3;
41563d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertex[0];
41663d025dbSVijay Mahadevan         jacobian[3] += dNi_dxi[i] * vertex[1];
41763d025dbSVijay Mahadevan         jacobian[6] += dNi_dxi[i] * vertex[2];
41863d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertex[0];
41963d025dbSVijay Mahadevan         jacobian[4] += dNi_deta[i] * vertex[1];
42063d025dbSVijay Mahadevan         jacobian[7] += dNi_deta[i] * vertex[2];
42163d025dbSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
42263d025dbSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
42363d025dbSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
424181f196bSVijay Mahadevan         if (phypts) {
425181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertex[0];
426181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertex[1];
427181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertex[2];
428181f196bSVijay Mahadevan         }
42963d025dbSVijay Mahadevan       }
43063d025dbSVijay Mahadevan 
43163d025dbSVijay Mahadevan       /* invert the jacobian */
4329566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
43308401ef6SPierre Jolivet       PetscCheck(*volume >= 1e-8, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Hexahedral element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
43463d025dbSVijay Mahadevan 
435a86ed394SVijay Mahadevan       if (jxw) jxw[j] *= (*volume);
43663d025dbSVijay Mahadevan 
43763d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
43863d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
439a86ed394SVijay Mahadevan         for (k = 0; k < 3; ++k) {
44063d025dbSVijay Mahadevan           if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0 * 3 + k];
44163d025dbSVijay Mahadevan           if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[1 * 3 + k];
44263d025dbSVijay Mahadevan           if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k];
44363d025dbSVijay Mahadevan         }
44463d025dbSVijay Mahadevan       }
44563d025dbSVijay Mahadevan     }
4462da392ccSBarry Smith   } else if (nverts == 4) { /* Linear Tetrahedra */
4472da392ccSBarry Smith     PetscReal       Dx[4] = {0, 0, 0, 0}, Dy[4] = {0, 0, 0, 0}, Dz[4] = {0, 0, 0, 0};
4482da392ccSBarry Smith     const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2];
44963d025dbSVijay Mahadevan 
4509566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(jacobian, 9));
4519566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ijacobian, 9));
45263d025dbSVijay Mahadevan 
453181f196bSVijay Mahadevan     /* compute the jacobian */
4549371c9d4SSatish Balay     jacobian[0] = coords[1 * 3 + 0] - x0;
4559371c9d4SSatish Balay     jacobian[1] = coords[2 * 3 + 0] - x0;
4569371c9d4SSatish Balay     jacobian[2] = coords[3 * 3 + 0] - x0;
4579371c9d4SSatish Balay     jacobian[3] = coords[1 * 3 + 1] - y0;
4589371c9d4SSatish Balay     jacobian[4] = coords[2 * 3 + 1] - y0;
4599371c9d4SSatish Balay     jacobian[5] = coords[3 * 3 + 1] - y0;
4609371c9d4SSatish Balay     jacobian[6] = coords[1 * 3 + 2] - z0;
4619371c9d4SSatish Balay     jacobian[7] = coords[2 * 3 + 2] - z0;
4629371c9d4SSatish Balay     jacobian[8] = coords[3 * 3 + 2] - z0;
46363d025dbSVijay Mahadevan 
46463d025dbSVijay Mahadevan     /* invert the jacobian */
4659566063dSJacob Faibussowitsch     PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
46608401ef6SPierre Jolivet     PetscCheck(*volume >= 1e-8, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral element has zero volume: %g. Degenerate element or invalid connectivity", (double)*volume);
46763d025dbSVijay Mahadevan 
468181f196bSVijay Mahadevan     if (dphidx) {
4699371c9d4SSatish Balay       Dx[0] = (coords[1 + 2 * 3] * (coords[2 + 1 * 3] - coords[2 + 3 * 3]) - coords[1 + 1 * 3] * (coords[2 + 2 * 3] - coords[2 + 3 * 3]) - coords[1 + 3 * 3] * (coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume;
4709371c9d4SSatish Balay       Dx[1] = -(coords[1 + 2 * 3] * (coords[2 + 0 * 3] - coords[2 + 3 * 3]) - coords[1 + 0 * 3] * (coords[2 + 2 * 3] - coords[2 + 3 * 3]) - coords[1 + 3 * 3] * (coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume;
4719371c9d4SSatish Balay       Dx[2] = (coords[1 + 1 * 3] * (coords[2 + 0 * 3] - coords[2 + 3 * 3]) - coords[1 + 0 * 3] * (coords[2 + 1 * 3] - coords[2 + 3 * 3]) - coords[1 + 3 * 3] * (coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume;
472181f196bSVijay Mahadevan       Dx[3] = -(Dx[0] + Dx[1] + Dx[2]);
4739371c9d4SSatish Balay       Dy[0] = (coords[0 + 1 * 3] * (coords[2 + 2 * 3] - coords[2 + 3 * 3]) - coords[0 + 2 * 3] * (coords[2 + 1 * 3] - coords[2 + 3 * 3]) + coords[0 + 3 * 3] * (coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume;
4749371c9d4SSatish Balay       Dy[1] = -(coords[0 + 0 * 3] * (coords[2 + 2 * 3] - coords[2 + 3 * 3]) - coords[0 + 2 * 3] * (coords[2 + 0 * 3] - coords[2 + 3 * 3]) + coords[0 + 3 * 3] * (coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume;
4759371c9d4SSatish Balay       Dy[2] = (coords[0 + 0 * 3] * (coords[2 + 1 * 3] - coords[2 + 3 * 3]) - coords[0 + 1 * 3] * (coords[2 + 0 * 3] - coords[2 + 3 * 3]) + coords[0 + 3 * 3] * (coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume;
476181f196bSVijay Mahadevan       Dy[3] = -(Dy[0] + Dy[1] + Dy[2]);
4779371c9d4SSatish Balay       Dz[0] = (coords[0 + 1 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3]) - coords[0 + 2 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3]) + coords[0 + 3 * 3] * (coords[1 + 2 * 3] - coords[1 + 1 * 3])) / *volume;
4789371c9d4SSatish Balay       Dz[1] = -(coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3]) + coords[0 + 2 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3]) - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 2 * 3])) / *volume;
4799371c9d4SSatish Balay       Dz[2] = (coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3]) + coords[0 + 1 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3]) - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 1 * 3])) / *volume;
480181f196bSVijay Mahadevan       Dz[3] = -(Dz[0] + Dz[1] + Dz[2]);
481181f196bSVijay Mahadevan     }
48263d025dbSVijay Mahadevan 
4832da392ccSBarry Smith     for (j = 0; j < npts; j++) {
484a86ed394SVijay Mahadevan       const PetscInt   offset = j * nverts;
485181f196bSVijay Mahadevan       const PetscReal &r      = quad[j * 3 + 0];
486181f196bSVijay Mahadevan       const PetscReal &s      = quad[j * 3 + 1];
487181f196bSVijay Mahadevan       const PetscReal &t      = quad[j * 3 + 2];
48863d025dbSVijay Mahadevan 
489181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
49063d025dbSVijay Mahadevan 
49163d025dbSVijay Mahadevan       phi[offset + 0] = 1.0 - r - s - t;
49263d025dbSVijay Mahadevan       phi[offset + 1] = r;
49363d025dbSVijay Mahadevan       phi[offset + 2] = s;
49463d025dbSVijay Mahadevan       phi[offset + 3] = t;
49563d025dbSVijay Mahadevan 
496181f196bSVijay Mahadevan       if (phypts) {
497181f196bSVijay Mahadevan         for (i = 0; i < nverts; ++i) {
498181f196bSVijay Mahadevan           const PetscScalar *vertices = coords + i * 3;
499181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
500181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
501181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
502181f196bSVijay Mahadevan         }
503181f196bSVijay Mahadevan       }
504181f196bSVijay Mahadevan 
505181f196bSVijay Mahadevan       /* Now set the derivatives */
50663d025dbSVijay Mahadevan       if (dphidx) {
507181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
508181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
509181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
510181f196bSVijay Mahadevan         dphidx[3 + offset] = Dx[3];
51163d025dbSVijay Mahadevan 
512181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
513181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
514181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
515181f196bSVijay Mahadevan         dphidy[3 + offset] = Dy[3];
51663d025dbSVijay Mahadevan 
517181f196bSVijay Mahadevan         dphidz[0 + offset] = Dz[0];
518181f196bSVijay Mahadevan         dphidz[1 + offset] = Dz[1];
519181f196bSVijay Mahadevan         dphidz[2 + offset] = Dz[2];
520181f196bSVijay Mahadevan         dphidz[3 + offset] = Dz[3];
52163d025dbSVijay Mahadevan       }
52263d025dbSVijay Mahadevan 
52363d025dbSVijay Mahadevan     } /* Tetrahedra -- ends */
52463a3b9bcSJacob Faibussowitsch   } else SETERRQ(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 : %" PetscInt_FMT, nverts);
52563d025dbSVijay Mahadevan   PetscFunctionReturn(0);
52663d025dbSVijay Mahadevan }
52763d025dbSVijay Mahadevan 
528cab5ea25SPierre Jolivet /*@C
52997b73a88SSatish Balay   DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
53063d025dbSVijay Mahadevan 
53163d025dbSVijay Mahadevan   The routine takes the coordinates of the vertices of an element and computes basis functions associated with
53263d025dbSVijay Mahadevan   each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
53363d025dbSVijay Mahadevan 
534d8d19677SJose E. Roman   Input Parameters:
53567b8a455SSatish Balay +  PetscInt  nverts -           the number of element vertices
53667b8a455SSatish Balay .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
53767b8a455SSatish Balay .  PetscInt  npts -             the number of evaluation points (quadrature points)
53867b8a455SSatish Balay -  PetscReal quad[3*npts] -     the evaluation points (quadrature points) in the reference space
53963d025dbSVijay Mahadevan 
540d8d19677SJose E. Roman   Output Parameters:
54167b8a455SSatish Balay +  PetscReal phypts[3*npts] -   the evaluation points (quadrature points) transformed to the physical space
54267b8a455SSatish Balay .  PetscReal jxw[npts] -        the jacobian determinant * quadrature weight necessary for assembling discrete contributions
54367b8a455SSatish Balay .  PetscReal fe_basis[npts] -   the bases values evaluated at the specified quadrature points
54467b8a455SSatish Balay -  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
54563d025dbSVijay Mahadevan 
546edc382c3SSatish Balay   Level: advanced
547edc382c3SSatish Balay 
54863d025dbSVijay Mahadevan @*/
5499371c9d4SSatish Balay PetscErrorCode DMMoabFEMComputeBasis(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature, PetscReal *phypts, PetscReal *jacobian_quadrature_weight_product, PetscReal *fe_basis, PetscReal **fe_basis_derivatives) {
550b21a1e88SVijay Mahadevan   PetscInt         npoints, idim;
55163d025dbSVijay Mahadevan   bool             compute_der;
55263d025dbSVijay Mahadevan   const PetscReal *quadpts, *quadwts;
553181f196bSVijay Mahadevan   PetscReal        jacobian[9], ijacobian[9], volume;
55463d025dbSVijay Mahadevan 
55563d025dbSVijay Mahadevan   PetscFunctionBegin;
55663d025dbSVijay Mahadevan   PetscValidPointer(coordinates, 3);
55778dc7ee3SMatthew G. Knepley   PetscValidHeaderSpecific(quadrature, PETSCQUADRATURE_CLASSID, 4);
55863d025dbSVijay Mahadevan   PetscValidPointer(fe_basis, 7);
55963d025dbSVijay Mahadevan   compute_der = (fe_basis_derivatives != NULL);
56063d025dbSVijay Mahadevan 
56163d025dbSVijay Mahadevan   /* Get the quadrature points and weights for the given quadrature rule */
5629566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts));
56363a3b9bcSJacob Faibussowitsch   PetscCheck(idim == dim, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%" PetscInt_FMT ") vs quadrature (%" PetscInt_FMT ")", idim, dim);
564*48a46eb9SPierre Jolivet   if (jacobian_quadrature_weight_product) PetscCall(PetscArraycpy(jacobian_quadrature_weight_product, quadwts, npoints));
56563d025dbSVijay Mahadevan 
56663d025dbSVijay Mahadevan   switch (dim) {
5679371c9d4SSatish Balay   case 1: PetscCall(Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, npoints, quadpts, phypts, jacobian_quadrature_weight_product, fe_basis, (compute_der ? fe_basis_derivatives[0] : NULL), jacobian, ijacobian, &volume)); break;
56863d025dbSVijay Mahadevan   case 2:
5699371c9d4SSatish Balay     PetscCall(Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, npoints, quadpts, phypts, jacobian_quadrature_weight_product, fe_basis, (compute_der ? fe_basis_derivatives[0] : NULL), (compute_der ? fe_basis_derivatives[1] : NULL), jacobian, ijacobian, &volume));
57063d025dbSVijay Mahadevan     break;
57163d025dbSVijay Mahadevan   case 3:
5729371c9d4SSatish Balay     PetscCall(Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, npoints, quadpts, phypts, jacobian_quadrature_weight_product, fe_basis, (compute_der ? fe_basis_derivatives[0] : NULL), (compute_der ? fe_basis_derivatives[1] : NULL), (compute_der ? fe_basis_derivatives[2] : NULL), jacobian, ijacobian, &volume));
57363d025dbSVijay Mahadevan     break;
5749371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
57563d025dbSVijay Mahadevan   }
57663d025dbSVijay Mahadevan   PetscFunctionReturn(0);
57763d025dbSVijay Mahadevan }
57863d025dbSVijay Mahadevan 
579cab5ea25SPierre Jolivet /*@C
58097b73a88SSatish Balay   DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given
58163d025dbSVijay Mahadevan   dimension and polynomial order (deciphered from number of element vertices).
58263d025dbSVijay Mahadevan 
583d8d19677SJose E. Roman   Input Parameters:
58463d025dbSVijay Mahadevan 
585a2b725a8SWilliam Gropp +  PetscInt  dim   -   the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
586a2b725a8SWilliam Gropp -  PetscInt nverts -   the number of vertices in the physical element
58763d025dbSVijay Mahadevan 
58863d025dbSVijay Mahadevan   Output Parameter:
589a2b725a8SWilliam Gropp .  PetscQuadrature quadrature -  the quadrature object with default settings to integrate polynomials defined over the element
59063d025dbSVijay Mahadevan 
591edc382c3SSatish Balay   Level: advanced
592edc382c3SSatish Balay 
59363d025dbSVijay Mahadevan @*/
5949371c9d4SSatish Balay PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature) {
59563d025dbSVijay Mahadevan   PetscReal *w, *x;
596b21a1e88SVijay Mahadevan   PetscInt   nc = 1;
59763d025dbSVijay Mahadevan 
59863d025dbSVijay Mahadevan   PetscFunctionBegin;
59963d025dbSVijay Mahadevan   /* Create an appropriate quadrature rule to sample basis */
6009371c9d4SSatish Balay   switch (dim) {
60163d025dbSVijay Mahadevan   case 1:
60263d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
6039566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(1, nc, nverts, 0, 1.0, quadrature));
60463d025dbSVijay Mahadevan     break;
60563d025dbSVijay Mahadevan   case 2:
60663d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
60763d025dbSVijay Mahadevan     if (nverts == 3) {
608a86ed394SVijay Mahadevan       const PetscInt order   = 2;
609a86ed394SVijay Mahadevan       const PetscInt npoints = (order == 2 ? 3 : 6);
6109566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(npoints * 2, &x, npoints, &w));
611181f196bSVijay Mahadevan       if (npoints == 3) {
61263d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0;
61363d025dbSVijay Mahadevan         x[3] = x[4] = 2.0 / 3.0;
61463d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 1.0 / 3.0;
6152da392ccSBarry Smith       } else if (npoints == 6) {
61663d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = 0.44594849091597;
61763d025dbSVijay Mahadevan         x[3] = x[4] = 0.10810301816807;
61863d025dbSVijay Mahadevan         x[5]        = 0.44594849091597;
61963d025dbSVijay Mahadevan         x[6] = x[7] = x[8] = 0.09157621350977;
62063d025dbSVijay Mahadevan         x[9] = x[10] = 0.81684757298046;
62163d025dbSVijay Mahadevan         x[11]        = 0.09157621350977;
62263d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 0.22338158967801;
623181f196bSVijay Mahadevan         w[3] = w[4] = w[5] = 0.10995174365532;
62463a3b9bcSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %" PetscInt_FMT, npoints);
6259566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
6269566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetOrder(*quadrature, order));
6279566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
6289566063dSJacob Faibussowitsch       /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
6291baa6e33SBarry Smith     } else PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
63063d025dbSVijay Mahadevan     break;
63163d025dbSVijay Mahadevan   case 3:
63263d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
63363d025dbSVijay Mahadevan     if (nverts == 4) {
634a86ed394SVijay Mahadevan       PetscInt       order;
635a86ed394SVijay Mahadevan       const PetscInt npoints = 4; // Choose between 4 and 10
6369566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(npoints * 3, &x, npoints, &w));
637181f196bSVijay Mahadevan       if (npoints == 4) { /*  KEAST1, K1,  N=4, O=4 */
6389371c9d4SSatish Balay         const PetscReal x_4[12] = {0.5854101966249685, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105,
6399371c9d4SSatish Balay                                    0.1381966011250105, 0.1381966011250105, 0.5854101966249685, 0.1381966011250105, 0.5854101966249685, 0.1381966011250105};
6409566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(x, x_4, 12));
641181f196bSVijay Mahadevan 
642181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0;
643181f196bSVijay Mahadevan         order                     = 4;
6442da392ccSBarry Smith       } else if (npoints == 10) { /*  KEAST3, K3  N=10, O=10 */
6459371c9d4SSatish Balay         const PetscReal x_10[30] = {0.5684305841968444, 0.1438564719343852, 0.1438564719343852, 0.1438564719343852, 0.1438564719343852, 0.1438564719343852, 0.1438564719343852, 0.1438564719343852, 0.5684305841968444, 0.1438564719343852,
6469371c9d4SSatish Balay                                     0.5684305841968444, 0.1438564719343852, 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 0.5000000000000000, 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 0.5000000000000000,
6479371c9d4SSatish Balay                                     0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.5000000000000000};
6489566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(x, x_10, 30));
649181f196bSVijay Mahadevan 
650181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 0.2177650698804054;
651181f196bSVijay Mahadevan         w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631;
652181f196bSVijay Mahadevan         order                                   = 10;
65363a3b9bcSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %" PetscInt_FMT, npoints);
6549566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
6559566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetOrder(*quadrature, order));
6569566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
6579566063dSJacob Faibussowitsch       /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
6581baa6e33SBarry Smith     } else PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
65963d025dbSVijay Mahadevan     break;
6609371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
66163d025dbSVijay Mahadevan   }
66263d025dbSVijay Mahadevan   PetscFunctionReturn(0);
66363d025dbSVijay Mahadevan }
66463d025dbSVijay Mahadevan 
665181f196bSVijay Mahadevan /* Compute Jacobians */
6669371c9d4SSatish Balay PetscErrorCode ComputeJacobian_Internal(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts, PetscReal *jacobian, PetscReal *ijacobian, PetscReal *dvolume) {
667a86ed394SVijay Mahadevan   PetscInt  i;
6682417220eSVijay Mahadevan   PetscReal volume = 1.0;
669181f196bSVijay Mahadevan 
670181f196bSVijay Mahadevan   PetscFunctionBegin;
671181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
672181f196bSVijay Mahadevan   PetscValidPointer(quad, 4);
673181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 5);
6749566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(jacobian, dim * dim));
675*48a46eb9SPierre Jolivet   if (ijacobian) PetscCall(PetscArrayzero(ijacobian, dim * dim));
676*48a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, /*npts=1 * */ 3));
677181f196bSVijay Mahadevan 
678181f196bSVijay Mahadevan   if (dim == 1) {
679181f196bSVijay Mahadevan     const PetscReal &r = quad[0];
680181f196bSVijay Mahadevan     if (nverts == 2) { /* Linear Edge */
681181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2] = {-1.0, 1.0};
682181f196bSVijay Mahadevan 
683181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
684181f196bSVijay Mahadevan         const PetscReal *vertices = coordinates + i * 3;
685181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
686181f196bSVijay Mahadevan       }
6872da392ccSBarry Smith     } else if (nverts == 3) { /* Quadratic Edge */
688181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3] = {4 * r - 3.0, 4 * (1.0 - 2.0 * r), 4.0 * r - 1.0};
689181f196bSVijay Mahadevan 
690181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
691181f196bSVijay Mahadevan         const PetscReal *vertices = coordinates + i * 3;
692181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
693181f196bSVijay Mahadevan       }
6941baa6e33SBarry Smith     } else SETERRQ(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 : %" PetscInt_FMT, nverts);
695181f196bSVijay Mahadevan 
696181f196bSVijay Mahadevan     if (ijacobian) {
697181f196bSVijay Mahadevan       /* invert the jacobian */
698181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
699181f196bSVijay Mahadevan     }
7002da392ccSBarry Smith   } else if (dim == 2) {
701181f196bSVijay Mahadevan     if (nverts == 4) { /* Linear Quadrangle */
702181f196bSVijay Mahadevan       const PetscReal &r = quad[0];
703181f196bSVijay Mahadevan       const PetscReal &s = quad[1];
704181f196bSVijay Mahadevan 
705181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = {-1.0 + s, 1.0 - s, s, -s};
706181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = {-1.0 + r, -r, r, 1.0 - r};
707181f196bSVijay Mahadevan 
708181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
709181f196bSVijay Mahadevan         const PetscReal *vertices = coordinates + i * 3;
710181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
711181f196bSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
712181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
713181f196bSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
714181f196bSVijay Mahadevan       }
7152da392ccSBarry Smith     } else if (nverts == 3) { /* Linear triangle */
716181f196bSVijay Mahadevan       const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1];
717181f196bSVijay Mahadevan 
718181f196bSVijay Mahadevan       /* Jacobian is constant */
7199371c9d4SSatish Balay       jacobian[0] = (coordinates[0 * 3 + 0] - x2);
7209371c9d4SSatish Balay       jacobian[1] = (coordinates[1 * 3 + 0] - x2);
7219371c9d4SSatish Balay       jacobian[2] = (coordinates[0 * 3 + 1] - y2);
7229371c9d4SSatish Balay       jacobian[3] = (coordinates[1 * 3 + 1] - y2);
72363a3b9bcSJacob Faibussowitsch     } else SETERRQ(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 : %" PetscInt_FMT, nverts);
724181f196bSVijay Mahadevan 
725181f196bSVijay Mahadevan     /* invert the jacobian */
726*48a46eb9SPierre Jolivet     if (ijacobian) PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume));
7272da392ccSBarry Smith   } else {
728181f196bSVijay Mahadevan     if (nverts == 8) { /* Linear Hexahedra */
729181f196bSVijay Mahadevan       const PetscReal &r          = quad[0];
730181f196bSVijay Mahadevan       const PetscReal &s          = quad[1];
731181f196bSVijay Mahadevan       const PetscReal &t          = quad[2];
7329371c9d4SSatish Balay       const PetscReal  dNi_dxi[8] = {-(1.0 - s) * (1.0 - t), (1.0 - s) * (1.0 - t), (s) * (1.0 - t), -(s) * (1.0 - t), -(1.0 - s) * (t), (1.0 - s) * (t), (s) * (t), -(s) * (t)};
733181f196bSVijay Mahadevan 
7349371c9d4SSatish Balay       const PetscReal dNi_deta[8] = {-(1.0 - r) * (1.0 - t), -(r) * (1.0 - t), (r) * (1.0 - t), (1.0 - r) * (1.0 - t), -(1.0 - r) * (t), -(r) * (t), (r) * (t), (1.0 - r) * (t)};
735181f196bSVijay Mahadevan 
7369371c9d4SSatish Balay       const PetscReal dNi_dzeta[8] = {-(1.0 - r) * (1.0 - s), -(r) * (1.0 - s), -(r) * (s), -(1.0 - r) * (s), (1.0 - r) * (1.0 - s), (r) * (1.0 - s), (r) * (s), (1.0 - r) * (s)};
737a86ed394SVijay Mahadevan 
738181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
739181f196bSVijay Mahadevan         const PetscReal *vertex = coordinates + i * 3;
740181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertex[0];
741181f196bSVijay Mahadevan         jacobian[3] += dNi_dxi[i] * vertex[1];
742181f196bSVijay Mahadevan         jacobian[6] += dNi_dxi[i] * vertex[2];
743181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertex[0];
744181f196bSVijay Mahadevan         jacobian[4] += dNi_deta[i] * vertex[1];
745181f196bSVijay Mahadevan         jacobian[7] += dNi_deta[i] * vertex[2];
746181f196bSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
747181f196bSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
748181f196bSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
749181f196bSVijay Mahadevan       }
7502da392ccSBarry Smith     } else if (nverts == 4) { /* Linear Tetrahedra */
751181f196bSVijay Mahadevan       const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2];
752181f196bSVijay Mahadevan 
753181f196bSVijay Mahadevan       /* compute the jacobian */
7549371c9d4SSatish Balay       jacobian[0] = coordinates[1 * 3 + 0] - x0;
7559371c9d4SSatish Balay       jacobian[1] = coordinates[2 * 3 + 0] - x0;
7569371c9d4SSatish Balay       jacobian[2] = coordinates[3 * 3 + 0] - x0;
7579371c9d4SSatish Balay       jacobian[3] = coordinates[1 * 3 + 1] - y0;
7589371c9d4SSatish Balay       jacobian[4] = coordinates[2 * 3 + 1] - y0;
7599371c9d4SSatish Balay       jacobian[5] = coordinates[3 * 3 + 1] - y0;
7609371c9d4SSatish Balay       jacobian[6] = coordinates[1 * 3 + 2] - z0;
7619371c9d4SSatish Balay       jacobian[7] = coordinates[2 * 3 + 2] - z0;
7629371c9d4SSatish Balay       jacobian[8] = coordinates[3 * 3 + 2] - z0;
76363a3b9bcSJacob Faibussowitsch     } else SETERRQ(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 : %" PetscInt_FMT, nverts);
764181f196bSVijay Mahadevan 
765181f196bSVijay Mahadevan     if (ijacobian) {
766181f196bSVijay Mahadevan       /* invert the jacobian */
7679566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume));
768181f196bSVijay Mahadevan     }
769181f196bSVijay Mahadevan   }
77008401ef6SPierre Jolivet   PetscCheck(volume >= 1e-12, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Element has zero volume: %g. Degenerate element or invalid connectivity", volume);
771a86ed394SVijay Mahadevan   if (dvolume) *dvolume = volume;
772181f196bSVijay Mahadevan   PetscFunctionReturn(0);
773181f196bSVijay Mahadevan }
774181f196bSVijay Mahadevan 
7759371c9d4SSatish Balay PetscErrorCode FEMComputeBasis_JandF(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quadrature, PetscReal *phypts, PetscReal *phibasis, PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) {
776181f196bSVijay Mahadevan   PetscFunctionBegin;
777181f196bSVijay Mahadevan   switch (dim) {
7789371c9d4SSatish Balay   case 1: PetscCall(Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, jacobian, ijacobian, volume)); break;
7799371c9d4SSatish Balay   case 2: PetscCall(Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume)); break;
7809371c9d4SSatish Balay   case 3: PetscCall(Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume)); break;
7819371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
782181f196bSVijay Mahadevan   }
783181f196bSVijay Mahadevan   PetscFunctionReturn(0);
784181f196bSVijay Mahadevan }
785181f196bSVijay Mahadevan 
786cab5ea25SPierre Jolivet /*@C
78797b73a88SSatish Balay   DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the
788a86ed394SVijay Mahadevan   canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration,
789a86ed394SVijay Mahadevan   the basis function at the parametric point is also evaluated optionally.
790a86ed394SVijay Mahadevan 
791d8d19677SJose E. Roman   Input Parameters:
792a2b725a8SWilliam Gropp +  PetscInt  dim -         the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
793a2b725a8SWilliam Gropp .  PetscInt nverts -       the number of vertices in the physical element
794a2b725a8SWilliam Gropp .  PetscReal coordinates - the coordinates of vertices in the physical element
795a2b725a8SWilliam Gropp -  PetscReal[3] xphy -     the coordinates of physical point for which natural coordinates (in reference frame) are sought
796a86ed394SVijay Mahadevan 
797d8d19677SJose E. Roman   Output Parameters:
798a2b725a8SWilliam Gropp +  PetscReal[3] natparam - the natural coordinates (in reference frame) corresponding to xphy
799a2b725a8SWilliam Gropp -  PetscReal[nverts] phi - the basis functions evaluated at the natural coordinates (natparam)
800a86ed394SVijay Mahadevan 
801edc382c3SSatish Balay   Level: advanced
802edc382c3SSatish Balay 
803a86ed394SVijay Mahadevan @*/
8049371c9d4SSatish Balay PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *xphy, PetscReal *natparam, PetscReal *phi) {
805a86ed394SVijay Mahadevan   /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */
806181f196bSVijay Mahadevan   const PetscReal tol            = 1.0e-10;
807181f196bSVijay Mahadevan   const PetscInt  max_iterations = 10;
808181f196bSVijay Mahadevan   const PetscReal error_tol_sqr  = tol * tol;
809181f196bSVijay Mahadevan   PetscReal       phibasis[8], jacobian[9], ijacobian[9], volume;
810181f196bSVijay Mahadevan   PetscReal       phypts[3] = {0.0, 0.0, 0.0};
811181f196bSVijay Mahadevan   PetscReal       delta[3]  = {0.0, 0.0, 0.0};
812181f196bSVijay Mahadevan   PetscInt        iters     = 0;
813181f196bSVijay Mahadevan   PetscReal       error     = 1.0;
814181f196bSVijay Mahadevan 
815181f196bSVijay Mahadevan   PetscFunctionBegin;
816181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
817181f196bSVijay Mahadevan   PetscValidPointer(xphy, 4);
818181f196bSVijay Mahadevan   PetscValidPointer(natparam, 5);
819181f196bSVijay Mahadevan 
8209566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(jacobian, dim * dim));
8219566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(ijacobian, dim * dim));
8229566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phibasis, nverts));
823181f196bSVijay Mahadevan 
824a86ed394SVijay Mahadevan   /* zero initial guess */
825a86ed394SVijay Mahadevan   natparam[0] = natparam[1] = natparam[2] = 0.0;
826181f196bSVijay Mahadevan 
827a86ed394SVijay Mahadevan   /* Compute delta = evaluate( xi) - x */
8289566063dSJacob Faibussowitsch   PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
829a86ed394SVijay Mahadevan 
830a86ed394SVijay Mahadevan   error = 0.0;
831a86ed394SVijay Mahadevan   switch (dim) {
8329371c9d4SSatish Balay   case 3: delta[2] = phypts[2] - xphy[2]; error += (delta[2] * delta[2]);
8339371c9d4SSatish Balay   case 2: delta[1] = phypts[1] - xphy[1]; error += (delta[1] * delta[1]);
834a86ed394SVijay Mahadevan   case 1:
835a86ed394SVijay Mahadevan     delta[0] = phypts[0] - xphy[0];
836a86ed394SVijay Mahadevan     error += (delta[0] * delta[0]);
837a86ed394SVijay Mahadevan     break;
838a86ed394SVijay Mahadevan   }
839a86ed394SVijay Mahadevan 
840181f196bSVijay Mahadevan   while (error > error_tol_sqr) {
8417a46b595SBarry Smith     PetscCheck(++iters <= max_iterations, 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]);
842181f196bSVijay Mahadevan 
843181f196bSVijay Mahadevan     /* Compute natparam -= J.inverse() * delta */
844181f196bSVijay Mahadevan     switch (dim) {
8459371c9d4SSatish Balay     case 1: natparam[0] -= ijacobian[0] * delta[0]; break;
846181f196bSVijay Mahadevan     case 2:
847181f196bSVijay Mahadevan       natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1];
848181f196bSVijay Mahadevan       natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1];
849181f196bSVijay Mahadevan       break;
850181f196bSVijay Mahadevan     case 3:
851181f196bSVijay Mahadevan       natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2];
852181f196bSVijay Mahadevan       natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2];
853181f196bSVijay Mahadevan       natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2];
854181f196bSVijay Mahadevan       break;
855181f196bSVijay Mahadevan     }
856181f196bSVijay Mahadevan 
857181f196bSVijay Mahadevan     /* Compute delta = evaluate( xi) - x */
8589566063dSJacob Faibussowitsch     PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
859181f196bSVijay Mahadevan 
860a86ed394SVijay Mahadevan     error = 0.0;
861a86ed394SVijay Mahadevan     switch (dim) {
8629371c9d4SSatish Balay     case 3: delta[2] = phypts[2] - xphy[2]; error += (delta[2] * delta[2]);
8639371c9d4SSatish Balay     case 2: delta[1] = phypts[1] - xphy[1]; error += (delta[1] * delta[1]);
864a86ed394SVijay Mahadevan     case 1:
865a86ed394SVijay Mahadevan       delta[0] = phypts[0] - xphy[0];
866a86ed394SVijay Mahadevan       error += (delta[0] * delta[0]);
867a86ed394SVijay Mahadevan       break;
868a86ed394SVijay Mahadevan     }
869181f196bSVijay Mahadevan   }
870*48a46eb9SPierre Jolivet   if (phi) PetscCall(PetscArraycpy(phi, phibasis, nverts));
871181f196bSVijay Mahadevan   PetscFunctionReturn(0);
872181f196bSVijay Mahadevan }
873