xref: /petsc/src/dm/interface/dmgenerate.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
1c0517cd5SMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I      "petscdm.h"          I*/
2c0517cd5SMatthew G. Knepley 
3c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMIsForest(DM, PetscBool *);
4c0517cd5SMatthew G. Knepley 
5c0517cd5SMatthew G. Knepley DMGeneratorFunctionList DMGenerateList              = NULL;
6c0517cd5SMatthew G. Knepley PetscBool               DMGenerateRegisterAllCalled = PETSC_FALSE;
7c0517cd5SMatthew G. Knepley 
8c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_TRIANGLE)
9c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM, PetscBool, DM *);
10c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM, double *, DM *);
11c0517cd5SMatthew G. Knepley #endif
12c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_TETGEN)
13c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM, PetscBool, DM *);
14c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM, double *, DM *);
15c0517cd5SMatthew G. Knepley #endif
16c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
17c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM, PetscBool, DM *);
18c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM, double *, DM *);
19c0517cd5SMatthew G. Knepley #endif
20c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_PRAGMATIC)
219fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric_Pragmatic_Plex(DM, Vec, DMLabel, DMLabel, DM *);
22c0517cd5SMatthew G. Knepley #endif
23c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_MMG)
249fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM, Vec, DMLabel, DMLabel, DM *);
25c0517cd5SMatthew G. Knepley #endif
26c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_PARMMG)
279fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM, Vec, DMLabel, DMLabel, DM *);
28c0517cd5SMatthew G. Knepley #endif
299fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMPlexTransformAdaptLabel(DM, Vec, DMLabel, DMLabel, DM *);
309fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptLabel_Forest(DM, Vec, DMLabel, DMLabel, DM *);
31c0517cd5SMatthew G. Knepley 
32c0517cd5SMatthew G. Knepley /*@C
3316a05f60SBarry Smith   DMGenerateRegisterAll - Registers all of the mesh generation methods in the `DM` package.
34c0517cd5SMatthew G. Knepley 
35c0517cd5SMatthew G. Knepley   Not Collective
36c0517cd5SMatthew G. Knepley 
37c0517cd5SMatthew G. Knepley   Level: advanced
38c0517cd5SMatthew G. Knepley 
3916a05f60SBarry Smith .seealso: `DM`, `DMGenerateRegisterDestroy()`
40c0517cd5SMatthew G. Knepley @*/
41d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGenerateRegisterAll(void)
42d71ae5a4SJacob Faibussowitsch {
43c0517cd5SMatthew G. Knepley   PetscFunctionBegin;
443ba16761SJacob Faibussowitsch   if (DMGenerateRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
45c0517cd5SMatthew G. Knepley   DMGenerateRegisterAllCalled = PETSC_TRUE;
46c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_TRIANGLE)
479566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("triangle", DMPlexGenerate_Triangle, DMPlexRefine_Triangle, NULL, 1));
48c0517cd5SMatthew G. Knepley #endif
49c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
509566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("ctetgen", DMPlexGenerate_CTetgen, DMPlexRefine_CTetgen, NULL, 2));
51c0517cd5SMatthew G. Knepley #endif
52c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_TETGEN)
539566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("tetgen", DMPlexGenerate_Tetgen, DMPlexRefine_Tetgen, NULL, 2));
54c0517cd5SMatthew G. Knepley #endif
55c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_PRAGMATIC)
569566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("pragmatic", NULL, NULL, DMAdaptMetric_Pragmatic_Plex, -1));
57c0517cd5SMatthew G. Knepley #endif
58c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_MMG)
599566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("mmg", NULL, NULL, DMAdaptMetric_Mmg_Plex, -1));
60c0517cd5SMatthew G. Knepley #endif
61c0517cd5SMatthew G. Knepley #if defined(PETSC_HAVE_PARMMG)
629566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("parmmg", NULL, NULL, DMAdaptMetric_ParMmg_Plex, -1));
63c0517cd5SMatthew G. Knepley #endif
649566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("cellrefiner", NULL, NULL, DMPlexTransformAdaptLabel, -1));
659566063dSJacob Faibussowitsch   PetscCall(DMGenerateRegister("forest", NULL, NULL, DMAdaptLabel_Forest, -1));
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67c0517cd5SMatthew G. Knepley }
68c0517cd5SMatthew G. Knepley 
69c0517cd5SMatthew G. Knepley /*@C
7016a05f60SBarry Smith   DMGenerateRegister -  Adds a grid generator to `DM`
71c0517cd5SMatthew G. Knepley 
72c0517cd5SMatthew G. Knepley   Not Collective
73c0517cd5SMatthew G. Knepley 
74c0517cd5SMatthew G. Knepley   Input Parameters:
7520f4b53cSBarry Smith + sname - name of a new user-defined grid generator
76c0517cd5SMatthew G. Knepley . fnc   - generator function
77c0517cd5SMatthew G. Knepley . rfnc  - refinement function
78c0517cd5SMatthew G. Knepley . alfnc - adapt by label function
79c0517cd5SMatthew G. Knepley - dim   - dimension of boundary of domain
80c0517cd5SMatthew G. Knepley 
8160225df5SJacob Faibussowitsch   Example Usage:
82c0517cd5SMatthew G. Knepley .vb
83c0517cd5SMatthew G. Knepley    DMGenerateRegister("my_generator", MyGeneratorCreate, MyGeneratorRefiner, MyGeneratorAdaptor, dim);
84c0517cd5SMatthew G. Knepley .ve
85c0517cd5SMatthew G. Knepley 
86c0517cd5SMatthew G. Knepley   Then, your generator can be chosen with the procedural interface via
87c0517cd5SMatthew G. Knepley $     DMGenerate(dm, "my_generator",...)
88c0517cd5SMatthew G. Knepley   or at runtime via the option
89c0517cd5SMatthew G. Knepley $     -dm_generator my_generator
90c0517cd5SMatthew G. Knepley 
91c0517cd5SMatthew G. Knepley   Level: advanced
92c0517cd5SMatthew G. Knepley 
9320f4b53cSBarry Smith   Note:
9420f4b53cSBarry Smith   `DMGenerateRegister()` may be called multiple times to add several user-defined generators
95c0517cd5SMatthew G. Knepley 
9620f4b53cSBarry Smith .seealso: `DM`, `DMGenerateRegisterAll()`, `DMPlexGenerate()`, `DMGenerateRegisterDestroy()`
97c0517cd5SMatthew G. Knepley @*/
98d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGenerateRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscBool, DM *), PetscErrorCode (*rfnc)(DM, PetscReal *, DM *), PetscErrorCode (*alfnc)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt dim)
99d71ae5a4SJacob Faibussowitsch {
100c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList entry;
101c0517cd5SMatthew G. Knepley 
102c0517cd5SMatthew G. Knepley   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(PetscNew(&entry));
1049566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(sname, &entry->name));
105c0517cd5SMatthew G. Knepley   entry->generate = fnc;
106c0517cd5SMatthew G. Knepley   entry->refine   = rfnc;
107c0517cd5SMatthew G. Knepley   entry->adapt    = alfnc;
108c0517cd5SMatthew G. Knepley   entry->dim      = dim;
109c0517cd5SMatthew G. Knepley   entry->next     = NULL;
110c0517cd5SMatthew G. Knepley   if (!DMGenerateList) DMGenerateList = entry;
111c0517cd5SMatthew G. Knepley   else {
112c0517cd5SMatthew G. Knepley     DMGeneratorFunctionList fl = DMGenerateList;
113c0517cd5SMatthew G. Knepley     while (fl->next) fl = fl->next;
114c0517cd5SMatthew G. Knepley     fl->next = entry;
115c0517cd5SMatthew G. Knepley   }
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117c0517cd5SMatthew G. Knepley }
118c0517cd5SMatthew G. Knepley 
119c0517cd5SMatthew G. Knepley extern PetscBool DMGenerateRegisterAllCalled;
120c0517cd5SMatthew G. Knepley 
121d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGenerateRegisterDestroy(void)
122d71ae5a4SJacob Faibussowitsch {
123c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList next, fl;
124c0517cd5SMatthew G. Knepley 
125c0517cd5SMatthew G. Knepley   PetscFunctionBegin;
126c0517cd5SMatthew G. Knepley   next = fl = DMGenerateList;
127c0517cd5SMatthew G. Knepley   while (next) {
128c0517cd5SMatthew G. Knepley     next = fl ? fl->next : NULL;
1299566063dSJacob Faibussowitsch     if (fl) PetscCall(PetscFree(fl->name));
1309566063dSJacob Faibussowitsch     PetscCall(PetscFree(fl));
131c0517cd5SMatthew G. Knepley     fl = next;
132c0517cd5SMatthew G. Knepley   }
133c0517cd5SMatthew G. Knepley   DMGenerateList              = NULL;
134c0517cd5SMatthew G. Knepley   DMGenerateRegisterAllCalled = PETSC_FALSE;
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136c0517cd5SMatthew G. Knepley }
137c0517cd5SMatthew G. Knepley 
138c0517cd5SMatthew G. Knepley /*@C
13920f4b53cSBarry Smith   DMAdaptLabel - Adapt a `DM` based on a `DMLabel` with values interpreted as coarsening and refining flags.  Specific implementations of `DM` maybe have
14020f4b53cSBarry Smith   specialized flags, but all implementations should accept flag values `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and,
14120f4b53cSBarry Smith   `DM_ADAPT_COARSEN`.
142c0517cd5SMatthew G. Knepley 
14320f4b53cSBarry Smith   Collective
144c0517cd5SMatthew G. Knepley 
14560225df5SJacob Faibussowitsch   Input Parameters:
14620f4b53cSBarry Smith + dm    - the pre-adaptation `DM` object
147c0517cd5SMatthew G. Knepley - label - label with the flags
148c0517cd5SMatthew G. Knepley 
14960225df5SJacob Faibussowitsch   Output Parameters:
15020f4b53cSBarry Smith . dmAdapt - the adapted `DM` object: may be `NULL` if an adapted `DM` could not be produced.
151c0517cd5SMatthew G. Knepley 
152c0517cd5SMatthew G. Knepley   Level: intermediate
153c0517cd5SMatthew G. Knepley 
15420f4b53cSBarry Smith .seealso: `DM`, `DMAdaptMetric()`, `DMCoarsen()`, `DMRefine()`
155c0517cd5SMatthew G. Knepley @*/
156d71ae5a4SJacob Faibussowitsch PetscErrorCode DMAdaptLabel(DM dm, DMLabel label, DM *dmAdapt)
157d71ae5a4SJacob Faibussowitsch {
158c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList fl;
159c0517cd5SMatthew G. Knepley   char                    adaptname[PETSC_MAX_PATH_LEN];
160c0517cd5SMatthew G. Knepley   const char             *name;
161c0517cd5SMatthew G. Knepley   PetscInt                dim;
162c0517cd5SMatthew G. Knepley   PetscBool               flg, isForest, found = PETSC_FALSE;
163c0517cd5SMatthew G. Knepley 
164c0517cd5SMatthew G. Knepley   PetscFunctionBegin;
165c0517cd5SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
166292bffcbSToby Isaac   if (label) PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 2);
1674f572ea9SToby Isaac   PetscAssertPointer(dmAdapt, 3);
168c0517cd5SMatthew G. Knepley   *dmAdapt = NULL;
1699566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1709566063dSJacob Faibussowitsch   PetscCall(DMIsForest(dm, &isForest));
171c0517cd5SMatthew G. Knepley   name = isForest ? "forest" : "cellrefiner";
1729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", adaptname, sizeof(adaptname), &flg));
173c0517cd5SMatthew G. Knepley   if (flg) name = adaptname;
174c0517cd5SMatthew G. Knepley 
175c0517cd5SMatthew G. Knepley   fl = DMGenerateList;
176c0517cd5SMatthew G. Knepley   while (fl) {
1779566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(fl->name, name, &flg));
178c0517cd5SMatthew G. Knepley     if (flg) {
1799566063dSJacob Faibussowitsch       PetscCall((*fl->adapt)(dm, NULL, label, NULL, dmAdapt));
180c0517cd5SMatthew G. Knepley       found = PETSC_TRUE;
181c0517cd5SMatthew G. Knepley     }
182c0517cd5SMatthew G. Knepley     fl = fl->next;
183c0517cd5SMatthew G. Knepley   }
1847a8be351SBarry Smith   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);
185c0517cd5SMatthew G. Knepley   if (*dmAdapt) {
186c0517cd5SMatthew G. Knepley     (*dmAdapt)->prealloc_only = dm->prealloc_only; /* maybe this should go .... */
1879566063dSJacob Faibussowitsch     PetscCall(PetscFree((*dmAdapt)->vectype));
1889566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*dmAdapt)->vectype));
1899566063dSJacob Faibussowitsch     PetscCall(PetscFree((*dmAdapt)->mattype));
1909566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*dmAdapt)->mattype));
191c0517cd5SMatthew G. Knepley   }
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
193c0517cd5SMatthew G. Knepley }
194c0517cd5SMatthew G. Knepley 
195c0517cd5SMatthew G. Knepley /*@C
196c0517cd5SMatthew G. Knepley   DMAdaptMetric - Generates a mesh adapted to the specified metric field.
197c0517cd5SMatthew G. Knepley 
198c0517cd5SMatthew G. Knepley   Input Parameters:
199c0517cd5SMatthew G. Knepley + dm      - The DM object
200c0517cd5SMatthew G. Knepley . metric  - The metric to which the mesh is adapted, defined vertex-wise.
2019fe9e680SJoe Wallwork . 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_".
2029fe9e680SJoe Wallwork - 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_".
203c0517cd5SMatthew G. Knepley 
204c0517cd5SMatthew G. Knepley   Output Parameter:
205c0517cd5SMatthew G. Knepley . dmAdapt - Pointer to the DM object containing the adapted mesh
206c0517cd5SMatthew G. Knepley 
207*a4e35b19SJacob Faibussowitsch   Note:
208*a4e35b19SJacob Faibussowitsch   The label in the adapted mesh will be registered under the name of the input `DMLabel` object
209c0517cd5SMatthew G. Knepley 
210c0517cd5SMatthew G. Knepley   Level: advanced
211c0517cd5SMatthew G. Knepley 
212db781477SPatrick Sanan .seealso: `DMAdaptLabel()`, `DMCoarsen()`, `DMRefine()`
213c0517cd5SMatthew G. Knepley @*/
214d71ae5a4SJacob Faibussowitsch PetscErrorCode DMAdaptMetric(DM dm, Vec metric, DMLabel bdLabel, DMLabel rgLabel, DM *dmAdapt)
215d71ae5a4SJacob Faibussowitsch {
216c0517cd5SMatthew G. Knepley   DMGeneratorFunctionList fl;
217c0517cd5SMatthew G. Knepley   char                    adaptname[PETSC_MAX_PATH_LEN];
21806716a2aSJoe Wallwork   const char             *name;
21906716a2aSJoe Wallwork   const char *const       adaptors[3] = {"pragmatic", "mmg", "parmmg"};
220c0517cd5SMatthew G. Knepley   PetscInt                dim;
221c0517cd5SMatthew G. Knepley   PetscBool               flg, found = PETSC_FALSE;
222c0517cd5SMatthew G. Knepley 
223c0517cd5SMatthew G. Knepley   PetscFunctionBegin;
224c0517cd5SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
225c0517cd5SMatthew G. Knepley   PetscValidHeaderSpecific(metric, VEC_CLASSID, 2);
226292bffcbSToby Isaac   if (bdLabel) PetscValidHeaderSpecific(bdLabel, DMLABEL_CLASSID, 3);
227292bffcbSToby Isaac   if (rgLabel) PetscValidHeaderSpecific(rgLabel, DMLABEL_CLASSID, 4);
2284f572ea9SToby Isaac   PetscAssertPointer(dmAdapt, 5);
229c0517cd5SMatthew G. Knepley   *dmAdapt = NULL;
2309566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", adaptname, sizeof(adaptname), &flg));
23206716a2aSJoe Wallwork 
23306716a2aSJoe Wallwork   /* Default to Mmg in serial and ParMmg in parallel */
234c0517cd5SMatthew G. Knepley   if (flg) name = adaptname;
23506716a2aSJoe Wallwork   else {
23606716a2aSJoe Wallwork     MPI_Comm    comm;
23706716a2aSJoe Wallwork     PetscMPIInt size;
23806716a2aSJoe Wallwork 
2399566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
24106716a2aSJoe Wallwork     if (size == 1) name = adaptors[1];
24206716a2aSJoe Wallwork     else name = adaptors[2];
24306716a2aSJoe Wallwork   }
244c0517cd5SMatthew G. Knepley 
245c0517cd5SMatthew G. Knepley   fl = DMGenerateList;
246c0517cd5SMatthew G. Knepley   while (fl) {
2479566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(fl->name, name, &flg));
248c0517cd5SMatthew G. Knepley     if (flg) {
2499566063dSJacob Faibussowitsch       PetscCall((*fl->adapt)(dm, metric, bdLabel, rgLabel, dmAdapt));
250c0517cd5SMatthew G. Knepley       found = PETSC_TRUE;
251c0517cd5SMatthew G. Knepley     }
252c0517cd5SMatthew G. Knepley     fl = fl->next;
253c0517cd5SMatthew G. Knepley   }
2547a8be351SBarry Smith   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);
255c0517cd5SMatthew G. Knepley   if (*dmAdapt) {
256c0517cd5SMatthew G. Knepley     (*dmAdapt)->prealloc_only = dm->prealloc_only; /* maybe this should go .... */
2579566063dSJacob Faibussowitsch     PetscCall(PetscFree((*dmAdapt)->vectype));
2589566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*dmAdapt)->vectype));
2599566063dSJacob Faibussowitsch     PetscCall(PetscFree((*dmAdapt)->mattype));
2609566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*dmAdapt)->mattype));
261c0517cd5SMatthew G. Knepley   }
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263c0517cd5SMatthew G. Knepley }
264