xref: /petsc/include/petsc/private/dmadaptorimpl.h (revision ad781fe3223b515d45fb60062f8062326362e786)
1 #pragma once
2 
3 #include <petscdmadaptor.h>
4 #include <petsc/private/petscimpl.h>
5 
6 #define MAXDMADAPTORMONITORS 16
7 
8 typedef struct _DMAdaptorOps *DMAdaptorOps;
9 struct _DMAdaptorOps {
10   PetscErrorCode (*setfromoptions)(DMAdaptor);
11   PetscErrorCode (*setup)(DMAdaptor);
12   PetscErrorCode (*view)(DMAdaptor, PetscViewer);
13   PetscErrorCode (*destroy)(DMAdaptor);
14   PetscErrorCode (*transfersolution)(DMAdaptor, DM, Vec, DM, Vec, void *);
15   PetscErrorCode (*mixedsetup)(DMAdaptor, DM);
16   PetscErrorCode (*computeerrorindicator)(DMAdaptor, Vec, Vec);
17   PetscErrorCode (*computecellerrorindicator)(DMAdaptor, PetscInt, PetscInt, const PetscScalar *, const PetscScalar *, const PetscFVCellGeom *, PetscReal *, void *);
18 };
19 
20 struct _p_DMAdaptor {
21   PETSCHEADER(struct _DMAdaptorOps);
22   void *data;
23 
24   /* Inputs */
25   DM        idm;                   /* Initial grid */
26   SNES      snes;                  /* Solver */
27   VecTagger refineTag, coarsenTag; /* Criteria for adaptivity */
28   /*   control */
29   DMAdaptationCriterion adaptCriterion;
30   PetscBool             femType;
31   PetscInt              numSeq;           /* Number of sequential adaptations */
32   PetscInt              Nadapt;           /* Target number of vertices */
33   PetscReal             refinementFactor; /* N_adapt = r^dim N_orig */
34   /*   FVM support */
35   PetscBool          computeGradient;
36   DM                 cellDM, gradDM;
37   Vec                cellGeom, faceGeom, cellGrad; /* Local vectors */
38   const PetscScalar *cellGeomArray, *cellGradArray;
39   // Monitors
40   PetscErrorCode (*monitor[MAXDMADAPTORMONITORS])(DMAdaptor, PetscInt, DM, DM, PetscInt, PetscReal[], Vec, void *);
41   PetscCtxDestroyFn *monitordestroy[MAXDMADAPTORMONITORS];
42   void              *monitorcontext[MAXDMADAPTORMONITORS];
43   PetscInt           numbermonitors;
44   /* Auxiliary objects */
45   PetscLimiter limiter;
46   PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
47   void **exactCtx;
48 };
49