Actual source code: dmadapt.c

  1: #include <petscdmadaptor.h>
  2: #include <petscdmplex.h>
  3: #include <petscdmforest.h>
  4: #include <petscds.h>
  5: #include <petscblaslapack.h>
  6: #include <petscsnes.h>
  7: #include <petscdraw.h>

  9: #include <petsc/private/dmadaptorimpl.h>
 10: #include <petsc/private/dmpleximpl.h>
 11: #include <petsc/private/petscfeimpl.h>

 13: PetscClassId DMADAPTOR_CLASSID;

 15: PetscFunctionList DMAdaptorList              = NULL;
 16: PetscBool         DMAdaptorRegisterAllCalled = PETSC_FALSE;

 18: PetscFunctionList DMAdaptorMonitorList              = NULL;
 19: PetscFunctionList DMAdaptorMonitorCreateList        = NULL;
 20: PetscFunctionList DMAdaptorMonitorDestroyList       = NULL;
 21: PetscBool         DMAdaptorMonitorRegisterAllCalled = PETSC_FALSE;

 23: const char *const DMAdaptationCriteria[] = {"NONE", "REFINE", "LABEL", "METRIC", "DMAdaptationCriterion", "DM_ADAPTATION_", NULL};

 25: /*@C
 26:   DMAdaptorRegister - Adds a new adaptor component implementation

 28:   Not Collective

 30:   Input Parameters:
 31: + name        - The name of a new user-defined creation routine
 32: - create_func - The creation routine

 34:   Example Usage:
 35: .vb
 36:   DMAdaptorRegister("my_adaptor", MyAdaptorCreate);
 37: .ve

 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

 49:   Level: advanced

 51:   Note:
 52:   `DMAdaptorRegister()` may be called multiple times to add several user-defined adaptors

 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: }

 64: PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Gradient(DMAdaptor);
 65: PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Flux(DMAdaptor);

 67: /*@C
 68:   DMAdaptorRegisterAll - Registers all of the adaptor components in the `DM` package.

 70:   Not Collective

 72:   Level: advanced

 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;

 82:   PetscCall(DMAdaptorRegister(DMADAPTORGRADIENT, DMAdaptorCreate_Gradient));
 83:   PetscCall(DMAdaptorRegister(DMADAPTORFLUX, DMAdaptorCreate_Flux));
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }

 87: /*@C
 88:   DMAdaptorRegisterDestroy - This function destroys the registered `DMAdaptorType`. It is called from `PetscFinalize()`.

 90:   Not collective

 92:   Level: developer

 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: }

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: }

115: /*@C
116:   DMAdaptorMonitorRegister -  Registers a mesh adaptation monitor routine that may be accessed with `DMAdaptorMonitorSetFromOptions()`

118:   Not Collective

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`

128:   Level: advanced

130:   Note:
131:   `DMAdaptorMonitorRegister()` may be called multiple times to add several user-defined monitors.

133:   Example Usage:
134: .vb
135:   DMAdaptorMonitorRegister("my_monitor", PETSCVIEWERASCII, PETSC_VIEWER_ASCII_INFO_DETAIL, MyMonitor, NULL, NULL);
136: .ve

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`

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];

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: }

159: /*@C
160:   DMAdaptorMonitorRegisterDestroy - This function destroys the registered monitors for `DMAdaptor`. It is called from `PetscFinalize()`.

162:   Not collective

164:   Level: developer

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: }

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`.

181:   Collective

183:   Input Parameter:
184: . comm - The communicator for the `DMAdaptor` object

186:   Output Parameter:
187: . adaptor - The `DMAdaptor` object

189:   Level: beginner

191: .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorDestroy()`, `DMAdaptorAdapt()`, `PetscConvEst`, `PetscConvEstCreate()`
192: @*/
193: PetscErrorCode DMAdaptorCreate(MPI_Comm comm, DMAdaptor *adaptor)
194: {
195:   VecTaggerBox refineBox, coarsenBox;

197:   PetscFunctionBegin;
198:   PetscAssertPointer(adaptor, 2);
199:   PetscCall(PetscSysInitializePackage());

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: }

219: /*@
220:   DMAdaptorDestroy - Destroys a `DMAdaptor` object

222:   Collective

224:   Input Parameter:
225: . adaptor - The `DMAdaptor` object

227:   Level: beginner

229: .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
230: @*/
231: PetscErrorCode DMAdaptorDestroy(DMAdaptor *adaptor)
232: {
233:   PetscFunctionBegin;
234:   if (!*adaptor) PetscFunctionReturn(PETSC_SUCCESS);
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: }

248: /*@
249:   DMAdaptorSetType - Sets the particular implementation for a adaptor.

251:   Collective

253:   Input Parameters:
254: + adaptor - The `DMAdaptor`
255: - method  - The name of the adaptor type

257:   Options Database Key:
258: . -adaptor_type <type> - Sets the adaptor type; see `DMAdaptorType`

260:   Level: intermediate

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;

269:   PetscFunctionBegin;
271:   PetscCall(PetscObjectTypeCompare((PetscObject)adaptor, method, &match));
272:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

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);

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: }

285: /*@
286:   DMAdaptorGetType - Gets the type name (as a string) from the adaptor.

288:   Not Collective

290:   Input Parameter:
291: . adaptor - The `DMAdaptor`

293:   Output Parameter:
294: . type - The `DMAdaptorType` name

296:   Level: intermediate

298: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptor`, `DMAdaptorType`, `DMAdaptorSetType()`, `DMAdaptorCreate()`
299: @*/
300: PetscErrorCode DMAdaptorGetType(DMAdaptor adaptor, DMAdaptorType *type)
301: {
302:   PetscFunctionBegin;
304:   PetscAssertPointer(type, 2);
305:   PetscCall(DMAdaptorRegisterAll());
306:   *type = ((PetscObject)adaptor)->type_name;
307:   PetscFunctionReturn(PETSC_SUCCESS);
308: }

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: }

318: /*@C
319:   DMAdaptorMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
320:   the error etc.

322:   Logically Collective

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`), see `PetscCtxDestroyFn` for its calling sequence

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()`

340:   Options Database Keys:
341: + -adaptor_monitor_size                - sets `DMAdaptorMonitorSize()`
342: . -adaptor_monitor_error               - sets `DMAdaptorMonitorError()`
343: . -adaptor_monitor_error draw          - sets `DMAdaptorMonitorErrorDraw()` and plots error
344: . -adaptor_monitor_error draw::draw_lg - sets `DMAdaptorMonitorErrorDrawLG()` and plots error
345: - -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.

347:   Level: beginner

349: .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptor`, `PetscCtxDestroyFn`
350: @*/
351: PetscErrorCode DMAdaptorMonitorSet(DMAdaptor adaptor, PetscErrorCode (*monitor)(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, void *ctx), void *ctx, PetscCtxDestroyFn *monitordestroy)
352: {
353:   PetscBool identical;

355:   PetscFunctionBegin;
357:   for (PetscInt i = 0; i < adaptor->numbermonitors; i++) {
358:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))monitor, ctx, monitordestroy, (PetscErrorCode (*)(void))adaptor->monitor[i], adaptor->monitorcontext[i], adaptor->monitordestroy[i], &identical));
359:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
360:   }
361:   PetscCheck(adaptor->numbermonitors < MAXDMADAPTORMONITORS, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_OUTOFRANGE, "Too many DMAdaptor monitors set");
362:   adaptor->monitor[adaptor->numbermonitors]          = monitor;
363:   adaptor->monitordestroy[adaptor->numbermonitors]   = monitordestroy;
364:   adaptor->monitorcontext[adaptor->numbermonitors++] = ctx;
365:   PetscFunctionReturn(PETSC_SUCCESS);
366: }

368: /*@
369:   DMAdaptorMonitorCancel - Clears all monitors for a `DMAdaptor` object.

371:   Logically Collective

373:   Input Parameter:
374: . adaptor - the `DMAdaptor`

376:   Options Database Key:
377: . -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.

379:   Level: intermediate

381: .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptorMonitorSet()`, `DMAdaptor`
382: @*/
383: PetscErrorCode DMAdaptorMonitorCancel(DMAdaptor adaptor)
384: {
385:   PetscFunctionBegin;
387:   for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) {
388:     if (adaptor->monitordestroy[i]) PetscCall((*adaptor->monitordestroy[i])(&adaptor->monitorcontext[i]));
389:   }
390:   adaptor->numbermonitors = 0;
391:   PetscFunctionReturn(PETSC_SUCCESS);
392: }

394: /*@C
395:   DMAdaptorMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user in the options database

397:   Collective

399:   Input Parameters:
400: + adaptor - `DMadaptor` object you wish to monitor
401: . opt     - the command line option for this monitor
402: . name    - the monitor type one is seeking
403: - ctx     - An optional user context for the monitor, or `NULL`

405:   Level: developer

407: .seealso: [](ch_snes), `DMAdaptorMonitorRegister()`, `DMAdaptorMonitorSet()`, `PetscOptionsGetViewer()`
408: @*/
409: PetscErrorCode DMAdaptorMonitorSetFromOptions(DMAdaptor adaptor, const char opt[], const char name[], void *ctx)
410: {
411:   PetscErrorCode (*mfunc)(DMAdaptor, PetscInt, DM, DM, PetscInt, PetscReal[], Vec, void *);
412:   PetscErrorCode (*cfunc)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
413:   PetscErrorCode (*dfunc)(PetscViewerAndFormat **);
414:   PetscViewerAndFormat *vf;
415:   PetscViewer           viewer;
416:   PetscViewerFormat     format;
417:   PetscViewerType       vtype;
418:   char                  key[PETSC_MAX_PATH_LEN];
419:   PetscBool             flg;
420:   const char           *prefix = NULL;

422:   PetscFunctionBegin;
423:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
424:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)adaptor), ((PetscObject)adaptor)->options, prefix, opt, &viewer, &format, &flg));
425:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

427:   PetscCall(PetscViewerGetType(viewer, &vtype));
428:   PetscCall(DMAdaptorMonitorMakeKey_Internal(name, vtype, format, key));
429:   PetscCall(PetscFunctionListFind(DMAdaptorMonitorList, key, &mfunc));
430:   PetscCall(PetscFunctionListFind(DMAdaptorMonitorCreateList, key, &cfunc));
431:   PetscCall(PetscFunctionListFind(DMAdaptorMonitorDestroyList, key, &dfunc));
432:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
433:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

435:   PetscCall((*cfunc)(viewer, format, ctx, &vf));
436:   PetscCall(PetscViewerDestroy(&viewer));
437:   PetscCall(DMAdaptorMonitorSet(adaptor, mfunc, vf, (PetscCtxDestroyFn *)dfunc));
438:   PetscFunctionReturn(PETSC_SUCCESS);
439: }

441: /*@
442:   DMAdaptorSetOptionsPrefix - Sets the prefix used for searching for all `DMAdaptor` options in the database.

444:   Logically Collective

446:   Input Parameters:
447: + adaptor - the `DMAdaptor`
448: - prefix  - the prefix to prepend to all option names

450:   Level: advanced

452:   Note:
453:   A hyphen (-) must NOT be given at the beginning of the prefix name.
454:   The first character of all runtime options is AUTOMATICALLY the hyphen.

456: .seealso: [](ch_snes), `DMAdaptor`, `SNESSetOptionsPrefix()`, `DMAdaptorSetFromOptions()`
457: @*/
458: PetscErrorCode DMAdaptorSetOptionsPrefix(DMAdaptor adaptor, const char prefix[])
459: {
460:   PetscFunctionBegin;
462:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor, prefix));
463:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor->refineTag, prefix));
464:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)adaptor->refineTag, "refine_"));
465:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor->coarsenTag, prefix));
466:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)adaptor->coarsenTag, "coarsen_"));
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: /*@
471:   DMAdaptorSetFromOptions - Sets properties of a `DMAdaptor` object from values in the options database

473:   Collective

475:   Input Parameter:
476: . adaptor - The `DMAdaptor` object

478:   Options Database Keys:
479: + -adaptor_monitor_size                - Monitor the mesh size
480: . -adaptor_monitor_error               - Monitor the solution error
481: . -adaptor_sequence_num <num>          - Number of adaptations to generate an optimal grid
482: . -adaptor_target_num <num>            - Set the target number of vertices N_adapt, -1 for automatic determination
483: . -adaptor_refinement_factor <r>       - Set r such that N_adapt = r^dim N_orig
484: - -adaptor_mixed_setup_function <func> - Set the function func that sets up the mixed problem

486:   Level: beginner

488: .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
489: @*/
490: PetscErrorCode DMAdaptorSetFromOptions(DMAdaptor adaptor)
491: {
492:   char                  typeName[PETSC_MAX_PATH_LEN];
493:   const char           *defName = DMADAPTORGRADIENT;
494:   char                  funcname[PETSC_MAX_PATH_LEN];
495:   DMAdaptationCriterion criterion = DM_ADAPTATION_NONE;
496:   PetscBool             flg;

498:   PetscFunctionBegin;
499:   PetscObjectOptionsBegin((PetscObject)adaptor);
500:   PetscCall(PetscOptionsFList("-adaptor_type", "DMAdaptor", "DMAdaptorSetType", DMAdaptorList, defName, typeName, 1024, &flg));
501:   if (flg) PetscCall(DMAdaptorSetType(adaptor, typeName));
502:   else if (!((PetscObject)adaptor)->type_name) PetscCall(DMAdaptorSetType(adaptor, defName));
503:   PetscCall(PetscOptionsEnum("-adaptor_criterion", "Criterion used to drive adaptation", "", DMAdaptationCriteria, (PetscEnum)criterion, (PetscEnum *)&criterion, &flg));
504:   if (flg) PetscCall(DMAdaptorSetCriterion(adaptor, criterion));
505:   PetscCall(PetscOptionsInt("-adaptor_sequence_num", "Number of adaptations to generate an optimal grid", "DMAdaptorSetSequenceLength", adaptor->numSeq, &adaptor->numSeq, NULL));
506:   PetscCall(PetscOptionsInt("-adaptor_target_num", "Set the target number of vertices N_adapt, -1 for automatic determination", "DMAdaptor", adaptor->Nadapt, &adaptor->Nadapt, NULL));
507:   PetscCall(PetscOptionsReal("-adaptor_refinement_factor", "Set r such that N_adapt = r^dim N_orig", "DMAdaptor", adaptor->refinementFactor, &adaptor->refinementFactor, NULL));
508:   PetscCall(PetscOptionsString("-adaptor_mixed_setup_function", "Function to setup the mixed problem", "DMAdaptorSetMixedSetupFunction", funcname, funcname, sizeof(funcname), &flg));
509:   if (flg) {
510:     PetscErrorCode (*setupFunc)(DMAdaptor, DM);

512:     PetscCall(PetscDLSym(NULL, funcname, (void **)&setupFunc));
513:     PetscCheck(setupFunc, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
514:     PetscCall(DMAdaptorSetMixedSetupFunction(adaptor, setupFunc));
515:   }
516:   PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_size", "size", adaptor));
517:   PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_error", "error", adaptor));
518:   PetscOptionsEnd();
519:   PetscCall(VecTaggerSetFromOptions(adaptor->refineTag));
520:   PetscCall(VecTaggerSetFromOptions(adaptor->coarsenTag));
521:   PetscFunctionReturn(PETSC_SUCCESS);
522: }

524: /*@
525:   DMAdaptorView - Views a `DMAdaptor` object

527:   Collective

529:   Input Parameters:
530: + adaptor - The `DMAdaptor` object
531: - viewer  - The `PetscViewer` object

533:   Level: beginner

535: .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
536: @*/
537: PetscErrorCode DMAdaptorView(DMAdaptor adaptor, PetscViewer viewer)
538: {
539:   PetscFunctionBegin;
540:   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adaptor, viewer));
541:   PetscCall(PetscViewerASCIIPrintf(viewer, "DM Adaptor\n"));
542:   PetscCall(PetscViewerASCIIPrintf(viewer, "  sequence length: %" PetscInt_FMT "\n", adaptor->numSeq));
543:   PetscCall(VecTaggerView(adaptor->refineTag, viewer));
544:   PetscCall(VecTaggerView(adaptor->coarsenTag, viewer));
545:   PetscFunctionReturn(PETSC_SUCCESS);
546: }

548: /*@
549:   DMAdaptorGetSolver - Gets the solver used to produce discrete solutions

551:   Not Collective

553:   Input Parameter:
554: . adaptor - The `DMAdaptor` object

556:   Output Parameter:
557: . snes - The solver

559:   Level: intermediate

561: .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
562: @*/
563: PetscErrorCode DMAdaptorGetSolver(DMAdaptor adaptor, SNES *snes)
564: {
565:   PetscFunctionBegin;
567:   PetscAssertPointer(snes, 2);
568:   *snes = adaptor->snes;
569:   PetscFunctionReturn(PETSC_SUCCESS);
570: }

572: /*@
573:   DMAdaptorSetSolver - Sets the solver used to produce discrete solutions

575:   Not Collective

577:   Input Parameters:
578: + adaptor - The `DMAdaptor` object
579: - snes    - The solver, this MUST have an attached `DM`/`PetscDS`, so that the exact solution can be computed

581:   Level: intermediate

583: .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorGetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
584: @*/
585: PetscErrorCode DMAdaptorSetSolver(DMAdaptor adaptor, SNES snes)
586: {
587:   PetscFunctionBegin;
590:   adaptor->snes = snes;
591:   PetscCall(SNESGetDM(adaptor->snes, &adaptor->idm));
592:   PetscFunctionReturn(PETSC_SUCCESS);
593: }

595: /*@
596:   DMAdaptorGetSequenceLength - Gets the number of sequential adaptations used by an adapter

598:   Not Collective

600:   Input Parameter:
601: . adaptor - The `DMAdaptor` object

603:   Output Parameter:
604: . num - The number of adaptations

606:   Level: intermediate

608: .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorSetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
609: @*/
610: PetscErrorCode DMAdaptorGetSequenceLength(DMAdaptor adaptor, PetscInt *num)
611: {
612:   PetscFunctionBegin;
614:   PetscAssertPointer(num, 2);
615:   *num = adaptor->numSeq;
616:   PetscFunctionReturn(PETSC_SUCCESS);
617: }

619: /*@
620:   DMAdaptorSetSequenceLength - Sets the number of sequential adaptations

622:   Not Collective

624:   Input Parameters:
625: + adaptor - The `DMAdaptor` object
626: - num     - The number of adaptations

628:   Level: intermediate

630: .seealso: [](ch_dmbase), `DMAdaptorGetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
631: @*/
632: PetscErrorCode DMAdaptorSetSequenceLength(DMAdaptor adaptor, PetscInt num)
633: {
634:   PetscFunctionBegin;
636:   adaptor->numSeq = num;
637:   PetscFunctionReturn(PETSC_SUCCESS);
638: }

640: static PetscErrorCode DMAdaptorTransferSolution_Exact_Private(DMAdaptor adaptor, DM dm, Vec u, DM adm, Vec au, void *ctx)
641: {
642:   PetscFunctionBeginUser;
643:   PetscCall(DMProjectFunction(adm, 0.0, adaptor->exactSol, adaptor->exactCtx, INSERT_ALL_VALUES, au));
644:   PetscFunctionReturn(PETSC_SUCCESS);
645: }

647: /*@
648:   DMAdaptorSetUp - After the solver is specified, creates data structures for controlling adaptivity

650:   Collective

652:   Input Parameter:
653: . adaptor - The `DMAdaptor` object

655:   Level: beginner

657: .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
658: @*/
659: PetscErrorCode DMAdaptorSetUp(DMAdaptor adaptor)
660: {
661:   PetscDS  prob;
662:   PetscInt Nf, f;

664:   PetscFunctionBegin;
665:   PetscCall(DMGetDS(adaptor->idm, &prob));
666:   PetscCall(VecTaggerSetUp(adaptor->refineTag));
667:   PetscCall(VecTaggerSetUp(adaptor->coarsenTag));
668:   PetscCall(PetscDSGetNumFields(prob, &Nf));
669:   PetscCall(PetscMalloc2(Nf, &adaptor->exactSol, Nf, &adaptor->exactCtx));
670:   for (f = 0; f < Nf; ++f) {
671:     PetscCall(PetscDSGetExactSolution(prob, f, &adaptor->exactSol[f], &adaptor->exactCtx[f]));
672:     /* TODO Have a flag that forces projection rather than using the exact solution */
673:     if (adaptor->exactSol[0]) PetscCall(DMAdaptorSetTransferFunction(adaptor, DMAdaptorTransferSolution_Exact_Private));
674:   }
675:   PetscTryTypeMethod(adaptor, setup);
676:   PetscFunctionReturn(PETSC_SUCCESS);
677: }

679: PetscErrorCode DMAdaptorGetTransferFunction(DMAdaptor adaptor, PetscErrorCode (**tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
680: {
681:   PetscFunctionBegin;
682:   *tfunc = adaptor->ops->transfersolution;
683:   PetscFunctionReturn(PETSC_SUCCESS);
684: }

686: PetscErrorCode DMAdaptorSetTransferFunction(DMAdaptor adaptor, PetscErrorCode (*tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
687: {
688:   PetscFunctionBegin;
689:   adaptor->ops->transfersolution = tfunc;
690:   PetscFunctionReturn(PETSC_SUCCESS);
691: }

693: static PetscErrorCode DMAdaptorPreAdapt(DMAdaptor adaptor, Vec locX)
694: {
695:   DM           plex;
696:   PetscDS      prob;
697:   PetscObject  obj;
698:   PetscClassId id;
699:   PetscBool    isForest;

701:   PetscFunctionBegin;
702:   PetscCall(DMConvert(adaptor->idm, DMPLEX, &plex));
703:   PetscCall(DMGetDS(adaptor->idm, &prob));
704:   PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
705:   PetscCall(PetscObjectGetClassId(obj, &id));
706:   PetscCall(DMIsForest(adaptor->idm, &isForest));
707:   if (adaptor->adaptCriterion == DM_ADAPTATION_NONE) {
708:     if (isForest) {
709:       adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
710:     }
711: #if defined(PETSC_HAVE_PRAGMATIC)
712:     else {
713:       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
714:     }
715: #elif defined(PETSC_HAVE_MMG)
716:     else {
717:       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
718:     }
719: #elif defined(PETSC_HAVE_PARMMG)
720:     else {
721:       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
722:     }
723: #else
724:     else {
725:       adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
726:     }
727: #endif
728:   }
729:   if (id == PETSCFV_CLASSID) {
730:     adaptor->femType = PETSC_FALSE;
731:   } else {
732:     adaptor->femType = PETSC_TRUE;
733:   }
734:   if (adaptor->femType) {
735:     /* Compute local solution bc */
736:     PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
737:   } else {
738:     PetscFV      fvm = (PetscFV)obj;
739:     PetscLimiter noneLimiter;
740:     Vec          grad;

742:     PetscCall(PetscFVGetComputeGradients(fvm, &adaptor->computeGradient));
743:     PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
744:     /* Use no limiting when reconstructing gradients for adaptivity */
745:     PetscCall(PetscFVGetLimiter(fvm, &adaptor->limiter));
746:     PetscCall(PetscObjectReference((PetscObject)adaptor->limiter));
747:     PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
748:     PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
749:     PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
750:     /* Get FVM data */
751:     PetscCall(DMPlexGetDataFVM(plex, fvm, &adaptor->cellGeom, &adaptor->faceGeom, &adaptor->gradDM));
752:     PetscCall(VecGetDM(adaptor->cellGeom, &adaptor->cellDM));
753:     PetscCall(VecGetArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
754:     /* Compute local solution bc */
755:     PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
756:     /* Compute gradients */
757:     PetscCall(DMCreateGlobalVector(adaptor->gradDM, &grad));
758:     PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
759:     PetscCall(DMGetLocalVector(adaptor->gradDM, &adaptor->cellGrad));
760:     PetscCall(DMGlobalToLocalBegin(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
761:     PetscCall(DMGlobalToLocalEnd(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
762:     PetscCall(VecDestroy(&grad));
763:     PetscCall(VecGetArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
764:   }
765:   PetscCall(DMDestroy(&plex));
766:   PetscFunctionReturn(PETSC_SUCCESS);
767: }

769: static PetscErrorCode DMAdaptorTransferSolution(DMAdaptor adaptor, DM dm, Vec x, DM adm, Vec ax)
770: {
771:   PetscReal time = 0.0;
772:   Mat       interp;
773:   void     *ctx;

775:   PetscFunctionBegin;
776:   PetscCall(DMGetApplicationContext(dm, &ctx));
777:   if (adaptor->ops->transfersolution) PetscUseTypeMethod(adaptor, transfersolution, dm, x, adm, ax, ctx);
778:   else {
779:     switch (adaptor->adaptCriterion) {
780:     case DM_ADAPTATION_LABEL:
781:       PetscCall(DMForestTransferVec(dm, x, adm, ax, PETSC_TRUE, time));
782:       break;
783:     case DM_ADAPTATION_REFINE:
784:     case DM_ADAPTATION_METRIC:
785:       PetscCall(DMCreateInterpolation(dm, adm, &interp, NULL));
786:       PetscCall(MatInterpolate(interp, x, ax));
787:       PetscCall(DMInterpolate(dm, interp, adm));
788:       PetscCall(MatDestroy(&interp));
789:       break;
790:     default:
791:       SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_SUP, "No built-in projection for this adaptation criterion: %d", adaptor->adaptCriterion);
792:     }
793:   }
794:   PetscFunctionReturn(PETSC_SUCCESS);
795: }

797: static PetscErrorCode DMAdaptorPostAdapt(DMAdaptor adaptor)
798: {
799:   PetscDS      prob;
800:   PetscObject  obj;
801:   PetscClassId id;

803:   PetscFunctionBegin;
804:   PetscCall(DMGetDS(adaptor->idm, &prob));
805:   PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
806:   PetscCall(PetscObjectGetClassId(obj, &id));
807:   if (id == PETSCFV_CLASSID) {
808:     PetscFV fvm = (PetscFV)obj;

810:     PetscCall(PetscFVSetComputeGradients(fvm, adaptor->computeGradient));
811:     /* Restore original limiter */
812:     PetscCall(PetscFVSetLimiter(fvm, adaptor->limiter));

814:     PetscCall(VecRestoreArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
815:     PetscCall(VecRestoreArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
816:     PetscCall(DMRestoreLocalVector(adaptor->gradDM, &adaptor->cellGrad));
817:   }
818:   PetscFunctionReturn(PETSC_SUCCESS);
819: }

821: /*
822:   DMAdaptorComputeCellErrorIndicator_Gradient - Use the integrated gradient as an error indicator in the `DMAdaptor`

824:   Input Parameters:
825: + adaptor  - The `DMAdaptor` object
826: . dim      - The topological dimension
827: . cell     - The cell
828: . field    - The field integrated over the cell
829: . gradient - The gradient integrated over the cell
830: . cg       - A `PetscFVCellGeom` struct
831: - ctx      - A user context

833:   Output Parameter:
834: . errInd   - The error indicator

836:   Developer Note:
837:   Some of the input arguments are absurdly specialized to special situations, it is not clear this is a good general API

839: .seealso: [](ch_dmbase), `DMAdaptor`
840: */
841: static PetscErrorCode DMAdaptorComputeCellErrorIndicator_Gradient(DMAdaptor adaptor, PetscInt dim, PetscInt Nc, const PetscScalar *field, const PetscScalar *gradient, const PetscFVCellGeom *cg, PetscReal *errInd, void *ctx)
842: {
843:   PetscReal err = 0.;
844:   PetscInt  c, d;

846:   PetscFunctionBeginHot;
847:   for (c = 0; c < Nc; c++) {
848:     for (d = 0; d < dim; ++d) err += PetscSqr(PetscRealPart(gradient[c * dim + d]));
849:   }
850:   *errInd = cg->volume * err;
851:   PetscFunctionReturn(PETSC_SUCCESS);
852: }

854: static PetscErrorCode DMAdaptorComputeErrorIndicator_Gradient(DMAdaptor adaptor, Vec locX, Vec errVec)
855: {
856:   DM              dm, plex, edm, eplex;
857:   PetscDS         ds;
858:   PetscObject     obj;
859:   PetscClassId    id;
860:   void           *ctx;
861:   PetscQuadrature quad;
862:   PetscScalar    *earray;
863:   PetscReal       minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2];
864:   PetscInt        dim, cdim, cStart, cEnd, Nf, Nc;

866:   PetscFunctionBegin;
867:   PetscCall(VecGetDM(locX, &dm));
868:   PetscCall(DMConvert(dm, DMPLEX, &plex));
869:   PetscCall(VecGetDM(errVec, &edm));
870:   PetscCall(DMConvert(edm, DMPLEX, &eplex));
871:   PetscCall(DMGetDimension(plex, &dim));
872:   PetscCall(DMGetCoordinateDim(plex, &cdim));
873:   PetscCall(DMGetApplicationContext(plex, &ctx));
874:   PetscCall(DMGetDS(plex, &ds));
875:   PetscCall(PetscDSGetNumFields(ds, &Nf));
876:   PetscCall(PetscDSGetDiscretization(ds, 0, &obj));
877:   PetscCall(PetscObjectGetClassId(obj, &id));

879:   PetscCall(VecGetArray(errVec, &earray));
880:   PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
881:   for (PetscInt cell = cStart; cell < cEnd; ++cell) {
882:     PetscScalar *eval;
883:     PetscReal    errInd = 0.;

885:     if (id == PETSCFV_CLASSID) {
886:       PetscFV            fv = (PetscFV)obj;
887:       const PetscScalar *pointSols;
888:       const PetscScalar *pointSol;
889:       const PetscScalar *pointGrad;
890:       PetscFVCellGeom   *cg;

892:       PetscCall(PetscFVGetNumComponents(fv, &Nc));
893:       PetscCall(VecGetArrayRead(locX, &pointSols));
894:       PetscCall(DMPlexPointLocalRead(plex, cell, pointSols, (void *)&pointSol));
895:       PetscCall(DMPlexPointLocalRead(adaptor->gradDM, cell, adaptor->cellGradArray, (void *)&pointGrad));
896:       PetscCall(DMPlexPointLocalRead(adaptor->cellDM, cell, adaptor->cellGeomArray, &cg));
897:       PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, pointSol, pointGrad, cg, &errInd, ctx);
898:       PetscCall(VecRestoreArrayRead(locX, &pointSols));
899:     } else {
900:       PetscFE          fe = (PetscFE)obj;
901:       PetscScalar     *x  = NULL, *field, *gradient, *interpolant, *interpolantGrad;
902:       PetscFVCellGeom  cg;
903:       PetscFEGeom      fegeom;
904:       const PetscReal *quadWeights;
905:       PetscReal       *coords;
906:       PetscInt         Nb, Nq, qNc;

908:       fegeom.dim      = dim;
909:       fegeom.dimEmbed = cdim;
910:       PetscCall(PetscFEGetNumComponents(fe, &Nc));
911:       PetscCall(PetscFEGetQuadrature(fe, &quad));
912:       PetscCall(PetscFEGetDimension(fe, &Nb));
913:       PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
914:       PetscCall(PetscCalloc6(Nc, &field, cdim * Nc, &gradient, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ));
915:       PetscCall(PetscMalloc2(Nc, &interpolant, cdim * Nc, &interpolantGrad));
916:       PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
917:       PetscCall(DMPlexComputeCellGeometryFVM(plex, cell, &cg.volume, NULL, NULL));
918:       PetscCall(PetscArrayzero(gradient, cdim * Nc));
919:       PetscCall(DMPlexVecGetClosure(plex, NULL, locX, cell, NULL, &x));
920:       for (PetscInt f = 0; f < Nf; ++f) {
921:         PetscInt qc = 0;

923:         PetscCall(PetscDSGetDiscretization(ds, f, &obj));
924:         PetscCall(PetscArrayzero(interpolant, Nc));
925:         PetscCall(PetscArrayzero(interpolantGrad, cdim * Nc));
926:         for (PetscInt q = 0; q < Nq; ++q) {
927:           PetscCall(PetscFEInterpolateFieldAndGradient_Static((PetscFE)obj, 1, x, &fegeom, q, interpolant, interpolantGrad));
928:           for (PetscInt fc = 0; fc < Nc; ++fc) {
929:             const PetscReal wt = quadWeights[q * qNc + qc + fc];

931:             field[fc] += interpolant[fc] * wt * fegeom.detJ[q];
932:             for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] += interpolantGrad[fc * dim + d] * wt * fegeom.detJ[q];
933:           }
934:         }
935:         qc += Nc;
936:       }
937:       PetscCall(PetscFree2(interpolant, interpolantGrad));
938:       PetscCall(DMPlexVecRestoreClosure(plex, NULL, locX, cell, NULL, &x));
939:       for (PetscInt fc = 0; fc < Nc; ++fc) {
940:         field[fc] /= cg.volume;
941:         for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] /= cg.volume;
942:       }
943:       PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, field, gradient, &cg, &errInd, ctx);
944:       PetscCall(PetscFree6(field, gradient, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
945:     }
946:     PetscCall(DMPlexPointGlobalRef(eplex, cell, earray, (void *)&eval));
947:     eval[0]      = errInd;
948:     minMaxInd[0] = PetscMin(minMaxInd[0], errInd);
949:     minMaxInd[1] = PetscMax(minMaxInd[1], errInd);
950:   }
951:   PetscCall(VecRestoreArray(errVec, &earray));
952:   PetscCall(DMDestroy(&plex));
953:   PetscCall(DMDestroy(&eplex));
954:   PetscCall(PetscGlobalMinMaxReal(PetscObjectComm((PetscObject)adaptor), minMaxInd, minMaxIndGlobal));
955:   PetscCall(PetscInfo(adaptor, "DMAdaptor: error indicator range (%g, %g)\n", (double)minMaxIndGlobal[0], (double)minMaxIndGlobal[1]));
956:   PetscFunctionReturn(PETSC_SUCCESS);
957: }

959: static PetscErrorCode DMAdaptorComputeErrorIndicator_Flux(DMAdaptor adaptor, Vec lu, Vec errVec)
960: {
961:   DM          dm, mdm;
962:   SNES        msnes;
963:   Vec         mu, lmu;
964:   void       *ctx;
965:   const char *prefix;

967:   PetscFunctionBegin;
968:   PetscCall(VecGetDM(lu, &dm));

970:   // Set up and solve mixed problem
971:   PetscCall(DMClone(dm, &mdm));
972:   PetscCall(SNESCreate(PetscObjectComm((PetscObject)mdm), &msnes));
973:   PetscCall(SNESSetDM(msnes, mdm));
974:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
975:   PetscCall(SNESSetOptionsPrefix(msnes, prefix));
976:   PetscCall(SNESAppendOptionsPrefix(msnes, "mixed_"));

978:   PetscTryTypeMethod(adaptor, mixedsetup, mdm);
979:   PetscCall(DMGetApplicationContext(dm, &ctx));
980:   PetscCall(DMPlexSetSNESLocalFEM(mdm, PETSC_FALSE, ctx));
981:   PetscCall(SNESSetFromOptions(msnes));

983:   PetscCall(DMCreateGlobalVector(mdm, &mu));
984:   PetscCall(PetscObjectSetName((PetscObject)mu, "Mixed Solution"));
985:   PetscCall(VecSet(mu, 0.0));
986:   PetscCall(SNESSolve(msnes, NULL, mu));
987:   PetscCall(VecViewFromOptions(mu, (PetscObject)adaptor, "-adapt_mixed_sol_vec_view"));

989:   PetscCall(DMGetLocalVector(mdm, &lmu));
990:   PetscCall(DMGlobalToLocal(mdm, mu, INSERT_VALUES, lmu));
991:   PetscCall(DMPlexInsertBoundaryValues(mdm, PETSC_TRUE, lmu, 0.0, NULL, NULL, NULL));
992:   PetscCall(DMPlexComputeL2FluxDiffVecLocal(lu, 0, lmu, 0, errVec));
993:   PetscCall(DMRestoreLocalVector(mdm, &lmu));
994:   PetscCall(VecDestroy(&mu));
995:   PetscCall(SNESDestroy(&msnes));
996:   PetscCall(DMDestroy(&mdm));
997:   PetscFunctionReturn(PETSC_SUCCESS);
998: }

1000: /*@
1001:   DMAdaptorMonitor - runs the user provided monitor routines, if they exist

1003:   Collective

1005:   Input Parameters:
1006: + adaptor - the `DMAdaptor`
1007: . it      - iteration number
1008: . odm     - the original `DM`
1009: . adm     - the adapted `DM`
1010: . Nf      - the number of fields
1011: . enorms  - the 2-norm error values for each field
1012: - error   - `Vec` of cellwise errors

1014:   Level: developer

1016:   Note:
1017:   This routine is called by the `DMAdaptor` implementations.
1018:   It does not typically need to be called by the user.

1020: .seealso: [](ch_snes), `DMAdaptorMonitorSet()`
1021: @*/
1022: PetscErrorCode DMAdaptorMonitor(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error)
1023: {
1024:   PetscFunctionBegin;
1025:   for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) PetscCall((*adaptor->monitor[i])(adaptor, it, odm, adm, Nf, enorms, error, adaptor->monitorcontext[i]));
1026:   PetscFunctionReturn(PETSC_SUCCESS);
1027: }

1029: /*@C
1030:   DMAdaptorMonitorSize - Prints the mesh sizes at each iteration of an adaptation loop.

1032:   Collective

1034:   Input Parameters:
1035: + adaptor - the `DMAdaptor`
1036: . n       - iteration number
1037: . odm     - the original `DM`
1038: . adm     - the adapted `DM`
1039: . Nf      - number of fields
1040: . enorms  - 2-norm error values for each field (may be estimated).
1041: . error   - `Vec` of cellwise errors
1042: - vf      - The viewer context

1044:   Options Database Key:
1045: . -adaptor_monitor_size - Activates `DMAdaptorMonitorSize()`

1047:   Level: intermediate

1049:   Note:
1050:   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1051:   to be used during the adaptation loop.

1053: .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorError()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
1054: @*/
1055: PetscErrorCode DMAdaptorMonitorSize(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1056: {
1057:   PetscViewer       viewer = vf->viewer;
1058:   PetscViewerFormat format = vf->format;
1059:   PetscInt          tablevel, cStart, cEnd, acStart, acEnd;
1060:   const char       *prefix;
1061:   PetscMPIInt       rank;

1063:   PetscFunctionBegin;
1065:   PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
1066:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
1067:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)adaptor), &rank));
1068:   PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
1069:   PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));

1071:   PetscCall(PetscViewerPushFormat(viewer, format));
1072:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1073:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Sizes for %s adaptation.\n", prefix));
1074:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor rank %d N_orig: %" PetscInt_FMT " N_adapt: %" PetscInt_FMT "\n", n, rank, cEnd - cStart, acEnd - acStart));
1075:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1076:   PetscCall(PetscViewerPopFormat(viewer));
1077:   PetscFunctionReturn(PETSC_SUCCESS);
1078: }

1080: /*@C
1081:   DMAdaptorMonitorError - Prints the error norm at each iteration of an adaptation loop.

1083:   Collective

1085:   Input Parameters:
1086: + adaptor - the `DMAdaptor`
1087: . n       - iteration number
1088: . odm     - the original `DM`
1089: . adm     - the adapted `DM`
1090: . Nf      - number of fields
1091: . enorms  - 2-norm error values for each field (may be estimated).
1092: . error   - `Vec` of cellwise errors
1093: - vf      - The viewer context

1095:   Options Database Key:
1096: . -adaptor_monitor_error - Activates `DMAdaptorMonitorError()`

1098:   Level: intermediate

1100:   Note:
1101:   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1102:   to be used during the adaptation loop.

1104: .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
1105: @*/
1106: PetscErrorCode DMAdaptorMonitorError(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1107: {
1108:   PetscViewer       viewer = vf->viewer;
1109:   PetscViewerFormat format = vf->format;
1110:   PetscInt          tablevel, cStart, cEnd, acStart, acEnd;
1111:   const char       *prefix;

1113:   PetscFunctionBegin;
1115:   PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
1116:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));

1118:   PetscCall(PetscViewerPushFormat(viewer, format));
1119:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1120:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Error norms for %s adaptation.\n", prefix));
1121:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor Error norm %s", n, Nf > 1 ? "[" : ""));
1122:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1123:   for (PetscInt f = 0; f < Nf; ++f) {
1124:     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1125:     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)enorms[f]));
1126:   }
1127:   PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
1128:   PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
1129:   PetscCall(PetscViewerASCIIPrintf(viewer, " N: %" PetscInt_FMT " Nadapt: %" PetscInt_FMT "\n", cEnd - cStart, acEnd - acStart));
1130:   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1131:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1132:   PetscCall(PetscViewerPopFormat(viewer));
1133:   PetscFunctionReturn(PETSC_SUCCESS);
1134: }

1136: /*@C
1137:   DMAdaptorMonitorErrorDraw - Plots the error at each iteration of an iterative solver.

1139:   Collective

1141:   Input Parameters:
1142: + adaptor - the `DMAdaptor`
1143: . n       - iteration number
1144: . odm     - the original `DM`
1145: . adm     - the adapted `DM`
1146: . Nf      - number of fields
1147: . enorms  - 2-norm error values for each field (may be estimated).
1148: . error   - `Vec` of cellwise errors
1149: - vf      - The viewer context

1151:   Options Database Key:
1152: . -adaptor_monitor_error draw - Activates `DMAdaptorMonitorErrorDraw()`

1154:   Level: intermediate

1156:   Note:
1157:   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1158:   to be used during the adaptation loop.

1160: .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
1161: @*/
1162: PetscErrorCode DMAdaptorMonitorErrorDraw(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1163: {
1164:   PetscViewer       viewer = vf->viewer;
1165:   PetscViewerFormat format = vf->format;

1167:   PetscFunctionBegin;
1169:   PetscCall(PetscViewerPushFormat(viewer, format));
1170:   PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
1171:   PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", (PetscObject)adaptor));
1172:   PetscCall(VecView(error, viewer));
1173:   PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", NULL));
1174:   PetscCall(PetscViewerPopFormat(viewer));
1175:   PetscFunctionReturn(PETSC_SUCCESS);
1176: }

1178: /*@C
1179:   DMAdaptorMonitorErrorDrawLGCreate - Creates the context for the error plotter `DMAdaptorMonitorErrorDrawLG()`

1181:   Collective

1183:   Input Parameters:
1184: + viewer - The `PetscViewer`
1185: . format - The viewer format
1186: - ctx    - An optional user context

1188:   Output Parameter:
1189: . vf - The viewer context

1191:   Level: intermediate

1193: .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
1194: @*/
1195: PetscErrorCode DMAdaptorMonitorErrorDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
1196: {
1197:   DMAdaptor adaptor = (DMAdaptor)ctx;
1198:   char    **names;
1199:   PetscInt  Nf;

1201:   PetscFunctionBegin;
1202:   PetscCall(DMGetNumFields(adaptor->idm, &Nf));
1203:   PetscCall(PetscMalloc1(Nf + 1, &names));
1204:   for (PetscInt f = 0; f < Nf; ++f) {
1205:     PetscObject disc;
1206:     const char *fname;
1207:     char        lname[PETSC_MAX_PATH_LEN];

1209:     PetscCall(DMGetField(adaptor->idm, f, NULL, &disc));
1210:     PetscCall(PetscObjectGetName(disc, &fname));
1211:     PetscCall(PetscStrncpy(lname, fname, PETSC_MAX_PATH_LEN));
1212:     PetscCall(PetscStrlcat(lname, " Error", PETSC_MAX_PATH_LEN));
1213:     PetscCall(PetscStrallocpy(lname, &names[f]));
1214:   }
1215:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
1216:   (*vf)->data = ctx;
1217:   PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Error Norm", Nf, (const char **)names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
1218:   for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(names[f]));
1219:   PetscCall(PetscFree(names));
1220:   PetscFunctionReturn(PETSC_SUCCESS);
1221: }

1223: /*@C
1224:   DMAdaptorMonitorErrorDrawLG - Plots the error norm at each iteration of an adaptive loop.

1226:   Collective

1228:   Input Parameters:
1229: + adaptor - the `DMAdaptor`
1230: . n       - iteration number
1231: . odm     - the original `DM`
1232: . adm     - the adapted `DM`
1233: . Nf      - number of fields
1234: . enorms  - 2-norm error values for each field (may be estimated).
1235: . error   - `Vec` of cellwise errors
1236: - vf      - The viewer context

1238:   Options Database Key:
1239: . -adaptor_error draw::draw_lg - Activates `DMAdaptorMonitorErrorDrawLG()`

1241:   Level: intermediate

1243:   Notes:
1244:   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1245:   to be used during the adaptation loop.

1247:   Call `DMAdaptorMonitorErrorDrawLGCreate()` to create the context needed for this monitor

1249: .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdadptorMonitorErrorDraw()`, `DMAdadptorMonitorError()`,
1250:           `DMAdaptorMonitorTrueResidualDrawLGCreate()`
1251: @*/
1252: PetscErrorCode DMAdaptorMonitorErrorDrawLG(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1253: {
1254:   PetscViewer       viewer = vf->viewer;
1255:   PetscViewerFormat format = vf->format;
1256:   PetscDrawLG       lg     = vf->lg;
1257:   PetscReal        *x, *e;

1259:   PetscFunctionBegin;
1262:   PetscCall(PetscCalloc2(Nf, &x, Nf, &e));
1263:   PetscCall(PetscViewerPushFormat(viewer, format));
1264:   if (!n) PetscCall(PetscDrawLGReset(lg));
1265:   for (PetscInt f = 0; f < Nf; ++f) {
1266:     x[f] = (PetscReal)n;
1267:     e[f] = enorms[f] > 0.0 ? PetscLog10Real(enorms[f]) : -15.;
1268:   }
1269:   PetscCall(PetscDrawLGAddPoint(lg, x, e));
1270:   PetscCall(PetscDrawLGDraw(lg));
1271:   PetscCall(PetscDrawLGSave(lg));
1272:   PetscCall(PetscViewerPopFormat(viewer));
1273:   PetscCall(PetscFree2(x, e));
1274:   PetscFunctionReturn(PETSC_SUCCESS);
1275: }

1277: /*@C
1278:   DMAdaptorMonitorRegisterAll - Registers all of the mesh adaptation monitors in the `SNES` package.

1280:   Not Collective

1282:   Level: advanced

1284: .seealso: [](ch_snes), `SNES`, `DM`, `DMAdaptorMonitorRegister()`, `DMAdaptorRegister()`
1285: @*/
1286: PetscErrorCode DMAdaptorMonitorRegisterAll(void)
1287: {
1288:   PetscFunctionBegin;
1289:   if (DMAdaptorMonitorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1290:   DMAdaptorMonitorRegisterAllCalled = PETSC_TRUE;

1292:   PetscCall(DMAdaptorMonitorRegister("size", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorSize, NULL, NULL));
1293:   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorError, NULL, NULL));
1294:   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorErrorDraw, NULL, NULL));
1295:   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DRAW_LG, DMAdaptorMonitorErrorDrawLG, DMAdaptorMonitorErrorDrawLGCreate, NULL));
1296:   PetscFunctionReturn(PETSC_SUCCESS);
1297: }

1299: 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[])
1300: {
1301:   const PetscInt Nc = uOff[1] - uOff[0];

1303:   for (PetscInt i = 0; i < Nc; ++i) f[i] = u[i];
1304: }

1306: 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[])
1307: {
1308:   for (PetscInt i = 0; i < dim; ++i) {
1309:     for (PetscInt j = 0; j < dim; ++j) f[i + dim * j] = u[i + dim * j];
1310:   }
1311: }

1313: static PetscErrorCode DMAdaptorAdapt_Sequence_Private(DMAdaptor adaptor, Vec inx, PetscBool doSolve, DM *adm, Vec *ax)
1314: {
1315:   PetscDS   ds;
1316:   PetscReal errorNorm = 0.;
1317:   PetscInt  numAdapt  = adaptor->numSeq, adaptIter;
1318:   PetscInt  dim, coordDim, Nf;
1319:   void     *ctx;
1320:   MPI_Comm  comm;

1322:   PetscFunctionBegin;
1323:   PetscCall(DMViewFromOptions(adaptor->idm, NULL, "-dm_adapt_pre_view"));
1324:   PetscCall(VecViewFromOptions(inx, NULL, "-sol_adapt_pre_view"));
1325:   PetscCall(PetscObjectGetComm((PetscObject)adaptor, &comm));
1326:   PetscCall(DMGetDimension(adaptor->idm, &dim));
1327:   PetscCall(DMGetCoordinateDim(adaptor->idm, &coordDim));
1328:   PetscCall(DMGetApplicationContext(adaptor->idm, &ctx));
1329:   PetscCall(DMGetDS(adaptor->idm, &ds));
1330:   PetscCall(PetscDSGetNumFields(ds, &Nf));
1331:   PetscCheck(Nf != 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot refine with no fields present!");

1333:   /* Adapt until nothing changes */
1334:   /* Adapt for a specified number of iterates */
1335:   for (adaptIter = 0; adaptIter < numAdapt - 1; ++adaptIter) PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(comm)));
1336:   for (adaptIter = 0; adaptIter < numAdapt; ++adaptIter) {
1337:     PetscBool adapted = PETSC_FALSE;
1338:     DM        dm      = adaptIter ? *adm : adaptor->idm, odm;
1339:     Vec       x       = adaptIter ? *ax : inx, locX, ox;
1340:     Vec       error   = NULL;

1342:     PetscCall(DMGetLocalVector(dm, &locX));
1343:     PetscCall(DMAdaptorPreAdapt(adaptor, locX));
1344:     if (doSolve) {
1345:       SNES snes;

1347:       PetscCall(DMAdaptorGetSolver(adaptor, &snes));
1348:       PetscCall(SNESSolve(snes, NULL, adaptIter ? *ax : x));
1349:     }
1350:     PetscCall(DMGlobalToLocalBegin(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
1351:     PetscCall(DMGlobalToLocalEnd(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
1352:     PetscCall(VecViewFromOptions(adaptIter ? *ax : x, (PetscObject)adaptor, "-adapt_primal_sol_vec_view"));
1353:     switch (adaptor->adaptCriterion) {
1354:     case DM_ADAPTATION_REFINE:
1355:       PetscCall(DMRefine(dm, comm, &odm));
1356:       PetscCheck(odm, comm, PETSC_ERR_ARG_INCOMP, "DMRefine() did not perform any refinement, cannot continue grid sequencing");
1357:       adapted = PETSC_TRUE;
1358:       PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, NULL));
1359:       break;
1360:     case DM_ADAPTATION_LABEL: {
1361:       /* Adapt DM
1362:            Create local solution
1363:            Reconstruct gradients (FVM) or solve adjoint equation (FEM)
1364:            Produce cellwise error indicator */
1365:       DM             edm, plex;
1366:       PetscDS        ds;
1367:       PetscFE        efe;
1368:       DMLabel        adaptLabel;
1369:       IS             refineIS, coarsenIS;
1370:       DMPolytopeType ct;
1371:       PetscScalar    errorVal;
1372:       PetscInt       nRefine, nCoarsen, cStart;

1374:       PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));

1376:       // TODO Move this creation to PreAdapt
1377:       PetscCall(DMClone(dm, &edm));
1378:       PetscCall(DMConvert(edm, DMPLEX, &plex));
1379:       PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, NULL));
1380:       PetscCall(DMPlexGetCellType(plex, cStart, &ct));
1381:       PetscCall(DMDestroy(&plex));
1382:       PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 0, PETSC_DEFAULT, &efe));
1383:       PetscCall(PetscObjectSetName((PetscObject)efe, "Error"));
1384:       PetscCall(DMSetField(edm, 0, NULL, (PetscObject)efe));
1385:       PetscCall(PetscFEDestroy(&efe));
1386:       PetscCall(DMCreateDS(edm));
1387:       PetscCall(DMGetGlobalVector(edm, &error));
1388:       PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));

1390:       PetscUseTypeMethod(adaptor, computeerrorindicator, locX, error);
1391:       PetscCall(VecViewFromOptions(error, (PetscObject)adaptor, "-adapt_error_vec_view"));
1392:       PetscCall(DMGetDS(edm, &ds));
1393:       PetscCall(PetscDSSetObjective(ds, 0, identity));
1394:       PetscCall(DMPlexComputeIntegralFEM(edm, error, &errorVal, NULL));
1395:       errorNorm = PetscRealPart(errorVal);

1397:       // Compute IS from VecTagger
1398:       PetscCall(VecTaggerComputeIS(adaptor->refineTag, error, &refineIS, NULL));
1399:       PetscCall(VecTaggerComputeIS(adaptor->coarsenTag, error, &coarsenIS, NULL));
1400:       PetscCall(ISViewFromOptions(refineIS, (PetscObject)adaptor->refineTag, "-is_view"));
1401:       PetscCall(ISViewFromOptions(coarsenIS, (PetscObject)adaptor->coarsenTag, "-is_view"));
1402:       PetscCall(ISGetSize(refineIS, &nRefine));
1403:       PetscCall(ISGetSize(coarsenIS, &nCoarsen));
1404:       PetscCall(PetscInfo(adaptor, "DMAdaptor: numRefine %" PetscInt_FMT ", numCoarsen %" PetscInt_FMT "\n", nRefine, nCoarsen));
1405:       if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
1406:       if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
1407:       PetscCall(ISDestroy(&coarsenIS));
1408:       PetscCall(ISDestroy(&refineIS));
1409:       // Adapt DM from label
1410:       if (nRefine || nCoarsen) {
1411:         char        oprefix[PETSC_MAX_PATH_LEN];
1412:         const char *p;
1413:         PetscBool   flg;

1415:         PetscCall(PetscOptionsHasName(NULL, adaptor->hdr.prefix, "-adapt_vec_view", &flg));
1416:         if (flg) {
1417:           Vec ref;

1419:           PetscCall(DMPlexCreateLabelField(dm, adaptLabel, &ref));
1420:           PetscCall(VecViewFromOptions(ref, (PetscObject)adaptor, "-adapt_vec_view"));
1421:           PetscCall(VecDestroy(&ref));
1422:         }

1424:         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &p));
1425:         PetscCall(PetscStrncpy(oprefix, p, PETSC_MAX_PATH_LEN));
1426:         PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)dm, "adapt_"));
1427:         PetscCall(DMAdaptLabel(dm, adaptLabel, &odm));
1428:         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, oprefix));
1429:         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)odm, oprefix));
1430:         PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, error));
1431:         adapted = PETSC_TRUE;
1432:       } else {
1433:         PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, error));
1434:       }
1435:       PetscCall(DMLabelDestroy(&adaptLabel));
1436:       PetscCall(DMRestoreGlobalVector(edm, &error));
1437:       PetscCall(DMDestroy(&edm));
1438:     } break;
1439:     case DM_ADAPTATION_METRIC: {
1440:       DM        dmGrad, dmHess, dmMetric, dmDet;
1441:       Vec       xGrad, xHess, metric, determinant;
1442:       PetscReal N;
1443:       DMLabel   bdLabel = NULL, rgLabel = NULL;
1444:       PetscBool higherOrder = PETSC_FALSE;
1445:       PetscInt  Nd          = coordDim * coordDim, f, vStart, vEnd;
1446:       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[]);

1448:       PetscCall(PetscMalloc(1, &funcs));
1449:       funcs[0] = identityFunc;

1451:       /*     Setup finite element spaces */
1452:       PetscCall(DMClone(dm, &dmGrad));
1453:       PetscCall(DMClone(dm, &dmHess));
1454:       PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple fields not yet considered"); // TODO
1455:       for (f = 0; f < Nf; ++f) {
1456:         PetscFE         fe, feGrad, feHess;
1457:         PetscDualSpace  Q;
1458:         PetscSpace      space;
1459:         DM              K;
1460:         PetscQuadrature q;
1461:         PetscInt        Nc, qorder, p;
1462:         const char     *prefix;

1464:         PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1465:         PetscCall(PetscFEGetNumComponents(fe, &Nc));
1466:         PetscCheck(Nc <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple components not yet considered"); // TODO
1467:         PetscCall(PetscFEGetBasisSpace(fe, &space));
1468:         PetscCall(PetscSpaceGetDegree(space, NULL, &p));
1469:         if (p > 1) higherOrder = PETSC_TRUE;
1470:         PetscCall(PetscFEGetDualSpace(fe, &Q));
1471:         PetscCall(PetscDualSpaceGetDM(Q, &K));
1472:         PetscCall(DMPlexGetDepthStratum(K, 0, &vStart, &vEnd));
1473:         PetscCall(PetscFEGetQuadrature(fe, &q));
1474:         PetscCall(PetscQuadratureGetOrder(q, &qorder));
1475:         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix));
1476:         PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmGrad), dim, Nc * coordDim, PETSC_TRUE, prefix, qorder, &feGrad));
1477:         PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmHess), dim, Nc * Nd, PETSC_TRUE, prefix, qorder, &feHess));
1478:         PetscCall(DMSetField(dmGrad, f, NULL, (PetscObject)feGrad));
1479:         PetscCall(DMSetField(dmHess, f, NULL, (PetscObject)feHess));
1480:         PetscCall(DMCreateDS(dmGrad));
1481:         PetscCall(DMCreateDS(dmHess));
1482:         PetscCall(PetscFEDestroy(&feGrad));
1483:         PetscCall(PetscFEDestroy(&feHess));
1484:       }
1485:       /*     Compute vertexwise gradients from cellwise gradients */
1486:       PetscCall(DMCreateLocalVector(dmGrad, &xGrad));
1487:       PetscCall(VecViewFromOptions(locX, NULL, "-sol_adapt_loc_pre_view"));
1488:       PetscCall(DMPlexComputeGradientClementInterpolant(dm, locX, xGrad));
1489:       PetscCall(VecViewFromOptions(xGrad, NULL, "-adapt_gradient_view"));
1490:       /*     Compute vertexwise Hessians from cellwise Hessians */
1491:       PetscCall(DMCreateLocalVector(dmHess, &xHess));
1492:       PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, xGrad, xHess));
1493:       PetscCall(VecViewFromOptions(xHess, NULL, "-adapt_hessian_view"));
1494:       PetscCall(VecDestroy(&xGrad));
1495:       PetscCall(DMDestroy(&dmGrad));
1496:       /*     Compute L-p normalized metric */
1497:       PetscCall(DMClone(dm, &dmMetric));
1498:       N = adaptor->Nadapt >= 0 ? adaptor->Nadapt : PetscPowRealInt(adaptor->refinementFactor, dim) * ((PetscReal)(vEnd - vStart));
1499:       // TODO This was where the old monitor was, figure out how to show metric and target N
1500:       PetscCall(DMPlexMetricSetTargetComplexity(dmMetric, N));
1501:       if (higherOrder) {
1502:         /*   Project Hessian into P1 space, if required */
1503:         PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1504:         PetscCall(DMProjectFieldLocal(dmMetric, 0.0, xHess, funcs, INSERT_ALL_VALUES, metric));
1505:         PetscCall(VecDestroy(&xHess));
1506:         xHess = metric;
1507:       }
1508:       PetscCall(PetscFree(funcs));
1509:       PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1510:       PetscCall(DMPlexMetricDeterminantCreate(dmMetric, 0, &determinant, &dmDet));
1511:       PetscCall(DMPlexMetricNormalize(dmMetric, xHess, PETSC_TRUE, PETSC_TRUE, metric, determinant));
1512:       PetscCall(VecDestroy(&determinant));
1513:       PetscCall(DMDestroy(&dmDet));
1514:       PetscCall(VecDestroy(&xHess));
1515:       PetscCall(DMDestroy(&dmHess));
1516:       /*     Adapt DM from metric */
1517:       PetscCall(DMGetLabel(dm, "marker", &bdLabel));
1518:       PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &odm));
1519:       PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, NULL));
1520:       adapted = PETSC_TRUE;
1521:       /*     Cleanup */
1522:       PetscCall(VecDestroy(&metric));
1523:       PetscCall(DMDestroy(&dmMetric));
1524:     } break;
1525:     default:
1526:       SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid adaptation type: %d", adaptor->adaptCriterion);
1527:     }
1528:     PetscCall(DMAdaptorPostAdapt(adaptor));
1529:     PetscCall(DMRestoreLocalVector(dm, &locX));
1530:     /* If DM was adapted, replace objects and recreate solution */
1531:     if (adapted) {
1532:       const char *name;

1534:       PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1535:       PetscCall(PetscObjectSetName((PetscObject)odm, name));
1536:       /* Reconfigure solver */
1537:       PetscCall(SNESReset(adaptor->snes));
1538:       PetscCall(SNESSetDM(adaptor->snes, odm));
1539:       PetscCall(DMAdaptorSetSolver(adaptor, adaptor->snes));
1540:       PetscCall(DMPlexSetSNESLocalFEM(odm, PETSC_FALSE, ctx));
1541:       PetscCall(SNESSetFromOptions(adaptor->snes));
1542:       /* Transfer system */
1543:       PetscCall(DMCopyDisc(dm, odm));
1544:       /* Transfer solution */
1545:       PetscCall(DMCreateGlobalVector(odm, &ox));
1546:       PetscCall(PetscObjectGetName((PetscObject)x, &name));
1547:       PetscCall(PetscObjectSetName((PetscObject)ox, name));
1548:       PetscCall(DMAdaptorTransferSolution(adaptor, dm, x, odm, ox));
1549:       /* Cleanup adaptivity info */
1550:       if (adaptIter > 0) PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(comm)));
1551:       PetscCall(DMForestSetAdaptivityForest(dm, NULL)); /* clear internal references to the previous dm */
1552:       PetscCall(DMDestroy(&dm));
1553:       PetscCall(VecDestroy(&x));
1554:       *adm = odm;
1555:       *ax  = ox;
1556:     } else {
1557:       *adm      = dm;
1558:       *ax       = x;
1559:       adaptIter = numAdapt;
1560:     }
1561:     if (adaptIter < numAdapt - 1) {
1562:       PetscCall(DMViewFromOptions(odm, NULL, "-dm_adapt_iter_view"));
1563:       PetscCall(VecViewFromOptions(ox, NULL, "-sol_adapt_iter_view"));
1564:     }
1565:   }
1566:   PetscCall(DMViewFromOptions(*adm, NULL, "-dm_adapt_view"));
1567:   PetscCall(VecViewFromOptions(*ax, NULL, "-sol_adapt_view"));
1568:   PetscFunctionReturn(PETSC_SUCCESS);
1569: }

1571: /*@
1572:   DMAdaptorAdapt - Creates a new `DM` that is adapted to the problem

1574:   Not Collective

1576:   Input Parameters:
1577: + adaptor  - The `DMAdaptor` object
1578: . x        - The global approximate solution
1579: - strategy - The adaptation strategy, see `DMAdaptationStrategy`

1581:   Output Parameters:
1582: + adm - The adapted `DM`
1583: - ax  - The adapted solution

1585:   Options Database Keys:
1586: + -snes_adapt <strategy> - initial, sequential, multigrid
1587: . -adapt_gradient_view   - View the Clement interpolant of the solution gradient
1588: . -adapt_hessian_view    - View the Clement interpolant of the solution Hessian
1589: - -adapt_metric_view     - View the metric tensor for adaptive mesh refinement

1591:   Level: intermediate

1593: .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`
1594: @*/
1595: PetscErrorCode DMAdaptorAdapt(DMAdaptor adaptor, Vec x, DMAdaptationStrategy strategy, DM *adm, Vec *ax)
1596: {
1597:   PetscFunctionBegin;
1598:   switch (strategy) {
1599:   case DM_ADAPTATION_INITIAL:
1600:     PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_FALSE, adm, ax));
1601:     break;
1602:   case DM_ADAPTATION_SEQUENTIAL:
1603:     PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_TRUE, adm, ax));
1604:     break;
1605:   default:
1606:     SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Unrecognized adaptation strategy %d", strategy);
1607:   }
1608:   PetscFunctionReturn(PETSC_SUCCESS);
1609: }

1611: /*@C
1612:   DMAdaptorGetMixedSetupFunction - Get the function setting up the mixed problem, if it exists

1614:   Not Collective

1616:   Input Parameter:
1617: . adaptor - the `DMAdaptor`

1619:   Output Parameter:
1620: . setupFunc - the function setting up the mixed problem, or `NULL`

1622:   Level: advanced

1624: .seealso: `DMAdaptor`, `DMAdaptorSetMixedSetupFunction()`, `DMAdaptorAdapt()`
1625: @*/
1626: PetscErrorCode DMAdaptorGetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (**setupFunc)(DMAdaptor, DM))
1627: {
1628:   PetscFunctionBegin;
1630:   PetscAssertPointer(setupFunc, 2);
1631:   *setupFunc = adaptor->ops->mixedsetup;
1632:   PetscFunctionReturn(PETSC_SUCCESS);
1633: }

1635: /*@C
1636:   DMAdaptorSetMixedSetupFunction - Set the function setting up the mixed problem

1638:   Not Collective

1640:   Input Parameters:
1641: + adaptor   - the `DMAdaptor`
1642: - setupFunc - the function setting up the mixed problem

1644:   Level: advanced

1646: .seealso: `DMAdaptor`, `DMAdaptorGetMixedSetupFunction()`, `DMAdaptorAdapt()`
1647: @*/
1648: PetscErrorCode DMAdaptorSetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (*setupFunc)(DMAdaptor, DM))
1649: {
1650:   PetscFunctionBegin;
1653:   adaptor->ops->mixedsetup = setupFunc;
1654:   PetscFunctionReturn(PETSC_SUCCESS);
1655: }

1657: /*@C
1658:   DMAdaptorGetCriterion - Get the adaptation criterion

1660:   Not Collective

1662:   Input Parameter:
1663: . adaptor - the `DMAdaptor`

1665:   Output Parameter:
1666: . criterion - the criterion for adaptation

1668:   Level: advanced

1670: .seealso: `DMAdaptor`, `DMAdaptorSetCrierion()`, `DMAdaptationCriterion`
1671: @*/
1672: PetscErrorCode DMAdaptorGetCriterion(DMAdaptor adaptor, DMAdaptationCriterion *criterion)
1673: {
1674:   PetscFunctionBegin;
1676:   PetscAssertPointer(criterion, 2);
1677:   *criterion = adaptor->adaptCriterion;
1678:   PetscFunctionReturn(PETSC_SUCCESS);
1679: }

1681: /*@C
1682:   DMAdaptorSetCriterion - Set the adaptation criterion

1684:   Not Collective

1686:   Input Parameters:
1687: + adaptor   - the `DMAdaptor`
1688: - criterion - the adaptation criterion

1690:   Level: advanced

1692: .seealso: `DMAdaptor`, `DMAdaptorGetCriterion()`, `DMAdaptationCriterion`
1693: @*/
1694: PetscErrorCode DMAdaptorSetCriterion(DMAdaptor adaptor, DMAdaptationCriterion criterion)
1695: {
1696:   PetscFunctionBegin;
1698:   adaptor->adaptCriterion = criterion;
1699:   PetscFunctionReturn(PETSC_SUCCESS);
1700: }

1702: static PetscErrorCode DMAdaptorInitialize_Gradient(DMAdaptor adaptor)
1703: {
1704:   PetscFunctionBegin;
1705:   adaptor->ops->computeerrorindicator     = DMAdaptorComputeErrorIndicator_Gradient;
1706:   adaptor->ops->computecellerrorindicator = DMAdaptorComputeCellErrorIndicator_Gradient;
1707:   PetscFunctionReturn(PETSC_SUCCESS);
1708: }

1710: PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Gradient(DMAdaptor adaptor)
1711: {
1712:   PetscFunctionBegin;
1714:   adaptor->data = NULL;

1716:   PetscCall(DMAdaptorInitialize_Gradient(adaptor));
1717:   PetscFunctionReturn(PETSC_SUCCESS);
1718: }

1720: static PetscErrorCode DMAdaptorInitialize_Flux(DMAdaptor adaptor)
1721: {
1722:   PetscFunctionBegin;
1723:   adaptor->ops->computeerrorindicator = DMAdaptorComputeErrorIndicator_Flux;
1724:   PetscFunctionReturn(PETSC_SUCCESS);
1725: }

1727: PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Flux(DMAdaptor adaptor)
1728: {
1729:   PetscFunctionBegin;
1731:   adaptor->data = NULL;

1733:   PetscCall(DMAdaptorInitialize_Flux(adaptor));
1734:   PetscFunctionReturn(PETSC_SUCCESS);
1735: }