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, PetscCtx 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, PetscCtx ctx), PetscCtx ctx, PetscCtxDestroyFn *monitordestroy)
352: {
353: PetscFunctionBegin;
355: for (PetscInt i = 0; i < adaptor->numbermonitors; i++) {
356: PetscBool identical;
358: PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)monitor, ctx, monitordestroy, (PetscErrorCode (*)(void))(PetscVoidFn *)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[], PetscCtx 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, PetscCtx 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;
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 (PetscInt 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 adaptor, DM dm, Vec xin, DM newdm, Vec xout, PetscCtx ctx))
680: {
681: PetscFunctionBegin;
682: *tfunc = adaptor->ops->transfersolution;
683: PetscFunctionReturn(PETSC_SUCCESS);
684: }
686: PetscErrorCode DMAdaptorSetTransferFunction(DMAdaptor adaptor, PetscErrorCode (*tfunc)(DMAdaptor adaptor, DM dm, Vec xin, DM newdm, Vec xout, PetscCtx ctx))
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) adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
709: #if defined(PETSC_HAVE_PRAGMATIC)
710: else {
711: adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
712: }
713: #elif defined(PETSC_HAVE_MMG)
714: else {
715: adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
716: }
717: #elif defined(PETSC_HAVE_PARMMG)
718: else {
719: adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
720: }
721: #else
722: else {
723: adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
724: }
725: #endif
726: }
727: if (id == PETSCFV_CLASSID) {
728: adaptor->femType = PETSC_FALSE;
729: } else {
730: adaptor->femType = PETSC_TRUE;
731: }
732: if (adaptor->femType) {
733: /* Compute local solution bc */
734: PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
735: } else {
736: PetscFV fvm = (PetscFV)obj;
737: PetscLimiter noneLimiter;
738: Vec grad;
740: PetscCall(PetscFVGetComputeGradients(fvm, &adaptor->computeGradient));
741: PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
742: /* Use no limiting when reconstructing gradients for adaptivity */
743: PetscCall(PetscFVGetLimiter(fvm, &adaptor->limiter));
744: PetscCall(PetscObjectReference((PetscObject)adaptor->limiter));
745: PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
746: PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
747: PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
748: /* Get FVM data */
749: PetscCall(DMPlexGetDataFVM(plex, fvm, &adaptor->cellGeom, &adaptor->faceGeom, &adaptor->gradDM));
750: PetscCall(VecGetDM(adaptor->cellGeom, &adaptor->cellDM));
751: PetscCall(VecGetArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
752: /* Compute local solution bc */
753: PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
754: /* Compute gradients */
755: PetscCall(DMCreateGlobalVector(adaptor->gradDM, &grad));
756: PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
757: PetscCall(DMGetLocalVector(adaptor->gradDM, &adaptor->cellGrad));
758: PetscCall(DMGlobalToLocalBegin(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
759: PetscCall(DMGlobalToLocalEnd(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
760: PetscCall(VecDestroy(&grad));
761: PetscCall(VecGetArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
762: }
763: PetscCall(DMDestroy(&plex));
764: PetscFunctionReturn(PETSC_SUCCESS);
765: }
767: static PetscErrorCode DMAdaptorTransferSolution(DMAdaptor adaptor, DM dm, Vec x, DM adm, Vec ax)
768: {
769: PetscReal time = 0.0;
770: Mat interp;
771: void *ctx;
773: PetscFunctionBegin;
774: PetscCall(DMGetApplicationContext(dm, &ctx));
775: if (adaptor->ops->transfersolution) PetscUseTypeMethod(adaptor, transfersolution, dm, x, adm, ax, ctx);
776: else {
777: switch (adaptor->adaptCriterion) {
778: case DM_ADAPTATION_LABEL:
779: PetscCall(DMForestTransferVec(dm, x, adm, ax, PETSC_TRUE, time));
780: break;
781: case DM_ADAPTATION_REFINE:
782: case DM_ADAPTATION_METRIC:
783: PetscCall(DMCreateInterpolation(dm, adm, &interp, NULL));
784: PetscCall(MatInterpolate(interp, x, ax));
785: PetscCall(DMInterpolate(dm, interp, adm));
786: PetscCall(MatDestroy(&interp));
787: break;
788: default:
789: SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_SUP, "No built-in projection for this adaptation criterion: %d", adaptor->adaptCriterion);
790: }
791: }
792: PetscFunctionReturn(PETSC_SUCCESS);
793: }
795: static PetscErrorCode DMAdaptorPostAdapt(DMAdaptor adaptor)
796: {
797: PetscDS prob;
798: PetscObject obj;
799: PetscClassId id;
801: PetscFunctionBegin;
802: PetscCall(DMGetDS(adaptor->idm, &prob));
803: PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
804: PetscCall(PetscObjectGetClassId(obj, &id));
805: if (id == PETSCFV_CLASSID) {
806: PetscFV fvm = (PetscFV)obj;
808: PetscCall(PetscFVSetComputeGradients(fvm, adaptor->computeGradient));
809: /* Restore original limiter */
810: PetscCall(PetscFVSetLimiter(fvm, adaptor->limiter));
812: PetscCall(VecRestoreArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
813: PetscCall(VecRestoreArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
814: PetscCall(DMRestoreLocalVector(adaptor->gradDM, &adaptor->cellGrad));
815: }
816: PetscFunctionReturn(PETSC_SUCCESS);
817: }
819: /*
820: DMAdaptorComputeCellErrorIndicator_Gradient - Use the integrated gradient as an error indicator in the `DMAdaptor`
822: Input Parameters:
823: + adaptor - The `DMAdaptor` object
824: . dim - The topological dimension
825: . cell - The cell
826: . field - The field integrated over the cell
827: . gradient - The gradient integrated over the cell
828: . cg - A `PetscFVCellGeom` struct
829: - ctx - A user context
831: Output Parameter:
832: . errInd - The error indicator
834: Developer Note:
835: Some of the input arguments are absurdly specialized to special situations, it is not clear this is a good general API
837: .seealso: [](ch_dmbase), `DMAdaptor`
838: */
839: static PetscErrorCode DMAdaptorComputeCellErrorIndicator_Gradient(DMAdaptor adaptor, PetscInt dim, PetscInt Nc, const PetscScalar *field, const PetscScalar *gradient, const PetscFVCellGeom *cg, PetscReal *errInd, PetscCtx ctx)
840: {
841: PetscReal err = 0.;
843: PetscFunctionBeginHot;
844: for (PetscInt c = 0; c < Nc; c++) {
845: for (PetscInt d = 0; d < dim; ++d) err += PetscSqr(PetscRealPart(gradient[c * dim + d]));
846: }
847: *errInd = cg->volume * err;
848: PetscFunctionReturn(PETSC_SUCCESS);
849: }
851: static PetscErrorCode DMAdaptorComputeErrorIndicator_Gradient(DMAdaptor adaptor, Vec locX, Vec errVec)
852: {
853: DM dm, plex, edm, eplex;
854: PetscDS ds;
855: PetscObject obj;
856: PetscClassId id;
857: void *ctx;
858: PetscQuadrature quad;
859: PetscScalar *earray;
860: PetscReal minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL};
861: PetscInt dim, cdim, cStart, cEnd, Nf, Nc;
863: PetscFunctionBegin;
864: PetscCall(VecGetDM(locX, &dm));
865: PetscCall(DMConvert(dm, DMPLEX, &plex));
866: PetscCall(VecGetDM(errVec, &edm));
867: PetscCall(DMConvert(edm, DMPLEX, &eplex));
868: PetscCall(DMGetDimension(plex, &dim));
869: PetscCall(DMGetCoordinateDim(plex, &cdim));
870: PetscCall(DMGetApplicationContext(plex, &ctx));
871: PetscCall(DMGetDS(plex, &ds));
872: PetscCall(PetscDSGetNumFields(ds, &Nf));
873: PetscCall(PetscDSGetDiscretization(ds, 0, &obj));
874: PetscCall(PetscObjectGetClassId(obj, &id));
876: PetscCall(VecGetArray(errVec, &earray));
877: PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
878: for (PetscInt cell = cStart; cell < cEnd; ++cell) {
879: PetscScalar *eval;
880: PetscReal errInd = 0.;
882: if (id == PETSCFV_CLASSID) {
883: PetscFV fv = (PetscFV)obj;
884: const PetscScalar *pointSols;
885: const PetscScalar *pointSol;
886: const PetscScalar *pointGrad;
887: PetscFVCellGeom *cg;
889: PetscCall(PetscFVGetNumComponents(fv, &Nc));
890: PetscCall(VecGetArrayRead(locX, &pointSols));
891: PetscCall(DMPlexPointLocalRead(plex, cell, pointSols, (void *)&pointSol));
892: PetscCall(DMPlexPointLocalRead(adaptor->gradDM, cell, adaptor->cellGradArray, (void *)&pointGrad));
893: PetscCall(DMPlexPointLocalRead(adaptor->cellDM, cell, adaptor->cellGeomArray, &cg));
894: PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, pointSol, pointGrad, cg, &errInd, ctx);
895: PetscCall(VecRestoreArrayRead(locX, &pointSols));
896: } else {
897: PetscFE fe = (PetscFE)obj;
898: PetscScalar *x = NULL, *field, *gradient, *interpolant, *interpolantGrad;
899: PetscFVCellGeom cg;
900: PetscFEGeom fegeom;
901: const PetscReal *quadWeights;
902: PetscReal *coords;
903: PetscInt Nb, Nq, qNc;
905: fegeom.dim = dim;
906: fegeom.dimEmbed = cdim;
907: PetscCall(PetscFEGetNumComponents(fe, &Nc));
908: PetscCall(PetscFEGetQuadrature(fe, &quad));
909: PetscCall(PetscFEGetDimension(fe, &Nb));
910: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
911: PetscCall(PetscCalloc6(Nc, &field, cdim * Nc, &gradient, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ));
912: PetscCall(PetscMalloc2(Nc, &interpolant, cdim * Nc, &interpolantGrad));
913: PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
914: PetscCall(DMPlexComputeCellGeometryFVM(plex, cell, &cg.volume, NULL, NULL));
915: PetscCall(PetscArrayzero(gradient, cdim * Nc));
916: PetscCall(DMPlexVecGetClosure(plex, NULL, locX, cell, NULL, &x));
917: for (PetscInt f = 0; f < Nf; ++f) {
918: PetscInt qc = 0;
920: PetscCall(PetscDSGetDiscretization(ds, f, &obj));
921: PetscCall(PetscArrayzero(interpolant, Nc));
922: PetscCall(PetscArrayzero(interpolantGrad, cdim * Nc));
923: for (PetscInt q = 0; q < Nq; ++q) {
924: PetscCall(PetscFEInterpolateFieldAndGradient_Static((PetscFE)obj, 1, x, &fegeom, q, interpolant, interpolantGrad));
925: for (PetscInt fc = 0; fc < Nc; ++fc) {
926: const PetscReal wt = quadWeights[q * qNc + qc + fc];
928: field[fc] += interpolant[fc] * wt * fegeom.detJ[q];
929: for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] += interpolantGrad[fc * dim + d] * wt * fegeom.detJ[q];
930: }
931: }
932: qc += Nc;
933: }
934: PetscCall(PetscFree2(interpolant, interpolantGrad));
935: PetscCall(DMPlexVecRestoreClosure(plex, NULL, locX, cell, NULL, &x));
936: for (PetscInt fc = 0; fc < Nc; ++fc) {
937: field[fc] /= cg.volume;
938: for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] /= cg.volume;
939: }
940: PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, field, gradient, &cg, &errInd, ctx);
941: PetscCall(PetscFree6(field, gradient, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
942: }
943: PetscCall(DMPlexPointGlobalRef(eplex, cell, earray, (void *)&eval));
944: eval[0] = errInd;
945: minMaxInd[0] = PetscMin(minMaxInd[0], errInd);
946: minMaxInd[1] = PetscMax(minMaxInd[1], errInd);
947: }
948: PetscCall(VecRestoreArray(errVec, &earray));
949: PetscCall(DMDestroy(&plex));
950: PetscCall(DMDestroy(&eplex));
951: PetscCall(PetscGlobalMinMaxReal(PetscObjectComm((PetscObject)adaptor), minMaxInd, minMaxInd));
952: PetscCall(PetscInfo(adaptor, "DMAdaptor: error indicator range (%g, %g)\n", (double)minMaxInd[0], (double)minMaxInd[1]));
953: PetscFunctionReturn(PETSC_SUCCESS);
954: }
956: static PetscErrorCode DMAdaptorComputeErrorIndicator_Flux(DMAdaptor adaptor, Vec lu, Vec errVec)
957: {
958: DM dm, mdm;
959: SNES msnes;
960: Vec mu, lmu;
961: void *ctx;
962: const char *prefix;
964: PetscFunctionBegin;
965: PetscCall(VecGetDM(lu, &dm));
967: // Set up and solve mixed problem
968: PetscCall(DMClone(dm, &mdm));
969: PetscCall(SNESCreate(PetscObjectComm((PetscObject)mdm), &msnes));
970: PetscCall(SNESSetDM(msnes, mdm));
971: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
972: PetscCall(SNESSetOptionsPrefix(msnes, prefix));
973: PetscCall(SNESAppendOptionsPrefix(msnes, "mixed_"));
975: PetscTryTypeMethod(adaptor, mixedsetup, mdm);
976: PetscCall(DMGetApplicationContext(dm, &ctx));
977: PetscCall(DMPlexSetSNESLocalFEM(mdm, PETSC_FALSE, ctx));
978: PetscCall(SNESSetFromOptions(msnes));
980: PetscCall(DMCreateGlobalVector(mdm, &mu));
981: PetscCall(PetscObjectSetName((PetscObject)mu, "Mixed Solution"));
982: PetscCall(SNESSolve(msnes, NULL, mu));
983: PetscCall(VecViewFromOptions(mu, (PetscObject)adaptor, "-adapt_mixed_sol_vec_view"));
985: PetscCall(DMGetLocalVector(mdm, &lmu));
986: PetscCall(DMGlobalToLocal(mdm, mu, INSERT_VALUES, lmu));
987: PetscCall(DMPlexInsertBoundaryValues(mdm, PETSC_TRUE, lmu, 0.0, NULL, NULL, NULL));
988: PetscCall(DMPlexComputeL2FluxDiffVecLocal(lu, 0, lmu, 0, errVec));
989: PetscCall(DMRestoreLocalVector(mdm, &lmu));
990: PetscCall(VecDestroy(&mu));
991: PetscCall(SNESDestroy(&msnes));
992: PetscCall(DMDestroy(&mdm));
993: PetscFunctionReturn(PETSC_SUCCESS);
994: }
996: /*@
997: DMAdaptorMonitor - runs the user provided monitor routines, if they exist
999: Collective
1001: Input Parameters:
1002: + adaptor - the `DMAdaptor`
1003: . it - iteration number
1004: . odm - the original `DM`
1005: . adm - the adapted `DM`
1006: . Nf - the number of fields
1007: . enorms - the 2-norm error values for each field
1008: - error - `Vec` of cellwise errors
1010: Level: developer
1012: Note:
1013: This routine is called by the `DMAdaptor` implementations.
1014: It does not typically need to be called by the user.
1016: .seealso: [](ch_snes), `DMAdaptorMonitorSet()`
1017: @*/
1018: PetscErrorCode DMAdaptorMonitor(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error)
1019: {
1020: PetscFunctionBegin;
1021: for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) PetscCall((*adaptor->monitor[i])(adaptor, it, odm, adm, Nf, enorms, error, adaptor->monitorcontext[i]));
1022: PetscFunctionReturn(PETSC_SUCCESS);
1023: }
1025: /*@C
1026: DMAdaptorMonitorSize - Prints the mesh sizes at each iteration of an adaptation loop.
1028: Collective
1030: Input Parameters:
1031: + adaptor - the `DMAdaptor`
1032: . n - iteration number
1033: . odm - the original `DM`
1034: . adm - the adapted `DM`
1035: . Nf - number of fields
1036: . enorms - 2-norm error values for each field (may be estimated).
1037: . error - `Vec` of cellwise errors
1038: - vf - The viewer context
1040: Options Database Key:
1041: . -adaptor_monitor_size - Activates `DMAdaptorMonitorSize()`
1043: Level: intermediate
1045: Note:
1046: This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1047: to be used during the adaptation loop.
1049: .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorError()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
1050: @*/
1051: PetscErrorCode DMAdaptorMonitorSize(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1052: {
1053: PetscViewer viewer = vf->viewer;
1054: PetscViewerFormat format = vf->format;
1055: PetscInt tablevel, cStart, cEnd, acStart, acEnd;
1056: const char *prefix;
1057: PetscMPIInt rank;
1059: PetscFunctionBegin;
1061: PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
1062: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
1063: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)adaptor), &rank));
1064: PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
1065: PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
1067: PetscCall(PetscViewerPushFormat(viewer, format));
1068: PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1069: if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Sizes for %s adaptation.\n", prefix));
1070: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor rank %d N_orig: %" PetscInt_FMT " N_adapt: %" PetscInt_FMT "\n", n, rank, cEnd - cStart, acEnd - acStart));
1071: PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1072: PetscCall(PetscViewerPopFormat(viewer));
1073: PetscFunctionReturn(PETSC_SUCCESS);
1074: }
1076: /*@C
1077: DMAdaptorMonitorError - Prints the error norm at each iteration of an adaptation loop.
1079: Collective
1081: Input Parameters:
1082: + adaptor - the `DMAdaptor`
1083: . n - iteration number
1084: . odm - the original `DM`
1085: . adm - the adapted `DM`
1086: . Nf - number of fields
1087: . enorms - 2-norm error values for each field (may be estimated).
1088: . error - `Vec` of cellwise errors
1089: - vf - The viewer context
1091: Options Database Key:
1092: . -adaptor_monitor_error - Activates `DMAdaptorMonitorError()`
1094: Level: intermediate
1096: Note:
1097: This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1098: to be used during the adaptation loop.
1100: .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
1101: @*/
1102: PetscErrorCode DMAdaptorMonitorError(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1103: {
1104: PetscViewer viewer = vf->viewer;
1105: PetscViewerFormat format = vf->format;
1106: PetscInt tablevel, cStart, cEnd, acStart, acEnd;
1107: const char *prefix;
1109: PetscFunctionBegin;
1111: PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
1112: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
1114: PetscCall(PetscViewerPushFormat(viewer, format));
1115: PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
1116: if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Error norms for %s adaptation.\n", prefix));
1117: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor Error norm %s", n, Nf > 1 ? "[" : ""));
1118: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1119: for (PetscInt f = 0; f < Nf; ++f) {
1120: if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1121: PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)enorms[f]));
1122: }
1123: PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
1124: PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
1125: PetscCall(PetscViewerASCIIPrintf(viewer, " N: %" PetscInt_FMT " Nadapt: %" PetscInt_FMT "\n", cEnd - cStart, acEnd - acStart));
1126: PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1127: PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
1128: PetscCall(PetscViewerPopFormat(viewer));
1129: PetscFunctionReturn(PETSC_SUCCESS);
1130: }
1132: /*@C
1133: DMAdaptorMonitorErrorDraw - Plots the error at each iteration of an iterative solver.
1135: Collective
1137: Input Parameters:
1138: + adaptor - the `DMAdaptor`
1139: . n - iteration number
1140: . odm - the original `DM`
1141: . adm - the adapted `DM`
1142: . Nf - number of fields
1143: . enorms - 2-norm error values for each field (may be estimated).
1144: . error - `Vec` of cellwise errors
1145: - vf - The viewer context
1147: Options Database Key:
1148: . -adaptor_monitor_error draw - Activates `DMAdaptorMonitorErrorDraw()`
1150: Level: intermediate
1152: Note:
1153: This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1154: to be used during the adaptation loop.
1156: .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
1157: @*/
1158: PetscErrorCode DMAdaptorMonitorErrorDraw(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1159: {
1160: PetscViewer viewer = vf->viewer;
1161: PetscViewerFormat format = vf->format;
1163: PetscFunctionBegin;
1165: PetscCall(PetscViewerPushFormat(viewer, format));
1166: PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
1167: PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", (PetscObject)adaptor));
1168: PetscCall(VecView(error, viewer));
1169: PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", NULL));
1170: PetscCall(PetscViewerPopFormat(viewer));
1171: PetscFunctionReturn(PETSC_SUCCESS);
1172: }
1174: /*@C
1175: DMAdaptorMonitorErrorDrawLGCreate - Creates the context for the error plotter `DMAdaptorMonitorErrorDrawLG()`
1177: Collective
1179: Input Parameters:
1180: + viewer - The `PetscViewer`
1181: . format - The viewer format
1182: - ctx - An optional user context
1184: Output Parameter:
1185: . vf - The viewer context
1187: Level: intermediate
1189: .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `PetscViewerMonitorGLSetUp()`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
1190: @*/
1191: PetscErrorCode DMAdaptorMonitorErrorDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
1192: {
1193: DMAdaptor adaptor = (DMAdaptor)ctx;
1194: char **names;
1195: PetscInt Nf;
1197: PetscFunctionBegin;
1198: PetscCall(DMGetNumFields(adaptor->idm, &Nf));
1199: PetscCall(PetscMalloc1(Nf + 1, &names));
1200: for (PetscInt f = 0; f < Nf; ++f) {
1201: PetscObject disc;
1202: const char *fname;
1203: char lname[PETSC_MAX_PATH_LEN];
1205: PetscCall(DMGetField(adaptor->idm, f, NULL, &disc));
1206: PetscCall(PetscObjectGetName(disc, &fname));
1207: PetscCall(PetscStrncpy(lname, fname, PETSC_MAX_PATH_LEN));
1208: PetscCall(PetscStrlcat(lname, " Error", PETSC_MAX_PATH_LEN));
1209: PetscCall(PetscStrallocpy(lname, &names[f]));
1210: }
1211: PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
1212: (*vf)->data = ctx;
1213: PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Error Norm", Nf, (const char **)names, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
1214: for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(names[f]));
1215: PetscCall(PetscFree(names));
1216: PetscFunctionReturn(PETSC_SUCCESS);
1217: }
1219: /*@C
1220: DMAdaptorMonitorErrorDrawLG - Plots the error norm at each iteration of an adaptive loop.
1222: Collective
1224: Input Parameters:
1225: + adaptor - the `DMAdaptor`
1226: . n - iteration number
1227: . odm - the original `DM`
1228: . adm - the adapted `DM`
1229: . Nf - number of fields
1230: . enorms - 2-norm error values for each field (may be estimated).
1231: . error - `Vec` of cellwise errors
1232: - vf - The viewer context, obtained via `DMAdaptorMonitorErrorDrawLGCreate()`
1234: Options Database Key:
1235: . -adaptor_error draw::draw_lg - Activates `DMAdaptorMonitorErrorDrawLG()`
1237: Level: intermediate
1239: Notes:
1240: This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
1241: to be used during the adaptation loop.
1243: Call `DMAdaptorMonitorErrorDrawLGCreate()` to create the context needed for this monitor
1245: .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorError()`,
1246: `DMAdaptorMonitorTrueResidualDrawLGCreate()`
1247: @*/
1248: PetscErrorCode DMAdaptorMonitorErrorDrawLG(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
1249: {
1250: PetscViewer viewer = vf->viewer;
1251: PetscViewerFormat format = vf->format;
1252: PetscDrawLG lg;
1253: PetscReal *x, *e;
1255: PetscFunctionBegin;
1257: PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
1258: PetscCall(PetscCalloc2(Nf, &x, Nf, &e));
1259: PetscCall(PetscViewerPushFormat(viewer, format));
1260: if (!n) PetscCall(PetscDrawLGReset(lg));
1261: for (PetscInt f = 0; f < Nf; ++f) {
1262: x[f] = (PetscReal)n;
1263: e[f] = enorms[f] > 0.0 ? PetscLog10Real(enorms[f]) : -15.;
1264: }
1265: PetscCall(PetscDrawLGAddPoint(lg, x, e));
1266: PetscCall(PetscDrawLGDraw(lg));
1267: PetscCall(PetscDrawLGSave(lg));
1268: PetscCall(PetscViewerPopFormat(viewer));
1269: PetscCall(PetscFree2(x, e));
1270: PetscFunctionReturn(PETSC_SUCCESS);
1271: }
1273: /*@C
1274: DMAdaptorMonitorRegisterAll - Registers all of the mesh adaptation monitors in the `SNES` package.
1276: Not Collective
1278: Level: advanced
1280: .seealso: [](ch_snes), `SNES`, `DM`, `DMAdaptorMonitorRegister()`, `DMAdaptorRegister()`
1281: @*/
1282: PetscErrorCode DMAdaptorMonitorRegisterAll(void)
1283: {
1284: PetscFunctionBegin;
1285: if (DMAdaptorMonitorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1286: DMAdaptorMonitorRegisterAllCalled = PETSC_TRUE;
1288: PetscCall(DMAdaptorMonitorRegister("size", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorSize, NULL, NULL));
1289: PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorError, NULL, NULL));
1290: PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorErrorDraw, NULL, NULL));
1291: PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DRAW_LG, DMAdaptorMonitorErrorDrawLG, DMAdaptorMonitorErrorDrawLGCreate, NULL));
1292: PetscFunctionReturn(PETSC_SUCCESS);
1293: }
1295: 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[])
1296: {
1297: const PetscInt Nc = uOff[1] - uOff[0];
1299: for (PetscInt i = 0; i < Nc; ++i) f[i] = u[i];
1300: }
1302: 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[])
1303: {
1304: for (PetscInt i = 0; i < dim; ++i) {
1305: for (PetscInt j = 0; j < dim; ++j) f[i + dim * j] = u[i + dim * j];
1306: }
1307: }
1309: static PetscErrorCode DMAdaptorAdapt_Sequence_Private(DMAdaptor adaptor, Vec inx, PetscBool doSolve, DM *adm, Vec *ax)
1310: {
1311: PetscDS ds;
1312: PetscReal errorNorm = 0.;
1313: PetscInt numAdapt = adaptor->numSeq, adaptIter;
1314: PetscInt dim, coordDim, Nf;
1315: void *ctx;
1316: MPI_Comm comm;
1318: PetscFunctionBegin;
1319: PetscCall(DMViewFromOptions(adaptor->idm, NULL, "-dm_adapt_pre_view"));
1320: PetscCall(VecViewFromOptions(inx, NULL, "-sol_adapt_pre_view"));
1321: PetscCall(PetscObjectGetComm((PetscObject)adaptor, &comm));
1322: PetscCall(DMGetDimension(adaptor->idm, &dim));
1323: PetscCall(DMGetCoordinateDim(adaptor->idm, &coordDim));
1324: PetscCall(DMGetApplicationContext(adaptor->idm, &ctx));
1325: PetscCall(DMGetDS(adaptor->idm, &ds));
1326: PetscCall(PetscDSGetNumFields(ds, &Nf));
1327: PetscCheck(Nf != 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot refine with no fields present!");
1329: /* Adapt until nothing changes */
1330: /* Adapt for a specified number of iterates */
1331: for (adaptIter = 0; adaptIter < numAdapt - 1; ++adaptIter) PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(comm)));
1332: for (adaptIter = 0; adaptIter < numAdapt; ++adaptIter) {
1333: PetscBool adapted = PETSC_FALSE;
1334: DM dm = adaptIter ? *adm : adaptor->idm, odm;
1335: Vec x = adaptIter ? *ax : inx, locX, ox;
1336: Vec error = NULL;
1338: PetscCall(DMGetLocalVector(dm, &locX));
1339: PetscCall(DMAdaptorPreAdapt(adaptor, locX));
1340: if (doSolve) {
1341: SNES snes;
1343: PetscCall(DMAdaptorGetSolver(adaptor, &snes));
1344: PetscCall(SNESSolve(snes, NULL, adaptIter ? *ax : x));
1345: }
1346: PetscCall(DMGlobalToLocalBegin(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
1347: PetscCall(DMGlobalToLocalEnd(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
1348: PetscCall(VecViewFromOptions(adaptIter ? *ax : x, (PetscObject)adaptor, "-adapt_primal_sol_vec_view"));
1349: switch (adaptor->adaptCriterion) {
1350: case DM_ADAPTATION_REFINE:
1351: PetscCall(DMRefine(dm, comm, &odm));
1352: PetscCheck(odm, comm, PETSC_ERR_ARG_INCOMP, "DMRefine() did not perform any refinement, cannot continue grid sequencing");
1353: adapted = PETSC_TRUE;
1354: PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, NULL));
1355: break;
1356: case DM_ADAPTATION_LABEL: {
1357: /* Adapt DM
1358: Create local solution
1359: Reconstruct gradients (FVM) or solve adjoint equation (FEM)
1360: Produce cellwise error indicator */
1361: DM edm, plex;
1362: PetscDS ds;
1363: PetscFE efe;
1364: DMLabel adaptLabel;
1365: IS refineIS, coarsenIS;
1366: DMPolytopeType ct;
1367: PetscScalar errorVal;
1368: PetscInt nRefine, nCoarsen, cStart;
1370: PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
1372: // TODO Move this creation to PreAdapt
1373: PetscCall(DMClone(dm, &edm));
1374: PetscCall(DMConvert(edm, DMPLEX, &plex));
1375: PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, NULL));
1376: PetscCall(DMPlexGetCellType(plex, cStart, &ct));
1377: PetscCall(DMDestroy(&plex));
1378: PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 0, PETSC_DEFAULT, &efe));
1379: PetscCall(PetscObjectSetName((PetscObject)efe, "Error"));
1380: PetscCall(DMSetField(edm, 0, NULL, (PetscObject)efe));
1381: PetscCall(PetscFEDestroy(&efe));
1382: PetscCall(DMCreateDS(edm));
1383: PetscCall(DMGetGlobalVector(edm, &error));
1384: PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
1386: PetscUseTypeMethod(adaptor, computeerrorindicator, locX, error);
1387: PetscCall(VecViewFromOptions(error, (PetscObject)adaptor, "-adapt_error_vec_view"));
1388: PetscCall(DMGetDS(edm, &ds));
1389: PetscCall(PetscDSSetObjective(ds, 0, identity));
1390: PetscCall(DMPlexComputeIntegralFEM(edm, error, &errorVal, NULL));
1391: errorNorm = PetscRealPart(errorVal);
1393: // Compute IS from VecTagger
1394: PetscCall(VecTaggerComputeIS(adaptor->refineTag, error, &refineIS, NULL));
1395: PetscCall(VecTaggerComputeIS(adaptor->coarsenTag, error, &coarsenIS, NULL));
1396: PetscCall(ISViewFromOptions(refineIS, (PetscObject)adaptor->refineTag, "-is_view"));
1397: PetscCall(ISViewFromOptions(coarsenIS, (PetscObject)adaptor->coarsenTag, "-is_view"));
1398: PetscCall(ISGetSize(refineIS, &nRefine));
1399: PetscCall(ISGetSize(coarsenIS, &nCoarsen));
1400: PetscCall(PetscInfo(adaptor, "DMAdaptor: numRefine %" PetscInt_FMT ", numCoarsen %" PetscInt_FMT "\n", nRefine, nCoarsen));
1401: if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
1402: if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
1403: PetscCall(ISDestroy(&coarsenIS));
1404: PetscCall(ISDestroy(&refineIS));
1405: // Adapt DM from label
1406: if (nRefine || nCoarsen) {
1407: char oprefix[PETSC_MAX_PATH_LEN];
1408: const char *p;
1409: PetscBool flg;
1411: PetscCall(PetscOptionsHasName(NULL, adaptor->hdr.prefix, "-adapt_vec_view", &flg));
1412: if (flg) {
1413: Vec ref;
1415: PetscCall(DMPlexCreateLabelField(dm, adaptLabel, &ref));
1416: PetscCall(VecViewFromOptions(ref, (PetscObject)adaptor, "-adapt_vec_view"));
1417: PetscCall(VecDestroy(&ref));
1418: }
1420: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &p));
1421: PetscCall(PetscStrncpy(oprefix, p, PETSC_MAX_PATH_LEN));
1422: PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)dm, "adapt_"));
1423: PetscCall(DMAdaptLabel(dm, adaptLabel, &odm));
1424: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, oprefix));
1425: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)odm, oprefix));
1426: PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, error));
1427: adapted = PETSC_TRUE;
1428: } else {
1429: PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, error));
1430: }
1431: PetscCall(DMLabelDestroy(&adaptLabel));
1432: PetscCall(DMRestoreGlobalVector(edm, &error));
1433: PetscCall(DMDestroy(&edm));
1434: } break;
1435: case DM_ADAPTATION_METRIC: {
1436: DM dmGrad, dmHess, dmMetric, dmDet;
1437: Vec xGrad, xHess, metric, determinant;
1438: PetscReal N;
1439: DMLabel bdLabel = NULL, rgLabel = NULL;
1440: PetscBool higherOrder = PETSC_FALSE;
1441: PetscInt Nd = coordDim * coordDim, f, vStart, vEnd;
1442: 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[]);
1444: PetscCall(PetscMalloc(1, &funcs));
1445: funcs[0] = identityFunc;
1447: /* Setup finite element spaces */
1448: PetscCall(DMClone(dm, &dmGrad));
1449: PetscCall(DMClone(dm, &dmHess));
1450: PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple fields not yet considered"); // TODO
1451: for (f = 0; f < Nf; ++f) {
1452: PetscFE fe, feGrad, feHess;
1453: PetscDualSpace Q;
1454: PetscSpace space;
1455: DM K;
1456: PetscQuadrature q;
1457: PetscInt Nc, qorder, p;
1458: const char *prefix;
1460: PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1461: PetscCall(PetscFEGetNumComponents(fe, &Nc));
1462: PetscCheck(Nc <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple components not yet considered"); // TODO
1463: PetscCall(PetscFEGetBasisSpace(fe, &space));
1464: PetscCall(PetscSpaceGetDegree(space, NULL, &p));
1465: if (p > 1) higherOrder = PETSC_TRUE;
1466: PetscCall(PetscFEGetDualSpace(fe, &Q));
1467: PetscCall(PetscDualSpaceGetDM(Q, &K));
1468: PetscCall(DMPlexGetDepthStratum(K, 0, &vStart, &vEnd));
1469: PetscCall(PetscFEGetQuadrature(fe, &q));
1470: PetscCall(PetscQuadratureGetOrder(q, &qorder));
1471: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix));
1472: PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmGrad), dim, Nc * coordDim, PETSC_TRUE, prefix, qorder, &feGrad));
1473: PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmHess), dim, Nc * Nd, PETSC_TRUE, prefix, qorder, &feHess));
1474: PetscCall(DMSetField(dmGrad, f, NULL, (PetscObject)feGrad));
1475: PetscCall(DMSetField(dmHess, f, NULL, (PetscObject)feHess));
1476: PetscCall(DMCreateDS(dmGrad));
1477: PetscCall(DMCreateDS(dmHess));
1478: PetscCall(PetscFEDestroy(&feGrad));
1479: PetscCall(PetscFEDestroy(&feHess));
1480: }
1481: /* Compute vertexwise gradients from cellwise gradients */
1482: PetscCall(DMCreateLocalVector(dmGrad, &xGrad));
1483: PetscCall(VecViewFromOptions(locX, NULL, "-sol_adapt_loc_pre_view"));
1484: PetscCall(DMPlexComputeGradientClementInterpolant(dm, locX, xGrad));
1485: PetscCall(VecViewFromOptions(xGrad, NULL, "-adapt_gradient_view"));
1486: /* Compute vertexwise Hessians from cellwise Hessians */
1487: PetscCall(DMCreateLocalVector(dmHess, &xHess));
1488: PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, xGrad, xHess));
1489: PetscCall(VecViewFromOptions(xHess, NULL, "-adapt_hessian_view"));
1490: PetscCall(VecDestroy(&xGrad));
1491: PetscCall(DMDestroy(&dmGrad));
1492: /* Compute L-p normalized metric */
1493: PetscCall(DMClone(dm, &dmMetric));
1494: N = adaptor->Nadapt >= 0 ? adaptor->Nadapt : PetscPowRealInt(adaptor->refinementFactor, dim) * ((PetscReal)(vEnd - vStart));
1495: // TODO This was where the old monitor was, figure out how to show metric and target N
1496: PetscCall(DMPlexMetricSetTargetComplexity(dmMetric, N));
1497: if (higherOrder) {
1498: /* Project Hessian into P1 space, if required */
1499: PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1500: PetscCall(DMProjectFieldLocal(dmMetric, 0.0, xHess, funcs, INSERT_ALL_VALUES, metric));
1501: PetscCall(VecDestroy(&xHess));
1502: xHess = metric;
1503: }
1504: PetscCall(PetscFree(funcs));
1505: PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1506: PetscCall(DMPlexMetricDeterminantCreate(dmMetric, 0, &determinant, &dmDet));
1507: PetscCall(DMPlexMetricNormalize(dmMetric, xHess, PETSC_TRUE, PETSC_TRUE, metric, determinant));
1508: PetscCall(VecDestroy(&determinant));
1509: PetscCall(DMDestroy(&dmDet));
1510: PetscCall(VecDestroy(&xHess));
1511: PetscCall(DMDestroy(&dmHess));
1512: /* Adapt DM from metric */
1513: PetscCall(DMGetLabel(dm, "marker", &bdLabel));
1514: PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &odm));
1515: PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, NULL));
1516: adapted = PETSC_TRUE;
1517: /* Cleanup */
1518: PetscCall(VecDestroy(&metric));
1519: PetscCall(DMDestroy(&dmMetric));
1520: } break;
1521: default:
1522: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid adaptation type: %d", adaptor->adaptCriterion);
1523: }
1524: PetscCall(DMAdaptorPostAdapt(adaptor));
1525: PetscCall(DMRestoreLocalVector(dm, &locX));
1526: /* If DM was adapted, replace objects and recreate solution */
1527: if (adapted) {
1528: const char *name;
1530: PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1531: PetscCall(PetscObjectSetName((PetscObject)odm, name));
1532: /* Reconfigure solver */
1533: PetscCall(SNESReset(adaptor->snes));
1534: PetscCall(SNESSetDM(adaptor->snes, odm));
1535: PetscCall(DMAdaptorSetSolver(adaptor, adaptor->snes));
1536: PetscCall(DMPlexSetSNESLocalFEM(odm, PETSC_FALSE, ctx));
1537: PetscCall(SNESSetFromOptions(adaptor->snes));
1538: /* Transfer system */
1539: PetscCall(DMCopyDisc(dm, odm));
1540: /* Transfer solution */
1541: PetscCall(DMCreateGlobalVector(odm, &ox));
1542: PetscCall(PetscObjectGetName((PetscObject)x, &name));
1543: PetscCall(PetscObjectSetName((PetscObject)ox, name));
1544: PetscCall(DMAdaptorTransferSolution(adaptor, dm, x, odm, ox));
1545: /* Cleanup adaptivity info */
1546: if (adaptIter > 0) PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(comm)));
1547: PetscCall(DMForestSetAdaptivityForest(dm, NULL)); /* clear internal references to the previous dm */
1548: PetscCall(DMDestroy(&dm));
1549: PetscCall(VecDestroy(&x));
1550: *adm = odm;
1551: *ax = ox;
1552: } else {
1553: *adm = dm;
1554: *ax = x;
1555: adaptIter = numAdapt;
1556: }
1557: if (adaptIter < numAdapt - 1) {
1558: PetscCall(DMViewFromOptions(odm, NULL, "-dm_adapt_iter_view"));
1559: PetscCall(VecViewFromOptions(ox, NULL, "-sol_adapt_iter_view"));
1560: }
1561: }
1562: PetscCall(DMViewFromOptions(*adm, NULL, "-dm_adapt_view"));
1563: PetscCall(VecViewFromOptions(*ax, NULL, "-sol_adapt_view"));
1564: PetscFunctionReturn(PETSC_SUCCESS);
1565: }
1567: /*@
1568: DMAdaptorAdapt - Creates a new `DM` that is adapted to the problem
1570: Not Collective
1572: Input Parameters:
1573: + adaptor - The `DMAdaptor` object
1574: . x - The global approximate solution
1575: - strategy - The adaptation strategy, see `DMAdaptationStrategy`
1577: Output Parameters:
1578: + adm - The adapted `DM`
1579: - ax - The adapted solution
1581: Options Database Keys:
1582: + -snes_adapt (initial|sequential|multigrid) - adaption strategy, see `DMAdaptationStrategy`
1583: . -adapt_gradient_view - View the Clement interpolant of the solution gradient
1584: . -adapt_hessian_view - View the Clement interpolant of the solution Hessian
1585: - -adapt_metric_view - View the metric tensor for adaptive mesh refinement
1587: Level: intermediate
1589: .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`
1590: @*/
1591: PetscErrorCode DMAdaptorAdapt(DMAdaptor adaptor, Vec x, DMAdaptationStrategy strategy, DM *adm, Vec *ax)
1592: {
1593: PetscFunctionBegin;
1594: switch (strategy) {
1595: case DM_ADAPTATION_INITIAL:
1596: PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_FALSE, adm, ax));
1597: break;
1598: case DM_ADAPTATION_SEQUENTIAL:
1599: PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_TRUE, adm, ax));
1600: break;
1601: default:
1602: SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Unrecognized adaptation strategy %d", strategy);
1603: }
1604: PetscFunctionReturn(PETSC_SUCCESS);
1605: }
1607: /*@C
1608: DMAdaptorGetMixedSetupFunction - Get the function setting up the mixed problem, if it exists
1610: Not Collective
1612: Input Parameter:
1613: . adaptor - the `DMAdaptor`
1615: Output Parameter:
1616: . setupFunc - the function setting up the mixed problem, or `NULL`
1618: Level: advanced
1620: .seealso: `DMAdaptor`, `DMAdaptorSetMixedSetupFunction()`, `DMAdaptorAdapt()`
1621: @*/
1622: PetscErrorCode DMAdaptorGetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (**setupFunc)(DMAdaptor, DM))
1623: {
1624: PetscFunctionBegin;
1626: PetscAssertPointer(setupFunc, 2);
1627: *setupFunc = adaptor->ops->mixedsetup;
1628: PetscFunctionReturn(PETSC_SUCCESS);
1629: }
1631: /*@C
1632: DMAdaptorSetMixedSetupFunction - Set the function setting up the mixed problem
1634: Not Collective
1636: Input Parameters:
1637: + adaptor - the `DMAdaptor`
1638: - setupFunc - the function setting up the mixed problem
1640: Calling sequence of setupFunc:
1641: + adaptor - the `DMAdaptor`
1642: - dm - the `DM`
1644: Level: advanced
1646: .seealso: `DMAdaptor`, `DMAdaptorGetMixedSetupFunction()`, `DMAdaptorAdapt()`
1647: @*/
1648: PetscErrorCode DMAdaptorSetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (*setupFunc)(DMAdaptor adaptor, DM dm))
1649: {
1650: PetscFunctionBegin;
1653: adaptor->ops->mixedsetup = setupFunc;
1654: PetscFunctionReturn(PETSC_SUCCESS);
1655: }
1657: /*@
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`, `DMAdaptorSetCriterion()`, `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: /*@
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: }