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 12181f196bSVijay Mahadevan #undef __FUNCT__ 13181f196bSVijay Mahadevan #define __FUNCT__ "DMatrix_Invert_2x2_Internal" 1463d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_2x2_Internal (const PetscReal *inmat, PetscReal *outmat, PetscReal *determinant) 1563d025dbSVijay Mahadevan { 1663d025dbSVijay Mahadevan if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 2x2 inversion."); 1763d025dbSVijay Mahadevan PetscReal det = DMatrix_Determinant_2x2_Internal(inmat); 1863d025dbSVijay Mahadevan if (outmat) { 1963d025dbSVijay Mahadevan outmat[0] = inmat[3] / det; 2063d025dbSVijay Mahadevan outmat[1] = -inmat[1] / det; 2163d025dbSVijay Mahadevan outmat[2] = -inmat[2] / det; 2263d025dbSVijay Mahadevan outmat[3] = inmat[0] / det; 2363d025dbSVijay Mahadevan } 2463d025dbSVijay Mahadevan if (determinant) *determinant = det; 2563d025dbSVijay Mahadevan PetscFunctionReturn(0); 2663d025dbSVijay Mahadevan } 2763d025dbSVijay Mahadevan 28181f196bSVijay Mahadevan static inline PetscReal DMatrix_Determinant_3x3_Internal ( const PetscReal inmat[3 * 3] ) 2963d025dbSVijay Mahadevan { 3063d025dbSVijay Mahadevan return inmat[0] * (inmat[8] * inmat[4] - inmat[7] * inmat[5]) 3163d025dbSVijay Mahadevan - inmat[3] * (inmat[8] * inmat[1] - inmat[7] * inmat[2]) 3263d025dbSVijay Mahadevan + inmat[6] * (inmat[5] * inmat[1] - inmat[4] * inmat[2]); 3363d025dbSVijay Mahadevan } 3463d025dbSVijay Mahadevan 35181f196bSVijay Mahadevan #undef __FUNCT__ 36181f196bSVijay Mahadevan #define __FUNCT__ "DMatrix_Invert_3x3_Internal" 3763d025dbSVijay Mahadevan static inline PetscErrorCode DMatrix_Invert_3x3_Internal (const PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant) 3863d025dbSVijay Mahadevan { 3963d025dbSVijay Mahadevan if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 3x3 inversion."); 40181f196bSVijay Mahadevan PetscReal det = DMatrix_Determinant_3x3_Internal(inmat); 4163d025dbSVijay Mahadevan if (outmat) { 4263d025dbSVijay Mahadevan outmat[0] = (inmat[8] * inmat[4] - inmat[7] * inmat[5]) / det; 4363d025dbSVijay Mahadevan outmat[1] = -(inmat[8] * inmat[1] - inmat[7] * inmat[2]) / det; 4463d025dbSVijay Mahadevan outmat[2] = (inmat[5] * inmat[1] - inmat[4] * inmat[2]) / det; 4563d025dbSVijay Mahadevan outmat[3] = -(inmat[8] * inmat[3] - inmat[6] * inmat[5]) / det; 4663d025dbSVijay Mahadevan outmat[4] = (inmat[8] * inmat[0] - inmat[6] * inmat[2]) / det; 4763d025dbSVijay Mahadevan outmat[5] = -(inmat[5] * inmat[0] - inmat[3] * inmat[2]) / det; 4863d025dbSVijay Mahadevan outmat[6] = (inmat[7] * inmat[3] - inmat[6] * inmat[4]) / det; 4963d025dbSVijay Mahadevan outmat[7] = -(inmat[7] * inmat[0] - inmat[6] * inmat[1]) / det; 5063d025dbSVijay Mahadevan outmat[8] = (inmat[4] * inmat[0] - inmat[3] * inmat[1]) / det; 5163d025dbSVijay Mahadevan } 5263d025dbSVijay Mahadevan if (determinant) *determinant = det; 5363d025dbSVijay Mahadevan PetscFunctionReturn(0); 5463d025dbSVijay Mahadevan } 5563d025dbSVijay Mahadevan 56181f196bSVijay Mahadevan inline PetscReal DMatrix_Determinant_4x4_Internal ( PetscReal inmat[4 * 4] ) 57181f196bSVijay Mahadevan { 58181f196bSVijay Mahadevan return 59181f196bSVijay Mahadevan inmat[0 + 0 * 4] * ( 60181f196bSVijay Mahadevan inmat[1 + 1 * 4] * ( inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4] ) 61181f196bSVijay Mahadevan - inmat[1 + 2 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4] ) 62181f196bSVijay Mahadevan + inmat[1 + 3 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4] ) ) 63181f196bSVijay Mahadevan - inmat[0 + 1 * 4] * ( 64181f196bSVijay Mahadevan inmat[1 + 0 * 4] * ( inmat[2 + 2 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 2 * 4] ) 65181f196bSVijay Mahadevan - inmat[1 + 2 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4] ) 66181f196bSVijay Mahadevan + inmat[1 + 3 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4] ) ) 67181f196bSVijay Mahadevan + inmat[0 + 2 * 4] * ( 68181f196bSVijay Mahadevan inmat[1 + 0 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 1 * 4] ) 69181f196bSVijay Mahadevan - inmat[1 + 1 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 3 * 4] - inmat[2 + 3 * 4] * inmat[3 + 0 * 4] ) 70181f196bSVijay Mahadevan + inmat[1 + 3 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4] ) ) 71181f196bSVijay Mahadevan - inmat[0 + 3 * 4] * ( 72181f196bSVijay Mahadevan inmat[1 + 0 * 4] * ( inmat[2 + 1 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 1 * 4] ) 73181f196bSVijay Mahadevan - inmat[1 + 1 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 2 * 4] - inmat[2 + 2 * 4] * inmat[3 + 0 * 4] ) 74181f196bSVijay Mahadevan + inmat[1 + 2 * 4] * ( inmat[2 + 0 * 4] * inmat[3 + 1 * 4] - inmat[2 + 1 * 4] * inmat[3 + 0 * 4] ) ); 75181f196bSVijay Mahadevan } 76181f196bSVijay Mahadevan 77181f196bSVijay Mahadevan #undef __FUNCT__ 78181f196bSVijay Mahadevan #define __FUNCT__ "DMatrix_Invert_4x4_Internal" 79181f196bSVijay Mahadevan inline PetscErrorCode DMatrix_Invert_4x4_Internal (PetscReal *inmat, PetscReal *outmat, PetscScalar *determinant) 80181f196bSVijay Mahadevan { 81181f196bSVijay Mahadevan if (!inmat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_POINTER, "Invalid input matrix specified for 4x4 inversion."); 82181f196bSVijay Mahadevan PetscReal det = DMatrix_Determinant_4x4_Internal(inmat); 83181f196bSVijay Mahadevan if (outmat) { 84181f196bSVijay 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; 85181f196bSVijay 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; 86181f196bSVijay 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; 87181f196bSVijay 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; 88181f196bSVijay 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; 89181f196bSVijay 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; 90181f196bSVijay 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; 91181f196bSVijay 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; 92181f196bSVijay 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; 93181f196bSVijay 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; 94181f196bSVijay 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; 95181f196bSVijay 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; 96181f196bSVijay 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; 97181f196bSVijay 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; 98181f196bSVijay 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; 99181f196bSVijay 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; 100181f196bSVijay Mahadevan } 101181f196bSVijay Mahadevan if (determinant) *determinant = det; 102181f196bSVijay Mahadevan PetscFunctionReturn(0); 103181f196bSVijay Mahadevan } 104181f196bSVijay Mahadevan 10563d025dbSVijay Mahadevan 10663d025dbSVijay Mahadevan /*@ 10763d025dbSVijay Mahadevan Compute_Lagrange_Basis_1D_Internal: Evaluate bases and derivatives at quadrature points for a EDGE2 or EDGE3 element. 10863d025dbSVijay Mahadevan 10963d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a linear or quadratic edge element. 11063d025dbSVijay Mahadevan 11163d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided, 11263d025dbSVijay Mahadevan and their derivatives with respect to X. 11363d025dbSVijay Mahadevan 11463d025dbSVijay Mahadevan Notes: 11563d025dbSVijay Mahadevan 11663d025dbSVijay Mahadevan Example Physical Element 11763d025dbSVijay Mahadevan 11863d025dbSVijay Mahadevan 1-------2 1----3----2 11963d025dbSVijay Mahadevan EDGE2 EDGE3 12063d025dbSVijay Mahadevan 12163d025dbSVijay Mahadevan Input Parameter: 12263d025dbSVijay Mahadevan 12363d025dbSVijay Mahadevan . PetscInt nverts, the number of element vertices 12463d025dbSVijay Mahadevan . PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering) 12563d025dbSVijay Mahadevan . PetscInt npts, the number of evaluation points (quadrature points) 12663d025dbSVijay Mahadevan . PetscReal quad[3*npts], the evaluation points (quadrature points) in the reference space 12763d025dbSVijay Mahadevan 12863d025dbSVijay Mahadevan Output Parameter: 12963d025dbSVijay Mahadevan . PetscReal phypts[3*npts], the evaluation points (quadrature points) transformed to the physical space 13063d025dbSVijay Mahadevan . PetscReal jxw[npts], the jacobian determinant * quadrature weight necessary for assembling discrete contributions 13163d025dbSVijay Mahadevan . PetscReal phi[npts], the bases evaluated at the specified quadrature points 13263d025dbSVijay Mahadevan . PetscReal dphidx[npts], the derivative of the bases wrt X-direction evaluated at the specified quadrature points 13363d025dbSVijay Mahadevan 13463d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 1-D 13563d025dbSVijay Mahadevan @*/ 13663d025dbSVijay Mahadevan #undef __FUNCT__ 13763d025dbSVijay Mahadevan #define __FUNCT__ "Compute_Lagrange_Basis_1D_Internal" 13863d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_1D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, 139181f196bSVijay Mahadevan PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, 140181f196bSVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) 14163d025dbSVijay Mahadevan { 14263d025dbSVijay Mahadevan int i, j; 14363d025dbSVijay Mahadevan PetscErrorCode ierr; 14463d025dbSVijay Mahadevan 145181f196bSVijay Mahadevan PetscFunctionBegin; 146181f196bSVijay Mahadevan PetscValidPointer(jacobian, 9); 147181f196bSVijay Mahadevan PetscValidPointer(ijacobian, 10); 148181f196bSVijay Mahadevan PetscValidPointer(volume, 11); 149181f196bSVijay Mahadevan if (phypts) { 15063d025dbSVijay Mahadevan ierr = PetscMemzero(phypts, npts * 3 * sizeof(PetscReal));CHKERRQ(ierr); 151181f196bSVijay Mahadevan } 15263d025dbSVijay Mahadevan if (dphidx) { /* Reset arrays. */ 15363d025dbSVijay Mahadevan ierr = PetscMemzero(dphidx, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr); 15463d025dbSVijay Mahadevan } 15563d025dbSVijay Mahadevan if (nverts == 2) { /* Linear Edge */ 15663d025dbSVijay Mahadevan 15763d025dbSVijay Mahadevan for (j = 0; j < npts; j++) 15863d025dbSVijay Mahadevan { 159*a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 160181f196bSVijay Mahadevan const PetscReal r = quad[j]; 16163d025dbSVijay Mahadevan 162181f196bSVijay Mahadevan phi[0 + offset] = ( 1.0 - r ); 163181f196bSVijay Mahadevan phi[1 + offset] = ( r ); 16463d025dbSVijay Mahadevan 165181f196bSVijay Mahadevan const PetscReal dNi_dxi[2] = { -1.0, 1.0 }; 16663d025dbSVijay Mahadevan 167181f196bSVijay Mahadevan jacobian[0] = ijacobian[0] = volume[0] = 0.0; 16863d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 169181f196bSVijay Mahadevan const PetscReal* vertices = coords + i * 3; 170181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 171181f196bSVijay Mahadevan if (phypts) { 172181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 173181f196bSVijay Mahadevan } 17463d025dbSVijay Mahadevan } 17563d025dbSVijay Mahadevan 17663d025dbSVijay Mahadevan /* invert the jacobian */ 177181f196bSVijay Mahadevan *volume = jacobian[0]; 178181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0]; 179181f196bSVijay Mahadevan jxw[j] *= *volume; 18063d025dbSVijay Mahadevan 18163d025dbSVijay Mahadevan /* Divide by element jacobian. */ 18263d025dbSVijay Mahadevan for ( i = 0; i < nverts; i++ ) { 183181f196bSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0]; 18463d025dbSVijay Mahadevan } 18563d025dbSVijay Mahadevan 18663d025dbSVijay Mahadevan } 18763d025dbSVijay Mahadevan } 18863d025dbSVijay Mahadevan else if (nverts == 3) { /* Quadratic Edge */ 18963d025dbSVijay Mahadevan 19063d025dbSVijay Mahadevan for (j = 0; j < npts; j++) 19163d025dbSVijay Mahadevan { 192*a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 193181f196bSVijay Mahadevan const PetscReal r = quad[j]; 19463d025dbSVijay Mahadevan 195181f196bSVijay Mahadevan phi[0 + offset] = 1.0 + r * ( 2.0 * r - 3.0 ); 196181f196bSVijay Mahadevan phi[1 + offset] = 4.0 * r * ( 1.0 - r ); 197181f196bSVijay Mahadevan phi[2 + offset] = r * ( 2.0 * r - 1.0 ); 19863d025dbSVijay Mahadevan 199181f196bSVijay Mahadevan const PetscReal dNi_dxi[3] = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r ), 4.0 * r - 1.0}; 20063d025dbSVijay Mahadevan 201181f196bSVijay Mahadevan jacobian[0] = ijacobian[0] = volume[0] = 0.0; 20263d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 203181f196bSVijay Mahadevan const PetscReal* vertices = coords + i * 3; 204181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 205181f196bSVijay Mahadevan if (phypts) { 206181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 207181f196bSVijay Mahadevan } 20863d025dbSVijay Mahadevan } 20963d025dbSVijay Mahadevan 21063d025dbSVijay Mahadevan /* invert the jacobian */ 211181f196bSVijay Mahadevan *volume = jacobian[0]; 212181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0]; 213181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume; 21463d025dbSVijay Mahadevan 21563d025dbSVijay Mahadevan /* Divide by element jacobian. */ 21663d025dbSVijay Mahadevan for ( i = 0; i < nverts; i++ ) { 217181f196bSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0]; 21863d025dbSVijay Mahadevan } 21963d025dbSVijay Mahadevan } 22063d025dbSVijay Mahadevan } 22163d025dbSVijay Mahadevan else { 22263d025dbSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support EDGE2 and EDGE3 basis evaluations in 1-D : %D", nverts); 22363d025dbSVijay Mahadevan } 22463d025dbSVijay Mahadevan #if 0 22563d025dbSVijay Mahadevan /* verify if the computed basis functions are consistent */ 22663d025dbSVijay Mahadevan for ( j = 0; j < npts; j++ ) { 22763d025dbSVijay Mahadevan PetscScalar phisum = 0, dphixsum = 0; 22863d025dbSVijay Mahadevan for ( i = 0; i < nverts; i++ ) { 22963d025dbSVijay Mahadevan phisum += phi[i + j * nverts]; 23063d025dbSVijay Mahadevan if (dphidx) dphixsum += dphidx[i + j * nverts]; 23163d025dbSVijay Mahadevan } 23263d025dbSVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D = %g, %g, %g\n", j, phisum, dphixsum); 23363d025dbSVijay Mahadevan } 23463d025dbSVijay Mahadevan #endif 23563d025dbSVijay Mahadevan PetscFunctionReturn(0); 23663d025dbSVijay Mahadevan } 23763d025dbSVijay Mahadevan 23863d025dbSVijay Mahadevan 23963d025dbSVijay Mahadevan /*@ 24063d025dbSVijay Mahadevan Compute_Lagrange_Basis_2D_Internal: Evaluate bases and derivatives at quadrature points for a QUAD4 or TRI3 element. 24163d025dbSVijay Mahadevan 24263d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a quadrangle or triangle. 24363d025dbSVijay Mahadevan 24463d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided, 24563d025dbSVijay Mahadevan and their derivatives with respect to X and Y. 24663d025dbSVijay Mahadevan 24763d025dbSVijay Mahadevan Notes: 24863d025dbSVijay Mahadevan 24963d025dbSVijay Mahadevan Example Physical Element (QUAD4) 25063d025dbSVijay Mahadevan 25163d025dbSVijay Mahadevan 4------3 s 25263d025dbSVijay Mahadevan | | | 25363d025dbSVijay Mahadevan | | | 25463d025dbSVijay Mahadevan | | | 25563d025dbSVijay Mahadevan 1------2 0-------r 25663d025dbSVijay Mahadevan 25763d025dbSVijay Mahadevan Input Parameter: 25863d025dbSVijay Mahadevan 25963d025dbSVijay Mahadevan . PetscInt nverts, the number of element vertices 26063d025dbSVijay Mahadevan . PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering) 26163d025dbSVijay Mahadevan . PetscInt npts, the number of evaluation points (quadrature points) 26263d025dbSVijay Mahadevan . PetscReal quad[3*npts], the evaluation points (quadrature points) in the reference space 26363d025dbSVijay Mahadevan 26463d025dbSVijay Mahadevan Output Parameter: 26563d025dbSVijay Mahadevan . PetscReal phypts[3*npts], the evaluation points (quadrature points) transformed to the physical space 26663d025dbSVijay Mahadevan . PetscReal jxw[npts], the jacobian determinant * quadrature weight necessary for assembling discrete contributions 26763d025dbSVijay Mahadevan . PetscReal phi[npts], the bases evaluated at the specified quadrature points 26863d025dbSVijay Mahadevan . PetscReal dphidx[npts], the derivative of the bases wrt X-direction evaluated at the specified quadrature points 26963d025dbSVijay Mahadevan . PetscReal dphidy[npts], the derivative of the bases wrt Y-direction evaluated at the specified quadrature points 27063d025dbSVijay Mahadevan 27163d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 2-D 27263d025dbSVijay Mahadevan @*/ 27363d025dbSVijay Mahadevan #undef __FUNCT__ 27463d025dbSVijay Mahadevan #define __FUNCT__ "Compute_Lagrange_Basis_2D_Internal" 27563d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_2D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, 276181f196bSVijay Mahadevan PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, 277181f196bSVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) 27863d025dbSVijay Mahadevan { 279*a86ed394SVijay Mahadevan PetscInt i, j, k; 28063d025dbSVijay Mahadevan PetscErrorCode ierr; 28163d025dbSVijay Mahadevan 28263d025dbSVijay Mahadevan PetscFunctionBegin; 283181f196bSVijay Mahadevan PetscValidPointer(jacobian, 10); 284181f196bSVijay Mahadevan PetscValidPointer(ijacobian, 11); 285181f196bSVijay Mahadevan PetscValidPointer(volume, 12); 28663d025dbSVijay Mahadevan ierr = PetscMemzero(phi, npts * sizeof(PetscReal));CHKERRQ(ierr); 287181f196bSVijay Mahadevan if (phypts) { 28863d025dbSVijay Mahadevan ierr = PetscMemzero(phypts, npts * 3 * sizeof(PetscReal));CHKERRQ(ierr); 289181f196bSVijay Mahadevan } 29063d025dbSVijay Mahadevan if (dphidx) { /* Reset arrays. */ 29163d025dbSVijay Mahadevan ierr = PetscMemzero(dphidx, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr); 29263d025dbSVijay Mahadevan ierr = PetscMemzero(dphidy, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr); 29363d025dbSVijay Mahadevan } 29463d025dbSVijay Mahadevan if (nverts == 4) { /* Linear Quadrangle */ 29563d025dbSVijay Mahadevan 29663d025dbSVijay Mahadevan for (j = 0; j < npts; j++) 29763d025dbSVijay Mahadevan { 298*a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 299181f196bSVijay Mahadevan const PetscReal r = quad[0 + j * 2]; 300181f196bSVijay Mahadevan const PetscReal s = quad[1 + j * 2]; 30163d025dbSVijay Mahadevan 30263d025dbSVijay Mahadevan phi[0 + offset] = ( 1.0 - r ) * ( 1.0 - s ); 30363d025dbSVijay Mahadevan phi[1 + offset] = r * ( 1.0 - s ); 30463d025dbSVijay Mahadevan phi[2 + offset] = r * s; 30563d025dbSVijay Mahadevan phi[3 + offset] = ( 1.0 - r ) * s; 30663d025dbSVijay Mahadevan 307181f196bSVijay Mahadevan const PetscReal dNi_dxi[4] = { -1.0 + s, 1.0 - s, s, -s }; 308181f196bSVijay Mahadevan const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r }; 30963d025dbSVijay Mahadevan 31063d025dbSVijay Mahadevan ierr = PetscMemzero(jacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr); 31163d025dbSVijay Mahadevan ierr = PetscMemzero(ijacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr); 31263d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 313181f196bSVijay Mahadevan const PetscReal* vertices = coords + i * 3; 31463d025dbSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 31563d025dbSVijay Mahadevan jacobian[2] += dNi_dxi[i] * vertices[1]; 31663d025dbSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertices[0]; 31763d025dbSVijay Mahadevan jacobian[3] += dNi_deta[i] * vertices[1]; 318181f196bSVijay Mahadevan if (phypts) { 319181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 320181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertices[1]; 321181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertices[2]; 322181f196bSVijay Mahadevan } 32363d025dbSVijay Mahadevan } 32463d025dbSVijay Mahadevan 32563d025dbSVijay Mahadevan /* invert the jacobian */ 326181f196bSVijay Mahadevan ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 327181f196bSVijay 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); 32863d025dbSVijay Mahadevan 329181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume; 33063d025dbSVijay Mahadevan 331181f196bSVijay Mahadevan /* Let us compute the bases derivatives by scaling with inverse jacobian. */ 332181f196bSVijay Mahadevan if (dphidx) { 33363d025dbSVijay Mahadevan for ( i = 0; i < nverts; i++ ) { 334*a86ed394SVijay Mahadevan for ( k = 0; k < 2; ++k) { 33563d025dbSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[k * 2 + 0]; 33663d025dbSVijay Mahadevan if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[k * 2 + 1]; 337181f196bSVijay Mahadevan } /* for k=[0..2) */ 338181f196bSVijay Mahadevan } /* for i=[0..nverts) */ 339181f196bSVijay Mahadevan } /* if (dphidx) */ 34063d025dbSVijay Mahadevan } 34163d025dbSVijay Mahadevan } 34263d025dbSVijay Mahadevan else if (nverts == 3) { /* Linear triangle */ 34363d025dbSVijay Mahadevan 34463d025dbSVijay Mahadevan ierr = PetscMemzero(jacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr); 34563d025dbSVijay Mahadevan ierr = PetscMemzero(ijacobian, 4 * sizeof(PetscReal));CHKERRQ(ierr); 34663d025dbSVijay Mahadevan 347181f196bSVijay Mahadevan const PetscReal x2 = coords[2 * 3 + 0], y2 = coords[2 * 3 + 1]; 348181f196bSVijay Mahadevan 34963d025dbSVijay Mahadevan /* Jacobian is constant */ 350181f196bSVijay Mahadevan jacobian[0] = (coords[0 * 3 + 0] - x2); jacobian[1] = (coords[1 * 3 + 0] - x2); 351181f196bSVijay Mahadevan jacobian[2] = (coords[0 * 3 + 1] - y2); jacobian[3] = (coords[1 * 3 + 1] - y2); 35263d025dbSVijay Mahadevan 35363d025dbSVijay Mahadevan /* invert the jacobian */ 354181f196bSVijay Mahadevan ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 355181f196bSVijay Mahadevan if ( *volume < 1e-12 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangular element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume); 356181f196bSVijay Mahadevan 357181f196bSVijay Mahadevan const PetscReal Dx[3] = { ijacobian[0], ijacobian[2], - ijacobian[0] - ijacobian[2] }; 358181f196bSVijay Mahadevan const PetscReal Dy[3] = { ijacobian[1], ijacobian[3], - ijacobian[1] - ijacobian[3] }; 35963d025dbSVijay Mahadevan 36063d025dbSVijay Mahadevan for (j = 0; j < npts; j++) 36163d025dbSVijay Mahadevan { 362*a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 363181f196bSVijay Mahadevan const PetscReal r = quad[0 + j * 2]; 364181f196bSVijay Mahadevan const PetscReal s = quad[1 + j * 2]; 36563d025dbSVijay Mahadevan 366181f196bSVijay Mahadevan if (jxw) jxw[j] *= 0.5; 36763d025dbSVijay Mahadevan 368181f196bSVijay Mahadevan /* const PetscReal Ni[3] = { r, s, 1.0 - r - s }; */ 369181f196bSVijay Mahadevan const PetscReal phipts_x = coords[2 * 3 + 0] + jacobian[0] * r + jacobian[1] * s; 370181f196bSVijay Mahadevan const PetscReal phipts_y = coords[2 * 3 + 1] + jacobian[2] * r + jacobian[3] * s; 371181f196bSVijay Mahadevan if (phypts) { 372181f196bSVijay Mahadevan phypts[offset + 0] = phipts_x; 373181f196bSVijay Mahadevan phypts[offset + 1] = phipts_y; 374181f196bSVijay Mahadevan } 37563d025dbSVijay Mahadevan 376181f196bSVijay Mahadevan /* \phi_0 = (b.y − c.y) x + (b.x − c.x) y + c.x b.y − b.x c.y */ 377181f196bSVijay Mahadevan phi[0 + offset] = ( ijacobian[0] * (phipts_x - x2) + ijacobian[1] * (phipts_y - y2) ); 378181f196bSVijay Mahadevan /* \phi_1 = (c.y − a.y) x + (a.x − c.x) y + c.x a.y − a.x c.y */ 379181f196bSVijay Mahadevan phi[1 + offset] = ( ijacobian[2] * (phipts_x - x2) + ijacobian[3] * (phipts_y - y2) ); 38063d025dbSVijay Mahadevan phi[2 + offset] = 1.0 - phi[0 + offset] - phi[1 + offset]; 38163d025dbSVijay Mahadevan 38263d025dbSVijay Mahadevan if (dphidx) { 383181f196bSVijay Mahadevan dphidx[0 + offset] = Dx[0]; 384181f196bSVijay Mahadevan dphidx[1 + offset] = Dx[1]; 385181f196bSVijay Mahadevan dphidx[2 + offset] = Dx[2]; 38663d025dbSVijay Mahadevan } 38763d025dbSVijay Mahadevan 38863d025dbSVijay Mahadevan if (dphidy) { 389181f196bSVijay Mahadevan dphidy[0 + offset] = Dy[0]; 390181f196bSVijay Mahadevan dphidy[1 + offset] = Dy[1]; 391181f196bSVijay Mahadevan dphidy[2 + offset] = Dy[2]; 39263d025dbSVijay Mahadevan } 39363d025dbSVijay Mahadevan 39463d025dbSVijay Mahadevan } 39563d025dbSVijay Mahadevan } 39663d025dbSVijay Mahadevan else { 39763d025dbSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support QUAD4 and TRI3 basis evaluations in 2-D : %D", nverts); 39863d025dbSVijay Mahadevan } 39963d025dbSVijay Mahadevan #if 0 40063d025dbSVijay Mahadevan /* verify if the computed basis functions are consistent */ 40163d025dbSVijay Mahadevan for ( j = 0; j < npts; j++ ) { 40263d025dbSVijay Mahadevan PetscScalar phisum = 0, dphixsum = 0, dphiysum = 0; 40363d025dbSVijay Mahadevan for ( i = 0; i < nverts; i++ ) { 40463d025dbSVijay Mahadevan phisum += phi[i + j * nverts]; 40563d025dbSVijay Mahadevan if (dphidx) dphixsum += dphidx[i + j * nverts]; 40663d025dbSVijay Mahadevan if (dphidy) dphiysum += dphidy[i + j * nverts]; 40763d025dbSVijay Mahadevan } 40863d025dbSVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D = %g, %g, %g\n", j, phisum, dphixsum, dphiysum); 40963d025dbSVijay Mahadevan } 41063d025dbSVijay Mahadevan #endif 41163d025dbSVijay Mahadevan PetscFunctionReturn(0); 41263d025dbSVijay Mahadevan } 41363d025dbSVijay Mahadevan 41463d025dbSVijay Mahadevan 41563d025dbSVijay Mahadevan /*@ 41663d025dbSVijay Mahadevan Compute_Lagrange_Basis_3D_Internal: Evaluate bases and derivatives at quadrature points for a HEX8 or TET4 element. 41763d025dbSVijay Mahadevan 41863d025dbSVijay Mahadevan The routine is given the coordinates of the vertices of a hexahedra or tetrahedra. 41963d025dbSVijay Mahadevan 42063d025dbSVijay Mahadevan The routine evaluates the basis functions associated with each quadrature point provided, 42163d025dbSVijay Mahadevan and their derivatives with respect to X, Y and Z. 42263d025dbSVijay Mahadevan 42363d025dbSVijay Mahadevan Notes: 42463d025dbSVijay Mahadevan 42563d025dbSVijay Mahadevan Example Physical Element (HEX8) 42663d025dbSVijay Mahadevan 42763d025dbSVijay Mahadevan 8------7 42863d025dbSVijay Mahadevan /| /| t s 42963d025dbSVijay Mahadevan 5------6 | | / 43063d025dbSVijay Mahadevan | | | | |/ 43163d025dbSVijay Mahadevan | 4----|-3 0-------r 43263d025dbSVijay Mahadevan |/ |/ 43363d025dbSVijay Mahadevan 1------2 43463d025dbSVijay Mahadevan 43563d025dbSVijay Mahadevan Input Parameter: 43663d025dbSVijay Mahadevan 43763d025dbSVijay Mahadevan . PetscInt nverts, the number of element vertices 43863d025dbSVijay Mahadevan . PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering) 43963d025dbSVijay Mahadevan . PetscInt npts, the number of evaluation points (quadrature points) 44063d025dbSVijay Mahadevan . PetscReal quad[3*npts], the evaluation points (quadrature points) in the reference space 44163d025dbSVijay Mahadevan 44263d025dbSVijay Mahadevan Output Parameter: 44363d025dbSVijay Mahadevan . PetscReal phypts[3*npts], the evaluation points (quadrature points) transformed to the physical space 44463d025dbSVijay Mahadevan . PetscReal jxw[npts], the jacobian determinant * quadrature weight necessary for assembling discrete contributions 44563d025dbSVijay Mahadevan . PetscReal phi[npts], the bases evaluated at the specified quadrature points 44663d025dbSVijay Mahadevan . PetscReal dphidx[npts], the derivative of the bases wrt X-direction evaluated at the specified quadrature points 44763d025dbSVijay Mahadevan . PetscReal dphidy[npts], the derivative of the bases wrt Y-direction evaluated at the specified quadrature points 44863d025dbSVijay Mahadevan . PetscReal dphidz[npts], the derivative of the bases wrt Z-direction evaluated at the specified quadrature points 44963d025dbSVijay Mahadevan 45063d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 3-D 45163d025dbSVijay Mahadevan @*/ 45263d025dbSVijay Mahadevan #undef __FUNCT__ 45363d025dbSVijay Mahadevan #define __FUNCT__ "Compute_Lagrange_Basis_3D_Internal" 45463d025dbSVijay Mahadevan PetscErrorCode Compute_Lagrange_Basis_3D_Internal ( const PetscInt nverts, const PetscReal *coords, const PetscInt npts, const PetscReal *quad, PetscReal *phypts, 455181f196bSVijay Mahadevan PetscReal *jxw, PetscReal *phi, PetscReal *dphidx, PetscReal *dphidy, PetscReal *dphidz, 456181f196bSVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal *volume) 45763d025dbSVijay Mahadevan { 458*a86ed394SVijay Mahadevan PetscInt i, j, k; 45963d025dbSVijay Mahadevan PetscErrorCode ierr; 46063d025dbSVijay Mahadevan 46163d025dbSVijay Mahadevan PetscFunctionBegin; 462181f196bSVijay Mahadevan PetscValidPointer(jacobian, 11); 463181f196bSVijay Mahadevan PetscValidPointer(ijacobian, 12); 464181f196bSVijay Mahadevan PetscValidPointer(volume, 13); 46563d025dbSVijay Mahadevan /* Reset arrays. */ 46663d025dbSVijay Mahadevan ierr = PetscMemzero(phi, npts * sizeof(PetscReal));CHKERRQ(ierr); 467181f196bSVijay Mahadevan if (phypts) { 46863d025dbSVijay Mahadevan ierr = PetscMemzero(phypts, npts * 3 * sizeof(PetscReal));CHKERRQ(ierr); 469181f196bSVijay Mahadevan } 47063d025dbSVijay Mahadevan if (dphidx) { 47163d025dbSVijay Mahadevan ierr = PetscMemzero(dphidx, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr); 47263d025dbSVijay Mahadevan ierr = PetscMemzero(dphidy, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr); 47363d025dbSVijay Mahadevan ierr = PetscMemzero(dphidz, npts * nverts * sizeof(PetscReal));CHKERRQ(ierr); 47463d025dbSVijay Mahadevan } 47563d025dbSVijay Mahadevan 47663d025dbSVijay Mahadevan if (nverts == 8) { /* Linear Hexahedra */ 47763d025dbSVijay Mahadevan 47863d025dbSVijay Mahadevan for (j = 0; j < npts; j++) 47963d025dbSVijay Mahadevan { 480*a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 481181f196bSVijay Mahadevan const PetscReal& r = quad[j * 3 + 0]; 482181f196bSVijay Mahadevan const PetscReal& s = quad[j * 3 + 1]; 483181f196bSVijay Mahadevan const PetscReal& t = quad[j * 3 + 2]; 48463d025dbSVijay Mahadevan 485*a86ed394SVijay Mahadevan phi[offset + 0] = ( 1.0 - r ) * ( 1.0 - s ) * ( 1.0 - t ); /* 0,0,0 */ 486*a86ed394SVijay Mahadevan phi[offset + 1] = ( r ) * ( 1.0 - s ) * ( 1.0 - t ); /* 1,0,0 */ 487*a86ed394SVijay Mahadevan phi[offset + 2] = ( r ) * ( s ) * ( 1.0 - t ); /* 1,1,0 */ 488*a86ed394SVijay Mahadevan phi[offset + 3] = ( 1.0 - r ) * ( s ) * ( 1.0 - t ); /* 0,1,0 */ 489*a86ed394SVijay Mahadevan phi[offset + 4] = ( 1.0 - r ) * ( 1.0 - s ) * ( t ); /* 0,0,1 */ 490*a86ed394SVijay Mahadevan phi[offset + 5] = ( r ) * ( 1.0 - s ) * ( t ); /* 1,0,1 */ 491*a86ed394SVijay Mahadevan phi[offset + 6] = ( r ) * ( s ) * ( t ); /* 1,1,1 */ 492*a86ed394SVijay Mahadevan phi[offset + 7] = ( 1.0 - r ) * ( s ) * ( t ); /* 0,1,1 */ 49363d025dbSVijay Mahadevan 494181f196bSVijay Mahadevan const PetscReal dNi_dxi[8] = {- ( 1.0 - s ) * ( 1.0 - t ), 49563d025dbSVijay Mahadevan ( 1.0 - s ) * ( 1.0 - t ), 496*a86ed394SVijay Mahadevan ( s ) * ( 1.0 - t ), 497*a86ed394SVijay Mahadevan - ( s ) * ( 1.0 - t ), 498*a86ed394SVijay Mahadevan - ( 1.0 - s ) * ( t ), 499*a86ed394SVijay Mahadevan ( 1.0 - s ) * ( t ), 500*a86ed394SVijay Mahadevan ( s ) * ( t ), 501*a86ed394SVijay Mahadevan - ( s ) * ( t ) 50263d025dbSVijay Mahadevan }; 50363d025dbSVijay Mahadevan 504181f196bSVijay Mahadevan const PetscReal dNi_deta[8] = { - ( 1.0 - r ) * ( 1.0 - t ), 505*a86ed394SVijay Mahadevan - ( r ) * ( 1.0 - t ), 506*a86ed394SVijay Mahadevan ( r ) * ( 1.0 - t ), 50763d025dbSVijay Mahadevan ( 1.0 - r ) * ( 1.0 - t ), 508*a86ed394SVijay Mahadevan - ( 1.0 - r ) * ( t ), 509*a86ed394SVijay Mahadevan - ( r ) * ( t ), 510*a86ed394SVijay Mahadevan ( r ) * ( t ), 511*a86ed394SVijay Mahadevan ( 1.0 - r ) * ( t ) 51263d025dbSVijay Mahadevan }; 51363d025dbSVijay Mahadevan 514181f196bSVijay Mahadevan const PetscReal dNi_dzeta[8] = { - ( 1.0 - r ) * ( 1.0 - s ), 515*a86ed394SVijay Mahadevan - ( r ) * ( 1.0 - s ), 516*a86ed394SVijay Mahadevan - ( r ) * ( s ), 517*a86ed394SVijay Mahadevan - ( 1.0 - r ) * ( s ), 51863d025dbSVijay Mahadevan ( 1.0 - r ) * ( 1.0 - s ), 519*a86ed394SVijay Mahadevan ( r ) * ( 1.0 - s ), 520*a86ed394SVijay Mahadevan ( r ) * ( s ), 521*a86ed394SVijay Mahadevan ( 1.0 - r ) * ( s ) 52263d025dbSVijay Mahadevan }; 52363d025dbSVijay Mahadevan 52463d025dbSVijay Mahadevan ierr = PetscMemzero(jacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr); 52563d025dbSVijay Mahadevan ierr = PetscMemzero(ijacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr); 52663d025dbSVijay Mahadevan for (i = 0; i < nverts; ++i) { 527181f196bSVijay Mahadevan const PetscReal* vertex = coords + i * 3; 52863d025dbSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertex[0]; 52963d025dbSVijay Mahadevan jacobian[3] += dNi_dxi[i] * vertex[1]; 53063d025dbSVijay Mahadevan jacobian[6] += dNi_dxi[i] * vertex[2]; 53163d025dbSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertex[0]; 53263d025dbSVijay Mahadevan jacobian[4] += dNi_deta[i] * vertex[1]; 53363d025dbSVijay Mahadevan jacobian[7] += dNi_deta[i] * vertex[2]; 53463d025dbSVijay Mahadevan jacobian[2] += dNi_dzeta[i] * vertex[0]; 53563d025dbSVijay Mahadevan jacobian[5] += dNi_dzeta[i] * vertex[1]; 53663d025dbSVijay Mahadevan jacobian[8] += dNi_dzeta[i] * vertex[2]; 537181f196bSVijay Mahadevan if (phypts) { 538181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertex[0]; 539181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertex[1]; 540181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertex[2]; 541181f196bSVijay Mahadevan } 54263d025dbSVijay Mahadevan } 54363d025dbSVijay Mahadevan 54463d025dbSVijay Mahadevan /* invert the jacobian */ 545181f196bSVijay Mahadevan ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 546*a86ed394SVijay 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); 54763d025dbSVijay Mahadevan 548*a86ed394SVijay Mahadevan if (jxw) jxw[j] *= (*volume); 54963d025dbSVijay Mahadevan 55063d025dbSVijay Mahadevan /* Divide by element jacobian. */ 55163d025dbSVijay Mahadevan for ( i = 0; i < nverts; ++i ) { 552*a86ed394SVijay Mahadevan for ( k = 0; k < 3; ++k) { 55363d025dbSVijay Mahadevan if (dphidx) dphidx[i + offset] += dNi_dxi[i] * ijacobian[0 * 3 + k]; 55463d025dbSVijay Mahadevan if (dphidy) dphidy[i + offset] += dNi_deta[i] * ijacobian[1 * 3 + k]; 55563d025dbSVijay Mahadevan if (dphidz) dphidz[i + offset] += dNi_dzeta[i] * ijacobian[2 * 3 + k]; 55663d025dbSVijay Mahadevan } 55763d025dbSVijay Mahadevan } 55863d025dbSVijay Mahadevan } 55963d025dbSVijay Mahadevan } 56063d025dbSVijay Mahadevan else if (nverts == 4) { /* Linear Tetrahedra */ 56163d025dbSVijay Mahadevan 56263d025dbSVijay Mahadevan ierr = PetscMemzero(jacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr); 56363d025dbSVijay Mahadevan ierr = PetscMemzero(ijacobian, 9 * sizeof(PetscReal));CHKERRQ(ierr); 56463d025dbSVijay Mahadevan 565181f196bSVijay Mahadevan const PetscReal x0 = coords[/*0 * 3 +*/ 0], y0 = coords[/*0 * 3 +*/ 1], z0 = coords[/*0 * 3 +*/ 2]; 566181f196bSVijay Mahadevan 567181f196bSVijay Mahadevan /* compute the jacobian */ 568181f196bSVijay Mahadevan jacobian[0] = coords[1 * 3 + 0] - x0; jacobian[1] = coords[2 * 3 + 0] - x0; jacobian[2] = coords[3 * 3 + 0] - x0; 569181f196bSVijay Mahadevan jacobian[3] = coords[1 * 3 + 1] - y0; jacobian[4] = coords[2 * 3 + 1] - y0; jacobian[5] = coords[3 * 3 + 1] - y0; 570181f196bSVijay Mahadevan jacobian[6] = coords[1 * 3 + 2] - z0; jacobian[7] = coords[2 * 3 + 2] - z0; jacobian[8] = coords[3 * 3 + 2] - z0; 57163d025dbSVijay Mahadevan 57263d025dbSVijay Mahadevan /* invert the jacobian */ 573181f196bSVijay Mahadevan ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, volume);CHKERRQ(ierr); 574*a86ed394SVijay Mahadevan if ( *volume < 1e-8 ) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral element has zero volume: %g. Degenerate element or invalid connectivity\n", *volume); 57563d025dbSVijay Mahadevan 576181f196bSVijay Mahadevan PetscReal Dx[4], Dy[4], Dz[4]; 577181f196bSVijay Mahadevan if (dphidx) { 578181f196bSVijay Mahadevan Dx[0] = ( coords[1 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] ) 579181f196bSVijay Mahadevan - coords[1 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] ) 580181f196bSVijay Mahadevan - coords[1 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3] ) 581181f196bSVijay Mahadevan ) / *volume; 582181f196bSVijay Mahadevan Dx[1] = - ( coords[1 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] ) 583181f196bSVijay Mahadevan - coords[1 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] ) 584181f196bSVijay Mahadevan - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3] ) 585181f196bSVijay Mahadevan ) / *volume; 586181f196bSVijay Mahadevan Dx[2] = ( coords[1 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] ) 587181f196bSVijay Mahadevan - coords[1 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] ) 588181f196bSVijay Mahadevan - coords[1 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3] ) 589181f196bSVijay Mahadevan ) / *volume; 590181f196bSVijay Mahadevan Dx[3] = - ( Dx[0] + Dx[1] + Dx[2] ); 591181f196bSVijay Mahadevan Dy[0] = ( coords[0 + 1 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] ) 592181f196bSVijay Mahadevan - coords[0 + 2 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] ) 593181f196bSVijay Mahadevan + coords[0 + 3 * 3] * ( coords[2 + 1 * 3] - coords[2 + 2 * 3] ) 594181f196bSVijay Mahadevan ) / *volume; 595181f196bSVijay Mahadevan Dy[1] = - ( coords[0 + 0 * 3] * ( coords[2 + 2 * 3] - coords[2 + 3 * 3] ) 596181f196bSVijay Mahadevan - coords[0 + 2 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] ) 597181f196bSVijay Mahadevan + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 2 * 3] ) 598181f196bSVijay Mahadevan ) / *volume; 599181f196bSVijay Mahadevan Dy[2] = ( coords[0 + 0 * 3] * ( coords[2 + 1 * 3] - coords[2 + 3 * 3] ) 600181f196bSVijay Mahadevan - coords[0 + 1 * 3] * ( coords[2 + 0 * 3] - coords[2 + 3 * 3] ) 601181f196bSVijay Mahadevan + coords[0 + 3 * 3] * ( coords[2 + 0 * 3] - coords[2 + 1 * 3] ) 602181f196bSVijay Mahadevan ) / *volume; 603181f196bSVijay Mahadevan Dy[3] = - ( Dy[0] + Dy[1] + Dy[2] ); 604181f196bSVijay Mahadevan Dz[0] = ( coords[0 + 1 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3] ) 605181f196bSVijay Mahadevan - coords[0 + 2 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3] ) 606181f196bSVijay Mahadevan + coords[0 + 3 * 3] * (coords[1 + 2 * 3] - coords[1 + 1 * 3] ) 607181f196bSVijay Mahadevan ) / *volume; 608181f196bSVijay Mahadevan Dz[1] = - ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 2 * 3] ) 609181f196bSVijay Mahadevan + coords[0 + 2 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3] ) 610181f196bSVijay Mahadevan - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 2 * 3] ) 611181f196bSVijay Mahadevan ) / *volume; 612181f196bSVijay Mahadevan Dz[2] = ( coords[0 + 0 * 3] * (coords[1 + 3 * 3] - coords[1 + 1 * 3] ) 613181f196bSVijay Mahadevan + coords[0 + 1 * 3] * (coords[1 + 0 * 3] - coords[1 + 3 * 3] ) 614181f196bSVijay Mahadevan - coords[0 + 3 * 3] * (coords[1 + 0 * 3] - coords[1 + 1 * 3] ) 615181f196bSVijay Mahadevan ) / *volume; 616181f196bSVijay Mahadevan Dz[3] = - ( Dz[0] + Dz[1] + Dz[2] ); 617181f196bSVijay Mahadevan } 61863d025dbSVijay Mahadevan 61963d025dbSVijay Mahadevan for ( j = 0; j < npts; j++ ) 62063d025dbSVijay Mahadevan { 621*a86ed394SVijay Mahadevan const PetscInt offset = j * nverts; 622181f196bSVijay Mahadevan const PetscReal& r = quad[j * 3 + 0]; 623181f196bSVijay Mahadevan const PetscReal& s = quad[j * 3 + 1]; 624181f196bSVijay Mahadevan const PetscReal& t = quad[j * 3 + 2]; 62563d025dbSVijay Mahadevan 626181f196bSVijay Mahadevan if (jxw) jxw[j] *= *volume; 62763d025dbSVijay Mahadevan 62863d025dbSVijay Mahadevan phi[offset + 0] = 1.0 - r - s - t; 62963d025dbSVijay Mahadevan phi[offset + 1] = r; 63063d025dbSVijay Mahadevan phi[offset + 2] = s; 63163d025dbSVijay Mahadevan phi[offset + 3] = t; 63263d025dbSVijay Mahadevan 633181f196bSVijay Mahadevan if (phypts) { 634181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 635181f196bSVijay Mahadevan const PetscScalar* vertices = coords + i * 3; 636181f196bSVijay Mahadevan phypts[3 * j + 0] += phi[i + offset] * vertices[0]; 637181f196bSVijay Mahadevan phypts[3 * j + 1] += phi[i + offset] * vertices[1]; 638181f196bSVijay Mahadevan phypts[3 * j + 2] += phi[i + offset] * vertices[2]; 639181f196bSVijay Mahadevan } 640181f196bSVijay Mahadevan } 641181f196bSVijay Mahadevan 642181f196bSVijay Mahadevan /* Now set the derivatives */ 64363d025dbSVijay Mahadevan if (dphidx) { 644181f196bSVijay Mahadevan dphidx[0 + offset] = Dx[0]; 645181f196bSVijay Mahadevan dphidx[1 + offset] = Dx[1]; 646181f196bSVijay Mahadevan dphidx[2 + offset] = Dx[2]; 647181f196bSVijay Mahadevan dphidx[3 + offset] = Dx[3]; 64863d025dbSVijay Mahadevan } 64963d025dbSVijay Mahadevan 65063d025dbSVijay Mahadevan if (dphidy) { 651181f196bSVijay Mahadevan dphidy[0 + offset] = Dy[0]; 652181f196bSVijay Mahadevan dphidy[1 + offset] = Dy[1]; 653181f196bSVijay Mahadevan dphidy[2 + offset] = Dy[2]; 654181f196bSVijay Mahadevan dphidy[3 + offset] = Dy[3]; 65563d025dbSVijay Mahadevan } 65663d025dbSVijay Mahadevan 65763d025dbSVijay Mahadevan if (dphidz) { 658181f196bSVijay Mahadevan dphidz[0 + offset] = Dz[0]; 659181f196bSVijay Mahadevan dphidz[1 + offset] = Dz[1]; 660181f196bSVijay Mahadevan dphidz[2 + offset] = Dz[2]; 661181f196bSVijay Mahadevan dphidz[3 + offset] = Dz[3]; 66263d025dbSVijay Mahadevan } 66363d025dbSVijay Mahadevan 66463d025dbSVijay Mahadevan } /* Tetrahedra -- ends */ 66563d025dbSVijay Mahadevan } 66663d025dbSVijay Mahadevan else 66763d025dbSVijay Mahadevan { 66863d025dbSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of entity vertices are invalid. Currently only support HEX8 and TET4 basis evaluations in 3-D : %D", nverts); 66963d025dbSVijay Mahadevan } 67063d025dbSVijay Mahadevan #if 0 67163d025dbSVijay Mahadevan /* verify if the computed basis functions are consistent */ 67263d025dbSVijay Mahadevan for ( j = 0; j < npts; j++ ) { 67363d025dbSVijay Mahadevan PetscScalar phisum = 0, dphixsum = 0, dphiysum = 0, dphizsum = 0; 67463d025dbSVijay Mahadevan const int offset = j * nverts; 67563d025dbSVijay Mahadevan for ( i = 0; i < nverts; i++ ) { 67663d025dbSVijay Mahadevan phisum += phi[i + offset]; 67763d025dbSVijay Mahadevan if (dphidx) dphixsum += dphidx[i + offset]; 67863d025dbSVijay Mahadevan if (dphidy) dphiysum += dphidy[i + offset]; 67963d025dbSVijay Mahadevan if (dphidz) dphizsum += dphidz[i + offset]; 68063d025dbSVijay Mahadevan if (dphidx) PetscPrintf(PETSC_COMM_WORLD, "\t Values [%d]: [JxW] [phi, dphidx, dphidy, dphidz] = %g, %g, %g, %g, %g\n", j, jxw[j], phi[i + offset], dphidx[i + offset], dphidy[i + offset], dphidz[i + offset]); 68163d025dbSVijay Mahadevan } 68263d025dbSVijay Mahadevan if (dphidx) PetscPrintf(PETSC_COMM_WORLD, "Sum of basis at quadrature point %D (%g, %g, %g) = %g, %g, %g, %g\n", j, quad[3 * j + 0], quad[3 * j + 1], quad[3 * j + 2], phisum, dphixsum, dphiysum, dphizsum); 68363d025dbSVijay Mahadevan } 68463d025dbSVijay Mahadevan #endif 68563d025dbSVijay Mahadevan PetscFunctionReturn(0); 68663d025dbSVijay Mahadevan } 68763d025dbSVijay Mahadevan 68863d025dbSVijay Mahadevan 68963d025dbSVijay Mahadevan 69063d025dbSVijay Mahadevan /*@ 69163d025dbSVijay Mahadevan DMMoabFEMComputeBasis: Evaluate bases and derivatives at quadrature points for a linear EDGE/QUAD/TRI/HEX/TET element. 69263d025dbSVijay Mahadevan 69363d025dbSVijay Mahadevan The routine takes the coordinates of the vertices of an element and computes basis functions associated with 69463d025dbSVijay Mahadevan each quadrature point provided, and their derivatives with respect to X, Y and Z as appropriate. 69563d025dbSVijay Mahadevan 69663d025dbSVijay Mahadevan Input Parameter: 69763d025dbSVijay Mahadevan 69863d025dbSVijay Mahadevan . PetscInt nverts, the number of element vertices 69963d025dbSVijay Mahadevan . PetscReal coords[3*nverts], the physical coordinates of the vertices (in canonical numbering) 70063d025dbSVijay Mahadevan . PetscInt npts, the number of evaluation points (quadrature points) 70163d025dbSVijay Mahadevan . PetscReal quad[3*npts], the evaluation points (quadrature points) in the reference space 70263d025dbSVijay Mahadevan 70363d025dbSVijay Mahadevan Output Parameter: 70463d025dbSVijay Mahadevan . PetscReal phypts[3*npts], the evaluation points (quadrature points) transformed to the physical space 70563d025dbSVijay Mahadevan . PetscReal jxw[npts], the jacobian determinant * quadrature weight necessary for assembling discrete contributions 70663d025dbSVijay Mahadevan . PetscReal fe_basis[npts], the bases values evaluated at the specified quadrature points 70763d025dbSVijay Mahadevan . 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 70863d025dbSVijay Mahadevan 70963d025dbSVijay Mahadevan .keywords: DMMoab, FEM, 3-D 71063d025dbSVijay Mahadevan @*/ 71163d025dbSVijay Mahadevan #undef __FUNCT__ 71263d025dbSVijay Mahadevan #define __FUNCT__ "DMMoabFEMComputeBasis" 713181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMComputeBasis ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscQuadrature quadrature, 714181f196bSVijay Mahadevan PetscReal *phypts, PetscReal *jacobian_quadrature_weight_product, 715181f196bSVijay Mahadevan PetscReal *fe_basis, PetscReal **fe_basis_derivatives) 71663d025dbSVijay Mahadevan { 71763d025dbSVijay Mahadevan PetscErrorCode ierr; 71863d025dbSVijay Mahadevan PetscInt npoints; 71963d025dbSVijay Mahadevan bool compute_der; 72063d025dbSVijay Mahadevan const PetscReal *quadpts, *quadwts; 721181f196bSVijay Mahadevan PetscReal jacobian[9], ijacobian[9], volume; 72263d025dbSVijay Mahadevan 72363d025dbSVijay Mahadevan PetscFunctionBegin; 72463d025dbSVijay Mahadevan PetscValidPointer(coordinates, 3); 72563d025dbSVijay Mahadevan PetscValidHeaderSpecific(quadrature, PETSC_OBJECT_CLASSID, 4); 72663d025dbSVijay Mahadevan PetscValidPointer(fe_basis, 7); 72763d025dbSVijay Mahadevan compute_der = (fe_basis_derivatives != NULL); 72863d025dbSVijay Mahadevan 72963d025dbSVijay Mahadevan /* Get the quadrature points and weights for the given quadrature rule */ 73063d025dbSVijay Mahadevan ierr = PetscQuadratureGetData(quadrature, NULL, &npoints, &quadpts, &quadwts);CHKERRQ(ierr); 731181f196bSVijay Mahadevan if (jacobian_quadrature_weight_product) { 73263d025dbSVijay Mahadevan ierr = PetscMemcpy(jacobian_quadrature_weight_product, quadwts, npoints * sizeof(PetscReal));CHKERRQ(ierr); 733181f196bSVijay Mahadevan } 73463d025dbSVijay Mahadevan 73563d025dbSVijay Mahadevan switch (dim) { 73663d025dbSVijay Mahadevan case 1: 73763d025dbSVijay Mahadevan ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, npoints, quadpts, phypts, 73863d025dbSVijay Mahadevan jacobian_quadrature_weight_product, fe_basis, 739181f196bSVijay Mahadevan (compute_der ? fe_basis_derivatives[0] : NULL), 740181f196bSVijay Mahadevan jacobian, ijacobian, &volume);CHKERRQ(ierr); 74163d025dbSVijay Mahadevan break; 74263d025dbSVijay Mahadevan case 2: 74363d025dbSVijay Mahadevan ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, npoints, quadpts, phypts, 74463d025dbSVijay Mahadevan jacobian_quadrature_weight_product, fe_basis, 74563d025dbSVijay Mahadevan (compute_der ? fe_basis_derivatives[0] : NULL), 746181f196bSVijay Mahadevan (compute_der ? fe_basis_derivatives[1] : NULL), 747181f196bSVijay Mahadevan jacobian, ijacobian, &volume);CHKERRQ(ierr); 74863d025dbSVijay Mahadevan break; 74963d025dbSVijay Mahadevan case 3: 75063d025dbSVijay Mahadevan ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, npoints, quadpts, phypts, 75163d025dbSVijay Mahadevan jacobian_quadrature_weight_product, fe_basis, 75263d025dbSVijay Mahadevan (compute_der ? fe_basis_derivatives[0] : NULL), 75363d025dbSVijay Mahadevan (compute_der ? fe_basis_derivatives[1] : NULL), 754181f196bSVijay Mahadevan (compute_der ? fe_basis_derivatives[2] : NULL), 755181f196bSVijay Mahadevan jacobian, ijacobian, &volume);CHKERRQ(ierr); 75663d025dbSVijay Mahadevan break; 75763d025dbSVijay Mahadevan default: 75863d025dbSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim); 75963d025dbSVijay Mahadevan } 76063d025dbSVijay Mahadevan PetscFunctionReturn(0); 76163d025dbSVijay Mahadevan } 76263d025dbSVijay Mahadevan 76363d025dbSVijay Mahadevan 76463d025dbSVijay Mahadevan 76563d025dbSVijay Mahadevan /*@ 76663d025dbSVijay Mahadevan DMMoabFEMCreateQuadratureDefault: Create default quadrature rules for integration over an element with a given 76763d025dbSVijay Mahadevan dimension and polynomial order (deciphered from number of element vertices). 76863d025dbSVijay Mahadevan 76963d025dbSVijay Mahadevan Input Parameter: 77063d025dbSVijay Mahadevan 77163d025dbSVijay Mahadevan . PetscInt dim, the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET) 77263d025dbSVijay Mahadevan . PetscInt nverts, the number of vertices in the physical element 77363d025dbSVijay Mahadevan 77463d025dbSVijay Mahadevan Output Parameter: 77563d025dbSVijay Mahadevan . PetscQuadrature quadrature, the quadrature object with default settings to integrate polynomials defined over the element 77663d025dbSVijay Mahadevan 77763d025dbSVijay Mahadevan .keywords: DMMoab, Quadrature, PetscDT 77863d025dbSVijay Mahadevan @*/ 77963d025dbSVijay Mahadevan #undef __FUNCT__ 78063d025dbSVijay Mahadevan #define __FUNCT__ "DMMoabFEMCreateQuadratureDefault" 781181f196bSVijay Mahadevan PetscErrorCode DMMoabFEMCreateQuadratureDefault ( const PetscInt dim, const PetscInt nverts, PetscQuadrature *quadrature ) 78263d025dbSVijay Mahadevan { 78363d025dbSVijay Mahadevan PetscReal *w, *x; 78463d025dbSVijay Mahadevan PetscErrorCode ierr; 78563d025dbSVijay Mahadevan 78663d025dbSVijay Mahadevan PetscFunctionBegin; 78763d025dbSVijay Mahadevan /* Create an appropriate quadrature rule to sample basis */ 78863d025dbSVijay Mahadevan switch (dim) 78963d025dbSVijay Mahadevan { 79063d025dbSVijay Mahadevan case 1: 79163d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */ 792181f196bSVijay Mahadevan ierr = PetscDTGaussJacobiQuadrature(1, nverts, 0, 1.0, quadrature);CHKERRQ(ierr); 79363d025dbSVijay Mahadevan break; 79463d025dbSVijay Mahadevan case 2: 79563d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */ 79663d025dbSVijay Mahadevan if (nverts == 3) { 797*a86ed394SVijay Mahadevan const PetscInt order = 2; 798*a86ed394SVijay Mahadevan const PetscInt npoints = (order == 2 ? 3 : 6); 79963d025dbSVijay Mahadevan ierr = PetscMalloc2(npoints * 2, &x, npoints, &w);CHKERRQ(ierr); 800181f196bSVijay Mahadevan if (npoints == 3) { 80163d025dbSVijay Mahadevan x[0] = x[1] = x[2] = x[5] = 1.0 / 6.0; 80263d025dbSVijay Mahadevan x[3] = x[4] = 2.0 / 3.0; 80363d025dbSVijay Mahadevan w[0] = w[1] = w[2] = 1.0 / 3.0; 804181f196bSVijay Mahadevan } 805181f196bSVijay Mahadevan else if (npoints == 6) { 80663d025dbSVijay Mahadevan x[0] = x[1] = x[2] = 0.44594849091597; 80763d025dbSVijay Mahadevan x[3] = x[4] = 0.10810301816807; 80863d025dbSVijay Mahadevan x[5] = 0.44594849091597; 80963d025dbSVijay Mahadevan x[6] = x[7] = x[8] = 0.09157621350977; 81063d025dbSVijay Mahadevan x[9] = x[10] = 0.81684757298046; 81163d025dbSVijay Mahadevan x[11] = 0.09157621350977; 81263d025dbSVijay Mahadevan w[0] = w[1] = w[2] = 0.22338158967801; 813181f196bSVijay Mahadevan w[3] = w[4] = w[5] = 0.10995174365532; 814181f196bSVijay Mahadevan } 815181f196bSVijay Mahadevan else { 816181f196bSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Triangle quadrature rules for points 3 and 6 supported; npoints : %D", npoints); 81763d025dbSVijay Mahadevan } 81863d025dbSVijay Mahadevan ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr); 81963d025dbSVijay Mahadevan ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr); 820181f196bSVijay Mahadevan ierr = PetscQuadratureSetData(*quadrature, dim, npoints, x, w);CHKERRQ(ierr); 821181f196bSVijay Mahadevan /* ierr = PetscDTGaussJacobiQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */ 82263d025dbSVijay Mahadevan } 82363d025dbSVijay Mahadevan else { 82463d025dbSVijay Mahadevan ierr = PetscDTGaussTensorQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); 82563d025dbSVijay Mahadevan } 82663d025dbSVijay Mahadevan break; 82763d025dbSVijay Mahadevan case 3: 82863d025dbSVijay Mahadevan /* Create Gauss quadrature rules with <order = nverts> in the span [-1, 1] */ 82963d025dbSVijay Mahadevan if (nverts == 4) { 830*a86ed394SVijay Mahadevan PetscInt order; 831*a86ed394SVijay Mahadevan const PetscInt npoints = 4; // Choose between 4 and 10 832181f196bSVijay Mahadevan ierr = PetscMalloc2(npoints * 3, &x, npoints, &w);CHKERRQ(ierr); 833181f196bSVijay Mahadevan if (npoints == 4) { /* KEAST1, K1, N=4, O=4 */ 834181f196bSVijay Mahadevan const PetscReal x_4[12] = { 0.5854101966249685, 0.1381966011250105, 0.1381966011250105, 835181f196bSVijay Mahadevan 0.1381966011250105, 0.1381966011250105, 0.1381966011250105, 836181f196bSVijay Mahadevan 0.1381966011250105, 0.1381966011250105, 0.5854101966249685, 837181f196bSVijay Mahadevan 0.1381966011250105, 0.5854101966249685, 0.1381966011250105 838181f196bSVijay Mahadevan }; 839181f196bSVijay Mahadevan ierr = PetscMemcpy(x, x_4, 12 * sizeof(PetscReal));CHKERRQ(ierr); 840181f196bSVijay Mahadevan 841181f196bSVijay Mahadevan w[0] = w[1] = w[2] = w[3] = 1.0 / 24.0; 842181f196bSVijay Mahadevan order = 4; 843181f196bSVijay Mahadevan } 844181f196bSVijay Mahadevan else if (npoints == 10) { /* KEAST3, K3 N=10, O=10 */ 845181f196bSVijay Mahadevan const PetscReal x_10[30] = { 0.5684305841968444, 0.1438564719343852, 0.1438564719343852, 846181f196bSVijay Mahadevan 0.1438564719343852, 0.1438564719343852, 0.1438564719343852, 847181f196bSVijay Mahadevan 0.1438564719343852, 0.1438564719343852, 0.5684305841968444, 848181f196bSVijay Mahadevan 0.1438564719343852, 0.5684305841968444, 0.1438564719343852, 849181f196bSVijay Mahadevan 0.0000000000000000, 0.5000000000000000, 0.5000000000000000, 850181f196bSVijay Mahadevan 0.5000000000000000, 0.0000000000000000, 0.5000000000000000, 851181f196bSVijay Mahadevan 0.5000000000000000, 0.5000000000000000, 0.0000000000000000, 852181f196bSVijay Mahadevan 0.5000000000000000, 0.0000000000000000, 0.0000000000000000, 853181f196bSVijay Mahadevan 0.0000000000000000, 0.5000000000000000, 0.0000000000000000, 854181f196bSVijay Mahadevan 0.0000000000000000, 0.0000000000000000, 0.5000000000000000 855181f196bSVijay Mahadevan }; 856181f196bSVijay Mahadevan ierr = PetscMemcpy(x, x_10, 30 * sizeof(PetscReal));CHKERRQ(ierr); 857181f196bSVijay Mahadevan 858181f196bSVijay Mahadevan w[0] = w[1] = w[2] = w[3] = 0.2177650698804054; 859181f196bSVijay Mahadevan w[4] = w[5] = w[6] = w[7] = w[8] = w[9] = 0.0214899534130631; 860181f196bSVijay Mahadevan order = 10; 861181f196bSVijay Mahadevan } 862181f196bSVijay Mahadevan else { 863181f196bSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Tetrahedral quadrature rules for points 4 and 10 supported; npoints : %D", npoints); 864181f196bSVijay Mahadevan } 865181f196bSVijay Mahadevan ierr = PetscQuadratureCreate(PETSC_COMM_SELF, quadrature);CHKERRQ(ierr); 866181f196bSVijay Mahadevan ierr = PetscQuadratureSetOrder(*quadrature, order);CHKERRQ(ierr); 867181f196bSVijay Mahadevan ierr = PetscQuadratureSetData(*quadrature, dim, npoints, x, w);CHKERRQ(ierr); 868181f196bSVijay Mahadevan /* ierr = PetscDTGaussJacobiQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); */ 86963d025dbSVijay Mahadevan } 87063d025dbSVijay Mahadevan else { 87163d025dbSVijay Mahadevan ierr = PetscDTGaussTensorQuadrature(dim, nverts, 0.0, 1.0, quadrature);CHKERRQ(ierr); 87263d025dbSVijay Mahadevan } 87363d025dbSVijay Mahadevan break; 87463d025dbSVijay Mahadevan default: 87563d025dbSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim); 87663d025dbSVijay Mahadevan } 87763d025dbSVijay Mahadevan PetscFunctionReturn(0); 87863d025dbSVijay Mahadevan } 87963d025dbSVijay Mahadevan 880181f196bSVijay Mahadevan /* Compute Jacobians */ 881181f196bSVijay Mahadevan #undef __FUNCT__ 882181f196bSVijay Mahadevan #define __FUNCT__ "ComputeJacobian_Internal" 883181f196bSVijay Mahadevan PetscErrorCode ComputeJacobian_Internal ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quad, PetscReal *phypts, 884*a86ed394SVijay Mahadevan PetscReal *jacobian, PetscReal *ijacobian, PetscReal* dvolume) 885181f196bSVijay Mahadevan { 886*a86ed394SVijay Mahadevan PetscInt i; 887*a86ed394SVijay Mahadevan PetscReal volume; 888181f196bSVijay Mahadevan PetscErrorCode ierr; 889181f196bSVijay Mahadevan 890181f196bSVijay Mahadevan PetscFunctionBegin; 891181f196bSVijay Mahadevan PetscValidPointer(coordinates, 3); 892181f196bSVijay Mahadevan PetscValidPointer(quad, 4); 893181f196bSVijay Mahadevan PetscValidPointer(jacobian, 5); 894181f196bSVijay Mahadevan ierr = PetscMemzero(jacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr); 895181f196bSVijay Mahadevan if (ijacobian) { 896181f196bSVijay Mahadevan ierr = PetscMemzero(ijacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr); 897181f196bSVijay Mahadevan } 898181f196bSVijay Mahadevan if (phypts) { 899181f196bSVijay Mahadevan ierr = PetscMemzero(phypts, /*npts=1 * */ 3 * sizeof(PetscReal));CHKERRQ(ierr); 900181f196bSVijay Mahadevan } 901181f196bSVijay Mahadevan 902181f196bSVijay Mahadevan if (dim == 1) { 903181f196bSVijay Mahadevan 904181f196bSVijay Mahadevan const PetscReal& r = quad[0]; 905181f196bSVijay Mahadevan if (nverts == 2) { /* Linear Edge */ 906181f196bSVijay Mahadevan const PetscReal dNi_dxi[2] = { -1.0, 1.0 }; 907181f196bSVijay Mahadevan 908181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 909181f196bSVijay Mahadevan const PetscReal* vertices = coordinates + i * 3; 910181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 911181f196bSVijay Mahadevan } 912181f196bSVijay Mahadevan } 913181f196bSVijay Mahadevan else if (nverts == 3) { /* Quadratic Edge */ 914181f196bSVijay Mahadevan const PetscReal dNi_dxi[3] = { 4 * r - 3.0, 4 * ( 1.0 - 2.0 * r ), 4.0 * r - 1.0}; 915181f196bSVijay Mahadevan 916181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 917181f196bSVijay Mahadevan const PetscReal* vertices = coordinates + i * 3; 918181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 919181f196bSVijay Mahadevan } 920181f196bSVijay Mahadevan } 921181f196bSVijay Mahadevan else { 922181f196bSVijay 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); 923181f196bSVijay Mahadevan } 924181f196bSVijay Mahadevan 925181f196bSVijay Mahadevan if (ijacobian) { 926181f196bSVijay Mahadevan /* invert the jacobian */ 927181f196bSVijay Mahadevan ijacobian[0] = 1.0 / jacobian[0]; 928181f196bSVijay Mahadevan } 929181f196bSVijay Mahadevan 930181f196bSVijay Mahadevan } 931181f196bSVijay Mahadevan else if (dim == 2) { 932181f196bSVijay Mahadevan 933181f196bSVijay Mahadevan if (nverts == 4) { /* Linear Quadrangle */ 934181f196bSVijay Mahadevan const PetscReal& r = quad[0]; 935181f196bSVijay Mahadevan const PetscReal& s = quad[1]; 936181f196bSVijay Mahadevan 937181f196bSVijay Mahadevan const PetscReal dNi_dxi[4] = { -1.0 + s, 1.0 - s, s, -s }; 938181f196bSVijay Mahadevan const PetscReal dNi_deta[4] = { -1.0 + r, -r, r, 1.0 - r }; 939181f196bSVijay Mahadevan 940181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 941181f196bSVijay Mahadevan const PetscReal* vertices = coordinates + i * 3; 942181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertices[0]; 943181f196bSVijay Mahadevan jacobian[2] += dNi_dxi[i] * vertices[1]; 944181f196bSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertices[0]; 945181f196bSVijay Mahadevan jacobian[3] += dNi_deta[i] * vertices[1]; 946181f196bSVijay Mahadevan } 947181f196bSVijay Mahadevan } 948181f196bSVijay Mahadevan else if (nverts == 3) { /* Linear triangle */ 949181f196bSVijay Mahadevan const PetscReal x2 = coordinates[2 * 3 + 0], y2 = coordinates[2 * 3 + 1]; 950181f196bSVijay Mahadevan 951181f196bSVijay Mahadevan /* Jacobian is constant */ 952181f196bSVijay Mahadevan jacobian[0] = (coordinates[0 * 3 + 0] - x2); jacobian[1] = (coordinates[1 * 3 + 0] - x2); 953181f196bSVijay Mahadevan jacobian[2] = (coordinates[0 * 3 + 1] - y2); jacobian[3] = (coordinates[1 * 3 + 1] - y2); 954181f196bSVijay Mahadevan } 955181f196bSVijay Mahadevan else { 956181f196bSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 2-D entity vertices are invalid. Currently only support QUAD4 and TRI3 basis evaluations in 2-D : %D", nverts); 957181f196bSVijay Mahadevan } 958181f196bSVijay Mahadevan 959181f196bSVijay Mahadevan /* invert the jacobian */ 960181f196bSVijay Mahadevan if (ijacobian) { 961*a86ed394SVijay Mahadevan ierr = DMatrix_Invert_2x2_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr); 962181f196bSVijay Mahadevan } 963181f196bSVijay Mahadevan 964181f196bSVijay Mahadevan } 965181f196bSVijay Mahadevan else { 966181f196bSVijay Mahadevan 967181f196bSVijay Mahadevan if (nverts == 8) { /* Linear Hexahedra */ 968181f196bSVijay Mahadevan const PetscReal& r = quad[0]; 969181f196bSVijay Mahadevan const PetscReal& s = quad[1]; 970181f196bSVijay Mahadevan const PetscReal& t = quad[2]; 971181f196bSVijay Mahadevan const PetscReal dNi_dxi[8] = {- ( 1.0 - s ) * ( 1.0 - t ), 972181f196bSVijay Mahadevan ( 1.0 - s ) * ( 1.0 - t ), 973*a86ed394SVijay Mahadevan ( s ) * ( 1.0 - t ), 974*a86ed394SVijay Mahadevan - ( s ) * ( 1.0 - t ), 975*a86ed394SVijay Mahadevan - ( 1.0 - s ) * ( t ), 976*a86ed394SVijay Mahadevan ( 1.0 - s ) * ( t ), 977*a86ed394SVijay Mahadevan ( s ) * ( t ), 978*a86ed394SVijay Mahadevan - ( s ) * ( t ) 979181f196bSVijay Mahadevan }; 980181f196bSVijay Mahadevan 981181f196bSVijay Mahadevan const PetscReal dNi_deta[8] = { - ( 1.0 - r ) * ( 1.0 - t ), 982*a86ed394SVijay Mahadevan - ( r ) * ( 1.0 - t ), 983*a86ed394SVijay Mahadevan ( r ) * ( 1.0 - t ), 984181f196bSVijay Mahadevan ( 1.0 - r ) * ( 1.0 - t ), 985*a86ed394SVijay Mahadevan - ( 1.0 - r ) * ( t ), 986*a86ed394SVijay Mahadevan - ( r ) * ( t ), 987*a86ed394SVijay Mahadevan ( r ) * ( t ), 988*a86ed394SVijay Mahadevan ( 1.0 - r ) * ( t ) 989181f196bSVijay Mahadevan }; 990181f196bSVijay Mahadevan 991181f196bSVijay Mahadevan const PetscReal dNi_dzeta[8] = { - ( 1.0 - r ) * ( 1.0 - s ), 992*a86ed394SVijay Mahadevan - ( r ) * ( 1.0 - s ), 993*a86ed394SVijay Mahadevan - ( r ) * ( s ), 994*a86ed394SVijay Mahadevan - ( 1.0 - r ) * ( s ), 995181f196bSVijay Mahadevan ( 1.0 - r ) * ( 1.0 - s ), 996*a86ed394SVijay Mahadevan ( r ) * ( 1.0 - s ), 997*a86ed394SVijay Mahadevan ( r ) * ( s ), 998*a86ed394SVijay Mahadevan ( 1.0 - r ) * ( s ) 999181f196bSVijay Mahadevan }; 1000*a86ed394SVijay Mahadevan 1001181f196bSVijay Mahadevan for (i = 0; i < nverts; ++i) { 1002181f196bSVijay Mahadevan const PetscReal* vertex = coordinates + i * 3; 1003181f196bSVijay Mahadevan jacobian[0] += dNi_dxi[i] * vertex[0]; 1004181f196bSVijay Mahadevan jacobian[3] += dNi_dxi[i] * vertex[1]; 1005181f196bSVijay Mahadevan jacobian[6] += dNi_dxi[i] * vertex[2]; 1006181f196bSVijay Mahadevan jacobian[1] += dNi_deta[i] * vertex[0]; 1007181f196bSVijay Mahadevan jacobian[4] += dNi_deta[i] * vertex[1]; 1008181f196bSVijay Mahadevan jacobian[7] += dNi_deta[i] * vertex[2]; 1009181f196bSVijay Mahadevan jacobian[2] += dNi_dzeta[i] * vertex[0]; 1010181f196bSVijay Mahadevan jacobian[5] += dNi_dzeta[i] * vertex[1]; 1011181f196bSVijay Mahadevan jacobian[8] += dNi_dzeta[i] * vertex[2]; 1012181f196bSVijay Mahadevan } 1013181f196bSVijay Mahadevan 1014181f196bSVijay Mahadevan } 1015181f196bSVijay Mahadevan else if (nverts == 4) { /* Linear Tetrahedra */ 1016181f196bSVijay Mahadevan const PetscReal x0 = coordinates[/*0 * 3 +*/ 0], y0 = coordinates[/*0 * 3 +*/ 1], z0 = coordinates[/*0 * 3 +*/ 2]; 1017181f196bSVijay Mahadevan 1018181f196bSVijay Mahadevan /* compute the jacobian */ 1019181f196bSVijay Mahadevan jacobian[0] = coordinates[1 * 3 + 0] - x0; jacobian[1] = coordinates[2 * 3 + 0] - x0; jacobian[2] = coordinates[3 * 3 + 0] - x0; 1020181f196bSVijay Mahadevan jacobian[3] = coordinates[1 * 3 + 1] - y0; jacobian[4] = coordinates[2 * 3 + 1] - y0; jacobian[5] = coordinates[3 * 3 + 1] - y0; 1021181f196bSVijay Mahadevan jacobian[6] = coordinates[1 * 3 + 2] - z0; jacobian[7] = coordinates[2 * 3 + 2] - z0; jacobian[8] = coordinates[3 * 3 + 2] - z0; 1022181f196bSVijay Mahadevan } /* Tetrahedra -- ends */ 1023181f196bSVijay Mahadevan else 1024181f196bSVijay Mahadevan { 1025181f196bSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of 3-D entity vertices are invalid. Currently only support HEX8 and TET4 basis evaluations in 3-D : %D", nverts); 1026181f196bSVijay Mahadevan } 1027181f196bSVijay Mahadevan 1028181f196bSVijay Mahadevan if (ijacobian) { 1029181f196bSVijay Mahadevan /* invert the jacobian */ 1030*a86ed394SVijay Mahadevan ierr = DMatrix_Invert_3x3_Internal(jacobian, ijacobian, &volume);CHKERRQ(ierr); 1031181f196bSVijay Mahadevan } 1032181f196bSVijay Mahadevan 1033181f196bSVijay Mahadevan } 1034*a86ed394SVijay 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); 1035*a86ed394SVijay Mahadevan if (dvolume) *dvolume = volume; 1036181f196bSVijay Mahadevan PetscFunctionReturn(0); 1037181f196bSVijay Mahadevan } 1038181f196bSVijay Mahadevan 1039181f196bSVijay Mahadevan 1040181f196bSVijay Mahadevan #undef __FUNCT__ 1041181f196bSVijay Mahadevan #define __FUNCT__ "FEMComputeBasis_JandF" 1042181f196bSVijay Mahadevan PetscErrorCode FEMComputeBasis_JandF ( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal *quadrature, PetscReal *phypts, 1043181f196bSVijay Mahadevan PetscReal *phibasis, PetscReal *jacobian, PetscReal *ijacobian, PetscReal* volume ) 1044181f196bSVijay Mahadevan { 1045181f196bSVijay Mahadevan PetscErrorCode ierr; 1046181f196bSVijay Mahadevan 1047181f196bSVijay Mahadevan PetscFunctionBegin; 1048181f196bSVijay Mahadevan switch (dim) { 1049181f196bSVijay Mahadevan case 1: 1050181f196bSVijay Mahadevan ierr = Compute_Lagrange_Basis_1D_Internal(nverts, coordinates, 1, quadrature, phypts, 1051181f196bSVijay Mahadevan NULL, phibasis, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr); 1052181f196bSVijay Mahadevan break; 1053181f196bSVijay Mahadevan case 2: 1054181f196bSVijay Mahadevan ierr = Compute_Lagrange_Basis_2D_Internal(nverts, coordinates, 1, quadrature, phypts, 1055181f196bSVijay Mahadevan NULL, phibasis, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr); 1056181f196bSVijay Mahadevan break; 1057181f196bSVijay Mahadevan case 3: 1058181f196bSVijay Mahadevan ierr = Compute_Lagrange_Basis_3D_Internal(nverts, coordinates, 1, quadrature, phypts, 1059181f196bSVijay Mahadevan NULL, phibasis, NULL, NULL, NULL, jacobian, ijacobian, volume);CHKERRQ(ierr); 1060181f196bSVijay Mahadevan break; 1061181f196bSVijay Mahadevan default: 1062181f196bSVijay Mahadevan SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension; should be in [1,3] : %D", dim); 1063181f196bSVijay Mahadevan } 1064181f196bSVijay Mahadevan PetscFunctionReturn(0); 1065181f196bSVijay Mahadevan } 1066181f196bSVijay Mahadevan 1067181f196bSVijay Mahadevan 1068181f196bSVijay Mahadevan 1069*a86ed394SVijay Mahadevan /*@ 1070*a86ed394SVijay Mahadevan DMMoabPToRMapping: Compute the mapping from the physical coordinate system for a given element to the 1071*a86ed394SVijay Mahadevan canonical reference element. In addition to finding the inverse mapping evaluation through Newton iteration, 1072*a86ed394SVijay Mahadevan the basis function at the parametric point is also evaluated optionally. 1073*a86ed394SVijay Mahadevan 1074*a86ed394SVijay Mahadevan Input Parameter: 1075*a86ed394SVijay Mahadevan 1076*a86ed394SVijay Mahadevan . PetscInt dim, the element dimension (1=EDGE, 2=QUAD/TRI, 3=HEX/TET) 1077*a86ed394SVijay Mahadevan . PetscInt nverts, the number of vertices in the physical element 1078*a86ed394SVijay Mahadevan . PetscReal coordinates, the coordinates of vertices in the physical element 1079*a86ed394SVijay Mahadevan . PetscReal[3] xphy, the coordinates of physical point for which natural coordinates (in reference frame) are sought 1080*a86ed394SVijay Mahadevan 1081*a86ed394SVijay Mahadevan Output Parameter: 1082*a86ed394SVijay Mahadevan . PetscReal[3] natparam, the natural coordinates (in reference frame) corresponding to xphy 1083*a86ed394SVijay Mahadevan . PetscReal[nverts] phi, the basis functions evaluated at the natural coordinates (natparam) 1084*a86ed394SVijay Mahadevan 1085*a86ed394SVijay Mahadevan .keywords: DMMoab, Mapping, FEM 1086*a86ed394SVijay Mahadevan @*/ 1087181f196bSVijay Mahadevan #undef __FUNCT__ 1088181f196bSVijay Mahadevan #define __FUNCT__ "DMMoabPToRMapping" 1089181f196bSVijay Mahadevan PetscErrorCode DMMoabPToRMapping( const PetscInt dim, const PetscInt nverts, const PetscReal *coordinates, const PetscReal* xphy, PetscReal* natparam, PetscReal* phi) 1090181f196bSVijay Mahadevan { 1091*a86ed394SVijay Mahadevan /* Perform inverse evaluation for the mapping with use of Newton Raphson iteration */ 1092181f196bSVijay Mahadevan const PetscReal tol = 1.0e-10; 1093181f196bSVijay Mahadevan const PetscInt max_iterations = 10; 1094181f196bSVijay Mahadevan const PetscReal error_tol_sqr = tol*tol; 1095181f196bSVijay Mahadevan PetscReal phibasis[8], jacobian[9], ijacobian[9], volume; 1096181f196bSVijay Mahadevan PetscReal phypts[3] = {0.0, 0.0, 0.0}; 1097181f196bSVijay Mahadevan PetscReal delta[3] = {0.0, 0.0, 0.0}; 1098181f196bSVijay Mahadevan PetscErrorCode ierr; 1099181f196bSVijay Mahadevan PetscInt iters=0; 1100181f196bSVijay Mahadevan PetscReal error=1.0; 1101181f196bSVijay Mahadevan 1102181f196bSVijay Mahadevan PetscFunctionBegin; 1103181f196bSVijay Mahadevan PetscValidPointer(coordinates, 3); 1104181f196bSVijay Mahadevan PetscValidPointer(xphy, 4); 1105181f196bSVijay Mahadevan PetscValidPointer(natparam, 5); 1106181f196bSVijay Mahadevan 1107181f196bSVijay Mahadevan ierr = PetscMemzero(jacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr); 1108181f196bSVijay Mahadevan ierr = PetscMemzero(ijacobian, dim * dim * sizeof(PetscReal));CHKERRQ(ierr); 1109181f196bSVijay Mahadevan ierr = PetscMemzero(phibasis, nverts * sizeof(PetscReal));CHKERRQ(ierr); 1110181f196bSVijay Mahadevan 1111*a86ed394SVijay Mahadevan /* zero initial guess */ 1112*a86ed394SVijay Mahadevan natparam[0] = natparam[1] = natparam[2] = 0.0; 1113181f196bSVijay Mahadevan 1114*a86ed394SVijay Mahadevan /* Compute delta = evaluate( xi ) - x */ 1115*a86ed394SVijay Mahadevan ierr = FEMComputeBasis_JandF ( dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume );CHKERRQ(ierr); 1116*a86ed394SVijay Mahadevan 1117*a86ed394SVijay Mahadevan error = 0.0; 1118*a86ed394SVijay Mahadevan switch(dim) { 1119*a86ed394SVijay Mahadevan case 3: 1120181f196bSVijay Mahadevan delta[2] = phypts[2] - xphy[2]; 1121*a86ed394SVijay Mahadevan error += (delta[2]*delta[2]); 1122*a86ed394SVijay Mahadevan case 2: 1123*a86ed394SVijay Mahadevan delta[1] = phypts[1] - xphy[1]; 1124*a86ed394SVijay Mahadevan error += (delta[1]*delta[1]); 1125*a86ed394SVijay Mahadevan case 1: 1126*a86ed394SVijay Mahadevan delta[0] = phypts[0] - xphy[0]; 1127*a86ed394SVijay Mahadevan error += (delta[0]*delta[0]); 1128*a86ed394SVijay Mahadevan break; 1129*a86ed394SVijay Mahadevan } 1130*a86ed394SVijay Mahadevan 1131*a86ed394SVijay Mahadevan #if 0 1132*a86ed394SVijay Mahadevan PetscInfo4(NULL, "Current point in physical space : (%g, %g, %g) and error : %g\n", xphy[0], xphy[1], xphy[2], error); 1133*a86ed394SVijay Mahadevan #endif 1134181f196bSVijay Mahadevan 1135181f196bSVijay Mahadevan while (error > error_tol_sqr) { 1136181f196bSVijay Mahadevan 1137181f196bSVijay Mahadevan if(++iters > max_iterations) 1138181f196bSVijay 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]); 1139181f196bSVijay Mahadevan 1140181f196bSVijay Mahadevan /* Compute natparam -= J.inverse() * delta */ 1141181f196bSVijay Mahadevan switch(dim) { 1142181f196bSVijay Mahadevan case 1: 1143181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0]; 1144181f196bSVijay Mahadevan break; 1145181f196bSVijay Mahadevan case 2: 1146181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1]; 1147181f196bSVijay Mahadevan natparam[1] -= ijacobian[2] * delta[0] + ijacobian[3] * delta[1]; 1148181f196bSVijay Mahadevan break; 1149181f196bSVijay Mahadevan case 3: 1150181f196bSVijay Mahadevan natparam[0] -= ijacobian[0] * delta[0] + ijacobian[1] * delta[1] + ijacobian[2] * delta[2]; 1151181f196bSVijay Mahadevan natparam[1] -= ijacobian[3] * delta[0] + ijacobian[4] * delta[1] + ijacobian[5] * delta[2]; 1152181f196bSVijay Mahadevan natparam[2] -= ijacobian[6] * delta[0] + ijacobian[7] * delta[1] + ijacobian[8] * delta[2]; 1153181f196bSVijay Mahadevan break; 1154181f196bSVijay Mahadevan } 1155181f196bSVijay Mahadevan 1156181f196bSVijay Mahadevan /* Compute delta = evaluate( xi ) - x */ 1157*a86ed394SVijay Mahadevan ierr = FEMComputeBasis_JandF ( dim, nverts, coordinates, natparam, phypts, phibasis, jacobian, ijacobian, &volume );CHKERRQ(ierr); 1158181f196bSVijay Mahadevan 1159*a86ed394SVijay Mahadevan error = 0.0; 1160*a86ed394SVijay Mahadevan switch(dim) { 1161*a86ed394SVijay Mahadevan case 3: 1162181f196bSVijay Mahadevan delta[2] = phypts[2] - xphy[2]; 1163*a86ed394SVijay Mahadevan error += (delta[2]*delta[2]); 1164*a86ed394SVijay Mahadevan case 2: 1165*a86ed394SVijay Mahadevan delta[1] = phypts[1] - xphy[1]; 1166*a86ed394SVijay Mahadevan error += (delta[1]*delta[1]); 1167*a86ed394SVijay Mahadevan case 1: 1168*a86ed394SVijay Mahadevan delta[0] = phypts[0] - xphy[0]; 1169*a86ed394SVijay Mahadevan error += (delta[0]*delta[0]); 1170*a86ed394SVijay Mahadevan break; 1171*a86ed394SVijay Mahadevan } 1172*a86ed394SVijay Mahadevan #if 0 1173*a86ed394SVijay Mahadevan PetscInfo5(NULL, "Newton [%D] : Current point in reference space : (%g, %g, %g) and error : %g\n", iters, natparam[0], natparam[1], natparam[2], error); 1174*a86ed394SVijay Mahadevan #endif 1175181f196bSVijay Mahadevan } 1176181f196bSVijay Mahadevan if (phi) { 1177181f196bSVijay Mahadevan ierr = PetscMemcpy(phi, phibasis, nverts * sizeof(PetscReal));CHKERRQ(ierr); 1178181f196bSVijay Mahadevan } 1179181f196bSVijay Mahadevan PetscFunctionReturn(0); 1180181f196bSVijay Mahadevan } 1181181f196bSVijay Mahadevan 1182