xref: /petsc/src/snes/tutorials/ex13.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Poisson Problem in 2d and 3d with finite elements.\n\
2c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4c4762a1bSJed Brown This example supports automatic convergence estimation\n\
5c4762a1bSJed Brown and eventually adaptivity.\n\n\n";
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscdmplex.h>
8157a3229SMatthew G. Knepley #include <petscdmceed.h>
9c4762a1bSJed Brown #include <petscsnes.h>
10c4762a1bSJed Brown #include <petscds.h>
11c4762a1bSJed Brown #include <petscconvest.h>
12c4762a1bSJed Brown 
13c4762a1bSJed Brown typedef struct {
14c4762a1bSJed Brown   /* Domain and mesh definition */
15c4762a1bSJed Brown   PetscBool spectral;    /* Look at the spectrum along planes in the solution */
16c4762a1bSJed Brown   PetscBool shear;       /* Shear the domain */
17c4762a1bSJed Brown   PetscBool adjoint;     /* Solve the adjoint problem */
18da81f932SPierre Jolivet   PetscBool homogeneous; /* Use homogeneous boundary conditions */
19ad1a7c93SMatthew Knepley   PetscBool viewError;   /* Output the solution error */
20c4762a1bSJed Brown } AppCtx;
21c4762a1bSJed Brown 
zero(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)22*2a8381b2SBarry Smith static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
23d71ae5a4SJacob Faibussowitsch {
24c4762a1bSJed Brown   *u = 0.0;
253ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
26c4762a1bSJed Brown }
27c4762a1bSJed Brown 
trig_inhomogeneous_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)28*2a8381b2SBarry Smith static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
29d71ae5a4SJacob Faibussowitsch {
30c4762a1bSJed Brown   PetscInt d;
31c4762a1bSJed Brown   *u = 0.0;
32c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
333ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
trig_homogeneous_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)36*2a8381b2SBarry Smith static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
37d71ae5a4SJacob Faibussowitsch {
38b724ec41SToby Isaac   PetscInt d;
39b724ec41SToby Isaac   *u = 1.0;
40b724ec41SToby Isaac   for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0 * PETSC_PI * x[d]);
413ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
42b724ec41SToby Isaac }
43b724ec41SToby Isaac 
44c4762a1bSJed Brown /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */
obj_error_u(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 obj[])45d71ae5a4SJacob Faibussowitsch static void obj_error_u(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 obj[])
46d71ae5a4SJacob Faibussowitsch {
47c4762a1bSJed Brown   obj[0] = a[aOff[0]] * (u[0] - a[aOff[1]]);
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
f0_trig_inhomogeneous_u(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 f0[])50d71ae5a4SJacob Faibussowitsch static void f0_trig_inhomogeneous_u(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 f0[])
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown   PetscInt d;
53c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
f0_trig_homogeneous_u(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 f0[])56d71ae5a4SJacob Faibussowitsch static void f0_trig_homogeneous_u(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 f0[])
57d71ae5a4SJacob Faibussowitsch {
58b724ec41SToby Isaac   PetscInt d;
59b724ec41SToby Isaac   for (d = 0; d < dim; ++d) {
60b724ec41SToby Isaac     PetscScalar v = 1.;
61b724ec41SToby Isaac     for (PetscInt e = 0; e < dim; e++) {
62b724ec41SToby Isaac       if (e == d) {
63b724ec41SToby Isaac         v *= -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
64b724ec41SToby Isaac       } else {
65b724ec41SToby Isaac         v *= PetscSinReal(2.0 * PETSC_PI * x[d]);
66b724ec41SToby Isaac       }
67b724ec41SToby Isaac     }
68b724ec41SToby Isaac     f0[0] += v;
69b724ec41SToby Isaac   }
70b724ec41SToby Isaac }
71b724ec41SToby Isaac 
f0_unity_u(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 f0[])72d71ae5a4SJacob Faibussowitsch static void f0_unity_u(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 f0[])
73d71ae5a4SJacob Faibussowitsch {
74c4762a1bSJed Brown   f0[0] = 1.0;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
f0_identityaux_u(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 f0[])77d71ae5a4SJacob Faibussowitsch static void f0_identityaux_u(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 f0[])
78d71ae5a4SJacob Faibussowitsch {
79c4762a1bSJed Brown   f0[0] = a[0];
80c4762a1bSJed Brown }
81c4762a1bSJed Brown 
f1_u(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 f1[])82d71ae5a4SJacob Faibussowitsch static void f1_u(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 f1[])
83d71ae5a4SJacob Faibussowitsch {
84c4762a1bSJed Brown   PetscInt d;
85c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
g3_uu(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,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g3[])88d71ae5a4SJacob Faibussowitsch static void g3_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
89d71ae5a4SJacob Faibussowitsch {
90c4762a1bSJed Brown   PetscInt d;
91c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
PLEXFE_QFUNCTION(Laplace,f0_trig_inhomogeneous_u,f1_u)94157a3229SMatthew G. Knepley PLEXFE_QFUNCTION(Laplace, f0_trig_inhomogeneous_u, f1_u)
95157a3229SMatthew G. Knepley 
96d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
97d71ae5a4SJacob Faibussowitsch {
98c4762a1bSJed Brown   PetscFunctionBeginUser;
99c4762a1bSJed Brown   options->shear       = PETSC_FALSE;
100c4762a1bSJed Brown   options->spectral    = PETSC_FALSE;
101c4762a1bSJed Brown   options->adjoint     = PETSC_FALSE;
102b724ec41SToby Isaac   options->homogeneous = PETSC_FALSE;
103ad1a7c93SMatthew Knepley   options->viewError   = PETSC_FALSE;
104c4762a1bSJed Brown 
105d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL));
1089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL));
1099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL));
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL));
111d0609cedSBarry Smith   PetscOptionsEnd();
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113c4762a1bSJed Brown }
114c4762a1bSJed Brown 
CreateSpectralPlanes(DM dm,PetscInt numPlanes,const PetscInt planeDir[],const PetscReal planeCoord[],AppCtx * user)115d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
116d71ae5a4SJacob Faibussowitsch {
117c4762a1bSJed Brown   PetscSection       coordSection;
118c4762a1bSJed Brown   Vec                coordinates;
119c4762a1bSJed Brown   const PetscScalar *coords;
120c4762a1bSJed Brown   PetscInt           dim, p, vStart, vEnd, v;
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   PetscFunctionBeginUser;
1239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dim));
1249566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1259566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1269566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
1279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
128c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
129c4762a1bSJed Brown     DMLabel label;
130c4762a1bSJed Brown     char    name[PETSC_MAX_PATH_LEN];
131c4762a1bSJed Brown 
13263a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p));
1339566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, name));
1349566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, name, &label));
1359566063dSJacob Faibussowitsch     PetscCall(DMLabelAddStratum(label, 1));
136c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
137c4762a1bSJed Brown       PetscInt off;
138c4762a1bSJed Brown 
1399566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
14048a46eb9SPierre Jolivet       if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off + planeDir[p]])) < PETSC_SMALL) PetscCall(DMLabelSetValue(label, v, 1));
141c4762a1bSJed Brown     }
142c4762a1bSJed Brown   }
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
1443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)147d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
148d71ae5a4SJacob Faibussowitsch {
149c4762a1bSJed Brown   PetscFunctionBeginUser;
1509566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1519566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1529566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1539566063dSJacob Faibussowitsch   if (user->shear) PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL));
1549566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
1559566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
156c4762a1bSJed Brown   if (user->spectral) {
157c4762a1bSJed Brown     PetscInt  planeDir[2]   = {0, 1};
158c4762a1bSJed Brown     PetscReal planeCoord[2] = {0., 1.};
159c4762a1bSJed Brown 
1609566063dSJacob Faibussowitsch     PetscCall(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user));
161c4762a1bSJed Brown   }
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
SetupPrimalProblem(DM dm,AppCtx * user)165d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
166d71ae5a4SJacob Faibussowitsch {
16745480ffeSMatthew G. Knepley   PetscDS        ds;
16845480ffeSMatthew G. Knepley   DMLabel        label;
169c4762a1bSJed Brown   const PetscInt id                                                                             = 1;
1702192575eSBarry Smith   PetscPointFn  *f0                                                                             = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u;
171b724ec41SToby Isaac   PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u;
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   PetscFunctionBeginUser;
1749566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1759566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, f0, f1_u));
1769566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
1779566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, ex, user));
1789566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
17957d50842SBarry Smith   if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)ex, NULL, user, NULL));
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
SetupAdjointProblem(DM dm,AppCtx * user)183d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user)
184d71ae5a4SJacob Faibussowitsch {
18545480ffeSMatthew G. Knepley   PetscDS        ds;
18645480ffeSMatthew G. Knepley   DMLabel        label;
187c4762a1bSJed Brown   const PetscInt id = 1;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   PetscFunctionBeginUser;
1909566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1919566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u));
1929566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
1939566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(ds, 0, obj_error_u));
1949566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
19557d50842SBarry Smith   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)zero, NULL, user, NULL));
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
SetupErrorProblem(DM dm,AppCtx * user)199d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user)
200d71ae5a4SJacob Faibussowitsch {
201c4762a1bSJed Brown   PetscDS prob;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   PetscFunctionBeginUser;
2049566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
SetupDiscretization(DM dm,const char name[],PetscErrorCode (* setup)(DM,AppCtx *),AppCtx * user)208d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
209d71ae5a4SJacob Faibussowitsch {
210c4762a1bSJed Brown   DM             cdm = dm;
211c4762a1bSJed Brown   PetscFE        fe;
2125be41b18SMatthew G. Knepley   DMPolytopeType ct;
2135be41b18SMatthew G. Knepley   PetscBool      simplex;
2145be41b18SMatthew G. Knepley   PetscInt       dim, cStart;
215c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   PetscFunctionBeginUser;
2189566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
2209566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
2215be41b18SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
222c4762a1bSJed Brown   /* Create finite element */
2239566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
2249566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
226c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
2279566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
2289566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
2299566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
230c4762a1bSJed Brown   while (cdm) {
2319566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
232c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
2339566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
234c4762a1bSJed Brown   }
2359566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
236157a3229SMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
237157a3229SMatthew G. Knepley   PetscBool useCeed;
238157a3229SMatthew G. Knepley   PetscCall(DMPlexGetUseCeed(dm, &useCeed));
239157a3229SMatthew G. Knepley   if (useCeed) PetscCall(DMCeedCreate(dm, PETSC_TRUE, PlexQFunctionLaplace, PlexQFunctionLaplace_loc));
240157a3229SMatthew G. Knepley #endif
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
ComputeSpectral(Vec u,PetscInt numPlanes,const PetscInt planeDir[],const PetscReal planeCoord[],AppCtx * user)244e8e188d2SZach Atkins static PetscErrorCode ComputeSpectral(Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
245d71ae5a4SJacob Faibussowitsch {
246c4762a1bSJed Brown   MPI_Comm           comm;
247e8e188d2SZach Atkins   DM                 dm;
248c4762a1bSJed Brown   PetscSection       coordSection, section;
249c4762a1bSJed Brown   Vec                coordinates, uloc;
250c4762a1bSJed Brown   const PetscScalar *coords, *array;
251c4762a1bSJed Brown   PetscInt           p;
252c4762a1bSJed Brown   PetscMPIInt        size, rank;
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   PetscFunctionBeginUser;
255e8e188d2SZach Atkins   if (!user->spectral) PetscFunctionReturn(PETSC_SUCCESS);
256e8e188d2SZach Atkins   PetscCall(VecGetDM(u, &dm));
2579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2609566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &uloc));
2619566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc));
2629566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc));
2639566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL));
2649566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(uloc, NULL, "-sol_view"));
2659566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
2669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(uloc, &array));
2679566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2689566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
2699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
270c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
271c4762a1bSJed Brown     DMLabel         label;
272c4762a1bSJed Brown     char            name[PETSC_MAX_PATH_LEN];
273c4762a1bSJed Brown     Mat             F;
274c4762a1bSJed Brown     Vec             x, y;
275c4762a1bSJed Brown     IS              stratum;
276c4762a1bSJed Brown     PetscReal      *ray, *gray;
277c4762a1bSJed Brown     PetscScalar    *rvals, *svals, *gsvals;
278c4762a1bSJed Brown     PetscInt       *perm, *nperm;
279c4762a1bSJed Brown     PetscInt        n, N, i, j, off, offu;
2806497c311SBarry Smith     PetscMPIInt     in;
281c4762a1bSJed Brown     const PetscInt *points;
282c4762a1bSJed Brown 
28363a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p));
2849566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, name, &label));
2859566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, 1, &stratum));
2869566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(stratum, &n));
2876497c311SBarry Smith     PetscCall(PetscMPIIntCast(n, &in));
2889566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(stratum, &points));
2899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n, &ray, n, &svals));
290c4762a1bSJed Brown     for (i = 0; i < n; ++i) {
2919566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, points[i], &off));
2929566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, points[i], &offu));
293c4762a1bSJed Brown       ray[i]   = PetscRealPart(coords[off + ((planeDir[p] + 1) % 2)]);
294c4762a1bSJed Brown       svals[i] = array[offu];
295c4762a1bSJed Brown     }
296c4762a1bSJed Brown     /* Gather the ray data to proc 0 */
297c4762a1bSJed Brown     if (size > 1) {
298c4762a1bSJed Brown       PetscMPIInt *cnt, *displs, p;
299c4762a1bSJed Brown 
3009566063dSJacob Faibussowitsch       PetscCall(PetscCalloc2(size, &cnt, size, &displs));
3019566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm));
302c4762a1bSJed Brown       for (p = 1; p < size; ++p) displs[p] = displs[p - 1] + cnt[p - 1];
303c4762a1bSJed Brown       N = displs[size - 1] + cnt[size - 1];
3049566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(N, &gray, N, &gsvals));
3056497c311SBarry Smith       PetscCallMPI(MPI_Gatherv(ray, in, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm));
3066497c311SBarry Smith       PetscCallMPI(MPI_Gatherv(svals, in, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm));
3079566063dSJacob Faibussowitsch       PetscCall(PetscFree2(cnt, displs));
308c4762a1bSJed Brown     } else {
309c4762a1bSJed Brown       N      = n;
310c4762a1bSJed Brown       gray   = ray;
311c4762a1bSJed Brown       gsvals = svals;
312c4762a1bSJed Brown     }
313dd400576SPatrick Sanan     if (rank == 0) {
314c4762a1bSJed Brown       /* Sort point along ray */
3159566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(N, &perm, N, &nperm));
316ad540459SPierre Jolivet       for (i = 0; i < N; ++i) perm[i] = i;
3179566063dSJacob Faibussowitsch       PetscCall(PetscSortRealWithPermutation(N, gray, perm));
318c4762a1bSJed Brown       /* Count duplicates and squish mapping */
319c4762a1bSJed Brown       nperm[0] = perm[0];
320c4762a1bSJed Brown       for (i = 1, j = 1; i < N; ++i) {
321c4762a1bSJed Brown         if (PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) > PETSC_SMALL) nperm[j++] = perm[i];
322c4762a1bSJed Brown       }
323c4762a1bSJed Brown       /* Create FFT structs */
3249566063dSJacob Faibussowitsch       PetscCall(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F));
3259566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(F, &x, &y));
3269566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)y, name));
3279566063dSJacob Faibussowitsch       PetscCall(VecGetArray(x, &rvals));
328c4762a1bSJed Brown       for (i = 0, j = 0; i < N; ++i) {
329c4762a1bSJed Brown         if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) < PETSC_SMALL) continue;
330c4762a1bSJed Brown         rvals[j] = gsvals[nperm[j]];
331c4762a1bSJed Brown         ++j;
332c4762a1bSJed Brown       }
3339566063dSJacob Faibussowitsch       PetscCall(PetscFree2(perm, nperm));
3349566063dSJacob Faibussowitsch       if (size > 1) PetscCall(PetscFree2(gray, gsvals));
3359566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(x, &rvals));
336c4762a1bSJed Brown       /* Do FFT along the ray */
3379566063dSJacob Faibussowitsch       PetscCall(MatMult(F, x, y));
338c4762a1bSJed Brown       /* Chop FFT */
33993ec0da9SPierre Jolivet       PetscCall(VecFilter(y, PETSC_SMALL));
3409566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(x, NULL, "-real_view"));
3419566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(y, NULL, "-fft_view"));
3429566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&x));
3439566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&y));
3449566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&F));
345c4762a1bSJed Brown     }
3469566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(stratum, &points));
3479566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&stratum));
3489566063dSJacob Faibussowitsch     PetscCall(PetscFree2(ray, svals));
349c4762a1bSJed Brown   }
3509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
3519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(uloc, &array));
3529566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &uloc));
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354c4762a1bSJed Brown }
355c4762a1bSJed Brown 
ComputeAdjoint(Vec u,AppCtx * user)356e8e188d2SZach Atkins static PetscErrorCode ComputeAdjoint(Vec u, AppCtx *user)
357d71ae5a4SJacob Faibussowitsch {
358e8e188d2SZach Atkins   PetscFunctionBegin;
359e8e188d2SZach Atkins   if (!user->adjoint) PetscFunctionReturn(PETSC_SUCCESS);
360e8e188d2SZach Atkins   DM   dm, dmAdj;
361c4762a1bSJed Brown   SNES snesAdj;
362c4762a1bSJed Brown   Vec  uAdj;
363c4762a1bSJed Brown 
364e8e188d2SZach Atkins   PetscCall(VecGetDM(u, &dm));
3659566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snesAdj));
3669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)snesAdj, "adjoint_"));
3679566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmAdj));
3689566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snesAdj, dmAdj));
369e8e188d2SZach Atkins   PetscCall(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, user));
3709566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dmAdj, &uAdj));
3719566063dSJacob Faibussowitsch   PetscCall(VecSet(uAdj, 0.0));
3729566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)uAdj, "adjoint"));
3736493148fSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(dmAdj, PETSC_FALSE, &user));
3749566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snesAdj));
3759566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snesAdj, NULL, uAdj));
3769566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snesAdj, &uAdj));
3779566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(uAdj, NULL, "-adjoint_view"));
378c4762a1bSJed Brown   /* Error representation */
379c4762a1bSJed Brown   {
380c4762a1bSJed Brown     DM        dmErr, dmErrAux, dms[2];
381c4762a1bSJed Brown     Vec       errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj;
382c4762a1bSJed Brown     IS       *subis;
383c4762a1bSJed Brown     PetscReal errorEstTot, errorL2Norm, errorL2Tot;
384c4762a1bSJed Brown     PetscInt  N, i;
385e8e188d2SZach Atkins     PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u};
3869371c9d4SSatish Balay     void (*identity[1])(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[]) = {f0_identityaux_u};
387*2a8381b2SBarry Smith     PetscCtx ctxs[1] = {0};
388c4762a1bSJed Brown 
389e8e188d2SZach Atkins     ctxs[0] = user;
3909566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmErr));
391e8e188d2SZach Atkins     PetscCall(SetupDiscretization(dmErr, "error", SetupErrorProblem, user));
3929566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmErr, &errorEst));
3939566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmErr, &errorL2));
394c4762a1bSJed Brown     /*   Compute auxiliary data (solution and projection of adjoint solution) */
3959566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dmAdj, &uAdjLoc));
3969566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
3979566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
3989566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &uAdjProj));
3999566063dSJacob Faibussowitsch     PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc));
4009566063dSJacob Faibussowitsch     PetscCall(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj));
4019566063dSJacob Faibussowitsch     PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
4029566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dmAdj, &uAdjLoc));
403c4762a1bSJed Brown     /*   Attach auxiliary data */
4049371c9d4SSatish Balay     dms[0] = dm;
4059371c9d4SSatish Balay     dms[1] = dm;
4069566063dSJacob Faibussowitsch     PetscCall(DMCreateSuperDM(dms, 2, &subis, &dmErrAux));
407c4762a1bSJed Brown     if (0) {
408c4762a1bSJed Brown       PetscSection sec;
409c4762a1bSJed Brown 
4109566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dms[0], &sec));
4119566063dSJacob Faibussowitsch       PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
4129566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dms[1], &sec));
4139566063dSJacob Faibussowitsch       PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
4149566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dmErrAux, &sec));
4159566063dSJacob Faibussowitsch       PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
416c4762a1bSJed Brown     }
4179566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view"));
4189566063dSJacob Faibussowitsch     PetscCall(ISViewFromOptions(subis[0], NULL, "-super_is_view"));
4199566063dSJacob Faibussowitsch     PetscCall(ISViewFromOptions(subis[1], NULL, "-super_is_view"));
4209566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmErrAux, &uErr));
4219566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-map_vec_view"));
4229566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view"));
4239566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(uErr, NULL, "-map_vec_view"));
4249566063dSJacob Faibussowitsch     PetscCall(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u));
4259566063dSJacob Faibussowitsch     PetscCall(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj));
4269566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &uAdjProj));
427157a3229SMatthew G. Knepley     for (i = 0; i < 2; ++i) PetscCall(ISDestroy(&subis[i]));
4289566063dSJacob Faibussowitsch     PetscCall(PetscFree(subis));
4299566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dmErrAux, &uErrLoc));
4309566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc));
4319566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc));
4329566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmErrAux, &uErr));
4339566063dSJacob Faibussowitsch     PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc));
434c4762a1bSJed Brown     /*   Compute cellwise error estimate */
4359566063dSJacob Faibussowitsch     PetscCall(VecSet(errorEst, 0.0));
436e8e188d2SZach Atkins     PetscCall(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, user));
4379566063dSJacob Faibussowitsch     PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL));
4389566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dmErrAux, &uErrLoc));
4399566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmErrAux));
440c4762a1bSJed Brown     /*   Plot cellwise error vector */
4419566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(errorEst, NULL, "-error_view"));
442c4762a1bSJed Brown     /*   Compute ratio of estimate (sum over cells) with actual L_2 error */
4439566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm));
4449566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2));
4459566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(errorL2, NULL, "-l2_error_view"));
4469566063dSJacob Faibussowitsch     PetscCall(VecNorm(errorL2, NORM_INFINITY, &errorL2Tot));
4479566063dSJacob Faibussowitsch     PetscCall(VecNorm(errorEst, NORM_INFINITY, &errorEstTot));
4489566063dSJacob Faibussowitsch     PetscCall(VecGetSize(errorEst, &N));
4499566063dSJacob Faibussowitsch     PetscCall(VecPointwiseDivide(errorEst, errorEst, errorL2));
4509566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)errorEst, "Error ratio"));
4519566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(errorEst, NULL, "-error_ratio_view"));
45263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g Error Ratio: %g/%g = %g\n", N, (double)errorL2Norm, (double)errorEstTot, (double)PetscSqrtReal(errorL2Tot), (double)(errorEstTot / PetscSqrtReal(errorL2Tot))));
4539566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmErr, &errorEst));
4549566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmErr, &errorL2));
4559566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmErr));
456c4762a1bSJed Brown   }
4579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmAdj));
4589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&uAdj));
4599566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snesAdj));
460e8e188d2SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
461c4762a1bSJed Brown }
462e8e188d2SZach Atkins 
ErrorView(Vec u,AppCtx * user)463e8e188d2SZach Atkins static PetscErrorCode ErrorView(Vec u, AppCtx *user)
464e8e188d2SZach Atkins {
465e8e188d2SZach Atkins   PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
466e8e188d2SZach Atkins   void     *ctx;
467e8e188d2SZach Atkins   DM        dm;
468e8e188d2SZach Atkins   PetscDS   ds;
469e8e188d2SZach Atkins   PetscReal error;
470e8e188d2SZach Atkins   PetscInt  N;
471e8e188d2SZach Atkins 
472e8e188d2SZach Atkins   PetscFunctionBegin;
473e8e188d2SZach Atkins   if (!user->viewError) PetscFunctionReturn(PETSC_SUCCESS);
474e8e188d2SZach Atkins   PetscCall(VecGetDM(u, &dm));
475e8e188d2SZach Atkins   PetscCall(DMGetDS(dm, &ds));
476e8e188d2SZach Atkins   PetscCall(PetscDSGetExactSolution(ds, 0, &sol, &ctx));
477e8e188d2SZach Atkins   PetscCall(VecGetSize(u, &N));
478e8e188d2SZach Atkins   PetscCall(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error));
479e8e188d2SZach Atkins   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g\n", N, (double)error));
480e8e188d2SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
481e8e188d2SZach Atkins }
482e8e188d2SZach Atkins 
main(int argc,char ** argv)483e8e188d2SZach Atkins int main(int argc, char **argv)
484e8e188d2SZach Atkins {
485e8e188d2SZach Atkins   DM        dm;   /* Problem specification */
486e8e188d2SZach Atkins   SNES      snes; /* Nonlinear solver */
487e8e188d2SZach Atkins   Vec       u;    /* Solutions */
488e8e188d2SZach Atkins   AppCtx    user; /* User-defined work context */
489e8e188d2SZach Atkins   PetscInt  planeDir[2]   = {0, 1};
490e8e188d2SZach Atkins   PetscReal planeCoord[2] = {0., 1.};
491e8e188d2SZach Atkins 
492e8e188d2SZach Atkins   PetscFunctionBeginUser;
493e8e188d2SZach Atkins   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
494e8e188d2SZach Atkins   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
495e8e188d2SZach Atkins   /* Primal system */
496e8e188d2SZach Atkins   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
497e8e188d2SZach Atkins   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
498e8e188d2SZach Atkins   PetscCall(SNESSetDM(snes, dm));
499e8e188d2SZach Atkins   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
500e8e188d2SZach Atkins   PetscCall(DMCreateGlobalVector(dm, &u));
501e8e188d2SZach Atkins   PetscCall(VecSet(u, 0.0));
502e8e188d2SZach Atkins   PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
5036493148fSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
504e8e188d2SZach Atkins   PetscCall(SNESSetFromOptions(snes));
505e8e188d2SZach Atkins   PetscCall(SNESSolve(snes, NULL, u));
506e8e188d2SZach Atkins   PetscCall(SNESGetSolution(snes, &u));
507e8e188d2SZach Atkins   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
508e8e188d2SZach Atkins   PetscCall(ErrorView(u, &user));
509e8e188d2SZach Atkins   PetscCall(ComputeSpectral(u, 2, planeDir, planeCoord, &user));
510e8e188d2SZach Atkins   PetscCall(ComputeAdjoint(u, &user));
511c4762a1bSJed Brown   /* Cleanup */
5129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
5139566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
5149566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
5159566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
516b122ec5aSJacob Faibussowitsch   return 0;
517c4762a1bSJed Brown }
518c4762a1bSJed Brown 
519c4762a1bSJed Brown /*TEST
520c4762a1bSJed Brown 
521c4762a1bSJed Brown   test:
5225be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5235be41b18SMatthew G. Knepley     suffix: 2d_p1_conv
524c4762a1bSJed Brown     requires: triangle
5255be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5265be41b18SMatthew G. Knepley   test:
5275be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5285be41b18SMatthew G. Knepley     suffix: 2d_p2_conv
5295be41b18SMatthew G. Knepley     requires: triangle
5305be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5315be41b18SMatthew G. Knepley   test:
5325be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5335be41b18SMatthew G. Knepley     suffix: 2d_p3_conv
5345be41b18SMatthew G. Knepley     requires: triangle
5355be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5365be41b18SMatthew G. Knepley   test:
5375be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5385be41b18SMatthew G. Knepley     suffix: 2d_q1_conv
53930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5405be41b18SMatthew G. Knepley   test:
5415be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5425be41b18SMatthew G. Knepley     suffix: 2d_q2_conv
54330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5445be41b18SMatthew G. Knepley   test:
5455be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5465be41b18SMatthew G. Knepley     suffix: 2d_q3_conv
54730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5485be41b18SMatthew G. Knepley   test:
5495be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
550157a3229SMatthew G. Knepley     suffix: 2d_q1_ceed_conv
551157a3229SMatthew G. Knepley     requires: libceed
552157a3229SMatthew G. Knepley     args: -dm_plex_use_ceed -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
553157a3229SMatthew G. Knepley   test:
554157a3229SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
555157a3229SMatthew G. Knepley     suffix: 2d_q2_ceed_conv
556157a3229SMatthew G. Knepley     requires: libceed
557157a3229SMatthew G. Knepley     args: -dm_plex_use_ceed -dm_plex_simplex 0 -potential_petscspace_degree 2 -cdm_default_quadrature_order 2 \
558157a3229SMatthew G. Knepley           -snes_convergence_estimate -convest_num_refine 2
559157a3229SMatthew G. Knepley   test:
560157a3229SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
561157a3229SMatthew G. Knepley     suffix: 2d_q3_ceed_conv
562157a3229SMatthew G. Knepley     requires: libceed
563157a3229SMatthew G. Knepley     args: -dm_plex_use_ceed -dm_plex_simplex 0 -potential_petscspace_degree 3 -cdm_default_quadrature_order 3 \
564157a3229SMatthew G. Knepley           -snes_convergence_estimate -convest_num_refine 2
565157a3229SMatthew G. Knepley   test:
566157a3229SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5675be41b18SMatthew G. Knepley     suffix: 2d_q1_shear_conv
56830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5695be41b18SMatthew G. Knepley   test:
5705be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5715be41b18SMatthew G. Knepley     suffix: 2d_q2_shear_conv
57230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5735be41b18SMatthew G. Knepley   test:
5745be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5755be41b18SMatthew G. Knepley     suffix: 2d_q3_shear_conv
57630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5775be41b18SMatthew G. Knepley   test:
5780fdc7489SMatthew Knepley     # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7
5795be41b18SMatthew G. Knepley     suffix: 3d_p1_conv
5805be41b18SMatthew G. Knepley     requires: ctetgen
58130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
5825be41b18SMatthew G. Knepley   test:
5835be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5845be41b18SMatthew G. Knepley     suffix: 3d_p2_conv
5855be41b18SMatthew G. Knepley     requires: ctetgen
58630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
5875be41b18SMatthew G. Knepley   test:
5885be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0
5895be41b18SMatthew G. Knepley     suffix: 3d_p3_conv
5905be41b18SMatthew G. Knepley     requires: ctetgen
59130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
5925be41b18SMatthew G. Knepley   test:
5935be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8
5945be41b18SMatthew G. Knepley     suffix: 3d_q1_conv
59530602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
5965be41b18SMatthew G. Knepley   test:
5975be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5985be41b18SMatthew G. Knepley     suffix: 3d_q2_conv
59930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
6005be41b18SMatthew G. Knepley   test:
6015be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8
6025be41b18SMatthew G. Knepley     suffix: 3d_q3_conv
60330602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
6049139b27eSToby Isaac   test:
6059139b27eSToby Isaac     suffix: 2d_p1_fas_full
6069139b27eSToby Isaac     requires: triangle
607e600fa54SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
6089139b27eSToby Isaac       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \
6099139b27eSToby Isaac         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
6109139b27eSToby Isaac         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
6119139b27eSToby Isaac           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
61273f7197eSJed Brown             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10
6139139b27eSToby Isaac   test:
6149139b27eSToby Isaac     suffix: 2d_p1_fas_full_homogeneous
6159139b27eSToby Isaac     requires: triangle
616e600fa54SMatthew G. Knepley     args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
6179139b27eSToby Isaac       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \
6189139b27eSToby Isaac         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
6199139b27eSToby Isaac         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
6209139b27eSToby Isaac           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
62173f7197eSJed Brown             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10
6225be41b18SMatthew G. Knepley 
623c4762a1bSJed Brown   test:
624c4762a1bSJed Brown     suffix: 2d_p1_scalable
6255be41b18SMatthew G. Knepley     requires: triangle
6265be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 3 \
627c4762a1bSJed Brown       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \
62873f7197eSJed Brown       -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \
629c4762a1bSJed Brown         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
630c4762a1bSJed Brown         -pc_gamg_coarse_eq_limit 1000 \
631c4762a1bSJed Brown         -pc_gamg_threshold 0.05 \
632c4762a1bSJed Brown         -pc_gamg_threshold_scale .0 \
633c4762a1bSJed Brown         -mg_levels_ksp_type chebyshev \
634c4762a1bSJed Brown         -mg_levels_ksp_max_it 1 \
635c4762a1bSJed Brown         -mg_levels_pc_type jacobi \
636c4762a1bSJed Brown       -matptap_via scalable
6373886731fSPierre Jolivet     output_file: output/empty.out
638c4762a1bSJed Brown   test:
639c4762a1bSJed Brown     suffix: 2d_p1_gmg_vcycle
640c4762a1bSJed Brown     requires: triangle
6413886731fSPierre Jolivet     output_file: output/empty.out
6425be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
643c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg \
644c4762a1bSJed Brown             -mg_levels_ksp_max_it 1 \
645c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
646c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
64773f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
648c4762a1bSJed Brown             -mg_levels_pc_type jacobi
64956245cf9SMatthew G. Knepley   # Run with -dm_refine_hierarchy 3 to get a better idea of the solver
6509a1f099bSMatthew G. Knepley   testset:
65156245cf9SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 2 \
652c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \
653c4762a1bSJed Brown             -mg_levels_ksp_max_it 2 \
654c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
655c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
65673f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
657c4762a1bSJed Brown             -mg_levels_pc_type jacobi
6583886731fSPierre Jolivet     output_file: output/empty.out
659c4762a1bSJed Brown     test:
6609a1f099bSMatthew G. Knepley       suffix: 2d_p1_gmg_fcycle
6619a1f099bSMatthew G. Knepley       requires: triangle
6629a1f099bSMatthew G. Knepley       args: -dm_plex_box_faces 2,2
6639a1f099bSMatthew G. Knepley     test:
6649a1f099bSMatthew G. Knepley       suffix: 2d_q1_gmg_fcycle
6659a1f099bSMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2
6669a1f099bSMatthew G. Knepley     test:
6679a1f099bSMatthew G. Knepley       suffix: 3d_p1_gmg_fcycle
6689a1f099bSMatthew G. Knepley       requires: ctetgen
66956245cf9SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,1
6709a1f099bSMatthew G. Knepley     test:
6719a1f099bSMatthew G. Knepley       suffix: 3d_q1_gmg_fcycle
67256245cf9SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,1
6739a1f099bSMatthew G. Knepley   test:
674d6837840SMatthew G. Knepley     suffix: 2d_p1_gmg_vcycle_adapt
6752b3cbbdaSStefano Zampini     requires: triangle
6763886731fSPierre Jolivet     output_file: output/empty.out
677908b9b43SStefano Zampini     args: -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
6782b3cbbdaSStefano Zampini           -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \
679d6837840SMatthew G. Knepley             -mg_levels_ksp_max_it 1 \
680d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_type cg \
681d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_max_it 10 \
68273f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
683d6837840SMatthew G. Knepley             -mg_levels_pc_type jacobi
684d6837840SMatthew G. Knepley   test:
685c4762a1bSJed Brown     suffix: 2d_p1_spectral_0
686c4762a1bSJed Brown     requires: triangle fftw !complex
6875be41b18SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view
688c4762a1bSJed Brown   test:
689c4762a1bSJed Brown     suffix: 2d_p1_spectral_1
690c4762a1bSJed Brown     requires: triangle fftw !complex
691c4762a1bSJed Brown     nsize: 2
692e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view
693c4762a1bSJed Brown   test:
694c4762a1bSJed Brown     suffix: 2d_p1_adj_0
695c4762a1bSJed Brown     requires: triangle
6965be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0
697b1ad28e3SJunchao Zhang   test:
698b1ad28e3SJunchao Zhang     nsize: 2
699dcfd994dSJunchao Zhang     requires: kokkos_kernels
700b1ad28e3SJunchao Zhang     suffix: kokkos
701e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,6 -petscpartitioner_type simple -dm_plex_simplex 0 -potential_petscspace_degree 1 \
702b1ad28e3SJunchao Zhang          -dm_refine 0 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_coarse_eq_limit 1000 -pc_gamg_threshold 0.0 \
70373f7197eSJed Brown          -pc_gamg_threshold_scale .5 -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 2 -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \
70473f7197eSJed Brown          -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos
705c4762a1bSJed Brown 
706c4762a1bSJed Brown TEST*/
707