Actual source code: plexegads.c
1: #include "petscsys.h"
2: #include <petsc/private/dmpleximpl.h>
3: #include <petsc/private/hashmapi.h>
5: /* We need to understand how to natively parse STEP files. There seems to be only one open-source implementation of
6: the STEP parser contained in the OpenCASCADE package. It is enough to make a strong man weep:
8: https://github.com/tpaviot/oce/tree/master/src/STEPControl
10: The STEP, and inner EXPRESS, formats are ISO standards, so they are documented
12: https://stackoverflow.com/questions/26774037/documentation-or-specification-for-step-and-stp-files
13: http://stepmod.sourceforge.net/express_model_spec/
15: but again it seems that there has been a deliberate effort at obfuscation, probably to raise the bar for entrants.
16: */
18: #ifdef PETSC_HAVE_EGADS
19: #include <petscdmplexegads.h>
21: PETSC_INTERN PetscErrorCode DMSnapToGeomModel_EGADS_Internal(DM, PetscInt, ego, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[], PetscBool);
22: PETSC_INTERN PetscErrorCode DMPlex_Geom_EDGE_XYZtoUV_Internal(const PetscScalar[], ego, const PetscScalar[], const PetscInt, const PetscInt, PetscScalar[], PetscBool);
23: PETSC_INTERN PetscErrorCode DMPlex_Geom_FACE_XYZtoUV_Internal(const PetscScalar[], ego, const PetscScalar[], const PetscInt, const PetscInt, PetscScalar[], PetscBool);
25: PetscErrorCode DMPlex_EGADS_GeomDecode_Internal(const PetscInt geomClass, const PetscInt geomType, char **retClass, char **retType)
26: {
27: PetscFunctionBeginHot;
28: /* EGADS Object Type */
29: if (geomClass == CONTXT) *retClass = (char *)"CONTEXT";
30: if (geomClass == TRANSFORM) *retClass = (char *)"TRANSFORM";
31: if (geomClass == TESSELLATION) *retClass = (char *)"TESSELLATION";
32: if (geomClass == NIL) *retClass = (char *)"NIL";
33: if (geomClass == EMPTY) *retClass = (char *)"EMPTY";
34: if (geomClass == REFERENCE) *retClass = (char *)"REFERENCE";
35: if (geomClass == PCURVE) *retClass = (char *)"PCURVE";
36: if (geomClass == CURVE) *retClass = (char *)"CURVE";
37: if (geomClass == SURFACE) *retClass = (char *)"SURFACE";
38: if (geomClass == NODE) *retClass = (char *)"NODE";
39: if (geomClass == EDGE) *retClass = (char *)"EDGE";
40: if (geomClass == LOOP) *retClass = (char *)"LOOP";
41: if (geomClass == FACE) *retClass = (char *)"FACE";
42: if (geomClass == SHELL) *retClass = (char *)"SHELL";
43: if (geomClass == BODY) *retClass = (char *)"BODY";
44: if (geomClass == MODEL) *retClass = (char *)"MODEL";
46: /* PCURVES & CURVES */
47: if (geomClass == PCURVE || geomClass == CURVE) {
48: if (geomType == LINE) *retType = (char *)"LINE";
49: if (geomType == CIRCLE) *retType = (char *)"CIRCLE";
50: if (geomType == ELLIPSE) *retType = (char *)"ELLIPSE";
51: if (geomType == PARABOLA) *retType = (char *)"PARABOLA";
52: if (geomType == HYPERBOLA) *retType = (char *)"HYPERBOLA";
53: if (geomType == TRIMMED) *retType = (char *)"TRIMMED";
54: if (geomType == BEZIER) *retType = (char *)"BEZIER";
55: if (geomType == BSPLINE) *retType = (char *)"BSPLINE";
56: if (geomType == OFFSET) *retType = (char *)"OFFSET";
57: }
59: /* SURFACE */
60: if (geomClass == SURFACE) {
61: if (geomType == PLANE) *retType = (char *)"PLANE";
62: if (geomType == SPHERICAL) *retType = (char *)"SPHERICAL";
63: if (geomType == CYLINDRICAL) *retType = (char *)"CYLINDRICAL";
64: if (geomType == REVOLUTION) *retType = (char *)"REVOLUTION";
65: if (geomType == TOROIDAL) *retType = (char *)"TOROIDAL";
66: if (geomType == CONICAL) *retType = (char *)"CONICAL";
67: if (geomType == EXTRUSION) *retType = (char *)"EXTRUSION";
68: if (geomType == BEZIER) *retType = (char *)"BEZIER";
69: if (geomType == BSPLINE) *retType = (char *)"BSPLINE";
70: }
72: /* TOPOLOGY */
73: if (geomClass == NODE || geomClass == EDGE || geomClass == LOOP || geomClass == FACE || geomClass == SHELL || geomClass == BODY || geomClass == MODEL) {
74: if (geomType == SREVERSE) *retType = (char *)"SREVERSE";
75: if (geomType == NOMTYPE) *retType = (char *)"NOMTYPE";
76: if (geomType == SFORWARD && geomClass == FACE) *retType = (char *)"SFORWARD";
77: if (geomType == ONENODE && geomClass == EDGE) *retType = (char *)"ONENODE";
78: if (geomType == TWONODE) *retType = (char *)"TWONODE";
79: if (geomType == OPEN) *retType = (char *)"OPEN";
80: if (geomType == CLOSED) *retType = (char *)"CLOSED";
81: if (geomType == DEGENERATE) *retType = (char *)"DEGENERATE";
82: if (geomType == WIREBODY) *retType = (char *)"WIREBODY";
83: if (geomType == FACEBODY) *retType = (char *)"FACEBODY";
84: if (geomType == SHEETBODY) *retType = (char *)"SHEETBODY";
85: if (geomType == SOLIDBODY) *retType = (char *)"SOLIDBODY";
86: }
87: PetscFunctionReturn(PETSC_SUCCESS);
88: }
90: PetscErrorCode DMPlex_EGADS_EDGE_XYZtoUV_Internal(const PetscScalar coords[], ego obj, const PetscScalar range[], const PetscInt v, const PetscInt dE, PetscScalar paramsV[])
91: {
92: //
93: //
94: // Depreciated. Changed all references to DMPlex_Geom_FACE_XYZtoUV_Internal()
95: //
96: //
98: PetscInt loopCntr = 0;
99: PetscScalar dx, dy, dz, lambda, tolr, obj_old, obj_tmp, target;
100: PetscScalar delta, A, b;
101: PetscScalar ts[2], tt[2], eval[18], data[18];
103: PetscFunctionBeginHot;
104: /* Initialize Levenberg-Marquardt parameters */
105: lambda = 1.0;
106: tolr = 1.0;
107: target = 1.0E-20;
108: ts[0] = (range[0] + range[1]) / 2.;
110: while (tolr >= target) {
111: PetscCall(EG_evaluate(obj, ts, eval));
112: dx = coords[v * dE + 0] - eval[0];
113: dy = coords[v * dE + 1] - eval[1];
114: dz = coords[v * dE + 2] - eval[2];
115: obj_old = dx * dx + dy * dy + dz * dz;
117: if (obj_old < target) {
118: tolr = obj_old;
119: break;
120: }
122: A = (eval[3] * eval[3] + eval[4] * eval[4] + eval[5] * eval[5]) * (1.0 + lambda);
123: if (A == 0.0) {
124: PetscCall(PetscPrintf(PETSC_COMM_SELF, "A = 0.0 \n"));
125: break;
126: }
127: b = eval[3] * dx + eval[4] * dy + eval[5] * dz;
129: /* Solve A*delta = b */
130: delta = b / A;
132: /* Find a temp (u,v) and associated objective function */
133: tt[0] = ts[0] + delta;
134: if (tt[0] < range[0]) {
135: tt[0] = range[0];
136: delta = tt[0] - ts[0];
137: }
138: if (tt[0] > range[1]) {
139: tt[0] = range[1];
140: delta = tt[0] - ts[0];
141: }
143: PetscCall(EG_evaluate(obj, tt, data));
145: obj_tmp = (coords[v * dE + 0] - data[0]) * (coords[v * dE + 0] - data[0]) + (coords[v * dE + 1] - data[1]) * (coords[v * dE + 1] - data[1]) + (coords[v * dE + 2] - data[2]) * (coords[v * dE + 2] - data[2]);
147: /* If step is better, accept it and halve lambda (making it more Newton-like) */
148: if (obj_tmp < obj_old) {
149: obj_old = obj_tmp;
150: ts[0] = tt[0];
151: for (int jj = 0; jj < 18; ++jj) eval[jj] = data[jj];
152: lambda /= 2.0;
153: if (lambda < 1.0E-14) lambda = 1.0E-14;
154: if (obj_old < target) {
155: tolr = obj_old;
156: break;
157: }
158: } else {
159: /* Otherwise reject it and double lambda (making it more gradient-descent like) */
160: lambda *= 2.0;
161: }
163: if ((tt[0] == range[0]) || (tt[0] == range[1])) break;
164: if (fabs(delta) < target) {
165: tolr = obj_old;
166: break;
167: }
169: tolr = obj_old;
171: loopCntr += 1;
172: if (loopCntr > 100) break;
173: }
174: paramsV[v * 3 + 0] = ts[0];
175: paramsV[v * 3 + 1] = 0.;
176: PetscFunctionReturn(PETSC_SUCCESS);
177: }
179: PetscErrorCode DMPlex_Geom_EDGE_XYZtoUV_Internal(const PetscScalar coords[], ego obj, const PetscScalar range[], const PetscInt v, const PetscInt dE, PetscScalar paramsV[], PetscBool islite)
180: {
181: PetscInt loopCntr = 0;
182: PetscScalar dx, dy, dz, lambda, tolr, obj_old, obj_tmp, target;
183: PetscScalar delta, A, b;
184: PetscScalar ts[2], tt[2], eval[18], data[18];
186: PetscFunctionBeginHot;
187: /* Initialize Levenberg-Marquardt parameters */
188: lambda = 1.0;
189: tolr = 1.0;
190: target = 1.0E-20;
191: ts[0] = (range[0] + range[1]) / 2.;
193: while (tolr >= target) {
194: if (islite) {
195: PetscCall(EGlite_evaluate(obj, ts, eval));
196: } else {
197: PetscCall(EG_evaluate(obj, ts, eval));
198: }
200: dx = coords[v * dE + 0] - eval[0];
201: dy = coords[v * dE + 1] - eval[1];
202: dz = coords[v * dE + 2] - eval[2];
203: obj_old = dx * dx + dy * dy + dz * dz;
205: if (obj_old < target) {
206: tolr = obj_old;
207: break;
208: }
210: A = (eval[3] * eval[3] + eval[4] * eval[4] + eval[5] * eval[5]) * (1.0 + lambda);
211: if (A == 0.0) {
212: PetscCall(PetscPrintf(PETSC_COMM_SELF, "A = 0.0 \n"));
213: break;
214: }
215: b = eval[3] * dx + eval[4] * dy + eval[5] * dz;
217: /* Solve A*delta = b */
218: delta = b / A;
220: /* Find a temp (u,v) and associated objective function */
221: tt[0] = ts[0] + delta;
222: if (tt[0] < range[0]) {
223: tt[0] = range[0];
224: delta = tt[0] - ts[0];
225: }
226: if (tt[0] > range[1]) {
227: tt[0] = range[1];
228: delta = tt[0] - ts[0];
229: }
231: if (islite) {
232: PetscCall(EGlite_evaluate(obj, tt, data));
233: } else {
234: PetscCall(EG_evaluate(obj, tt, data));
235: }
237: obj_tmp = (coords[v * dE + 0] - data[0]) * (coords[v * dE + 0] - data[0]) + (coords[v * dE + 1] - data[1]) * (coords[v * dE + 1] - data[1]) + (coords[v * dE + 2] - data[2]) * (coords[v * dE + 2] - data[2]);
239: /* If step is better, accept it and halve lambda (making it more Newton-like) */
240: if (obj_tmp < obj_old) {
241: obj_old = obj_tmp;
242: ts[0] = tt[0];
243: for (int jj = 0; jj < 18; ++jj) eval[jj] = data[jj];
244: lambda /= 2.0;
245: if (lambda < 1.0E-14) lambda = 1.0E-14;
246: if (obj_old < target) {
247: tolr = obj_old;
248: break;
249: }
250: } else {
251: /* Otherwise reject it and double lambda (making it more gradient-descent like) */
252: lambda *= 2.0;
253: }
255: if ((tt[0] == range[0]) || (tt[0] == range[1])) break;
256: if (fabs(delta) < target) {
257: tolr = obj_old;
258: break;
259: }
261: tolr = obj_old;
263: loopCntr += 1;
264: if (loopCntr > 100) break;
265: }
266: paramsV[v * 3 + 0] = ts[0];
267: paramsV[v * 3 + 1] = 0.;
268: PetscFunctionReturn(PETSC_SUCCESS);
269: }
271: PetscErrorCode DMPlex_EGADS_FACE_XYZtoUV_Internal(const PetscScalar coords[], ego obj, const PetscScalar range[], const PetscInt v, const PetscInt dE, PetscScalar paramsV[])
272: {
273: //
274: //
275: // Depreciated. Changed all references to DMPlex_Geom_FACE_XYZtoUV_Internal()
276: //
277: //
279: PetscInt loopCntr = 0;
280: PetscScalar dx, dy, dz, lambda, tolr, denom, obj_old, obj_tmp, target;
281: PetscScalar uvs[2], uvt[2], delta[2], A[4], b[2], eval[18], data[18];
283: PetscFunctionBeginHot;
284: /* Initialize Levenberg-Marquardt parameters */
285: lambda = 1.0;
286: tolr = 1.0;
287: target = 1.0E-20;
288: uvs[0] = (range[0] + range[1]) / 2.;
289: uvs[1] = (range[2] + range[3]) / 2.;
291: while (tolr >= target) {
292: PetscCall(EG_evaluate(obj, uvs, eval));
293: dx = coords[v * dE + 0] - eval[0];
294: dy = coords[v * dE + 1] - eval[1];
295: dz = coords[v * dE + 2] - eval[2];
296: obj_old = dx * dx + dy * dy + dz * dz;
298: if (obj_old < target) {
299: tolr = obj_old;
300: break;
301: }
303: A[0] = (eval[3] * eval[3] + eval[4] * eval[4] + eval[5] * eval[5]) * (1.0 + lambda);
304: A[1] = eval[3] * eval[6] + eval[4] * eval[7] + eval[5] * eval[8];
305: A[2] = A[1];
306: A[3] = (eval[6] * eval[6] + eval[7] * eval[7] + eval[8] * eval[8]) * (1.0 + lambda);
308: b[0] = eval[3] * dx + eval[4] * dy + eval[5] * dz;
309: b[1] = eval[6] * dx + eval[7] * dy + eval[8] * dz;
311: /* Solve A*delta = b using Cramer's Rule */
312: denom = A[0] * A[3] - A[2] * A[1];
313: if (denom == 0.0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "denom = 0.0 \n"));
314: delta[0] = (b[0] * A[3] - b[1] * A[1]) / denom;
315: delta[1] = (A[0] * b[1] - A[2] * b[0]) / denom;
317: /* Find a temp (u,v) and associated objective function */
318: uvt[0] = uvs[0] + delta[0];
319: uvt[1] = uvs[1] + delta[1];
321: if (uvt[0] < range[0]) {
322: uvt[0] = range[0];
323: delta[0] = uvt[0] - uvs[0];
324: }
325: if (uvt[0] > range[1]) {
326: uvt[0] = range[1];
327: delta[0] = uvt[0] - uvs[0];
328: }
329: if (uvt[1] < range[2]) {
330: uvt[1] = range[2];
331: delta[1] = uvt[1] - uvs[1];
332: }
333: if (uvt[1] > range[3]) {
334: uvt[1] = range[3];
335: delta[1] = uvt[1] - uvs[1];
336: }
338: PetscCall(EG_evaluate(obj, uvt, data));
340: obj_tmp = (coords[v * dE + 0] - data[0]) * (coords[v * dE + 0] - data[0]) + (coords[v * dE + 1] - data[1]) * (coords[v * dE + 1] - data[1]) + (coords[v * dE + 2] - data[2]) * (coords[v * dE + 2] - data[2]);
342: /* If step is better, accept it and halve lambda (making it more Newton-like) */
343: if (obj_tmp < obj_old) {
344: obj_old = obj_tmp;
345: uvs[0] = uvt[0];
346: uvs[1] = uvt[1];
347: for (int jj = 0; jj < 18; ++jj) eval[jj] = data[jj];
348: lambda /= 2.0;
349: if (lambda < 1.0E-14) lambda = 1.0E-14;
350: if (obj_old < target) {
351: tolr = obj_old;
352: break;
353: }
354: } else {
355: /* Otherwise reject it and double lambda (making it more gradient-descent like) */
356: lambda *= 2.0;
357: }
359: if (sqrt(delta[0] * delta[0] + delta[1] * delta[1]) < target) {
360: tolr = obj_old;
361: break;
362: }
364: tolr = obj_old;
366: loopCntr += 1;
367: if (loopCntr > 100) break;
368: }
369: paramsV[v * 3 + 0] = uvs[0];
370: paramsV[v * 3 + 1] = uvs[1];
371: PetscFunctionReturn(PETSC_SUCCESS);
372: }
374: PetscErrorCode DMPlex_Geom_FACE_XYZtoUV_Internal(const PetscScalar coords[], ego obj, const PetscScalar range[], const PetscInt v, const PetscInt dE, PetscScalar paramsV[], PetscBool islite)
375: {
376: PetscInt loopCntr = 0;
377: PetscScalar dx, dy, dz, lambda, tolr, denom, obj_old, obj_tmp, target;
378: PetscScalar uvs[2], uvt[2], delta[2], A[4], b[2], eval[18], data[18];
380: PetscFunctionBeginHot;
381: /* Initialize Levenberg-Marquardt parameters */
382: lambda = 1.0;
383: tolr = 1.0;
384: target = 1.0E-20;
385: uvs[0] = (range[0] + range[1]) / 2.;
386: uvs[1] = (range[2] + range[3]) / 2.;
388: while (tolr >= target) {
389: if (islite) {
390: PetscCallEGADS(EGlite_evaluate, (obj, uvs, eval));
391: } else {
392: PetscCallEGADS(EG_evaluate, (obj, uvs, eval));
393: }
395: dx = coords[v * dE + 0] - eval[0];
396: dy = coords[v * dE + 1] - eval[1];
397: dz = coords[v * dE + 2] - eval[2];
398: obj_old = dx * dx + dy * dy + dz * dz;
400: if (obj_old < target) {
401: tolr = obj_old;
402: break;
403: }
405: A[0] = (eval[3] * eval[3] + eval[4] * eval[4] + eval[5] * eval[5]) * (1.0 + lambda);
406: A[1] = eval[3] * eval[6] + eval[4] * eval[7] + eval[5] * eval[8];
407: A[2] = A[1];
408: A[3] = (eval[6] * eval[6] + eval[7] * eval[7] + eval[8] * eval[8]) * (1.0 + lambda);
410: b[0] = eval[3] * dx + eval[4] * dy + eval[5] * dz;
411: b[1] = eval[6] * dx + eval[7] * dy + eval[8] * dz;
413: /* Solve A*delta = b using Cramer's Rule */
414: denom = A[0] * A[3] - A[2] * A[1];
415: if (denom == 0.0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "denom = 0.0 \n"));
416: delta[0] = (b[0] * A[3] - b[1] * A[1]) / denom;
417: delta[1] = (A[0] * b[1] - A[2] * b[0]) / denom;
419: /* Find a temp (u,v) and associated objective function */
420: uvt[0] = uvs[0] + delta[0];
421: uvt[1] = uvs[1] + delta[1];
423: if (uvt[0] < range[0]) {
424: uvt[0] = range[0];
425: delta[0] = uvt[0] - uvs[0];
426: }
427: if (uvt[0] > range[1]) {
428: uvt[0] = range[1];
429: delta[0] = uvt[0] - uvs[0];
430: }
431: if (uvt[1] < range[2]) {
432: uvt[1] = range[2];
433: delta[1] = uvt[1] - uvs[1];
434: }
435: if (uvt[1] > range[3]) {
436: uvt[1] = range[3];
437: delta[1] = uvt[1] - uvs[1];
438: }
440: if (islite) {
441: PetscCall(EGlite_evaluate(obj, uvt, data));
442: } else {
443: PetscCall(EG_evaluate(obj, uvt, data));
444: }
446: obj_tmp = (coords[v * dE + 0] - data[0]) * (coords[v * dE + 0] - data[0]) + (coords[v * dE + 1] - data[1]) * (coords[v * dE + 1] - data[1]) + (coords[v * dE + 2] - data[2]) * (coords[v * dE + 2] - data[2]);
448: /* If step is better, accept it and halve lambda (making it more Newton-like) */
449: if (obj_tmp < obj_old) {
450: obj_old = obj_tmp;
451: uvs[0] = uvt[0];
452: uvs[1] = uvt[1];
453: for (int jj = 0; jj < 18; ++jj) eval[jj] = data[jj];
454: lambda /= 2.0;
455: if (lambda < 1.0E-14) lambda = 1.0E-14;
456: if (obj_old < target) {
457: tolr = obj_old;
458: break;
459: }
460: } else {
461: /* Otherwise reject it and double lambda (making it more gradient-descent like) */
462: lambda *= 2.0;
463: }
465: if (sqrt(delta[0] * delta[0] + delta[1] * delta[1]) < target) {
466: tolr = obj_old;
467: break;
468: }
470: tolr = obj_old;
472: loopCntr += 1;
473: if (loopCntr > 100) break;
474: }
475: paramsV[v * 3 + 0] = uvs[0];
476: paramsV[v * 3 + 1] = uvs[1];
477: PetscFunctionReturn(PETSC_SUCCESS);
478: }
480: PetscErrorCode DMSnapToGeomModel_EGADS_Internal(DM dm, PetscInt p, ego model, PetscInt bodyID, PetscInt faceID, PetscInt edgeID, const PetscScalar mcoords[], PetscScalar gcoords[], PetscBool islite)
481: {
482: /* PETSc Variables */
483: DM cdm;
484: ego *bodies;
485: ego geom, body, obj;
486: /* result has to hold derivatives, along with the value */
487: double params[3], result[18], paramsV[16 * 3], range[4];
488: int Nb, oclass, mtype, *senses, peri;
489: Vec coordinatesLocal;
490: PetscScalar *coords = NULL;
491: PetscInt Nv, v, Np = 0, pm;
492: PetscInt dE, d;
493: PetscReal pTolr = 1.0e-14;
495: PetscFunctionBeginHot;
496: PetscCall(DMGetCoordinateDM(dm, &cdm));
497: PetscCall(DMGetCoordinateDim(dm, &dE));
498: PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal));
500: if (islite) {
501: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
502: } else {
503: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
504: }
506: PetscCheck(bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", bodyID, Nb);
507: body = bodies[bodyID];
509: if (edgeID >= 0) {
510: if (islite) {
511: PetscCall(EGlite_objectBodyTopo(body, EDGE, edgeID, &obj));
512: Np = 1;
513: } else {
514: PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &obj));
515: Np = 1;
516: }
517: } else if (faceID >= 0) {
518: if (islite) {
519: PetscCall(EGlite_objectBodyTopo(body, FACE, faceID, &obj));
520: Np = 2;
521: } else {
522: PetscCall(EG_objectBodyTopo(body, FACE, faceID, &obj));
523: Np = 2;
524: }
525: } else {
526: for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
527: PetscFunctionReturn(PETSC_SUCCESS);
528: }
530: /* Calculate parameters (t or u,v) for vertices */
531: PetscCall(DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords));
532: Nv /= dE;
533: if (Nv == 1) {
534: PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords));
535: for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
536: PetscFunctionReturn(PETSC_SUCCESS);
537: }
538: PetscCheck(Nv <= 16, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " coordinates associated to point %" PetscInt_FMT, Nv, p);
540: /* Correct EGADS/EGADSlite 2pi bug when calculating nearest point on Periodic Surfaces */
541: if (islite) {
542: PetscCall(EGlite_getRange(obj, range, &peri));
543: } else {
544: PetscCall(EG_getRange(obj, range, &peri));
545: }
547: for (v = 0; v < Nv; ++v) {
548: if (edgeID > 0) {
549: PetscCall(DMPlex_Geom_EDGE_XYZtoUV_Internal(coords, obj, range, v, dE, paramsV, islite));
550: } else {
551: PetscCall(DMPlex_Geom_FACE_XYZtoUV_Internal(coords, obj, range, v, dE, paramsV, islite));
552: }
553: }
554: PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords));
556: /* Calculate parameters (t or u,v) for new vertex at edge midpoint */
557: for (pm = 0; pm < Np; ++pm) {
558: params[pm] = 0.;
559: for (v = 0; v < Nv; ++v) params[pm] += paramsV[v * 3 + pm];
560: params[pm] /= Nv;
561: }
562: PetscCheck((params[0] + pTolr >= range[0]) || (params[0] - pTolr <= range[1]), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " had bad interpolation", p);
563: PetscCheck(Np < 2 || ((params[1] + pTolr >= range[2]) || (params[1] - pTolr <= range[3])), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d had bad interpolation on v", p);
565: /* Put coordinates for new vertex in result[] */
566: if (islite) {
567: PetscCall(EGlite_evaluate(obj, params, result));
568: } else {
569: PetscCall(EG_evaluate(obj, params, result));
570: }
572: for (d = 0; d < dE; ++d) gcoords[d] = result[d];
573: PetscFunctionReturn(PETSC_SUCCESS);
574: }
575: #endif
577: PetscErrorCode DMSnapToGeomModel_EGADS(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
578: {
579: PetscFunctionBeginHot;
580: #ifdef PETSC_HAVE_EGADS
581: DMLabel bodyLabel, faceLabel, edgeLabel;
582: PetscInt bodyID, faceID, edgeID;
583: PetscContainer modelObj;
584: ego model;
585: PetscBool islite = PETSC_FALSE;
587: // FIXME: Change -dm_plex_refine_without_snap_to_geom to DM to shut off snapping
588: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
589: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
590: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
591: PetscCheck(bodyLabel && faceLabel && edgeLabel, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "EGADS meshes must have body, face, and edge labels defined");
592: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
593: if (!modelObj) {
594: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
595: islite = PETSC_TRUE;
596: }
597: PetscCheck(modelObj, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "EGADS mesh missing model object");
599: PetscCall(PetscContainerGetPointer(modelObj, &model));
600: PetscCall(DMLabelGetValue(bodyLabel, p, &bodyID));
601: PetscCall(DMLabelGetValue(faceLabel, p, &faceID));
602: PetscCall(DMLabelGetValue(edgeLabel, p, &edgeID));
603: /* Allows for "Connective" Plex Edges present in models with multiple non-touching Entities */
604: if (bodyID < 0) {
605: for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
606: PetscFunctionReturn(PETSC_SUCCESS);
607: }
608: PetscCall(DMSnapToGeomModel_EGADS_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords, islite));
609: #endif
610: PetscFunctionReturn(PETSC_SUCCESS);
611: }
613: #if defined(PETSC_HAVE_EGADS)
614: PetscErrorCode DMPlexGeomPrintModel_Internal(ego model, PetscBool islite)
615: {
616: /* PETSc Variables */
617: ego geom, *bodies, *nobjs, *mobjs, *lobjs, *shobjs, *fobjs, *eobjs;
618: int oclass, mtype, *senses, *shsenses, *fsenses, *lsenses, *esenses;
619: int Nb, b;
621: PetscFunctionBeginUser;
622: /* test bodyTopo functions */
623: if (islite) {
624: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
625: } else {
626: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
627: }
628: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %" PetscInt_FMT " \n", Nb));
630: for (b = 0; b < Nb; ++b) {
631: ego body = bodies[b];
632: int id, sh, Nsh, f, Nf, l, Nl, e, Ne, v, Nv;
634: /* List Topology of Bodies */
635: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
636: PetscCall(PetscPrintf(PETSC_COMM_SELF, " BODY %d TOPOLOGY SUMMARY \n", b));
638: /* Output Basic Model Topology */
639: if (islite) {
640: PetscCall(EGlite_getBodyTopos(body, NULL, SHELL, &Nsh, &shobjs));
641: } else {
642: PetscCall(EG_getBodyTopos(body, NULL, SHELL, &Nsh, &shobjs));
643: }
644: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of SHELLS: %d \n", Nsh));
646: if (islite) {
647: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
648: } else {
649: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
650: }
651: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of FACES: %d \n", Nf));
653: if (islite) {
654: PetscCall(EGlite_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
655: } else {
656: PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
657: }
658: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of LOOPS: %d \n", Nl));
660: if (islite) {
661: PetscCall(EGlite_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
662: } else {
663: PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
664: }
665: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of EDGES: %d \n", Ne));
667: if (islite) {
668: PetscCall(EGlite_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
669: } else {
670: PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
671: }
672: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of NODES: %d \n", Nv));
674: if (islite) {
675: EGlite_free(shobjs);
676: EGlite_free(fobjs);
677: EGlite_free(lobjs);
678: EGlite_free(eobjs);
679: EGlite_free(nobjs);
680: } else {
681: EG_free(shobjs);
682: EG_free(fobjs);
683: EG_free(lobjs);
684: EG_free(eobjs);
685: EG_free(nobjs);
686: }
688: /* List Topology of Bodies */
689: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
690: PetscCall(PetscPrintf(PETSC_COMM_SELF, " BODY %d TOPOLOGY DETAILS \n", b));
692: /* Get SHELL info which associated with the current BODY */
693: if (islite) {
694: PetscCall(EGlite_getTopology(body, &geom, &oclass, &mtype, NULL, &Nsh, &shobjs, &shsenses));
695: } else {
696: PetscCall(EG_getTopology(body, &geom, &oclass, &mtype, NULL, &Nsh, &shobjs, &shsenses));
697: }
699: for (sh = 0; sh < Nsh; ++sh) {
700: ego shell = shobjs[sh];
701: int shsense = shsenses[sh];
703: if (islite) {
704: id = EGlite_indexBodyTopo(body, shell);
705: } else {
706: id = EG_indexBodyTopo(body, shell);
707: }
708: PetscCall(PetscPrintf(PETSC_COMM_SELF, " SHELL ID: %d :: sense = %d\n", id, shsense));
710: /* Get FACE information associated with current SHELL */
711: if (islite) {
712: PetscCall(EGlite_getTopology(shell, &geom, &oclass, &mtype, NULL, &Nf, &fobjs, &fsenses));
713: } else {
714: PetscCall(EG_getTopology(shell, &geom, &oclass, &mtype, NULL, &Nf, &fobjs, &fsenses));
715: }
717: for (f = 0; f < Nf; ++f) {
718: ego face = fobjs[f];
719: ego gRef, gPrev, gNext;
720: int goclass, gmtype, *gpinfo;
721: double *gprv;
722: char *gClass = (char *)"", *gType = (char *)"";
723: double fdata[4];
724: ego fRef, fPrev, fNext;
725: int foclass, fmtype;
727: if (islite) {
728: id = EGlite_indexBodyTopo(body, face);
729: } else {
730: id = EG_indexBodyTopo(body, face);
731: }
733: /* Get LOOP info associated with current FACE */
734: if (islite) {
735: PetscCall(EGlite_getTopology(face, &geom, &oclass, &mtype, fdata, &Nl, &lobjs, &lsenses));
736: PetscCall(EGlite_getInfo(face, &foclass, &fmtype, &fRef, &fPrev, &fNext));
737: PetscCall(EGlite_getGeometry(geom, &goclass, &gmtype, &gRef, &gpinfo, &gprv));
738: PetscCall(EGlite_getInfo(geom, &goclass, &gmtype, &gRef, &gPrev, &gNext));
739: } else {
740: PetscCall(EG_getTopology(face, &geom, &oclass, &mtype, fdata, &Nl, &lobjs, &lsenses));
741: PetscCall(EG_getInfo(face, &foclass, &fmtype, &fRef, &fPrev, &fNext));
742: PetscCall(EG_getGeometry(geom, &goclass, &gmtype, &gRef, &gpinfo, &gprv));
743: PetscCall(EG_getInfo(geom, &goclass, &gmtype, &gRef, &gPrev, &gNext));
744: }
746: PetscCall(DMPlex_EGADS_GeomDecode_Internal(goclass, gmtype, &gClass, &gType));
747: PetscCall(PetscPrintf(PETSC_COMM_SELF, " FACE ID: %d :: sense = %d \n", id, fmtype));
748: PetscCall(PetscPrintf(PETSC_COMM_SELF, " GEOMETRY CLASS: %s \n", gClass));
749: PetscCall(PetscPrintf(PETSC_COMM_SELF, " GEOMETRY TYPE: %s \n\n", gType));
750: PetscCall(PetscPrintf(PETSC_COMM_SELF, " RANGE (umin, umax) = (%f, %f) \n", fdata[0], fdata[1]));
751: PetscCall(PetscPrintf(PETSC_COMM_SELF, " (vmin, vmax) = (%f, %f) \n\n", fdata[2], fdata[3]));
753: for (l = 0; l < Nl; ++l) {
754: ego loop = lobjs[l];
755: int lsense = lsenses[l];
757: if (islite) {
758: id = EGlite_indexBodyTopo(body, loop);
759: } else {
760: id = EG_indexBodyTopo(body, loop);
761: }
763: PetscCall(PetscPrintf(PETSC_COMM_SELF, " LOOP ID: %d :: sense = %d\n", id, lsense));
765: /* Get EDGE info associated with the current LOOP */
766: if (islite) {
767: PetscCall(EGlite_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &esenses));
768: } else {
769: PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &esenses));
770: }
772: for (e = 0; e < Ne; ++e) {
773: ego edge = eobjs[e];
774: ego topRef, prev, next;
775: int esense = esenses[e];
776: double range[4] = {0., 0., 0., 0.};
777: int peri;
778: ego gEref, gEprev, gEnext;
779: int gEoclass, gEmtype;
780: char *gEclass = (char *)"", *gEtype = (char *)"";
782: if (islite) {
783: PetscCall(EGlite_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
784: if (mtype != DEGENERATE) PetscCall(EGlite_getInfo(geom, &gEoclass, &gEmtype, &gEref, &gEprev, &gEnext));
785: } else {
786: PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
787: PetscCall(EG_getInfo(geom, &gEoclass, &gEmtype, &gEref, &gEprev, &gEnext));
788: }
790: if (mtype != DEGENERATE) PetscCall(DMPlex_EGADS_GeomDecode_Internal(gEoclass, gEmtype, &gEclass, &gEtype));
792: if (islite) {
793: id = EGlite_indexBodyTopo(body, edge);
794: PetscCall(EGlite_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
795: } else {
796: id = EG_indexBodyTopo(body, edge);
797: PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
798: }
800: PetscCall(PetscPrintf(PETSC_COMM_SELF, " EDGE ID: %d :: sense = %d\n", id, esense));
801: if (mtype != DEGENERATE) {
802: PetscCall(PetscPrintf(PETSC_COMM_SELF, " GEOMETRY CLASS: %s \n", gEclass));
803: PetscCall(PetscPrintf(PETSC_COMM_SELF, " GEOMETRY TYPE: %s \n", gEtype));
804: }
806: if (mtype == DEGENERATE) PetscCall(PetscPrintf(PETSC_COMM_SELF, " EDGE %d is DEGENERATE \n", id));
808: if (islite) {
809: PetscCall(EGlite_getRange(edge, range, &peri));
810: } else {
811: PetscCall(EG_getRange(edge, range, &peri));
812: }
814: PetscCall(PetscPrintf(PETSC_COMM_SELF, " Peri = %d :: Range = %lf, %lf, %lf, %lf \n", peri, range[0], range[1], range[2], range[3]));
816: /* Get NODE info associated with the current EDGE */
817: if (islite) {
818: PetscCall(EGlite_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
819: } else {
820: PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
821: }
823: for (v = 0; v < Nv; ++v) {
824: ego vertex = nobjs[v];
825: double limits[4];
826: int unused;
828: if (islite) {
829: PetscCall(EGlite_getTopology(vertex, &geom, &oclass, &mtype, limits, &unused, &mobjs, &senses));
830: id = EGlite_indexBodyTopo(body, vertex);
831: } else {
832: PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &unused, &mobjs, &senses));
833: id = EG_indexBodyTopo(body, vertex);
834: }
835: PetscCall(PetscPrintf(PETSC_COMM_SELF, " NODE ID: %d \n", id));
836: PetscCall(PetscPrintf(PETSC_COMM_SELF, " (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2]));
837: }
838: }
839: }
840: }
841: }
842: }
843: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
844: PetscFunctionReturn(PETSC_SUCCESS);
845: }
847: static PetscErrorCode DMPlexEGADSDestroy_Private(PetscCtxRt context)
848: {
849: if (*(void **)context) EG_deleteObject((ego) * (void **)context);
850: return PETSC_SUCCESS;
851: }
853: static PetscErrorCode DMPlexEGADSClose_Private(PetscCtxRt context)
854: {
855: if (*(void **)context) EG_close((ego) * (void **)context);
856: return PETSC_SUCCESS;
857: }
859: PetscErrorCode DMPlexEGADSliteDestroy_Private(PetscCtxRt context)
860: {
861: if (*(void **)context) EGlite_deleteObject((ego) * (void **)context);
862: return PETSC_SUCCESS;
863: }
865: PetscErrorCode DMPlexEGADSliteClose_Private(PetscCtxRt context)
866: {
867: if (*(void **)context) EGlite_close((ego) * (void **)context);
868: return PETSC_SUCCESS;
869: }
871: PetscErrorCode DMPlexCreateGeom_Internal(MPI_Comm comm, ego context, ego model, DM *newdm, PetscBool islite)
872: {
873: /* EGADS variables */
874: ego geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
875: int oclass, mtype, nbodies, *senses;
876: int b;
877: /* PETSc variables */
878: DM dm;
879: DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel;
880: PetscHMapI edgeMap = NULL;
881: PetscInt cStart, cEnd, c;
882: PetscInt dim = -1, cdim = -1, numCorners = 0, maxCorners = 0, numVertices = 0, newVertices = 0, numEdges = 0, numCells = 0, newCells = 0, numQuads = 0, cOff = 0, fOff = 0;
883: PetscInt *cells = NULL, *cone = NULL;
884: PetscReal *coords = NULL;
885: PetscMPIInt rank;
887: PetscFunctionBegin;
888: PetscCallMPI(MPI_Comm_rank(comm, &rank));
889: if (rank == 0) {
890: const PetscInt debug = 0;
892: /*
893: Generate PETSc DMPlex
894: Get all Nodes in model, record coordinates in a correctly formatted array
895: Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array
896: We need to uniformly refine the initial geometry to guarantee a valid mesh
897: */
899: /* Calculate cell and vertex sizes */
900: if (islite) {
901: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
902: } else {
903: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
904: }
906: PetscCall(PetscHMapICreate(&edgeMap));
907: numEdges = 0;
908: for (b = 0; b < nbodies; ++b) {
909: ego body = bodies[b];
910: int id, Nl, l, Nv, v;
912: if (islite) {
913: PetscCall(EGlite_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
914: } else {
915: PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
916: }
918: for (l = 0; l < Nl; ++l) {
919: ego loop = lobjs[l];
920: int Ner = 0, Ne, e, Nc;
922: if (islite) {
923: PetscCall(EGlite_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
924: } else {
925: PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
926: }
928: for (e = 0; e < Ne; ++e) {
929: ego edge = objs[e];
930: int Nv, id;
931: PetscHashIter iter;
932: PetscBool found;
934: if (islite) {
935: PetscCall(EGlite_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
936: } else {
937: PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
938: }
940: if (mtype == DEGENERATE) continue;
942: if (islite) {
943: id = EGlite_indexBodyTopo(body, edge);
944: } else {
945: id = EG_indexBodyTopo(body, edge);
946: }
948: PetscCall(PetscHMapIFind(edgeMap, id - 1, &iter, &found));
949: if (!found) PetscCall(PetscHMapISet(edgeMap, id - 1, numEdges++));
950: ++Ner;
951: }
952: if (Ner == 2) {
953: Nc = 2;
954: } else if (Ner == 3) {
955: Nc = 4;
956: } else if (Ner == 4) {
957: Nc = 8;
958: ++numQuads;
959: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot support loop with %d edges", Ner);
960: numCells += Nc;
961: newCells += Nc - 1;
962: maxCorners = PetscMax(Ner * 2 + 1, maxCorners);
963: }
964: if (islite) {
965: PetscCall(EGlite_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
966: } else {
967: PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
968: }
970: for (v = 0; v < Nv; ++v) {
971: ego vertex = nobjs[v];
973: if (islite) {
974: id = EGlite_indexBodyTopo(body, vertex);
975: } else {
976: id = EG_indexBodyTopo(body, vertex);
977: }
978: /* TODO: Instead of assuming contiguous ids, we could use a hash table */
979: numVertices = PetscMax(id, numVertices);
980: }
981: if (islite) {
982: EGlite_free(lobjs);
983: EGlite_free(nobjs);
984: } else {
985: EG_free(lobjs);
986: EG_free(nobjs);
987: }
988: }
989: PetscCall(PetscHMapIGetSize(edgeMap, &numEdges));
990: newVertices = numEdges + numQuads;
991: numVertices += newVertices;
993: dim = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
994: cdim = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
995: numCorners = 3; /* Split cells into triangles */
996: PetscCall(PetscMalloc3(numVertices * cdim, &coords, numCells * numCorners, &cells, maxCorners, &cone));
998: /* Get vertex coordinates */
999: for (b = 0; b < nbodies; ++b) {
1000: ego body = bodies[b];
1001: int id, Nv, v;
1003: if (islite) {
1004: PetscCall(EGlite_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
1005: } else {
1006: PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
1007: }
1009: for (v = 0; v < Nv; ++v) {
1010: ego vertex = nobjs[v];
1011: double limits[4];
1012: int unused;
1014: if (islite) {
1015: PetscCall(EGlite_getTopology(vertex, &geom, &oclass, &mtype, limits, &unused, &mobjs, &senses));
1016: id = EGlite_indexBodyTopo(body, vertex);
1017: } else {
1018: PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &unused, &mobjs, &senses));
1019: id = EG_indexBodyTopo(body, vertex);
1020: }
1022: coords[(id - 1) * cdim + 0] = limits[0];
1023: coords[(id - 1) * cdim + 1] = limits[1];
1024: coords[(id - 1) * cdim + 2] = limits[2];
1025: }
1026: if (islite) {
1027: EGlite_free(nobjs);
1028: } else {
1029: EG_free(nobjs);
1030: }
1031: }
1032: PetscCall(PetscHMapIClear(edgeMap));
1033: fOff = numVertices - newVertices + numEdges;
1034: numEdges = 0;
1035: numQuads = 0;
1036: for (b = 0; b < nbodies; ++b) {
1037: ego body = bodies[b];
1038: int Nl, l;
1040: if (islite) {
1041: PetscCall(EGlite_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
1042: } else {
1043: PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
1044: }
1046: for (l = 0; l < Nl; ++l) {
1047: ego loop = lobjs[l];
1048: int lid, Ner = 0, Ne, e;
1050: if (islite) {
1051: lid = EGlite_indexBodyTopo(body, loop);
1052: PetscCall(EGlite_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
1053: } else {
1054: lid = EG_indexBodyTopo(body, loop);
1055: PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
1056: }
1058: for (e = 0; e < Ne; ++e) {
1059: ego edge = objs[e];
1060: int eid, Nv;
1061: PetscHashIter iter;
1062: PetscBool found;
1064: if (islite) {
1065: PetscCall(EGlite_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1066: } else {
1067: PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1068: }
1070: if (mtype == DEGENERATE) continue;
1071: ++Ner;
1073: if (islite) {
1074: eid = EGlite_indexBodyTopo(body, edge);
1075: } else {
1076: eid = EG_indexBodyTopo(body, edge);
1077: }
1079: PetscCall(PetscHMapIFind(edgeMap, eid - 1, &iter, &found));
1080: if (!found) {
1081: PetscInt v = numVertices - newVertices + numEdges;
1082: double range[4], params[3] = {0., 0., 0.}, result[18];
1083: int periodic[2];
1085: PetscCall(PetscHMapISet(edgeMap, eid - 1, numEdges++));
1087: if (islite) {
1088: PetscCall(EGlite_getRange(edge, range, periodic));
1089: } else {
1090: PetscCall(EG_getRange(edge, range, periodic));
1091: }
1093: params[0] = 0.5 * (range[0] + range[1]);
1094: if (islite) {
1095: PetscCall(EGlite_evaluate(edge, params, result));
1096: } else {
1097: PetscCall(EG_evaluate(edge, params, result));
1098: }
1099: coords[v * cdim + 0] = result[0];
1100: coords[v * cdim + 1] = result[1];
1101: coords[v * cdim + 2] = result[2];
1102: }
1103: }
1104: if (Ner == 4) {
1105: PetscInt v = fOff + numQuads++;
1106: ego *fobjs, face;
1107: double range[4], params[3] = {0., 0., 0.}, result[18];
1108: int Nf, fid, periodic[2];
1110: if (islite) {
1111: PetscCall(EGlite_getBodyTopos(body, loop, FACE, &Nf, &fobjs));
1112: } else {
1113: PetscCall(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs));
1114: }
1115: face = fobjs[0];
1117: if (islite) {
1118: fid = EGlite_indexBodyTopo(body, face);
1119: } else {
1120: fid = EG_indexBodyTopo(body, face);
1121: }
1123: PetscCheck(Nf != 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Loop %d has %" PetscInt_FMT " faces, instead of 1 (%" PetscInt_FMT ")", lid - 1, Nf, fid);
1124: if (islite) {
1125: PetscCall(EGlite_getRange(face, range, periodic));
1126: } else {
1127: PetscCall(EG_getRange(face, range, periodic));
1128: }
1129: params[0] = 0.5 * (range[0] + range[1]);
1130: params[1] = 0.5 * (range[2] + range[3]);
1131: if (islite) {
1132: PetscCall(EGlite_evaluate(face, params, result));
1133: } else {
1134: PetscCall(EG_evaluate(face, params, result));
1135: }
1136: coords[v * cdim + 0] = result[0];
1137: coords[v * cdim + 1] = result[1];
1138: coords[v * cdim + 2] = result[2];
1139: }
1140: }
1141: }
1142: PetscCheck(numEdges + numQuads == newVertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of new vertices %d != %d previous count", numEdges + numQuads, newVertices);
1144: /* Get cell vertices by traversing loops */
1145: numQuads = 0;
1146: cOff = 0;
1147: for (b = 0; b < nbodies; ++b) {
1148: ego body = bodies[b];
1149: int id, Nl, l;
1151: if (islite) {
1152: PetscCall(EGlite_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
1153: } else {
1154: PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
1155: }
1156: for (l = 0; l < Nl; ++l) {
1157: ego loop = lobjs[l];
1158: int lid, Ner = 0, Ne, e, nc = 0, c, Nt, t;
1160: if (islite) {
1161: lid = EGlite_indexBodyTopo(body, loop);
1162: PetscCall(EGlite_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
1163: } else {
1164: lid = EG_indexBodyTopo(body, loop);
1165: PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
1166: }
1168: for (e = 0; e < Ne; ++e) {
1169: ego edge = objs[e];
1170: int points[3];
1171: int eid, Nv, v, tmp;
1173: if (islite) {
1174: eid = EGlite_indexBodyTopo(body, edge);
1175: PetscCall(EGlite_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1176: } else {
1177: eid = EG_indexBodyTopo(body, edge);
1178: PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1179: }
1181: if (mtype == DEGENERATE) continue;
1182: else ++Ner;
1183: PetscCheck(Nv == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %" PetscInt_FMT " has %" PetscInt_FMT " vertices != 2", eid, Nv);
1185: for (v = 0; v < Nv; ++v) {
1186: ego vertex = nobjs[v];
1188: if (islite) {
1189: id = EGlite_indexBodyTopo(body, vertex);
1190: } else {
1191: id = EG_indexBodyTopo(body, vertex);
1192: }
1193: points[v * 2] = id - 1;
1194: }
1195: {
1196: PetscInt edgeNum;
1198: PetscCall(PetscHMapIGet(edgeMap, eid - 1, &edgeNum));
1199: points[1] = numVertices - newVertices + edgeNum;
1200: }
1201: /* EGADS loops are not oriented, but seem to be in order, so we must piece them together */
1202: if (!nc) {
1203: for (v = 0; v < Nv + 1; ++v) cone[nc++] = points[v];
1204: } else {
1205: if (cone[nc - 1] == points[0]) {
1206: cone[nc++] = points[1];
1207: if (cone[0] != points[2]) cone[nc++] = points[2];
1208: } else if (cone[nc - 1] == points[2]) {
1209: cone[nc++] = points[1];
1210: if (cone[0] != points[0]) cone[nc++] = points[0];
1211: } else if (cone[nc - 3] == points[0]) {
1212: tmp = cone[nc - 3];
1213: cone[nc - 3] = cone[nc - 1];
1214: cone[nc - 1] = tmp;
1215: cone[nc++] = points[1];
1216: if (cone[0] != points[2]) cone[nc++] = points[2];
1217: } else if (cone[nc - 3] == points[2]) {
1218: tmp = cone[nc - 3];
1219: cone[nc - 3] = cone[nc - 1];
1220: cone[nc - 1] = tmp;
1221: cone[nc++] = points[1];
1222: if (cone[0] != points[0]) cone[nc++] = points[0];
1223: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %d does not match its predecessor", eid);
1224: }
1225: }
1226: PetscCheck(nc == 2 * Ner, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %" PetscInt_FMT " != %" PetscInt_FMT, nc, 2 * Ner);
1227: if (Ner == 4) cone[nc++] = numVertices - newVertices + numEdges + numQuads++;
1228: PetscCheck(nc <= maxCorners, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %" PetscInt_FMT " > %" PetscInt_FMT " max", nc, maxCorners);
1229: /* Triangulate the loop */
1230: switch (Ner) {
1231: case 2: /* Bi-Segment -> 2 triangles */
1232: Nt = 2;
1233: cells[cOff * numCorners + 0] = cone[0];
1234: cells[cOff * numCorners + 1] = cone[1];
1235: cells[cOff * numCorners + 2] = cone[2];
1236: ++cOff;
1237: cells[cOff * numCorners + 0] = cone[0];
1238: cells[cOff * numCorners + 1] = cone[2];
1239: cells[cOff * numCorners + 2] = cone[3];
1240: ++cOff;
1241: break;
1242: case 3: /* Triangle -> 4 triangles */
1243: Nt = 4;
1244: cells[cOff * numCorners + 0] = cone[0];
1245: cells[cOff * numCorners + 1] = cone[1];
1246: cells[cOff * numCorners + 2] = cone[5];
1247: ++cOff;
1248: cells[cOff * numCorners + 0] = cone[1];
1249: cells[cOff * numCorners + 1] = cone[2];
1250: cells[cOff * numCorners + 2] = cone[3];
1251: ++cOff;
1252: cells[cOff * numCorners + 0] = cone[5];
1253: cells[cOff * numCorners + 1] = cone[3];
1254: cells[cOff * numCorners + 2] = cone[4];
1255: ++cOff;
1256: cells[cOff * numCorners + 0] = cone[1];
1257: cells[cOff * numCorners + 1] = cone[3];
1258: cells[cOff * numCorners + 2] = cone[5];
1259: ++cOff;
1260: break;
1261: case 4: /* Quad -> 8 triangles */
1262: Nt = 8;
1263: cells[cOff * numCorners + 0] = cone[0];
1264: cells[cOff * numCorners + 1] = cone[1];
1265: cells[cOff * numCorners + 2] = cone[7];
1266: ++cOff;
1267: cells[cOff * numCorners + 0] = cone[1];
1268: cells[cOff * numCorners + 1] = cone[2];
1269: cells[cOff * numCorners + 2] = cone[3];
1270: ++cOff;
1271: cells[cOff * numCorners + 0] = cone[3];
1272: cells[cOff * numCorners + 1] = cone[4];
1273: cells[cOff * numCorners + 2] = cone[5];
1274: ++cOff;
1275: cells[cOff * numCorners + 0] = cone[5];
1276: cells[cOff * numCorners + 1] = cone[6];
1277: cells[cOff * numCorners + 2] = cone[7];
1278: ++cOff;
1279: cells[cOff * numCorners + 0] = cone[8];
1280: cells[cOff * numCorners + 1] = cone[1];
1281: cells[cOff * numCorners + 2] = cone[3];
1282: ++cOff;
1283: cells[cOff * numCorners + 0] = cone[8];
1284: cells[cOff * numCorners + 1] = cone[3];
1285: cells[cOff * numCorners + 2] = cone[5];
1286: ++cOff;
1287: cells[cOff * numCorners + 0] = cone[8];
1288: cells[cOff * numCorners + 1] = cone[5];
1289: cells[cOff * numCorners + 2] = cone[7];
1290: ++cOff;
1291: cells[cOff * numCorners + 0] = cone[8];
1292: cells[cOff * numCorners + 1] = cone[7];
1293: cells[cOff * numCorners + 2] = cone[1];
1294: ++cOff;
1295: break;
1296: default:
1297: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d edges, which we do not support", lid, Ner);
1298: }
1299: if (debug) {
1300: for (t = 0; t < Nt; ++t) {
1301: PetscCall(PetscPrintf(PETSC_COMM_SELF, " LOOP Corner NODEs Triangle %d (", t));
1302: for (c = 0; c < numCorners; ++c) {
1303: if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
1304: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%d", cells[(cOff - Nt + t) * numCorners + c]));
1305: }
1306: PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
1307: }
1308: }
1309: }
1310: if (islite) {
1311: EGlite_free(lobjs);
1312: } else {
1313: EG_free(lobjs);
1314: }
1315: }
1316: }
1317: PetscCheck(cOff == numCells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count of total cells %d != %d previous count", cOff, numCells);
1318: PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm));
1319: PetscCall(PetscFree3(coords, cells, cone));
1320: PetscCall(PetscInfo(dm, " Total Number of Unique Cells = %d (%d)\n", numCells, newCells));
1321: PetscCall(PetscInfo(dm, " Total Number of Unique Vertices = %d (%d)\n", numVertices, newVertices));
1322: /* Embed EGADS model in DM */
1323: {
1324: PetscContainer modelObj, contextObj;
1326: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj));
1327: PetscCall(PetscContainerSetPointer(modelObj, model));
1328: PetscCall(PetscContainerSetCtxDestroy(modelObj, (PetscCtxDestroyFn *)DMPlexEGADSDestroy_Private));
1329: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj));
1330: PetscCall(PetscContainerDestroy(&modelObj));
1332: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj));
1333: PetscCall(PetscContainerSetPointer(contextObj, context));
1334: PetscCall(PetscContainerSetCtxDestroy(contextObj, (PetscCtxDestroyFn *)DMPlexEGADSClose_Private));
1335: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj));
1336: PetscCall(PetscContainerDestroy(&contextObj));
1337: }
1338: /* Label points */
1339: PetscCall(DMCreateLabel(dm, "EGADS Body ID"));
1340: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
1341: PetscCall(DMCreateLabel(dm, "EGADS Face ID"));
1342: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
1343: PetscCall(DMCreateLabel(dm, "EGADS Edge ID"));
1344: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
1345: PetscCall(DMCreateLabel(dm, "EGADS Vertex ID"));
1346: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
1347: cOff = 0;
1348: for (b = 0; b < nbodies; ++b) {
1349: ego body = bodies[b];
1350: int id, Nl, l;
1352: if (islite) {
1353: PetscCall(EGlite_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
1354: } else {
1355: PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
1356: }
1357: for (l = 0; l < Nl; ++l) {
1358: ego loop = lobjs[l];
1359: ego *fobjs;
1360: int lid, Nf, fid, Ner = 0, Ne, e, Nt = 0, t;
1362: if (islite) {
1363: lid = EGlite_indexBodyTopo(body, loop);
1364: PetscCall(EGlite_getBodyTopos(body, loop, FACE, &Nf, &fobjs));
1365: } else {
1366: lid = EG_indexBodyTopo(body, loop);
1367: PetscCall(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs));
1368: }
1370: PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf);
1371: if (islite) {
1372: fid = EGlite_indexBodyTopo(body, fobjs[0]);
1373: EGlite_free(fobjs);
1374: PetscCall(EGlite_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
1375: } else {
1376: fid = EG_indexBodyTopo(body, fobjs[0]);
1377: EG_free(fobjs);
1378: PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
1379: }
1381: for (e = 0; e < Ne; ++e) {
1382: ego edge = objs[e];
1383: int eid, Nv, v;
1384: PetscInt points[3], support[2], numEdges, edgeNum;
1385: const PetscInt *edges;
1387: if (islite) {
1388: eid = EGlite_indexBodyTopo(body, edge);
1389: PetscCall(EGlite_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1390: } else {
1391: eid = EG_indexBodyTopo(body, edge);
1392: PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1393: }
1395: if (mtype == DEGENERATE) continue;
1396: else ++Ner;
1397: for (v = 0; v < Nv; ++v) {
1398: ego vertex = nobjs[v];
1400: if (islite) {
1401: id = EGlite_indexBodyTopo(body, vertex);
1402: } else {
1403: id = EG_indexBodyTopo(body, vertex);
1404: }
1406: PetscCall(DMLabelSetValue(edgeLabel, numCells + id - 1, eid));
1407: points[v * 2] = numCells + id - 1;
1408: }
1409: PetscCall(PetscHMapIGet(edgeMap, eid - 1, &edgeNum));
1410: points[1] = numCells + numVertices - newVertices + edgeNum;
1412: PetscCall(DMLabelSetValue(edgeLabel, points[1], eid));
1413: support[0] = points[0];
1414: support[1] = points[1];
1415: PetscCall(DMPlexGetJoin(dm, 2, support, &numEdges, &edges));
1416: PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%d, %d) should only bound 1 edge, not %d", support[0], support[1], numEdges);
1417: PetscCall(DMLabelSetValue(edgeLabel, edges[0], eid));
1418: PetscCall(DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges));
1419: support[0] = points[1];
1420: support[1] = points[2];
1421: PetscCall(DMPlexGetJoin(dm, 2, support, &numEdges, &edges));
1422: PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%d, %d) should only bound 1 edge, not %d", support[0], support[1], numEdges);
1423: PetscCall(DMLabelSetValue(edgeLabel, edges[0], eid));
1424: PetscCall(DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges));
1425: }
1426: switch (Ner) {
1427: case 2:
1428: Nt = 2;
1429: break;
1430: case 3:
1431: Nt = 4;
1432: break;
1433: case 4:
1434: Nt = 8;
1435: break;
1436: default:
1437: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Loop with %d edges is unsupported", Ner);
1438: }
1439: for (t = 0; t < Nt; ++t) {
1440: PetscCall(DMLabelSetValue(bodyLabel, cOff + t, b));
1441: PetscCall(DMLabelSetValue(faceLabel, cOff + t, fid));
1442: }
1443: cOff += Nt;
1444: }
1445: if (islite) {
1446: EGlite_free(lobjs);
1447: } else {
1448: EG_free(lobjs);
1449: }
1450: }
1451: PetscCall(PetscHMapIDestroy(&edgeMap));
1452: PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1453: for (c = cStart; c < cEnd; ++c) {
1454: PetscInt *closure = NULL;
1455: PetscInt clSize, cl, bval, fval;
1457: PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure));
1458: PetscCall(DMLabelGetValue(bodyLabel, c, &bval));
1459: PetscCall(DMLabelGetValue(faceLabel, c, &fval));
1460: for (cl = 0; cl < clSize * 2; cl += 2) {
1461: PetscCall(DMLabelSetValue(bodyLabel, closure[cl], bval));
1462: PetscCall(DMLabelSetValue(faceLabel, closure[cl], fval));
1463: }
1464: PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure));
1465: }
1466: *newdm = dm;
1467: PetscFunctionReturn(PETSC_SUCCESS);
1468: }
1470: PetscErrorCode DMPlexCreateGeom(MPI_Comm comm, ego context, ego model, DM *newdm, PetscBool islite)
1471: {
1472: // EGADS variables
1473: ego geom, *bodies, *mobjs, *fobjs, *lobjs, *eobjs, *nobjs;
1474: ego topRef, prev, next;
1475: int oclass, mtype, nbodies, *senses, *lSenses, *eSenses;
1476: int b;
1477: // PETSc variables
1478: DM dm;
1479: DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel;
1480: PetscHMapI edgeMap = NULL, bodyIndexMap = NULL, bodyVertexMap = NULL, bodyEdgeMap = NULL, bodyFaceMap = NULL, bodyEdgeGlobalMap = NULL;
1481: PetscInt dim = -1, cdim = -1, numCorners = 0, numVertices = 0, numEdges = 0, numFaces = 0, numCells = 0, edgeCntr = 0;
1482: PetscInt cellCntr = 0, numPoints = 0;
1483: PetscInt *cells = NULL;
1484: const PetscInt *cone = NULL;
1485: PetscReal *coords = NULL;
1486: PetscMPIInt rank;
1488: PetscFunctionBeginUser;
1489: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1490: if (rank == 0) {
1491: // Generate PETSc DMPlex
1492: // Get all Nodes in model, record coordinates in a correctly formatted array
1493: // Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array
1494: // We need to uniformly refine the initial geometry to guarantee a valid mesh
1496: // Calculate cell and vertex sizes
1497: if (islite) {
1498: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
1499: } else {
1500: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
1501: }
1502: PetscCall(PetscHMapICreate(&edgeMap));
1503: PetscCall(PetscHMapICreate(&bodyIndexMap));
1504: PetscCall(PetscHMapICreate(&bodyVertexMap));
1505: PetscCall(PetscHMapICreate(&bodyEdgeMap));
1506: PetscCall(PetscHMapICreate(&bodyEdgeGlobalMap));
1507: PetscCall(PetscHMapICreate(&bodyFaceMap));
1509: for (b = 0; b < nbodies; ++b) {
1510: ego body = bodies[b];
1511: int Nf, Ne, Nv;
1512: PetscHashIter BIiter, BViter, BEiter, BEGiter, BFiter, EMiter;
1513: PetscBool BIfound, BVfound, BEfound, BEGfound, BFfound, EMfound;
1515: PetscCall(PetscHMapIFind(bodyIndexMap, b, &BIiter, &BIfound));
1516: PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound));
1517: PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound));
1518: PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound));
1519: PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound));
1521: if (!BIfound) PetscCall(PetscHMapISet(bodyIndexMap, b, numFaces + numEdges + numVertices));
1522: if (!BVfound) PetscCall(PetscHMapISet(bodyVertexMap, b, numVertices));
1523: if (!BEfound) PetscCall(PetscHMapISet(bodyEdgeMap, b, numEdges));
1524: if (!BEGfound) PetscCall(PetscHMapISet(bodyEdgeGlobalMap, b, edgeCntr));
1525: if (!BFfound) PetscCall(PetscHMapISet(bodyFaceMap, b, numFaces));
1527: if (islite) {
1528: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1529: PetscCall(EGlite_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
1530: PetscCall(EGlite_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
1531: EGlite_free(fobjs);
1532: EGlite_free(eobjs);
1533: EGlite_free(nobjs);
1534: } else {
1535: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1536: PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
1537: PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
1538: EG_free(fobjs);
1539: EG_free(eobjs);
1540: EG_free(nobjs);
1541: }
1543: // Remove DEGENERATE EDGES from Edge count
1544: if (islite) {
1545: PetscCall(EGlite_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
1546: } else {
1547: PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
1548: }
1550: int Netemp = 0;
1551: for (int e = 0; e < Ne; ++e) {
1552: ego edge = eobjs[e];
1553: int eid;
1555: if (islite) {
1556: PetscCall(EGlite_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1557: eid = EGlite_indexBodyTopo(body, edge);
1558: } else {
1559: PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1560: eid = EG_indexBodyTopo(body, edge);
1561: }
1563: PetscCall(PetscHMapIFind(edgeMap, edgeCntr + eid - 1, &EMiter, &EMfound));
1564: if (mtype == DEGENERATE) {
1565: if (!EMfound) PetscCall(PetscHMapISet(edgeMap, edgeCntr + eid - 1, -1));
1566: } else {
1567: ++Netemp;
1568: if (!EMfound) PetscCall(PetscHMapISet(edgeMap, edgeCntr + eid - 1, Netemp));
1569: }
1570: }
1571: if (islite) {
1572: EGlite_free(eobjs);
1573: } else {
1574: EG_free(eobjs);
1575: }
1577: // Determine Number of Cells
1578: if (islite) {
1579: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1580: } else {
1581: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1582: }
1584: for (int f = 0; f < Nf; ++f) {
1585: ego face = fobjs[f];
1586: int edgeTemp = 0;
1588: if (islite) {
1589: PetscCall(EGlite_getBodyTopos(body, face, EDGE, &Ne, &eobjs));
1590: } else {
1591: PetscCall(EG_getBodyTopos(body, face, EDGE, &Ne, &eobjs));
1592: }
1594: for (int e = 0; e < Ne; ++e) {
1595: ego edge = eobjs[e];
1597: if (islite) {
1598: PetscCall(EGlite_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1599: } else {
1600: PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1601: }
1602: if (mtype != DEGENERATE) ++edgeTemp;
1603: }
1604: numCells += (2 * edgeTemp);
1605: if (islite) {
1606: EGlite_free(eobjs);
1607: } else {
1608: EG_free(eobjs);
1609: }
1610: }
1611: if (islite) {
1612: EGlite_free(fobjs);
1613: } else {
1614: EG_free(fobjs);
1615: }
1617: numFaces += Nf;
1618: numEdges += Netemp;
1619: numVertices += Nv;
1620: edgeCntr += Ne;
1621: }
1623: // Set up basic DMPlex parameters
1624: dim = 2; // Assumes 3D Models :: Need to handle 2D models in the future
1625: cdim = 3; // Assumes 3D Models :: Need to update to handle 2D models in future
1626: numCorners = 3; // Split Faces into triangles
1627: numPoints = numVertices + numEdges + numFaces; // total number of coordinate points
1629: PetscCall(PetscMalloc2(numPoints * cdim, &coords, numCells * numCorners, &cells));
1631: // Get Vertex Coordinates and Set up Cells
1632: for (b = 0; b < nbodies; ++b) {
1633: ego body = bodies[b];
1634: int Nf, Ne, Nv;
1635: PetscInt bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart;
1636: PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter;
1637: PetscBool BVfound, BEfound, BEGfound, BFfound, EMfound;
1639: // Vertices on Current Body
1640: if (islite) {
1641: PetscCall(EGlite_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
1642: } else {
1643: PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
1644: }
1646: PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound));
1647: PetscCheck(BVfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %" PetscInt_FMT " not found in bodyVertexMap", b);
1648: PetscCall(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart));
1650: for (int v = 0; v < Nv; ++v) {
1651: ego vertex = nobjs[v];
1652: double limits[4];
1653: int id, unused;
1655: if (islite) {
1656: PetscCall(EGlite_getTopology(vertex, &geom, &oclass, &mtype, limits, &unused, &mobjs, &senses));
1657: id = EGlite_indexBodyTopo(body, vertex);
1658: } else {
1659: PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &unused, &mobjs, &senses));
1660: id = EG_indexBodyTopo(body, vertex);
1661: }
1663: coords[(bodyVertexIndexStart + id - 1) * cdim + 0] = limits[0];
1664: coords[(bodyVertexIndexStart + id - 1) * cdim + 1] = limits[1];
1665: coords[(bodyVertexIndexStart + id - 1) * cdim + 2] = limits[2];
1666: }
1667: if (islite) {
1668: EGlite_free(nobjs);
1669: } else {
1670: EG_free(nobjs);
1671: }
1673: // Edge Midpoint Vertices on Current Body
1674: if (islite) {
1675: PetscCall(EGlite_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
1676: } else {
1677: PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
1678: }
1680: PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound));
1681: PetscCheck(BEfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %" PetscInt_FMT " not found in bodyEdgeMap", b);
1682: PetscCall(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart));
1684: PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound));
1685: PetscCheck(BEGfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %" PetscInt_FMT " not found in bodyEdgeGlobalMap", b);
1686: PetscCall(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart));
1688: for (int e = 0; e < Ne; ++e) {
1689: ego edge = eobjs[e];
1690: double range[2], avgt[1], cntrPnt[9];
1691: int eid, eOffset;
1692: int periodic;
1694: if (islite) {
1695: PetscCall(EGlite_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1696: } else {
1697: PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1698: }
1699: if (mtype == DEGENERATE) continue;
1701: if (islite) {
1702: eid = EGlite_indexBodyTopo(body, edge);
1703: } else {
1704: eid = EG_indexBodyTopo(body, edge);
1705: }
1706: // get relative offset from globalEdgeID Vector
1707: PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound));
1708: PetscCheck(EMfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %" PetscInt_FMT " not found in edgeMap", bodyEdgeGlobalIndexStart + eid - 1);
1709: PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset));
1711: if (islite) {
1712: PetscCall(EGlite_getRange(edge, range, &periodic));
1713: } else {
1714: PetscCall(EG_getRange(edge, range, &periodic));
1715: }
1716: avgt[0] = (range[0] + range[1]) / 2.;
1718: if (islite) {
1719: PetscCall(EGlite_evaluate(edge, avgt, cntrPnt));
1720: } else {
1721: PetscCall(EG_evaluate(edge, avgt, cntrPnt));
1722: }
1723: coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 0] = cntrPnt[0];
1724: coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 1] = cntrPnt[1];
1725: coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 2] = cntrPnt[2];
1726: }
1727: if (islite) {
1728: EGlite_free(eobjs);
1729: } else {
1730: EG_free(eobjs);
1731: }
1732: // Face Midpoint Vertices on Current Body
1733: if (islite) {
1734: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1735: } else {
1736: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1737: }
1738: PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound));
1739: PetscCheck(BFfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b);
1740: PetscCall(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart));
1742: for (int f = 0; f < Nf; ++f) {
1743: ego face = fobjs[f];
1744: double range[4], avgUV[2], cntrPnt[18];
1745: int peri, id;
1747: if (islite) {
1748: id = EGlite_indexBodyTopo(body, face);
1749: PetscCall(EGlite_getRange(face, range, &peri));
1750: } else {
1751: id = EG_indexBodyTopo(body, face);
1752: PetscCall(EG_getRange(face, range, &peri));
1753: }
1755: avgUV[0] = (range[0] + range[1]) / 2.;
1756: avgUV[1] = (range[2] + range[3]) / 2.;
1758: if (islite) {
1759: PetscCall(EGlite_evaluate(face, avgUV, cntrPnt));
1760: } else {
1761: PetscCall(EG_evaluate(face, avgUV, cntrPnt));
1762: }
1764: coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 0] = cntrPnt[0];
1765: coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 1] = cntrPnt[1];
1766: coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 2] = cntrPnt[2];
1767: }
1768: if (islite) {
1769: EGlite_free(fobjs);
1770: } else {
1771: EG_free(fobjs);
1772: }
1774: // Define Cells :: Note - This could be incorporated in the Face Midpoint Vertices Loop but was kept separate for clarity
1775: if (islite) {
1776: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1777: } else {
1778: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1779: }
1780: for (int f = 0; f < Nf; ++f) {
1781: ego face = fobjs[f];
1782: int fID, midFaceID, midPntID, startID, endID, Nl;
1784: if (islite) {
1785: fID = EGlite_indexBodyTopo(body, face);
1786: } else {
1787: fID = EG_indexBodyTopo(body, face);
1788: }
1790: midFaceID = numVertices + numEdges + bodyFaceIndexStart + fID - 1;
1791: // Must Traverse Loop to ensure we have all necessary information like the sense (+/- 1) of the edges.
1792: // TODO :: Only handles single loop faces (No holes). The choices for handling multiloop faces are:
1793: // 1) Use the DMPlexCreateGeomFromFile() with the -dm_plex_geom_with_tess = 1 option.
1794: // This will use a default EGADS tessellation as an initial surface mesh.
1795: // 2) Create the initial surface mesh via a 2D mesher :: Currently not available (?future?)
1796: // May I suggest the XXXX as a starting point?
1798: if (islite) {
1799: PetscCall(EGlite_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lSenses));
1800: } else {
1801: PetscCall(EG_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lSenses));
1802: }
1804: PetscCheck(Nl == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face has %" PetscInt_FMT " Loops. Can only handle Faces with 1 Loop. Please use --dm_plex_geom_with_tess = 1 Option", Nl);
1805: for (int l = 0; l < Nl; ++l) {
1806: ego loop = lobjs[l];
1808: if (islite) {
1809: PetscCall(EGlite_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses));
1810: } else {
1811: PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses));
1812: }
1814: for (int e = 0; e < Ne; ++e) {
1815: ego edge = eobjs[e];
1816: int eid, eOffset;
1818: if (islite) {
1819: PetscCall(EGlite_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1820: eid = EGlite_indexBodyTopo(body, edge);
1821: } else {
1822: PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1823: eid = EG_indexBodyTopo(body, edge);
1824: }
1825: if (mtype == DEGENERATE) continue;
1827: // get relative offset from globalEdgeID Vector
1828: PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound));
1829: PetscCheck(EMfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %" PetscInt_FMT " of Body %" PetscInt_FMT " not found in edgeMap. Global Edge ID :: %" PetscInt_FMT, eid, b, bodyEdgeGlobalIndexStart + eid - 1);
1830: PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset));
1832: midPntID = numVertices + bodyEdgeIndexStart + eOffset - 1;
1834: if (islite) {
1835: PetscCall(EGlite_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1836: } else {
1837: PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
1838: }
1840: if (eSenses[e] > 0) {
1841: if (islite) {
1842: startID = EGlite_indexBodyTopo(body, nobjs[0]);
1843: endID = EGlite_indexBodyTopo(body, nobjs[1]);
1844: } else {
1845: startID = EG_indexBodyTopo(body, nobjs[0]);
1846: endID = EG_indexBodyTopo(body, nobjs[1]);
1847: }
1848: } else {
1849: if (islite) {
1850: startID = EGlite_indexBodyTopo(body, nobjs[1]);
1851: endID = EGlite_indexBodyTopo(body, nobjs[0]);
1852: } else {
1853: startID = EG_indexBodyTopo(body, nobjs[1]);
1854: endID = EG_indexBodyTopo(body, nobjs[0]);
1855: }
1856: }
1858: // Define 2 Cells per Edge with correct orientation
1859: cells[cellCntr * numCorners + 0] = midFaceID;
1860: cells[cellCntr * numCorners + 1] = bodyVertexIndexStart + startID - 1;
1861: cells[cellCntr * numCorners + 2] = midPntID;
1863: cells[cellCntr * numCorners + 3] = midFaceID;
1864: cells[cellCntr * numCorners + 4] = midPntID;
1865: cells[cellCntr * numCorners + 5] = bodyVertexIndexStart + endID - 1;
1867: cellCntr = cellCntr + 2;
1868: }
1869: }
1870: }
1871: if (islite) {
1872: EGlite_free(fobjs);
1873: } else {
1874: EG_free(fobjs);
1875: }
1876: }
1877: }
1879: // Generate DMPlex
1880: PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm));
1881: PetscCall(PetscFree2(coords, cells));
1882: PetscCall(PetscInfo(dm, " Total Number of Unique Cells = %" PetscInt_FMT " \n", numCells));
1883: PetscCall(PetscInfo(dm, " Total Number of Unique Vertices = %" PetscInt_FMT " \n", numVertices));
1885: // Embed EGADS model in DM
1886: {
1887: PetscContainer modelObj, contextObj;
1889: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj));
1890: PetscCall(PetscContainerSetPointer(modelObj, model));
1891: if (islite) {
1892: PetscCall(PetscContainerSetCtxDestroy(modelObj, DMPlexEGADSliteDestroy_Private));
1893: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADSlite Model", (PetscObject)modelObj));
1894: } else {
1895: PetscCall(PetscContainerSetCtxDestroy(modelObj, DMPlexEGADSDestroy_Private));
1896: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj));
1897: }
1898: PetscCall(PetscContainerDestroy(&modelObj));
1900: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj));
1901: PetscCall(PetscContainerSetPointer(contextObj, context));
1903: if (islite) {
1904: PetscCall(PetscContainerSetCtxDestroy(contextObj, DMPlexEGADSliteClose_Private));
1905: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADSlite Context", (PetscObject)contextObj));
1906: } else {
1907: PetscCall(PetscContainerSetCtxDestroy(contextObj, DMPlexEGADSClose_Private));
1908: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj));
1909: }
1910: PetscCall(PetscContainerDestroy(&contextObj));
1911: }
1912: // Label points
1913: PetscInt nStart, nEnd;
1915: PetscCall(DMCreateLabel(dm, "EGADS Body ID"));
1916: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
1917: PetscCall(DMCreateLabel(dm, "EGADS Face ID"));
1918: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
1919: PetscCall(DMCreateLabel(dm, "EGADS Edge ID"));
1920: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
1921: PetscCall(DMCreateLabel(dm, "EGADS Vertex ID"));
1922: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
1924: PetscCall(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd));
1926: cellCntr = 0;
1927: for (b = 0; b < nbodies; ++b) {
1928: ego body = bodies[b];
1929: int Nv, Ne, Nf;
1930: PetscInt bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart;
1931: PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter;
1932: PetscBool BVfound, BEfound, BEGfound, BFfound, EMfound;
1934: PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound));
1935: PetscCheck(BVfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b);
1936: PetscCall(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart));
1938: PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound));
1939: PetscCheck(BEfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b);
1940: PetscCall(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart));
1942: PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound));
1943: PetscCheck(BFfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b);
1944: PetscCall(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart));
1946: PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound));
1947: PetscCheck(BEGfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b);
1948: PetscCall(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart));
1950: if (islite) {
1951: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1952: } else {
1953: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1954: }
1956: for (int f = 0; f < Nf; ++f) {
1957: ego face = fobjs[f];
1958: int fID, Nl;
1960: if (islite) {
1961: fID = EGlite_indexBodyTopo(body, face);
1962: PetscCall(EGlite_getBodyTopos(body, face, LOOP, &Nl, &lobjs));
1963: } else {
1964: fID = EG_indexBodyTopo(body, face);
1965: PetscCall(EG_getBodyTopos(body, face, LOOP, &Nl, &lobjs));
1966: }
1968: for (int l = 0; l < Nl; ++l) {
1969: ego loop = lobjs[l];
1970: int lid;
1972: if (islite) {
1973: lid = EGlite_indexBodyTopo(body, loop);
1974: } else {
1975: lid = EG_indexBodyTopo(body, loop);
1976: }
1978: PetscCheck(Nl == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %" PetscInt_FMT " has %" PetscInt_FMT " > 1 faces, which is not supported", lid, Nf);
1980: if (islite) {
1981: PetscCall(EGlite_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses));
1982: } else {
1983: PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses));
1984: }
1986: for (int e = 0; e < Ne; ++e) {
1987: ego edge = eobjs[e];
1988: int eid, eOffset;
1990: // Skip DEGENERATE Edges
1991: if (islite) {
1992: PetscCall(EGlite_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1993: } else {
1994: PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1995: }
1997: if (mtype == DEGENERATE) continue;
1999: if (islite) {
2000: eid = EGlite_indexBodyTopo(body, edge);
2001: } else {
2002: eid = EG_indexBodyTopo(body, edge);
2003: }
2005: // get relative offset from globalEdgeID Vector
2006: PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound));
2007: PetscCheck(EMfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %" PetscInt_FMT " of Body %" PetscInt_FMT " not found in edgeMap. Global Edge ID :: %" PetscInt_FMT, eid, b, bodyEdgeGlobalIndexStart + eid - 1);
2008: PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset));
2010: if (islite) {
2011: PetscCall(EGlite_getBodyTopos(body, edge, NODE, &Nv, &nobjs));
2012: } else {
2013: PetscCall(EG_getBodyTopos(body, edge, NODE, &Nv, &nobjs));
2014: }
2016: for (int v = 0; v < Nv; ++v) {
2017: ego vertex = nobjs[v];
2018: int vID;
2020: if (islite) {
2021: vID = EGlite_indexBodyTopo(body, vertex);
2022: } else {
2023: vID = EG_indexBodyTopo(body, vertex);
2024: }
2026: PetscCall(DMLabelSetValue(bodyLabel, nStart + bodyVertexIndexStart + vID - 1, b));
2027: PetscCall(DMLabelSetValue(vertexLabel, nStart + bodyVertexIndexStart + vID - 1, vID));
2028: }
2029: if (islite) {
2030: EGlite_free(nobjs);
2031: } else {
2032: EG_free(nobjs);
2033: }
2035: PetscCall(DMLabelSetValue(bodyLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, b));
2036: PetscCall(DMLabelSetValue(edgeLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, eid));
2038: // Define Cell faces
2039: for (int jj = 0; jj < 2; ++jj) {
2040: PetscCall(DMLabelSetValue(bodyLabel, cellCntr, b));
2041: PetscCall(DMLabelSetValue(faceLabel, cellCntr, fID));
2042: PetscCall(DMPlexGetCone(dm, cellCntr, &cone));
2044: PetscCall(DMLabelSetValue(bodyLabel, cone[0], b));
2045: PetscCall(DMLabelSetValue(faceLabel, cone[0], fID));
2047: PetscCall(DMLabelSetValue(bodyLabel, cone[1], b));
2048: PetscCall(DMLabelSetValue(edgeLabel, cone[1], eid));
2050: PetscCall(DMLabelSetValue(bodyLabel, cone[2], b));
2051: PetscCall(DMLabelSetValue(faceLabel, cone[2], fID));
2053: cellCntr = cellCntr + 1;
2054: }
2055: }
2056: }
2057: if (islite) {
2058: EGlite_free(lobjs);
2059: } else {
2060: EG_free(lobjs);
2061: }
2063: PetscCall(DMLabelSetValue(bodyLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, b));
2064: PetscCall(DMLabelSetValue(faceLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, fID));
2065: }
2066: if (islite) {
2067: EGlite_free(fobjs);
2068: } else {
2069: EG_free(fobjs);
2070: }
2071: }
2073: PetscCall(PetscHMapIDestroy(&edgeMap));
2074: PetscCall(PetscHMapIDestroy(&bodyIndexMap));
2075: PetscCall(PetscHMapIDestroy(&bodyVertexMap));
2076: PetscCall(PetscHMapIDestroy(&bodyEdgeMap));
2077: PetscCall(PetscHMapIDestroy(&bodyEdgeGlobalMap));
2078: PetscCall(PetscHMapIDestroy(&bodyFaceMap));
2080: *newdm = dm;
2081: PetscFunctionReturn(PETSC_SUCCESS);
2082: }
2084: PetscErrorCode DMPlexCreateGeom_Tess_Internal(MPI_Comm comm, ego context, ego model, DM *newdm, PetscBool islite)
2085: {
2086: /* EGADSlite variables */
2087: ego geom, *bodies, *fobjs;
2088: int b, oclass, mtype, nbodies, *senses;
2089: int totalNumTris = 0, totalNumPoints = 0;
2090: double boundBox[6] = {0., 0., 0., 0., 0., 0.}, tessSize;
2091: /* PETSc variables */
2092: DM dm;
2093: DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel;
2094: PetscHMapI pointIndexStartMap = NULL, triIndexStartMap = NULL, pTypeLabelMap = NULL, pIndexLabelMap = NULL;
2095: PetscHMapI pBodyIndexLabelMap = NULL, triFaceIDLabelMap = NULL, triBodyIDLabelMap = NULL;
2096: PetscInt dim = -1, cdim = -1, numCorners = 0, counter = 0;
2097: PetscInt *cells = NULL;
2098: const PetscInt *cone = NULL;
2099: PetscReal *coords = NULL;
2100: PetscMPIInt rank;
2102: PetscFunctionBeginUser;
2103: PetscCallMPI(MPI_Comm_rank(comm, &rank));
2104: if (rank == 0) {
2105: // Generate PETSc DMPlex from EGADSlite created Tessellation of geometry
2107: // Calculate cell and vertex sizes
2108: if (islite) {
2109: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
2110: } else {
2111: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
2112: }
2114: PetscCall(PetscHMapICreate(&pointIndexStartMap));
2115: PetscCall(PetscHMapICreate(&triIndexStartMap));
2116: PetscCall(PetscHMapICreate(&pTypeLabelMap));
2117: PetscCall(PetscHMapICreate(&pIndexLabelMap));
2118: PetscCall(PetscHMapICreate(&pBodyIndexLabelMap));
2119: PetscCall(PetscHMapICreate(&triFaceIDLabelMap));
2120: PetscCall(PetscHMapICreate(&triBodyIDLabelMap));
2122: /* Create Tessellation of Bodies */
2123: ego *tessArray;
2125: PetscCall(PetscMalloc1(nbodies, &tessArray));
2126: for (b = 0; b < nbodies; ++b) {
2127: ego body = bodies[b];
2128: double params[3] = {0.0, 0.0, 0.0}; // Parameters for Tessellation
2129: int Nf, bodyNumPoints = 0, bodyNumTris = 0;
2130: PetscHashIter PISiter, TISiter;
2131: PetscBool PISfound, TISfound;
2133: /* Store Start Index for each Body's Point and Tris */
2134: PetscCall(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound));
2135: PetscCall(PetscHMapIFind(triIndexStartMap, b, &TISiter, &TISfound));
2137: if (!PISfound) PetscCall(PetscHMapISet(pointIndexStartMap, b, totalNumPoints));
2138: if (!TISfound) PetscCall(PetscHMapISet(triIndexStartMap, b, totalNumTris));
2140: /* Calculate Tessellation parameters based on Bounding Box */
2141: /* Get Bounding Box Dimensions of the BODY */
2142: if (islite) {
2143: PetscCall(EGlite_getBoundingBox(body, boundBox));
2144: } else {
2145: PetscCall(EG_getBoundingBox(body, boundBox));
2146: }
2148: tessSize = boundBox[3] - boundBox[0];
2149: if (tessSize < boundBox[4] - boundBox[1]) tessSize = boundBox[4] - boundBox[1];
2150: if (tessSize < boundBox[5] - boundBox[2]) tessSize = boundBox[5] - boundBox[2];
2152: // TODO :: May want to give users tessellation parameter options //
2153: params[0] = 0.0250 * tessSize;
2154: params[1] = 0.0075 * tessSize;
2155: params[2] = 15.0;
2157: if (islite) {
2158: PetscCall(EGlite_makeTessBody(body, params, &tessArray[b]));
2159: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
2160: } else {
2161: PetscCall(EG_makeTessBody(body, params, &tessArray[b]));
2162: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
2163: }
2165: for (int f = 0; f < Nf; ++f) {
2166: ego face = fobjs[f];
2167: int len, fID, ntris;
2168: const int *ptype, *pindex, *ptris, *ptric;
2169: const double *pxyz, *puv;
2171: // Get Face ID //
2172: if (islite) {
2173: fID = EGlite_indexBodyTopo(body, face);
2174: } else {
2175: fID = EG_indexBodyTopo(body, face);
2176: }
2178: // Checkout the Surface Tessellation //
2179: if (islite) {
2180: PetscCall(EGlite_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric));
2181: } else {
2182: PetscCall(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric));
2183: }
2185: // Determine total number of triangle cells in the tessellation //
2186: bodyNumTris += (int)ntris;
2188: // Check out the point index and coordinate //
2189: for (int p = 0; p < len; ++p) {
2190: int global;
2192: if (islite) {
2193: PetscCall(EGlite_localToGlobal(tessArray[b], fID, p + 1, &global));
2194: } else {
2195: PetscCall(EG_localToGlobal(tessArray[b], fID, p + 1, &global));
2196: }
2198: // Determine the total number of points in the tessellation //
2199: bodyNumPoints = PetscMax(bodyNumPoints, global);
2200: }
2201: }
2202: if (islite) {
2203: EGlite_free(fobjs);
2204: } else {
2205: EG_free(fobjs);
2206: }
2208: totalNumPoints += bodyNumPoints;
2209: totalNumTris += bodyNumTris;
2210: }
2212: dim = 2;
2213: cdim = 3;
2214: numCorners = 3;
2216: /* NEED TO DEFINE MATRICES/VECTORS TO STORE GEOM REFERENCE DATA */
2217: /* Fill in below and use to define DMLabels after DMPlex creation */
2218: PetscCall(PetscMalloc2(totalNumPoints * cdim, &coords, totalNumTris * numCorners, &cells));
2220: for (b = 0; b < nbodies; ++b) {
2221: ego body = bodies[b];
2222: int Nf;
2223: PetscInt pointIndexStart;
2224: PetscHashIter PISiter;
2225: PetscBool PISfound;
2227: PetscCall(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound));
2228: PetscCheck(PISfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %" PetscInt_FMT " not found in pointIndexStartMap", b);
2229: PetscCall(PetscHMapIGet(pointIndexStartMap, b, &pointIndexStart));
2231: if (islite) {
2232: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
2233: } else {
2234: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
2235: }
2237: for (int f = 0; f < Nf; ++f) {
2238: /* Get Face Object */
2239: ego face = fobjs[f];
2240: int len, fID, ntris;
2241: const int *ptype, *pindex, *ptris, *ptric;
2242: const double *pxyz, *puv;
2244: /* Get Face ID */
2245: if (islite) {
2246: fID = EGlite_indexBodyTopo(body, face);
2247: } else {
2248: fID = EG_indexBodyTopo(body, face);
2249: }
2251: /* Checkout the Surface Tessellation */
2252: if (islite) {
2253: PetscCall(EGlite_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric));
2254: } else {
2255: PetscCall(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric));
2256: }
2258: /* Check out the point index and coordinate */
2259: for (int p = 0; p < len; ++p) {
2260: int global;
2261: PetscHashIter PTLiter, PILiter, PBLiter;
2262: PetscBool PTLfound, PILfound, PBLfound;
2264: if (islite) {
2265: PetscCall(EGlite_localToGlobal(tessArray[b], fID, p + 1, &global));
2266: } else {
2267: PetscCall(EG_localToGlobal(tessArray[b], fID, p + 1, &global));
2268: }
2270: /* Set the coordinates array for DAG */
2271: coords[((global - 1 + pointIndexStart) * 3) + 0] = pxyz[(p * 3) + 0];
2272: coords[((global - 1 + pointIndexStart) * 3) + 1] = pxyz[(p * 3) + 1];
2273: coords[((global - 1 + pointIndexStart) * 3) + 2] = pxyz[(p * 3) + 2];
2275: /* Store Geometry Label Information for DMLabel assignment later */
2276: PetscCall(PetscHMapIFind(pTypeLabelMap, global - 1 + pointIndexStart, &PTLiter, &PTLfound));
2277: PetscCall(PetscHMapIFind(pIndexLabelMap, global - 1 + pointIndexStart, &PILiter, &PILfound));
2278: PetscCall(PetscHMapIFind(pBodyIndexLabelMap, global - 1 + pointIndexStart, &PBLiter, &PBLfound));
2280: if (!PTLfound) PetscCall(PetscHMapISet(pTypeLabelMap, global - 1 + pointIndexStart, ptype[p]));
2281: if (!PILfound) PetscCall(PetscHMapISet(pIndexLabelMap, global - 1 + pointIndexStart, pindex[p]));
2282: if (!PBLfound) PetscCall(PetscHMapISet(pBodyIndexLabelMap, global - 1 + pointIndexStart, b));
2284: if (ptype[p] < 0) PetscCall(PetscHMapISet(pIndexLabelMap, global - 1 + pointIndexStart, fID));
2285: }
2287: for (int t = 0; t < (int)ntris; ++t) {
2288: int global, globalA, globalB;
2289: PetscHashIter TFLiter, TBLiter;
2290: PetscBool TFLfound, TBLfound;
2292: if (islite) {
2293: PetscCall(EGlite_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 0], &global));
2294: } else {
2295: PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 0], &global));
2296: }
2297: cells[(counter * 3) + 0] = global - 1 + pointIndexStart;
2299: if (islite) {
2300: PetscCall(EGlite_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 1], &globalA));
2301: } else {
2302: PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 1], &globalA));
2303: }
2304: cells[(counter * 3) + 1] = globalA - 1 + pointIndexStart;
2306: if (islite) {
2307: PetscCall(EGlite_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 2], &globalB));
2308: } else {
2309: PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 2], &globalB));
2310: }
2311: cells[(counter * 3) + 2] = globalB - 1 + pointIndexStart;
2313: PetscCall(PetscHMapIFind(triFaceIDLabelMap, counter, &TFLiter, &TFLfound));
2314: PetscCall(PetscHMapIFind(triBodyIDLabelMap, counter, &TBLiter, &TBLfound));
2316: if (!TFLfound) PetscCall(PetscHMapISet(triFaceIDLabelMap, counter, fID));
2317: if (!TBLfound) PetscCall(PetscHMapISet(triBodyIDLabelMap, counter, b));
2319: counter += 1;
2320: }
2321: }
2322: if (islite) {
2323: EGlite_free(fobjs);
2324: } else {
2325: EG_free(fobjs);
2326: }
2327: }
2328: PetscCall(PetscFree(tessArray));
2329: }
2331: //Build DMPlex
2332: PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, totalNumTris, totalNumPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm));
2333: PetscCall(PetscFree2(coords, cells));
2335: // Embed EGADS model in DM
2336: {
2337: PetscContainer modelObj, contextObj;
2339: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj));
2340: PetscCall(PetscContainerSetPointer(modelObj, model));
2341: if (islite) {
2342: PetscCall(PetscContainerSetCtxDestroy(modelObj, DMPlexEGADSliteDestroy_Private));
2343: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADSlite Model", (PetscObject)modelObj));
2344: } else {
2345: PetscCall(PetscContainerSetCtxDestroy(modelObj, DMPlexEGADSDestroy_Private));
2346: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj));
2347: }
2348: PetscCall(PetscContainerDestroy(&modelObj));
2350: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj));
2351: PetscCall(PetscContainerSetPointer(contextObj, context));
2353: if (islite) {
2354: PetscCall(PetscContainerSetCtxDestroy(contextObj, DMPlexEGADSliteClose_Private));
2355: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADSlite Context", (PetscObject)contextObj));
2356: } else {
2357: PetscCall(PetscContainerSetCtxDestroy(contextObj, DMPlexEGADSClose_Private));
2358: PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj));
2359: }
2360: PetscCall(PetscContainerDestroy(&contextObj));
2361: }
2363: // Label Points
2364: PetscCall(DMCreateLabel(dm, "EGADS Body ID"));
2365: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
2366: PetscCall(DMCreateLabel(dm, "EGADS Face ID"));
2367: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
2368: PetscCall(DMCreateLabel(dm, "EGADS Edge ID"));
2369: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
2370: PetscCall(DMCreateLabel(dm, "EGADS Vertex ID"));
2371: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
2373: /* Get Number of DAG Nodes at each level */
2374: int fStart, fEnd, eStart, eEnd, nStart, nEnd;
2376: PetscCall(DMPlexGetHeightStratum(dm, 0, &fStart, &fEnd));
2377: PetscCall(DMPlexGetHeightStratum(dm, 1, &eStart, &eEnd));
2378: PetscCall(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd));
2380: /* Set DMLabels for NODES */
2381: for (int n = nStart; n < nEnd; ++n) {
2382: int pTypeVal, pIndexVal, pBodyVal;
2383: PetscHashIter PTLiter, PILiter, PBLiter;
2384: PetscBool PTLfound, PILfound, PBLfound;
2386: //Converted to Hash Tables
2387: PetscCall(PetscHMapIFind(pTypeLabelMap, n - nStart, &PTLiter, &PTLfound));
2388: PetscCheck(PTLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %" PetscInt_FMT " not found in pTypeLabelMap", n);
2389: PetscCall(PetscHMapIGet(pTypeLabelMap, n - nStart, &pTypeVal));
2391: PetscCall(PetscHMapIFind(pIndexLabelMap, n - nStart, &PILiter, &PILfound));
2392: PetscCheck(PILfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %" PetscInt_FMT " not found in pIndexLabelMap", n);
2393: PetscCall(PetscHMapIGet(pIndexLabelMap, n - nStart, &pIndexVal));
2395: PetscCall(PetscHMapIFind(pBodyIndexLabelMap, n - nStart, &PBLiter, &PBLfound));
2396: PetscCheck(PBLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %" PetscInt_FMT " not found in pBodyLabelMap", n);
2397: PetscCall(PetscHMapIGet(pBodyIndexLabelMap, n - nStart, &pBodyVal));
2399: PetscCall(DMLabelSetValue(bodyLabel, n, pBodyVal));
2400: if (pTypeVal == 0) PetscCall(DMLabelSetValue(vertexLabel, n, pIndexVal));
2401: if (pTypeVal > 0) PetscCall(DMLabelSetValue(edgeLabel, n, pIndexVal));
2402: if (pTypeVal < 0) PetscCall(DMLabelSetValue(faceLabel, n, pIndexVal));
2403: }
2405: /* Set DMLabels for Edges - Based on the DMLabels of the EDGE's NODES */
2406: for (int e = eStart; e < eEnd; ++e) {
2407: int bodyID_0, vertexID_0, vertexID_1, edgeID_0, edgeID_1, faceID_0, faceID_1;
2409: PetscCall(DMPlexGetCone(dm, e, &cone));
2410: PetscCall(DMLabelGetValue(bodyLabel, cone[0], &bodyID_0)); // Do I need to check the other end?
2411: PetscCall(DMLabelGetValue(vertexLabel, cone[0], &vertexID_0));
2412: PetscCall(DMLabelGetValue(vertexLabel, cone[1], &vertexID_1));
2413: PetscCall(DMLabelGetValue(edgeLabel, cone[0], &edgeID_0));
2414: PetscCall(DMLabelGetValue(edgeLabel, cone[1], &edgeID_1));
2415: PetscCall(DMLabelGetValue(faceLabel, cone[0], &faceID_0));
2416: PetscCall(DMLabelGetValue(faceLabel, cone[1], &faceID_1));
2418: PetscCall(DMLabelSetValue(bodyLabel, e, bodyID_0));
2420: if (edgeID_0 == edgeID_1) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_0));
2421: else if (vertexID_0 > 0 && edgeID_1 > 0) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_1));
2422: else if (vertexID_1 > 0 && edgeID_0 > 0) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_0));
2423: else { /* Do Nothing */ }
2424: }
2426: /* Set DMLabels for Cells */
2427: for (int f = fStart; f < fEnd; ++f) {
2428: int edgeID_0;
2429: PetscInt triBodyVal, triFaceVal;
2430: PetscHashIter TFLiter, TBLiter;
2431: PetscBool TFLfound, TBLfound;
2433: // Convert to Hash Table
2434: PetscCall(PetscHMapIFind(triFaceIDLabelMap, f - fStart, &TFLiter, &TFLfound));
2435: PetscCheck(TFLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %" PetscInt_FMT " not found in triFaceIDLabelMap", f);
2436: PetscCall(PetscHMapIGet(triFaceIDLabelMap, f - fStart, &triFaceVal));
2438: PetscCall(PetscHMapIFind(triBodyIDLabelMap, f - fStart, &TBLiter, &TBLfound));
2439: PetscCheck(TBLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %" PetscInt_FMT " not found in triBodyIDLabelMap", f);
2440: PetscCall(PetscHMapIGet(triBodyIDLabelMap, f - fStart, &triBodyVal));
2442: PetscCall(DMLabelSetValue(bodyLabel, f, triBodyVal));
2443: PetscCall(DMLabelSetValue(faceLabel, f, triFaceVal));
2445: /* Finish Labeling previously unlabeled DMPlex Edges - Assumes Triangular Cell (3 Edges Max) */
2446: PetscCall(DMPlexGetCone(dm, f, &cone));
2448: for (int jj = 0; jj < 3; ++jj) {
2449: PetscCall(DMLabelGetValue(edgeLabel, cone[jj], &edgeID_0));
2451: if (edgeID_0 < 0) {
2452: PetscCall(DMLabelSetValue(bodyLabel, cone[jj], triBodyVal));
2453: PetscCall(DMLabelSetValue(faceLabel, cone[jj], triFaceVal));
2454: }
2455: }
2456: }
2458: *newdm = dm;
2459: PetscFunctionReturn(PETSC_SUCCESS);
2460: }
2461: #endif
2463: /*@C
2464: DMPlexInflateToGeomModelUseXYZ - Snaps the vertex coordinates of a `DMPLEX` object representing the mesh to its geometry if some vertices depart from the model. This usually happens with non-conforming refinement.
2466: Collective
2468: Input Parameter:
2469: . dm - The uninflated `DM` object representing the mesh
2471: Level: intermediate
2473: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()`, `DMPlexCreateEGADS()`
2474: @*/
2475: PetscErrorCode DMPlexInflateToGeomModelUseXYZ(DM dm) PeNS
2476: {
2477: // please don't fucking write code like this with #ifdef all of the place!
2478: #if defined(PETSC_HAVE_EGADS)
2479: /* EGADS Variables */
2480: ego model, geom, body, face, edge, vertex;
2481: ego *bodies;
2482: int Nb, oclass, mtype, *senses;
2483: double result[4];
2484: /* PETSc Variables */
2485: DM cdm;
2486: PetscContainer modelObj;
2487: DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel;
2488: Vec coordinates;
2489: PetscScalar *coords;
2490: PetscInt bodyID, faceID, edgeID, vertexID;
2491: PetscInt cdim, d, vStart, vEnd, v;
2492: PetscBool islite = PETSC_FALSE;
2493: #endif
2495: PetscFunctionBegin;
2496: #if defined(PETSC_HAVE_EGADS)
2497: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
2498: if (!modelObj) {
2499: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
2500: islite = PETSC_TRUE;
2501: }
2502: if (!modelObj) PetscFunctionReturn(PETSC_SUCCESS);
2503: PetscCall(DMGetCoordinateDim(dm, &cdim));
2504: PetscCall(DMGetCoordinateDM(dm, &cdm));
2505: PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2506: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
2507: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
2508: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
2509: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
2511: PetscCall(PetscContainerGetPointer(modelObj, &model));
2513: if (islite) {
2514: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2515: } else {
2516: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2517: }
2519: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2520: PetscCall(VecGetArrayWrite(coordinates, &coords));
2521: for (v = vStart; v < vEnd; ++v) {
2522: PetscScalar *vcoords;
2524: PetscCall(DMLabelGetValue(bodyLabel, v, &bodyID));
2525: PetscCall(DMLabelGetValue(faceLabel, v, &faceID));
2526: PetscCall(DMLabelGetValue(edgeLabel, v, &edgeID));
2527: PetscCall(DMLabelGetValue(vertexLabel, v, &vertexID));
2529: PetscCheck(bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", bodyID, Nb);
2530: body = bodies[bodyID];
2532: PetscCall(DMPlexPointLocalRef(cdm, v, coords, (void *)&vcoords));
2533: if (vertexID > 0) {
2534: if (islite) {
2535: PetscCall(EGlite_objectBodyTopo(body, NODE, vertexID, &vertex));
2536: PetscCall(EGlite_evaluate(vertex, NULL, result));
2537: } else {
2538: PetscCall(EG_objectBodyTopo(body, NODE, vertexID, &vertex));
2539: PetscCall(EG_evaluate(vertex, NULL, result));
2540: }
2541: for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
2542: } else if (edgeID > 0) {
2543: /* Snap to EDGE at nearest location */
2544: double params[1];
2545: if (islite) {
2546: PetscCall(EGlite_objectBodyTopo(body, EDGE, edgeID, &edge));
2547: PetscCall(EGlite_invEvaluate(edge, vcoords, params, result));
2548: } // Get (x,y,z) of nearest point on EDGE
2549: else {
2550: PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &edge));
2551: PetscCall(EG_invEvaluate(edge, vcoords, params, result));
2552: }
2553: for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
2554: } else if (faceID > 0) {
2555: /* Snap to FACE at nearest location */
2556: double params[2];
2557: if (islite) {
2558: PetscCall(EGlite_objectBodyTopo(body, FACE, faceID, &face));
2559: PetscCall(EGlite_invEvaluate(face, vcoords, params, result));
2560: } // Get (x,y,z) of nearest point on FACE
2561: else {
2562: PetscCall(EG_objectBodyTopo(body, FACE, faceID, &face));
2563: PetscCall(EG_invEvaluate(face, vcoords, params, result));
2564: }
2565: for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
2566: }
2567: }
2568: PetscCall(VecRestoreArrayWrite(coordinates, &coords));
2569: /* Clear out global coordinates */
2570: PetscCall(VecDestroy(&dm->coordinates[0].x));
2571: #endif
2572: PetscFunctionReturn(PETSC_SUCCESS);
2573: }
2575: #if defined(PETSC_HAVE_EGADS)
2576: // This replaces the model in-place
2577: PetscErrorCode ConvertGeomModelToAllBSplines(PetscBool islite, ego *model) PeNS
2578: {
2579: /* EGADS/EGADSlite Variables */
2580: ego context = NULL, geom, *bodies, *fobjs;
2581: int oclass, mtype;
2582: int *senses;
2583: int Nb, Nf;
2585: PetscFunctionBegin;
2586: // Get the number of bodies and body objects in the model
2587: if (islite) PetscCallEGADS(EGlite_getTopology, (*model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2588: else PetscCallEGADS(EG_getTopology, (*model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2590: // Get all Faces on the body <-- Only working with 1 body at the moment.
2591: ego body = bodies[0];
2592: if (islite) PetscCallEGADS(EGlite_getBodyTopos, (body, NULL, FACE, &Nf, &fobjs));
2593: else PetscCallEGADS(EG_getBodyTopos, (body, NULL, FACE, &Nf, &fobjs));
2594: ego newGeom[Nf];
2595: ego newFaces[Nf];
2597: // Convert the 1st Face to a BSpline Geometry
2598: for (int ii = 0; ii < Nf; ++ii) {
2599: ego face = fobjs[ii];
2600: ego gRef, gPrev, gNext, *lobjs;
2601: int goclass, gmtype, *gpinfo;
2602: int Nl, *lsenses;
2603: double *gprv;
2604: char *gClass = (char *)"", *gType = (char *)"";
2606: /* Shape Optimization is NOT available for EGADSlite geometry files. */
2607: /* Note :: islite options are left below in case future versions of EGADSlite includes this capability */
2608: PetscCheck(!islite, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot convert geometric entities to all BSplines for geometries defined by EGADSlite (.egadslite)! Please use another geometry file format STEP, IGES, EGADS or BRep");
2610: if (islite) {
2611: PetscCallEGADS(EGlite_getTopology, (face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lsenses)); // Get FACES Geometry object (geom_
2612: PetscCallEGADS(EGlite_getGeometry, (geom, &goclass, &gmtype, &gRef, &gpinfo, &gprv)); // Get geometry object info
2613: PetscCallEGADS(EGlite_getInfo, (geom, &goclass, &gmtype, &gRef, &gPrev, &gNext));
2614: } // Get geometry info
2615: else {
2616: PetscCallEGADS(EG_getTopology, (face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lsenses)); // Get FACES Geometry object (geom_
2617: PetscCallEGADS(EG_getGeometry, (geom, &goclass, &gmtype, &gRef, &gpinfo, &gprv)); // Get geometry object info
2618: PetscCallEGADS(EG_getInfo, (geom, &goclass, &gmtype, &gRef, &gPrev, &gNext));
2619: } // Get geometry info
2621: PetscCall(DMPlex_EGADS_GeomDecode_Internal(goclass, gmtype, &gClass, &gType)); // Decode Geometry integers
2623: // Convert current FACE to a BSpline Surface
2624: ego bspline;
2625: ego bRef, bPrev, bNext;
2626: int boclass, bmtype, *bpinfo;
2627: double *bprv;
2628: char *bClass = (char *)"", *bType = (char *)"";
2630: PetscCallEGADS(EG_convertToBSpline, (face, &bspline)); // Does not have an EGlite_ version
2632: if (islite) {
2633: PetscCallEGADS(EGlite_getGeometry, (bspline, &boclass, &bmtype, &bRef, &bpinfo, &bprv)); // Get geometry object info
2634: PetscCallEGADS(EGlite_getInfo, (bspline, &boclass, &bmtype, &bRef, &bPrev, &bNext));
2635: } // Get geometry info
2636: else {
2637: PetscCallEGADS(EG_getGeometry, (bspline, &boclass, &bmtype, &bRef, &bpinfo, &bprv)); // Get geometry object info
2638: PetscCallEGADS(EG_getInfo, (bspline, &boclass, &bmtype, &bRef, &bPrev, &bNext));
2639: } // Get geometry info
2641: PetscCall(DMPlex_EGADS_GeomDecode_Internal(boclass, bmtype, &bClass, &bType)); // Decode Geometry integers
2643: // Get Context from FACE
2644: context = NULL;
2645: PetscCallEGADS(EG_getContext, (face, &context)); // Does not have an EGlite_ version
2647: // Silence WARNING Regarding OPENCASCADE 7.5
2648: if (islite) PetscCallEGADS(EGlite_setOutLevel, (context, 0));
2649: else PetscCallEGADS(EG_setOutLevel, (context, 0));
2651: ego newgeom;
2652: PetscCallEGADS(EG_makeGeometry, (context, SURFACE, BSPLINE, NULL, bpinfo, bprv, &newgeom)); // Does not have an EGlite_ version
2654: PetscCallEGADS(EG_deleteObject, (bspline));
2656: // Create new FACE based on new SURFACE geometry
2657: double data[4];
2658: int periodic;
2659: if (islite) PetscCallEGADS(EGlite_getRange, (newgeom, data, &periodic));
2660: else PetscCallEGADS(EG_getRange, (newgeom, data, &periodic));
2662: ego newface;
2663: PetscCallEGADS(EG_makeFace, (newgeom, SFORWARD, data, &newface)); // Does not have an EGlite_ version
2664: //PetscCallEGADS(EG_deleteObject, (newgeom));
2665: //PetscCallEGADS(EG_deleteObject, (newface));
2666: newFaces[ii] = newface;
2667: newGeom[ii] = newgeom;
2669: // Reinstate WARNING Regarding OPENCASCADE 7.5
2670: if (islite) PetscCallEGADS(EGlite_setOutLevel, (context, 1));
2671: else PetscCallEGADS(EG_setOutLevel, (context, 1));
2672: }
2674: // Sew New Faces together to get a new model
2675: ego newmodel;
2676: PetscCallEGADS(EG_sewFaces, (Nf, newFaces, 0.0, 0, &newmodel)); // Does not have an EGlite_ version
2677: for (int ii = 0; ii < Nf; ++ii) {
2678: PetscCallEGADS(EG_deleteObject, (newFaces[ii]));
2679: PetscCallEGADS(EG_deleteObject, (newGeom[ii]));
2680: }
2681: PetscCallEGADS(EG_deleteObject, (*model));
2682: *model = newmodel;
2683: PetscFunctionReturn(PETSC_SUCCESS);
2684: }
2685: #endif
2687: /*@C
2688: DMPlexCreateGeomFromFile - Create a `DMPLEX` mesh from an EGADS, IGES, or STEP file.
2690: Collective
2692: Input Parameters:
2693: + comm - The MPI communicator
2694: . filename - The name of the EGADS, IGES, or STEP file
2695: - islite - Flag for EGADSlite support
2697: Output Parameter:
2698: . dm - The `DM` object representing the mesh
2700: Level: beginner
2702: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()`, `DMPlexCreateEGADS()`, `DMPlexCreateEGADSliteFromFile()`
2703: @*/
2704: PetscErrorCode DMPlexCreateGeomFromFile(MPI_Comm comm, const char filename[], DM *dm, PetscBool islite) PeNS
2705: {
2706: /* PETSc Variables */
2707: PetscMPIInt rank;
2708: PetscBool printModel = PETSC_FALSE, tessModel = PETSC_FALSE, newModel = PETSC_FALSE;
2709: PetscBool shapeOpt = PETSC_FALSE;
2711: #if defined(PETSC_HAVE_EGADS)
2712: ego context = NULL, model = NULL;
2713: #endif
2715: PetscFunctionBegin;
2716: PetscAssertPointer(filename, 2);
2717: PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_geom_print_model", &printModel, NULL));
2718: PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_geom_tess_model", &tessModel, NULL));
2719: PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_geom_new_model", &newModel, NULL));
2720: PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_geom_shape_opt", &shapeOpt, NULL));
2721: PetscCallMPI(MPI_Comm_rank(comm, &rank));
2722: #if defined(PETSC_HAVE_EGADS)
2723: if (rank == 0) {
2724: /* EGADSlite files cannot be used for Shape Optimization Work. It lacks the ability to make new geometry. */
2725: /* Must use EGADS, STEP, IGES or BRep files to perform this work. */
2726: if (islite) {
2727: PetscCallEGADS(EGlite_open, (&context));
2728: PetscCallEGADS(EGlite_loadModel, (context, 0, filename, &model));
2729: if (shapeOpt) PetscCall(ConvertGeomModelToAllBSplines(islite, &model));
2730: if (printModel) PetscCall(DMPlexGeomPrintModel_Internal(model, islite));
2731: } else {
2732: PetscCallEGADS(EG_open, (&context));
2733: PetscCallEGADS(EG_loadModel, (context, 0, filename, &model));
2734: if (shapeOpt) PetscCall(ConvertGeomModelToAllBSplines(islite, &model));
2735: if (printModel) PetscCall(DMPlexGeomPrintModel_Internal(model, islite));
2736: }
2737: }
2738: if (tessModel) PetscCall(DMPlexCreateGeom_Tess_Internal(comm, context, model, dm, islite));
2739: else if (newModel) PetscCall(DMPlexCreateGeom_Internal(comm, context, model, dm, islite));
2740: else PetscCall(DMPlexCreateGeom(comm, context, model, dm, islite));
2741: PetscFunctionReturn(PETSC_SUCCESS);
2742: #else
2743: SETERRQ(comm, PETSC_ERR_SUP, "This method requires EGADS support. Reconfigure using --download-egads");
2744: #endif
2745: }
2747: #if defined(PETSC_HAVE_EGADS)
2748: /*@C
2749: DMPlex_Surface_Grad - Exposes the Geometry's Control Points and Weights and Calculates the Mesh Topology Boundary Nodes Gradient
2750: with respect the associated geometry's Control Points and Weights.
2752: // ----- Depreciated ---- See DMPlexGeomDataAndGrads ------ //
2754: Collective
2756: Input Parameters:
2757: . dm - The DM object representing the mesh with PetscContainer containing an EGADS geometry model
2759: Output Parameter:
2760: . dm - The DM object representing the mesh with PetscContainers containing the EGADS geometry model, Array-Hash Table Geometry Control Point Pair, Array-Hash Table Geometry Weights Pair and Matrix-Hash Table Surface Gradient Pair
2762: Level: intermediate
2764: .seealso:
2765: @*/
2766: PetscErrorCode DMPlex_Surface_Grad(DM dm)
2767: {
2768: ego model, geom, *bodies, *fobjs;
2769: PetscContainer modelObj;
2770: int oclass, mtype, *senses;
2771: int Nb, Nf;
2772: PetscHMapI faceCntrlPtRow_Start = NULL, faceCPWeightsRow_Start = NULL;
2773: PetscHMapI pointSurfGradRow_Start = NULL;
2774: Mat pointSurfGrad;
2775: IS faceLabelValues, edgeLabelValues, vertexLabelValues;
2776: PetscInt faceLabelSize, edgeLabelSize, vertexLabelSize;
2777: PetscBool islite = PETSC_FALSE;
2779: PetscFunctionBegin;
2780: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
2781: if (!modelObj) {
2782: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
2783: islite = PETSC_TRUE;
2784: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, " Cannot provide geometric data or associated calculated gradients for geometries defined by EGADSlite (.egadslite)! \n Please use another geometry file format STEP, IGES, EGADS or BRep");
2785: }
2787: // Get attached EGADS model (pointer)
2788: PetscCall(PetscContainerGetPointer(modelObj, &model));
2790: // Get the bodies in the model
2791: if (islite) {
2792: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2793: } else {
2794: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
2795: }
2797: ego body = bodies[0]; // Only operate on 1st body. Model should only have 1 body.
2799: // Get the total number of FACEs in the model
2800: if (islite) {
2801: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
2802: } else {
2803: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
2804: }
2806: // Get the total number of points and IDs in the DMPlex with a "EGADS Face Label"
2807: // This will provide the total number of DMPlex points on the boundary of the geometry
2808: PetscCall(DMGetLabelIdIS(dm, "EGADS Face ID", &faceLabelValues));
2809: PetscCall(DMGetLabelSize(dm, "EGADS Face ID", &faceLabelSize));
2811: PetscCall(DMGetLabelIdIS(dm, "EGADS Edge ID", &edgeLabelValues));
2812: PetscCall(DMGetLabelSize(dm, "EGADS Edge ID", &edgeLabelSize));
2814: PetscCall(DMGetLabelIdIS(dm, "EGADS Vertex ID", &vertexLabelValues));
2815: PetscCall(DMGetLabelSize(dm, "EGADS Vertex ID", &vertexLabelSize));
2817: const PetscInt *faceIndices, *edgeIndices, *vertexIndices;
2818: PetscCall(ISGetIndices(faceLabelValues, &faceIndices));
2819: PetscCall(ISGetIndices(edgeLabelValues, &edgeIndices));
2820: PetscCall(ISGetIndices(vertexLabelValues, &vertexIndices));
2822: // Get the points associated with each FACE, EDGE and VERTEX label in the DM
2823: PetscInt totalNumPoints = 0;
2824: for (int ii = 0; ii < faceLabelSize; ++ii) {
2825: // Cycle through FACE labels
2826: PetscInt size;
2827: PetscCall(DMGetStratumSize(dm, "EGADS Face ID", faceIndices[ii], &size));
2828: totalNumPoints += size;
2829: }
2830: PetscCall(ISRestoreIndices(faceLabelValues, &faceIndices));
2831: PetscCall(ISDestroy(&faceLabelValues));
2833: for (int ii = 0; ii < edgeLabelSize; ++ii) {
2834: // Cycle Through EDGE Labels
2835: PetscInt size;
2836: PetscCall(DMGetStratumSize(dm, "EGADS Edge ID", edgeIndices[ii], &size));
2837: totalNumPoints += size;
2838: }
2839: PetscCall(ISRestoreIndices(edgeLabelValues, &edgeIndices));
2840: PetscCall(ISDestroy(&edgeLabelValues));
2842: for (int ii = 0; ii < vertexLabelSize; ++ii) {
2843: // Cycle Through VERTEX Labels
2844: PetscInt size;
2845: PetscCall(DMGetStratumSize(dm, "EGADS Vertex ID", vertexIndices[ii], &size));
2846: totalNumPoints += size;
2847: }
2848: PetscCall(ISRestoreIndices(vertexLabelValues, &vertexIndices));
2849: PetscCall(ISDestroy(&vertexLabelValues));
2851: int maxNumCPs = 0;
2852: int totalNumCPs = 0;
2853: ego bRef, bPrev, bNext, fgeom, *lobjs;
2854: int id, boclass, bmtype, *bpinfo;
2855: int foclass, fmtype, Nl, *lsenses;
2856: double *bprv;
2857: double fdata[4];
2859: // Create Hash Tables
2860: PetscInt cntr = 0, wcntr = 0;
2861: PetscCall(PetscHMapICreate(&faceCntrlPtRow_Start));
2862: PetscCall(PetscHMapICreate(&faceCPWeightsRow_Start));
2864: for (int ii = 0; ii < Nf; ++ii) {
2865: // Need to get the maximum number of Control Points defining the FACEs
2866: ego face = fobjs[ii];
2867: int maxNumCPs_temp;
2869: if (islite) {
2870: id = EGlite_indexBodyTopo(body, face);
2871: PetscCall(EGlite_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
2872: PetscCall(EGlite_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
2873: PetscCall(EGlite_getInfo(fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
2874: } else {
2875: id = EG_indexBodyTopo(body, face);
2876: PetscCall(EG_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
2877: PetscCall(EG_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
2878: PetscCall(EG_getInfo(fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
2879: }
2881: maxNumCPs_temp = bpinfo[2] * bpinfo[5];
2882: totalNumCPs += bpinfo[2] * bpinfo[5];
2884: if (maxNumCPs_temp > maxNumCPs) maxNumCPs = maxNumCPs_temp;
2885: }
2887: PetscInt *cpCoordDataLengthPtr, *wDataLengthPtr;
2888: PetscInt cpCoordDataLength = 3 * totalNumCPs;
2889: PetscInt wDataLength = totalNumCPs;
2890: cpCoordDataLengthPtr = &cpCoordDataLength;
2891: wDataLengthPtr = &wDataLength;
2892: PetscScalar *cntrlPtCoords, *cntrlPtWeights;
2893: PetscMalloc1(cpCoordDataLength, &cntrlPtCoords);
2894: PetscMalloc1(wDataLength, &cntrlPtWeights);
2895: for (int ii = 0; ii < Nf; ++ii) {
2896: // Need to Populate Control Point Coordinates and Weight Vectors
2897: ego face = fobjs[ii];
2898: PetscHashIter hashKeyIter, wHashKeyIter;
2899: PetscBool hashKeyFound, wHashKeyFound;
2901: if (islite) {
2902: id = EGlite_indexBodyTopo(body, face);
2903: PetscCall(EGlite_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
2904: PetscCall(EGlite_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
2905: PetscCall(EGlite_getInfo(fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
2906: } else {
2907: id = EG_indexBodyTopo(body, face);
2908: PetscCall(EG_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
2909: PetscCall(EG_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
2910: PetscCall(EG_getInfo(fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
2911: }
2913: // Store Face ID to 1st Row of Control Point Vector
2914: PetscCall(PetscHMapIFind(faceCntrlPtRow_Start, id, &hashKeyIter, &hashKeyFound));
2916: if (!hashKeyFound) PetscCall(PetscHMapISet(faceCntrlPtRow_Start, id, cntr));
2918: int offsetCoord = bpinfo[3] + bpinfo[6];
2919: for (int jj = 0; jj < 3 * bpinfo[2] * bpinfo[5]; ++jj) {
2920: cntrlPtCoords[cntr] = bprv[offsetCoord + jj];
2921: cntr += 1;
2922: }
2924: // Store Face ID to 1st Row of Control Point Weight Vector
2925: PetscCall(PetscHMapIFind(faceCPWeightsRow_Start, id, &wHashKeyIter, &wHashKeyFound));
2927: if (!wHashKeyFound) PetscCall(PetscHMapISet(faceCPWeightsRow_Start, id, wcntr));
2929: int offsetWeight = bpinfo[3] + bpinfo[6] + (3 * bpinfo[2] * bpinfo[5]);
2930: for (int jj = 0; jj < bpinfo[2] * bpinfo[5]; ++jj) {
2931: cntrlPtWeights[wcntr] = bprv[offsetWeight + jj];
2932: wcntr += 1;
2933: }
2934: }
2936: // Attach Control Point and Weight Data to DM
2937: {
2938: PetscContainer cpOrgObj, cpCoordObj, cpCoordLengthObj;
2939: PetscContainer wOrgObj, wValObj, wDataLengthObj;
2941: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cpOrgObj));
2942: PetscCall(PetscContainerSetPointer(cpOrgObj, faceCntrlPtRow_Start));
2943: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Hash Table", (PetscObject)cpOrgObj));
2944: PetscCall(PetscContainerDestroy(&cpOrgObj));
2946: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cpCoordObj));
2947: PetscCall(PetscContainerSetPointer(cpCoordObj, cntrlPtCoords));
2948: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Coordinates", (PetscObject)cpCoordObj));
2949: PetscCall(PetscContainerDestroy(&cpCoordObj));
2951: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cpCoordLengthObj));
2952: PetscCall(PetscContainerSetPointer(cpCoordLengthObj, cpCoordDataLengthPtr));
2953: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Coordinate Data Length", (PetscObject)cpCoordLengthObj));
2954: PetscCall(PetscContainerDestroy(&cpCoordLengthObj));
2956: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &wOrgObj));
2957: PetscCall(PetscContainerSetPointer(wOrgObj, faceCPWeightsRow_Start));
2958: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weights Hash Table", (PetscObject)wOrgObj));
2959: PetscCall(PetscContainerDestroy(&wOrgObj));
2961: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &wValObj));
2962: PetscCall(PetscContainerSetPointer(wValObj, cntrlPtWeights));
2963: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weight Data", (PetscObject)wValObj));
2964: PetscCall(PetscContainerDestroy(&wValObj));
2966: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &wDataLengthObj));
2967: PetscCall(PetscContainerSetPointer(wDataLengthObj, wDataLengthPtr));
2968: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weight Data Length", (PetscObject)wDataLengthObj));
2969: PetscCall(PetscContainerDestroy(&wDataLengthObj));
2970: }
2972: // Define Matrix to store Surface Gradient information dx_i/dCPj_i
2973: PetscInt gcntr = 0;
2974: const PetscInt rowSize = 3 * maxNumCPs * totalNumPoints;
2975: const PetscInt colSize = 4 * Nf;
2977: // Create Point Surface Gradient Matrix
2978: MatCreate(PETSC_COMM_WORLD, &pointSurfGrad);
2979: MatSetSizes(pointSurfGrad, PETSC_DECIDE, PETSC_DECIDE, rowSize, colSize);
2980: MatSetType(pointSurfGrad, MATAIJ);
2981: MatSetUp(pointSurfGrad);
2983: // Create Hash Table to store Point's stare row in surfaceGrad[][]
2984: PetscCall(PetscHMapICreate(&pointSurfGradRow_Start));
2986: // Get Coordinates for the DMPlex point
2987: DM cdm;
2988: PetscInt dE, Nv;
2989: Vec coordinatesLocal;
2990: PetscScalar *coords = NULL;
2991: PetscCall(DMGetCoordinateDM(dm, &cdm));
2992: PetscCall(DMGetCoordinateDim(dm, &dE));
2993: PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal));
2995: // CYCLE THROUGH FACEs
2996: for (int ii = 0; ii < Nf; ++ii) {
2997: ego face = fobjs[ii];
2998: ego *eobjs, *nobjs;
2999: PetscInt fid, Ne, Nn;
3000: DMLabel faceLabel, edgeLabel, nodeLabel;
3001: PetscHMapI currFaceUniquePoints = NULL;
3002: IS facePoints, edgePoints, nodePoints;
3003: const PetscInt *fIndices, *eIndices, *nIndices;
3004: PetscInt fSize, eSize, nSize;
3005: PetscHashIter fHashKeyIter, eHashKeyIter, nHashKeyIter, pHashKeyIter;
3006: PetscBool fHashKeyFound, eHashKeyFound, nHashKeyFound, pHashKeyFound;
3007: PetscInt cfCntr = 0;
3009: // Get Geometry Object for the Current FACE
3010: if (islite) {
3011: PetscCall(EGlite_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3012: PetscCall(EGlite_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3013: } else {
3014: PetscCall(EG_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3015: PetscCall(EG_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3016: }
3018: // Get all EDGE and NODE objects attached to the current FACE
3019: if (islite) {
3020: PetscCall(EGlite_getBodyTopos(body, face, EDGE, &Ne, &eobjs));
3021: PetscCall(EGlite_getBodyTopos(body, face, NODE, &Nn, &nobjs));
3022: } else {
3023: PetscCall(EG_getBodyTopos(body, face, EDGE, &Ne, &eobjs));
3024: PetscCall(EG_getBodyTopos(body, face, NODE, &Nn, &nobjs));
3025: }
3027: // Get all DMPlex Points that have DMLabel "EGADS Face ID" and store them in a Hash Table for later use
3028: if (islite) {
3029: fid = EGlite_indexBodyTopo(body, face);
3030: } else {
3031: fid = EG_indexBodyTopo(body, face);
3032: }
3034: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
3035: PetscCall(DMLabelGetStratumIS(faceLabel, fid, &facePoints));
3036: PetscCall(ISGetIndices(facePoints, &fIndices));
3037: PetscCall(ISGetSize(facePoints, &fSize));
3039: PetscCall(PetscHMapICreate(&currFaceUniquePoints));
3041: for (int jj = 0; jj < fSize; ++jj) {
3042: PetscCall(PetscHMapIFind(currFaceUniquePoints, fIndices[jj], &fHashKeyIter, &fHashKeyFound));
3044: if (!fHashKeyFound) {
3045: PetscCall(PetscHMapISet(currFaceUniquePoints, fIndices[jj], cfCntr));
3046: cfCntr += 1;
3047: }
3049: PetscCall(PetscHMapIFind(pointSurfGradRow_Start, fIndices[jj], &pHashKeyIter, &pHashKeyFound));
3051: if (!pHashKeyFound) {
3052: PetscCall(PetscHMapISet(pointSurfGradRow_Start, fIndices[jj], gcntr));
3053: gcntr += 3 * maxNumCPs;
3054: }
3055: }
3056: PetscCall(ISRestoreIndices(facePoints, &fIndices));
3057: PetscCall(ISDestroy(&facePoints));
3059: // Get all DMPlex Points that have DMLable "EGADS Edge ID" attached to the current FACE and store them in a Hash Table for later use.
3060: for (int jj = 0; jj < Ne; ++jj) {
3061: ego edge = eobjs[jj];
3062: PetscBool containLabelValue;
3064: if (islite) {
3065: id = EGlite_indexBodyTopo(body, edge);
3066: } else {
3067: id = EG_indexBodyTopo(body, edge);
3068: }
3070: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
3071: PetscCall(DMLabelHasValue(edgeLabel, id, &containLabelValue));
3073: if (containLabelValue) {
3074: PetscCall(DMLabelGetStratumIS(edgeLabel, id, &edgePoints));
3075: PetscCall(ISGetIndices(edgePoints, &eIndices));
3076: PetscCall(ISGetSize(edgePoints, &eSize));
3078: for (int kk = 0; kk < eSize; ++kk) {
3079: PetscCall(PetscHMapIFind(currFaceUniquePoints, eIndices[kk], &eHashKeyIter, &eHashKeyFound));
3081: if (!eHashKeyFound) {
3082: PetscCall(PetscHMapISet(currFaceUniquePoints, eIndices[kk], cfCntr));
3083: cfCntr += 1;
3084: }
3086: PetscCall(PetscHMapIFind(pointSurfGradRow_Start, eIndices[kk], &pHashKeyIter, &pHashKeyFound));
3088: if (!pHashKeyFound) {
3089: PetscCall(PetscHMapISet(pointSurfGradRow_Start, eIndices[kk], gcntr));
3090: gcntr += 3 * maxNumCPs;
3091: }
3092: }
3093: PetscCall(ISRestoreIndices(edgePoints, &eIndices));
3094: PetscCall(ISDestroy(&edgePoints));
3095: }
3096: }
3098: // Get all DMPlex Points that have DMLabel "EGADS Vertex ID" attached to the current FACE and store them in a Hash Table for later use.
3099: for (int jj = 0; jj < Nn; ++jj) {
3100: ego node = nobjs[jj];
3102: if (islite) {
3103: id = EGlite_indexBodyTopo(body, node);
3104: } else {
3105: id = EG_indexBodyTopo(body, node);
3106: }
3108: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &nodeLabel));
3109: PetscCall(DMLabelGetStratumIS(nodeLabel, id, &nodePoints));
3110: PetscCall(ISGetIndices(nodePoints, &nIndices));
3111: PetscCall(ISGetSize(nodePoints, &nSize));
3113: for (int kk = 0; kk < nSize; ++kk) {
3114: PetscCall(PetscHMapIFind(currFaceUniquePoints, nIndices[kk], &nHashKeyIter, &nHashKeyFound));
3116: if (!nHashKeyFound) {
3117: PetscCall(PetscHMapISet(currFaceUniquePoints, nIndices[kk], cfCntr));
3118: cfCntr += 1;
3119: }
3121: PetscCall(PetscHMapIFind(pointSurfGradRow_Start, nIndices[kk], &pHashKeyIter, &pHashKeyFound));
3122: if (!pHashKeyFound) {
3123: PetscCall(PetscHMapISet(pointSurfGradRow_Start, nIndices[kk], gcntr));
3124: gcntr += 3 * maxNumCPs;
3125: }
3126: }
3127: PetscCall(ISRestoreIndices(nodePoints, &nIndices));
3128: PetscCall(ISDestroy(&nodePoints));
3129: }
3131: // Get the Total Number of entries in the Hash Table
3132: PetscInt currFaceUPSize;
3133: PetscCall(PetscHMapIGetSize(currFaceUniquePoints, &currFaceUPSize));
3135: // Get Keys
3136: PetscInt currFaceUPKeys[currFaceUPSize], off = 0;
3137: PetscCall(PetscHMapIGetKeys(currFaceUniquePoints, &off, currFaceUPKeys));
3139: // Cycle through all points on the current FACE
3140: for (int jj = 0; jj < currFaceUPSize; ++jj) {
3141: PetscInt currPointID = currFaceUPKeys[jj];
3142: PetscCall(DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, currPointID, &Nv, &coords));
3144: // Get UV position of FACE
3145: double params[2], range[4], eval[18];
3146: int peri;
3148: if (islite) {
3149: PetscCall(EGlite_getRange(face, range, &peri));
3150: } else {
3151: PetscCall(EG_getRange(face, range, &peri));
3152: }
3154: PetscCall(DMPlex_Geom_FACE_XYZtoUV_Internal(coords, face, range, 0, dE, params, islite));
3156: if (islite) {
3157: PetscCall(EGlite_evaluate(face, params, eval));
3158: } else {
3159: PetscCall(EG_evaluate(face, params, eval));
3160: }
3162: // Make a new SURFACE Geometry by changing the location of the Control Points
3163: int prvSize = bpinfo[3] + bpinfo[6] + (4 * bpinfo[2] * bpinfo[5]);
3164: double nbprv[prvSize];
3166: // Cycle through each Control Point
3167: double deltaCoord = 1.0E-4;
3168: int offset = bpinfo[3] + bpinfo[6];
3169: int wOffset = offset + (3 * bpinfo[2] * bpinfo[5]);
3170: for (int ii = 0; ii < bpinfo[2] * bpinfo[5]; ++ii) {
3171: // Cycle through each direction (x, then y, then z)
3172: for (int kk = 0; kk < 4; ++kk) {
3173: // Reinitialize nbprv[] values because we only want to change one value at a time
3174: for (int mm = 0; mm < prvSize; ++mm) nbprv[mm] = bprv[mm];
3176: if (kk == 0) { //X
3177: nbprv[offset + 0] = bprv[offset + 0] + deltaCoord;
3178: nbprv[offset + 1] = bprv[offset + 1];
3179: nbprv[offset + 2] = bprv[offset + 2];
3180: } else if (kk == 1) { //Y
3181: nbprv[offset + 0] = bprv[offset + 0];
3182: nbprv[offset + 1] = bprv[offset + 1] + deltaCoord;
3183: nbprv[offset + 2] = bprv[offset + 2];
3184: } else if (kk == 2) { //Z
3185: nbprv[offset + 0] = bprv[offset + 0];
3186: nbprv[offset + 1] = bprv[offset + 1];
3187: nbprv[offset + 2] = bprv[offset + 2] + deltaCoord;
3188: } else if (kk == 3) { // Weights
3189: nbprv[wOffset + ii] = bprv[wOffset + ii] + deltaCoord;
3190: } else {
3191: // currently do nothing
3192: }
3194: // Create New Surface Based on New Control Points or Weights
3195: ego newgeom, context;
3196: if (islite) {
3197: PetscCall(EGlite_open(&context));
3198: PetscCall(EGlite_setOutLevel(context, 0));
3199: } else {
3200: PetscCall(EG_open(&context));
3201: PetscCall(EG_setOutLevel(context, 0));
3202: }
3204: PetscCall(EG_makeGeometry(context, SURFACE, BSPLINE, NULL, bpinfo, nbprv, &newgeom)); // Does not have an EGlite_ version KNOWN_ISSUE
3206: if (islite) {
3207: PetscCall(EGlite_setOutLevel(context, 1));
3208: } else {
3209: PetscCall(EG_setOutLevel(context, 1));
3210: }
3212: // Evaluate new (x, y, z) Point Position based on new Surface Definition
3213: double newCoords[18];
3214: if (islite) {
3215: PetscCall(EGlite_getRange(newgeom, range, &peri));
3216: } else {
3217: PetscCall(EG_getRange(newgeom, range, &peri));
3218: }
3220: PetscCall(DMPlex_Geom_FACE_XYZtoUV_Internal(coords, newgeom, range, 0, dE, params, islite));
3222: if (islite) {
3223: PetscCall(EGlite_evaluate(newgeom, params, newCoords));
3224: } else {
3225: PetscCall(EG_evaluate(newgeom, params, newCoords));
3226: }
3228: // Now Calculate the Surface Gradient for the change in x-component Control Point
3229: PetscScalar dxdCx = (newCoords[0] - coords[0]) / deltaCoord;
3230: PetscScalar dxdCy = (newCoords[1] - coords[1]) / deltaCoord;
3231: PetscScalar dxdCz = (newCoords[2] - coords[2]) / deltaCoord;
3233: // Store Gradient Information in surfaceGrad[][] Matrix
3234: PetscInt startRow;
3235: PetscCall(PetscHMapIGet(pointSurfGradRow_Start, currPointID, &startRow));
3237: // Store Results in PETSc Mat
3238: PetscCall(MatSetValue(pointSurfGrad, startRow + (ii * 3) + 0, ((fid - 1) * 4) + kk, dxdCx, INSERT_VALUES));
3239: PetscCall(MatSetValue(pointSurfGrad, startRow + (ii * 3) + 1, ((fid - 1) * 4) + kk, dxdCy, INSERT_VALUES));
3240: PetscCall(MatSetValue(pointSurfGrad, startRow + (ii * 3) + 2, ((fid - 1) * 4) + kk, dxdCz, INSERT_VALUES));
3241: }
3242: offset += 3;
3243: }
3244: PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, currPointID, &Nv, &coords));
3245: }
3246: }
3248: // Assemble Point Surface Grad Matrix
3249: MatAssemblyBegin(pointSurfGrad, MAT_FINAL_ASSEMBLY);
3250: MatAssemblyEnd(pointSurfGrad, MAT_FINAL_ASSEMBLY);
3252: // Attach Surface Gradient Hash Table and Matrix to DM
3253: {
3254: PetscContainer surfGradOrgObj, surfGradObj;
3256: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &surfGradOrgObj));
3257: PetscCall(PetscContainerSetPointer(surfGradOrgObj, pointSurfGradRow_Start));
3258: PetscCall(PetscObjectCompose((PetscObject)dm, "Surface Gradient Hash Table", (PetscObject)surfGradOrgObj));
3259: PetscCall(PetscContainerDestroy(&surfGradOrgObj));
3261: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &surfGradObj));
3262: PetscCall(PetscContainerSetPointer(surfGradObj, pointSurfGrad));
3263: PetscCall(PetscObjectCompose((PetscObject)dm, "Surface Gradient Matrix", (PetscObject)surfGradObj));
3264: PetscCall(PetscContainerDestroy(&surfGradObj));
3265: }
3266: if (islite) EGlite_free(fobjs);
3267: else EG_free(fobjs);
3268: PetscFunctionReturn(PETSC_SUCCESS);
3269: }
3271: static PetscErrorCode DestroyHashMap(PetscCtxRt p)
3272: {
3273: PetscFunctionBegin;
3274: PetscCall(PetscHMapIDestroy((PetscHMapI *)p));
3275: PetscFunctionReturn(PETSC_SUCCESS);
3276: }
3277: #endif
3279: /*@C
3280: DMPlexGeomDataAndGrads - Exposes Control Points and Control Point Weights defining the underlying geometry allowing user manipulation of the geometry.
3282: Collective
3284: Input Parameters:
3285: + dm - The DM object representing the mesh with PetscContainer containing an EGADS geometry model
3286: - fullGeomGrad - PetscBool flag. Determines how the Surface Area and Volume Gradients wrt to Control Points and Control Point Weights are calculated.
3287: PETSC_FALSE :: Surface Area Gradient wrt Control Points and Control Point Weights are calculated using the change in the local
3288: FACE changes (not the entire body). Volume Gradients are not calculated. Faster computations.
3289: PETSC_TRUE :: Surface Area Gradietn wrt to Control Points and Control Point Weights are calculated using the change observed in
3290: the entire solid body. Volume Gradients are calculated. Slower computation due to the need to generate a new solid
3291: body geometry for every Control Point and Control Point Weight change.
3293: Output Parameter:
3294: . dm - The updated DM object representing the mesh with PetscContainers containing the Control Point, Control Point Weight and Gradient Data.
3296: Level: intermediate
3298: Note:
3299: Calculates the DM Point location, surface area and volume gradients wrt to Control Point and Control Point Weights using Finite Difference (small perturbation of Control Point coordinates or Control Point Weight value).
3301: .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateGeom()`, `DMPlexModifyEGADSGeomModel()`
3302: @*/
3303: PetscErrorCode DMPlexGeomDataAndGrads(DM dm, PetscBool fullGeomGrad) PeNS
3304: {
3305: #if defined(PETSC_HAVE_EGADS)
3306: /* PETSc Variables */
3307: PetscContainer modelObj;
3308: PetscHMapI faceCntrlPtRow_Start = NULL, faceCPWeightsRow_Start = NULL;
3309: PetscHMapI pointSurfGradRow_Start = NULL;
3310: Mat pointSurfGrad, cpEquiv;
3311: IS faceLabelValues, edgeLabelValues, vertexLabelValues;
3312: PetscInt faceLabelSize, edgeLabelSize, vertexLabelSize;
3313: PetscBool islite = PETSC_FALSE;
3314: /* EGADS Variables */
3315: ego model, geom, *bodies, *fobjs = NULL;
3316: int oclass, mtype, *senses;
3317: int Nb, Nf;
3318: #endif
3320: PetscFunctionBegin;
3321: #if defined(PETSC_HAVE_EGADS)
3323: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
3324: if (!modelObj) {
3325: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
3326: PetscCheck(modelObj, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Input DM must have attached EGADS Geometry Model");
3327: islite = PETSC_TRUE;
3328: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot provide geometric data or associated calculated gradients for geometries defined by EGADSlite (.egadslite)!\nPlease use another geometry file format STEP, IGES, EGADS or BRep");
3329: }
3331: // Get attached EGADS model (pointer)
3332: PetscCall(PetscContainerGetPointer(modelObj, &model));
3334: // Get the bodies in the model
3335: if (islite) {
3336: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
3337: } else {
3338: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
3339: }
3341: ego body = bodies[0]; // Only operate on 1st body. Model should only have 1 body.
3343: // Get the total number of FACEs in the model
3344: if (islite) {
3345: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
3346: } else {
3347: PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
3348: }
3350: // Get the total number of points and IDs in the DMPlex with a "EGADS Face Label"
3351: // This will provide the total number of DMPlex points on the boundary of the geometry
3352: PetscCall(DMGetLabelIdIS(dm, "EGADS Face ID", &faceLabelValues));
3353: PetscCall(DMGetLabelSize(dm, "EGADS Face ID", &faceLabelSize));
3355: PetscCall(DMGetLabelIdIS(dm, "EGADS Edge ID", &edgeLabelValues));
3356: PetscCall(DMGetLabelSize(dm, "EGADS Edge ID", &edgeLabelSize));
3358: PetscCall(DMGetLabelIdIS(dm, "EGADS Vertex ID", &vertexLabelValues));
3359: PetscCall(DMGetLabelSize(dm, "EGADS Vertex ID", &vertexLabelSize));
3361: const PetscInt *faceIndices, *edgeIndices, *vertexIndices;
3362: PetscCall(ISGetIndices(faceLabelValues, &faceIndices));
3363: PetscCall(ISGetIndices(edgeLabelValues, &edgeIndices));
3364: PetscCall(ISGetIndices(vertexLabelValues, &vertexIndices));
3366: // Get the points associated with each FACE, EDGE and VERTEX label in the DM
3367: PetscInt totalNumPoints = 0;
3368: for (int f = 0; f < faceLabelSize; ++f) {
3369: // Cycle through FACE labels
3370: PetscInt size;
3371: PetscCall(DMGetStratumSize(dm, "EGADS Face ID", faceIndices[f], &size));
3372: totalNumPoints += size;
3373: }
3374: PetscCall(ISRestoreIndices(faceLabelValues, &faceIndices));
3375: PetscCall(ISDestroy(&faceLabelValues));
3377: for (int e = 0; e < edgeLabelSize; ++e) {
3378: // Cycle Through EDGE Labels
3379: PetscInt size;
3380: PetscCall(DMGetStratumSize(dm, "EGADS Edge ID", edgeIndices[e], &size));
3381: totalNumPoints += size;
3382: }
3383: PetscCall(ISRestoreIndices(edgeLabelValues, &edgeIndices));
3384: PetscCall(ISDestroy(&edgeLabelValues));
3386: for (int ii = 0; ii < vertexLabelSize; ++ii) {
3387: // Cycle Through VERTEX Labels
3388: PetscInt size;
3389: PetscCall(DMGetStratumSize(dm, "EGADS Vertex ID", vertexIndices[ii], &size));
3390: totalNumPoints += size;
3391: }
3392: PetscCall(ISRestoreIndices(vertexLabelValues, &vertexIndices));
3393: PetscCall(ISDestroy(&vertexLabelValues));
3395: int maxNumCPs = 0;
3396: int totalNumCPs = 0;
3397: ego bRef, bPrev, bNext, fgeom, *lobjs;
3398: int id, boclass, bmtype, *bpinfo;
3399: int foclass, fmtype, Nl, *lsenses;
3400: double *bprv;
3401: double fdata[4];
3403: // Create Hash Tables
3404: PetscInt cntr = 0, wcntr = 0, vcntr = 0;
3405: PetscCall(PetscHMapICreate(&faceCntrlPtRow_Start));
3406: PetscCall(PetscHMapICreate(&faceCPWeightsRow_Start));
3408: for (int f = 0; f < Nf; ++f) {
3409: // Need to get the maximum number of Control Points defining the FACEs
3410: ego face = fobjs[f];
3411: int maxNumCPs_temp;
3413: if (islite) {
3414: id = EGlite_indexBodyTopo(body, face);
3415: PetscCall(EGlite_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3416: PetscCall(EGlite_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3417: PetscCall(EGlite_getInfo(fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
3418: } else {
3419: id = EG_indexBodyTopo(body, face);
3420: PetscCall(EG_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3421: PetscCall(EG_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3422: PetscCall(EG_getInfo(fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
3423: }
3424: maxNumCPs_temp = bpinfo[2] * bpinfo[5];
3425: totalNumCPs += bpinfo[2] * bpinfo[5];
3427: if (maxNumCPs_temp > maxNumCPs) maxNumCPs = maxNumCPs_temp;
3428: }
3430: PetscInt *cpCoordDataLengthPtr, *wDataLengthPtr;
3431: PetscInt cpCoordDataLength = 3 * totalNumCPs;
3432: PetscInt wDataLength = totalNumCPs;
3433: cpCoordDataLengthPtr = &cpCoordDataLength;
3434: wDataLengthPtr = &wDataLength;
3436: Vec cntrlPtCoordsVec, cntrlPtWeightsVec;
3437: PetscScalar *cntrlPtCoords, *cntrlPtWeights;
3438: PetscCall(VecCreateSeq(PETSC_COMM_SELF, cpCoordDataLength, &cntrlPtCoordsVec));
3439: PetscCall(VecCreateSeq(PETSC_COMM_SELF, wDataLength, &cntrlPtWeightsVec));
3441: // For dSA/dCPi
3442: Vec gradSACPVec, gradSAWVec, gradVCPVec, gradVWVec;
3443: PetscScalar *gradSACP, *gradSAW, *gradVCP, *gradVW;
3444: PetscCall(VecCreateSeq(PETSC_COMM_SELF, cpCoordDataLength, &gradSACPVec));
3445: PetscCall(VecCreateSeq(PETSC_COMM_SELF, wDataLength, &gradSAWVec));
3446: PetscCall(VecCreateSeq(PETSC_COMM_SELF, cpCoordDataLength, &gradVCPVec));
3447: PetscCall(VecCreateSeq(PETSC_COMM_SELF, wDataLength, &gradVWVec));
3449: // Control Point - Vertex/Edge/Face Relationship
3450: PetscInt *cp_vertex, *cp_edge, *cp_face;
3451: PetscInt *w_vertex, *w_edge, *w_face;
3452: PetscCall(PetscMalloc1(totalNumCPs, &cp_vertex));
3453: PetscCall(PetscMalloc1(totalNumCPs, &cp_edge));
3454: PetscCall(PetscMalloc1(totalNumCPs, &cp_face));
3455: PetscCall(PetscMalloc1(wDataLength, &w_vertex));
3456: PetscCall(PetscMalloc1(wDataLength, &w_edge));
3457: PetscCall(PetscMalloc1(wDataLength, &w_face));
3459: for (int f = 0; f < Nf; ++f) {
3460: // Need to Populate Control Point Coordinates and Weight Vectors
3461: ego face = fobjs[f];
3462: ego *vobjs, *eobjs;
3463: int offsetCoord, offsetWeight;
3464: PetscInt Nv, Ne, wRowStart = 0;
3465: PetscHashIter hashKeyIter, wHashKeyIter;
3466: PetscBool hashKeyFound, wHashKeyFound;
3468: if (islite) {
3469: id = EGlite_indexBodyTopo(body, face);
3470: PetscCallEGADS(EGlite_getTopology, (face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3471: PetscCallEGADS(EGlite_getGeometry, (fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3472: PetscCallEGADS(EGlite_getInfo, (fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
3473: PetscCallEGADS(EGlite_getBodyTopos, (body, face, NODE, &Nv, &vobjs));
3474: } else {
3475: id = EG_indexBodyTopo(body, face);
3476: PetscCallEGADS(EG_getTopology, (face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3477: PetscCallEGADS(EG_getGeometry, (fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3478: PetscCallEGADS(EG_getInfo, (fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
3479: PetscCallEGADS(EG_getBodyTopos, (body, face, NODE, &Nv, &vobjs));
3480: }
3482: // Store Face ID to 1st Row of Control Point Vector
3483: PetscCall(PetscHMapIFind(faceCntrlPtRow_Start, id, &hashKeyIter, &hashKeyFound));
3485: if (!hashKeyFound) PetscCall(PetscHMapISet(faceCntrlPtRow_Start, id, cntr));
3487: PetscCall(VecGetArrayWrite(cntrlPtCoordsVec, &cntrlPtCoords));
3488: offsetCoord = bpinfo[3] + bpinfo[6];
3489: for (int jj = 0; jj < 3 * bpinfo[2] * bpinfo[5]; ++jj) {
3490: cntrlPtCoords[cntr] = bprv[offsetCoord + jj];
3491: cntr += 1;
3492: }
3494: // Store Face ID to 1st Row of Control Point Weight Vector
3495: PetscCall(PetscHMapIFind(faceCPWeightsRow_Start, id, &wHashKeyIter, &wHashKeyFound));
3497: if (!wHashKeyFound) {
3498: PetscCall(PetscHMapISet(faceCPWeightsRow_Start, id, wcntr));
3499: wRowStart = wcntr;
3500: }
3502: PetscCall(VecGetArrayWrite(cntrlPtWeightsVec, &cntrlPtWeights));
3503: offsetWeight = bpinfo[3] + bpinfo[6] + (3 * bpinfo[2] * bpinfo[5]);
3504: for (int jj = 0; jj < bpinfo[2] * bpinfo[5]; ++jj) {
3505: cntrlPtWeights[wcntr] = bprv[offsetWeight + jj];
3506: cp_face[wcntr] = id;
3507: w_face[wcntr] = id;
3508: wcntr += 1;
3509: }
3510: PetscCall(VecRestoreArrayWrite(cntrlPtWeightsVec, &cntrlPtWeights));
3512: // Associate Control Points with Vertex IDs
3513: PetscScalar xcp, ycp, zcp;
3514: offsetCoord = bpinfo[3] + bpinfo[6];
3515: for (int jj = 0; jj < 3 * bpinfo[2] * bpinfo[5]; jj += 3) {
3516: xcp = bprv[offsetCoord + jj + 0];
3517: ycp = bprv[offsetCoord + jj + 1];
3518: zcp = bprv[offsetCoord + jj + 2];
3520: //Initialize Control Point and Weight to Vertex ID relationship to -1
3521: cp_vertex[vcntr] = -1;
3522: w_vertex[vcntr] = -1;
3523: cp_edge[vcntr] = -1;
3524: w_edge[vcntr] = -1;
3526: for (int kk = 0; kk < Nv; ++kk) {
3527: int vid;
3528: double vCoords[3];
3529: PetscScalar vDelta;
3530: ego vertex = vobjs[kk];
3532: if (islite) {
3533: vid = EGlite_indexBodyTopo(body, vertex);
3534: PetscCallEGADS(EGlite_evaluate, (vertex, NULL, vCoords));
3535: } else {
3536: vid = EG_indexBodyTopo(body, vertex);
3537: PetscCallEGADS(EG_evaluate, (vertex, NULL, vCoords));
3538: }
3539: vDelta = PetscSqrtReal(PetscSqr(vCoords[0] - xcp) + PetscSqr(vCoords[1] - ycp) + PetscSqr(vCoords[2] - zcp));
3541: if (vDelta < 1.0E-15) {
3542: cp_vertex[vcntr] = vid;
3543: w_vertex[vcntr] = vid;
3544: }
3545: }
3546: vcntr += 1;
3547: }
3548: // These two line could be replaced with DMPlexFreeGeomObject()
3549: if (islite) EGlite_free(vobjs);
3550: else EG_free(vobjs);
3552: // Associate Control Points with Edge IDs
3553: if (islite) PetscCallEGADS(EGlite_getBodyTopos, (body, face, EDGE, &Ne, &eobjs));
3554: else PetscCallEGADS(EG_getBodyTopos, (body, face, EDGE, &Ne, &eobjs));
3556: int cpV1, cpV2;
3557: int minID, maxID;
3559: // Along vmin axis
3560: minID = wRowStart;
3561: maxID = wRowStart + (bpinfo[2] - 1);
3562: cpV1 = cp_vertex[minID];
3563: cpV2 = cp_vertex[maxID];
3564: for (int jj = 0; jj < Ne; ++jj) {
3565: ego edge = eobjs[jj];
3566: ego egeom, *nobjs;
3567: int eoclass, emtype, Nn, *nsenses;
3568: int n1ID, n2ID, eid;
3570: if (islite) {
3571: eid = EGlite_indexBodyTopo(body, edge);
3572: PetscCallEGADS(EGlite_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3573: } else {
3574: eid = EG_indexBodyTopo(body, edge);
3575: PetscCallEGADS(EG_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3576: }
3578: if (emtype != DEGENERATE) {
3579: // Get IDs for current Edge's End Vertices
3580: if (islite) {
3581: n1ID = EGlite_indexBodyTopo(body, nobjs[0]);
3582: n2ID = EGlite_indexBodyTopo(body, nobjs[1]);
3583: } else {
3584: n1ID = EG_indexBodyTopo(body, nobjs[0]);
3585: n2ID = EG_indexBodyTopo(body, nobjs[1]);
3586: }
3588: if ((cpV1 == n1ID || cpV1 == n2ID) && (cpV2 == n1ID || cpV2 == n2ID)) {
3589: for (int kk = minID + 1; kk < maxID; ++kk) {
3590: cp_edge[kk] = eid;
3591: w_edge[kk] = eid;
3592: }
3593: }
3594: }
3595: }
3597: // Along vmax axis
3598: minID = wRowStart + (bpinfo[2] * (bpinfo[5] - 1));
3599: maxID = wRowStart + (bpinfo[2] * bpinfo[5] - 1);
3601: cpV1 = cp_vertex[minID];
3602: cpV2 = cp_vertex[maxID];
3603: for (int jj = 0; jj < Ne; ++jj) {
3604: ego edge = eobjs[jj];
3605: ego egeom, *nobjs;
3606: int eoclass, emtype, Nn, *nsenses;
3607: int n1ID, n2ID, eid;
3609: if (islite) {
3610: eid = EGlite_indexBodyTopo(body, edge);
3611: PetscCallEGADS(EGlite_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3612: } else {
3613: eid = EG_indexBodyTopo(body, edge);
3614: PetscCallEGADS(EG_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3615: }
3617: if (emtype != DEGENERATE) {
3618: // Get IDs for current Edge's End Vertices
3619: if (islite) {
3620: n1ID = EGlite_indexBodyTopo(body, nobjs[0]);
3621: n2ID = EGlite_indexBodyTopo(body, nobjs[1]);
3622: } else {
3623: n1ID = EG_indexBodyTopo(body, nobjs[0]);
3624: n2ID = EG_indexBodyTopo(body, nobjs[1]);
3625: }
3627: if ((cpV1 == n1ID || cpV1 == n2ID) && (cpV2 == n1ID || cpV2 == n2ID)) {
3628: for (int kk = minID + 1; kk < maxID - 1; ++kk) {
3629: cp_edge[kk] = eid;
3630: w_edge[kk] = eid;
3631: }
3632: }
3633: }
3634: }
3636: // Along umin axis
3637: minID = wRowStart;
3638: maxID = wRowStart + (bpinfo[2] * (bpinfo[5] - 1));
3640: cpV1 = cp_vertex[minID];
3641: cpV2 = cp_vertex[maxID];
3642: for (int jj = 0; jj < Ne; ++jj) {
3643: ego edge = eobjs[jj];
3644: ego egeom, *nobjs;
3645: int eoclass, emtype, Nn, *nsenses;
3646: int n1ID, n2ID, eid;
3648: if (islite) {
3649: eid = EGlite_indexBodyTopo(body, edge);
3650: PetscCallEGADS(EGlite_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3651: } else {
3652: eid = EG_indexBodyTopo(body, edge);
3653: PetscCallEGADS(EG_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3654: }
3656: if (emtype != DEGENERATE) {
3657: // Get IDs for current Edge's End Vertices
3658: if (islite) {
3659: n1ID = EGlite_indexBodyTopo(body, nobjs[0]);
3660: n2ID = EGlite_indexBodyTopo(body, nobjs[1]);
3661: } else {
3662: n1ID = EG_indexBodyTopo(body, nobjs[0]);
3663: n2ID = EG_indexBodyTopo(body, nobjs[1]);
3664: }
3666: if ((cpV1 == n1ID || cpV1 == n2ID) && (cpV2 == n1ID || cpV2 == n2ID)) {
3667: for (int kk = minID + bpinfo[2]; kk < maxID; kk += bpinfo[2]) {
3668: cp_edge[kk] = eid;
3669: w_edge[kk] = eid;
3670: }
3671: }
3672: }
3673: }
3675: // Along umax axis
3676: minID = wRowStart + (bpinfo[2] - 1);
3677: maxID = wRowStart + (bpinfo[2] * bpinfo[5]) - 1;
3678: cpV1 = cp_vertex[minID];
3679: cpV2 = cp_vertex[maxID];
3680: for (int jj = 0; jj < Ne; ++jj) {
3681: ego edge = eobjs[jj];
3682: ego egeom, *nobjs;
3683: int eoclass, emtype, Nn, *nsenses;
3684: int n1ID, n2ID, eid;
3686: if (islite) {
3687: eid = EGlite_indexBodyTopo(body, edge);
3688: PetscCallEGADS(EGlite_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3689: } else {
3690: eid = EG_indexBodyTopo(body, edge);
3691: PetscCallEGADS(EG_getTopology, (edge, &egeom, &eoclass, &emtype, NULL, &Nn, &nobjs, &nsenses));
3692: }
3694: if (emtype != DEGENERATE) {
3695: // Get IDs for current Edge's End Vertices
3696: if (islite) {
3697: n1ID = EGlite_indexBodyTopo(body, nobjs[0]);
3698: n2ID = EGlite_indexBodyTopo(body, nobjs[1]);
3699: } else {
3700: n1ID = EG_indexBodyTopo(body, nobjs[0]);
3701: n2ID = EG_indexBodyTopo(body, nobjs[1]);
3702: }
3704: if ((cpV1 == n1ID || cpV1 == n2ID) && (cpV2 == n1ID || cpV2 == n2ID)) {
3705: for (int kk = minID + bpinfo[2]; kk < maxID; kk += bpinfo[2]) {
3706: cp_edge[kk] = eid;
3707: w_edge[kk] = eid;
3708: }
3709: }
3710: }
3711: }
3712: // These two lines could be replaced with DMPlexFreeGeomObject()
3713: if (islite) EGlite_free(eobjs);
3714: else EG_free(eobjs);
3715: }
3717: // Determine Control Point Equivalence Matrix relating Control Points between Surfaces
3718: // Note: The Weights will also be tied together in the same manner
3719: // Also can use the Weight Hash Table for Row Start ID of each Face
3720: const PetscInt cpRowSize = totalNumCPs;
3721: const PetscInt cpColSize = cpRowSize;
3722: PetscInt *maxNumRelatePtr;
3723: PetscInt maxNumRelate = 0;
3725: // Create Point Surface Gradient Matrix
3726: PetscCall(MatCreate(PETSC_COMM_WORLD, &cpEquiv));
3727: PetscCall(MatSetSizes(cpEquiv, PETSC_DECIDE, PETSC_DECIDE, cpRowSize, cpColSize));
3728: PetscCall(MatSetType(cpEquiv, MATAIJ));
3729: PetscCall(MatSetUp(cpEquiv));
3731: for (int ii = 0; ii < totalNumCPs; ++ii) {
3732: PetscScalar x1, y1, z1;
3733: PetscInt maxRelateTemp = 0;
3735: x1 = cntrlPtCoords[(3 * ii) + 0];
3736: y1 = cntrlPtCoords[(3 * ii) + 1];
3737: z1 = cntrlPtCoords[(3 * ii) + 2];
3739: for (int jj = 0; jj < totalNumCPs; ++jj) {
3740: PetscScalar x2, y2, z2;
3741: PetscScalar cpDelta, eqFactor;
3742: x2 = cntrlPtCoords[(3 * jj) + 0];
3743: y2 = cntrlPtCoords[(3 * jj) + 1];
3744: z2 = cntrlPtCoords[(3 * jj) + 2];
3746: cpDelta = PetscSqrtReal(PetscSqr(x2 - x1) + PetscSqr(y2 - y1) + PetscSqr(z2 - z1));
3747: if (cpDelta < 1.0E-15) {
3748: eqFactor = 1.0;
3749: maxRelateTemp += 1;
3750: } else {
3751: eqFactor = 0.0;
3752: }
3754: // Store Results in PETSc Mat
3755: PetscCall(MatSetValue(cpEquiv, ii, jj, eqFactor, INSERT_VALUES));
3756: }
3757: if (maxRelateTemp > maxNumRelate) maxNumRelate = maxRelateTemp;
3758: }
3759: maxNumRelatePtr = &maxNumRelate;
3760: PetscCall(VecRestoreArrayWrite(cntrlPtCoordsVec, &cntrlPtCoords));
3762: // Assemble Point Surface Grad Matrix
3763: PetscCall(MatAssemblyBegin(cpEquiv, MAT_FINAL_ASSEMBLY));
3764: PetscCall(MatAssemblyEnd(cpEquiv, MAT_FINAL_ASSEMBLY));
3766: // Attach Control Point and Weight Data to DM
3767: {
3768: PetscContainer cpOrgObj, cpCoordLengthObj;
3769: PetscContainer wOrgObj, wDataLengthObj;
3770: PetscContainer cp_faceObj, cp_edgeObj, cp_vertexObj;
3771: PetscContainer w_faceObj, w_edgeObj, w_vertexObj;
3772: PetscContainer maxNumRelateObj;
3774: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Hash Table", (PetscObject *)&cpOrgObj));
3775: if (!cpOrgObj) {
3776: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cpOrgObj));
3777: PetscCall(PetscContainerSetPointer(cpOrgObj, faceCntrlPtRow_Start));
3778: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Hash Table", (PetscObject)cpOrgObj));
3779: PetscCall(PetscContainerDestroy(&cpOrgObj));
3780: } else {
3781: PetscCall(PetscContainerSetPointer(cpOrgObj, faceCntrlPtRow_Start));
3782: }
3784: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Coordinates", (PetscObject)cntrlPtCoordsVec));
3785: PetscCall(VecDestroy(&cntrlPtCoordsVec));
3787: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Coordinate Data Length", (PetscObject *)&cpCoordLengthObj));
3788: if (!cpCoordLengthObj) {
3789: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cpCoordLengthObj));
3790: PetscCall(PetscContainerSetPointer(cpCoordLengthObj, cpCoordDataLengthPtr));
3791: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Coordinate Data Length", (PetscObject)cpCoordLengthObj));
3792: PetscCall(PetscContainerDestroy(&cpCoordLengthObj));
3793: } else {
3794: PetscCall(PetscContainerSetPointer(cpCoordLengthObj, cpCoordDataLengthPtr));
3795: }
3797: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weights Hash Table", (PetscObject *)&wOrgObj));
3798: if (!wOrgObj) {
3799: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &wOrgObj));
3800: PetscCall(PetscContainerSetPointer(wOrgObj, faceCPWeightsRow_Start));
3801: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weights Hash Table", (PetscObject)wOrgObj));
3802: PetscCall(PetscContainerDestroy(&wOrgObj));
3803: } else {
3804: PetscCall(PetscContainerSetPointer(wOrgObj, faceCPWeightsRow_Start));
3805: }
3807: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weight Data", (PetscObject)cntrlPtWeightsVec));
3808: PetscCall(VecDestroy(&cntrlPtWeightsVec));
3810: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data Length", (PetscObject *)&wDataLengthObj));
3811: if (!wDataLengthObj) {
3812: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &wDataLengthObj));
3813: PetscCall(PetscContainerSetPointer(wDataLengthObj, wDataLengthPtr));
3814: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weight Data Length", (PetscObject)wDataLengthObj));
3815: PetscCall(PetscContainerDestroy(&wDataLengthObj));
3816: } else {
3817: PetscCall(PetscContainerSetPointer(wDataLengthObj, wDataLengthPtr));
3818: }
3820: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Equivalency Matrix", (PetscObject)cpEquiv));
3822: PetscCall(PetscObjectQuery((PetscObject)dm, "Maximum Number Control Point Equivalency", (PetscObject *)&maxNumRelateObj));
3823: if (!maxNumRelateObj) {
3824: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &maxNumRelateObj));
3825: PetscCall(PetscContainerSetPointer(maxNumRelateObj, maxNumRelatePtr));
3826: PetscCall(PetscObjectCompose((PetscObject)dm, "Maximum Number Control Point Equivalency", (PetscObject)maxNumRelateObj));
3827: PetscCall(PetscContainerDestroy(&maxNumRelateObj));
3828: } else {
3829: PetscCall(PetscContainerSetPointer(maxNumRelateObj, maxNumRelatePtr));
3830: }
3832: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point - Face Map", (PetscObject *)&cp_faceObj));
3833: if (!cp_faceObj) {
3834: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cp_faceObj));
3835: PetscCall(PetscContainerSetPointer(cp_faceObj, cp_face));
3836: PetscCall(PetscContainerSetCtxDestroy(cp_faceObj, PetscCtxDestroyDefault));
3837: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point - Face Map", (PetscObject)cp_faceObj));
3838: PetscCall(PetscContainerDestroy(&cp_faceObj));
3839: } else {
3840: void *tmp;
3842: PetscCall(PetscContainerGetPointer(cp_faceObj, &tmp));
3843: PetscCall(PetscFree(tmp));
3844: PetscCall(PetscContainerSetPointer(cp_faceObj, cp_face));
3845: }
3847: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight - Face Map", (PetscObject *)&w_faceObj));
3848: if (!w_faceObj) {
3849: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &w_faceObj));
3850: PetscCall(PetscContainerSetPointer(w_faceObj, w_face));
3851: PetscCall(PetscContainerSetCtxDestroy(w_faceObj, PetscCtxDestroyDefault));
3852: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weight - Face Map", (PetscObject)w_faceObj));
3853: PetscCall(PetscContainerDestroy(&w_faceObj));
3854: } else {
3855: void *tmp;
3857: PetscCall(PetscContainerGetPointer(w_faceObj, &tmp));
3858: PetscCall(PetscFree(tmp));
3859: PetscCall(PetscContainerSetPointer(w_faceObj, w_face));
3860: }
3862: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point - Edge Map", (PetscObject *)&cp_edgeObj));
3863: if (!cp_edgeObj) {
3864: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cp_edgeObj));
3865: PetscCall(PetscContainerSetPointer(cp_edgeObj, cp_edge));
3866: PetscCall(PetscContainerSetCtxDestroy(cp_edgeObj, PetscCtxDestroyDefault));
3867: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point - Edge Map", (PetscObject)cp_edgeObj));
3868: PetscCall(PetscContainerDestroy(&cp_edgeObj));
3869: } else {
3870: void *tmp;
3872: PetscCall(PetscContainerGetPointer(cp_edgeObj, &tmp));
3873: PetscCall(PetscFree(tmp));
3874: PetscCall(PetscContainerSetPointer(cp_edgeObj, cp_edge));
3875: }
3877: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight - Edge Map", (PetscObject *)&w_edgeObj));
3878: if (!w_edgeObj) {
3879: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &w_edgeObj));
3880: PetscCall(PetscContainerSetPointer(w_edgeObj, w_edge));
3881: PetscCall(PetscContainerSetCtxDestroy(w_edgeObj, PetscCtxDestroyDefault));
3882: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weight - Edge Map", (PetscObject)w_edgeObj));
3883: PetscCall(PetscContainerDestroy(&w_edgeObj));
3884: } else {
3885: void *tmp;
3887: PetscCall(PetscContainerGetPointer(w_edgeObj, &tmp));
3888: PetscCall(PetscFree(tmp));
3889: PetscCall(PetscContainerSetPointer(w_edgeObj, w_edge));
3890: }
3892: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point - Vertex Map", (PetscObject *)&cp_vertexObj));
3893: if (!cp_vertexObj) {
3894: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &cp_vertexObj));
3895: PetscCall(PetscContainerSetPointer(cp_vertexObj, cp_vertex));
3896: PetscCall(PetscContainerSetCtxDestroy(cp_vertexObj, PetscCtxDestroyDefault));
3897: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point - Vertex Map", (PetscObject)cp_vertexObj));
3898: PetscCall(PetscContainerDestroy(&cp_vertexObj));
3899: } else {
3900: void *tmp;
3902: PetscCall(PetscContainerGetPointer(cp_vertexObj, &tmp));
3903: PetscCall(PetscFree(tmp));
3904: PetscCall(PetscContainerSetPointer(cp_vertexObj, cp_vertex));
3905: }
3907: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight - Vertex Map", (PetscObject *)&w_vertexObj));
3908: if (!w_vertexObj) {
3909: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &w_vertexObj));
3910: PetscCall(PetscContainerSetPointer(w_vertexObj, w_vertex));
3911: PetscCall(PetscContainerSetCtxDestroy(w_vertexObj, PetscCtxDestroyDefault));
3912: PetscCall(PetscObjectCompose((PetscObject)dm, "Control Point Weight - Vertex Map", (PetscObject)w_vertexObj));
3913: PetscCall(PetscContainerDestroy(&w_vertexObj));
3914: } else {
3915: void *tmp;
3917: PetscCall(PetscContainerGetPointer(w_vertexObj, &tmp));
3918: PetscCall(PetscFree(tmp));
3919: PetscCall(PetscContainerSetPointer(w_vertexObj, w_vertex));
3920: }
3921: }
3923: // Define Matrix to store Geometry Gradient information dGeom_i/dCPj_i
3924: PetscInt gcntr = 0;
3925: const PetscInt rowSize = 3 * maxNumCPs * totalNumPoints;
3926: const PetscInt colSize = 4 * Nf;
3928: // Create Point Surface Gradient Matrix
3929: PetscCall(MatCreate(PETSC_COMM_WORLD, &pointSurfGrad));
3930: PetscCall(MatSetSizes(pointSurfGrad, PETSC_DECIDE, PETSC_DECIDE, rowSize, colSize));
3931: PetscCall(MatSetType(pointSurfGrad, MATAIJ));
3932: PetscCall(MatSetUp(pointSurfGrad));
3934: // Create Hash Table to store Point's stare row in surfaceGrad[][]
3935: PetscCall(PetscHMapICreate(&pointSurfGradRow_Start));
3937: // Get Coordinates for the DMPlex point
3938: DM cdm;
3939: PetscInt dE, Nv;
3940: Vec coordinatesLocal;
3941: PetscScalar *coords = NULL;
3943: PetscCall(DMGetCoordinateDM(dm, &cdm));
3944: PetscCall(DMGetCoordinateDim(dm, &dE));
3945: PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal));
3947: // CYCLE THROUGH FACEs
3948: PetscScalar maxGrad = 0.;
3949: PetscCall(VecGetArrayWrite(gradSACPVec, &gradSACP));
3950: PetscCall(VecGetArrayWrite(gradSAWVec, &gradSAW));
3951: PetscCall(VecGetArrayWrite(gradVCPVec, &gradVCP));
3952: PetscCall(VecGetArrayWrite(gradVWVec, &gradVW));
3953: for (int f = 0; f < Nf; ++f) {
3954: ego face = fobjs[f];
3955: ego *eobjs, *nobjs;
3956: PetscInt fid, Ne, Nn;
3957: DMLabel faceLabel, edgeLabel, nodeLabel;
3958: PetscHMapI currFaceUniquePoints = NULL;
3959: IS facePoints, edgePoints, nodePoints;
3960: const PetscInt *fIndices, *eIndices, *nIndices;
3961: PetscInt fSize, eSize, nSize;
3962: PetscHashIter fHashKeyIter, eHashKeyIter, nHashKeyIter, pHashKeyIter;
3963: PetscBool fHashKeyFound, eHashKeyFound, nHashKeyFound, pHashKeyFound;
3964: PetscInt cfCntr = 0;
3966: // Get Geometry Object for the Current FACE
3967: if (islite) {
3968: PetscCall(EGlite_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3969: PetscCall(EGlite_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3970: } else {
3971: PetscCall(EG_getTopology(face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
3972: PetscCall(EG_getGeometry(fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
3973: }
3975: // Get all EDGE and NODE objects attached to the current FACE
3976: if (islite) {
3977: PetscCall(EGlite_getBodyTopos(body, face, EDGE, &Ne, &eobjs));
3978: PetscCall(EGlite_getBodyTopos(body, face, NODE, &Nn, &nobjs));
3979: } else {
3980: PetscCall(EG_getBodyTopos(body, face, EDGE, &Ne, &eobjs));
3981: PetscCall(EG_getBodyTopos(body, face, NODE, &Nn, &nobjs));
3982: }
3984: // Get all DMPlex Points that have DMLabel "EGADS Face ID" and store them in a Hash Table for later use
3985: if (islite) {
3986: fid = EGlite_indexBodyTopo(body, face);
3987: } else {
3988: fid = EG_indexBodyTopo(body, face);
3989: }
3991: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
3992: PetscCall(DMLabelGetStratumIS(faceLabel, fid, &facePoints));
3993: PetscCall(ISGetIndices(facePoints, &fIndices));
3994: PetscCall(ISGetSize(facePoints, &fSize));
3996: PetscCall(PetscHMapICreate(&currFaceUniquePoints));
3998: for (int jj = 0; jj < fSize; ++jj) {
3999: PetscCall(PetscHMapIFind(currFaceUniquePoints, fIndices[jj], &fHashKeyIter, &fHashKeyFound));
4001: if (!fHashKeyFound) {
4002: PetscCall(PetscHMapISet(currFaceUniquePoints, fIndices[jj], cfCntr));
4003: cfCntr += 1;
4004: }
4006: PetscCall(PetscHMapIFind(pointSurfGradRow_Start, fIndices[jj], &pHashKeyIter, &pHashKeyFound));
4008: if (!pHashKeyFound) {
4009: PetscCall(PetscHMapISet(pointSurfGradRow_Start, fIndices[jj], gcntr));
4010: gcntr += 3 * maxNumCPs;
4011: }
4012: }
4013: PetscCall(ISRestoreIndices(facePoints, &fIndices));
4014: PetscCall(ISDestroy(&facePoints));
4016: // Get all DMPlex Points that have DMLable "EGADS Edge ID" attached to the current FACE and store them in a Hash Table for later use.
4017: for (int jj = 0; jj < Ne; ++jj) {
4018: ego edge = eobjs[jj];
4019: PetscBool containLabelValue;
4021: if (islite) {
4022: id = EGlite_indexBodyTopo(body, edge);
4023: } else {
4024: id = EG_indexBodyTopo(body, edge);
4025: }
4027: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
4028: PetscCall(DMLabelHasValue(edgeLabel, id, &containLabelValue));
4030: if (containLabelValue) {
4031: PetscCall(DMLabelGetStratumIS(edgeLabel, id, &edgePoints));
4032: PetscCall(ISGetIndices(edgePoints, &eIndices));
4033: PetscCall(ISGetSize(edgePoints, &eSize));
4035: for (int kk = 0; kk < eSize; ++kk) {
4036: PetscCall(PetscHMapIFind(currFaceUniquePoints, eIndices[kk], &eHashKeyIter, &eHashKeyFound));
4038: if (!eHashKeyFound) {
4039: PetscCall(PetscHMapISet(currFaceUniquePoints, eIndices[kk], cfCntr));
4040: cfCntr += 1;
4041: }
4043: PetscCall(PetscHMapIFind(pointSurfGradRow_Start, eIndices[kk], &pHashKeyIter, &pHashKeyFound));
4045: if (!pHashKeyFound) {
4046: PetscCall(PetscHMapISet(pointSurfGradRow_Start, eIndices[kk], gcntr));
4047: gcntr += 3 * maxNumCPs;
4048: }
4049: }
4050: PetscCall(ISRestoreIndices(edgePoints, &eIndices));
4051: PetscCall(ISDestroy(&edgePoints));
4052: }
4053: }
4055: // Get all DMPlex Points that have DMLabel "EGADS Vertex ID" attached to the current FACE and store them in a Hash Table for later use.
4056: for (int jj = 0; jj < Nn; ++jj) {
4057: ego node = nobjs[jj];
4059: if (islite) {
4060: id = EGlite_indexBodyTopo(body, node);
4061: } else {
4062: id = EG_indexBodyTopo(body, node);
4063: }
4065: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &nodeLabel));
4066: PetscCall(DMLabelGetStratumIS(nodeLabel, id, &nodePoints));
4067: PetscCall(ISGetIndices(nodePoints, &nIndices));
4068: PetscCall(ISGetSize(nodePoints, &nSize));
4070: for (int kk = 0; kk < nSize; ++kk) {
4071: PetscCall(PetscHMapIFind(currFaceUniquePoints, nIndices[kk], &nHashKeyIter, &nHashKeyFound));
4073: if (!nHashKeyFound) {
4074: PetscCall(PetscHMapISet(currFaceUniquePoints, nIndices[kk], cfCntr));
4075: cfCntr += 1;
4076: }
4078: PetscCall(PetscHMapIFind(pointSurfGradRow_Start, nIndices[kk], &pHashKeyIter, &pHashKeyFound));
4079: if (!pHashKeyFound) {
4080: PetscCall(PetscHMapISet(pointSurfGradRow_Start, nIndices[kk], gcntr));
4081: gcntr += 3 * maxNumCPs;
4082: }
4083: }
4084: PetscCall(ISRestoreIndices(nodePoints, &nIndices));
4085: PetscCall(ISDestroy(&nodePoints));
4086: }
4088: // Get the Total Number of entries in the Hash Table
4089: PetscInt currFaceUPSize;
4090: PetscCall(PetscHMapIGetSize(currFaceUniquePoints, &currFaceUPSize));
4092: // Get Keys
4093: PetscInt currFaceUPKeys[currFaceUPSize], off = 0;
4094: PetscCall(PetscHMapIGetKeys(currFaceUniquePoints, &off, currFaceUPKeys));
4095: PetscCall(PetscHMapIDestroy(&currFaceUniquePoints));
4097: // Get Current Face Surface Area
4098: PetscScalar fSA, faceData[14];
4099: PetscCall(EG_getMassProperties(face, faceData)); // This doesn't have a EGlite version. Will it work for EGADSlite files?? KNOWN_ISSUE
4100: fSA = faceData[1];
4102: // Get Start Row in cpEquiv Matrix
4103: PetscHashIter Witer;
4104: PetscBool Wfound;
4105: PetscInt faceWStartRow;
4106: PetscCall(PetscHMapIFind(faceCPWeightsRow_Start, fid, &Witer, &Wfound));
4107: PetscCheck(Wfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "FACE ID not found in Control Point Weights Hash Table");
4108: PetscCall(PetscHMapIGet(faceCPWeightsRow_Start, fid, &faceWStartRow));
4110: // Cycle through all points on the current FACE
4111: for (int jj = 0; jj < currFaceUPSize; ++jj) {
4112: PetscInt currPointID = currFaceUPKeys[jj];
4113: PetscCall(DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, currPointID, &Nv, &coords));
4115: // Get UV position of FACE
4116: double params[2], range[4], eval[18];
4117: int peri;
4119: if (islite) PetscCall(EGlite_getRange(face, range, &peri));
4120: else PetscCall(EG_getRange(face, range, &peri));
4122: PetscCall(DMPlex_Geom_FACE_XYZtoUV_Internal(coords, face, range, 0, dE, params, islite));
4124: if (islite) PetscCall(EGlite_evaluate(face, params, eval));
4125: else PetscCall(EG_evaluate(face, params, eval));
4127: // Make a new SURFACE Geometry by changing the location of the Control Points
4128: int prvSize = bpinfo[3] + bpinfo[6] + (4 * bpinfo[2] * bpinfo[5]);
4129: double nbprv[prvSize];
4131: // Cycle through each Control Point
4132: double denomNew, denomOld;
4133: double deltaCoord = 1.0E-4;
4134: int offset = bpinfo[3] + bpinfo[6];
4135: int wOffset = offset + (3 * bpinfo[2] * bpinfo[5]);
4136: for (int ii = 0; ii < bpinfo[2] * bpinfo[5]; ++ii) {
4137: PetscCheck(face->blind, PETSC_COMM_SELF, PETSC_ERR_LIB, "Face %d is corrupted: %d %d", f, jj, ii);
4138: #if 0
4139: // Cycle through each direction (x, then y, then z)
4140: if (jj == 0) {
4141: // Get the Number Control Points that are the same as the current points
4142: // We are looking for repeated Control Points
4143: PetscInt commonCPcntr = 0;
4144: for (int mm = 0; mm < bpinfo[2]*bpinfo[5]; ++mm) {
4145: PetscScalar matValue;
4146: PetscCall(MatGetValue(cpEquiv, faceWStartRow + ii, faceWStartRow + mm, &matValue));
4148: if (matValue > 0.0) commonCPcntr += 1;
4149: }
4150: }
4151: #endif
4153: for (int kk = 0; kk < 4; ++kk) {
4154: // Reinitialize nbprv[] values because we only want to change one value at a time
4155: for (int mm = 0; mm < prvSize; ++mm) nbprv[mm] = bprv[mm];
4156: PetscCheck(face->blind, PETSC_COMM_SELF, PETSC_ERR_LIB, "Face %d is corrupted: %d %d %d", f, jj, ii, kk);
4158: if (kk == 0) { //X
4159: nbprv[offset + 0] = bprv[offset + 0] + deltaCoord;
4160: nbprv[offset + 1] = bprv[offset + 1];
4161: nbprv[offset + 2] = bprv[offset + 2];
4162: denomNew = nbprv[offset + 0];
4163: denomOld = bprv[offset + 0];
4164: } else if (kk == 1) { //Y
4165: nbprv[offset + 0] = bprv[offset + 0];
4166: nbprv[offset + 1] = bprv[offset + 1] + deltaCoord;
4167: nbprv[offset + 2] = bprv[offset + 2];
4168: denomNew = nbprv[offset + 1];
4169: denomOld = bprv[offset + 1];
4170: } else if (kk == 2) { //Z
4171: nbprv[offset + 0] = bprv[offset + 0];
4172: nbprv[offset + 1] = bprv[offset + 1];
4173: nbprv[offset + 2] = bprv[offset + 2] + deltaCoord;
4174: denomNew = nbprv[offset + 2];
4175: denomOld = bprv[offset + 2];
4176: } else if (kk == 3) { // Weights
4177: nbprv[wOffset + ii] = bprv[wOffset + ii] + deltaCoord;
4178: denomNew = nbprv[wOffset + ii];
4179: denomOld = bprv[wOffset + ii];
4180: } else {
4181: // currently do nothing
4182: }
4184: // Create New Surface Based on New Control Points or Weights
4185: ego newgeom, context;
4186: PetscCallEGADS(EG_getContext, (face, &context)); // This does not have an EGlite_ version KNOWN_ISSUE
4187: PetscCallEGADS(EG_makeGeometry, (context, SURFACE, BSPLINE, NULL, bpinfo, nbprv, &newgeom)); // This does not have an EGlite_ version KNOWN_ISSUE
4188: PetscCheck(face->blind, PETSC_COMM_SELF, PETSC_ERR_LIB, "Face %d is corrupted: %d %d %d", f, jj, ii, kk);
4190: // Evaluate new (x, y, z) Point Position based on new Surface Definition
4191: double newCoords[18];
4192: if (islite) PetscCall(EGlite_getRange(newgeom, range, &peri));
4193: else PetscCall(EG_getRange(newgeom, range, &peri));
4195: PetscCall(DMPlex_Geom_FACE_XYZtoUV_Internal(coords, face, range, 0, dE, params, islite));
4196: PetscCheck(face->blind, PETSC_COMM_SELF, PETSC_ERR_LIB, "Face %d is corrupted: %d %d %d", f, jj, ii, kk);
4198: if (islite) PetscCall(EGlite_evaluate(newgeom, params, newCoords));
4199: else PetscCall(EG_evaluate(newgeom, params, newCoords));
4201: // Calculate Surface Area Gradients wrt Control Points and Weights using the local discrete FACE only
4202: // NOTE 1: Will not provide Volume Gradient wrt to Control Points and Weights.
4203: // NOTE 2: This is faster than below where an entire new solid geometry is created for each
4204: // Control Point and Weight gradient
4205: if (!fullGeomGrad) {
4206: // Create new FACE based on new SURFACE geometry
4207: if (jj == 0) { // only for 1st DMPlex Point because we only per CP or Weight
4208: double newFaceRange[4];
4209: int newFacePeri;
4210: if (islite) PetscCall(EGlite_getRange(newgeom, newFaceRange, &newFacePeri));
4211: else PetscCall(EG_getRange(newgeom, newFaceRange, &newFacePeri));
4213: ego newface;
4214: PetscCallEGADS(EG_makeFace, (newgeom, SFORWARD, newFaceRange, &newface)); // Does not have EGlite version KNOWN_ISSUE
4215: PetscCheck(face->blind, PETSC_COMM_SELF, PETSC_ERR_LIB, "Face %d is corrupted: %d %d %d", f, jj, ii, kk);
4217: // Get New Face Surface Area
4218: PetscScalar newfSA, newFaceData[14];
4219: PetscCall(EG_getMassProperties(newface, newFaceData)); // Does not have EGlite version KNOWN_ISSUE
4220: newfSA = newFaceData[1];
4221: PetscCallEGADS(EG_deleteObject, (newface));
4222: PetscCheck(face->blind, PETSC_COMM_SELF, PETSC_ERR_LIB, "Face %d is corrupted: %d %d %d", f, jj, ii, kk);
4224: // Update Control Points
4225: PetscHashIter CPiter, Witer;
4226: PetscBool CPfound, Wfound;
4227: PetscInt faceCPStartRow, faceWStartRow;
4229: PetscScalar dSAdCPi;
4230: dSAdCPi = (newfSA - fSA) / (denomNew - denomOld);
4232: if (kk < 3) {
4233: PetscCall(PetscHMapIFind(faceCntrlPtRow_Start, fid, &CPiter, &CPfound));
4234: PetscCheck(CPfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "FACE ID not found in Control Point Hash Table");
4235: PetscCall(PetscHMapIGet(faceCntrlPtRow_Start, fid, &faceCPStartRow));
4237: gradSACP[faceCPStartRow + (ii * 3) + kk] = dSAdCPi;
4239: if (PetscAbsReal(dSAdCPi) > maxGrad) maxGrad = PetscAbsReal(dSAdCPi);
4241: } else if (kk == 3) {
4242: PetscCall(PetscHMapIFind(faceCPWeightsRow_Start, fid, &Witer, &Wfound));
4243: PetscCheck(Wfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "FACE ID not found in Control Point Hash Table");
4244: PetscCall(PetscHMapIGet(faceCPWeightsRow_Start, fid, &faceWStartRow));
4246: gradSAW[faceWStartRow + ii] = dSAdCPi;
4248: } else {
4249: // Do Nothing
4250: }
4251: }
4252: }
4253: PetscCallEGADS(EG_deleteObject, (newgeom));
4255: // Now Calculate the Surface Gradient for the change in x-component Control Point
4256: PetscScalar dxdCx = (newCoords[0] - coords[0]) / deltaCoord;
4257: PetscScalar dxdCy = (newCoords[1] - coords[1]) / deltaCoord;
4258: PetscScalar dxdCz = (newCoords[2] - coords[2]) / deltaCoord;
4260: // Store Gradient Information in surfaceGrad[][] Matrix
4261: PetscInt startRow;
4262: PetscCall(PetscHMapIGet(pointSurfGradRow_Start, currPointID, &startRow));
4264: // Store Results in PETSc Mat
4265: PetscCall(MatSetValue(pointSurfGrad, startRow + (ii * 3) + 0, ((fid - 1) * 4) + kk, dxdCx, INSERT_VALUES));
4266: PetscCall(MatSetValue(pointSurfGrad, startRow + (ii * 3) + 1, ((fid - 1) * 4) + kk, dxdCy, INSERT_VALUES));
4267: PetscCall(MatSetValue(pointSurfGrad, startRow + (ii * 3) + 2, ((fid - 1) * 4) + kk, dxdCz, INSERT_VALUES));
4269: //PetscCallEGADS(EG_deleteObject, (newgeom));
4270: PetscCheck(face->blind, PETSC_COMM_SELF, PETSC_ERR_LIB, "Face is corrupted");
4271: }
4272: offset += 3;
4273: }
4274: PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, currPointID, &Nv, &coords));
4275: }
4276: }
4278: // Assemble Point Surface Grad Matrix
4279: PetscCall(MatAssemblyBegin(pointSurfGrad, MAT_FINAL_ASSEMBLY));
4280: PetscCall(MatAssemblyEnd(pointSurfGrad, MAT_FINAL_ASSEMBLY));
4282: if (fullGeomGrad) {
4283: // Calculate Surface Area and Volume Control Point and Control Point Weight Gradients
4284: // Note: This is much slower than above due to a new solid geometry being created for
4285: // each change in Control Point and Control Point Weight. However, this method
4286: // will provide the Volume Gradient.
4288: // Get Current Face Surface Area
4289: PetscScalar bodyVol, bodySA, bodyData[14];
4290: PetscCall(EG_getMassProperties(body, bodyData)); // Does not have an EGlite version KNOWN_ISSUE
4291: bodyVol = bodyData[0];
4292: bodySA = bodyData[1];
4294: // Cycle through Control Points
4295: for (int ii = 0; ii < totalNumCPs; ++ii) { // ii should also be the row in cpEquiv for the Control Point
4296: // Cycle through X, Y, Z, W changes
4297: for (int jj = 0; jj < 4; ++jj) {
4298: // Cycle Through Faces
4299: double denomNew = 0.0, denomOld = 0.0;
4300: double deltaCoord = 1.0E-4;
4301: ego newGeom[Nf];
4302: ego newFaces[Nf];
4303: for (int kk = 0; kk < Nf; ++kk) {
4304: ego face;
4305: PetscInt currFID = kk + 1;
4307: if (islite) {
4308: // Get Current FACE
4309: PetscCallEGADS(EGlite_objectBodyTopo, (body, FACE, currFID, &face));
4311: // Get Geometry Object for the Current FACE
4312: PetscCallEGADS(EGlite_getTopology, (face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
4313: PetscCallEGADS(EGlite_getGeometry, (fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
4314: } else {
4315: // Get Current FACE
4316: PetscCallEGADS(EG_objectBodyTopo, (body, FACE, currFID, &face));
4318: // Get Geometry Object for the Current FACE
4319: PetscCallEGADS(EG_getTopology, (face, &fgeom, &foclass, &fmtype, fdata, &Nl, &lobjs, &lsenses));
4320: PetscCallEGADS(EG_getGeometry, (fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
4321: }
4323: // Make a new SURFACE Geometry by changing the location of the Control Points
4324: int prvSize = bpinfo[3] + bpinfo[6] + (4 * bpinfo[2] * bpinfo[5]);
4325: double nbprv[prvSize];
4327: // Reinitialize nbprv[] values because we only want to change one value at a time
4328: for (int mm = 0; mm < prvSize; ++mm) nbprv[mm] = bprv[mm];
4330: // Get Control Point Row and Column Start for cpEquiv
4331: PetscHashIter Witer;
4332: PetscBool Wfound;
4333: PetscInt faceWStartRow;
4334: PetscCall(PetscHMapIFind(faceCPWeightsRow_Start, currFID, &Witer, &Wfound));
4335: PetscCheck(Wfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "FACE ID not found in Control Point Weights Hash Table");
4336: PetscCall(PetscHMapIGet(faceCPWeightsRow_Start, currFID, &faceWStartRow));
4338: // Modify the Current Control Point on this FACE and All Other FACES
4339: // IMPORTANT!!! If you do not move all identical Control Points on other FACES
4340: // you will not generate a solid body. You will generate a set of
4341: // disconnected surfaces that have gap(s) between them.
4342: int offset = bpinfo[3] + bpinfo[6];
4343: int wOffset = offset + (3 * bpinfo[2] * bpinfo[5]);
4344: for (int mm = 0; mm < bpinfo[2] * bpinfo[5]; ++mm) {
4345: PetscScalar matValue;
4346: PetscCall(MatGetValue(cpEquiv, ii, faceWStartRow + mm, &matValue));
4348: if (matValue > 0.0) {
4349: if (jj == 0) { //X
4350: nbprv[offset + (3 * mm) + 0] = bprv[offset + (3 * mm) + 0] + deltaCoord;
4351: nbprv[offset + (3 * mm) + 1] = bprv[offset + (3 * mm) + 1];
4352: nbprv[offset + (3 * mm) + 2] = bprv[offset + (3 * mm) + 2];
4353: denomNew = nbprv[offset + (3 * mm) + 0];
4354: denomOld = bprv[offset + (3 * mm) + 0];
4355: } else if (jj == 1) { //Y
4356: nbprv[offset + (3 * mm) + 0] = bprv[offset + (3 * mm) + 0];
4357: nbprv[offset + (3 * mm) + 1] = bprv[offset + (3 * mm) + 1] + deltaCoord;
4358: nbprv[offset + (3 * mm) + 2] = bprv[offset + (3 * mm) + 2];
4359: denomNew = nbprv[offset + (3 * mm) + 1];
4360: denomOld = bprv[offset + (3 * mm) + 1];
4361: } else if (jj == 2) { //Z
4362: nbprv[offset + (3 * mm) + 0] = bprv[offset + (3 * mm) + 0];
4363: nbprv[offset + (3 * mm) + 1] = bprv[offset + (3 * mm) + 1];
4364: nbprv[offset + (3 * mm) + 2] = bprv[offset + (3 * mm) + 2] + deltaCoord;
4365: denomNew = nbprv[offset + (3 * mm) + 2];
4366: denomOld = bprv[offset + (3 * mm) + 2];
4367: } else if (jj == 3) { // Weights
4368: nbprv[wOffset + mm] = bprv[wOffset + mm] + deltaCoord;
4369: denomNew = nbprv[wOffset + mm];
4370: denomOld = bprv[wOffset + mm];
4371: } else {
4372: // currently do nothing
4373: }
4374: }
4375: }
4377: // Create New Surface Based on New Control Points or Weights
4378: ego newgeom, context;
4379: PetscCallEGADS(EG_getContext, (face, &context)); // Does not have an EGlite_ versions KNOWN_ISSUE
4380: PetscCallEGADS(EG_makeGeometry, (context, SURFACE, BSPLINE, NULL, bpinfo, nbprv, &newgeom)); // Does not have an EGlite_ version KNOWN_ISSUE
4382: // Create New FACE based on modified geometry
4383: double newFaceRange[4];
4384: int newFacePeri;
4385: if (islite) PetscCallEGADS(EGlite_getRange, (newgeom, newFaceRange, &newFacePeri));
4386: else PetscCallEGADS(EG_getRange, (newgeom, newFaceRange, &newFacePeri));
4388: ego newface;
4389: PetscCallEGADS(EG_makeFace, (newgeom, SFORWARD, newFaceRange, &newface)); // Does not have an EGlite_ version KNOWN_ISSUE
4391: // store new face for later assembly
4392: newGeom[kk] = newgeom;
4393: newFaces[kk] = newface;
4394: }
4396: // X-WANT TO BUILD THE NEW GEOMETRY, X-GET NEW SA AND PERFORM dSA/dCPi CALCS HERE <---
4397: // Sew New Faces together to get a new model
4398: ego newmodel;
4399: PetscCall(EG_sewFaces(Nf, newFaces, 0.0, 0, &newmodel)); // Does not have an EGlite_ version KNOWN_ISSUE
4401: // Get Surface Area and Volume of New/Updated Solid Body
4402: PetscScalar newData[14];
4403: if (islite) PetscCallEGADS(EGlite_getTopology, (newmodel, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4404: else PetscCallEGADS(EG_getTopology, (newmodel, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4406: ego nbody = bodies[0];
4407: PetscCall(EG_getMassProperties(nbody, newData)); // Does not have an EGlite_ version KNOWN_ISSUE
4409: PetscScalar dSAdCPi, dVdCPi;
4410: PetscScalar nbodyVol = newData[0], nbodySA = newData[1];
4412: // Calculate Gradients wrt to Control Points and Control Points Weights depending on jj value
4413: dSAdCPi = (nbodySA - bodySA) / (denomNew - denomOld);
4414: dVdCPi = (nbodyVol - bodyVol) / (denomNew - denomOld);
4416: if (jj < 3) {
4417: // Gradienst wrt to Control Points
4418: gradSACP[(ii * 3) + jj] = dSAdCPi;
4419: gradVCP[(ii * 3) + jj] = dVdCPi;
4420: } else if (jj == 3) {
4421: // Gradients wrt to Control Point Weights
4422: gradSAW[ii] = dSAdCPi;
4423: gradVW[ii] = dVdCPi;
4424: } else {
4425: // Do Nothing
4426: }
4427: PetscCallEGADS(EG_deleteObject, (newmodel));
4428: for (int kk = 0; kk < Nf; ++kk) {
4429: PetscCallEGADS(EG_deleteObject, (newFaces[kk]));
4430: PetscCallEGADS(EG_deleteObject, (newGeom[kk]));
4431: }
4432: }
4433: }
4434: }
4435: PetscCall(VecRestoreArrayWrite(gradSACPVec, &gradSACP));
4436: PetscCall(VecRestoreArrayWrite(gradSAWVec, &gradSAW));
4437: PetscCall(VecRestoreArrayWrite(gradVCPVec, &gradVCP));
4438: PetscCall(VecRestoreArrayWrite(gradVWVec, &gradVW));
4439: PetscCall(MatDestroy(&cpEquiv));
4441: // Attach Surface Gradient Hash Table and Matrix to DM
4442: {
4443: PetscContainer surfGradOrgObj;
4445: PetscCall(PetscObjectQuery((PetscObject)dm, "Surface Gradient Hash Table", (PetscObject *)&surfGradOrgObj));
4446: if (!surfGradOrgObj) {
4447: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &surfGradOrgObj));
4448: PetscCall(PetscContainerSetPointer(surfGradOrgObj, pointSurfGradRow_Start));
4449: PetscCall(PetscContainerSetCtxDestroy(surfGradOrgObj, DestroyHashMap));
4450: PetscCall(PetscObjectCompose((PetscObject)dm, "Surface Gradient Hash Table", (PetscObject)surfGradOrgObj));
4451: PetscCall(PetscContainerDestroy(&surfGradOrgObj));
4452: } else {
4453: PetscCall(PetscContainerSetPointer(surfGradOrgObj, pointSurfGradRow_Start));
4454: }
4456: PetscCall(PetscObjectCompose((PetscObject)dm, "Surface Gradient Matrix", (PetscObject)pointSurfGrad));
4457: PetscCall(MatDestroy(&pointSurfGrad));
4459: PetscCall(PetscObjectCompose((PetscObject)dm, "Surface Area Control Point Gradient", (PetscObject)gradSACPVec));
4460: PetscCall(VecDestroy(&gradSACPVec));
4462: PetscCall(PetscObjectCompose((PetscObject)dm, "Surface Area Weights Gradient", (PetscObject)gradSAWVec));
4463: PetscCall(VecDestroy(&gradSAWVec));
4465: if (fullGeomGrad) {
4466: PetscCall(PetscObjectCompose((PetscObject)dm, "Volume Control Point Gradient", (PetscObject)gradVCPVec));
4467: PetscCall(PetscObjectCompose((PetscObject)dm, "Volume Weights Gradient", (PetscObject)gradVWVec));
4468: }
4469: PetscCall(VecDestroy(&gradVCPVec));
4470: PetscCall(VecDestroy(&gradVWVec));
4471: }
4473: // Could be replaced with DMPlexFreeGeomObject()
4474: if (islite) EGlite_free(fobjs);
4475: else EG_free(fobjs);
4476: PetscFunctionReturn(PETSC_SUCCESS);
4477: #else
4478: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This method requires EGADS support. Reconfigure using --download-egads");
4479: #endif
4480: }
4482: /*@C
4483: DMPlexModifyGeomModel - Generates a new EGADS geometry model based in user provided Control Points and Control Points Weights. Optionally, the function will inflate the DM to the new geometry and save the new geometry to a file.
4485: Collective
4487: Input Parameters:
4488: + dm - The DM object representing the mesh with PetscContainer containing an EGADS geometry model
4489: . comm - MPI_Comm object
4490: . newCP - C Array of [x, y, z] New/Updated Control Point Coordinates defining the geometry (See DMPlexGeomDataAndGrads() for format)
4491: . newW - C Array of New/Updated Control Point Weights associated with the Control Points defining the new geometry (See DMPlexGemGrads() for format)
4492: . autoInflate - PetscBool Flag denoting if the user would like to inflate the DM points to the new geometry.
4493: . saveGeom - PetscBool Flag denoting if the user would iike to save the new geometry to a file.
4494: - stpName - Char Array indicating the name of the file to save the new geometry to. Extension must be included and will denote type of file written.
4495: *.stp or *.step = STEP File
4496: *.igs or *.iges = IGES File
4497: *.egads = EGADS File
4498: *.brep = BRep File (OpenCASCADE File)
4500: Output Parameter:
4501: . dm - The updated DM object representing the mesh with PetscContainers containing the updated/modified geometry
4503: Level: intermediate
4505: Note:
4506: Functionality not available for DMPlexes with attached EGADSlite geometry files (.egadslite).
4508: .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateGeom()`, `DMPlexGeomDataAndGrads()`
4509: @*/
4510: PetscErrorCode DMPlexModifyGeomModel(DM dm, MPI_Comm comm, PetscScalar newCP[], PetscScalar newW[], PetscBool autoInflate, PetscBool saveGeom, const char *stpName) PeNS
4511: {
4512: #if defined(PETSC_HAVE_EGADS)
4513: /* EGADS/EGADSlite variables */
4514: ego context, model, geom, *bodies, *lobjs, *fobjs;
4515: int oclass, mtype, *senses, *lsenses;
4516: int Nb, Nf, Nl, id;
4517: /* PETSc variables */
4518: DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel;
4519: PetscContainer modelObj, cpHashTableObj, wHashTableObj;
4520: PetscHMapI cpHashTable = NULL, wHashTable = NULL;
4521: PetscBool islite = PETSC_FALSE;
4522: #endif
4524: #if defined(PETSC_HAVE_EGADS)
4525: PetscFunctionBegin;
4526: // Look to see if DM has a Container with either a EGADS or EGADSlite Model
4527: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
4528: if (!modelObj) {
4529: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
4530: islite = PETSC_TRUE;
4531: }
4532: PetscCheck(!islite, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot modify geometries defined by EGADSlite (.egadslite)! Please use another geometry file format STEP, IGES, EGADS or BRep");
4533: PetscCheck(modelObj, PETSC_COMM_SELF, PETSC_ERR_SUP, "DM does not have a EGADS Geometry Model attached to it!");
4535: // Get attached EGADS model (pointer)
4536: PetscCall(PetscContainerGetPointer(modelObj, &model));
4538: // Look to see if DM has Container for Geometry Control Point Data
4539: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Hash Table", (PetscObject *)&cpHashTableObj));
4540: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weights Hash Table", (PetscObject *)&wHashTableObj));
4542: PetscCheck(cpHashTableObj && wHashTableObj, PETSC_COMM_SELF, PETSC_ERR_SUP, "DM does not have required Geometry Data attached! Please run DMPlexGeomDataAndGrads() Function first.");
4544: // Get attached EGADS model Control Point and Weights Hash Tables and Data Arrays (pointer)
4545: PetscCall(PetscContainerGetPointer(cpHashTableObj, &cpHashTable));
4546: PetscCall(PetscContainerGetPointer(wHashTableObj, &wHashTable));
4548: // Get the number of bodies and body objects in the model
4549: if (islite) PetscCallEGADS(EGlite_getTopology, (model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4550: else PetscCallEGADS(EG_getTopology, (model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4552: // Get all Faces on the body
4553: ego body = bodies[0];
4554: if (islite) PetscCallEGADS(EGlite_getBodyTopos, (body, NULL, FACE, &Nf, &fobjs));
4555: else PetscCallEGADS(EG_getBodyTopos, (body, NULL, FACE, &Nf, &fobjs));
4557: ego newGeom[Nf];
4558: ego newFaces[Nf];
4560: // Update Control Point and Weight definitions for each surface
4561: for (int jj = 0; jj < Nf; ++jj) {
4562: ego face = fobjs[jj];
4563: ego bRef, bPrev, bNext;
4564: ego fgeom;
4565: int offset;
4566: int boclass, bmtype, *bpinfo;
4567: double *bprv;
4569: // Get FACE ID and other Geometry Data
4570: if (islite) {
4571: id = EGlite_indexBodyTopo(body, face);
4572: PetscCallEGADS(EGlite_getTopology, (face, &fgeom, &oclass, &mtype, NULL, &Nl, &lobjs, &lsenses));
4573: PetscCallEGADS(EGlite_getGeometry, (fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
4574: PetscCallEGADS(EGlite_getInfo, (fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
4575: } else {
4576: id = EG_indexBodyTopo(body, face);
4577: PetscCallEGADS(EG_getTopology, (face, &fgeom, &oclass, &mtype, NULL, &Nl, &lobjs, &lsenses));
4578: PetscCallEGADS(EG_getGeometry, (fgeom, &boclass, &bmtype, &bRef, &bpinfo, &bprv));
4579: PetscCallEGADS(EG_getInfo, (fgeom, &boclass, &bmtype, &bRef, &bPrev, &bNext));
4580: }
4582: // Update Control Points
4583: PetscHashIter CPiter, Witer;
4584: PetscBool CPfound, Wfound;
4585: PetscInt faceCPStartRow, faceWStartRow;
4587: PetscCall(PetscHMapIFind(cpHashTable, id, &CPiter, &CPfound));
4588: PetscCheck(CPfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "FACE ID not found in Control Point Hash Table");
4589: PetscCall(PetscHMapIGet(cpHashTable, id, &faceCPStartRow));
4591: PetscCall(PetscHMapIFind(wHashTable, id, &Witer, &Wfound));
4592: PetscCheck(Wfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "FACE ID not found in Control Point Weights Hash Table");
4593: PetscCall(PetscHMapIGet(wHashTable, id, &faceWStartRow));
4595: // UPDATE CONTROL POINTS Locations
4596: offset = bpinfo[3] + bpinfo[6];
4597: for (int ii = 0; ii < 3 * bpinfo[2] * bpinfo[5]; ++ii) bprv[offset + ii] = newCP[faceCPStartRow + ii];
4599: // UPDATE CONTROL POINT WEIGHTS
4600: offset = bpinfo[3] + bpinfo[6] + 3 * bpinfo[2] * bpinfo[5];
4601: for (int ii = 0; ii < bpinfo[2] * bpinfo[5]; ++ii) bprv[offset + ii] = newW[faceWStartRow + ii];
4603: // Get Context from FACE
4604: context = NULL;
4605: PetscCallEGADS(EG_getContext, (face, &context)); // Does not have an EGlite_ version KNOWN_ISSUE
4607: // Create New Surface
4608: ego newgeom;
4609: PetscCallEGADS(EG_makeGeometry, (context, SURFACE, BSPLINE, NULL, bpinfo, bprv, &newgeom)); // Does not have an EGlite_ version KNOWN_ISSUE
4611: // Create new FACE based on new SURFACE geometry
4612: double data[4];
4613: int periodic;
4614: if (islite) PetscCallEGADS(EGlite_getRange, (newgeom, data, &periodic));
4615: else PetscCallEGADS(EG_getRange, (newgeom, data, &periodic));
4617: ego newface;
4618: PetscCallEGADS(EG_makeFace, (newgeom, SFORWARD, data, &newface)); // Does not have an EGlite_ version KNOWN_ISSUE
4619: newGeom[jj] = newgeom;
4620: newFaces[jj] = newface;
4621: }
4622: // Could be replaced by DMPlexFreeGeomObject
4623: if (islite) EGlite_free(fobjs);
4624: else EG_free(fobjs);
4626: // Sew New Faces together to get a new model
4627: ego newmodel;
4628: PetscCall(EG_sewFaces(Nf, newFaces, 0.0, 0, &newmodel)); // Does not have an EGlite_ version KNOWN_ISSUE
4629: for (PetscInt f = 0; f < Nf; ++f) {
4630: PetscCallEGADS(EG_deleteObject, (newFaces[f]));
4631: PetscCallEGADS(EG_deleteObject, (newGeom[f]));
4632: }
4634: // Get the total number of NODEs on the original geometry. (This will be the same for the new geometry)
4635: int totalNumNode;
4636: ego *nobjTotal;
4637: if (islite) {
4638: PetscCallEGADS(EGlite_getBodyTopos, (body, NULL, NODE, &totalNumNode, &nobjTotal));
4639: EGlite_free(nobjTotal);
4640: } else {
4641: PetscCallEGADS(EG_getBodyTopos, (body, NULL, NODE, &totalNumNode, &nobjTotal));
4642: EG_free(nobjTotal);
4643: } // Could be replaced with DMPlexFreeGeomObject
4645: // Initialize vector to store equivalent NODE indices between the 2 geometries
4646: // FORMAT :: vector index is the Original Geometry's NODE ID, the vector Value is the New Geometry's NODE ID
4647: int nodeIDEquiv[totalNumNode + 1];
4649: // Now we need to Map the NODE and EDGE IDs from each Model
4650: if (islite) PetscCallEGADS(EGlite_getBodyTopos, (body, NULL, FACE, &Nf, &fobjs));
4651: else PetscCallEGADS(EG_getBodyTopos, (body, NULL, FACE, &Nf, &fobjs));
4653: // New CAD
4654: ego *newbodies, newgeomtest, *nfobjs;
4655: int nNf, newNb, newoclass, newmtype, *newsenses;
4656: if (islite) PetscCallEGADS(EGlite_getTopology, (newmodel, &newgeomtest, &newoclass, &newmtype, NULL, &newNb, &newbodies, &newsenses));
4657: else PetscCallEGADS(EG_getTopology, (newmodel, &newgeomtest, &newoclass, &newmtype, NULL, &newNb, &newbodies, &newsenses));
4659: ego newbody = newbodies[0];
4660: if (islite) PetscCallEGADS(EGlite_getBodyTopos, (newbody, NULL, FACE, &nNf, &nfobjs));
4661: else PetscCallEGADS(EG_getBodyTopos, (newbody, NULL, FACE, &nNf, &nfobjs));
4663: PetscCheck(newNb == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ERROR :: newNb > 1 || newNb = %d", newNb);
4665: // Find Equivalent Nodes
4666: for (int ii = 0; ii < Nf; ++ii) {
4667: double fdata[4];
4668: int peri;
4670: // Get Current FACE [u, v] Ranges
4671: if (islite) PetscCallEGADS(EGlite_getRange, (fobjs[ii], fdata, &peri));
4672: else PetscCallEGADS(EG_getRange, (fobjs[ii], fdata, &peri));
4674: // Equate NODE IDs between 2 FACEs by working through (u, v) limits of FACE
4675: for (int jj = 0; jj < 2; ++jj) {
4676: for (int kk = 2; kk < 4; ++kk) {
4677: double params[2] = {fdata[jj], fdata[kk]};
4678: double eval[18];
4679: if (islite) PetscCallEGADS(EGlite_evaluate, (fobjs[ii], params, eval));
4680: else PetscCallEGADS(EG_evaluate, (fobjs[ii], params, eval));
4682: // Original Body
4683: ego *nobjsOrigFace;
4684: int origNn;
4685: if (islite) PetscCallEGADS(EGlite_getBodyTopos, (body, fobjs[ii], NODE, &origNn, &nobjsOrigFace));
4686: else PetscCallEGADS(EG_getBodyTopos, (body, fobjs[ii], NODE, &origNn, &nobjsOrigFace));
4688: double minVal = 1.0E10;
4689: double evalCheck[18];
4690: int equivOrigNodeID = -1;
4691: for (int mm = 0; mm < origNn; ++mm) {
4692: double delta = 1.0E10;
4693: if (islite) PetscCallEGADS(EGlite_evaluate, (nobjsOrigFace[mm], NULL, evalCheck));
4694: else PetscCallEGADS(EG_evaluate, (nobjsOrigFace[mm], NULL, evalCheck));
4696: delta = PetscSqrtReal(PetscSqr(evalCheck[0] - eval[0]) + PetscSqr(evalCheck[1] - eval[1]) + PetscSqr(evalCheck[2] - eval[2]));
4698: if (delta < minVal) {
4699: if (islite) equivOrigNodeID = EGlite_indexBodyTopo(body, nobjsOrigFace[mm]);
4700: else equivOrigNodeID = EG_indexBodyTopo(body, nobjsOrigFace[mm]);
4702: minVal = delta;
4703: }
4704: }
4705: // Could be replaced with DMPlexFreeGeomObject
4706: if (islite) EGlite_free(nobjsOrigFace);
4707: else EG_free(nobjsOrigFace);
4709: // New Body
4710: ego *nobjsNewFace;
4711: int newNn;
4712: if (islite) PetscCallEGADS(EGlite_getBodyTopos, (newbody, nfobjs[ii], NODE, &newNn, &nobjsNewFace));
4713: else PetscCallEGADS(EG_getBodyTopos, (newbody, nfobjs[ii], NODE, &newNn, &nobjsNewFace));
4715: minVal = 1.0E10;
4716: int equivNewNodeID = -1;
4717: for (int mm = 0; mm < newNn; ++mm) {
4718: double delta = 1.0E10;
4719: if (islite) PetscCallEGADS(EGlite_evaluate, (nobjsNewFace[mm], NULL, evalCheck));
4720: else PetscCallEGADS(EG_evaluate, (nobjsNewFace[mm], NULL, evalCheck));
4722: delta = PetscSqrtReal(PetscSqr(evalCheck[0] - eval[0]) + PetscSqr(evalCheck[1] - eval[1]) + PetscSqr(evalCheck[2] - eval[2]));
4724: if (delta < minVal) {
4725: if (islite) equivNewNodeID = EGlite_indexBodyTopo(newbody, nobjsNewFace[mm]);
4726: else equivNewNodeID = EG_indexBodyTopo(newbody, nobjsNewFace[mm]);
4728: minVal = delta;
4729: }
4730: }
4731: if (islite) EGlite_free(nobjsNewFace);
4732: else EG_free(nobjsNewFace);
4734: // Store equivalent NODE IDs
4735: nodeIDEquiv[equivOrigNodeID] = equivNewNodeID;
4736: }
4737: }
4738: }
4740: // Find Equivalent EDGEs
4741: // Get total number of EDGEs on Original Geometry
4742: int totalNumEdge;
4743: ego *eobjsOrig;
4744: if (islite) {
4745: PetscCallEGADS(EGlite_getBodyTopos, (body, NULL, EDGE, &totalNumEdge, &eobjsOrig));
4746: EGlite_free(eobjsOrig);
4747: } else {
4748: PetscCallEGADS(EG_getBodyTopos, (body, NULL, EDGE, &totalNumEdge, &eobjsOrig));
4749: EG_free(eobjsOrig);
4750: }
4752: // Get total number of EDGEs on New Geometry
4753: int totalNumEdgeNew;
4754: ego *eobjsNew;
4755: if (islite) {
4756: PetscCallEGADS(EGlite_getBodyTopos, (newbody, NULL, EDGE, &totalNumEdgeNew, &eobjsNew));
4757: EGlite_free(eobjsNew);
4758: } else {
4759: PetscCallEGADS(EG_getBodyTopos, (newbody, NULL, EDGE, &totalNumEdgeNew, &eobjsNew));
4760: EG_free(eobjsNew);
4761: }
4763: // Initialize EDGE ID equivalent vector
4764: // FORMAT :: vector index is the Original Geometry's EDGE ID, the vector Value is the New Geometry's EDGE ID
4765: int edgeIDEquiv[totalNumEdge + 1];
4767: // Find Equivalent EDGEs
4768: for (int ii = 0; ii < Nf; ++ii) {
4769: // Get Original Geometry EDGE's NODEs
4770: int numOrigEdge, numNewEdge;
4771: if (islite) {
4772: PetscCallEGADS(EGlite_getBodyTopos, (body, fobjs[ii], EDGE, &numOrigEdge, &eobjsOrig));
4773: PetscCallEGADS(EGlite_getBodyTopos, (newbody, nfobjs[ii], EDGE, &numNewEdge, &eobjsNew));
4774: } else {
4775: PetscCallEGADS(EG_getBodyTopos, (body, fobjs[ii], EDGE, &numOrigEdge, &eobjsOrig));
4776: PetscCallEGADS(EG_getBodyTopos, (newbody, nfobjs[ii], EDGE, &numNewEdge, &eobjsNew));
4777: }
4779: // new loop below
4780: for (int nn = 0; nn < numOrigEdge; ++nn) {
4781: ego origEdge = eobjsOrig[nn];
4782: ego geomEdgeOrig, *nobjsOrig;
4783: int oclassEdgeOrig, mtypeEdgeOrig;
4784: int NnOrig, *nsensesEdgeOrig;
4786: if (islite) PetscCallEGADS(EGlite_getTopology, (origEdge, &geomEdgeOrig, &oclassEdgeOrig, &mtypeEdgeOrig, NULL, &NnOrig, &nobjsOrig, &nsensesEdgeOrig));
4787: else PetscCallEGADS(EG_getTopology, (origEdge, &geomEdgeOrig, &oclassEdgeOrig, &mtypeEdgeOrig, NULL, &NnOrig, &nobjsOrig, &nsensesEdgeOrig));
4789: PetscBool isSame = PETSC_FALSE;
4790: for (int jj = 0; jj < numNewEdge; ++jj) {
4791: ego newEdge = eobjsNew[jj];
4792: ego geomEdgeNew, *nobjsNew;
4793: int oclassEdgeNew, mtypeEdgeNew;
4794: int NnNew, *nsensesEdgeNew;
4796: if (islite) PetscCallEGADS(EGlite_getTopology, (newEdge, &geomEdgeNew, &oclassEdgeNew, &mtypeEdgeNew, NULL, &NnNew, &nobjsNew, &nsensesEdgeNew));
4797: else PetscCallEGADS(EG_getTopology, (newEdge, &geomEdgeNew, &oclassEdgeNew, &mtypeEdgeNew, NULL, &NnNew, &nobjsNew, &nsensesEdgeNew));
4799: if (mtypeEdgeOrig == mtypeEdgeNew) {
4800: // Only operate if the EDGE types are the same
4801: for (int kk = 0; kk < NnNew; ++kk) {
4802: int nodeIDOrigGeom, nodeIDNewGeom;
4803: if (islite) {
4804: nodeIDOrigGeom = EGlite_indexBodyTopo(body, nobjsOrig[kk]);
4805: nodeIDNewGeom = EGlite_indexBodyTopo(newbody, nobjsNew[kk]);
4806: } else {
4807: nodeIDOrigGeom = EG_indexBodyTopo(body, nobjsOrig[kk]);
4808: nodeIDNewGeom = EG_indexBodyTopo(newbody, nobjsNew[kk]);
4809: }
4811: if (nodeIDNewGeom == nodeIDEquiv[nodeIDOrigGeom]) {
4812: isSame = PETSC_TRUE;
4813: } else {
4814: isSame = PETSC_FALSE;
4815: kk = NnNew; // skip ahead because first NODE failed test and order is important
4816: }
4817: }
4819: if (isSame == PETSC_TRUE) {
4820: int edgeIDOrig, edgeIDNew;
4821: if (islite) {
4822: edgeIDOrig = EGlite_indexBodyTopo(body, origEdge);
4823: edgeIDNew = EGlite_indexBodyTopo(newbody, newEdge);
4824: } else {
4825: edgeIDOrig = EG_indexBodyTopo(body, origEdge);
4826: edgeIDNew = EG_indexBodyTopo(newbody, newEdge);
4827: }
4829: edgeIDEquiv[edgeIDOrig] = edgeIDNew;
4830: jj = numNewEdge;
4831: }
4832: }
4833: }
4834: }
4835: if (islite) {
4836: EGlite_free(eobjsOrig);
4837: EGlite_free(eobjsNew);
4838: } else {
4839: EG_free(eobjsOrig);
4840: EG_free(eobjsNew);
4841: }
4842: }
4843: if (islite) {
4844: EGlite_free(fobjs);
4845: EGlite_free(nfobjs);
4846: } else {
4847: EG_free(fobjs);
4848: EG_free(nfobjs);
4849: }
4851: // Modify labels to point to the IDs on the new Geometry
4852: IS isNodeID, isEdgeID;
4854: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
4855: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
4856: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
4857: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
4859: PetscCall(ISCreateGeneral(comm, totalNumNode + 1, nodeIDEquiv, PETSC_COPY_VALUES, &isNodeID));
4860: PetscCall(ISCreateGeneral(comm, totalNumEdge + 1, edgeIDEquiv, PETSC_COPY_VALUES, &isEdgeID));
4861: /* Do not perform check. Np may != Nv due to Degenerate Geometry which is not stored in labels. */
4862: /* We do not know in advance which IDs have been omitted. This may also change due to geometry modifications. */
4863: PetscCall(DMLabelRewriteValues(vertexLabel, isNodeID));
4864: PetscCall(DMLabelRewriteValues(edgeLabel, isEdgeID));
4865: PetscCall(ISDestroy(&isNodeID));
4866: PetscCall(ISDestroy(&isEdgeID));
4868: // Attempt to point to the new geometry
4869: PetscCallEGADS(EG_deleteObject, (model));
4870: PetscCall(PetscContainerSetPointer(modelObj, newmodel));
4872: // save updated model to file
4873: if (saveGeom == PETSC_TRUE && stpName != NULL) PetscCall(EG_saveModel(newmodel, stpName));
4875: // Inflate Mesh to EGADS Model
4876: if (autoInflate == PETSC_TRUE) PetscCall(DMPlexInflateToGeomModel(dm, PETSC_TRUE));
4877: PetscFunctionReturn(PETSC_SUCCESS);
4878: #else
4879: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This method requires EGADS support. Reconfigure using --download-egads");
4880: #endif
4881: }
4883: /*@C
4884: DMPlexGetGeomModelTUV - Gets the [t] (EDGES) and [u, v] (FACES) geometry parameters of DM points that are associated geometry relationships. Requires a DM with a EGADS model attached.
4886: Collective
4888: Input Parameter:
4889: . dm - The DM object representing the mesh with PetscContainer containing an EGADS geometry model
4891: Level: intermediate
4893: .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateGeom()`, `DMPlexGeomDataAndGrads()`
4894: @*/
4895: PetscErrorCode DMPlexGetGeomModelTUV(DM dm) PeNS
4896: {
4897: #if defined(PETSC_HAVE_EGADS)
4898: /* EGADS Variables */
4899: ego model, geom, body, face, edge;
4900: ego *bodies;
4901: int Nb, oclass, mtype, *senses;
4902: double result[4];
4903: /* PETSc Variables */
4904: DM cdm;
4905: PetscContainer modelObj;
4906: DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel;
4907: Vec coordinates;
4908: PetscScalar *coords;
4909: PetscInt bodyID, faceID, edgeID, vertexID;
4910: PetscInt cdim, vStart, vEnd, v;
4911: PetscBool islite = PETSC_FALSE;
4912: #endif
4914: PetscFunctionBegin;
4915: #if defined(PETSC_HAVE_EGADS)
4916: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
4917: if (!modelObj) {
4918: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
4919: islite = PETSC_TRUE;
4920: }
4921: if (!modelObj) PetscFunctionReturn(PETSC_SUCCESS);
4923: PetscCall(DMGetCoordinateDim(dm, &cdim));
4924: PetscCall(DMGetCoordinateDM(dm, &cdm));
4925: PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
4926: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
4927: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
4928: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
4929: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
4931: PetscCall(PetscContainerGetPointer(modelObj, &model));
4933: if (islite) PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4934: else PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4936: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4937: PetscCall(VecGetArrayWrite(coordinates, &coords));
4939: // Define t, u, v arrays to be stored in a PetscContainer after populated
4940: PetscScalar *t_point, *u_point, *v_point;
4941: PetscCall(PetscMalloc1(vEnd - vStart, &t_point));
4942: PetscCall(PetscMalloc1(vEnd - vStart, &u_point));
4943: PetscCall(PetscMalloc1(vEnd - vStart, &v_point));
4945: for (v = vStart; v < vEnd; ++v) {
4946: PetscScalar *vcoords;
4948: PetscCall(DMLabelGetValue(bodyLabel, v, &bodyID));
4949: PetscCall(DMLabelGetValue(faceLabel, v, &faceID));
4950: PetscCall(DMLabelGetValue(edgeLabel, v, &edgeID));
4951: PetscCall(DMLabelGetValue(vertexLabel, v, &vertexID));
4953: // TODO Figure out why this is unknown sometimes
4954: if (bodyID < 0 && Nb == 1) bodyID = 0;
4955: PetscCheck(bodyID >= 0 && bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " for vertex %" PetscInt_FMT " is not in [0, %d)", bodyID, v, Nb);
4956: body = bodies[bodyID];
4958: PetscCall(DMPlexPointLocalRef(cdm, v, coords, (void *)&vcoords));
4959: if (edgeID > 0) {
4960: /* Snap to EDGE at nearest location */
4961: double params[1];
4963: if (islite) {
4964: PetscCall(EGlite_objectBodyTopo(body, EDGE, edgeID, &edge));
4965: PetscCall(EGlite_invEvaluate(edge, vcoords, params, result));
4966: } // Get (t) of nearest point on EDGE
4967: else {
4968: PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &edge));
4969: PetscCall(EG_invEvaluate(edge, vcoords, params, result));
4970: } // Get (t) of nearest point on EDGE
4972: t_point[v - vStart] = params[0];
4973: u_point[v - vStart] = 0.0;
4974: v_point[v - vStart] = 0.0;
4975: } else if (faceID > 0) {
4976: /* Snap to FACE at nearest location */
4977: double params[2];
4979: if (islite) {
4980: PetscCall(EGlite_objectBodyTopo(body, FACE, faceID, &face));
4981: PetscCall(EGlite_invEvaluate(face, vcoords, params, result));
4982: } // Get (x,y,z) of nearest point on FACE
4983: else {
4984: PetscCall(EG_objectBodyTopo(body, FACE, faceID, &face));
4985: PetscCall(EG_invEvaluate(face, vcoords, params, result));
4986: } // Get (x,y,z) of nearest point on FACE
4988: t_point[v - vStart] = 0.0;
4989: u_point[v - vStart] = params[0];
4990: v_point[v - vStart] = params[1];
4991: } else {
4992: t_point[v - vStart] = 0.0;
4993: u_point[v - vStart] = 0.0;
4994: v_point[v - vStart] = 0.0;
4995: }
4996: }
4997: PetscCall(VecRestoreArrayWrite(coordinates, &coords));
4998: /* Clear out global coordinates */
4999: PetscCall(VecDestroy(&dm->coordinates[0].x));
5001: /* Store in PetscContainters */
5002: {
5003: PetscContainer t_pointObj, u_pointObj, v_pointObj;
5005: PetscCall(PetscObjectQuery((PetscObject)dm, "Point - Edge t Parameter", (PetscObject *)&t_pointObj));
5006: if (!t_pointObj) {
5007: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &t_pointObj));
5008: PetscCall(PetscContainerSetPointer(t_pointObj, t_point));
5009: PetscCall(PetscObjectCompose((PetscObject)dm, "Point - Edge t Parameter", (PetscObject)t_pointObj));
5010: PetscCall(PetscContainerSetCtxDestroy(t_pointObj, PetscCtxDestroyDefault));
5011: PetscCall(PetscContainerDestroy(&t_pointObj));
5012: } else {
5013: void *old;
5015: PetscCall(PetscContainerGetPointer(t_pointObj, &old));
5016: PetscCall(PetscFree(old));
5017: PetscCall(PetscContainerSetPointer(t_pointObj, t_point));
5018: }
5020: PetscCall(PetscObjectQuery((PetscObject)dm, "Point - Face u Parameter", (PetscObject *)&u_pointObj));
5021: if (!u_pointObj) {
5022: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &u_pointObj));
5023: PetscCall(PetscContainerSetPointer(u_pointObj, u_point));
5024: PetscCall(PetscObjectCompose((PetscObject)dm, "Point - Face u Parameter", (PetscObject)u_pointObj));
5025: PetscCall(PetscContainerSetCtxDestroy(u_pointObj, PetscCtxDestroyDefault));
5026: PetscCall(PetscContainerDestroy(&u_pointObj));
5027: } else {
5028: void *old;
5030: PetscCall(PetscContainerGetPointer(u_pointObj, &old));
5031: PetscCall(PetscFree(old));
5032: PetscCall(PetscContainerSetPointer(u_pointObj, u_point));
5033: }
5035: PetscCall(PetscObjectQuery((PetscObject)dm, "Point - Face v Parameter", (PetscObject *)&v_pointObj));
5036: if (!v_pointObj) {
5037: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &v_pointObj));
5038: PetscCall(PetscContainerSetPointer(v_pointObj, v_point));
5039: PetscCall(PetscObjectCompose((PetscObject)dm, "Point - Face v Parameter", (PetscObject)v_pointObj));
5040: PetscCall(PetscContainerSetCtxDestroy(v_pointObj, PetscCtxDestroyDefault));
5041: PetscCall(PetscContainerDestroy(&v_pointObj));
5042: } else {
5043: void *old;
5045: PetscCall(PetscContainerGetPointer(v_pointObj, &old));
5046: PetscCall(PetscFree(old));
5047: PetscCall(PetscContainerSetPointer(v_pointObj, v_point));
5048: }
5049: }
5050: #endif
5051: PetscFunctionReturn(PETSC_SUCCESS);
5052: }
5054: /*@C
5055: DMPlexInflateToGeomModelUseTUV - Inflates the DM to the associated underlying geometry using the [t] {EDGES) and [u, v] (FACES} associated parameters. Requires a DM with an EGADS model attached and a previous call to DMPlexGetGeomModelTUV().
5057: Collective
5059: Input Parameter:
5060: . dm - The DM object representing the mesh with PetscContainer containing an EGADS geometry model
5062: Level: intermediate
5064: Note:
5065: The updated DM object inflated to the associated underlying geometry. This updates the [x, y, z] coordinates of DM points associated with geometry.
5067: .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateGeom()`, `DMPlexGeomDataAndGrads()`, `DMPlexGetGeomModelTUV()`
5068: @*/
5069: PetscErrorCode DMPlexInflateToGeomModelUseTUV(DM dm) PeNS
5070: {
5071: #if defined(PETSC_HAVE_EGADS)
5072: /* EGADS Variables */
5073: ego model, geom, body, face, edge, vertex;
5074: ego *bodies;
5075: int Nb, oclass, mtype, *senses;
5076: double result[18], params[2];
5077: /* PETSc Variables */
5078: DM cdm;
5079: PetscContainer modelObj;
5080: PetscContainer t_pointObj, u_pointObj, v_pointObj;
5081: DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel;
5082: Vec coordinates;
5083: PetscScalar *coords;
5084: PetscScalar *t_point, *u_point, *v_point;
5085: PetscInt bodyID, faceID, edgeID, vertexID;
5086: PetscInt cdim, d, vStart, vEnd, v;
5087: PetscBool islite = PETSC_FALSE;
5088: #endif
5090: PetscFunctionBegin;
5091: #if defined(PETSC_HAVE_EGADS)
5092: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5093: if (!modelObj) {
5094: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5095: islite = PETSC_TRUE;
5096: }
5098: PetscCall(PetscObjectQuery((PetscObject)dm, "Point - Edge t Parameter", (PetscObject *)&t_pointObj));
5099: PetscCall(PetscObjectQuery((PetscObject)dm, "Point - Face u Parameter", (PetscObject *)&u_pointObj));
5100: PetscCall(PetscObjectQuery((PetscObject)dm, "Point - Face v Parameter", (PetscObject *)&v_pointObj));
5102: if (!modelObj) PetscFunctionReturn(PETSC_SUCCESS);
5103: if (!t_pointObj) PetscFunctionReturn(PETSC_SUCCESS);
5104: if (!u_pointObj) PetscFunctionReturn(PETSC_SUCCESS);
5105: if (!v_pointObj) PetscFunctionReturn(PETSC_SUCCESS);
5107: PetscCall(DMGetCoordinateDim(dm, &cdim));
5108: PetscCall(DMGetCoordinateDM(dm, &cdm));
5109: PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
5110: PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
5111: PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
5112: PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
5113: PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
5115: PetscCall(PetscContainerGetPointer(t_pointObj, &t_point));
5116: PetscCall(PetscContainerGetPointer(u_pointObj, &u_point));
5117: PetscCall(PetscContainerGetPointer(v_pointObj, &v_point));
5119: PetscCall(PetscContainerGetPointer(modelObj, &model));
5121: if (islite) {
5122: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
5123: } else {
5124: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
5125: }
5127: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
5128: PetscCall(VecGetArrayWrite(coordinates, &coords));
5130: for (v = vStart; v < vEnd; ++v) {
5131: PetscScalar *vcoords;
5133: PetscCall(DMLabelGetValue(bodyLabel, v, &bodyID));
5134: PetscCall(DMLabelGetValue(faceLabel, v, &faceID));
5135: PetscCall(DMLabelGetValue(edgeLabel, v, &edgeID));
5136: PetscCall(DMLabelGetValue(vertexLabel, v, &vertexID));
5138: // TODO Figure out why this is unknown sometimes
5139: if (bodyID < 0 && Nb == 1) bodyID = 0;
5140: PetscCheck(bodyID >= 0 && bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " for vertex %" PetscInt_FMT " is not in [0, %d)", bodyID, v, Nb);
5141: body = bodies[bodyID];
5143: PetscCall(DMPlexPointLocalRef(cdm, v, coords, (void *)&vcoords));
5144: if (vertexID > 0) {
5145: /* Snap to Vertices */
5146: if (islite) {
5147: PetscCall(EGlite_objectBodyTopo(body, NODE, vertexID, &vertex));
5148: PetscCall(EGlite_evaluate(vertex, NULL, result));
5149: } else {
5150: PetscCall(EG_objectBodyTopo(body, NODE, vertexID, &vertex));
5151: PetscCall(EG_evaluate(vertex, NULL, result));
5152: }
5153: for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
5154: } else if (edgeID > 0) {
5155: /* Snap to EDGE */
5156: params[0] = t_point[v - vStart];
5157: if (islite) {
5158: PetscCall(EGlite_objectBodyTopo(body, EDGE, edgeID, &edge));
5159: PetscCall(EGlite_evaluate(edge, params, result));
5160: } else {
5161: PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &edge));
5162: PetscCall(EG_evaluate(edge, params, result));
5163: }
5164: for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
5165: } else if (faceID > 0) {
5166: /* Snap to FACE */
5167: params[0] = u_point[v - vStart];
5168: params[1] = v_point[v - vStart];
5169: if (islite) {
5170: PetscCall(EGlite_objectBodyTopo(body, FACE, faceID, &face));
5171: PetscCall(EGlite_evaluate(face, params, result));
5172: } else {
5173: PetscCall(EG_objectBodyTopo(body, FACE, faceID, &face));
5174: PetscCall(EG_evaluate(face, params, result));
5175: }
5176: for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
5177: }
5178: }
5179: PetscCall(VecRestoreArrayWrite(coordinates, &coords));
5180: /* Clear out global coordinates */
5181: PetscCall(VecDestroy(&dm->coordinates[0].x));
5182: #endif
5183: PetscFunctionReturn(PETSC_SUCCESS);
5184: }
5186: /*@
5187: DMPlexInflateToGeomModel - Wrapper function allowing two methods for inflating refined meshes to the underlying geometric domain.
5189: Collective
5191: Input Parameters:
5192: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5193: - useTUV - PetscBool indicating if the user would like to inflate the DMPlex to the underlying geometry
5194: using (t) for nodes on EDGEs and (u, v) for nodes on FACEs or using the nodes (x, y, z) coordinates
5195: and shortest distance routine.
5196: If useTUV = PETSC_TRUE, use the (t) or (u, v) parameters to inflate the DMPlex to the CAD geometry.
5197: If useTUV = PETSC_FALSE, use the nodes (x, y, z) coordinates and the shortest disctance routine.
5199: Notes:
5200: DM with nodal coordinates modified so that they lie on the EDGEs and FACEs of the underlying geometry.
5202: (t) and (u, v) parameters for all DMPlex nodes on EDGEs and FACEs are stored in arrays within PetscContainers attached to the DM.
5203: The containers have names "Point - Edge t Parameter", "Point - Face u Parameter", and "Point - Face v Parameter".
5204: The arrays are organized by Point 0-based ID (i.e. [v-vstart] as defined in the DMPlex.
5206: Level: intermediate
5208: .seealso: `DMPlexGetGeomModelTUV()`, `DMPlexInflateToGeomModelUseTUV()`, `DMPlexInflateToGeomModelUseXYZ()`
5209: @*/
5210: PetscErrorCode DMPlexInflateToGeomModel(DM dm, PetscBool useTUV) PeNS
5211: {
5212: PetscFunctionBeginHot;
5213: if (useTUV) {
5214: PetscCall(DMPlexGetGeomModelTUV(dm));
5215: PetscCall(DMPlexInflateToGeomModelUseTUV(dm));
5216: } else {
5217: PetscCall(DMPlexInflateToGeomModelUseXYZ(dm));
5218: }
5219: PetscFunctionReturn(PETSC_SUCCESS);
5220: }
5222: #ifdef PETSC_HAVE_EGADS
5223: /*@C
5224: DMPlexGetGeomModelBodies - Returns an array of `PetscGeom` BODY objects attached to the referenced geometric model entity as well as the number of BODYs.
5226: Collective
5228: Input Parameter:
5229: . dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5231: Output Parameters:
5232: + bodies - Array of PetscGeom BODY objects referenced by the geometric model.
5233: - numBodies - Number of BODYs referenced by the geometric model. Also the size of **bodies array.
5235: Level: intermediate
5237: .seealso:
5238: @*/
5239: PetscErrorCode DMPlexGetGeomModelBodies(DM dm, PetscGeom **bodies, PetscInt *numBodies) PeNS
5240: {
5241: PetscFunctionBeginHot;
5242: PetscContainer modelObj;
5243: PetscBool islite = PETSC_FALSE;
5244: ego model, geom;
5245: int oclass, mtype;
5246: int *senses;
5248: /* Determine which type of EGADS model is attached to the DM */
5249: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5250: if (!modelObj) {
5251: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5252: islite = PETSC_TRUE;
5253: }
5255: // Get attached EGADS or EGADSlite model (pointer)
5256: PetscCall(PetscContainerGetPointer(modelObj, &model));
5258: if (islite) {
5259: PetscCall(EGlite_getTopology(model, &geom, &oclass, &mtype, NULL, numBodies, bodies, &senses));
5260: } else {
5261: PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, numBodies, bodies, &senses));
5262: }
5263: PetscFunctionReturn(PETSC_SUCCESS);
5264: }
5266: /*@C
5267: DMPlexGetGeomModelBodyShells - Returns an array of `PetscGeom` SHELL objects attached to the referenced BODY geometric entity as well as the number of SHELLs.
5269: Collective
5271: Input Parameters:
5272: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5273: - body - PetscGeom BODY object containing the SHELL objects of interest.
5275: Output Parameters:
5276: + shells - Array of PetscGeom SHELL objects referenced by the PetscGeom BODY object
5277: - numShells - Number of SHELLs referenced by the PetscGeom BODY object. Also the size of **shells array.
5279: Level: intermediate
5281: .seealso:
5282: @*/
5283: PetscErrorCode DMPlexGetGeomModelBodyShells(DM dm, PetscGeom body, PetscGeom **shells, PetscInt *numShells) PeNS
5284: {
5285: PetscFunctionBeginHot;
5286: #ifdef PETSC_HAVE_EGADS
5287: PetscContainer modelObj;
5288: PetscBool islite = PETSC_FALSE;
5290: /* Determine which type of EGADS model is attached to the DM */
5291: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5292: if (!modelObj) {
5293: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5294: islite = PETSC_TRUE;
5295: }
5297: if (islite) {
5298: PetscCall(EGlite_getBodyTopos(body, NULL, SHELL, numShells, shells));
5299: } else {
5300: PetscCall(EG_getBodyTopos(body, NULL, SHELL, numShells, shells));
5301: }
5302: #endif
5303: PetscFunctionReturn(PETSC_SUCCESS);
5304: }
5306: /*@C
5307: DMPlexGetGeomModelBodyFaces - Returns an array of `PetscGeom` FACE objects attached to the referenced BODY geometric entity as well as the number of FACEs.
5309: Collective
5311: Input Parameters:
5312: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5313: - body - PetscGeom BODY object containing the FACE objects of interest.
5315: Output Parameters:
5316: + faces - Array of PetscGeom FACE objects referenced by the PetscGeom BODY object
5317: - numFaces - Number of FACEs referenced by the PetscGeom BODY object. Also the size of **faces array.
5319: Level: intermediate
5321: .seealso:
5322: @*/
5323: PetscErrorCode DMPlexGetGeomModelBodyFaces(DM dm, PetscGeom body, PetscGeom **faces, PetscInt *numFaces) PeNS
5324: {
5325: PetscFunctionBeginHot;
5326: #ifdef PETSC_HAVE_EGADS
5327: PetscContainer modelObj;
5328: PetscBool islite = PETSC_FALSE;
5330: /* Determine which type of EGADS model is attached to the DM */
5331: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5332: if (!modelObj) {
5333: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5334: islite = PETSC_TRUE;
5335: }
5337: if (islite) {
5338: PetscCall(EGlite_getBodyTopos(body, NULL, FACE, numFaces, faces));
5339: } else {
5340: PetscCall(EG_getBodyTopos(body, NULL, FACE, numFaces, faces));
5341: }
5342: #endif
5343: PetscFunctionReturn(PETSC_SUCCESS);
5344: }
5346: /*@C
5347: DMPlexGetGeomModelBodyLoops - Returns an array of `PetscGeom` Loop objects attached to the referenced BODY geometric entity as well as the number of LOOPs.
5349: Collective
5351: Input Parameters:
5352: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5353: - body - PetscGeom BODY object containing the LOOP objects of interest.
5355: Output Parameters:
5356: + loops - Array of PetscGeom FACE objects referenced by the PetscGeom SHELL object
5357: - numLoops - Number of LOOPs referenced by the PetscGeom BODY object. Also the size of **loops array.
5359: Level: intermediate
5361: .seealso:
5362: @*/
5363: PetscErrorCode DMPlexGetGeomModelBodyLoops(DM dm, PetscGeom body, PetscGeom **loops, PetscInt *numLoops) PeNS
5364: {
5365: PetscFunctionBeginHot;
5366: #ifdef PETSC_HAVE_EGADS
5367: PetscContainer modelObj;
5368: PetscBool islite = PETSC_FALSE;
5370: /* Determine which type of EGADS model is attached to the DM */
5371: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5372: if (!modelObj) {
5373: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5374: islite = PETSC_TRUE;
5375: }
5377: if (islite) {
5378: PetscCall(EGlite_getBodyTopos(body, NULL, LOOP, numLoops, loops));
5379: } else {
5380: PetscCall(EG_getBodyTopos(body, NULL, LOOP, numLoops, loops));
5381: }
5382: #endif
5383: PetscFunctionReturn(PETSC_SUCCESS);
5384: }
5386: /*@C
5387: DMPlexGetGeomModelShellFaces - Returns an array of `PetscGeom` FACE objects attached to the referenced SHELL geometric entity as well as the number of FACEs.
5389: Collective
5391: Input Parameters:
5392: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5393: . body - PetscGeom BODY object containing the FACE objects of interest.
5394: - shell - PetscGeom SHELL object with FACEs of interest.
5396: Output Parameters:
5397: + faces - Array of PetscGeom FACE objects referenced by the PetscGeom SHELL object
5398: - numFaces - Number of FACEs referenced by the PetscGeom SHELL object. Also the size of **faces array.
5400: Level: intermediate
5402: .seealso:
5403: @*/
5404: PetscErrorCode DMPlexGetGeomModelShellFaces(DM dm, PetscGeom body, PetscGeom shell, PetscGeom **faces, PetscInt *numFaces) PeNS
5405: {
5406: PetscFunctionBeginHot;
5407: #ifdef PETSC_HAVE_EGADS
5408: PetscContainer modelObj;
5409: PetscBool islite = PETSC_FALSE;
5411: /* Determine which type of EGADS model is attached to the DM */
5412: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5413: if (!modelObj) {
5414: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5415: islite = PETSC_TRUE;
5416: }
5418: if (islite) {
5419: PetscCall(EGlite_getBodyTopos(body, shell, FACE, numFaces, faces));
5420: } else {
5421: PetscCall(EG_getBodyTopos(body, shell, FACE, numFaces, faces));
5422: }
5423: #endif
5424: PetscFunctionReturn(PETSC_SUCCESS);
5425: }
5427: /*@C
5428: DMPlexGetGeomModelFaceLoops - Returns an array of `PetscGeom` LOOP objects attached to the referenced FACE geometric entity as well as the number of LOOPs.
5430: Collective
5432: Input Parameters:
5433: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5434: . body - PetscGeom BODY object containing the LOOP objects of interest.
5435: - face - PetscGeom FACE object with LOOPs of interest.
5437: Output Parameters:
5438: + loops - Array of PetscGeom LOOP objects referenced by the PetscGeom FACE object
5439: - numLoops - Number of LOOPs referenced by the PetscGeom FACE object. Also the size of **loops array.
5441: Level: intermediate
5443: .seealso:
5444: @*/
5445: PetscErrorCode DMPlexGetGeomModelFaceLoops(DM dm, PetscGeom body, PetscGeom face, PetscGeom **loops, PetscInt *numLoops) PeNS
5446: {
5447: PetscFunctionBeginHot;
5448: #ifdef PETSC_HAVE_EGADS
5449: PetscContainer modelObj;
5450: PetscBool islite = PETSC_FALSE;
5452: /* Determine which type of EGADS model is attached to the DM */
5453: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5454: if (!modelObj) {
5455: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5456: islite = PETSC_TRUE;
5457: }
5459: if (islite) {
5460: PetscCall(EGlite_getBodyTopos(body, face, LOOP, numLoops, loops));
5461: } else {
5462: PetscCall(EG_getBodyTopos(body, face, LOOP, numLoops, loops));
5463: }
5464: #endif
5465: PetscFunctionReturn(PETSC_SUCCESS);
5466: }
5468: /*@C
5469: DMPlexGetGeomModelFaceEdges - Returns an array of `PetscGeom` EDGE objects attached to the referenced FACE geometric entity as well as the number of EDGEs.
5471: Collective
5473: Input Parameters:
5474: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5475: . body - PetscGeom Body object containing the EDGE objects of interest.
5476: - face - PetscGeom FACE object with EDGEs of interest.
5478: Output Parameters:
5479: + edges - Array of PetscGeom EDGE objects referenced by the PetscGeom FACE object
5480: - numEdges - Number of EDGEs referenced by the PetscGeom FACE object. Also the size of **edges array.
5482: Level: intermediate
5484: .seealso:
5485: @*/
5486: PetscErrorCode DMPlexGetGeomModelFaceEdges(DM dm, PetscGeom body, PetscGeom face, PetscGeom **edges, PetscInt *numEdges) PeNS
5487: {
5488: PetscFunctionBeginHot;
5489: #ifdef PETSC_HAVE_EGADS
5490: PetscContainer modelObj;
5491: PetscBool islite = PETSC_FALSE;
5493: /* Determine which type of EGADS model is attached to the DM */
5494: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5495: if (!modelObj) {
5496: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5497: islite = PETSC_TRUE;
5498: }
5500: if (islite) {
5501: PetscCall(EGlite_getBodyTopos(body, face, EDGE, numEdges, edges));
5502: } else {
5503: PetscCall(EG_getBodyTopos(body, face, EDGE, numEdges, edges));
5504: }
5505: #endif
5506: PetscFunctionReturn(PETSC_SUCCESS);
5507: }
5509: /*@C
5510: DMPlexGetGeomModelBodyEdges - Returns an array of `PetscGeom` EDGE objects attached to the referenced BODY geometric entity as well as the number of EDGEs.
5512: Collective
5514: Input Parameters:
5515: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5516: - body - PetscGeom body object of interest.
5518: Output Parameters:
5519: + edges - Array of PetscGeom EDGE objects referenced by the PetscGeom BODY object
5520: - numEdges - Number of EDGEs referenced by the PetscGeom BODY object. Also the size of **edges array.
5522: Level: intermediate
5524: .seealso:
5525: @*/
5526: PetscErrorCode DMPlexGetGeomModelBodyEdges(DM dm, PetscGeom body, PetscGeom **edges, PetscInt *numEdges) PeNS
5527: {
5528: PetscFunctionBeginHot;
5529: #ifdef PETSC_HAVE_EGADS
5530: PetscContainer modelObj;
5531: PetscBool islite = PETSC_FALSE;
5533: /* Determine which type of EGADS model is attached to the DM */
5534: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5535: if (!modelObj) {
5536: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5537: islite = PETSC_TRUE;
5538: }
5540: if (islite) {
5541: PetscCall(EGlite_getBodyTopos(body, NULL, EDGE, numEdges, edges));
5542: } else {
5543: PetscCall(EG_getBodyTopos(body, NULL, EDGE, numEdges, edges));
5544: }
5545: #endif
5546: PetscFunctionReturn(PETSC_SUCCESS);
5547: }
5549: /*@C
5550: DMPlexGetGeomModelBodyNodes - Returns an array of `PetscGeom` NODE objects attached to the referenced BODY geometric entity as well as the number of NODES.
5552: Collective
5554: Input Parameters:
5555: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5556: - body - PetscGeom body object of interest.
5558: Output Parameters:
5559: + nodes - Array of PetscGeom NODE objects referenced by the PetscGeom BODY object
5560: - numNodes - Number of NODEs referenced by the PetscGeom BODY object. Also the size of **nodes array.
5562: Level: intermediate
5564: .seealso:
5565: @*/
5566: PetscErrorCode DMPlexGetGeomModelBodyNodes(DM dm, PetscGeom body, PetscGeom **nodes, PetscInt *numNodes) PeNS
5567: {
5568: PetscFunctionBeginHot;
5569: #ifdef PETSC_HAVE_EGADS
5570: PetscContainer modelObj;
5571: PetscBool islite = PETSC_FALSE;
5573: /* Determine which type of EGADS model is attached to the DM */
5574: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5575: if (!modelObj) {
5576: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5577: islite = PETSC_TRUE;
5578: }
5580: if (islite) {
5581: PetscCall(EGlite_getBodyTopos(body, NULL, NODE, numNodes, nodes));
5582: } else {
5583: PetscCall(EG_getBodyTopos(body, NULL, NODE, numNodes, nodes));
5584: }
5585: #endif
5586: PetscFunctionReturn(PETSC_SUCCESS);
5587: }
5589: /*@C
5590: DMPlexGetGeomModelEdgeNodes - Returns an array of `PetscGeom` NODE objects attached to the referenced EDGE geometric entity as well as the number of NODES.
5592: Collective
5594: Input Parameters:
5595: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5596: . body - PetscGeom body object containing the EDGE object of interest.
5597: - edge - PetscGeom EDGE object with NODEs of interest.
5599: Output Parameters:
5600: + nodes - Array of PetscGeom NODE objects referenced by the PetscGeom EDGE object
5601: - numNodes - Number of Nodes referenced by the PetscGeom EDGE object. Also the size of **nodes array.
5603: Level: intermediate
5605: .seealso:
5606: @*/
5607: PetscErrorCode DMPlexGetGeomModelEdgeNodes(DM dm, PetscGeom body, PetscGeom edge, PetscGeom **nodes, PetscInt *numNodes) PeNS
5608: {
5609: PetscFunctionBeginHot;
5610: #ifdef PETSC_HAVE_EGADS
5611: PetscContainer modelObj;
5612: PetscBool islite = PETSC_FALSE;
5614: /* Determine which type of EGADS model is attached to the DM */
5615: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5616: if (!modelObj) {
5617: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5618: islite = PETSC_TRUE;
5619: }
5621: if (islite) {
5622: PetscCall(EGlite_getBodyTopos(body, edge, NODE, numNodes, nodes));
5623: } else {
5624: PetscCall(EG_getBodyTopos(body, edge, NODE, numNodes, nodes));
5625: }
5626: #endif
5627: PetscFunctionReturn(PETSC_SUCCESS);
5628: }
5630: /*@C
5631: DMPlexGetGeomID - Returns ID number of the entity in the geometric (CAD) model
5633: Collective
5635: Input Parameters:
5636: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5637: . body - PetscGeom body object containing the lower level entity the ID number is being requested.
5638: - topoObj - PetscGeom SHELL, FACE, LOOP, EDGE, or NODE object for which ID number is being requested.
5640: Output Parameter:
5641: . id - ID number of the entity
5643: Level: intermediate
5645: .seealso:
5646: @*/
5647: PetscErrorCode DMPlexGetGeomID(DM dm, PetscGeom body, PetscGeom topoObj, PetscInt *id) PeNS
5648: {
5649: PetscFunctionBeginHot;
5650: #ifdef PETSC_HAVE_EGADS
5651: PetscContainer modelObj;
5652: PetscBool islite = PETSC_FALSE;
5653: int topoID;
5655: /* Determine which type of EGADS model is attached to the DM */
5656: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5657: if (!modelObj) {
5658: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5659: islite = PETSC_TRUE;
5660: }
5662: // Get Topology Object's ID
5663: if (islite) {
5664: topoID = EGlite_indexBodyTopo(body, topoObj);
5665: } else {
5666: topoID = EG_indexBodyTopo(body, topoObj);
5667: }
5669: *id = topoID;
5670: #endif
5671: PetscFunctionReturn(PETSC_SUCCESS);
5672: }
5674: /*@C
5675: DMPlexGetGeomObject - Returns Geometry Object using the objects ID in the geometric (CAD) model
5677: Collective
5679: Input Parameters:
5680: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5681: . body - PetscGeom body object containing the lower level entity the referenced by the ID.
5682: . geomType - Keyword SHELL, FACE, LOOP, EDGE, or NODE of the geometry type for which ID number is being requested.
5683: - geomID - ID number of the geometry entity being requested.
5685: Output Parameter:
5686: . geomObj - Geometry Object referenced by the ID number requested.
5688: Level: intermediate
5690: .seealso:
5691: @*/
5692: PetscErrorCode DMPlexGetGeomObject(DM dm, PetscGeom body, PetscInt geomType, PetscInt geomID, PetscGeom *geomObj) PeNS
5693: {
5694: PetscFunctionBeginHot;
5695: #ifdef PETSC_HAVE_EGADS
5696: PetscContainer modelObj;
5697: PetscBool islite = PETSC_FALSE;
5699: /* Determine which type of EGADS model is attached to the DM */
5700: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5701: if (!modelObj) {
5702: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5703: islite = PETSC_TRUE;
5704: }
5706: // Get Topology Object's ID
5707: if (islite) {
5708: PetscCall(EGlite_objectBodyTopo(body, geomType, geomID, geomObj));
5709: } else {
5710: PetscCall(EG_objectBodyTopo(body, geomType, geomID, geomObj));
5711: }
5712: #endif
5713: PetscFunctionReturn(PETSC_SUCCESS);
5714: }
5716: /*@C
5717: DMPlexGetGeomFaceNumOfControlPoints - Returns the total number of Control Points (and associated Weights) defining a FACE of a Geometry
5719: Not collective
5721: Input Parameters:
5722: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5723: - face - PetscGeom FACE object
5725: Output Parameter:
5726: . numCntrlPnts - Number of Control Points (and Weights) defining the FACE
5728: Level: intermediate
5730: .seealso:
5731: @*/
5732: PetscErrorCode DMPlexGetGeomFaceNumOfControlPoints(DM dm, PetscGeom face, PetscInt *numCntrlPnts) PeNS
5733: {
5734: PetscFunctionBeginHot;
5735: #ifdef PETSC_HAVE_EGADS
5736: PetscContainer modelObj;
5737: PetscBool islite = PETSC_FALSE;
5738: PetscGeom geom, gRef;
5739: PetscGeom *lobjs;
5740: int Nl, oclass, mtype, goclass, gmtype;
5741: int *lsenses, *gpinfo;
5742: double *gprv;
5744: /* Determine which type of EGADS model is attached to the DM */
5745: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5746: if (!modelObj) {
5747: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5748: islite = PETSC_TRUE;
5749: }
5751: // Get Total Number of Control Points on FACE
5752: if (islite) {
5753: PetscCall(EGlite_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lsenses));
5754: PetscCall(EGlite_getGeometry(geom, &goclass, &gmtype, &gRef, &gpinfo, &gprv));
5755: } else {
5756: PetscCall(EG_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lsenses));
5757: PetscCall(EG_getGeometry(geom, &goclass, &gmtype, &gRef, &gpinfo, &gprv));
5758: }
5760: *numCntrlPnts = gpinfo[2] * gpinfo[5];
5761: #endif
5762: PetscFunctionReturn(PETSC_SUCCESS);
5763: }
5765: /*@C
5766: DMPlexGetGeomBodyMassProperties - Returns the Volume, Surface Area, Center of Gravity, and Inertia about the Body's Center of Gravity
5768: Not collective
5770: Input Parameters:
5771: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5772: - body - PetscGeom BODY object
5774: Output Parameters:
5775: + volume - Volume of the CAD Body attached to the DM Plex
5776: . surfArea - Surface Area of the CAD Body attached to the DM Plex
5777: . centerOfGravity - Array with the Center of Gravity coordinates of the CAD Body attached to the DM Plex [x, y, z]
5778: . COGszie - Size of centerOfGravity[] Array
5779: . inertiaMatrixCOG - Array containing the Inertia about the Body's Center of Gravity [Ixx, Ixy, Ixz, Iyx, Iyy, Iyz, Izx, Izy, Izz]
5780: - IMCOGsize - Size of inertiaMatrixCOG[] Array
5782: Level: intermediate
5784: .seealso:
5785: @*/
5786: PetscErrorCode DMPlexGetGeomBodyMassProperties(DM dm, PetscGeom body, PetscScalar *volume, PetscScalar *surfArea, PetscScalar **centerOfGravity, PetscInt *COGsize, PetscScalar **inertiaMatrixCOG, PetscInt *IMCOGsize) PeNS
5787: {
5788: PetscFunctionBeginHot;
5789: #ifdef PETSC_HAVE_EGADS
5790: PetscContainer modelObj;
5791: PetscBool islite = PETSC_FALSE;
5792: PetscScalar geomData[14];
5794: /* Determine which type of EGADS model is attached to the DM */
5795: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5796: if (!modelObj) {
5797: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5798: islite = PETSC_TRUE;
5799: PetscCheck(modelObj, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot provide geometric mass properties for geometries defined by EGADSlite (.egadslite)! Please use another geometry file format STEP, IGES, EGADS or BRep");
5800: }
5802: if (islite) {
5803: PetscCall(PetscPrintf(PETSC_COMM_SELF, " WARNING!! This functionality is not supported for EGADSlite files. \n"));
5804: PetscCall(PetscPrintf(PETSC_COMM_SELF, " All returned values are equal to 0 \n"));
5805: } else {
5806: PetscCall(EG_getMassProperties(body, geomData));
5807: }
5809: PetscCall(PetscMalloc2(3, centerOfGravity, 9, inertiaMatrixCOG));
5811: if (!islite) {
5812: *volume = geomData[0];
5813: *surfArea = geomData[1];
5814: for (int ii = 2; ii < 5; ++ii) (*centerOfGravity)[ii - 2] = geomData[ii];
5815: *COGsize = 3;
5816: for (int ii = 5; ii < 14; ++ii) (*inertiaMatrixCOG)[ii - 5] = geomData[ii];
5817: *IMCOGsize = 9;
5818: } else {
5819: *volume = 0.;
5820: *surfArea = 0.;
5821: for (int ii = 2; ii < 5; ++ii) (*centerOfGravity)[ii - 2] = 0.;
5822: *COGsize = 0;
5823: for (int ii = 5; ii < 14; ++ii) (*inertiaMatrixCOG)[ii - 5] = 0.;
5824: *IMCOGsize = 0;
5825: }
5826: #endif
5827: PetscFunctionReturn(PETSC_SUCCESS);
5828: }
5830: PetscErrorCode DMPlexRestoreGeomBodyMassProperties(DM dm, PetscGeom body, PetscScalar *volume, PetscScalar *surfArea, PetscScalar **centerOfGravity, PetscInt *COGsize, PetscScalar **inertiaMatrixCOG, PetscInt *IMCOGsize) PeNS
5831: {
5832: PetscFunctionBegin;
5833: PetscCall(PetscFree2(*centerOfGravity, *inertiaMatrixCOG));
5834: PetscFunctionReturn(PETSC_SUCCESS);
5835: }
5837: /*@C
5838: DMPlexFreeGeomObject - Frees PetscGeom Objects
5840: Not collective
5842: Input Parameters:
5843: + dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5844: - geomObj - PetscGeom object
5846: Level: intermediate
5848: .seealso:
5849: @*/
5850: PetscErrorCode DMPlexFreeGeomObject(DM dm, PetscGeom *geomObj) PeNS
5851: {
5852: PetscFunctionBeginHot;
5853: #ifdef PETSC_HAVE_EGADS
5854: PetscContainer modelObj;
5855: PetscBool islite = PETSC_FALSE;
5857: /* Determine which type of EGADS model is attached to the DM */
5858: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5859: if (!modelObj) {
5860: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5861: islite = PETSC_TRUE;
5862: }
5864: if (islite) {
5865: EGlite_free(geomObj);
5866: } else {
5867: EG_free(geomObj);
5868: }
5869: #endif
5870: PetscFunctionReturn(PETSC_SUCCESS);
5871: }
5873: /*@C
5874: DMPlexGetGeomCntrlPntAndWeightData - Gets Control Point and Associated Weight Data for the Geometry attached to the DMPlex
5876: Not collective
5878: Input Parameter:
5879: . dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5881: Output Parameters:
5882: + cpHashTable - Hash Table containing the relationship between FACE ID and Control Point IDs.
5883: . cpCoordDataLength - Length of cpCoordData Array.
5884: . cpCoordData - Array holding the Geometry Control Point Coordinate Data.
5885: . maxNumEquiv - Maximum Number of Equivalent Control Points (Control Points with the same coordinates but different IDs).
5886: . cpEquiv - Matrix with a size(Number of Control Points, Number or Control Points) which stores a value of 1.0 in locations where Control Points with different IDS (row or column) have the same coordinates
5887: . wHashTable - Hash Table containing the relationship between FACE ID and Control Point Weight.
5888: . wDataLength - Length of wData Array.
5889: - wData - Array holding the Weight for an associated Geometry Control Point.
5891: Note:
5892: Must Call DMPLexGeomDataAndGrads() before calling this function.
5894: Level: intermediate
5896: .seealso:
5897: @*/
5898: PetscErrorCode DMPlexGetGeomCntrlPntAndWeightData(DM dm, PetscHMapI *cpHashTable, PetscInt *cpCoordDataLength, PetscScalar **cpCoordData, PetscInt *maxNumEquiv, Mat *cpEquiv, PetscHMapI *wHashTable, PetscInt *wDataLength, PetscScalar **wData) PeNS
5899: {
5900: PetscContainer modelObj, cpHashTableObj, wHashTableObj, cpCoordDataLengthObj, wDataLengthObj, maxNumRelateObj;
5901: Vec cntrlPtCoordsVec, cntrlPtWeightsVec;
5902: PetscInt *cpCoordDataLengthPtr, *wDataLengthPtr, *maxNumEquivPtr;
5903: PetscHMapI cpHashTableTemp, wHashTableTemp;
5905: PetscFunctionBeginHot;
5906: /* Determine which type of EGADS model is attached to the DM */
5907: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5908: if (!modelObj) PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5910: if (!modelObj) PetscFunctionReturn(PETSC_SUCCESS);
5912: // Look to see if DM has Container for Geometry Control Point Data
5913: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Hash Table", (PetscObject *)&cpHashTableObj));
5914: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Coordinates", (PetscObject *)&cntrlPtCoordsVec));
5915: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Coordinate Data Length", (PetscObject *)&cpCoordDataLengthObj));
5916: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weights Hash Table", (PetscObject *)&wHashTableObj));
5917: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data", (PetscObject *)&cntrlPtWeightsVec));
5918: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data Length", (PetscObject *)&wDataLengthObj));
5919: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Equivalency Matrix", (PetscObject *)cpEquiv));
5920: PetscCall(PetscObjectQuery((PetscObject)dm, "Maximum Number Control Point Equivalency", (PetscObject *)&maxNumRelateObj));
5922: // Get attached EGADS model Control Point and Weights Hash Tables and Data Arrays (pointer)
5923: PetscCall(PetscContainerGetPointer(cpHashTableObj, &cpHashTableTemp));
5924: PetscCall(PetscContainerGetPointer(cpCoordDataLengthObj, &cpCoordDataLengthPtr));
5925: PetscCall(PetscContainerGetPointer(wHashTableObj, &wHashTableTemp));
5926: PetscCall(PetscContainerGetPointer(wDataLengthObj, &wDataLengthPtr));
5927: PetscCall(PetscContainerGetPointer(maxNumRelateObj, &maxNumEquivPtr));
5929: *cpCoordDataLength = *cpCoordDataLengthPtr;
5930: *wDataLength = *wDataLengthPtr;
5931: *maxNumEquiv = *maxNumEquivPtr;
5932: *cpHashTable = cpHashTableTemp;
5933: *wHashTable = wHashTableTemp;
5934: PetscCall(VecGetArrayWrite(cntrlPtCoordsVec, cpCoordData));
5935: PetscCall(VecGetArrayWrite(cntrlPtWeightsVec, wData));
5936: PetscFunctionReturn(PETSC_SUCCESS);
5937: }
5939: PetscErrorCode DMPlexRestoreGeomCntrlPntAndWeightData(DM dm, PetscHMapI *cpHashTable, PetscInt *cpCoordDataLength, PetscScalar **cpCoordData, PetscInt *maxNumEquiv, Mat *cpEquiv, PetscHMapI *wHashTable, PetscInt *wDataLength, PetscScalar **wData)
5940: {
5941: Vec cntrlPtCoordsVec, cntrlPtWeightsVec;
5943: PetscFunctionBeginHot;
5944: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Coordinates", (PetscObject *)&cntrlPtCoordsVec));
5945: PetscCall(VecRestoreArrayWrite(cntrlPtCoordsVec, cpCoordData));
5946: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data", (PetscObject *)&cntrlPtWeightsVec));
5947: PetscCall(VecRestoreArrayWrite(cntrlPtWeightsVec, wData));
5948: PetscFunctionReturn(PETSC_SUCCESS);
5949: }
5951: /*@C
5952: DMPlexGetGeomGradData - Gets Point, Surface and Volume Gradients with respect to changes in Control Points and their associated Weights for the Geometry attached to the DMPlex .
5954: Not collective
5956: Input Parameter:
5957: . dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
5959: Output Parameters:
5960: + cpSurfGradHashTable - Hash Table Relating the Control Point ID to the the Row in the cpSurfGrad Matrix
5961: . cpSurfGrad - Matrix containing the Surface Gradient with respect to the Control Point Data. Data is ranged where the Row corresponds to Control Point ID and the Columns are associated with the Geometric FACE.
5962: . cpArraySize - The size of arrays gradSACP and gradVolCP and is equal to 3 * total number of Control Points in the Geometry
5963: . gradSACP - Array containing the Surface Area Gradient with respect to Control Point Data. Data is arranged by Control Point ID * 3 where 3 is for the coordinate dimension.
5964: . gradVolCP - Array containing the Volume Gradient with respect to Control Point Data. Data is arranged by Control Point ID * 3 where 3 is for the coordinate dimension.
5965: . wArraySize - The size of arrayws gradSAW and gradVolW and is equal to the total number of Control Points in the Geometry.
5966: . gradSAW - Array containing the Surface Area Gradient with respect to Control Point Weight. Data is arranged by Control Point ID.
5967: - gradVolW - Array containing the Volume Gradient with respect to Control Point Weight. Data is arranged by Control Point ID.
5969: Notes:
5970: Must Call DMPLexGeomDataAndGrads() before calling this function.
5972: gradVolCP and gradVolW are only available when DMPlexGeomDataAndGrads() is called with fullGeomGrad = PETSC_TRUE.
5974: Level: intermediate
5976: .seealso: DMPlexGeomDataAndGrads
5977: @*/
5978: PetscErrorCode DMPlexGetGeomGradData(DM dm, PetscHMapI *cpSurfGradHashTable, Mat *cpSurfGrad, PetscInt *cpArraySize, PetscScalar **gradSACP, PetscScalar **gradVolCP, PetscInt *wArraySize, PetscScalar **gradSAW, PetscScalar **gradVolW)
5979: {
5980: PetscContainer modelObj, cpSurfGradHashTableObj, cpArraySizeObj, wArraySizeObj;
5981: Vec gradSACPVec, gradVolCPVec, gradSAWVec, gradVolWVec;
5982: PetscInt *cpArraySizePtr, *wArraySizePtr;
5983: PetscHMapI cpSurfGradHashTableTemp;
5985: PetscFunctionBeginHot;
5986: /* Determine which type of EGADS model is attached to the DM */
5987: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
5988: if (!modelObj) PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
5990: if (!modelObj) PetscFunctionReturn(PETSC_SUCCESS);
5992: // Look to see if DM has Container for Geometry Control Point Data
5993: PetscCall(PetscObjectQuery((PetscObject)dm, "Surface Gradient Hash Table", (PetscObject *)&cpSurfGradHashTableObj));
5994: PetscCall(PetscObjectQuery((PetscObject)dm, "Surface Gradient Matrix", (PetscObject *)cpSurfGrad));
5995: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Coordinate Data Length", (PetscObject *)&cpArraySizeObj));
5996: PetscCall(PetscObjectQuery((PetscObject)dm, "Surface Area Control Point Gradient", (PetscObject *)&gradSACPVec));
5997: PetscCall(PetscObjectQuery((PetscObject)dm, "Volume Control Point Gradient", (PetscObject *)&gradVolCPVec));
5998: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data Length", (PetscObject *)&wArraySizeObj));
5999: PetscCall(PetscObjectQuery((PetscObject)dm, "Surface Area Weights Gradient", (PetscObject *)&gradSAWVec));
6000: PetscCall(PetscObjectQuery((PetscObject)dm, "Volume Weights Gradient", (PetscObject *)&gradVolWVec));
6002: // Get attached EGADS model Control Point and Weights Hash Tables and Data Arrays (pointer)
6003: if (cpSurfGradHashTableObj) {
6004: PetscCall(PetscContainerGetPointer(cpSurfGradHashTableObj, &cpSurfGradHashTableTemp));
6005: *cpSurfGradHashTable = cpSurfGradHashTableTemp;
6006: }
6008: if (cpArraySizeObj) {
6009: PetscCall(PetscContainerGetPointer(cpArraySizeObj, &cpArraySizePtr));
6010: *cpArraySize = *cpArraySizePtr;
6011: }
6013: if (gradSACPVec) PetscCall(VecGetArrayWrite(gradSACPVec, gradSACP));
6014: if (gradVolCPVec) PetscCall(VecGetArrayWrite(gradVolCPVec, gradVolCP));
6015: if (gradSAWVec) PetscCall(VecGetArrayWrite(gradSAWVec, gradSAW));
6016: if (gradVolWVec) PetscCall(VecGetArrayWrite(gradVolWVec, gradVolW));
6018: if (wArraySizeObj) {
6019: PetscCall(PetscContainerGetPointer(wArraySizeObj, &wArraySizePtr));
6020: *wArraySize = *wArraySizePtr;
6021: }
6022: PetscFunctionReturn(PETSC_SUCCESS);
6023: }
6025: PetscErrorCode DMPlexRestoreGeomGradData(DM dm, PetscHMapI *cpSurfGradHashTable, Mat *cpSurfGrad, PetscInt *cpArraySize, PetscScalar **gradSACP, PetscScalar **gradVolCP, PetscInt *wArraySize, PetscScalar **gradSAW, PetscScalar **gradVolW)
6026: {
6027: Vec gradSACPVec, gradVolCPVec, gradSAWVec, gradVolWVec;
6029: PetscFunctionBegin;
6030: PetscCall(PetscObjectQuery((PetscObject)dm, "Surface Area Control Point Gradient", (PetscObject *)&gradSACPVec));
6031: PetscCall(PetscObjectQuery((PetscObject)dm, "Volume Control Point Gradient", (PetscObject *)&gradVolCPVec));
6032: PetscCall(PetscObjectQuery((PetscObject)dm, "Surface Area Weights Gradient", (PetscObject *)&gradSAWVec));
6033: PetscCall(PetscObjectQuery((PetscObject)dm, "Volume Weights Gradient", (PetscObject *)&gradVolWVec));
6035: if (gradSACPVec) PetscCall(VecRestoreArrayWrite(gradSACPVec, gradSACP));
6036: if (gradVolCPVec) PetscCall(VecRestoreArrayWrite(gradVolCPVec, gradVolCP));
6037: if (gradSAWVec) PetscCall(VecRestoreArrayWrite(gradSAWVec, gradSAW));
6038: if (gradVolWVec) PetscCall(VecRestoreArrayWrite(gradVolWVec, gradVolW));
6039: PetscFunctionReturn(PETSC_SUCCESS);
6040: }
6042: /*@C
6043: DMPlexGetGeomCntrlPntMaps - Gets arrays which maps Control Point IDs to their associated Geometry FACE, EDGE, and VERTEX.
6045: Not collective
6047: Input Parameter:
6048: . dm - The DMPlex object with an attached PetscContainer storing a CAD Geometry object
6050: Output Parameters:
6051: + numCntrlPnts - Number of Control Points defining the Geometry attached to the DMPlex
6052: . cntrlPntFaceMap - Array containing the FACE ID for the Control Point. Array index corresponds to Control Point ID.
6053: . cntrlPntWeightFaceMap - Array containing the FACE ID for the Control Point Weight. Array index corresponds to Control Point ID.
6054: . cntrlPntEdgeMap - Array containing the EDGE ID for the Control Point. Array index corresponds to Control Point ID.
6055: . cntrlPntWeightEdgeMap - Array containing the EDGE ID for the Control Point Weight. Array index corresponds to Control Point ID.
6056: . cntrlPntVertexMap - Array containing the VERTEX ID for the Control Point. Array index corresponds to Control Point ID.
6057: - cntrlPntWeightVertexMap - Array containing the VERTEX ID for the Control Point Weight. Array index corresponds to Control Point ID.
6059: Note:
6060: Arrays are initialized to -1. Array elements with a -1 value indicates that the Control Point or Control Point Weight not associated with the referenced Geometric entity in the array name.
6062: Level: intermediate
6064: .seealso: DMPlexGeomDataAndGrads
6065: @*/
6066: PetscErrorCode DMPlexGetGeomCntrlPntMaps(DM dm, PetscInt *numCntrlPnts, PetscInt **cntrlPntFaceMap, PetscInt **cntrlPntWeightFaceMap, PetscInt **cntrlPntEdgeMap, PetscInt **cntrlPntWeightEdgeMap, PetscInt **cntrlPntVertexMap, PetscInt **cntrlPntWeightVertexMap)
6067: {
6068: PetscFunctionBeginHot;
6069: #ifdef PETSC_HAVE_EGADS
6070: PetscContainer modelObj, numCntrlPntsObj, cntrlPntFaceMapObj, cntrlPntWeightFaceMapObj, cntrlPntEdgeMapObj, cntrlPntWeightEdgeMapObj, cntrlPntVertexMapObj, cntrlPntWeightVertexMapObj;
6071: PetscInt *numCntrlPntsPtr, *cntrlPntFaceMapPtr, *cntrlPntWeightFaceMapPtr, *cntrlPntEdgeMapPtr, *cntrlPntWeightEdgeMapPtr, *cntrlPntVertexMapPtr, *cntrlPntWeightVertexMapPtr;
6073: /* Determine which type of EGADS model is attached to the DM */
6074: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
6075: if (!modelObj) PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", (PetscObject *)&modelObj));
6077: if (!modelObj) PetscFunctionReturn(PETSC_SUCCESS);
6079: // Look to see if DM has Container for Geometry Control Point Data
6080: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight Data Length", (PetscObject *)&numCntrlPntsObj));
6081: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point - Face Map", (PetscObject *)&cntrlPntFaceMapObj));
6082: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight - Face Map", (PetscObject *)&cntrlPntWeightFaceMapObj));
6083: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point - Edge Map", (PetscObject *)&cntrlPntEdgeMapObj));
6084: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight - Edge Map", (PetscObject *)&cntrlPntWeightEdgeMapObj));
6085: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point - Vertex Map", (PetscObject *)&cntrlPntVertexMapObj));
6086: PetscCall(PetscObjectQuery((PetscObject)dm, "Control Point Weight - Vertex Map", (PetscObject *)&cntrlPntWeightVertexMapObj));
6088: // Get attached EGADS model Control Point and Weights Hash Tables and Data Arrays (pointer)
6089: if (numCntrlPntsObj) {
6090: PetscCall(PetscContainerGetPointer(numCntrlPntsObj, &numCntrlPntsPtr));
6091: *numCntrlPnts = *numCntrlPntsPtr;
6092: }
6094: if (cntrlPntFaceMapObj) {
6095: PetscCall(PetscContainerGetPointer(cntrlPntFaceMapObj, &cntrlPntFaceMapPtr));
6096: *cntrlPntFaceMap = cntrlPntFaceMapPtr;
6097: }
6099: if (cntrlPntWeightFaceMapObj) {
6100: PetscCall(PetscContainerGetPointer(cntrlPntWeightFaceMapObj, &cntrlPntWeightFaceMapPtr));
6101: *cntrlPntWeightFaceMap = cntrlPntWeightFaceMapPtr;
6102: }
6104: if (cntrlPntEdgeMapObj) {
6105: PetscCall(PetscContainerGetPointer(cntrlPntEdgeMapObj, &cntrlPntEdgeMapPtr));
6106: *cntrlPntEdgeMap = cntrlPntEdgeMapPtr;
6107: }
6109: if (cntrlPntWeightEdgeMapObj) {
6110: PetscCall(PetscContainerGetPointer(cntrlPntWeightEdgeMapObj, &cntrlPntWeightEdgeMapPtr));
6111: *cntrlPntWeightEdgeMap = cntrlPntWeightEdgeMapPtr;
6112: }
6114: if (cntrlPntVertexMapObj) {
6115: PetscCall(PetscContainerGetPointer(cntrlPntVertexMapObj, &cntrlPntVertexMapPtr));
6116: *cntrlPntVertexMap = cntrlPntVertexMapPtr;
6117: }
6119: if (cntrlPntWeightVertexMapObj) {
6120: PetscCall(PetscContainerGetPointer(cntrlPntWeightVertexMapObj, &cntrlPntWeightVertexMapPtr));
6121: *cntrlPntWeightVertexMap = cntrlPntWeightVertexMapPtr;
6122: }
6124: #endif
6125: PetscFunctionReturn(PETSC_SUCCESS);
6126: }
6128: #endif