163d025dbSVijay Mahadevan 263d025dbSVijay Mahadevan #include <petscconf.h> 363d025dbSVijay Mahadevan #include <petscdt.h> /*I "petscdt.h" I*/ 463d025dbSVijay Mahadevan #include <petsc/private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 563d025dbSVijay Mahadevan 663d025dbSVijay Mahadevan /* Utility functions */ 763d025dbSVijay Mahadevan static inline PetscReal DMatrix_Determinant_2x2_Internal (const PetscReal inmat[2 * 2]) 863d025dbSVijay Mahadevan { 963d025dbSVijay Mahadevan return inmat[0] * inmat[3] - inmat[1] * inmat[2]; 1063d025dbSVijay Mahadevan } 1163d025dbSVijay Mahadevan 1263d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_2x2_Internal(const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant) 1363d025dbSVijay Mahadevan { 1463d025dbSVijay Mahadevan PetscReal det = DMatrix_Determinant_2x2_Internal(inmat); 1563d025dbSVijay Mahadevan if (outmat) { 1663d025dbSVijay Mahadevan outmat[0] = inmat[3] / det; 1763d025dbSVijay Mahadevan outmat[1] = -inmat[1] / det; 1863d025dbSVijay Mahadevan outmat[2] = -inmat[2] / det; 1963d025dbSVijay Mahadevan outmat[3] = inmat[0] / det; 2063d025dbSVijay Mahadevan } 2163d025dbSVijay Mahadevan if (determinant) *determinant = det; 22362febeeSStefano Zampini return 0; 2363d025dbSVijay Mahadevan } 2463d025dbSVijay Mahadevan 25181f196bSVijay Mahadevan static inline PetscReal DMatrix_Determinant_3x3_Internal(const PetscReal inmat[3 * 3]) 2663d025dbSVijay Mahadevan { 2763d025dbSVijay Mahadevan return inmat[0] * (inmat[8] * inmat[4] - inmat[7] * inmat[5]) 2863d025dbSVijay Mahadevan - inmat[3] * (inmat[8] * inmat[1] - inmat[7] * inmat[2]) 2963d025dbSVijay Mahadevan + inmat[6] * (inmat[5] * inmat[1] - inmat[4] * inmat[2]); 3063d025dbSVijay Mahadevan } 3163d025dbSVijay Mahadevan 3263d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_3x3_Internal (const PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant) 3363d025dbSVijay Mahadevan { 34181f196bSVijay Mahadevan PetscReal det = DMatrix_Determinant_3x3_Internal(inmat); 3563d025dbSVijay Mahadevan if (outmat) { 3663d025dbSVijay Mahadevan outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det; 3763d025dbSVijay Mahadevan outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det; 3863d025dbSVijay Mahadevan outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det; 3963d025dbSVijay Mahadevan outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det; 4063d025dbSVijay Mahadevan outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det; 4163d025dbSVijay Mahadevan outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det; 4263d025dbSVijay Mahadevan outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det; 4363d025dbSVijay Mahadevan outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det; 4463d025dbSVijay Mahadevan outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det; 4563d025dbSVijay Mahadevan } 4663d025dbSVijay Mahadevan if (determinant) *determinant = det; 47362febeeSStefano Zampini return 0; 4863d025dbSVijay Mahadevan } 4963d025dbSVijay Mahadevan 50181f196bSVijay Mahadevan inline PetscReal DMatrix_Determinant_4x4_Internal(PetscReal inmat[4 * 4]) 51181f196bSVijay Mahadevan { 52181f196bSVijay Mahadevan return 53181f196bSVijay Mahadevan inmat[0 + 0 * 4] * ( 54181f196bSVijay Mahadevan inmat[1 + 1 * 4] * (inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4]) 55181f196bSVijay Mahadevan - inmat[1 + 2 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4]) 56181f196bSVijay Mahadevan + inmat[1 + 3 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4])) 57181f196bSVijay Mahadevan - inmat[0 + 1 * 4] * ( 58181f196bSVijay Mahadevan inmat[1 + 0 * 4] * (inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4]) 59181f196bSVijay Mahadevan - inmat[1 + 2 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4]) 60181f196bSVijay Mahadevan + inmat[1 + 3 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4])) 61181f196bSVijay Mahadevan + inmat[0 + 2 * 4] * ( 62181f196bSVijay Mahadevan inmat[1 + 0 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4]) 63181f196bSVijay Mahadevan - inmat[1 + 1 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4]) 64181f196bSVijay Mahadevan + inmat[1 + 3 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4])) 65181f196bSVijay Mahadevan - inmat[0 + 3 * 4] * ( 66181f196bSVijay Mahadevan inmat[1 + 0 * 4] * (inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4]) 67181f196bSVijay Mahadevan - inmat[1 + 1 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4]) 68181f196bSVijay Mahadevan + inmat[1 + 2 * 4] * (inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4])); 69181f196bSVijay Mahadevan } 70181f196bSVijay Mahadevan 71181f196bSVijay Mahadevan inline PetscErrorCode DMatrix_Invert_4x4_Internal (PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant) 72181f196bSVijay Mahadevan { 73181f196bSVijay Mahadevan PetscReal det = DMatrix_Determinant_4x4_Internal(inmat); 74181f196bSVijay Mahadevan if (outmat) { 75181f196bSVijay Mahadevan outmat[0] = (inmat[5] * inmat[10] * inmat[15] + inmat[6] * inmat[11] * inmat[13] + inmat[7] * inmat[9] * inmat[14] - inmat[5] * inmat[11] * inmat[14] - inmat[6] * inmat[9] * inmat[15] - inmat[7] * inmat[10] * inmat[13]) / det; 76181f196bSVijay Mahadevan outmat[1] = (inmat[1] * inmat[11] * inmat[14] + inmat[2] * inmat[9] * inmat[15] + inmat[3] * inmat[10] * inmat[13] - inmat[1] * inmat[10] * inmat[15] - inmat[2] * inmat[11] * inmat[13] - inmat[3] * inmat[9] * inmat[14]) / det; 77181f196bSVijay Mahadevan outmat[2] = (inmat[1] * inmat[6] * inmat[15] + inmat[2] * inmat[7] * inmat[13] + inmat[3] * inmat[5] * inmat[14] - inmat[1] * inmat[7] * inmat[14] - inmat[2] * inmat[5] * inmat[15] - inmat[3] * inmat[6] * inmat[13]) / det; 78181f196bSVijay Mahadevan outmat[3] = (inmat[1] * inmat[7] * inmat[10] + inmat[2] * inmat[5] * inmat[11] + inmat[3] * inmat[6] * inmat[9] - inmat[1] * inmat[6] * inmat[11] - inmat[2] * inmat[7] * inmat[9] - inmat[3] * inmat[5] * inmat[10]) / det; 79181f196bSVijay Mahadevan outmat[4] = (inmat[4] * inmat[11] * inmat[14] + inmat[6] * inmat[8] * inmat[15] + inmat[7] * inmat[10] * inmat[12] - inmat[4] * inmat[10] * inmat[15] - inmat[6] * inmat[11] * inmat[12] - inmat[7] * inmat[8] * inmat[14]) / det; 80181f196bSVijay Mahadevan outmat[5] = (inmat[0] * inmat[10] * inmat[15] + inmat[2] * inmat[11] * inmat[12] + inmat[3] * inmat[8] * inmat[14] - inmat[0] * inmat[11] * inmat[14] - inmat[2] * inmat[8] * inmat[15] - inmat[3] * inmat[10] * inmat[12]) / det; 81181f196bSVijay Mahadevan outmat[6] = (inmat[0] * inmat[7] * inmat[14] + inmat[2] * inmat[4] * inmat[15] + inmat[3] * inmat[6] * inmat[12] - inmat[0] * inmat[6] * inmat[15] - inmat[2] * inmat[7] * inmat[12] - inmat[3] * inmat[4] * inmat[14]) / det; 82181f196bSVijay Mahadevan outmat[7] = (inmat[0] * inmat[6] * inmat[11] + inmat[2] * inmat[7] * inmat[8] + inmat[3] * inmat[4] * inmat[10] - inmat[0] * inmat[7] * inmat[10] - inmat[2] * inmat[4] * inmat[11] - inmat[3] * inmat[6] * inmat[8]) / det; 83181f196bSVijay Mahadevan outmat[8] = (inmat[4] * inmat[9] * inmat[15] + inmat[5] * inmat[11] * inmat[12] + inmat[7] * inmat[8] * inmat[13] - inmat[4] * inmat[11] * inmat[13] - inmat[5] * inmat[8] * inmat[15] - inmat[7] * inmat[9] * inmat[12]) / det; 84181f196bSVijay Mahadevan outmat[9] = (inmat[0] * inmat[11] * inmat[13] + inmat[1] * inmat[8] * inmat[15] + inmat[3] * inmat[9] * inmat[12] - inmat[0] * inmat[9] * inmat[15] - inmat[1] * inmat[11] * inmat[12] - inmat[3] * inmat[8] * inmat[13]) / det; 85181f196bSVijay Mahadevan outmat[10] = (inmat[0] * inmat[5] * inmat[15] + inmat[1] * inmat[7] * inmat[12] + inmat[3] * inmat[4] * inmat[13] - inmat[0] * inmat[7] * inmat[13] - inmat[1] * inmat[4] * inmat[15] - inmat[3] * inmat[5] * inmat[12]) / det; 86181f196bSVijay Mahadevan outmat[11] = (inmat[0] * inmat[7] * inmat[9] + inmat[1] * inmat[4] * inmat[11] + inmat[3] * inmat[5] * inmat[8] - inmat[0] * inmat[5] * inmat[11] - inmat[1] * inmat[7] * inmat[8] - inmat[3] * inmat[4] * inmat[9]) / det; 87181f196bSVijay Mahadevan outmat[12] = (inmat[4] * inmat[10] * inmat[13] + inmat[5] * inmat[8] * inmat[14] + inmat[6] * inmat[9] * inmat[12] - inmat[4] * inmat[9] * inmat[14] - inmat[5] * inmat[10] * inmat[12] - inmat[6] * inmat[8] * inmat[13]) / det; 88181f196bSVijay Mahadevan outmat[13] = (inmat[0] * inmat[9] * inmat[14] + inmat[1] * inmat[10] * inmat[12] + inmat[2] * inmat[8] * inmat[13] - inmat[0] * inmat[10] * inmat[13] - inmat[1] * inmat[8] * inmat[14] - inmat[2] * inmat[9] * inmat[12]) / det; 89181f196bSVijay Mahadevan outmat[14] = (inmat[0] * inmat[6] * inmat[13] + inmat[1] * inmat[4] * inmat[14] + inmat[2] * inmat[5] * inmat[12] - inmat[0] * inmat[5] * inmat[14] - inmat[1] * inmat[6] * inmat[12] - inmat[2] * inmat[4] * inmat[13]) / det; 90181f196bSVijay Mahadevan outmat[15] = (inmat[0] * inmat[5] * inmat[10] + inmat[1] * inmat[6] * inmat[8] + inmat[2] * inmat[4] * inmat[9] - inmat[0] * inmat[6] * inmat[9] - inmat[1] * inmat[4] * inmat[10] - inmat[2] * inmat[5] * inmat[8]) / det; 91181f196bSVijay Mahadevan } 92181f196bSVijay Mahadevan if (determinant) *determinant = det; 93362febeeSStefano Zampini return 0; 94181f196bSVijay Mahadevan } 95181f196bSVijay Mahadevan 96cab5ea25SPierre Jolivet /*@C 9797b73a88SSatish Balay Compute_Lagrange_Basis_1D_Internal - Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element. 9863d025dbSVijay Mahadevan 9963d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a linear or quadratic edge element. 10063d025dbSVijay Mahadevan 10163d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided, 10263d025dbSVijay Mahadevan and their derivatives with respect to X. 10363d025dbSVijay Mahadevan 10463d025dbSVijay Mahadevan Notes: 10563d025dbSVijay Mahadevan 10663d025dbSVijay Mahadevan Example Physical Element 107a2b725a8SWilliam Gropp .vb 10863d025dbSVijay Mahadevan 1-------2 1----3----2 10963d025dbSVijay Mahadevan EDGE2 EDGE3 110a2b725a8SWilliam Gropp .ve 11163d025dbSVijay Mahadevan 112d8d19677SJose E. Roman Input Parameters: 113a2b725a8SWilliam Gropp + PetscInt nverts - the number of element vertices 114a2b725a8SWilliam Gropp . PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering) 115a2b725a8SWilliam Gropp . PetscInt npts - the number of evaluation points (quadrature points) 116a2b725a8SWilliam Gropp - PetscReal quad[3*npts] - the evaluation points (quadrature points) in the reference space 11763d025dbSVijay Mahadevan 118d8d19677SJose E. Roman Output Parameters: 119a2b725a8SWilliam Gropp + PetscReal phypts[3*npts] - the evaluation points (quadrature points) transformed to the physical space 120a2b725a8SWilliam Gropp . PetscReal jxw[npts] - the jacobian determinant * quadrature weight necessary for assembling discrete contributions 121a2b725a8SWilliam Gropp . PetscReal phi[npts] - the bases evaluated at the specified quadrature points 122*6b867d5aSJose E. Roman . PetscReal dphidx[npts] - the derivative of the bases wrt X-direction evaluated at the specified quadrature points 123*6b867d5aSJose E. Roman . jacobian - 124*6b867d5aSJose E. Roman . ijacobian - 125*6b867d5aSJose E. Roman - volume 12663d025dbSVijay Mahadevan 127edc382c3SSatish Balay Level: advanced 128edc382c3SSatish Balay 12963d025dbSVijay Mahadevan @*/ 13063d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_1D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, 131181f196bSVijay Mahadevan PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, 132181f196bSVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) 13363d025dbSVijay Mahadevan { 13463d025dbSVijay Mahadevan int i, j; 13563d025dbSVijay Mahadevan PetscErrorCode ierr; 13663d025dbSVijay Mahadevan 137181f196bSVijay Mahadevan PetscFunctionBegin; 138181f196bSVijay Mahadevan PetscValidPointer(jacobian, 9); 139181f196bSVijay Mahadevan PetscValidPointer(ijacobian, 10); 140181f196bSVijay Mahadevan PetscValidPointer(volume, 11); 141181f196bSVijay Mahadevan if (phypts) { 142580bdb30SBarry Smith ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr); 143181f196bSVijay Mahadevan } 14463d025dbSVijay Mahadevan if (dphidx) { /* Reset arrays. */ 145580bdb30SBarry Smith ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr); 14663d025dbSVijay Mahadevan } 14763d025dbSVijay Mahadevan if (nverts == 2) { /* Linear Edge */ 14863d025dbSVijay Mahadevan 1492da392ccSBarry Smith for (j = 0; j < npts; j++) { 150a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 151181f196bSVijay Mahadevan const PetscReal r = quad[j]; 15263d025dbSVijay Mahadevan 153181f196bSVijay Mahadevan phi[0 + offset] = ( 1.0 - r); 154181f196bSVijay Mahadevan phi[1 + offset] = ( r); 15563d025dbSVijay Mahadevan 156181f196bSVijay Mahadevan const PetscReal dNi_dxi[2] = { -1.0, 1.0 }; 15763d025dbSVijay Mahadevan 158181f196bSVijay Mahadevan jacobian[0] = ijacobian[0] = volume[0] = 0.0; 15963d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 160181f196bSVijay Mahadevan const PetscReal* vertices = coords + i * 3; 161181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 162181f196bSVijay Mahadevan if (phypts) { 163181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 164181f196bSVijay Mahadevan } 16563d025dbSVijay Mahadevan } 16663d025dbSVijay Mahadevan 16763d025dbSVijay Mahadevan /* invert the jacobian */ 168181f196bSVijay Mahadevan *volume = jacobian[0]; 169181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0]; 170181f196bSVijay Mahadevan jxw[j] *= *volume; 17163d025dbSVijay Mahadevan 17263d025dbSVijay Mahadevan /* Divide by element jacobian. */ 17363d025dbSVijay Mahadevan for (i = 0; i < nverts; i++) { 174181f196bSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0]; 17563d025dbSVijay Mahadevan } 17663d025dbSVijay Mahadevan } 1772da392ccSBarry Smith } else if (nverts == 3) { /* Quadratic Edge */ 17863d025dbSVijay Mahadevan 1792da392ccSBarry Smith for (j = 0; j < npts; j++) { 180a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 181181f196bSVijay Mahadevan const PetscReal r = quad[j]; 18263d025dbSVijay Mahadevan 183181f196bSVijay Mahadevan phi[0 + offset] = 1.0 + r * ( 2.0 * r - 3.0); 184181f196bSVijay Mahadevan phi[1 + offset] = 4.0 * r * ( 1.0 - r); 185181f196bSVijay Mahadevan phi[2 + offset] = r * ( 2.0 * r - 1.0); 18663d025dbSVijay Mahadevan 187181f196bSVijay Mahadevan const PetscReal dNi_dxi[3] = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r), 4.0 * r - 1.0}; 18863d025dbSVijay Mahadevan 189181f196bSVijay Mahadevan jacobian[0] = ijacobian[0] = volume[0] = 0.0; 19063d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 191181f196bSVijay Mahadevan const PetscReal* vertices = coords + i * 3; 192181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 193181f196bSVijay Mahadevan if (phypts) { 194181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 195181f196bSVijay Mahadevan } 19663d025dbSVijay Mahadevan } 19763d025dbSVijay Mahadevan 19863d025dbSVijay Mahadevan /* invert the jacobian */ 199181f196bSVijay Mahadevan *volume = jacobian[0]; 200181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0]; 201181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume; 20263d025dbSVijay Mahadevan 20363d025dbSVijay Mahadevan /* Divide by element jacobian. */ 20463d025dbSVijay Mahadevan for (i = 0; i < nverts; i++) { 205181f196bSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0]; 20663d025dbSVijay Mahadevan } 20763d025dbSVijay Mahadevan } 2082da392ccSBarry Smith } else 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); 20963d025dbSVijay Mahadevan PetscFunctionReturn(0); 21063d025dbSVijay Mahadevan } 21163d025dbSVijay Mahadevan 212cab5ea25SPierre Jolivet /*@C 21397b73a88SSatish Balay Compute_Lagrange_Basis_2D_Internal - Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element. 21463d025dbSVijay Mahadevan 21563d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a quadrangle or triangle. 21663d025dbSVijay Mahadevan 21763d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided, 21863d025dbSVijay Mahadevan and their derivatives with respect to X and Y. 21963d025dbSVijay Mahadevan 22063d025dbSVijay Mahadevan Notes: 22163d025dbSVijay Mahadevan 22263d025dbSVijay Mahadevan Example Physical Element (QUAD4) 223a2b725a8SWilliam Gropp .vb 22463d025dbSVijay Mahadevan 4------3 s 22563d025dbSVijay Mahadevan | | | 22663d025dbSVijay Mahadevan | | | 22763d025dbSVijay Mahadevan | | | 22863d025dbSVijay Mahadevan 1------2 0-------r 229a2b725a8SWilliam Gropp .ve 23063d025dbSVijay Mahadevan 231d8d19677SJose E. Roman Input Parameters: 232a2b725a8SWilliam Gropp + PetscInt nverts - the number of element vertices 233a2b725a8SWilliam Gropp . PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering) 234a2b725a8SWilliam Gropp . PetscInt npts - the number of evaluation points (quadrature points) 235a2b725a8SWilliam Gropp - PetscReal quad[3*npts] - the evaluation points (quadrature points) in the reference space 23663d025dbSVijay Mahadevan 237d8d19677SJose E. Roman Output Parameters: 238a2b725a8SWilliam Gropp + PetscReal phypts[3*npts] - the evaluation points (quadrature points) transformed to the physical space 239a2b725a8SWilliam Gropp . PetscReal jxw[npts] - the jacobian determinant * quadrature weight necessary for assembling discrete contributions 240a2b725a8SWilliam Gropp . PetscReal phi[npts] - the bases evaluated at the specified quadrature points 241a2b725a8SWilliam Gropp . PetscReal dphidx[npts] - the derivative of the bases wrt X-direction evaluated at the specified quadrature points 242*6b867d5aSJose E. Roman . PetscReal dphidy[npts] - the derivative of the bases wrt Y-direction evaluated at the specified quadrature points 243*6b867d5aSJose E. Roman . jacobian - 244*6b867d5aSJose E. Roman . ijacobian - 245*6b867d5aSJose E. Roman - volume 24663d025dbSVijay Mahadevan 247edc382c3SSatish Balay Level: advanced 248edc382c3SSatish Balay 24963d025dbSVijay Mahadevan @*/ 25063d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_2D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, 251181f196bSVijay Mahadevan PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, 252181f196bSVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) 25363d025dbSVijay Mahadevan { 254a86ed394SVijay Mahadevan PetscInt i, j, k; 25563d025dbSVijay Mahadevan PetscErrorCode ierr; 25663d025dbSVijay Mahadevan 25763d025dbSVijay Mahadevan PetscFunctionBegin; 258181f196bSVijay Mahadevan PetscValidPointer(jacobian, 10); 259181f196bSVijay Mahadevan PetscValidPointer(ijacobian, 11); 260181f196bSVijay Mahadevan PetscValidPointer(volume, 12); 261580bdb30SBarry Smith ierr = PetscArrayzero(phi, npts);CHKERRQ(ierr); 262181f196bSVijay Mahadevan if (phypts) { 263580bdb30SBarry Smith ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr); 264181f196bSVijay Mahadevan } 26563d025dbSVijay Mahadevan if (dphidx) { /* Reset arrays. */ 266580bdb30SBarry Smith ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr); 267580bdb30SBarry Smith ierr = PetscArrayzero(dphidy, npts * nverts);CHKERRQ(ierr); 26863d025dbSVijay Mahadevan } 26963d025dbSVijay Mahadevan if (nverts == 4) { /* Linear Quadrangle */ 27063d025dbSVijay Mahadevan 27163d025dbSVijay Mahadevan for (j = 0; j < npts; j++) 27263d025dbSVijay Mahadevan { 273a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 274181f196bSVijay Mahadevan const PetscReal r = quad[0 + j * 2]; 275181f196bSVijay Mahadevan const PetscReal s = quad[1 + j * 2]; 27663d025dbSVijay Mahadevan 27763d025dbSVijay Mahadevan phi[0 + offset] = ( 1.0 - r) * ( 1.0 - s); 27863d025dbSVijay Mahadevan phi[1 + offset] = r * ( 1.0 - s); 27963d025dbSVijay Mahadevan phi[2 + offset] = r * s; 28063d025dbSVijay Mahadevan phi[3 + offset] = ( 1.0 - r) * s; 28163d025dbSVijay Mahadevan 282181f196bSVijay Mahadevan const PetscReal dNi_dxi[4] = { -1.0 + s, 1.0 - s, s, -s }; 283181f196bSVijay Mahadevan const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r }; 28463d025dbSVijay Mahadevan 285580bdb30SBarry Smith ierr = PetscArrayzero(jacobian, 4);CHKERRQ(ierr); 286580bdb30SBarry Smith ierr = PetscArrayzero(ijacobian, 4);CHKERRQ(ierr); 28763d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 288181f196bSVijay Mahadevan const PetscReal* vertices = coords + i * 3; 28963d025dbSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 29063d025dbSVijay Mahadevan jacobian[2] += dNi_dxi[i] * vertices[1]; 29163d025dbSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertices[0]; 29263d025dbSVijay Mahadevan jacobian[3] += dNi_deta[i] * vertices[1]; 293181f196bSVijay Mahadevan if (phypts) { 294181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 295181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertices[1]; 296181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertices[2]; 297181f196bSVijay Mahadevan } 29863d025dbSVijay Mahadevan } 29963d025dbSVijay Mahadevan 30063d025dbSVijay Mahadevan /* invert the jacobian */ 301181f196bSVijay Mahadevan ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 302181f196bSVijay 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); 30363d025dbSVijay Mahadevan 304181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume; 30563d025dbSVijay Mahadevan 306181f196bSVijay Mahadevan /* Let us compute the bases derivatives by scaling with inverse jacobian. */ 307181f196bSVijay Mahadevan if (dphidx) { 30863d025dbSVijay Mahadevan for (i = 0; i < nverts; i++) { 309a86ed394SVijay Mahadevan for (k = 0; k < 2; ++k) { 31063d025dbSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0]; 31163d025dbSVijay Mahadevan if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1]; 312181f196bSVijay Mahadevan } /* for k=[0..2) */ 313181f196bSVijay Mahadevan } /* for i=[0..nverts) */ 314181f196bSVijay Mahadevan } /* if (dphidx) */ 31563d025dbSVijay Mahadevan } 3162da392ccSBarry Smith } else if (nverts == 3) { /* Linear triangle */ 3172da392ccSBarry Smith const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1]; 31863d025dbSVijay Mahadevan 319580bdb30SBarry Smith ierr = PetscArrayzero(jacobian, 4);CHKERRQ(ierr); 320580bdb30SBarry Smith ierr = PetscArrayzero(ijacobian, 4);CHKERRQ(ierr); 32163d025dbSVijay Mahadevan 32263d025dbSVijay Mahadevan /* Jacobian is constant */ 323181f196bSVijay Mahadevan jacobian[0] = (coords[0 * 3 + 0] - x2); jacobian[1] = (coords[1 * 3 + 0] - x2); 324181f196bSVijay Mahadevan jacobian[2] = (coords[0 * 3 + 1] - y2); jacobian[3] = (coords[1 * 3 + 1] - y2); 32563d025dbSVijay Mahadevan 32663d025dbSVijay Mahadevan /* invert the jacobian */ 327181f196bSVijay Mahadevan ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 3282da392ccSBarry Smith if (*volume < 1e-12) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangular element has zero volume: %g. Degenerate element or invalid connectivity\n", (double)*volume); 329181f196bSVijay Mahadevan 330181f196bSVijay Mahadevan const PetscReal Dx[3] = { ijacobian[0], ijacobian[2], - ijacobian[0] - ijacobian[2] }; 331181f196bSVijay Mahadevan const PetscReal Dy[3] = { ijacobian[1], ijacobian[3], - ijacobian[1] - ijacobian[3] }; 33263d025dbSVijay Mahadevan 3332da392ccSBarry Smith for (j = 0; j < npts; j++) { 334a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 335181f196bSVijay Mahadevan const PetscReal r = quad[0 + j * 2]; 336181f196bSVijay Mahadevan const PetscReal s = quad[1 + j * 2]; 33763d025dbSVijay Mahadevan 338181f196bSVijay Mahadevan if (jxw) jxw[j] *= 0.5; 33963d025dbSVijay Mahadevan 340181f196bSVijay Mahadevan /* const PetscReal Ni[3] = { r, s, 1.0 - r - s }; */ 341181f196bSVijay Mahadevan const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s; 342181f196bSVijay Mahadevan const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s; 343181f196bSVijay Mahadevan if (phypts) { 344181f196bSVijay Mahadevan phypts[offset + 0] = phipts_x; 345181f196bSVijay Mahadevan phypts[offset + 1] = phipts_y; 346181f196bSVijay Mahadevan } 34763d025dbSVijay Mahadevan 348110fc3b0SBarry Smith /* \phi_0 = (b.y - c.y) x + (b.x - c.x) y + c.x b.y - b.x c.y */ 349181f196bSVijay Mahadevan phi[0 + offset] = ( ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2)); 350110fc3b0SBarry Smith /* \phi_1 = (c.y - a.y) x + (a.x - c.x) y + c.x a.y - a.x c.y */ 351181f196bSVijay Mahadevan phi[1 + offset] = ( ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2)); 35263d025dbSVijay Mahadevan phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset]; 35363d025dbSVijay Mahadevan 35463d025dbSVijay Mahadevan if (dphidx) { 355181f196bSVijay Mahadevan dphidx[0 + offset] = Dx[0]; 356181f196bSVijay Mahadevan dphidx[1 + offset] = Dx[1]; 357181f196bSVijay Mahadevan dphidx[2 + offset] = Dx[2]; 35863d025dbSVijay Mahadevan } 35963d025dbSVijay Mahadevan 36063d025dbSVijay Mahadevan if (dphidy) { 361181f196bSVijay Mahadevan dphidy[0 + offset] = Dy[0]; 362181f196bSVijay Mahadevan dphidy[1 + offset] = Dy[1]; 363181f196bSVijay Mahadevan dphidy[2 + offset] = Dy[2]; 36463d025dbSVijay Mahadevan } 36563d025dbSVijay Mahadevan 36663d025dbSVijay Mahadevan } 3672da392ccSBarry Smith } else 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); 36863d025dbSVijay Mahadevan PetscFunctionReturn(0); 36963d025dbSVijay Mahadevan } 37063d025dbSVijay Mahadevan 371cab5ea25SPierre Jolivet /*@C 37297b73a88SSatish Balay Compute_Lagrange_Basis_3D_Internal - Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element. 37363d025dbSVijay Mahadevan 37463d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a hexahedra or tetrahedra. 37563d025dbSVijay Mahadevan 37663d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided, 37763d025dbSVijay Mahadevan and their derivatives with respect to X, Y and Z. 37863d025dbSVijay Mahadevan 37963d025dbSVijay Mahadevan Notes: 38063d025dbSVijay Mahadevan 38163d025dbSVijay Mahadevan Example Physical Element (HEX8) 382a2b725a8SWilliam Gropp .vb 38363d025dbSVijay Mahadevan 8------7 38463d025dbSVijay Mahadevan /| /| t s 38563d025dbSVijay Mahadevan 5------6 | | / 38663d025dbSVijay Mahadevan | | | | |/ 38763d025dbSVijay Mahadevan | 4----|-3 0-------r 38863d025dbSVijay Mahadevan |/ |/ 38963d025dbSVijay Mahadevan 1------2 390a2b725a8SWilliam Gropp .ve 39163d025dbSVijay Mahadevan 392d8d19677SJose E. Roman Input Parameters: 393a2b725a8SWilliam Gropp + PetscInt nverts - the number of element vertices 394a2b725a8SWilliam Gropp . PetscReal coords[3*nverts] - the physical coordinates of the vertices (in canonical numbering) 395a2b725a8SWilliam Gropp . PetscInt npts - the number of evaluation points (quadrature points) 396a2b725a8SWilliam Gropp - PetscReal quad[3*npts] - the evaluation points (quadrature points) in the reference space 39763d025dbSVijay Mahadevan 398d8d19677SJose E. Roman Output Parameters: 399a2b725a8SWilliam Gropp + PetscReal phypts[3*npts] - the evaluation points (quadrature points) transformed to the physical space 400a2b725a8SWilliam Gropp . PetscReal jxw[npts] - the jacobian determinant * quadrature weight necessary for assembling discrete contributions 401a2b725a8SWilliam Gropp . PetscReal phi[npts] - the bases evaluated at the specified quadrature points 402a2b725a8SWilliam Gropp . PetscReal dphidx[npts] - the derivative of the bases wrt X-direction evaluated at the specified quadrature points 403a2b725a8SWilliam Gropp . PetscReal dphidy[npts] - the derivative of the bases wrt Y-direction evaluated at the specified quadrature points 404*6b867d5aSJose E. Roman . PetscReal dphidz[npts] - the derivative of the bases wrt Z-direction evaluated at the specified quadrature points 405*6b867d5aSJose E. Roman . jacobian - 406*6b867d5aSJose E. Roman . ijacobian - 407*6b867d5aSJose E. Roman - volume 40863d025dbSVijay Mahadevan 409edc382c3SSatish Balay Level: advanced 410edc382c3SSatish Balay 41163d025dbSVijay Mahadevan @*/ 41263d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_3D_Internal(const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, 413181f196bSVijay Mahadevan PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz, 414181f196bSVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) 41563d025dbSVijay Mahadevan { 416a86ed394SVijay Mahadevan PetscInt i, j, k; 41763d025dbSVijay Mahadevan PetscErrorCode ierr; 41863d025dbSVijay Mahadevan 41963d025dbSVijay Mahadevan PetscFunctionBegin; 420181f196bSVijay Mahadevan PetscValidPointer(jacobian, 11); 421181f196bSVijay Mahadevan PetscValidPointer(ijacobian, 12); 422181f196bSVijay Mahadevan PetscValidPointer(volume, 13); 4232da392ccSBarry Smith 424580bdb30SBarry Smith ierr = PetscArrayzero(phi, npts);CHKERRQ(ierr); 425181f196bSVijay Mahadevan if (phypts) { 426580bdb30SBarry Smith ierr = PetscArrayzero(phypts, npts * 3);CHKERRQ(ierr); 427181f196bSVijay Mahadevan } 42863d025dbSVijay Mahadevan if (dphidx) { 429580bdb30SBarry Smith ierr = PetscArrayzero(dphidx, npts * nverts);CHKERRQ(ierr); 430580bdb30SBarry Smith ierr = PetscArrayzero(dphidy, npts * nverts);CHKERRQ(ierr); 431580bdb30SBarry Smith ierr = PetscArrayzero(dphidz, npts * nverts);CHKERRQ(ierr); 43263d025dbSVijay Mahadevan } 43363d025dbSVijay Mahadevan 43463d025dbSVijay Mahadevan if (nverts == 8) { /* Linear Hexahedra */ 43563d025dbSVijay Mahadevan 4362da392ccSBarry Smith for (j = 0; j < npts; j++) { 437a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 438181f196bSVijay Mahadevan const PetscReal& r = quad[j * 3 + 0]; 439181f196bSVijay Mahadevan const PetscReal& s = quad[j * 3 + 1]; 440181f196bSVijay Mahadevan const PetscReal& t = quad[j * 3 + 2]; 44163d025dbSVijay Mahadevan 442a86ed394SVijay Mahadevan phi[offset + 0] = ( 1.0 - r) * ( 1.0 - s) * ( 1.0 - t); /* 0,0,0 */ 443a86ed394SVijay Mahadevan phi[offset + 1] = ( r) * ( 1.0 - s) * ( 1.0 - t); /* 1,0,0 */ 444a86ed394SVijay Mahadevan phi[offset + 2] = ( r) * ( s) * ( 1.0 - t); /* 1,1,0 */ 445a86ed394SVijay Mahadevan phi[offset + 3] = ( 1.0 - r) * ( s) * ( 1.0 - t); /* 0,1,0 */ 446a86ed394SVijay Mahadevan phi[offset + 4] = ( 1.0 - r) * ( 1.0 - s) * ( t); /* 0,0,1 */ 447a86ed394SVijay Mahadevan phi[offset + 5] = ( r) * ( 1.0 - s) * ( t); /* 1,0,1 */ 448a86ed394SVijay Mahadevan phi[offset + 6] = ( r) * ( s) * ( t); /* 1,1,1 */ 449a86ed394SVijay Mahadevan phi[offset + 7] = ( 1.0 - r) * ( s) * ( t); /* 0,1,1 */ 45063d025dbSVijay Mahadevan 451181f196bSVijay Mahadevan const PetscReal dNi_dxi[8] = {- ( 1.0 - s) * ( 1.0 - t), 45263d025dbSVijay Mahadevan ( 1.0 - s) * ( 1.0 - t), 453a86ed394SVijay Mahadevan ( s) * ( 1.0 - t), 454a86ed394SVijay Mahadevan - ( s) * ( 1.0 - t), 455a86ed394SVijay Mahadevan - ( 1.0 - s) * ( t), 456a86ed394SVijay Mahadevan ( 1.0 - s) * ( t), 457a86ed394SVijay Mahadevan ( s) * ( t), 458a86ed394SVijay Mahadevan - ( s) * ( t) 45963d025dbSVijay Mahadevan }; 46063d025dbSVijay Mahadevan 461181f196bSVijay Mahadevan const PetscReal dNi_deta[8] = { - ( 1.0 - r) * ( 1.0 - t), 462a86ed394SVijay Mahadevan - ( r) * ( 1.0 - t), 463a86ed394SVijay Mahadevan ( r) * ( 1.0 - t), 46463d025dbSVijay Mahadevan ( 1.0 - r) * ( 1.0 - t), 465a86ed394SVijay Mahadevan - ( 1.0 - r) * ( t), 466a86ed394SVijay Mahadevan - ( r) * ( t), 467a86ed394SVijay Mahadevan ( r) * ( t), 468a86ed394SVijay Mahadevan ( 1.0 - r) * ( t) 46963d025dbSVijay Mahadevan }; 47063d025dbSVijay Mahadevan 471181f196bSVijay Mahadevan const PetscReal dNi_dzeta[8] = { - ( 1.0 - r) * ( 1.0 - s), 472a86ed394SVijay Mahadevan - ( r) * ( 1.0 - s), 473a86ed394SVijay Mahadevan - ( r) * ( s), 474a86ed394SVijay Mahadevan - ( 1.0 - r) * ( s), 47563d025dbSVijay Mahadevan ( 1.0 - r) * ( 1.0 - s), 476a86ed394SVijay Mahadevan ( r) * ( 1.0 - s), 477a86ed394SVijay Mahadevan ( r) * ( s), 478a86ed394SVijay Mahadevan ( 1.0 - r) * ( s) 47963d025dbSVijay Mahadevan }; 48063d025dbSVijay Mahadevan 481580bdb30SBarry Smith ierr = PetscArrayzero(jacobian, 9);CHKERRQ(ierr); 482580bdb30SBarry Smith ierr = PetscArrayzero(ijacobian, 9);CHKERRQ(ierr); 48363d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 484181f196bSVijay Mahadevan const PetscReal* vertex = coords + i * 3; 48563d025dbSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertex[0]; 48663d025dbSVijay Mahadevan jacobian[3] += dNi_dxi[i] * vertex[1]; 48763d025dbSVijay Mahadevan jacobian[6] += dNi_dxi[i] * vertex[2]; 48863d025dbSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertex[0]; 48963d025dbSVijay Mahadevan jacobian[4] += dNi_deta[i] * vertex[1]; 49063d025dbSVijay Mahadevan jacobian[7] += dNi_deta[i] * vertex[2]; 49163d025dbSVijay Mahadevan jacobian[2] += dNi_dzeta[i] * vertex[0]; 49263d025dbSVijay Mahadevan jacobian[5] += dNi_dzeta[i] * vertex[1]; 49363d025dbSVijay Mahadevan jacobian[8] += dNi_dzeta[i] * vertex[2]; 494181f196bSVijay Mahadevan if (phypts) { 495181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertex[0]; 496181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertex[1]; 497181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertex[2]; 498181f196bSVijay Mahadevan } 49963d025dbSVijay Mahadevan } 50063d025dbSVijay Mahadevan 50163d025dbSVijay Mahadevan /* invert the jacobian */ 502181f196bSVijay Mahadevan ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 503a86ed394SVijay 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); 50463d025dbSVijay Mahadevan 505a86ed394SVijay Mahadevan if (jxw) jxw[j] *= (*volume); 50663d025dbSVijay Mahadevan 50763d025dbSVijay Mahadevan /* Divide by element jacobian. */ 50863d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 509a86ed394SVijay Mahadevan for (k = 0; k < 3; ++k) { 51063d025dbSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0 * 3 + k]; 51163d025dbSVijay Mahadevan if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[1 * 3 + k]; 51263d025dbSVijay Mahadevan if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k]; 51363d025dbSVijay Mahadevan } 51463d025dbSVijay Mahadevan } 51563d025dbSVijay Mahadevan } 5162da392ccSBarry Smith } else if (nverts == 4) { /* Linear Tetrahedra */ 5172da392ccSBarry Smith PetscReal Dx[4]={0,0,0,0}, Dy[4]={0,0,0,0}, Dz[4]={0,0,0,0}; 5182da392ccSBarry Smith const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2]; 51963d025dbSVijay Mahadevan 520580bdb30SBarry Smith ierr = PetscArrayzero(jacobian, 9);CHKERRQ(ierr); 521580bdb30SBarry Smith ierr = PetscArrayzero(ijacobian, 9);CHKERRQ(ierr); 52263d025dbSVijay Mahadevan 523181f196bSVijay Mahadevan /* compute the jacobian */ 524181f196bSVijay Mahadevan jacobian[0] = coords[1 * 3 + 0] - x0; jacobian[1] = coords[2 * 3 + 0] - x0; jacobian[2] = coords[3 * 3 + 0] - x0; 525181f196bSVijay Mahadevan jacobian[3] = coords[1 * 3 + 1] - y0; jacobian[4] = coords[2 * 3 + 1] - y0; jacobian[5] = coords[3 * 3 + 1] - y0; 526181f196bSVijay Mahadevan jacobian[6] = coords[1 * 3 + 2] - z0; jacobian[7] = coords[2 * 3 + 2] - z0; jacobian[8] = coords[3 * 3 + 2] - z0; 52763d025dbSVijay Mahadevan 52863d025dbSVijay Mahadevan /* invert the jacobian */ 529181f196bSVijay Mahadevan ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 5302da392ccSBarry Smith if (*volume < 1e-8) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral element has zero volume: %g. Degenerate element or invalid connectivity\n", (double)*volume); 53163d025dbSVijay Mahadevan 532181f196bSVijay Mahadevan if (dphidx) { 533181f196bSVijay Mahadevan Dx[0] = ( coords[1 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3]) 534181f196bSVijay Mahadevan - coords[1 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3]) 5352da392ccSBarry Smith - coords[1 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume; 536181f196bSVijay Mahadevan Dx[1] = - ( coords[1 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3]) 537181f196bSVijay Mahadevan - coords[1 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3]) 5382da392ccSBarry Smith - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume; 539181f196bSVijay Mahadevan Dx[2] = ( coords[1 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3]) 540181f196bSVijay Mahadevan - coords[1 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3]) 5412da392ccSBarry Smith - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume; 542181f196bSVijay Mahadevan Dx[3] = - ( Dx[0] + Dx[1] + Dx[2]); 543181f196bSVijay Mahadevan Dy[0] = ( coords[0 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3]) 544181f196bSVijay Mahadevan - coords[0 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3]) 5452da392ccSBarry Smith + coords[0 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3])) / *volume; 546181f196bSVijay Mahadevan Dy[1] = - ( coords[0 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3]) 547181f196bSVijay Mahadevan - coords[0 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3]) 5482da392ccSBarry Smith + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3])) / *volume; 549181f196bSVijay Mahadevan Dy[2] = ( coords[0 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3]) 550181f196bSVijay Mahadevan - coords[0 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3]) 5512da392ccSBarry Smith + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3])) / *volume; 552181f196bSVijay Mahadevan Dy[3] = - ( Dy[0] + Dy[1] + Dy[2]); 553181f196bSVijay Mahadevan Dz[0] = ( coords[0 + 1 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3]) 554181f196bSVijay Mahadevan - coords[0 + 2 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3]) 5552da392ccSBarry Smith + coords[0 + 3 * 3] * (coords[1 + 2 * 3] - coords[1 + 1 * 3])) / *volume; 556181f196bSVijay Mahadevan Dz[1] = - ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3]) 557181f196bSVijay Mahadevan + coords[0 + 2 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3]) 5582da392ccSBarry Smith - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 2 * 3])) / *volume; 559181f196bSVijay Mahadevan Dz[2] = ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3]) 560181f196bSVijay Mahadevan + coords[0 + 1 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3]) 5612da392ccSBarry Smith - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 1 * 3])) / *volume; 562181f196bSVijay Mahadevan Dz[3] = - ( Dz[0] + Dz[1] + Dz[2]); 563181f196bSVijay Mahadevan } 56463d025dbSVijay Mahadevan 5652da392ccSBarry Smith for (j = 0; j < npts; j++) { 566a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 567181f196bSVijay Mahadevan const PetscReal& r = quad[j * 3 + 0]; 568181f196bSVijay Mahadevan const PetscReal& s = quad[j * 3 + 1]; 569181f196bSVijay Mahadevan const PetscReal& t = quad[j * 3 + 2]; 57063d025dbSVijay Mahadevan 571181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume; 57263d025dbSVijay Mahadevan 57363d025dbSVijay Mahadevan phi[offset + 0] = 1.0 - r - s - t; 57463d025dbSVijay Mahadevan phi[offset + 1] = r; 57563d025dbSVijay Mahadevan phi[offset + 2] = s; 57663d025dbSVijay Mahadevan phi[offset + 3] = t; 57763d025dbSVijay Mahadevan 578181f196bSVijay Mahadevan if (phypts) { 579181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 580181f196bSVijay Mahadevan const PetscScalar* vertices = coords + i * 3; 581181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 582181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertices[1]; 583181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertices[2]; 584181f196bSVijay Mahadevan } 585181f196bSVijay Mahadevan } 586181f196bSVijay Mahadevan 587181f196bSVijay Mahadevan /* Now set the derivatives */ 58863d025dbSVijay Mahadevan if (dphidx) { 589181f196bSVijay Mahadevan dphidx[0 + offset] = Dx[0]; 590181f196bSVijay Mahadevan dphidx[1 + offset] = Dx[1]; 591181f196bSVijay Mahadevan dphidx[2 + offset] = Dx[2]; 592181f196bSVijay Mahadevan dphidx[3 + offset] = Dx[3]; 59363d025dbSVijay Mahadevan 594181f196bSVijay Mahadevan dphidy[0 + offset] = Dy[0]; 595181f196bSVijay Mahadevan dphidy[1 + offset] = Dy[1]; 596181f196bSVijay Mahadevan dphidy[2 + offset] = Dy[2]; 597181f196bSVijay Mahadevan dphidy[3 + offset] = Dy[3]; 59863d025dbSVijay Mahadevan 599181f196bSVijay Mahadevan dphidz[0 + offset] = Dz[0]; 600181f196bSVijay Mahadevan dphidz[1 + offset] = Dz[1]; 601181f196bSVijay Mahadevan dphidz[2 + offset] = Dz[2]; 602181f196bSVijay Mahadevan dphidz[3 + offset] = Dz[3]; 60363d025dbSVijay Mahadevan } 60463d025dbSVijay Mahadevan 60563d025dbSVijay Mahadevan } /* Tetrahedra -- ends */ 6062da392ccSBarry Smith } else 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); 60763d025dbSVijay Mahadevan PetscFunctionReturn(0); 60863d025dbSVijay Mahadevan } 60963d025dbSVijay Mahadevan 610cab5ea25SPierre Jolivet /*@C 61197b73a88SSatish Balay DMMoabFEMComputeBasis - Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element. 61263d025dbSVijay Mahadevan 61363d025dbSVijay Mahadevan The routine takes the coordinates of the vertices of an element and computes basis functions associated with 61463d025dbSVijay Mahadevan each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate. 61563d025dbSVijay Mahadevan 616d8d19677SJose E. Roman Input Parameters: 617a2b725a8SWilliam Gropp + PetscInt nverts, the number of element vertices 61863d025dbSVijay Mahadevan . PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering) 61963d025dbSVijay Mahadevan . PetscInt npts, the number of evaluation points (quadrature points) 620a2b725a8SWilliam Gropp - PetscReal quad[3*npts], the evaluation points (quadrature points) in the reference space 62163d025dbSVijay Mahadevan 622d8d19677SJose E. Roman Output Parameters: 623a2b725a8SWilliam Gropp + PetscReal phypts[3*npts], the evaluation points (quadrature points) transformed to the physical space 62463d025dbSVijay Mahadevan . PetscReal jxw[npts], the jacobian determinant * quadrature weight necessary for assembling discrete contributions 62563d025dbSVijay Mahadevan . PetscReal fe_basis[npts], the bases values evaluated at the specified quadrature points 626a2b725a8SWilliam 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 62763d025dbSVijay Mahadevan 628edc382c3SSatish Balay Level: advanced 629edc382c3SSatish Balay 63063d025dbSVijay Mahadevan @*/ 631181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMComputeBasis(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature, 632181f196bSVijay Mahadevan PetscReal *phypts, PetscReal *jacobian_quadrature_weight_product, 633181f196bSVijay Mahadevan PetscReal *fe_basis, PetscReal **fe_basis_derivatives) 63463d025dbSVijay Mahadevan { 63563d025dbSVijay Mahadevan PetscErrorCode ierr; 636b21a1e88SVijay Mahadevan PetscInt npoints,idim; 63763d025dbSVijay Mahadevan bool compute_der; 63863d025dbSVijay Mahadevan const PetscReal *quadpts, *quadwts; 639181f196bSVijay Mahadevan PetscReal jacobian[9], ijacobian[9], volume; 64063d025dbSVijay Mahadevan 64163d025dbSVijay Mahadevan PetscFunctionBegin; 64263d025dbSVijay Mahadevan PetscValidPointer(coordinates, 3); 64378dc7ee3SMatthew G. Knepley PetscValidHeaderSpecific(quadrature, PETSCQUADRATURE_CLASSID, 4); 64463d025dbSVijay Mahadevan PetscValidPointer(fe_basis, 7); 64563d025dbSVijay Mahadevan compute_der = (fe_basis_derivatives != NULL); 64663d025dbSVijay Mahadevan 64763d025dbSVijay Mahadevan /* Get the quadrature points and weights for the given quadrature rule */ 648b21a1e88SVijay Mahadevan ierr = PetscQuadratureGetData(quadrature, &idim, NULL, &npoints, &quadpts, &quadwts);CHKERRQ(ierr); 649b21a1e88SVijay Mahadevan if (idim != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Dimension mismatch: provided (%D) vs quadrature (%D)\n",idim,dim); 650181f196bSVijay Mahadevan if (jacobian_quadrature_weight_product) { 651580bdb30SBarry Smith ierr = PetscArraycpy(jacobian_quadrature_weight_product, quadwts, npoints);CHKERRQ(ierr); 652181f196bSVijay Mahadevan } 65363d025dbSVijay Mahadevan 65463d025dbSVijay Mahadevan switch (dim) { 65563d025dbSVijay Mahadevan case 1: 65663d025dbSVijay Mahadevan ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, npoints, quadpts, phypts, 65763d025dbSVijay Mahadevan jacobian_quadrature_weight_product, fe_basis, 658181f196bSVijay Mahadevan (compute_der ? fe_basis_derivatives[0] : NULL), 659181f196bSVijay Mahadevan jacobian, ijacobian, &volume);CHKERRQ(ierr); 66063d025dbSVijay Mahadevan break; 66163d025dbSVijay Mahadevan case 2: 66263d025dbSVijay Mahadevan ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, npoints, quadpts, phypts, 66363d025dbSVijay Mahadevan jacobian_quadrature_weight_product, fe_basis, 66463d025dbSVijay Mahadevan (compute_der ? fe_basis_derivatives[0] : NULL), 665181f196bSVijay Mahadevan (compute_der ? fe_basis_derivatives[1] : NULL), 666181f196bSVijay Mahadevan jacobian, ijacobian, &volume);CHKERRQ(ierr); 66763d025dbSVijay Mahadevan break; 66863d025dbSVijay Mahadevan case 3: 66963d025dbSVijay Mahadevan ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, npoints, quadpts, phypts, 67063d025dbSVijay Mahadevan jacobian_quadrature_weight_product, fe_basis, 67163d025dbSVijay Mahadevan (compute_der ? fe_basis_derivatives[0] : NULL), 67263d025dbSVijay Mahadevan (compute_der ? fe_basis_derivatives[1] : NULL), 673181f196bSVijay Mahadevan (compute_der ? fe_basis_derivatives[2] : NULL), 674181f196bSVijay Mahadevan jacobian, ijacobian, &volume);CHKERRQ(ierr); 67563d025dbSVijay Mahadevan break; 67663d025dbSVijay Mahadevan default: 67763d025dbSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim); 67863d025dbSVijay Mahadevan } 67963d025dbSVijay Mahadevan PetscFunctionReturn(0); 68063d025dbSVijay Mahadevan } 68163d025dbSVijay Mahadevan 682cab5ea25SPierre Jolivet /*@C 68397b73a88SSatish Balay DMMoabFEMCreateQuadratureDefault - Create default quadrature rules for integration over an element with a given 68463d025dbSVijay Mahadevan dimension and polynomial order (deciphered from number of element vertices). 68563d025dbSVijay Mahadevan 686d8d19677SJose E. Roman Input Parameters: 68763d025dbSVijay Mahadevan 688a2b725a8SWilliam Gropp + PetscInt dim - the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET) 689a2b725a8SWilliam Gropp - PetscInt nverts - the number of vertices in the physical element 69063d025dbSVijay Mahadevan 69163d025dbSVijay Mahadevan Output Parameter: 692a2b725a8SWilliam Gropp . PetscQuadrature quadrature - the quadrature object with default settings to integrate polynomials defined over the element 69363d025dbSVijay Mahadevan 694edc382c3SSatish Balay Level: advanced 695edc382c3SSatish Balay 69663d025dbSVijay Mahadevan @*/ 697181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMCreateQuadratureDefault(const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature) 69863d025dbSVijay Mahadevan { 69963d025dbSVijay Mahadevan PetscReal *w, *x; 700b21a1e88SVijay Mahadevan PetscInt nc=1; 70163d025dbSVijay Mahadevan PetscErrorCode ierr; 70263d025dbSVijay Mahadevan 70363d025dbSVijay Mahadevan PetscFunctionBegin; 70463d025dbSVijay Mahadevan /* Create an appropriate quadrature rule to sample basis */ 70563d025dbSVijay Mahadevan switch (dim) 70663d025dbSVijay Mahadevan { 70763d025dbSVijay Mahadevan case 1: 70863d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */ 709e6a796c3SToby Isaac ierr = PetscDTStroudConicalQuadrature(1, nc, nverts, 0, 1.0, quadrature);CHKERRQ(ierr); 71063d025dbSVijay Mahadevan break; 71163d025dbSVijay Mahadevan case 2: 71263d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */ 71363d025dbSVijay Mahadevan if (nverts == 3) { 714a86ed394SVijay Mahadevan const PetscInt order = 2; 715a86ed394SVijay Mahadevan const PetscInt npoints = (order == 2 ? 3 : 6); 71663d025dbSVijay Mahadevan ierr = PetscMalloc2(npoints * 2, &x, npoints, &w);CHKERRQ(ierr); 717181f196bSVijay Mahadevan if (npoints == 3) { 71863d025dbSVijay Mahadevan x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0; 71963d025dbSVijay Mahadevan x[3] = x[4] = 2.0 / 3.0; 72063d025dbSVijay Mahadevan w[0] = w[1] = w[2] = 1.0 / 3.0; 7212da392ccSBarry Smith } else if (npoints == 6) { 72263d025dbSVijay Mahadevan x[0] = x[1] = x[2] = 0.44594849091597; 72363d025dbSVijay Mahadevan x[3] = x[4] = 0.10810301816807; 72463d025dbSVijay Mahadevan x[5] = 0.44594849091597; 72563d025dbSVijay Mahadevan x[6] = x[7] = x[8] = 0.09157621350977; 72663d025dbSVijay Mahadevan x[9] = x[10] = 0.81684757298046; 72763d025dbSVijay Mahadevan x[11] = 0.09157621350977; 72863d025dbSVijay Mahadevan w[0] = w[1] = w[2] = 0.22338158967801; 729181f196bSVijay Mahadevan w[3] = w[4] = w[5] = 0.10995174365532; 7302da392ccSBarry Smith } else SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %D", npoints); 73163d025dbSVijay Mahadevan ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr); 73263d025dbSVijay Mahadevan ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr); 733b21a1e88SVijay Mahadevan ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr); 734e6a796c3SToby Isaac /* ierr = PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */ 7352da392ccSBarry Smith } else { 736b21a1e88SVijay Mahadevan ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); 73763d025dbSVijay Mahadevan } 73863d025dbSVijay Mahadevan break; 73963d025dbSVijay Mahadevan case 3: 74063d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */ 74163d025dbSVijay Mahadevan if (nverts == 4) { 742a86ed394SVijay Mahadevan PetscInt order; 743a86ed394SVijay Mahadevan const PetscInt npoints = 4; // Choose between 4 and 10 744181f196bSVijay Mahadevan ierr = PetscMalloc2(npoints * 3, &x, npoints, &w);CHKERRQ(ierr); 745181f196bSVijay Mahadevan if (npoints == 4) { /* KEAST1, K1, N=4, O=4 */ 746181f196bSVijay Mahadevan const PetscReal x_4[12] = { 0.5854101966249685, 0.1381966011250105, 0.1381966011250105, 747181f196bSVijay Mahadevan 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 748181f196bSVijay Mahadevan 0.1381966011250105, 0.1381966011250105, 0.5854101966249685, 749181f196bSVijay Mahadevan 0.1381966011250105, 0.5854101966249685, 0.1381966011250105 750181f196bSVijay Mahadevan }; 751580bdb30SBarry Smith ierr = PetscArraycpy(x, x_4, 12);CHKERRQ(ierr); 752181f196bSVijay Mahadevan 753181f196bSVijay Mahadevan w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0; 754181f196bSVijay Mahadevan order = 4; 7552da392ccSBarry Smith } else if (npoints == 10) { /* KEAST3, K3 N=10, O=10 */ 756181f196bSVijay Mahadevan const PetscReal x_10[30] = { 0.5684305841968444, 0.1438564719343852, 0.1438564719343852, 757181f196bSVijay Mahadevan 0.1438564719343852, 0.1438564719343852, 0.1438564719343852, 758181f196bSVijay Mahadevan 0.1438564719343852, 0.1438564719343852, 0.5684305841968444, 759181f196bSVijay Mahadevan 0.1438564719343852, 0.5684305841968444, 0.1438564719343852, 760181f196bSVijay Mahadevan 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 761181f196bSVijay Mahadevan 0.5000000000000000, 0.0000000000000000, 0.5000000000000000, 762181f196bSVijay Mahadevan 0.5000000000000000, 0.5000000000000000, 0.0000000000000000, 763181f196bSVijay Mahadevan 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 764181f196bSVijay Mahadevan 0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 765181f196bSVijay Mahadevan 0.0000000000000000, 0.0000000000000000, 0.5000000000000000 766181f196bSVijay Mahadevan }; 767580bdb30SBarry Smith ierr = PetscArraycpy(x, x_10, 30);CHKERRQ(ierr); 768181f196bSVijay Mahadevan 769181f196bSVijay Mahadevan w[0] = w[1] = w[2] = w[3] = 0.2177650698804054; 770181f196bSVijay Mahadevan w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631; 771181f196bSVijay Mahadevan order = 10; 7722da392ccSBarry Smith } else SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %D", npoints); 773181f196bSVijay Mahadevan ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr); 774181f196bSVijay Mahadevan ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr); 775b21a1e88SVijay Mahadevan ierr = PetscQuadratureSetData(*quadrature, dim, nc, npoints, x, w);CHKERRQ(ierr); 776e6a796c3SToby Isaac /* ierr = PetscDTStroudConicalQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */ 7772da392ccSBarry Smith } else { 778b21a1e88SVijay Mahadevan ierr = PetscDTGaussTensorQuadrature(dim, nc, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); 77963d025dbSVijay Mahadevan } 78063d025dbSVijay Mahadevan break; 78163d025dbSVijay Mahadevan default: 78263d025dbSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim); 78363d025dbSVijay Mahadevan } 78463d025dbSVijay Mahadevan PetscFunctionReturn(0); 78563d025dbSVijay Mahadevan } 78663d025dbSVijay Mahadevan 787181f196bSVijay Mahadevan /* Compute Jacobians */ 788181f196bSVijay Mahadevan PetscErrorCode ComputeJacobian_Internal (const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts, 789a86ed394SVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal* dvolume) 790181f196bSVijay Mahadevan { 791a86ed394SVijay Mahadevan PetscInt i; 7922417220eSVijay Mahadevan PetscReal volume=1.0; 793181f196bSVijay Mahadevan PetscErrorCode ierr; 794181f196bSVijay Mahadevan 795181f196bSVijay Mahadevan PetscFunctionBegin; 796181f196bSVijay Mahadevan PetscValidPointer(coordinates, 3); 797181f196bSVijay Mahadevan PetscValidPointer(quad, 4); 798181f196bSVijay Mahadevan PetscValidPointer(jacobian, 5); 799580bdb30SBarry Smith ierr = PetscArrayzero(jacobian, dim * dim);CHKERRQ(ierr); 800181f196bSVijay Mahadevan if (ijacobian) { 801580bdb30SBarry Smith ierr = PetscArrayzero(ijacobian, dim * dim);CHKERRQ(ierr); 802181f196bSVijay Mahadevan } 803181f196bSVijay Mahadevan if (phypts) { 804580bdb30SBarry Smith ierr = PetscArrayzero(phypts, /*npts=1 * */ 3);CHKERRQ(ierr); 805181f196bSVijay Mahadevan } 806181f196bSVijay Mahadevan 807181f196bSVijay Mahadevan if (dim == 1) { 808181f196bSVijay Mahadevan const PetscReal& r = quad[0]; 809181f196bSVijay Mahadevan if (nverts == 2) { /* Linear Edge */ 810181f196bSVijay Mahadevan const PetscReal dNi_dxi[2] = { -1.0, 1.0 }; 811181f196bSVijay Mahadevan 812181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 813181f196bSVijay Mahadevan const PetscReal* vertices = coordinates + i * 3; 814181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 815181f196bSVijay Mahadevan } 8162da392ccSBarry Smith } else if (nverts == 3) { /* Quadratic Edge */ 817181f196bSVijay Mahadevan const PetscReal dNi_dxi[3] = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r), 4.0 * r - 1.0}; 818181f196bSVijay Mahadevan 819181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 820181f196bSVijay Mahadevan const PetscReal* vertices = coordinates + i * 3; 821181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 822181f196bSVijay Mahadevan } 8232da392ccSBarry Smith } else { 824181f196bSVijay 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); 825181f196bSVijay Mahadevan } 826181f196bSVijay Mahadevan 827181f196bSVijay Mahadevan if (ijacobian) { 828181f196bSVijay Mahadevan /* invert the jacobian */ 829181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0]; 830181f196bSVijay Mahadevan } 8312da392ccSBarry Smith } else if (dim == 2) { 832181f196bSVijay Mahadevan 833181f196bSVijay Mahadevan if (nverts == 4) { /* Linear Quadrangle */ 834181f196bSVijay Mahadevan const PetscReal& r = quad[0]; 835181f196bSVijay Mahadevan const PetscReal& s = quad[1]; 836181f196bSVijay Mahadevan 837181f196bSVijay Mahadevan const PetscReal dNi_dxi[4] = { -1.0 + s, 1.0 - s, s, -s }; 838181f196bSVijay Mahadevan const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r }; 839181f196bSVijay Mahadevan 840181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 841181f196bSVijay Mahadevan const PetscReal* vertices = coordinates + i * 3; 842181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 843181f196bSVijay Mahadevan jacobian[2] += dNi_dxi[i] * vertices[1]; 844181f196bSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertices[0]; 845181f196bSVijay Mahadevan jacobian[3] += dNi_deta[i] * vertices[1]; 846181f196bSVijay Mahadevan } 8472da392ccSBarry Smith } else if (nverts == 3) { /* Linear triangle */ 848181f196bSVijay Mahadevan const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1]; 849181f196bSVijay Mahadevan 850181f196bSVijay Mahadevan /* Jacobian is constant */ 851181f196bSVijay Mahadevan jacobian[0] = (coordinates[0 * 3 + 0] - x2); jacobian[1] = (coordinates[1 * 3 + 0] - x2); 852181f196bSVijay Mahadevan jacobian[2] = (coordinates[0 * 3 + 1] - y2); jacobian[3] = (coordinates[1 * 3 + 1] - y2); 8532da392ccSBarry Smith } else 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); 854181f196bSVijay Mahadevan 855181f196bSVijay Mahadevan /* invert the jacobian */ 856181f196bSVijay Mahadevan if (ijacobian) { 857a86ed394SVijay Mahadevan ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr); 858181f196bSVijay Mahadevan } 8592da392ccSBarry Smith } else { 860181f196bSVijay Mahadevan 861181f196bSVijay Mahadevan if (nverts == 8) { /* Linear Hexahedra */ 862181f196bSVijay Mahadevan const PetscReal &r = quad[0]; 863181f196bSVijay Mahadevan const PetscReal &s = quad[1]; 864181f196bSVijay Mahadevan const PetscReal &t = quad[2]; 865181f196bSVijay Mahadevan const PetscReal dNi_dxi[8] = {- ( 1.0 - s) * ( 1.0 - t), 866181f196bSVijay Mahadevan ( 1.0 - s) * ( 1.0 - t), 867a86ed394SVijay Mahadevan ( s) * ( 1.0 - t), 868a86ed394SVijay Mahadevan - ( s) * ( 1.0 - t), 869a86ed394SVijay Mahadevan - ( 1.0 - s) * ( t), 870a86ed394SVijay Mahadevan ( 1.0 - s) * ( t), 871a86ed394SVijay Mahadevan ( s) * ( t), 872a86ed394SVijay Mahadevan - ( s) * ( t) 873181f196bSVijay Mahadevan }; 874181f196bSVijay Mahadevan 875181f196bSVijay Mahadevan const PetscReal dNi_deta[8] = { - ( 1.0 - r) * ( 1.0 - t), 876a86ed394SVijay Mahadevan - ( r) * ( 1.0 - t), 877a86ed394SVijay Mahadevan ( r) * ( 1.0 - t), 878181f196bSVijay Mahadevan ( 1.0 - r) * ( 1.0 - t), 879a86ed394SVijay Mahadevan - ( 1.0 - r) * ( t), 880a86ed394SVijay Mahadevan - ( r) * ( t), 881a86ed394SVijay Mahadevan ( r) * ( t), 882a86ed394SVijay Mahadevan ( 1.0 - r) * ( t) 883181f196bSVijay Mahadevan }; 884181f196bSVijay Mahadevan 885181f196bSVijay Mahadevan const PetscReal dNi_dzeta[8] = { - ( 1.0 - r) * ( 1.0 - s), 886a86ed394SVijay Mahadevan - ( r) * ( 1.0 - s), 887a86ed394SVijay Mahadevan - ( r) * ( s), 888a86ed394SVijay Mahadevan - ( 1.0 - r) * ( s), 889181f196bSVijay Mahadevan ( 1.0 - r) * ( 1.0 - s), 890a86ed394SVijay Mahadevan ( r) * ( 1.0 - s), 891a86ed394SVijay Mahadevan ( r) * ( s), 892a86ed394SVijay Mahadevan ( 1.0 - r) * ( s) 893181f196bSVijay Mahadevan }; 894a86ed394SVijay Mahadevan 895181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 896181f196bSVijay Mahadevan const PetscReal* vertex = coordinates + i * 3; 897181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertex[0]; 898181f196bSVijay Mahadevan jacobian[3] += dNi_dxi[i] * vertex[1]; 899181f196bSVijay Mahadevan jacobian[6] += dNi_dxi[i] * vertex[2]; 900181f196bSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertex[0]; 901181f196bSVijay Mahadevan jacobian[4] += dNi_deta[i] * vertex[1]; 902181f196bSVijay Mahadevan jacobian[7] += dNi_deta[i] * vertex[2]; 903181f196bSVijay Mahadevan jacobian[2] += dNi_dzeta[i] * vertex[0]; 904181f196bSVijay Mahadevan jacobian[5] += dNi_dzeta[i] * vertex[1]; 905181f196bSVijay Mahadevan jacobian[8] += dNi_dzeta[i] * vertex[2]; 906181f196bSVijay Mahadevan } 9072da392ccSBarry Smith } else if (nverts == 4) { /* Linear Tetrahedra */ 908181f196bSVijay Mahadevan const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2]; 909181f196bSVijay Mahadevan 910181f196bSVijay Mahadevan /* compute the jacobian */ 911181f196bSVijay Mahadevan jacobian[0] = coordinates[1 * 3 + 0] - x0; jacobian[1] = coordinates[2 * 3 + 0] - x0; jacobian[2] = coordinates[3 * 3 + 0] - x0; 912181f196bSVijay Mahadevan jacobian[3] = coordinates[1 * 3 + 1] - y0; jacobian[4] = coordinates[2 * 3 + 1] - y0; jacobian[5] = coordinates[3 * 3 + 1] - y0; 913181f196bSVijay Mahadevan jacobian[6] = coordinates[1 * 3 + 2] - z0; jacobian[7] = coordinates[2 * 3 + 2] - z0; jacobian[8] = coordinates[3 * 3 + 2] - z0; 9142da392ccSBarry Smith } else 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); 915181f196bSVijay Mahadevan 916181f196bSVijay Mahadevan if (ijacobian) { 917181f196bSVijay Mahadevan /* invert the jacobian */ 918a86ed394SVijay Mahadevan ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr); 919181f196bSVijay Mahadevan } 920181f196bSVijay Mahadevan 921181f196bSVijay Mahadevan } 922a86ed394SVijay 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); 923a86ed394SVijay Mahadevan if (dvolume) *dvolume = volume; 924181f196bSVijay Mahadevan PetscFunctionReturn(0); 925181f196bSVijay Mahadevan } 926181f196bSVijay Mahadevan 927181f196bSVijay Mahadevan PetscErrorCode FEMComputeBasis_JandF(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quadrature, PetscReal *phypts, 928181f196bSVijay Mahadevan PetscReal *phibasis, PetscReal *jacobian, PetscReal *ijacobian, PetscReal* volume) 929181f196bSVijay Mahadevan { 930181f196bSVijay Mahadevan PetscErrorCode ierr; 931181f196bSVijay Mahadevan 932181f196bSVijay Mahadevan PetscFunctionBegin; 933181f196bSVijay Mahadevan switch (dim) { 934181f196bSVijay Mahadevan case 1: 935181f196bSVijay Mahadevan ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts, 936181f196bSVijay Mahadevan NULL, phibasis, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr); 937181f196bSVijay Mahadevan break; 938181f196bSVijay Mahadevan case 2: 939181f196bSVijay Mahadevan ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts, 940181f196bSVijay Mahadevan NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr); 941181f196bSVijay Mahadevan break; 942181f196bSVijay Mahadevan case 3: 943181f196bSVijay Mahadevan ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts, 944181f196bSVijay Mahadevan NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr); 945181f196bSVijay Mahadevan break; 946181f196bSVijay Mahadevan default: 947181f196bSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim); 948181f196bSVijay Mahadevan } 949181f196bSVijay Mahadevan PetscFunctionReturn(0); 950181f196bSVijay Mahadevan } 951181f196bSVijay Mahadevan 952cab5ea25SPierre Jolivet /*@C 95397b73a88SSatish Balay DMMoabPToRMapping - Compute the mapping from the physical coordinate system for a given element to the 954a86ed394SVijay Mahadevan canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration, 955a86ed394SVijay Mahadevan the basis function at the parametric point is also evaluated optionally. 956a86ed394SVijay Mahadevan 957d8d19677SJose E. Roman Input Parameters: 958a2b725a8SWilliam Gropp + PetscInt dim - the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET) 959a2b725a8SWilliam Gropp . PetscInt nverts - the number of vertices in the physical element 960a2b725a8SWilliam Gropp . PetscReal coordinates - the coordinates of vertices in the physical element 961a2b725a8SWilliam Gropp - PetscReal[3] xphy - the coordinates of physical point for which natural coordinates (in reference frame) are sought 962a86ed394SVijay Mahadevan 963d8d19677SJose E. Roman Output Parameters: 964a2b725a8SWilliam Gropp + PetscReal[3] natparam - the natural coordinates (in reference frame) corresponding to xphy 965a2b725a8SWilliam Gropp - PetscReal[nverts] phi - the basis functions evaluated at the natural coordinates (natparam) 966a86ed394SVijay Mahadevan 967edc382c3SSatish Balay Level: advanced 968edc382c3SSatish Balay 969a86ed394SVijay Mahadevan @*/ 970181f196bSVijay Mahadevan PetscErrorCode DMMoabPToRMapping(const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal* xphy, PetscReal* natparam, PetscReal* phi) 971181f196bSVijay Mahadevan { 972a86ed394SVijay Mahadevan /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */ 973181f196bSVijay Mahadevan const PetscReal tol = 1.0e-10; 974181f196bSVijay Mahadevan const PetscInt max_iterations = 10; 975181f196bSVijay Mahadevan const PetscReal error_tol_sqr = tol*tol; 976181f196bSVijay Mahadevan PetscReal phibasis[8], jacobian[9], ijacobian[9], volume; 977181f196bSVijay Mahadevan PetscReal phypts[3] = {0.0, 0.0, 0.0}; 978181f196bSVijay Mahadevan PetscReal delta[3] = {0.0, 0.0, 0.0}; 979181f196bSVijay Mahadevan PetscErrorCode ierr; 980181f196bSVijay Mahadevan PetscInt iters=0; 981181f196bSVijay Mahadevan PetscReal error=1.0; 982181f196bSVijay Mahadevan 983181f196bSVijay Mahadevan PetscFunctionBegin; 984181f196bSVijay Mahadevan PetscValidPointer(coordinates, 3); 985181f196bSVijay Mahadevan PetscValidPointer(xphy, 4); 986181f196bSVijay Mahadevan PetscValidPointer(natparam, 5); 987181f196bSVijay Mahadevan 988580bdb30SBarry Smith ierr = PetscArrayzero(jacobian, dim * dim);CHKERRQ(ierr); 989580bdb30SBarry Smith ierr = PetscArrayzero(ijacobian, dim * dim);CHKERRQ(ierr); 990580bdb30SBarry Smith ierr = PetscArrayzero(phibasis, nverts);CHKERRQ(ierr); 991181f196bSVijay Mahadevan 992a86ed394SVijay Mahadevan /* zero initial guess */ 993a86ed394SVijay Mahadevan natparam[0] = natparam[1] = natparam[2] = 0.0; 994181f196bSVijay Mahadevan 995a86ed394SVijay Mahadevan /* Compute delta = evaluate( xi) - x */ 996a86ed394SVijay Mahadevan ierr = FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume);CHKERRQ(ierr); 997a86ed394SVijay Mahadevan 998a86ed394SVijay Mahadevan error = 0.0; 999a86ed394SVijay Mahadevan switch(dim) { 1000a86ed394SVijay Mahadevan case 3: 1001181f196bSVijay Mahadevan delta[2] = phypts[2] - xphy[2]; 1002a86ed394SVijay Mahadevan error += (delta[2]*delta[2]); 1003a86ed394SVijay Mahadevan case 2: 1004a86ed394SVijay Mahadevan delta[1] = phypts[1] - xphy[1]; 1005a86ed394SVijay Mahadevan error += (delta[1]*delta[1]); 1006a86ed394SVijay Mahadevan case 1: 1007a86ed394SVijay Mahadevan delta[0] = phypts[0] - xphy[0]; 1008a86ed394SVijay Mahadevan error += (delta[0]*delta[0]); 1009a86ed394SVijay Mahadevan break; 1010a86ed394SVijay Mahadevan } 1011a86ed394SVijay Mahadevan 1012181f196bSVijay Mahadevan while (error > error_tol_sqr) { 1013181f196bSVijay Mahadevan 1014181f196bSVijay Mahadevan if (++iters > max_iterations) 1015181f196bSVijay 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]); 1016181f196bSVijay Mahadevan 1017181f196bSVijay Mahadevan /* Compute natparam -= J.inverse() * delta */ 1018181f196bSVijay Mahadevan switch(dim) { 1019181f196bSVijay Mahadevan case 1: 1020181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0]; 1021181f196bSVijay Mahadevan break; 1022181f196bSVijay Mahadevan case 2: 1023181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1]; 1024181f196bSVijay Mahadevan natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1]; 1025181f196bSVijay Mahadevan break; 1026181f196bSVijay Mahadevan case 3: 1027181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2]; 1028181f196bSVijay Mahadevan natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2]; 1029181f196bSVijay Mahadevan natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2]; 1030181f196bSVijay Mahadevan break; 1031181f196bSVijay Mahadevan } 1032181f196bSVijay Mahadevan 1033181f196bSVijay Mahadevan /* Compute delta = evaluate( xi) - x */ 1034a86ed394SVijay Mahadevan ierr = FEMComputeBasis_JandF(dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume);CHKERRQ(ierr); 1035181f196bSVijay Mahadevan 1036a86ed394SVijay Mahadevan error = 0.0; 1037a86ed394SVijay Mahadevan switch(dim) { 1038a86ed394SVijay Mahadevan case 3: 1039181f196bSVijay Mahadevan delta[2] = phypts[2] - xphy[2]; 1040a86ed394SVijay Mahadevan error += (delta[2]*delta[2]); 1041a86ed394SVijay Mahadevan case 2: 1042a86ed394SVijay Mahadevan delta[1] = phypts[1] - xphy[1]; 1043a86ed394SVijay Mahadevan error += (delta[1]*delta[1]); 1044a86ed394SVijay Mahadevan case 1: 1045a86ed394SVijay Mahadevan delta[0] = phypts[0] - xphy[0]; 1046a86ed394SVijay Mahadevan error += (delta[0]*delta[0]); 1047a86ed394SVijay Mahadevan break; 1048a86ed394SVijay Mahadevan } 1049181f196bSVijay Mahadevan } 1050181f196bSVijay Mahadevan if (phi) { 1051580bdb30SBarry Smith ierr = PetscArraycpy(phi, phibasis, nverts);CHKERRQ(ierr); 1052181f196bSVijay Mahadevan } 1053181f196bSVijay Mahadevan PetscFunctionReturn(0); 1054181f196bSVijay Mahadevan } 1055