xref: /petsc/src/dm/impls/moab/dmmbfem.cxx (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
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 */
7d71ae5a4SJacob Faibussowitsch static inline PetscReal DMatrix_Determinant_2x2_Internal(const PetscReal inmat[2 * 2])
8d71ae5a4SJacob Faibussowitsch {
963d025dbSVijay Mahadevan   return inmat[0] * inmat[3] - inmat[1] * inmat[2];
1063d025dbSVijay Mahadevan }
1163d025dbSVijay Mahadevan 
12d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMatrix_Invert_2x2_Internal(const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant)
13d71ae5a4SJacob Faibussowitsch {
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;
223ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
2363d025dbSVijay Mahadevan }
2463d025dbSVijay Mahadevan 
25d71ae5a4SJacob Faibussowitsch static inline PetscReal DMatrix_Determinant_3x3_Internal(const PetscReal inmat[3 * 3])
26d71ae5a4SJacob Faibussowitsch {
279371c9d4SSatish 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]);
2863d025dbSVijay Mahadevan }
2963d025dbSVijay Mahadevan 
30d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMatrix_Invert_3x3_Internal(const PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
31d71ae5a4SJacob Faibussowitsch {
32181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_3x3_Internal(inmat);
3363d025dbSVijay Mahadevan   if (outmat) {
3463d025dbSVijay Mahadevan     outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det;
3563d025dbSVijay Mahadevan     outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det;
3663d025dbSVijay Mahadevan     outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det;
3763d025dbSVijay Mahadevan     outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det;
3863d025dbSVijay Mahadevan     outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det;
3963d025dbSVijay Mahadevan     outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det;
4063d025dbSVijay Mahadevan     outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det;
4163d025dbSVijay Mahadevan     outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det;
4263d025dbSVijay Mahadevan     outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det;
4363d025dbSVijay Mahadevan   }
4463d025dbSVijay Mahadevan   if (determinant) *determinant = det;
453ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
4663d025dbSVijay Mahadevan }
4763d025dbSVijay Mahadevan 
48*a4e35b19SJacob Faibussowitsch static inline PetscReal DMatrix_Determinant_4x4_Internal(PetscReal inmat[4 * 4])
49d71ae5a4SJacob Faibussowitsch {
509371c9d4SSatish 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]));
51181f196bSVijay Mahadevan }
52181f196bSVijay Mahadevan 
53*a4e35b19SJacob Faibussowitsch static inline PetscErrorCode DMatrix_Invert_4x4_Internal(PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
54d71ae5a4SJacob Faibussowitsch {
55181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_4x4_Internal(inmat);
56181f196bSVijay Mahadevan   if (outmat) {
57181f196bSVijay 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;
58181f196bSVijay 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;
59181f196bSVijay 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;
60181f196bSVijay 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;
61181f196bSVijay 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;
62181f196bSVijay 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;
63181f196bSVijay 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;
64181f196bSVijay 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;
65181f196bSVijay 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;
66181f196bSVijay 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;
67181f196bSVijay 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;
68181f196bSVijay 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;
69181f196bSVijay 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;
70181f196bSVijay 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;
71181f196bSVijay 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;
72181f196bSVijay 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;
73181f196bSVijay Mahadevan   }
74181f196bSVijay Mahadevan   if (determinant) *determinant = det;
753ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
76181f196bSVijay Mahadevan }
77181f196bSVijay Mahadevan 
78*a4e35b19SJacob Faibussowitsch /*
7997b73a88SSatish Balay   Compute_Lagrange_Basis_1D_Internal - Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element.
8063d025dbSVijay Mahadevan 
8163d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a linear or quadratic edge element.
8263d025dbSVijay Mahadevan 
8363d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
8463d025dbSVijay Mahadevan   and their derivatives with respect to X.
8563d025dbSVijay Mahadevan 
8663d025dbSVijay Mahadevan   Notes:
8763d025dbSVijay Mahadevan 
8863d025dbSVijay Mahadevan   Example Physical Element
89a2b725a8SWilliam Gropp .vb
9063d025dbSVijay Mahadevan     1-------2        1----3----2
9163d025dbSVijay Mahadevan       EDGE2             EDGE3
92a2b725a8SWilliam Gropp .ve
9363d025dbSVijay Mahadevan 
94d8d19677SJose E. Roman   Input Parameters:
95a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
96a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
97a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
98a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
9963d025dbSVijay Mahadevan 
100d8d19677SJose E. Roman   Output Parameters:
101a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
102a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
103a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
1046b867d5aSJose E. Roman .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
10567b8a455SSatish Balay . jacobian  - jacobian
10667b8a455SSatish Balay . ijacobian - ijacobian
10767b8a455SSatish Balay - volume    - volume
10863d025dbSVijay Mahadevan 
109edc382c3SSatish Balay   Level: advanced
110*a4e35b19SJacob Faibussowitsch */
111*a4e35b19SJacob Faibussowitsch static 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)
112d71ae5a4SJacob Faibussowitsch {
11363d025dbSVijay Mahadevan   int i, j;
11463d025dbSVijay Mahadevan 
115181f196bSVijay Mahadevan   PetscFunctionBegin;
1164f572ea9SToby Isaac   PetscAssertPointer(jacobian, 9);
1174f572ea9SToby Isaac   PetscAssertPointer(ijacobian, 10);
1184f572ea9SToby Isaac   PetscAssertPointer(volume, 11);
11948a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
12063d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
1219566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
12263d025dbSVijay Mahadevan   }
12363d025dbSVijay Mahadevan   if (nverts == 2) { /* Linear Edge */
12463d025dbSVijay Mahadevan 
1252da392ccSBarry Smith     for (j = 0; j < npts; j++) {
126a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
127181f196bSVijay Mahadevan       const PetscReal r      = quad[j];
12863d025dbSVijay Mahadevan 
129181f196bSVijay Mahadevan       phi[0 + offset] = (1.0 - r);
130181f196bSVijay Mahadevan       phi[1 + offset] = (r);
13163d025dbSVijay Mahadevan 
132181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2] = {-1.0, 1.0};
13363d025dbSVijay Mahadevan 
134181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
13563d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
136181f196bSVijay Mahadevan         const PetscReal *vertices = coords + i * 3;
137181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
138ad540459SPierre Jolivet         if (phypts) phypts[3 * j + 0] += phi[i + offset] * vertices[0];
13963d025dbSVijay Mahadevan       }
14063d025dbSVijay Mahadevan 
14163d025dbSVijay Mahadevan       /* invert the jacobian */
142181f196bSVijay Mahadevan       *volume      = jacobian[0];
143181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
144181f196bSVijay Mahadevan       jxw[j] *= *volume;
14563d025dbSVijay Mahadevan 
14663d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
14763d025dbSVijay Mahadevan       for (i = 0; i < nverts; i++) {
148181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
14963d025dbSVijay Mahadevan       }
15063d025dbSVijay Mahadevan     }
1512da392ccSBarry Smith   } else if (nverts == 3) { /* Quadratic Edge */
15263d025dbSVijay Mahadevan 
1532da392ccSBarry Smith     for (j = 0; j < npts; j++) {
154a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
155181f196bSVijay Mahadevan       const PetscReal r      = quad[j];
15663d025dbSVijay Mahadevan 
157181f196bSVijay Mahadevan       phi[0 + offset] = 1.0 + r * (2.0 * r - 3.0);
158181f196bSVijay Mahadevan       phi[1 + offset] = 4.0 * r * (1.0 - r);
159181f196bSVijay Mahadevan       phi[2 + offset] = r * (2.0 * r - 1.0);
16063d025dbSVijay Mahadevan 
161181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3] = {4 * r - 3.0, 4 * (1.0 - 2.0 * r), 4.0 * r - 1.0};
16263d025dbSVijay Mahadevan 
163181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
16463d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
165181f196bSVijay Mahadevan         const PetscReal *vertices = coords + i * 3;
166181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
167ad540459SPierre Jolivet         if (phypts) phypts[3 * j + 0] += phi[i + offset] * vertices[0];
16863d025dbSVijay Mahadevan       }
16963d025dbSVijay Mahadevan 
17063d025dbSVijay Mahadevan       /* invert the jacobian */
171181f196bSVijay Mahadevan       *volume      = jacobian[0];
172181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
173181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
17463d025dbSVijay Mahadevan 
17563d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
17663d025dbSVijay Mahadevan       for (i = 0; i < nverts; i++) {
177181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
17863d025dbSVijay Mahadevan       }
17963d025dbSVijay Mahadevan     }
18063a3b9bcSJacob 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);
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18263d025dbSVijay Mahadevan }
18363d025dbSVijay Mahadevan 
184*a4e35b19SJacob Faibussowitsch /*
18597b73a88SSatish Balay   Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.
18663d025dbSVijay Mahadevan 
18763d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a quadrangle or triangle.
18863d025dbSVijay Mahadevan 
18963d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
19063d025dbSVijay Mahadevan   and their derivatives with respect to X and Y.
19163d025dbSVijay Mahadevan 
19263d025dbSVijay Mahadevan   Notes:
19363d025dbSVijay Mahadevan 
19463d025dbSVijay Mahadevan   Example Physical Element (QUAD4)
195a2b725a8SWilliam Gropp .vb
19663d025dbSVijay Mahadevan     4------3        s
19763d025dbSVijay Mahadevan     |      |        |
19863d025dbSVijay Mahadevan     |      |        |
19963d025dbSVijay Mahadevan     |      |        |
20063d025dbSVijay Mahadevan     1------2        0-------r
201a2b725a8SWilliam Gropp .ve
20263d025dbSVijay Mahadevan 
203d8d19677SJose E. Roman   Input Parameters:
204a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
205a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
206a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
207a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
20863d025dbSVijay Mahadevan 
209d8d19677SJose E. Roman   Output Parameters:
210a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
211a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
212a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
213a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
2146b867d5aSJose E. Roman .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
21567b8a455SSatish Balay . jacobian  - jacobian
21667b8a455SSatish Balay . ijacobian - ijacobian
21767b8a455SSatish Balay - volume    - volume
21863d025dbSVijay Mahadevan 
219edc382c3SSatish Balay   Level: advanced
220*a4e35b19SJacob Faibussowitsch */
221*a4e35b19SJacob Faibussowitsch static 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)
222d71ae5a4SJacob Faibussowitsch {
223a86ed394SVijay Mahadevan   PetscInt i, j, k;
22463d025dbSVijay Mahadevan 
22563d025dbSVijay Mahadevan   PetscFunctionBegin;
2264f572ea9SToby Isaac   PetscAssertPointer(jacobian, 10);
2274f572ea9SToby Isaac   PetscAssertPointer(ijacobian, 11);
2284f572ea9SToby Isaac   PetscAssertPointer(volume, 12);
2299566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phi, npts));
23048a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
23163d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
2329566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
2339566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidy, npts * nverts));
23463d025dbSVijay Mahadevan   }
23563d025dbSVijay Mahadevan   if (nverts == 4) { /* Linear Quadrangle */
23663d025dbSVijay Mahadevan 
2379371c9d4SSatish Balay     for (j = 0; j < npts; j++) {
238a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
239181f196bSVijay Mahadevan       const PetscReal r      = quad[0 + j * 2];
240181f196bSVijay Mahadevan       const PetscReal s      = quad[1 + j * 2];
24163d025dbSVijay Mahadevan 
24263d025dbSVijay Mahadevan       phi[0 + offset] = (1.0 - r) * (1.0 - s);
24363d025dbSVijay Mahadevan       phi[1 + offset] = r * (1.0 - s);
24463d025dbSVijay Mahadevan       phi[2 + offset] = r * s;
24563d025dbSVijay Mahadevan       phi[3 + offset] = (1.0 - r) * s;
24663d025dbSVijay Mahadevan 
247181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = {-1.0 + s, 1.0 - s, s, -s};
248181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = {-1.0 + r, -r, r, 1.0 - r};
24963d025dbSVijay Mahadevan 
2509566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(jacobian, 4));
2519566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ijacobian, 4));
25263d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
253181f196bSVijay Mahadevan         const PetscReal *vertices = coords + i * 3;
25463d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
25563d025dbSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
25663d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
25763d025dbSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
258181f196bSVijay Mahadevan         if (phypts) {
259181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
260181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
261181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
262181f196bSVijay Mahadevan         }
26363d025dbSVijay Mahadevan       }
26463d025dbSVijay Mahadevan 
26563d025dbSVijay Mahadevan       /* invert the jacobian */
2669566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
26708401ef6SPierre Jolivet       PetscCheck(*volume >= 1e-12, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Quadrangular element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
26863d025dbSVijay Mahadevan 
269181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
27063d025dbSVijay Mahadevan 
271181f196bSVijay Mahadevan       /*  Let us compute the bases derivatives by scaling with inverse jacobian. */
272181f196bSVijay Mahadevan       if (dphidx) {
27363d025dbSVijay Mahadevan         for (i = 0; i < nverts; i++) {
274a86ed394SVijay Mahadevan           for (k = 0; k < 2; ++k) {
27563d025dbSVijay Mahadevan             if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0];
27663d025dbSVijay Mahadevan             if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1];
277181f196bSVijay Mahadevan           } /* for k=[0..2) */
278181f196bSVijay Mahadevan         }   /* for i=[0..nverts) */
279181f196bSVijay Mahadevan       }     /* if (dphidx) */
28063d025dbSVijay Mahadevan     }
2812da392ccSBarry Smith   } else if (nverts == 3) { /* Linear triangle */
2822da392ccSBarry Smith     const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1];
28363d025dbSVijay Mahadevan 
2849566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(jacobian, 4));
2859566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ijacobian, 4));
28663d025dbSVijay Mahadevan 
28763d025dbSVijay Mahadevan     /* Jacobian is constant */
2889371c9d4SSatish Balay     jacobian[0] = (coords[0 * 3 + 0] - x2);
2899371c9d4SSatish Balay     jacobian[1] = (coords[1 * 3 + 0] - x2);
2909371c9d4SSatish Balay     jacobian[2] = (coords[0 * 3 + 1] - y2);
2919371c9d4SSatish Balay     jacobian[3] = (coords[1 * 3 + 1] - y2);
29263d025dbSVijay Mahadevan 
29363d025dbSVijay Mahadevan     /* invert the jacobian */
2949566063dSJacob Faibussowitsch     PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
29508401ef6SPierre 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);
296181f196bSVijay Mahadevan 
297181f196bSVijay Mahadevan     const PetscReal Dx[3] = {ijacobian[0], ijacobian[2], -ijacobian[0] - ijacobian[2]};
298181f196bSVijay Mahadevan     const PetscReal Dy[3] = {ijacobian[1], ijacobian[3], -ijacobian[1] - ijacobian[3]};
29963d025dbSVijay Mahadevan 
3002da392ccSBarry Smith     for (j = 0; j < npts; j++) {
301a86ed394SVijay Mahadevan       const PetscInt  offset = j * nverts;
302181f196bSVijay Mahadevan       const PetscReal r      = quad[0 + j * 2];
303181f196bSVijay Mahadevan       const PetscReal s      = quad[1 + j * 2];
30463d025dbSVijay Mahadevan 
305181f196bSVijay Mahadevan       if (jxw) jxw[j] *= 0.5;
30663d025dbSVijay Mahadevan 
307181f196bSVijay Mahadevan       /* const PetscReal Ni[3]  = { r, s, 1.0 - r - s }; */
308181f196bSVijay Mahadevan       const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s;
309181f196bSVijay Mahadevan       const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s;
310181f196bSVijay Mahadevan       if (phypts) {
311181f196bSVijay Mahadevan         phypts[offset + 0] = phipts_x;
312181f196bSVijay Mahadevan         phypts[offset + 1] = phipts_y;
313181f196bSVijay Mahadevan       }
31463d025dbSVijay Mahadevan 
315110fc3b0SBarry Smith       /* \phi_0 = (b.y - c.y) x + (b.x - c.x) y + c.x b.y - b.x c.y */
316181f196bSVijay Mahadevan       phi[0 + offset] = (ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2));
317110fc3b0SBarry Smith       /* \phi_1 = (c.y - a.y) x + (a.x - c.x) y + c.x a.y - a.x c.y */
318181f196bSVijay Mahadevan       phi[1 + offset] = (ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2));
31963d025dbSVijay Mahadevan       phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset];
32063d025dbSVijay Mahadevan 
32163d025dbSVijay Mahadevan       if (dphidx) {
322181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
323181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
324181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
32563d025dbSVijay Mahadevan       }
32663d025dbSVijay Mahadevan 
32763d025dbSVijay Mahadevan       if (dphidy) {
328181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
329181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
330181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
33163d025dbSVijay Mahadevan       }
33263d025dbSVijay Mahadevan     }
33363a3b9bcSJacob 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);
3343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33563d025dbSVijay Mahadevan }
33663d025dbSVijay Mahadevan 
337*a4e35b19SJacob Faibussowitsch /*
33897b73a88SSatish Balay   Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element.
33963d025dbSVijay Mahadevan 
34063d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a hexahedra or tetrahedra.
34163d025dbSVijay Mahadevan 
34263d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
34363d025dbSVijay Mahadevan   and their derivatives with respect to X, Y and Z.
34463d025dbSVijay Mahadevan 
34563d025dbSVijay Mahadevan   Notes:
34663d025dbSVijay Mahadevan 
34763d025dbSVijay Mahadevan   Example Physical Element (HEX8)
348a2b725a8SWilliam Gropp .vb
34963d025dbSVijay Mahadevan       8------7
35063d025dbSVijay Mahadevan      /|     /|        t  s
35163d025dbSVijay Mahadevan     5------6 |        | /
35263d025dbSVijay Mahadevan     | |    | |        |/
35363d025dbSVijay Mahadevan     | 4----|-3        0-------r
35463d025dbSVijay Mahadevan     |/     |/
35563d025dbSVijay Mahadevan     1------2
356a2b725a8SWilliam Gropp .ve
35763d025dbSVijay Mahadevan 
358d8d19677SJose E. Roman   Input Parameters:
359a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
360a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
361a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
362a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
36363d025dbSVijay Mahadevan 
364d8d19677SJose E. Roman   Output Parameters:
365a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
366a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
367a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
368a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
369a2b725a8SWilliam Gropp .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
3706b867d5aSJose E. Roman .  PetscReal dphidz[npts] -    the derivative of the bases wrt Z-direction evaluated at the specified quadrature points
37167b8a455SSatish Balay . jacobian  - jacobian
37267b8a455SSatish Balay . ijacobian - ijacobian
37367b8a455SSatish Balay - volume    - volume
37463d025dbSVijay Mahadevan 
375edc382c3SSatish Balay   Level: advanced
376*a4e35b19SJacob Faibussowitsch */
377*a4e35b19SJacob Faibussowitsch static 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)
378d71ae5a4SJacob Faibussowitsch {
379a86ed394SVijay Mahadevan   PetscInt i, j, k;
38063d025dbSVijay Mahadevan 
38163d025dbSVijay Mahadevan   PetscFunctionBegin;
3824f572ea9SToby Isaac   PetscAssertPointer(jacobian, 11);
3834f572ea9SToby Isaac   PetscAssertPointer(ijacobian, 12);
3844f572ea9SToby Isaac   PetscAssertPointer(volume, 13);
3852da392ccSBarry Smith 
3869566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phi, npts));
38748a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
38863d025dbSVijay Mahadevan   if (dphidx) {
3899566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidx, npts * nverts));
3909566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidy, npts * nverts));
3919566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(dphidz, npts * nverts));
39263d025dbSVijay Mahadevan   }
39363d025dbSVijay Mahadevan 
39463d025dbSVijay Mahadevan   if (nverts == 8) { /* Linear Hexahedra */
39563d025dbSVijay Mahadevan 
3962da392ccSBarry Smith     for (j = 0; j < npts; j++) {
397a86ed394SVijay Mahadevan       const PetscInt   offset = j * nverts;
398181f196bSVijay Mahadevan       const PetscReal &r      = quad[j * 3 + 0];
399181f196bSVijay Mahadevan       const PetscReal &s      = quad[j * 3 + 1];
400181f196bSVijay Mahadevan       const PetscReal &t      = quad[j * 3 + 2];
40163d025dbSVijay Mahadevan 
402a86ed394SVijay Mahadevan       phi[offset + 0] = (1.0 - r) * (1.0 - s) * (1.0 - t); /* 0,0,0 */
403a86ed394SVijay Mahadevan       phi[offset + 1] = (r) * (1.0 - s) * (1.0 - t);       /* 1,0,0 */
404a86ed394SVijay Mahadevan       phi[offset + 2] = (r) * (s) * (1.0 - t);             /* 1,1,0 */
405a86ed394SVijay Mahadevan       phi[offset + 3] = (1.0 - r) * (s) * (1.0 - t);       /* 0,1,0 */
406a86ed394SVijay Mahadevan       phi[offset + 4] = (1.0 - r) * (1.0 - s) * (t);       /* 0,0,1 */
407a86ed394SVijay Mahadevan       phi[offset + 5] = (r) * (1.0 - s) * (t);             /* 1,0,1 */
408a86ed394SVijay Mahadevan       phi[offset + 6] = (r) * (s) * (t);                   /* 1,1,1 */
409a86ed394SVijay Mahadevan       phi[offset + 7] = (1.0 - r) * (s) * (t);             /* 0,1,1 */
41063d025dbSVijay Mahadevan 
4119371c9d4SSatish 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)};
41263d025dbSVijay Mahadevan 
4139371c9d4SSatish 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)};
41463d025dbSVijay Mahadevan 
4159371c9d4SSatish 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)};
41663d025dbSVijay Mahadevan 
4179566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(jacobian, 9));
4189566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ijacobian, 9));
41963d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
420181f196bSVijay Mahadevan         const PetscReal *vertex = coords + i * 3;
42163d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertex[0];
42263d025dbSVijay Mahadevan         jacobian[3] += dNi_dxi[i] * vertex[1];
42363d025dbSVijay Mahadevan         jacobian[6] += dNi_dxi[i] * vertex[2];
42463d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertex[0];
42563d025dbSVijay Mahadevan         jacobian[4] += dNi_deta[i] * vertex[1];
42663d025dbSVijay Mahadevan         jacobian[7] += dNi_deta[i] * vertex[2];
42763d025dbSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
42863d025dbSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
42963d025dbSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
430181f196bSVijay Mahadevan         if (phypts) {
431181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertex[0];
432181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertex[1];
433181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertex[2];
434181f196bSVijay Mahadevan         }
43563d025dbSVijay Mahadevan       }
43663d025dbSVijay Mahadevan 
43763d025dbSVijay Mahadevan       /* invert the jacobian */
4389566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
43908401ef6SPierre Jolivet       PetscCheck(*volume >= 1e-8, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Hexahedral element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
44063d025dbSVijay Mahadevan 
441a86ed394SVijay Mahadevan       if (jxw) jxw[j] *= (*volume);
44263d025dbSVijay Mahadevan 
44363d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
44463d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
445a86ed394SVijay Mahadevan         for (k = 0; k < 3; ++k) {
44663d025dbSVijay Mahadevan           if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0 * 3 + k];
44763d025dbSVijay Mahadevan           if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[1 * 3 + k];
44863d025dbSVijay Mahadevan           if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k];
44963d025dbSVijay Mahadevan         }
45063d025dbSVijay Mahadevan       }
45163d025dbSVijay Mahadevan     }
4522da392ccSBarry Smith   } else if (nverts == 4) { /* Linear Tetrahedra */
4532da392ccSBarry Smith     PetscReal       Dx[4] = {0, 0, 0, 0}, Dy[4] = {0, 0, 0, 0}, Dz[4] = {0, 0, 0, 0};
4542da392ccSBarry Smith     const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2];
45563d025dbSVijay Mahadevan 
4569566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(jacobian, 9));
4579566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ijacobian, 9));
45863d025dbSVijay Mahadevan 
459181f196bSVijay Mahadevan     /* compute the jacobian */
4609371c9d4SSatish Balay     jacobian[0] = coords[1 * 3 + 0] - x0;
4619371c9d4SSatish Balay     jacobian[1] = coords[2 * 3 + 0] - x0;
4629371c9d4SSatish Balay     jacobian[2] = coords[3 * 3 + 0] - x0;
4639371c9d4SSatish Balay     jacobian[3] = coords[1 * 3 + 1] - y0;
4649371c9d4SSatish Balay     jacobian[4] = coords[2 * 3 + 1] - y0;
4659371c9d4SSatish Balay     jacobian[5] = coords[3 * 3 + 1] - y0;
4669371c9d4SSatish Balay     jacobian[6] = coords[1 * 3 + 2] - z0;
4679371c9d4SSatish Balay     jacobian[7] = coords[2 * 3 + 2] - z0;
4689371c9d4SSatish Balay     jacobian[8] = coords[3 * 3 + 2] - z0;
46963d025dbSVijay Mahadevan 
47063d025dbSVijay Mahadevan     /* invert the jacobian */
4719566063dSJacob Faibussowitsch     PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
47208401ef6SPierre 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);
47363d025dbSVijay Mahadevan 
474181f196bSVijay Mahadevan     if (dphidx) {
4759371c9d4SSatish 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;
4769371c9d4SSatish 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;
4779371c9d4SSatish 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;
478181f196bSVijay Mahadevan       Dx[3] = -(Dx[0] + Dx[1] + Dx[2]);
4799371c9d4SSatish 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;
4809371c9d4SSatish 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;
4819371c9d4SSatish 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;
482181f196bSVijay Mahadevan       Dy[3] = -(Dy[0] + Dy[1] + Dy[2]);
4839371c9d4SSatish 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;
4849371c9d4SSatish 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;
4859371c9d4SSatish 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;
486181f196bSVijay Mahadevan       Dz[3] = -(Dz[0] + Dz[1] + Dz[2]);
487181f196bSVijay Mahadevan     }
48863d025dbSVijay Mahadevan 
4892da392ccSBarry Smith     for (j = 0; j < npts; j++) {
490a86ed394SVijay Mahadevan       const PetscInt   offset = j * nverts;
491181f196bSVijay Mahadevan       const PetscReal &r      = quad[j * 3 + 0];
492181f196bSVijay Mahadevan       const PetscReal &s      = quad[j * 3 + 1];
493181f196bSVijay Mahadevan       const PetscReal &t      = quad[j * 3 + 2];
49463d025dbSVijay Mahadevan 
495181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
49663d025dbSVijay Mahadevan 
49763d025dbSVijay Mahadevan       phi[offset + 0] = 1.0 - r - s - t;
49863d025dbSVijay Mahadevan       phi[offset + 1] = r;
49963d025dbSVijay Mahadevan       phi[offset + 2] = s;
50063d025dbSVijay Mahadevan       phi[offset + 3] = t;
50163d025dbSVijay Mahadevan 
502181f196bSVijay Mahadevan       if (phypts) {
503181f196bSVijay Mahadevan         for (i = 0; i < nverts; ++i) {
504181f196bSVijay Mahadevan           const PetscScalar *vertices = coords + i * 3;
505181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
506181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
507181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
508181f196bSVijay Mahadevan         }
509181f196bSVijay Mahadevan       }
510181f196bSVijay Mahadevan 
511181f196bSVijay Mahadevan       /* Now set the derivatives */
51263d025dbSVijay Mahadevan       if (dphidx) {
513181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
514181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
515181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
516181f196bSVijay Mahadevan         dphidx[3 + offset] = Dx[3];
51763d025dbSVijay Mahadevan 
518181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
519181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
520181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
521181f196bSVijay Mahadevan         dphidy[3 + offset] = Dy[3];
52263d025dbSVijay Mahadevan 
523181f196bSVijay Mahadevan         dphidz[0 + offset] = Dz[0];
524181f196bSVijay Mahadevan         dphidz[1 + offset] = Dz[1];
525181f196bSVijay Mahadevan         dphidz[2 + offset] = Dz[2];
526181f196bSVijay Mahadevan         dphidz[3 + offset] = Dz[3];
52763d025dbSVijay Mahadevan       }
52863d025dbSVijay Mahadevan 
52963d025dbSVijay Mahadevan     } /* Tetrahedra -- ends */
53063a3b9bcSJacob 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);
5313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53263d025dbSVijay Mahadevan }
53363d025dbSVijay Mahadevan 
534cab5ea25SPierre Jolivet /*@C
53597b73a88SSatish Balay   DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
53663d025dbSVijay Mahadevan 
53763d025dbSVijay Mahadevan   The routine takes the coordinates of the vertices of an element and computes basis functions associated with
53863d025dbSVijay Mahadevan   each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
53963d025dbSVijay Mahadevan 
540d8d19677SJose E. Roman   Input Parameters:
541*a4e35b19SJacob Faibussowitsch + dim         - the dimension
542*a4e35b19SJacob Faibussowitsch . nverts      - the number of element vertices
543*a4e35b19SJacob Faibussowitsch . coordinates - the physical coordinates of the vertices (in canonical numbering)
544*a4e35b19SJacob Faibussowitsch - quadrature  - the evaluation points (quadrature points) in the reference space
54563d025dbSVijay Mahadevan 
546d8d19677SJose E. Roman   Output Parameters:
547*a4e35b19SJacob Faibussowitsch + phypts                             - the evaluation points (quadrature points) transformed to the physical space
548*a4e35b19SJacob Faibussowitsch . jacobian_quadrature_weight_product - the jacobian determinant * quadrature weight necessary for assembling discrete contributions
549*a4e35b19SJacob Faibussowitsch . fe_basis                           - the bases values evaluated at the specified quadrature points
550*a4e35b19SJacob Faibussowitsch - fe_basis_derivatives               - the derivative of the bases wrt (X,Y,Z)-directions (depending on the dimension) evaluated at the specified quadrature points
55163d025dbSVijay Mahadevan 
552edc382c3SSatish Balay   Level: advanced
553edc382c3SSatish Balay 
554*a4e35b19SJacob Faibussowitsch .seealso: `DMMoabCreate()`
55563d025dbSVijay Mahadevan @*/
556d71ae5a4SJacob Faibussowitsch 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)
557d71ae5a4SJacob Faibussowitsch {
558b21a1e88SVijay Mahadevan   PetscInt         npoints, idim;
55963d025dbSVijay Mahadevan   bool             compute_der;
56063d025dbSVijay Mahadevan   const PetscReal *quadpts, *quadwts;
561181f196bSVijay Mahadevan   PetscReal        jacobian[9], ijacobian[9], volume;
56263d025dbSVijay Mahadevan 
56363d025dbSVijay Mahadevan   PetscFunctionBegin;
5644f572ea9SToby Isaac   PetscAssertPointer(coordinates, 3);
56578dc7ee3SMatthew G. Knepley   PetscValidHeaderSpecific(quadrature, PETSCQUADRATURE_CLASSID, 4);
5664f572ea9SToby Isaac   PetscAssertPointer(fe_basis, 7);
56763d025dbSVijay Mahadevan   compute_der = (fe_basis_derivatives != NULL);
56863d025dbSVijay Mahadevan 
56963d025dbSVijay Mahadevan   /* Get the quadrature points and weights for the given quadrature rule */
5709566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts));
57163a3b9bcSJacob Faibussowitsch   PetscCheck(idim == dim, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%" PetscInt_FMT ") vs quadrature (%" PetscInt_FMT ")", idim, dim);
57248a46eb9SPierre Jolivet   if (jacobian_quadrature_weight_product) PetscCall(PetscArraycpy(jacobian_quadrature_weight_product, quadwts, npoints));
57363d025dbSVijay Mahadevan 
57463d025dbSVijay Mahadevan   switch (dim) {
575d71ae5a4SJacob Faibussowitsch   case 1:
576d71ae5a4SJacob Faibussowitsch     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));
577d71ae5a4SJacob Faibussowitsch     break;
57863d025dbSVijay Mahadevan   case 2:
5799371c9d4SSatish 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));
58063d025dbSVijay Mahadevan     break;
58163d025dbSVijay Mahadevan   case 3:
5829371c9d4SSatish 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));
58363d025dbSVijay Mahadevan     break;
584d71ae5a4SJacob Faibussowitsch   default:
585d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
58663d025dbSVijay Mahadevan   }
5873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58863d025dbSVijay Mahadevan }
58963d025dbSVijay Mahadevan 
590cab5ea25SPierre Jolivet /*@C
59197b73a88SSatish Balay   DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given
59263d025dbSVijay Mahadevan   dimension and polynomial order (deciphered from number of element vertices).
59363d025dbSVijay Mahadevan 
594d8d19677SJose E. Roman   Input Parameters:
595*a4e35b19SJacob Faibussowitsch + dim    - the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
596*a4e35b19SJacob Faibussowitsch - nverts - the number of vertices in the physical element
59763d025dbSVijay Mahadevan 
59863d025dbSVijay Mahadevan   Output Parameter:
599*a4e35b19SJacob Faibussowitsch . quadrature - the quadrature object with default settings to integrate polynomials defined over the element
60063d025dbSVijay Mahadevan 
601edc382c3SSatish Balay   Level: advanced
602edc382c3SSatish Balay 
603*a4e35b19SJacob Faibussowitsch .seealso: `DMMoabCreate()`
60463d025dbSVijay Mahadevan @*/
605d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature)
606d71ae5a4SJacob Faibussowitsch {
60763d025dbSVijay Mahadevan   PetscReal *w, *x;
608b21a1e88SVijay Mahadevan   PetscInt   nc = 1;
60963d025dbSVijay Mahadevan 
61063d025dbSVijay Mahadevan   PetscFunctionBegin;
61163d025dbSVijay Mahadevan   /* Create an appropriate quadrature rule to sample basis */
6129371c9d4SSatish Balay   switch (dim) {
61363d025dbSVijay Mahadevan   case 1:
61463d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
6159566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(1, nc, nverts, 0, 1.0, quadrature));
61663d025dbSVijay Mahadevan     break;
61763d025dbSVijay Mahadevan   case 2:
61863d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
61963d025dbSVijay Mahadevan     if (nverts == 3) {
620a86ed394SVijay Mahadevan       const PetscInt order   = 2;
621a86ed394SVijay Mahadevan       const PetscInt npoints = (order == 2 ? 3 : 6);
6229566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(npoints * 2, &x, npoints, &w));
623181f196bSVijay Mahadevan       if (npoints == 3) {
62463d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0;
62563d025dbSVijay Mahadevan         x[3] = x[4] = 2.0 / 3.0;
62663d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 1.0 / 3.0;
6272da392ccSBarry Smith       } else if (npoints == 6) {
62863d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = 0.44594849091597;
62963d025dbSVijay Mahadevan         x[3] = x[4] = 0.10810301816807;
63063d025dbSVijay Mahadevan         x[5]        = 0.44594849091597;
63163d025dbSVijay Mahadevan         x[6] = x[7] = x[8] = 0.09157621350977;
63263d025dbSVijay Mahadevan         x[9] = x[10] = 0.81684757298046;
63363d025dbSVijay Mahadevan         x[11]        = 0.09157621350977;
63463d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 0.22338158967801;
635181f196bSVijay Mahadevan         w[3] = w[4] = w[5] = 0.10995174365532;
63663a3b9bcSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %" PetscInt_FMT, npoints);
6379566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
6389566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetOrder(*quadrature, order));
6399566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
6409566063dSJacob Faibussowitsch       /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
6411baa6e33SBarry Smith     } else PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
64263d025dbSVijay Mahadevan     break;
64363d025dbSVijay Mahadevan   case 3:
64463d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
64563d025dbSVijay Mahadevan     if (nverts == 4) {
646a86ed394SVijay Mahadevan       PetscInt       order;
647a86ed394SVijay Mahadevan       const PetscInt npoints = 4; // Choose between 4 and 10
6489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(npoints * 3, &x, npoints, &w));
649181f196bSVijay Mahadevan       if (npoints == 4) { /*  KEAST1, K1,  N=4, O=4 */
6509371c9d4SSatish Balay         const PetscReal x_4[12] = {0.5854101966249685, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105,
6519371c9d4SSatish Balay                                    0.1381966011250105, 0.1381966011250105, 0.5854101966249685, 0.1381966011250105, 0.5854101966249685, 0.1381966011250105};
6529566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(x, x_4, 12));
653181f196bSVijay Mahadevan 
654181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0;
655181f196bSVijay Mahadevan         order                     = 4;
6562da392ccSBarry Smith       } else if (npoints == 10) { /*  KEAST3, K3  N=10, O=10 */
6579371c9d4SSatish 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,
6589371c9d4SSatish Balay                                     0.5684305841968444, 0.1438564719343852, 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 0.5000000000000000, 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 0.5000000000000000,
6599371c9d4SSatish Balay                                     0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.5000000000000000};
6609566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(x, x_10, 30));
661181f196bSVijay Mahadevan 
662181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 0.2177650698804054;
663181f196bSVijay Mahadevan         w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631;
664181f196bSVijay Mahadevan         order                                   = 10;
66563a3b9bcSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %" PetscInt_FMT, npoints);
6669566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
6679566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetOrder(*quadrature, order));
6689566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
6699566063dSJacob Faibussowitsch       /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
6701baa6e33SBarry Smith     } else PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
67163d025dbSVijay Mahadevan     break;
672d71ae5a4SJacob Faibussowitsch   default:
673d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
67463d025dbSVijay Mahadevan   }
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67663d025dbSVijay Mahadevan }
67763d025dbSVijay Mahadevan 
678181f196bSVijay Mahadevan /* Compute Jacobians */
679d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeJacobian_Internal(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts, PetscReal *jacobian, PetscReal *ijacobian, PetscReal *dvolume)
680d71ae5a4SJacob Faibussowitsch {
681a86ed394SVijay Mahadevan   PetscInt  i;
6822417220eSVijay Mahadevan   PetscReal volume = 1.0;
683181f196bSVijay Mahadevan 
684181f196bSVijay Mahadevan   PetscFunctionBegin;
6854f572ea9SToby Isaac   PetscAssertPointer(coordinates, 3);
6864f572ea9SToby Isaac   PetscAssertPointer(quad, 4);
6874f572ea9SToby Isaac   PetscAssertPointer(jacobian, 5);
6889566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(jacobian, dim * dim));
68948a46eb9SPierre Jolivet   if (ijacobian) PetscCall(PetscArrayzero(ijacobian, dim * dim));
69048a46eb9SPierre Jolivet   if (phypts) PetscCall(PetscArrayzero(phypts, /*npts=1 * */ 3));
691181f196bSVijay Mahadevan 
692181f196bSVijay Mahadevan   if (dim == 1) {
693181f196bSVijay Mahadevan     const PetscReal &r = quad[0];
694181f196bSVijay Mahadevan     if (nverts == 2) { /* Linear Edge */
695181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2] = {-1.0, 1.0};
696181f196bSVijay Mahadevan 
697181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
698181f196bSVijay Mahadevan         const PetscReal *vertices = coordinates + i * 3;
699181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
700181f196bSVijay Mahadevan       }
7012da392ccSBarry Smith     } else if (nverts == 3) { /* Quadratic Edge */
702181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3] = {4 * r - 3.0, 4 * (1.0 - 2.0 * r), 4.0 * r - 1.0};
703181f196bSVijay Mahadevan 
704181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
705181f196bSVijay Mahadevan         const PetscReal *vertices = coordinates + i * 3;
706181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
707181f196bSVijay Mahadevan       }
7081baa6e33SBarry 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);
709181f196bSVijay Mahadevan 
710181f196bSVijay Mahadevan     if (ijacobian) {
711181f196bSVijay Mahadevan       /* invert the jacobian */
712181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
713181f196bSVijay Mahadevan     }
7142da392ccSBarry Smith   } else if (dim == 2) {
715181f196bSVijay Mahadevan     if (nverts == 4) { /* Linear Quadrangle */
716181f196bSVijay Mahadevan       const PetscReal &r = quad[0];
717181f196bSVijay Mahadevan       const PetscReal &s = quad[1];
718181f196bSVijay Mahadevan 
719181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = {-1.0 + s, 1.0 - s, s, -s};
720181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = {-1.0 + r, -r, r, 1.0 - r};
721181f196bSVijay Mahadevan 
722181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
723181f196bSVijay Mahadevan         const PetscReal *vertices = coordinates + i * 3;
724181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
725181f196bSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
726181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
727181f196bSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
728181f196bSVijay Mahadevan       }
7292da392ccSBarry Smith     } else if (nverts == 3) { /* Linear triangle */
730181f196bSVijay Mahadevan       const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1];
731181f196bSVijay Mahadevan 
732181f196bSVijay Mahadevan       /* Jacobian is constant */
7339371c9d4SSatish Balay       jacobian[0] = (coordinates[0 * 3 + 0] - x2);
7349371c9d4SSatish Balay       jacobian[1] = (coordinates[1 * 3 + 0] - x2);
7359371c9d4SSatish Balay       jacobian[2] = (coordinates[0 * 3 + 1] - y2);
7369371c9d4SSatish Balay       jacobian[3] = (coordinates[1 * 3 + 1] - y2);
73763a3b9bcSJacob 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);
738181f196bSVijay Mahadevan 
739181f196bSVijay Mahadevan     /* invert the jacobian */
74048a46eb9SPierre Jolivet     if (ijacobian) PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume));
7412da392ccSBarry Smith   } else {
742181f196bSVijay Mahadevan     if (nverts == 8) { /* Linear Hexahedra */
743181f196bSVijay Mahadevan       const PetscReal &r          = quad[0];
744181f196bSVijay Mahadevan       const PetscReal &s          = quad[1];
745181f196bSVijay Mahadevan       const PetscReal &t          = quad[2];
7469371c9d4SSatish 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)};
747181f196bSVijay Mahadevan 
7489371c9d4SSatish 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)};
749181f196bSVijay Mahadevan 
7509371c9d4SSatish 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)};
751a86ed394SVijay Mahadevan 
752181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
753181f196bSVijay Mahadevan         const PetscReal *vertex = coordinates + i * 3;
754181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertex[0];
755181f196bSVijay Mahadevan         jacobian[3] += dNi_dxi[i] * vertex[1];
756181f196bSVijay Mahadevan         jacobian[6] += dNi_dxi[i] * vertex[2];
757181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertex[0];
758181f196bSVijay Mahadevan         jacobian[4] += dNi_deta[i] * vertex[1];
759181f196bSVijay Mahadevan         jacobian[7] += dNi_deta[i] * vertex[2];
760181f196bSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
761181f196bSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
762181f196bSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
763181f196bSVijay Mahadevan       }
7642da392ccSBarry Smith     } else if (nverts == 4) { /* Linear Tetrahedra */
765181f196bSVijay Mahadevan       const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2];
766181f196bSVijay Mahadevan 
767181f196bSVijay Mahadevan       /* compute the jacobian */
7689371c9d4SSatish Balay       jacobian[0] = coordinates[1 * 3 + 0] - x0;
7699371c9d4SSatish Balay       jacobian[1] = coordinates[2 * 3 + 0] - x0;
7709371c9d4SSatish Balay       jacobian[2] = coordinates[3 * 3 + 0] - x0;
7719371c9d4SSatish Balay       jacobian[3] = coordinates[1 * 3 + 1] - y0;
7729371c9d4SSatish Balay       jacobian[4] = coordinates[2 * 3 + 1] - y0;
7739371c9d4SSatish Balay       jacobian[5] = coordinates[3 * 3 + 1] - y0;
7749371c9d4SSatish Balay       jacobian[6] = coordinates[1 * 3 + 2] - z0;
7759371c9d4SSatish Balay       jacobian[7] = coordinates[2 * 3 + 2] - z0;
7769371c9d4SSatish Balay       jacobian[8] = coordinates[3 * 3 + 2] - z0;
77763a3b9bcSJacob 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);
778181f196bSVijay Mahadevan 
779181f196bSVijay Mahadevan     if (ijacobian) {
780181f196bSVijay Mahadevan       /* invert the jacobian */
7819566063dSJacob Faibussowitsch       PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume));
782181f196bSVijay Mahadevan     }
783181f196bSVijay Mahadevan   }
78408401ef6SPierre Jolivet   PetscCheck(volume >= 1e-12, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Element has zero volume: %g. Degenerate element or invalid connectivity", volume);
785a86ed394SVijay Mahadevan   if (dvolume) *dvolume = volume;
7863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
787181f196bSVijay Mahadevan }
788181f196bSVijay Mahadevan 
789d71ae5a4SJacob Faibussowitsch 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)
790d71ae5a4SJacob Faibussowitsch {
791181f196bSVijay Mahadevan   PetscFunctionBegin;
792181f196bSVijay Mahadevan   switch (dim) {
793d71ae5a4SJacob Faibussowitsch   case 1:
794d71ae5a4SJacob Faibussowitsch     PetscCall(Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, jacobian, ijacobian, volume));
795d71ae5a4SJacob Faibussowitsch     break;
796d71ae5a4SJacob Faibussowitsch   case 2:
797d71ae5a4SJacob Faibussowitsch     PetscCall(Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume));
798d71ae5a4SJacob Faibussowitsch     break;
799d71ae5a4SJacob Faibussowitsch   case 3:
800d71ae5a4SJacob Faibussowitsch     PetscCall(Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume));
801d71ae5a4SJacob Faibussowitsch     break;
802d71ae5a4SJacob Faibussowitsch   default:
803d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
804181f196bSVijay Mahadevan   }
8053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
806181f196bSVijay Mahadevan }
807181f196bSVijay Mahadevan 
808cab5ea25SPierre Jolivet /*@C
80997b73a88SSatish Balay   DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the
810*a4e35b19SJacob Faibussowitsch   canonical reference element.
811a86ed394SVijay Mahadevan 
812d8d19677SJose E. Roman   Input Parameters:
813*a4e35b19SJacob Faibussowitsch + dim         - the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
814*a4e35b19SJacob Faibussowitsch . nverts      - the number of vertices in the physical element
815*a4e35b19SJacob Faibussowitsch . coordinates - the coordinates of vertices in the physical element
816*a4e35b19SJacob Faibussowitsch - xphy        - the coordinates of physical point for which natural coordinates (in reference frame) are sought
817a86ed394SVijay Mahadevan 
818d8d19677SJose E. Roman   Output Parameters:
819*a4e35b19SJacob Faibussowitsch + natparam - the natural coordinates (in reference frame) corresponding to xphy
820*a4e35b19SJacob Faibussowitsch - phi      - the basis functions evaluated at the natural coordinates (natparam)
821a86ed394SVijay Mahadevan 
822edc382c3SSatish Balay   Level: advanced
823edc382c3SSatish Balay 
824*a4e35b19SJacob Faibussowitsch   Notes:
825*a4e35b19SJacob Faibussowitsch   In addition to finding the inverse mapping evaluation through Newton iteration, the basis
826*a4e35b19SJacob Faibussowitsch   function at the parametric point is also evaluated optionally.
827*a4e35b19SJacob Faibussowitsch 
828*a4e35b19SJacob Faibussowitsch .seealso: `DMMoabCreate()`
829a86ed394SVijay Mahadevan @*/
830d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *xphy, PetscReal *natparam, PetscReal *phi)
831d71ae5a4SJacob Faibussowitsch {
832a86ed394SVijay Mahadevan   /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */
833181f196bSVijay Mahadevan   const PetscReal tol            = 1.0e-10;
834181f196bSVijay Mahadevan   const PetscInt  max_iterations = 10;
835181f196bSVijay Mahadevan   const PetscReal error_tol_sqr  = tol * tol;
836181f196bSVijay Mahadevan   PetscReal       phibasis[8], jacobian[9], ijacobian[9], volume;
837181f196bSVijay Mahadevan   PetscReal       phypts[3] = {0.0, 0.0, 0.0};
838181f196bSVijay Mahadevan   PetscReal       delta[3]  = {0.0, 0.0, 0.0};
839181f196bSVijay Mahadevan   PetscInt        iters     = 0;
840181f196bSVijay Mahadevan   PetscReal       error     = 1.0;
841181f196bSVijay Mahadevan 
842181f196bSVijay Mahadevan   PetscFunctionBegin;
8434f572ea9SToby Isaac   PetscAssertPointer(coordinates, 3);
8444f572ea9SToby Isaac   PetscAssertPointer(xphy, 4);
8454f572ea9SToby Isaac   PetscAssertPointer(natparam, 5);
846181f196bSVijay Mahadevan 
8479566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(jacobian, dim * dim));
8489566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(ijacobian, dim * dim));
8499566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(phibasis, nverts));
850181f196bSVijay Mahadevan 
851a86ed394SVijay Mahadevan   /* zero initial guess */
852a86ed394SVijay Mahadevan   natparam[0] = natparam[1] = natparam[2] = 0.0;
853181f196bSVijay Mahadevan 
854a86ed394SVijay Mahadevan   /* Compute delta = evaluate( xi) - x */
8559566063dSJacob Faibussowitsch   PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
856a86ed394SVijay Mahadevan 
857a86ed394SVijay Mahadevan   error = 0.0;
858a86ed394SVijay Mahadevan   switch (dim) {
859d71ae5a4SJacob Faibussowitsch   case 3:
860d71ae5a4SJacob Faibussowitsch     delta[2] = phypts[2] - xphy[2];
861d71ae5a4SJacob Faibussowitsch     error += (delta[2] * delta[2]);
862d71ae5a4SJacob Faibussowitsch   case 2:
863d71ae5a4SJacob Faibussowitsch     delta[1] = phypts[1] - xphy[1];
864d71ae5a4SJacob Faibussowitsch     error += (delta[1] * delta[1]);
865a86ed394SVijay Mahadevan   case 1:
866a86ed394SVijay Mahadevan     delta[0] = phypts[0] - xphy[0];
867a86ed394SVijay Mahadevan     error += (delta[0] * delta[0]);
868a86ed394SVijay Mahadevan     break;
869a86ed394SVijay Mahadevan   }
870a86ed394SVijay Mahadevan 
871181f196bSVijay Mahadevan   while (error > error_tol_sqr) {
8727a46b595SBarry 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]);
873181f196bSVijay Mahadevan 
874181f196bSVijay Mahadevan     /* Compute natparam -= J.inverse() * delta */
875181f196bSVijay Mahadevan     switch (dim) {
876d71ae5a4SJacob Faibussowitsch     case 1:
877d71ae5a4SJacob Faibussowitsch       natparam[0] -= ijacobian[0] * delta[0];
878d71ae5a4SJacob Faibussowitsch       break;
879181f196bSVijay Mahadevan     case 2:
880181f196bSVijay Mahadevan       natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1];
881181f196bSVijay Mahadevan       natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1];
882181f196bSVijay Mahadevan       break;
883181f196bSVijay Mahadevan     case 3:
884181f196bSVijay Mahadevan       natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2];
885181f196bSVijay Mahadevan       natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2];
886181f196bSVijay Mahadevan       natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2];
887181f196bSVijay Mahadevan       break;
888181f196bSVijay Mahadevan     }
889181f196bSVijay Mahadevan 
890181f196bSVijay Mahadevan     /* Compute delta = evaluate( xi) - x */
8919566063dSJacob Faibussowitsch     PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
892181f196bSVijay Mahadevan 
893a86ed394SVijay Mahadevan     error = 0.0;
894a86ed394SVijay Mahadevan     switch (dim) {
895d71ae5a4SJacob Faibussowitsch     case 3:
896d71ae5a4SJacob Faibussowitsch       delta[2] = phypts[2] - xphy[2];
897d71ae5a4SJacob Faibussowitsch       error += (delta[2] * delta[2]);
898d71ae5a4SJacob Faibussowitsch     case 2:
899d71ae5a4SJacob Faibussowitsch       delta[1] = phypts[1] - xphy[1];
900d71ae5a4SJacob Faibussowitsch       error += (delta[1] * delta[1]);
901a86ed394SVijay Mahadevan     case 1:
902a86ed394SVijay Mahadevan       delta[0] = phypts[0] - xphy[0];
903a86ed394SVijay Mahadevan       error += (delta[0] * delta[0]);
904a86ed394SVijay Mahadevan       break;
905a86ed394SVijay Mahadevan     }
906181f196bSVijay Mahadevan   }
90748a46eb9SPierre Jolivet   if (phi) PetscCall(PetscArraycpy(phi, phibasis, nverts));
9083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
909181f196bSVijay Mahadevan }
910