xref: /petsc/src/dm/impls/moab/dmmbfem.cxx (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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;
22362febeeSStefano 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;
47362febeeSStefano 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;
93362febeeSStefano 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 
112d8d19677SJose E. Roman   Input Parameters:
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 
118d8d19677SJose E. Roman   Output Parameters:
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
1226b867d5aSJose E. Roman .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
12367b8a455SSatish Balay .  jacobian -                  jacobian
12467b8a455SSatish Balay .  ijacobian -                 ijacobian
12567b8a455SSatish Balay -  volume -                    volume
12663d025dbSVijay Mahadevan 
127edc382c3SSatish Balay   Level: advanced
128edc382c3SSatish Balay 
12963d025dbSVijay Mahadevan @*/
13063d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_1D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
131181f196bSVijay Mahadevan                                                   PetscReal *jxw, PetscReal *phi, PetscReal *dphidx,
132181f196bSVijay Mahadevan                                                   PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
13363d025dbSVijay Mahadevan {
13463d025dbSVijay Mahadevan   int             i, j;
13563d025dbSVijay Mahadevan 
136181f196bSVijay Mahadevan   PetscFunctionBegin;
137181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 9);
138181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 10);
139181f196bSVijay Mahadevan   PetscValidPointer(volume, 11);
140181f196bSVijay Mahadevan   if (phypts) {
1419566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(phypts, npts * 3));
142181f196bSVijay Mahadevan   }
14363d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
1449566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
14563d025dbSVijay Mahadevan   }
14663d025dbSVijay Mahadevan   if (nverts == 2) { /* Linear Edge */
14763d025dbSVijay Mahadevan 
1482da392ccSBarry Smith     for (j = 0; j < npts; j++) {
149a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
150181f196bSVijay Mahadevan       const PetscReal r = quad[j];
15163d025dbSVijay Mahadevan 
152181f196bSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r);
153181f196bSVijay Mahadevan       phi[1 + offset] = (       r);
15463d025dbSVijay Mahadevan 
155181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 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];
161181f196bSVijay Mahadevan         if (phypts) {
162181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
163181f196bSVijay Mahadevan         }
16463d025dbSVijay Mahadevan       }
16563d025dbSVijay Mahadevan 
16663d025dbSVijay Mahadevan       /* invert the jacobian */
167181f196bSVijay Mahadevan       *volume = jacobian[0];
168181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
169181f196bSVijay Mahadevan       jxw[j] *= *volume;
17063d025dbSVijay Mahadevan 
17163d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
17263d025dbSVijay Mahadevan       for (i = 0; i < nverts; i++) {
173181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
17463d025dbSVijay Mahadevan       }
17563d025dbSVijay Mahadevan     }
1762da392ccSBarry Smith   } else if (nverts == 3) { /* Quadratic Edge */
17763d025dbSVijay Mahadevan 
1782da392ccSBarry Smith     for (j = 0; j < npts; j++) {
179a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
180181f196bSVijay Mahadevan       const PetscReal r = quad[j];
18163d025dbSVijay Mahadevan 
182181f196bSVijay Mahadevan       phi[0 + offset] = 1.0 + r * ( 2.0 * r - 3.0);
183181f196bSVijay Mahadevan       phi[1 + offset] = 4.0 * r * ( 1.0 - r);
184181f196bSVijay Mahadevan       phi[2 + offset] = r * ( 2.0 * r - 1.0);
18563d025dbSVijay Mahadevan 
186181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r), 4.0 * r - 1.0};
18763d025dbSVijay Mahadevan 
188181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
18963d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
190181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
191181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
192181f196bSVijay Mahadevan         if (phypts) {
193181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
194181f196bSVijay Mahadevan         }
19563d025dbSVijay Mahadevan       }
19663d025dbSVijay Mahadevan 
19763d025dbSVijay Mahadevan       /* invert the jacobian */
198181f196bSVijay Mahadevan       *volume = jacobian[0];
199181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
200181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
20163d025dbSVijay Mahadevan 
20263d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
20363d025dbSVijay Mahadevan       for (i = 0; i < nverts; i++) {
204181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
20563d025dbSVijay Mahadevan       }
20663d025dbSVijay Mahadevan     }
20798921bdaSJacob 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 : %D", nverts);
20863d025dbSVijay Mahadevan   PetscFunctionReturn(0);
20963d025dbSVijay Mahadevan }
21063d025dbSVijay Mahadevan 
211cab5ea25SPierre Jolivet /*@C
21297b73a88SSatish Balay   Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.
21363d025dbSVijay Mahadevan 
21463d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a quadrangle or triangle.
21563d025dbSVijay Mahadevan 
21663d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
21763d025dbSVijay Mahadevan   and their derivatives with respect to X and Y.
21863d025dbSVijay Mahadevan 
21963d025dbSVijay Mahadevan   Notes:
22063d025dbSVijay Mahadevan 
22163d025dbSVijay Mahadevan   Example Physical Element (QUAD4)
222a2b725a8SWilliam Gropp .vb
22363d025dbSVijay Mahadevan     4------3        s
22463d025dbSVijay Mahadevan     |      |        |
22563d025dbSVijay Mahadevan     |      |        |
22663d025dbSVijay Mahadevan     |      |        |
22763d025dbSVijay Mahadevan     1------2        0-------r
228a2b725a8SWilliam Gropp .ve
22963d025dbSVijay Mahadevan 
230d8d19677SJose E. Roman   Input Parameters:
231a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
232a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
233a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
234a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
23563d025dbSVijay Mahadevan 
236d8d19677SJose E. Roman   Output Parameters:
237a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
238a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
239a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
240a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
2416b867d5aSJose E. Roman .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
24267b8a455SSatish Balay .  jacobian -                  jacobian
24367b8a455SSatish Balay .  ijacobian -                 ijacobian
24467b8a455SSatish Balay -  volume -                    volume
24563d025dbSVijay Mahadevan 
246edc382c3SSatish Balay   Level: advanced
247edc382c3SSatish Balay 
24863d025dbSVijay Mahadevan @*/
24963d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_2D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
250181f196bSVijay Mahadevan                                                   PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy,
251181f196bSVijay Mahadevan                                                   PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
25263d025dbSVijay Mahadevan {
253a86ed394SVijay Mahadevan   PetscInt       i, j, k;
25463d025dbSVijay Mahadevan 
25563d025dbSVijay Mahadevan   PetscFunctionBegin;
256181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 10);
257181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 11);
258181f196bSVijay Mahadevan   PetscValidPointer(volume, 12);
2599566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phi, npts));
260181f196bSVijay Mahadevan   if (phypts) {
2619566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(phypts, npts * 3));
262181f196bSVijay Mahadevan   }
26363d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
2649566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
2659566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidy, npts * nverts));
26663d025dbSVijay Mahadevan   }
26763d025dbSVijay Mahadevan   if (nverts == 4) { /* Linear Quadrangle */
26863d025dbSVijay Mahadevan 
26963d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
27063d025dbSVijay Mahadevan     {
271a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
272181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
273181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
27463d025dbSVijay Mahadevan 
27563d025dbSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r) * ( 1.0 - s);
27663d025dbSVijay Mahadevan       phi[1 + offset] =         r   * ( 1.0 - s);
27763d025dbSVijay Mahadevan       phi[2 + offset] =         r   *         s;
27863d025dbSVijay Mahadevan       phi[3 + offset] = ( 1.0 - r) *         s;
27963d025dbSVijay Mahadevan 
280181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
281181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
28263d025dbSVijay Mahadevan 
2839566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(jacobian, 4));
2849566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ijacobian, 4));
28563d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
286181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
28763d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
28863d025dbSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
28963d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
29063d025dbSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
291181f196bSVijay Mahadevan         if (phypts) {
292181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
293181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
294181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
295181f196bSVijay Mahadevan         }
29663d025dbSVijay Mahadevan       }
29763d025dbSVijay Mahadevan 
29863d025dbSVijay Mahadevan       /* invert the jacobian */
2999566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
300*08401ef6SPierre Jolivet       PetscCheck(*volume >= 1e-12,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Quadrangular element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
30163d025dbSVijay Mahadevan 
302181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
30363d025dbSVijay Mahadevan 
304181f196bSVijay Mahadevan       /*  Let us compute the bases derivatives by scaling with inverse jacobian. */
305181f196bSVijay Mahadevan       if (dphidx) {
30663d025dbSVijay Mahadevan         for (i = 0; i < nverts; i++) {
307a86ed394SVijay Mahadevan           for (k = 0; k < 2; ++k) {
30863d025dbSVijay Mahadevan             if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0];
30963d025dbSVijay Mahadevan             if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1];
310181f196bSVijay Mahadevan           } /* for k=[0..2) */
311181f196bSVijay Mahadevan         } /* for i=[0..nverts) */
312181f196bSVijay Mahadevan       } /* if (dphidx) */
31363d025dbSVijay Mahadevan     }
3142da392ccSBarry Smith   } else if (nverts == 3) { /* Linear triangle */
3152da392ccSBarry Smith     const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1];
31663d025dbSVijay Mahadevan 
3179566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(jacobian, 4));
3189566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ijacobian, 4));
31963d025dbSVijay Mahadevan 
32063d025dbSVijay Mahadevan     /* Jacobian is constant */
321181f196bSVijay Mahadevan     jacobian[0] = (coords[0 * 3 + 0] - x2); jacobian[1] = (coords[1 * 3 + 0] - x2);
322181f196bSVijay Mahadevan     jacobian[2] = (coords[0 * 3 + 1] - y2); jacobian[3] = (coords[1 * 3 + 1] - y2);
32363d025dbSVijay Mahadevan 
32463d025dbSVijay Mahadevan     /* invert the jacobian */
3259566063dSJacob Faibussowitsch     PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
326*08401ef6SPierre 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);
327181f196bSVijay Mahadevan 
328181f196bSVijay Mahadevan     const PetscReal Dx[3] = { ijacobian[0], ijacobian[2], - ijacobian[0] - ijacobian[2] };
329181f196bSVijay Mahadevan     const PetscReal Dy[3] = { ijacobian[1], ijacobian[3], - ijacobian[1] - ijacobian[3] };
33063d025dbSVijay Mahadevan 
3312da392ccSBarry Smith     for (j = 0; j < npts; j++) {
332a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
333181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
334181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
33563d025dbSVijay Mahadevan 
336181f196bSVijay Mahadevan       if (jxw) jxw[j] *= 0.5;
33763d025dbSVijay Mahadevan 
338181f196bSVijay Mahadevan       /* const PetscReal Ni[3]  = { r, s, 1.0 - r - s }; */
339181f196bSVijay Mahadevan       const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s;
340181f196bSVijay Mahadevan       const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s;
341181f196bSVijay Mahadevan       if (phypts) {
342181f196bSVijay Mahadevan         phypts[offset + 0] = phipts_x;
343181f196bSVijay Mahadevan         phypts[offset + 1] = phipts_y;
344181f196bSVijay Mahadevan       }
34563d025dbSVijay Mahadevan 
346110fc3b0SBarry Smith       /* \phi_0 = (b.y - c.y) x + (b.x - c.x) y + c.x b.y - b.x c.y */
347181f196bSVijay Mahadevan       phi[0 + offset] = (  ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2));
348110fc3b0SBarry Smith       /* \phi_1 = (c.y - a.y) x + (a.x - c.x) y + c.x a.y - a.x c.y */
349181f196bSVijay Mahadevan       phi[1 + offset] = (  ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2));
35063d025dbSVijay Mahadevan       phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset];
35163d025dbSVijay Mahadevan 
35263d025dbSVijay Mahadevan       if (dphidx) {
353181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
354181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
355181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
35663d025dbSVijay Mahadevan       }
35763d025dbSVijay Mahadevan 
35863d025dbSVijay Mahadevan       if (dphidy) {
359181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
360181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
361181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
36263d025dbSVijay Mahadevan       }
36363d025dbSVijay Mahadevan 
36463d025dbSVijay Mahadevan     }
36598921bdaSJacob 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 : %D", nverts);
36663d025dbSVijay Mahadevan   PetscFunctionReturn(0);
36763d025dbSVijay Mahadevan }
36863d025dbSVijay Mahadevan 
369cab5ea25SPierre Jolivet /*@C
37097b73a88SSatish Balay   Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element.
37163d025dbSVijay Mahadevan 
37263d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a hexahedra or tetrahedra.
37363d025dbSVijay Mahadevan 
37463d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
37563d025dbSVijay Mahadevan   and their derivatives with respect to X, Y and Z.
37663d025dbSVijay Mahadevan 
37763d025dbSVijay Mahadevan   Notes:
37863d025dbSVijay Mahadevan 
37963d025dbSVijay Mahadevan   Example Physical Element (HEX8)
380a2b725a8SWilliam Gropp .vb
38163d025dbSVijay Mahadevan       8------7
38263d025dbSVijay Mahadevan      /|     /|        t  s
38363d025dbSVijay Mahadevan     5------6 |        | /
38463d025dbSVijay Mahadevan     | |    | |        |/
38563d025dbSVijay Mahadevan     | 4----|-3        0-------r
38663d025dbSVijay Mahadevan     |/     |/
38763d025dbSVijay Mahadevan     1------2
388a2b725a8SWilliam Gropp .ve
38963d025dbSVijay Mahadevan 
390d8d19677SJose E. Roman   Input Parameters:
391a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
392a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
393a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
394a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
39563d025dbSVijay Mahadevan 
396d8d19677SJose E. Roman   Output Parameters:
397a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
398a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
399a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
400a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
401a2b725a8SWilliam Gropp .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
4026b867d5aSJose E. Roman .  PetscReal dphidz[npts] -    the derivative of the bases wrt Z-direction evaluated at the specified quadrature points
40367b8a455SSatish Balay .  jacobian -                  jacobian
40467b8a455SSatish Balay .  ijacobian -                 ijacobian
40567b8a455SSatish Balay -  volume -                    volume
40663d025dbSVijay Mahadevan 
407edc382c3SSatish Balay   Level: advanced
408edc382c3SSatish Balay 
40963d025dbSVijay Mahadevan @*/
41063d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_3D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
411181f196bSVijay Mahadevan                                                   PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz,
412181f196bSVijay Mahadevan                                                   PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
41363d025dbSVijay Mahadevan {
414a86ed394SVijay Mahadevan   PetscInt       i, j, k;
41563d025dbSVijay Mahadevan 
41663d025dbSVijay Mahadevan   PetscFunctionBegin;
417181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 11);
418181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 12);
419181f196bSVijay Mahadevan   PetscValidPointer(volume, 13);
4202da392ccSBarry Smith 
4219566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phi, npts));
422181f196bSVijay Mahadevan   if (phypts) {
4239566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(phypts, npts * 3));
424181f196bSVijay Mahadevan   }
42563d025dbSVijay Mahadevan   if (dphidx) {
4269566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
4279566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidy, npts * nverts));
4289566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidz, npts * nverts));
42963d025dbSVijay Mahadevan   }
43063d025dbSVijay Mahadevan 
43163d025dbSVijay Mahadevan   if (nverts == 8) { /* Linear Hexahedra */
43263d025dbSVijay Mahadevan 
4332da392ccSBarry Smith     for (j = 0; j < npts; j++) {
434a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
435181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
436181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
437181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
43863d025dbSVijay Mahadevan 
439a86ed394SVijay Mahadevan       phi[offset + 0] = ( 1.0 - r) * ( 1.0 - s) * ( 1.0 - t); /* 0,0,0 */
440a86ed394SVijay Mahadevan       phi[offset + 1] = (       r) * ( 1.0 - s) * ( 1.0 - t); /* 1,0,0 */
441a86ed394SVijay Mahadevan       phi[offset + 2] = (       r) * (       s) * ( 1.0 - t); /* 1,1,0 */
442a86ed394SVijay Mahadevan       phi[offset + 3] = ( 1.0 - r) * (       s) * ( 1.0 - t); /* 0,1,0 */
443a86ed394SVijay Mahadevan       phi[offset + 4] = ( 1.0 - r) * ( 1.0 - s) * (       t); /* 0,0,1 */
444a86ed394SVijay Mahadevan       phi[offset + 5] = (       r) * ( 1.0 - s) * (       t); /* 1,0,1 */
445a86ed394SVijay Mahadevan       phi[offset + 6] = (       r) * (       s) * (       t); /* 1,1,1 */
446a86ed394SVijay Mahadevan       phi[offset + 7] = ( 1.0 - r) * (       s) * (       t); /* 0,1,1 */
44763d025dbSVijay Mahadevan 
448181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s) * ( 1.0 - t),
44963d025dbSVijay Mahadevan                                        ( 1.0 - s) * ( 1.0 - t),
450a86ed394SVijay Mahadevan                                        (       s) * ( 1.0 - t),
451a86ed394SVijay Mahadevan                                      - (       s) * ( 1.0 - t),
452a86ed394SVijay Mahadevan                                      - ( 1.0 - s) * (       t),
453a86ed394SVijay Mahadevan                                        ( 1.0 - s) * (       t),
454a86ed394SVijay Mahadevan                                        (       s) * (       t),
455a86ed394SVijay Mahadevan                                      - (       s) * (       t)
45663d025dbSVijay Mahadevan                                     };
45763d025dbSVijay Mahadevan 
458181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r) * ( 1.0 - t),
459a86ed394SVijay Mahadevan                                        - (       r) * ( 1.0 - t),
460a86ed394SVijay Mahadevan                                          (       r) * ( 1.0 - t),
46163d025dbSVijay Mahadevan                                          ( 1.0 - r) * ( 1.0 - t),
462a86ed394SVijay Mahadevan                                        - ( 1.0 - r) * (       t),
463a86ed394SVijay Mahadevan                                        - (       r) * (       t),
464a86ed394SVijay Mahadevan                                          (       r) * (       t),
465a86ed394SVijay Mahadevan                                          ( 1.0 - r) * (       t)
46663d025dbSVijay Mahadevan                                       };
46763d025dbSVijay Mahadevan 
468181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r) * ( 1.0 - s),
469a86ed394SVijay Mahadevan                                         - (       r) * ( 1.0 - s),
470a86ed394SVijay Mahadevan                                         - (       r) * (       s),
471a86ed394SVijay Mahadevan                                         - ( 1.0 - r) * (       s),
47263d025dbSVijay Mahadevan                                           ( 1.0 - r) * ( 1.0 - s),
473a86ed394SVijay Mahadevan                                           (       r) * ( 1.0 - s),
474a86ed394SVijay Mahadevan                                           (       r) * (       s),
475a86ed394SVijay Mahadevan                                           ( 1.0 - r) * (       s)
47663d025dbSVijay Mahadevan                                      };
47763d025dbSVijay Mahadevan 
4789566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(jacobian, 9));
4799566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ijacobian, 9));
48063d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
481181f196bSVijay Mahadevan         const PetscReal* vertex = coords + i * 3;
48263d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
48363d025dbSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
48463d025dbSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
48563d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
48663d025dbSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
48763d025dbSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
48863d025dbSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
48963d025dbSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
49063d025dbSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
491181f196bSVijay Mahadevan         if (phypts) {
492181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertex[0];
493181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertex[1];
494181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertex[2];
495181f196bSVijay Mahadevan         }
49663d025dbSVijay Mahadevan       }
49763d025dbSVijay Mahadevan 
49863d025dbSVijay Mahadevan       /* invert the jacobian */
4999566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
500*08401ef6SPierre Jolivet       PetscCheck(*volume >= 1e-8,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Hexahedral element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
50163d025dbSVijay Mahadevan 
502a86ed394SVijay Mahadevan       if (jxw) jxw[j] *= (*volume);
50363d025dbSVijay Mahadevan 
50463d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
50563d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
506a86ed394SVijay Mahadevan         for (k = 0; k < 3; ++k) {
50763d025dbSVijay Mahadevan           if (dphidx) dphidx[i + offset] += dNi_dxi[i]   * ijacobian[0 * 3 + k];
50863d025dbSVijay Mahadevan           if (dphidy) dphidy[i + offset] += dNi_deta[i]  * ijacobian[1 * 3 + k];
50963d025dbSVijay Mahadevan           if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k];
51063d025dbSVijay Mahadevan         }
51163d025dbSVijay Mahadevan       }
51263d025dbSVijay Mahadevan     }
5132da392ccSBarry Smith   } else if (nverts == 4) { /* Linear Tetrahedra */
5142da392ccSBarry Smith     PetscReal       Dx[4]={0,0,0,0}, Dy[4]={0,0,0,0}, Dz[4]={0,0,0,0};
5152da392ccSBarry Smith     const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2];
51663d025dbSVijay Mahadevan 
5179566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(jacobian, 9));
5189566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ijacobian, 9));
51963d025dbSVijay Mahadevan 
520181f196bSVijay Mahadevan     /* compute the jacobian */
521181f196bSVijay Mahadevan     jacobian[0] = coords[1 * 3 + 0] - x0;  jacobian[1] = coords[2 * 3 + 0] - x0; jacobian[2] = coords[3 * 3 + 0] - x0;
522181f196bSVijay Mahadevan     jacobian[3] = coords[1 * 3 + 1] - y0;  jacobian[4] = coords[2 * 3 + 1] - y0; jacobian[5] = coords[3 * 3 + 1] - y0;
523181f196bSVijay Mahadevan     jacobian[6] = coords[1 * 3 + 2] - z0;  jacobian[7] = coords[2 * 3 + 2] - z0; jacobian[8] = coords[3 * 3 + 2] - z0;
52463d025dbSVijay Mahadevan 
52563d025dbSVijay Mahadevan     /* invert the jacobian */
5269566063dSJacob Faibussowitsch     PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
527*08401ef6SPierre 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);
52863d025dbSVijay Mahadevan 
529181f196bSVijay Mahadevan     if (dphidx) {
530181f196bSVijay Mahadevan       Dx[0] =   ( coords[1 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
531181f196bSVijay Mahadevan                  - coords[1 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
5322da392ccSBarry Smith                  - coords[1 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume;
533181f196bSVijay Mahadevan       Dx[1] = - ( coords[1 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
534181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
5352da392ccSBarry Smith                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume;
536181f196bSVijay Mahadevan       Dx[2] =   ( coords[1 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
537181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
5382da392ccSBarry Smith                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume;
539181f196bSVijay Mahadevan       Dx[3] =  - ( Dx[0] + Dx[1] + Dx[2]);
540181f196bSVijay Mahadevan       Dy[0] =   ( coords[0 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
541181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
5422da392ccSBarry Smith                  + coords[0 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume;
543181f196bSVijay Mahadevan       Dy[1] = - ( coords[0 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3])
544181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
5452da392ccSBarry Smith                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume;
546181f196bSVijay Mahadevan       Dy[2] =   ( coords[0 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3])
547181f196bSVijay Mahadevan                  - coords[0 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3])
5482da392ccSBarry Smith                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume;
549181f196bSVijay Mahadevan       Dy[3] =  - ( Dy[0] + Dy[1] + Dy[2]);
550181f196bSVijay Mahadevan       Dz[0] =   ( coords[0 + 1 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3])
551181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3])
5522da392ccSBarry Smith                  + coords[0 + 3 * 3] * (coords[1 + 2 * 3] - coords[1 + 1 * 3])) / *volume;
553181f196bSVijay Mahadevan       Dz[1] = - ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3])
554181f196bSVijay Mahadevan                   + coords[0 + 2 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3])
5552da392ccSBarry Smith                   - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 2 * 3])) / *volume;
556181f196bSVijay Mahadevan       Dz[2] =   ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3])
557181f196bSVijay Mahadevan                  + coords[0 + 1 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3])
5582da392ccSBarry Smith                  - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 1 * 3])) / *volume;
559181f196bSVijay Mahadevan       Dz[3] =  - ( Dz[0] + Dz[1] + Dz[2]);
560181f196bSVijay Mahadevan     }
56163d025dbSVijay Mahadevan 
5622da392ccSBarry Smith     for (j = 0; j < npts; j++) {
563a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
564181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
565181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
566181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
56763d025dbSVijay Mahadevan 
568181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
56963d025dbSVijay Mahadevan 
57063d025dbSVijay Mahadevan       phi[offset + 0] = 1.0 - r - s - t;
57163d025dbSVijay Mahadevan       phi[offset + 1] = r;
57263d025dbSVijay Mahadevan       phi[offset + 2] = s;
57363d025dbSVijay Mahadevan       phi[offset + 3] = t;
57463d025dbSVijay Mahadevan 
575181f196bSVijay Mahadevan       if (phypts) {
576181f196bSVijay Mahadevan         for (i = 0; i < nverts; ++i) {
577181f196bSVijay Mahadevan           const PetscScalar* vertices = coords + i * 3;
578181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
579181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
580181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
581181f196bSVijay Mahadevan         }
582181f196bSVijay Mahadevan       }
583181f196bSVijay Mahadevan 
584181f196bSVijay Mahadevan       /* Now set the derivatives */
58563d025dbSVijay Mahadevan       if (dphidx) {
586181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
587181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
588181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
589181f196bSVijay Mahadevan         dphidx[3 + offset] = Dx[3];
59063d025dbSVijay Mahadevan 
591181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
592181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
593181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
594181f196bSVijay Mahadevan         dphidy[3 + offset] = Dy[3];
59563d025dbSVijay Mahadevan 
596181f196bSVijay Mahadevan         dphidz[0 + offset] = Dz[0];
597181f196bSVijay Mahadevan         dphidz[1 + offset] = Dz[1];
598181f196bSVijay Mahadevan         dphidz[2 + offset] = Dz[2];
599181f196bSVijay Mahadevan         dphidz[3 + offset] = Dz[3];
60063d025dbSVijay Mahadevan       }
60163d025dbSVijay Mahadevan 
60263d025dbSVijay Mahadevan     } /* Tetrahedra -- ends */
60398921bdaSJacob 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 : %D", nverts);
60463d025dbSVijay Mahadevan   PetscFunctionReturn(0);
60563d025dbSVijay Mahadevan }
60663d025dbSVijay Mahadevan 
607cab5ea25SPierre Jolivet /*@C
60897b73a88SSatish Balay   DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
60963d025dbSVijay Mahadevan 
61063d025dbSVijay Mahadevan   The routine takes the coordinates of the vertices of an element and computes basis functions associated with
61163d025dbSVijay Mahadevan   each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
61263d025dbSVijay Mahadevan 
613d8d19677SJose E. Roman   Input Parameters:
61467b8a455SSatish Balay +  PetscInt  nverts -           the number of element vertices
61567b8a455SSatish Balay .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
61667b8a455SSatish Balay .  PetscInt  npts -             the number of evaluation points (quadrature points)
61767b8a455SSatish Balay -  PetscReal quad[3*npts] -     the evaluation points (quadrature points) in the reference space
61863d025dbSVijay Mahadevan 
619d8d19677SJose E. Roman   Output Parameters:
62067b8a455SSatish Balay +  PetscReal phypts[3*npts] -   the evaluation points (quadrature points) transformed to the physical space
62167b8a455SSatish Balay .  PetscReal jxw[npts] -        the jacobian determinant * quadrature weight necessary for assembling discrete contributions
62267b8a455SSatish Balay .  PetscReal fe_basis[npts] -   the bases values evaluated at the specified quadrature points
62367b8a455SSatish 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
62463d025dbSVijay Mahadevan 
625edc382c3SSatish Balay   Level: advanced
626edc382c3SSatish Balay 
62763d025dbSVijay Mahadevan @*/
628181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMComputeBasis(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature,
629181f196bSVijay Mahadevan                                      PetscReal *phypts, PetscReal *jacobian_quadrature_weight_product,
630181f196bSVijay Mahadevan                                      PetscReal *fe_basis, PetscReal **fe_basis_derivatives)
63163d025dbSVijay Mahadevan {
632b21a1e88SVijay Mahadevan   PetscInt        npoints,idim;
63363d025dbSVijay Mahadevan   bool            compute_der;
63463d025dbSVijay Mahadevan   const PetscReal *quadpts, *quadwts;
635181f196bSVijay Mahadevan   PetscReal       jacobian[9], ijacobian[9], volume;
63663d025dbSVijay Mahadevan 
63763d025dbSVijay Mahadevan   PetscFunctionBegin;
63863d025dbSVijay Mahadevan   PetscValidPointer(coordinates, 3);
63978dc7ee3SMatthew G. Knepley   PetscValidHeaderSpecific(quadrature, PETSCQUADRATURE_CLASSID, 4);
64063d025dbSVijay Mahadevan   PetscValidPointer(fe_basis, 7);
64163d025dbSVijay Mahadevan   compute_der = (fe_basis_derivatives != NULL);
64263d025dbSVijay Mahadevan 
64363d025dbSVijay Mahadevan   /* Get the quadrature points and weights for the given quadrature rule */
6449566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts));
645*08401ef6SPierre Jolivet   PetscCheck(idim == dim,PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%D) vs quadrature (%D)",idim,dim);
646181f196bSVijay Mahadevan   if (jacobian_quadrature_weight_product) {
6479566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jacobian_quadrature_weight_product, quadwts, npoints));
648181f196bSVijay Mahadevan   }
64963d025dbSVijay Mahadevan 
65063d025dbSVijay Mahadevan   switch (dim) {
65163d025dbSVijay Mahadevan   case 1:
6529566063dSJacob Faibussowitsch     PetscCall(Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, npoints, quadpts, phypts,
65363d025dbSVijay Mahadevan                                                jacobian_quadrature_weight_product, fe_basis,
654181f196bSVijay Mahadevan                                                (compute_der ? fe_basis_derivatives[0] : NULL),
6555f80ce2aSJacob Faibussowitsch                                                jacobian, ijacobian, &volume));
65663d025dbSVijay Mahadevan     break;
65763d025dbSVijay Mahadevan   case 2:
6589566063dSJacob Faibussowitsch     PetscCall(Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, npoints, quadpts, phypts,
65963d025dbSVijay Mahadevan                                                jacobian_quadrature_weight_product, fe_basis,
66063d025dbSVijay Mahadevan                                                (compute_der ? fe_basis_derivatives[0] : NULL),
661181f196bSVijay Mahadevan                                                (compute_der ? fe_basis_derivatives[1] : NULL),
6625f80ce2aSJacob Faibussowitsch                                                jacobian, ijacobian, &volume));
66363d025dbSVijay Mahadevan     break;
66463d025dbSVijay Mahadevan   case 3:
6659566063dSJacob Faibussowitsch     PetscCall(Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, npoints, quadpts, phypts,
66663d025dbSVijay Mahadevan                                                jacobian_quadrature_weight_product, fe_basis,
66763d025dbSVijay Mahadevan                                                (compute_der ? fe_basis_derivatives[0] : NULL),
66863d025dbSVijay Mahadevan                                                (compute_der ? fe_basis_derivatives[1] : NULL),
669181f196bSVijay Mahadevan                                                (compute_der ? fe_basis_derivatives[2] : NULL),
6705f80ce2aSJacob Faibussowitsch                                                jacobian, ijacobian, &volume));
67163d025dbSVijay Mahadevan     break;
67263d025dbSVijay Mahadevan   default:
67398921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
67463d025dbSVijay Mahadevan   }
67563d025dbSVijay Mahadevan   PetscFunctionReturn(0);
67663d025dbSVijay Mahadevan }
67763d025dbSVijay Mahadevan 
678cab5ea25SPierre Jolivet /*@C
67997b73a88SSatish Balay   DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given
68063d025dbSVijay Mahadevan   dimension and polynomial order (deciphered from number of element vertices).
68163d025dbSVijay Mahadevan 
682d8d19677SJose E. Roman   Input Parameters:
68363d025dbSVijay Mahadevan 
684a2b725a8SWilliam Gropp +  PetscInt  dim   -   the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
685a2b725a8SWilliam Gropp -  PetscInt nverts -   the number of vertices in the physical element
68663d025dbSVijay Mahadevan 
68763d025dbSVijay Mahadevan   Output Parameter:
688a2b725a8SWilliam Gropp .  PetscQuadrature quadrature -  the quadrature object with default settings to integrate polynomials defined over the element
68963d025dbSVijay Mahadevan 
690edc382c3SSatish Balay   Level: advanced
691edc382c3SSatish Balay 
69263d025dbSVijay Mahadevan @*/
693181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature)
69463d025dbSVijay Mahadevan {
69563d025dbSVijay Mahadevan   PetscReal       *w, *x;
696b21a1e88SVijay Mahadevan   PetscInt        nc=1;
69763d025dbSVijay Mahadevan 
69863d025dbSVijay Mahadevan   PetscFunctionBegin;
69963d025dbSVijay Mahadevan   /* Create an appropriate quadrature rule to sample basis */
70063d025dbSVijay Mahadevan   switch (dim)
70163d025dbSVijay Mahadevan   {
70263d025dbSVijay Mahadevan   case 1:
70363d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
7049566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(1, nc, nverts, 0, 1.0, quadrature));
70563d025dbSVijay Mahadevan     break;
70663d025dbSVijay Mahadevan   case 2:
70763d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
70863d025dbSVijay Mahadevan     if (nverts == 3) {
709a86ed394SVijay Mahadevan       const PetscInt order = 2;
710a86ed394SVijay Mahadevan       const PetscInt npoints = (order == 2 ? 3 : 6);
7119566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(npoints * 2, &x, npoints, &w));
712181f196bSVijay Mahadevan       if (npoints == 3) {
71363d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0;
71463d025dbSVijay Mahadevan         x[3] = x[4] = 2.0 / 3.0;
71563d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 1.0 / 3.0;
7162da392ccSBarry Smith       } else if (npoints == 6) {
71763d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = 0.44594849091597;
71863d025dbSVijay Mahadevan         x[3] = x[4] = 0.10810301816807;
71963d025dbSVijay Mahadevan         x[5] = 0.44594849091597;
72063d025dbSVijay Mahadevan         x[6] = x[7] = x[8] = 0.09157621350977;
72163d025dbSVijay Mahadevan         x[9] = x[10] = 0.81684757298046;
72263d025dbSVijay Mahadevan         x[11] = 0.09157621350977;
72363d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 0.22338158967801;
724181f196bSVijay Mahadevan         w[3] = w[4] = w[5] = 0.10995174365532;
72598921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %D", npoints);
7269566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
7279566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetOrder(*quadrature, order));
7289566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
7299566063dSJacob Faibussowitsch       /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
7302da392ccSBarry Smith     } else {
7319566063dSJacob Faibussowitsch       PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
73263d025dbSVijay Mahadevan     }
73363d025dbSVijay Mahadevan     break;
73463d025dbSVijay Mahadevan   case 3:
73563d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
73663d025dbSVijay Mahadevan     if (nverts == 4) {
737a86ed394SVijay Mahadevan       PetscInt order;
738a86ed394SVijay Mahadevan       const PetscInt npoints = 4; // Choose between 4 and 10
7399566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(npoints * 3, &x, npoints, &w));
740181f196bSVijay Mahadevan       if (npoints == 4) { /*  KEAST1, K1,  N=4, O=4 */
741181f196bSVijay Mahadevan         const PetscReal x_4[12] = { 0.5854101966249685, 0.1381966011250105, 0.1381966011250105,
742181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.1381966011250105,
743181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.5854101966249685,
744181f196bSVijay Mahadevan                                     0.1381966011250105, 0.5854101966249685, 0.1381966011250105
745181f196bSVijay Mahadevan                                   };
7469566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(x, x_4, 12));
747181f196bSVijay Mahadevan 
748181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0;
749181f196bSVijay Mahadevan         order = 4;
7502da392ccSBarry Smith       } else if (npoints == 10) { /*  KEAST3, K3  N=10, O=10 */
751181f196bSVijay Mahadevan         const PetscReal x_10[30] = { 0.5684305841968444, 0.1438564719343852, 0.1438564719343852,
752181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.1438564719343852,
753181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.5684305841968444,
754181f196bSVijay Mahadevan                                      0.1438564719343852, 0.5684305841968444, 0.1438564719343852,
755181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.5000000000000000,
756181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.5000000000000000,
757181f196bSVijay Mahadevan                                      0.5000000000000000, 0.5000000000000000, 0.0000000000000000,
758181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.0000000000000000,
759181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.0000000000000000,
760181f196bSVijay Mahadevan                                      0.0000000000000000, 0.0000000000000000, 0.5000000000000000
761181f196bSVijay Mahadevan                                    };
7629566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(x, x_10, 30));
763181f196bSVijay Mahadevan 
764181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 0.2177650698804054;
765181f196bSVijay Mahadevan         w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631;
766181f196bSVijay Mahadevan         order = 10;
76798921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %D", npoints);
7689566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
7699566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetOrder(*quadrature, order));
7709566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
7719566063dSJacob Faibussowitsch       /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
7722da392ccSBarry Smith     } else {
7739566063dSJacob Faibussowitsch       PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
77463d025dbSVijay Mahadevan     }
77563d025dbSVijay Mahadevan     break;
77663d025dbSVijay Mahadevan   default:
77798921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
77863d025dbSVijay Mahadevan   }
77963d025dbSVijay Mahadevan   PetscFunctionReturn(0);
78063d025dbSVijay Mahadevan }
78163d025dbSVijay Mahadevan 
782181f196bSVijay Mahadevan /* Compute Jacobians */
783181f196bSVijay Mahadevan PetscErrorCode ComputeJacobian_Internal (const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts,
784a86ed394SVijay Mahadevan   PetscReal *jacobian, PetscReal *ijacobian, PetscReal* dvolume)
785181f196bSVijay Mahadevan {
786a86ed394SVijay Mahadevan   PetscInt       i;
7872417220eSVijay Mahadevan   PetscReal      volume=1.0;
788181f196bSVijay Mahadevan 
789181f196bSVijay Mahadevan   PetscFunctionBegin;
790181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
791181f196bSVijay Mahadevan   PetscValidPointer(quad, 4);
792181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 5);
7939566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(jacobian, dim * dim));
794181f196bSVijay Mahadevan   if (ijacobian) {
7959566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ijacobian, dim * dim));
796181f196bSVijay Mahadevan   }
797181f196bSVijay Mahadevan   if (phypts) {
7989566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(phypts, /*npts=1 * */ 3));
799181f196bSVijay Mahadevan   }
800181f196bSVijay Mahadevan 
801181f196bSVijay Mahadevan   if (dim == 1) {
802181f196bSVijay Mahadevan     const PetscReal& r = quad[0];
803181f196bSVijay Mahadevan     if (nverts == 2) { /* Linear Edge */
804181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 1.0 };
805181f196bSVijay Mahadevan 
806181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
807181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
808181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
809181f196bSVijay Mahadevan       }
8102da392ccSBarry Smith     } else if (nverts == 3) { /* Quadratic Edge */
811181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r), 4.0 * r - 1.0};
812181f196bSVijay Mahadevan 
813181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
814181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
815181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
816181f196bSVijay Mahadevan       }
8172da392ccSBarry Smith     } else {
81898921bdaSJacob Faibussowitsch       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 : %D", nverts);
819181f196bSVijay Mahadevan     }
820181f196bSVijay Mahadevan 
821181f196bSVijay Mahadevan     if (ijacobian) {
822181f196bSVijay Mahadevan       /* invert the jacobian */
823181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
824181f196bSVijay Mahadevan     }
8252da392ccSBarry Smith   } else if (dim == 2) {
826181f196bSVijay Mahadevan 
827181f196bSVijay Mahadevan     if (nverts == 4) { /* Linear Quadrangle */
828181f196bSVijay Mahadevan       const PetscReal& r = quad[0];
829181f196bSVijay Mahadevan       const PetscReal& s = quad[1];
830181f196bSVijay Mahadevan 
831181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
832181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
833181f196bSVijay Mahadevan 
834181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
835181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
836181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]  * vertices[0];
837181f196bSVijay Mahadevan         jacobian[2] += dNi_dxi[i]  * vertices[1];
838181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
839181f196bSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
840181f196bSVijay Mahadevan       }
8412da392ccSBarry Smith     } else if (nverts == 3) { /* Linear triangle */
842181f196bSVijay Mahadevan       const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1];
843181f196bSVijay Mahadevan 
844181f196bSVijay Mahadevan       /* Jacobian is constant */
845181f196bSVijay Mahadevan       jacobian[0] = (coordinates[0 * 3 + 0] - x2); jacobian[1] = (coordinates[1 * 3 + 0] - x2);
846181f196bSVijay Mahadevan       jacobian[2] = (coordinates[0 * 3 + 1] - y2); jacobian[3] = (coordinates[1 * 3 + 1] - y2);
84798921bdaSJacob 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 : %D", nverts);
848181f196bSVijay Mahadevan 
849181f196bSVijay Mahadevan     /* invert the jacobian */
850181f196bSVijay Mahadevan     if (ijacobian) {
8519566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume));
852181f196bSVijay Mahadevan     }
8532da392ccSBarry Smith   } else {
854181f196bSVijay Mahadevan 
855181f196bSVijay Mahadevan     if (nverts == 8) { /* Linear Hexahedra */
856181f196bSVijay Mahadevan       const PetscReal &r = quad[0];
857181f196bSVijay Mahadevan       const PetscReal &s = quad[1];
858181f196bSVijay Mahadevan       const PetscReal &t = quad[2];
859181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s) * ( 1.0 - t),
860181f196bSVijay Mahadevan                                        ( 1.0 - s) * ( 1.0 - t),
861a86ed394SVijay Mahadevan                                        (       s) * ( 1.0 - t),
862a86ed394SVijay Mahadevan                                      - (       s) * ( 1.0 - t),
863a86ed394SVijay Mahadevan                                      - ( 1.0 - s) * (       t),
864a86ed394SVijay Mahadevan                                        ( 1.0 - s) * (       t),
865a86ed394SVijay Mahadevan                                        (       s) * (       t),
866a86ed394SVijay Mahadevan                                      - (       s) * (       t)
867181f196bSVijay Mahadevan                                     };
868181f196bSVijay Mahadevan 
869181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r) * ( 1.0 - t),
870a86ed394SVijay Mahadevan                                        - (       r) * ( 1.0 - t),
871a86ed394SVijay Mahadevan                                          (       r) * ( 1.0 - t),
872181f196bSVijay Mahadevan                                          ( 1.0 - r) * ( 1.0 - t),
873a86ed394SVijay Mahadevan                                        - ( 1.0 - r) * (       t),
874a86ed394SVijay Mahadevan                                        - (       r) * (       t),
875a86ed394SVijay Mahadevan                                          (       r) * (       t),
876a86ed394SVijay Mahadevan                                          ( 1.0 - r) * (       t)
877181f196bSVijay Mahadevan                                       };
878181f196bSVijay Mahadevan 
879181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r) * ( 1.0 - s),
880a86ed394SVijay Mahadevan                                         - (       r) * ( 1.0 - s),
881a86ed394SVijay Mahadevan                                         - (       r) * (       s),
882a86ed394SVijay Mahadevan                                         - ( 1.0 - r) * (       s),
883181f196bSVijay Mahadevan                                           ( 1.0 - r) * ( 1.0 - s),
884a86ed394SVijay Mahadevan                                           (       r) * ( 1.0 - s),
885a86ed394SVijay Mahadevan                                           (       r) * (       s),
886a86ed394SVijay Mahadevan                                           ( 1.0 - r) * (       s)
887181f196bSVijay Mahadevan                                      };
888a86ed394SVijay Mahadevan 
889181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
890181f196bSVijay Mahadevan         const PetscReal* vertex = coordinates + i * 3;
891181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
892181f196bSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
893181f196bSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
894181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
895181f196bSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
896181f196bSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
897181f196bSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
898181f196bSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
899181f196bSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
900181f196bSVijay Mahadevan       }
9012da392ccSBarry Smith     } else if (nverts == 4) { /* Linear Tetrahedra */
902181f196bSVijay Mahadevan       const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2];
903181f196bSVijay Mahadevan 
904181f196bSVijay Mahadevan       /* compute the jacobian */
905181f196bSVijay Mahadevan       jacobian[0] = coordinates[1 * 3 + 0] - x0;  jacobian[1] = coordinates[2 * 3 + 0] - x0; jacobian[2] = coordinates[3 * 3 + 0] - x0;
906181f196bSVijay Mahadevan       jacobian[3] = coordinates[1 * 3 + 1] - y0;  jacobian[4] = coordinates[2 * 3 + 1] - y0; jacobian[5] = coordinates[3 * 3 + 1] - y0;
907181f196bSVijay Mahadevan       jacobian[6] = coordinates[1 * 3 + 2] - z0;  jacobian[7] = coordinates[2 * 3 + 2] - z0; jacobian[8] = coordinates[3 * 3 + 2] - z0;
90898921bdaSJacob 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 : %D", nverts);
909181f196bSVijay Mahadevan 
910181f196bSVijay Mahadevan     if (ijacobian) {
911181f196bSVijay Mahadevan       /* invert the jacobian */
9129566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume));
913181f196bSVijay Mahadevan     }
914181f196bSVijay Mahadevan 
915181f196bSVijay Mahadevan   }
916*08401ef6SPierre Jolivet   PetscCheck(volume >= 1e-12,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Element has zero volume: %g. Degenerate element or invalid connectivity", volume);
917a86ed394SVijay Mahadevan   if (dvolume) *dvolume = volume;
918181f196bSVijay Mahadevan   PetscFunctionReturn(0);
919181f196bSVijay Mahadevan }
920181f196bSVijay Mahadevan 
921181f196bSVijay Mahadevan PetscErrorCode FEMComputeBasis_JandF(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quadrature, PetscReal *phypts,
922181f196bSVijay Mahadevan                                      PetscReal *phibasis, PetscReal *jacobian, PetscReal *ijacobian, PetscReal* volume)
923181f196bSVijay Mahadevan {
924181f196bSVijay Mahadevan   PetscFunctionBegin;
925181f196bSVijay Mahadevan   switch (dim) {
926181f196bSVijay Mahadevan     case 1:
9279566063dSJacob Faibussowitsch       PetscCall(Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, jacobian, ijacobian, volume));
928181f196bSVijay Mahadevan       break;
929181f196bSVijay Mahadevan     case 2:
9309566063dSJacob Faibussowitsch       PetscCall(Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume));
931181f196bSVijay Mahadevan       break;
932181f196bSVijay Mahadevan     case 3:
9339566063dSJacob Faibussowitsch       PetscCall(Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume));
934181f196bSVijay Mahadevan       break;
935181f196bSVijay Mahadevan     default:
93698921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
937181f196bSVijay Mahadevan   }
938181f196bSVijay Mahadevan   PetscFunctionReturn(0);
939181f196bSVijay Mahadevan }
940181f196bSVijay Mahadevan 
941cab5ea25SPierre Jolivet /*@C
94297b73a88SSatish Balay   DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the
943a86ed394SVijay Mahadevan   canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration,
944a86ed394SVijay Mahadevan   the basis function at the parametric point is also evaluated optionally.
945a86ed394SVijay Mahadevan 
946d8d19677SJose E. Roman   Input Parameters:
947a2b725a8SWilliam Gropp +  PetscInt  dim -         the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
948a2b725a8SWilliam Gropp .  PetscInt nverts -       the number of vertices in the physical element
949a2b725a8SWilliam Gropp .  PetscReal coordinates - the coordinates of vertices in the physical element
950a2b725a8SWilliam Gropp -  PetscReal[3] xphy -     the coordinates of physical point for which natural coordinates (in reference frame) are sought
951a86ed394SVijay Mahadevan 
952d8d19677SJose E. Roman   Output Parameters:
953a2b725a8SWilliam Gropp +  PetscReal[3] natparam - the natural coordinates (in reference frame) corresponding to xphy
954a2b725a8SWilliam Gropp -  PetscReal[nverts] phi - the basis functions evaluated at the natural coordinates (natparam)
955a86ed394SVijay Mahadevan 
956edc382c3SSatish Balay   Level: advanced
957edc382c3SSatish Balay 
958a86ed394SVijay Mahadevan @*/
959181f196bSVijay Mahadevan PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal* xphy, PetscReal* natparam, PetscReal* phi)
960181f196bSVijay Mahadevan {
961a86ed394SVijay Mahadevan   /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */
962181f196bSVijay Mahadevan   const PetscReal tol = 1.0e-10;
963181f196bSVijay Mahadevan   const PetscInt  max_iterations = 10;
964181f196bSVijay Mahadevan   const PetscReal error_tol_sqr = tol*tol;
965181f196bSVijay Mahadevan   PetscReal       phibasis[8], jacobian[9], ijacobian[9], volume;
966181f196bSVijay Mahadevan   PetscReal       phypts[3] = {0.0, 0.0, 0.0};
967181f196bSVijay Mahadevan   PetscReal       delta[3] = {0.0, 0.0, 0.0};
968181f196bSVijay Mahadevan   PetscInt        iters=0;
969181f196bSVijay Mahadevan   PetscReal       error=1.0;
970181f196bSVijay Mahadevan 
971181f196bSVijay Mahadevan   PetscFunctionBegin;
972181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
973181f196bSVijay Mahadevan   PetscValidPointer(xphy, 4);
974181f196bSVijay Mahadevan   PetscValidPointer(natparam, 5);
975181f196bSVijay Mahadevan 
9769566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(jacobian, dim * dim));
9779566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(ijacobian, dim * dim));
9789566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phibasis, nverts));
979181f196bSVijay Mahadevan 
980a86ed394SVijay Mahadevan   /* zero initial guess */
981a86ed394SVijay Mahadevan   natparam[0] = natparam[1] = natparam[2] = 0.0;
982181f196bSVijay Mahadevan 
983a86ed394SVijay Mahadevan   /* Compute delta = evaluate( xi) - x */
9849566063dSJacob Faibussowitsch   PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
985a86ed394SVijay Mahadevan 
986a86ed394SVijay Mahadevan   error = 0.0;
987a86ed394SVijay Mahadevan   switch(dim) {
988a86ed394SVijay Mahadevan     case 3:
989181f196bSVijay Mahadevan       delta[2] = phypts[2] - xphy[2];
990a86ed394SVijay Mahadevan       error += (delta[2]*delta[2]);
991a86ed394SVijay Mahadevan     case 2:
992a86ed394SVijay Mahadevan       delta[1] = phypts[1] - xphy[1];
993a86ed394SVijay Mahadevan       error += (delta[1]*delta[1]);
994a86ed394SVijay Mahadevan     case 1:
995a86ed394SVijay Mahadevan       delta[0] = phypts[0] - xphy[0];
996a86ed394SVijay Mahadevan       error += (delta[0]*delta[0]);
997a86ed394SVijay Mahadevan       break;
998a86ed394SVijay Mahadevan   }
999a86ed394SVijay Mahadevan 
1000181f196bSVijay Mahadevan   while (error > error_tol_sqr) {
1001181f196bSVijay Mahadevan 
1002181f196bSVijay Mahadevan     if (++iters > max_iterations)
100398921bdaSJacob Faibussowitsch       SETERRQ(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]);
1004181f196bSVijay Mahadevan 
1005181f196bSVijay Mahadevan     /* Compute natparam -= J.inverse() * delta */
1006181f196bSVijay Mahadevan     switch(dim) {
1007181f196bSVijay Mahadevan       case 1:
1008181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0];
1009181f196bSVijay Mahadevan         break;
1010181f196bSVijay Mahadevan       case 2:
1011181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1];
1012181f196bSVijay Mahadevan         natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1];
1013181f196bSVijay Mahadevan         break;
1014181f196bSVijay Mahadevan       case 3:
1015181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2];
1016181f196bSVijay Mahadevan         natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2];
1017181f196bSVijay Mahadevan         natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2];
1018181f196bSVijay Mahadevan         break;
1019181f196bSVijay Mahadevan     }
1020181f196bSVijay Mahadevan 
1021181f196bSVijay Mahadevan     /* Compute delta = evaluate( xi) - x */
10229566063dSJacob Faibussowitsch     PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
1023181f196bSVijay Mahadevan 
1024a86ed394SVijay Mahadevan     error = 0.0;
1025a86ed394SVijay Mahadevan     switch(dim) {
1026a86ed394SVijay Mahadevan       case 3:
1027181f196bSVijay Mahadevan         delta[2] = phypts[2] - xphy[2];
1028a86ed394SVijay Mahadevan         error += (delta[2]*delta[2]);
1029a86ed394SVijay Mahadevan       case 2:
1030a86ed394SVijay Mahadevan         delta[1] = phypts[1] - xphy[1];
1031a86ed394SVijay Mahadevan         error += (delta[1]*delta[1]);
1032a86ed394SVijay Mahadevan       case 1:
1033a86ed394SVijay Mahadevan         delta[0] = phypts[0] - xphy[0];
1034a86ed394SVijay Mahadevan         error += (delta[0]*delta[0]);
1035a86ed394SVijay Mahadevan         break;
1036a86ed394SVijay Mahadevan     }
1037181f196bSVijay Mahadevan   }
1038181f196bSVijay Mahadevan   if (phi) {
10399566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(phi, phibasis, nverts));
1040181f196bSVijay Mahadevan   }
1041181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1042181f196bSVijay Mahadevan }
1043