1a8ededdfSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2a8ededdfSMatthew G. Knepley 3a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS 4a8ededdfSMatthew G. Knepley #include <egads.h> 5a8ededdfSMatthew G. Knepley #endif 6a8ededdfSMatthew G. Knepley 7a8ededdfSMatthew G. Knepley /* We need to understand how to natively parse STEP files. There seems to be only one open source implementation of 8a8ededdfSMatthew G. Knepley the STEP parser contained in the OpenCASCADE package. It is enough to make a strong man weep: 9a8ededdfSMatthew G. Knepley 10a8ededdfSMatthew G. Knepley https://github.com/tpaviot/oce/tree/master/src/STEPControl 11a8ededdfSMatthew G. Knepley 12a8ededdfSMatthew G. Knepley The STEP, and inner EXPRESS, formats are ISO standards, so they are documented 13a8ededdfSMatthew G. Knepley 14a8ededdfSMatthew G. Knepley https://stackoverflow.com/questions/26774037/documentation-or-specification-for-step-and-stp-files 15a8ededdfSMatthew G. Knepley http://stepmod.sourceforge.net/express_model_spec/ 16a8ededdfSMatthew G. Knepley 17a8ededdfSMatthew G. Knepley but again it seems that there has been a deliberate effort at obfuscation, probably to raise the bar for entrants. 18a8ededdfSMatthew G. Knepley */ 19a8ededdfSMatthew G. Knepley 20a8ededdfSMatthew G. Knepley 21a8ededdfSMatthew G. Knepley /*@ 22a8ededdfSMatthew 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. 23a8ededdfSMatthew G. Knepley 24a8ededdfSMatthew G. Knepley Not collective 25a8ededdfSMatthew G. Knepley 26a8ededdfSMatthew G. Knepley Input Parameters: 27a8ededdfSMatthew G. Knepley + dm - The DMPlex object 28a8ededdfSMatthew G. Knepley . p - The mesh point 29a8ededdfSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point 30a8ededdfSMatthew G. Knepley 31a8ededdfSMatthew G. Knepley Output Parameter: 32a8ededdfSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point 33a8ededdfSMatthew G. Knepley 34a8ededdfSMatthew G. Knepley Note: Returns the original coordinates if no geometry model is found. Right now the only supported geometry model is EGADS. 35a8ededdfSMatthew G. Knepley 36a8ededdfSMatthew G. Knepley Level: intermediate 37a8ededdfSMatthew G. Knepley 38a8ededdfSMatthew G. Knepley .seealso: DMRefine(), DMPlexCreate(), DMPlexSetRefinementUniform() 39a8ededdfSMatthew G. Knepley @*/ 40a8ededdfSMatthew G. Knepley PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, const PetscScalar mcoords[], PetscScalar gcoords[]) 41a8ededdfSMatthew G. Knepley { 42a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS 43*f0fcf0adSMatthew G. Knepley DM cdm; 44a8ededdfSMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel; 45a8ededdfSMatthew G. Knepley PetscContainer modelObj; 46a8ededdfSMatthew G. Knepley PetscInt bodyID, faceID, edgeID; 47a8ededdfSMatthew G. Knepley ego *bodies; 48*f0fcf0adSMatthew G. Knepley ego model, geom, body, obj; 49*f0fcf0adSMatthew G. Knepley /* result has to hold derviatives, along with the value */ 50*f0fcf0adSMatthew G. Knepley double params[3], result[18], paramsV[16*3], resultV[16*3]; 51a8ededdfSMatthew G. Knepley int Nb, oclass, mtype, *senses; 52*f0fcf0adSMatthew G. Knepley Vec coordinatesLocal; 53*f0fcf0adSMatthew G. Knepley PetscScalar *coords = NULL; 54*f0fcf0adSMatthew G. Knepley PetscInt Nv, v, Np = 0, pm; 55a8ededdfSMatthew G. Knepley #endif 56*f0fcf0adSMatthew G. Knepley PetscInt dE, d; 57a8ededdfSMatthew G. Knepley PetscErrorCode ierr; 58a8ededdfSMatthew G. Knepley 59a8ededdfSMatthew G. Knepley PetscFunctionBegin; 60*f0fcf0adSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr); 61a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS 62a8ededdfSMatthew G. Knepley ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr); 63a8ededdfSMatthew G. Knepley ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr); 64a8ededdfSMatthew G. Knepley ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr); 65a8ededdfSMatthew G. Knepley if (!bodyLabel || !faceLabel || !edgeLabel) { 66*f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 67a8ededdfSMatthew G. Knepley PetscFunctionReturn(0); 68a8ededdfSMatthew G. Knepley } 69*f0fcf0adSMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 70*f0fcf0adSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinatesLocal);CHKERRQ(ierr); 71a8ededdfSMatthew G. Knepley ierr = DMLabelGetValue(bodyLabel, p, &bodyID);CHKERRQ(ierr); 72a8ededdfSMatthew G. Knepley ierr = DMLabelGetValue(faceLabel, p, &faceID);CHKERRQ(ierr); 73a8ededdfSMatthew G. Knepley ierr = DMLabelGetValue(edgeLabel, p, &edgeID);CHKERRQ(ierr); 74a8ededdfSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr); 75a8ededdfSMatthew G. Knepley ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr); 76a8ededdfSMatthew G. Knepley ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr); 77*f0fcf0adSMatthew G. Knepley if (bodyID < 0 || bodyID >= Nb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %D is not in [0, %d)", bodyID, Nb); 78a8ededdfSMatthew G. Knepley body = bodies[bodyID]; 79*f0fcf0adSMatthew G. Knepley 80*f0fcf0adSMatthew G. Knepley if (edgeID >= 0) {ierr = EG_objectBodyTopo(body, EDGE, edgeID, &obj);CHKERRQ(ierr); Np = 1;} 81*f0fcf0adSMatthew G. Knepley else if (faceID >= 0) {ierr = EG_objectBodyTopo(body, FACE, faceID, &obj);CHKERRQ(ierr); Np = 2;} 82*f0fcf0adSMatthew G. Knepley else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D is not in edge or face label for EGADS", p); 83*f0fcf0adSMatthew G. Knepley /* Calculate parameters (t or u,v) for vertices */ 84*f0fcf0adSMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr); 85*f0fcf0adSMatthew G. Knepley Nv /= dE; 86*f0fcf0adSMatthew G. Knepley if (Nv == 1) { 87*f0fcf0adSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr); 88*f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 89*f0fcf0adSMatthew G. Knepley PetscFunctionReturn(0); 90a8ededdfSMatthew G. Knepley } 91*f0fcf0adSMatthew G. Knepley if (Nv > 16) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %D coordinates associated to point %D", Nv, p); 92*f0fcf0adSMatthew G. Knepley for (v = 0; v < Nv; ++v) {ierr = EG_invEvaluate(obj, &coords[v*dE], ¶msV[v*3], &resultV[v*3]);CHKERRQ(ierr);} 93*f0fcf0adSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr); 94*f0fcf0adSMatthew G. Knepley /* Calculate parameters (t or u,v) for new vertex at edge midpoint */ 95*f0fcf0adSMatthew G. Knepley for (pm = 0; pm < Np; ++pm) { 96*f0fcf0adSMatthew G. Knepley params[pm] = 0.; 97*f0fcf0adSMatthew G. Knepley for (v = 0; v < Nv; ++v) {params[pm] += paramsV[v*3+pm];} 98*f0fcf0adSMatthew G. Knepley params[pm] /= Nv; 99*f0fcf0adSMatthew G. Knepley } 100*f0fcf0adSMatthew G. Knepley /* TODO Check 101*f0fcf0adSMatthew G. Knepley double range[4]; // [umin, umax, vmin, vmax] 102*f0fcf0adSMatthew G. Knepley int peri; 103*f0fcf0adSMatthew G. Knepley ierr = EG_getRange(face, range, &peri); CHKERRQ(ierr); 104*f0fcf0adSMatthew G. Knepley if ((paramsNew[0] < range[0]) || (paramsNew[0] > range[1]) || (paramsNew[1] < range[2]) || (paramsNew[1] > range[3])) SETERRQ(); 105*f0fcf0adSMatthew G. Knepley */ 106*f0fcf0adSMatthew G. Knepley /* Put coordinates for new vertex in result[] */ 107*f0fcf0adSMatthew G. Knepley ierr = EG_evaluate(obj, params, result);CHKERRQ(ierr); 108*f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = result[d]; 109a8ededdfSMatthew G. Knepley #else 110*f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 111a8ededdfSMatthew G. Knepley #endif 112a8ededdfSMatthew G. Knepley PetscFunctionReturn(0); 113a8ededdfSMatthew G. Knepley } 114