xref: /petsc/src/dm/interface/dmgenerate.c (revision 4d5500e8ab0f8e7b1fe2746469a330fafa8aaa90)
1 #include <petsc/private/dmimpl.h> /*I      "petscdm.h"          I*/
2 
3 PETSC_EXTERN PetscErrorCode DMIsForest(DM, PetscBool *);
4 
5 DMGeneratorFunctionList DMGenerateList              = NULL;
6 PetscBool               DMGenerateRegisterAllCalled = PETSC_FALSE;
7 
8 #if defined(PETSC_HAVE_TRIANGLE)
9 PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM, PetscBool, DM *);
10 PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM, double *, DM *);
11 #endif
12 #if defined(PETSC_HAVE_TETGEN)
13 PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM, PetscBool, DM *);
14 PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM, double *, DM *);
15 #endif
16 #if defined(PETSC_HAVE_CTETGEN)
17 PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM, PetscBool, DM *);
18 PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM, double *, DM *);
19 #endif
20 #if defined(PETSC_HAVE_PRAGMATIC)
21 PETSC_EXTERN PetscErrorCode DMAdaptMetric_Pragmatic_Plex(DM, Vec, DMLabel, DMLabel, DM *);
22 #endif
23 #if defined(PETSC_HAVE_MMG)
24 PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM, Vec, DMLabel, DMLabel, DM *);
25 #endif
26 #if defined(PETSC_HAVE_PARMMG)
27 PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM, Vec, DMLabel, DMLabel, DM *);
28 #endif
29 PETSC_EXTERN PetscErrorCode DMPlexTransformAdaptLabel(DM, Vec, DMLabel, DMLabel, DM *);
30 PETSC_EXTERN PetscErrorCode DMAdaptLabel_Forest(DM, Vec, DMLabel, DMLabel, DM *);
31 
32 /*@C
33   DMGenerateRegisterAll - Registers all of the mesh generation methods in the `DM` package.
34 
35   Not Collective
36 
37   Level: advanced
38 
39 .seealso: `DM`, `DMGenerateRegisterDestroy()`
40 @*/
41 PetscErrorCode DMGenerateRegisterAll(void)
42 {
43   PetscFunctionBegin;
44   if (DMGenerateRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
45   DMGenerateRegisterAllCalled = PETSC_TRUE;
46 #if defined(PETSC_HAVE_TRIANGLE)
47   PetscCall(DMGenerateRegister("triangle", DMPlexGenerate_Triangle, DMPlexRefine_Triangle, NULL, 1));
48 #endif
49 #if defined(PETSC_HAVE_CTETGEN)
50   PetscCall(DMGenerateRegister("ctetgen", DMPlexGenerate_CTetgen, DMPlexRefine_CTetgen, NULL, 2));
51 #endif
52 #if defined(PETSC_HAVE_TETGEN)
53   PetscCall(DMGenerateRegister("tetgen", DMPlexGenerate_Tetgen, DMPlexRefine_Tetgen, NULL, 2));
54 #endif
55 #if defined(PETSC_HAVE_PRAGMATIC)
56   PetscCall(DMGenerateRegister("pragmatic", NULL, NULL, DMAdaptMetric_Pragmatic_Plex, -1));
57 #endif
58 #if defined(PETSC_HAVE_MMG)
59   PetscCall(DMGenerateRegister("mmg", NULL, NULL, DMAdaptMetric_Mmg_Plex, -1));
60 #endif
61 #if defined(PETSC_HAVE_PARMMG)
62   PetscCall(DMGenerateRegister("parmmg", NULL, NULL, DMAdaptMetric_ParMmg_Plex, -1));
63 #endif
64   PetscCall(DMGenerateRegister("cellrefiner", NULL, NULL, DMPlexTransformAdaptLabel, -1));
65   PetscCall(DMGenerateRegister("forest", NULL, NULL, DMAdaptLabel_Forest, -1));
66   PetscFunctionReturn(PETSC_SUCCESS);
67 }
68 
69 /*@C
70   DMGenerateRegister -  Adds a grid generator to `DM`
71 
72   Not Collective, No Fortran Support
73 
74   Input Parameters:
75 + sname - name of a new user-defined grid generator
76 . fnc   - generator function
77 . rfnc  - refinement function
78 . alfnc - adapt by label function
79 - dim   - dimension of boundary of domain
80 
81   Example Usage:
82 .vb
83    DMGenerateRegister("my_generator", MyGeneratorCreate, MyGeneratorRefiner, MyGeneratorAdaptor, dim);
84 .ve
85 
86   Then, your generator can be chosen with the procedural interface via
87 $     DMGenerate(dm, "my_generator",...)
88   or at runtime via the option
89 $     -dm_generator my_generator
90 
91   Level: advanced
92 
93   Note:
94   `DMGenerateRegister()` may be called multiple times to add several user-defined generators
95 
96 .seealso: `DM`, `DMGenerateRegisterAll()`, `DMPlexGenerate()`, `DMGenerateRegisterDestroy()`
97 @*/
98 PetscErrorCode DMGenerateRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscBool, DM *), PetscErrorCode (*rfnc)(DM, PetscReal *, DM *), PetscErrorCode (*alfnc)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt dim)
99 {
100   DMGeneratorFunctionList entry;
101 
102   PetscFunctionBegin;
103   PetscCall(PetscNew(&entry));
104   PetscCall(PetscStrallocpy(sname, &entry->name));
105   entry->generate = fnc;
106   entry->refine   = rfnc;
107   entry->adapt    = alfnc;
108   entry->dim      = dim;
109   entry->next     = NULL;
110   if (!DMGenerateList) DMGenerateList = entry;
111   else {
112     DMGeneratorFunctionList fl = DMGenerateList;
113     while (fl->next) fl = fl->next;
114     fl->next = entry;
115   }
116   PetscFunctionReturn(PETSC_SUCCESS);
117 }
118 
119 extern PetscBool DMGenerateRegisterAllCalled;
120 
121 PetscErrorCode DMGenerateRegisterDestroy(void)
122 {
123   DMGeneratorFunctionList next, fl;
124 
125   PetscFunctionBegin;
126   next = fl = DMGenerateList;
127   while (next) {
128     next = fl ? fl->next : NULL;
129     if (fl) PetscCall(PetscFree(fl->name));
130     PetscCall(PetscFree(fl));
131     fl = next;
132   }
133   DMGenerateList              = NULL;
134   DMGenerateRegisterAllCalled = PETSC_FALSE;
135   PetscFunctionReturn(PETSC_SUCCESS);
136 }
137 
138 /*@
139   DMAdaptLabel - Adapt a `DM` based on a `DMLabel` with values interpreted as coarsening and refining flags.  Specific implementations of `DM` maybe have
140   specialized flags, but all implementations should accept flag values `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and,
141   `DM_ADAPT_COARSEN`.
142 
143   Collective
144 
145   Input Parameters:
146 + dm    - the pre-adaptation `DM` object
147 - label - label with the flags
148 
149   Output Parameters:
150 . dmAdapt - the adapted `DM` object: may be `NULL` if an adapted `DM` could not be produced.
151 
152   Level: intermediate
153 
154 .seealso: `DM`, `DMAdaptMetric()`, `DMCoarsen()`, `DMRefine()`
155 @*/
156 PetscErrorCode DMAdaptLabel(DM dm, DMLabel label, DM *dmAdapt)
157 {
158   DMGeneratorFunctionList fl;
159   char                    adaptname[PETSC_MAX_PATH_LEN];
160   const char             *name;
161   PetscInt                dim;
162   PetscBool               flg, isForest, found = PETSC_FALSE;
163 
164   PetscFunctionBegin;
165   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
166   if (label) PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 2);
167   PetscAssertPointer(dmAdapt, 3);
168   *dmAdapt = NULL;
169   PetscCall(DMGetDimension(dm, &dim));
170   PetscCall(DMIsForest(dm, &isForest));
171   name = isForest ? "forest" : "cellrefiner";
172   PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", adaptname, sizeof(adaptname), &flg));
173   if (flg) name = adaptname;
174 
175   fl = DMGenerateList;
176   while (fl) {
177     PetscCall(PetscStrcmp(fl->name, name, &flg));
178     if (flg) {
179       PetscCall((*fl->adapt)(dm, NULL, label, NULL, dmAdapt));
180       found = PETSC_TRUE;
181     }
182     fl = fl->next;
183   }
184   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Grid adaptor %s not registered; you may need to add --download-%s to your ./configure options", name, name);
185   if (*dmAdapt) {
186     void *ctx;
187 
188     (*dmAdapt)->prealloc_only = dm->prealloc_only; /* maybe this should go .... */
189     PetscCall(PetscFree((*dmAdapt)->vectype));
190     PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*dmAdapt)->vectype));
191     PetscCall(PetscFree((*dmAdapt)->mattype));
192     PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*dmAdapt)->mattype));
193     PetscCall(DMGetApplicationContext(dm, &ctx));
194     PetscCall(DMSetApplicationContext(*dmAdapt, ctx));
195   }
196   PetscFunctionReturn(PETSC_SUCCESS);
197 }
198 
199 /*@
200   DMAdaptMetric - Generates a mesh adapted to the specified metric field.
201 
202   Input Parameters:
203 + dm      - The DM object
204 . metric  - The metric to which the mesh is adapted, defined vertex-wise.
205 . bdLabel - Label for boundary tags, which will be preserved in the output mesh. `bdLabel` should be `NULL` if there is no such label, and should be different from "_boundary_".
206 - rgLabel - Label for cell tags, which will be preserved in the output mesh. `rgLabel` should be `NULL` if there is no such label, and should be different from "_regions_".
207 
208   Output Parameter:
209 . dmAdapt - Pointer to the `DM` object containing the adapted mesh
210 
211   Note:
212   The label in the adapted mesh will be registered under the name of the input `DMLabel` object
213 
214   Level: advanced
215 
216 .seealso: `DMAdaptLabel()`, `DMCoarsen()`, `DMRefine()`
217 @*/
218 PetscErrorCode DMAdaptMetric(DM dm, Vec metric, DMLabel bdLabel, DMLabel rgLabel, DM *dmAdapt)
219 {
220   DMGeneratorFunctionList fl;
221   char                    adaptname[PETSC_MAX_PATH_LEN];
222   const char             *name;
223   const char *const       adaptors[3] = {"pragmatic", "mmg", "parmmg"};
224   PetscInt                dim;
225   PetscBool               flg, found = PETSC_FALSE;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
229   PetscValidHeaderSpecific(metric, VEC_CLASSID, 2);
230   if (bdLabel) PetscValidHeaderSpecific(bdLabel, DMLABEL_CLASSID, 3);
231   if (rgLabel) PetscValidHeaderSpecific(rgLabel, DMLABEL_CLASSID, 4);
232   PetscAssertPointer(dmAdapt, 5);
233   *dmAdapt = NULL;
234   PetscCall(DMGetDimension(dm, &dim));
235   PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", adaptname, sizeof(adaptname), &flg));
236 
237   /* Default to Mmg in serial and ParMmg in parallel */
238   if (flg) name = adaptname;
239   else {
240     MPI_Comm    comm;
241     PetscMPIInt size;
242 
243     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
244     PetscCallMPI(MPI_Comm_size(comm, &size));
245     if (size == 1) name = adaptors[1];
246     else name = adaptors[2];
247   }
248 
249   fl = DMGenerateList;
250   while (fl) {
251     PetscCall(PetscStrcmp(fl->name, name, &flg));
252     if (flg) {
253       PetscCall((*fl->adapt)(dm, metric, bdLabel, rgLabel, dmAdapt));
254       found = PETSC_TRUE;
255     }
256     fl = fl->next;
257   }
258   PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Grid adaptor %s not registered; you may need to add --download-%s to your ./configure options", name, name);
259   if (*dmAdapt) {
260     (*dmAdapt)->prealloc_only = dm->prealloc_only; /* maybe this should go .... */
261     PetscCall(PetscFree((*dmAdapt)->vectype));
262     PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*dmAdapt)->vectype));
263     PetscCall(PetscFree((*dmAdapt)->mattype));
264     PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*dmAdapt)->mattype));
265   }
266   PetscFunctionReturn(PETSC_SUCCESS);
267 }
268