xref: /petsc/src/snes/utils/dm/dmadapt.c (revision 8a7d4057d9226490dba4e1a062f54f84e7d90861)
1 #include <petscdmadaptor.h> /*I "petscdmadaptor.h" I*/
2 #include <petscdmplex.h>
3 #include <petscdmforest.h>
4 #include <petscds.h>
5 #include <petscblaslapack.h>
6 #include <petscsnes.h>
7 #include <petscdraw.h>
8 
9 #include <petsc/private/dmadaptorimpl.h>
10 #include <petsc/private/dmpleximpl.h>
11 #include <petsc/private/petscfeimpl.h>
12 
13 PetscClassId DMADAPTOR_CLASSID;
14 
15 PetscFunctionList DMAdaptorList              = NULL;
16 PetscBool         DMAdaptorRegisterAllCalled = PETSC_FALSE;
17 
18 PetscFunctionList DMAdaptorMonitorList              = NULL;
19 PetscFunctionList DMAdaptorMonitorCreateList        = NULL;
20 PetscFunctionList DMAdaptorMonitorDestroyList       = NULL;
21 PetscBool         DMAdaptorMonitorRegisterAllCalled = PETSC_FALSE;
22 
23 const char *const DMAdaptationCriteria[] = {"NONE", "REFINE", "LABEL", "METRIC", "DMAdaptationCriterion", "DM_ADAPTATION_", NULL};
24 
25 /*@C
26   DMAdaptorRegister - Adds a new adaptor component implementation
27 
28   Not Collective
29 
30   Input Parameters:
31 + name        - The name of a new user-defined creation routine
32 - create_func - The creation routine
33 
34   Example Usage:
35 .vb
36   DMAdaptorRegister("my_adaptor", MyAdaptorCreate);
37 .ve
38 
39   Then, your adaptor type can be chosen with the procedural interface via
40 .vb
41   DMAdaptorCreate(MPI_Comm, DMAdaptor *);
42   DMAdaptorSetType(DMAdaptor, "my_adaptor");
43 .ve
44   or at runtime via the option
45 .vb
46   -adaptor_type my_adaptor
47 .ve
48 
49   Level: advanced
50 
51   Note:
52   `DMAdaptorRegister()` may be called multiple times to add several user-defined adaptors
53 
54 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptor`, `DMAdaptorRegisterAll()`, `DMAdaptorRegisterDestroy()`
55 @*/
56 PetscErrorCode DMAdaptorRegister(const char name[], PetscErrorCode (*create_func)(DMAdaptor))
57 {
58   PetscFunctionBegin;
59   PetscCall(DMInitializePackage());
60   PetscCall(PetscFunctionListAdd(&DMAdaptorList, name, create_func));
61   PetscFunctionReturn(PETSC_SUCCESS);
62 }
63 
64 PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Gradient(DMAdaptor);
65 PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Flux(DMAdaptor);
66 
67 /*@C
68   DMAdaptorRegisterAll - Registers all of the adaptor components in the `DM` package.
69 
70   Not Collective
71 
72   Level: advanced
73 
74 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptorType`, `DMRegisterAll()`, `DMAdaptorRegisterDestroy()`
75 @*/
76 PetscErrorCode DMAdaptorRegisterAll(void)
77 {
78   PetscFunctionBegin;
79   if (DMAdaptorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
80   DMAdaptorRegisterAllCalled = PETSC_TRUE;
81 
82   PetscCall(DMAdaptorRegister(DMADAPTORGRADIENT, DMAdaptorCreate_Gradient));
83   PetscCall(DMAdaptorRegister(DMADAPTORFLUX, DMAdaptorCreate_Flux));
84   PetscFunctionReturn(PETSC_SUCCESS);
85 }
86 
87 /*@C
88   DMAdaptorRegisterDestroy - This function destroys the registered `DMAdaptorType`. It is called from `PetscFinalize()`.
89 
90   Not collective
91 
92   Level: developer
93 
94 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptorRegisterAll()`, `DMAdaptorType`, `PetscFinalize()`
95 @*/
96 PetscErrorCode DMAdaptorRegisterDestroy(void)
97 {
98   PetscFunctionBegin;
99   PetscCall(PetscFunctionListDestroy(&DMAdaptorList));
100   DMAdaptorRegisterAllCalled = PETSC_FALSE;
101   PetscFunctionReturn(PETSC_SUCCESS);
102 }
103 
104 static PetscErrorCode DMAdaptorMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[])
105 {
106   PetscFunctionBegin;
107   PetscCall(PetscStrncpy(key, name, PETSC_MAX_PATH_LEN));
108   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
109   PetscCall(PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN));
110   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
111   PetscCall(PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN));
112   PetscFunctionReturn(PETSC_SUCCESS);
113 }
114 
115 /*@C
116   DMAdaptorMonitorRegister -  Registers a mesh adaptation monitor routine that may be accessed with `DMAdaptorMonitorSetFromOptions()`
117 
118   Not Collective
119 
120   Input Parameters:
121 + name    - name of a new monitor routine
122 . vtype   - A `PetscViewerType` for the output
123 . format  - A `PetscViewerFormat` for the output
124 . monitor - Monitor routine
125 . create  - Creation routine, or `NULL`
126 - destroy - Destruction routine, or `NULL`
127 
128   Level: advanced
129 
130   Note:
131   `DMAdaptorMonitorRegister()` may be called multiple times to add several user-defined monitors.
132 
133   Example Usage:
134 .vb
135   DMAdaptorMonitorRegister("my_monitor", PETSCVIEWERASCII, PETSC_VIEWER_ASCII_INFO_DETAIL, MyMonitor, NULL, NULL);
136 .ve
137 
138   Then, your monitor can be chosen with the procedural interface via
139 .vb
140   DMAdaptorMonitorSetFromOptions(ksp, "-adaptor_monitor_my_monitor", "my_monitor", NULL)
141 .ve
142   or at runtime via the option `-adaptor_monitor_my_monitor`
143 
144 .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorRegisterAll()`, `DMAdaptorMonitorSetFromOptions()`
145 @*/
146 PetscErrorCode DMAdaptorMonitorRegister(const char name[], PetscViewerType vtype, PetscViewerFormat format, PetscErrorCode (*monitor)(DMAdaptor, PetscInt, DM, DM, PetscInt, PetscReal[], Vec, PetscViewerAndFormat *), PetscErrorCode (*create)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*destroy)(PetscViewerAndFormat **))
147 {
148   char key[PETSC_MAX_PATH_LEN];
149 
150   PetscFunctionBegin;
151   PetscCall(SNESInitializePackage());
152   PetscCall(DMAdaptorMonitorMakeKey_Internal(name, vtype, format, key));
153   PetscCall(PetscFunctionListAdd(&DMAdaptorMonitorList, key, monitor));
154   if (create) PetscCall(PetscFunctionListAdd(&DMAdaptorMonitorCreateList, key, create));
155   if (destroy) PetscCall(PetscFunctionListAdd(&DMAdaptorMonitorDestroyList, key, destroy));
156   PetscFunctionReturn(PETSC_SUCCESS);
157 }
158 
159 /*@C
160   DMAdaptorMonitorRegisterDestroy - This function destroys the registered monitors for `DMAdaptor`. It is called from `PetscFinalize()`.
161 
162   Not collective
163 
164   Level: developer
165 
166 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptorMonitorRegisterAll()`, `DMAdaptor`, `PetscFinalize()`
167 @*/
168 PetscErrorCode DMAdaptorMonitorRegisterDestroy(void)
169 {
170   PetscFunctionBegin;
171   PetscCall(PetscFunctionListDestroy(&DMAdaptorMonitorList));
172   PetscCall(PetscFunctionListDestroy(&DMAdaptorMonitorCreateList));
173   PetscCall(PetscFunctionListDestroy(&DMAdaptorMonitorDestroyList));
174   DMAdaptorMonitorRegisterAllCalled = PETSC_FALSE;
175   PetscFunctionReturn(PETSC_SUCCESS);
176 }
177 
178 /*@
179   DMAdaptorCreate - Create a `DMAdaptor` object. Its purpose is to construct a adaptation `DMLabel` or metric `Vec` that can be used to modify the `DM`.
180 
181   Collective
182 
183   Input Parameter:
184 . comm - The communicator for the `DMAdaptor` object
185 
186   Output Parameter:
187 . adaptor - The `DMAdaptor` object
188 
189   Level: beginner
190 
191 .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorDestroy()`, `DMAdaptorAdapt()`, `PetscConvEst`, `PetscConvEstCreate()`
192 @*/
193 PetscErrorCode DMAdaptorCreate(MPI_Comm comm, DMAdaptor *adaptor)
194 {
195   VecTaggerBox refineBox, coarsenBox;
196 
197   PetscFunctionBegin;
198   PetscAssertPointer(adaptor, 2);
199   PetscCall(PetscSysInitializePackage());
200 
201   PetscCall(PetscHeaderCreate(*adaptor, DMADAPTOR_CLASSID, "DMAdaptor", "DM Adaptor", "DMAdaptor", comm, DMAdaptorDestroy, DMAdaptorView));
202   (*adaptor)->adaptCriterion   = DM_ADAPTATION_NONE;
203   (*adaptor)->numSeq           = 1;
204   (*adaptor)->Nadapt           = -1;
205   (*adaptor)->refinementFactor = 2.0;
206   refineBox.min = refineBox.max = PETSC_MAX_REAL;
207   PetscCall(VecTaggerCreate(PetscObjectComm((PetscObject)*adaptor), &(*adaptor)->refineTag));
208   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(*adaptor)->refineTag, "refine_"));
209   PetscCall(VecTaggerSetType((*adaptor)->refineTag, VECTAGGERABSOLUTE));
210   PetscCall(VecTaggerAbsoluteSetBox((*adaptor)->refineTag, &refineBox));
211   coarsenBox.min = coarsenBox.max = PETSC_MAX_REAL;
212   PetscCall(VecTaggerCreate(PetscObjectComm((PetscObject)*adaptor), &(*adaptor)->coarsenTag));
213   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(*adaptor)->coarsenTag, "coarsen_"));
214   PetscCall(VecTaggerSetType((*adaptor)->coarsenTag, VECTAGGERABSOLUTE));
215   PetscCall(VecTaggerAbsoluteSetBox((*adaptor)->coarsenTag, &coarsenBox));
216   PetscFunctionReturn(PETSC_SUCCESS);
217 }
218 
219 /*@
220   DMAdaptorDestroy - Destroys a `DMAdaptor` object
221 
222   Collective
223 
224   Input Parameter:
225 . adaptor - The `DMAdaptor` object
226 
227   Level: beginner
228 
229 .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
230 @*/
231 PetscErrorCode DMAdaptorDestroy(DMAdaptor *adaptor)
232 {
233   PetscFunctionBegin;
234   if (!*adaptor) PetscFunctionReturn(PETSC_SUCCESS);
235   PetscValidHeaderSpecific(*adaptor, DMADAPTOR_CLASSID, 1);
236   if (--((PetscObject)*adaptor)->refct > 0) {
237     *adaptor = NULL;
238     PetscFunctionReturn(PETSC_SUCCESS);
239   }
240   PetscCall(VecTaggerDestroy(&(*adaptor)->refineTag));
241   PetscCall(VecTaggerDestroy(&(*adaptor)->coarsenTag));
242   PetscCall(PetscFree2((*adaptor)->exactSol, (*adaptor)->exactCtx));
243   PetscCall(DMAdaptorMonitorCancel(*adaptor));
244   PetscCall(PetscHeaderDestroy(adaptor));
245   PetscFunctionReturn(PETSC_SUCCESS);
246 }
247 
248 /*@
249   DMAdaptorSetType - Sets the particular implementation for a adaptor.
250 
251   Collective
252 
253   Input Parameters:
254 + adaptor - The `DMAdaptor`
255 - method  - The name of the adaptor type
256 
257   Options Database Key:
258 . -adaptor_type <type> - Sets the adaptor type; see `DMAdaptorType`
259 
260   Level: intermediate
261 
262 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptor`, `DMAdaptorType`, `DMAdaptorGetType()`, `DMAdaptorCreate()`
263 @*/
264 PetscErrorCode DMAdaptorSetType(DMAdaptor adaptor, DMAdaptorType method)
265 {
266   PetscErrorCode (*r)(DMAdaptor);
267   PetscBool match;
268 
269   PetscFunctionBegin;
270   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
271   PetscCall(PetscObjectTypeCompare((PetscObject)adaptor, method, &match));
272   if (match) PetscFunctionReturn(PETSC_SUCCESS);
273 
274   PetscCall(DMAdaptorRegisterAll());
275   PetscCall(PetscFunctionListFind(DMAdaptorList, method, &r));
276   PetscCheck(r, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DMAdaptor type: %s", method);
277 
278   PetscTryTypeMethod(adaptor, destroy);
279   PetscCall(PetscMemzero(adaptor->ops, sizeof(*adaptor->ops)));
280   PetscCall(PetscObjectChangeTypeName((PetscObject)adaptor, method));
281   PetscCall((*r)(adaptor));
282   PetscFunctionReturn(PETSC_SUCCESS);
283 }
284 
285 /*@
286   DMAdaptorGetType - Gets the type name (as a string) from the adaptor.
287 
288   Not Collective
289 
290   Input Parameter:
291 . adaptor - The `DMAdaptor`
292 
293   Output Parameter:
294 . type - The `DMAdaptorType` name
295 
296   Level: intermediate
297 
298 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptor`, `DMAdaptorType`, `DMAdaptorSetType()`, `DMAdaptorCreate()`
299 @*/
300 PetscErrorCode DMAdaptorGetType(DMAdaptor adaptor, DMAdaptorType *type)
301 {
302   PetscFunctionBegin;
303   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
304   PetscAssertPointer(type, 2);
305   PetscCall(DMAdaptorRegisterAll());
306   *type = ((PetscObject)adaptor)->type_name;
307   PetscFunctionReturn(PETSC_SUCCESS);
308 }
309 
310 static PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
311 {
312   PetscFunctionBegin;
313   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
314   (*vf)->data = ctx;
315   PetscFunctionReturn(PETSC_SUCCESS);
316 }
317 
318 /*@C
319   DMAdaptorMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
320   the error etc.
321 
322   Logically Collective
323 
324   Input Parameters:
325 + adaptor        - the `DMAdaptor`
326 . monitor        - pointer to function (if this is `NULL`, it turns off monitoring
327 . ctx            - [optional] context for private data for the monitor routine (use `NULL` if no context is needed)
328 - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
329 
330   Calling sequence of `monitor`:
331 + adaptor - the `DMAdaptor`
332 . it      - iteration number
333 . odm     - the original `DM`
334 . adm     - the adapted `DM`
335 . Nf      - number of fields
336 . enorms  - (estimated) 2-norm of the error for each field
337 . error   - `Vec` of cellwise errors
338 - ctx     - optional monitoring context, as set by `DMAdaptorMonitorSet()`
339 
340   Calling sequence of `monitordestroy`:
341 . ctx - optional monitoring context, as set by `DMAdaptorMonitorSet()`
342 
343   Options Database Keys:
344 + -adaptor_monitor_size                - sets `DMAdaptorMonitorSize()`
345 . -adaptor_monitor_error               - sets `DMAdaptorMonitorError()`
346 . -adaptor_monitor_error draw          - sets `DMAdaptorMonitorErrorDraw()` and plots error
347 . -adaptor_monitor_error draw::draw_lg - sets `DMAdaptorMonitorErrorDrawLG()` and plots error
348 - -dm_adaptor_monitor_cancel           - Cancels all monitors that have been hardwired into a code by calls to `DMAdaptorMonitorSet()`, but does not cancel those set via the options database.
349 
350   Level: beginner
351 
352 .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptor`
353 @*/
354 PetscErrorCode DMAdaptorMonitorSet(DMAdaptor adaptor, PetscErrorCode (*monitor)(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, void *ctx), void *ctx, PetscErrorCode (*monitordestroy)(void **ctx))
355 {
356   PetscBool identical;
357 
358   PetscFunctionBegin;
359   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
360   for (PetscInt i = 0; i < adaptor->numbermonitors; i++) {
361     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))monitor, ctx, monitordestroy, (PetscErrorCode(*)(void))adaptor->monitor[i], adaptor->monitorcontext[i], adaptor->monitordestroy[i], &identical));
362     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
363   }
364   PetscCheck(adaptor->numbermonitors < MAXDMADAPTORMONITORS, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_OUTOFRANGE, "Too many DMAdaptor monitors set");
365   adaptor->monitor[adaptor->numbermonitors]          = monitor;
366   adaptor->monitordestroy[adaptor->numbermonitors]   = monitordestroy;
367   adaptor->monitorcontext[adaptor->numbermonitors++] = (void *)ctx;
368   PetscFunctionReturn(PETSC_SUCCESS);
369 }
370 
371 /*@
372   DMAdaptorMonitorCancel - Clears all monitors for a `DMAdaptor` object.
373 
374   Logically Collective
375 
376   Input Parameter:
377 . adaptor - the `DMAdaptor`
378 
379   Options Database Key:
380 . -dm_adaptor_monitor_cancel - Cancels all monitors that have been hardwired into a code by calls to `DMAdaptorMonitorSet()`, but does not cancel those set via the options database.
381 
382   Level: intermediate
383 
384 .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptorMonitorSet()`, `DMAdaptor`
385 @*/
386 PetscErrorCode DMAdaptorMonitorCancel(DMAdaptor adaptor)
387 {
388   PetscFunctionBegin;
389   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
390   for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) {
391     if (adaptor->monitordestroy[i]) PetscCall((*adaptor->monitordestroy[i])(&adaptor->monitorcontext[i]));
392   }
393   adaptor->numbermonitors = 0;
394   PetscFunctionReturn(PETSC_SUCCESS);
395 }
396 
397 /*@C
398   DMAdaptorMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user in the options database
399 
400   Collective
401 
402   Input Parameters:
403 + adaptor - `DMadaptor` object you wish to monitor
404 . opt     - the command line option for this monitor
405 . name    - the monitor type one is seeking
406 - ctx     - An optional user context for the monitor, or `NULL`
407 
408   Level: developer
409 
410 .seealso: [](ch_snes), `DMAdaptorMonitorRegister()`, `DMAdaptorMonitorSet()`, `PetscOptionsGetViewer()`
411 @*/
412 PetscErrorCode DMAdaptorMonitorSetFromOptions(DMAdaptor adaptor, const char opt[], const char name[], void *ctx)
413 {
414   PetscErrorCode (*mfunc)(DMAdaptor, PetscInt, DM, DM, PetscInt, PetscReal[], Vec, void *);
415   PetscErrorCode (*cfunc)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
416   PetscErrorCode (*dfunc)(PetscViewerAndFormat **);
417   PetscViewerAndFormat *vf;
418   PetscViewer           viewer;
419   PetscViewerFormat     format;
420   PetscViewerType       vtype;
421   char                  key[PETSC_MAX_PATH_LEN];
422   PetscBool             flg;
423   const char           *prefix = NULL;
424 
425   PetscFunctionBegin;
426   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
427   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)adaptor), ((PetscObject)adaptor)->options, prefix, opt, &viewer, &format, &flg));
428   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
429 
430   PetscCall(PetscViewerGetType(viewer, &vtype));
431   PetscCall(DMAdaptorMonitorMakeKey_Internal(name, vtype, format, key));
432   PetscCall(PetscFunctionListFind(DMAdaptorMonitorList, key, &mfunc));
433   PetscCall(PetscFunctionListFind(DMAdaptorMonitorCreateList, key, &cfunc));
434   PetscCall(PetscFunctionListFind(DMAdaptorMonitorDestroyList, key, &dfunc));
435   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
436   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
437 
438   PetscCall((*cfunc)(viewer, format, ctx, &vf));
439   PetscCall(PetscViewerDestroy(&viewer));
440   PetscCall(DMAdaptorMonitorSet(adaptor, mfunc, vf, (PetscErrorCode(*)(void **))dfunc));
441   PetscFunctionReturn(PETSC_SUCCESS);
442 }
443 
444 /*@
445   DMAdaptorSetFromOptions - Sets properties of a `DMAdaptor` object from values in the options database
446 
447   Collective
448 
449   Input Parameter:
450 . adaptor - The `DMAdaptor` object
451 
452   Options Database Keys:
453 + -adaptor_monitor_size                - Monitor the mesh size
454 . -adaptor_monitor_error               - Monitor the solution error
455 . -adaptor_sequence_num <num>          - Number of adaptations to generate an optimal grid
456 . -adaptor_target_num <num>            - Set the target number of vertices N_adapt, -1 for automatic determination
457 . -adaptor_refinement_factor <r>       - Set r such that N_adapt = r^dim N_orig
458 - -adaptor_mixed_setup_function <func> - Set the fnction func that sets up the mixed problem
459 
460   Level: beginner
461 
462 .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
463 @*/
464 PetscErrorCode DMAdaptorSetFromOptions(DMAdaptor adaptor)
465 {
466   char                  typeName[PETSC_MAX_PATH_LEN];
467   const char           *defName = DMADAPTORGRADIENT;
468   char                  funcname[PETSC_MAX_PATH_LEN];
469   DMAdaptationCriterion criterion = DM_ADAPTATION_NONE;
470   PetscBool             flg;
471 
472   PetscFunctionBegin;
473   PetscObjectOptionsBegin((PetscObject)adaptor);
474   PetscCall(PetscOptionsFList("-adaptor_type", "DMAdaptor", "DMAdaptorSetType", DMAdaptorList, defName, typeName, 1024, &flg));
475   if (flg) PetscCall(DMAdaptorSetType(adaptor, typeName));
476   else if (!((PetscObject)adaptor)->type_name) PetscCall(DMAdaptorSetType(adaptor, defName));
477   PetscCall(PetscOptionsEnum("-adaptor_criterion", "Criterion used to drive adaptation", "", DMAdaptationCriteria, (PetscEnum)criterion, (PetscEnum *)&criterion, &flg));
478   if (flg) PetscCall(DMAdaptorSetCriterion(adaptor, criterion));
479   PetscCall(PetscOptionsInt("-adaptor_sequence_num", "Number of adaptations to generate an optimal grid", "DMAdaptorSetSequenceLength", adaptor->numSeq, &adaptor->numSeq, NULL));
480   PetscCall(PetscOptionsInt("-adaptor_target_num", "Set the target number of vertices N_adapt, -1 for automatic determination", "DMAdaptor", adaptor->Nadapt, &adaptor->Nadapt, NULL));
481   PetscCall(PetscOptionsReal("-adaptor_refinement_factor", "Set r such that N_adapt = r^dim N_orig", "DMAdaptor", adaptor->refinementFactor, &adaptor->refinementFactor, NULL));
482   PetscCall(PetscOptionsString("-adaptor_mixed_setup_function", "Function to setup the mixed problem", "DMAdaptorSetMixedSetupFunction", funcname, funcname, sizeof(funcname), &flg));
483   if (flg) {
484     PetscErrorCode (*setupFunc)(DMAdaptor, DM);
485 
486     PetscCall(PetscDLSym(NULL, funcname, (void **)&setupFunc));
487     PetscCheck(setupFunc, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
488     PetscCall(DMAdaptorSetMixedSetupFunction(adaptor, setupFunc));
489   }
490   PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_size", "size", adaptor));
491   PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_error", "error", adaptor));
492   PetscOptionsEnd();
493   PetscCall(VecTaggerSetFromOptions(adaptor->refineTag));
494   PetscCall(VecTaggerSetFromOptions(adaptor->coarsenTag));
495   PetscFunctionReturn(PETSC_SUCCESS);
496 }
497 
498 /*@
499   DMAdaptorView - Views a `DMAdaptor` object
500 
501   Collective
502 
503   Input Parameters:
504 + adaptor - The `DMAdaptor` object
505 - viewer  - The `PetscViewer` object
506 
507   Level: beginner
508 
509 .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
510 @*/
511 PetscErrorCode DMAdaptorView(DMAdaptor adaptor, PetscViewer viewer)
512 {
513   PetscFunctionBegin;
514   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adaptor, viewer));
515   PetscCall(PetscViewerASCIIPrintf(viewer, "DM Adaptor\n"));
516   PetscCall(PetscViewerASCIIPrintf(viewer, "  sequence length: %" PetscInt_FMT "\n", adaptor->numSeq));
517   PetscCall(VecTaggerView(adaptor->refineTag, viewer));
518   PetscCall(VecTaggerView(adaptor->coarsenTag, viewer));
519   PetscFunctionReturn(PETSC_SUCCESS);
520 }
521 
522 /*@
523   DMAdaptorGetSolver - Gets the solver used to produce discrete solutions
524 
525   Not Collective
526 
527   Input Parameter:
528 . adaptor - The `DMAdaptor` object
529 
530   Output Parameter:
531 . snes - The solver
532 
533   Level: intermediate
534 
535 .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
536 @*/
537 PetscErrorCode DMAdaptorGetSolver(DMAdaptor adaptor, SNES *snes)
538 {
539   PetscFunctionBegin;
540   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
541   PetscAssertPointer(snes, 2);
542   *snes = adaptor->snes;
543   PetscFunctionReturn(PETSC_SUCCESS);
544 }
545 
546 /*@
547   DMAdaptorSetSolver - Sets the solver used to produce discrete solutions
548 
549   Not Collective
550 
551   Input Parameters:
552 + adaptor - The `DMAdaptor` object
553 - snes    - The solver, this MUST have an attached `DM`/`PetscDS`, so that the exact solution can be computed
554 
555   Level: intermediate
556 
557 .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorGetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
558 @*/
559 PetscErrorCode DMAdaptorSetSolver(DMAdaptor adaptor, SNES snes)
560 {
561   PetscFunctionBegin;
562   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
563   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
564   adaptor->snes = snes;
565   PetscCall(SNESGetDM(adaptor->snes, &adaptor->idm));
566   PetscFunctionReturn(PETSC_SUCCESS);
567 }
568 
569 /*@
570   DMAdaptorGetSequenceLength - Gets the number of sequential adaptations used by an adapter
571 
572   Not Collective
573 
574   Input Parameter:
575 . adaptor - The `DMAdaptor` object
576 
577   Output Parameter:
578 . num - The number of adaptations
579 
580   Level: intermediate
581 
582 .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorSetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
583 @*/
584 PetscErrorCode DMAdaptorGetSequenceLength(DMAdaptor adaptor, PetscInt *num)
585 {
586   PetscFunctionBegin;
587   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
588   PetscAssertPointer(num, 2);
589   *num = adaptor->numSeq;
590   PetscFunctionReturn(PETSC_SUCCESS);
591 }
592 
593 /*@
594   DMAdaptorSetSequenceLength - Sets the number of sequential adaptations
595 
596   Not Collective
597 
598   Input Parameters:
599 + adaptor - The `DMAdaptor` object
600 - num     - The number of adaptations
601 
602   Level: intermediate
603 
604 .seealso: [](ch_dmbase), `DMAdaptorGetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
605 @*/
606 PetscErrorCode DMAdaptorSetSequenceLength(DMAdaptor adaptor, PetscInt num)
607 {
608   PetscFunctionBegin;
609   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
610   adaptor->numSeq = num;
611   PetscFunctionReturn(PETSC_SUCCESS);
612 }
613 
614 static PetscErrorCode DMAdaptorTransferSolution_Exact_Private(DMAdaptor adaptor, DM dm, Vec u, DM adm, Vec au, void *ctx)
615 {
616   PetscFunctionBeginUser;
617   PetscCall(DMProjectFunction(adm, 0.0, adaptor->exactSol, adaptor->exactCtx, INSERT_ALL_VALUES, au));
618   PetscFunctionReturn(PETSC_SUCCESS);
619 }
620 
621 /*@
622   DMAdaptorSetUp - After the solver is specified, creates data structures for controlling adaptivity
623 
624   Collective
625 
626   Input Parameter:
627 . adaptor - The `DMAdaptor` object
628 
629   Level: beginner
630 
631 .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
632 @*/
633 PetscErrorCode DMAdaptorSetUp(DMAdaptor adaptor)
634 {
635   PetscDS  prob;
636   PetscInt Nf, f;
637 
638   PetscFunctionBegin;
639   PetscCall(DMGetDS(adaptor->idm, &prob));
640   PetscCall(VecTaggerSetUp(adaptor->refineTag));
641   PetscCall(VecTaggerSetUp(adaptor->coarsenTag));
642   PetscCall(PetscDSGetNumFields(prob, &Nf));
643   PetscCall(PetscMalloc2(Nf, &adaptor->exactSol, Nf, &adaptor->exactCtx));
644   for (f = 0; f < Nf; ++f) {
645     PetscCall(PetscDSGetExactSolution(prob, f, &adaptor->exactSol[f], &adaptor->exactCtx[f]));
646     /* TODO Have a flag that forces projection rather than using the exact solution */
647     if (adaptor->exactSol[0]) PetscCall(DMAdaptorSetTransferFunction(adaptor, DMAdaptorTransferSolution_Exact_Private));
648   }
649   PetscTryTypeMethod(adaptor, setup);
650   PetscFunctionReturn(PETSC_SUCCESS);
651 }
652 
653 PetscErrorCode DMAdaptorGetTransferFunction(DMAdaptor adaptor, PetscErrorCode (**tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
654 {
655   PetscFunctionBegin;
656   *tfunc = adaptor->ops->transfersolution;
657   PetscFunctionReturn(PETSC_SUCCESS);
658 }
659 
660 PetscErrorCode DMAdaptorSetTransferFunction(DMAdaptor adaptor, PetscErrorCode (*tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
661 {
662   PetscFunctionBegin;
663   adaptor->ops->transfersolution = tfunc;
664   PetscFunctionReturn(PETSC_SUCCESS);
665 }
666 
667 static PetscErrorCode DMAdaptorPreAdapt(DMAdaptor adaptor, Vec locX)
668 {
669   DM           plex;
670   PetscDS      prob;
671   PetscObject  obj;
672   PetscClassId id;
673   PetscBool    isForest;
674 
675   PetscFunctionBegin;
676   PetscCall(DMConvert(adaptor->idm, DMPLEX, &plex));
677   PetscCall(DMGetDS(adaptor->idm, &prob));
678   PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
679   PetscCall(PetscObjectGetClassId(obj, &id));
680   PetscCall(DMIsForest(adaptor->idm, &isForest));
681   if (adaptor->adaptCriterion == DM_ADAPTATION_NONE) {
682     if (isForest) {
683       adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
684     }
685 #if defined(PETSC_HAVE_PRAGMATIC)
686     else {
687       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
688     }
689 #elif defined(PETSC_HAVE_MMG)
690     else {
691       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
692     }
693 #elif defined(PETSC_HAVE_PARMMG)
694     else {
695       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
696     }
697 #else
698     else {
699       adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
700     }
701 #endif
702   }
703   if (id == PETSCFV_CLASSID) {
704     adaptor->femType = PETSC_FALSE;
705   } else {
706     adaptor->femType = PETSC_TRUE;
707   }
708   if (adaptor->femType) {
709     /* Compute local solution bc */
710     PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
711   } else {
712     PetscFV      fvm = (PetscFV)obj;
713     PetscLimiter noneLimiter;
714     Vec          grad;
715 
716     PetscCall(PetscFVGetComputeGradients(fvm, &adaptor->computeGradient));
717     PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
718     /* Use no limiting when reconstructing gradients for adaptivity */
719     PetscCall(PetscFVGetLimiter(fvm, &adaptor->limiter));
720     PetscCall(PetscObjectReference((PetscObject)adaptor->limiter));
721     PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
722     PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
723     PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
724     /* Get FVM data */
725     PetscCall(DMPlexGetDataFVM(plex, fvm, &adaptor->cellGeom, &adaptor->faceGeom, &adaptor->gradDM));
726     PetscCall(VecGetDM(adaptor->cellGeom, &adaptor->cellDM));
727     PetscCall(VecGetArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
728     /* Compute local solution bc */
729     PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
730     /* Compute gradients */
731     PetscCall(DMCreateGlobalVector(adaptor->gradDM, &grad));
732     PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
733     PetscCall(DMGetLocalVector(adaptor->gradDM, &adaptor->cellGrad));
734     PetscCall(DMGlobalToLocalBegin(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
735     PetscCall(DMGlobalToLocalEnd(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
736     PetscCall(VecDestroy(&grad));
737     PetscCall(VecGetArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
738   }
739   PetscCall(DMDestroy(&plex));
740   PetscFunctionReturn(PETSC_SUCCESS);
741 }
742 
743 static PetscErrorCode DMAdaptorTransferSolution(DMAdaptor adaptor, DM dm, Vec x, DM adm, Vec ax)
744 {
745   PetscReal time = 0.0;
746   Mat       interp;
747   void     *ctx;
748 
749   PetscFunctionBegin;
750   PetscCall(DMGetApplicationContext(dm, &ctx));
751   if (adaptor->ops->transfersolution) PetscUseTypeMethod(adaptor, transfersolution, dm, x, adm, ax, ctx);
752   else {
753     switch (adaptor->adaptCriterion) {
754     case DM_ADAPTATION_LABEL:
755       PetscCall(DMForestTransferVec(dm, x, adm, ax, PETSC_TRUE, time));
756       break;
757     case DM_ADAPTATION_REFINE:
758     case DM_ADAPTATION_METRIC:
759       PetscCall(DMCreateInterpolation(dm, adm, &interp, NULL));
760       PetscCall(MatInterpolate(interp, x, ax));
761       PetscCall(DMInterpolate(dm, interp, adm));
762       PetscCall(MatDestroy(&interp));
763       break;
764     default:
765       SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_SUP, "No built-in projection for this adaptation criterion: %d", adaptor->adaptCriterion);
766     }
767   }
768   PetscFunctionReturn(PETSC_SUCCESS);
769 }
770 
771 static PetscErrorCode DMAdaptorPostAdapt(DMAdaptor adaptor)
772 {
773   PetscDS      prob;
774   PetscObject  obj;
775   PetscClassId id;
776 
777   PetscFunctionBegin;
778   PetscCall(DMGetDS(adaptor->idm, &prob));
779   PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
780   PetscCall(PetscObjectGetClassId(obj, &id));
781   if (id == PETSCFV_CLASSID) {
782     PetscFV fvm = (PetscFV)obj;
783 
784     PetscCall(PetscFVSetComputeGradients(fvm, adaptor->computeGradient));
785     /* Restore original limiter */
786     PetscCall(PetscFVSetLimiter(fvm, adaptor->limiter));
787 
788     PetscCall(VecRestoreArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
789     PetscCall(VecRestoreArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
790     PetscCall(DMRestoreLocalVector(adaptor->gradDM, &adaptor->cellGrad));
791   }
792   PetscFunctionReturn(PETSC_SUCCESS);
793 }
794 
795 /*
796   DMAdaptorComputeCellErrorIndicator_Gradient - Use the integrated gradient as an error indicator in the `DMAdaptor`
797 
798   Input Parameters:
799 + adaptor  - The `DMAdaptor` object
800 . dim      - The topological dimension
801 . cell     - The cell
802 . field    - The field integrated over the cell
803 . gradient - The gradient integrated over the cell
804 . cg       - A `PetscFVCellGeom` struct
805 - ctx      - A user context
806 
807   Output Parameter:
808 . errInd   - The error indicator
809 
810   Developer Note:
811   Some of the input arguments are absurdly specialized to special situations, it is not clear this is a good general API
812 
813 .seealso: [](ch_dmbase), `DMAdaptor`
814 */
815 static PetscErrorCode DMAdaptorComputeCellErrorIndicator_Gradient(DMAdaptor adaptor, PetscInt dim, PetscInt Nc, const PetscScalar *field, const PetscScalar *gradient, const PetscFVCellGeom *cg, PetscReal *errInd, void *ctx)
816 {
817   PetscReal err = 0.;
818   PetscInt  c, d;
819 
820   PetscFunctionBeginHot;
821   for (c = 0; c < Nc; c++) {
822     for (d = 0; d < dim; ++d) err += PetscSqr(PetscRealPart(gradient[c * dim + d]));
823   }
824   *errInd = cg->volume * err;
825   PetscFunctionReturn(PETSC_SUCCESS);
826 }
827 
828 static PetscErrorCode DMAdaptorComputeErrorIndicator_Gradient(DMAdaptor adaptor, Vec locX, Vec errVec)
829 {
830   DM              dm, plex, edm, eplex;
831   PetscDS         ds;
832   PetscObject     obj;
833   PetscClassId    id;
834   void           *ctx;
835   PetscQuadrature quad;
836   PetscScalar    *earray;
837   PetscReal       minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2];
838   PetscInt        dim, cdim, cStart, cEnd, Nf, Nc;
839 
840   PetscFunctionBegin;
841   PetscCall(VecGetDM(locX, &dm));
842   PetscCall(DMConvert(dm, DMPLEX, &plex));
843   PetscCall(VecGetDM(errVec, &edm));
844   PetscCall(DMConvert(edm, DMPLEX, &eplex));
845   PetscCall(DMGetDimension(plex, &dim));
846   PetscCall(DMGetCoordinateDim(plex, &cdim));
847   PetscCall(DMGetApplicationContext(plex, &ctx));
848   PetscCall(DMGetDS(plex, &ds));
849   PetscCall(PetscDSGetNumFields(ds, &Nf));
850   PetscCall(PetscDSGetDiscretization(ds, 0, &obj));
851   PetscCall(PetscObjectGetClassId(obj, &id));
852 
853   PetscCall(VecGetArray(errVec, &earray));
854   PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
855   for (PetscInt cell = cStart; cell < cEnd; ++cell) {
856     PetscScalar *eval;
857     PetscReal    errInd = 0.;
858 
859     if (id == PETSCFV_CLASSID) {
860       PetscFV            fv = (PetscFV)obj;
861       const PetscScalar *pointSols;
862       const PetscScalar *pointSol;
863       const PetscScalar *pointGrad;
864       PetscFVCellGeom   *cg;
865 
866       PetscCall(PetscFVGetNumComponents(fv, &Nc));
867       PetscCall(VecGetArrayRead(locX, &pointSols));
868       PetscCall(DMPlexPointLocalRead(plex, cell, pointSols, (void *)&pointSol));
869       PetscCall(DMPlexPointLocalRead(adaptor->gradDM, cell, adaptor->cellGradArray, (void *)&pointGrad));
870       PetscCall(DMPlexPointLocalRead(adaptor->cellDM, cell, adaptor->cellGeomArray, &cg));
871       PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, pointSol, pointGrad, cg, &errInd, ctx);
872       PetscCall(VecRestoreArrayRead(locX, &pointSols));
873     } else {
874       PetscFE          fe = (PetscFE)obj;
875       PetscScalar     *x  = NULL, *field, *gradient, *interpolant, *interpolantGrad;
876       PetscFVCellGeom  cg;
877       PetscFEGeom      fegeom;
878       const PetscReal *quadWeights;
879       PetscReal       *coords;
880       PetscInt         Nb, Nq, qNc;
881 
882       fegeom.dim      = dim;
883       fegeom.dimEmbed = cdim;
884       PetscCall(PetscFEGetNumComponents(fe, &Nc));
885       PetscCall(PetscFEGetQuadrature(fe, &quad));
886       PetscCall(PetscFEGetDimension(fe, &Nb));
887       PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
888       PetscCall(PetscCalloc6(Nc, &field, cdim * Nc, &gradient, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ));
889       PetscCall(PetscMalloc2(Nc, &interpolant, cdim * Nc, &interpolantGrad));
890       PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
891       PetscCall(DMPlexComputeCellGeometryFVM(plex, cell, &cg.volume, NULL, NULL));
892       PetscCall(PetscArrayzero(gradient, cdim * Nc));
893       PetscCall(DMPlexVecGetClosure(plex, NULL, locX, cell, NULL, &x));
894       for (PetscInt f = 0; f < Nf; ++f) {
895         PetscInt qc = 0;
896 
897         PetscCall(PetscDSGetDiscretization(ds, f, &obj));
898         PetscCall(PetscArrayzero(interpolant, Nc));
899         PetscCall(PetscArrayzero(interpolantGrad, cdim * Nc));
900         for (PetscInt q = 0; q < Nq; ++q) {
901           PetscCall(PetscFEInterpolateFieldAndGradient_Static((PetscFE)obj, 1, x, &fegeom, q, interpolant, interpolantGrad));
902           for (PetscInt fc = 0; fc < Nc; ++fc) {
903             const PetscReal wt = quadWeights[q * qNc + qc + fc];
904 
905             field[fc] += interpolant[fc] * wt * fegeom.detJ[q];
906             for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] += interpolantGrad[fc * dim + d] * wt * fegeom.detJ[q];
907           }
908         }
909         qc += Nc;
910       }
911       PetscCall(PetscFree2(interpolant, interpolantGrad));
912       PetscCall(DMPlexVecRestoreClosure(plex, NULL, locX, cell, NULL, &x));
913       for (PetscInt fc = 0; fc < Nc; ++fc) {
914         field[fc] /= cg.volume;
915         for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] /= cg.volume;
916       }
917       PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, field, gradient, &cg, &errInd, ctx);
918       PetscCall(PetscFree6(field, gradient, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
919     }
920     PetscCall(DMPlexPointGlobalRef(eplex, cell, earray, (void *)&eval));
921     eval[0]      = errInd;
922     minMaxInd[0] = PetscMin(minMaxInd[0], errInd);
923     minMaxInd[1] = PetscMax(minMaxInd[1], errInd);
924   }
925   PetscCall(VecRestoreArray(errVec, &earray));
926   PetscCall(DMDestroy(&plex));
927   PetscCall(DMDestroy(&eplex));
928   PetscCall(PetscGlobalMinMaxReal(PetscObjectComm((PetscObject)adaptor), minMaxInd, minMaxIndGlobal));
929   PetscCall(PetscInfo(adaptor, "DMAdaptor: error indicator range (%g, %g)\n", (double)minMaxIndGlobal[0], (double)minMaxIndGlobal[1]));
930   PetscFunctionReturn(PETSC_SUCCESS);
931 }
932 
933 static PetscErrorCode DMAdaptorComputeErrorIndicator_Flux(DMAdaptor adaptor, Vec lu, Vec errVec)
934 {
935   DM    dm, mdm;
936   SNES  msnes;
937   Vec   mu, lmu;
938   void *ctx;
939 
940   PetscFunctionBegin;
941   PetscCall(VecGetDM(lu, &dm));
942 
943   // Set up and solve mixed problem
944   PetscCall(DMClone(dm, &mdm));
945   PetscCall(SNESCreate(PetscObjectComm((PetscObject)mdm), &msnes));
946   PetscCall(SNESSetDM(msnes, mdm));
947   PetscCall(SNESAppendOptionsPrefix(msnes, "mixed_"));
948 
949   PetscTryTypeMethod(adaptor, mixedsetup, mdm);
950   PetscCall(DMGetApplicationContext(dm, &ctx));
951   PetscCall(DMPlexSetSNESLocalFEM(mdm, PETSC_FALSE, ctx));
952   PetscCall(SNESSetFromOptions(msnes));
953 
954   PetscCall(DMCreateGlobalVector(mdm, &mu));
955   PetscCall(PetscObjectSetName((PetscObject)mu, "Mixed Solution"));
956   PetscCall(VecSet(mu, 0.0));
957   PetscCall(SNESSolve(msnes, NULL, mu));
958   PetscCall(VecViewFromOptions(mu, (PetscObject)adaptor, "-adapt_mixed_sol_vec_view"));
959 
960   PetscCall(DMGetLocalVector(mdm, &lmu));
961   PetscCall(DMGlobalToLocal(mdm, mu, INSERT_VALUES, lmu));
962   PetscCall(DMPlexInsertBoundaryValues(mdm, PETSC_TRUE, lmu, 0.0, NULL, NULL, NULL));
963   PetscCall(DMPlexComputeL2FluxDiffVecLocal(lu, 0, lmu, 0, errVec));
964   PetscCall(DMRestoreLocalVector(mdm, &lmu));
965   PetscCall(VecDestroy(&mu));
966   PetscCall(SNESDestroy(&msnes));
967   PetscCall(DMDestroy(&mdm));
968   PetscFunctionReturn(PETSC_SUCCESS);
969 }
970 
971 /*@
972   DMAdaptorMonitor - runs the user provided monitor routines, if they exist
973 
974   Collective
975 
976   Input Parameters:
977 + adaptor - the `DMAdaptor`
978 . it      - iteration number
979 . odm     - the original `DM`
980 . adm     - the adapted `DM`
981 . Nf      - the number of fields
982 . enorms  - the 2-norm error values for each field
983 - error   - `Vec` of cellwise errors
984 
985   Level: developer
986 
987   Note:
988   This routine is called by the `DMAdaptor` implementations.
989   It does not typically need to be called by the user.
990 
991 .seealso: [](ch_snes), `DMAdaptorMonitorSet()`
992 @*/
993 PetscErrorCode DMAdaptorMonitor(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error)
994 {
995   PetscFunctionBegin;
996   for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) PetscCall((*adaptor->monitor[i])(adaptor, it, odm, adm, Nf, enorms, error, adaptor->monitorcontext[i]));
997   PetscFunctionReturn(PETSC_SUCCESS);
998 }
999 
1000 /*@C
1001   DMAdaptorMonitorSize - Prints the mesh sizes at each iteration of an adaptation loop.
1002 
1003   Collective
1004 
1005   Input Parameters:
1006 + adaptor - the `DMAdaptor`
1007 . n       - iteration number
1008 . odm     - the original `DM`
1009 . adm     - the adapted `DM`
1010 . Nf      - number of fields
1011 . enorms  - 2-norm error values for each field (may be estimated).
1012 . error   - `Vec` of cellwise errors
1013 - vf      - The viewer context
1014 
1015   Options Database Key:
1016 . -adaptor_monitor_size - Activates `DMAdaptorMonitorSize()`
1017 
1018   Level: intermediate
1019 
1020   Note:
1021   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1022   to be used during the adaptation loop.
1023 
1024 .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorError()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
1025 @*/
1026 PetscErrorCode DMAdaptorMonitorSize(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1027 {
1028   PetscViewer       viewer = vf->viewer;
1029   PetscViewerFormat format = vf->format;
1030   PetscInt          tablevel, cStart, cEnd, acStart, acEnd;
1031   const char       *prefix;
1032   PetscMPIInt       rank;
1033 
1034   PetscFunctionBegin;
1035   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
1036   PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
1037   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
1038   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)adaptor), &rank));
1039   PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
1040   PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
1041 
1042   PetscCall(PetscViewerPushFormat(viewer, format));
1043   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1044   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Sizes for %s adaptation.\n", prefix));
1045   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor rank %d N_orig: %" PetscInt_FMT " N_adapt: %" PetscInt_FMT "\n", n, rank, cEnd - cStart, acEnd - acStart));
1046   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1047   PetscCall(PetscViewerPopFormat(viewer));
1048   PetscFunctionReturn(PETSC_SUCCESS);
1049 }
1050 
1051 /*@C
1052   DMAdaptorMonitorError - Prints the error norm at each iteration of an adaptation loop.
1053 
1054   Collective
1055 
1056   Input Parameters:
1057 + adaptor - the `DMAdaptor`
1058 . n       - iteration number
1059 . odm     - the original `DM`
1060 . adm     - the adapted `DM`
1061 . Nf      - number of fields
1062 . enorms  - 2-norm error values for each field (may be estimated).
1063 . error   - `Vec` of cellwise errors
1064 - vf      - The viewer context
1065 
1066   Options Database Key:
1067 . -adaptor_monitor_error - Activates `DMAdaptorMonitorError()`
1068 
1069   Level: intermediate
1070 
1071   Note:
1072   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1073   to be used during the adaptation loop.
1074 
1075 .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
1076 @*/
1077 PetscErrorCode DMAdaptorMonitorError(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1078 {
1079   PetscViewer       viewer = vf->viewer;
1080   PetscViewerFormat format = vf->format;
1081   PetscInt          tablevel, cStart, cEnd, acStart, acEnd;
1082   const char       *prefix;
1083 
1084   PetscFunctionBegin;
1085   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
1086   PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
1087   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
1088 
1089   PetscCall(PetscViewerPushFormat(viewer, format));
1090   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1091   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Error norms for %s adaptation.\n", prefix));
1092   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor Error norm %s", n, Nf > 1 ? "[" : ""));
1093   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1094   for (PetscInt f = 0; f < Nf; ++f) {
1095     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1096     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)enorms[f]));
1097   }
1098   PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
1099   PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
1100   PetscCall(PetscViewerASCIIPrintf(viewer, " N: %" PetscInt_FMT " Nadapt: %" PetscInt_FMT "\n", cEnd - cStart, acEnd - acStart));
1101   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1102   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1103   PetscCall(PetscViewerPopFormat(viewer));
1104   PetscFunctionReturn(PETSC_SUCCESS);
1105 }
1106 
1107 /*@C
1108   DMAdaptorMonitorErrorDraw - Plots the error at each iteration of an iterative solver.
1109 
1110   Collective
1111 
1112   Input Parameters:
1113 + adaptor - the `DMAdaptor`
1114 . n       - iteration number
1115 . odm     - the original `DM`
1116 . adm     - the adapted `DM`
1117 . Nf      - number of fields
1118 . enorms  - 2-norm error values for each field (may be estimated).
1119 . error   - `Vec` of cellwise errors
1120 - vf      - The viewer context
1121 
1122   Options Database Key:
1123 . -adaptor_monitor_error draw - Activates `DMAdaptorMonitorErrorDraw()`
1124 
1125   Level: intermediate
1126 
1127   Note:
1128   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1129   to be used during the adaptation loop.
1130 
1131 .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
1132 @*/
1133 PetscErrorCode DMAdaptorMonitorErrorDraw(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1134 {
1135   PetscViewer       viewer = vf->viewer;
1136   PetscViewerFormat format = vf->format;
1137 
1138   PetscFunctionBegin;
1139   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
1140   PetscCall(PetscViewerPushFormat(viewer, format));
1141   PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
1142   PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", (PetscObject)adaptor));
1143   PetscCall(VecView(error, viewer));
1144   PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", NULL));
1145   PetscCall(PetscViewerPopFormat(viewer));
1146   PetscFunctionReturn(PETSC_SUCCESS);
1147 }
1148 
1149 /*@C
1150   DMAdaptorMonitorErrorDrawLGCreate - Creates the context for the erro plotter `DMAdaptorMonitorErrorDrawLG()`
1151 
1152   Collective
1153 
1154   Input Parameters:
1155 + viewer - The `PetscViewer`
1156 . format - The viewer format
1157 - ctx    - An optional user context
1158 
1159   Output Parameter:
1160 . vf - The viewer context
1161 
1162   Level: intermediate
1163 
1164 .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
1165 @*/
1166 PetscErrorCode DMAdaptorMonitorErrorDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
1167 {
1168   DMAdaptor adaptor = (DMAdaptor)ctx;
1169   char    **names;
1170   PetscInt  Nf;
1171 
1172   PetscFunctionBegin;
1173   PetscCall(DMGetNumFields(adaptor->idm, &Nf));
1174   PetscCall(PetscMalloc1(Nf + 1, &names));
1175   for (PetscInt f = 0; f < Nf; ++f) {
1176     PetscObject disc;
1177     const char *fname;
1178     char        lname[PETSC_MAX_PATH_LEN];
1179 
1180     PetscCall(DMGetField(adaptor->idm, f, NULL, &disc));
1181     PetscCall(PetscObjectGetName(disc, &fname));
1182     PetscCall(PetscStrncpy(lname, fname, PETSC_MAX_PATH_LEN));
1183     PetscCall(PetscStrlcat(lname, " Error", PETSC_MAX_PATH_LEN));
1184     PetscCall(PetscStrallocpy(lname, &names[f]));
1185   }
1186   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
1187   (*vf)->data = ctx;
1188   PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Error Norm", Nf, (const char **)names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
1189   for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(names[f]));
1190   PetscCall(PetscFree(names));
1191   PetscFunctionReturn(PETSC_SUCCESS);
1192 }
1193 
1194 /*@C
1195   DMAdaptorMonitorErrorDrawLG - Plots the error norm at each iteration of an adaptive loop.
1196 
1197   Collective
1198 
1199   Input Parameters:
1200 + adaptor - the `DMAdaptor`
1201 . n       - iteration number
1202 . odm     - the original `DM`
1203 . adm     - the adapted `DM`
1204 . Nf      - number of fields
1205 . enorms  - 2-norm error values for each field (may be estimated).
1206 . error   - `Vec` of cellwise errors
1207 - vf      - The viewer context
1208 
1209   Options Database Key:
1210 . -adaptor_error draw::draw_lg - Activates `DMAdaptorMonitorErrorDrawLG()`
1211 
1212   Level: intermediate
1213 
1214   Notes:
1215   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1216   to be used during the adaptation loop.
1217 
1218   Call `DMAdaptorMonitorErrorDrawLGCreate()` to create the context needed for this monitor
1219 
1220 .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdadptorMonitorErrorDraw()`, `DMAdadptorMonitorError()`,
1221           `DMAdaptorMonitorTrueResidualDrawLGCreate()`
1222 @*/
1223 PetscErrorCode DMAdaptorMonitorErrorDrawLG(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1224 {
1225   PetscViewer       viewer = vf->viewer;
1226   PetscViewerFormat format = vf->format;
1227   PetscDrawLG       lg     = vf->lg;
1228   PetscReal        *x, *e;
1229 
1230   PetscFunctionBegin;
1231   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
1232   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 8);
1233   PetscCall(PetscCalloc2(Nf, &x, Nf, &e));
1234   PetscCall(PetscViewerPushFormat(viewer, format));
1235   if (!n) PetscCall(PetscDrawLGReset(lg));
1236   for (PetscInt f = 0; f < Nf; ++f) {
1237     x[f] = (PetscReal)n;
1238     e[f] = enorms[f] > 0.0 ? PetscLog10Real(enorms[f]) : -15.;
1239   }
1240   PetscCall(PetscDrawLGAddPoint(lg, x, e));
1241   PetscCall(PetscDrawLGDraw(lg));
1242   PetscCall(PetscDrawLGSave(lg));
1243   PetscCall(PetscViewerPopFormat(viewer));
1244   PetscCall(PetscFree2(x, e));
1245   PetscFunctionReturn(PETSC_SUCCESS);
1246 }
1247 
1248 /*@C
1249   DMAdaptorMonitorRegisterAll - Registers all of the mesh adaptation monitors in the `SNES` package.
1250 
1251   Not Collective
1252 
1253   Level: advanced
1254 
1255 .seealso: [](ch_snes), `SNES`, `DM`, `DMAdaptorMonitorRegister()`, `DMAdaptorRegister()`
1256 @*/
1257 PetscErrorCode DMAdaptorMonitorRegisterAll(void)
1258 {
1259   PetscFunctionBegin;
1260   if (DMAdaptorMonitorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1261   DMAdaptorMonitorRegisterAllCalled = PETSC_TRUE;
1262 
1263   PetscCall(DMAdaptorMonitorRegister("size", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorSize, NULL, NULL));
1264   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorError, NULL, NULL));
1265   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorErrorDraw, NULL, NULL));
1266   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DRAW_LG, DMAdaptorMonitorErrorDrawLG, DMAdaptorMonitorErrorDrawLGCreate, NULL));
1267   PetscFunctionReturn(PETSC_SUCCESS);
1268 }
1269 
1270 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
1271 {
1272   const PetscInt Nc = uOff[1] - uOff[0];
1273 
1274   for (PetscInt i = 0; i < Nc; ++i) f[i] = u[i];
1275 }
1276 
1277 static void identityFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
1278 {
1279   for (PetscInt i = 0; i < dim; ++i) {
1280     for (PetscInt j = 0; j < dim; ++j) f[i + dim * j] = u[i + dim * j];
1281   }
1282 }
1283 
1284 static PetscErrorCode DMAdaptorAdapt_Sequence_Private(DMAdaptor adaptor, Vec inx, PetscBool doSolve, DM *adm, Vec *ax)
1285 {
1286   PetscDS   ds;
1287   PetscReal errorNorm = 0.;
1288   PetscInt  numAdapt  = adaptor->numSeq, adaptIter;
1289   PetscInt  dim, coordDim, Nf;
1290   void     *ctx;
1291   MPI_Comm  comm;
1292 
1293   PetscFunctionBegin;
1294   PetscCall(DMViewFromOptions(adaptor->idm, NULL, "-dm_adapt_pre_view"));
1295   PetscCall(VecViewFromOptions(inx, NULL, "-sol_adapt_pre_view"));
1296   PetscCall(PetscObjectGetComm((PetscObject)adaptor, &comm));
1297   PetscCall(DMGetDimension(adaptor->idm, &dim));
1298   PetscCall(DMGetCoordinateDim(adaptor->idm, &coordDim));
1299   PetscCall(DMGetApplicationContext(adaptor->idm, &ctx));
1300   PetscCall(DMGetDS(adaptor->idm, &ds));
1301   PetscCall(PetscDSGetNumFields(ds, &Nf));
1302   PetscCheck(Nf != 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot refine with no fields present!");
1303 
1304   /* Adapt until nothing changes */
1305   /* Adapt for a specified number of iterates */
1306   for (adaptIter = 0; adaptIter < numAdapt - 1; ++adaptIter) PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(comm)));
1307   for (adaptIter = 0; adaptIter < numAdapt; ++adaptIter) {
1308     PetscBool adapted = PETSC_FALSE;
1309     DM        dm      = adaptIter ? *adm : adaptor->idm, odm;
1310     Vec       x       = adaptIter ? *ax : inx, locX, ox;
1311     Vec       error   = NULL;
1312 
1313     PetscCall(DMGetLocalVector(dm, &locX));
1314     PetscCall(DMAdaptorPreAdapt(adaptor, locX));
1315     if (doSolve) {
1316       SNES snes;
1317 
1318       PetscCall(DMAdaptorGetSolver(adaptor, &snes));
1319       PetscCall(SNESSolve(snes, NULL, adaptIter ? *ax : x));
1320     }
1321     PetscCall(DMGlobalToLocalBegin(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
1322     PetscCall(DMGlobalToLocalEnd(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
1323     PetscCall(VecViewFromOptions(adaptIter ? *ax : x, (PetscObject)adaptor, "-adapt_primal_sol_vec_view"));
1324     switch (adaptor->adaptCriterion) {
1325     case DM_ADAPTATION_REFINE:
1326       PetscCall(DMRefine(dm, comm, &odm));
1327       PetscCheck(odm, comm, PETSC_ERR_ARG_INCOMP, "DMRefine() did not perform any refinement, cannot continue grid sequencing");
1328       adapted = PETSC_TRUE;
1329       PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, NULL));
1330       break;
1331     case DM_ADAPTATION_LABEL: {
1332       /* Adapt DM
1333            Create local solution
1334            Reconstruct gradients (FVM) or solve adjoint equation (FEM)
1335            Produce cellwise error indicator */
1336       DM             edm, plex;
1337       PetscDS        ds;
1338       PetscFE        efe;
1339       DMLabel        adaptLabel;
1340       IS             refineIS, coarsenIS;
1341       DMPolytopeType ct;
1342       PetscScalar    errorVal;
1343       PetscInt       nRefine, nCoarsen, cStart;
1344 
1345       PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
1346 
1347       // TODO Move this creation to PreAdapt
1348       PetscCall(DMClone(dm, &edm));
1349       PetscCall(DMConvert(edm, DMPLEX, &plex));
1350       PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, NULL));
1351       PetscCall(DMPlexGetCellType(plex, cStart, &ct));
1352       PetscCall(DMDestroy(&plex));
1353       PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 0, PETSC_DEFAULT, &efe));
1354       PetscCall(PetscObjectSetName((PetscObject)efe, "Error"));
1355       PetscCall(DMSetField(edm, 0, NULL, (PetscObject)efe));
1356       PetscCall(PetscFEDestroy(&efe));
1357       PetscCall(DMCreateDS(edm));
1358       PetscCall(DMGetGlobalVector(edm, &error));
1359       PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
1360 
1361       PetscUseTypeMethod(adaptor, computeerrorindicator, locX, error);
1362       PetscCall(VecViewFromOptions(error, (PetscObject)adaptor, "-adapt_error_vec_view"));
1363       PetscCall(DMGetDS(edm, &ds));
1364       PetscCall(PetscDSSetObjective(ds, 0, identity));
1365       PetscCall(DMPlexComputeIntegralFEM(edm, error, &errorVal, NULL));
1366       errorNorm = PetscRealPart(errorVal);
1367 
1368       // Compute IS from VecTagger
1369       PetscCall(VecTaggerComputeIS(adaptor->refineTag, error, &refineIS, NULL));
1370       PetscCall(VecTaggerComputeIS(adaptor->coarsenTag, error, &coarsenIS, NULL));
1371       PetscCall(ISViewFromOptions(refineIS, (PetscObject)adaptor->refineTag, "-is_view"));
1372       PetscCall(ISViewFromOptions(coarsenIS, (PetscObject)adaptor->coarsenTag, "-is_view"));
1373       PetscCall(ISGetSize(refineIS, &nRefine));
1374       PetscCall(ISGetSize(coarsenIS, &nCoarsen));
1375       PetscCall(PetscInfo(adaptor, "DMAdaptor: numRefine %" PetscInt_FMT ", numCoarsen %" PetscInt_FMT "\n", nRefine, nCoarsen));
1376       if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
1377       if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
1378       PetscCall(ISDestroy(&coarsenIS));
1379       PetscCall(ISDestroy(&refineIS));
1380       // Adapt DM from label
1381       if (nRefine || nCoarsen) {
1382         char        oprefix[PETSC_MAX_PATH_LEN];
1383         const char *p;
1384         PetscBool   flg;
1385 
1386         PetscCall(PetscOptionsHasName(NULL, adaptor->hdr.prefix, "-adapt_vec_view", &flg));
1387         if (flg) {
1388           Vec ref;
1389 
1390           PetscCall(DMPlexCreateLabelField(dm, adaptLabel, &ref));
1391           PetscCall(VecViewFromOptions(ref, (PetscObject)adaptor, "-adapt_vec_view"));
1392           PetscCall(VecDestroy(&ref));
1393         }
1394 
1395         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &p));
1396         PetscCall(PetscStrncpy(oprefix, p, PETSC_MAX_PATH_LEN));
1397         PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)dm, "adapt_"));
1398         PetscCall(DMAdaptLabel(dm, adaptLabel, &odm));
1399         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, oprefix));
1400         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)odm, oprefix));
1401         PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, error));
1402         adapted = PETSC_TRUE;
1403       } else {
1404         PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, error));
1405       }
1406       PetscCall(DMLabelDestroy(&adaptLabel));
1407       PetscCall(DMRestoreGlobalVector(edm, &error));
1408       PetscCall(DMDestroy(&edm));
1409     } break;
1410     case DM_ADAPTATION_METRIC: {
1411       DM        dmGrad, dmHess, dmMetric, dmDet;
1412       Vec       xGrad, xHess, metric, determinant;
1413       PetscReal N;
1414       DMLabel   bdLabel = NULL, rgLabel = NULL;
1415       PetscBool higherOrder = PETSC_FALSE;
1416       PetscInt  Nd          = coordDim * coordDim, f, vStart, vEnd;
1417       void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
1418 
1419       PetscCall(PetscMalloc(1, &funcs));
1420       funcs[0] = identityFunc;
1421 
1422       /*     Setup finite element spaces */
1423       PetscCall(DMClone(dm, &dmGrad));
1424       PetscCall(DMClone(dm, &dmHess));
1425       PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple fields not yet considered"); // TODO
1426       for (f = 0; f < Nf; ++f) {
1427         PetscFE         fe, feGrad, feHess;
1428         PetscDualSpace  Q;
1429         PetscSpace      space;
1430         DM              K;
1431         PetscQuadrature q;
1432         PetscInt        Nc, qorder, p;
1433         const char     *prefix;
1434 
1435         PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1436         PetscCall(PetscFEGetNumComponents(fe, &Nc));
1437         PetscCheck(Nc <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple components not yet considered"); // TODO
1438         PetscCall(PetscFEGetBasisSpace(fe, &space));
1439         PetscCall(PetscSpaceGetDegree(space, NULL, &p));
1440         if (p > 1) higherOrder = PETSC_TRUE;
1441         PetscCall(PetscFEGetDualSpace(fe, &Q));
1442         PetscCall(PetscDualSpaceGetDM(Q, &K));
1443         PetscCall(DMPlexGetDepthStratum(K, 0, &vStart, &vEnd));
1444         PetscCall(PetscFEGetQuadrature(fe, &q));
1445         PetscCall(PetscQuadratureGetOrder(q, &qorder));
1446         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix));
1447         PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmGrad), dim, Nc * coordDim, PETSC_TRUE, prefix, qorder, &feGrad));
1448         PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmHess), dim, Nc * Nd, PETSC_TRUE, prefix, qorder, &feHess));
1449         PetscCall(DMSetField(dmGrad, f, NULL, (PetscObject)feGrad));
1450         PetscCall(DMSetField(dmHess, f, NULL, (PetscObject)feHess));
1451         PetscCall(DMCreateDS(dmGrad));
1452         PetscCall(DMCreateDS(dmHess));
1453         PetscCall(PetscFEDestroy(&feGrad));
1454         PetscCall(PetscFEDestroy(&feHess));
1455       }
1456       /*     Compute vertexwise gradients from cellwise gradients */
1457       PetscCall(DMCreateLocalVector(dmGrad, &xGrad));
1458       PetscCall(VecViewFromOptions(locX, NULL, "-sol_adapt_loc_pre_view"));
1459       PetscCall(DMPlexComputeGradientClementInterpolant(dm, locX, xGrad));
1460       PetscCall(VecViewFromOptions(xGrad, NULL, "-adapt_gradient_view"));
1461       /*     Compute vertexwise Hessians from cellwise Hessians */
1462       PetscCall(DMCreateLocalVector(dmHess, &xHess));
1463       PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, xGrad, xHess));
1464       PetscCall(VecViewFromOptions(xHess, NULL, "-adapt_hessian_view"));
1465       PetscCall(VecDestroy(&xGrad));
1466       PetscCall(DMDestroy(&dmGrad));
1467       /*     Compute L-p normalized metric */
1468       PetscCall(DMClone(dm, &dmMetric));
1469       N = adaptor->Nadapt >= 0 ? adaptor->Nadapt : PetscPowRealInt(adaptor->refinementFactor, dim) * ((PetscReal)(vEnd - vStart));
1470       // TODO This was where the old monitor was, figure out how to show metric and target N
1471       PetscCall(DMPlexMetricSetTargetComplexity(dmMetric, (PetscReal)N));
1472       if (higherOrder) {
1473         /*   Project Hessian into P1 space, if required */
1474         PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1475         PetscCall(DMProjectFieldLocal(dmMetric, 0.0, xHess, funcs, INSERT_ALL_VALUES, metric));
1476         PetscCall(VecDestroy(&xHess));
1477         xHess = metric;
1478       }
1479       PetscCall(PetscFree(funcs));
1480       PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1481       PetscCall(DMPlexMetricDeterminantCreate(dmMetric, 0, &determinant, &dmDet));
1482       PetscCall(DMPlexMetricNormalize(dmMetric, xHess, PETSC_TRUE, PETSC_TRUE, metric, determinant));
1483       PetscCall(VecDestroy(&determinant));
1484       PetscCall(DMDestroy(&dmDet));
1485       PetscCall(VecDestroy(&xHess));
1486       PetscCall(DMDestroy(&dmHess));
1487       /*     Adapt DM from metric */
1488       PetscCall(DMGetLabel(dm, "marker", &bdLabel));
1489       PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &odm));
1490       PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, NULL));
1491       adapted = PETSC_TRUE;
1492       /*     Cleanup */
1493       PetscCall(VecDestroy(&metric));
1494       PetscCall(DMDestroy(&dmMetric));
1495     } break;
1496     default:
1497       SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid adaptation type: %d", adaptor->adaptCriterion);
1498     }
1499     PetscCall(DMAdaptorPostAdapt(adaptor));
1500     PetscCall(DMRestoreLocalVector(dm, &locX));
1501     /* If DM was adapted, replace objects and recreate solution */
1502     if (adapted) {
1503       const char *name;
1504 
1505       PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1506       PetscCall(PetscObjectSetName((PetscObject)odm, name));
1507       /* Reconfigure solver */
1508       PetscCall(SNESReset(adaptor->snes));
1509       PetscCall(SNESSetDM(adaptor->snes, odm));
1510       PetscCall(DMAdaptorSetSolver(adaptor, adaptor->snes));
1511       PetscCall(DMPlexSetSNESLocalFEM(odm, PETSC_FALSE, ctx));
1512       PetscCall(SNESSetFromOptions(adaptor->snes));
1513       /* Transfer system */
1514       PetscCall(DMCopyDisc(dm, odm));
1515       /* Transfer solution */
1516       PetscCall(DMCreateGlobalVector(odm, &ox));
1517       PetscCall(PetscObjectGetName((PetscObject)x, &name));
1518       PetscCall(PetscObjectSetName((PetscObject)ox, name));
1519       PetscCall(DMAdaptorTransferSolution(adaptor, dm, x, odm, ox));
1520       /* Cleanup adaptivity info */
1521       if (adaptIter > 0) PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(comm)));
1522       PetscCall(DMForestSetAdaptivityForest(dm, NULL)); /* clear internal references to the previous dm */
1523       PetscCall(DMDestroy(&dm));
1524       PetscCall(VecDestroy(&x));
1525       *adm = odm;
1526       *ax  = ox;
1527     } else {
1528       *adm      = dm;
1529       *ax       = x;
1530       adaptIter = numAdapt;
1531     }
1532     if (adaptIter < numAdapt - 1) {
1533       PetscCall(DMViewFromOptions(odm, NULL, "-dm_adapt_iter_view"));
1534       PetscCall(VecViewFromOptions(ox, NULL, "-sol_adapt_iter_view"));
1535     }
1536   }
1537   PetscCall(DMViewFromOptions(*adm, NULL, "-dm_adapt_view"));
1538   PetscCall(VecViewFromOptions(*ax, NULL, "-sol_adapt_view"));
1539   PetscFunctionReturn(PETSC_SUCCESS);
1540 }
1541 
1542 /*@
1543   DMAdaptorAdapt - Creates a new `DM` that is adapted to the problem
1544 
1545   Not Collective
1546 
1547   Input Parameters:
1548 + adaptor  - The `DMAdaptor` object
1549 . x        - The global approximate solution
1550 - strategy - The adaptation strategy, see `DMAdaptationStrategy`
1551 
1552   Output Parameters:
1553 + adm - The adapted `DM`
1554 - ax  - The adapted solution
1555 
1556   Options Database Keys:
1557 + -snes_adapt <strategy> - initial, sequential, multigrid
1558 . -adapt_gradient_view   - View the Clement interpolant of the solution gradient
1559 . -adapt_hessian_view    - View the Clement interpolant of the solution Hessian
1560 - -adapt_metric_view     - View the metric tensor for adaptive mesh refinement
1561 
1562   Level: intermediate
1563 
1564 .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`
1565 @*/
1566 PetscErrorCode DMAdaptorAdapt(DMAdaptor adaptor, Vec x, DMAdaptationStrategy strategy, DM *adm, Vec *ax)
1567 {
1568   PetscFunctionBegin;
1569   switch (strategy) {
1570   case DM_ADAPTATION_INITIAL:
1571     PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_FALSE, adm, ax));
1572     break;
1573   case DM_ADAPTATION_SEQUENTIAL:
1574     PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_TRUE, adm, ax));
1575     break;
1576   default:
1577     SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Unrecognized adaptation strategy %d", strategy);
1578   }
1579   PetscFunctionReturn(PETSC_SUCCESS);
1580 }
1581 
1582 /*@C
1583   DMAdaptorGetMixedSetupFunction - Get the function setting up the mixed problem, if it exists
1584 
1585   Not Collective
1586 
1587   Input Parameter:
1588 . adaptor - the `DMAdaptor`
1589 
1590   Output Parameter:
1591 . setupFunc - the function setting up the mixed problem, or `NULL`
1592 
1593   Level: advanced
1594 
1595 .seealso: `DMAdaptor`, `DMAdaptorSetMixedSetupFunction()`, `DMAdaptorAdapt()`
1596 @*/
1597 PetscErrorCode DMAdaptorGetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (**setupFunc)(DMAdaptor, DM))
1598 {
1599   PetscFunctionBegin;
1600   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
1601   PetscAssertPointer(setupFunc, 2);
1602   *setupFunc = adaptor->ops->mixedsetup;
1603   PetscFunctionReturn(PETSC_SUCCESS);
1604 }
1605 
1606 /*@C
1607   DMAdaptorSetMixedSetupFunction - Set the function setting up the mixed problem
1608 
1609   Not Collective
1610 
1611   Input Parameters:
1612 + adaptor   - the `DMAdaptor`
1613 - setupFunc - the function setting up the mixed problem
1614 
1615   Level: advanced
1616 
1617 .seealso: `DMAdaptor`, `DMAdaptorGetMixedSetupFunction()`, `DMAdaptorAdapt()`
1618 @*/
1619 PetscErrorCode DMAdaptorSetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (*setupFunc)(DMAdaptor, DM))
1620 {
1621   PetscFunctionBegin;
1622   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
1623   PetscValidFunction(setupFunc, 2);
1624   adaptor->ops->mixedsetup = setupFunc;
1625   PetscFunctionReturn(PETSC_SUCCESS);
1626 }
1627 
1628 /*@C
1629   DMAdaptorGetCriterion - Get the adaptation criterion
1630 
1631   Not Collective
1632 
1633   Input Parameter:
1634 . adaptor - the `DMAdaptor`
1635 
1636   Output Parameter:
1637 . criterion - the criterion for adaptation
1638 
1639   Level: advanced
1640 
1641 .seealso: `DMAdaptor`, `DMAdaptorSetCrierion()`, `DMAdaptationCriterion`
1642 @*/
1643 PetscErrorCode DMAdaptorGetCriterion(DMAdaptor adaptor, DMAdaptationCriterion *criterion)
1644 {
1645   PetscFunctionBegin;
1646   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
1647   PetscAssertPointer(criterion, 2);
1648   *criterion = adaptor->adaptCriterion;
1649   PetscFunctionReturn(PETSC_SUCCESS);
1650 }
1651 
1652 /*@C
1653   DMAdaptorSetCriterion - Set the adaptation criterion
1654 
1655   Not Collective
1656 
1657   Input Parameters:
1658 + adaptor   - the `DMAdaptor`
1659 - criterion - the adaptation criterion
1660 
1661   Level: advanced
1662 
1663 .seealso: `DMAdaptor`, `DMAdaptorGetCriterion()`, `DMAdaptationCriterion`
1664 @*/
1665 PetscErrorCode DMAdaptorSetCriterion(DMAdaptor adaptor, DMAdaptationCriterion criterion)
1666 {
1667   PetscFunctionBegin;
1668   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
1669   adaptor->adaptCriterion = criterion;
1670   PetscFunctionReturn(PETSC_SUCCESS);
1671 }
1672 
1673 static PetscErrorCode DMAdaptorInitialize_Gradient(DMAdaptor adaptor)
1674 {
1675   PetscFunctionBegin;
1676   adaptor->ops->computeerrorindicator     = DMAdaptorComputeErrorIndicator_Gradient;
1677   adaptor->ops->computecellerrorindicator = DMAdaptorComputeCellErrorIndicator_Gradient;
1678   PetscFunctionReturn(PETSC_SUCCESS);
1679 }
1680 
1681 PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Gradient(DMAdaptor adaptor)
1682 {
1683   PetscFunctionBegin;
1684   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
1685   adaptor->data = NULL;
1686 
1687   PetscCall(DMAdaptorInitialize_Gradient(adaptor));
1688   PetscFunctionReturn(PETSC_SUCCESS);
1689 }
1690 
1691 static PetscErrorCode DMAdaptorInitialize_Flux(DMAdaptor adaptor)
1692 {
1693   PetscFunctionBegin;
1694   adaptor->ops->computeerrorindicator = DMAdaptorComputeErrorIndicator_Flux;
1695   PetscFunctionReturn(PETSC_SUCCESS);
1696 }
1697 
1698 PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Flux(DMAdaptor adaptor)
1699 {
1700   PetscFunctionBegin;
1701   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
1702   adaptor->data = NULL;
1703 
1704   PetscCall(DMAdaptorInitialize_Flux(adaptor));
1705   PetscFunctionReturn(PETSC_SUCCESS);
1706 }
1707