xref: /petsc/src/dm/interface/dmgenerate.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
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: `DMGenerateRegisterDestroy()`
40 @*/
41 PetscErrorCode DMGenerateRegisterAll(void)
42 {
43   PetscFunctionBegin;
44   if (DMGenerateRegisterAllCalled) PetscFunctionReturn(0);
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(0);
67 }
68 
69 /*@C
70   DMGenerateRegister -  Adds a grid generator to DM
71 
72    Not Collective
73 
74    Input Parameters:
75 +  name_solver - 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    Notes:
82    DMGenerateRegister() may be called multiple times to add several user-defined solvers.
83 
84    Sample usage:
85 .vb
86    DMGenerateRegister("my_generator",MyGeneratorCreate,MyGeneratorRefiner,MyGeneratorAdaptor,dim);
87 .ve
88 
89    Then, your generator can be chosen with the procedural interface via
90 $     DMGenerate(dm,"my_generator",...)
91    or at runtime via the option
92 $     -dm_generator my_generator
93 
94    Level: advanced
95 
96 .seealso: `DMGenerateRegisterAll()`, `DMPlexGenerate()`, `DMGenerateRegisterDestroy()`
97 
98 @*/
99 PetscErrorCode DMGenerateRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscBool, DM*), PetscErrorCode (*rfnc)(DM, PetscReal*, DM*), PetscErrorCode (*alfnc)(DM, Vec, DMLabel, DMLabel, DM*), PetscInt dim)
100 {
101   DMGeneratorFunctionList entry;
102 
103   PetscFunctionBegin;
104   PetscCall(PetscNew(&entry));
105   PetscCall(PetscStrallocpy(sname,&entry->name));
106   entry->generate = fnc;
107   entry->refine   = rfnc;
108   entry->adapt    = alfnc;
109   entry->dim      = dim;
110   entry->next     = NULL;
111   if (!DMGenerateList) DMGenerateList = entry;
112   else {
113     DMGeneratorFunctionList fl = DMGenerateList;
114     while (fl->next) fl = fl->next;
115     fl->next = entry;
116   }
117   PetscFunctionReturn(0);
118 }
119 
120 extern PetscBool DMGenerateRegisterAllCalled;
121 
122 PetscErrorCode DMGenerateRegisterDestroy(void)
123 {
124   DMGeneratorFunctionList next, fl;
125 
126   PetscFunctionBegin;
127   next = fl = DMGenerateList;
128   while (next) {
129     next = fl ? fl->next : NULL;
130     if (fl) PetscCall(PetscFree(fl->name));
131     PetscCall(PetscFree(fl));
132     fl = next;
133   }
134   DMGenerateList              = NULL;
135   DMGenerateRegisterAllCalled = PETSC_FALSE;
136   PetscFunctionReturn(0);
137 }
138 
139 /*@C
140   DMAdaptLabel - Adapt a dm based on a label with values interpreted as coarsening and refining flags.  Specific implementations of DM maybe have
141                  specialized flags, but all implementations should accept flag values DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN.
142 
143   Collective on dm
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: `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) PetscValidPointer(label, 2);
167   PetscValidPointer(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     (*dmAdapt)->prealloc_only = dm->prealloc_only;  /* maybe this should go .... */
187     PetscCall(PetscFree((*dmAdapt)->vectype));
188     PetscCall(PetscStrallocpy(dm->vectype,(char**)&(*dmAdapt)->vectype));
189     PetscCall(PetscFree((*dmAdapt)->mattype));
190     PetscCall(PetscStrallocpy(dm->mattype,(char**)&(*dmAdapt)->mattype));
191   }
192   PetscFunctionReturn(0);
193 }
194 
195 /*@C
196   DMAdaptMetric - Generates a mesh adapted to the specified metric field.
197 
198   Input Parameters:
199 + dm - The DM object
200 . metric - The metric to which the mesh is adapted, defined vertex-wise.
201 . 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_".
202 - 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_".
203 
204   Output Parameter:
205 . dmAdapt  - Pointer to the DM object containing the adapted mesh
206 
207   Note: The label in the adapted mesh will be registered under the name of the input DMLabel object
208 
209   Level: advanced
210 
211 .seealso: `DMAdaptLabel()`, `DMCoarsen()`, `DMRefine()`
212 @*/
213 PetscErrorCode DMAdaptMetric(DM dm, Vec metric, DMLabel bdLabel, DMLabel rgLabel, DM *dmAdapt)
214 {
215   DMGeneratorFunctionList fl;
216   char                    adaptname[PETSC_MAX_PATH_LEN];
217   const char             *name;
218   const char * const      adaptors[3] = {"pragmatic", "mmg", "parmmg"};
219   PetscInt                dim;
220   PetscBool               flg, found = PETSC_FALSE;
221 
222   PetscFunctionBegin;
223   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
224   PetscValidHeaderSpecific(metric, VEC_CLASSID, 2);
225   if (bdLabel) PetscValidPointer(bdLabel, 3);
226   if (rgLabel) PetscValidPointer(rgLabel, 4);
227   PetscValidPointer(dmAdapt, 5);
228   *dmAdapt = NULL;
229   PetscCall(DMGetDimension(dm, &dim));
230   PetscCall(PetscOptionsGetString(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_adaptor", adaptname, sizeof(adaptname), &flg));
231 
232   /* Default to Mmg in serial and ParMmg in parallel */
233   if (flg) name = adaptname;
234   else {
235     MPI_Comm                comm;
236     PetscMPIInt             size;
237 
238     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
239     PetscCallMPI(MPI_Comm_size(comm, &size));
240     if (size == 1) name = adaptors[1];
241     else           name = adaptors[2];
242   }
243 
244   fl = DMGenerateList;
245   while (fl) {
246     PetscCall(PetscStrcmp(fl->name, name, &flg));
247     if (flg) {
248       PetscCall((*fl->adapt)(dm, metric, bdLabel, rgLabel, dmAdapt));
249       found = PETSC_TRUE;
250     }
251     fl = fl->next;
252   }
253   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);
254   if (*dmAdapt) {
255     (*dmAdapt)->prealloc_only = dm->prealloc_only;  /* maybe this should go .... */
256     PetscCall(PetscFree((*dmAdapt)->vectype));
257     PetscCall(PetscStrallocpy(dm->vectype,(char**)&(*dmAdapt)->vectype));
258     PetscCall(PetscFree((*dmAdapt)->mattype));
259     PetscCall(PetscStrallocpy(dm->mattype,(char**)&(*dmAdapt)->mattype));
260   }
261   PetscFunctionReturn(0);
262 }
263