xref: /petsc/src/dm/impls/moab/dmmbfem.cxx (revision cab5ea258f2d0a024a16bc9456058de7d6803488)
163d025dbSVijay Mahadevan 
263d025dbSVijay Mahadevan #include <petscconf.h>
363d025dbSVijay Mahadevan #include <petscdt.h>            /*I "petscdt.h" I*/
463d025dbSVijay Mahadevan #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
563d025dbSVijay Mahadevan 
663d025dbSVijay Mahadevan /* Utility functions */
763d025dbSVijay Mahadevan static inline PetscReal DMatrix_Determinant_2x2_Internal ( const PetscReal inmat[2 * 2] )
863d025dbSVijay Mahadevan {
963d025dbSVijay Mahadevan   return  inmat[0] * inmat[3] - inmat[1] * inmat[2];
1063d025dbSVijay Mahadevan }
1163d025dbSVijay Mahadevan 
1263d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_2x2_Internal (const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant)
1363d025dbSVijay Mahadevan {
1463d025dbSVijay Mahadevan   if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 2x2 inversion.");
1563d025dbSVijay Mahadevan   PetscReal det = DMatrix_Determinant_2x2_Internal(inmat);
1663d025dbSVijay Mahadevan   if (outmat) {
1763d025dbSVijay Mahadevan     outmat[0] = inmat[3] / det;
1863d025dbSVijay Mahadevan     outmat[1] = -inmat[1] / det;
1963d025dbSVijay Mahadevan     outmat[2] = -inmat[2] / det;
2063d025dbSVijay Mahadevan     outmat[3] = inmat[0] / det;
2163d025dbSVijay Mahadevan   }
2263d025dbSVijay Mahadevan   if (determinant) *determinant = det;
2363d025dbSVijay Mahadevan   PetscFunctionReturn(0);
2463d025dbSVijay Mahadevan }
2563d025dbSVijay Mahadevan 
26181f196bSVijay Mahadevan static inline PetscReal DMatrix_Determinant_3x3_Internal ( const PetscReal inmat[3 * 3] )
2763d025dbSVijay Mahadevan {
2863d025dbSVijay Mahadevan   return   inmat[0] * (inmat[8] * inmat[4] - inmat[7] * inmat[5])
2963d025dbSVijay Mahadevan            - inmat[3] * (inmat[8] * inmat[1] - inmat[7] * inmat[2])
3063d025dbSVijay Mahadevan            + inmat[6] * (inmat[5] * inmat[1] - inmat[4] * inmat[2]);
3163d025dbSVijay Mahadevan }
3263d025dbSVijay Mahadevan 
3363d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_3x3_Internal (const PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
3463d025dbSVijay Mahadevan {
3563d025dbSVijay Mahadevan   if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 3x3 inversion.");
36181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_3x3_Internal(inmat);
3763d025dbSVijay Mahadevan   if (outmat) {
3863d025dbSVijay Mahadevan     outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det;
3963d025dbSVijay Mahadevan     outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det;
4063d025dbSVijay Mahadevan     outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det;
4163d025dbSVijay Mahadevan     outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det;
4263d025dbSVijay Mahadevan     outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det;
4363d025dbSVijay Mahadevan     outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det;
4463d025dbSVijay Mahadevan     outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det;
4563d025dbSVijay Mahadevan     outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det;
4663d025dbSVijay Mahadevan     outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det;
4763d025dbSVijay Mahadevan   }
4863d025dbSVijay Mahadevan   if (determinant) *determinant = det;
4963d025dbSVijay Mahadevan   PetscFunctionReturn(0);
5063d025dbSVijay Mahadevan }
5163d025dbSVijay Mahadevan 
52181f196bSVijay Mahadevan inline PetscReal DMatrix_Determinant_4x4_Internal ( PetscReal inmat[4 * 4] )
53181f196bSVijay Mahadevan {
54181f196bSVijay Mahadevan   return
55181f196bSVijay Mahadevan     inmat[0 + 0 * 4] * (
56181f196bSVijay Mahadevan       inmat[1 + 1 * 4] * ( inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4] )
57181f196bSVijay Mahadevan       - inmat[1 + 2 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4] )
58181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4] ) )
59181f196bSVijay Mahadevan     - inmat[0 + 1 * 4] * (
60181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * ( inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4] )
61181f196bSVijay Mahadevan       - inmat[1 + 2 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4] )
62181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4] ) )
63181f196bSVijay Mahadevan     + inmat[0 + 2 * 4] * (
64181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4] )
65181f196bSVijay Mahadevan       - inmat[1 + 1 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4] )
66181f196bSVijay Mahadevan       + inmat[1 + 3 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4] ) )
67181f196bSVijay Mahadevan     - inmat[0 + 3 * 4] * (
68181f196bSVijay Mahadevan       inmat[1 + 0 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4] )
69181f196bSVijay Mahadevan       - inmat[1 + 1 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4] )
70181f196bSVijay Mahadevan       + inmat[1 + 2 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4] ) );
71181f196bSVijay Mahadevan }
72181f196bSVijay Mahadevan 
73181f196bSVijay Mahadevan inline PetscErrorCode DMatrix_Invert_4x4_Internal (PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
74181f196bSVijay Mahadevan {
75181f196bSVijay Mahadevan   if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 4x4 inversion.");
76181f196bSVijay Mahadevan   PetscReal det = DMatrix_Determinant_4x4_Internal(inmat);
77181f196bSVijay Mahadevan   if (outmat) {
78181f196bSVijay Mahadevan     outmat[0] =  (inmat[5] * inmat[10] * inmat[15] + inmat[6] * inmat[11] * inmat[13] + inmat[7] * inmat[9] * inmat[14] - inmat[5] * inmat[11] * inmat[14] - inmat[6] * inmat[9] * inmat[15] - inmat[7] * inmat[10] * inmat[13]) / det;
79181f196bSVijay Mahadevan     outmat[1] =  (inmat[1] * inmat[11] * inmat[14] + inmat[2] * inmat[9] * inmat[15] + inmat[3] * inmat[10] * inmat[13] - inmat[1] * inmat[10] * inmat[15] - inmat[2] * inmat[11] * inmat[13] - inmat[3] * inmat[9] * inmat[14]) / det;
80181f196bSVijay Mahadevan     outmat[2] =  (inmat[1] * inmat[6] * inmat[15] + inmat[2] * inmat[7] * inmat[13] + inmat[3] * inmat[5] * inmat[14] - inmat[1] * inmat[7] * inmat[14] - inmat[2] * inmat[5] * inmat[15] - inmat[3] * inmat[6] * inmat[13]) / det;
81181f196bSVijay Mahadevan     outmat[3] =  (inmat[1] * inmat[7] * inmat[10] + inmat[2] * inmat[5] * inmat[11] + inmat[3] * inmat[6] * inmat[9] - inmat[1] * inmat[6] * inmat[11] - inmat[2] * inmat[7] * inmat[9] - inmat[3] * inmat[5] * inmat[10]) / det;
82181f196bSVijay Mahadevan     outmat[4] =  (inmat[4] * inmat[11] * inmat[14] + inmat[6] * inmat[8] * inmat[15] + inmat[7] * inmat[10] * inmat[12] - inmat[4] * inmat[10] * inmat[15] - inmat[6] * inmat[11] * inmat[12] - inmat[7] * inmat[8] * inmat[14]) / det;
83181f196bSVijay Mahadevan     outmat[5] =  (inmat[0] * inmat[10] * inmat[15] + inmat[2] * inmat[11] * inmat[12] + inmat[3] * inmat[8] * inmat[14] - inmat[0] * inmat[11] * inmat[14] - inmat[2] * inmat[8] * inmat[15] - inmat[3] * inmat[10] * inmat[12]) / det;
84181f196bSVijay Mahadevan     outmat[6] =  (inmat[0] * inmat[7] * inmat[14] + inmat[2] * inmat[4] * inmat[15] + inmat[3] * inmat[6] * inmat[12] - inmat[0] * inmat[6] * inmat[15] - inmat[2] * inmat[7] * inmat[12] - inmat[3] * inmat[4] * inmat[14]) / det;
85181f196bSVijay Mahadevan     outmat[7] =  (inmat[0] * inmat[6] * inmat[11] + inmat[2] * inmat[7] * inmat[8] + inmat[3] * inmat[4] * inmat[10] - inmat[0] * inmat[7] * inmat[10] - inmat[2] * inmat[4] * inmat[11] - inmat[3] * inmat[6] * inmat[8]) / det;
86181f196bSVijay Mahadevan     outmat[8] =  (inmat[4] * inmat[9] * inmat[15] + inmat[5] * inmat[11] * inmat[12] + inmat[7] * inmat[8] * inmat[13] - inmat[4] * inmat[11] * inmat[13] - inmat[5] * inmat[8] * inmat[15] - inmat[7] * inmat[9] * inmat[12]) / det;
87181f196bSVijay Mahadevan     outmat[9] =  (inmat[0] * inmat[11] * inmat[13] + inmat[1] * inmat[8] * inmat[15] + inmat[3] * inmat[9] * inmat[12] - inmat[0] * inmat[9] * inmat[15] - inmat[1] * inmat[11] * inmat[12] - inmat[3] * inmat[8] * inmat[13]) / det;
88181f196bSVijay Mahadevan     outmat[10] = (inmat[0] * inmat[5] * inmat[15] + inmat[1] * inmat[7] * inmat[12] + inmat[3] * inmat[4] * inmat[13] - inmat[0] * inmat[7] * inmat[13] - inmat[1] * inmat[4] * inmat[15] - inmat[3] * inmat[5] * inmat[12]) / det;
89181f196bSVijay Mahadevan     outmat[11] = (inmat[0] * inmat[7] * inmat[9] + inmat[1] * inmat[4] * inmat[11] + inmat[3] * inmat[5] * inmat[8] - inmat[0] * inmat[5] * inmat[11] - inmat[1] * inmat[7] * inmat[8] - inmat[3] * inmat[4] * inmat[9]) / det;
90181f196bSVijay Mahadevan     outmat[12] = (inmat[4] * inmat[10] * inmat[13] + inmat[5] * inmat[8] * inmat[14] + inmat[6] * inmat[9] * inmat[12] - inmat[4] * inmat[9] * inmat[14] - inmat[5] * inmat[10] * inmat[12] - inmat[6] * inmat[8] * inmat[13]) / det;
91181f196bSVijay Mahadevan     outmat[13] = (inmat[0] * inmat[9] * inmat[14] + inmat[1] * inmat[10] * inmat[12] + inmat[2] * inmat[8] * inmat[13] - inmat[0] * inmat[10] * inmat[13] - inmat[1] * inmat[8] * inmat[14] - inmat[2] * inmat[9] * inmat[12]) / det;
92181f196bSVijay Mahadevan     outmat[14] = (inmat[0] * inmat[6] * inmat[13] + inmat[1] * inmat[4] * inmat[14] + inmat[2] * inmat[5] * inmat[12] - inmat[0] * inmat[5] * inmat[14] - inmat[1] * inmat[6] * inmat[12] - inmat[2] * inmat[4] * inmat[13]) / det;
93181f196bSVijay Mahadevan     outmat[15] = (inmat[0] * inmat[5] * inmat[10] + inmat[1] * inmat[6] * inmat[8] + inmat[2] * inmat[4] * inmat[9] - inmat[0] * inmat[6] * inmat[9] - inmat[1] * inmat[4] * inmat[10] - inmat[2] * inmat[5] * inmat[8]) / det;
94181f196bSVijay Mahadevan   }
95181f196bSVijay Mahadevan   if (determinant) *determinant = det;
96181f196bSVijay Mahadevan   PetscFunctionReturn(0);
97181f196bSVijay Mahadevan }
98181f196bSVijay Mahadevan 
9963d025dbSVijay Mahadevan 
100*cab5ea25SPierre Jolivet /*@C
10197b73a88SSatish Balay   Compute_Lagrange_Basis_1D_Internal - Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element.
10263d025dbSVijay Mahadevan 
10363d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a linear or quadratic edge element.
10463d025dbSVijay Mahadevan 
10563d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
10663d025dbSVijay Mahadevan   and their derivatives with respect to X.
10763d025dbSVijay Mahadevan 
10863d025dbSVijay Mahadevan   Notes:
10963d025dbSVijay Mahadevan 
11063d025dbSVijay Mahadevan   Example Physical Element
111a2b725a8SWilliam Gropp .vb
11263d025dbSVijay Mahadevan     1-------2        1----3----2
11363d025dbSVijay Mahadevan       EDGE2             EDGE3
114a2b725a8SWilliam Gropp .ve
11563d025dbSVijay Mahadevan 
11663d025dbSVijay Mahadevan   Input Parameter:
117a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
118a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
119a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
120a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
12163d025dbSVijay Mahadevan 
12263d025dbSVijay Mahadevan   Output Parameter:
123a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
124a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
125a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
126a2b725a8SWilliam Gropp -  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
12763d025dbSVijay Mahadevan 
128edc382c3SSatish Balay   Level: advanced
129edc382c3SSatish Balay 
13063d025dbSVijay Mahadevan @*/
13163d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_1D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
132181f196bSVijay Mahadevan     PetscReal *jxw, PetscReal *phi, PetscReal *dphidx,
133181f196bSVijay Mahadevan     PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
13463d025dbSVijay Mahadevan {
13563d025dbSVijay Mahadevan   int i, j;
13663d025dbSVijay Mahadevan   PetscErrorCode  ierr;
13763d025dbSVijay Mahadevan 
138181f196bSVijay Mahadevan   PetscFunctionBegin;
139181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 9);
140181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 10);
141181f196bSVijay Mahadevan   PetscValidPointer(volume, 11);
142181f196bSVijay Mahadevan   if (phypts) {
143580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr);
144181f196bSVijay Mahadevan   }
14563d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
146580bdb30SBarry Smith     ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr);
14763d025dbSVijay Mahadevan   }
14863d025dbSVijay Mahadevan   if (nverts == 2) { /* Linear Edge */
14963d025dbSVijay Mahadevan 
15063d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
15163d025dbSVijay Mahadevan     {
152a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
153181f196bSVijay Mahadevan       const PetscReal r = quad[j];
15463d025dbSVijay Mahadevan 
155181f196bSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r );
156181f196bSVijay Mahadevan       phi[1 + offset] = (       r );
15763d025dbSVijay Mahadevan 
158181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 1.0 };
15963d025dbSVijay Mahadevan 
160181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
16163d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
162181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
163181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
164181f196bSVijay Mahadevan         if (phypts) {
165181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
166181f196bSVijay Mahadevan         }
16763d025dbSVijay Mahadevan       }
16863d025dbSVijay Mahadevan 
16963d025dbSVijay Mahadevan       /* invert the jacobian */
170181f196bSVijay Mahadevan       *volume = jacobian[0];
171181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
172181f196bSVijay Mahadevan       jxw[j] *= *volume;
17363d025dbSVijay Mahadevan 
17463d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
17563d025dbSVijay Mahadevan       for ( i = 0; i < nverts; i++ ) {
176181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
17763d025dbSVijay Mahadevan       }
17863d025dbSVijay Mahadevan 
17963d025dbSVijay Mahadevan     }
18063d025dbSVijay Mahadevan   }
18163d025dbSVijay Mahadevan   else if (nverts == 3) { /* Quadratic Edge */
18263d025dbSVijay Mahadevan 
18363d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
18463d025dbSVijay Mahadevan     {
185a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
186181f196bSVijay Mahadevan       const PetscReal r = quad[j];
18763d025dbSVijay Mahadevan 
188181f196bSVijay Mahadevan       phi[0 + offset] = 1.0 + r * ( 2.0 * r - 3.0 );
189181f196bSVijay Mahadevan       phi[1 + offset] = 4.0 * r * ( 1.0 - r );
190181f196bSVijay Mahadevan       phi[2 + offset] = r * ( 2.0 * r - 1.0 );
19163d025dbSVijay Mahadevan 
192181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r ), 4.0 * r - 1.0};
19363d025dbSVijay Mahadevan 
194181f196bSVijay Mahadevan       jacobian[0] = ijacobian[0] = volume[0] = 0.0;
19563d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
196181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
197181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
198181f196bSVijay Mahadevan         if (phypts) {
199181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
200181f196bSVijay Mahadevan         }
20163d025dbSVijay Mahadevan       }
20263d025dbSVijay Mahadevan 
20363d025dbSVijay Mahadevan       /* invert the jacobian */
204181f196bSVijay Mahadevan       *volume = jacobian[0];
205181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
206181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
20763d025dbSVijay Mahadevan 
20863d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
20963d025dbSVijay Mahadevan       for ( i = 0; i < nverts; i++ ) {
210181f196bSVijay Mahadevan         if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
21163d025dbSVijay Mahadevan       }
21263d025dbSVijay Mahadevan     }
21363d025dbSVijay Mahadevan   }
21463d025dbSVijay Mahadevan   else {
21563d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support EDGE2 and EDGE3 basis evaluations in 1-D : %D", nverts);
21663d025dbSVijay Mahadevan   }
21763d025dbSVijay Mahadevan #if 0
21863d025dbSVijay Mahadevan   /* verify if the computed basis functions are consistent */
21963d025dbSVijay Mahadevan   for ( j = 0; j < npts; j++ ) {
22063d025dbSVijay Mahadevan     PetscScalar phisum = 0, dphixsum = 0;
22163d025dbSVijay Mahadevan     for ( i = 0; i < nverts; i++ ) {
22263d025dbSVijay Mahadevan       phisum += phi[i + j * nverts];
22363d025dbSVijay Mahadevan       if (dphidx) dphixsum += dphidx[i + j * nverts];
22463d025dbSVijay Mahadevan     }
22563d025dbSVijay Mahadevan     PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D = %g, %g, %g\n", j, phisum, dphixsum);
22663d025dbSVijay Mahadevan   }
22763d025dbSVijay Mahadevan #endif
22863d025dbSVijay Mahadevan   PetscFunctionReturn(0);
22963d025dbSVijay Mahadevan }
23063d025dbSVijay Mahadevan 
23163d025dbSVijay Mahadevan 
232*cab5ea25SPierre Jolivet /*@C
23397b73a88SSatish Balay   Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.
23463d025dbSVijay Mahadevan 
23563d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a quadrangle or triangle.
23663d025dbSVijay Mahadevan 
23763d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
23863d025dbSVijay Mahadevan   and their derivatives with respect to X and Y.
23963d025dbSVijay Mahadevan 
24063d025dbSVijay Mahadevan   Notes:
24163d025dbSVijay Mahadevan 
24263d025dbSVijay Mahadevan   Example Physical Element (QUAD4)
243a2b725a8SWilliam Gropp .vb
24463d025dbSVijay Mahadevan     4------3        s
24563d025dbSVijay Mahadevan     |      |        |
24663d025dbSVijay Mahadevan     |      |        |
24763d025dbSVijay Mahadevan     |      |        |
24863d025dbSVijay Mahadevan     1------2        0-------r
249a2b725a8SWilliam Gropp .ve
25063d025dbSVijay Mahadevan 
25163d025dbSVijay Mahadevan   Input Parameter:
252a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
253a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
254a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
255a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
25663d025dbSVijay Mahadevan 
25763d025dbSVijay Mahadevan   Output Parameter:
258a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
259a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
260a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
261a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
262a2b725a8SWilliam Gropp -  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
26363d025dbSVijay Mahadevan 
264edc382c3SSatish Balay   Level: advanced
265edc382c3SSatish Balay 
26663d025dbSVijay Mahadevan @*/
26763d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_2D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
268181f196bSVijay Mahadevan     PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy,
269181f196bSVijay Mahadevan     PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
27063d025dbSVijay Mahadevan {
271a86ed394SVijay Mahadevan   PetscInt i, j, k;
27263d025dbSVijay Mahadevan   PetscErrorCode   ierr;
27363d025dbSVijay Mahadevan 
27463d025dbSVijay Mahadevan   PetscFunctionBegin;
275181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 10);
276181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 11);
277181f196bSVijay Mahadevan   PetscValidPointer(volume, 12);
278580bdb30SBarry Smith   ierr = PetscArrayzero(phi, npts);CHKERRQ(ierr);
279181f196bSVijay Mahadevan   if (phypts) {
280580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr);
281181f196bSVijay Mahadevan   }
28263d025dbSVijay Mahadevan   if (dphidx) { /* Reset arrays. */
283580bdb30SBarry Smith     ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr);
284580bdb30SBarry Smith     ierr = PetscArrayzero(dphidy, npts * nverts);CHKERRQ(ierr);
28563d025dbSVijay Mahadevan   }
28663d025dbSVijay Mahadevan   if (nverts == 4) { /* Linear Quadrangle */
28763d025dbSVijay Mahadevan 
28863d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
28963d025dbSVijay Mahadevan     {
290a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
291181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
292181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
29363d025dbSVijay Mahadevan 
29463d025dbSVijay Mahadevan       phi[0 + offset] = ( 1.0 - r ) * ( 1.0 - s );
29563d025dbSVijay Mahadevan       phi[1 + offset] =         r   * ( 1.0 - s );
29663d025dbSVijay Mahadevan       phi[2 + offset] =         r   *         s;
29763d025dbSVijay Mahadevan       phi[3 + offset] = ( 1.0 - r ) *         s;
29863d025dbSVijay Mahadevan 
299181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
300181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
30163d025dbSVijay Mahadevan 
302580bdb30SBarry Smith       ierr = PetscArrayzero(jacobian, 4);CHKERRQ(ierr);
303580bdb30SBarry Smith       ierr = PetscArrayzero(ijacobian, 4);CHKERRQ(ierr);
30463d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
305181f196bSVijay Mahadevan         const PetscReal* vertices = coords + i * 3;
30663d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
30763d025dbSVijay Mahadevan         jacobian[2] += dNi_dxi[i] * vertices[1];
30863d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
30963d025dbSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
310181f196bSVijay Mahadevan         if (phypts) {
311181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
312181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
313181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
314181f196bSVijay Mahadevan         }
31563d025dbSVijay Mahadevan       }
31663d025dbSVijay Mahadevan 
31763d025dbSVijay Mahadevan       /* invert the jacobian */
318181f196bSVijay Mahadevan       ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
319181f196bSVijay Mahadevan       if ( *volume < 1e-12 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Quadrangular element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
32063d025dbSVijay Mahadevan 
321181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
32263d025dbSVijay Mahadevan 
323181f196bSVijay Mahadevan       /*  Let us compute the bases derivatives by scaling with inverse jacobian. */
324181f196bSVijay Mahadevan       if (dphidx) {
32563d025dbSVijay Mahadevan         for ( i = 0; i < nverts; i++ ) {
326a86ed394SVijay Mahadevan           for ( k = 0; k < 2; ++k) {
32763d025dbSVijay Mahadevan             if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0];
32863d025dbSVijay Mahadevan             if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1];
329181f196bSVijay Mahadevan           } /* for k=[0..2) */
330181f196bSVijay Mahadevan         } /* for i=[0..nverts) */
331181f196bSVijay Mahadevan       } /* if (dphidx) */
33263d025dbSVijay Mahadevan     }
33363d025dbSVijay Mahadevan   }
33463d025dbSVijay Mahadevan   else if (nverts == 3) { /* Linear triangle */
33563d025dbSVijay Mahadevan 
336580bdb30SBarry Smith     ierr = PetscArrayzero(jacobian, 4);CHKERRQ(ierr);
337580bdb30SBarry Smith     ierr = PetscArrayzero(ijacobian, 4);CHKERRQ(ierr);
33863d025dbSVijay Mahadevan 
339181f196bSVijay Mahadevan     const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1];
340181f196bSVijay Mahadevan 
34163d025dbSVijay Mahadevan     /* Jacobian is constant */
342181f196bSVijay Mahadevan     jacobian[0] = (coords[0 * 3 + 0] - x2); jacobian[1] = (coords[1 * 3 + 0] - x2);
343181f196bSVijay Mahadevan     jacobian[2] = (coords[0 * 3 + 1] - y2); jacobian[3] = (coords[1 * 3 + 1] - y2);
34463d025dbSVijay Mahadevan 
34563d025dbSVijay Mahadevan     /* invert the jacobian */
346181f196bSVijay Mahadevan     ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
347181f196bSVijay Mahadevan     if ( *volume < 1e-12 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangular element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
348181f196bSVijay Mahadevan 
349181f196bSVijay Mahadevan     const PetscReal Dx[3] = { ijacobian[0], ijacobian[2], - ijacobian[0] - ijacobian[2] };
350181f196bSVijay Mahadevan     const PetscReal Dy[3] = { ijacobian[1], ijacobian[3], - ijacobian[1] - ijacobian[3] };
35163d025dbSVijay Mahadevan 
35263d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
35363d025dbSVijay Mahadevan     {
354a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
355181f196bSVijay Mahadevan       const PetscReal r = quad[0 + j * 2];
356181f196bSVijay Mahadevan       const PetscReal s = quad[1 + j * 2];
35763d025dbSVijay Mahadevan 
358181f196bSVijay Mahadevan       if (jxw) jxw[j] *= 0.5;
35963d025dbSVijay Mahadevan 
360181f196bSVijay Mahadevan       /* const PetscReal Ni[3]  = { r, s, 1.0 - r - s }; */
361181f196bSVijay Mahadevan       const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s;
362181f196bSVijay Mahadevan       const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s;
363181f196bSVijay Mahadevan       if (phypts) {
364181f196bSVijay Mahadevan         phypts[offset + 0] = phipts_x;
365181f196bSVijay Mahadevan         phypts[offset + 1] = phipts_y;
366181f196bSVijay Mahadevan       }
36763d025dbSVijay Mahadevan 
368181f196bSVijay Mahadevan       /* \phi_0 = (b.y − c.y) x + (b.x − c.x) y + c.x b.y − b.x c.y */
369181f196bSVijay Mahadevan       phi[0 + offset] = (  ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2) );
370181f196bSVijay Mahadevan       /* \phi_1 = (c.y − a.y) x + (a.x − c.x) y + c.x a.y − a.x c.y */
371181f196bSVijay Mahadevan       phi[1 + offset] = (  ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2) );
37263d025dbSVijay Mahadevan       phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset];
37363d025dbSVijay Mahadevan 
37463d025dbSVijay Mahadevan       if (dphidx) {
375181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
376181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
377181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
37863d025dbSVijay Mahadevan       }
37963d025dbSVijay Mahadevan 
38063d025dbSVijay Mahadevan       if (dphidy) {
381181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
382181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
383181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
38463d025dbSVijay Mahadevan       }
38563d025dbSVijay Mahadevan 
38663d025dbSVijay Mahadevan     }
38763d025dbSVijay Mahadevan   }
38863d025dbSVijay Mahadevan   else {
38963d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support QUAD4 and TRI3 basis evaluations in 2-D : %D", nverts);
39063d025dbSVijay Mahadevan   }
39163d025dbSVijay Mahadevan #if 0
39263d025dbSVijay Mahadevan   /* verify if the computed basis functions are consistent */
39363d025dbSVijay Mahadevan   for ( j = 0; j < npts; j++ ) {
39463d025dbSVijay Mahadevan     PetscScalar phisum = 0, dphixsum = 0, dphiysum = 0;
39563d025dbSVijay Mahadevan     for ( i = 0; i < nverts; i++ ) {
39663d025dbSVijay Mahadevan       phisum += phi[i + j * nverts];
39763d025dbSVijay Mahadevan       if (dphidx) dphixsum += dphidx[i + j * nverts];
39863d025dbSVijay Mahadevan       if (dphidy) dphiysum += dphidy[i + j * nverts];
39963d025dbSVijay Mahadevan     }
40063d025dbSVijay Mahadevan     PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D = %g, %g, %g\n", j, phisum, dphixsum, dphiysum);
40163d025dbSVijay Mahadevan   }
40263d025dbSVijay Mahadevan #endif
40363d025dbSVijay Mahadevan   PetscFunctionReturn(0);
40463d025dbSVijay Mahadevan }
40563d025dbSVijay Mahadevan 
40663d025dbSVijay Mahadevan 
407*cab5ea25SPierre Jolivet /*@C
40897b73a88SSatish Balay   Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element.
40963d025dbSVijay Mahadevan 
41063d025dbSVijay Mahadevan   The routine is given the coordinates of the vertices of a hexahedra or tetrahedra.
41163d025dbSVijay Mahadevan 
41263d025dbSVijay Mahadevan   The routine evaluates the basis functions associated with each quadrature point provided,
41363d025dbSVijay Mahadevan   and their derivatives with respect to X, Y and Z.
41463d025dbSVijay Mahadevan 
41563d025dbSVijay Mahadevan   Notes:
41663d025dbSVijay Mahadevan 
41763d025dbSVijay Mahadevan   Example Physical Element (HEX8)
418a2b725a8SWilliam Gropp .vb
41963d025dbSVijay Mahadevan       8------7
42063d025dbSVijay Mahadevan      /|     /|        t  s
42163d025dbSVijay Mahadevan     5------6 |        | /
42263d025dbSVijay Mahadevan     | |    | |        |/
42363d025dbSVijay Mahadevan     | 4----|-3        0-------r
42463d025dbSVijay Mahadevan     |/     |/
42563d025dbSVijay Mahadevan     1------2
426a2b725a8SWilliam Gropp .ve
42763d025dbSVijay Mahadevan 
42863d025dbSVijay Mahadevan   Input Parameter:
429a2b725a8SWilliam Gropp +  PetscInt  nverts -          the number of element vertices
430a2b725a8SWilliam Gropp .  PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
431a2b725a8SWilliam Gropp .  PetscInt  npts -            the number of evaluation points (quadrature points)
432a2b725a8SWilliam Gropp -  PetscReal quad[3*npts] -    the evaluation points (quadrature points) in the reference space
43363d025dbSVijay Mahadevan 
43463d025dbSVijay Mahadevan   Output Parameter:
435a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts] -  the evaluation points (quadrature points) transformed to the physical space
436a2b725a8SWilliam Gropp .  PetscReal jxw[npts] -       the jacobian determinant * quadrature weight necessary for assembling discrete contributions
437a2b725a8SWilliam Gropp .  PetscReal phi[npts] -       the bases evaluated at the specified quadrature points
438a2b725a8SWilliam Gropp .  PetscReal dphidx[npts] -    the derivative of the bases wrt X-direction evaluated at the specified quadrature points
439a2b725a8SWilliam Gropp .  PetscReal dphidy[npts] -    the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
440a2b725a8SWilliam Gropp -  PetscReal dphidz[npts] -    the derivative of the bases wrt Z-direction evaluated at the specified quadrature points
44163d025dbSVijay Mahadevan 
442edc382c3SSatish Balay   Level: advanced
443edc382c3SSatish Balay 
44463d025dbSVijay Mahadevan @*/
44563d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_3D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts,
446181f196bSVijay Mahadevan     PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz,
447181f196bSVijay Mahadevan     PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume)
44863d025dbSVijay Mahadevan {
449a86ed394SVijay Mahadevan   PetscInt i, j, k;
45063d025dbSVijay Mahadevan   PetscErrorCode ierr;
45163d025dbSVijay Mahadevan 
45263d025dbSVijay Mahadevan   PetscFunctionBegin;
453181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 11);
454181f196bSVijay Mahadevan   PetscValidPointer(ijacobian, 12);
455181f196bSVijay Mahadevan   PetscValidPointer(volume, 13);
45663d025dbSVijay Mahadevan   /* Reset arrays. */
457580bdb30SBarry Smith   ierr = PetscArrayzero(phi, npts);CHKERRQ(ierr);
458181f196bSVijay Mahadevan   if (phypts) {
459580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr);
460181f196bSVijay Mahadevan   }
46163d025dbSVijay Mahadevan   if (dphidx) {
462580bdb30SBarry Smith     ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr);
463580bdb30SBarry Smith     ierr = PetscArrayzero(dphidy, npts * nverts);CHKERRQ(ierr);
464580bdb30SBarry Smith     ierr = PetscArrayzero(dphidz, npts * nverts);CHKERRQ(ierr);
46563d025dbSVijay Mahadevan   }
46663d025dbSVijay Mahadevan 
46763d025dbSVijay Mahadevan   if (nverts == 8) { /* Linear Hexahedra */
46863d025dbSVijay Mahadevan 
46963d025dbSVijay Mahadevan     for (j = 0; j < npts; j++)
47063d025dbSVijay Mahadevan     {
471a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
472181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
473181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
474181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
47563d025dbSVijay Mahadevan 
476a86ed394SVijay Mahadevan       phi[offset + 0] = ( 1.0 - r ) * ( 1.0 - s ) * ( 1.0 - t ); /* 0,0,0 */
477a86ed394SVijay Mahadevan       phi[offset + 1] = (       r ) * ( 1.0 - s ) * ( 1.0 - t ); /* 1,0,0 */
478a86ed394SVijay Mahadevan       phi[offset + 2] = (       r ) * (       s ) * ( 1.0 - t ); /* 1,1,0 */
479a86ed394SVijay Mahadevan       phi[offset + 3] = ( 1.0 - r ) * (       s ) * ( 1.0 - t ); /* 0,1,0 */
480a86ed394SVijay Mahadevan       phi[offset + 4] = ( 1.0 - r ) * ( 1.0 - s ) * (       t ); /* 0,0,1 */
481a86ed394SVijay Mahadevan       phi[offset + 5] = (       r ) * ( 1.0 - s ) * (       t ); /* 1,0,1 */
482a86ed394SVijay Mahadevan       phi[offset + 6] = (       r ) * (       s ) * (       t ); /* 1,1,1 */
483a86ed394SVijay Mahadevan       phi[offset + 7] = ( 1.0 - r ) * (       s ) * (       t ); /* 0,1,1 */
48463d025dbSVijay Mahadevan 
485181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s ) * ( 1.0 - t ),
48663d025dbSVijay Mahadevan                                        ( 1.0 - s ) * ( 1.0 - t ),
487a86ed394SVijay Mahadevan                                        (       s ) * ( 1.0 - t ),
488a86ed394SVijay Mahadevan                                      - (       s ) * ( 1.0 - t ),
489a86ed394SVijay Mahadevan                                      - ( 1.0 - s ) * (       t ),
490a86ed394SVijay Mahadevan                                        ( 1.0 - s ) * (       t ),
491a86ed394SVijay Mahadevan                                        (       s ) * (       t ),
492a86ed394SVijay Mahadevan                                      - (       s ) * (       t )
49363d025dbSVijay Mahadevan                                     };
49463d025dbSVijay Mahadevan 
495181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r ) * ( 1.0 - t ),
496a86ed394SVijay Mahadevan                                        - (       r ) * ( 1.0 - t ),
497a86ed394SVijay Mahadevan                                          (       r ) * ( 1.0 - t ),
49863d025dbSVijay Mahadevan                                          ( 1.0 - r ) * ( 1.0 - t ),
499a86ed394SVijay Mahadevan                                        - ( 1.0 - r ) * (       t ),
500a86ed394SVijay Mahadevan                                        - (       r ) * (       t ),
501a86ed394SVijay Mahadevan                                          (       r ) * (       t ),
502a86ed394SVijay Mahadevan                                          ( 1.0 - r ) * (       t )
50363d025dbSVijay Mahadevan                                       };
50463d025dbSVijay Mahadevan 
505181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r ) * ( 1.0 - s ),
506a86ed394SVijay Mahadevan                                         - (       r ) * ( 1.0 - s ),
507a86ed394SVijay Mahadevan                                         - (       r ) * (       s ),
508a86ed394SVijay Mahadevan                                         - ( 1.0 - r ) * (       s ),
50963d025dbSVijay Mahadevan                                           ( 1.0 - r ) * ( 1.0 - s ),
510a86ed394SVijay Mahadevan                                           (       r ) * ( 1.0 - s ),
511a86ed394SVijay Mahadevan                                           (       r ) * (       s ),
512a86ed394SVijay Mahadevan                                           ( 1.0 - r ) * (       s )
51363d025dbSVijay Mahadevan                                      };
51463d025dbSVijay Mahadevan 
515580bdb30SBarry Smith       ierr = PetscArrayzero(jacobian, 9);CHKERRQ(ierr);
516580bdb30SBarry Smith       ierr = PetscArrayzero(ijacobian, 9);CHKERRQ(ierr);
51763d025dbSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
518181f196bSVijay Mahadevan         const PetscReal* vertex = coords + i * 3;
51963d025dbSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
52063d025dbSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
52163d025dbSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
52263d025dbSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
52363d025dbSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
52463d025dbSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
52563d025dbSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
52663d025dbSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
52763d025dbSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
528181f196bSVijay Mahadevan         if (phypts) {
529181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertex[0];
530181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertex[1];
531181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertex[2];
532181f196bSVijay Mahadevan         }
53363d025dbSVijay Mahadevan       }
53463d025dbSVijay Mahadevan 
53563d025dbSVijay Mahadevan       /* invert the jacobian */
536181f196bSVijay Mahadevan       ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
537a86ed394SVijay Mahadevan       if ( *volume < 1e-8 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Hexahedral element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
53863d025dbSVijay Mahadevan 
539a86ed394SVijay Mahadevan       if (jxw) jxw[j] *= (*volume);
54063d025dbSVijay Mahadevan 
54163d025dbSVijay Mahadevan       /*  Divide by element jacobian. */
54263d025dbSVijay Mahadevan       for ( i = 0; i < nverts; ++i ) {
543a86ed394SVijay Mahadevan         for ( k = 0; k < 3; ++k) {
54463d025dbSVijay Mahadevan           if (dphidx) dphidx[i + offset] += dNi_dxi[i]   * ijacobian[0 * 3 + k];
54563d025dbSVijay Mahadevan           if (dphidy) dphidy[i + offset] += dNi_deta[i]  * ijacobian[1 * 3 + k];
54663d025dbSVijay Mahadevan           if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k];
54763d025dbSVijay Mahadevan         }
54863d025dbSVijay Mahadevan       }
54963d025dbSVijay Mahadevan     }
55063d025dbSVijay Mahadevan   }
55163d025dbSVijay Mahadevan   else if (nverts == 4) { /* Linear Tetrahedra */
55263d025dbSVijay Mahadevan 
553580bdb30SBarry Smith     ierr = PetscArrayzero(jacobian, 9);CHKERRQ(ierr);
554580bdb30SBarry Smith     ierr = PetscArrayzero(ijacobian, 9);CHKERRQ(ierr);
55563d025dbSVijay Mahadevan 
556181f196bSVijay Mahadevan     const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2];
557181f196bSVijay Mahadevan 
558181f196bSVijay Mahadevan     /* compute the jacobian */
559181f196bSVijay Mahadevan     jacobian[0] = coords[1 * 3 + 0] - x0;  jacobian[1] = coords[2 * 3 + 0] - x0; jacobian[2] = coords[3 * 3 + 0] - x0;
560181f196bSVijay Mahadevan     jacobian[3] = coords[1 * 3 + 1] - y0;  jacobian[4] = coords[2 * 3 + 1] - y0; jacobian[5] = coords[3 * 3 + 1] - y0;
561181f196bSVijay Mahadevan     jacobian[6] = coords[1 * 3 + 2] - z0;  jacobian[7] = coords[2 * 3 + 2] - z0; jacobian[8] = coords[3 * 3 + 2] - z0;
56263d025dbSVijay Mahadevan 
56363d025dbSVijay Mahadevan     /* invert the jacobian */
564181f196bSVijay Mahadevan     ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr);
565a86ed394SVijay Mahadevan     if ( *volume < 1e-8 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume);
56663d025dbSVijay Mahadevan 
5676ea892caSVijay Mahadevan     PetscReal Dx[4]={0,0,0,0}, Dy[4]={0,0,0,0}, Dz[4]={0,0,0,0};
568181f196bSVijay Mahadevan     if (dphidx) {
569181f196bSVijay Mahadevan       Dx[0] =   ( coords[1 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
570181f196bSVijay Mahadevan                  - coords[1 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
571181f196bSVijay Mahadevan                  - coords[1 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3] )
572181f196bSVijay Mahadevan                 ) / *volume;
573181f196bSVijay Mahadevan       Dx[1] = - ( coords[1 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
574181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
575181f196bSVijay Mahadevan                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3] )
576181f196bSVijay Mahadevan                 ) / *volume;
577181f196bSVijay Mahadevan       Dx[2] =   ( coords[1 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
578181f196bSVijay Mahadevan                  - coords[1 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
579181f196bSVijay Mahadevan                  - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3] )
580181f196bSVijay Mahadevan                 ) / *volume;
581181f196bSVijay Mahadevan       Dx[3] =  - ( Dx[0] + Dx[1] + Dx[2] );
582181f196bSVijay Mahadevan       Dy[0] =   ( coords[0 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
583181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
584181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3] )
585181f196bSVijay Mahadevan                 ) / *volume;
586181f196bSVijay Mahadevan       Dy[1] = - ( coords[0 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] )
587181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
588181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3] )
589181f196bSVijay Mahadevan                 ) / *volume;
590181f196bSVijay Mahadevan       Dy[2] =   ( coords[0 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] )
591181f196bSVijay Mahadevan                  - coords[0 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] )
592181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3] )
593181f196bSVijay Mahadevan                 ) / *volume;
594181f196bSVijay Mahadevan       Dy[3] =  - ( Dy[0] + Dy[1] + Dy[2] );
595181f196bSVijay Mahadevan       Dz[0] =   ( coords[0 + 1 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3] )
596181f196bSVijay Mahadevan                  - coords[0 + 2 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3] )
597181f196bSVijay Mahadevan                  + coords[0 + 3 * 3] * (coords[1 + 2 * 3] - coords[1 + 1 * 3] )
598181f196bSVijay Mahadevan                 ) / *volume;
599181f196bSVijay Mahadevan       Dz[1] = - ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3] )
600181f196bSVijay Mahadevan                   + coords[0 + 2 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3] )
601181f196bSVijay Mahadevan                   - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 2 * 3] )
602181f196bSVijay Mahadevan                 ) / *volume;
603181f196bSVijay Mahadevan       Dz[2] =   ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3] )
604181f196bSVijay Mahadevan                  + coords[0 + 1 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3] )
605181f196bSVijay Mahadevan                  - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 1 * 3] )
606181f196bSVijay Mahadevan                 ) / *volume;
607181f196bSVijay Mahadevan       Dz[3] =  - ( Dz[0] + Dz[1] + Dz[2] );
608181f196bSVijay Mahadevan     }
60963d025dbSVijay Mahadevan 
61063d025dbSVijay Mahadevan     for ( j = 0; j < npts; j++ )
61163d025dbSVijay Mahadevan     {
612a86ed394SVijay Mahadevan       const PetscInt offset = j * nverts;
613181f196bSVijay Mahadevan       const PetscReal& r = quad[j * 3 + 0];
614181f196bSVijay Mahadevan       const PetscReal& s = quad[j * 3 + 1];
615181f196bSVijay Mahadevan       const PetscReal& t = quad[j * 3 + 2];
61663d025dbSVijay Mahadevan 
617181f196bSVijay Mahadevan       if (jxw) jxw[j] *= *volume;
61863d025dbSVijay Mahadevan 
61963d025dbSVijay Mahadevan       phi[offset + 0] = 1.0 - r - s - t;
62063d025dbSVijay Mahadevan       phi[offset + 1] = r;
62163d025dbSVijay Mahadevan       phi[offset + 2] = s;
62263d025dbSVijay Mahadevan       phi[offset + 3] = t;
62363d025dbSVijay Mahadevan 
624181f196bSVijay Mahadevan       if (phypts) {
625181f196bSVijay Mahadevan         for (i = 0; i < nverts; ++i) {
626181f196bSVijay Mahadevan           const PetscScalar* vertices = coords + i * 3;
627181f196bSVijay Mahadevan           phypts[3 * j + 0] += phi[i + offset] * vertices[0];
628181f196bSVijay Mahadevan           phypts[3 * j + 1] += phi[i + offset] * vertices[1];
629181f196bSVijay Mahadevan           phypts[3 * j + 2] += phi[i + offset] * vertices[2];
630181f196bSVijay Mahadevan         }
631181f196bSVijay Mahadevan       }
632181f196bSVijay Mahadevan 
633181f196bSVijay Mahadevan       /* Now set the derivatives */
63463d025dbSVijay Mahadevan       if (dphidx) {
635181f196bSVijay Mahadevan         dphidx[0 + offset] = Dx[0];
636181f196bSVijay Mahadevan         dphidx[1 + offset] = Dx[1];
637181f196bSVijay Mahadevan         dphidx[2 + offset] = Dx[2];
638181f196bSVijay Mahadevan         dphidx[3 + offset] = Dx[3];
63963d025dbSVijay Mahadevan 
640181f196bSVijay Mahadevan         dphidy[0 + offset] = Dy[0];
641181f196bSVijay Mahadevan         dphidy[1 + offset] = Dy[1];
642181f196bSVijay Mahadevan         dphidy[2 + offset] = Dy[2];
643181f196bSVijay Mahadevan         dphidy[3 + offset] = Dy[3];
64463d025dbSVijay Mahadevan 
645181f196bSVijay Mahadevan         dphidz[0 + offset] = Dz[0];
646181f196bSVijay Mahadevan         dphidz[1 + offset] = Dz[1];
647181f196bSVijay Mahadevan         dphidz[2 + offset] = Dz[2];
648181f196bSVijay Mahadevan         dphidz[3 + offset] = Dz[3];
64963d025dbSVijay Mahadevan       }
65063d025dbSVijay Mahadevan 
65163d025dbSVijay Mahadevan     } /* Tetrahedra -- ends */
65263d025dbSVijay Mahadevan   }
65363d025dbSVijay Mahadevan   else
65463d025dbSVijay Mahadevan   {
65563d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support HEX8 and TET4 basis evaluations in 3-D : %D", nverts);
65663d025dbSVijay Mahadevan   }
65763d025dbSVijay Mahadevan #if 0
65863d025dbSVijay Mahadevan   /* verify if the computed basis functions are consistent */
65963d025dbSVijay Mahadevan   for ( j = 0; j < npts; j++ ) {
66063d025dbSVijay Mahadevan     PetscScalar phisum = 0, dphixsum = 0, dphiysum = 0, dphizsum = 0;
66163d025dbSVijay Mahadevan     const int offset = j * nverts;
66263d025dbSVijay Mahadevan     for ( i = 0; i < nverts; i++ ) {
66363d025dbSVijay Mahadevan       phisum += phi[i + offset];
66463d025dbSVijay Mahadevan       if (dphidx) dphixsum += dphidx[i + offset];
66563d025dbSVijay Mahadevan       if (dphidy) dphiysum += dphidy[i + offset];
66663d025dbSVijay Mahadevan       if (dphidz) dphizsum += dphidz[i + offset];
66763d025dbSVijay Mahadevan       if (dphidx) PetscPrintf(PETSC_COMM_WORLD, "\t Values [%d]: [JxW] [phi, dphidx, dphidy, dphidz] = %g, %g, %g, %g, %g\n", j, jxw[j], phi[i + offset], dphidx[i + offset], dphidy[i + offset], dphidz[i + offset]);
66863d025dbSVijay Mahadevan     }
66963d025dbSVijay Mahadevan     if (dphidx) PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D (%g, %g, %g) = %g, %g, %g, %g\n", j, quad[3 * j + 0], quad[3 * j + 1], quad[3 * j + 2], phisum, dphixsum, dphiysum, dphizsum);
67063d025dbSVijay Mahadevan   }
67163d025dbSVijay Mahadevan #endif
67263d025dbSVijay Mahadevan   PetscFunctionReturn(0);
67363d025dbSVijay Mahadevan }
67463d025dbSVijay Mahadevan 
67563d025dbSVijay Mahadevan 
67663d025dbSVijay Mahadevan 
677*cab5ea25SPierre Jolivet /*@C
67897b73a88SSatish Balay   DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
67963d025dbSVijay Mahadevan 
68063d025dbSVijay Mahadevan   The routine takes the coordinates of the vertices of an element and computes basis functions associated with
68163d025dbSVijay Mahadevan   each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
68263d025dbSVijay Mahadevan 
68363d025dbSVijay Mahadevan   Input Parameter:
684a2b725a8SWilliam Gropp +  PetscInt  nverts,           the number of element vertices
68563d025dbSVijay Mahadevan .  PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering)
68663d025dbSVijay Mahadevan .  PetscInt  npts,             the number of evaluation points (quadrature points)
687a2b725a8SWilliam Gropp -  PetscReal quad[3*npts],     the evaluation points (quadrature points) in the reference space
68863d025dbSVijay Mahadevan 
68963d025dbSVijay Mahadevan   Output Parameter:
690a2b725a8SWilliam Gropp +  PetscReal phypts[3*npts],   the evaluation points (quadrature points) transformed to the physical space
69163d025dbSVijay Mahadevan .  PetscReal jxw[npts],        the jacobian determinant * quadrature weight necessary for assembling discrete contributions
69263d025dbSVijay Mahadevan .  PetscReal fe_basis[npts],        the bases values evaluated at the specified quadrature points
693a2b725a8SWilliam Gropp -  PetscReal fe_basis_derivatives[dim][npts],  the derivative of the bases wrt (X,Y,Z)-directions (depending on the dimension) evaluated at the specified quadrature points
69463d025dbSVijay Mahadevan 
695edc382c3SSatish Balay   Level: advanced
696edc382c3SSatish Balay 
69763d025dbSVijay Mahadevan @*/
698181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMComputeBasis ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature,
699181f196bSVijay Mahadevan                                        PetscReal *phypts, PetscReal *jacobian_quadrature_weight_product,
700181f196bSVijay Mahadevan                                        PetscReal *fe_basis, PetscReal **fe_basis_derivatives)
70163d025dbSVijay Mahadevan {
70263d025dbSVijay Mahadevan   PetscErrorCode  ierr;
703b21a1e88SVijay Mahadevan   PetscInt        npoints,idim;
70463d025dbSVijay Mahadevan   bool            compute_der;
70563d025dbSVijay Mahadevan   const PetscReal *quadpts, *quadwts;
706181f196bSVijay Mahadevan   PetscReal       jacobian[9], ijacobian[9], volume;
70763d025dbSVijay Mahadevan 
70863d025dbSVijay Mahadevan   PetscFunctionBegin;
70963d025dbSVijay Mahadevan   PetscValidPointer(coordinates, 3);
71078dc7ee3SMatthew G. Knepley   PetscValidHeaderSpecific(quadrature, PETSCQUADRATURE_CLASSID, 4);
71163d025dbSVijay Mahadevan   PetscValidPointer(fe_basis, 7);
71263d025dbSVijay Mahadevan   compute_der = (fe_basis_derivatives != NULL);
71363d025dbSVijay Mahadevan 
71463d025dbSVijay Mahadevan   /* Get the quadrature points and weights for the given quadrature rule */
715b21a1e88SVijay Mahadevan   ierr = PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts);CHKERRQ(ierr);
716b21a1e88SVijay Mahadevan   if (idim != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%D) vs quadrature (%D)\n",idim,dim);
717181f196bSVijay Mahadevan   if (jacobian_quadrature_weight_product) {
718580bdb30SBarry Smith     ierr = PetscArraycpy(jacobian_quadrature_weight_product, quadwts, npoints);CHKERRQ(ierr);
719181f196bSVijay Mahadevan   }
72063d025dbSVijay Mahadevan 
72163d025dbSVijay Mahadevan   switch (dim) {
72263d025dbSVijay Mahadevan   case 1:
72363d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, npoints, quadpts, phypts,
72463d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
725181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
726181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
72763d025dbSVijay Mahadevan     break;
72863d025dbSVijay Mahadevan   case 2:
72963d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, npoints, quadpts, phypts,
73063d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
73163d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
732181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[1] : NULL),
733181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
73463d025dbSVijay Mahadevan     break;
73563d025dbSVijay Mahadevan   case 3:
73663d025dbSVijay Mahadevan     ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, npoints, quadpts, phypts,
73763d025dbSVijay Mahadevan            jacobian_quadrature_weight_product, fe_basis,
73863d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[0] : NULL),
73963d025dbSVijay Mahadevan            (compute_der ? fe_basis_derivatives[1] : NULL),
740181f196bSVijay Mahadevan            (compute_der ? fe_basis_derivatives[2] : NULL),
741181f196bSVijay Mahadevan            jacobian, ijacobian, &volume);CHKERRQ(ierr);
74263d025dbSVijay Mahadevan     break;
74363d025dbSVijay Mahadevan   default:
74463d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
74563d025dbSVijay Mahadevan   }
74663d025dbSVijay Mahadevan   PetscFunctionReturn(0);
74763d025dbSVijay Mahadevan }
74863d025dbSVijay Mahadevan 
74963d025dbSVijay Mahadevan 
75063d025dbSVijay Mahadevan 
751*cab5ea25SPierre Jolivet /*@C
75297b73a88SSatish Balay   DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given
75363d025dbSVijay Mahadevan   dimension and polynomial order (deciphered from number of element vertices).
75463d025dbSVijay Mahadevan 
75563d025dbSVijay Mahadevan   Input Parameter:
75663d025dbSVijay Mahadevan 
757a2b725a8SWilliam Gropp +  PetscInt  dim   -   the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
758a2b725a8SWilliam Gropp -  PetscInt nverts -   the number of vertices in the physical element
75963d025dbSVijay Mahadevan 
76063d025dbSVijay Mahadevan   Output Parameter:
761a2b725a8SWilliam Gropp .  PetscQuadrature quadrature -  the quadrature object with default settings to integrate polynomials defined over the element
76263d025dbSVijay Mahadevan 
763edc382c3SSatish Balay   Level: advanced
764edc382c3SSatish Balay 
76563d025dbSVijay Mahadevan @*/
766181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMCreateQuadratureDefault ( const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature )
76763d025dbSVijay Mahadevan {
76863d025dbSVijay Mahadevan   PetscReal *w, *x;
769b21a1e88SVijay Mahadevan   PetscInt nc=1;
77063d025dbSVijay Mahadevan   PetscErrorCode  ierr;
77163d025dbSVijay Mahadevan 
77263d025dbSVijay Mahadevan   PetscFunctionBegin;
77363d025dbSVijay Mahadevan   /* Create an appropriate quadrature rule to sample basis */
77463d025dbSVijay Mahadevan   switch (dim)
77563d025dbSVijay Mahadevan   {
77663d025dbSVijay Mahadevan   case 1:
77763d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
778b21a1e88SVijay Mahadevan     ierr = PetscDTGaussJacobiQuadrature(1, nc, nverts, 0, 1.0, quadrature);CHKERRQ(ierr);
77963d025dbSVijay Mahadevan     break;
78063d025dbSVijay Mahadevan   case 2:
78163d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
78263d025dbSVijay Mahadevan     if (nverts == 3) {
783a86ed394SVijay Mahadevan       const PetscInt order = 2;
784a86ed394SVijay Mahadevan       const PetscInt npoints = (order == 2 ? 3 : 6);
78563d025dbSVijay Mahadevan       ierr = PetscMalloc2(npoints * 2, &x, npoints, &w);CHKERRQ(ierr);
786181f196bSVijay Mahadevan       if (npoints == 3) {
78763d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0;
78863d025dbSVijay Mahadevan         x[3] = x[4] = 2.0 / 3.0;
78963d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 1.0 / 3.0;
790181f196bSVijay Mahadevan       }
791181f196bSVijay Mahadevan       else if (npoints == 6) {
79263d025dbSVijay Mahadevan         x[0] = x[1] = x[2] = 0.44594849091597;
79363d025dbSVijay Mahadevan         x[3] = x[4] = 0.10810301816807;
79463d025dbSVijay Mahadevan         x[5] = 0.44594849091597;
79563d025dbSVijay Mahadevan         x[6] = x[7] = x[8] = 0.09157621350977;
79663d025dbSVijay Mahadevan         x[9] = x[10] = 0.81684757298046;
79763d025dbSVijay Mahadevan         x[11] = 0.09157621350977;
79863d025dbSVijay Mahadevan         w[0] = w[1] = w[2] = 0.22338158967801;
799181f196bSVijay Mahadevan         w[3] = w[4] = w[5] = 0.10995174365532;
800181f196bSVijay Mahadevan       }
801181f196bSVijay Mahadevan       else {
802181f196bSVijay Mahadevan         SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %D", npoints);
80363d025dbSVijay Mahadevan       }
80463d025dbSVijay Mahadevan       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr);
80563d025dbSVijay Mahadevan       ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr);
806b21a1e88SVijay Mahadevan       ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr);
807181f196bSVijay Mahadevan       /* ierr = PetscDTGaussJacobiQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */
80863d025dbSVijay Mahadevan     }
80963d025dbSVijay Mahadevan     else {
810b21a1e88SVijay Mahadevan       ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr);
81163d025dbSVijay Mahadevan     }
81263d025dbSVijay Mahadevan     break;
81363d025dbSVijay Mahadevan   case 3:
81463d025dbSVijay Mahadevan     /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
81563d025dbSVijay Mahadevan     if (nverts == 4) {
816a86ed394SVijay Mahadevan       PetscInt order;
817a86ed394SVijay Mahadevan       const PetscInt npoints = 4; // Choose between 4 and 10
818181f196bSVijay Mahadevan       ierr = PetscMalloc2(npoints * 3, &x, npoints, &w);CHKERRQ(ierr);
819181f196bSVijay Mahadevan       if (npoints == 4) { /*  KEAST1, K1,  N=4, O=4 */
820181f196bSVijay Mahadevan         const PetscReal x_4[12] = { 0.5854101966249685, 0.1381966011250105, 0.1381966011250105,
821181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.1381966011250105,
822181f196bSVijay Mahadevan                                     0.1381966011250105, 0.1381966011250105, 0.5854101966249685,
823181f196bSVijay Mahadevan                                     0.1381966011250105, 0.5854101966249685, 0.1381966011250105
824181f196bSVijay Mahadevan                                   };
825580bdb30SBarry Smith         ierr = PetscArraycpy(x, x_4, 12);CHKERRQ(ierr);
826181f196bSVijay Mahadevan 
827181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0;
828181f196bSVijay Mahadevan         order = 4;
829181f196bSVijay Mahadevan       }
830181f196bSVijay Mahadevan       else if (npoints == 10) { /*  KEAST3, K3  N=10, O=10 */
831181f196bSVijay Mahadevan         const PetscReal x_10[30] = { 0.5684305841968444, 0.1438564719343852, 0.1438564719343852,
832181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.1438564719343852,
833181f196bSVijay Mahadevan                                      0.1438564719343852, 0.1438564719343852, 0.5684305841968444,
834181f196bSVijay Mahadevan                                      0.1438564719343852, 0.5684305841968444, 0.1438564719343852,
835181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.5000000000000000,
836181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.5000000000000000,
837181f196bSVijay Mahadevan                                      0.5000000000000000, 0.5000000000000000, 0.0000000000000000,
838181f196bSVijay Mahadevan                                      0.5000000000000000, 0.0000000000000000, 0.0000000000000000,
839181f196bSVijay Mahadevan                                      0.0000000000000000, 0.5000000000000000, 0.0000000000000000,
840181f196bSVijay Mahadevan                                      0.0000000000000000, 0.0000000000000000, 0.5000000000000000
841181f196bSVijay Mahadevan                                    };
842580bdb30SBarry Smith         ierr = PetscArraycpy(x, x_10, 30);CHKERRQ(ierr);
843181f196bSVijay Mahadevan 
844181f196bSVijay Mahadevan         w[0] = w[1] = w[2] = w[3] = 0.2177650698804054;
845181f196bSVijay Mahadevan         w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631;
846181f196bSVijay Mahadevan         order = 10;
847181f196bSVijay Mahadevan       }
848181f196bSVijay Mahadevan       else {
849181f196bSVijay Mahadevan         SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %D", npoints);
850181f196bSVijay Mahadevan       }
851181f196bSVijay Mahadevan       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr);
852181f196bSVijay Mahadevan       ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr);
853b21a1e88SVijay Mahadevan       ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr);
854181f196bSVijay Mahadevan       /* ierr = PetscDTGaussJacobiQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */
85563d025dbSVijay Mahadevan     }
85663d025dbSVijay Mahadevan     else {
857b21a1e88SVijay Mahadevan       ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr);
85863d025dbSVijay Mahadevan     }
85963d025dbSVijay Mahadevan     break;
86063d025dbSVijay Mahadevan   default:
86163d025dbSVijay Mahadevan     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
86263d025dbSVijay Mahadevan   }
86363d025dbSVijay Mahadevan   PetscFunctionReturn(0);
86463d025dbSVijay Mahadevan }
86563d025dbSVijay Mahadevan 
866181f196bSVijay Mahadevan /* Compute Jacobians */
867181f196bSVijay Mahadevan PetscErrorCode ComputeJacobian_Internal ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts,
868a86ed394SVijay Mahadevan   PetscReal *jacobian, PetscReal *ijacobian, PetscReal* dvolume)
869181f196bSVijay Mahadevan {
870a86ed394SVijay Mahadevan   PetscInt i;
8712417220eSVijay Mahadevan   PetscReal volume=1.0;
872181f196bSVijay Mahadevan   PetscErrorCode ierr;
873181f196bSVijay Mahadevan 
874181f196bSVijay Mahadevan   PetscFunctionBegin;
875181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
876181f196bSVijay Mahadevan   PetscValidPointer(quad, 4);
877181f196bSVijay Mahadevan   PetscValidPointer(jacobian, 5);
878580bdb30SBarry Smith   ierr = PetscArrayzero(jacobian, dim * dim);CHKERRQ(ierr);
879181f196bSVijay Mahadevan   if (ijacobian) {
880580bdb30SBarry Smith     ierr = PetscArrayzero(ijacobian, dim * dim);CHKERRQ(ierr);
881181f196bSVijay Mahadevan   }
882181f196bSVijay Mahadevan   if (phypts) {
883580bdb30SBarry Smith     ierr = PetscArrayzero(phypts, /*npts=1 * */ 3);CHKERRQ(ierr);
884181f196bSVijay Mahadevan   }
885181f196bSVijay Mahadevan 
886181f196bSVijay Mahadevan   if (dim == 1) {
887181f196bSVijay Mahadevan 
888181f196bSVijay Mahadevan     const PetscReal& r = quad[0];
889181f196bSVijay Mahadevan     if (nverts == 2) { /* Linear Edge */
890181f196bSVijay Mahadevan       const PetscReal dNi_dxi[2]  = { -1.0, 1.0 };
891181f196bSVijay Mahadevan 
892181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
893181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
894181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
895181f196bSVijay Mahadevan       }
896181f196bSVijay Mahadevan     }
897181f196bSVijay Mahadevan     else if (nverts == 3) { /* Quadratic Edge */
898181f196bSVijay Mahadevan       const PetscReal dNi_dxi[3]  = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r ), 4.0 * r - 1.0};
899181f196bSVijay Mahadevan 
900181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
901181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
902181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i] * vertices[0];
903181f196bSVijay Mahadevan       }
904181f196bSVijay Mahadevan     }
905181f196bSVijay Mahadevan     else {
906181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 1-D entity vertices are invalid. Currently only support EDGE2 and EDGE3 basis evaluations in 1-D : %D", nverts);
907181f196bSVijay Mahadevan     }
908181f196bSVijay Mahadevan 
909181f196bSVijay Mahadevan     if (ijacobian) {
910181f196bSVijay Mahadevan       /* invert the jacobian */
911181f196bSVijay Mahadevan       ijacobian[0] = 1.0 / jacobian[0];
912181f196bSVijay Mahadevan     }
913181f196bSVijay Mahadevan 
914181f196bSVijay Mahadevan   }
915181f196bSVijay Mahadevan   else if (dim == 2) {
916181f196bSVijay Mahadevan 
917181f196bSVijay Mahadevan     if (nverts == 4) { /* Linear Quadrangle */
918181f196bSVijay Mahadevan       const PetscReal& r = quad[0];
919181f196bSVijay Mahadevan       const PetscReal& s = quad[1];
920181f196bSVijay Mahadevan 
921181f196bSVijay Mahadevan       const PetscReal dNi_dxi[4]  = { -1.0 + s, 1.0 - s, s, -s };
922181f196bSVijay Mahadevan       const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r };
923181f196bSVijay Mahadevan 
924181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
925181f196bSVijay Mahadevan         const PetscReal* vertices = coordinates + i * 3;
926181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]  * vertices[0];
927181f196bSVijay Mahadevan         jacobian[2] += dNi_dxi[i]  * vertices[1];
928181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i] * vertices[0];
929181f196bSVijay Mahadevan         jacobian[3] += dNi_deta[i] * vertices[1];
930181f196bSVijay Mahadevan       }
931181f196bSVijay Mahadevan     }
932181f196bSVijay Mahadevan     else if (nverts == 3) { /* Linear triangle */
933181f196bSVijay Mahadevan       const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1];
934181f196bSVijay Mahadevan 
935181f196bSVijay Mahadevan       /* Jacobian is constant */
936181f196bSVijay Mahadevan       jacobian[0] = (coordinates[0 * 3 + 0] - x2); jacobian[1] = (coordinates[1 * 3 + 0] - x2);
937181f196bSVijay Mahadevan       jacobian[2] = (coordinates[0 * 3 + 1] - y2); jacobian[3] = (coordinates[1 * 3 + 1] - y2);
938181f196bSVijay Mahadevan     }
939181f196bSVijay Mahadevan     else {
940181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 2-D entity vertices are invalid. Currently only support QUAD4 and TRI3 basis evaluations in 2-D : %D", nverts);
941181f196bSVijay Mahadevan     }
942181f196bSVijay Mahadevan 
943181f196bSVijay Mahadevan     /* invert the jacobian */
944181f196bSVijay Mahadevan     if (ijacobian) {
945a86ed394SVijay Mahadevan       ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr);
946181f196bSVijay Mahadevan     }
947181f196bSVijay Mahadevan 
948181f196bSVijay Mahadevan   }
949181f196bSVijay Mahadevan   else {
950181f196bSVijay Mahadevan 
951181f196bSVijay Mahadevan     if (nverts == 8) { /* Linear Hexahedra */
952181f196bSVijay Mahadevan       const PetscReal& r = quad[0];
953181f196bSVijay Mahadevan       const PetscReal& s = quad[1];
954181f196bSVijay Mahadevan       const PetscReal& t = quad[2];
955181f196bSVijay Mahadevan       const PetscReal dNi_dxi[8]  = {- ( 1.0 - s ) * ( 1.0 - t ),
956181f196bSVijay Mahadevan                                        ( 1.0 - s ) * ( 1.0 - t ),
957a86ed394SVijay Mahadevan                                        (       s ) * ( 1.0 - t ),
958a86ed394SVijay Mahadevan                                      - (       s ) * ( 1.0 - t ),
959a86ed394SVijay Mahadevan                                      - ( 1.0 - s ) * (       t ),
960a86ed394SVijay Mahadevan                                        ( 1.0 - s ) * (       t ),
961a86ed394SVijay Mahadevan                                        (       s ) * (       t ),
962a86ed394SVijay Mahadevan                                      - (       s ) * (       t )
963181f196bSVijay Mahadevan                                     };
964181f196bSVijay Mahadevan 
965181f196bSVijay Mahadevan       const PetscReal dNi_deta[8]  = { - ( 1.0 - r ) * ( 1.0 - t ),
966a86ed394SVijay Mahadevan                                        - (       r ) * ( 1.0 - t ),
967a86ed394SVijay Mahadevan                                          (       r ) * ( 1.0 - t ),
968181f196bSVijay Mahadevan                                          ( 1.0 - r ) * ( 1.0 - t ),
969a86ed394SVijay Mahadevan                                        - ( 1.0 - r ) * (       t ),
970a86ed394SVijay Mahadevan                                        - (       r ) * (       t ),
971a86ed394SVijay Mahadevan                                          (       r ) * (       t ),
972a86ed394SVijay Mahadevan                                          ( 1.0 - r ) * (       t )
973181f196bSVijay Mahadevan                                       };
974181f196bSVijay Mahadevan 
975181f196bSVijay Mahadevan       const PetscReal dNi_dzeta[8]  = { - ( 1.0 - r ) * ( 1.0 - s ),
976a86ed394SVijay Mahadevan                                         - (       r ) * ( 1.0 - s ),
977a86ed394SVijay Mahadevan                                         - (       r ) * (       s ),
978a86ed394SVijay Mahadevan                                         - ( 1.0 - r ) * (       s ),
979181f196bSVijay Mahadevan                                           ( 1.0 - r ) * ( 1.0 - s ),
980a86ed394SVijay Mahadevan                                           (       r ) * ( 1.0 - s ),
981a86ed394SVijay Mahadevan                                           (       r ) * (       s ),
982a86ed394SVijay Mahadevan                                           ( 1.0 - r ) * (       s )
983181f196bSVijay Mahadevan                                      };
984a86ed394SVijay Mahadevan 
985181f196bSVijay Mahadevan       for (i = 0; i < nverts; ++i) {
986181f196bSVijay Mahadevan         const PetscReal* vertex = coordinates + i * 3;
987181f196bSVijay Mahadevan         jacobian[0] += dNi_dxi[i]   * vertex[0];
988181f196bSVijay Mahadevan         jacobian[3] += dNi_dxi[i]   * vertex[1];
989181f196bSVijay Mahadevan         jacobian[6] += dNi_dxi[i]   * vertex[2];
990181f196bSVijay Mahadevan         jacobian[1] += dNi_deta[i]  * vertex[0];
991181f196bSVijay Mahadevan         jacobian[4] += dNi_deta[i]  * vertex[1];
992181f196bSVijay Mahadevan         jacobian[7] += dNi_deta[i]  * vertex[2];
993181f196bSVijay Mahadevan         jacobian[2] += dNi_dzeta[i] * vertex[0];
994181f196bSVijay Mahadevan         jacobian[5] += dNi_dzeta[i] * vertex[1];
995181f196bSVijay Mahadevan         jacobian[8] += dNi_dzeta[i] * vertex[2];
996181f196bSVijay Mahadevan       }
997181f196bSVijay Mahadevan 
998181f196bSVijay Mahadevan     }
999181f196bSVijay Mahadevan     else if (nverts == 4) { /* Linear Tetrahedra */
1000181f196bSVijay Mahadevan       const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2];
1001181f196bSVijay Mahadevan 
1002181f196bSVijay Mahadevan       /* compute the jacobian */
1003181f196bSVijay Mahadevan       jacobian[0] = coordinates[1 * 3 + 0] - x0;  jacobian[1] = coordinates[2 * 3 + 0] - x0; jacobian[2] = coordinates[3 * 3 + 0] - x0;
1004181f196bSVijay Mahadevan       jacobian[3] = coordinates[1 * 3 + 1] - y0;  jacobian[4] = coordinates[2 * 3 + 1] - y0; jacobian[5] = coordinates[3 * 3 + 1] - y0;
1005181f196bSVijay Mahadevan       jacobian[6] = coordinates[1 * 3 + 2] - z0;  jacobian[7] = coordinates[2 * 3 + 2] - z0; jacobian[8] = coordinates[3 * 3 + 2] - z0;
1006181f196bSVijay Mahadevan     } /* Tetrahedra -- ends */
1007181f196bSVijay Mahadevan     else
1008181f196bSVijay Mahadevan     {
1009181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 3-D entity vertices are invalid. Currently only support HEX8 and TET4 basis evaluations in 3-D : %D", nverts);
1010181f196bSVijay Mahadevan     }
1011181f196bSVijay Mahadevan 
1012181f196bSVijay Mahadevan     if (ijacobian) {
1013181f196bSVijay Mahadevan       /* invert the jacobian */
1014a86ed394SVijay Mahadevan       ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr);
1015181f196bSVijay Mahadevan     }
1016181f196bSVijay Mahadevan 
1017181f196bSVijay Mahadevan   }
1018a86ed394SVijay Mahadevan   if ( volume < 1e-12 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Element has zero volume: %g. Degenerate element or invalid connectivity\n", volume);
1019a86ed394SVijay Mahadevan   if (dvolume) *dvolume = volume;
1020181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1021181f196bSVijay Mahadevan }
1022181f196bSVijay Mahadevan 
1023181f196bSVijay Mahadevan 
1024181f196bSVijay Mahadevan PetscErrorCode FEMComputeBasis_JandF ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quadrature, PetscReal *phypts,
1025181f196bSVijay Mahadevan                                        PetscReal *phibasis, PetscReal *jacobian, PetscReal *ijacobian, PetscReal* volume  )
1026181f196bSVijay Mahadevan {
1027181f196bSVijay Mahadevan   PetscErrorCode  ierr;
1028181f196bSVijay Mahadevan 
1029181f196bSVijay Mahadevan   PetscFunctionBegin;
1030181f196bSVijay Mahadevan   switch (dim) {
1031181f196bSVijay Mahadevan     case 1:
1032181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts,
1033181f196bSVijay Mahadevan             NULL, phibasis, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
1034181f196bSVijay Mahadevan       break;
1035181f196bSVijay Mahadevan     case 2:
1036181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts,
1037181f196bSVijay Mahadevan             NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
1038181f196bSVijay Mahadevan       break;
1039181f196bSVijay Mahadevan     case 3:
1040181f196bSVijay Mahadevan       ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts,
1041181f196bSVijay Mahadevan             NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr);
1042181f196bSVijay Mahadevan       break;
1043181f196bSVijay Mahadevan     default:
1044181f196bSVijay Mahadevan       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim);
1045181f196bSVijay Mahadevan   }
1046181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1047181f196bSVijay Mahadevan }
1048181f196bSVijay Mahadevan 
1049181f196bSVijay Mahadevan 
1050181f196bSVijay Mahadevan 
1051*cab5ea25SPierre Jolivet /*@C
105297b73a88SSatish Balay   DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the
1053a86ed394SVijay Mahadevan   canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration,
1054a86ed394SVijay Mahadevan   the basis function at the parametric point is also evaluated optionally.
1055a86ed394SVijay Mahadevan 
1056a86ed394SVijay Mahadevan   Input Parameter:
1057a2b725a8SWilliam Gropp +  PetscInt  dim -         the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
1058a2b725a8SWilliam Gropp .  PetscInt nverts -       the number of vertices in the physical element
1059a2b725a8SWilliam Gropp .  PetscReal coordinates - the coordinates of vertices in the physical element
1060a2b725a8SWilliam Gropp -  PetscReal[3] xphy -     the coordinates of physical point for which natural coordinates (in reference frame) are sought
1061a86ed394SVijay Mahadevan 
1062a86ed394SVijay Mahadevan   Output Parameter:
1063a2b725a8SWilliam Gropp +  PetscReal[3] natparam - the natural coordinates (in reference frame) corresponding to xphy
1064a2b725a8SWilliam Gropp -  PetscReal[nverts] phi - the basis functions evaluated at the natural coordinates (natparam)
1065a86ed394SVijay Mahadevan 
1066edc382c3SSatish Balay   Level: advanced
1067edc382c3SSatish Balay 
1068a86ed394SVijay Mahadevan @*/
1069181f196bSVijay Mahadevan PetscErrorCode DMMoabPToRMapping( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal* xphy, PetscReal* natparam, PetscReal* phi)
1070181f196bSVijay Mahadevan {
1071a86ed394SVijay Mahadevan   /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */
1072181f196bSVijay Mahadevan   const PetscReal tol = 1.0e-10;
1073181f196bSVijay Mahadevan   const PetscInt max_iterations = 10;
1074181f196bSVijay Mahadevan   const PetscReal error_tol_sqr = tol*tol;
1075181f196bSVijay Mahadevan   PetscReal phibasis[8], jacobian[9], ijacobian[9], volume;
1076181f196bSVijay Mahadevan   PetscReal phypts[3] = {0.0, 0.0, 0.0};
1077181f196bSVijay Mahadevan   PetscReal delta[3] = {0.0, 0.0, 0.0};
1078181f196bSVijay Mahadevan   PetscErrorCode  ierr;
1079181f196bSVijay Mahadevan   PetscInt iters=0;
1080181f196bSVijay Mahadevan   PetscReal error=1.0;
1081181f196bSVijay Mahadevan 
1082181f196bSVijay Mahadevan   PetscFunctionBegin;
1083181f196bSVijay Mahadevan   PetscValidPointer(coordinates, 3);
1084181f196bSVijay Mahadevan   PetscValidPointer(xphy, 4);
1085181f196bSVijay Mahadevan   PetscValidPointer(natparam, 5);
1086181f196bSVijay Mahadevan 
1087580bdb30SBarry Smith   ierr = PetscArrayzero(jacobian, dim * dim);CHKERRQ(ierr);
1088580bdb30SBarry Smith   ierr = PetscArrayzero(ijacobian, dim * dim);CHKERRQ(ierr);
1089580bdb30SBarry Smith   ierr = PetscArrayzero(phibasis, nverts);CHKERRQ(ierr);
1090181f196bSVijay Mahadevan 
1091a86ed394SVijay Mahadevan   /* zero initial guess */
1092a86ed394SVijay Mahadevan   natparam[0] = natparam[1] = natparam[2] = 0.0;
1093181f196bSVijay Mahadevan 
1094a86ed394SVijay Mahadevan   /* Compute delta = evaluate( xi ) - x */
1095a86ed394SVijay Mahadevan   ierr = FEMComputeBasis_JandF ( dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume );CHKERRQ(ierr);
1096a86ed394SVijay Mahadevan 
1097a86ed394SVijay Mahadevan   error = 0.0;
1098a86ed394SVijay Mahadevan   switch(dim) {
1099a86ed394SVijay Mahadevan     case 3:
1100181f196bSVijay Mahadevan       delta[2] = phypts[2] - xphy[2];
1101a86ed394SVijay Mahadevan       error += (delta[2]*delta[2]);
1102a86ed394SVijay Mahadevan     case 2:
1103a86ed394SVijay Mahadevan       delta[1] = phypts[1] - xphy[1];
1104a86ed394SVijay Mahadevan       error += (delta[1]*delta[1]);
1105a86ed394SVijay Mahadevan     case 1:
1106a86ed394SVijay Mahadevan       delta[0] = phypts[0] - xphy[0];
1107a86ed394SVijay Mahadevan       error += (delta[0]*delta[0]);
1108a86ed394SVijay Mahadevan       break;
1109a86ed394SVijay Mahadevan   }
1110a86ed394SVijay Mahadevan 
1111a86ed394SVijay Mahadevan #if 0
1112a86ed394SVijay Mahadevan   PetscInfo4(NULL, "Current point in physical space : (%g, %g, %g) and error : %g\n", xphy[0], xphy[1], xphy[2], error);
1113a86ed394SVijay Mahadevan #endif
1114181f196bSVijay Mahadevan 
1115181f196bSVijay Mahadevan   while (error > error_tol_sqr) {
1116181f196bSVijay Mahadevan 
1117181f196bSVijay Mahadevan     if(++iters > max_iterations)
1118181f196bSVijay Mahadevan       SETERRQ3(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Maximum Newton iterations (10) reached. Current point in reference space : (%g, %g, %g)", natparam[0], natparam[1], natparam[2]);
1119181f196bSVijay Mahadevan 
1120181f196bSVijay Mahadevan     /* Compute natparam -= J.inverse() * delta */
1121181f196bSVijay Mahadevan     switch(dim) {
1122181f196bSVijay Mahadevan       case 1:
1123181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0];
1124181f196bSVijay Mahadevan         break;
1125181f196bSVijay Mahadevan       case 2:
1126181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1];
1127181f196bSVijay Mahadevan         natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1];
1128181f196bSVijay Mahadevan         break;
1129181f196bSVijay Mahadevan       case 3:
1130181f196bSVijay Mahadevan         natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2];
1131181f196bSVijay Mahadevan         natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2];
1132181f196bSVijay Mahadevan         natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2];
1133181f196bSVijay Mahadevan         break;
1134181f196bSVijay Mahadevan     }
1135181f196bSVijay Mahadevan 
1136181f196bSVijay Mahadevan     /* Compute delta = evaluate( xi ) - x */
1137a86ed394SVijay Mahadevan     ierr = FEMComputeBasis_JandF ( dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume );CHKERRQ(ierr);
1138181f196bSVijay Mahadevan 
1139a86ed394SVijay Mahadevan     error = 0.0;
1140a86ed394SVijay Mahadevan     switch(dim) {
1141a86ed394SVijay Mahadevan       case 3:
1142181f196bSVijay Mahadevan         delta[2] = phypts[2] - xphy[2];
1143a86ed394SVijay Mahadevan         error += (delta[2]*delta[2]);
1144a86ed394SVijay Mahadevan       case 2:
1145a86ed394SVijay Mahadevan         delta[1] = phypts[1] - xphy[1];
1146a86ed394SVijay Mahadevan         error += (delta[1]*delta[1]);
1147a86ed394SVijay Mahadevan       case 1:
1148a86ed394SVijay Mahadevan         delta[0] = phypts[0] - xphy[0];
1149a86ed394SVijay Mahadevan         error += (delta[0]*delta[0]);
1150a86ed394SVijay Mahadevan         break;
1151a86ed394SVijay Mahadevan     }
1152a86ed394SVijay Mahadevan #if 0
1153a86ed394SVijay Mahadevan     PetscInfo5(NULL, "Newton [%D] : Current point in reference space : (%g, %g, %g) and error : %g\n", iters, natparam[0], natparam[1], natparam[2], error);
1154a86ed394SVijay Mahadevan #endif
1155181f196bSVijay Mahadevan   }
1156181f196bSVijay Mahadevan   if (phi) {
1157580bdb30SBarry Smith     ierr = PetscArraycpy(phi, phibasis, nverts);CHKERRQ(ierr);
1158181f196bSVijay Mahadevan   }
1159181f196bSVijay Mahadevan   PetscFunctionReturn(0);
1160181f196bSVijay Mahadevan }
1161