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