xref: /petsc/src/dm/impls/plex/plexegads.c (revision a8ededdf47aa46be0940f1dfc7fe543fa88f0313)
1*a8ededdfSMatthew G. Knepley #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*a8ededdfSMatthew G. Knepley 
3*a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
4*a8ededdfSMatthew G. Knepley #include <egads.h>
5*a8ededdfSMatthew G. Knepley #endif
6*a8ededdfSMatthew G. Knepley 
7*a8ededdfSMatthew G. Knepley /* We need to understand how to natively parse STEP files. There seems to be only one open source implementation of
8*a8ededdfSMatthew G. Knepley    the STEP parser contained in the OpenCASCADE package. It is enough to make a strong man weep:
9*a8ededdfSMatthew G. Knepley 
10*a8ededdfSMatthew G. Knepley      https://github.com/tpaviot/oce/tree/master/src/STEPControl
11*a8ededdfSMatthew G. Knepley 
12*a8ededdfSMatthew G. Knepley    The STEP, and inner EXPRESS, formats are ISO standards, so they are documented
13*a8ededdfSMatthew G. Knepley 
14*a8ededdfSMatthew G. Knepley      https://stackoverflow.com/questions/26774037/documentation-or-specification-for-step-and-stp-files
15*a8ededdfSMatthew G. Knepley      http://stepmod.sourceforge.net/express_model_spec/
16*a8ededdfSMatthew G. Knepley 
17*a8ededdfSMatthew G. Knepley    but again it seems that there has been a deliberate effort at obfuscation, probably to raise the bar for entrants.
18*a8ededdfSMatthew G. Knepley */
19*a8ededdfSMatthew G. Knepley 
20*a8ededdfSMatthew G. Knepley 
21*a8ededdfSMatthew G. Knepley /*@
22*a8ededdfSMatthew G. Knepley   DMPlexSnapToGeomModel - Given a coordinate point 'mcoords' on the mesh point 'p', return the closest coordinate point 'gcoords' on the geometry model associated with that point.
23*a8ededdfSMatthew G. Knepley 
24*a8ededdfSMatthew G. Knepley   Not collective
25*a8ededdfSMatthew G. Knepley 
26*a8ededdfSMatthew G. Knepley   Input Parameters:
27*a8ededdfSMatthew G. Knepley + dm      - The DMPlex object
28*a8ededdfSMatthew G. Knepley . p       - The mesh point
29*a8ededdfSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point
30*a8ededdfSMatthew G. Knepley 
31*a8ededdfSMatthew G. Knepley   Output Parameter:
32*a8ededdfSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point
33*a8ededdfSMatthew G. Knepley 
34*a8ededdfSMatthew G. Knepley   Note: Returns the original coordinates if no geometry model is found. Right now the only supported geometry model is EGADS.
35*a8ededdfSMatthew G. Knepley 
36*a8ededdfSMatthew G. Knepley   Level: intermediate
37*a8ededdfSMatthew G. Knepley 
38*a8ededdfSMatthew G. Knepley .seealso: DMRefine(), DMPlexCreate(), DMPlexSetRefinementUniform()
39*a8ededdfSMatthew G. Knepley @*/
40*a8ededdfSMatthew G. Knepley PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, const PetscScalar mcoords[], PetscScalar gcoords[])
41*a8ededdfSMatthew G. Knepley {
42*a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
43*a8ededdfSMatthew G. Knepley   DMLabel        bodyLabel, faceLabel, edgeLabel;
44*a8ededdfSMatthew G. Knepley   PetscContainer modelObj;
45*a8ededdfSMatthew G. Knepley   PetscInt       bodyID, faceID, edgeID;
46*a8ededdfSMatthew G. Knepley   ego           *bodies;
47*a8ededdfSMatthew G. Knepley   ego            model, geom, body, face, edge;
48*a8ededdfSMatthew G. Knepley   double         point[3], params[3], result[3];
49*a8ededdfSMatthew G. Knepley   int            Nb, oclass, mtype, *senses;
50*a8ededdfSMatthew G. Knepley #endif
51*a8ededdfSMatthew G. Knepley   PetscInt       cdim, d;
52*a8ededdfSMatthew G. Knepley   PetscErrorCode ierr;
53*a8ededdfSMatthew G. Knepley 
54*a8ededdfSMatthew G. Knepley   PetscFunctionBegin;
55*a8ededdfSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
56*a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
57*a8ededdfSMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
58*a8ededdfSMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
59*a8ededdfSMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
60*a8ededdfSMatthew G. Knepley   if (!bodyLabel || !faceLabel || !edgeLabel) {
61*a8ededdfSMatthew G. Knepley     for (d = 0; d < cdim; ++d) gcoords[d] = mcoords[d];
62*a8ededdfSMatthew G. Knepley     PetscFunctionReturn(0);
63*a8ededdfSMatthew G. Knepley   }
64*a8ededdfSMatthew G. Knepley   ierr = DMLabelGetValue(bodyLabel, p, &bodyID);CHKERRQ(ierr);
65*a8ededdfSMatthew G. Knepley   ierr = DMLabelGetValue(faceLabel, p, &faceID);CHKERRQ(ierr);
66*a8ededdfSMatthew G. Knepley   ierr = DMLabelGetValue(edgeLabel, p, &edgeID);CHKERRQ(ierr);
67*a8ededdfSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
68*a8ededdfSMatthew G. Knepley   ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
69*a8ededdfSMatthew G. Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
70*a8ededdfSMatthew G. Knepley   if (bodyID >= Nb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %D is not in [0, %d)", bodyID, Nb);
71*a8ededdfSMatthew G. Knepley   body = bodies[bodyID];
72*a8ededdfSMatthew G. Knepley   for (d = 0; d < cdim; ++d) point[d] = mcoords[d];
73*a8ededdfSMatthew G. Knepley   if (edgeID >= 0) {
74*a8ededdfSMatthew G. Knepley     ierr = EG_objectBodyTopo(body, EDGE, edgeID, &edge);CHKERRQ(ierr);
75*a8ededdfSMatthew G. Knepley     ierr = EG_invEvaluate(edge, point, params, result);
76*a8ededdfSMatthew G. Knepley   } else {
77*a8ededdfSMatthew G. Knepley     ierr = EG_objectBodyTopo(body, FACE, faceID, &face);CHKERRQ(ierr);
78*a8ededdfSMatthew G. Knepley     ierr = EG_invEvaluate(face, point, params, result);
79*a8ededdfSMatthew G. Knepley   }
80*a8ededdfSMatthew G. Knepley   for (d = 0; d < cdim; ++d) gcoords[d] = result[d];
81*a8ededdfSMatthew G. Knepley #else
82*a8ededdfSMatthew G. Knepley   for (d = 0; d < cdim; ++d) gcoords[d] = mcoords[d];
83*a8ededdfSMatthew G. Knepley #endif
84*a8ededdfSMatthew G. Knepley   PetscFunctionReturn(0);
85*a8ededdfSMatthew G. Knepley }
86