120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac 32b99622eSMatthew G. Knepley /*@C 42b99622eSMatthew G. Knepley PetscFEGeomCreate - Create a PetscFEGeom object to manage geometry for a group of cells 52b99622eSMatthew G. Knepley 62b99622eSMatthew G. Knepley Input Parameters: 72b99622eSMatthew G. Knepley + quad - A PetscQuadrature determining the tabulation 82b99622eSMatthew G. Knepley . numCells - The number of cells in the group 92b99622eSMatthew G. Knepley . dimEmbed - The coordinate dimension 102b99622eSMatthew G. Knepley - faceData - Flag to construct geometry data for the faces 112b99622eSMatthew G. Knepley 122b99622eSMatthew G. Knepley Output Parameter: 132b99622eSMatthew G. Knepley . geom - The PetscFEGeom object 142b99622eSMatthew G. Knepley 152b99622eSMatthew G. Knepley Level: beginner 162b99622eSMatthew G. Knepley 17db781477SPatrick Sanan .seealso: `PetscFEGeomDestroy()`, `PetscFEGeomComplete()` 182b99622eSMatthew G. Knepley @*/ 19*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom) 20*d71ae5a4SJacob Faibussowitsch { 2120cf1dd8SToby Isaac PetscFEGeom *g; 227489efa5SMatthew G. Knepley PetscInt dim, Nq, N; 2320cf1dd8SToby Isaac const PetscReal *p; 2420cf1dd8SToby Isaac 2520cf1dd8SToby Isaac PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &dim, NULL, &Nq, &p, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscNew(&g)); 2820cf1dd8SToby Isaac g->xi = p; 2920cf1dd8SToby Isaac g->numCells = numCells; 307489efa5SMatthew G. Knepley g->numPoints = Nq; 3120cf1dd8SToby Isaac g->dim = dim; 3220cf1dd8SToby Isaac g->dimEmbed = dimEmbed; 335fedec97SMatthew G. Knepley g->isCohesive = PETSC_FALSE; 347489efa5SMatthew G. Knepley N = numCells * Nq; 359566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ)); 3620cf1dd8SToby Isaac if (faceData) { 379566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numCells, &g->face, N * dimEmbed, &g->n)); 389371c9d4SSatish Balay PetscCall(PetscCalloc6(N * dimEmbed * dimEmbed, &(g->suppJ[0]), N * dimEmbed * dimEmbed, &(g->suppJ[1]), N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]), N, &(g->suppDetJ[0]), N, &(g->suppDetJ[1]))); 3920cf1dd8SToby Isaac } 409566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ)); 4120cf1dd8SToby Isaac *geom = g; 4220cf1dd8SToby Isaac PetscFunctionReturn(0); 4320cf1dd8SToby Isaac } 4420cf1dd8SToby Isaac 452b99622eSMatthew G. Knepley /*@C 462b99622eSMatthew G. Knepley PetscFEGeomDestroy - Destroy a PetscFEGeom object 472b99622eSMatthew G. Knepley 482b99622eSMatthew G. Knepley Input Parameter: 492b99622eSMatthew G. Knepley . geom - PetscFEGeom object 502b99622eSMatthew G. Knepley 512b99622eSMatthew G. Knepley Level: beginner 522b99622eSMatthew G. Knepley 53db781477SPatrick Sanan .seealso: `PetscFEGeomCreate()` 542b99622eSMatthew G. Knepley @*/ 55*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom) 56*d71ae5a4SJacob Faibussowitsch { 5720cf1dd8SToby Isaac PetscFunctionBegin; 5820cf1dd8SToby Isaac if (!*geom) PetscFunctionReturn(0); 599566063dSJacob Faibussowitsch PetscCall(PetscFree3((*geom)->v, (*geom)->J, (*geom)->detJ)); 609566063dSJacob Faibussowitsch PetscCall(PetscFree((*geom)->invJ)); 619566063dSJacob Faibussowitsch PetscCall(PetscFree2((*geom)->face, (*geom)->n)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree6((*geom)->suppJ[0], (*geom)->suppJ[1], (*geom)->suppInvJ[0], (*geom)->suppInvJ[1], (*geom)->suppDetJ[0], (*geom)->suppDetJ[1])); 639566063dSJacob Faibussowitsch PetscCall(PetscFree(*geom)); 6420cf1dd8SToby Isaac PetscFunctionReturn(0); 6520cf1dd8SToby Isaac } 6620cf1dd8SToby Isaac 672b99622eSMatthew G. Knepley /*@C 682b99622eSMatthew G. Knepley PetscFEGeomGetChunk - Get a chunk of cells in the group as a PetscFEGeom 692b99622eSMatthew G. Knepley 702b99622eSMatthew G. Knepley Input Parameters: 712b99622eSMatthew G. Knepley + geom - PetscFEGeom object 722b99622eSMatthew G. Knepley . cStart - The first cell in the chunk 732b99622eSMatthew G. Knepley - cEnd - The first cell not in the chunk 742b99622eSMatthew G. Knepley 752b99622eSMatthew G. Knepley Output Parameter: 762b99622eSMatthew G. Knepley . chunkGeom - The chunk of cells 772b99622eSMatthew G. Knepley 782b99622eSMatthew G. Knepley Level: intermediate 792b99622eSMatthew G. Knepley 80db781477SPatrick Sanan .seealso: `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()` 812b99622eSMatthew G. Knepley @*/ 82*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 83*d71ae5a4SJacob Faibussowitsch { 8420cf1dd8SToby Isaac PetscInt Nq; 8520cf1dd8SToby Isaac PetscInt dE; 8620cf1dd8SToby Isaac 8720cf1dd8SToby Isaac PetscFunctionBegin; 8820cf1dd8SToby Isaac PetscValidPointer(geom, 1); 89064a246eSJacob Faibussowitsch PetscValidPointer(chunkGeom, 4); 9048a46eb9SPierre Jolivet if (!(*chunkGeom)) PetscCall(PetscNew(chunkGeom)); 9120cf1dd8SToby Isaac Nq = geom->numPoints; 9220cf1dd8SToby Isaac dE = geom->dimEmbed; 9320cf1dd8SToby Isaac (*chunkGeom)->dim = geom->dim; 9420cf1dd8SToby Isaac (*chunkGeom)->dimEmbed = geom->dimEmbed; 9520cf1dd8SToby Isaac (*chunkGeom)->numPoints = geom->numPoints; 9620cf1dd8SToby Isaac (*chunkGeom)->numCells = cEnd - cStart; 9720cf1dd8SToby Isaac (*chunkGeom)->xi = geom->xi; 9820cf1dd8SToby Isaac (*chunkGeom)->v = &geom->v[Nq * dE * cStart]; 9920cf1dd8SToby Isaac (*chunkGeom)->J = &geom->J[Nq * dE * dE * cStart]; 10020cf1dd8SToby Isaac (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq * dE * dE * cStart] : NULL; 10120cf1dd8SToby Isaac (*chunkGeom)->detJ = &geom->detJ[Nq * cStart]; 10220cf1dd8SToby Isaac (*chunkGeom)->n = geom->n ? &geom->n[Nq * dE * cStart] : NULL; 10320cf1dd8SToby Isaac (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL; 10441418987SMatthew G. Knepley (*chunkGeom)->suppJ[0] = geom->suppJ[0] ? &geom->suppJ[0][Nq * dE * dE * cStart] : NULL; 10541418987SMatthew G. Knepley (*chunkGeom)->suppJ[1] = geom->suppJ[1] ? &geom->suppJ[1][Nq * dE * dE * cStart] : NULL; 10620cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq * dE * dE * cStart] : NULL; 10720cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq * dE * dE * cStart] : NULL; 10841418987SMatthew G. Knepley (*chunkGeom)->suppDetJ[0] = geom->suppDetJ[0] ? &geom->suppDetJ[0][Nq * cStart] : NULL; 10941418987SMatthew G. Knepley (*chunkGeom)->suppDetJ[1] = geom->suppDetJ[1] ? &geom->suppDetJ[1][Nq * cStart] : NULL; 11020cf1dd8SToby Isaac (*chunkGeom)->isAffine = geom->isAffine; 11120cf1dd8SToby Isaac PetscFunctionReturn(0); 11220cf1dd8SToby Isaac } 11320cf1dd8SToby Isaac 1142b99622eSMatthew G. Knepley /*@C 1152b99622eSMatthew G. Knepley PetscFEGeomRestoreChunk - Restore the chunk 1162b99622eSMatthew G. Knepley 1172b99622eSMatthew G. Knepley Input Parameters: 1182b99622eSMatthew G. Knepley + geom - PetscFEGeom object 1192b99622eSMatthew G. Knepley . cStart - The first cell in the chunk 1202b99622eSMatthew G. Knepley . cEnd - The first cell not in the chunk 1212b99622eSMatthew G. Knepley - chunkGeom - The chunk of cells 1222b99622eSMatthew G. Knepley 1232b99622eSMatthew G. Knepley Level: intermediate 1242b99622eSMatthew G. Knepley 125db781477SPatrick Sanan .seealso: `PetscFEGeomGetChunk()`, `PetscFEGeomCreate()` 1262b99622eSMatthew G. Knepley @*/ 127*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 128*d71ae5a4SJacob Faibussowitsch { 12920cf1dd8SToby Isaac PetscFunctionBegin; 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(*chunkGeom)); 13120cf1dd8SToby Isaac PetscFunctionReturn(0); 13220cf1dd8SToby Isaac } 13320cf1dd8SToby Isaac 1346587ee25SMatthew G. Knepley /*@C 1356587ee25SMatthew G. Knepley PetscFEGeomGetPoint - Get the geometry for cell c at point p as a PetscFEGeom 1366587ee25SMatthew G. Knepley 1376587ee25SMatthew G. Knepley Input Parameters: 1386587ee25SMatthew G. Knepley + geom - PetscFEGeom object 1396587ee25SMatthew G. Knepley . c - The cell 1406587ee25SMatthew G. Knepley . p - The point 1416587ee25SMatthew G. Knepley - pcoords - The reference coordinates of point p, or NULL 1426587ee25SMatthew G. Knepley 1436587ee25SMatthew G. Knepley Output Parameter: 1446587ee25SMatthew G. Knepley . pgeom - The geometry of cell c at point p 1456587ee25SMatthew G. Knepley 1466587ee25SMatthew G. Knepley Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom, 1476587ee25SMatthew G. Knepley nothing needs to be done with it afterwards. 1486587ee25SMatthew G. Knepley 1496587ee25SMatthew G. Knepley In the affine case, pgeom must have storage for the integration point coordinates in pgeom->v if pcoords is passed in. 1506587ee25SMatthew G. Knepley 1516587ee25SMatthew G. Knepley Level: intermediate 1526587ee25SMatthew G. Knepley 153db781477SPatrick Sanan .seealso: `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()` 1546587ee25SMatthew G. Knepley @*/ 155*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, const PetscReal pcoords[], PetscFEGeom *pgeom) 156*d71ae5a4SJacob Faibussowitsch { 1576587ee25SMatthew G. Knepley const PetscInt dim = geom->dim; 1586587ee25SMatthew G. Knepley const PetscInt dE = geom->dimEmbed; 1596587ee25SMatthew G. Knepley const PetscInt Np = geom->numPoints; 1606587ee25SMatthew G. Knepley 1616587ee25SMatthew G. Knepley PetscFunctionBeginHot; 1626587ee25SMatthew G. Knepley pgeom->dim = dim; 1636587ee25SMatthew G. Knepley pgeom->dimEmbed = dE; 1646587ee25SMatthew G. Knepley //pgeom->isAffine = geom->isAffine; 1656587ee25SMatthew G. Knepley if (geom->isAffine) { 1666587ee25SMatthew G. Knepley if (!p) { 1676587ee25SMatthew G. Knepley pgeom->xi = geom->xi; 1686587ee25SMatthew G. Knepley pgeom->J = &geom->J[c * Np * dE * dE]; 1696587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[c * Np * dE * dE]; 1706587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np]; 1716587ee25SMatthew G. Knepley pgeom->n = geom->n ? &geom->n[c * Np * dE] : NULL; 1726587ee25SMatthew G. Knepley } 173ad540459SPierre Jolivet if (pcoords) CoordinatesRefToReal(dE, dim, pgeom->xi, &geom->v[c * Np * dE], pgeom->J, pcoords, pgeom->v); 1746587ee25SMatthew G. Knepley } else { 1756587ee25SMatthew G. Knepley pgeom->v = &geom->v[(c * Np + p) * dE]; 1766587ee25SMatthew G. Knepley pgeom->J = &geom->J[(c * Np + p) * dE * dE]; 1776587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[(c * Np + p) * dE * dE]; 1786587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np + p]; 1796587ee25SMatthew G. Knepley pgeom->n = geom->n ? &geom->n[(c * Np + p) * dE] : NULL; 1806587ee25SMatthew G. Knepley } 1816587ee25SMatthew G. Knepley PetscFunctionReturn(0); 1826587ee25SMatthew G. Knepley } 1836587ee25SMatthew G. Knepley 1846587ee25SMatthew G. Knepley /*@C 1856587ee25SMatthew G. Knepley PetscFEGeomGetCellPoint - Get the cell geometry for face f at point p as a PetscFEGeom 1866587ee25SMatthew G. Knepley 1876587ee25SMatthew G. Knepley Input Parameters: 1886587ee25SMatthew G. Knepley + geom - PetscFEGeom object 1896587ee25SMatthew G. Knepley . f - The face 1906587ee25SMatthew G. Knepley - p - The point 1916587ee25SMatthew G. Knepley 1926587ee25SMatthew G. Knepley Output Parameter: 1936587ee25SMatthew G. Knepley . pgeom - The cell geometry of face f at point p 1946587ee25SMatthew G. Knepley 1956587ee25SMatthew G. Knepley Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom, 1966587ee25SMatthew G. Knepley nothing needs to be done with it afterwards. 1976587ee25SMatthew G. Knepley 1986587ee25SMatthew G. Knepley Level: intermediate 1996587ee25SMatthew G. Knepley 200db781477SPatrick Sanan .seealso: `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()` 2016587ee25SMatthew G. Knepley @*/ 202*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, PetscFEGeom *pgeom) 203*d71ae5a4SJacob Faibussowitsch { 2045fedec97SMatthew G. Knepley const PetscBool bd = geom->dimEmbed > geom->dim && !geom->isCohesive ? PETSC_TRUE : PETSC_FALSE; 2056587ee25SMatthew G. Knepley const PetscInt dim = bd ? geom->dimEmbed : geom->dim; 2066587ee25SMatthew G. Knepley const PetscInt dE = geom->dimEmbed; 2076587ee25SMatthew G. Knepley const PetscInt Np = geom->numPoints; 2086587ee25SMatthew G. Knepley 2096587ee25SMatthew G. Knepley PetscFunctionBeginHot; 2106587ee25SMatthew G. Knepley pgeom->dim = dim; 2116587ee25SMatthew G. Knepley pgeom->dimEmbed = dE; 2126587ee25SMatthew G. Knepley //pgeom->isAffine = geom->isAffine; 2136587ee25SMatthew G. Knepley if (geom->isAffine) { 2146587ee25SMatthew G. Knepley if (!p) { 2156587ee25SMatthew G. Knepley if (bd) { 2166587ee25SMatthew G. Knepley pgeom->J = &geom->suppJ[0][c * Np * dE * dE]; 2176587ee25SMatthew G. Knepley pgeom->invJ = &geom->suppInvJ[0][c * Np * dE * dE]; 2186587ee25SMatthew G. Knepley pgeom->detJ = &geom->suppDetJ[0][c * Np]; 2196587ee25SMatthew G. Knepley } else { 2206587ee25SMatthew G. Knepley pgeom->J = &geom->J[c * Np * dE * dE]; 2216587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[c * Np * dE * dE]; 2226587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np]; 2236587ee25SMatthew G. Knepley } 2246587ee25SMatthew G. Knepley } 2256587ee25SMatthew G. Knepley } else { 2266587ee25SMatthew G. Knepley if (bd) { 2276587ee25SMatthew G. Knepley pgeom->J = &geom->suppJ[0][(c * Np + p) * dE * dE]; 2286587ee25SMatthew G. Knepley pgeom->invJ = &geom->suppInvJ[0][(c * Np + p) * dE * dE]; 2296587ee25SMatthew G. Knepley pgeom->detJ = &geom->suppDetJ[0][c * Np + p]; 2306587ee25SMatthew G. Knepley } else { 2316587ee25SMatthew G. Knepley pgeom->J = &geom->J[(c * Np + p) * dE * dE]; 2326587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[(c * Np + p) * dE * dE]; 2336587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np + p]; 2346587ee25SMatthew G. Knepley } 2356587ee25SMatthew G. Knepley } 2366587ee25SMatthew G. Knepley PetscFunctionReturn(0); 2376587ee25SMatthew G. Knepley } 2386587ee25SMatthew G. Knepley 2392b99622eSMatthew G. Knepley /*@ 2402b99622eSMatthew G. Knepley PetscFEGeomComplete - Calculate derived quntites from base geometry specification 2412b99622eSMatthew G. Knepley 2422b99622eSMatthew G. Knepley Input Parameter: 2432b99622eSMatthew G. Knepley . geom - PetscFEGeom object 2442b99622eSMatthew G. Knepley 2452b99622eSMatthew G. Knepley Level: intermediate 2462b99622eSMatthew G. Knepley 247db781477SPatrick Sanan .seealso: `PetscFEGeomCreate()` 2482b99622eSMatthew G. Knepley @*/ 249*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom) 250*d71ae5a4SJacob Faibussowitsch { 25120cf1dd8SToby Isaac PetscInt i, j, N, dE; 25220cf1dd8SToby Isaac 25320cf1dd8SToby Isaac PetscFunctionBeginHot; 25420cf1dd8SToby Isaac N = geom->numPoints * geom->numCells; 25520cf1dd8SToby Isaac dE = geom->dimEmbed; 25620cf1dd8SToby Isaac switch (dE) { 25720cf1dd8SToby Isaac case 3: 25820cf1dd8SToby Isaac for (i = 0; i < N; i++) { 25920cf1dd8SToby Isaac DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE * dE * i]); 26020cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE * dE * i], &geom->J[dE * dE * i], geom->detJ[i]); 26120cf1dd8SToby Isaac } 26220cf1dd8SToby Isaac break; 26320cf1dd8SToby Isaac case 2: 26420cf1dd8SToby Isaac for (i = 0; i < N; i++) { 26520cf1dd8SToby Isaac DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE * dE * i]); 26620cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE * dE * i], &geom->J[dE * dE * i], geom->detJ[i]); 26720cf1dd8SToby Isaac } 26820cf1dd8SToby Isaac break; 26920cf1dd8SToby Isaac case 1: 27020cf1dd8SToby Isaac for (i = 0; i < N; i++) { 27120cf1dd8SToby Isaac geom->detJ[i] = PetscAbsReal(geom->J[i]); 27220cf1dd8SToby Isaac if (geom->invJ) geom->invJ[i] = 1. / geom->J[i]; 27320cf1dd8SToby Isaac } 27420cf1dd8SToby Isaac break; 27520cf1dd8SToby Isaac } 27620cf1dd8SToby Isaac if (geom->n) { 27720cf1dd8SToby Isaac for (i = 0; i < N; i++) { 278ad540459SPierre Jolivet for (j = 0; j < dE; j++) geom->n[dE * i + j] = geom->J[dE * dE * i + dE * j + dE - 1] * ((dE == 2) ? -1. : 1.); 27920cf1dd8SToby Isaac } 28020cf1dd8SToby Isaac } 28120cf1dd8SToby Isaac PetscFunctionReturn(0); 28220cf1dd8SToby Isaac } 283