xref: /petsc/src/dm/impls/plex/tutorials/ex19.c (revision 5552b385b77b214b234683fbe1b434751e6350f0)
1*5552b385SBrandon static const char help[] = "Test of PETSc CAD Shape Optimization & Mesh Modification Technology";
2*5552b385SBrandon 
3*5552b385SBrandon #include <petscdmplexegads.h>
4*5552b385SBrandon #include <petsc/private/hashmapi.h>
5*5552b385SBrandon 
6*5552b385SBrandon typedef struct {
7*5552b385SBrandon   char     filename[PETSC_MAX_PATH_LEN];
8*5552b385SBrandon   PetscInt saMaxIter; // Maximum number of iterations of shape optimization loop
9*5552b385SBrandon } AppCtx;
10*5552b385SBrandon 
11*5552b385SBrandon static PetscErrorCode surfArea(DM dm)
12*5552b385SBrandon {
13*5552b385SBrandon   DMLabel     bodyLabel, faceLabel;
14*5552b385SBrandon   double      surfaceArea = 0., volume = 0.;
15*5552b385SBrandon   PetscReal   vol, centroid[3], normal[3];
16*5552b385SBrandon   PetscInt    dim, cStart, cEnd, fStart, fEnd;
17*5552b385SBrandon   PetscInt    bodyID, faceID;
18*5552b385SBrandon   MPI_Comm    comm;
19*5552b385SBrandon   const char *name;
20*5552b385SBrandon 
21*5552b385SBrandon   PetscFunctionBeginUser;
22*5552b385SBrandon   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
23*5552b385SBrandon   PetscCall(PetscObjectGetName((PetscObject)dm, &name));
24*5552b385SBrandon   PetscCall(DMGetDimension(dm, &dim));
25*5552b385SBrandon   PetscCall(PetscPrintf(comm, "    dim = %" PetscInt_FMT "\n", dim));
26*5552b385SBrandon   PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
27*5552b385SBrandon   PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
28*5552b385SBrandon 
29*5552b385SBrandon   if (dim == 2) {
30*5552b385SBrandon     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
31*5552b385SBrandon     for (PetscInt ii = cStart; ii < cEnd; ++ii) {
32*5552b385SBrandon       PetscCall(DMLabelGetValue(faceLabel, ii, &faceID));
33*5552b385SBrandon       if (faceID >= 0) {
34*5552b385SBrandon         PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal));
35*5552b385SBrandon         surfaceArea += vol;
36*5552b385SBrandon       }
37*5552b385SBrandon     }
38*5552b385SBrandon   }
39*5552b385SBrandon 
40*5552b385SBrandon   if (dim == 3) {
41*5552b385SBrandon     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
42*5552b385SBrandon     for (PetscInt ii = fStart; ii < fEnd; ++ii) {
43*5552b385SBrandon       PetscCall(DMLabelGetValue(faceLabel, ii, &faceID));
44*5552b385SBrandon       if (faceID >= 0) {
45*5552b385SBrandon         PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal));
46*5552b385SBrandon         surfaceArea += vol;
47*5552b385SBrandon       }
48*5552b385SBrandon     }
49*5552b385SBrandon 
50*5552b385SBrandon     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
51*5552b385SBrandon     for (PetscInt ii = cStart; ii < cEnd; ++ii) {
52*5552b385SBrandon       PetscCall(DMLabelGetValue(bodyLabel, ii, &bodyID));
53*5552b385SBrandon       if (bodyID >= 0) {
54*5552b385SBrandon         PetscCall(DMPlexComputeCellGeometryFVM(dm, ii, &vol, centroid, normal));
55*5552b385SBrandon         volume += vol;
56*5552b385SBrandon       }
57*5552b385SBrandon     }
58*5552b385SBrandon   }
59*5552b385SBrandon 
60*5552b385SBrandon   if (dim == 2) {
61*5552b385SBrandon     PetscCall(PetscPrintf(comm, "%s Surface Area = %.6e \n\n", name, (double)surfaceArea));
62*5552b385SBrandon   } else if (dim == 3) {
63*5552b385SBrandon     PetscCall(PetscPrintf(comm, "%s Volume = %.6e \n", name, (double)volume));
64*5552b385SBrandon     PetscCall(PetscPrintf(comm, "%s Surface Area = %.6e \n\n", name, (double)surfaceArea));
65*5552b385SBrandon   }
66*5552b385SBrandon   PetscFunctionReturn(PETSC_SUCCESS);
67*5552b385SBrandon }
68*5552b385SBrandon 
69*5552b385SBrandon static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
70*5552b385SBrandon {
71*5552b385SBrandon   PetscFunctionBeginUser;
72*5552b385SBrandon   options->filename[0] = '\0';
73*5552b385SBrandon   options->saMaxIter   = 200;
74*5552b385SBrandon 
75*5552b385SBrandon   PetscOptionsBegin(comm, "", "DMPlex w/CAD Options", "DMPlex w/CAD");
76*5552b385SBrandon   PetscCall(PetscOptionsString("-filename", "The CAD/Geometry file", __FILE__, options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL));
77*5552b385SBrandon   PetscCall(PetscOptionsBoundedInt("-sa_max_iter", "The maximum number of iterates for the shape optimization loop", __FILE__, options->saMaxIter, &options->saMaxIter, NULL, 0));
78*5552b385SBrandon   PetscOptionsEnd();
79*5552b385SBrandon   PetscFunctionReturn(0);
80*5552b385SBrandon }
81*5552b385SBrandon 
82*5552b385SBrandon int main(int argc, char *argv[])
83*5552b385SBrandon {
84*5552b385SBrandon   /* EGADS/EGADSlite variables */
85*5552b385SBrandon   PetscGeom model, *bodies, *fobjs;
86*5552b385SBrandon   int       Nb, Nf, id;
87*5552b385SBrandon   /* PETSc variables */
88*5552b385SBrandon   DM          dmNozzle = NULL;
89*5552b385SBrandon   MPI_Comm    comm;
90*5552b385SBrandon   AppCtx      ctx;
91*5552b385SBrandon   PetscScalar equivR     = 0.0;
92*5552b385SBrandon   const char  baseName[] = "Nozzle_Mesh";
93*5552b385SBrandon 
94*5552b385SBrandon   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
95*5552b385SBrandon   comm = PETSC_COMM_WORLD;
96*5552b385SBrandon   PetscCall(ProcessOptions(comm, &ctx));
97*5552b385SBrandon 
98*5552b385SBrandon   PetscCall(DMPlexCreateFromFile(comm, ctx.filename, "EGADS", PETSC_TRUE, &dmNozzle));
99*5552b385SBrandon   PetscCall(PetscObjectSetName((PetscObject)dmNozzle, baseName));
100*5552b385SBrandon   //PetscCall(DMCreate(PETSC_COMM_WORLD, &dmNozzle));
101*5552b385SBrandon   //PetscCall(DMSetType(dmNozzle, DMPLEX));
102*5552b385SBrandon   //PetscCall(DMPlexDistributeSetDefault(dmNozzle, PETSC_FALSE));
103*5552b385SBrandon   //PetscCall(DMSetFromOptions(dmNozzle));
104*5552b385SBrandon 
105*5552b385SBrandon   // Get Common Viewer to store all Mesh Results
106*5552b385SBrandon   PetscViewer       viewer;
107*5552b385SBrandon   PetscViewerFormat format;
108*5552b385SBrandon   PetscBool         flg;
109*5552b385SBrandon   PetscInt          num = 0;
110*5552b385SBrandon   PetscViewerType   viewType;
111*5552b385SBrandon   PetscViewerFormat viewFormat;
112*5552b385SBrandon   PetscCall(PetscOptionsCreateViewer(PETSC_COMM_WORLD, NULL, NULL, "-dm_view_test", &viewer, &format, &flg));
113*5552b385SBrandon   //PetscCall(PetscViewerPushFormat(viewer, format));
114*5552b385SBrandon   if (flg) {
115*5552b385SBrandon     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  flg = TRUE \n"));
116*5552b385SBrandon     PetscCall(PetscViewerGetType(viewer, &viewType));
117*5552b385SBrandon     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.
118*5552b385SBrandon     PetscCall(PetscViewerGetFormat(viewer, &viewFormat));
119*5552b385SBrandon     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  viewer type = %s \n", viewType));
120*5552b385SBrandon     PetscCall(PetscPrintf(PETSC_COMM_SELF, "  viewer format = %d \n", viewFormat));
121*5552b385SBrandon   }
122*5552b385SBrandon 
123*5552b385SBrandon   // Refines Surface Mesh per option -dm_refine
124*5552b385SBrandon   PetscCall(DMSetFromOptions(dmNozzle));
125*5552b385SBrandon   PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view"));
126*5552b385SBrandon   //PetscCall(DMGetOutputSequenceNumber(dmNozzle, &num, NULL));
127*5552b385SBrandon   //PetscCall(DMSetOutputSequenceNumber(dmNozzle, num, -1));
128*5552b385SBrandon   num += 1;
129*5552b385SBrandon   PetscCall(DMView(dmNozzle, viewer));
130*5552b385SBrandon   PetscCall(surfArea(dmNozzle));
131*5552b385SBrandon 
132*5552b385SBrandon   for (PetscInt saloop = 0, seqNum = 0; saloop < ctx.saMaxIter; ++saloop) {
133*5552b385SBrandon     PetscContainer modelObj;
134*5552b385SBrandon     //PetscContainer gradSACPObj, gradSAWObj;
135*5552b385SBrandon     PetscScalar *cpCoordData, *wData, *gradSACP, *gradSAW;
136*5552b385SBrandon     PetscInt     cpCoordDataLength, wDataLength, maxNumEquiv;
137*5552b385SBrandon     PetscHMapI   cpHashTable, wHashTable;
138*5552b385SBrandon     Mat          cpEquiv;
139*5552b385SBrandon     char         stpName[PETSC_MAX_PATH_LEN];
140*5552b385SBrandon     char         meshName[PETSC_MAX_PATH_LEN];
141*5552b385SBrandon 
142*5552b385SBrandon     PetscCall(PetscSNPrintf(meshName, PETSC_MAX_PATH_LEN, "%s_Step_%" PetscInt_FMT, baseName, saloop));
143*5552b385SBrandon     PetscCall(PetscObjectSetName((PetscObject)dmNozzle, meshName));
144*5552b385SBrandon     // Save Step File of Updated Geometry at designated iterations
145*5552b385SBrandon     if (saloop == 1 || saloop == 2 || saloop == 5 || saloop == 20 || saloop == 50 || saloop == 100 || saloop == 150 || saloop == 200 || saloop == 300 || saloop == 400 || saloop == 500) {
146*5552b385SBrandon       PetscCall(PetscSNPrintf(stpName, sizeof(stpName), "newGeom_clean_%d.stp", saloop));
147*5552b385SBrandon     }
148*5552b385SBrandon 
149*5552b385SBrandon     if (saloop == 0) PetscCall(DMSetFromOptions(dmNozzle));
150*5552b385SBrandon 
151*5552b385SBrandon     // Expose Geometry Definition Data and Calculate Surface Gradients
152*5552b385SBrandon     PetscCall(DMPlexGeomDataAndGrads(dmNozzle, PETSC_FALSE));
153*5552b385SBrandon 
154*5552b385SBrandon     PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "EGADS Model", (PetscObject *)&modelObj));
155*5552b385SBrandon     if (!modelObj) PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "EGADSlite Model", (PetscObject *)&modelObj));
156*5552b385SBrandon 
157*5552b385SBrandon     // Get attached EGADS model (pointer)
158*5552b385SBrandon     PetscCall(PetscContainerGetPointer(modelObj, (void **)&model));
159*5552b385SBrandon 
160*5552b385SBrandon     // Look to see if DM has Container for Geometry Control Point Data
161*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Hash Table", (PetscObject *)&cpHashTableObj));
162*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Coordinates", (PetscObject *)&cpCoordDataObj));
163*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Coordinate Data Length", (PetscObject *)&cpCoordDataLengthObj));
164*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Weights Hash Table", (PetscObject *)&wHashTableObj));
165*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Weight Data", (PetscObject *)&wDataObj));
166*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Weight Data Length", (PetscObject *)&wDataLengthObj));
167*5552b385SBrandon     ////PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Surface Area Control Point Gradient", (PetscObject *)&gradSACPObj));
168*5552b385SBrandon     ////PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Surface Area Weights Gradient", (PetscObject *)&gradSAWObj));
169*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Control Point Equivalancy Matrix", (PetscObject *)&cpEquivObj));
170*5552b385SBrandon     //PetscCall(PetscObjectQuery((PetscObject)dmNozzle, "Maximum Number Control Point Equivalency", (PetscObject *)&maxNumRelateObj));
171*5552b385SBrandon 
172*5552b385SBrandon     // Get attached EGADS model Control Point and Weights Hash Tables and Data Arrays (pointer)
173*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(cpHashTableObj, (void **)&cpHashTable));
174*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(cpCoordDataObj, (void **)&cpCoordData));
175*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(cpCoordDataLengthObj, (void **)&cpCoordDataLengthPtr));
176*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(wHashTableObj, (void **)&wHashTable));
177*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(wDataObj, (void **)&wData));
178*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(wDataLengthObj, (void **)&wDataLengthPtr));
179*5552b385SBrandon     ////PetscCall(PetscContainerGetPointer(gradSACPObj, (void **)&gradSACP));
180*5552b385SBrandon     ////PetscCall(PetscContainerGetPointer(gradSAWObj, (void **)&gradSAW));
181*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(cpEquivObj, (void **)&cpEquiv));
182*5552b385SBrandon     //PetscCall(PetscContainerGetPointer(maxNumRelateObj, (void **)&maxNumRelatePtr));
183*5552b385SBrandon 
184*5552b385SBrandon     // Trying out new Function
185*5552b385SBrandon     PetscInt     cpArraySize, wArraySize;
186*5552b385SBrandon     PetscHMapI   cpSurfGradHT;
187*5552b385SBrandon     Mat          cpSurfGrad;
188*5552b385SBrandon     PetscScalar *gradVolCP, *gradVolW;
189*5552b385SBrandon 
190*5552b385SBrandon     PetscCall(DMPlexGetGeomCntrlPntAndWeightData(dmNozzle, &cpHashTable, &cpCoordDataLength, &cpCoordData, &maxNumEquiv, &cpEquiv, &wHashTable, &wDataLength, &wData));
191*5552b385SBrandon     PetscCall(DMPlexGetGeomGradData(dmNozzle, &cpSurfGradHT, &cpSurfGrad, &cpArraySize, &gradSACP, &gradVolCP, &wArraySize, &gradSAW, &gradVolW));
192*5552b385SBrandon 
193*5552b385SBrandon     // Get the number of bodies and body objects in the model
194*5552b385SBrandon     PetscCall(DMPlexGetGeomModelBodies(dmNozzle, &bodies, &Nb));
195*5552b385SBrandon 
196*5552b385SBrandon     // Get all FACES of the Body
197*5552b385SBrandon     PetscGeom body = bodies[0];
198*5552b385SBrandon     PetscCall(DMPlexGetGeomModelBodyFaces(dmNozzle, body, &fobjs, &Nf));
199*5552b385SBrandon 
200*5552b385SBrandon     // Print out Surface Area and Volume using EGADS' Function
201*5552b385SBrandon     PetscScalar  volume, surfaceArea;
202*5552b385SBrandon     PetscInt     COGsize, IMCOGsize;
203*5552b385SBrandon     PetscScalar *centerOfGravity, *inertiaMatrixCOG;
204*5552b385SBrandon     PetscCall(DMPlexGetGeomBodyMassProperties(dmNozzle, body, &volume, &surfaceArea, &centerOfGravity, &COGsize, &inertiaMatrixCOG, &IMCOGsize));
205*5552b385SBrandon 
206*5552b385SBrandon     // Calculate Radius of Sphere for Equivalent Volume
207*5552b385SBrandon     if (saloop == 0) equivR = PetscPowReal(volume * (3. / 4.) / PETSC_PI, 1. / 3.);
208*5552b385SBrandon 
209*5552b385SBrandon     // Get Size of Control Point Equivalancy Matrix
210*5552b385SBrandon     PetscInt numRows, numCols;
211*5552b385SBrandon     PetscCall(MatGetSize(cpEquiv, &numRows, &numCols));
212*5552b385SBrandon 
213*5552b385SBrandon     // Get max number of relations
214*5552b385SBrandon     PetscInt maxRelateNew = 0;
215*5552b385SBrandon     for (PetscInt ii = 0; ii < numRows; ++ii) {
216*5552b385SBrandon       PetscInt maxRelateNewTemp = 0;
217*5552b385SBrandon       for (PetscInt jj = 0; jj < numCols; ++jj) {
218*5552b385SBrandon         PetscScalar matValue;
219*5552b385SBrandon         PetscCall(MatGetValue(cpEquiv, ii, jj, &matValue));
220*5552b385SBrandon 
221*5552b385SBrandon         if (matValue > 0.0) maxRelateNewTemp += 1;
222*5552b385SBrandon       }
223*5552b385SBrandon       if (maxRelateNewTemp > maxRelateNew) maxRelateNew = maxRelateNewTemp;
224*5552b385SBrandon     }
225*5552b385SBrandon 
226*5552b385SBrandon     // Update Control Point and Weight definitions for each surface
227*5552b385SBrandon     for (PetscInt jj = 0; jj < Nf; ++jj) {
228*5552b385SBrandon       PetscGeom face         = fobjs[jj];
229*5552b385SBrandon       PetscInt  numCntrlPnts = 0;
230*5552b385SBrandon       PetscCall(DMPlexGetGeomID(dmNozzle, body, face, &id));
231*5552b385SBrandon       PetscCall(DMPlexGetGeomFaceNumOfControlPoints(dmNozzle, face, &numCntrlPnts));
232*5552b385SBrandon 
233*5552b385SBrandon       // Update Control Points
234*5552b385SBrandon       PetscHashIter CPiter, Witer;
235*5552b385SBrandon       PetscBool     CPfound, Wfound;
236*5552b385SBrandon       PetscInt      faceCPStartRow, faceWStartRow;
237*5552b385SBrandon 
238*5552b385SBrandon       PetscCall(PetscHMapIFind(cpHashTable, id, &CPiter, &CPfound));
239*5552b385SBrandon       PetscCheck(CPfound, comm, PETSC_ERR_SUP, "FACE ID not found in Control Point Hash Table");
240*5552b385SBrandon       PetscCall(PetscHMapIGet(cpHashTable, id, &faceCPStartRow));
241*5552b385SBrandon 
242*5552b385SBrandon       PetscCall(PetscHMapIFind(wHashTable, id, &Witer, &Wfound));
243*5552b385SBrandon       PetscCheck(Wfound, comm, PETSC_ERR_SUP, "FACE ID not found in Control Point Weights Hash Table");
244*5552b385SBrandon       PetscCall(PetscHMapIGet(wHashTable, id, &faceWStartRow));
245*5552b385SBrandon 
246*5552b385SBrandon       for (PetscInt ii = 0; ii < numCntrlPnts; ++ii) {
247*5552b385SBrandon         PetscReal xold, yold, zold, rold, phi, theta;
248*5552b385SBrandon 
249*5552b385SBrandon         // Update Control Points - Original Code
250*5552b385SBrandon         xold  = cpCoordData[faceCPStartRow + (3 * ii) + 0];
251*5552b385SBrandon         yold  = cpCoordData[faceCPStartRow + (3 * ii) + 1];
252*5552b385SBrandon         zold  = cpCoordData[faceCPStartRow + (3 * ii) + 2];
253*5552b385SBrandon         rold  = PetscSqrtReal(xold * xold + yold * yold + zold * zold);
254*5552b385SBrandon         phi   = PetscAtan2Real(yold, xold);
255*5552b385SBrandon         theta = PetscAtan2Real(PetscSqrtReal(xold * xold + yold * yold), zold);
256*5552b385SBrandon 
257*5552b385SBrandon         // Account for Different Weights for Control Points on Degenerate Edges (multiple control points have same location and different weights
258*5552b385SBrandon         // only use the largest weight
259*5552b385SBrandon         PetscScalar maxWeight = 0.0;
260*5552b385SBrandon         //PetscInt    wCntr = 0;
261*5552b385SBrandon         for (PetscInt kk = faceWStartRow; kk < faceWStartRow + numCntrlPnts; ++kk) {
262*5552b385SBrandon           PetscScalar matValue;
263*5552b385SBrandon           PetscCall(MatGetValue(cpEquiv, faceWStartRow + ii, kk, &matValue));
264*5552b385SBrandon 
265*5552b385SBrandon           PetscScalar weight = 0.0;
266*5552b385SBrandon           if (matValue > 0.0) {
267*5552b385SBrandon             weight = wData[kk];
268*5552b385SBrandon 
269*5552b385SBrandon             if (weight > maxWeight) { maxWeight = weight; }
270*5552b385SBrandon             //wCntr += 1;
271*5552b385SBrandon           }
272*5552b385SBrandon         }
273*5552b385SBrandon 
274*5552b385SBrandon         //Reduce to Constant R = 0.0254
275*5552b385SBrandon         PetscScalar deltaR;
276*5552b385SBrandon         PetscScalar localTargetR;
277*5552b385SBrandon         localTargetR = equivR / maxWeight;
278*5552b385SBrandon         deltaR       = rold - localTargetR;
279*5552b385SBrandon 
280*5552b385SBrandon         cpCoordData[faceCPStartRow + (3 * ii) + 0] += -0.05 * deltaR * PetscCosReal(phi) * PetscSinReal(theta);
281*5552b385SBrandon         cpCoordData[faceCPStartRow + (3 * ii) + 1] += -0.05 * deltaR * PetscSinReal(phi) * PetscSinReal(theta);
282*5552b385SBrandon         cpCoordData[faceCPStartRow + (3 * ii) + 2] += -0.05 * deltaR * PetscCosReal(theta);
283*5552b385SBrandon       }
284*5552b385SBrandon     }
285*5552b385SBrandon     PetscCall(DMPlexFreeGeomObject(dmNozzle, fobjs));
286*5552b385SBrandon     PetscBool writeFile = PETSC_FALSE;
287*5552b385SBrandon 
288*5552b385SBrandon     // Modify Control Points of Geometry
289*5552b385SBrandon     PetscCall(DMPlexModifyGeomModel(dmNozzle, comm, cpCoordData, wData, PETSC_FALSE, writeFile, stpName));
290*5552b385SBrandon     writeFile = PETSC_FALSE;
291*5552b385SBrandon 
292*5552b385SBrandon     // Get attached EGADS model (pointer)
293*5552b385SBrandon     PetscGeom newmodel;
294*5552b385SBrandon     PetscCall(PetscContainerGetPointer(modelObj, (void **)&newmodel));
295*5552b385SBrandon 
296*5552b385SBrandon     // Get the number of bodies and body objects in the model
297*5552b385SBrandon     PetscCall(DMPlexGetGeomModelBodies(dmNozzle, &bodies, &Nb));
298*5552b385SBrandon 
299*5552b385SBrandon     // Get all FACES of the Body
300*5552b385SBrandon     PetscGeom newbody = bodies[0];
301*5552b385SBrandon     PetscCall(DMPlexGetGeomModelBodyFaces(dmNozzle, newbody, &fobjs, &Nf));
302*5552b385SBrandon 
303*5552b385SBrandon     // Save Step File of Updated Geometry at designated iterations
304*5552b385SBrandon     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; }
305*5552b385SBrandon 
306*5552b385SBrandon     // Modify Geometry and Inflate Mesh to New Geoemetry
307*5552b385SBrandon     PetscCall(DMPlexModifyGeomModel(dmNozzle, comm, cpCoordData, wData, PETSC_FALSE, writeFile, stpName));
308*5552b385SBrandon     PetscCall(DMPlexInflateToGeomModel(dmNozzle, PETSC_TRUE));
309*5552b385SBrandon 
310*5552b385SBrandon     // Periodically Refine and Write Mesh to hdf5 file
311*5552b385SBrandon     if (saloop == 0) {
312*5552b385SBrandon       PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view7"));
313*5552b385SBrandon       //PetscCall(DMGetOutputSequenceNumber(dmNozzle, &num, NULL));
314*5552b385SBrandon       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  num = %d \n", num));
315*5552b385SBrandon       //PetscCall(DMGetOutputSequenceNumber(dmNozzle, NULL, &num));
316*5552b385SBrandon       //PetscCall(DMSetOutputSequenceNumber(dmNozzle, num, saloop));
317*5552b385SBrandon       num += 1;
318*5552b385SBrandon       PetscCall(PetscObjectSetName((PetscObject)dmNozzle, "nozzle_meshes_1"));
319*5552b385SBrandon       //PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_APPEND));
320*5552b385SBrandon       PetscCall(DMView(dmNozzle, viewer));
321*5552b385SBrandon     }
322*5552b385SBrandon     if (saloop == 1 || saloop == 5 || saloop == 20 || saloop == 50 || saloop == 100 || saloop == 150 || saloop == 200 || saloop == 300 || saloop == 400 || saloop == 500) {
323*5552b385SBrandon       PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view"));
324*5552b385SBrandon       PetscCall(DMSetOutputSequenceNumber(dmNozzle, seqNum++, saloop));
325*5552b385SBrandon       PetscCall(DMView(dmNozzle, viewer));
326*5552b385SBrandon       if (saloop == 200 || saloop == 500) {
327*5552b385SBrandon         PetscCall(DMSetFromOptions(dmNozzle));
328*5552b385SBrandon         PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view"));
329*5552b385SBrandon         PetscCall(DMSetOutputSequenceNumber(dmNozzle, seqNum++, saloop));
330*5552b385SBrandon         PetscCall(DMView(dmNozzle, viewer));
331*5552b385SBrandon 
332*5552b385SBrandon         PetscCall(DMSetFromOptions(dmNozzle));
333*5552b385SBrandon         PetscCall(DMViewFromOptions(dmNozzle, NULL, "-dm_view"));
334*5552b385SBrandon         PetscCall(DMSetOutputSequenceNumber(dmNozzle, seqNum++, saloop));
335*5552b385SBrandon         PetscCall(DMView(dmNozzle, viewer));
336*5552b385SBrandon       }
337*5552b385SBrandon     }
338*5552b385SBrandon     PetscCall(DMPlexRestoreGeomBodyMassProperties(dmNozzle, body, &volume, &surfaceArea, &centerOfGravity, &COGsize, &inertiaMatrixCOG, &IMCOGsize));
339*5552b385SBrandon     PetscCall(DMPlexRestoreGeomCntrlPntAndWeightData(dmNozzle, &cpHashTable, &cpCoordDataLength, &cpCoordData, &maxNumEquiv, &cpEquiv, &wHashTable, &wDataLength, &wData));
340*5552b385SBrandon     PetscCall(DMPlexRestoreGeomGradData(dmNozzle, &cpSurfGradHT, &cpSurfGrad, &cpArraySize, &gradSACP, &gradVolCP, &wArraySize, &gradSAW, &gradVolW));
341*5552b385SBrandon     PetscCall(DMPlexFreeGeomObject(dmNozzle, fobjs));
342*5552b385SBrandon   }
343*5552b385SBrandon   PetscCall(DMDestroy(&dmNozzle));
344*5552b385SBrandon   PetscCall(PetscFinalize());
345*5552b385SBrandon }
346*5552b385SBrandon 
347*5552b385SBrandon /*TEST
348*5552b385SBrandon   build:
349*5552b385SBrandon     requires: egads
350*5552b385SBrandon 
351*5552b385SBrandon   # To view mesh use -dm_plex_view_hdf5_storage_version 3.1.0 -dm_view hdf5:mesh_minSA_abstract.h5:hdf5_petsc:append
352*5552b385SBrandon   test:
353*5552b385SBrandon     suffix: minSA
354*5552b385SBrandon     requires: datafilespath
355*5552b385SBrandon     temporaries: newGeom_clean_1.stp newGeom_clean_2.stp newGeom_clean_5.stp
356*5552b385SBrandon     args: -filename ${DATAFILESPATH}/meshes/cad/sphere_example.stp \
357*5552b385SBrandon           -dm_refine 1 -dm_plex_geom_print_model 1 -dm_plex_geom_shape_opt 1 \
358*5552b385SBrandon           -sa_max_iter 2
359*5552b385SBrandon 
360*5552b385SBrandon TEST*/
361