120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac 32b99622eSMatthew G. Knepley /*@C 4*dce8aebaSBarry Smith PetscFEGeomCreate - Create a `PetscFEGeom` object to manage geometry for a group of cells 52b99622eSMatthew G. Knepley 62b99622eSMatthew G. Knepley Input Parameters: 7*dce8aebaSBarry Smith + 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: 13*dce8aebaSBarry Smith . geom - The `PetscFEGeom` object 142b99622eSMatthew G. Knepley 152b99622eSMatthew G. Knepley Level: beginner 162b99622eSMatthew G. Knepley 17*dce8aebaSBarry Smith .seealso: `PetscFEGeom`, `PetscQuadrature`, `PetscFEGeom`, `PetscFEGeomDestroy()`, `PetscFEGeomComplete()` 182b99622eSMatthew G. Knepley @*/ 19d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom) 20d71ae5a4SJacob 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 46*dce8aebaSBarry Smith PetscFEGeomDestroy - Destroy a `PetscFEGeom` object 472b99622eSMatthew G. Knepley 482b99622eSMatthew G. Knepley Input Parameter: 49*dce8aebaSBarry Smith . geom - `PetscFEGeom` object 502b99622eSMatthew G. Knepley 512b99622eSMatthew G. Knepley Level: beginner 522b99622eSMatthew G. Knepley 53*dce8aebaSBarry Smith .seealso: `PetscFEGeom`, `PetscFEGeomCreate()` 542b99622eSMatthew G. Knepley @*/ 55d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom) 56d71ae5a4SJacob 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 68*dce8aebaSBarry Smith PetscFEGeomGetChunk - Get a chunk of cells in the group as a `PetscFEGeom` 692b99622eSMatthew G. Knepley 702b99622eSMatthew G. Knepley Input Parameters: 71*dce8aebaSBarry Smith + 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 80*dce8aebaSBarry Smith Note: 81*dce8aebaSBarry Smith Use `PetscFEGeomRestoreChunk()` to return the result 82*dce8aebaSBarry Smith 83*dce8aebaSBarry Smith .seealso: `PetscFEGeom`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()` 842b99622eSMatthew G. Knepley @*/ 85d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 86d71ae5a4SJacob Faibussowitsch { 8720cf1dd8SToby Isaac PetscInt Nq; 8820cf1dd8SToby Isaac PetscInt dE; 8920cf1dd8SToby Isaac 9020cf1dd8SToby Isaac PetscFunctionBegin; 9120cf1dd8SToby Isaac PetscValidPointer(geom, 1); 92064a246eSJacob Faibussowitsch PetscValidPointer(chunkGeom, 4); 9348a46eb9SPierre Jolivet if (!(*chunkGeom)) PetscCall(PetscNew(chunkGeom)); 9420cf1dd8SToby Isaac Nq = geom->numPoints; 9520cf1dd8SToby Isaac dE = geom->dimEmbed; 9620cf1dd8SToby Isaac (*chunkGeom)->dim = geom->dim; 9720cf1dd8SToby Isaac (*chunkGeom)->dimEmbed = geom->dimEmbed; 9820cf1dd8SToby Isaac (*chunkGeom)->numPoints = geom->numPoints; 9920cf1dd8SToby Isaac (*chunkGeom)->numCells = cEnd - cStart; 10020cf1dd8SToby Isaac (*chunkGeom)->xi = geom->xi; 10120cf1dd8SToby Isaac (*chunkGeom)->v = &geom->v[Nq * dE * cStart]; 10220cf1dd8SToby Isaac (*chunkGeom)->J = &geom->J[Nq * dE * dE * cStart]; 10320cf1dd8SToby Isaac (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq * dE * dE * cStart] : NULL; 10420cf1dd8SToby Isaac (*chunkGeom)->detJ = &geom->detJ[Nq * cStart]; 10520cf1dd8SToby Isaac (*chunkGeom)->n = geom->n ? &geom->n[Nq * dE * cStart] : NULL; 10620cf1dd8SToby Isaac (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL; 10741418987SMatthew G. Knepley (*chunkGeom)->suppJ[0] = geom->suppJ[0] ? &geom->suppJ[0][Nq * dE * dE * cStart] : NULL; 10841418987SMatthew G. Knepley (*chunkGeom)->suppJ[1] = geom->suppJ[1] ? &geom->suppJ[1][Nq * dE * dE * cStart] : NULL; 10920cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq * dE * dE * cStart] : NULL; 11020cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq * dE * dE * cStart] : NULL; 11141418987SMatthew G. Knepley (*chunkGeom)->suppDetJ[0] = geom->suppDetJ[0] ? &geom->suppDetJ[0][Nq * cStart] : NULL; 11241418987SMatthew G. Knepley (*chunkGeom)->suppDetJ[1] = geom->suppDetJ[1] ? &geom->suppDetJ[1][Nq * cStart] : NULL; 11320cf1dd8SToby Isaac (*chunkGeom)->isAffine = geom->isAffine; 11420cf1dd8SToby Isaac PetscFunctionReturn(0); 11520cf1dd8SToby Isaac } 11620cf1dd8SToby Isaac 1172b99622eSMatthew G. Knepley /*@C 118*dce8aebaSBarry Smith PetscFEGeomRestoreChunk - Restore the chunk obtained with `PetscFEGeomCreateChunk()` 1192b99622eSMatthew G. Knepley 1202b99622eSMatthew G. Knepley Input Parameters: 121*dce8aebaSBarry Smith + geom - `PetscFEGeom` object 1222b99622eSMatthew G. Knepley . cStart - The first cell in the chunk 1232b99622eSMatthew G. Knepley . cEnd - The first cell not in the chunk 1242b99622eSMatthew G. Knepley - chunkGeom - The chunk of cells 1252b99622eSMatthew G. Knepley 1262b99622eSMatthew G. Knepley Level: intermediate 1272b99622eSMatthew G. Knepley 128*dce8aebaSBarry Smith .seealso: `PetscFEGeom`, `PetscFEGeomGetChunk()`, `PetscFEGeomCreate()` 1292b99622eSMatthew G. Knepley @*/ 130d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 131d71ae5a4SJacob Faibussowitsch { 13220cf1dd8SToby Isaac PetscFunctionBegin; 1339566063dSJacob Faibussowitsch PetscCall(PetscFree(*chunkGeom)); 13420cf1dd8SToby Isaac PetscFunctionReturn(0); 13520cf1dd8SToby Isaac } 13620cf1dd8SToby Isaac 1376587ee25SMatthew G. Knepley /*@C 138*dce8aebaSBarry Smith PetscFEGeomGetPoint - Get the geometry for cell c at point p as a `PetscFEGeom` 1396587ee25SMatthew G. Knepley 1406587ee25SMatthew G. Knepley Input Parameters: 141*dce8aebaSBarry Smith + geom - `PetscFEGeom` object 1426587ee25SMatthew G. Knepley . c - The cell 1436587ee25SMatthew G. Knepley . p - The point 1446587ee25SMatthew G. Knepley - pcoords - The reference coordinates of point p, or NULL 1456587ee25SMatthew G. Knepley 1466587ee25SMatthew G. Knepley Output Parameter: 1476587ee25SMatthew G. Knepley . pgeom - The geometry of cell c at point p 1486587ee25SMatthew G. Knepley 149*dce8aebaSBarry Smith Level: intermediate 150*dce8aebaSBarry Smith 151*dce8aebaSBarry Smith Notes: 152*dce8aebaSBarry Smith For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom, 1536587ee25SMatthew G. Knepley nothing needs to be done with it afterwards. 1546587ee25SMatthew G. Knepley 1556587ee25SMatthew G. Knepley In the affine case, pgeom must have storage for the integration point coordinates in pgeom->v if pcoords is passed in. 1566587ee25SMatthew G. Knepley 157*dce8aebaSBarry Smith .seealso: `PetscFEGeom`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()` 1586587ee25SMatthew G. Knepley @*/ 159d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, const PetscReal pcoords[], PetscFEGeom *pgeom) 160d71ae5a4SJacob Faibussowitsch { 1616587ee25SMatthew G. Knepley const PetscInt dim = geom->dim; 1626587ee25SMatthew G. Knepley const PetscInt dE = geom->dimEmbed; 1636587ee25SMatthew G. Knepley const PetscInt Np = geom->numPoints; 1646587ee25SMatthew G. Knepley 1656587ee25SMatthew G. Knepley PetscFunctionBeginHot; 1666587ee25SMatthew G. Knepley pgeom->dim = dim; 1676587ee25SMatthew G. Knepley pgeom->dimEmbed = dE; 1686587ee25SMatthew G. Knepley //pgeom->isAffine = geom->isAffine; 1696587ee25SMatthew G. Knepley if (geom->isAffine) { 1706587ee25SMatthew G. Knepley if (!p) { 1716587ee25SMatthew G. Knepley pgeom->xi = geom->xi; 1726587ee25SMatthew G. Knepley pgeom->J = &geom->J[c * Np * dE * dE]; 1736587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[c * Np * dE * dE]; 1746587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np]; 1756587ee25SMatthew G. Knepley pgeom->n = geom->n ? &geom->n[c * Np * dE] : NULL; 1766587ee25SMatthew G. Knepley } 177ad540459SPierre Jolivet if (pcoords) CoordinatesRefToReal(dE, dim, pgeom->xi, &geom->v[c * Np * dE], pgeom->J, pcoords, pgeom->v); 1786587ee25SMatthew G. Knepley } else { 1796587ee25SMatthew G. Knepley pgeom->v = &geom->v[(c * Np + p) * dE]; 1806587ee25SMatthew G. Knepley pgeom->J = &geom->J[(c * Np + p) * dE * dE]; 1816587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[(c * Np + p) * dE * dE]; 1826587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np + p]; 1836587ee25SMatthew G. Knepley pgeom->n = geom->n ? &geom->n[(c * Np + p) * dE] : NULL; 1846587ee25SMatthew G. Knepley } 1856587ee25SMatthew G. Knepley PetscFunctionReturn(0); 1866587ee25SMatthew G. Knepley } 1876587ee25SMatthew G. Knepley 1886587ee25SMatthew G. Knepley /*@C 189*dce8aebaSBarry Smith PetscFEGeomGetCellPoint - Get the cell geometry for face f at point p as a `PetscFEGeom` 1906587ee25SMatthew G. Knepley 1916587ee25SMatthew G. Knepley Input Parameters: 192*dce8aebaSBarry Smith + geom - `PetscFEGeom` object 1936587ee25SMatthew G. Knepley . f - The face 1946587ee25SMatthew G. Knepley - p - The point 1956587ee25SMatthew G. Knepley 1966587ee25SMatthew G. Knepley Output Parameter: 1976587ee25SMatthew G. Knepley . pgeom - The cell geometry of face f at point p 1986587ee25SMatthew G. Knepley 1996587ee25SMatthew G. Knepley Level: intermediate 2006587ee25SMatthew G. Knepley 201*dce8aebaSBarry Smith Note: 202*dce8aebaSBarry Smith For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom, 203*dce8aebaSBarry Smith nothing needs to be done with it afterwards. 204*dce8aebaSBarry Smith 205*dce8aebaSBarry Smith .seealso: `PetscFEGeom()`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()` 2066587ee25SMatthew G. Knepley @*/ 207d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, PetscFEGeom *pgeom) 208d71ae5a4SJacob Faibussowitsch { 2095fedec97SMatthew G. Knepley const PetscBool bd = geom->dimEmbed > geom->dim && !geom->isCohesive ? PETSC_TRUE : PETSC_FALSE; 2106587ee25SMatthew G. Knepley const PetscInt dim = bd ? geom->dimEmbed : geom->dim; 2116587ee25SMatthew G. Knepley const PetscInt dE = geom->dimEmbed; 2126587ee25SMatthew G. Knepley const PetscInt Np = geom->numPoints; 2136587ee25SMatthew G. Knepley 2146587ee25SMatthew G. Knepley PetscFunctionBeginHot; 2156587ee25SMatthew G. Knepley pgeom->dim = dim; 2166587ee25SMatthew G. Knepley pgeom->dimEmbed = dE; 2176587ee25SMatthew G. Knepley //pgeom->isAffine = geom->isAffine; 2186587ee25SMatthew G. Knepley if (geom->isAffine) { 2196587ee25SMatthew G. Knepley if (!p) { 2206587ee25SMatthew G. Knepley if (bd) { 2216587ee25SMatthew G. Knepley pgeom->J = &geom->suppJ[0][c * Np * dE * dE]; 2226587ee25SMatthew G. Knepley pgeom->invJ = &geom->suppInvJ[0][c * Np * dE * dE]; 2236587ee25SMatthew G. Knepley pgeom->detJ = &geom->suppDetJ[0][c * Np]; 2246587ee25SMatthew G. Knepley } else { 2256587ee25SMatthew G. Knepley pgeom->J = &geom->J[c * Np * dE * dE]; 2266587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[c * Np * dE * dE]; 2276587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np]; 2286587ee25SMatthew G. Knepley } 2296587ee25SMatthew G. Knepley } 2306587ee25SMatthew G. Knepley } else { 2316587ee25SMatthew G. Knepley if (bd) { 2326587ee25SMatthew G. Knepley pgeom->J = &geom->suppJ[0][(c * Np + p) * dE * dE]; 2336587ee25SMatthew G. Knepley pgeom->invJ = &geom->suppInvJ[0][(c * Np + p) * dE * dE]; 2346587ee25SMatthew G. Knepley pgeom->detJ = &geom->suppDetJ[0][c * Np + p]; 2356587ee25SMatthew G. Knepley } else { 2366587ee25SMatthew G. Knepley pgeom->J = &geom->J[(c * Np + p) * dE * dE]; 2376587ee25SMatthew G. Knepley pgeom->invJ = &geom->invJ[(c * Np + p) * dE * dE]; 2386587ee25SMatthew G. Knepley pgeom->detJ = &geom->detJ[c * Np + p]; 2396587ee25SMatthew G. Knepley } 2406587ee25SMatthew G. Knepley } 2416587ee25SMatthew G. Knepley PetscFunctionReturn(0); 2426587ee25SMatthew G. Knepley } 2436587ee25SMatthew G. Knepley 2442b99622eSMatthew G. Knepley /*@ 245*dce8aebaSBarry Smith PetscFEGeomComplete - Calculate derived quantities from base geometry specification 2462b99622eSMatthew G. Knepley 2472b99622eSMatthew G. Knepley Input Parameter: 248*dce8aebaSBarry Smith . geom - `PetscFEGeom` object 2492b99622eSMatthew G. Knepley 2502b99622eSMatthew G. Knepley Level: intermediate 2512b99622eSMatthew G. Knepley 252*dce8aebaSBarry Smith .seealso: `PetscFEGeom`, `PetscFEGeomCreate()` 2532b99622eSMatthew G. Knepley @*/ 254d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom) 255d71ae5a4SJacob Faibussowitsch { 25620cf1dd8SToby Isaac PetscInt i, j, N, dE; 25720cf1dd8SToby Isaac 25820cf1dd8SToby Isaac PetscFunctionBeginHot; 25920cf1dd8SToby Isaac N = geom->numPoints * geom->numCells; 26020cf1dd8SToby Isaac dE = geom->dimEmbed; 26120cf1dd8SToby Isaac switch (dE) { 26220cf1dd8SToby Isaac case 3: 26320cf1dd8SToby Isaac for (i = 0; i < N; i++) { 26420cf1dd8SToby Isaac DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE * dE * i]); 26520cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE * dE * i], &geom->J[dE * dE * i], geom->detJ[i]); 26620cf1dd8SToby Isaac } 26720cf1dd8SToby Isaac break; 26820cf1dd8SToby Isaac case 2: 26920cf1dd8SToby Isaac for (i = 0; i < N; i++) { 27020cf1dd8SToby Isaac DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE * dE * i]); 27120cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE * dE * i], &geom->J[dE * dE * i], geom->detJ[i]); 27220cf1dd8SToby Isaac } 27320cf1dd8SToby Isaac break; 27420cf1dd8SToby Isaac case 1: 27520cf1dd8SToby Isaac for (i = 0; i < N; i++) { 27620cf1dd8SToby Isaac geom->detJ[i] = PetscAbsReal(geom->J[i]); 27720cf1dd8SToby Isaac if (geom->invJ) geom->invJ[i] = 1. / geom->J[i]; 27820cf1dd8SToby Isaac } 27920cf1dd8SToby Isaac break; 28020cf1dd8SToby Isaac } 28120cf1dd8SToby Isaac if (geom->n) { 28220cf1dd8SToby Isaac for (i = 0; i < N; i++) { 283ad540459SPierre 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.); 28420cf1dd8SToby Isaac } 28520cf1dd8SToby Isaac } 28620cf1dd8SToby Isaac PetscFunctionReturn(0); 28720cf1dd8SToby Isaac } 288