xref: /petsc/src/snes/tutorials/ex23.c (revision d94ba0934425393760e91090bbd98ca594b056f0)
1*d94ba093SMatthew G. Knepley static char help[] = "Poisson Problem with a split domain.\n\
2*d94ba093SMatthew G. Knepley We solve the Poisson problem in two halves of a domain.\n\
3*d94ba093SMatthew G. Knepley In one half, we include an additional field.\n\n\n";
4*d94ba093SMatthew G. Knepley 
5*d94ba093SMatthew G. Knepley #include <petscdmplex.h>
6*d94ba093SMatthew G. Knepley #include <petscsnes.h>
7*d94ba093SMatthew G. Knepley #include <petscds.h>
8*d94ba093SMatthew G. Knepley #include <petscconvest.h>
9*d94ba093SMatthew G. Knepley 
10*d94ba093SMatthew G. Knepley typedef struct {
11*d94ba093SMatthew G. Knepley   /* Domain and mesh definition */
12*d94ba093SMatthew G. Knepley   PetscInt  dim;               /* The topological mesh dimension */
13*d94ba093SMatthew G. Knepley   PetscBool simplex;           /* Simplicial mesh */
14*d94ba093SMatthew G. Knepley } AppCtx;
15*d94ba093SMatthew G. Knepley 
16*d94ba093SMatthew G. Knepley static PetscErrorCode quad_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
17*d94ba093SMatthew G. Knepley {
18*d94ba093SMatthew G. Knepley   u[0] = x[0]*x[0] + x[1]*x[1];
19*d94ba093SMatthew G. Knepley   return 0;
20*d94ba093SMatthew G. Knepley }
21*d94ba093SMatthew G. Knepley 
22*d94ba093SMatthew G. Knepley static PetscErrorCode quad_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
23*d94ba093SMatthew G. Knepley {
24*d94ba093SMatthew G. Knepley   u[0] = 2.0;
25*d94ba093SMatthew G. Knepley   return 0;
26*d94ba093SMatthew G. Knepley }
27*d94ba093SMatthew G. Knepley 
28*d94ba093SMatthew G. Knepley static void f0_quad_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
29*d94ba093SMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
30*d94ba093SMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
31*d94ba093SMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
32*d94ba093SMatthew G. Knepley {
33*d94ba093SMatthew G. Knepley   PetscInt d;
34*d94ba093SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += 2.0;
35*d94ba093SMatthew G. Knepley }
36*d94ba093SMatthew G. Knepley 
37*d94ba093SMatthew G. Knepley static void f0_quad_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
38*d94ba093SMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
39*d94ba093SMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
40*d94ba093SMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
41*d94ba093SMatthew G. Knepley {
42*d94ba093SMatthew G. Knepley   f0[0] = u[uOff[1]] - 2.0;
43*d94ba093SMatthew G. Knepley }
44*d94ba093SMatthew G. Knepley 
45*d94ba093SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
46*d94ba093SMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
47*d94ba093SMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
48*d94ba093SMatthew G. Knepley                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
49*d94ba093SMatthew G. Knepley {
50*d94ba093SMatthew G. Knepley   PetscInt d;
51*d94ba093SMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
52*d94ba093SMatthew G. Knepley }
53*d94ba093SMatthew G. Knepley 
54*d94ba093SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
55*d94ba093SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
56*d94ba093SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
57*d94ba093SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
58*d94ba093SMatthew G. Knepley {
59*d94ba093SMatthew G. Knepley   PetscInt d;
60*d94ba093SMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
61*d94ba093SMatthew G. Knepley }
62*d94ba093SMatthew G. Knepley 
63*d94ba093SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
64*d94ba093SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
65*d94ba093SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
66*d94ba093SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
67*d94ba093SMatthew G. Knepley {
68*d94ba093SMatthew G. Knepley   g0[0] = 1.0;
69*d94ba093SMatthew G. Knepley }
70*d94ba093SMatthew G. Knepley 
71*d94ba093SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
72*d94ba093SMatthew G. Knepley {
73*d94ba093SMatthew G. Knepley   PetscErrorCode ierr;
74*d94ba093SMatthew G. Knepley 
75*d94ba093SMatthew G. Knepley   PetscFunctionBeginUser;
76*d94ba093SMatthew G. Knepley   options->dim     = 2;
77*d94ba093SMatthew G. Knepley   options->simplex = PETSC_TRUE;
78*d94ba093SMatthew G. Knepley 
79*d94ba093SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
80*d94ba093SMatthew G. Knepley   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex23.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
81*d94ba093SMatthew G. Knepley   ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex23.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
82*d94ba093SMatthew G. Knepley   ierr = PetscOptionsEnd();
83*d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
84*d94ba093SMatthew G. Knepley }
85*d94ba093SMatthew G. Knepley 
86*d94ba093SMatthew G. Knepley static PetscErrorCode DivideDomain(DM dm, AppCtx *user)
87*d94ba093SMatthew G. Knepley {
88*d94ba093SMatthew G. Knepley   DMLabel        top, bottom;
89*d94ba093SMatthew G. Knepley   PetscReal      low[3], high[3], midy;
90*d94ba093SMatthew G. Knepley   PetscInt       cStart, cEnd, c;
91*d94ba093SMatthew G. Knepley   PetscErrorCode ierr;
92*d94ba093SMatthew G. Knepley 
93*d94ba093SMatthew G. Knepley   PetscFunctionBeginUser;
94*d94ba093SMatthew G. Knepley   ierr = DMCreateLabel(dm, "top");CHKERRQ(ierr);
95*d94ba093SMatthew G. Knepley   ierr = DMCreateLabel(dm, "bottom");CHKERRQ(ierr);
96*d94ba093SMatthew G. Knepley   ierr = DMGetLabel(dm, "top", &top);CHKERRQ(ierr);
97*d94ba093SMatthew G. Knepley   ierr = DMGetLabel(dm, "bottom", &bottom);CHKERRQ(ierr);
98*d94ba093SMatthew G. Knepley   ierr = DMGetBoundingBox(dm, low, high);CHKERRQ(ierr);
99*d94ba093SMatthew G. Knepley   midy = 0.5*(high[1] - low[1]);
100*d94ba093SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
101*d94ba093SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
102*d94ba093SMatthew G. Knepley     PetscReal centroid[3];
103*d94ba093SMatthew G. Knepley 
104*d94ba093SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
105*d94ba093SMatthew G. Knepley     if (centroid[1] > midy) {ierr = DMLabelSetValue(top, c, 1);CHKERRQ(ierr);}
106*d94ba093SMatthew G. Knepley     else                    {ierr = DMLabelSetValue(bottom, c, 1);CHKERRQ(ierr);}
107*d94ba093SMatthew G. Knepley   }
108*d94ba093SMatthew G. Knepley   ierr = DMPlexLabelComplete(dm, top);CHKERRQ(ierr);
109*d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
110*d94ba093SMatthew G. Knepley }
111*d94ba093SMatthew G. Knepley 
112*d94ba093SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
113*d94ba093SMatthew G. Knepley {
114*d94ba093SMatthew G. Knepley   PetscErrorCode ierr;
115*d94ba093SMatthew G. Knepley 
116*d94ba093SMatthew G. Knepley   PetscFunctionBeginUser;
117*d94ba093SMatthew G. Knepley   /* Create box mesh */
118*d94ba093SMatthew G. Knepley   ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
119*d94ba093SMatthew G. Knepley   /* Distribute mesh over processes */
120*d94ba093SMatthew G. Knepley   {
121*d94ba093SMatthew G. Knepley     DM               dmDist = NULL;
122*d94ba093SMatthew G. Knepley     PetscPartitioner part;
123*d94ba093SMatthew G. Knepley 
124*d94ba093SMatthew G. Knepley     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
125*d94ba093SMatthew G. Knepley     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
126*d94ba093SMatthew G. Knepley     ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr);
127*d94ba093SMatthew G. Knepley     if (dmDist) {
128*d94ba093SMatthew G. Knepley       ierr = DMDestroy(dm);CHKERRQ(ierr);
129*d94ba093SMatthew G. Knepley       *dm  = dmDist;
130*d94ba093SMatthew G. Knepley     }
131*d94ba093SMatthew G. Knepley   }
132*d94ba093SMatthew G. Knepley   /* TODO: This should be pulled into the library */
133*d94ba093SMatthew G. Knepley   {
134*d94ba093SMatthew G. Knepley     char      convType[256];
135*d94ba093SMatthew G. Knepley     PetscBool flg;
136*d94ba093SMatthew G. Knepley 
137*d94ba093SMatthew G. Knepley     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
138*d94ba093SMatthew G. Knepley     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
139*d94ba093SMatthew G. Knepley     ierr = PetscOptionsEnd();
140*d94ba093SMatthew G. Knepley     if (flg) {
141*d94ba093SMatthew G. Knepley       DM dmConv;
142*d94ba093SMatthew G. Knepley 
143*d94ba093SMatthew G. Knepley       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
144*d94ba093SMatthew G. Knepley       if (dmConv) {
145*d94ba093SMatthew G. Knepley         ierr = DMDestroy(dm);CHKERRQ(ierr);
146*d94ba093SMatthew G. Knepley         *dm  = dmConv;
147*d94ba093SMatthew G. Knepley       }
148*d94ba093SMatthew G. Knepley     }
149*d94ba093SMatthew G. Knepley   }
150*d94ba093SMatthew G. Knepley   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
151*d94ba093SMatthew G. Knepley 
152*d94ba093SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
153*d94ba093SMatthew G. Knepley   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
154*d94ba093SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
155*d94ba093SMatthew G. Knepley   ierr = DivideDomain(*dm, user);CHKERRQ(ierr);
156*d94ba093SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
157*d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
158*d94ba093SMatthew G. Knepley }
159*d94ba093SMatthew G. Knepley 
160*d94ba093SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
161*d94ba093SMatthew G. Knepley {
162*d94ba093SMatthew G. Knepley   PetscDS        ds;
163*d94ba093SMatthew G. Knepley   const PetscInt id = 1;
164*d94ba093SMatthew G. Knepley   PetscErrorCode ierr;
165*d94ba093SMatthew G. Knepley 
166*d94ba093SMatthew G. Knepley   PetscFunctionBeginUser;
167*d94ba093SMatthew G. Knepley   ierr = DMGetRegionNumDS(dm, 0, NULL, NULL, &ds);CHKERRQ(ierr);
168*d94ba093SMatthew G. Knepley   ierr = PetscDSSetResidual(ds, 0, f0_quad_u, f1_u);CHKERRQ(ierr);
169*d94ba093SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
170*d94ba093SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr);
171*d94ba093SMatthew G. Knepley   ierr = DMGetRegionNumDS(dm, 1, NULL, NULL, &ds);CHKERRQ(ierr);
172*d94ba093SMatthew G. Knepley   ierr = PetscDSSetResidual(ds, 0, f0_quad_u, f1_u);CHKERRQ(ierr);
173*d94ba093SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
174*d94ba093SMatthew G. Knepley   ierr = PetscDSSetResidual(ds, 1, f0_quad_p, NULL);CHKERRQ(ierr);
175*d94ba093SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 1, 1, g0_pp, NULL, NULL, NULL);CHKERRQ(ierr);
176*d94ba093SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 0, quad_u, user);CHKERRQ(ierr);
177*d94ba093SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 1, quad_p, user);CHKERRQ(ierr);
178*d94ba093SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) quad_u, 1, &id, user);CHKERRQ(ierr);
179*d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
180*d94ba093SMatthew G. Knepley }
181*d94ba093SMatthew G. Knepley 
182*d94ba093SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
183*d94ba093SMatthew G. Knepley {
184*d94ba093SMatthew G. Knepley   DM              cdm = dm;
185*d94ba093SMatthew G. Knepley   DMLabel         top;
186*d94ba093SMatthew G. Knepley   PetscFE         fe, feTop;
187*d94ba093SMatthew G. Knepley   PetscQuadrature q;
188*d94ba093SMatthew G. Knepley   const char     *nameTop = "pressure";
189*d94ba093SMatthew G. Knepley   char            prefix[PETSC_MAX_PATH_LEN];
190*d94ba093SMatthew G. Knepley   PetscErrorCode  ierr;
191*d94ba093SMatthew G. Knepley 
192*d94ba093SMatthew G. Knepley   PetscFunctionBeginUser;
193*d94ba093SMatthew G. Knepley   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
194*d94ba093SMatthew G. Knepley   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
195*d94ba093SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
196*d94ba093SMatthew G. Knepley   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
197*d94ba093SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
198*d94ba093SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
199*d94ba093SMatthew G. Knepley   ierr = DMGetLabel(dm, "top", &top);CHKERRQ(ierr);
200*d94ba093SMatthew G. Knepley   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_top_", nameTop);CHKERRQ(ierr);
201*d94ba093SMatthew G. Knepley   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), user->dim, 1, user->simplex, name ? prefix : NULL, -1, &feTop);CHKERRQ(ierr);
202*d94ba093SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) feTop, nameTop);CHKERRQ(ierr);
203*d94ba093SMatthew G. Knepley   ierr = PetscFESetQuadrature(feTop, q);CHKERRQ(ierr);
204*d94ba093SMatthew G. Knepley   ierr = DMSetField(dm, 1, top, (PetscObject) feTop);CHKERRQ(ierr);
205*d94ba093SMatthew G. Knepley   ierr = PetscFEDestroy(&feTop);CHKERRQ(ierr);
206*d94ba093SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
207*d94ba093SMatthew G. Knepley   ierr = (*setup)(dm, user);CHKERRQ(ierr);
208*d94ba093SMatthew G. Knepley   while (cdm) {
209*d94ba093SMatthew G. Knepley     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
210*d94ba093SMatthew G. Knepley     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
211*d94ba093SMatthew G. Knepley   }
212*d94ba093SMatthew G. Knepley   PetscFunctionReturn(0);
213*d94ba093SMatthew G. Knepley }
214*d94ba093SMatthew G. Knepley 
215*d94ba093SMatthew G. Knepley int main(int argc, char **argv)
216*d94ba093SMatthew G. Knepley {
217*d94ba093SMatthew G. Knepley   DM             dm;   /* Problem specification */
218*d94ba093SMatthew G. Knepley   SNES           snes; /* Nonlinear solver */
219*d94ba093SMatthew G. Knepley   Vec            u;    /* Solutions */
220*d94ba093SMatthew G. Knepley   AppCtx         user; /* User-defined work context */
221*d94ba093SMatthew G. Knepley   PetscErrorCode ierr;
222*d94ba093SMatthew G. Knepley 
223*d94ba093SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
224*d94ba093SMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
225*d94ba093SMatthew G. Knepley   /* Primal system */
226*d94ba093SMatthew G. Knepley   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
227*d94ba093SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
228*d94ba093SMatthew G. Knepley   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
229*d94ba093SMatthew G. Knepley   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
230*d94ba093SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
231*d94ba093SMatthew G. Knepley   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
232*d94ba093SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr);
233*d94ba093SMatthew G. Knepley   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
234*d94ba093SMatthew G. Knepley   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
235*d94ba093SMatthew G. Knepley   ierr = DMSNESCheckFromOptions(snes, u, NULL, NULL);CHKERRQ(ierr);
236*d94ba093SMatthew G. Knepley   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
237*d94ba093SMatthew G. Knepley   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
238*d94ba093SMatthew G. Knepley   ierr = VecViewFromOptions(u, NULL, "-sol_view");CHKERRQ(ierr);
239*d94ba093SMatthew G. Knepley   /* Cleanup */
240*d94ba093SMatthew G. Knepley   ierr = VecDestroy(&u);CHKERRQ(ierr);
241*d94ba093SMatthew G. Knepley   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
242*d94ba093SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
243*d94ba093SMatthew G. Knepley   ierr = PetscFinalize();
244*d94ba093SMatthew G. Knepley   return ierr;
245*d94ba093SMatthew G. Knepley }
246*d94ba093SMatthew G. Knepley 
247*d94ba093SMatthew G. Knepley /*TEST
248*d94ba093SMatthew G. Knepley 
249*d94ba093SMatthew G. Knepley   test:
250*d94ba093SMatthew G. Knepley     suffix: 2d_p1_0
251*d94ba093SMatthew G. Knepley     requires: triangle
252*d94ba093SMatthew G. Knepley     args: -potential_petscspace_degree 2 -pressure_top_petscspace_degree 0 -dm_refine 0 -dmsnes_check
253*d94ba093SMatthew G. Knepley 
254*d94ba093SMatthew G. Knepley   test:
255*d94ba093SMatthew G. Knepley     suffix: 2d_p1_1
256*d94ba093SMatthew G. Knepley     requires: triangle
257*d94ba093SMatthew G. Knepley     args: -potential_petscspace_degree 1 -pressure_top_petscspace_degree 0 -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate
258*d94ba093SMatthew G. Knepley 
259*d94ba093SMatthew G. Knepley TEST*/
260