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 @*/
DMAdaptorRegister(const char name[],PetscErrorCode (* create_func)(DMAdaptor))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 @*/
DMAdaptorRegisterAll(void)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 @*/
DMAdaptorRegisterDestroy(void)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
DMAdaptorMonitorMakeKey_Internal(const char name[],PetscViewerType vtype,PetscViewerFormat format,char key[])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 @*/
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 **))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 @*/
DMAdaptorMonitorRegisterDestroy(void)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 @*/
DMAdaptorCreate(MPI_Comm comm,DMAdaptor * adaptor)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 @*/
DMAdaptorDestroy(DMAdaptor * adaptor)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 @*/
DMAdaptorSetType(DMAdaptor adaptor,DMAdaptorType method)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 @*/
DMAdaptorGetType(DMAdaptor adaptor,DMAdaptorType * type)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
PetscViewerAndFormatCreate_Internal(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat ** vf)310*2a8381b2SBarry Smith static PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer viewer, PetscViewerFormat format, PetscCtx 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)
32849abdd8aSBarry Smith - monitordestroy - [optional] routine that frees monitor context (may be `NULL`), see `PetscCtxDestroyFn` for its calling sequence
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 Options Database Keys:
3418b724c91SMatthew G. Knepley + -adaptor_monitor_size - sets `DMAdaptorMonitorSize()`
3428b724c91SMatthew G. Knepley . -adaptor_monitor_error - sets `DMAdaptorMonitorError()`
3438b724c91SMatthew G. Knepley . -adaptor_monitor_error draw - sets `DMAdaptorMonitorErrorDraw()` and plots error
3448b724c91SMatthew G. Knepley . -adaptor_monitor_error draw::draw_lg - sets `DMAdaptorMonitorErrorDrawLG()` and plots error
3458b724c91SMatthew 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.
3468b724c91SMatthew G. Knepley
3478b724c91SMatthew G. Knepley Level: beginner
3488b724c91SMatthew G. Knepley
34949abdd8aSBarry Smith .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptor`, `PetscCtxDestroyFn`
3508b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorSet(DMAdaptor adaptor,PetscErrorCode (* monitor)(DMAdaptor adaptor,PetscInt it,DM odm,DM adm,PetscInt Nf,PetscReal enorms[],Vec error,PetscCtx ctx),PetscCtx ctx,PetscCtxDestroyFn * monitordestroy)351*2a8381b2SBarry Smith PetscErrorCode DMAdaptorMonitorSet(DMAdaptor adaptor, PetscErrorCode (*monitor)(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscCtx ctx), PetscCtx ctx, PetscCtxDestroyFn *monitordestroy)
3528b724c91SMatthew G. Knepley {
3538b724c91SMatthew G. Knepley PetscFunctionBegin;
3548b724c91SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
3558b724c91SMatthew G. Knepley for (PetscInt i = 0; i < adaptor->numbermonitors; i++) {
356453a69bbSBarry Smith PetscBool identical;
357453a69bbSBarry Smith
358453a69bbSBarry Smith PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)monitor, ctx, monitordestroy, (PetscErrorCode (*)(void))(PetscVoidFn *)adaptor->monitor[i], adaptor->monitorcontext[i], adaptor->monitordestroy[i], &identical));
3598b724c91SMatthew G. Knepley if (identical) PetscFunctionReturn(PETSC_SUCCESS);
3608b724c91SMatthew G. Knepley }
3618b724c91SMatthew G. Knepley PetscCheck(adaptor->numbermonitors < MAXDMADAPTORMONITORS, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_OUTOFRANGE, "Too many DMAdaptor monitors set");
3628b724c91SMatthew G. Knepley adaptor->monitor[adaptor->numbermonitors] = monitor;
3638b724c91SMatthew G. Knepley adaptor->monitordestroy[adaptor->numbermonitors] = monitordestroy;
364835f2295SStefano Zampini adaptor->monitorcontext[adaptor->numbermonitors++] = ctx;
3658b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
3668b724c91SMatthew G. Knepley }
3678b724c91SMatthew G. Knepley
3688b724c91SMatthew G. Knepley /*@
3698b724c91SMatthew G. Knepley DMAdaptorMonitorCancel - Clears all monitors for a `DMAdaptor` object.
3708b724c91SMatthew G. Knepley
3718b724c91SMatthew G. Knepley Logically Collective
3728b724c91SMatthew G. Knepley
3738b724c91SMatthew G. Knepley Input Parameter:
3748b724c91SMatthew G. Knepley . adaptor - the `DMAdaptor`
3758b724c91SMatthew G. Knepley
3768b724c91SMatthew G. Knepley Options Database Key:
3778b724c91SMatthew 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.
3788b724c91SMatthew G. Knepley
3798b724c91SMatthew G. Knepley Level: intermediate
3808b724c91SMatthew G. Knepley
3818b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptorMonitorError()`, `DMAdaptorMonitorSet()`, `DMAdaptor`
3828b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorCancel(DMAdaptor adaptor)3838b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorCancel(DMAdaptor adaptor)
3848b724c91SMatthew G. Knepley {
3858b724c91SMatthew G. Knepley PetscFunctionBegin;
3868b724c91SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
3878b724c91SMatthew G. Knepley for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) {
3888b724c91SMatthew G. Knepley if (adaptor->monitordestroy[i]) PetscCall((*adaptor->monitordestroy[i])(&adaptor->monitorcontext[i]));
3898b724c91SMatthew G. Knepley }
3908b724c91SMatthew G. Knepley adaptor->numbermonitors = 0;
3918b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
3928b724c91SMatthew G. Knepley }
3938b724c91SMatthew G. Knepley
3948b724c91SMatthew G. Knepley /*@C
3958b724c91SMatthew G. Knepley DMAdaptorMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user in the options database
3968b724c91SMatthew G. Knepley
3978b724c91SMatthew G. Knepley Collective
3988b724c91SMatthew G. Knepley
3998b724c91SMatthew G. Knepley Input Parameters:
4008b724c91SMatthew G. Knepley + adaptor - `DMadaptor` object you wish to monitor
4018b724c91SMatthew G. Knepley . opt - the command line option for this monitor
4028b724c91SMatthew G. Knepley . name - the monitor type one is seeking
4038b724c91SMatthew G. Knepley - ctx - An optional user context for the monitor, or `NULL`
4048b724c91SMatthew G. Knepley
4058b724c91SMatthew G. Knepley Level: developer
4068b724c91SMatthew G. Knepley
4078b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptorMonitorRegister()`, `DMAdaptorMonitorSet()`, `PetscOptionsGetViewer()`
4088b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorSetFromOptions(DMAdaptor adaptor,const char opt[],const char name[],PetscCtx ctx)409*2a8381b2SBarry Smith PetscErrorCode DMAdaptorMonitorSetFromOptions(DMAdaptor adaptor, const char opt[], const char name[], PetscCtx ctx)
4108b724c91SMatthew G. Knepley {
4118b724c91SMatthew G. Knepley PetscErrorCode (*mfunc)(DMAdaptor, PetscInt, DM, DM, PetscInt, PetscReal[], Vec, void *);
4128b724c91SMatthew G. Knepley PetscErrorCode (*cfunc)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
4138b724c91SMatthew G. Knepley PetscErrorCode (*dfunc)(PetscViewerAndFormat **);
4148b724c91SMatthew G. Knepley PetscViewerAndFormat *vf;
4158b724c91SMatthew G. Knepley PetscViewer viewer;
4168b724c91SMatthew G. Knepley PetscViewerFormat format;
4178b724c91SMatthew G. Knepley PetscViewerType vtype;
4188b724c91SMatthew G. Knepley char key[PETSC_MAX_PATH_LEN];
4198b724c91SMatthew G. Knepley PetscBool flg;
4208b724c91SMatthew G. Knepley const char *prefix = NULL;
4218b724c91SMatthew G. Knepley
4228b724c91SMatthew G. Knepley PetscFunctionBegin;
4238b724c91SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
4248b724c91SMatthew G. Knepley PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)adaptor), ((PetscObject)adaptor)->options, prefix, opt, &viewer, &format, &flg));
4258b724c91SMatthew G. Knepley if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
4268b724c91SMatthew G. Knepley
4278b724c91SMatthew G. Knepley PetscCall(PetscViewerGetType(viewer, &vtype));
4288b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitorMakeKey_Internal(name, vtype, format, key));
4298b724c91SMatthew G. Knepley PetscCall(PetscFunctionListFind(DMAdaptorMonitorList, key, &mfunc));
4308b724c91SMatthew G. Knepley PetscCall(PetscFunctionListFind(DMAdaptorMonitorCreateList, key, &cfunc));
4318b724c91SMatthew G. Knepley PetscCall(PetscFunctionListFind(DMAdaptorMonitorDestroyList, key, &dfunc));
4328b724c91SMatthew G. Knepley if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
4338b724c91SMatthew G. Knepley if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
4348b724c91SMatthew G. Knepley
4358b724c91SMatthew G. Knepley PetscCall((*cfunc)(viewer, format, ctx, &vf));
4368b724c91SMatthew G. Knepley PetscCall(PetscViewerDestroy(&viewer));
43749abdd8aSBarry Smith PetscCall(DMAdaptorMonitorSet(adaptor, mfunc, vf, (PetscCtxDestroyFn *)dfunc));
4388b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
4398b724c91SMatthew G. Knepley }
4408b724c91SMatthew G. Knepley
4413a336bb1SMatthew G. Knepley /*@
442e03fd340SMatthew G. Knepley DMAdaptorSetOptionsPrefix - Sets the prefix used for searching for all `DMAdaptor` options in the database.
443e03fd340SMatthew G. Knepley
444e03fd340SMatthew G. Knepley Logically Collective
445e03fd340SMatthew G. Knepley
446e03fd340SMatthew G. Knepley Input Parameters:
447e03fd340SMatthew G. Knepley + adaptor - the `DMAdaptor`
448e03fd340SMatthew G. Knepley - prefix - the prefix to prepend to all option names
449e03fd340SMatthew G. Knepley
450e03fd340SMatthew G. Knepley Level: advanced
451e03fd340SMatthew G. Knepley
452e03fd340SMatthew G. Knepley Note:
453e03fd340SMatthew G. Knepley A hyphen (-) must NOT be given at the beginning of the prefix name.
454e03fd340SMatthew G. Knepley The first character of all runtime options is AUTOMATICALLY the hyphen.
455e03fd340SMatthew G. Knepley
456e03fd340SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptor`, `SNESSetOptionsPrefix()`, `DMAdaptorSetFromOptions()`
457e03fd340SMatthew G. Knepley @*/
DMAdaptorSetOptionsPrefix(DMAdaptor adaptor,const char prefix[])458e03fd340SMatthew G. Knepley PetscErrorCode DMAdaptorSetOptionsPrefix(DMAdaptor adaptor, const char prefix[])
459e03fd340SMatthew G. Knepley {
460e03fd340SMatthew G. Knepley PetscFunctionBegin;
461e03fd340SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
462e03fd340SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor, prefix));
463e03fd340SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor->refineTag, prefix));
464e03fd340SMatthew G. Knepley PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)adaptor->refineTag, "refine_"));
465e03fd340SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adaptor->coarsenTag, prefix));
466e03fd340SMatthew G. Knepley PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)adaptor->coarsenTag, "coarsen_"));
467e03fd340SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
468e03fd340SMatthew G. Knepley }
469e03fd340SMatthew G. Knepley
470e03fd340SMatthew G. Knepley /*@
471c6f61ee2SBarry Smith DMAdaptorSetFromOptions - Sets properties of a `DMAdaptor` object from values in the options database
472c6f61ee2SBarry Smith
473c6f61ee2SBarry Smith Collective
474c6f61ee2SBarry Smith
475c6f61ee2SBarry Smith Input Parameter:
476c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
477c6f61ee2SBarry Smith
478c6f61ee2SBarry Smith Options Database Keys:
4798b724c91SMatthew G. Knepley + -adaptor_monitor_size - Monitor the mesh size
4808b724c91SMatthew G. Knepley . -adaptor_monitor_error - Monitor the solution error
481c6f61ee2SBarry Smith . -adaptor_sequence_num <num> - Number of adaptations to generate an optimal grid
482c6f61ee2SBarry Smith . -adaptor_target_num <num> - Set the target number of vertices N_adapt, -1 for automatic determination
4833a336bb1SMatthew G. Knepley . -adaptor_refinement_factor <r> - Set r such that N_adapt = r^dim N_orig
484d7c1f440SPierre Jolivet - -adaptor_mixed_setup_function <func> - Set the function func that sets up the mixed problem
485c6f61ee2SBarry Smith
486c6f61ee2SBarry Smith Level: beginner
487c6f61ee2SBarry Smith
488df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
489c6f61ee2SBarry Smith @*/
DMAdaptorSetFromOptions(DMAdaptor adaptor)490c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetFromOptions(DMAdaptor adaptor)
491c6f61ee2SBarry Smith {
4923a336bb1SMatthew G. Knepley char typeName[PETSC_MAX_PATH_LEN];
4933a336bb1SMatthew G. Knepley const char *defName = DMADAPTORGRADIENT;
4943a336bb1SMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN];
4958b724c91SMatthew G. Knepley DMAdaptationCriterion criterion = DM_ADAPTATION_NONE;
4963a336bb1SMatthew G. Knepley PetscBool flg;
4973a336bb1SMatthew G. Knepley
498c6f61ee2SBarry Smith PetscFunctionBegin;
4993a336bb1SMatthew G. Knepley PetscObjectOptionsBegin((PetscObject)adaptor);
5003a336bb1SMatthew G. Knepley PetscCall(PetscOptionsFList("-adaptor_type", "DMAdaptor", "DMAdaptorSetType", DMAdaptorList, defName, typeName, 1024, &flg));
5013a336bb1SMatthew G. Knepley if (flg) PetscCall(DMAdaptorSetType(adaptor, typeName));
5023a336bb1SMatthew G. Knepley else if (!((PetscObject)adaptor)->type_name) PetscCall(DMAdaptorSetType(adaptor, defName));
5038b724c91SMatthew G. Knepley PetscCall(PetscOptionsEnum("-adaptor_criterion", "Criterion used to drive adaptation", "", DMAdaptationCriteria, (PetscEnum)criterion, (PetscEnum *)&criterion, &flg));
5048b724c91SMatthew G. Knepley if (flg) PetscCall(DMAdaptorSetCriterion(adaptor, criterion));
505c6f61ee2SBarry Smith PetscCall(PetscOptionsInt("-adaptor_sequence_num", "Number of adaptations to generate an optimal grid", "DMAdaptorSetSequenceLength", adaptor->numSeq, &adaptor->numSeq, NULL));
506c6f61ee2SBarry Smith PetscCall(PetscOptionsInt("-adaptor_target_num", "Set the target number of vertices N_adapt, -1 for automatic determination", "DMAdaptor", adaptor->Nadapt, &adaptor->Nadapt, NULL));
507c6f61ee2SBarry Smith PetscCall(PetscOptionsReal("-adaptor_refinement_factor", "Set r such that N_adapt = r^dim N_orig", "DMAdaptor", adaptor->refinementFactor, &adaptor->refinementFactor, NULL));
5083a336bb1SMatthew G. Knepley PetscCall(PetscOptionsString("-adaptor_mixed_setup_function", "Function to setup the mixed problem", "DMAdaptorSetMixedSetupFunction", funcname, funcname, sizeof(funcname), &flg));
5093a336bb1SMatthew G. Knepley if (flg) {
5103a336bb1SMatthew G. Knepley PetscErrorCode (*setupFunc)(DMAdaptor, DM);
5113a336bb1SMatthew G. Knepley
5123a336bb1SMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&setupFunc));
5133a336bb1SMatthew G. Knepley PetscCheck(setupFunc, PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname);
5143a336bb1SMatthew G. Knepley PetscCall(DMAdaptorSetMixedSetupFunction(adaptor, setupFunc));
5153a336bb1SMatthew G. Knepley }
5168b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_size", "size", adaptor));
5178b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitorSetFromOptions(adaptor, "-adaptor_monitor_error", "error", adaptor));
518c6f61ee2SBarry Smith PetscOptionsEnd();
519c6f61ee2SBarry Smith PetscCall(VecTaggerSetFromOptions(adaptor->refineTag));
520c6f61ee2SBarry Smith PetscCall(VecTaggerSetFromOptions(adaptor->coarsenTag));
521c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
522c6f61ee2SBarry Smith }
523c6f61ee2SBarry Smith
524c6f61ee2SBarry Smith /*@
525c6f61ee2SBarry Smith DMAdaptorView - Views a `DMAdaptor` object
526c6f61ee2SBarry Smith
527c6f61ee2SBarry Smith Collective
528c6f61ee2SBarry Smith
529c6f61ee2SBarry Smith Input Parameters:
530c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
531c6f61ee2SBarry Smith - viewer - The `PetscViewer` object
532c6f61ee2SBarry Smith
533c6f61ee2SBarry Smith Level: beginner
534c6f61ee2SBarry Smith
535df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
536c6f61ee2SBarry Smith @*/
DMAdaptorView(DMAdaptor adaptor,PetscViewer viewer)537c6f61ee2SBarry Smith PetscErrorCode DMAdaptorView(DMAdaptor adaptor, PetscViewer viewer)
538c6f61ee2SBarry Smith {
539c6f61ee2SBarry Smith PetscFunctionBegin;
540c6f61ee2SBarry Smith PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adaptor, viewer));
541c6f61ee2SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "DM Adaptor\n"));
542c6f61ee2SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " sequence length: %" PetscInt_FMT "\n", adaptor->numSeq));
543c6f61ee2SBarry Smith PetscCall(VecTaggerView(adaptor->refineTag, viewer));
544c6f61ee2SBarry Smith PetscCall(VecTaggerView(adaptor->coarsenTag, viewer));
545c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
546c6f61ee2SBarry Smith }
547c6f61ee2SBarry Smith
548c6f61ee2SBarry Smith /*@
549c6f61ee2SBarry Smith DMAdaptorGetSolver - Gets the solver used to produce discrete solutions
550c6f61ee2SBarry Smith
551c6f61ee2SBarry Smith Not Collective
552c6f61ee2SBarry Smith
553c6f61ee2SBarry Smith Input Parameter:
554c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
555c6f61ee2SBarry Smith
556c6f61ee2SBarry Smith Output Parameter:
557c6f61ee2SBarry Smith . snes - The solver
558c6f61ee2SBarry Smith
559c6f61ee2SBarry Smith Level: intermediate
560c6f61ee2SBarry Smith
561df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
562c6f61ee2SBarry Smith @*/
DMAdaptorGetSolver(DMAdaptor adaptor,SNES * snes)563c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetSolver(DMAdaptor adaptor, SNES *snes)
564c6f61ee2SBarry Smith {
565c6f61ee2SBarry Smith PetscFunctionBegin;
5663a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
567c6f61ee2SBarry Smith PetscAssertPointer(snes, 2);
568c6f61ee2SBarry Smith *snes = adaptor->snes;
569c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
570c6f61ee2SBarry Smith }
571c6f61ee2SBarry Smith
572c6f61ee2SBarry Smith /*@
573c6f61ee2SBarry Smith DMAdaptorSetSolver - Sets the solver used to produce discrete solutions
574c6f61ee2SBarry Smith
575c6f61ee2SBarry Smith Not Collective
576c6f61ee2SBarry Smith
577c6f61ee2SBarry Smith Input Parameters:
578c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
579c6f61ee2SBarry Smith - snes - The solver, this MUST have an attached `DM`/`PetscDS`, so that the exact solution can be computed
580c6f61ee2SBarry Smith
581c6f61ee2SBarry Smith Level: intermediate
582c6f61ee2SBarry Smith
583df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorGetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
584c6f61ee2SBarry Smith @*/
DMAdaptorSetSolver(DMAdaptor adaptor,SNES snes)585c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetSolver(DMAdaptor adaptor, SNES snes)
586c6f61ee2SBarry Smith {
587c6f61ee2SBarry Smith PetscFunctionBegin;
5883a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
589c6f61ee2SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
590c6f61ee2SBarry Smith adaptor->snes = snes;
591c6f61ee2SBarry Smith PetscCall(SNESGetDM(adaptor->snes, &adaptor->idm));
592c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
593c6f61ee2SBarry Smith }
594c6f61ee2SBarry Smith
595c6f61ee2SBarry Smith /*@
596c6f61ee2SBarry Smith DMAdaptorGetSequenceLength - Gets the number of sequential adaptations used by an adapter
597c6f61ee2SBarry Smith
598c6f61ee2SBarry Smith Not Collective
599c6f61ee2SBarry Smith
600c6f61ee2SBarry Smith Input Parameter:
601c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
602c6f61ee2SBarry Smith
603c6f61ee2SBarry Smith Output Parameter:
604c6f61ee2SBarry Smith . num - The number of adaptations
605c6f61ee2SBarry Smith
606c6f61ee2SBarry Smith Level: intermediate
607c6f61ee2SBarry Smith
608df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorSetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
609c6f61ee2SBarry Smith @*/
DMAdaptorGetSequenceLength(DMAdaptor adaptor,PetscInt * num)610c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetSequenceLength(DMAdaptor adaptor, PetscInt *num)
611c6f61ee2SBarry Smith {
612c6f61ee2SBarry Smith PetscFunctionBegin;
6133a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
614c6f61ee2SBarry Smith PetscAssertPointer(num, 2);
615c6f61ee2SBarry Smith *num = adaptor->numSeq;
616c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
617c6f61ee2SBarry Smith }
618c6f61ee2SBarry Smith
619c6f61ee2SBarry Smith /*@
620c6f61ee2SBarry Smith DMAdaptorSetSequenceLength - Sets the number of sequential adaptations
621c6f61ee2SBarry Smith
622c6f61ee2SBarry Smith Not Collective
623c6f61ee2SBarry Smith
624c6f61ee2SBarry Smith Input Parameters:
625c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
626c6f61ee2SBarry Smith - num - The number of adaptations
627c6f61ee2SBarry Smith
628c6f61ee2SBarry Smith Level: intermediate
629c6f61ee2SBarry Smith
630df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptorGetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
631c6f61ee2SBarry Smith @*/
DMAdaptorSetSequenceLength(DMAdaptor adaptor,PetscInt num)632c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetSequenceLength(DMAdaptor adaptor, PetscInt num)
633c6f61ee2SBarry Smith {
634c6f61ee2SBarry Smith PetscFunctionBegin;
6353a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
636c6f61ee2SBarry Smith adaptor->numSeq = num;
637c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
638c6f61ee2SBarry Smith }
639c6f61ee2SBarry Smith
DMAdaptorTransferSolution_Exact_Private(DMAdaptor adaptor,DM dm,Vec u,DM adm,Vec au,PetscCtx ctx)640*2a8381b2SBarry Smith static PetscErrorCode DMAdaptorTransferSolution_Exact_Private(DMAdaptor adaptor, DM dm, Vec u, DM adm, Vec au, PetscCtx ctx)
6413a336bb1SMatthew G. Knepley {
6423a336bb1SMatthew G. Knepley PetscFunctionBeginUser;
6433a336bb1SMatthew G. Knepley PetscCall(DMProjectFunction(adm, 0.0, adaptor->exactSol, adaptor->exactCtx, INSERT_ALL_VALUES, au));
6443a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
6453a336bb1SMatthew G. Knepley }
6463a336bb1SMatthew G. Knepley
647c6f61ee2SBarry Smith /*@
648c6f61ee2SBarry Smith DMAdaptorSetUp - After the solver is specified, creates data structures for controlling adaptivity
649c6f61ee2SBarry Smith
650c6f61ee2SBarry Smith Collective
651c6f61ee2SBarry Smith
652c6f61ee2SBarry Smith Input Parameter:
653c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object
654c6f61ee2SBarry Smith
655c6f61ee2SBarry Smith Level: beginner
656c6f61ee2SBarry Smith
657df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()`
658c6f61ee2SBarry Smith @*/
DMAdaptorSetUp(DMAdaptor adaptor)659c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetUp(DMAdaptor adaptor)
660c6f61ee2SBarry Smith {
661c6f61ee2SBarry Smith PetscDS prob;
662c6f61ee2SBarry Smith PetscInt Nf, f;
663c6f61ee2SBarry Smith
664c6f61ee2SBarry Smith PetscFunctionBegin;
665c6f61ee2SBarry Smith PetscCall(DMGetDS(adaptor->idm, &prob));
666c6f61ee2SBarry Smith PetscCall(VecTaggerSetUp(adaptor->refineTag));
667c6f61ee2SBarry Smith PetscCall(VecTaggerSetUp(adaptor->coarsenTag));
668c6f61ee2SBarry Smith PetscCall(PetscDSGetNumFields(prob, &Nf));
669c6f61ee2SBarry Smith PetscCall(PetscMalloc2(Nf, &adaptor->exactSol, Nf, &adaptor->exactCtx));
670c6f61ee2SBarry Smith for (f = 0; f < Nf; ++f) {
671c6f61ee2SBarry Smith PetscCall(PetscDSGetExactSolution(prob, f, &adaptor->exactSol[f], &adaptor->exactCtx[f]));
672c6f61ee2SBarry Smith /* TODO Have a flag that forces projection rather than using the exact solution */
673c6f61ee2SBarry Smith if (adaptor->exactSol[0]) PetscCall(DMAdaptorSetTransferFunction(adaptor, DMAdaptorTransferSolution_Exact_Private));
674c6f61ee2SBarry Smith }
6753a336bb1SMatthew G. Knepley PetscTryTypeMethod(adaptor, setup);
676c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
677c6f61ee2SBarry Smith }
678c6f61ee2SBarry Smith
DMAdaptorGetTransferFunction(DMAdaptor adaptor,PetscErrorCode (** tfunc)(DMAdaptor,DM,Vec,DM,Vec,void *))679c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetTransferFunction(DMAdaptor adaptor, PetscErrorCode (**tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
680c6f61ee2SBarry Smith {
681c6f61ee2SBarry Smith PetscFunctionBegin;
682c6f61ee2SBarry Smith *tfunc = adaptor->ops->transfersolution;
683c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
684c6f61ee2SBarry Smith }
685c6f61ee2SBarry Smith
DMAdaptorSetTransferFunction(DMAdaptor adaptor,PetscErrorCode (* tfunc)(DMAdaptor,DM,Vec,DM,Vec,void *))686c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetTransferFunction(DMAdaptor adaptor, PetscErrorCode (*tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *))
687c6f61ee2SBarry Smith {
688c6f61ee2SBarry Smith PetscFunctionBegin;
689c6f61ee2SBarry Smith adaptor->ops->transfersolution = tfunc;
690c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
691c6f61ee2SBarry Smith }
692c6f61ee2SBarry Smith
DMAdaptorPreAdapt(DMAdaptor adaptor,Vec locX)693c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorPreAdapt(DMAdaptor adaptor, Vec locX)
694c6f61ee2SBarry Smith {
695c6f61ee2SBarry Smith DM plex;
696c6f61ee2SBarry Smith PetscDS prob;
697c6f61ee2SBarry Smith PetscObject obj;
698c6f61ee2SBarry Smith PetscClassId id;
699c6f61ee2SBarry Smith PetscBool isForest;
700c6f61ee2SBarry Smith
701c6f61ee2SBarry Smith PetscFunctionBegin;
702c6f61ee2SBarry Smith PetscCall(DMConvert(adaptor->idm, DMPLEX, &plex));
703c6f61ee2SBarry Smith PetscCall(DMGetDS(adaptor->idm, &prob));
704c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
705c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id));
706c6f61ee2SBarry Smith PetscCall(DMIsForest(adaptor->idm, &isForest));
707c6f61ee2SBarry Smith if (adaptor->adaptCriterion == DM_ADAPTATION_NONE) {
7083a7d0413SPierre Jolivet if (isForest) adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
709c6f61ee2SBarry Smith #if defined(PETSC_HAVE_PRAGMATIC)
710c6f61ee2SBarry Smith else {
711c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
712c6f61ee2SBarry Smith }
713c6f61ee2SBarry Smith #elif defined(PETSC_HAVE_MMG)
714c6f61ee2SBarry Smith else {
715c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
716c6f61ee2SBarry Smith }
717c6f61ee2SBarry Smith #elif defined(PETSC_HAVE_PARMMG)
718c6f61ee2SBarry Smith else {
719c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_METRIC;
720c6f61ee2SBarry Smith }
721c6f61ee2SBarry Smith #else
722c6f61ee2SBarry Smith else {
7233a336bb1SMatthew G. Knepley adaptor->adaptCriterion = DM_ADAPTATION_LABEL;
724c6f61ee2SBarry Smith }
725c6f61ee2SBarry Smith #endif
726c6f61ee2SBarry Smith }
727c6f61ee2SBarry Smith if (id == PETSCFV_CLASSID) {
728c6f61ee2SBarry Smith adaptor->femType = PETSC_FALSE;
729c6f61ee2SBarry Smith } else {
730c6f61ee2SBarry Smith adaptor->femType = PETSC_TRUE;
731c6f61ee2SBarry Smith }
732c6f61ee2SBarry Smith if (adaptor->femType) {
733c6f61ee2SBarry Smith /* Compute local solution bc */
734c6f61ee2SBarry Smith PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
735c6f61ee2SBarry Smith } else {
736c6f61ee2SBarry Smith PetscFV fvm = (PetscFV)obj;
737c6f61ee2SBarry Smith PetscLimiter noneLimiter;
738c6f61ee2SBarry Smith Vec grad;
739c6f61ee2SBarry Smith
740c6f61ee2SBarry Smith PetscCall(PetscFVGetComputeGradients(fvm, &adaptor->computeGradient));
741c6f61ee2SBarry Smith PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
742c6f61ee2SBarry Smith /* Use no limiting when reconstructing gradients for adaptivity */
743c6f61ee2SBarry Smith PetscCall(PetscFVGetLimiter(fvm, &adaptor->limiter));
744c6f61ee2SBarry Smith PetscCall(PetscObjectReference((PetscObject)adaptor->limiter));
745c6f61ee2SBarry Smith PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter));
746c6f61ee2SBarry Smith PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE));
747c6f61ee2SBarry Smith PetscCall(PetscFVSetLimiter(fvm, noneLimiter));
748c6f61ee2SBarry Smith /* Get FVM data */
749c6f61ee2SBarry Smith PetscCall(DMPlexGetDataFVM(plex, fvm, &adaptor->cellGeom, &adaptor->faceGeom, &adaptor->gradDM));
750c6f61ee2SBarry Smith PetscCall(VecGetDM(adaptor->cellGeom, &adaptor->cellDM));
751c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
752c6f61ee2SBarry Smith /* Compute local solution bc */
753c6f61ee2SBarry Smith PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL));
754c6f61ee2SBarry Smith /* Compute gradients */
755c6f61ee2SBarry Smith PetscCall(DMCreateGlobalVector(adaptor->gradDM, &grad));
756c6f61ee2SBarry Smith PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad));
757c6f61ee2SBarry Smith PetscCall(DMGetLocalVector(adaptor->gradDM, &adaptor->cellGrad));
758c6f61ee2SBarry Smith PetscCall(DMGlobalToLocalBegin(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
759c6f61ee2SBarry Smith PetscCall(DMGlobalToLocalEnd(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad));
760c6f61ee2SBarry Smith PetscCall(VecDestroy(&grad));
761c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
762c6f61ee2SBarry Smith }
763c6f61ee2SBarry Smith PetscCall(DMDestroy(&plex));
764c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
765c6f61ee2SBarry Smith }
766c6f61ee2SBarry Smith
DMAdaptorTransferSolution(DMAdaptor adaptor,DM dm,Vec x,DM adm,Vec ax)767c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorTransferSolution(DMAdaptor adaptor, DM dm, Vec x, DM adm, Vec ax)
768c6f61ee2SBarry Smith {
769c6f61ee2SBarry Smith PetscReal time = 0.0;
770c6f61ee2SBarry Smith Mat interp;
771c6f61ee2SBarry Smith void *ctx;
772c6f61ee2SBarry Smith
773c6f61ee2SBarry Smith PetscFunctionBegin;
774c6f61ee2SBarry Smith PetscCall(DMGetApplicationContext(dm, &ctx));
775c6f61ee2SBarry Smith if (adaptor->ops->transfersolution) PetscUseTypeMethod(adaptor, transfersolution, dm, x, adm, ax, ctx);
776c6f61ee2SBarry Smith else {
777c6f61ee2SBarry Smith switch (adaptor->adaptCriterion) {
778c6f61ee2SBarry Smith case DM_ADAPTATION_LABEL:
779c6f61ee2SBarry Smith PetscCall(DMForestTransferVec(dm, x, adm, ax, PETSC_TRUE, time));
780c6f61ee2SBarry Smith break;
781c6f61ee2SBarry Smith case DM_ADAPTATION_REFINE:
782c6f61ee2SBarry Smith case DM_ADAPTATION_METRIC:
783c6f61ee2SBarry Smith PetscCall(DMCreateInterpolation(dm, adm, &interp, NULL));
784c6f61ee2SBarry Smith PetscCall(MatInterpolate(interp, x, ax));
785c6f61ee2SBarry Smith PetscCall(DMInterpolate(dm, interp, adm));
786c6f61ee2SBarry Smith PetscCall(MatDestroy(&interp));
787c6f61ee2SBarry Smith break;
788c6f61ee2SBarry Smith default:
789c6f61ee2SBarry Smith SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_SUP, "No built-in projection for this adaptation criterion: %d", adaptor->adaptCriterion);
790c6f61ee2SBarry Smith }
791c6f61ee2SBarry Smith }
792c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
793c6f61ee2SBarry Smith }
794c6f61ee2SBarry Smith
DMAdaptorPostAdapt(DMAdaptor adaptor)795c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorPostAdapt(DMAdaptor adaptor)
796c6f61ee2SBarry Smith {
797c6f61ee2SBarry Smith PetscDS prob;
798c6f61ee2SBarry Smith PetscObject obj;
799c6f61ee2SBarry Smith PetscClassId id;
800c6f61ee2SBarry Smith
801c6f61ee2SBarry Smith PetscFunctionBegin;
802c6f61ee2SBarry Smith PetscCall(DMGetDS(adaptor->idm, &prob));
803c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(prob, 0, &obj));
804c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id));
805c6f61ee2SBarry Smith if (id == PETSCFV_CLASSID) {
806c6f61ee2SBarry Smith PetscFV fvm = (PetscFV)obj;
807c6f61ee2SBarry Smith
808c6f61ee2SBarry Smith PetscCall(PetscFVSetComputeGradients(fvm, adaptor->computeGradient));
809c6f61ee2SBarry Smith /* Restore original limiter */
810c6f61ee2SBarry Smith PetscCall(PetscFVSetLimiter(fvm, adaptor->limiter));
811c6f61ee2SBarry Smith
812c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray));
813c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(adaptor->cellGrad, &adaptor->cellGradArray));
814c6f61ee2SBarry Smith PetscCall(DMRestoreLocalVector(adaptor->gradDM, &adaptor->cellGrad));
815c6f61ee2SBarry Smith }
816c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
817c6f61ee2SBarry Smith }
818c6f61ee2SBarry Smith
819c6f61ee2SBarry Smith /*
8203a336bb1SMatthew G. Knepley DMAdaptorComputeCellErrorIndicator_Gradient - Use the integrated gradient as an error indicator in the `DMAdaptor`
821c6f61ee2SBarry Smith
822c6f61ee2SBarry Smith Input Parameters:
823c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
824c6f61ee2SBarry Smith . dim - The topological dimension
825c6f61ee2SBarry Smith . cell - The cell
826c6f61ee2SBarry Smith . field - The field integrated over the cell
827c6f61ee2SBarry Smith . gradient - The gradient integrated over the cell
828c6f61ee2SBarry Smith . cg - A `PetscFVCellGeom` struct
829c6f61ee2SBarry Smith - ctx - A user context
830c6f61ee2SBarry Smith
831c6f61ee2SBarry Smith Output Parameter:
832c6f61ee2SBarry Smith . errInd - The error indicator
833c6f61ee2SBarry Smith
834c6f61ee2SBarry Smith Developer Note:
835c6f61ee2SBarry Smith Some of the input arguments are absurdly specialized to special situations, it is not clear this is a good general API
836c6f61ee2SBarry Smith
8373a336bb1SMatthew G. Knepley .seealso: [](ch_dmbase), `DMAdaptor`
838c6f61ee2SBarry Smith */
DMAdaptorComputeCellErrorIndicator_Gradient(DMAdaptor adaptor,PetscInt dim,PetscInt Nc,const PetscScalar * field,const PetscScalar * gradient,const PetscFVCellGeom * cg,PetscReal * errInd,PetscCtx ctx)839*2a8381b2SBarry Smith static PetscErrorCode DMAdaptorComputeCellErrorIndicator_Gradient(DMAdaptor adaptor, PetscInt dim, PetscInt Nc, const PetscScalar *field, const PetscScalar *gradient, const PetscFVCellGeom *cg, PetscReal *errInd, PetscCtx ctx)
840c6f61ee2SBarry Smith {
841c6f61ee2SBarry Smith PetscReal err = 0.;
842c6f61ee2SBarry Smith PetscInt c, d;
843c6f61ee2SBarry Smith
844c6f61ee2SBarry Smith PetscFunctionBeginHot;
845c6f61ee2SBarry Smith for (c = 0; c < Nc; c++) {
846c6f61ee2SBarry Smith for (d = 0; d < dim; ++d) err += PetscSqr(PetscRealPart(gradient[c * dim + d]));
847c6f61ee2SBarry Smith }
848c6f61ee2SBarry Smith *errInd = cg->volume * err;
849c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
850c6f61ee2SBarry Smith }
851c6f61ee2SBarry Smith
DMAdaptorComputeErrorIndicator_Gradient(DMAdaptor adaptor,Vec locX,Vec errVec)8523a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorComputeErrorIndicator_Gradient(DMAdaptor adaptor, Vec locX, Vec errVec)
853c6f61ee2SBarry Smith {
8543a336bb1SMatthew G. Knepley DM dm, plex, edm, eplex;
8553a336bb1SMatthew G. Knepley PetscDS ds;
856c6f61ee2SBarry Smith PetscObject obj;
857c6f61ee2SBarry Smith PetscClassId id;
858c6f61ee2SBarry Smith void *ctx;
859c6f61ee2SBarry Smith PetscQuadrature quad;
8603a336bb1SMatthew G. Knepley PetscScalar *earray;
8613a336bb1SMatthew G. Knepley PetscReal minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2];
8623a336bb1SMatthew G. Knepley PetscInt dim, cdim, cStart, cEnd, Nf, Nc;
863c6f61ee2SBarry Smith
864c6f61ee2SBarry Smith PetscFunctionBegin;
8653a336bb1SMatthew G. Knepley PetscCall(VecGetDM(locX, &dm));
8663a336bb1SMatthew G. Knepley PetscCall(DMConvert(dm, DMPLEX, &plex));
8673a336bb1SMatthew G. Knepley PetscCall(VecGetDM(errVec, &edm));
8683a336bb1SMatthew G. Knepley PetscCall(DMConvert(edm, DMPLEX, &eplex));
869c6f61ee2SBarry Smith PetscCall(DMGetDimension(plex, &dim));
870c6f61ee2SBarry Smith PetscCall(DMGetCoordinateDim(plex, &cdim));
871c6f61ee2SBarry Smith PetscCall(DMGetApplicationContext(plex, &ctx));
8723a336bb1SMatthew G. Knepley PetscCall(DMGetDS(plex, &ds));
8733a336bb1SMatthew G. Knepley PetscCall(PetscDSGetNumFields(ds, &Nf));
8743a336bb1SMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, 0, &obj));
875c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id));
8763a336bb1SMatthew G. Knepley
8773a336bb1SMatthew G. Knepley PetscCall(VecGetArray(errVec, &earray));
8783a336bb1SMatthew G. Knepley PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd));
8793a336bb1SMatthew G. Knepley for (PetscInt cell = cStart; cell < cEnd; ++cell) {
8803a336bb1SMatthew G. Knepley PetscScalar *eval;
8813a336bb1SMatthew G. Knepley PetscReal errInd = 0.;
8823a336bb1SMatthew G. Knepley
883c6f61ee2SBarry Smith if (id == PETSCFV_CLASSID) {
8843a336bb1SMatthew G. Knepley PetscFV fv = (PetscFV)obj;
885c6f61ee2SBarry Smith const PetscScalar *pointSols;
886c6f61ee2SBarry Smith const PetscScalar *pointSol;
887c6f61ee2SBarry Smith const PetscScalar *pointGrad;
888c6f61ee2SBarry Smith PetscFVCellGeom *cg;
889c6f61ee2SBarry Smith
8903a336bb1SMatthew G. Knepley PetscCall(PetscFVGetNumComponents(fv, &Nc));
891c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(locX, &pointSols));
892c6f61ee2SBarry Smith PetscCall(DMPlexPointLocalRead(plex, cell, pointSols, (void *)&pointSol));
893c6f61ee2SBarry Smith PetscCall(DMPlexPointLocalRead(adaptor->gradDM, cell, adaptor->cellGradArray, (void *)&pointGrad));
894c6f61ee2SBarry Smith PetscCall(DMPlexPointLocalRead(adaptor->cellDM, cell, adaptor->cellGeomArray, &cg));
8953a336bb1SMatthew G. Knepley PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, pointSol, pointGrad, cg, &errInd, ctx);
896c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(locX, &pointSols));
897c6f61ee2SBarry Smith } else {
8983a336bb1SMatthew G. Knepley PetscFE fe = (PetscFE)obj;
899c6f61ee2SBarry Smith PetscScalar *x = NULL, *field, *gradient, *interpolant, *interpolantGrad;
900c6f61ee2SBarry Smith PetscFVCellGeom cg;
901c6f61ee2SBarry Smith PetscFEGeom fegeom;
902c6f61ee2SBarry Smith const PetscReal *quadWeights;
903c6f61ee2SBarry Smith PetscReal *coords;
9043a336bb1SMatthew G. Knepley PetscInt Nb, Nq, qNc;
905c6f61ee2SBarry Smith
906c6f61ee2SBarry Smith fegeom.dim = dim;
907c6f61ee2SBarry Smith fegeom.dimEmbed = cdim;
9083a336bb1SMatthew G. Knepley PetscCall(PetscFEGetNumComponents(fe, &Nc));
9093a336bb1SMatthew G. Knepley PetscCall(PetscFEGetQuadrature(fe, &quad));
9103a336bb1SMatthew G. Knepley PetscCall(PetscFEGetDimension(fe, &Nb));
911c6f61ee2SBarry Smith PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
912c6f61ee2SBarry Smith PetscCall(PetscCalloc6(Nc, &field, cdim * Nc, &gradient, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ));
913c6f61ee2SBarry Smith PetscCall(PetscMalloc2(Nc, &interpolant, cdim * Nc, &interpolantGrad));
914c6f61ee2SBarry Smith PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
915c6f61ee2SBarry Smith PetscCall(DMPlexComputeCellGeometryFVM(plex, cell, &cg.volume, NULL, NULL));
916c6f61ee2SBarry Smith PetscCall(PetscArrayzero(gradient, cdim * Nc));
9173a336bb1SMatthew G. Knepley PetscCall(DMPlexVecGetClosure(plex, NULL, locX, cell, NULL, &x));
9183a336bb1SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) {
9193a336bb1SMatthew G. Knepley PetscInt qc = 0;
920c6f61ee2SBarry Smith
9213a336bb1SMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, f, &obj));
922c6f61ee2SBarry Smith PetscCall(PetscArrayzero(interpolant, Nc));
923c6f61ee2SBarry Smith PetscCall(PetscArrayzero(interpolantGrad, cdim * Nc));
9243a336bb1SMatthew G. Knepley for (PetscInt q = 0; q < Nq; ++q) {
925c6f61ee2SBarry Smith PetscCall(PetscFEInterpolateFieldAndGradient_Static((PetscFE)obj, 1, x, &fegeom, q, interpolant, interpolantGrad));
9263a336bb1SMatthew G. Knepley for (PetscInt fc = 0; fc < Nc; ++fc) {
927c6f61ee2SBarry Smith const PetscReal wt = quadWeights[q * qNc + qc + fc];
928c6f61ee2SBarry Smith
929c6f61ee2SBarry Smith field[fc] += interpolant[fc] * wt * fegeom.detJ[q];
9303a336bb1SMatthew G. Knepley for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] += interpolantGrad[fc * dim + d] * wt * fegeom.detJ[q];
931c6f61ee2SBarry Smith }
932c6f61ee2SBarry Smith }
933c6f61ee2SBarry Smith qc += Nc;
934c6f61ee2SBarry Smith }
935c6f61ee2SBarry Smith PetscCall(PetscFree2(interpolant, interpolantGrad));
936c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(plex, NULL, locX, cell, NULL, &x));
9373a336bb1SMatthew G. Knepley for (PetscInt fc = 0; fc < Nc; ++fc) {
938c6f61ee2SBarry Smith field[fc] /= cg.volume;
9393a336bb1SMatthew G. Knepley for (PetscInt d = 0; d < cdim; ++d) gradient[fc * cdim + d] /= cg.volume;
940c6f61ee2SBarry Smith }
9413a336bb1SMatthew G. Knepley PetscUseTypeMethod(adaptor, computecellerrorindicator, dim, Nc, field, gradient, &cg, &errInd, ctx);
942c6f61ee2SBarry Smith PetscCall(PetscFree6(field, gradient, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
943c6f61ee2SBarry Smith }
9443a336bb1SMatthew G. Knepley PetscCall(DMPlexPointGlobalRef(eplex, cell, earray, (void *)&eval));
9453a336bb1SMatthew G. Knepley eval[0] = errInd;
9463a336bb1SMatthew G. Knepley minMaxInd[0] = PetscMin(minMaxInd[0], errInd);
9473a336bb1SMatthew G. Knepley minMaxInd[1] = PetscMax(minMaxInd[1], errInd);
9483a336bb1SMatthew G. Knepley }
9493a336bb1SMatthew G. Knepley PetscCall(VecRestoreArray(errVec, &earray));
9503a336bb1SMatthew G. Knepley PetscCall(DMDestroy(&plex));
9513a336bb1SMatthew G. Knepley PetscCall(DMDestroy(&eplex));
9523a336bb1SMatthew G. Knepley PetscCall(PetscGlobalMinMaxReal(PetscObjectComm((PetscObject)adaptor), minMaxInd, minMaxIndGlobal));
9533a336bb1SMatthew G. Knepley PetscCall(PetscInfo(adaptor, "DMAdaptor: error indicator range (%g, %g)\n", (double)minMaxIndGlobal[0], (double)minMaxIndGlobal[1]));
9543a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
9553a336bb1SMatthew G. Knepley }
9563a336bb1SMatthew G. Knepley
DMAdaptorComputeErrorIndicator_Flux(DMAdaptor adaptor,Vec lu,Vec errVec)9573a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorComputeErrorIndicator_Flux(DMAdaptor adaptor, Vec lu, Vec errVec)
9583a336bb1SMatthew G. Knepley {
9593a336bb1SMatthew G. Knepley DM dm, mdm;
9603a336bb1SMatthew G. Knepley SNES msnes;
9613a336bb1SMatthew G. Knepley Vec mu, lmu;
9623a336bb1SMatthew G. Knepley void *ctx;
963e03fd340SMatthew G. Knepley const char *prefix;
9643a336bb1SMatthew G. Knepley
9653a336bb1SMatthew G. Knepley PetscFunctionBegin;
9663a336bb1SMatthew G. Knepley PetscCall(VecGetDM(lu, &dm));
9673a336bb1SMatthew G. Knepley
9683a336bb1SMatthew G. Knepley // Set up and solve mixed problem
9693a336bb1SMatthew G. Knepley PetscCall(DMClone(dm, &mdm));
9703a336bb1SMatthew G. Knepley PetscCall(SNESCreate(PetscObjectComm((PetscObject)mdm), &msnes));
9713a336bb1SMatthew G. Knepley PetscCall(SNESSetDM(msnes, mdm));
972e03fd340SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
973e03fd340SMatthew G. Knepley PetscCall(SNESSetOptionsPrefix(msnes, prefix));
9743a336bb1SMatthew G. Knepley PetscCall(SNESAppendOptionsPrefix(msnes, "mixed_"));
9753a336bb1SMatthew G. Knepley
9763a336bb1SMatthew G. Knepley PetscTryTypeMethod(adaptor, mixedsetup, mdm);
9773a336bb1SMatthew G. Knepley PetscCall(DMGetApplicationContext(dm, &ctx));
9783a336bb1SMatthew G. Knepley PetscCall(DMPlexSetSNESLocalFEM(mdm, PETSC_FALSE, ctx));
9793a336bb1SMatthew G. Knepley PetscCall(SNESSetFromOptions(msnes));
9803a336bb1SMatthew G. Knepley
9813a336bb1SMatthew G. Knepley PetscCall(DMCreateGlobalVector(mdm, &mu));
9823a336bb1SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)mu, "Mixed Solution"));
9833a336bb1SMatthew G. Knepley PetscCall(VecSet(mu, 0.0));
9843a336bb1SMatthew G. Knepley PetscCall(SNESSolve(msnes, NULL, mu));
9853a336bb1SMatthew G. Knepley PetscCall(VecViewFromOptions(mu, (PetscObject)adaptor, "-adapt_mixed_sol_vec_view"));
9863a336bb1SMatthew G. Knepley
9873a336bb1SMatthew G. Knepley PetscCall(DMGetLocalVector(mdm, &lmu));
9883a336bb1SMatthew G. Knepley PetscCall(DMGlobalToLocal(mdm, mu, INSERT_VALUES, lmu));
9893a336bb1SMatthew G. Knepley PetscCall(DMPlexInsertBoundaryValues(mdm, PETSC_TRUE, lmu, 0.0, NULL, NULL, NULL));
9903a336bb1SMatthew G. Knepley PetscCall(DMPlexComputeL2FluxDiffVecLocal(lu, 0, lmu, 0, errVec));
9913a336bb1SMatthew G. Knepley PetscCall(DMRestoreLocalVector(mdm, &lmu));
9923a336bb1SMatthew G. Knepley PetscCall(VecDestroy(&mu));
9933a336bb1SMatthew G. Knepley PetscCall(SNESDestroy(&msnes));
9943a336bb1SMatthew G. Knepley PetscCall(DMDestroy(&mdm));
995c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
996c6f61ee2SBarry Smith }
997c6f61ee2SBarry Smith
9988b724c91SMatthew G. Knepley /*@
9998b724c91SMatthew G. Knepley DMAdaptorMonitor - runs the user provided monitor routines, if they exist
10008b724c91SMatthew G. Knepley
10018b724c91SMatthew G. Knepley Collective
10028b724c91SMatthew G. Knepley
10038b724c91SMatthew G. Knepley Input Parameters:
10048b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
10058b724c91SMatthew G. Knepley . it - iteration number
10068b724c91SMatthew G. Knepley . odm - the original `DM`
10078b724c91SMatthew G. Knepley . adm - the adapted `DM`
10088b724c91SMatthew G. Knepley . Nf - the number of fields
10098b724c91SMatthew G. Knepley . enorms - the 2-norm error values for each field
10108b724c91SMatthew G. Knepley - error - `Vec` of cellwise errors
10118b724c91SMatthew G. Knepley
10128b724c91SMatthew G. Knepley Level: developer
10138b724c91SMatthew G. Knepley
10148b724c91SMatthew G. Knepley Note:
10158b724c91SMatthew G. Knepley This routine is called by the `DMAdaptor` implementations.
10168b724c91SMatthew G. Knepley It does not typically need to be called by the user.
10178b724c91SMatthew G. Knepley
10188b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptorMonitorSet()`
10198b724c91SMatthew G. Knepley @*/
DMAdaptorMonitor(DMAdaptor adaptor,PetscInt it,DM odm,DM adm,PetscInt Nf,PetscReal enorms[],Vec error)10208b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitor(DMAdaptor adaptor, PetscInt it, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error)
10218b724c91SMatthew G. Knepley {
10228b724c91SMatthew G. Knepley PetscFunctionBegin;
10238b724c91SMatthew G. Knepley for (PetscInt i = 0; i < adaptor->numbermonitors; ++i) PetscCall((*adaptor->monitor[i])(adaptor, it, odm, adm, Nf, enorms, error, adaptor->monitorcontext[i]));
10248b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
10258b724c91SMatthew G. Knepley }
10268b724c91SMatthew G. Knepley
10278b724c91SMatthew G. Knepley /*@C
10288b724c91SMatthew G. Knepley DMAdaptorMonitorSize - Prints the mesh sizes at each iteration of an adaptation loop.
10298b724c91SMatthew G. Knepley
10308b724c91SMatthew G. Knepley Collective
10318b724c91SMatthew G. Knepley
10328b724c91SMatthew G. Knepley Input Parameters:
10338b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
10348b724c91SMatthew G. Knepley . n - iteration number
10358b724c91SMatthew G. Knepley . odm - the original `DM`
10368b724c91SMatthew G. Knepley . adm - the adapted `DM`
10378b724c91SMatthew G. Knepley . Nf - number of fields
10388b724c91SMatthew G. Knepley . enorms - 2-norm error values for each field (may be estimated).
10398b724c91SMatthew G. Knepley . error - `Vec` of cellwise errors
10408b724c91SMatthew G. Knepley - vf - The viewer context
10418b724c91SMatthew G. Knepley
10428b724c91SMatthew G. Knepley Options Database Key:
10438b724c91SMatthew G. Knepley . -adaptor_monitor_size - Activates `DMAdaptorMonitorSize()`
10448b724c91SMatthew G. Knepley
10458b724c91SMatthew G. Knepley Level: intermediate
10468b724c91SMatthew G. Knepley
10478b724c91SMatthew G. Knepley Note:
10488b724c91SMatthew G. Knepley This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
10498b724c91SMatthew G. Knepley to be used during the adaptation loop.
10508b724c91SMatthew G. Knepley
10518b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorError()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
10528b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorSize(DMAdaptor adaptor,PetscInt n,DM odm,DM adm,PetscInt Nf,PetscReal enorms[],Vec error,PetscViewerAndFormat * vf)10538b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorSize(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
10548b724c91SMatthew G. Knepley {
10558b724c91SMatthew G. Knepley PetscViewer viewer = vf->viewer;
10568b724c91SMatthew G. Knepley PetscViewerFormat format = vf->format;
10578b724c91SMatthew G. Knepley PetscInt tablevel, cStart, cEnd, acStart, acEnd;
10588b724c91SMatthew G. Knepley const char *prefix;
10598b724c91SMatthew G. Knepley PetscMPIInt rank;
10608b724c91SMatthew G. Knepley
10618b724c91SMatthew G. Knepley PetscFunctionBegin;
10628b724c91SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
10638b724c91SMatthew G. Knepley PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
10648b724c91SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
10658b724c91SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)adaptor), &rank));
10668b724c91SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
10678b724c91SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
10688b724c91SMatthew G. Knepley
10698b724c91SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format));
10708b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
10718b724c91SMatthew G. Knepley if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Sizes for %s adaptation.\n", prefix));
10728b724c91SMatthew 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));
10738b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
10748b724c91SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer));
10758b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
10768b724c91SMatthew G. Knepley }
10778b724c91SMatthew G. Knepley
10788b724c91SMatthew G. Knepley /*@C
10798b724c91SMatthew G. Knepley DMAdaptorMonitorError - Prints the error norm at each iteration of an adaptation loop.
10808b724c91SMatthew G. Knepley
10818b724c91SMatthew G. Knepley Collective
10828b724c91SMatthew G. Knepley
10838b724c91SMatthew G. Knepley Input Parameters:
10848b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
10858b724c91SMatthew G. Knepley . n - iteration number
10868b724c91SMatthew G. Knepley . odm - the original `DM`
10878b724c91SMatthew G. Knepley . adm - the adapted `DM`
10888b724c91SMatthew G. Knepley . Nf - number of fields
10898b724c91SMatthew G. Knepley . enorms - 2-norm error values for each field (may be estimated).
10908b724c91SMatthew G. Knepley . error - `Vec` of cellwise errors
10918b724c91SMatthew G. Knepley - vf - The viewer context
10928b724c91SMatthew G. Knepley
10938b724c91SMatthew G. Knepley Options Database Key:
10948b724c91SMatthew G. Knepley . -adaptor_monitor_error - Activates `DMAdaptorMonitorError()`
10958b724c91SMatthew G. Knepley
10968b724c91SMatthew G. Knepley Level: intermediate
10978b724c91SMatthew G. Knepley
10988b724c91SMatthew G. Knepley Note:
10998b724c91SMatthew G. Knepley This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
11008b724c91SMatthew G. Knepley to be used during the adaptation loop.
11018b724c91SMatthew G. Knepley
11028b724c91SMatthew G. Knepley .seealso: [](ch_snes), `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorErrorDrawLG()`
11038b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorError(DMAdaptor adaptor,PetscInt n,DM odm,DM adm,PetscInt Nf,PetscReal enorms[],Vec error,PetscViewerAndFormat * vf)11048b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorError(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
11058b724c91SMatthew G. Knepley {
11068b724c91SMatthew G. Knepley PetscViewer viewer = vf->viewer;
11078b724c91SMatthew G. Knepley PetscViewerFormat format = vf->format;
11088b724c91SMatthew G. Knepley PetscInt tablevel, cStart, cEnd, acStart, acEnd;
11098b724c91SMatthew G. Knepley const char *prefix;
11108b724c91SMatthew G. Knepley
11118b724c91SMatthew G. Knepley PetscFunctionBegin;
11128b724c91SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
11138b724c91SMatthew G. Knepley PetscCall(PetscObjectGetTabLevel((PetscObject)adaptor, &tablevel));
11148b724c91SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)adaptor, &prefix));
11158b724c91SMatthew G. Knepley
11168b724c91SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format));
11178b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
11188b724c91SMatthew G. Knepley if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, " Error norms for %s adaptation.\n", prefix));
11198b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " DMAdaptor Error norm %s", n, Nf > 1 ? "[" : ""));
11208b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
11218b724c91SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) {
11228b724c91SMatthew G. Knepley if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
11238b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)enorms[f]));
11248b724c91SMatthew G. Knepley }
11258b724c91SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
11268b724c91SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(adm, 0, &acStart, &acEnd));
11278b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(viewer, " N: %" PetscInt_FMT " Nadapt: %" PetscInt_FMT "\n", cEnd - cStart, acEnd - acStart));
11288b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
11298b724c91SMatthew G. Knepley PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
11308b724c91SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer));
11318b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
11328b724c91SMatthew G. Knepley }
11338b724c91SMatthew G. Knepley
11348b724c91SMatthew G. Knepley /*@C
11358b724c91SMatthew G. Knepley DMAdaptorMonitorErrorDraw - Plots the error at each iteration of an iterative solver.
11368b724c91SMatthew G. Knepley
11378b724c91SMatthew G. Knepley Collective
11388b724c91SMatthew G. Knepley
11398b724c91SMatthew G. Knepley Input Parameters:
11408b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
11418b724c91SMatthew G. Knepley . n - iteration number
11428b724c91SMatthew G. Knepley . odm - the original `DM`
11438b724c91SMatthew G. Knepley . adm - the adapted `DM`
11448b724c91SMatthew G. Knepley . Nf - number of fields
11458b724c91SMatthew G. Knepley . enorms - 2-norm error values for each field (may be estimated).
11468b724c91SMatthew G. Knepley . error - `Vec` of cellwise errors
11478b724c91SMatthew G. Knepley - vf - The viewer context
11488b724c91SMatthew G. Knepley
11498b724c91SMatthew G. Knepley Options Database Key:
11508b724c91SMatthew G. Knepley . -adaptor_monitor_error draw - Activates `DMAdaptorMonitorErrorDraw()`
11518b724c91SMatthew G. Knepley
11528b724c91SMatthew G. Knepley Level: intermediate
11538b724c91SMatthew G. Knepley
11548b724c91SMatthew G. Knepley Note:
11558b724c91SMatthew G. Knepley This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
11568b724c91SMatthew G. Knepley to be used during the adaptation loop.
11578b724c91SMatthew G. Knepley
11588b724c91SMatthew G. Knepley .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
11598b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorErrorDraw(DMAdaptor adaptor,PetscInt n,DM odm,DM adm,PetscInt Nf,PetscReal enorms[],Vec error,PetscViewerAndFormat * vf)11608b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorErrorDraw(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
11618b724c91SMatthew G. Knepley {
11628b724c91SMatthew G. Knepley PetscViewer viewer = vf->viewer;
11638b724c91SMatthew G. Knepley PetscViewerFormat format = vf->format;
11648b724c91SMatthew G. Knepley
11658b724c91SMatthew G. Knepley PetscFunctionBegin;
11668b724c91SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
11678b724c91SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format));
11688b724c91SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
11698b724c91SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", (PetscObject)adaptor));
11708b724c91SMatthew G. Knepley PetscCall(VecView(error, viewer));
11718b724c91SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)error, "__Vec_bc_zero__", NULL));
11728b724c91SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer));
11738b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
11748b724c91SMatthew G. Knepley }
11758b724c91SMatthew G. Knepley
11768b724c91SMatthew G. Knepley /*@C
1177d7c1f440SPierre Jolivet DMAdaptorMonitorErrorDrawLGCreate - Creates the context for the error plotter `DMAdaptorMonitorErrorDrawLG()`
11788b724c91SMatthew G. Knepley
11798b724c91SMatthew G. Knepley Collective
11808b724c91SMatthew G. Knepley
11818b724c91SMatthew G. Knepley Input Parameters:
11828b724c91SMatthew G. Knepley + viewer - The `PetscViewer`
11838b724c91SMatthew G. Knepley . format - The viewer format
11848b724c91SMatthew G. Knepley - ctx - An optional user context
11858b724c91SMatthew G. Knepley
11868b724c91SMatthew G. Knepley Output Parameter:
11878b724c91SMatthew G. Knepley . vf - The viewer context
11888b724c91SMatthew G. Knepley
11898b724c91SMatthew G. Knepley Level: intermediate
11908b724c91SMatthew G. Knepley
1191ce78bad3SBarry Smith .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `PetscViewerMonitorGLSetUp()`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDrawLG()`
11928b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorErrorDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat ** vf)1193*2a8381b2SBarry Smith PetscErrorCode DMAdaptorMonitorErrorDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
11948b724c91SMatthew G. Knepley {
11958b724c91SMatthew G. Knepley DMAdaptor adaptor = (DMAdaptor)ctx;
11968b724c91SMatthew G. Knepley char **names;
11978b724c91SMatthew G. Knepley PetscInt Nf;
11988b724c91SMatthew G. Knepley
11998b724c91SMatthew G. Knepley PetscFunctionBegin;
12008b724c91SMatthew G. Knepley PetscCall(DMGetNumFields(adaptor->idm, &Nf));
12018b724c91SMatthew G. Knepley PetscCall(PetscMalloc1(Nf + 1, &names));
12028b724c91SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) {
12038b724c91SMatthew G. Knepley PetscObject disc;
12048b724c91SMatthew G. Knepley const char *fname;
12058b724c91SMatthew G. Knepley char lname[PETSC_MAX_PATH_LEN];
12068b724c91SMatthew G. Knepley
12078b724c91SMatthew G. Knepley PetscCall(DMGetField(adaptor->idm, f, NULL, &disc));
12088b724c91SMatthew G. Knepley PetscCall(PetscObjectGetName(disc, &fname));
12098b724c91SMatthew G. Knepley PetscCall(PetscStrncpy(lname, fname, PETSC_MAX_PATH_LEN));
12108b724c91SMatthew G. Knepley PetscCall(PetscStrlcat(lname, " Error", PETSC_MAX_PATH_LEN));
12118b724c91SMatthew G. Knepley PetscCall(PetscStrallocpy(lname, &names[f]));
12128b724c91SMatthew G. Knepley }
12138b724c91SMatthew G. Knepley PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
12148b724c91SMatthew G. Knepley (*vf)->data = ctx;
1215ce78bad3SBarry Smith PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Error Norm", Nf, (const char **)names, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
12168b724c91SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(names[f]));
12178b724c91SMatthew G. Knepley PetscCall(PetscFree(names));
12188b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
12198b724c91SMatthew G. Knepley }
12208b724c91SMatthew G. Knepley
12218b724c91SMatthew G. Knepley /*@C
12228b724c91SMatthew G. Knepley DMAdaptorMonitorErrorDrawLG - Plots the error norm at each iteration of an adaptive loop.
12238b724c91SMatthew G. Knepley
12248b724c91SMatthew G. Knepley Collective
12258b724c91SMatthew G. Knepley
12268b724c91SMatthew G. Knepley Input Parameters:
12278b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
12288b724c91SMatthew G. Knepley . n - iteration number
12298b724c91SMatthew G. Knepley . odm - the original `DM`
12308b724c91SMatthew G. Knepley . adm - the adapted `DM`
12318b724c91SMatthew G. Knepley . Nf - number of fields
12328b724c91SMatthew G. Knepley . enorms - 2-norm error values for each field (may be estimated).
12338b724c91SMatthew G. Knepley . error - `Vec` of cellwise errors
1234ce78bad3SBarry Smith - vf - The viewer context, obtained via `DMAdaptorMonitorErrorDrawLGCreate()`
12358b724c91SMatthew G. Knepley
12368b724c91SMatthew G. Knepley Options Database Key:
12378b724c91SMatthew G. Knepley . -adaptor_error draw::draw_lg - Activates `DMAdaptorMonitorErrorDrawLG()`
12388b724c91SMatthew G. Knepley
12398b724c91SMatthew G. Knepley Level: intermediate
12408b724c91SMatthew G. Knepley
12418b724c91SMatthew G. Knepley Notes:
12428b724c91SMatthew G. Knepley This is not called directly by users, rather one calls `DMAdaptorMonitorSet()`, with this function as an argument, to cause the monitor
12438b724c91SMatthew G. Knepley to be used during the adaptation loop.
12448b724c91SMatthew G. Knepley
12458b724c91SMatthew G. Knepley Call `DMAdaptorMonitorErrorDrawLGCreate()` to create the context needed for this monitor
12468b724c91SMatthew G. Knepley
1247f8d70eaaSPierre Jolivet .seealso: [](ch_snes), `PETSCVIEWERDRAW`, `DMAdaptor`, `DMAdaptorMonitorSet()`, `DMAdaptorMonitorErrorDraw()`, `DMAdaptorMonitorError()`,
12488b724c91SMatthew G. Knepley `DMAdaptorMonitorTrueResidualDrawLGCreate()`
12498b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorErrorDrawLG(DMAdaptor adaptor,PetscInt n,DM odm,DM adm,PetscInt Nf,PetscReal enorms[],Vec error,PetscViewerAndFormat * vf)12508b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorErrorDrawLG(DMAdaptor adaptor, PetscInt n, DM odm, DM adm, PetscInt Nf, PetscReal enorms[], Vec error, PetscViewerAndFormat *vf)
12518b724c91SMatthew G. Knepley {
12528b724c91SMatthew G. Knepley PetscViewer viewer = vf->viewer;
12538b724c91SMatthew G. Knepley PetscViewerFormat format = vf->format;
1254ce78bad3SBarry Smith PetscDrawLG lg;
12558b724c91SMatthew G. Knepley PetscReal *x, *e;
12568b724c91SMatthew G. Knepley
12578b724c91SMatthew G. Knepley PetscFunctionBegin;
12588b724c91SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
1259ce78bad3SBarry Smith PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
12608b724c91SMatthew G. Knepley PetscCall(PetscCalloc2(Nf, &x, Nf, &e));
12618b724c91SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format));
12628b724c91SMatthew G. Knepley if (!n) PetscCall(PetscDrawLGReset(lg));
12638b724c91SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) {
12648b724c91SMatthew G. Knepley x[f] = (PetscReal)n;
12658b724c91SMatthew G. Knepley e[f] = enorms[f] > 0.0 ? PetscLog10Real(enorms[f]) : -15.;
12668b724c91SMatthew G. Knepley }
12678b724c91SMatthew G. Knepley PetscCall(PetscDrawLGAddPoint(lg, x, e));
12688b724c91SMatthew G. Knepley PetscCall(PetscDrawLGDraw(lg));
12698b724c91SMatthew G. Knepley PetscCall(PetscDrawLGSave(lg));
12708b724c91SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer));
12718b724c91SMatthew G. Knepley PetscCall(PetscFree2(x, e));
12728b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
12738b724c91SMatthew G. Knepley }
12748b724c91SMatthew G. Knepley
12758b724c91SMatthew G. Knepley /*@C
12768b724c91SMatthew G. Knepley DMAdaptorMonitorRegisterAll - Registers all of the mesh adaptation monitors in the `SNES` package.
12778b724c91SMatthew G. Knepley
12788b724c91SMatthew G. Knepley Not Collective
12798b724c91SMatthew G. Knepley
12808b724c91SMatthew G. Knepley Level: advanced
12818b724c91SMatthew G. Knepley
12828b724c91SMatthew G. Knepley .seealso: [](ch_snes), `SNES`, `DM`, `DMAdaptorMonitorRegister()`, `DMAdaptorRegister()`
12838b724c91SMatthew G. Knepley @*/
DMAdaptorMonitorRegisterAll(void)12848b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorMonitorRegisterAll(void)
12858b724c91SMatthew G. Knepley {
12868b724c91SMatthew G. Knepley PetscFunctionBegin;
12878b724c91SMatthew G. Knepley if (DMAdaptorMonitorRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
12888b724c91SMatthew G. Knepley DMAdaptorMonitorRegisterAllCalled = PETSC_TRUE;
12898b724c91SMatthew G. Knepley
12908b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitorRegister("size", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorSize, NULL, NULL));
12918b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERASCII, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorError, NULL, NULL));
12928b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DEFAULT, DMAdaptorMonitorErrorDraw, NULL, NULL));
12938b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitorRegister("error", PETSCVIEWERDRAW, PETSC_VIEWER_DRAW_LG, DMAdaptorMonitorErrorDrawLG, DMAdaptorMonitorErrorDrawLGCreate, NULL));
12948b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
12958b724c91SMatthew G. Knepley }
12968b724c91SMatthew G. Knepley
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[])12978b724c91SMatthew 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[])
12988b724c91SMatthew G. Knepley {
12998b724c91SMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
13008b724c91SMatthew G. Knepley
13018b724c91SMatthew G. Knepley for (PetscInt i = 0; i < Nc; ++i) f[i] = u[i];
13028b724c91SMatthew G. Knepley }
13038b724c91SMatthew G. Knepley
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[])1304c6f61ee2SBarry 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[])
1305c6f61ee2SBarry Smith {
13068b724c91SMatthew G. Knepley for (PetscInt i = 0; i < dim; ++i) {
13078b724c91SMatthew G. Knepley for (PetscInt j = 0; j < dim; ++j) f[i + dim * j] = u[i + dim * j];
1308c6f61ee2SBarry Smith }
1309c6f61ee2SBarry Smith }
1310c6f61ee2SBarry Smith
DMAdaptorAdapt_Sequence_Private(DMAdaptor adaptor,Vec inx,PetscBool doSolve,DM * adm,Vec * ax)1311c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorAdapt_Sequence_Private(DMAdaptor adaptor, Vec inx, PetscBool doSolve, DM *adm, Vec *ax)
1312c6f61ee2SBarry Smith {
13133a336bb1SMatthew G. Knepley PetscDS ds;
13148b724c91SMatthew G. Knepley PetscReal errorNorm = 0.;
1315c6f61ee2SBarry Smith PetscInt numAdapt = adaptor->numSeq, adaptIter;
13163a336bb1SMatthew G. Knepley PetscInt dim, coordDim, Nf;
13178b724c91SMatthew G. Knepley void *ctx;
13188b724c91SMatthew G. Knepley MPI_Comm comm;
1319c6f61ee2SBarry Smith
1320c6f61ee2SBarry Smith PetscFunctionBegin;
1321c6f61ee2SBarry Smith PetscCall(DMViewFromOptions(adaptor->idm, NULL, "-dm_adapt_pre_view"));
1322c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(inx, NULL, "-sol_adapt_pre_view"));
1323c6f61ee2SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)adaptor, &comm));
1324c6f61ee2SBarry Smith PetscCall(DMGetDimension(adaptor->idm, &dim));
1325c6f61ee2SBarry Smith PetscCall(DMGetCoordinateDim(adaptor->idm, &coordDim));
1326c6f61ee2SBarry Smith PetscCall(DMGetApplicationContext(adaptor->idm, &ctx));
13273a336bb1SMatthew G. Knepley PetscCall(DMGetDS(adaptor->idm, &ds));
13283a336bb1SMatthew G. Knepley PetscCall(PetscDSGetNumFields(ds, &Nf));
13293a336bb1SMatthew G. Knepley PetscCheck(Nf != 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot refine with no fields present!");
1330c6f61ee2SBarry Smith
1331c6f61ee2SBarry Smith /* Adapt until nothing changes */
1332c6f61ee2SBarry Smith /* Adapt for a specified number of iterates */
1333c6f61ee2SBarry Smith for (adaptIter = 0; adaptIter < numAdapt - 1; ++adaptIter) PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(comm)));
1334c6f61ee2SBarry Smith for (adaptIter = 0; adaptIter < numAdapt; ++adaptIter) {
1335c6f61ee2SBarry Smith PetscBool adapted = PETSC_FALSE;
1336c6f61ee2SBarry Smith DM dm = adaptIter ? *adm : adaptor->idm, odm;
1337c6f61ee2SBarry Smith Vec x = adaptIter ? *ax : inx, locX, ox;
13388b724c91SMatthew G. Knepley Vec error = NULL;
1339c6f61ee2SBarry Smith
1340c6f61ee2SBarry Smith PetscCall(DMGetLocalVector(dm, &locX));
1341c6f61ee2SBarry Smith PetscCall(DMAdaptorPreAdapt(adaptor, locX));
1342c6f61ee2SBarry Smith if (doSolve) {
1343c6f61ee2SBarry Smith SNES snes;
1344c6f61ee2SBarry Smith
1345c6f61ee2SBarry Smith PetscCall(DMAdaptorGetSolver(adaptor, &snes));
1346c6f61ee2SBarry Smith PetscCall(SNESSolve(snes, NULL, adaptIter ? *ax : x));
1347c6f61ee2SBarry Smith }
13483a336bb1SMatthew G. Knepley PetscCall(DMGlobalToLocalBegin(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
13493a336bb1SMatthew G. Knepley PetscCall(DMGlobalToLocalEnd(dm, adaptIter ? *ax : x, INSERT_VALUES, locX));
13503a336bb1SMatthew G. Knepley PetscCall(VecViewFromOptions(adaptIter ? *ax : x, (PetscObject)adaptor, "-adapt_primal_sol_vec_view"));
1351c6f61ee2SBarry Smith switch (adaptor->adaptCriterion) {
1352c6f61ee2SBarry Smith case DM_ADAPTATION_REFINE:
1353c6f61ee2SBarry Smith PetscCall(DMRefine(dm, comm, &odm));
1354c6f61ee2SBarry Smith PetscCheck(odm, comm, PETSC_ERR_ARG_INCOMP, "DMRefine() did not perform any refinement, cannot continue grid sequencing");
1355c6f61ee2SBarry Smith adapted = PETSC_TRUE;
13568b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, NULL));
1357c6f61ee2SBarry Smith break;
1358c6f61ee2SBarry Smith case DM_ADAPTATION_LABEL: {
1359c6f61ee2SBarry Smith /* Adapt DM
1360c6f61ee2SBarry Smith Create local solution
1361c6f61ee2SBarry Smith Reconstruct gradients (FVM) or solve adjoint equation (FEM)
1362c6f61ee2SBarry Smith Produce cellwise error indicator */
13633a336bb1SMatthew G. Knepley DM edm, plex;
13648b724c91SMatthew G. Knepley PetscDS ds;
13653a336bb1SMatthew G. Knepley PetscFE efe;
1366c6f61ee2SBarry Smith DMLabel adaptLabel;
1367c6f61ee2SBarry Smith IS refineIS, coarsenIS;
13683a336bb1SMatthew G. Knepley DMPolytopeType ct;
13698b724c91SMatthew G. Knepley PetscScalar errorVal;
13703a336bb1SMatthew G. Knepley PetscInt nRefine, nCoarsen, cStart;
1371c6f61ee2SBarry Smith
1372c6f61ee2SBarry Smith PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
1373c6f61ee2SBarry Smith
13743a336bb1SMatthew G. Knepley // TODO Move this creation to PreAdapt
13753a336bb1SMatthew G. Knepley PetscCall(DMClone(dm, &edm));
13763a336bb1SMatthew G. Knepley PetscCall(DMConvert(edm, DMPLEX, &plex));
13773a336bb1SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, NULL));
13783a336bb1SMatthew G. Knepley PetscCall(DMPlexGetCellType(plex, cStart, &ct));
13793a336bb1SMatthew G. Knepley PetscCall(DMDestroy(&plex));
13803a336bb1SMatthew G. Knepley PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 0, PETSC_DEFAULT, &efe));
13813a336bb1SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)efe, "Error"));
13823a336bb1SMatthew G. Knepley PetscCall(DMSetField(edm, 0, NULL, (PetscObject)efe));
13833a336bb1SMatthew G. Knepley PetscCall(PetscFEDestroy(&efe));
13843a336bb1SMatthew G. Knepley PetscCall(DMCreateDS(edm));
13858b724c91SMatthew G. Knepley PetscCall(DMGetGlobalVector(edm, &error));
13868b724c91SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)error, "Error Estimator"));
1387c6f61ee2SBarry Smith
13888b724c91SMatthew G. Knepley PetscUseTypeMethod(adaptor, computeerrorindicator, locX, error);
13898b724c91SMatthew G. Knepley PetscCall(VecViewFromOptions(error, (PetscObject)adaptor, "-adapt_error_vec_view"));
13908b724c91SMatthew G. Knepley PetscCall(DMGetDS(edm, &ds));
13918b724c91SMatthew G. Knepley PetscCall(PetscDSSetObjective(ds, 0, identity));
13928b724c91SMatthew G. Knepley PetscCall(DMPlexComputeIntegralFEM(edm, error, &errorVal, NULL));
13938b724c91SMatthew G. Knepley errorNorm = PetscRealPart(errorVal);
13943a336bb1SMatthew G. Knepley
13953a336bb1SMatthew G. Knepley // Compute IS from VecTagger
13968b724c91SMatthew G. Knepley PetscCall(VecTaggerComputeIS(adaptor->refineTag, error, &refineIS, NULL));
13978b724c91SMatthew G. Knepley PetscCall(VecTaggerComputeIS(adaptor->coarsenTag, error, &coarsenIS, NULL));
13983a336bb1SMatthew G. Knepley PetscCall(ISViewFromOptions(refineIS, (PetscObject)adaptor->refineTag, "-is_view"));
13993a336bb1SMatthew G. Knepley PetscCall(ISViewFromOptions(coarsenIS, (PetscObject)adaptor->coarsenTag, "-is_view"));
1400c6f61ee2SBarry Smith PetscCall(ISGetSize(refineIS, &nRefine));
1401c6f61ee2SBarry Smith PetscCall(ISGetSize(coarsenIS, &nCoarsen));
1402c6f61ee2SBarry Smith PetscCall(PetscInfo(adaptor, "DMAdaptor: numRefine %" PetscInt_FMT ", numCoarsen %" PetscInt_FMT "\n", nRefine, nCoarsen));
1403c6f61ee2SBarry Smith if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
1404c6f61ee2SBarry Smith if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS));
1405c6f61ee2SBarry Smith PetscCall(ISDestroy(&coarsenIS));
1406c6f61ee2SBarry Smith PetscCall(ISDestroy(&refineIS));
14073a336bb1SMatthew G. Knepley // Adapt DM from label
1408c6f61ee2SBarry Smith if (nRefine || nCoarsen) {
14093a336bb1SMatthew G. Knepley char oprefix[PETSC_MAX_PATH_LEN];
14103a336bb1SMatthew G. Knepley const char *p;
14113a336bb1SMatthew G. Knepley PetscBool flg;
14123a336bb1SMatthew G. Knepley
14133a336bb1SMatthew G. Knepley PetscCall(PetscOptionsHasName(NULL, adaptor->hdr.prefix, "-adapt_vec_view", &flg));
14143a336bb1SMatthew G. Knepley if (flg) {
14153a336bb1SMatthew G. Knepley Vec ref;
14163a336bb1SMatthew G. Knepley
14173a336bb1SMatthew G. Knepley PetscCall(DMPlexCreateLabelField(dm, adaptLabel, &ref));
14183a336bb1SMatthew G. Knepley PetscCall(VecViewFromOptions(ref, (PetscObject)adaptor, "-adapt_vec_view"));
14193a336bb1SMatthew G. Knepley PetscCall(VecDestroy(&ref));
14203a336bb1SMatthew G. Knepley }
14213a336bb1SMatthew G. Knepley
14223a336bb1SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &p));
14233a336bb1SMatthew G. Knepley PetscCall(PetscStrncpy(oprefix, p, PETSC_MAX_PATH_LEN));
14243a336bb1SMatthew G. Knepley PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)dm, "adapt_"));
1425c6f61ee2SBarry Smith PetscCall(DMAdaptLabel(dm, adaptLabel, &odm));
14263a336bb1SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, oprefix));
14273a336bb1SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)odm, oprefix));
14288b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, error));
1429c6f61ee2SBarry Smith adapted = PETSC_TRUE;
14308b724c91SMatthew G. Knepley } else {
14318b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, dm, 1, &errorNorm, error));
1432c6f61ee2SBarry Smith }
1433c6f61ee2SBarry Smith PetscCall(DMLabelDestroy(&adaptLabel));
14348b724c91SMatthew G. Knepley PetscCall(DMRestoreGlobalVector(edm, &error));
14358b724c91SMatthew G. Knepley PetscCall(DMDestroy(&edm));
1436c6f61ee2SBarry Smith } break;
1437c6f61ee2SBarry Smith case DM_ADAPTATION_METRIC: {
1438c6f61ee2SBarry Smith DM dmGrad, dmHess, dmMetric, dmDet;
1439c6f61ee2SBarry Smith Vec xGrad, xHess, metric, determinant;
1440c6f61ee2SBarry Smith PetscReal N;
1441c6f61ee2SBarry Smith DMLabel bdLabel = NULL, rgLabel = NULL;
1442c6f61ee2SBarry Smith PetscBool higherOrder = PETSC_FALSE;
1443c6f61ee2SBarry Smith PetscInt Nd = coordDim * coordDim, f, vStart, vEnd;
1444c6f61ee2SBarry 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[]);
1445c6f61ee2SBarry Smith
1446c6f61ee2SBarry Smith PetscCall(PetscMalloc(1, &funcs));
1447c6f61ee2SBarry Smith funcs[0] = identityFunc;
1448c6f61ee2SBarry Smith
1449c6f61ee2SBarry Smith /* Setup finite element spaces */
1450c6f61ee2SBarry Smith PetscCall(DMClone(dm, &dmGrad));
1451c6f61ee2SBarry Smith PetscCall(DMClone(dm, &dmHess));
14523a336bb1SMatthew G. Knepley PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple fields not yet considered"); // TODO
14533a336bb1SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
1454c6f61ee2SBarry Smith PetscFE fe, feGrad, feHess;
1455c6f61ee2SBarry Smith PetscDualSpace Q;
1456c6f61ee2SBarry Smith PetscSpace space;
1457c6f61ee2SBarry Smith DM K;
1458c6f61ee2SBarry Smith PetscQuadrature q;
1459c6f61ee2SBarry Smith PetscInt Nc, qorder, p;
1460c6f61ee2SBarry Smith const char *prefix;
1461c6f61ee2SBarry Smith
14623a336bb1SMatthew G. Knepley PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
1463c6f61ee2SBarry Smith PetscCall(PetscFEGetNumComponents(fe, &Nc));
1464c6f61ee2SBarry Smith PetscCheck(Nc <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple components not yet considered"); // TODO
1465c6f61ee2SBarry Smith PetscCall(PetscFEGetBasisSpace(fe, &space));
1466c6f61ee2SBarry Smith PetscCall(PetscSpaceGetDegree(space, NULL, &p));
1467c6f61ee2SBarry Smith if (p > 1) higherOrder = PETSC_TRUE;
1468c6f61ee2SBarry Smith PetscCall(PetscFEGetDualSpace(fe, &Q));
1469c6f61ee2SBarry Smith PetscCall(PetscDualSpaceGetDM(Q, &K));
1470c6f61ee2SBarry Smith PetscCall(DMPlexGetDepthStratum(K, 0, &vStart, &vEnd));
1471c6f61ee2SBarry Smith PetscCall(PetscFEGetQuadrature(fe, &q));
1472c6f61ee2SBarry Smith PetscCall(PetscQuadratureGetOrder(q, &qorder));
1473c6f61ee2SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix));
1474c6f61ee2SBarry Smith PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmGrad), dim, Nc * coordDim, PETSC_TRUE, prefix, qorder, &feGrad));
1475c6f61ee2SBarry Smith PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmHess), dim, Nc * Nd, PETSC_TRUE, prefix, qorder, &feHess));
1476c6f61ee2SBarry Smith PetscCall(DMSetField(dmGrad, f, NULL, (PetscObject)feGrad));
1477c6f61ee2SBarry Smith PetscCall(DMSetField(dmHess, f, NULL, (PetscObject)feHess));
1478c6f61ee2SBarry Smith PetscCall(DMCreateDS(dmGrad));
1479c6f61ee2SBarry Smith PetscCall(DMCreateDS(dmHess));
1480c6f61ee2SBarry Smith PetscCall(PetscFEDestroy(&feGrad));
1481c6f61ee2SBarry Smith PetscCall(PetscFEDestroy(&feHess));
1482c6f61ee2SBarry Smith }
1483c6f61ee2SBarry Smith /* Compute vertexwise gradients from cellwise gradients */
1484c6f61ee2SBarry Smith PetscCall(DMCreateLocalVector(dmGrad, &xGrad));
1485c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(locX, NULL, "-sol_adapt_loc_pre_view"));
1486c6f61ee2SBarry Smith PetscCall(DMPlexComputeGradientClementInterpolant(dm, locX, xGrad));
1487c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(xGrad, NULL, "-adapt_gradient_view"));
1488c6f61ee2SBarry Smith /* Compute vertexwise Hessians from cellwise Hessians */
1489c6f61ee2SBarry Smith PetscCall(DMCreateLocalVector(dmHess, &xHess));
1490c6f61ee2SBarry Smith PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, xGrad, xHess));
1491c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(xHess, NULL, "-adapt_hessian_view"));
1492c6f61ee2SBarry Smith PetscCall(VecDestroy(&xGrad));
1493c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmGrad));
1494c6f61ee2SBarry Smith /* Compute L-p normalized metric */
1495c6f61ee2SBarry Smith PetscCall(DMClone(dm, &dmMetric));
1496c6f61ee2SBarry Smith N = adaptor->Nadapt >= 0 ? adaptor->Nadapt : PetscPowRealInt(adaptor->refinementFactor, dim) * ((PetscReal)(vEnd - vStart));
14978b724c91SMatthew G. Knepley // TODO This was where the old monitor was, figure out how to show metric and target N
1498835f2295SStefano Zampini PetscCall(DMPlexMetricSetTargetComplexity(dmMetric, N));
1499c6f61ee2SBarry Smith if (higherOrder) {
1500c6f61ee2SBarry Smith /* Project Hessian into P1 space, if required */
1501c6f61ee2SBarry Smith PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1502c6f61ee2SBarry Smith PetscCall(DMProjectFieldLocal(dmMetric, 0.0, xHess, funcs, INSERT_ALL_VALUES, metric));
1503c6f61ee2SBarry Smith PetscCall(VecDestroy(&xHess));
1504c6f61ee2SBarry Smith xHess = metric;
1505c6f61ee2SBarry Smith }
1506c6f61ee2SBarry Smith PetscCall(PetscFree(funcs));
1507c6f61ee2SBarry Smith PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric));
1508c6f61ee2SBarry Smith PetscCall(DMPlexMetricDeterminantCreate(dmMetric, 0, &determinant, &dmDet));
1509c6f61ee2SBarry Smith PetscCall(DMPlexMetricNormalize(dmMetric, xHess, PETSC_TRUE, PETSC_TRUE, metric, determinant));
1510c6f61ee2SBarry Smith PetscCall(VecDestroy(&determinant));
1511c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmDet));
1512c6f61ee2SBarry Smith PetscCall(VecDestroy(&xHess));
1513c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmHess));
1514c6f61ee2SBarry Smith /* Adapt DM from metric */
1515c6f61ee2SBarry Smith PetscCall(DMGetLabel(dm, "marker", &bdLabel));
1516c6f61ee2SBarry Smith PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &odm));
15178b724c91SMatthew G. Knepley PetscCall(DMAdaptorMonitor(adaptor, adaptIter, dm, odm, 1, &errorNorm, NULL));
1518c6f61ee2SBarry Smith adapted = PETSC_TRUE;
1519c6f61ee2SBarry Smith /* Cleanup */
1520c6f61ee2SBarry Smith PetscCall(VecDestroy(&metric));
1521c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmMetric));
1522c6f61ee2SBarry Smith } break;
1523c6f61ee2SBarry Smith default:
1524c6f61ee2SBarry Smith SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid adaptation type: %d", adaptor->adaptCriterion);
1525c6f61ee2SBarry Smith }
1526c6f61ee2SBarry Smith PetscCall(DMAdaptorPostAdapt(adaptor));
1527c6f61ee2SBarry Smith PetscCall(DMRestoreLocalVector(dm, &locX));
1528c6f61ee2SBarry Smith /* If DM was adapted, replace objects and recreate solution */
1529c6f61ee2SBarry Smith if (adapted) {
1530c6f61ee2SBarry Smith const char *name;
1531c6f61ee2SBarry Smith
1532c6f61ee2SBarry Smith PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1533c6f61ee2SBarry Smith PetscCall(PetscObjectSetName((PetscObject)odm, name));
1534c6f61ee2SBarry Smith /* Reconfigure solver */
1535c6f61ee2SBarry Smith PetscCall(SNESReset(adaptor->snes));
1536c6f61ee2SBarry Smith PetscCall(SNESSetDM(adaptor->snes, odm));
1537c6f61ee2SBarry Smith PetscCall(DMAdaptorSetSolver(adaptor, adaptor->snes));
1538c6f61ee2SBarry Smith PetscCall(DMPlexSetSNESLocalFEM(odm, PETSC_FALSE, ctx));
1539c6f61ee2SBarry Smith PetscCall(SNESSetFromOptions(adaptor->snes));
1540c6f61ee2SBarry Smith /* Transfer system */
1541c6f61ee2SBarry Smith PetscCall(DMCopyDisc(dm, odm));
1542c6f61ee2SBarry Smith /* Transfer solution */
1543c6f61ee2SBarry Smith PetscCall(DMCreateGlobalVector(odm, &ox));
1544c6f61ee2SBarry Smith PetscCall(PetscObjectGetName((PetscObject)x, &name));
1545c6f61ee2SBarry Smith PetscCall(PetscObjectSetName((PetscObject)ox, name));
1546c6f61ee2SBarry Smith PetscCall(DMAdaptorTransferSolution(adaptor, dm, x, odm, ox));
1547c6f61ee2SBarry Smith /* Cleanup adaptivity info */
1548c6f61ee2SBarry Smith if (adaptIter > 0) PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(comm)));
1549c6f61ee2SBarry Smith PetscCall(DMForestSetAdaptivityForest(dm, NULL)); /* clear internal references to the previous dm */
1550c6f61ee2SBarry Smith PetscCall(DMDestroy(&dm));
1551c6f61ee2SBarry Smith PetscCall(VecDestroy(&x));
1552c6f61ee2SBarry Smith *adm = odm;
1553c6f61ee2SBarry Smith *ax = ox;
1554c6f61ee2SBarry Smith } else {
1555c6f61ee2SBarry Smith *adm = dm;
1556c6f61ee2SBarry Smith *ax = x;
1557c6f61ee2SBarry Smith adaptIter = numAdapt;
1558c6f61ee2SBarry Smith }
1559c6f61ee2SBarry Smith if (adaptIter < numAdapt - 1) {
1560c6f61ee2SBarry Smith PetscCall(DMViewFromOptions(odm, NULL, "-dm_adapt_iter_view"));
1561c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(ox, NULL, "-sol_adapt_iter_view"));
1562c6f61ee2SBarry Smith }
1563c6f61ee2SBarry Smith }
1564c6f61ee2SBarry Smith PetscCall(DMViewFromOptions(*adm, NULL, "-dm_adapt_view"));
1565c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(*ax, NULL, "-sol_adapt_view"));
1566c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
1567c6f61ee2SBarry Smith }
1568c6f61ee2SBarry Smith
1569c6f61ee2SBarry Smith /*@
1570c6f61ee2SBarry Smith DMAdaptorAdapt - Creates a new `DM` that is adapted to the problem
1571c6f61ee2SBarry Smith
1572c6f61ee2SBarry Smith Not Collective
1573c6f61ee2SBarry Smith
1574c6f61ee2SBarry Smith Input Parameters:
1575c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object
1576c6f61ee2SBarry Smith . x - The global approximate solution
1577c6f61ee2SBarry Smith - strategy - The adaptation strategy, see `DMAdaptationStrategy`
1578c6f61ee2SBarry Smith
1579c6f61ee2SBarry Smith Output Parameters:
1580c6f61ee2SBarry Smith + adm - The adapted `DM`
1581c6f61ee2SBarry Smith - ax - The adapted solution
1582c6f61ee2SBarry Smith
1583c6f61ee2SBarry Smith Options Database Keys:
1584c6f61ee2SBarry Smith + -snes_adapt <strategy> - initial, sequential, multigrid
1585c6f61ee2SBarry Smith . -adapt_gradient_view - View the Clement interpolant of the solution gradient
1586c6f61ee2SBarry Smith . -adapt_hessian_view - View the Clement interpolant of the solution Hessian
1587c6f61ee2SBarry Smith - -adapt_metric_view - View the metric tensor for adaptive mesh refinement
1588c6f61ee2SBarry Smith
1589c6f61ee2SBarry Smith Level: intermediate
1590c6f61ee2SBarry Smith
1591df514481SBarry Smith .seealso: [](ch_dmbase), `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`
1592c6f61ee2SBarry Smith @*/
DMAdaptorAdapt(DMAdaptor adaptor,Vec x,DMAdaptationStrategy strategy,DM * adm,Vec * ax)1593c6f61ee2SBarry Smith PetscErrorCode DMAdaptorAdapt(DMAdaptor adaptor, Vec x, DMAdaptationStrategy strategy, DM *adm, Vec *ax)
1594c6f61ee2SBarry Smith {
1595c6f61ee2SBarry Smith PetscFunctionBegin;
1596c6f61ee2SBarry Smith switch (strategy) {
1597c6f61ee2SBarry Smith case DM_ADAPTATION_INITIAL:
1598c6f61ee2SBarry Smith PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_FALSE, adm, ax));
1599c6f61ee2SBarry Smith break;
1600c6f61ee2SBarry Smith case DM_ADAPTATION_SEQUENTIAL:
1601c6f61ee2SBarry Smith PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_TRUE, adm, ax));
1602c6f61ee2SBarry Smith break;
1603c6f61ee2SBarry Smith default:
1604c6f61ee2SBarry Smith SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Unrecognized adaptation strategy %d", strategy);
1605c6f61ee2SBarry Smith }
1606c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
1607c6f61ee2SBarry Smith }
16083a336bb1SMatthew G. Knepley
16093a336bb1SMatthew G. Knepley /*@C
16103a336bb1SMatthew G. Knepley DMAdaptorGetMixedSetupFunction - Get the function setting up the mixed problem, if it exists
16113a336bb1SMatthew G. Knepley
16123a336bb1SMatthew G. Knepley Not Collective
16133a336bb1SMatthew G. Knepley
16143a336bb1SMatthew G. Knepley Input Parameter:
16153a336bb1SMatthew G. Knepley . adaptor - the `DMAdaptor`
16163a336bb1SMatthew G. Knepley
16173a336bb1SMatthew G. Knepley Output Parameter:
16183a336bb1SMatthew G. Knepley . setupFunc - the function setting up the mixed problem, or `NULL`
16193a336bb1SMatthew G. Knepley
16203a336bb1SMatthew G. Knepley Level: advanced
16213a336bb1SMatthew G. Knepley
16223a336bb1SMatthew G. Knepley .seealso: `DMAdaptor`, `DMAdaptorSetMixedSetupFunction()`, `DMAdaptorAdapt()`
16233a336bb1SMatthew G. Knepley @*/
DMAdaptorGetMixedSetupFunction(DMAdaptor adaptor,PetscErrorCode (** setupFunc)(DMAdaptor,DM))16243a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorGetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (**setupFunc)(DMAdaptor, DM))
16253a336bb1SMatthew G. Knepley {
16263a336bb1SMatthew G. Knepley PetscFunctionBegin;
16273a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
16283a336bb1SMatthew G. Knepley PetscAssertPointer(setupFunc, 2);
16293a336bb1SMatthew G. Knepley *setupFunc = adaptor->ops->mixedsetup;
16303a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
16313a336bb1SMatthew G. Knepley }
16323a336bb1SMatthew G. Knepley
16333a336bb1SMatthew G. Knepley /*@C
16343a336bb1SMatthew G. Knepley DMAdaptorSetMixedSetupFunction - Set the function setting up the mixed problem
16353a336bb1SMatthew G. Knepley
16363a336bb1SMatthew G. Knepley Not Collective
16373a336bb1SMatthew G. Knepley
16383a336bb1SMatthew G. Knepley Input Parameters:
16393a336bb1SMatthew G. Knepley + adaptor - the `DMAdaptor`
16403a336bb1SMatthew G. Knepley - setupFunc - the function setting up the mixed problem
16413a336bb1SMatthew G. Knepley
16423a336bb1SMatthew G. Knepley Level: advanced
16433a336bb1SMatthew G. Knepley
16443a336bb1SMatthew G. Knepley .seealso: `DMAdaptor`, `DMAdaptorGetMixedSetupFunction()`, `DMAdaptorAdapt()`
16453a336bb1SMatthew G. Knepley @*/
DMAdaptorSetMixedSetupFunction(DMAdaptor adaptor,PetscErrorCode (* setupFunc)(DMAdaptor,DM))16463a336bb1SMatthew G. Knepley PetscErrorCode DMAdaptorSetMixedSetupFunction(DMAdaptor adaptor, PetscErrorCode (*setupFunc)(DMAdaptor, DM))
16473a336bb1SMatthew G. Knepley {
16483a336bb1SMatthew G. Knepley PetscFunctionBegin;
16493a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
16503a336bb1SMatthew G. Knepley PetscValidFunction(setupFunc, 2);
16513a336bb1SMatthew G. Knepley adaptor->ops->mixedsetup = setupFunc;
16523a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
16533a336bb1SMatthew G. Knepley }
16543a336bb1SMatthew G. Knepley
1655ce78bad3SBarry Smith /*@
16568b724c91SMatthew G. Knepley DMAdaptorGetCriterion - Get the adaptation criterion
16578b724c91SMatthew G. Knepley
16588b724c91SMatthew G. Knepley Not Collective
16598b724c91SMatthew G. Knepley
16608b724c91SMatthew G. Knepley Input Parameter:
16618b724c91SMatthew G. Knepley . adaptor - the `DMAdaptor`
16628b724c91SMatthew G. Knepley
16638b724c91SMatthew G. Knepley Output Parameter:
16648b724c91SMatthew G. Knepley . criterion - the criterion for adaptation
16658b724c91SMatthew G. Knepley
16668b724c91SMatthew G. Knepley Level: advanced
16678b724c91SMatthew G. Knepley
1668f8d70eaaSPierre Jolivet .seealso: `DMAdaptor`, `DMAdaptorSetCriterion()`, `DMAdaptationCriterion`
16698b724c91SMatthew G. Knepley @*/
DMAdaptorGetCriterion(DMAdaptor adaptor,DMAdaptationCriterion * criterion)16708b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorGetCriterion(DMAdaptor adaptor, DMAdaptationCriterion *criterion)
16718b724c91SMatthew G. Knepley {
16728b724c91SMatthew G. Knepley PetscFunctionBegin;
16738b724c91SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
16748b724c91SMatthew G. Knepley PetscAssertPointer(criterion, 2);
16758b724c91SMatthew G. Knepley *criterion = adaptor->adaptCriterion;
16768b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
16778b724c91SMatthew G. Knepley }
16788b724c91SMatthew G. Knepley
1679ce78bad3SBarry Smith /*@
16808b724c91SMatthew G. Knepley DMAdaptorSetCriterion - Set the adaptation criterion
16818b724c91SMatthew G. Knepley
16828b724c91SMatthew G. Knepley Not Collective
16838b724c91SMatthew G. Knepley
16848b724c91SMatthew G. Knepley Input Parameters:
16858b724c91SMatthew G. Knepley + adaptor - the `DMAdaptor`
16868b724c91SMatthew G. Knepley - criterion - the adaptation criterion
16878b724c91SMatthew G. Knepley
16888b724c91SMatthew G. Knepley Level: advanced
16898b724c91SMatthew G. Knepley
16908b724c91SMatthew G. Knepley .seealso: `DMAdaptor`, `DMAdaptorGetCriterion()`, `DMAdaptationCriterion`
16918b724c91SMatthew G. Knepley @*/
DMAdaptorSetCriterion(DMAdaptor adaptor,DMAdaptationCriterion criterion)16928b724c91SMatthew G. Knepley PetscErrorCode DMAdaptorSetCriterion(DMAdaptor adaptor, DMAdaptationCriterion criterion)
16938b724c91SMatthew G. Knepley {
16948b724c91SMatthew G. Knepley PetscFunctionBegin;
16958b724c91SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
16968b724c91SMatthew G. Knepley adaptor->adaptCriterion = criterion;
16978b724c91SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
16988b724c91SMatthew G. Knepley }
16998b724c91SMatthew G. Knepley
DMAdaptorInitialize_Gradient(DMAdaptor adaptor)17003a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorInitialize_Gradient(DMAdaptor adaptor)
17013a336bb1SMatthew G. Knepley {
17023a336bb1SMatthew G. Knepley PetscFunctionBegin;
17033a336bb1SMatthew G. Knepley adaptor->ops->computeerrorindicator = DMAdaptorComputeErrorIndicator_Gradient;
17043a336bb1SMatthew G. Knepley adaptor->ops->computecellerrorindicator = DMAdaptorComputeCellErrorIndicator_Gradient;
17053a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
17063a336bb1SMatthew G. Knepley }
17073a336bb1SMatthew G. Knepley
DMAdaptorCreate_Gradient(DMAdaptor adaptor)17083a336bb1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Gradient(DMAdaptor adaptor)
17093a336bb1SMatthew G. Knepley {
17103a336bb1SMatthew G. Knepley PetscFunctionBegin;
17113a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
17123a336bb1SMatthew G. Knepley adaptor->data = NULL;
17133a336bb1SMatthew G. Knepley
17143a336bb1SMatthew G. Knepley PetscCall(DMAdaptorInitialize_Gradient(adaptor));
17153a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
17163a336bb1SMatthew G. Knepley }
17173a336bb1SMatthew G. Knepley
DMAdaptorInitialize_Flux(DMAdaptor adaptor)17183a336bb1SMatthew G. Knepley static PetscErrorCode DMAdaptorInitialize_Flux(DMAdaptor adaptor)
17193a336bb1SMatthew G. Knepley {
17203a336bb1SMatthew G. Knepley PetscFunctionBegin;
17213a336bb1SMatthew G. Knepley adaptor->ops->computeerrorindicator = DMAdaptorComputeErrorIndicator_Flux;
17223a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
17233a336bb1SMatthew G. Knepley }
17243a336bb1SMatthew G. Knepley
DMAdaptorCreate_Flux(DMAdaptor adaptor)17253a336bb1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptorCreate_Flux(DMAdaptor adaptor)
17263a336bb1SMatthew G. Knepley {
17273a336bb1SMatthew G. Knepley PetscFunctionBegin;
17283a336bb1SMatthew G. Knepley PetscValidHeaderSpecific(adaptor, DMADAPTOR_CLASSID, 1);
17293a336bb1SMatthew G. Knepley adaptor->data = NULL;
17303a336bb1SMatthew G. Knepley
17313a336bb1SMatthew G. Knepley PetscCall(DMAdaptorInitialize_Flux(adaptor));
17323a336bb1SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
17333a336bb1SMatthew G. Knepley }
1734