xref: /petsc/src/dm/interface/dmgeommodel.c (revision d7c1f4409a34685d8dcd545a97d161d483d89f66)
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)
77625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADS(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
87625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADSLite(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
97625649eSMatthew G. Knepley #endif
107625649eSMatthew G. Knepley 
117625649eSMatthew G. Knepley static PetscErrorCode DMSnapToGeomModelBall(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
127625649eSMatthew G. Knepley {
137625649eSMatthew G. Knepley   PetscInt val;
147625649eSMatthew G. Knepley 
157625649eSMatthew G. Knepley   PetscFunctionBeginUser;
167625649eSMatthew G. Knepley   PetscCall(DMGetLabelValue(dm, "marker", p, &val));
177625649eSMatthew G. Knepley   if (val >= 0) {
187625649eSMatthew G. Knepley     PetscReal norm = 0.;
197625649eSMatthew G. Knepley 
207625649eSMatthew G. Knepley     for (PetscInt d = 0; d < dE; ++d) norm += PetscSqr(PetscRealPart(mcoords[d]));
217625649eSMatthew G. Knepley     norm = PetscSqrtReal(norm);
227625649eSMatthew G. Knepley     for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d] / norm;
237625649eSMatthew G. Knepley   } else {
247625649eSMatthew G. Knepley     for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
257625649eSMatthew G. Knepley   }
267625649eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
277625649eSMatthew G. Knepley }
287625649eSMatthew G. Knepley 
297625649eSMatthew G. Knepley static PetscErrorCode DMSnapToGeomModelCylinder(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
307625649eSMatthew G. Knepley {
317625649eSMatthew G. Knepley   PetscReal gmin[3], gmax[3];
327625649eSMatthew G. Knepley   PetscInt  val;
337625649eSMatthew G. Knepley 
347625649eSMatthew G. Knepley   PetscFunctionBeginUser;
357625649eSMatthew G. Knepley   PetscCall(DMGetBoundingBox(dm, gmin, gmax));
367625649eSMatthew G. Knepley   PetscCall(DMGetLabelValue(dm, "generatrix", p, &val));
377625649eSMatthew G. Knepley   if (val >= 0) {
387625649eSMatthew G. Knepley     PetscReal norm = 0.;
397625649eSMatthew G. Knepley 
407625649eSMatthew G. Knepley     for (PetscInt d = 0; d < dE - 1; ++d) norm += PetscSqr(PetscRealPart(mcoords[d]));
417625649eSMatthew G. Knepley     norm = PetscSqrtReal(norm);
427625649eSMatthew G. Knepley     for (PetscInt d = 0; d < dE - 1; ++d) gcoords[d] = gmax[0] * mcoords[d] / norm;
437625649eSMatthew G. Knepley     gcoords[dE - 1] = mcoords[dE - 1];
447625649eSMatthew G. Knepley   } else {
457625649eSMatthew G. Knepley     for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
467625649eSMatthew G. Knepley   }
477625649eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
487625649eSMatthew G. Knepley }
497625649eSMatthew G. Knepley 
507625649eSMatthew G. Knepley /*@C
517625649eSMatthew G. Knepley   DMGeomModelRegisterAll - Registers all of the geometry model methods in the `DM` package.
527625649eSMatthew G. Knepley 
537625649eSMatthew G. Knepley   Not Collective
547625649eSMatthew G. Knepley 
557625649eSMatthew G. Knepley   Level: advanced
567625649eSMatthew G. Knepley 
577625649eSMatthew G. Knepley .seealso: `DM`, `DMGeomModelRegisterDestroy()`
587625649eSMatthew G. Knepley @*/
597625649eSMatthew G. Knepley PetscErrorCode DMGeomModelRegisterAll(void)
607625649eSMatthew G. Knepley {
617625649eSMatthew G. Knepley   PetscFunctionBegin;
627625649eSMatthew G. Knepley   if (DMGeomModelRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
637625649eSMatthew G. Knepley   DMGeomModelRegisterAllCalled = PETSC_TRUE;
647625649eSMatthew G. Knepley   PetscCall(DMGeomModelRegister("ball", DMSnapToGeomModelBall));
657625649eSMatthew G. Knepley   PetscCall(DMGeomModelRegister("cylinder", DMSnapToGeomModelCylinder));
667625649eSMatthew G. Knepley #if defined(PETSC_HAVE_EGADS)
677625649eSMatthew G. Knepley   PetscCall(DMGeomModelRegister("egads", DMSnapToGeomModel_EGADS));
687625649eSMatthew G. Knepley   PetscCall(DMGeomModelRegister("egadslite", DMSnapToGeomModel_EGADSLite));
697625649eSMatthew G. Knepley #endif
707625649eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
717625649eSMatthew G. Knepley }
727625649eSMatthew G. Knepley 
737625649eSMatthew G. Knepley /*@C
747625649eSMatthew G. Knepley   DMGeomModelRegister -  Adds a geometry model to `DM`
757625649eSMatthew G. Knepley 
767625649eSMatthew G. Knepley   Not Collective, No Fortran Support
777625649eSMatthew G. Knepley 
787625649eSMatthew G. Knepley   Input Parameters:
79*d7c1f440SPierre Jolivet + sname - name of a new user-defined geometry model
807625649eSMatthew G. Knepley - fnc   - geometry model function
817625649eSMatthew G. Knepley 
827625649eSMatthew G. Knepley   Example Usage:
837625649eSMatthew G. Knepley .vb
847625649eSMatthew G. Knepley    DMGeomModelRegister("my_geom_model", MySnapToGeomModel);
857625649eSMatthew G. Knepley .ve
867625649eSMatthew G. Knepley 
877625649eSMatthew G. Knepley   Then, your generator can be chosen with the procedural interface via
887625649eSMatthew G. Knepley $     DMSetGeomModel(dm, "my_geom_model",...)
897625649eSMatthew G. Knepley   or at runtime via the option
907625649eSMatthew G. Knepley $     -dm_geom_model my_geom_model
917625649eSMatthew G. Knepley 
927625649eSMatthew G. Knepley   Level: advanced
937625649eSMatthew G. Knepley 
947625649eSMatthew G. Knepley   Note:
957625649eSMatthew G. Knepley   `DMGeomModelRegister()` may be called multiple times to add several user-defined generators
967625649eSMatthew G. Knepley 
977625649eSMatthew G. Knepley .seealso: `DM`, `DMGeomModelRegisterAll()`, `DMPlexGeomModel()`, `DMGeomModelRegisterDestroy()`
987625649eSMatthew G. Knepley @*/
997625649eSMatthew G. Knepley PetscErrorCode DMGeomModelRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]))
1007625649eSMatthew G. Knepley {
1017625649eSMatthew G. Knepley   PetscFunctionBegin;
1027625649eSMatthew G. Knepley   PetscCall(PetscFunctionListAdd(&DMGeomModelList, sname, (PetscVoidFn *)fnc));
1037625649eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1047625649eSMatthew G. Knepley }
1057625649eSMatthew G. Knepley 
1067625649eSMatthew G. Knepley extern PetscBool DMGeomModelRegisterAllCalled;
1077625649eSMatthew G. Knepley 
1087625649eSMatthew G. Knepley PetscErrorCode DMGeomModelRegisterDestroy(void)
1097625649eSMatthew G. Knepley {
1107625649eSMatthew G. Knepley   PetscFunctionBegin;
1117625649eSMatthew G. Knepley   PetscCall(PetscFunctionListDestroy(&DMGeomModelList));
1127625649eSMatthew G. Knepley   DMGeomModelRegisterAllCalled = PETSC_FALSE;
1137625649eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1147625649eSMatthew G. Knepley }
1157625649eSMatthew G. Knepley 
1167625649eSMatthew G. Knepley /*@
1177625649eSMatthew G. Knepley   DMSetSnapToGeomModel - Choose a geometry model for this `DM`.
1187625649eSMatthew G. Knepley 
1197625649eSMatthew G. Knepley   Not Collective
1207625649eSMatthew G. Knepley 
1217625649eSMatthew G. Knepley   Input Parameters:
1227625649eSMatthew G. Knepley + dm   - The `DM` object
1237625649eSMatthew G. Knepley - name - A geometry model name, or `NULL` for the default
1247625649eSMatthew G. Knepley 
1257625649eSMatthew G. Knepley   Level: intermediate
1267625649eSMatthew G. Knepley 
1277625649eSMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexCreate()`, `DMSnapToGeomModel()`
1287625649eSMatthew G. Knepley @*/
1297625649eSMatthew G. Knepley PetscErrorCode DMSetSnapToGeomModel(DM dm, const char name[])
1307625649eSMatthew G. Knepley {
1317625649eSMatthew G. Knepley   char      geomname[PETSC_MAX_PATH_LEN];
1327625649eSMatthew G. Knepley   PetscBool flg;
1337625649eSMatthew G. Knepley 
1347625649eSMatthew G. Knepley   PetscFunctionBegin;
1357625649eSMatthew G. Knepley   if (!name && dm->ops->snaptogeommodel) PetscFunctionReturn(PETSC_SUCCESS);
1367625649eSMatthew G. Knepley   PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_geom_model", geomname, sizeof(geomname), &flg));
1377625649eSMatthew G. Knepley   if (flg) name = geomname;
1387625649eSMatthew G. Knepley   if (!name) {
1397625649eSMatthew G. Knepley     PetscObject modelObj;
1407625649eSMatthew G. Knepley 
1417625649eSMatthew G. Knepley     PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
1427625649eSMatthew G. Knepley     if (modelObj) name = "egads";
1437625649eSMatthew G. Knepley     else {
1447625649eSMatthew G. Knepley       PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSLite Model", (PetscObject *)&modelObj));
1457625649eSMatthew G. Knepley       if (modelObj) name = "egadslite";
1467625649eSMatthew G. Knepley     }
1477625649eSMatthew G. Knepley   }
1487625649eSMatthew G. Knepley   if (!name) PetscFunctionReturn(PETSC_SUCCESS);
1497625649eSMatthew G. Knepley 
1507625649eSMatthew G. Knepley   PetscCall(PetscFunctionListFind(DMGeomModelList, name, &dm->ops->snaptogeommodel));
1517625649eSMatthew 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);
1527625649eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1537625649eSMatthew G. Knepley }
1547625649eSMatthew G. Knepley 
1557625649eSMatthew G. Knepley /*@
1567625649eSMatthew 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.
1577625649eSMatthew G. Knepley 
1587625649eSMatthew G. Knepley   Not Collective
1597625649eSMatthew G. Knepley 
1607625649eSMatthew G. Knepley   Input Parameters:
1617625649eSMatthew G. Knepley + dm      - The `DMPLEX` object
1627625649eSMatthew G. Knepley . p       - The mesh point
1637625649eSMatthew G. Knepley . dE      - The coordinate dimension
1647625649eSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point
1657625649eSMatthew G. Knepley 
1667625649eSMatthew G. Knepley   Output Parameter:
1677625649eSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point
1687625649eSMatthew G. Knepley 
1697625649eSMatthew G. Knepley   Level: intermediate
1707625649eSMatthew G. Knepley 
1717625649eSMatthew G. Knepley   Note:
1727625649eSMatthew G. Knepley   Returns the original coordinates if no geometry model is found.
1737625649eSMatthew G. Knepley 
1747625649eSMatthew G. Knepley   The coordinate dimension may be different from the coordinate dimension of the `dm`, for example if the transformation is extrusion.
1757625649eSMatthew G. Knepley 
1767625649eSMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexCreate()`, `DMPlexSetRefinementUniform()`
1777625649eSMatthew G. Knepley @*/
1787625649eSMatthew G. Knepley PetscErrorCode DMSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
1797625649eSMatthew G. Knepley {
1807625649eSMatthew G. Knepley   PetscFunctionBegin;
1817625649eSMatthew G. Knepley   if (!dm->ops->snaptogeommodel)
1827625649eSMatthew G. Knepley     for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
1837625649eSMatthew G. Knepley   else PetscUseTypeMethod(dm, snaptogeommodel, p, dE, mcoords, gcoords);
1847625649eSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1857625649eSMatthew G. Knepley }
186