163d025dbSVijay Mahadevan #include <petscconf.h>
263d025dbSVijay Mahadevan #include <petscdt.h> /*I "petscdt.h" I*/
363d025dbSVijay Mahadevan #include <petsc/private/dmmbimpl.h> /*I "petscdmmoab.h" I*/
463d025dbSVijay Mahadevan
563d025dbSVijay Mahadevan /* Utility functions */
DMatrix_Determinant_2x2_Internal(const PetscReal inmat[2* 2])6d71ae5a4SJacob Faibussowitsch static inline PetscReal DMatrix_Determinant_2x2_Internal(const PetscReal inmat[2 * 2])
7d71ae5a4SJacob Faibussowitsch {
863d025dbSVijay Mahadevan return inmat[0] * inmat[3] - inmat[1] * inmat[2];
963d025dbSVijay Mahadevan }
1063d025dbSVijay Mahadevan
DMatrix_Invert_2x2_Internal(const PetscReal * inmat,PetscReal * outmat,PetscReal * determinant)11d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMatrix_Invert_2x2_Internal(const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant)
12d71ae5a4SJacob Faibussowitsch {
1363d025dbSVijay Mahadevan PetscReal det = DMatrix_Determinant_2x2_Internal(inmat);
1463d025dbSVijay Mahadevan if (outmat) {
1563d025dbSVijay Mahadevan outmat[0] = inmat[3] / det;
1663d025dbSVijay Mahadevan outmat[1] = -inmat[1] / det;
1763d025dbSVijay Mahadevan outmat[2] = -inmat[2] / det;
1863d025dbSVijay Mahadevan outmat[3] = inmat[0] / det;
1963d025dbSVijay Mahadevan }
2063d025dbSVijay Mahadevan if (determinant) *determinant = det;
213ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
2263d025dbSVijay Mahadevan }
2363d025dbSVijay Mahadevan
DMatrix_Determinant_3x3_Internal(const PetscReal inmat[3* 3])24d71ae5a4SJacob Faibussowitsch static inline PetscReal DMatrix_Determinant_3x3_Internal(const PetscReal inmat[3 * 3])
25d71ae5a4SJacob Faibussowitsch {
269371c9d4SSatish 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]);
2763d025dbSVijay Mahadevan }
2863d025dbSVijay Mahadevan
DMatrix_Invert_3x3_Internal(const PetscReal * inmat,PetscReal * outmat,PetscReal * determinant)29ca4445c7SIlya Fursov static inline PetscErrorCode DMatrix_Invert_3x3_Internal(const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant)
30d71ae5a4SJacob Faibussowitsch {
31181f196bSVijay Mahadevan PetscReal det = DMatrix_Determinant_3x3_Internal(inmat);
3263d025dbSVijay Mahadevan if (outmat) {
3363d025dbSVijay Mahadevan outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det;
3463d025dbSVijay Mahadevan outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det;
3563d025dbSVijay Mahadevan outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det;
3663d025dbSVijay Mahadevan outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det;
3763d025dbSVijay Mahadevan outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det;
3863d025dbSVijay Mahadevan outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det;
3963d025dbSVijay Mahadevan outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det;
4063d025dbSVijay Mahadevan outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det;
4163d025dbSVijay Mahadevan outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det;
4263d025dbSVijay Mahadevan }
4363d025dbSVijay Mahadevan if (determinant) *determinant = det;
443ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
4563d025dbSVijay Mahadevan }
4663d025dbSVijay Mahadevan
DMatrix_Determinant_4x4_Internal(PetscReal inmat[4* 4])47a4e35b19SJacob Faibussowitsch static inline PetscReal DMatrix_Determinant_4x4_Internal(PetscReal inmat[4 * 4])
48d71ae5a4SJacob Faibussowitsch {
499371c9d4SSatish 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]));
50181f196bSVijay Mahadevan }
51181f196bSVijay Mahadevan
DMatrix_Invert_4x4_Internal(PetscReal * inmat,PetscReal * outmat,PetscScalar * determinant)52621c8970SIulian Grindeanu static inline PETSC_UNUSED PetscErrorCode DMatrix_Invert_4x4_Internal(PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant)
53d71ae5a4SJacob Faibussowitsch {
54181f196bSVijay Mahadevan PetscReal det = DMatrix_Determinant_4x4_Internal(inmat);
55181f196bSVijay Mahadevan if (outmat) {
56181f196bSVijay 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;
57181f196bSVijay 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;
58181f196bSVijay 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;
59181f196bSVijay 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;
60181f196bSVijay 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;
61181f196bSVijay 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;
62181f196bSVijay 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;
63181f196bSVijay 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;
64181f196bSVijay 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;
65181f196bSVijay 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;
66181f196bSVijay 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;
67181f196bSVijay 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;
68181f196bSVijay 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;
69181f196bSVijay 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;
70181f196bSVijay 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;
71181f196bSVijay 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;
72181f196bSVijay Mahadevan }
73181f196bSVijay Mahadevan if (determinant) *determinant = det;
743ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
75181f196bSVijay Mahadevan }
76181f196bSVijay Mahadevan
77a4e35b19SJacob Faibussowitsch /*
7897b73a88SSatish Balay Compute_Lagrange_Basis_1D_Internal - Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element.
7963d025dbSVijay Mahadevan
8063d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a linear or quadratic edge element.
8163d025dbSVijay Mahadevan
8263d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided,
8363d025dbSVijay Mahadevan and their derivatives with respect to X.
8463d025dbSVijay Mahadevan
8563d025dbSVijay Mahadevan Notes:
8663d025dbSVijay Mahadevan
8763d025dbSVijay Mahadevan Example Physical Element
88a2b725a8SWilliam Gropp .vb
8963d025dbSVijay Mahadevan 1-------2 1----3----2
9063d025dbSVijay Mahadevan EDGE2 EDGE3
91a2b725a8SWilliam Gropp .ve
9263d025dbSVijay Mahadevan
93d8d19677SJose E. Roman Input Parameters:
94a2b725a8SWilliam Gropp + PetscInt nverts - the number of element vertices
95a2b725a8SWilliam Gropp . PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
96a2b725a8SWilliam Gropp . PetscInt npts - the number of evaluation points (quadrature points)
97a2b725a8SWilliam Gropp - PetscReal quad[3*npts] - the evaluation points (quadrature points) in the reference space
9863d025dbSVijay Mahadevan
99d8d19677SJose E. Roman Output Parameters:
100a2b725a8SWilliam Gropp + PetscReal phypts[3*npts] - the evaluation points (quadrature points) transformed to the physical space
101a2b725a8SWilliam Gropp . PetscReal jxw[npts] - the jacobian determinant * quadrature weight necessary for assembling discrete contributions
102a2b725a8SWilliam Gropp . PetscReal phi[npts] - the bases evaluated at the specified quadrature points
1036b867d5aSJose E. Roman . PetscReal dphidx[npts] - the derivative of the bases wrt X-direction evaluated at the specified quadrature points
10467b8a455SSatish Balay . jacobian - jacobian
10567b8a455SSatish Balay . ijacobian - ijacobian
10667b8a455SSatish Balay - volume - volume
10763d025dbSVijay Mahadevan
108edc382c3SSatish Balay Level: advanced
109a4e35b19SJacob Faibussowitsch */
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)110a4e35b19SJacob 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)
111d71ae5a4SJacob Faibussowitsch {
11263d025dbSVijay Mahadevan int i, j;
11363d025dbSVijay Mahadevan
114181f196bSVijay Mahadevan PetscFunctionBegin;
1154f572ea9SToby Isaac PetscAssertPointer(jacobian, 9);
1164f572ea9SToby Isaac PetscAssertPointer(ijacobian, 10);
1174f572ea9SToby Isaac PetscAssertPointer(volume, 11);
11848a46eb9SPierre Jolivet if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
11963d025dbSVijay Mahadevan if (dphidx) { /* Reset arrays. */
1209566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(dphidx, npts * nverts));
12163d025dbSVijay Mahadevan }
12263d025dbSVijay Mahadevan if (nverts == 2) { /* Linear Edge */
12363d025dbSVijay Mahadevan
1242da392ccSBarry Smith for (j = 0; j < npts; j++) {
125a86ed394SVijay Mahadevan const PetscInt offset = j * nverts;
126181f196bSVijay Mahadevan const PetscReal r = quad[j];
12763d025dbSVijay Mahadevan
128181f196bSVijay Mahadevan phi[0 + offset] = (1.0 - r);
129181f196bSVijay Mahadevan phi[1 + offset] = (r);
13063d025dbSVijay Mahadevan
131181f196bSVijay Mahadevan const PetscReal dNi_dxi[2] = {-1.0, 1.0};
13263d025dbSVijay Mahadevan
133181f196bSVijay Mahadevan jacobian[0] = ijacobian[0] = volume[0] = 0.0;
13463d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) {
135181f196bSVijay Mahadevan const PetscReal *vertices = coords + i * 3;
136181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0];
137ad540459SPierre Jolivet if (phypts) phypts[3 * j + 0] += phi[i + offset] * vertices[0];
13863d025dbSVijay Mahadevan }
13963d025dbSVijay Mahadevan
14063d025dbSVijay Mahadevan /* invert the jacobian */
141181f196bSVijay Mahadevan *volume = jacobian[0];
142181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0];
143181f196bSVijay Mahadevan jxw[j] *= *volume;
14463d025dbSVijay Mahadevan
14563d025dbSVijay Mahadevan /* Divide by element jacobian. */
14663d025dbSVijay Mahadevan for (i = 0; i < nverts; i++) {
147181f196bSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0];
14863d025dbSVijay Mahadevan }
14963d025dbSVijay Mahadevan }
1502da392ccSBarry Smith } else if (nverts == 3) { /* Quadratic Edge */
15163d025dbSVijay Mahadevan
1522da392ccSBarry Smith for (j = 0; j < npts; j++) {
153a86ed394SVijay Mahadevan const PetscInt offset = j * nverts;
154181f196bSVijay Mahadevan const PetscReal r = quad[j];
15563d025dbSVijay Mahadevan
156181f196bSVijay Mahadevan phi[0 + offset] = 1.0 + r * (2.0 * r - 3.0);
157181f196bSVijay Mahadevan phi[1 + offset] = 4.0 * r * (1.0 - r);
158181f196bSVijay Mahadevan phi[2 + offset] = r * (2.0 * r - 1.0);
15963d025dbSVijay Mahadevan
160181f196bSVijay Mahadevan const PetscReal dNi_dxi[3] = {4 * r - 3.0, 4 * (1.0 - 2.0 * r), 4.0 * r - 1.0};
16163d025dbSVijay Mahadevan
162181f196bSVijay Mahadevan jacobian[0] = ijacobian[0] = volume[0] = 0.0;
16363d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) {
164181f196bSVijay Mahadevan const PetscReal *vertices = coords + i * 3;
165181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0];
166ad540459SPierre Jolivet if (phypts) phypts[3 * j + 0] += phi[i + offset] * vertices[0];
16763d025dbSVijay Mahadevan }
16863d025dbSVijay Mahadevan
16963d025dbSVijay Mahadevan /* invert the jacobian */
170181f196bSVijay Mahadevan *volume = jacobian[0];
171181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0];
172181f196bSVijay Mahadevan if (jxw) 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 }
17963a3b9bcSJacob 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);
1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
18163d025dbSVijay Mahadevan }
18263d025dbSVijay Mahadevan
183a4e35b19SJacob Faibussowitsch /*
18497b73a88SSatish Balay Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element.
18563d025dbSVijay Mahadevan
18663d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a quadrangle or triangle.
18763d025dbSVijay Mahadevan
18863d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided,
18963d025dbSVijay Mahadevan and their derivatives with respect to X and Y.
19063d025dbSVijay Mahadevan
19163d025dbSVijay Mahadevan Notes:
19263d025dbSVijay Mahadevan
19363d025dbSVijay Mahadevan Example Physical Element (QUAD4)
194a2b725a8SWilliam Gropp .vb
19563d025dbSVijay Mahadevan 4------3 s
19663d025dbSVijay Mahadevan | | |
19763d025dbSVijay Mahadevan | | |
19863d025dbSVijay Mahadevan | | |
19963d025dbSVijay Mahadevan 1------2 0-------r
200a2b725a8SWilliam Gropp .ve
20163d025dbSVijay Mahadevan
202d8d19677SJose E. Roman Input Parameters:
203a2b725a8SWilliam Gropp + PetscInt nverts - the number of element vertices
204a2b725a8SWilliam Gropp . PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
205a2b725a8SWilliam Gropp . PetscInt npts - the number of evaluation points (quadrature points)
206a2b725a8SWilliam Gropp - PetscReal quad[3*npts] - the evaluation points (quadrature points) in the reference space
20763d025dbSVijay Mahadevan
208d8d19677SJose E. Roman Output Parameters:
209a2b725a8SWilliam Gropp + PetscReal phypts[3*npts] - the evaluation points (quadrature points) transformed to the physical space
210a2b725a8SWilliam Gropp . PetscReal jxw[npts] - the jacobian determinant * quadrature weight necessary for assembling discrete contributions
211a2b725a8SWilliam Gropp . PetscReal phi[npts] - the bases evaluated at the specified quadrature points
212a2b725a8SWilliam Gropp . PetscReal dphidx[npts] - the derivative of the bases wrt X-direction evaluated at the specified quadrature points
2136b867d5aSJose E. Roman . PetscReal dphidy[npts] - the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
21467b8a455SSatish Balay . jacobian - jacobian
21567b8a455SSatish Balay . ijacobian - ijacobian
21667b8a455SSatish Balay - volume - volume
21763d025dbSVijay Mahadevan
218edc382c3SSatish Balay Level: advanced
219a4e35b19SJacob Faibussowitsch */
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)220a4e35b19SJacob 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)
221d71ae5a4SJacob Faibussowitsch {
222a86ed394SVijay Mahadevan PetscInt i, j, k;
22363d025dbSVijay Mahadevan
22463d025dbSVijay Mahadevan PetscFunctionBegin;
2254f572ea9SToby Isaac PetscAssertPointer(jacobian, 10);
2264f572ea9SToby Isaac PetscAssertPointer(ijacobian, 11);
2274f572ea9SToby Isaac PetscAssertPointer(volume, 12);
2289566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(phi, npts));
22948a46eb9SPierre Jolivet if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
23063d025dbSVijay Mahadevan if (dphidx) { /* Reset arrays. */
2319566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(dphidx, npts * nverts));
2329566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(dphidy, npts * nverts));
23363d025dbSVijay Mahadevan }
23463d025dbSVijay Mahadevan if (nverts == 4) { /* Linear Quadrangle */
23563d025dbSVijay Mahadevan
2369371c9d4SSatish Balay for (j = 0; j < npts; j++) {
237a86ed394SVijay Mahadevan const PetscInt offset = j * nverts;
238181f196bSVijay Mahadevan const PetscReal r = quad[0 + j * 2];
239181f196bSVijay Mahadevan const PetscReal s = quad[1 + j * 2];
24063d025dbSVijay Mahadevan
24163d025dbSVijay Mahadevan phi[0 + offset] = (1.0 - r) * (1.0 - s);
24263d025dbSVijay Mahadevan phi[1 + offset] = r * (1.0 - s);
24363d025dbSVijay Mahadevan phi[2 + offset] = r * s;
24463d025dbSVijay Mahadevan phi[3 + offset] = (1.0 - r) * s;
24563d025dbSVijay Mahadevan
246181f196bSVijay Mahadevan const PetscReal dNi_dxi[4] = {-1.0 + s, 1.0 - s, s, -s};
247181f196bSVijay Mahadevan const PetscReal dNi_deta[4] = {-1.0 + r, -r, r, 1.0 - r};
24863d025dbSVijay Mahadevan
2499566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(jacobian, 4));
2509566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ijacobian, 4));
25163d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) {
252181f196bSVijay Mahadevan const PetscReal *vertices = coords + i * 3;
25363d025dbSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0];
25463d025dbSVijay Mahadevan jacobian[2] += dNi_dxi[i] * vertices[1];
25563d025dbSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertices[0];
25663d025dbSVijay Mahadevan jacobian[3] += dNi_deta[i] * vertices[1];
257181f196bSVijay Mahadevan if (phypts) {
258181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0];
259181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertices[1];
260181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertices[2];
261181f196bSVijay Mahadevan }
26263d025dbSVijay Mahadevan }
26363d025dbSVijay Mahadevan
26463d025dbSVijay Mahadevan /* invert the jacobian */
2659566063dSJacob Faibussowitsch PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
26608401ef6SPierre Jolivet PetscCheck(*volume >= 1e-12, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Quadrangular element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
26763d025dbSVijay Mahadevan
268181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume;
26963d025dbSVijay Mahadevan
270181f196bSVijay Mahadevan /* Let us compute the bases derivatives by scaling with inverse jacobian. */
271181f196bSVijay Mahadevan if (dphidx) {
27263d025dbSVijay Mahadevan for (i = 0; i < nverts; i++) {
273a86ed394SVijay Mahadevan for (k = 0; k < 2; ++k) {
27463d025dbSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0];
27563d025dbSVijay Mahadevan if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1];
276181f196bSVijay Mahadevan } /* for k=[0..2) */
277181f196bSVijay Mahadevan } /* for i=[0..nverts) */
278181f196bSVijay Mahadevan } /* if (dphidx) */
27963d025dbSVijay Mahadevan }
2802da392ccSBarry Smith } else if (nverts == 3) { /* Linear triangle */
2812da392ccSBarry Smith const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1];
28263d025dbSVijay Mahadevan
2839566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(jacobian, 4));
2849566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ijacobian, 4));
28563d025dbSVijay Mahadevan
28663d025dbSVijay Mahadevan /* Jacobian is constant */
2879371c9d4SSatish Balay jacobian[0] = (coords[0 * 3 + 0] - x2);
2889371c9d4SSatish Balay jacobian[1] = (coords[1 * 3 + 0] - x2);
2899371c9d4SSatish Balay jacobian[2] = (coords[0 * 3 + 1] - y2);
2909371c9d4SSatish Balay jacobian[3] = (coords[1 * 3 + 1] - y2);
29163d025dbSVijay Mahadevan
29263d025dbSVijay Mahadevan /* invert the jacobian */
2939566063dSJacob Faibussowitsch PetscCall(DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume));
29408401ef6SPierre 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);
295181f196bSVijay Mahadevan
296181f196bSVijay Mahadevan const PetscReal Dx[3] = {ijacobian[0], ijacobian[2], -ijacobian[0] - ijacobian[2]};
297181f196bSVijay Mahadevan const PetscReal Dy[3] = {ijacobian[1], ijacobian[3], -ijacobian[1] - ijacobian[3]};
29863d025dbSVijay Mahadevan
2992da392ccSBarry Smith for (j = 0; j < npts; j++) {
300a86ed394SVijay Mahadevan const PetscInt offset = j * nverts;
301181f196bSVijay Mahadevan const PetscReal r = quad[0 + j * 2];
302181f196bSVijay Mahadevan const PetscReal s = quad[1 + j * 2];
30363d025dbSVijay Mahadevan
304181f196bSVijay Mahadevan if (jxw) jxw[j] *= 0.5;
30563d025dbSVijay Mahadevan
306181f196bSVijay Mahadevan /* const PetscReal Ni[3] = { r, s, 1.0 - r - s }; */
307181f196bSVijay Mahadevan const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s;
308181f196bSVijay Mahadevan const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s;
309181f196bSVijay Mahadevan if (phypts) {
310181f196bSVijay Mahadevan phypts[offset + 0] = phipts_x;
311181f196bSVijay Mahadevan phypts[offset + 1] = phipts_y;
312181f196bSVijay Mahadevan }
31363d025dbSVijay Mahadevan
314110fc3b0SBarry Smith /* \phi_0 = (b.y - c.y) x + (b.x - c.x) y + c.x b.y - b.x c.y */
315181f196bSVijay Mahadevan phi[0 + offset] = (ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2));
316110fc3b0SBarry Smith /* \phi_1 = (c.y - a.y) x + (a.x - c.x) y + c.x a.y - a.x c.y */
317181f196bSVijay Mahadevan phi[1 + offset] = (ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2));
31863d025dbSVijay Mahadevan phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset];
31963d025dbSVijay Mahadevan
32063d025dbSVijay Mahadevan if (dphidx) {
321181f196bSVijay Mahadevan dphidx[0 + offset] = Dx[0];
322181f196bSVijay Mahadevan dphidx[1 + offset] = Dx[1];
323181f196bSVijay Mahadevan dphidx[2 + offset] = Dx[2];
32463d025dbSVijay Mahadevan }
32563d025dbSVijay Mahadevan
32663d025dbSVijay Mahadevan if (dphidy) {
327181f196bSVijay Mahadevan dphidy[0 + offset] = Dy[0];
328181f196bSVijay Mahadevan dphidy[1 + offset] = Dy[1];
329181f196bSVijay Mahadevan dphidy[2 + offset] = Dy[2];
33063d025dbSVijay Mahadevan }
33163d025dbSVijay Mahadevan }
33263a3b9bcSJacob 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);
3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
33463d025dbSVijay Mahadevan }
33563d025dbSVijay Mahadevan
336a4e35b19SJacob Faibussowitsch /*
33797b73a88SSatish Balay Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element.
33863d025dbSVijay Mahadevan
33963d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a hexahedra or tetrahedra.
34063d025dbSVijay Mahadevan
34163d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided,
34263d025dbSVijay Mahadevan and their derivatives with respect to X, Y and Z.
34363d025dbSVijay Mahadevan
34463d025dbSVijay Mahadevan Notes:
34563d025dbSVijay Mahadevan
34663d025dbSVijay Mahadevan Example Physical Element (HEX8)
347a2b725a8SWilliam Gropp .vb
34863d025dbSVijay Mahadevan 8------7
34963d025dbSVijay Mahadevan /| /| t s
35063d025dbSVijay Mahadevan 5------6 | | /
35163d025dbSVijay Mahadevan | | | | |/
35263d025dbSVijay Mahadevan | 4----|-3 0-------r
35363d025dbSVijay Mahadevan |/ |/
35463d025dbSVijay Mahadevan 1------2
355a2b725a8SWilliam Gropp .ve
35663d025dbSVijay Mahadevan
357d8d19677SJose E. Roman Input Parameters:
358a2b725a8SWilliam Gropp + PetscInt nverts - the number of element vertices
359a2b725a8SWilliam Gropp . PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering)
360a2b725a8SWilliam Gropp . PetscInt npts - the number of evaluation points (quadrature points)
361a2b725a8SWilliam Gropp - PetscReal quad[3*npts] - the evaluation points (quadrature points) in the reference space
36263d025dbSVijay Mahadevan
363d8d19677SJose E. Roman Output Parameters:
364a2b725a8SWilliam Gropp + PetscReal phypts[3*npts] - the evaluation points (quadrature points) transformed to the physical space
365a2b725a8SWilliam Gropp . PetscReal jxw[npts] - the jacobian determinant * quadrature weight necessary for assembling discrete contributions
366a2b725a8SWilliam Gropp . PetscReal phi[npts] - the bases evaluated at the specified quadrature points
367a2b725a8SWilliam Gropp . PetscReal dphidx[npts] - the derivative of the bases wrt X-direction evaluated at the specified quadrature points
368a2b725a8SWilliam Gropp . PetscReal dphidy[npts] - the derivative of the bases wrt Y-direction evaluated at the specified quadrature points
3696b867d5aSJose E. Roman . PetscReal dphidz[npts] - the derivative of the bases wrt Z-direction evaluated at the specified quadrature points
37067b8a455SSatish Balay . jacobian - jacobian
37167b8a455SSatish Balay . ijacobian - ijacobian
37267b8a455SSatish Balay - volume - volume
37363d025dbSVijay Mahadevan
374edc382c3SSatish Balay Level: advanced
375a4e35b19SJacob Faibussowitsch */
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)376a4e35b19SJacob 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)
377d71ae5a4SJacob Faibussowitsch {
378a86ed394SVijay Mahadevan PetscInt i, j, k;
37963d025dbSVijay Mahadevan
38063d025dbSVijay Mahadevan PetscFunctionBegin;
3814f572ea9SToby Isaac PetscAssertPointer(jacobian, 11);
3824f572ea9SToby Isaac PetscAssertPointer(ijacobian, 12);
3834f572ea9SToby Isaac PetscAssertPointer(volume, 13);
3842da392ccSBarry Smith
3859566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(phi, npts));
38648a46eb9SPierre Jolivet if (phypts) PetscCall(PetscArrayzero(phypts, npts * 3));
38763d025dbSVijay Mahadevan if (dphidx) {
3889566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(dphidx, npts * nverts));
3899566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(dphidy, npts * nverts));
3909566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(dphidz, npts * nverts));
39163d025dbSVijay Mahadevan }
39263d025dbSVijay Mahadevan
39363d025dbSVijay Mahadevan if (nverts == 8) { /* Linear Hexahedra */
39463d025dbSVijay Mahadevan
3952da392ccSBarry Smith for (j = 0; j < npts; j++) {
396a86ed394SVijay Mahadevan const PetscInt offset = j * nverts;
397181f196bSVijay Mahadevan const PetscReal &r = quad[j * 3 + 0];
398181f196bSVijay Mahadevan const PetscReal &s = quad[j * 3 + 1];
399181f196bSVijay Mahadevan const PetscReal &t = quad[j * 3 + 2];
40063d025dbSVijay Mahadevan
401a86ed394SVijay Mahadevan phi[offset + 0] = (1.0 - r) * (1.0 - s) * (1.0 - t); /* 0,0,0 */
402a86ed394SVijay Mahadevan phi[offset + 1] = (r) * (1.0 - s) * (1.0 - t); /* 1,0,0 */
403a86ed394SVijay Mahadevan phi[offset + 2] = (r) * (s) * (1.0 - t); /* 1,1,0 */
404a86ed394SVijay Mahadevan phi[offset + 3] = (1.0 - r) * (s) * (1.0 - t); /* 0,1,0 */
405a86ed394SVijay Mahadevan phi[offset + 4] = (1.0 - r) * (1.0 - s) * (t); /* 0,0,1 */
406a86ed394SVijay Mahadevan phi[offset + 5] = (r) * (1.0 - s) * (t); /* 1,0,1 */
407a86ed394SVijay Mahadevan phi[offset + 6] = (r) * (s) * (t); /* 1,1,1 */
408a86ed394SVijay Mahadevan phi[offset + 7] = (1.0 - r) * (s) * (t); /* 0,1,1 */
40963d025dbSVijay Mahadevan
4109371c9d4SSatish 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)};
41163d025dbSVijay Mahadevan
4129371c9d4SSatish 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)};
41363d025dbSVijay Mahadevan
4149371c9d4SSatish 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)};
41563d025dbSVijay Mahadevan
4169566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(jacobian, 9));
4179566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ijacobian, 9));
41863d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) {
419181f196bSVijay Mahadevan const PetscReal *vertex = coords + i * 3;
42063d025dbSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertex[0];
42163d025dbSVijay Mahadevan jacobian[3] += dNi_dxi[i] * vertex[1];
42263d025dbSVijay Mahadevan jacobian[6] += dNi_dxi[i] * vertex[2];
42363d025dbSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertex[0];
42463d025dbSVijay Mahadevan jacobian[4] += dNi_deta[i] * vertex[1];
42563d025dbSVijay Mahadevan jacobian[7] += dNi_deta[i] * vertex[2];
42663d025dbSVijay Mahadevan jacobian[2] += dNi_dzeta[i] * vertex[0];
42763d025dbSVijay Mahadevan jacobian[5] += dNi_dzeta[i] * vertex[1];
42863d025dbSVijay Mahadevan jacobian[8] += dNi_dzeta[i] * vertex[2];
429181f196bSVijay Mahadevan if (phypts) {
430181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertex[0];
431181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertex[1];
432181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertex[2];
433181f196bSVijay Mahadevan }
43463d025dbSVijay Mahadevan }
43563d025dbSVijay Mahadevan
43663d025dbSVijay Mahadevan /* invert the jacobian */
4379566063dSJacob Faibussowitsch PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
43808401ef6SPierre Jolivet PetscCheck(*volume >= 1e-8, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Hexahedral element has zero volume: %g. Degenerate element or invalid connectivity", *volume);
43963d025dbSVijay Mahadevan
440a86ed394SVijay Mahadevan if (jxw) jxw[j] *= (*volume);
44163d025dbSVijay Mahadevan
44263d025dbSVijay Mahadevan /* Divide by element jacobian. */
44363d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) {
444a86ed394SVijay Mahadevan for (k = 0; k < 3; ++k) {
44563d025dbSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0 * 3 + k];
44663d025dbSVijay Mahadevan if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[1 * 3 + k];
44763d025dbSVijay Mahadevan if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k];
44863d025dbSVijay Mahadevan }
44963d025dbSVijay Mahadevan }
45063d025dbSVijay Mahadevan }
4512da392ccSBarry Smith } else if (nverts == 4) { /* Linear Tetrahedra */
4522da392ccSBarry Smith PetscReal Dx[4] = {0, 0, 0, 0}, Dy[4] = {0, 0, 0, 0}, Dz[4] = {0, 0, 0, 0};
4532da392ccSBarry Smith const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2];
45463d025dbSVijay Mahadevan
4559566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(jacobian, 9));
4569566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ijacobian, 9));
45763d025dbSVijay Mahadevan
458181f196bSVijay Mahadevan /* compute the jacobian */
4599371c9d4SSatish Balay jacobian[0] = coords[1 * 3 + 0] - x0;
4609371c9d4SSatish Balay jacobian[1] = coords[2 * 3 + 0] - x0;
4619371c9d4SSatish Balay jacobian[2] = coords[3 * 3 + 0] - x0;
4629371c9d4SSatish Balay jacobian[3] = coords[1 * 3 + 1] - y0;
4639371c9d4SSatish Balay jacobian[4] = coords[2 * 3 + 1] - y0;
4649371c9d4SSatish Balay jacobian[5] = coords[3 * 3 + 1] - y0;
4659371c9d4SSatish Balay jacobian[6] = coords[1 * 3 + 2] - z0;
4669371c9d4SSatish Balay jacobian[7] = coords[2 * 3 + 2] - z0;
4679371c9d4SSatish Balay jacobian[8] = coords[3 * 3 + 2] - z0;
46863d025dbSVijay Mahadevan
46963d025dbSVijay Mahadevan /* invert the jacobian */
4709566063dSJacob Faibussowitsch PetscCall(DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume));
47108401ef6SPierre 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);
47263d025dbSVijay Mahadevan
473181f196bSVijay Mahadevan if (dphidx) {
4749371c9d4SSatish 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;
4759371c9d4SSatish 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;
4769371c9d4SSatish 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;
477181f196bSVijay Mahadevan Dx[3] = -(Dx[0] + Dx[1] + Dx[2]);
4789371c9d4SSatish 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;
4799371c9d4SSatish 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;
4809371c9d4SSatish 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;
481181f196bSVijay Mahadevan Dy[3] = -(Dy[0] + Dy[1] + Dy[2]);
4829371c9d4SSatish 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;
4839371c9d4SSatish 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;
4849371c9d4SSatish 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;
485181f196bSVijay Mahadevan Dz[3] = -(Dz[0] + Dz[1] + Dz[2]);
486181f196bSVijay Mahadevan }
48763d025dbSVijay Mahadevan
4882da392ccSBarry Smith for (j = 0; j < npts; j++) {
489a86ed394SVijay Mahadevan const PetscInt offset = j * nverts;
490181f196bSVijay Mahadevan const PetscReal &r = quad[j * 3 + 0];
491181f196bSVijay Mahadevan const PetscReal &s = quad[j * 3 + 1];
492181f196bSVijay Mahadevan const PetscReal &t = quad[j * 3 + 2];
49363d025dbSVijay Mahadevan
494181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume;
49563d025dbSVijay Mahadevan
49663d025dbSVijay Mahadevan phi[offset + 0] = 1.0 - r - s - t;
49763d025dbSVijay Mahadevan phi[offset + 1] = r;
49863d025dbSVijay Mahadevan phi[offset + 2] = s;
49963d025dbSVijay Mahadevan phi[offset + 3] = t;
50063d025dbSVijay Mahadevan
501181f196bSVijay Mahadevan if (phypts) {
502181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) {
503ca4445c7SIlya Fursov const PetscReal *vertices = coords + i * 3;
504181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0];
505181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertices[1];
506181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertices[2];
507181f196bSVijay Mahadevan }
508181f196bSVijay Mahadevan }
509181f196bSVijay Mahadevan
510181f196bSVijay Mahadevan /* Now set the derivatives */
51163d025dbSVijay Mahadevan if (dphidx) {
512181f196bSVijay Mahadevan dphidx[0 + offset] = Dx[0];
513181f196bSVijay Mahadevan dphidx[1 + offset] = Dx[1];
514181f196bSVijay Mahadevan dphidx[2 + offset] = Dx[2];
515181f196bSVijay Mahadevan dphidx[3 + offset] = Dx[3];
51663d025dbSVijay Mahadevan
517181f196bSVijay Mahadevan dphidy[0 + offset] = Dy[0];
518181f196bSVijay Mahadevan dphidy[1 + offset] = Dy[1];
519181f196bSVijay Mahadevan dphidy[2 + offset] = Dy[2];
520181f196bSVijay Mahadevan dphidy[3 + offset] = Dy[3];
52163d025dbSVijay Mahadevan
522181f196bSVijay Mahadevan dphidz[0 + offset] = Dz[0];
523181f196bSVijay Mahadevan dphidz[1 + offset] = Dz[1];
524181f196bSVijay Mahadevan dphidz[2 + offset] = Dz[2];
525181f196bSVijay Mahadevan dphidz[3 + offset] = Dz[3];
52663d025dbSVijay Mahadevan }
52763d025dbSVijay Mahadevan
52863d025dbSVijay Mahadevan } /* Tetrahedra -- ends */
52963a3b9bcSJacob 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);
5303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
53163d025dbSVijay Mahadevan }
53263d025dbSVijay Mahadevan
533cab5ea25SPierre Jolivet /*@C
53497b73a88SSatish Balay DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element.
53563d025dbSVijay Mahadevan
53663d025dbSVijay Mahadevan The routine takes the coordinates of the vertices of an element and computes basis functions associated with
53763d025dbSVijay Mahadevan each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate.
53863d025dbSVijay Mahadevan
539d8d19677SJose E. Roman Input Parameters:
540a4e35b19SJacob Faibussowitsch + dim - the dimension
541a4e35b19SJacob Faibussowitsch . nverts - the number of element vertices
542a4e35b19SJacob Faibussowitsch . coordinates - the physical coordinates of the vertices (in canonical numbering)
543a4e35b19SJacob Faibussowitsch - quadrature - the evaluation points (quadrature points) in the reference space
54463d025dbSVijay Mahadevan
545d8d19677SJose E. Roman Output Parameters:
546a4e35b19SJacob Faibussowitsch + phypts - the evaluation points (quadrature points) transformed to the physical space
547a4e35b19SJacob Faibussowitsch . jacobian_quadrature_weight_product - the jacobian determinant * quadrature weight necessary for assembling discrete contributions
548a4e35b19SJacob Faibussowitsch . fe_basis - the bases values evaluated at the specified quadrature points
549a4e35b19SJacob Faibussowitsch - fe_basis_derivatives - the derivative of the bases wrt (X,Y,Z)-directions (depending on the dimension) evaluated at the specified quadrature points
55063d025dbSVijay Mahadevan
551edc382c3SSatish Balay Level: advanced
552edc382c3SSatish Balay
553a4e35b19SJacob Faibussowitsch .seealso: `DMMoabCreate()`
55463d025dbSVijay Mahadevan @*/
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)555d71ae5a4SJacob 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)
556d71ae5a4SJacob Faibussowitsch {
557b21a1e88SVijay Mahadevan PetscInt npoints, idim;
55863d025dbSVijay Mahadevan bool compute_der;
55963d025dbSVijay Mahadevan const PetscReal *quadpts, *quadwts;
560181f196bSVijay Mahadevan PetscReal jacobian[9], ijacobian[9], volume;
56163d025dbSVijay Mahadevan
56263d025dbSVijay Mahadevan PetscFunctionBegin;
5634f572ea9SToby Isaac PetscAssertPointer(coordinates, 3);
56478dc7ee3SMatthew G. Knepley PetscValidHeaderSpecific(quadrature, PETSCQUADRATURE_CLASSID, 4);
5654f572ea9SToby Isaac PetscAssertPointer(fe_basis, 7);
56663d025dbSVijay Mahadevan compute_der = (fe_basis_derivatives != NULL);
56763d025dbSVijay Mahadevan
56863d025dbSVijay Mahadevan /* Get the quadrature points and weights for the given quadrature rule */
5699566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts));
57063a3b9bcSJacob Faibussowitsch PetscCheck(idim == dim, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%" PetscInt_FMT ") vs quadrature (%" PetscInt_FMT ")", idim, dim);
57148a46eb9SPierre Jolivet if (jacobian_quadrature_weight_product) PetscCall(PetscArraycpy(jacobian_quadrature_weight_product, quadwts, npoints));
57263d025dbSVijay Mahadevan
57363d025dbSVijay Mahadevan switch (dim) {
574d71ae5a4SJacob Faibussowitsch case 1:
575*57508eceSPierre Jolivet 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));
576d71ae5a4SJacob Faibussowitsch break;
57763d025dbSVijay Mahadevan case 2:
578*57508eceSPierre Jolivet 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));
57963d025dbSVijay Mahadevan break;
58063d025dbSVijay Mahadevan case 3:
581*57508eceSPierre Jolivet 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));
58263d025dbSVijay Mahadevan break;
583d71ae5a4SJacob Faibussowitsch default:
584d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
58563d025dbSVijay Mahadevan }
5863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
58763d025dbSVijay Mahadevan }
58863d025dbSVijay Mahadevan
589cab5ea25SPierre Jolivet /*@C
59097b73a88SSatish Balay DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given
59163d025dbSVijay Mahadevan dimension and polynomial order (deciphered from number of element vertices).
59263d025dbSVijay Mahadevan
593d8d19677SJose E. Roman Input Parameters:
594a4e35b19SJacob Faibussowitsch + dim - the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
595a4e35b19SJacob Faibussowitsch - nverts - the number of vertices in the physical element
59663d025dbSVijay Mahadevan
59763d025dbSVijay Mahadevan Output Parameter:
598a4e35b19SJacob Faibussowitsch . quadrature - the quadrature object with default settings to integrate polynomials defined over the element
59963d025dbSVijay Mahadevan
600edc382c3SSatish Balay Level: advanced
601edc382c3SSatish Balay
602a4e35b19SJacob Faibussowitsch .seealso: `DMMoabCreate()`
60363d025dbSVijay Mahadevan @*/
DMMoabFEMCreateQuadratureDefault(const PetscInt dim,const PetscInt nverts,PetscQuadrature * quadrature)604d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature)
605d71ae5a4SJacob Faibussowitsch {
60663d025dbSVijay Mahadevan PetscReal *w, *x;
607b21a1e88SVijay Mahadevan PetscInt nc = 1;
60863d025dbSVijay Mahadevan
60963d025dbSVijay Mahadevan PetscFunctionBegin;
61063d025dbSVijay Mahadevan /* Create an appropriate quadrature rule to sample basis */
6119371c9d4SSatish Balay switch (dim) {
61263d025dbSVijay Mahadevan case 1:
61363d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
6149566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(1, nc, nverts, 0, 1.0, quadrature));
61563d025dbSVijay Mahadevan break;
61663d025dbSVijay Mahadevan case 2:
61763d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
61863d025dbSVijay Mahadevan if (nverts == 3) {
619a86ed394SVijay Mahadevan const PetscInt order = 2;
620a86ed394SVijay Mahadevan const PetscInt npoints = (order == 2 ? 3 : 6);
6219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(npoints * 2, &x, npoints, &w));
622181f196bSVijay Mahadevan if (npoints == 3) {
62363d025dbSVijay Mahadevan x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0;
62463d025dbSVijay Mahadevan x[3] = x[4] = 2.0 / 3.0;
62563d025dbSVijay Mahadevan w[0] = w[1] = w[2] = 1.0 / 3.0;
6262da392ccSBarry Smith } else if (npoints == 6) {
62763d025dbSVijay Mahadevan x[0] = x[1] = x[2] = 0.44594849091597;
62863d025dbSVijay Mahadevan x[3] = x[4] = 0.10810301816807;
62963d025dbSVijay Mahadevan x[5] = 0.44594849091597;
63063d025dbSVijay Mahadevan x[6] = x[7] = x[8] = 0.09157621350977;
63163d025dbSVijay Mahadevan x[9] = x[10] = 0.81684757298046;
63263d025dbSVijay Mahadevan x[11] = 0.09157621350977;
63363d025dbSVijay Mahadevan w[0] = w[1] = w[2] = 0.22338158967801;
634181f196bSVijay Mahadevan w[3] = w[4] = w[5] = 0.10995174365532;
63563a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %" PetscInt_FMT, npoints);
6369566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
6379566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetOrder(*quadrature, order));
6389566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
6399566063dSJacob Faibussowitsch /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
6401baa6e33SBarry Smith } else PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
64163d025dbSVijay Mahadevan break;
64263d025dbSVijay Mahadevan case 3:
64363d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */
64463d025dbSVijay Mahadevan if (nverts == 4) {
645a86ed394SVijay Mahadevan PetscInt order;
646a86ed394SVijay Mahadevan const PetscInt npoints = 4; // Choose between 4 and 10
6479566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(npoints * 3, &x, npoints, &w));
648181f196bSVijay Mahadevan if (npoints == 4) { /* KEAST1, K1, N=4, O=4 */
6499371c9d4SSatish Balay const PetscReal x_4[12] = {0.5854101966249685, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 0.1381966011250105,
6509371c9d4SSatish Balay 0.1381966011250105, 0.1381966011250105, 0.5854101966249685, 0.1381966011250105, 0.5854101966249685, 0.1381966011250105};
6519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(x, x_4, 12));
652181f196bSVijay Mahadevan
653181f196bSVijay Mahadevan w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0;
654181f196bSVijay Mahadevan order = 4;
6552da392ccSBarry Smith } else if (npoints == 10) { /* KEAST3, K3 N=10, O=10 */
6569371c9d4SSatish 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,
6579371c9d4SSatish Balay 0.5684305841968444, 0.1438564719343852, 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 0.5000000000000000, 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 0.5000000000000000,
6589371c9d4SSatish Balay 0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 0.0000000000000000, 0.5000000000000000};
6599566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(x, x_10, 30));
660181f196bSVijay Mahadevan
661181f196bSVijay Mahadevan w[0] = w[1] = w[2] = w[3] = 0.2177650698804054;
662181f196bSVijay Mahadevan w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631;
663181f196bSVijay Mahadevan order = 10;
66463a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %" PetscInt_FMT, npoints);
6659566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, quadrature));
6669566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetOrder(*quadrature, order));
6679566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w));
6689566063dSJacob Faibussowitsch /* PetscCall(PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature)); */
6691baa6e33SBarry Smith } else PetscCall(PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature));
67063d025dbSVijay Mahadevan break;
671d71ae5a4SJacob Faibussowitsch default:
672d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
67363d025dbSVijay Mahadevan }
6743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
67563d025dbSVijay Mahadevan }
67663d025dbSVijay Mahadevan
677181f196bSVijay Mahadevan /* Compute Jacobians */
FEMComputeBasis_JandF(const PetscInt dim,const PetscInt nverts,const PetscReal * coordinates,const PetscReal * quadrature,PetscReal * phypts,PetscReal * phibasis,PetscReal * jacobian,PetscReal * ijacobian,PetscReal * volume)678ba38deedSJacob Faibussowitsch static 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)
679d71ae5a4SJacob Faibussowitsch {
680181f196bSVijay Mahadevan PetscFunctionBegin;
681181f196bSVijay Mahadevan switch (dim) {
682d71ae5a4SJacob Faibussowitsch case 1:
683d71ae5a4SJacob Faibussowitsch PetscCall(Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, jacobian, ijacobian, volume));
684d71ae5a4SJacob Faibussowitsch break;
685d71ae5a4SJacob Faibussowitsch case 2:
686d71ae5a4SJacob Faibussowitsch PetscCall(Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume));
687d71ae5a4SJacob Faibussowitsch break;
688d71ae5a4SJacob Faibussowitsch case 3:
689d71ae5a4SJacob Faibussowitsch PetscCall(Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts, NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume));
690d71ae5a4SJacob Faibussowitsch break;
691d71ae5a4SJacob Faibussowitsch default:
692d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %" PetscInt_FMT, dim);
693181f196bSVijay Mahadevan }
6943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
695181f196bSVijay Mahadevan }
696181f196bSVijay Mahadevan
697cab5ea25SPierre Jolivet /*@C
69897b73a88SSatish Balay DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the
699a4e35b19SJacob Faibussowitsch canonical reference element.
700a86ed394SVijay Mahadevan
701d8d19677SJose E. Roman Input Parameters:
702a4e35b19SJacob Faibussowitsch + dim - the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET)
703a4e35b19SJacob Faibussowitsch . nverts - the number of vertices in the physical element
704a4e35b19SJacob Faibussowitsch . coordinates - the coordinates of vertices in the physical element
705a4e35b19SJacob Faibussowitsch - xphy - the coordinates of physical point for which natural coordinates (in reference frame) are sought
706a86ed394SVijay Mahadevan
707d8d19677SJose E. Roman Output Parameters:
708a4e35b19SJacob Faibussowitsch + natparam - the natural coordinates (in reference frame) corresponding to xphy
709a4e35b19SJacob Faibussowitsch - phi - the basis functions evaluated at the natural coordinates (natparam)
710a86ed394SVijay Mahadevan
711edc382c3SSatish Balay Level: advanced
712edc382c3SSatish Balay
713a4e35b19SJacob Faibussowitsch Notes:
714a4e35b19SJacob Faibussowitsch In addition to finding the inverse mapping evaluation through Newton iteration, the basis
715a4e35b19SJacob Faibussowitsch function at the parametric point is also evaluated optionally.
716a4e35b19SJacob Faibussowitsch
717a4e35b19SJacob Faibussowitsch .seealso: `DMMoabCreate()`
718a86ed394SVijay Mahadevan @*/
DMMoabPToRMapping(const PetscInt dim,const PetscInt nverts,const PetscReal * coordinates,const PetscReal * xphy,PetscReal * natparam,PetscReal * phi)719d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *xphy, PetscReal *natparam, PetscReal *phi)
720d71ae5a4SJacob Faibussowitsch {
721a86ed394SVijay Mahadevan /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */
722181f196bSVijay Mahadevan const PetscReal tol = 1.0e-10;
723181f196bSVijay Mahadevan const PetscInt max_iterations = 10;
724181f196bSVijay Mahadevan const PetscReal error_tol_sqr = tol * tol;
725181f196bSVijay Mahadevan PetscReal phibasis[8], jacobian[9], ijacobian[9], volume;
726181f196bSVijay Mahadevan PetscReal phypts[3] = {0.0, 0.0, 0.0};
727181f196bSVijay Mahadevan PetscReal delta[3] = {0.0, 0.0, 0.0};
728181f196bSVijay Mahadevan PetscInt iters = 0;
729181f196bSVijay Mahadevan PetscReal error = 1.0;
730181f196bSVijay Mahadevan
731181f196bSVijay Mahadevan PetscFunctionBegin;
7324f572ea9SToby Isaac PetscAssertPointer(coordinates, 3);
7334f572ea9SToby Isaac PetscAssertPointer(xphy, 4);
7344f572ea9SToby Isaac PetscAssertPointer(natparam, 5);
735181f196bSVijay Mahadevan
7369566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(jacobian, dim * dim));
7379566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ijacobian, dim * dim));
7389566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(phibasis, nverts));
739181f196bSVijay Mahadevan
740a86ed394SVijay Mahadevan /* zero initial guess */
741a86ed394SVijay Mahadevan natparam[0] = natparam[1] = natparam[2] = 0.0;
742181f196bSVijay Mahadevan
743a86ed394SVijay Mahadevan /* Compute delta = evaluate( xi) - x */
7449566063dSJacob Faibussowitsch PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
745a86ed394SVijay Mahadevan
746a86ed394SVijay Mahadevan error = 0.0;
747a86ed394SVijay Mahadevan switch (dim) {
748d71ae5a4SJacob Faibussowitsch case 3:
749d71ae5a4SJacob Faibussowitsch delta[2] = phypts[2] - xphy[2];
750d71ae5a4SJacob Faibussowitsch error += (delta[2] * delta[2]);
751d71ae5a4SJacob Faibussowitsch case 2:
752d71ae5a4SJacob Faibussowitsch delta[1] = phypts[1] - xphy[1];
753d71ae5a4SJacob Faibussowitsch error += (delta[1] * delta[1]);
754a86ed394SVijay Mahadevan case 1:
755a86ed394SVijay Mahadevan delta[0] = phypts[0] - xphy[0];
756a86ed394SVijay Mahadevan error += (delta[0] * delta[0]);
757a86ed394SVijay Mahadevan break;
758a86ed394SVijay Mahadevan }
759a86ed394SVijay Mahadevan
760181f196bSVijay Mahadevan while (error > error_tol_sqr) {
7617a46b595SBarry 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]);
762181f196bSVijay Mahadevan
763181f196bSVijay Mahadevan /* Compute natparam -= J.inverse() * delta */
764181f196bSVijay Mahadevan switch (dim) {
765d71ae5a4SJacob Faibussowitsch case 1:
766d71ae5a4SJacob Faibussowitsch natparam[0] -= ijacobian[0] * delta[0];
767d71ae5a4SJacob Faibussowitsch break;
768181f196bSVijay Mahadevan case 2:
769181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1];
770181f196bSVijay Mahadevan natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1];
771181f196bSVijay Mahadevan break;
772181f196bSVijay Mahadevan case 3:
773181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2];
774181f196bSVijay Mahadevan natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2];
775181f196bSVijay Mahadevan natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2];
776181f196bSVijay Mahadevan break;
777181f196bSVijay Mahadevan }
778181f196bSVijay Mahadevan
779181f196bSVijay Mahadevan /* Compute delta = evaluate( xi) - x */
7809566063dSJacob Faibussowitsch PetscCall(FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume));
781181f196bSVijay Mahadevan
782a86ed394SVijay Mahadevan error = 0.0;
783a86ed394SVijay Mahadevan switch (dim) {
784d71ae5a4SJacob Faibussowitsch case 3:
785d71ae5a4SJacob Faibussowitsch delta[2] = phypts[2] - xphy[2];
786d71ae5a4SJacob Faibussowitsch error += (delta[2] * delta[2]);
787d71ae5a4SJacob Faibussowitsch case 2:
788d71ae5a4SJacob Faibussowitsch delta[1] = phypts[1] - xphy[1];
789d71ae5a4SJacob Faibussowitsch error += (delta[1] * delta[1]);
790a86ed394SVijay Mahadevan case 1:
791a86ed394SVijay Mahadevan delta[0] = phypts[0] - xphy[0];
792a86ed394SVijay Mahadevan error += (delta[0] * delta[0]);
793a86ed394SVijay Mahadevan break;
794a86ed394SVijay Mahadevan }
795181f196bSVijay Mahadevan }
79648a46eb9SPierre Jolivet if (phi) PetscCall(PetscArraycpy(phi, phibasis, nverts));
7973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
798181f196bSVijay Mahadevan }
799