17625649eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
27625649eSMatthew G. Knepley
37625649eSMatthew G. Knepley PetscFunctionList DMGeomModelList = NULL;
47625649eSMatthew G. Knepley PetscBool DMGeomModelRegisterAllCalled = PETSC_FALSE;
57625649eSMatthew G. Knepley
67625649eSMatthew G. Knepley #if defined(PETSC_HAVE_EGADS)
7d6acfc2dSPierre Jolivet PETSC_INTERN PetscErrorCode DMSnapToGeomModel_EGADS(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
87625649eSMatthew G. Knepley #endif
97625649eSMatthew G. Knepley
DMSnapToGeomModelBall(DM dm,PetscInt p,PetscInt dE,const PetscScalar mcoords[],PetscScalar gcoords[])107625649eSMatthew G. Knepley static PetscErrorCode DMSnapToGeomModelBall(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
117625649eSMatthew G. Knepley {
127625649eSMatthew G. Knepley PetscInt val;
137625649eSMatthew G. Knepley
147625649eSMatthew G. Knepley PetscFunctionBeginUser;
157625649eSMatthew G. Knepley PetscCall(DMGetLabelValue(dm, "marker", p, &val));
167625649eSMatthew G. Knepley if (val >= 0) {
177625649eSMatthew G. Knepley PetscReal norm = 0.;
187625649eSMatthew G. Knepley
197625649eSMatthew G. Knepley for (PetscInt d = 0; d < dE; ++d) norm += PetscSqr(PetscRealPart(mcoords[d]));
207625649eSMatthew G. Knepley norm = PetscSqrtReal(norm);
217625649eSMatthew G. Knepley for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d] / norm;
227625649eSMatthew G. Knepley } else {
237625649eSMatthew G. Knepley for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
247625649eSMatthew G. Knepley }
257625649eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
267625649eSMatthew G. Knepley }
277625649eSMatthew G. Knepley
DMSnapToGeomModelCylinder(DM dm,PetscInt p,PetscInt dE,const PetscScalar mcoords[],PetscScalar gcoords[])287625649eSMatthew G. Knepley static PetscErrorCode DMSnapToGeomModelCylinder(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
297625649eSMatthew G. Knepley {
307625649eSMatthew G. Knepley PetscReal gmin[3], gmax[3];
317625649eSMatthew G. Knepley PetscInt val;
327625649eSMatthew G. Knepley
337625649eSMatthew G. Knepley PetscFunctionBeginUser;
347625649eSMatthew G. Knepley PetscCall(DMGetBoundingBox(dm, gmin, gmax));
357625649eSMatthew G. Knepley PetscCall(DMGetLabelValue(dm, "generatrix", p, &val));
367625649eSMatthew G. Knepley if (val >= 0) {
377625649eSMatthew G. Knepley PetscReal norm = 0.;
387625649eSMatthew G. Knepley
397625649eSMatthew G. Knepley for (PetscInt d = 0; d < dE - 1; ++d) norm += PetscSqr(PetscRealPart(mcoords[d]));
407625649eSMatthew G. Knepley norm = PetscSqrtReal(norm);
41883424caSPierre Jolivet for (PetscInt d = 0; d < dE - 1; ++d) gcoords[d] = mcoords[d] * gmax[0] / norm;
427625649eSMatthew G. Knepley gcoords[dE - 1] = mcoords[dE - 1];
437625649eSMatthew G. Knepley } else {
447625649eSMatthew G. Knepley for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
457625649eSMatthew G. Knepley }
467625649eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
477625649eSMatthew G. Knepley }
487625649eSMatthew G. Knepley
497625649eSMatthew G. Knepley /*@C
507625649eSMatthew G. Knepley DMGeomModelRegisterAll - Registers all of the geometry model methods in the `DM` package.
517625649eSMatthew G. Knepley
527625649eSMatthew G. Knepley Not Collective
537625649eSMatthew G. Knepley
547625649eSMatthew G. Knepley Level: advanced
557625649eSMatthew G. Knepley
567625649eSMatthew G. Knepley .seealso: `DM`, `DMGeomModelRegisterDestroy()`
577625649eSMatthew G. Knepley @*/
DMGeomModelRegisterAll(void)587625649eSMatthew G. Knepley PetscErrorCode DMGeomModelRegisterAll(void)
597625649eSMatthew G. Knepley {
607625649eSMatthew G. Knepley PetscFunctionBegin;
617625649eSMatthew G. Knepley if (DMGeomModelRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
627625649eSMatthew G. Knepley DMGeomModelRegisterAllCalled = PETSC_TRUE;
637625649eSMatthew G. Knepley PetscCall(DMGeomModelRegister("ball", DMSnapToGeomModelBall));
647625649eSMatthew G. Knepley PetscCall(DMGeomModelRegister("cylinder", DMSnapToGeomModelCylinder));
657625649eSMatthew G. Knepley #if defined(PETSC_HAVE_EGADS)
665552b385SBrandon // FIXME: Brandon uses DMPlexSnapToGeomModel() here instead
677625649eSMatthew G. Knepley PetscCall(DMGeomModelRegister("egads", DMSnapToGeomModel_EGADS));
687625649eSMatthew G. Knepley #endif
697625649eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
707625649eSMatthew G. Knepley }
717625649eSMatthew G. Knepley
727625649eSMatthew G. Knepley /*@C
737625649eSMatthew G. Knepley DMGeomModelRegister - Adds a geometry model to `DM`
747625649eSMatthew G. Knepley
757625649eSMatthew G. Knepley Not Collective, No Fortran Support
767625649eSMatthew G. Knepley
777625649eSMatthew G. Knepley Input Parameters:
78d7c1f440SPierre Jolivet + sname - name of a new user-defined geometry model
797625649eSMatthew G. Knepley - fnc - geometry model function
807625649eSMatthew G. Knepley
817625649eSMatthew G. Knepley Example Usage:
827625649eSMatthew G. Knepley .vb
837625649eSMatthew G. Knepley DMGeomModelRegister("my_geom_model", MySnapToGeomModel);
847625649eSMatthew G. Knepley .ve
857625649eSMatthew G. Knepley
867625649eSMatthew G. Knepley Then, your generator can be chosen with the procedural interface via
87*b44f4de4SBarry Smith .vb
88*b44f4de4SBarry Smith DMSetGeomModel(dm, "my_geom_model",...)
89*b44f4de4SBarry Smith .ve
907625649eSMatthew G. Knepley or at runtime via the option
91*b44f4de4SBarry Smith .vb
92*b44f4de4SBarry Smith -dm_geom_model my_geom_model
93*b44f4de4SBarry Smith .ve
947625649eSMatthew G. Knepley
957625649eSMatthew G. Knepley Level: advanced
967625649eSMatthew G. Knepley
977625649eSMatthew G. Knepley Note:
987625649eSMatthew G. Knepley `DMGeomModelRegister()` may be called multiple times to add several user-defined generators
997625649eSMatthew G. Knepley
1007625649eSMatthew G. Knepley .seealso: `DM`, `DMGeomModelRegisterAll()`, `DMPlexGeomModel()`, `DMGeomModelRegisterDestroy()`
1017625649eSMatthew G. Knepley @*/
DMGeomModelRegister(const char sname[],PetscErrorCode (* fnc)(DM,PetscInt,PetscInt,const PetscScalar[],PetscScalar[]))1027625649eSMatthew G. Knepley PetscErrorCode DMGeomModelRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]))
1037625649eSMatthew G. Knepley {
1047625649eSMatthew G. Knepley PetscFunctionBegin;
105835f2295SStefano Zampini PetscCall(PetscFunctionListAdd(&DMGeomModelList, sname, fnc));
1067625649eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1077625649eSMatthew G. Knepley }
1087625649eSMatthew G. Knepley
DMGeomModelRegisterDestroy(void)1097625649eSMatthew G. Knepley PetscErrorCode DMGeomModelRegisterDestroy(void)
1107625649eSMatthew G. Knepley {
1117625649eSMatthew G. Knepley PetscFunctionBegin;
1127625649eSMatthew G. Knepley PetscCall(PetscFunctionListDestroy(&DMGeomModelList));
1137625649eSMatthew G. Knepley DMGeomModelRegisterAllCalled = PETSC_FALSE;
1147625649eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1157625649eSMatthew G. Knepley }
1167625649eSMatthew G. Knepley
1177625649eSMatthew G. Knepley /*@
1187625649eSMatthew G. Knepley DMSetSnapToGeomModel - Choose a geometry model for this `DM`.
1197625649eSMatthew G. Knepley
1207625649eSMatthew G. Knepley Not Collective
1217625649eSMatthew G. Knepley
1227625649eSMatthew G. Knepley Input Parameters:
1237625649eSMatthew G. Knepley + dm - The `DM` object
1247625649eSMatthew G. Knepley - name - A geometry model name, or `NULL` for the default
1257625649eSMatthew G. Knepley
1267625649eSMatthew G. Knepley Level: intermediate
1277625649eSMatthew G. Knepley
1287625649eSMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexCreate()`, `DMSnapToGeomModel()`
1297625649eSMatthew G. Knepley @*/
DMSetSnapToGeomModel(DM dm,const char name[])1307625649eSMatthew G. Knepley PetscErrorCode DMSetSnapToGeomModel(DM dm, const char name[])
1317625649eSMatthew G. Knepley {
1327625649eSMatthew G. Knepley char geomname[PETSC_MAX_PATH_LEN];
1337625649eSMatthew G. Knepley PetscBool flg;
1347625649eSMatthew G. Knepley
1357625649eSMatthew G. Knepley PetscFunctionBegin;
1367625649eSMatthew G. Knepley if (!name && dm->ops->snaptogeommodel) PetscFunctionReturn(PETSC_SUCCESS);
1377625649eSMatthew G. Knepley PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_geom_model", geomname, sizeof(geomname), &flg));
1387625649eSMatthew G. Knepley if (flg) name = geomname;
1397625649eSMatthew G. Knepley if (!name) {
1407625649eSMatthew G. Knepley PetscObject modelObj;
1417625649eSMatthew G. Knepley
142835f2295SStefano Zampini PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", &modelObj));
1437625649eSMatthew G. Knepley if (modelObj) name = "egads";
1447625649eSMatthew G. Knepley else {
1455552b385SBrandon PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSlite Model", &modelObj));
1465552b385SBrandon if (modelObj) name = "egads";
1477625649eSMatthew G. Knepley }
1487625649eSMatthew G. Knepley }
1497625649eSMatthew G. Knepley if (!name) PetscFunctionReturn(PETSC_SUCCESS);
1507625649eSMatthew G. Knepley
1517625649eSMatthew G. Knepley PetscCall(PetscFunctionListFind(DMGeomModelList, name, &dm->ops->snaptogeommodel));
1527625649eSMatthew G. Knepley PetscCheck(dm->ops->snaptogeommodel, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Geometry model %s not registered; you may need to add --download-%s to your ./configure options", name, name);
1537625649eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1547625649eSMatthew G. Knepley }
1557625649eSMatthew G. Knepley
1567625649eSMatthew G. Knepley /*@
1577625649eSMatthew G. Knepley DMSnapToGeomModel - Given a coordinate point 'mcoords' on the mesh point 'p', return the closest coordinate point 'gcoords' on the geometry model associated with that point.
1587625649eSMatthew G. Knepley
1597625649eSMatthew G. Knepley Not Collective
1607625649eSMatthew G. Knepley
1617625649eSMatthew G. Knepley Input Parameters:
1627625649eSMatthew G. Knepley + dm - The `DMPLEX` object
1637625649eSMatthew G. Knepley . p - The mesh point
1647625649eSMatthew G. Knepley . dE - The coordinate dimension
1657625649eSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point
1667625649eSMatthew G. Knepley
1677625649eSMatthew G. Knepley Output Parameter:
1687625649eSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point
1697625649eSMatthew G. Knepley
1707625649eSMatthew G. Knepley Level: intermediate
1717625649eSMatthew G. Knepley
1727625649eSMatthew G. Knepley Note:
1737625649eSMatthew G. Knepley Returns the original coordinates if no geometry model is found.
1747625649eSMatthew G. Knepley
1757625649eSMatthew G. Knepley The coordinate dimension may be different from the coordinate dimension of the `dm`, for example if the transformation is extrusion.
1767625649eSMatthew G. Knepley
1777625649eSMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexCreate()`, `DMPlexSetRefinementUniform()`
1787625649eSMatthew G. Knepley @*/
DMSnapToGeomModel(DM dm,PetscInt p,PetscInt dE,const PetscScalar mcoords[],PetscScalar gcoords[])1797625649eSMatthew G. Knepley PetscErrorCode DMSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
1807625649eSMatthew G. Knepley {
1817625649eSMatthew G. Knepley PetscFunctionBegin;
1827625649eSMatthew G. Knepley if (!dm->ops->snaptogeommodel)
1837625649eSMatthew G. Knepley for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
1847625649eSMatthew G. Knepley else PetscUseTypeMethod(dm, snaptogeommodel, p, dE, mcoords, gcoords);
1857625649eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
1867625649eSMatthew G. Knepley }
187