15552b385SBrandon static const char help[] = "Test of PETSc CAD Shape Optimization & Mesh Modification Technology"; 25552b385SBrandon 35552b385SBrandon #include <petscdmplexegads.h> 45552b385SBrandon #include <petsc/private/hashmapi.h> 55552b385SBrandon 65552b385SBrandon typedef struct { 75552b385SBrandon char filename[PETSC_MAX_PATH_LEN]; 85552b385SBrandon PetscInt saMaxIter; // Maximum number of iterations of shape optimization loop 95552b385SBrandon } AppCtx; 105552b385SBrandon 115552b385SBrandon static PetscErrorCode surfArea(DM dm) 125552b385SBrandon { 135552b385SBrandon DMLabel bodyLabel, faceLabel; 145552b385SBrandon double surfaceArea = 0., volume = 0.; 155552b385SBrandon PetscReal vol, centroid[3], normal[3]; 165552b385SBrandon PetscInt dim, cStart, cEnd, fStart, fEnd; 175552b385SBrandon PetscInt bodyID, faceID; 185552b385SBrandon MPI_Comm comm; 195552b385SBrandon const char *name; 205552b385SBrandon 215552b385SBrandon PetscFunctionBeginUser; 225552b385SBrandon PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 235552b385SBrandon PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 245552b385SBrandon PetscCall(DMGetDimension(dm, &dim)); 255552b385SBrandon PetscCall(PetscPrintf(comm, " dim = %" PetscInt_FMT "\n", dim)); 265552b385SBrandon PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 275552b385SBrandon PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 285552b385SBrandon 295552b385SBrandon if (dim == 2) { 305552b385SBrandon PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 315552b385SBrandon for (PetscInt ii = cStart; ii < cEnd; ++ii) { 325552b385SBrandon PetscCall(DMLabelGetValue(faceLabel, ii, &faceID)); 335552b385SBrandon if (faceID >= 0) { 345552b385SBrandon PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal)); 355552b385SBrandon surfaceArea += vol; 365552b385SBrandon } 375552b385SBrandon } 385552b385SBrandon } 395552b385SBrandon 405552b385SBrandon if (dim == 3) { 415552b385SBrandon PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 425552b385SBrandon for (PetscInt ii = fStart; ii < fEnd; ++ii) { 435552b385SBrandon PetscCall(DMLabelGetValue(faceLabel, ii, &faceID)); 445552b385SBrandon if (faceID >= 0) { 455552b385SBrandon PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal)); 465552b385SBrandon surfaceArea += vol; 475552b385SBrandon } 485552b385SBrandon } 495552b385SBrandon 505552b385SBrandon PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 515552b385SBrandon for (PetscInt ii = cStart; ii < cEnd; ++ii) { 525552b385SBrandon PetscCall(DMLabelGetValue(bodyLabel, ii, &bodyID)); 535552b385SBrandon if (bodyID >= 0) { 545552b385SBrandon PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal)); 555552b385SBrandon volume += vol; 565552b385SBrandon } 575552b385SBrandon } 585552b385SBrandon } 595552b385SBrandon 605552b385SBrandon if (dim == 2) { 615552b385SBrandon PetscCall(PetscPrintf(comm, "%s Surface Area = %.6e \n\n", name, (double)surfaceArea)); 625552b385SBrandon } else if (dim == 3) { 635552b385SBrandon PetscCall(PetscPrintf(comm, "%s Volume = %.6e \n", name, (double)volume)); 645552b385SBrandon PetscCall(PetscPrintf(comm, "%s Surface Area = %.6e \n\n", name, (double)surfaceArea)); 655552b385SBrandon } 665552b385SBrandon PetscFunctionReturn(PETSC_SUCCESS); 675552b385SBrandon } 685552b385SBrandon 695552b385SBrandon static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 705552b385SBrandon { 715552b385SBrandon PetscFunctionBeginUser; 725552b385SBrandon options->filename[0] = '\0'; 735552b385SBrandon options->saMaxIter = 200; 745552b385SBrandon 755552b385SBrandon PetscOptionsBegin(comm, "", "DMPlex w/CAD Options", "DMPlex w/CAD"); 765552b385SBrandon PetscCall(PetscOptionsString("-filename", "The CAD/Geometry file", __FILE__, options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL)); 775552b385SBrandon PetscCall(PetscOptionsBoundedInt("-sa_max_iter", "The maximum number of iterates for the shape optimization loop", __FILE__, options->saMaxIter, &options->saMaxIter, NULL, 0)); 785552b385SBrandon PetscOptionsEnd(); 795552b385SBrandon PetscFunctionReturn(0); 805552b385SBrandon } 815552b385SBrandon 825552b385SBrandon int main(int argc, char *argv[]) 835552b385SBrandon { 845552b385SBrandon /* EGADS/EGADSlite variables */ 855552b385SBrandon PetscGeom model, *bodies, *fobjs; 865552b385SBrandon int Nb, Nf, id; 875552b385SBrandon /* PETSc variables */ 885552b385SBrandon DM dmNozzle = NULL; 895552b385SBrandon MPI_Comm comm; 905552b385SBrandon AppCtx ctx; 915552b385SBrandon PetscScalar equivR = 0.0; 925552b385SBrandon const char baseName[] = "Nozzle_Mesh"; 935552b385SBrandon 945552b385SBrandon PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 955552b385SBrandon comm = PETSC_COMM_WORLD; 965552b385SBrandon PetscCall(ProcessOptions(comm, &ctx)); 975552b385SBrandon 985552b385SBrandon PetscCall(DMPlexCreateFromFile(comm, ctx.filename, "EGADS", PETSC_TRUE, &dmNozzle)); 995552b385SBrandon PetscCall(PetscObjectSetName((PetscObject)dmNozzle, baseName)); 1005552b385SBrandon //PetscCall(DMCreate(PETSC_COMM_WORLD, &dmNozzle)); 1015552b385SBrandon //PetscCall(DMSetType(dmNozzle, DMPLEX)); 1025552b385SBrandon //PetscCall(DMPlexDistributeSetDefault(dmNozzle, PETSC_FALSE)); 1035552b385SBrandon //PetscCall(DMSetFromOptions(dmNozzle)); 1045552b385SBrandon 1055552b385SBrandon // Get Common Viewer to store all Mesh Results 1065552b385SBrandon PetscViewer viewer; 1075552b385SBrandon PetscViewerFormat format; 1085552b385SBrandon PetscBool flg; 1095552b385SBrandon PetscInt num = 0; 1105552b385SBrandon PetscViewerType viewType; 1115552b385SBrandon PetscViewerFormat viewFormat; 1125552b385SBrandon PetscCall(PetscOptionsCreateViewer(PETSC_COMM_WORLD, NULL, NULL, "-dm_view_test", &viewer, &format, &flg)); 1135552b385SBrandon //PetscCall(PetscViewerPushFormat(viewer, format)); 1145552b385SBrandon if (flg) { 1155552b385SBrandon PetscCall(PetscPrintf(PETSC_COMM_SELF, " flg = TRUE \n")); 1165552b385SBrandon PetscCall(PetscViewerGetType(viewer, &viewType)); 1175552b385SBrandon PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_PETSC)); // PetscOptionsGetViewer returns &format as PETSC_VIEWER_DEFAULT need PETSC_VIEWER_HDF5_PETSC to save multiple DMPlexes in a single .h5 file. 1185552b385SBrandon PetscCall(PetscViewerGetFormat(viewer, &viewFormat)); 1195552b385SBrandon PetscCall(PetscPrintf(PETSC_COMM_SELF, " viewer type = %s \n", viewType)); 1205552b385SBrandon PetscCall(PetscPrintf(PETSC_COMM_SELF, " viewer format = %d \n", viewFormat)); 1215552b385SBrandon } 1225552b385SBrandon 1235552b385SBrandon // Refines Surface Mesh per option -dm_refine 1245552b385SBrandon PetscCall(DMSetFromOptions(dmNozzle)); 1255552b385SBrandon PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view")); 1265552b385SBrandon //PetscCall(DMGetOutputSequenceNumber(dmNozzle, &num, NULL)); 1275552b385SBrandon //PetscCall(DMSetOutputSequenceNumber(dmNozzle, num, -1)); 1285552b385SBrandon num += 1; 1295552b385SBrandon PetscCall(DMView(dmNozzle, viewer)); 1305552b385SBrandon PetscCall(surfArea(dmNozzle)); 1315552b385SBrandon 1325552b385SBrandon for (PetscInt saloop = 0, seqNum = 0; saloop < ctx.saMaxIter; ++saloop) { 1335552b385SBrandon PetscContainer modelObj; 1345552b385SBrandon //PetscContainer gradSACPObj, gradSAWObj; 1355552b385SBrandon PetscScalar *cpCoordData, *wData, *gradSACP, *gradSAW; 1365552b385SBrandon PetscInt cpCoordDataLength, wDataLength, maxNumEquiv; 1375552b385SBrandon PetscHMapI cpHashTable, wHashTable; 1385552b385SBrandon Mat cpEquiv; 1395552b385SBrandon char stpName[PETSC_MAX_PATH_LEN]; 1405552b385SBrandon char meshName[PETSC_MAX_PATH_LEN]; 1415552b385SBrandon 1425552b385SBrandon PetscCall(PetscSNPrintf(meshName, PETSC_MAX_PATH_LEN, "%s_Step_%" PetscInt_FMT, baseName, saloop)); 1435552b385SBrandon PetscCall(PetscObjectSetName((PetscObject)dmNozzle, meshName)); 1445552b385SBrandon // Save Step File of Updated Geometry at designated iterations 1455552b385SBrandon if (saloop == 1 || saloop == 2 || saloop == 5 || saloop == 20 || saloop == 50 || saloop == 100 || saloop == 150 || saloop == 200 || saloop == 300 || saloop == 400 || saloop == 500) { 1465552b385SBrandon PetscCall(PetscSNPrintf(stpName, sizeof(stpName), "newGeom_clean_%d.stp", saloop)); 1475552b385SBrandon } 1485552b385SBrandon 1495552b385SBrandon if (saloop == 0) PetscCall(DMSetFromOptions(dmNozzle)); 1505552b385SBrandon 1515552b385SBrandon // Expose Geometry Definition Data and Calculate Surface Gradients 1525552b385SBrandon PetscCall(DMPlexGeomDataAndGrads(dmNozzle, PETSC_FALSE)); 1535552b385SBrandon 1545552b385SBrandon PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "EGADS Model", (PetscObject *)&modelObj)); 1555552b385SBrandon if (!modelObj) PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "EGADSlite Model", (PetscObject *)&modelObj)); 1565552b385SBrandon 1575552b385SBrandon // Get attached EGADS model (pointer) 1585552b385SBrandon PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 1595552b385SBrandon 1605552b385SBrandon // Look to see if DM has Container for Geometry Control Point Data 1615552b385SBrandon //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Hash Table", (PetscObject *)&cpHashTableObj)); 1625552b385SBrandon //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Coordinates", (PetscObject *)&cpCoordDataObj)); 1635552b385SBrandon //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Coordinate Data Length", (PetscObject *)&cpCoordDataLengthObj)); 1645552b385SBrandon //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Weights Hash Table", (PetscObject *)&wHashTableObj)); 1655552b385SBrandon //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Weight Data", (PetscObject *)&wDataObj)); 1665552b385SBrandon //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Weight Data Length", (PetscObject *)&wDataLengthObj)); 1675552b385SBrandon ////PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Surface Area Control Point Gradient", (PetscObject *)&gradSACPObj)); 1685552b385SBrandon ////PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Surface Area Weights Gradient", (PetscObject *)&gradSAWObj)); 169*bfe80ac4SPierre Jolivet //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Equivalency Matrix", (PetscObject *)&cpEquivObj)); 1705552b385SBrandon //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Maximum Number Control Point Equivalency", (PetscObject *)&maxNumRelateObj)); 1715552b385SBrandon 1725552b385SBrandon // Get attached EGADS model Control Point and Weights Hash Tables and Data Arrays (pointer) 1735552b385SBrandon //PetscCall(PetscContainerGetPointer(cpHashTableObj, (void **)&cpHashTable)); 1745552b385SBrandon //PetscCall(PetscContainerGetPointer(cpCoordDataObj, (void **)&cpCoordData)); 1755552b385SBrandon //PetscCall(PetscContainerGetPointer(cpCoordDataLengthObj, (void **)&cpCoordDataLengthPtr)); 1765552b385SBrandon //PetscCall(PetscContainerGetPointer(wHashTableObj, (void **)&wHashTable)); 1775552b385SBrandon //PetscCall(PetscContainerGetPointer(wDataObj, (void **)&wData)); 1785552b385SBrandon //PetscCall(PetscContainerGetPointer(wDataLengthObj, (void **)&wDataLengthPtr)); 1795552b385SBrandon ////PetscCall(PetscContainerGetPointer(gradSACPObj, (void **)&gradSACP)); 1805552b385SBrandon ////PetscCall(PetscContainerGetPointer(gradSAWObj, (void **)&gradSAW)); 1815552b385SBrandon //PetscCall(PetscContainerGetPointer(cpEquivObj, (void **)&cpEquiv)); 1825552b385SBrandon //PetscCall(PetscContainerGetPointer(maxNumRelateObj, (void **)&maxNumRelatePtr)); 1835552b385SBrandon 1845552b385SBrandon // Trying out new Function 1855552b385SBrandon PetscInt cpArraySize, wArraySize; 1865552b385SBrandon PetscHMapI cpSurfGradHT; 1875552b385SBrandon Mat cpSurfGrad; 1885552b385SBrandon PetscScalar *gradVolCP, *gradVolW; 1895552b385SBrandon 1905552b385SBrandon PetscCall(DMPlexGetGeomCntrlPntAndWeightData(dmNozzle, &cpHashTable, &cpCoordDataLength, &cpCoordData, &maxNumEquiv, &cpEquiv, &wHashTable, &wDataLength, &wData)); 1915552b385SBrandon PetscCall(DMPlexGetGeomGradData(dmNozzle, &cpSurfGradHT, &cpSurfGrad, &cpArraySize, &gradSACP, &gradVolCP, &wArraySize, &gradSAW, &gradVolW)); 1925552b385SBrandon 1935552b385SBrandon // Get the number of bodies and body objects in the model 1945552b385SBrandon PetscCall(DMPlexGetGeomModelBodies(dmNozzle, &bodies, &Nb)); 1955552b385SBrandon 1965552b385SBrandon // Get all FACES of the Body 1975552b385SBrandon PetscGeom body = bodies[0]; 1985552b385SBrandon PetscCall(DMPlexGetGeomModelBodyFaces(dmNozzle, body, &fobjs, &Nf)); 1995552b385SBrandon 2005552b385SBrandon // Print out Surface Area and Volume using EGADS' Function 2015552b385SBrandon PetscScalar volume, surfaceArea; 2025552b385SBrandon PetscInt COGsize, IMCOGsize; 2035552b385SBrandon PetscScalar *centerOfGravity, *inertiaMatrixCOG; 2045552b385SBrandon PetscCall(DMPlexGetGeomBodyMassProperties(dmNozzle, body, &volume, &surfaceArea, ¢erOfGravity, &COGsize, &inertiaMatrixCOG, &IMCOGsize)); 2055552b385SBrandon 2065552b385SBrandon // Calculate Radius of Sphere for Equivalent Volume 2075552b385SBrandon if (saloop == 0) equivR = PetscPowReal(volume * (3. / 4.) / PETSC_PI, 1. / 3.); 2085552b385SBrandon 209*bfe80ac4SPierre Jolivet // Get Size of Control Point Equivalency Matrix 2105552b385SBrandon PetscInt numRows, numCols; 2115552b385SBrandon PetscCall(MatGetSize(cpEquiv, &numRows, &numCols)); 2125552b385SBrandon 2135552b385SBrandon // Get max number of relations 2145552b385SBrandon PetscInt maxRelateNew = 0; 2155552b385SBrandon for (PetscInt ii = 0; ii < numRows; ++ii) { 2165552b385SBrandon PetscInt maxRelateNewTemp = 0; 2175552b385SBrandon for (PetscInt jj = 0; jj < numCols; ++jj) { 2185552b385SBrandon PetscScalar matValue; 2195552b385SBrandon PetscCall(MatGetValue(cpEquiv, ii, jj, &matValue)); 2205552b385SBrandon 2215552b385SBrandon if (matValue > 0.0) maxRelateNewTemp += 1; 2225552b385SBrandon } 2235552b385SBrandon if (maxRelateNewTemp > maxRelateNew) maxRelateNew = maxRelateNewTemp; 2245552b385SBrandon } 2255552b385SBrandon 2265552b385SBrandon // Update Control Point and Weight definitions for each surface 2275552b385SBrandon for (PetscInt jj = 0; jj < Nf; ++jj) { 2285552b385SBrandon PetscGeom face = fobjs[jj]; 2295552b385SBrandon PetscInt numCntrlPnts = 0; 2305552b385SBrandon PetscCall(DMPlexGetGeomID(dmNozzle, body, face, &id)); 2315552b385SBrandon PetscCall(DMPlexGetGeomFaceNumOfControlPoints(dmNozzle, face, &numCntrlPnts)); 2325552b385SBrandon 2335552b385SBrandon // Update Control Points 2345552b385SBrandon PetscHashIter CPiter, Witer; 2355552b385SBrandon PetscBool CPfound, Wfound; 2365552b385SBrandon PetscInt faceCPStartRow, faceWStartRow; 2375552b385SBrandon 2385552b385SBrandon PetscCall(PetscHMapIFind(cpHashTable, id, &CPiter, &CPfound)); 2395552b385SBrandon PetscCheck(CPfound, comm, PETSC_ERR_SUP, "FACE ID not found in Control Point Hash Table"); 2405552b385SBrandon PetscCall(PetscHMapIGet(cpHashTable, id, &faceCPStartRow)); 2415552b385SBrandon 2425552b385SBrandon PetscCall(PetscHMapIFind(wHashTable, id, &Witer, &Wfound)); 2435552b385SBrandon PetscCheck(Wfound, comm, PETSC_ERR_SUP, "FACE ID not found in Control Point Weights Hash Table"); 2445552b385SBrandon PetscCall(PetscHMapIGet(wHashTable, id, &faceWStartRow)); 2455552b385SBrandon 2465552b385SBrandon for (PetscInt ii = 0; ii < numCntrlPnts; ++ii) { 2475552b385SBrandon PetscReal xold, yold, zold, rold, phi, theta; 2485552b385SBrandon 2495552b385SBrandon // Update Control Points - Original Code 2505552b385SBrandon xold = cpCoordData[faceCPStartRow + (3 * ii) + 0]; 2515552b385SBrandon yold = cpCoordData[faceCPStartRow + (3 * ii) + 1]; 2525552b385SBrandon zold = cpCoordData[faceCPStartRow + (3 * ii) + 2]; 2535552b385SBrandon rold = PetscSqrtReal(xold * xold + yold * yold + zold * zold); 2545552b385SBrandon phi = PetscAtan2Real(yold, xold); 2555552b385SBrandon theta = PetscAtan2Real(PetscSqrtReal(xold * xold + yold * yold), zold); 2565552b385SBrandon 2575552b385SBrandon // Account for Different Weights for Control Points on Degenerate Edges (multiple control points have same location and different weights 2585552b385SBrandon // only use the largest weight 2595552b385SBrandon PetscScalar maxWeight = 0.0; 2605552b385SBrandon //PetscInt wCntr = 0; 2615552b385SBrandon for (PetscInt kk = faceWStartRow; kk < faceWStartRow + numCntrlPnts; ++kk) { 2625552b385SBrandon PetscScalar matValue; 2635552b385SBrandon PetscCall(MatGetValue(cpEquiv, faceWStartRow + ii, kk, &matValue)); 2645552b385SBrandon 2655552b385SBrandon PetscScalar weight = 0.0; 2665552b385SBrandon if (matValue > 0.0) { 2675552b385SBrandon weight = wData[kk]; 2685552b385SBrandon 2695552b385SBrandon if (weight > maxWeight) { maxWeight = weight; } 2705552b385SBrandon //wCntr += 1; 2715552b385SBrandon } 2725552b385SBrandon } 2735552b385SBrandon 2745552b385SBrandon //Reduce to Constant R = 0.0254 2755552b385SBrandon PetscScalar deltaR; 2765552b385SBrandon PetscScalar localTargetR; 2775552b385SBrandon localTargetR = equivR / maxWeight; 2785552b385SBrandon deltaR = rold - localTargetR; 2795552b385SBrandon 2805552b385SBrandon cpCoordData[faceCPStartRow + (3 * ii) + 0] += -0.05 * deltaR * PetscCosReal(phi) * PetscSinReal(theta); 2815552b385SBrandon cpCoordData[faceCPStartRow + (3 * ii) + 1] += -0.05 * deltaR * PetscSinReal(phi) * PetscSinReal(theta); 2825552b385SBrandon cpCoordData[faceCPStartRow + (3 * ii) + 2] += -0.05 * deltaR * PetscCosReal(theta); 2835552b385SBrandon } 2845552b385SBrandon } 2855552b385SBrandon PetscCall(DMPlexFreeGeomObject(dmNozzle, fobjs)); 2865552b385SBrandon PetscBool writeFile = PETSC_FALSE; 2875552b385SBrandon 2885552b385SBrandon // Modify Control Points of Geometry 2895552b385SBrandon PetscCall(DMPlexModifyGeomModel(dmNozzle, comm, cpCoordData, wData, PETSC_FALSE, writeFile, stpName)); 2905552b385SBrandon writeFile = PETSC_FALSE; 2915552b385SBrandon 2925552b385SBrandon // Get attached EGADS model (pointer) 2935552b385SBrandon PetscGeom newmodel; 2945552b385SBrandon PetscCall(PetscContainerGetPointer(modelObj, (void **)&newmodel)); 2955552b385SBrandon 2965552b385SBrandon // Get the number of bodies and body objects in the model 2975552b385SBrandon PetscCall(DMPlexGetGeomModelBodies(dmNozzle, &bodies, &Nb)); 2985552b385SBrandon 2995552b385SBrandon // Get all FACES of the Body 3005552b385SBrandon PetscGeom newbody = bodies[0]; 3015552b385SBrandon PetscCall(DMPlexGetGeomModelBodyFaces(dmNozzle, newbody, &fobjs, &Nf)); 3025552b385SBrandon 3035552b385SBrandon // Save Step File of Updated Geometry at designated iterations 3045552b385SBrandon if (saloop == 1 || saloop == 2 || saloop == 5 || saloop == 20 || saloop == 50 || saloop == 100 || saloop == 150 || saloop == 200 || saloop == 300 || saloop == 400 || saloop == 500) { writeFile = PETSC_TRUE; } 3055552b385SBrandon 3065552b385SBrandon // Modify Geometry and Inflate Mesh to New Geoemetry 3075552b385SBrandon PetscCall(DMPlexModifyGeomModel(dmNozzle, comm, cpCoordData, wData, PETSC_FALSE, writeFile, stpName)); 3085552b385SBrandon PetscCall(DMPlexInflateToGeomModel(dmNozzle, PETSC_TRUE)); 3095552b385SBrandon 3105552b385SBrandon // Periodically Refine and Write Mesh to hdf5 file 3115552b385SBrandon if (saloop == 0) { 3125552b385SBrandon PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view7")); 3135552b385SBrandon //PetscCall(DMGetOutputSequenceNumber(dmNozzle, &num, NULL)); 3145552b385SBrandon PetscCall(PetscPrintf(PETSC_COMM_SELF, " num = %d \n", num)); 3155552b385SBrandon //PetscCall(DMGetOutputSequenceNumber(dmNozzle, NULL, &num)); 3165552b385SBrandon //PetscCall(DMSetOutputSequenceNumber(dmNozzle, num, saloop)); 3175552b385SBrandon num += 1; 3185552b385SBrandon PetscCall(PetscObjectSetName((PetscObject)dmNozzle, "nozzle_meshes_1")); 3195552b385SBrandon //PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_APPEND)); 3205552b385SBrandon PetscCall(DMView(dmNozzle, viewer)); 3215552b385SBrandon } 3225552b385SBrandon if (saloop == 1 || saloop == 5 || saloop == 20 || saloop == 50 || saloop == 100 || saloop == 150 || saloop == 200 || saloop == 300 || saloop == 400 || saloop == 500) { 3235552b385SBrandon PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view")); 3245552b385SBrandon PetscCall(DMSetOutputSequenceNumber(dmNozzle, seqNum++, saloop)); 3255552b385SBrandon PetscCall(DMView(dmNozzle, viewer)); 3265552b385SBrandon if (saloop == 200 || saloop == 500) { 3275552b385SBrandon PetscCall(DMSetFromOptions(dmNozzle)); 3285552b385SBrandon PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view")); 3295552b385SBrandon PetscCall(DMSetOutputSequenceNumber(dmNozzle, seqNum++, saloop)); 3305552b385SBrandon PetscCall(DMView(dmNozzle, viewer)); 3315552b385SBrandon 3325552b385SBrandon PetscCall(DMSetFromOptions(dmNozzle)); 3335552b385SBrandon PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view")); 3345552b385SBrandon PetscCall(DMSetOutputSequenceNumber(dmNozzle, seqNum++, saloop)); 3355552b385SBrandon PetscCall(DMView(dmNozzle, viewer)); 3365552b385SBrandon } 3375552b385SBrandon } 3385552b385SBrandon PetscCall(DMPlexRestoreGeomBodyMassProperties(dmNozzle, body, &volume, &surfaceArea, ¢erOfGravity, &COGsize, &inertiaMatrixCOG, &IMCOGsize)); 3395552b385SBrandon PetscCall(DMPlexRestoreGeomCntrlPntAndWeightData(dmNozzle, &cpHashTable, &cpCoordDataLength, &cpCoordData, &maxNumEquiv, &cpEquiv, &wHashTable, &wDataLength, &wData)); 3405552b385SBrandon PetscCall(DMPlexRestoreGeomGradData(dmNozzle, &cpSurfGradHT, &cpSurfGrad, &cpArraySize, &gradSACP, &gradVolCP, &wArraySize, &gradSAW, &gradVolW)); 3415552b385SBrandon PetscCall(DMPlexFreeGeomObject(dmNozzle, fobjs)); 3425552b385SBrandon } 3435552b385SBrandon PetscCall(DMDestroy(&dmNozzle)); 3445552b385SBrandon PetscCall(PetscFinalize()); 3455552b385SBrandon } 3465552b385SBrandon 3475552b385SBrandon /*TEST 3485552b385SBrandon build: 3495552b385SBrandon requires: egads 3505552b385SBrandon 3515552b385SBrandon # To view mesh use -dm_plex_view_hdf5_storage_version 3.1.0 -dm_view hdf5:mesh_minSA_abstract.h5:hdf5_petsc:append 3525552b385SBrandon test: 3535552b385SBrandon suffix: minSA 3545552b385SBrandon requires: datafilespath 3555552b385SBrandon temporaries: newGeom_clean_1.stp newGeom_clean_2.stp newGeom_clean_5.stp 3565552b385SBrandon args: -filename ${DATAFILESPATH}/meshes/cad/sphere_example.stp \ 3575552b385SBrandon -dm_refine 1 -dm_plex_geom_print_model 1 -dm_plex_geom_shape_opt 1 \ 3585552b385SBrandon -sa_max_iter 2 3595552b385SBrandon 3605552b385SBrandon TEST*/ 361