Actual source code: dmgeommodel.c
1: #include <petsc/private/dmimpl.h>
3: PetscFunctionList DMGeomModelList = NULL;
4: PetscBool DMGeomModelRegisterAllCalled = PETSC_FALSE;
6: #if defined(PETSC_HAVE_EGADS)
7: PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADS(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
8: PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADSLite(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
9: #endif
11: static PetscErrorCode DMSnapToGeomModelBall(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
12: {
13: PetscInt val;
15: PetscFunctionBeginUser;
16: PetscCall(DMGetLabelValue(dm, "marker", p, &val));
17: if (val >= 0) {
18: PetscReal norm = 0.;
20: for (PetscInt d = 0; d < dE; ++d) norm += PetscSqr(PetscRealPart(mcoords[d]));
21: norm = PetscSqrtReal(norm);
22: for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d] / norm;
23: } else {
24: for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
25: }
26: PetscFunctionReturn(PETSC_SUCCESS);
27: }
29: static PetscErrorCode DMSnapToGeomModelCylinder(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
30: {
31: PetscReal gmin[3], gmax[3];
32: PetscInt val;
34: PetscFunctionBeginUser;
35: PetscCall(DMGetBoundingBox(dm, gmin, gmax));
36: PetscCall(DMGetLabelValue(dm, "generatrix", p, &val));
37: if (val >= 0) {
38: PetscReal norm = 0.;
40: for (PetscInt d = 0; d < dE - 1; ++d) norm += PetscSqr(PetscRealPart(mcoords[d]));
41: norm = PetscSqrtReal(norm);
42: for (PetscInt d = 0; d < dE - 1; ++d) gcoords[d] = gmax[0] * mcoords[d] / norm;
43: gcoords[dE - 1] = mcoords[dE - 1];
44: } else {
45: for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
46: }
47: PetscFunctionReturn(PETSC_SUCCESS);
48: }
50: /*@C
51: DMGeomModelRegisterAll - Registers all of the geometry model methods in the `DM` package.
53: Not Collective
55: Level: advanced
57: .seealso: `DM`, `DMGeomModelRegisterDestroy()`
58: @*/
59: PetscErrorCode DMGeomModelRegisterAll(void)
60: {
61: PetscFunctionBegin;
62: if (DMGeomModelRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
63: DMGeomModelRegisterAllCalled = PETSC_TRUE;
64: PetscCall(DMGeomModelRegister("ball", DMSnapToGeomModelBall));
65: PetscCall(DMGeomModelRegister("cylinder", DMSnapToGeomModelCylinder));
66: #if defined(PETSC_HAVE_EGADS)
67: PetscCall(DMGeomModelRegister("egads", DMSnapToGeomModel_EGADS));
68: PetscCall(DMGeomModelRegister("egadslite", DMSnapToGeomModel_EGADSLite));
69: #endif
70: PetscFunctionReturn(PETSC_SUCCESS);
71: }
73: /*@C
74: DMGeomModelRegister - Adds a geometry model to `DM`
76: Not Collective, No Fortran Support
78: Input Parameters:
79: + sname - name of a new user-defined geometry model
80: - fnc - geometry model function
82: Example Usage:
83: .vb
84: DMGeomModelRegister("my_geom_model", MySnapToGeomModel);
85: .ve
87: Then, your generator can be chosen with the procedural interface via
88: $ DMSetGeomModel(dm, "my_geom_model",...)
89: or at runtime via the option
90: $ -dm_geom_model my_geom_model
92: Level: advanced
94: Note:
95: `DMGeomModelRegister()` may be called multiple times to add several user-defined generators
97: .seealso: `DM`, `DMGeomModelRegisterAll()`, `DMPlexGeomModel()`, `DMGeomModelRegisterDestroy()`
98: @*/
99: PetscErrorCode DMGeomModelRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]))
100: {
101: PetscFunctionBegin;
102: PetscCall(PetscFunctionListAdd(&DMGeomModelList, sname, (PetscVoidFn *)fnc));
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
106: extern PetscBool DMGeomModelRegisterAllCalled;
108: PetscErrorCode DMGeomModelRegisterDestroy(void)
109: {
110: PetscFunctionBegin;
111: PetscCall(PetscFunctionListDestroy(&DMGeomModelList));
112: DMGeomModelRegisterAllCalled = PETSC_FALSE;
113: PetscFunctionReturn(PETSC_SUCCESS);
114: }
116: /*@
117: DMSetSnapToGeomModel - Choose a geometry model for this `DM`.
119: Not Collective
121: Input Parameters:
122: + dm - The `DM` object
123: - name - A geometry model name, or `NULL` for the default
125: Level: intermediate
127: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexCreate()`, `DMSnapToGeomModel()`
128: @*/
129: PetscErrorCode DMSetSnapToGeomModel(DM dm, const char name[])
130: {
131: char geomname[PETSC_MAX_PATH_LEN];
132: PetscBool flg;
134: PetscFunctionBegin;
135: if (!name && dm->ops->snaptogeommodel) PetscFunctionReturn(PETSC_SUCCESS);
136: PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_geom_model", geomname, sizeof(geomname), &flg));
137: if (flg) name = geomname;
138: if (!name) {
139: PetscObject modelObj;
141: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
142: if (modelObj) name = "egads";
143: else {
144: PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSLite Model", (PetscObject *)&modelObj));
145: if (modelObj) name = "egadslite";
146: }
147: }
148: if (!name) PetscFunctionReturn(PETSC_SUCCESS);
150: PetscCall(PetscFunctionListFind(DMGeomModelList, name, &dm->ops->snaptogeommodel));
151: 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);
152: PetscFunctionReturn(PETSC_SUCCESS);
153: }
155: /*@
156: 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.
158: Not Collective
160: Input Parameters:
161: + dm - The `DMPLEX` object
162: . p - The mesh point
163: . dE - The coordinate dimension
164: - mcoords - A coordinate point lying on the mesh point
166: Output Parameter:
167: . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point
169: Level: intermediate
171: Note:
172: Returns the original coordinates if no geometry model is found.
174: The coordinate dimension may be different from the coordinate dimension of the `dm`, for example if the transformation is extrusion.
176: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexCreate()`, `DMPlexSetRefinementUniform()`
177: @*/
178: PetscErrorCode DMSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
179: {
180: PetscFunctionBegin;
181: if (!dm->ops->snaptogeommodel)
182: for (PetscInt d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
183: else PetscUseTypeMethod(dm, snaptogeommodel, p, dE, mcoords, gcoords);
184: PetscFunctionReturn(PETSC_SUCCESS);
185: }