xref: /petsc/src/snes/utils/dm/dmadapt.c (revision f8d70eaab67b53d355abe83957ffb5a3120bcf52)
1c6f61ee2SBarry Smith #include <petscdmadaptor.h> /*I "petscdmadaptor.h" I*/
2c6f61ee2SBarry Smith #include <petscdmplex.h>
3c6f61ee2SBarry Smith #include <petscdmforest.h>
4c6f61ee2SBarry Smith #include <petscds.h>
5c6f61ee2SBarry Smith #include <petscblaslapack.h>
63a336bb1SMatthew G. Knepley #include <petscsnes.h>
78b724c91SMatthew G. Knepley #include <petscdraw.h>
8c6f61ee2SBarry Smith 
9c6f61ee2SBarry Smith #include <petsc/private/dmadaptorimpl.h>
10c6f61ee2SBarry Smith #include <petsc/private/dmpleximpl.h>
11c6f61ee2SBarry Smith #include <petsc/private/petscfeimpl.h>
12c6f61ee2SBarry Smith 
133a336bb1SMatthew G. Knepley PetscClassId DMADAPTOR_CLASSID;
14c6f61ee2SBarry Smith 
153a336bb1SMatthew G. Knepley PetscFunctionList DMAdaptorList              = NULL;
163a336bb1SMatthew G. Knepley PetscBool         DMAdaptorRegisterAllCalled = PETSC_FALSE;
173a336bb1SMatthew G. Knepley 
188b724c91SMatthew G. Knepley PetscFunctionList DMAdaptorMonitorList              = NULL;
198b724c91SMatthew G. Knepley PetscFunctionList DMAdaptorMonitorCreateList        = NULL;
208b724c91SMatthew G. Knepley PetscFunctionList DMAdaptorMonitorDestroyList       = NULL;
218b724c91SMatthew G. Knepley PetscBool         DMAdaptorMonitorRegisterAllCalled = PETSC_FALSE;
228b724c91SMatthew G. Knepley 
238b724c91SMatthew G. Knepley const char *const DMAdaptationCriteria[] = {"NONE", "REFINE", "LABEL", "METRIC", "DMAdaptationCriterion", "DM_ADAPTATION_", NULL};
248b724c91SMatthew G. Knepley 
253a336bb1SMatthew G. Knepley /*@C
263a336bb1SMatthew G. Knepley   DMAdaptorRegister - Adds a new adaptor component implementation
273a336bb1SMatthew G. Knepley 
283a336bb1SMatthew G. Knepley   Not Collective
293a336bb1SMatthew G. Knepley 
303a336bb1SMatthew G. Knepley   Input Parameters:
313a336bb1SMatthew G. Knepley + name        - The name of a new user-defined creation routine
323a336bb1SMatthew G. Knepley - create_func - The creation routine
333a336bb1SMatthew G. Knepley 
343a336bb1SMatthew G. Knepley   Example Usage:
353a336bb1SMatthew G. Knepley .vb
363a336bb1SMatthew G. Knepley   DMAdaptorRegister("my_adaptor", MyAdaptorCreate);
373a336bb1SMatthew G. Knepley .ve
383a336bb1SMatthew G. Knepley 
393a336bb1SMatthew G. Knepley   Then, your adaptor type can be chosen with the procedural interface via
403a336bb1SMatthew G. Knepley .vb
413a336bb1SMatthew G. Knepley   DMAdaptorCreate(MPI_Comm, DMAdaptor *);
423a336bb1SMatthew G. Knepley   DMAdaptorSetType(DMAdaptor, "my_adaptor");
433a336bb1SMatthew G. Knepley .ve
443a336bb1SMatthew G. Knepley   or at runtime via the option
453a336bb1SMatthew G. Knepley .vb
463a336bb1SMatthew G. Knepley   -adaptor_type my_adaptor
473a336bb1SMatthew G. Knepley .ve
483a336bb1SMatthew G. Knepley 
493a336bb1SMatthew G. Knepley   Level: advanced
503a336bb1SMatthew G. Knepley 
513a336bb1SMatthew G. Knepley   Note:
523a336bb1SMatthew G. Knepley   `DMAdaptorRegister()` may be called multiple times to add several user-defined adaptors
533a336bb1SMatthew G. Knepley 
543a336bb1SMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptor`, `DMAdaptorRegisterAll()`, `DMAdaptorRegisterDestroy()`
553a336bb1SMatthew G. Knepley @*/
563a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorRegister(const char name[], PetscErrorCode (*create_func)(DMAdaptor))
57c6f61ee2SBarry Smith {
583a336bb1SMatthew G. Knepley   PetscFunctionBegin;
593a336bb1SMatthew G. Knepley   PetscCall(DMInitializePackage());
603a336bb1SMatthew G. Knepley   PetscCall(PetscFunctionListAdd(&DMAdaptorList, name, create_func));
613a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
623a336bb1SMatthew G. Knepley }
633a336bb1SMatthew G. Knepley 
643a336bb1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Gradient(DMAdaptor);
653a336bb1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Flux(DMAdaptor);
663a336bb1SMatthew G. Knepley 
673a336bb1SMatthew G. Knepley /*@C
683a336bb1SMatthew G. Knepley   DMAdaptorRegisterAll - Registers all of the adaptor components in the `DM` package.
693a336bb1SMatthew G. Knepley 
703a336bb1SMatthew G. Knepley   Not Collective
713a336bb1SMatthew G. Knepley 
723a336bb1SMatthew G. Knepley   Level: advanced
733a336bb1SMatthew G. Knepley 
743a336bb1SMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptorType`, `DMRegisterAll()`, `DMAdaptorRegisterDestroy()`
753a336bb1SMatthew G. Knepley @*/
763a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorRegisterAll(void)
773a336bb1SMatthew G. Knepley {
783a336bb1SMatthew G. Knepley   PetscFunctionBegin;
793a336bb1SMatthew G. Knepley   if (DMAdaptorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
803a336bb1SMatthew G. Knepley   DMAdaptorRegisterAllCalled = PETSC_TRUE;
813a336bb1SMatthew G. Knepley 
823a336bb1SMatthew G. Knepley   PetscCall(DMAdaptorRegister(DMADAPTORGRADIENT, DMAdaptorCreate_Gradient));
833a336bb1SMatthew G. Knepley   PetscCall(DMAdaptorRegister(DMADAPTORFLUX, DMAdaptorCreate_Flux));
843a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
853a336bb1SMatthew G. Knepley }
863a336bb1SMatthew G. Knepley 
873a336bb1SMatthew G. Knepley /*@C
883a336bb1SMatthew G. Knepley   DMAdaptorRegisterDestroy - This function destroys the registered `DMAdaptorType`. It is called from `PetscFinalize()`.
893a336bb1SMatthew G. Knepley 
903a336bb1SMatthew G. Knepley   Not collective
913a336bb1SMatthew G. Knepley 
923a336bb1SMatthew G. Knepley   Level: developer
933a336bb1SMatthew G. Knepley 
948b724c91SMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptorRegisterAll()`, `DMAdaptorType`, `PetscFinalize()`
953a336bb1SMatthew G. Knepley @*/
963a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorRegisterDestroy(void)
973a336bb1SMatthew G. Knepley {
983a336bb1SMatthew G. Knepley   PetscFunctionBegin;
993a336bb1SMatthew G. Knepley   PetscCall(PetscFunctionListDestroy(&DMAdaptorList));
1003a336bb1SMatthew G. Knepley   DMAdaptorRegisterAllCalled = PETSC_FALSE;
101c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
102c6f61ee2SBarry Smith }
103c6f61ee2SBarry Smith 
1048b724c91SMatthew G. Knepley static PetscErrorCode DMAdaptorMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[])
1058b724c91SMatthew G. Knepley {
1068b724c91SMatthew G. Knepley   PetscFunctionBegin;
1078b724c91SMatthew G. Knepley   PetscCall(PetscStrncpy(key, name, PETSC_MAX_PATH_LEN));
1088b724c91SMatthew G. Knepley   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
1098b724c91SMatthew G. Knepley   PetscCall(PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN));
1108b724c91SMatthew G. Knepley   PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN));
1118b724c91SMatthew G. Knepley   PetscCall(PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN));
1128b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1138b724c91SMatthew G. Knepley }
1148b724c91SMatthew G. Knepley 
1158b724c91SMatthew G. Knepley /*@C
1168b724c91SMatthew G. Knepley   DMAdaptorMonitorRegister -  Registers a mesh adaptation monitor routine that may be accessed with `DMAdaptorMonitorSetFromOptions()`
1178b724c91SMatthew G. Knepley 
1188b724c91SMatthew G. Knepley   Not Collective
1198b724c91SMatthew G. Knepley 
1208b724c91SMatthew G. Knepley   Input Parameters:
1218b724c91SMatthew G. Knepley + name    - name of a new monitor routine
1228b724c91SMatthew G. Knepley . vtype   - A `PetscViewerType` for the output
1238b724c91SMatthew G. Knepley . format  - A `PetscViewerFormat` for the output
1248b724c91SMatthew G. Knepley . monitor - Monitor routine
1258b724c91SMatthew G. Knepley . create  - Creation routine, or `NULL`
1268b724c91SMatthew G. Knepley - destroy - Destruction routine, or `NULL`
1278b724c91SMatthew G. Knepley 
1288b724c91SMatthew G. Knepley   Level: advanced
1298b724c91SMatthew G. Knepley 
1308b724c91SMatthew G. Knepley   Note:
1318b724c91SMatthew G. Knepley   `DMAdaptorMonitorRegister()` may be called multiple times to add several user-defined monitors.
1328b724c91SMatthew G. Knepley 
1338b724c91SMatthew G. Knepley   Example Usage:
1348b724c91SMatthew G. Knepley .vb
1358b724c91SMatthew G. Knepley   DMAdaptorMonitorRegister("my_monitor", PETSCVIEWERASCII, PETSC_VIEWER_ASCII_INFO_DETAIL, MyMonitor, NULL, NULL);
1368b724c91SMatthew G. Knepley .ve
1378b724c91SMatthew G. Knepley 
1388b724c91SMatthew G. Knepley   Then, your monitor can be chosen with the procedural interface via
1398b724c91SMatthew G. Knepley .vb
1408b724c91SMatthew G. Knepley   DMAdaptorMonitorSetFromOptions(ksp, "-adaptor_monitor_my_monitor", "my_monitor", NULL)
1418b724c91SMatthew G. Knepley .ve
1428b724c91SMatthew G. Knepley   or at runtime via the option `-adaptor_monitor_my_monitor`
1438b724c91SMatthew G. Knepley 
1448b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorRegisterAll()`, `DMAdaptorMonitorSetFromOptions()`
1458b724c91SMatthew G. Knepley @*/
1468b724c91SMatthew G. Knepley 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 **))
1478b724c91SMatthew G. Knepley {
1488b724c91SMatthew G. Knepley   char key[PETSC_MAX_PATH_LEN];
1498b724c91SMatthew G. Knepley 
1508b724c91SMatthew G. Knepley   PetscFunctionBegin;
1518b724c91SMatthew G. Knepley   PetscCall(SNESInitializePackage());
1528b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorMakeKey_Internal(name, vtype, format, key));
1538b724c91SMatthew G. Knepley   PetscCall(PetscFunctionListAdd(&DMAdaptorMonitorList, key, monitor));
1548b724c91SMatthew G. Knepley   if (create) PetscCall(PetscFunctionListAdd(&DMAdaptorMonitorCreateList, key, create));
1558b724c91SMatthew G. Knepley   if (destroy) PetscCall(PetscFunctionListAdd(&DMAdaptorMonitorDestroyList, key, destroy));
1568b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1578b724c91SMatthew G. Knepley }
1588b724c91SMatthew G. Knepley 
1598b724c91SMatthew G. Knepley /*@C
1608b724c91SMatthew G. Knepley   DMAdaptorMonitorRegisterDestroy - This function destroys the registered monitors for `DMAdaptor`. It is called from `PetscFinalize()`.
1618b724c91SMatthew G. Knepley 
1628b724c91SMatthew G. Knepley   Not collective
1638b724c91SMatthew G. Knepley 
1648b724c91SMatthew G. Knepley   Level: developer
1658b724c91SMatthew G. Knepley 
1668b724c91SMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptorMonitorRegisterAll()`, `DMAdaptor`, `PetscFinalize()`
1678b724c91SMatthew G. Knepley @*/
1688b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorRegisterDestroy(void)
1698b724c91SMatthew G. Knepley {
1708b724c91SMatthew G. Knepley   PetscFunctionBegin;
1718b724c91SMatthew G. Knepley   PetscCall(PetscFunctionListDestroy(&DMAdaptorMonitorList));
1728b724c91SMatthew G. Knepley   PetscCall(PetscFunctionListDestroy(&DMAdaptorMonitorCreateList));
1738b724c91SMatthew G. Knepley   PetscCall(PetscFunctionListDestroy(&DMAdaptorMonitorDestroyList));
1748b724c91SMatthew G. Knepley   DMAdaptorMonitorRegisterAllCalled = PETSC_FALSE;
1758b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1768b724c91SMatthew G. Knepley }
1778b724c91SMatthew G. Knepley 
178c6f61ee2SBarry Smith /*@
179c6f61ee2SBarry Smith   DMAdaptorCreate - Create a `DMAdaptor` object. Its purpose is to construct a adaptation `DMLabel` or metric `Vec` that can be used to modify the `DM`.
180c6f61ee2SBarry Smith 
181c6f61ee2SBarry Smith   Collective
182c6f61ee2SBarry Smith 
183c6f61ee2SBarry Smith   Input Parameter:
184c6f61ee2SBarry Smith . comm - The communicator for the `DMAdaptor` object
185c6f61ee2SBarry Smith 
186c6f61ee2SBarry Smith   Output Parameter:
187c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
188c6f61ee2SBarry Smith 
189c6f61ee2SBarry Smith   Level: beginner
190c6f61ee2SBarry Smith 
191df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorDestroy()`, `DMAdaptorAdapt()`, `PetscConvEst`, `PetscConvEstCreate()`
192c6f61ee2SBarry Smith @*/
193c6f61ee2SBarry Smith PetscErrorCode DMAdaptorCreate(MPI_Comm comm, DMAdaptor *adaptor)
194c6f61ee2SBarry Smith {
195c6f61ee2SBarry Smith   VecTaggerBox refineBox, coarsenBox;
196c6f61ee2SBarry Smith 
197c6f61ee2SBarry Smith   PetscFunctionBegin;
198c6f61ee2SBarry Smith   PetscAssertPointer(adaptor, 2);
199c6f61ee2SBarry Smith   PetscCall(PetscSysInitializePackage());
200c6f61ee2SBarry Smith 
2013a336bb1SMatthew G. Knepley   PetscCall(PetscHeaderCreate(*adaptor, DMADAPTOR_CLASSID, "DMAdaptor", "DM Adaptor", "DMAdaptor", comm, DMAdaptorDestroy, DMAdaptorView));
202c6f61ee2SBarry Smith   (*adaptor)->adaptCriterion   = DM_ADAPTATION_NONE;
203c6f61ee2SBarry Smith   (*adaptor)->numSeq           = 1;
204c6f61ee2SBarry Smith   (*adaptor)->Nadapt           = -1;
205c6f61ee2SBarry Smith   (*adaptor)->refinementFactor = 2.0;
206c6f61ee2SBarry Smith   refineBox.min = refineBox.max = PETSC_MAX_REAL;
207c6f61ee2SBarry Smith   PetscCall(VecTaggerCreate(PetscObjectComm((PetscObject)*adaptor), &(*adaptor)->refineTag));
208c6f61ee2SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(*adaptor)->refineTag, "refine_"));
209c6f61ee2SBarry Smith   PetscCall(VecTaggerSetType((*adaptor)->refineTag, VECTAGGERABSOLUTE));
210c6f61ee2SBarry Smith   PetscCall(VecTaggerAbsoluteSetBox((*adaptor)->refineTag, &refineBox));
211c6f61ee2SBarry Smith   coarsenBox.min = coarsenBox.max = PETSC_MAX_REAL;
212c6f61ee2SBarry Smith   PetscCall(VecTaggerCreate(PetscObjectComm((PetscObject)*adaptor), &(*adaptor)->coarsenTag));
213c6f61ee2SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(*adaptor)->coarsenTag, "coarsen_"));
214c6f61ee2SBarry Smith   PetscCall(VecTaggerSetType((*adaptor)->coarsenTag, VECTAGGERABSOLUTE));
215c6f61ee2SBarry Smith   PetscCall(VecTaggerAbsoluteSetBox((*adaptor)->coarsenTag, &coarsenBox));
216c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
217c6f61ee2SBarry Smith }
218c6f61ee2SBarry Smith 
219c6f61ee2SBarry Smith /*@
220c6f61ee2SBarry Smith   DMAdaptorDestroy - Destroys a `DMAdaptor` object
221c6f61ee2SBarry Smith 
222c6f61ee2SBarry Smith   Collective
223c6f61ee2SBarry Smith 
224c6f61ee2SBarry Smith   Input Parameter:
225c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
226c6f61ee2SBarry Smith 
227c6f61ee2SBarry Smith   Level: beginner
228c6f61ee2SBarry Smith 
229df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
230c6f61ee2SBarry Smith @*/
231c6f61ee2SBarry Smith PetscErrorCode DMAdaptorDestroy(DMAdaptor *adaptor)
232c6f61ee2SBarry Smith {
233c6f61ee2SBarry Smith   PetscFunctionBegin;
234c6f61ee2SBarry Smith   if (!*adaptor) PetscFunctionReturn(PETSC_SUCCESS);
2353a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(*adaptor, DMADAPTOR_CLASSID, 1);
236f4f49eeaSPierre Jolivet   if (--((PetscObject)*adaptor)->refct > 0) {
237c6f61ee2SBarry Smith     *adaptor = NULL;
238c6f61ee2SBarry Smith     PetscFunctionReturn(PETSC_SUCCESS);
239c6f61ee2SBarry Smith   }
240c6f61ee2SBarry Smith   PetscCall(VecTaggerDestroy(&(*adaptor)->refineTag));
241c6f61ee2SBarry Smith   PetscCall(VecTaggerDestroy(&(*adaptor)->coarsenTag));
242c6f61ee2SBarry Smith   PetscCall(PetscFree2((*adaptor)->exactSol, (*adaptor)->exactCtx));
2438b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorCancel(*adaptor));
244c6f61ee2SBarry Smith   PetscCall(PetscHeaderDestroy(adaptor));
245c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
246c6f61ee2SBarry Smith }
247c6f61ee2SBarry Smith 
248c6f61ee2SBarry Smith /*@
2493a336bb1SMatthew G. Knepley   DMAdaptorSetType - Sets the particular implementation for a adaptor.
2503a336bb1SMatthew G. Knepley 
2513a336bb1SMatthew G. Knepley   Collective
2523a336bb1SMatthew G. Knepley 
2533a336bb1SMatthew G. Knepley   Input Parameters:
2543a336bb1SMatthew G. Knepley + adaptor - The `DMAdaptor`
2553a336bb1SMatthew G. Knepley - method  - The name of the adaptor type
2563a336bb1SMatthew G. Knepley 
2573a336bb1SMatthew G. Knepley   Options Database Key:
2583a336bb1SMatthew G. Knepley . -adaptor_type <type> - Sets the adaptor type; see `DMAdaptorType`
2593a336bb1SMatthew G. Knepley 
2603a336bb1SMatthew G. Knepley   Level: intermediate
2613a336bb1SMatthew G. Knepley 
2623a336bb1SMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptor`, `DMAdaptorType`, `DMAdaptorGetType()`, `DMAdaptorCreate()`
2633a336bb1SMatthew G. Knepley @*/
2643a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorSetType(DMAdaptor adaptor, DMAdaptorType method)
2653a336bb1SMatthew G. Knepley {
2663a336bb1SMatthew G. Knepley   PetscErrorCode (*r)(DMAdaptor);
2673a336bb1SMatthew G. Knepley   PetscBool match;
2683a336bb1SMatthew G. Knepley 
2693a336bb1SMatthew G. Knepley   PetscFunctionBegin;
2703a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
2713a336bb1SMatthew G. Knepley   PetscCall(PetscObjectTypeCompare((PetscObject)adaptor, method, &match));
2723a336bb1SMatthew G. Knepley   if (match) PetscFunctionReturn(PETSC_SUCCESS);
2733a336bb1SMatthew G. Knepley 
2743a336bb1SMatthew G. Knepley   PetscCall(DMAdaptorRegisterAll());
2753a336bb1SMatthew G. Knepley   PetscCall(PetscFunctionListFind(DMAdaptorList, method, &r));
2763a336bb1SMatthew G. Knepley   PetscCheck(r, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DMAdaptor type: %s", method);
2773a336bb1SMatthew G. Knepley 
2783a336bb1SMatthew G. Knepley   PetscTryTypeMethod(adaptor, destroy);
2793a336bb1SMatthew G. Knepley   PetscCall(PetscMemzero(adaptor->ops, sizeof(*adaptor->ops)));
2803a336bb1SMatthew G. Knepley   PetscCall(PetscObjectChangeTypeName((PetscObject)adaptor, method));
2813a336bb1SMatthew G. Knepley   PetscCall((*r)(adaptor));
2823a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2833a336bb1SMatthew G. Knepley }
2843a336bb1SMatthew G. Knepley 
2853a336bb1SMatthew G. Knepley /*@
2863a336bb1SMatthew G. Knepley   DMAdaptorGetType - Gets the type name (as a string) from the adaptor.
2873a336bb1SMatthew G. Knepley 
2883a336bb1SMatthew G. Knepley   Not Collective
2893a336bb1SMatthew G. Knepley 
2903a336bb1SMatthew G. Knepley   Input Parameter:
2913a336bb1SMatthew G. Knepley . adaptor - The `DMAdaptor`
2923a336bb1SMatthew G. Knepley 
2933a336bb1SMatthew G. Knepley   Output Parameter:
2943a336bb1SMatthew G. Knepley . type - The `DMAdaptorType` name
2953a336bb1SMatthew G. Knepley 
2963a336bb1SMatthew G. Knepley   Level: intermediate
2973a336bb1SMatthew G. Knepley 
2983a336bb1SMatthew G. Knepley .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMAdaptor`, `DMAdaptorType`, `DMAdaptorSetType()`, `DMAdaptorCreate()`
2993a336bb1SMatthew G. Knepley @*/
3003a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorGetType(DMAdaptor adaptor, DMAdaptorType *type)
3013a336bb1SMatthew G. Knepley {
3023a336bb1SMatthew G. Knepley   PetscFunctionBegin;
3033a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
3043a336bb1SMatthew G. Knepley   PetscAssertPointer(type, 2);
3053a336bb1SMatthew G. Knepley   PetscCall(DMAdaptorRegisterAll());
3063a336bb1SMatthew G. Knepley   *type = ((PetscObject)adaptor)->type_name;
3073a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3083a336bb1SMatthew G. Knepley }
3093a336bb1SMatthew G. Knepley 
3108b724c91SMatthew G. Knepley static PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
3118b724c91SMatthew G. Knepley {
3128b724c91SMatthew G. Knepley   PetscFunctionBegin;
3138b724c91SMatthew G. Knepley   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
3148b724c91SMatthew G. Knepley   (*vf)->data = ctx;
3158b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3168b724c91SMatthew G. Knepley }
3178b724c91SMatthew G. Knepley 
3188b724c91SMatthew G. Knepley /*@C
3198b724c91SMatthew G. Knepley   DMAdaptorMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
3208b724c91SMatthew G. Knepley   the error etc.
3218b724c91SMatthew G. Knepley 
3228b724c91SMatthew G. Knepley   Logically Collective
3238b724c91SMatthew G. Knepley 
3248b724c91SMatthew G. Knepley   Input Parameters:
3258b724c91SMatthew G. Knepley + adaptor        - the `DMAdaptor`
3268b724c91SMatthew G. Knepley . monitor        - pointer to function (if this is `NULL`, it turns off monitoring
3278b724c91SMatthew G. Knepley . ctx            - [optional] context for private data for the monitor routine (use `NULL` if no context is needed)
3288b724c91SMatthew G. Knepley - monitordestroy - [optional] routine that frees monitor context (may be `NULL`)
3298b724c91SMatthew G. Knepley 
3308b724c91SMatthew G. Knepley   Calling sequence of `monitor`:
3318b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
3328b724c91SMatthew G. Knepley . it      - iteration number
3338b724c91SMatthew G. Knepley . odm     - the original `DM`
3348b724c91SMatthew G. Knepley . adm     - the adapted `DM`
3358b724c91SMatthew G. Knepley . Nf      - number of fields
3368b724c91SMatthew G. Knepley . enorms  - (estimated) 2-norm of the error for each field
3378b724c91SMatthew G. Knepley . error   - `Vec` of cellwise errors
3388b724c91SMatthew G. Knepley - ctx     - optional monitoring context, as set by `DMAdaptorMonitorSet()`
3398b724c91SMatthew G. Knepley 
3408b724c91SMatthew G. Knepley   Calling sequence of `monitordestroy`:
3418b724c91SMatthew G. Knepley . ctx - optional monitoring context, as set by `DMAdaptorMonitorSet()`
3428b724c91SMatthew G. Knepley 
3438b724c91SMatthew G. Knepley   Options Database Keys:
3448b724c91SMatthew G. Knepley + -adaptor_monitor_size                - sets `DMAdaptorMonitorSize()`
3458b724c91SMatthew G. Knepley . -adaptor_monitor_error               - sets `DMAdaptorMonitorError()`
3468b724c91SMatthew G. Knepley . -adaptor_monitor_error draw          - sets `DMAdaptorMonitorErrorDraw()` and plots error
3478b724c91SMatthew G. Knepley . -adaptor_monitor_error draw::draw_lg - sets `DMAdaptorMonitorErrorDrawLG()` and plots error
3488b724c91SMatthew G. Knepley - -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.
3498b724c91SMatthew G. Knepley 
3508b724c91SMatthew G. Knepley   Level: beginner
3518b724c91SMatthew G. Knepley 
3528b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptor`
3538b724c91SMatthew G. Knepley @*/
3548b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorSet(DMAdaptor adaptor, PetscErrorCode (*monitor)(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, void *ctx), void *ctx, PetscErrorCode (*monitordestroy)(void **ctx))
3558b724c91SMatthew G. Knepley {
3568b724c91SMatthew G. Knepley   PetscBool identical;
3578b724c91SMatthew G. Knepley 
3588b724c91SMatthew G. Knepley   PetscFunctionBegin;
3598b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
3608b724c91SMatthew G. Knepley   for (PetscInt i = 0; i < adaptor->numbermonitors; i++) {
3618b724c91SMatthew G. Knepley     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))monitor, ctx, monitordestroy, (PetscErrorCode (*)(void))adaptor->monitor[i], adaptor->monitorcontext[i], adaptor->monitordestroy[i], &identical));
3628b724c91SMatthew G. Knepley     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
3638b724c91SMatthew G. Knepley   }
3648b724c91SMatthew G. Knepley   PetscCheck(adaptor->numbermonitors < MAXDMADAPTORMONITORS, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_OUTOFRANGE, "Too many DMAdaptor monitors set");
3658b724c91SMatthew G. Knepley   adaptor->monitor[adaptor->numbermonitors]          = monitor;
3668b724c91SMatthew G. Knepley   adaptor->monitordestroy[adaptor->numbermonitors]   = monitordestroy;
3678b724c91SMatthew G. Knepley   adaptor->monitorcontext[adaptor->numbermonitors++] = (void *)ctx;
3688b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3698b724c91SMatthew G. Knepley }
3708b724c91SMatthew G. Knepley 
3718b724c91SMatthew G. Knepley /*@
3728b724c91SMatthew G. Knepley   DMAdaptorMonitorCancel - Clears all monitors for a `DMAdaptor` object.
3738b724c91SMatthew G. Knepley 
3748b724c91SMatthew G. Knepley   Logically Collective
3758b724c91SMatthew G. Knepley 
3768b724c91SMatthew G. Knepley   Input Parameter:
3778b724c91SMatthew G. Knepley . adaptor - the `DMAdaptor`
3788b724c91SMatthew G. Knepley 
3798b724c91SMatthew G. Knepley   Options Database Key:
3808b724c91SMatthew G. Knepley . -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.
3818b724c91SMatthew G. Knepley 
3828b724c91SMatthew G. Knepley   Level: intermediate
3838b724c91SMatthew G. Knepley 
3848b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptorMonitorSet()`, `DMAdaptor`
3858b724c91SMatthew G. Knepley @*/
3868b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorCancel(DMAdaptor adaptor)
3878b724c91SMatthew G. Knepley {
3888b724c91SMatthew G. Knepley   PetscFunctionBegin;
3898b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
3908b724c91SMatthew G. Knepley   for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) {
3918b724c91SMatthew G. Knepley     if (adaptor->monitordestroy[i]) PetscCall((*adaptor->monitordestroy[i])(&adaptor->monitorcontext[i]));
3928b724c91SMatthew G. Knepley   }
3938b724c91SMatthew G. Knepley   adaptor->numbermonitors = 0;
3948b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3958b724c91SMatthew G. Knepley }
3968b724c91SMatthew G. Knepley 
3978b724c91SMatthew G. Knepley /*@C
3988b724c91SMatthew G. Knepley   DMAdaptorMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user in the options database
3998b724c91SMatthew G. Knepley 
4008b724c91SMatthew G. Knepley   Collective
4018b724c91SMatthew G. Knepley 
4028b724c91SMatthew G. Knepley   Input Parameters:
4038b724c91SMatthew G. Knepley + adaptor - `DMadaptor` object you wish to monitor
4048b724c91SMatthew G. Knepley . opt     - the command line option for this monitor
4058b724c91SMatthew G. Knepley . name    - the monitor type one is seeking
4068b724c91SMatthew G. Knepley - ctx     - An optional user context for the monitor, or `NULL`
4078b724c91SMatthew G. Knepley 
4088b724c91SMatthew G. Knepley   Level: developer
4098b724c91SMatthew G. Knepley 
4108b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptorMonitorRegister()`, `DMAdaptorMonitorSet()`, `PetscOptionsGetViewer()`
4118b724c91SMatthew G. Knepley @*/
4128b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorSetFromOptions(DMAdaptor adaptor, const char opt[], const char name[], void *ctx)
4138b724c91SMatthew G. Knepley {
4148b724c91SMatthew G. Knepley   PetscErrorCode (*mfunc)(DMAdaptor, PetscInt, DM, DM, PetscInt, PetscReal[], Vec, void *);
4158b724c91SMatthew G. Knepley   PetscErrorCode (*cfunc)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
4168b724c91SMatthew G. Knepley   PetscErrorCode (*dfunc)(PetscViewerAndFormat **);
4178b724c91SMatthew G. Knepley   PetscViewerAndFormat *vf;
4188b724c91SMatthew G. Knepley   PetscViewer           viewer;
4198b724c91SMatthew G. Knepley   PetscViewerFormat     format;
4208b724c91SMatthew G. Knepley   PetscViewerType       vtype;
4218b724c91SMatthew G. Knepley   char                  key[PETSC_MAX_PATH_LEN];
4228b724c91SMatthew G. Knepley   PetscBool             flg;
4238b724c91SMatthew G. Knepley   const char           *prefix = NULL;
4248b724c91SMatthew G. Knepley 
4258b724c91SMatthew G. Knepley   PetscFunctionBegin;
4268b724c91SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
4278b724c91SMatthew G. Knepley   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)adaptor), ((PetscObject)adaptor)->options, prefix, opt, &viewer, &format, &flg));
4288b724c91SMatthew G. Knepley   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
4298b724c91SMatthew G. Knepley 
4308b724c91SMatthew G. Knepley   PetscCall(PetscViewerGetType(viewer, &vtype));
4318b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorMakeKey_Internal(name, vtype, format, key));
4328b724c91SMatthew G. Knepley   PetscCall(PetscFunctionListFind(DMAdaptorMonitorList, key, &mfunc));
4338b724c91SMatthew G. Knepley   PetscCall(PetscFunctionListFind(DMAdaptorMonitorCreateList, key, &cfunc));
4348b724c91SMatthew G. Knepley   PetscCall(PetscFunctionListFind(DMAdaptorMonitorDestroyList, key, &dfunc));
4358b724c91SMatthew G. Knepley   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
4368b724c91SMatthew G. Knepley   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
4378b724c91SMatthew G. Knepley 
4388b724c91SMatthew G. Knepley   PetscCall((*cfunc)(viewer, format, ctx, &vf));
4398b724c91SMatthew G. Knepley   PetscCall(PetscViewerDestroy(&viewer));
4408b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorSet(adaptor, mfunc, vf, (PetscErrorCode (*)(void **))dfunc));
4418b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
4428b724c91SMatthew G. Knepley }
4438b724c91SMatthew G. Knepley 
4443a336bb1SMatthew G. Knepley /*@
445e03fd340SMatthew G. Knepley   DMAdaptorSetOptionsPrefix - Sets the prefix used for searching for all `DMAdaptor` options in the database.
446e03fd340SMatthew G. Knepley 
447e03fd340SMatthew G. Knepley   Logically Collective
448e03fd340SMatthew G. Knepley 
449e03fd340SMatthew G. Knepley   Input Parameters:
450e03fd340SMatthew G. Knepley + adaptor - the `DMAdaptor`
451e03fd340SMatthew G. Knepley - prefix  - the prefix to prepend to all option names
452e03fd340SMatthew G. Knepley 
453e03fd340SMatthew G. Knepley   Level: advanced
454e03fd340SMatthew G. Knepley 
455e03fd340SMatthew G. Knepley   Note:
456e03fd340SMatthew G. Knepley   A hyphen (-) must NOT be given at the beginning of the prefix name.
457e03fd340SMatthew G. Knepley   The first character of all runtime options is AUTOMATICALLY the hyphen.
458e03fd340SMatthew G. Knepley 
459e03fd340SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptor`, `SNESSetOptionsPrefix()`, `DMAdaptorSetFromOptions()`
460e03fd340SMatthew G. Knepley @*/
461e03fd340SMatthew G. Knepley PetscErrorCode DMAdaptorSetOptionsPrefix(DMAdaptor adaptor, const char prefix[])
462e03fd340SMatthew G. Knepley {
463e03fd340SMatthew G. Knepley   PetscFunctionBegin;
464e03fd340SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
465e03fd340SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor, prefix));
466e03fd340SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor->refineTag, prefix));
467e03fd340SMatthew G. Knepley   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)adaptor->refineTag, "refine_"));
468e03fd340SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor->coarsenTag, prefix));
469e03fd340SMatthew G. Knepley   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)adaptor->coarsenTag, "coarsen_"));
470e03fd340SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
471e03fd340SMatthew G. Knepley }
472e03fd340SMatthew G. Knepley 
473e03fd340SMatthew G. Knepley /*@
474c6f61ee2SBarry Smith   DMAdaptorSetFromOptions - Sets properties of a `DMAdaptor` object from values in the options database
475c6f61ee2SBarry Smith 
476c6f61ee2SBarry Smith   Collective
477c6f61ee2SBarry Smith 
478c6f61ee2SBarry Smith   Input Parameter:
479c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
480c6f61ee2SBarry Smith 
481c6f61ee2SBarry Smith   Options Database Keys:
4828b724c91SMatthew G. Knepley + -adaptor_monitor_size                - Monitor the mesh size
4838b724c91SMatthew G. Knepley . -adaptor_monitor_error               - Monitor the solution error
484c6f61ee2SBarry Smith . -adaptor_sequence_num <num>          - Number of adaptations to generate an optimal grid
485c6f61ee2SBarry Smith . -adaptor_target_num <num>            - Set the target number of vertices N_adapt, -1 for automatic determination
4863a336bb1SMatthew G. Knepley . -adaptor_refinement_factor <r>       - Set r such that N_adapt = r^dim N_orig
487d7c1f440SPierre Jolivet - -adaptor_mixed_setup_function <func> - Set the function func that sets up the mixed problem
488c6f61ee2SBarry Smith 
489c6f61ee2SBarry Smith   Level: beginner
490c6f61ee2SBarry Smith 
491df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
492c6f61ee2SBarry Smith @*/
493c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetFromOptions(DMAdaptor adaptor)
494c6f61ee2SBarry Smith {
4953a336bb1SMatthew G. Knepley   char                  typeName[PETSC_MAX_PATH_LEN];
4963a336bb1SMatthew G. Knepley   const char           *defName = DMADAPTORGRADIENT;
4973a336bb1SMatthew G. Knepley   char                  funcname[PETSC_MAX_PATH_LEN];
4988b724c91SMatthew G. Knepley   DMAdaptationCriterion criterion = DM_ADAPTATION_NONE;
4993a336bb1SMatthew G. Knepley   PetscBool             flg;
5003a336bb1SMatthew G. Knepley 
501c6f61ee2SBarry Smith   PetscFunctionBegin;
5023a336bb1SMatthew G. Knepley   PetscObjectOptionsBegin((PetscObject)adaptor);
5033a336bb1SMatthew G. Knepley   PetscCall(PetscOptionsFList("-adaptor_type", "DMAdaptor", "DMAdaptorSetType", DMAdaptorList, defName, typeName, 1024, &flg));
5043a336bb1SMatthew G. Knepley   if (flg) PetscCall(DMAdaptorSetType(adaptor, typeName));
5053a336bb1SMatthew G. Knepley   else if (!((PetscObject)adaptor)->type_name) PetscCall(DMAdaptorSetType(adaptor, defName));
5068b724c91SMatthew G. Knepley   PetscCall(PetscOptionsEnum("-adaptor_criterion", "Criterion used to drive adaptation", "", DMAdaptationCriteria, (PetscEnum)criterion, (PetscEnum *)&criterion, &flg));
5078b724c91SMatthew G. Knepley   if (flg) PetscCall(DMAdaptorSetCriterion(adaptor, criterion));
508c6f61ee2SBarry Smith   PetscCall(PetscOptionsInt("-adaptor_sequence_num", "Number of adaptations to generate an optimal grid", "DMAdaptorSetSequenceLength", adaptor->numSeq, &adaptor->numSeq, NULL));
509c6f61ee2SBarry Smith   PetscCall(PetscOptionsInt("-adaptor_target_num", "Set the target number of vertices N_adapt, -1 for automatic determination", "DMAdaptor", adaptor->Nadapt, &adaptor->Nadapt, NULL));
510c6f61ee2SBarry Smith   PetscCall(PetscOptionsReal("-adaptor_refinement_factor", "Set r such that N_adapt = r^dim N_orig", "DMAdaptor", adaptor->refinementFactor, &adaptor->refinementFactor, NULL));
5113a336bb1SMatthew G. Knepley   PetscCall(PetscOptionsString("-adaptor_mixed_setup_function", "Function to setup the mixed problem", "DMAdaptorSetMixedSetupFunction", funcname, funcname, sizeof(funcname), &flg));
5123a336bb1SMatthew G. Knepley   if (flg) {
5133a336bb1SMatthew G. Knepley     PetscErrorCode (*setupFunc)(DMAdaptor, DM);
5143a336bb1SMatthew G. Knepley 
5153a336bb1SMatthew G. Knepley     PetscCall(PetscDLSym(NULL, funcname, (void **)&setupFunc));
5163a336bb1SMatthew G. Knepley     PetscCheck(setupFunc, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
5173a336bb1SMatthew G. Knepley     PetscCall(DMAdaptorSetMixedSetupFunction(adaptor, setupFunc));
5183a336bb1SMatthew G. Knepley   }
5198b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_size", "size", adaptor));
5208b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_error", "error", adaptor));
521c6f61ee2SBarry Smith   PetscOptionsEnd();
522c6f61ee2SBarry Smith   PetscCall(VecTaggerSetFromOptions(adaptor->refineTag));
523c6f61ee2SBarry Smith   PetscCall(VecTaggerSetFromOptions(adaptor->coarsenTag));
524c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
525c6f61ee2SBarry Smith }
526c6f61ee2SBarry Smith 
527c6f61ee2SBarry Smith /*@
528c6f61ee2SBarry Smith   DMAdaptorView - Views a `DMAdaptor` object
529c6f61ee2SBarry Smith 
530c6f61ee2SBarry Smith   Collective
531c6f61ee2SBarry Smith 
532c6f61ee2SBarry Smith   Input Parameters:
533c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
534c6f61ee2SBarry Smith - viewer  - The `PetscViewer` object
535c6f61ee2SBarry Smith 
536c6f61ee2SBarry Smith   Level: beginner
537c6f61ee2SBarry Smith 
538df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
539c6f61ee2SBarry Smith @*/
540c6f61ee2SBarry Smith PetscErrorCode DMAdaptorView(DMAdaptor adaptor, PetscViewer viewer)
541c6f61ee2SBarry Smith {
542c6f61ee2SBarry Smith   PetscFunctionBegin;
543c6f61ee2SBarry Smith   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adaptor, viewer));
544c6f61ee2SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "DM Adaptor\n"));
545c6f61ee2SBarry Smith   PetscCall(PetscViewerASCIIPrintf(viewer, "  sequence length: %" PetscInt_FMT "\n", adaptor->numSeq));
546c6f61ee2SBarry Smith   PetscCall(VecTaggerView(adaptor->refineTag, viewer));
547c6f61ee2SBarry Smith   PetscCall(VecTaggerView(adaptor->coarsenTag, viewer));
548c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
549c6f61ee2SBarry Smith }
550c6f61ee2SBarry Smith 
551c6f61ee2SBarry Smith /*@
552c6f61ee2SBarry Smith   DMAdaptorGetSolver - Gets the solver used to produce discrete solutions
553c6f61ee2SBarry Smith 
554c6f61ee2SBarry Smith   Not Collective
555c6f61ee2SBarry Smith 
556c6f61ee2SBarry Smith   Input Parameter:
557c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
558c6f61ee2SBarry Smith 
559c6f61ee2SBarry Smith   Output Parameter:
560c6f61ee2SBarry Smith . snes - The solver
561c6f61ee2SBarry Smith 
562c6f61ee2SBarry Smith   Level: intermediate
563c6f61ee2SBarry Smith 
564df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
565c6f61ee2SBarry Smith @*/
566c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetSolver(DMAdaptor adaptor, SNES *snes)
567c6f61ee2SBarry Smith {
568c6f61ee2SBarry Smith   PetscFunctionBegin;
5693a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
570c6f61ee2SBarry Smith   PetscAssertPointer(snes, 2);
571c6f61ee2SBarry Smith   *snes = adaptor->snes;
572c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
573c6f61ee2SBarry Smith }
574c6f61ee2SBarry Smith 
575c6f61ee2SBarry Smith /*@
576c6f61ee2SBarry Smith   DMAdaptorSetSolver - Sets the solver used to produce discrete solutions
577c6f61ee2SBarry Smith 
578c6f61ee2SBarry Smith   Not Collective
579c6f61ee2SBarry Smith 
580c6f61ee2SBarry Smith   Input Parameters:
581c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
582c6f61ee2SBarry Smith - snes    - The solver, this MUST have an attached `DM`/`PetscDS`, so that the exact solution can be computed
583c6f61ee2SBarry Smith 
584c6f61ee2SBarry Smith   Level: intermediate
585c6f61ee2SBarry Smith 
586df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorGetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
587c6f61ee2SBarry Smith @*/
588c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetSolver(DMAdaptor adaptor, SNES snes)
589c6f61ee2SBarry Smith {
590c6f61ee2SBarry Smith   PetscFunctionBegin;
5913a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
592c6f61ee2SBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
593c6f61ee2SBarry Smith   adaptor->snes = snes;
594c6f61ee2SBarry Smith   PetscCall(SNESGetDM(adaptor->snes, &adaptor->idm));
595c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
596c6f61ee2SBarry Smith }
597c6f61ee2SBarry Smith 
598c6f61ee2SBarry Smith /*@
599c6f61ee2SBarry Smith   DMAdaptorGetSequenceLength - Gets the number of sequential adaptations used by an adapter
600c6f61ee2SBarry Smith 
601c6f61ee2SBarry Smith   Not Collective
602c6f61ee2SBarry Smith 
603c6f61ee2SBarry Smith   Input Parameter:
604c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
605c6f61ee2SBarry Smith 
606c6f61ee2SBarry Smith   Output Parameter:
607c6f61ee2SBarry Smith . num - The number of adaptations
608c6f61ee2SBarry Smith 
609c6f61ee2SBarry Smith   Level: intermediate
610c6f61ee2SBarry Smith 
611df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorSetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
612c6f61ee2SBarry Smith @*/
613c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetSequenceLength(DMAdaptor adaptor, PetscInt *num)
614c6f61ee2SBarry Smith {
615c6f61ee2SBarry Smith   PetscFunctionBegin;
6163a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
617c6f61ee2SBarry Smith   PetscAssertPointer(num, 2);
618c6f61ee2SBarry Smith   *num = adaptor->numSeq;
619c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
620c6f61ee2SBarry Smith }
621c6f61ee2SBarry Smith 
622c6f61ee2SBarry Smith /*@
623c6f61ee2SBarry Smith   DMAdaptorSetSequenceLength - Sets the number of sequential adaptations
624c6f61ee2SBarry Smith 
625c6f61ee2SBarry Smith   Not Collective
626c6f61ee2SBarry Smith 
627c6f61ee2SBarry Smith   Input Parameters:
628c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
629c6f61ee2SBarry Smith - num     - The number of adaptations
630c6f61ee2SBarry Smith 
631c6f61ee2SBarry Smith   Level: intermediate
632c6f61ee2SBarry Smith 
633df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptorGetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
634c6f61ee2SBarry Smith @*/
635c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetSequenceLength(DMAdaptor adaptor, PetscInt num)
636c6f61ee2SBarry Smith {
637c6f61ee2SBarry Smith   PetscFunctionBegin;
6383a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
639c6f61ee2SBarry Smith   adaptor->numSeq = num;
640c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
641c6f61ee2SBarry Smith }
642c6f61ee2SBarry Smith 
6433a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorTransferSolution_Exact_Private(DMAdaptor adaptor, DM dm, Vec u, DM adm, Vec au, void *ctx)
6443a336bb1SMatthew G. Knepley {
6453a336bb1SMatthew G. Knepley   PetscFunctionBeginUser;
6463a336bb1SMatthew G. Knepley   PetscCall(DMProjectFunction(adm, 0.0, adaptor->exactSol, adaptor->exactCtx, INSERT_ALL_VALUES, au));
6473a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
6483a336bb1SMatthew G. Knepley }
6493a336bb1SMatthew G. Knepley 
650c6f61ee2SBarry Smith /*@
651c6f61ee2SBarry Smith   DMAdaptorSetUp - After the solver is specified, creates data structures for controlling adaptivity
652c6f61ee2SBarry Smith 
653c6f61ee2SBarry Smith   Collective
654c6f61ee2SBarry Smith 
655c6f61ee2SBarry Smith   Input Parameter:
656c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
657c6f61ee2SBarry Smith 
658c6f61ee2SBarry Smith   Level: beginner
659c6f61ee2SBarry Smith 
660df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
661c6f61ee2SBarry Smith @*/
662c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetUp(DMAdaptor adaptor)
663c6f61ee2SBarry Smith {
664c6f61ee2SBarry Smith   PetscDS  prob;
665c6f61ee2SBarry Smith   PetscInt Nf, f;
666c6f61ee2SBarry Smith 
667c6f61ee2SBarry Smith   PetscFunctionBegin;
668c6f61ee2SBarry Smith   PetscCall(DMGetDS(adaptor->idm, &prob));
669c6f61ee2SBarry Smith   PetscCall(VecTaggerSetUp(adaptor->refineTag));
670c6f61ee2SBarry Smith   PetscCall(VecTaggerSetUp(adaptor->coarsenTag));
671c6f61ee2SBarry Smith   PetscCall(PetscDSGetNumFields(prob, &Nf));
672c6f61ee2SBarry Smith   PetscCall(PetscMalloc2(Nf, &adaptor->exactSol, Nf, &adaptor->exactCtx));
673c6f61ee2SBarry Smith   for (f = 0; f < Nf; ++f) {
674c6f61ee2SBarry Smith     PetscCall(PetscDSGetExactSolution(prob, f, &adaptor->exactSol[f], &adaptor->exactCtx[f]));
675c6f61ee2SBarry Smith     /* TODO Have a flag that forces projection rather than using the exact solution */
676c6f61ee2SBarry Smith     if (adaptor->exactSol[0]) PetscCall(DMAdaptorSetTransferFunction(adaptor, DMAdaptorTransferSolution_Exact_Private));
677c6f61ee2SBarry Smith   }
6783a336bb1SMatthew G. Knepley   PetscTryTypeMethod(adaptor, setup);
679c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
680c6f61ee2SBarry Smith }
681c6f61ee2SBarry Smith 
682c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetTransferFunction(DMAdaptor adaptor, PetscErrorCode (**tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
683c6f61ee2SBarry Smith {
684c6f61ee2SBarry Smith   PetscFunctionBegin;
685c6f61ee2SBarry Smith   *tfunc = adaptor->ops->transfersolution;
686c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
687c6f61ee2SBarry Smith }
688c6f61ee2SBarry Smith 
689c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetTransferFunction(DMAdaptor adaptor, PetscErrorCode (*tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
690c6f61ee2SBarry Smith {
691c6f61ee2SBarry Smith   PetscFunctionBegin;
692c6f61ee2SBarry Smith   adaptor->ops->transfersolution = tfunc;
693c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
694c6f61ee2SBarry Smith }
695c6f61ee2SBarry Smith 
696c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorPreAdapt(DMAdaptor adaptor, Vec locX)
697c6f61ee2SBarry Smith {
698c6f61ee2SBarry Smith   DM           plex;
699c6f61ee2SBarry Smith   PetscDS      prob;
700c6f61ee2SBarry Smith   PetscObject  obj;
701c6f61ee2SBarry Smith   PetscClassId id;
702c6f61ee2SBarry Smith   PetscBool    isForest;
703c6f61ee2SBarry Smith 
704c6f61ee2SBarry Smith   PetscFunctionBegin;
705c6f61ee2SBarry Smith   PetscCall(DMConvert(adaptor->idm, DMPLEX, &plex));
706c6f61ee2SBarry Smith   PetscCall(DMGetDS(adaptor->idm, &prob));
707c6f61ee2SBarry Smith   PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
708c6f61ee2SBarry Smith   PetscCall(PetscObjectGetClassId(obj, &id));
709c6f61ee2SBarry Smith   PetscCall(DMIsForest(adaptor->idm, &isForest));
710c6f61ee2SBarry Smith   if (adaptor->adaptCriterion == DM_ADAPTATION_NONE) {
711c6f61ee2SBarry Smith     if (isForest) {
712c6f61ee2SBarry Smith       adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
713c6f61ee2SBarry Smith     }
714c6f61ee2SBarry Smith #if defined(PETSC_HAVE_PRAGMATIC)
715c6f61ee2SBarry Smith     else {
716c6f61ee2SBarry Smith       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
717c6f61ee2SBarry Smith     }
718c6f61ee2SBarry Smith #elif defined(PETSC_HAVE_MMG)
719c6f61ee2SBarry Smith     else {
720c6f61ee2SBarry Smith       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
721c6f61ee2SBarry Smith     }
722c6f61ee2SBarry Smith #elif defined(PETSC_HAVE_PARMMG)
723c6f61ee2SBarry Smith     else {
724c6f61ee2SBarry Smith       adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
725c6f61ee2SBarry Smith     }
726c6f61ee2SBarry Smith #else
727c6f61ee2SBarry Smith     else {
7283a336bb1SMatthew G. Knepley       adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
729c6f61ee2SBarry Smith     }
730c6f61ee2SBarry Smith #endif
731c6f61ee2SBarry Smith   }
732c6f61ee2SBarry Smith   if (id == PETSCFV_CLASSID) {
733c6f61ee2SBarry Smith     adaptor->femType = PETSC_FALSE;
734c6f61ee2SBarry Smith   } else {
735c6f61ee2SBarry Smith     adaptor->femType = PETSC_TRUE;
736c6f61ee2SBarry Smith   }
737c6f61ee2SBarry Smith   if (adaptor->femType) {
738c6f61ee2SBarry Smith     /* Compute local solution bc */
739c6f61ee2SBarry Smith     PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
740c6f61ee2SBarry Smith   } else {
741c6f61ee2SBarry Smith     PetscFV      fvm = (PetscFV)obj;
742c6f61ee2SBarry Smith     PetscLimiter noneLimiter;
743c6f61ee2SBarry Smith     Vec          grad;
744c6f61ee2SBarry Smith 
745c6f61ee2SBarry Smith     PetscCall(PetscFVGetComputeGradients(fvm, &adaptor->computeGradient));
746c6f61ee2SBarry Smith     PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
747c6f61ee2SBarry Smith     /* Use no limiting when reconstructing gradients for adaptivity */
748c6f61ee2SBarry Smith     PetscCall(PetscFVGetLimiter(fvm, &adaptor->limiter));
749c6f61ee2SBarry Smith     PetscCall(PetscObjectReference((PetscObject)adaptor->limiter));
750c6f61ee2SBarry Smith     PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
751c6f61ee2SBarry Smith     PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
752c6f61ee2SBarry Smith     PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
753c6f61ee2SBarry Smith     /* Get FVM data */
754c6f61ee2SBarry Smith     PetscCall(DMPlexGetDataFVM(plex, fvm, &adaptor->cellGeom, &adaptor->faceGeom, &adaptor->gradDM));
755c6f61ee2SBarry Smith     PetscCall(VecGetDM(adaptor->cellGeom, &adaptor->cellDM));
756c6f61ee2SBarry Smith     PetscCall(VecGetArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
757c6f61ee2SBarry Smith     /* Compute local solution bc */
758c6f61ee2SBarry Smith     PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
759c6f61ee2SBarry Smith     /* Compute gradients */
760c6f61ee2SBarry Smith     PetscCall(DMCreateGlobalVector(adaptor->gradDM, &grad));
761c6f61ee2SBarry Smith     PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
762c6f61ee2SBarry Smith     PetscCall(DMGetLocalVector(adaptor->gradDM, &adaptor->cellGrad));
763c6f61ee2SBarry Smith     PetscCall(DMGlobalToLocalBegin(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
764c6f61ee2SBarry Smith     PetscCall(DMGlobalToLocalEnd(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
765c6f61ee2SBarry Smith     PetscCall(VecDestroy(&grad));
766c6f61ee2SBarry Smith     PetscCall(VecGetArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
767c6f61ee2SBarry Smith   }
768c6f61ee2SBarry Smith   PetscCall(DMDestroy(&plex));
769c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
770c6f61ee2SBarry Smith }
771c6f61ee2SBarry Smith 
772c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorTransferSolution(DMAdaptor adaptor, DM dm, Vec x, DM adm, Vec ax)
773c6f61ee2SBarry Smith {
774c6f61ee2SBarry Smith   PetscReal time = 0.0;
775c6f61ee2SBarry Smith   Mat       interp;
776c6f61ee2SBarry Smith   void     *ctx;
777c6f61ee2SBarry Smith 
778c6f61ee2SBarry Smith   PetscFunctionBegin;
779c6f61ee2SBarry Smith   PetscCall(DMGetApplicationContext(dm, &ctx));
780c6f61ee2SBarry Smith   if (adaptor->ops->transfersolution) PetscUseTypeMethod(adaptor, transfersolution, dm, x, adm, ax, ctx);
781c6f61ee2SBarry Smith   else {
782c6f61ee2SBarry Smith     switch (adaptor->adaptCriterion) {
783c6f61ee2SBarry Smith     case DM_ADAPTATION_LABEL:
784c6f61ee2SBarry Smith       PetscCall(DMForestTransferVec(dm, x, adm, ax, PETSC_TRUE, time));
785c6f61ee2SBarry Smith       break;
786c6f61ee2SBarry Smith     case DM_ADAPTATION_REFINE:
787c6f61ee2SBarry Smith     case DM_ADAPTATION_METRIC:
788c6f61ee2SBarry Smith       PetscCall(DMCreateInterpolation(dm, adm, &interp, NULL));
789c6f61ee2SBarry Smith       PetscCall(MatInterpolate(interp, x, ax));
790c6f61ee2SBarry Smith       PetscCall(DMInterpolate(dm, interp, adm));
791c6f61ee2SBarry Smith       PetscCall(MatDestroy(&interp));
792c6f61ee2SBarry Smith       break;
793c6f61ee2SBarry Smith     default:
794c6f61ee2SBarry Smith       SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_SUP, "No built-in projection for this adaptation criterion: %d", adaptor->adaptCriterion);
795c6f61ee2SBarry Smith     }
796c6f61ee2SBarry Smith   }
797c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
798c6f61ee2SBarry Smith }
799c6f61ee2SBarry Smith 
800c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorPostAdapt(DMAdaptor adaptor)
801c6f61ee2SBarry Smith {
802c6f61ee2SBarry Smith   PetscDS      prob;
803c6f61ee2SBarry Smith   PetscObject  obj;
804c6f61ee2SBarry Smith   PetscClassId id;
805c6f61ee2SBarry Smith 
806c6f61ee2SBarry Smith   PetscFunctionBegin;
807c6f61ee2SBarry Smith   PetscCall(DMGetDS(adaptor->idm, &prob));
808c6f61ee2SBarry Smith   PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
809c6f61ee2SBarry Smith   PetscCall(PetscObjectGetClassId(obj, &id));
810c6f61ee2SBarry Smith   if (id == PETSCFV_CLASSID) {
811c6f61ee2SBarry Smith     PetscFV fvm = (PetscFV)obj;
812c6f61ee2SBarry Smith 
813c6f61ee2SBarry Smith     PetscCall(PetscFVSetComputeGradients(fvm, adaptor->computeGradient));
814c6f61ee2SBarry Smith     /* Restore original limiter */
815c6f61ee2SBarry Smith     PetscCall(PetscFVSetLimiter(fvm, adaptor->limiter));
816c6f61ee2SBarry Smith 
817c6f61ee2SBarry Smith     PetscCall(VecRestoreArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
818c6f61ee2SBarry Smith     PetscCall(VecRestoreArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
819c6f61ee2SBarry Smith     PetscCall(DMRestoreLocalVector(adaptor->gradDM, &adaptor->cellGrad));
820c6f61ee2SBarry Smith   }
821c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
822c6f61ee2SBarry Smith }
823c6f61ee2SBarry Smith 
824c6f61ee2SBarry Smith /*
8253a336bb1SMatthew G. Knepley   DMAdaptorComputeCellErrorIndicator_Gradient - Use the integrated gradient as an error indicator in the `DMAdaptor`
826c6f61ee2SBarry Smith 
827c6f61ee2SBarry Smith   Input Parameters:
828c6f61ee2SBarry Smith + adaptor  - The `DMAdaptor` object
829c6f61ee2SBarry Smith . dim      - The topological dimension
830c6f61ee2SBarry Smith . cell     - The cell
831c6f61ee2SBarry Smith . field    - The field integrated over the cell
832c6f61ee2SBarry Smith . gradient - The gradient integrated over the cell
833c6f61ee2SBarry Smith . cg       - A `PetscFVCellGeom` struct
834c6f61ee2SBarry Smith - ctx      - A user context
835c6f61ee2SBarry Smith 
836c6f61ee2SBarry Smith   Output Parameter:
837c6f61ee2SBarry Smith . errInd   - The error indicator
838c6f61ee2SBarry Smith 
839c6f61ee2SBarry Smith   Developer Note:
840c6f61ee2SBarry Smith   Some of the input arguments are absurdly specialized to special situations, it is not clear this is a good general API
841c6f61ee2SBarry Smith 
8423a336bb1SMatthew G. Knepley .seealso: [](ch_dmbase), `DMAdaptor`
843c6f61ee2SBarry Smith */
8443a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorComputeCellErrorIndicator_Gradient(DMAdaptor adaptor, PetscInt dim, PetscInt Nc, const PetscScalar *field, const PetscScalar *gradient, const PetscFVCellGeom *cg, PetscReal *errInd, void *ctx)
845c6f61ee2SBarry Smith {
846c6f61ee2SBarry Smith   PetscReal err = 0.;
847c6f61ee2SBarry Smith   PetscInt  c, d;
848c6f61ee2SBarry Smith 
849c6f61ee2SBarry Smith   PetscFunctionBeginHot;
850c6f61ee2SBarry Smith   for (c = 0; c < Nc; c++) {
851c6f61ee2SBarry Smith     for (d = 0; d < dim; ++d) err += PetscSqr(PetscRealPart(gradient[c * dim + d]));
852c6f61ee2SBarry Smith   }
853c6f61ee2SBarry Smith   *errInd = cg->volume * err;
854c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
855c6f61ee2SBarry Smith }
856c6f61ee2SBarry Smith 
8573a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorComputeErrorIndicator_Gradient(DMAdaptor adaptor, Vec locX, Vec errVec)
858c6f61ee2SBarry Smith {
8593a336bb1SMatthew G. Knepley   DM              dm, plex, edm, eplex;
8603a336bb1SMatthew G. Knepley   PetscDS         ds;
861c6f61ee2SBarry Smith   PetscObject     obj;
862c6f61ee2SBarry Smith   PetscClassId    id;
863c6f61ee2SBarry Smith   void           *ctx;
864c6f61ee2SBarry Smith   PetscQuadrature quad;
8653a336bb1SMatthew G. Knepley   PetscScalar    *earray;
8663a336bb1SMatthew G. Knepley   PetscReal       minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2];
8673a336bb1SMatthew G. Knepley   PetscInt        dim, cdim, cStart, cEnd, Nf, Nc;
868c6f61ee2SBarry Smith 
869c6f61ee2SBarry Smith   PetscFunctionBegin;
8703a336bb1SMatthew G. Knepley   PetscCall(VecGetDM(locX, &dm));
8713a336bb1SMatthew G. Knepley   PetscCall(DMConvert(dm, DMPLEX, &plex));
8723a336bb1SMatthew G. Knepley   PetscCall(VecGetDM(errVec, &edm));
8733a336bb1SMatthew G. Knepley   PetscCall(DMConvert(edm, DMPLEX, &eplex));
874c6f61ee2SBarry Smith   PetscCall(DMGetDimension(plex, &dim));
875c6f61ee2SBarry Smith   PetscCall(DMGetCoordinateDim(plex, &cdim));
876c6f61ee2SBarry Smith   PetscCall(DMGetApplicationContext(plex, &ctx));
8773a336bb1SMatthew G. Knepley   PetscCall(DMGetDS(plex, &ds));
8783a336bb1SMatthew G. Knepley   PetscCall(PetscDSGetNumFields(ds, &Nf));
8793a336bb1SMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(ds, 0, &obj));
880c6f61ee2SBarry Smith   PetscCall(PetscObjectGetClassId(obj, &id));
8813a336bb1SMatthew G. Knepley 
8823a336bb1SMatthew G. Knepley   PetscCall(VecGetArray(errVec, &earray));
8833a336bb1SMatthew G. Knepley   PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
8843a336bb1SMatthew G. Knepley   for (PetscInt cell = cStart; cell < cEnd; ++cell) {
8853a336bb1SMatthew G. Knepley     PetscScalar *eval;
8863a336bb1SMatthew G. Knepley     PetscReal    errInd = 0.;
8873a336bb1SMatthew G. Knepley 
888c6f61ee2SBarry Smith     if (id == PETSCFV_CLASSID) {
8893a336bb1SMatthew G. Knepley       PetscFV            fv = (PetscFV)obj;
890c6f61ee2SBarry Smith       const PetscScalar *pointSols;
891c6f61ee2SBarry Smith       const PetscScalar *pointSol;
892c6f61ee2SBarry Smith       const PetscScalar *pointGrad;
893c6f61ee2SBarry Smith       PetscFVCellGeom   *cg;
894c6f61ee2SBarry Smith 
8953a336bb1SMatthew G. Knepley       PetscCall(PetscFVGetNumComponents(fv, &Nc));
896c6f61ee2SBarry Smith       PetscCall(VecGetArrayRead(locX, &pointSols));
897c6f61ee2SBarry Smith       PetscCall(DMPlexPointLocalRead(plex, cell, pointSols, (void *)&pointSol));
898c6f61ee2SBarry Smith       PetscCall(DMPlexPointLocalRead(adaptor->gradDM, cell, adaptor->cellGradArray, (void *)&pointGrad));
899c6f61ee2SBarry Smith       PetscCall(DMPlexPointLocalRead(adaptor->cellDM, cell, adaptor->cellGeomArray, &cg));
9003a336bb1SMatthew G. Knepley       PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, pointSol, pointGrad, cg, &errInd, ctx);
901c6f61ee2SBarry Smith       PetscCall(VecRestoreArrayRead(locX, &pointSols));
902c6f61ee2SBarry Smith     } else {
9033a336bb1SMatthew G. Knepley       PetscFE          fe = (PetscFE)obj;
904c6f61ee2SBarry Smith       PetscScalar     *x  = NULL, *field, *gradient, *interpolant, *interpolantGrad;
905c6f61ee2SBarry Smith       PetscFVCellGeom  cg;
906c6f61ee2SBarry Smith       PetscFEGeom      fegeom;
907c6f61ee2SBarry Smith       const PetscReal *quadWeights;
908c6f61ee2SBarry Smith       PetscReal       *coords;
9093a336bb1SMatthew G. Knepley       PetscInt         Nb, Nq, qNc;
910c6f61ee2SBarry Smith 
911c6f61ee2SBarry Smith       fegeom.dim      = dim;
912c6f61ee2SBarry Smith       fegeom.dimEmbed = cdim;
9133a336bb1SMatthew G. Knepley       PetscCall(PetscFEGetNumComponents(fe, &Nc));
9143a336bb1SMatthew G. Knepley       PetscCall(PetscFEGetQuadrature(fe, &quad));
9153a336bb1SMatthew G. Knepley       PetscCall(PetscFEGetDimension(fe, &Nb));
916c6f61ee2SBarry Smith       PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
917c6f61ee2SBarry Smith       PetscCall(PetscCalloc6(Nc, &field, cdim * Nc, &gradient, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ));
918c6f61ee2SBarry Smith       PetscCall(PetscMalloc2(Nc, &interpolant, cdim * Nc, &interpolantGrad));
919c6f61ee2SBarry Smith       PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
920c6f61ee2SBarry Smith       PetscCall(DMPlexComputeCellGeometryFVM(plex, cell, &cg.volume, NULL, NULL));
921c6f61ee2SBarry Smith       PetscCall(PetscArrayzero(gradient, cdim * Nc));
9223a336bb1SMatthew G. Knepley       PetscCall(DMPlexVecGetClosure(plex, NULL, locX, cell, NULL, &x));
9233a336bb1SMatthew G. Knepley       for (PetscInt f = 0; f < Nf; ++f) {
9243a336bb1SMatthew G. Knepley         PetscInt qc = 0;
925c6f61ee2SBarry Smith 
9263a336bb1SMatthew G. Knepley         PetscCall(PetscDSGetDiscretization(ds, f, &obj));
927c6f61ee2SBarry Smith         PetscCall(PetscArrayzero(interpolant, Nc));
928c6f61ee2SBarry Smith         PetscCall(PetscArrayzero(interpolantGrad, cdim * Nc));
9293a336bb1SMatthew G. Knepley         for (PetscInt q = 0; q < Nq; ++q) {
930c6f61ee2SBarry Smith           PetscCall(PetscFEInterpolateFieldAndGradient_Static((PetscFE)obj, 1, x, &fegeom, q, interpolant, interpolantGrad));
9313a336bb1SMatthew G. Knepley           for (PetscInt fc = 0; fc < Nc; ++fc) {
932c6f61ee2SBarry Smith             const PetscReal wt = quadWeights[q * qNc + qc + fc];
933c6f61ee2SBarry Smith 
934c6f61ee2SBarry Smith             field[fc] += interpolant[fc] * wt * fegeom.detJ[q];
9353a336bb1SMatthew G. Knepley             for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] += interpolantGrad[fc * dim + d] * wt * fegeom.detJ[q];
936c6f61ee2SBarry Smith           }
937c6f61ee2SBarry Smith         }
938c6f61ee2SBarry Smith         qc += Nc;
939c6f61ee2SBarry Smith       }
940c6f61ee2SBarry Smith       PetscCall(PetscFree2(interpolant, interpolantGrad));
941c6f61ee2SBarry Smith       PetscCall(DMPlexVecRestoreClosure(plex, NULL, locX, cell, NULL, &x));
9423a336bb1SMatthew G. Knepley       for (PetscInt fc = 0; fc < Nc; ++fc) {
943c6f61ee2SBarry Smith         field[fc] /= cg.volume;
9443a336bb1SMatthew G. Knepley         for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] /= cg.volume;
945c6f61ee2SBarry Smith       }
9463a336bb1SMatthew G. Knepley       PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, field, gradient, &cg, &errInd, ctx);
947c6f61ee2SBarry Smith       PetscCall(PetscFree6(field, gradient, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
948c6f61ee2SBarry Smith     }
9493a336bb1SMatthew G. Knepley     PetscCall(DMPlexPointGlobalRef(eplex, cell, earray, (void *)&eval));
9503a336bb1SMatthew G. Knepley     eval[0]      = errInd;
9513a336bb1SMatthew G. Knepley     minMaxInd[0] = PetscMin(minMaxInd[0], errInd);
9523a336bb1SMatthew G. Knepley     minMaxInd[1] = PetscMax(minMaxInd[1], errInd);
9533a336bb1SMatthew G. Knepley   }
9543a336bb1SMatthew G. Knepley   PetscCall(VecRestoreArray(errVec, &earray));
9553a336bb1SMatthew G. Knepley   PetscCall(DMDestroy(&plex));
9563a336bb1SMatthew G. Knepley   PetscCall(DMDestroy(&eplex));
9573a336bb1SMatthew G. Knepley   PetscCall(PetscGlobalMinMaxReal(PetscObjectComm((PetscObject)adaptor), minMaxInd, minMaxIndGlobal));
9583a336bb1SMatthew G. Knepley   PetscCall(PetscInfo(adaptor, "DMAdaptor: error indicator range (%g, %g)\n", (double)minMaxIndGlobal[0], (double)minMaxIndGlobal[1]));
9593a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
9603a336bb1SMatthew G. Knepley }
9613a336bb1SMatthew G. Knepley 
9623a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorComputeErrorIndicator_Flux(DMAdaptor adaptor, Vec lu, Vec errVec)
9633a336bb1SMatthew G. Knepley {
9643a336bb1SMatthew G. Knepley   DM          dm, mdm;
9653a336bb1SMatthew G. Knepley   SNES        msnes;
9663a336bb1SMatthew G. Knepley   Vec         mu, lmu;
9673a336bb1SMatthew G. Knepley   void       *ctx;
968e03fd340SMatthew G. Knepley   const char *prefix;
9693a336bb1SMatthew G. Knepley 
9703a336bb1SMatthew G. Knepley   PetscFunctionBegin;
9713a336bb1SMatthew G. Knepley   PetscCall(VecGetDM(lu, &dm));
9723a336bb1SMatthew G. Knepley 
9733a336bb1SMatthew G. Knepley   // Set up and solve mixed problem
9743a336bb1SMatthew G. Knepley   PetscCall(DMClone(dm, &mdm));
9753a336bb1SMatthew G. Knepley   PetscCall(SNESCreate(PetscObjectComm((PetscObject)mdm), &msnes));
9763a336bb1SMatthew G. Knepley   PetscCall(SNESSetDM(msnes, mdm));
977e03fd340SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
978e03fd340SMatthew G. Knepley   PetscCall(SNESSetOptionsPrefix(msnes, prefix));
9793a336bb1SMatthew G. Knepley   PetscCall(SNESAppendOptionsPrefix(msnes, "mixed_"));
9803a336bb1SMatthew G. Knepley 
9813a336bb1SMatthew G. Knepley   PetscTryTypeMethod(adaptor, mixedsetup, mdm);
9823a336bb1SMatthew G. Knepley   PetscCall(DMGetApplicationContext(dm, &ctx));
9833a336bb1SMatthew G. Knepley   PetscCall(DMPlexSetSNESLocalFEM(mdm, PETSC_FALSE, ctx));
9843a336bb1SMatthew G. Knepley   PetscCall(SNESSetFromOptions(msnes));
9853a336bb1SMatthew G. Knepley 
9863a336bb1SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(mdm, &mu));
9873a336bb1SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)mu, "Mixed Solution"));
9883a336bb1SMatthew G. Knepley   PetscCall(VecSet(mu, 0.0));
9893a336bb1SMatthew G. Knepley   PetscCall(SNESSolve(msnes, NULL, mu));
9903a336bb1SMatthew G. Knepley   PetscCall(VecViewFromOptions(mu, (PetscObject)adaptor, "-adapt_mixed_sol_vec_view"));
9913a336bb1SMatthew G. Knepley 
9923a336bb1SMatthew G. Knepley   PetscCall(DMGetLocalVector(mdm, &lmu));
9933a336bb1SMatthew G. Knepley   PetscCall(DMGlobalToLocal(mdm, mu, INSERT_VALUES, lmu));
9943a336bb1SMatthew G. Knepley   PetscCall(DMPlexInsertBoundaryValues(mdm, PETSC_TRUE, lmu, 0.0, NULL, NULL, NULL));
9953a336bb1SMatthew G. Knepley   PetscCall(DMPlexComputeL2FluxDiffVecLocal(lu, 0, lmu, 0, errVec));
9963a336bb1SMatthew G. Knepley   PetscCall(DMRestoreLocalVector(mdm, &lmu));
9973a336bb1SMatthew G. Knepley   PetscCall(VecDestroy(&mu));
9983a336bb1SMatthew G. Knepley   PetscCall(SNESDestroy(&msnes));
9993a336bb1SMatthew G. Knepley   PetscCall(DMDestroy(&mdm));
1000c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1001c6f61ee2SBarry Smith }
1002c6f61ee2SBarry Smith 
10038b724c91SMatthew G. Knepley /*@
10048b724c91SMatthew G. Knepley   DMAdaptorMonitor - runs the user provided monitor routines, if they exist
10058b724c91SMatthew G. Knepley 
10068b724c91SMatthew G. Knepley   Collective
10078b724c91SMatthew G. Knepley 
10088b724c91SMatthew G. Knepley   Input Parameters:
10098b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
10108b724c91SMatthew G. Knepley . it      - iteration number
10118b724c91SMatthew G. Knepley . odm     - the original `DM`
10128b724c91SMatthew G. Knepley . adm     - the adapted `DM`
10138b724c91SMatthew G. Knepley . Nf      - the number of fields
10148b724c91SMatthew G. Knepley . enorms  - the 2-norm error values for each field
10158b724c91SMatthew G. Knepley - error   - `Vec` of cellwise errors
10168b724c91SMatthew G. Knepley 
10178b724c91SMatthew G. Knepley   Level: developer
10188b724c91SMatthew G. Knepley 
10198b724c91SMatthew G. Knepley   Note:
10208b724c91SMatthew G. Knepley   This routine is called by the `DMAdaptor` implementations.
10218b724c91SMatthew G. Knepley   It does not typically need to be called by the user.
10228b724c91SMatthew G. Knepley 
10238b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptorMonitorSet()`
10248b724c91SMatthew G. Knepley @*/
10258b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitor(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error)
10268b724c91SMatthew G. Knepley {
10278b724c91SMatthew G. Knepley   PetscFunctionBegin;
10288b724c91SMatthew G. Knepley   for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) PetscCall((*adaptor->monitor[i])(adaptor, it, odm, adm, Nf, enorms, error, adaptor->monitorcontext[i]));
10298b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
10308b724c91SMatthew G. Knepley }
10318b724c91SMatthew G. Knepley 
10328b724c91SMatthew G. Knepley /*@C
10338b724c91SMatthew G. Knepley   DMAdaptorMonitorSize - Prints the mesh sizes at each iteration of an adaptation loop.
10348b724c91SMatthew G. Knepley 
10358b724c91SMatthew G. Knepley   Collective
10368b724c91SMatthew G. Knepley 
10378b724c91SMatthew G. Knepley   Input Parameters:
10388b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
10398b724c91SMatthew G. Knepley . n       - iteration number
10408b724c91SMatthew G. Knepley . odm     - the original `DM`
10418b724c91SMatthew G. Knepley . adm     - the adapted `DM`
10428b724c91SMatthew G. Knepley . Nf      - number of fields
10438b724c91SMatthew G. Knepley . enorms  - 2-norm error values for each field (may be estimated).
10448b724c91SMatthew G. Knepley . error   - `Vec` of cellwise errors
10458b724c91SMatthew G. Knepley - vf      - The viewer context
10468b724c91SMatthew G. Knepley 
10478b724c91SMatthew G. Knepley   Options Database Key:
10488b724c91SMatthew G. Knepley . -adaptor_monitor_size - Activates `DMAdaptorMonitorSize()`
10498b724c91SMatthew G. Knepley 
10508b724c91SMatthew G. Knepley   Level: intermediate
10518b724c91SMatthew G. Knepley 
10528b724c91SMatthew G. Knepley   Note:
10538b724c91SMatthew G. Knepley   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
10548b724c91SMatthew G. Knepley   to be used during the adaptation loop.
10558b724c91SMatthew G. Knepley 
10568b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorError()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
10578b724c91SMatthew G. Knepley @*/
10588b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorSize(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
10598b724c91SMatthew G. Knepley {
10608b724c91SMatthew G. Knepley   PetscViewer       viewer = vf->viewer;
10618b724c91SMatthew G. Knepley   PetscViewerFormat format = vf->format;
10628b724c91SMatthew G. Knepley   PetscInt          tablevel, cStart, cEnd, acStart, acEnd;
10638b724c91SMatthew G. Knepley   const char       *prefix;
10648b724c91SMatthew G. Knepley   PetscMPIInt       rank;
10658b724c91SMatthew G. Knepley 
10668b724c91SMatthew G. Knepley   PetscFunctionBegin;
10678b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
10688b724c91SMatthew G. Knepley   PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
10698b724c91SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
10708b724c91SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)adaptor), &rank));
10718b724c91SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
10728b724c91SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
10738b724c91SMatthew G. Knepley 
10748b724c91SMatthew G. Knepley   PetscCall(PetscViewerPushFormat(viewer, format));
10758b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
10768b724c91SMatthew G. Knepley   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Sizes for %s adaptation.\n", prefix));
10778b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor rank %d N_orig: %" PetscInt_FMT " N_adapt: %" PetscInt_FMT "\n", n, rank, cEnd - cStart, acEnd - acStart));
10788b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
10798b724c91SMatthew G. Knepley   PetscCall(PetscViewerPopFormat(viewer));
10808b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
10818b724c91SMatthew G. Knepley }
10828b724c91SMatthew G. Knepley 
10838b724c91SMatthew G. Knepley /*@C
10848b724c91SMatthew G. Knepley   DMAdaptorMonitorError - Prints the error norm at each iteration of an adaptation loop.
10858b724c91SMatthew G. Knepley 
10868b724c91SMatthew G. Knepley   Collective
10878b724c91SMatthew G. Knepley 
10888b724c91SMatthew G. Knepley   Input Parameters:
10898b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
10908b724c91SMatthew G. Knepley . n       - iteration number
10918b724c91SMatthew G. Knepley . odm     - the original `DM`
10928b724c91SMatthew G. Knepley . adm     - the adapted `DM`
10938b724c91SMatthew G. Knepley . Nf      - number of fields
10948b724c91SMatthew G. Knepley . enorms  - 2-norm error values for each field (may be estimated).
10958b724c91SMatthew G. Knepley . error   - `Vec` of cellwise errors
10968b724c91SMatthew G. Knepley - vf      - The viewer context
10978b724c91SMatthew G. Knepley 
10988b724c91SMatthew G. Knepley   Options Database Key:
10998b724c91SMatthew G. Knepley . -adaptor_monitor_error - Activates `DMAdaptorMonitorError()`
11008b724c91SMatthew G. Knepley 
11018b724c91SMatthew G. Knepley   Level: intermediate
11028b724c91SMatthew G. Knepley 
11038b724c91SMatthew G. Knepley   Note:
11048b724c91SMatthew G. Knepley   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
11058b724c91SMatthew G. Knepley   to be used during the adaptation loop.
11068b724c91SMatthew G. Knepley 
11078b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
11088b724c91SMatthew G. Knepley @*/
11098b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorError(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
11108b724c91SMatthew G. Knepley {
11118b724c91SMatthew G. Knepley   PetscViewer       viewer = vf->viewer;
11128b724c91SMatthew G. Knepley   PetscViewerFormat format = vf->format;
11138b724c91SMatthew G. Knepley   PetscInt          tablevel, cStart, cEnd, acStart, acEnd;
11148b724c91SMatthew G. Knepley   const char       *prefix;
11158b724c91SMatthew G. Knepley 
11168b724c91SMatthew G. Knepley   PetscFunctionBegin;
11178b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
11188b724c91SMatthew G. Knepley   PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
11198b724c91SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
11208b724c91SMatthew G. Knepley 
11218b724c91SMatthew G. Knepley   PetscCall(PetscViewerPushFormat(viewer, format));
11228b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
11238b724c91SMatthew G. Knepley   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Error norms for %s adaptation.\n", prefix));
11248b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor Error norm %s", n, Nf > 1 ? "[" : ""));
11258b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
11268b724c91SMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
11278b724c91SMatthew G. Knepley     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
11288b724c91SMatthew G. Knepley     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)enorms[f]));
11298b724c91SMatthew G. Knepley   }
11308b724c91SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
11318b724c91SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
11328b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, " N: %" PetscInt_FMT " Nadapt: %" PetscInt_FMT "\n", cEnd - cStart, acEnd - acStart));
11338b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
11348b724c91SMatthew G. Knepley   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
11358b724c91SMatthew G. Knepley   PetscCall(PetscViewerPopFormat(viewer));
11368b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
11378b724c91SMatthew G. Knepley }
11388b724c91SMatthew G. Knepley 
11398b724c91SMatthew G. Knepley /*@C
11408b724c91SMatthew G. Knepley   DMAdaptorMonitorErrorDraw - Plots the error at each iteration of an iterative solver.
11418b724c91SMatthew G. Knepley 
11428b724c91SMatthew G. Knepley   Collective
11438b724c91SMatthew G. Knepley 
11448b724c91SMatthew G. Knepley   Input Parameters:
11458b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
11468b724c91SMatthew G. Knepley . n       - iteration number
11478b724c91SMatthew G. Knepley . odm     - the original `DM`
11488b724c91SMatthew G. Knepley . adm     - the adapted `DM`
11498b724c91SMatthew G. Knepley . Nf      - number of fields
11508b724c91SMatthew G. Knepley . enorms  - 2-norm error values for each field (may be estimated).
11518b724c91SMatthew G. Knepley . error   - `Vec` of cellwise errors
11528b724c91SMatthew G. Knepley - vf      - The viewer context
11538b724c91SMatthew G. Knepley 
11548b724c91SMatthew G. Knepley   Options Database Key:
11558b724c91SMatthew G. Knepley . -adaptor_monitor_error draw - Activates `DMAdaptorMonitorErrorDraw()`
11568b724c91SMatthew G. Knepley 
11578b724c91SMatthew G. Knepley   Level: intermediate
11588b724c91SMatthew G. Knepley 
11598b724c91SMatthew G. Knepley   Note:
11608b724c91SMatthew G. Knepley   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
11618b724c91SMatthew G. Knepley   to be used during the adaptation loop.
11628b724c91SMatthew G. Knepley 
11638b724c91SMatthew G. Knepley .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
11648b724c91SMatthew G. Knepley @*/
11658b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorErrorDraw(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
11668b724c91SMatthew G. Knepley {
11678b724c91SMatthew G. Knepley   PetscViewer       viewer = vf->viewer;
11688b724c91SMatthew G. Knepley   PetscViewerFormat format = vf->format;
11698b724c91SMatthew G. Knepley 
11708b724c91SMatthew G. Knepley   PetscFunctionBegin;
11718b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
11728b724c91SMatthew G. Knepley   PetscCall(PetscViewerPushFormat(viewer, format));
11738b724c91SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
11748b724c91SMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", (PetscObject)adaptor));
11758b724c91SMatthew G. Knepley   PetscCall(VecView(error, viewer));
11768b724c91SMatthew G. Knepley   PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", NULL));
11778b724c91SMatthew G. Knepley   PetscCall(PetscViewerPopFormat(viewer));
11788b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
11798b724c91SMatthew G. Knepley }
11808b724c91SMatthew G. Knepley 
11818b724c91SMatthew G. Knepley /*@C
1182d7c1f440SPierre Jolivet   DMAdaptorMonitorErrorDrawLGCreate - Creates the context for the error plotter `DMAdaptorMonitorErrorDrawLG()`
11838b724c91SMatthew G. Knepley 
11848b724c91SMatthew G. Knepley   Collective
11858b724c91SMatthew G. Knepley 
11868b724c91SMatthew G. Knepley   Input Parameters:
11878b724c91SMatthew G. Knepley + viewer - The `PetscViewer`
11888b724c91SMatthew G. Knepley . format - The viewer format
11898b724c91SMatthew G. Knepley - ctx    - An optional user context
11908b724c91SMatthew G. Knepley 
11918b724c91SMatthew G. Knepley   Output Parameter:
11928b724c91SMatthew G. Knepley . vf - The viewer context
11938b724c91SMatthew G. Knepley 
11948b724c91SMatthew G. Knepley   Level: intermediate
11958b724c91SMatthew G. Knepley 
11968b724c91SMatthew G. Knepley .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
11978b724c91SMatthew G. Knepley @*/
11988b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorErrorDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
11998b724c91SMatthew G. Knepley {
12008b724c91SMatthew G. Knepley   DMAdaptor adaptor = (DMAdaptor)ctx;
12018b724c91SMatthew G. Knepley   char    **names;
12028b724c91SMatthew G. Knepley   PetscInt  Nf;
12038b724c91SMatthew G. Knepley 
12048b724c91SMatthew G. Knepley   PetscFunctionBegin;
12058b724c91SMatthew G. Knepley   PetscCall(DMGetNumFields(adaptor->idm, &Nf));
12068b724c91SMatthew G. Knepley   PetscCall(PetscMalloc1(Nf + 1, &names));
12078b724c91SMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
12088b724c91SMatthew G. Knepley     PetscObject disc;
12098b724c91SMatthew G. Knepley     const char *fname;
12108b724c91SMatthew G. Knepley     char        lname[PETSC_MAX_PATH_LEN];
12118b724c91SMatthew G. Knepley 
12128b724c91SMatthew G. Knepley     PetscCall(DMGetField(adaptor->idm, f, NULL, &disc));
12138b724c91SMatthew G. Knepley     PetscCall(PetscObjectGetName(disc, &fname));
12148b724c91SMatthew G. Knepley     PetscCall(PetscStrncpy(lname, fname, PETSC_MAX_PATH_LEN));
12158b724c91SMatthew G. Knepley     PetscCall(PetscStrlcat(lname, " Error", PETSC_MAX_PATH_LEN));
12168b724c91SMatthew G. Knepley     PetscCall(PetscStrallocpy(lname, &names[f]));
12178b724c91SMatthew G. Knepley   }
12188b724c91SMatthew G. Knepley   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
12198b724c91SMatthew G. Knepley   (*vf)->data = ctx;
12208b724c91SMatthew G. Knepley   PetscCall(KSPMonitorLGCreate(PetscObjectComm((PetscObject)viewer), NULL, NULL, "Log Error Norm", Nf, (const char **)names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg));
12218b724c91SMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(names[f]));
12228b724c91SMatthew G. Knepley   PetscCall(PetscFree(names));
12238b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
12248b724c91SMatthew G. Knepley }
12258b724c91SMatthew G. Knepley 
12268b724c91SMatthew G. Knepley /*@C
12278b724c91SMatthew G. Knepley   DMAdaptorMonitorErrorDrawLG - Plots the error norm at each iteration of an adaptive loop.
12288b724c91SMatthew G. Knepley 
12298b724c91SMatthew G. Knepley   Collective
12308b724c91SMatthew G. Knepley 
12318b724c91SMatthew G. Knepley   Input Parameters:
12328b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
12338b724c91SMatthew G. Knepley . n       - iteration number
12348b724c91SMatthew G. Knepley . odm     - the original `DM`
12358b724c91SMatthew G. Knepley . adm     - the adapted `DM`
12368b724c91SMatthew G. Knepley . Nf      - number of fields
12378b724c91SMatthew G. Knepley . enorms  - 2-norm error values for each field (may be estimated).
12388b724c91SMatthew G. Knepley . error   - `Vec` of cellwise errors
12398b724c91SMatthew G. Knepley - vf      - The viewer context
12408b724c91SMatthew G. Knepley 
12418b724c91SMatthew G. Knepley   Options Database Key:
12428b724c91SMatthew G. Knepley . -adaptor_error draw::draw_lg - Activates `DMAdaptorMonitorErrorDrawLG()`
12438b724c91SMatthew G. Knepley 
12448b724c91SMatthew G. Knepley   Level: intermediate
12458b724c91SMatthew G. Knepley 
12468b724c91SMatthew G. Knepley   Notes:
12478b724c91SMatthew G. Knepley   This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
12488b724c91SMatthew G. Knepley   to be used during the adaptation loop.
12498b724c91SMatthew G. Knepley 
12508b724c91SMatthew G. Knepley   Call `DMAdaptorMonitorErrorDrawLGCreate()` to create the context needed for this monitor
12518b724c91SMatthew G. Knepley 
1252*f8d70eaaSPierre Jolivet .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorError()`,
12538b724c91SMatthew G. Knepley           `DMAdaptorMonitorTrueResidualDrawLGCreate()`
12548b724c91SMatthew G. Knepley @*/
12558b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorErrorDrawLG(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
12568b724c91SMatthew G. Knepley {
12578b724c91SMatthew G. Knepley   PetscViewer       viewer = vf->viewer;
12588b724c91SMatthew G. Knepley   PetscViewerFormat format = vf->format;
12598b724c91SMatthew G. Knepley   PetscDrawLG       lg     = vf->lg;
12608b724c91SMatthew G. Knepley   PetscReal        *x, *e;
12618b724c91SMatthew G. Knepley 
12628b724c91SMatthew G. Knepley   PetscFunctionBegin;
12638b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
12648b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(lg, PETSC_DRAWLG_CLASSID, 8);
12658b724c91SMatthew G. Knepley   PetscCall(PetscCalloc2(Nf, &x, Nf, &e));
12668b724c91SMatthew G. Knepley   PetscCall(PetscViewerPushFormat(viewer, format));
12678b724c91SMatthew G. Knepley   if (!n) PetscCall(PetscDrawLGReset(lg));
12688b724c91SMatthew G. Knepley   for (PetscInt f = 0; f < Nf; ++f) {
12698b724c91SMatthew G. Knepley     x[f] = (PetscReal)n;
12708b724c91SMatthew G. Knepley     e[f] = enorms[f] > 0.0 ? PetscLog10Real(enorms[f]) : -15.;
12718b724c91SMatthew G. Knepley   }
12728b724c91SMatthew G. Knepley   PetscCall(PetscDrawLGAddPoint(lg, x, e));
12738b724c91SMatthew G. Knepley   PetscCall(PetscDrawLGDraw(lg));
12748b724c91SMatthew G. Knepley   PetscCall(PetscDrawLGSave(lg));
12758b724c91SMatthew G. Knepley   PetscCall(PetscViewerPopFormat(viewer));
12768b724c91SMatthew G. Knepley   PetscCall(PetscFree2(x, e));
12778b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
12788b724c91SMatthew G. Knepley }
12798b724c91SMatthew G. Knepley 
12808b724c91SMatthew G. Knepley /*@C
12818b724c91SMatthew G. Knepley   DMAdaptorMonitorRegisterAll - Registers all of the mesh adaptation monitors in the `SNES` package.
12828b724c91SMatthew G. Knepley 
12838b724c91SMatthew G. Knepley   Not Collective
12848b724c91SMatthew G. Knepley 
12858b724c91SMatthew G. Knepley   Level: advanced
12868b724c91SMatthew G. Knepley 
12878b724c91SMatthew G. Knepley .seealso: [](ch_snes), `SNES`, `DM`, `DMAdaptorMonitorRegister()`, `DMAdaptorRegister()`
12888b724c91SMatthew G. Knepley @*/
12898b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorRegisterAll(void)
12908b724c91SMatthew G. Knepley {
12918b724c91SMatthew G. Knepley   PetscFunctionBegin;
12928b724c91SMatthew G. Knepley   if (DMAdaptorMonitorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
12938b724c91SMatthew G. Knepley   DMAdaptorMonitorRegisterAllCalled = PETSC_TRUE;
12948b724c91SMatthew G. Knepley 
12958b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorRegister("size", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorSize, NULL, NULL));
12968b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorError, NULL, NULL));
12978b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorErrorDraw, NULL, NULL));
12988b724c91SMatthew G. Knepley   PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DRAW_LG, DMAdaptorMonitorErrorDrawLG, DMAdaptorMonitorErrorDrawLGCreate, NULL));
12998b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
13008b724c91SMatthew G. Knepley }
13018b724c91SMatthew G. Knepley 
13028b724c91SMatthew G. Knepley 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[])
13038b724c91SMatthew G. Knepley {
13048b724c91SMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
13058b724c91SMatthew G. Knepley 
13068b724c91SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) f[i] = u[i];
13078b724c91SMatthew G. Knepley }
13088b724c91SMatthew G. Knepley 
1309c6f61ee2SBarry Smith 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[])
1310c6f61ee2SBarry Smith {
13118b724c91SMatthew G. Knepley   for (PetscInt i = 0; i < dim; ++i) {
13128b724c91SMatthew G. Knepley     for (PetscInt j = 0; j < dim; ++j) f[i + dim * j] = u[i + dim * j];
1313c6f61ee2SBarry Smith   }
1314c6f61ee2SBarry Smith }
1315c6f61ee2SBarry Smith 
1316c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorAdapt_Sequence_Private(DMAdaptor adaptor, Vec inx, PetscBool doSolve, DM *adm, Vec *ax)
1317c6f61ee2SBarry Smith {
13183a336bb1SMatthew G. Knepley   PetscDS   ds;
13198b724c91SMatthew G. Knepley   PetscReal errorNorm = 0.;
1320c6f61ee2SBarry Smith   PetscInt  numAdapt  = adaptor->numSeq, adaptIter;
13213a336bb1SMatthew G. Knepley   PetscInt  dim, coordDim, Nf;
13228b724c91SMatthew G. Knepley   void     *ctx;
13238b724c91SMatthew G. Knepley   MPI_Comm  comm;
1324c6f61ee2SBarry Smith 
1325c6f61ee2SBarry Smith   PetscFunctionBegin;
1326c6f61ee2SBarry Smith   PetscCall(DMViewFromOptions(adaptor->idm, NULL, "-dm_adapt_pre_view"));
1327c6f61ee2SBarry Smith   PetscCall(VecViewFromOptions(inx, NULL, "-sol_adapt_pre_view"));
1328c6f61ee2SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)adaptor, &comm));
1329c6f61ee2SBarry Smith   PetscCall(DMGetDimension(adaptor->idm, &dim));
1330c6f61ee2SBarry Smith   PetscCall(DMGetCoordinateDim(adaptor->idm, &coordDim));
1331c6f61ee2SBarry Smith   PetscCall(DMGetApplicationContext(adaptor->idm, &ctx));
13323a336bb1SMatthew G. Knepley   PetscCall(DMGetDS(adaptor->idm, &ds));
13333a336bb1SMatthew G. Knepley   PetscCall(PetscDSGetNumFields(ds, &Nf));
13343a336bb1SMatthew G. Knepley   PetscCheck(Nf != 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot refine with no fields present!");
1335c6f61ee2SBarry Smith 
1336c6f61ee2SBarry Smith   /* Adapt until nothing changes */
1337c6f61ee2SBarry Smith   /* Adapt for a specified number of iterates */
1338c6f61ee2SBarry Smith   for (adaptIter = 0; adaptIter < numAdapt - 1; ++adaptIter) PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(comm)));
1339c6f61ee2SBarry Smith   for (adaptIter = 0; adaptIter < numAdapt; ++adaptIter) {
1340c6f61ee2SBarry Smith     PetscBool adapted = PETSC_FALSE;
1341c6f61ee2SBarry Smith     DM        dm      = adaptIter ? *adm : adaptor->idm, odm;
1342c6f61ee2SBarry Smith     Vec       x       = adaptIter ? *ax : inx, locX, ox;
13438b724c91SMatthew G. Knepley     Vec       error   = NULL;
1344c6f61ee2SBarry Smith 
1345c6f61ee2SBarry Smith     PetscCall(DMGetLocalVector(dm, &locX));
1346c6f61ee2SBarry Smith     PetscCall(DMAdaptorPreAdapt(adaptor, locX));
1347c6f61ee2SBarry Smith     if (doSolve) {
1348c6f61ee2SBarry Smith       SNES snes;
1349c6f61ee2SBarry Smith 
1350c6f61ee2SBarry Smith       PetscCall(DMAdaptorGetSolver(adaptor, &snes));
1351c6f61ee2SBarry Smith       PetscCall(SNESSolve(snes, NULL, adaptIter ? *ax : x));
1352c6f61ee2SBarry Smith     }
13533a336bb1SMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
13543a336bb1SMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
13553a336bb1SMatthew G. Knepley     PetscCall(VecViewFromOptions(adaptIter ? *ax : x, (PetscObject)adaptor, "-adapt_primal_sol_vec_view"));
1356c6f61ee2SBarry Smith     switch (adaptor->adaptCriterion) {
1357c6f61ee2SBarry Smith     case DM_ADAPTATION_REFINE:
1358c6f61ee2SBarry Smith       PetscCall(DMRefine(dm, comm, &odm));
1359c6f61ee2SBarry Smith       PetscCheck(odm, comm, PETSC_ERR_ARG_INCOMP, "DMRefine() did not perform any refinement, cannot continue grid sequencing");
1360c6f61ee2SBarry Smith       adapted = PETSC_TRUE;
13618b724c91SMatthew G. Knepley       PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, NULL));
1362c6f61ee2SBarry Smith       break;
1363c6f61ee2SBarry Smith     case DM_ADAPTATION_LABEL: {
1364c6f61ee2SBarry Smith       /* Adapt DM
1365c6f61ee2SBarry Smith            Create local solution
1366c6f61ee2SBarry Smith            Reconstruct gradients (FVM) or solve adjoint equation (FEM)
1367c6f61ee2SBarry Smith            Produce cellwise error indicator */
13683a336bb1SMatthew G. Knepley       DM             edm, plex;
13698b724c91SMatthew G. Knepley       PetscDS        ds;
13703a336bb1SMatthew G. Knepley       PetscFE        efe;
1371c6f61ee2SBarry Smith       DMLabel        adaptLabel;
1372c6f61ee2SBarry Smith       IS             refineIS, coarsenIS;
13733a336bb1SMatthew G. Knepley       DMPolytopeType ct;
13748b724c91SMatthew G. Knepley       PetscScalar    errorVal;
13753a336bb1SMatthew G. Knepley       PetscInt       nRefine, nCoarsen, cStart;
1376c6f61ee2SBarry Smith 
1377c6f61ee2SBarry Smith       PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
1378c6f61ee2SBarry Smith 
13793a336bb1SMatthew G. Knepley       // TODO Move this creation to PreAdapt
13803a336bb1SMatthew G. Knepley       PetscCall(DMClone(dm, &edm));
13813a336bb1SMatthew G. Knepley       PetscCall(DMConvert(edm, DMPLEX, &plex));
13823a336bb1SMatthew G. Knepley       PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, NULL));
13833a336bb1SMatthew G. Knepley       PetscCall(DMPlexGetCellType(plex, cStart, &ct));
13843a336bb1SMatthew G. Knepley       PetscCall(DMDestroy(&plex));
13853a336bb1SMatthew G. Knepley       PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 0, PETSC_DEFAULT, &efe));
13863a336bb1SMatthew G. Knepley       PetscCall(PetscObjectSetName((PetscObject)efe, "Error"));
13873a336bb1SMatthew G. Knepley       PetscCall(DMSetField(edm, 0, NULL, (PetscObject)efe));
13883a336bb1SMatthew G. Knepley       PetscCall(PetscFEDestroy(&efe));
13893a336bb1SMatthew G. Knepley       PetscCall(DMCreateDS(edm));
13908b724c91SMatthew G. Knepley       PetscCall(DMGetGlobalVector(edm, &error));
13918b724c91SMatthew G. Knepley       PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
1392c6f61ee2SBarry Smith 
13938b724c91SMatthew G. Knepley       PetscUseTypeMethod(adaptor, computeerrorindicator, locX, error);
13948b724c91SMatthew G. Knepley       PetscCall(VecViewFromOptions(error, (PetscObject)adaptor, "-adapt_error_vec_view"));
13958b724c91SMatthew G. Knepley       PetscCall(DMGetDS(edm, &ds));
13968b724c91SMatthew G. Knepley       PetscCall(PetscDSSetObjective(ds, 0, identity));
13978b724c91SMatthew G. Knepley       PetscCall(DMPlexComputeIntegralFEM(edm, error, &errorVal, NULL));
13988b724c91SMatthew G. Knepley       errorNorm = PetscRealPart(errorVal);
13993a336bb1SMatthew G. Knepley 
14003a336bb1SMatthew G. Knepley       // Compute IS from VecTagger
14018b724c91SMatthew G. Knepley       PetscCall(VecTaggerComputeIS(adaptor->refineTag, error, &refineIS, NULL));
14028b724c91SMatthew G. Knepley       PetscCall(VecTaggerComputeIS(adaptor->coarsenTag, error, &coarsenIS, NULL));
14033a336bb1SMatthew G. Knepley       PetscCall(ISViewFromOptions(refineIS, (PetscObject)adaptor->refineTag, "-is_view"));
14043a336bb1SMatthew G. Knepley       PetscCall(ISViewFromOptions(coarsenIS, (PetscObject)adaptor->coarsenTag, "-is_view"));
1405c6f61ee2SBarry Smith       PetscCall(ISGetSize(refineIS, &nRefine));
1406c6f61ee2SBarry Smith       PetscCall(ISGetSize(coarsenIS, &nCoarsen));
1407c6f61ee2SBarry Smith       PetscCall(PetscInfo(adaptor, "DMAdaptor: numRefine %" PetscInt_FMT ", numCoarsen %" PetscInt_FMT "\n", nRefine, nCoarsen));
1408c6f61ee2SBarry Smith       if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
1409c6f61ee2SBarry Smith       if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
1410c6f61ee2SBarry Smith       PetscCall(ISDestroy(&coarsenIS));
1411c6f61ee2SBarry Smith       PetscCall(ISDestroy(&refineIS));
14123a336bb1SMatthew G. Knepley       // Adapt DM from label
1413c6f61ee2SBarry Smith       if (nRefine || nCoarsen) {
14143a336bb1SMatthew G. Knepley         char        oprefix[PETSC_MAX_PATH_LEN];
14153a336bb1SMatthew G. Knepley         const char *p;
14163a336bb1SMatthew G. Knepley         PetscBool   flg;
14173a336bb1SMatthew G. Knepley 
14183a336bb1SMatthew G. Knepley         PetscCall(PetscOptionsHasName(NULL, adaptor->hdr.prefix, "-adapt_vec_view", &flg));
14193a336bb1SMatthew G. Knepley         if (flg) {
14203a336bb1SMatthew G. Knepley           Vec ref;
14213a336bb1SMatthew G. Knepley 
14223a336bb1SMatthew G. Knepley           PetscCall(DMPlexCreateLabelField(dm, adaptLabel, &ref));
14233a336bb1SMatthew G. Knepley           PetscCall(VecViewFromOptions(ref, (PetscObject)adaptor, "-adapt_vec_view"));
14243a336bb1SMatthew G. Knepley           PetscCall(VecDestroy(&ref));
14253a336bb1SMatthew G. Knepley         }
14263a336bb1SMatthew G. Knepley 
14273a336bb1SMatthew G. Knepley         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &p));
14283a336bb1SMatthew G. Knepley         PetscCall(PetscStrncpy(oprefix, p, PETSC_MAX_PATH_LEN));
14293a336bb1SMatthew G. Knepley         PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)dm, "adapt_"));
1430c6f61ee2SBarry Smith         PetscCall(DMAdaptLabel(dm, adaptLabel, &odm));
14313a336bb1SMatthew G. Knepley         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, oprefix));
14323a336bb1SMatthew G. Knepley         PetscCall(PetscObjectSetOptionsPrefix((PetscObject)odm, oprefix));
14338b724c91SMatthew G. Knepley         PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, error));
1434c6f61ee2SBarry Smith         adapted = PETSC_TRUE;
14358b724c91SMatthew G. Knepley       } else {
14368b724c91SMatthew G. Knepley         PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, error));
1437c6f61ee2SBarry Smith       }
1438c6f61ee2SBarry Smith       PetscCall(DMLabelDestroy(&adaptLabel));
14398b724c91SMatthew G. Knepley       PetscCall(DMRestoreGlobalVector(edm, &error));
14408b724c91SMatthew G. Knepley       PetscCall(DMDestroy(&edm));
1441c6f61ee2SBarry Smith     } break;
1442c6f61ee2SBarry Smith     case DM_ADAPTATION_METRIC: {
1443c6f61ee2SBarry Smith       DM        dmGrad, dmHess, dmMetric, dmDet;
1444c6f61ee2SBarry Smith       Vec       xGrad, xHess, metric, determinant;
1445c6f61ee2SBarry Smith       PetscReal N;
1446c6f61ee2SBarry Smith       DMLabel   bdLabel = NULL, rgLabel = NULL;
1447c6f61ee2SBarry Smith       PetscBool higherOrder = PETSC_FALSE;
1448c6f61ee2SBarry Smith       PetscInt  Nd          = coordDim * coordDim, f, vStart, vEnd;
1449c6f61ee2SBarry Smith       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[]);
1450c6f61ee2SBarry Smith 
1451c6f61ee2SBarry Smith       PetscCall(PetscMalloc(1, &funcs));
1452c6f61ee2SBarry Smith       funcs[0] = identityFunc;
1453c6f61ee2SBarry Smith 
1454c6f61ee2SBarry Smith       /*     Setup finite element spaces */
1455c6f61ee2SBarry Smith       PetscCall(DMClone(dm, &dmGrad));
1456c6f61ee2SBarry Smith       PetscCall(DMClone(dm, &dmHess));
14573a336bb1SMatthew G. Knepley       PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple fields not yet considered"); // TODO
14583a336bb1SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1459c6f61ee2SBarry Smith         PetscFE         fe, feGrad, feHess;
1460c6f61ee2SBarry Smith         PetscDualSpace  Q;
1461c6f61ee2SBarry Smith         PetscSpace      space;
1462c6f61ee2SBarry Smith         DM              K;
1463c6f61ee2SBarry Smith         PetscQuadrature q;
1464c6f61ee2SBarry Smith         PetscInt        Nc, qorder, p;
1465c6f61ee2SBarry Smith         const char     *prefix;
1466c6f61ee2SBarry Smith 
14673a336bb1SMatthew G. Knepley         PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1468c6f61ee2SBarry Smith         PetscCall(PetscFEGetNumComponents(fe, &Nc));
1469c6f61ee2SBarry Smith         PetscCheck(Nc <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple components not yet considered"); // TODO
1470c6f61ee2SBarry Smith         PetscCall(PetscFEGetBasisSpace(fe, &space));
1471c6f61ee2SBarry Smith         PetscCall(PetscSpaceGetDegree(space, NULL, &p));
1472c6f61ee2SBarry Smith         if (p > 1) higherOrder = PETSC_TRUE;
1473c6f61ee2SBarry Smith         PetscCall(PetscFEGetDualSpace(fe, &Q));
1474c6f61ee2SBarry Smith         PetscCall(PetscDualSpaceGetDM(Q, &K));
1475c6f61ee2SBarry Smith         PetscCall(DMPlexGetDepthStratum(K, 0, &vStart, &vEnd));
1476c6f61ee2SBarry Smith         PetscCall(PetscFEGetQuadrature(fe, &q));
1477c6f61ee2SBarry Smith         PetscCall(PetscQuadratureGetOrder(q, &qorder));
1478c6f61ee2SBarry Smith         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix));
1479c6f61ee2SBarry Smith         PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmGrad), dim, Nc * coordDim, PETSC_TRUE, prefix, qorder, &feGrad));
1480c6f61ee2SBarry Smith         PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmHess), dim, Nc * Nd, PETSC_TRUE, prefix, qorder, &feHess));
1481c6f61ee2SBarry Smith         PetscCall(DMSetField(dmGrad, f, NULL, (PetscObject)feGrad));
1482c6f61ee2SBarry Smith         PetscCall(DMSetField(dmHess, f, NULL, (PetscObject)feHess));
1483c6f61ee2SBarry Smith         PetscCall(DMCreateDS(dmGrad));
1484c6f61ee2SBarry Smith         PetscCall(DMCreateDS(dmHess));
1485c6f61ee2SBarry Smith         PetscCall(PetscFEDestroy(&feGrad));
1486c6f61ee2SBarry Smith         PetscCall(PetscFEDestroy(&feHess));
1487c6f61ee2SBarry Smith       }
1488c6f61ee2SBarry Smith       /*     Compute vertexwise gradients from cellwise gradients */
1489c6f61ee2SBarry Smith       PetscCall(DMCreateLocalVector(dmGrad, &xGrad));
1490c6f61ee2SBarry Smith       PetscCall(VecViewFromOptions(locX, NULL, "-sol_adapt_loc_pre_view"));
1491c6f61ee2SBarry Smith       PetscCall(DMPlexComputeGradientClementInterpolant(dm, locX, xGrad));
1492c6f61ee2SBarry Smith       PetscCall(VecViewFromOptions(xGrad, NULL, "-adapt_gradient_view"));
1493c6f61ee2SBarry Smith       /*     Compute vertexwise Hessians from cellwise Hessians */
1494c6f61ee2SBarry Smith       PetscCall(DMCreateLocalVector(dmHess, &xHess));
1495c6f61ee2SBarry Smith       PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, xGrad, xHess));
1496c6f61ee2SBarry Smith       PetscCall(VecViewFromOptions(xHess, NULL, "-adapt_hessian_view"));
1497c6f61ee2SBarry Smith       PetscCall(VecDestroy(&xGrad));
1498c6f61ee2SBarry Smith       PetscCall(DMDestroy(&dmGrad));
1499c6f61ee2SBarry Smith       /*     Compute L-p normalized metric */
1500c6f61ee2SBarry Smith       PetscCall(DMClone(dm, &dmMetric));
1501c6f61ee2SBarry Smith       N = adaptor->Nadapt >= 0 ? adaptor->Nadapt : PetscPowRealInt(adaptor->refinementFactor, dim) * ((PetscReal)(vEnd - vStart));
15028b724c91SMatthew G. Knepley       // TODO This was where the old monitor was, figure out how to show metric and target N
1503c6f61ee2SBarry Smith       PetscCall(DMPlexMetricSetTargetComplexity(dmMetric, (PetscReal)N));
1504c6f61ee2SBarry Smith       if (higherOrder) {
1505c6f61ee2SBarry Smith         /*   Project Hessian into P1 space, if required */
1506c6f61ee2SBarry Smith         PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1507c6f61ee2SBarry Smith         PetscCall(DMProjectFieldLocal(dmMetric, 0.0, xHess, funcs, INSERT_ALL_VALUES, metric));
1508c6f61ee2SBarry Smith         PetscCall(VecDestroy(&xHess));
1509c6f61ee2SBarry Smith         xHess = metric;
1510c6f61ee2SBarry Smith       }
1511c6f61ee2SBarry Smith       PetscCall(PetscFree(funcs));
1512c6f61ee2SBarry Smith       PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1513c6f61ee2SBarry Smith       PetscCall(DMPlexMetricDeterminantCreate(dmMetric, 0, &determinant, &dmDet));
1514c6f61ee2SBarry Smith       PetscCall(DMPlexMetricNormalize(dmMetric, xHess, PETSC_TRUE, PETSC_TRUE, metric, determinant));
1515c6f61ee2SBarry Smith       PetscCall(VecDestroy(&determinant));
1516c6f61ee2SBarry Smith       PetscCall(DMDestroy(&dmDet));
1517c6f61ee2SBarry Smith       PetscCall(VecDestroy(&xHess));
1518c6f61ee2SBarry Smith       PetscCall(DMDestroy(&dmHess));
1519c6f61ee2SBarry Smith       /*     Adapt DM from metric */
1520c6f61ee2SBarry Smith       PetscCall(DMGetLabel(dm, "marker", &bdLabel));
1521c6f61ee2SBarry Smith       PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &odm));
15228b724c91SMatthew G. Knepley       PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, NULL));
1523c6f61ee2SBarry Smith       adapted = PETSC_TRUE;
1524c6f61ee2SBarry Smith       /*     Cleanup */
1525c6f61ee2SBarry Smith       PetscCall(VecDestroy(&metric));
1526c6f61ee2SBarry Smith       PetscCall(DMDestroy(&dmMetric));
1527c6f61ee2SBarry Smith     } break;
1528c6f61ee2SBarry Smith     default:
1529c6f61ee2SBarry Smith       SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid adaptation type: %d", adaptor->adaptCriterion);
1530c6f61ee2SBarry Smith     }
1531c6f61ee2SBarry Smith     PetscCall(DMAdaptorPostAdapt(adaptor));
1532c6f61ee2SBarry Smith     PetscCall(DMRestoreLocalVector(dm, &locX));
1533c6f61ee2SBarry Smith     /* If DM was adapted, replace objects and recreate solution */
1534c6f61ee2SBarry Smith     if (adapted) {
1535c6f61ee2SBarry Smith       const char *name;
1536c6f61ee2SBarry Smith 
1537c6f61ee2SBarry Smith       PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1538c6f61ee2SBarry Smith       PetscCall(PetscObjectSetName((PetscObject)odm, name));
1539c6f61ee2SBarry Smith       /* Reconfigure solver */
1540c6f61ee2SBarry Smith       PetscCall(SNESReset(adaptor->snes));
1541c6f61ee2SBarry Smith       PetscCall(SNESSetDM(adaptor->snes, odm));
1542c6f61ee2SBarry Smith       PetscCall(DMAdaptorSetSolver(adaptor, adaptor->snes));
1543c6f61ee2SBarry Smith       PetscCall(DMPlexSetSNESLocalFEM(odm, PETSC_FALSE, ctx));
1544c6f61ee2SBarry Smith       PetscCall(SNESSetFromOptions(adaptor->snes));
1545c6f61ee2SBarry Smith       /* Transfer system */
1546c6f61ee2SBarry Smith       PetscCall(DMCopyDisc(dm, odm));
1547c6f61ee2SBarry Smith       /* Transfer solution */
1548c6f61ee2SBarry Smith       PetscCall(DMCreateGlobalVector(odm, &ox));
1549c6f61ee2SBarry Smith       PetscCall(PetscObjectGetName((PetscObject)x, &name));
1550c6f61ee2SBarry Smith       PetscCall(PetscObjectSetName((PetscObject)ox, name));
1551c6f61ee2SBarry Smith       PetscCall(DMAdaptorTransferSolution(adaptor, dm, x, odm, ox));
1552c6f61ee2SBarry Smith       /* Cleanup adaptivity info */
1553c6f61ee2SBarry Smith       if (adaptIter > 0) PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(comm)));
1554c6f61ee2SBarry Smith       PetscCall(DMForestSetAdaptivityForest(dm, NULL)); /* clear internal references to the previous dm */
1555c6f61ee2SBarry Smith       PetscCall(DMDestroy(&dm));
1556c6f61ee2SBarry Smith       PetscCall(VecDestroy(&x));
1557c6f61ee2SBarry Smith       *adm = odm;
1558c6f61ee2SBarry Smith       *ax  = ox;
1559c6f61ee2SBarry Smith     } else {
1560c6f61ee2SBarry Smith       *adm      = dm;
1561c6f61ee2SBarry Smith       *ax       = x;
1562c6f61ee2SBarry Smith       adaptIter = numAdapt;
1563c6f61ee2SBarry Smith     }
1564c6f61ee2SBarry Smith     if (adaptIter < numAdapt - 1) {
1565c6f61ee2SBarry Smith       PetscCall(DMViewFromOptions(odm, NULL, "-dm_adapt_iter_view"));
1566c6f61ee2SBarry Smith       PetscCall(VecViewFromOptions(ox, NULL, "-sol_adapt_iter_view"));
1567c6f61ee2SBarry Smith     }
1568c6f61ee2SBarry Smith   }
1569c6f61ee2SBarry Smith   PetscCall(DMViewFromOptions(*adm, NULL, "-dm_adapt_view"));
1570c6f61ee2SBarry Smith   PetscCall(VecViewFromOptions(*ax, NULL, "-sol_adapt_view"));
1571c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1572c6f61ee2SBarry Smith }
1573c6f61ee2SBarry Smith 
1574c6f61ee2SBarry Smith /*@
1575c6f61ee2SBarry Smith   DMAdaptorAdapt - Creates a new `DM` that is adapted to the problem
1576c6f61ee2SBarry Smith 
1577c6f61ee2SBarry Smith   Not Collective
1578c6f61ee2SBarry Smith 
1579c6f61ee2SBarry Smith   Input Parameters:
1580c6f61ee2SBarry Smith + adaptor  - The `DMAdaptor` object
1581c6f61ee2SBarry Smith . x        - The global approximate solution
1582c6f61ee2SBarry Smith - strategy - The adaptation strategy, see `DMAdaptationStrategy`
1583c6f61ee2SBarry Smith 
1584c6f61ee2SBarry Smith   Output Parameters:
1585c6f61ee2SBarry Smith + adm - The adapted `DM`
1586c6f61ee2SBarry Smith - ax  - The adapted solution
1587c6f61ee2SBarry Smith 
1588c6f61ee2SBarry Smith   Options Database Keys:
1589c6f61ee2SBarry Smith + -snes_adapt <strategy> - initial, sequential, multigrid
1590c6f61ee2SBarry Smith . -adapt_gradient_view   - View the Clement interpolant of the solution gradient
1591c6f61ee2SBarry Smith . -adapt_hessian_view    - View the Clement interpolant of the solution Hessian
1592c6f61ee2SBarry Smith - -adapt_metric_view     - View the metric tensor for adaptive mesh refinement
1593c6f61ee2SBarry Smith 
1594c6f61ee2SBarry Smith   Level: intermediate
1595c6f61ee2SBarry Smith 
1596df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`
1597c6f61ee2SBarry Smith @*/
1598c6f61ee2SBarry Smith PetscErrorCode DMAdaptorAdapt(DMAdaptor adaptor, Vec x, DMAdaptationStrategy strategy, DM *adm, Vec *ax)
1599c6f61ee2SBarry Smith {
1600c6f61ee2SBarry Smith   PetscFunctionBegin;
1601c6f61ee2SBarry Smith   switch (strategy) {
1602c6f61ee2SBarry Smith   case DM_ADAPTATION_INITIAL:
1603c6f61ee2SBarry Smith     PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_FALSE, adm, ax));
1604c6f61ee2SBarry Smith     break;
1605c6f61ee2SBarry Smith   case DM_ADAPTATION_SEQUENTIAL:
1606c6f61ee2SBarry Smith     PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_TRUE, adm, ax));
1607c6f61ee2SBarry Smith     break;
1608c6f61ee2SBarry Smith   default:
1609c6f61ee2SBarry Smith     SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Unrecognized adaptation strategy %d", strategy);
1610c6f61ee2SBarry Smith   }
1611c6f61ee2SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1612c6f61ee2SBarry Smith }
16133a336bb1SMatthew G. Knepley 
16143a336bb1SMatthew G. Knepley /*@C
16153a336bb1SMatthew G. Knepley   DMAdaptorGetMixedSetupFunction - Get the function setting up the mixed problem, if it exists
16163a336bb1SMatthew G. Knepley 
16173a336bb1SMatthew G. Knepley   Not Collective
16183a336bb1SMatthew G. Knepley 
16193a336bb1SMatthew G. Knepley   Input Parameter:
16203a336bb1SMatthew G. Knepley . adaptor - the `DMAdaptor`
16213a336bb1SMatthew G. Knepley 
16223a336bb1SMatthew G. Knepley   Output Parameter:
16233a336bb1SMatthew G. Knepley . setupFunc - the function setting up the mixed problem, or `NULL`
16243a336bb1SMatthew G. Knepley 
16253a336bb1SMatthew G. Knepley   Level: advanced
16263a336bb1SMatthew G. Knepley 
16273a336bb1SMatthew G. Knepley .seealso: `DMAdaptor`, `DMAdaptorSetMixedSetupFunction()`, `DMAdaptorAdapt()`
16283a336bb1SMatthew G. Knepley @*/
16293a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorGetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (**setupFunc)(DMAdaptor, DM))
16303a336bb1SMatthew G. Knepley {
16313a336bb1SMatthew G. Knepley   PetscFunctionBegin;
16323a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
16333a336bb1SMatthew G. Knepley   PetscAssertPointer(setupFunc, 2);
16343a336bb1SMatthew G. Knepley   *setupFunc = adaptor->ops->mixedsetup;
16353a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
16363a336bb1SMatthew G. Knepley }
16373a336bb1SMatthew G. Knepley 
16383a336bb1SMatthew G. Knepley /*@C
16393a336bb1SMatthew G. Knepley   DMAdaptorSetMixedSetupFunction - Set the function setting up the mixed problem
16403a336bb1SMatthew G. Knepley 
16413a336bb1SMatthew G. Knepley   Not Collective
16423a336bb1SMatthew G. Knepley 
16433a336bb1SMatthew G. Knepley   Input Parameters:
16443a336bb1SMatthew G. Knepley + adaptor   - the `DMAdaptor`
16453a336bb1SMatthew G. Knepley - setupFunc - the function setting up the mixed problem
16463a336bb1SMatthew G. Knepley 
16473a336bb1SMatthew G. Knepley   Level: advanced
16483a336bb1SMatthew G. Knepley 
16493a336bb1SMatthew G. Knepley .seealso: `DMAdaptor`, `DMAdaptorGetMixedSetupFunction()`, `DMAdaptorAdapt()`
16503a336bb1SMatthew G. Knepley @*/
16513a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorSetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (*setupFunc)(DMAdaptor, DM))
16523a336bb1SMatthew G. Knepley {
16533a336bb1SMatthew G. Knepley   PetscFunctionBegin;
16543a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
16553a336bb1SMatthew G. Knepley   PetscValidFunction(setupFunc, 2);
16563a336bb1SMatthew G. Knepley   adaptor->ops->mixedsetup = setupFunc;
16573a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
16583a336bb1SMatthew G. Knepley }
16593a336bb1SMatthew G. Knepley 
16608b724c91SMatthew G. Knepley /*@C
16618b724c91SMatthew G. Knepley   DMAdaptorGetCriterion - Get the adaptation criterion
16628b724c91SMatthew G. Knepley 
16638b724c91SMatthew G. Knepley   Not Collective
16648b724c91SMatthew G. Knepley 
16658b724c91SMatthew G. Knepley   Input Parameter:
16668b724c91SMatthew G. Knepley . adaptor - the `DMAdaptor`
16678b724c91SMatthew G. Knepley 
16688b724c91SMatthew G. Knepley   Output Parameter:
16698b724c91SMatthew G. Knepley . criterion - the criterion for adaptation
16708b724c91SMatthew G. Knepley 
16718b724c91SMatthew G. Knepley   Level: advanced
16728b724c91SMatthew G. Knepley 
1673*f8d70eaaSPierre Jolivet .seealso: `DMAdaptor`, `DMAdaptorSetCriterion()`, `DMAdaptationCriterion`
16748b724c91SMatthew G. Knepley @*/
16758b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorGetCriterion(DMAdaptor adaptor, DMAdaptationCriterion *criterion)
16768b724c91SMatthew G. Knepley {
16778b724c91SMatthew G. Knepley   PetscFunctionBegin;
16788b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
16798b724c91SMatthew G. Knepley   PetscAssertPointer(criterion, 2);
16808b724c91SMatthew G. Knepley   *criterion = adaptor->adaptCriterion;
16818b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
16828b724c91SMatthew G. Knepley }
16838b724c91SMatthew G. Knepley 
16848b724c91SMatthew G. Knepley /*@C
16858b724c91SMatthew G. Knepley   DMAdaptorSetCriterion - Set the adaptation criterion
16868b724c91SMatthew G. Knepley 
16878b724c91SMatthew G. Knepley   Not Collective
16888b724c91SMatthew G. Knepley 
16898b724c91SMatthew G. Knepley   Input Parameters:
16908b724c91SMatthew G. Knepley + adaptor   - the `DMAdaptor`
16918b724c91SMatthew G. Knepley - criterion - the adaptation criterion
16928b724c91SMatthew G. Knepley 
16938b724c91SMatthew G. Knepley   Level: advanced
16948b724c91SMatthew G. Knepley 
16958b724c91SMatthew G. Knepley .seealso: `DMAdaptor`, `DMAdaptorGetCriterion()`, `DMAdaptationCriterion`
16968b724c91SMatthew G. Knepley @*/
16978b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorSetCriterion(DMAdaptor adaptor, DMAdaptationCriterion criterion)
16988b724c91SMatthew G. Knepley {
16998b724c91SMatthew G. Knepley   PetscFunctionBegin;
17008b724c91SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
17018b724c91SMatthew G. Knepley   adaptor->adaptCriterion = criterion;
17028b724c91SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
17038b724c91SMatthew G. Knepley }
17048b724c91SMatthew G. Knepley 
17053a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorInitialize_Gradient(DMAdaptor adaptor)
17063a336bb1SMatthew G. Knepley {
17073a336bb1SMatthew G. Knepley   PetscFunctionBegin;
17083a336bb1SMatthew G. Knepley   adaptor->ops->computeerrorindicator     = DMAdaptorComputeErrorIndicator_Gradient;
17093a336bb1SMatthew G. Knepley   adaptor->ops->computecellerrorindicator = DMAdaptorComputeCellErrorIndicator_Gradient;
17103a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
17113a336bb1SMatthew G. Knepley }
17123a336bb1SMatthew G. Knepley 
17133a336bb1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Gradient(DMAdaptor adaptor)
17143a336bb1SMatthew G. Knepley {
17153a336bb1SMatthew G. Knepley   PetscFunctionBegin;
17163a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
17173a336bb1SMatthew G. Knepley   adaptor->data = NULL;
17183a336bb1SMatthew G. Knepley 
17193a336bb1SMatthew G. Knepley   PetscCall(DMAdaptorInitialize_Gradient(adaptor));
17203a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
17213a336bb1SMatthew G. Knepley }
17223a336bb1SMatthew G. Knepley 
17233a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorInitialize_Flux(DMAdaptor adaptor)
17243a336bb1SMatthew G. Knepley {
17253a336bb1SMatthew G. Knepley   PetscFunctionBegin;
17263a336bb1SMatthew G. Knepley   adaptor->ops->computeerrorindicator = DMAdaptorComputeErrorIndicator_Flux;
17273a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
17283a336bb1SMatthew G. Knepley }
17293a336bb1SMatthew G. Knepley 
17303a336bb1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Flux(DMAdaptor adaptor)
17313a336bb1SMatthew G. Knepley {
17323a336bb1SMatthew G. Knepley   PetscFunctionBegin;
17333a336bb1SMatthew G. Knepley   PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
17343a336bb1SMatthew G. Knepley   adaptor->data = NULL;
17353a336bb1SMatthew G. Knepley 
17363a336bb1SMatthew G. Knepley   PetscCall(DMAdaptorInitialize_Flux(adaptor));
17373a336bb1SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
17383a336bb1SMatthew G. Knepley }
1739