xref: /petsc/src/dm/dt/fe/interface/fegeom.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 @*/
199371c9d4SSatish Balay PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom) {
2020cf1dd8SToby Isaac   PetscFEGeom     *g;
217489efa5SMatthew G. Knepley   PetscInt         dim, Nq, N;
2220cf1dd8SToby Isaac   const PetscReal *p;
2320cf1dd8SToby Isaac 
2420cf1dd8SToby Isaac   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &dim, NULL, &Nq, &p, NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscNew(&g));
2720cf1dd8SToby Isaac   g->xi         = p;
2820cf1dd8SToby Isaac   g->numCells   = numCells;
297489efa5SMatthew G. Knepley   g->numPoints  = Nq;
3020cf1dd8SToby Isaac   g->dim        = dim;
3120cf1dd8SToby Isaac   g->dimEmbed   = dimEmbed;
325fedec97SMatthew G. Knepley   g->isCohesive = PETSC_FALSE;
337489efa5SMatthew G. Knepley   N             = numCells * Nq;
349566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ));
3520cf1dd8SToby Isaac   if (faceData) {
369566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(numCells, &g->face, N * dimEmbed, &g->n));
379371c9d4SSatish 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])));
3820cf1dd8SToby Isaac   }
399566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ));
4020cf1dd8SToby Isaac   *geom = g;
4120cf1dd8SToby Isaac   PetscFunctionReturn(0);
4220cf1dd8SToby Isaac }
4320cf1dd8SToby Isaac 
442b99622eSMatthew G. Knepley /*@C
452b99622eSMatthew G. Knepley   PetscFEGeomDestroy - Destroy a PetscFEGeom object
462b99622eSMatthew G. Knepley 
472b99622eSMatthew G. Knepley   Input Parameter:
482b99622eSMatthew G. Knepley . geom - PetscFEGeom object
492b99622eSMatthew G. Knepley 
502b99622eSMatthew G. Knepley   Level: beginner
512b99622eSMatthew G. Knepley 
52db781477SPatrick Sanan .seealso: `PetscFEGeomCreate()`
532b99622eSMatthew G. Knepley @*/
549371c9d4SSatish Balay PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom) {
5520cf1dd8SToby Isaac   PetscFunctionBegin;
5620cf1dd8SToby Isaac   if (!*geom) PetscFunctionReturn(0);
579566063dSJacob Faibussowitsch   PetscCall(PetscFree3((*geom)->v, (*geom)->J, (*geom)->detJ));
589566063dSJacob Faibussowitsch   PetscCall(PetscFree((*geom)->invJ));
599566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*geom)->face, (*geom)->n));
609566063dSJacob Faibussowitsch   PetscCall(PetscFree6((*geom)->suppJ[0], (*geom)->suppJ[1], (*geom)->suppInvJ[0], (*geom)->suppInvJ[1], (*geom)->suppDetJ[0], (*geom)->suppDetJ[1]));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(*geom));
6220cf1dd8SToby Isaac   PetscFunctionReturn(0);
6320cf1dd8SToby Isaac }
6420cf1dd8SToby Isaac 
652b99622eSMatthew G. Knepley /*@C
662b99622eSMatthew G. Knepley   PetscFEGeomGetChunk - Get a chunk of cells in the group as a PetscFEGeom
672b99622eSMatthew G. Knepley 
682b99622eSMatthew G. Knepley   Input Parameters:
692b99622eSMatthew G. Knepley + geom   - PetscFEGeom object
702b99622eSMatthew G. Knepley . cStart - The first cell in the chunk
712b99622eSMatthew G. Knepley - cEnd   - The first cell not in the chunk
722b99622eSMatthew G. Knepley 
732b99622eSMatthew G. Knepley   Output Parameter:
742b99622eSMatthew G. Knepley . chunkGeom - The chunk of cells
752b99622eSMatthew G. Knepley 
762b99622eSMatthew G. Knepley   Level: intermediate
772b99622eSMatthew G. Knepley 
78db781477SPatrick Sanan .seealso: `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()`
792b99622eSMatthew G. Knepley @*/
809371c9d4SSatish Balay PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) {
8120cf1dd8SToby Isaac   PetscInt Nq;
8220cf1dd8SToby Isaac   PetscInt dE;
8320cf1dd8SToby Isaac 
8420cf1dd8SToby Isaac   PetscFunctionBegin;
8520cf1dd8SToby Isaac   PetscValidPointer(geom, 1);
86064a246eSJacob Faibussowitsch   PetscValidPointer(chunkGeom, 4);
87*48a46eb9SPierre Jolivet   if (!(*chunkGeom)) PetscCall(PetscNew(chunkGeom));
8820cf1dd8SToby Isaac   Nq                        = geom->numPoints;
8920cf1dd8SToby Isaac   dE                        = geom->dimEmbed;
9020cf1dd8SToby Isaac   (*chunkGeom)->dim         = geom->dim;
9120cf1dd8SToby Isaac   (*chunkGeom)->dimEmbed    = geom->dimEmbed;
9220cf1dd8SToby Isaac   (*chunkGeom)->numPoints   = geom->numPoints;
9320cf1dd8SToby Isaac   (*chunkGeom)->numCells    = cEnd - cStart;
9420cf1dd8SToby Isaac   (*chunkGeom)->xi          = geom->xi;
9520cf1dd8SToby Isaac   (*chunkGeom)->v           = &geom->v[Nq * dE * cStart];
9620cf1dd8SToby Isaac   (*chunkGeom)->J           = &geom->J[Nq * dE * dE * cStart];
9720cf1dd8SToby Isaac   (*chunkGeom)->invJ        = (geom->invJ) ? &geom->invJ[Nq * dE * dE * cStart] : NULL;
9820cf1dd8SToby Isaac   (*chunkGeom)->detJ        = &geom->detJ[Nq * cStart];
9920cf1dd8SToby Isaac   (*chunkGeom)->n           = geom->n ? &geom->n[Nq * dE * cStart] : NULL;
10020cf1dd8SToby Isaac   (*chunkGeom)->face        = geom->face ? &geom->face[cStart] : NULL;
10141418987SMatthew G. Knepley   (*chunkGeom)->suppJ[0]    = geom->suppJ[0] ? &geom->suppJ[0][Nq * dE * dE * cStart] : NULL;
10241418987SMatthew G. Knepley   (*chunkGeom)->suppJ[1]    = geom->suppJ[1] ? &geom->suppJ[1][Nq * dE * dE * cStart] : NULL;
10320cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq * dE * dE * cStart] : NULL;
10420cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq * dE * dE * cStart] : NULL;
10541418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[0] = geom->suppDetJ[0] ? &geom->suppDetJ[0][Nq * cStart] : NULL;
10641418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[1] = geom->suppDetJ[1] ? &geom->suppDetJ[1][Nq * cStart] : NULL;
10720cf1dd8SToby Isaac   (*chunkGeom)->isAffine    = geom->isAffine;
10820cf1dd8SToby Isaac   PetscFunctionReturn(0);
10920cf1dd8SToby Isaac }
11020cf1dd8SToby Isaac 
1112b99622eSMatthew G. Knepley /*@C
1122b99622eSMatthew G. Knepley   PetscFEGeomRestoreChunk - Restore the chunk
1132b99622eSMatthew G. Knepley 
1142b99622eSMatthew G. Knepley   Input Parameters:
1152b99622eSMatthew G. Knepley + geom      - PetscFEGeom object
1162b99622eSMatthew G. Knepley . cStart    - The first cell in the chunk
1172b99622eSMatthew G. Knepley . cEnd      - The first cell not in the chunk
1182b99622eSMatthew G. Knepley - chunkGeom - The chunk of cells
1192b99622eSMatthew G. Knepley 
1202b99622eSMatthew G. Knepley   Level: intermediate
1212b99622eSMatthew G. Knepley 
122db781477SPatrick Sanan .seealso: `PetscFEGeomGetChunk()`, `PetscFEGeomCreate()`
1232b99622eSMatthew G. Knepley @*/
1249371c9d4SSatish Balay PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) {
12520cf1dd8SToby Isaac   PetscFunctionBegin;
1269566063dSJacob Faibussowitsch   PetscCall(PetscFree(*chunkGeom));
12720cf1dd8SToby Isaac   PetscFunctionReturn(0);
12820cf1dd8SToby Isaac }
12920cf1dd8SToby Isaac 
1306587ee25SMatthew G. Knepley /*@C
1316587ee25SMatthew G. Knepley   PetscFEGeomGetPoint - Get the geometry for cell c at point p as a PetscFEGeom
1326587ee25SMatthew G. Knepley 
1336587ee25SMatthew G. Knepley   Input Parameters:
1346587ee25SMatthew G. Knepley + geom    - PetscFEGeom object
1356587ee25SMatthew G. Knepley . c       - The cell
1366587ee25SMatthew G. Knepley . p       - The point
1376587ee25SMatthew G. Knepley - pcoords - The reference coordinates of point p, or NULL
1386587ee25SMatthew G. Knepley 
1396587ee25SMatthew G. Knepley   Output Parameter:
1406587ee25SMatthew G. Knepley . pgeom - The geometry of cell c at point p
1416587ee25SMatthew G. Knepley 
1426587ee25SMatthew G. Knepley   Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom,
1436587ee25SMatthew G. Knepley   nothing needs to be done with it afterwards.
1446587ee25SMatthew G. Knepley 
1456587ee25SMatthew G. Knepley   In the affine case, pgeom must have storage for the integration point coordinates in pgeom->v if pcoords is passed in.
1466587ee25SMatthew G. Knepley 
1476587ee25SMatthew G. Knepley   Level: intermediate
1486587ee25SMatthew G. Knepley 
149db781477SPatrick Sanan .seealso: `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()`
1506587ee25SMatthew G. Knepley @*/
1519371c9d4SSatish Balay PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, const PetscReal pcoords[], PetscFEGeom *pgeom) {
1526587ee25SMatthew G. Knepley   const PetscInt dim = geom->dim;
1536587ee25SMatthew G. Knepley   const PetscInt dE  = geom->dimEmbed;
1546587ee25SMatthew G. Knepley   const PetscInt Np  = geom->numPoints;
1556587ee25SMatthew G. Knepley 
1566587ee25SMatthew G. Knepley   PetscFunctionBeginHot;
1576587ee25SMatthew G. Knepley   pgeom->dim      = dim;
1586587ee25SMatthew G. Knepley   pgeom->dimEmbed = dE;
1596587ee25SMatthew G. Knepley   //pgeom->isAffine = geom->isAffine;
1606587ee25SMatthew G. Knepley   if (geom->isAffine) {
1616587ee25SMatthew G. Knepley     if (!p) {
1626587ee25SMatthew G. Knepley       pgeom->xi   = geom->xi;
1636587ee25SMatthew G. Knepley       pgeom->J    = &geom->J[c * Np * dE * dE];
1646587ee25SMatthew G. Knepley       pgeom->invJ = &geom->invJ[c * Np * dE * dE];
1656587ee25SMatthew G. Knepley       pgeom->detJ = &geom->detJ[c * Np];
1666587ee25SMatthew G. Knepley       pgeom->n    = geom->n ? &geom->n[c * Np * dE] : NULL;
1676587ee25SMatthew G. Knepley     }
1686587ee25SMatthew G. Knepley     if (pcoords) { CoordinatesRefToReal(dE, dim, pgeom->xi, &geom->v[c * Np * dE], pgeom->J, pcoords, pgeom->v); }
1696587ee25SMatthew G. Knepley   } else {
1706587ee25SMatthew G. Knepley     pgeom->v    = &geom->v[(c * Np + p) * dE];
1716587ee25SMatthew G. Knepley     pgeom->J    = &geom->J[(c * Np + p) * dE * dE];
1726587ee25SMatthew G. Knepley     pgeom->invJ = &geom->invJ[(c * Np + p) * dE * dE];
1736587ee25SMatthew G. Knepley     pgeom->detJ = &geom->detJ[c * Np + p];
1746587ee25SMatthew G. Knepley     pgeom->n    = geom->n ? &geom->n[(c * Np + p) * dE] : NULL;
1756587ee25SMatthew G. Knepley   }
1766587ee25SMatthew G. Knepley   PetscFunctionReturn(0);
1776587ee25SMatthew G. Knepley }
1786587ee25SMatthew G. Knepley 
1796587ee25SMatthew G. Knepley /*@C
1806587ee25SMatthew G. Knepley   PetscFEGeomGetCellPoint - Get the cell geometry for face f at point p as a PetscFEGeom
1816587ee25SMatthew G. Knepley 
1826587ee25SMatthew G. Knepley   Input Parameters:
1836587ee25SMatthew G. Knepley + geom    - PetscFEGeom object
1846587ee25SMatthew G. Knepley . f       - The face
1856587ee25SMatthew G. Knepley - p       - The point
1866587ee25SMatthew G. Knepley 
1876587ee25SMatthew G. Knepley   Output Parameter:
1886587ee25SMatthew G. Knepley . pgeom - The cell geometry of face f at point p
1896587ee25SMatthew G. Knepley 
1906587ee25SMatthew G. Knepley   Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom,
1916587ee25SMatthew G. Knepley   nothing needs to be done with it afterwards.
1926587ee25SMatthew G. Knepley 
1936587ee25SMatthew G. Knepley   Level: intermediate
1946587ee25SMatthew G. Knepley 
195db781477SPatrick Sanan .seealso: `PetscFEGeomRestoreChunk()`, `PetscFEGeomCreate()`
1966587ee25SMatthew G. Knepley @*/
1979371c9d4SSatish Balay PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, PetscFEGeom *pgeom) {
1985fedec97SMatthew G. Knepley   const PetscBool bd  = geom->dimEmbed > geom->dim && !geom->isCohesive ? PETSC_TRUE : PETSC_FALSE;
1996587ee25SMatthew G. Knepley   const PetscInt  dim = bd ? geom->dimEmbed : geom->dim;
2006587ee25SMatthew G. Knepley   const PetscInt  dE  = geom->dimEmbed;
2016587ee25SMatthew G. Knepley   const PetscInt  Np  = geom->numPoints;
2026587ee25SMatthew G. Knepley 
2036587ee25SMatthew G. Knepley   PetscFunctionBeginHot;
2046587ee25SMatthew G. Knepley   pgeom->dim      = dim;
2056587ee25SMatthew G. Knepley   pgeom->dimEmbed = dE;
2066587ee25SMatthew G. Knepley   //pgeom->isAffine = geom->isAffine;
2076587ee25SMatthew G. Knepley   if (geom->isAffine) {
2086587ee25SMatthew G. Knepley     if (!p) {
2096587ee25SMatthew G. Knepley       if (bd) {
2106587ee25SMatthew G. Knepley         pgeom->J    = &geom->suppJ[0][c * Np * dE * dE];
2116587ee25SMatthew G. Knepley         pgeom->invJ = &geom->suppInvJ[0][c * Np * dE * dE];
2126587ee25SMatthew G. Knepley         pgeom->detJ = &geom->suppDetJ[0][c * Np];
2136587ee25SMatthew G. Knepley       } else {
2146587ee25SMatthew G. Knepley         pgeom->J    = &geom->J[c * Np * dE * dE];
2156587ee25SMatthew G. Knepley         pgeom->invJ = &geom->invJ[c * Np * dE * dE];
2166587ee25SMatthew G. Knepley         pgeom->detJ = &geom->detJ[c * Np];
2176587ee25SMatthew G. Knepley       }
2186587ee25SMatthew G. Knepley     }
2196587ee25SMatthew G. Knepley   } else {
2206587ee25SMatthew G. Knepley     if (bd) {
2216587ee25SMatthew G. Knepley       pgeom->J    = &geom->suppJ[0][(c * Np + p) * dE * dE];
2226587ee25SMatthew G. Knepley       pgeom->invJ = &geom->suppInvJ[0][(c * Np + p) * dE * dE];
2236587ee25SMatthew G. Knepley       pgeom->detJ = &geom->suppDetJ[0][c * Np + p];
2246587ee25SMatthew G. Knepley     } else {
2256587ee25SMatthew G. Knepley       pgeom->J    = &geom->J[(c * Np + p) * dE * dE];
2266587ee25SMatthew G. Knepley       pgeom->invJ = &geom->invJ[(c * Np + p) * dE * dE];
2276587ee25SMatthew G. Knepley       pgeom->detJ = &geom->detJ[c * Np + p];
2286587ee25SMatthew G. Knepley     }
2296587ee25SMatthew G. Knepley   }
2306587ee25SMatthew G. Knepley   PetscFunctionReturn(0);
2316587ee25SMatthew G. Knepley }
2326587ee25SMatthew G. Knepley 
2332b99622eSMatthew G. Knepley /*@
2342b99622eSMatthew G. Knepley   PetscFEGeomComplete - Calculate derived quntites from base geometry specification
2352b99622eSMatthew G. Knepley 
2362b99622eSMatthew G. Knepley   Input Parameter:
2372b99622eSMatthew G. Knepley . geom - PetscFEGeom object
2382b99622eSMatthew G. Knepley 
2392b99622eSMatthew G. Knepley   Level: intermediate
2402b99622eSMatthew G. Knepley 
241db781477SPatrick Sanan .seealso: `PetscFEGeomCreate()`
2422b99622eSMatthew G. Knepley @*/
2439371c9d4SSatish Balay PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom) {
24420cf1dd8SToby Isaac   PetscInt i, j, N, dE;
24520cf1dd8SToby Isaac 
24620cf1dd8SToby Isaac   PetscFunctionBeginHot;
24720cf1dd8SToby Isaac   N  = geom->numPoints * geom->numCells;
24820cf1dd8SToby Isaac   dE = geom->dimEmbed;
24920cf1dd8SToby Isaac   switch (dE) {
25020cf1dd8SToby Isaac   case 3:
25120cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
25220cf1dd8SToby Isaac       DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE * dE * i]);
25320cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE * dE * i], &geom->J[dE * dE * i], geom->detJ[i]);
25420cf1dd8SToby Isaac     }
25520cf1dd8SToby Isaac     break;
25620cf1dd8SToby Isaac   case 2:
25720cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
25820cf1dd8SToby Isaac       DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE * dE * i]);
25920cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE * dE * i], &geom->J[dE * dE * i], geom->detJ[i]);
26020cf1dd8SToby Isaac     }
26120cf1dd8SToby Isaac     break;
26220cf1dd8SToby Isaac   case 1:
26320cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
26420cf1dd8SToby Isaac       geom->detJ[i] = PetscAbsReal(geom->J[i]);
26520cf1dd8SToby Isaac       if (geom->invJ) geom->invJ[i] = 1. / geom->J[i];
26620cf1dd8SToby Isaac     }
26720cf1dd8SToby Isaac     break;
26820cf1dd8SToby Isaac   }
26920cf1dd8SToby Isaac   if (geom->n) {
27020cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
2719371c9d4SSatish Balay       for (j = 0; j < dE; j++) { geom->n[dE * i + j] = geom->J[dE * dE * i + dE * j + dE - 1] * ((dE == 2) ? -1. : 1.); }
27220cf1dd8SToby Isaac     }
27320cf1dd8SToby Isaac   }
27420cf1dd8SToby Isaac   PetscFunctionReturn(0);
27520cf1dd8SToby Isaac }
276