xref: /petsc/src/dm/impls/plex/tutorials/ex18.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
15552b385SBrandon #include "petscsys.h"
27766ff4bSPierre Jolivet static const char help[] = "Test of PETSc/CAD Shape Modification Technology";
35552b385SBrandon 
45552b385SBrandon #include <petscdmplexegads.h>
55552b385SBrandon #include <petsc/private/hashmapi.h>
65552b385SBrandon 
surfArea(DM dm)75552b385SBrandon static PetscErrorCode surfArea(DM dm)
85552b385SBrandon {
95552b385SBrandon   DMLabel     bodyLabel, faceLabel;
105552b385SBrandon   double      surfaceArea = 0., volume = 0.;
115552b385SBrandon   PetscReal   vol, centroid[3], normal[3];
125552b385SBrandon   PetscInt    dim, cStart, cEnd, fStart, fEnd;
135552b385SBrandon   PetscInt    bodyID, faceID;
145552b385SBrandon   MPI_Comm    comm;
155552b385SBrandon   const char *name;
165552b385SBrandon 
175552b385SBrandon   PetscFunctionBeginUser;
185552b385SBrandon   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
195552b385SBrandon   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
205552b385SBrandon   PetscCall(DMGetDimension(dm, &dim));
215552b385SBrandon   PetscCall(PetscPrintf(comm, "    dim = %" PetscInt_FMT "\n", dim));
225552b385SBrandon   PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
235552b385SBrandon   PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
245552b385SBrandon   PetscCall(DMGetCoordinatesLocalSetUp(dm));
255552b385SBrandon 
265552b385SBrandon   if (dim == 2) {
275552b385SBrandon     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
285552b385SBrandon     for (PetscInt ii = cStart; ii < cEnd; ++ii) {
295552b385SBrandon       PetscCall(DMLabelGetValue(faceLabel, ii, &faceID));
305552b385SBrandon       if (faceID >= 0) {
315552b385SBrandon         PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal));
325552b385SBrandon         surfaceArea += vol;
335552b385SBrandon       }
345552b385SBrandon     }
355552b385SBrandon   }
365552b385SBrandon 
375552b385SBrandon   if (dim == 3) {
385552b385SBrandon     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
395552b385SBrandon     for (PetscInt ii = fStart; ii < fEnd; ++ii) {
405552b385SBrandon       PetscCall(DMLabelGetValue(faceLabel, ii, &faceID));
415552b385SBrandon       if (faceID >= 0) {
425552b385SBrandon         PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal));
435552b385SBrandon         surfaceArea += vol;
445552b385SBrandon       }
455552b385SBrandon     }
465552b385SBrandon 
475552b385SBrandon     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
485552b385SBrandon     for (PetscInt ii = cStart; ii < cEnd; ++ii) {
495552b385SBrandon       PetscCall(DMLabelGetValue(bodyLabel, ii, &bodyID));
505552b385SBrandon       if (bodyID >= 0) {
515552b385SBrandon         PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal));
525552b385SBrandon         volume += vol;
535552b385SBrandon       }
545552b385SBrandon     }
555552b385SBrandon   }
565552b385SBrandon 
575552b385SBrandon   if (dim == 2) {
585552b385SBrandon     PetscCall(PetscPrintf(comm, "%s Surface Area = %.6e \n\n", name, (double)surfaceArea));
595552b385SBrandon   } else if (dim == 3) {
605552b385SBrandon     PetscCall(PetscPrintf(comm, "%s Volume = %.6e \n", name, (double)volume));
615552b385SBrandon     PetscCall(PetscPrintf(comm, "%s Surface Area = %.6e \n\n", name, (double)surfaceArea));
625552b385SBrandon   }
635552b385SBrandon   PetscFunctionReturn(PETSC_SUCCESS);
645552b385SBrandon }
655552b385SBrandon 
main(int argc,char * argv[])665552b385SBrandon int main(int argc, char *argv[])
675552b385SBrandon {
685552b385SBrandon   /* EGADSlite variables */
695552b385SBrandon   PetscGeom model, *bodies, *fobjs;
705552b385SBrandon   int       Nb, Nf, id;
715552b385SBrandon   /* PETSc variables */
725552b385SBrandon   DM             dm;
735552b385SBrandon   MPI_Comm       comm;
745552b385SBrandon   PetscContainer modelObj, cpHashTableObj, wHashTableObj, cpCoordDataLengthObj, wDataLengthObj;
755552b385SBrandon   Vec            cntrlPtCoordsVec, cntrlPtWeightsVec;
765552b385SBrandon   PetscScalar   *cpCoordData, *wData;
775552b385SBrandon   PetscInt       cpCoordDataLength, wDataLength;
785552b385SBrandon   PetscInt      *cpCoordDataLengthPtr, *wDataLengthPtr;
795552b385SBrandon   PetscHMapI     cpHashTable, wHashTable;
805552b385SBrandon   PetscInt       Nr = 2;
815552b385SBrandon 
825552b385SBrandon   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
835552b385SBrandon   comm = PETSC_COMM_WORLD;
845552b385SBrandon   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
855552b385SBrandon   PetscCall(DMSetType(dm, DMPLEX));
865552b385SBrandon   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
875552b385SBrandon   PetscCall(DMSetFromOptions(dm));
885552b385SBrandon 
895552b385SBrandon   PetscCall(PetscObjectSetName((PetscObject)dm, "Original Surface"));
905552b385SBrandon   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
915552b385SBrandon   PetscCall(surfArea(dm));
925552b385SBrandon 
935552b385SBrandon   // Expose Geometry Definition Data and Calculate Surface Gradients
945552b385SBrandon   PetscCall(DMPlexGeomDataAndGrads(dm, PETSC_FALSE));
955552b385SBrandon 
965552b385SBrandon   // Get Attached EGADS model
975552b385SBrandon   PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
985552b385SBrandon   if (!modelObj) PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
995552b385SBrandon 
1005552b385SBrandon   // Get attached EGADS model (pointer)
101*2a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(modelObj, &model));
1025552b385SBrandon 
1035552b385SBrandon   // Look to see if DM has Container for Geometry Control Point Data
1045552b385SBrandon   PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Hash Table", (PetscObject *)&cpHashTableObj));
1055552b385SBrandon   PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Coordinates", (PetscObject *)&cntrlPtCoordsVec));
1065552b385SBrandon   PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Coordinate Data Length", (PetscObject *)&cpCoordDataLengthObj));
1075552b385SBrandon   PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weights Hash Table", (PetscObject *)&wHashTableObj));
1085552b385SBrandon   PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data", (PetscObject *)&cntrlPtWeightsVec));
1095552b385SBrandon   PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data Length", (PetscObject *)&wDataLengthObj));
1105552b385SBrandon 
1115552b385SBrandon   // Get attached EGADS model Control Point and Weights Hash Tables and Data Arrays (pointer)
112*2a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(cpHashTableObj, &cpHashTable));
1135552b385SBrandon   PetscCall(VecGetArrayWrite(cntrlPtCoordsVec, &cpCoordData));
114*2a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(cpCoordDataLengthObj, &cpCoordDataLengthPtr));
115*2a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(wHashTableObj, &wHashTable));
1165552b385SBrandon   PetscCall(VecGetArrayWrite(cntrlPtWeightsVec, &wData));
117*2a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(wDataLengthObj, &wDataLengthPtr));
1185552b385SBrandon 
1195552b385SBrandon   cpCoordDataLength = *cpCoordDataLengthPtr;
1205552b385SBrandon   wDataLength       = *wDataLengthPtr;
1215552b385SBrandon   PetscCheck(cpCoordDataLength && wDataLength, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Data sizes must be positive");
1225552b385SBrandon 
1235552b385SBrandon   // Get the number of bodies and body objects in the model
1245552b385SBrandon   PetscCall(DMPlexGetGeomModelBodies(dm, &bodies, &Nb));
1255552b385SBrandon 
1265552b385SBrandon   // Get all FACES of the Body
1275552b385SBrandon   PetscGeom body = bodies[0];
1285552b385SBrandon   PetscCall(DMPlexGetGeomModelBodyFaces(dm, body, &fobjs, &Nf));
1295552b385SBrandon 
1305552b385SBrandon   // Update Control Point and Weight definitions for each surface
1315552b385SBrandon   for (PetscInt jj = 0; jj < Nf; ++jj) {
1325552b385SBrandon     PetscGeom face         = fobjs[jj];
1335552b385SBrandon     PetscInt  numCntrlPnts = 0;
1345552b385SBrandon 
1355552b385SBrandon     PetscCall(DMPlexGetGeomID(dm, body, face, &id));
1365552b385SBrandon     PetscCall(DMPlexGetGeomFaceNumOfControlPoints(dm, face, &numCntrlPnts));
1375552b385SBrandon 
1385552b385SBrandon     // Update Control Points
1395552b385SBrandon     PetscHashIter CPiter, Witer;
1405552b385SBrandon     PetscBool     CPfound, Wfound;
1415552b385SBrandon     PetscInt      faceCPStartRow, faceWStartRow;
1425552b385SBrandon 
1435552b385SBrandon     PetscCall(PetscHMapIFind(cpHashTable, id, &CPiter, &CPfound));
1445552b385SBrandon     PetscCheck(CPfound, comm, PETSC_ERR_SUP, "FACE ID not found in Control Point Hash Table");
1455552b385SBrandon     PetscCall(PetscHMapIGet(cpHashTable, id, &faceCPStartRow));
1465552b385SBrandon 
1475552b385SBrandon     PetscCall(PetscHMapIFind(wHashTable, id, &Witer, &Wfound));
1485552b385SBrandon     PetscCheck(Wfound, comm, PETSC_ERR_SUP, "FACE ID not found in Control Point Weights Hash Table");
1495552b385SBrandon     PetscCall(PetscHMapIGet(wHashTable, id, &faceWStartRow));
1505552b385SBrandon 
1515552b385SBrandon     for (PetscInt ii = 0; ii < numCntrlPnts; ++ii) {
1525552b385SBrandon       if (ii == 4) {
1535552b385SBrandon         // Update Control Points - Change the location of the center control point of the faces
1545552b385SBrandon         // Note :: Modification geometry requires knowledge of how the geometry is defined.
1555552b385SBrandon         cpCoordData[faceCPStartRow + (3 * ii) + 0] = 2.0 * cpCoordData[faceCPStartRow + (3 * ii) + 0];
1565552b385SBrandon         cpCoordData[faceCPStartRow + (3 * ii) + 1] = 2.0 * cpCoordData[faceCPStartRow + (3 * ii) + 1];
1575552b385SBrandon         cpCoordData[faceCPStartRow + (3 * ii) + 2] = 2.0 * cpCoordData[faceCPStartRow + (3 * ii) + 2];
1585552b385SBrandon       } else {
1595552b385SBrandon         // Do Nothing
1605552b385SBrandon         // Note: Could use section the change location of other face control points.
1615552b385SBrandon       }
1625552b385SBrandon     }
1635552b385SBrandon   }
1645552b385SBrandon   PetscCall(DMPlexFreeGeomObject(dm, fobjs));
1655552b385SBrandon 
1665552b385SBrandon   // Modify Control Points of Geometry
1675552b385SBrandon   PetscCall(PetscObjectSetName((PetscObject)dm, "Modified Surface"));
1685552b385SBrandon   // TODO Wrap EG_saveModel() in a Plex viewer to manage file access
1695552b385SBrandon   PetscCall(DMPlexModifyGeomModel(dm, comm, cpCoordData, wData, PETSC_FALSE, PETSC_TRUE, "newModel.stp"));
1705552b385SBrandon   PetscCall(VecRestoreArrayWrite(cntrlPtCoordsVec, &cpCoordData));
1715552b385SBrandon   PetscCall(VecRestoreArrayWrite(cntrlPtWeightsVec, &wData));
1725552b385SBrandon 
1735552b385SBrandon   // Inflate Mesh to Geometry
1745552b385SBrandon   PetscCall(DMPlexInflateToGeomModel(dm, PETSC_FALSE));
1755552b385SBrandon   PetscCall(surfArea(dm));
1765552b385SBrandon 
1775552b385SBrandon   // Output .hdf5 file
1785552b385SBrandon   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1795552b385SBrandon 
1805552b385SBrandon   for (PetscInt r = 0; r < Nr; ++r) {
1815552b385SBrandon     char name[PETSC_MAX_PATH_LEN];
1825552b385SBrandon 
1835552b385SBrandon     // Perform refinement on the Mesh attached to the new geometry
1845552b385SBrandon     PetscCall(DMSetFromOptions(dm));
1855552b385SBrandon     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Modified Surface refinement %" PetscInt_FMT, r));
1865552b385SBrandon     PetscCall(PetscObjectSetName((PetscObject)dm, name));
1875552b385SBrandon     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1885552b385SBrandon     PetscCall(surfArea(dm));
1895552b385SBrandon   }
1905552b385SBrandon 
1915552b385SBrandon   PetscCall(DMDestroy(&dm));
1925552b385SBrandon   PetscCall(PetscFinalize());
1935552b385SBrandon   return 0;
1945552b385SBrandon }
1955552b385SBrandon 
1965552b385SBrandon /*TEST
1975552b385SBrandon   build:
1985552b385SBrandon     requires: egads
1995552b385SBrandon 
2005552b385SBrandon   # Use -dm_view hdf5:mesh_shapeMod_sphere.h5 to view the mesh
2015552b385SBrandon   test:
2025552b385SBrandon     suffix: sphere_shapeMod
2035552b385SBrandon     requires: datafilespath
2045552b385SBrandon     temporaries: newModel.stp
2055552b385SBrandon     args: -dm_plex_filename ${DATAFILESPATH}/meshes/cad/sphere_example.stp \
2065552b385SBrandon           -dm_refine -dm_plex_geom_print_model 1 -dm_plex_geom_shape_opt 1
2075552b385SBrandon 
2085552b385SBrandon TEST*/
209