xref: /petsc/src/snes/tutorials/ex24.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Poisson Problem in mixed form with 2d and 3d with finite elements.\n\
2*c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\
3*c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4*c4762a1bSJed Brown This example supports automatic convergence estimation\n\
5*c4762a1bSJed Brown and Hdiv elements.\n\n\n";
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown #include <petscdmplex.h>
8*c4762a1bSJed Brown #include <petscsnes.h>
9*c4762a1bSJed Brown #include <petscds.h>
10*c4762a1bSJed Brown #include <petscconvest.h>
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown typedef enum {SOL_LINEAR, SOL_QUADRATIC, SOL_QUARTIC, SOL_UNKNOWN, NUM_SOLTYPE} SolType;
13*c4762a1bSJed Brown const char *SolTypeNames[NUM_SOLTYPE+3] = {"linear", "quadratic", "quartic", "unknown", "SolType", "SOL_", NULL};
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown typedef struct {
16*c4762a1bSJed Brown   /* Domain and mesh definition */
17*c4762a1bSJed Brown   PetscInt  dim;     /* The topological mesh dimension */
18*c4762a1bSJed Brown   PetscBool simplex; /* Simplicial mesh */
19*c4762a1bSJed Brown   SolType   solType; /* The type of exact solution */
20*c4762a1bSJed Brown } AppCtx;
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown /* 2D Dirichlet potential example
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   u = x
25*c4762a1bSJed Brown   q = <1, 0>
26*c4762a1bSJed Brown   f = 0
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   We will need a boundary integral of u over \Gamma.
29*c4762a1bSJed Brown */
30*c4762a1bSJed Brown static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
31*c4762a1bSJed Brown {
32*c4762a1bSJed Brown   u[0] = x[0];
33*c4762a1bSJed Brown   return 0;
34*c4762a1bSJed Brown }
35*c4762a1bSJed Brown static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
36*c4762a1bSJed Brown {
37*c4762a1bSJed Brown   PetscInt c;
38*c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0;
39*c4762a1bSJed Brown   return 0;
40*c4762a1bSJed Brown }
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown /* 2D Dirichlet potential example
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   u = x^2 + y^2
45*c4762a1bSJed Brown   q = <2x, 2y>
46*c4762a1bSJed Brown   f = 4
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   We will need a boundary integral of u over \Gamma.
49*c4762a1bSJed Brown */
50*c4762a1bSJed Brown static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51*c4762a1bSJed Brown {
52*c4762a1bSJed Brown   PetscInt d;
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown   u[0] = 0.0;
55*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[0] += x[d]*x[d];
56*c4762a1bSJed Brown   return 0;
57*c4762a1bSJed Brown }
58*c4762a1bSJed Brown static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
59*c4762a1bSJed Brown {
60*c4762a1bSJed Brown   PetscInt c;
61*c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = 2.0*x[c];
62*c4762a1bSJed Brown   return 0;
63*c4762a1bSJed Brown }
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown /* 2D Dirichlet potential example
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown   u = x (1-x) y (1-y)
68*c4762a1bSJed Brown   q = <(1-2x) y (1-y), x (1-x) (1-2y)>
69*c4762a1bSJed Brown   f = -y (1-y) - x (1-x)
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   u|_\Gamma = 0 so that the boundary integral vanishes
72*c4762a1bSJed Brown */
73*c4762a1bSJed Brown static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
74*c4762a1bSJed Brown {
75*c4762a1bSJed Brown   PetscInt d;
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   u[0] = 1.0;
78*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[0] *= x[d]*(1.0 - x[d]);
79*c4762a1bSJed Brown   return 0;
80*c4762a1bSJed Brown }
81*c4762a1bSJed Brown static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
82*c4762a1bSJed Brown {
83*c4762a1bSJed Brown   PetscInt c, d;
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
86*c4762a1bSJed Brown     u[c] = 1.0;
87*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
88*c4762a1bSJed Brown       if (c == d) u[c] *= 1 - 2.0*x[d];
89*c4762a1bSJed Brown       else        u[c] *= x[d]*(1.0 - x[d]);
90*c4762a1bSJed Brown     }
91*c4762a1bSJed Brown   }
92*c4762a1bSJed Brown   return 0;
93*c4762a1bSJed Brown }
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown /* <v, -\nabla\cdot q> + <v, f> */
96*c4762a1bSJed Brown static void f0_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
97*c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
98*c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
99*c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
100*c4762a1bSJed Brown {
101*c4762a1bSJed Brown   f0[0] = 0.0;
102*c4762a1bSJed Brown }
103*c4762a1bSJed Brown static void f0_bd_linear_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
104*c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
105*c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106*c4762a1bSJed Brown                            PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
107*c4762a1bSJed Brown {
108*c4762a1bSJed Brown   PetscScalar potential;
109*c4762a1bSJed Brown   PetscInt    d;
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown   linear_u(dim, t, x, dim, &potential, NULL);
112*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[d] = -potential*n[d];
113*c4762a1bSJed Brown }
114*c4762a1bSJed Brown static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
115*c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
116*c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
117*c4762a1bSJed Brown                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
118*c4762a1bSJed Brown {
119*c4762a1bSJed Brown   PetscInt d;
120*c4762a1bSJed Brown   f0[0] = 0.0;
121*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
122*c4762a1bSJed Brown     f0[0] -= u_x[uOff_x[0]+d*dim+d];
123*c4762a1bSJed Brown   }
124*c4762a1bSJed Brown   f0[0] += 4.0;
125*c4762a1bSJed Brown }
126*c4762a1bSJed Brown static void f0_bd_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
127*c4762a1bSJed Brown                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
128*c4762a1bSJed Brown                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
129*c4762a1bSJed Brown                               PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
130*c4762a1bSJed Brown {
131*c4762a1bSJed Brown   PetscScalar potential;
132*c4762a1bSJed Brown   PetscInt    d;
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown   quadratic_u(dim, t, x, dim, &potential, NULL);
135*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[d] = -potential*n[d];
136*c4762a1bSJed Brown }
137*c4762a1bSJed Brown static void f0_quartic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
138*c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
139*c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
140*c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
141*c4762a1bSJed Brown {
142*c4762a1bSJed Brown   PetscInt d;
143*c4762a1bSJed Brown   f0[0] = 0.0;
144*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
145*c4762a1bSJed Brown     f0[0] -= u_x[uOff_x[0]+d*dim+d];
146*c4762a1bSJed Brown     f0[0] += -2.0*x[d]*(1.0 - x[d]);
147*c4762a1bSJed Brown   }
148*c4762a1bSJed Brown }
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown /* <w, q> */
151*c4762a1bSJed Brown static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
152*c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
153*c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
154*c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
155*c4762a1bSJed Brown {
156*c4762a1bSJed Brown   PetscInt c;
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown   for (c = 0; c < dim; ++c) {
159*c4762a1bSJed Brown     f0[c] = u[uOff[0]+c];
160*c4762a1bSJed Brown   }
161*c4762a1bSJed Brown }
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */
164*c4762a1bSJed Brown static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
165*c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
166*c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
167*c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
168*c4762a1bSJed Brown {
169*c4762a1bSJed Brown   PetscInt c, d;
170*c4762a1bSJed Brown   for (c = 0; c < dim; ++c) {
171*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
172*c4762a1bSJed Brown       if (c == d) f1[c*dim+d] = u[uOff[1]];
173*c4762a1bSJed Brown     }
174*c4762a1bSJed Brown   }
175*c4762a1bSJed Brown }
176*c4762a1bSJed Brown 
177*c4762a1bSJed Brown /* <w, q> */
178*c4762a1bSJed Brown static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux,
179*c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
180*c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
181*c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
182*c4762a1bSJed Brown {
183*c4762a1bSJed Brown   PetscInt c;
184*c4762a1bSJed Brown   for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0;
185*c4762a1bSJed Brown }
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */
188*c4762a1bSJed Brown static void g2_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
189*c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
190*c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
191*c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
192*c4762a1bSJed Brown {
193*c4762a1bSJed Brown   PetscInt d;
194*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d*dim+d] = 1.0;
195*c4762a1bSJed Brown }
196*c4762a1bSJed Brown 
197*c4762a1bSJed Brown /* <v, -\nabla\cdot q> */
198*c4762a1bSJed Brown static void g1_uq(PetscInt dim, PetscInt Nf, PetscInt NfAux,
199*c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
200*c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
201*c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
202*c4762a1bSJed Brown {
203*c4762a1bSJed Brown   PetscInt d;
204*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0;
205*c4762a1bSJed Brown }
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
208*c4762a1bSJed Brown {
209*c4762a1bSJed Brown   PetscErrorCode ierr;
210*c4762a1bSJed Brown 
211*c4762a1bSJed Brown   PetscFunctionBeginUser;
212*c4762a1bSJed Brown   options->dim     = 2;
213*c4762a1bSJed Brown   options->simplex = PETSC_TRUE;
214*c4762a1bSJed Brown   options->solType = SOL_LINEAR;
215*c4762a1bSJed Brown 
216*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex24.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex24.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum) options->solType, (PetscEnum *) &options->solType, NULL);CHKERRQ(ierr);
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
222*c4762a1bSJed Brown   PetscFunctionReturn(0);
223*c4762a1bSJed Brown }
224*c4762a1bSJed Brown 
225*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
226*c4762a1bSJed Brown {
227*c4762a1bSJed Brown   PetscErrorCode ierr;
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   PetscFunctionBeginUser;
230*c4762a1bSJed Brown   if (0) {
231*c4762a1bSJed Brown     DMLabel     label;
232*c4762a1bSJed Brown     const char *name = "marker";
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown     ierr = DMPlexCreateReferenceCell(comm, user->dim, user->simplex, dm);CHKERRQ(ierr);
235*c4762a1bSJed Brown     ierr = DMCreateLabel(*dm, name);CHKERRQ(ierr);
236*c4762a1bSJed Brown     ierr = DMGetLabel(*dm, name, &label);CHKERRQ(ierr);
237*c4762a1bSJed Brown     ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr);
238*c4762a1bSJed Brown     ierr = DMPlexLabelComplete(*dm, label);CHKERRQ(ierr);
239*c4762a1bSJed Brown   } else {
240*c4762a1bSJed Brown     /* Create box mesh */
241*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
242*c4762a1bSJed Brown   }
243*c4762a1bSJed Brown   /* Distribute mesh over processes */
244*c4762a1bSJed Brown   {
245*c4762a1bSJed Brown     DM               dmDist = NULL;
246*c4762a1bSJed Brown     PetscPartitioner part;
247*c4762a1bSJed Brown 
248*c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
249*c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
250*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr);
251*c4762a1bSJed Brown     if (dmDist) {
252*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
253*c4762a1bSJed Brown       *dm  = dmDist;
254*c4762a1bSJed Brown     }
255*c4762a1bSJed Brown   }
256*c4762a1bSJed Brown   /* TODO: This should be pulled into the library */
257*c4762a1bSJed Brown   {
258*c4762a1bSJed Brown     char      convType[256];
259*c4762a1bSJed Brown     PetscBool flg;
260*c4762a1bSJed Brown 
261*c4762a1bSJed Brown     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
262*c4762a1bSJed Brown     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
263*c4762a1bSJed Brown     ierr = PetscOptionsEnd();
264*c4762a1bSJed Brown     if (flg) {
265*c4762a1bSJed Brown       DM dmConv;
266*c4762a1bSJed Brown 
267*c4762a1bSJed Brown       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
268*c4762a1bSJed Brown       if (dmConv) {
269*c4762a1bSJed Brown         ierr = DMDestroy(dm);CHKERRQ(ierr);
270*c4762a1bSJed Brown         *dm  = dmConv;
271*c4762a1bSJed Brown       }
272*c4762a1bSJed Brown     }
273*c4762a1bSJed Brown   }
274*c4762a1bSJed Brown   /* TODO: This should be pulled into the library */
275*c4762a1bSJed Brown   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
276*c4762a1bSJed Brown 
277*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
278*c4762a1bSJed Brown   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
279*c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
281*c4762a1bSJed Brown   PetscFunctionReturn(0);
282*c4762a1bSJed Brown }
283*c4762a1bSJed Brown 
284*c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
285*c4762a1bSJed Brown {
286*c4762a1bSJed Brown   PetscDS        prob;
287*c4762a1bSJed Brown   const PetscInt id = 1;
288*c4762a1bSJed Brown   PetscErrorCode ierr;
289*c4762a1bSJed Brown 
290*c4762a1bSJed Brown   PetscFunctionBeginUser;
291*c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
292*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 0, f0_q, f1_q);CHKERRQ(ierr);
293*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, g0_qq, NULL, NULL, NULL);CHKERRQ(ierr);
294*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_qu, NULL);CHKERRQ(ierr);
295*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_uq, NULL, NULL);CHKERRQ(ierr);
296*c4762a1bSJed Brown   switch (user->solType)
297*c4762a1bSJed Brown   {
298*c4762a1bSJed Brown     case SOL_LINEAR:
299*c4762a1bSJed Brown       ierr = PetscDSSetResidual(prob, 1, f0_linear_u, NULL);CHKERRQ(ierr);
300*c4762a1bSJed Brown       ierr = PetscDSSetBdResidual(prob, 0, f0_bd_linear_q, NULL);CHKERRQ(ierr);
301*c4762a1bSJed Brown       ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL, "Dirichlet Bd Integral", "marker", 0, 0, NULL, (void (*)(void)) NULL, 1, &id, user);CHKERRQ(ierr);
302*c4762a1bSJed Brown       ierr = PetscDSSetExactSolution(prob, 0, linear_q, user);CHKERRQ(ierr);
303*c4762a1bSJed Brown       ierr = PetscDSSetExactSolution(prob, 1, linear_u, user);CHKERRQ(ierr);
304*c4762a1bSJed Brown       break;
305*c4762a1bSJed Brown     case SOL_QUADRATIC:
306*c4762a1bSJed Brown       ierr = PetscDSSetResidual(prob, 1, f0_quadratic_u, NULL);CHKERRQ(ierr);
307*c4762a1bSJed Brown       ierr = PetscDSSetBdResidual(prob, 0, f0_bd_quadratic_q, NULL);CHKERRQ(ierr);
308*c4762a1bSJed Brown       ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL, "Dirichlet Bd Integral", "marker", 0, 0, NULL, (void (*)(void)) NULL, 1, &id, user);CHKERRQ(ierr);
309*c4762a1bSJed Brown       ierr = PetscDSSetExactSolution(prob, 0, quadratic_q, user);CHKERRQ(ierr);
310*c4762a1bSJed Brown       ierr = PetscDSSetExactSolution(prob, 1, quadratic_u, user);CHKERRQ(ierr);
311*c4762a1bSJed Brown       break;
312*c4762a1bSJed Brown     case SOL_QUARTIC:
313*c4762a1bSJed Brown       ierr = PetscDSSetResidual(prob, 1, f0_quartic_u, NULL);CHKERRQ(ierr);
314*c4762a1bSJed Brown       ierr = PetscDSSetExactSolution(prob, 0, quartic_q, user);CHKERRQ(ierr);
315*c4762a1bSJed Brown       ierr = PetscDSSetExactSolution(prob, 1, quartic_u, user);CHKERRQ(ierr);
316*c4762a1bSJed Brown       break;
317*c4762a1bSJed Brown     default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]);
318*c4762a1bSJed Brown   }
319*c4762a1bSJed Brown   PetscFunctionReturn(0);
320*c4762a1bSJed Brown }
321*c4762a1bSJed Brown 
322*c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
323*c4762a1bSJed Brown {
324*c4762a1bSJed Brown   DM              cdm = dm;
325*c4762a1bSJed Brown   PetscFE         feq, feu;
326*c4762a1bSJed Brown   const PetscInt  dim = user->dim;
327*c4762a1bSJed Brown   PetscErrorCode  ierr;
328*c4762a1bSJed Brown 
329*c4762a1bSJed Brown   PetscFunctionBeginUser;
330*c4762a1bSJed Brown   /* Create finite element */
331*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, user->simplex, "field_",     -1, &feq);CHKERRQ(ierr);
332*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) feq, "field");CHKERRQ(ierr);
333*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1,   user->simplex, "potential_", -1, &feu);CHKERRQ(ierr);
334*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) feu, "potential");CHKERRQ(ierr);
335*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(feq, feu);CHKERRQ(ierr);
336*c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
337*c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) feq);CHKERRQ(ierr);
338*c4762a1bSJed Brown   ierr = DMSetField(dm, 1, NULL, (PetscObject) feu);CHKERRQ(ierr);
339*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = (*setup)(dm, user);CHKERRQ(ierr);
341*c4762a1bSJed Brown   while (cdm) {
342*c4762a1bSJed Brown     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
343*c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
344*c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
345*c4762a1bSJed Brown   }
346*c4762a1bSJed Brown   ierr = PetscFEDestroy(&feq);CHKERRQ(ierr);
347*c4762a1bSJed Brown   ierr = PetscFEDestroy(&feu);CHKERRQ(ierr);
348*c4762a1bSJed Brown   PetscFunctionReturn(0);
349*c4762a1bSJed Brown }
350*c4762a1bSJed Brown 
351*c4762a1bSJed Brown int main(int argc, char **argv)
352*c4762a1bSJed Brown {
353*c4762a1bSJed Brown   DM             dm;   /* Problem specification */
354*c4762a1bSJed Brown   SNES           snes; /* Nonlinear solver */
355*c4762a1bSJed Brown   Vec            u;    /* Solutions */
356*c4762a1bSJed Brown   AppCtx         user; /* User-defined work context */
357*c4762a1bSJed Brown   PetscErrorCode ierr;
358*c4762a1bSJed Brown 
359*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
360*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
361*c4762a1bSJed Brown   /* Primal system */
362*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
363*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
364*c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
365*c4762a1bSJed Brown   ierr = SetupDiscretization(dm, SetupPrimalProblem, &user);CHKERRQ(ierr);
366*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
367*c4762a1bSJed Brown   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
368*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
369*c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
370*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
371*c4762a1bSJed Brown   ierr = DMSNESCheckFromOptions(snes, u, NULL, NULL);CHKERRQ(ierr);
372*c4762a1bSJed Brown   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
373*c4762a1bSJed Brown   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
374*c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
375*c4762a1bSJed Brown   /* Cleanup */
376*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
377*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
378*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
379*c4762a1bSJed Brown   ierr = PetscFinalize();
380*c4762a1bSJed Brown   return ierr;
381*c4762a1bSJed Brown }
382*c4762a1bSJed Brown 
383*c4762a1bSJed Brown /*TEST
384*c4762a1bSJed Brown 
385*c4762a1bSJed Brown   test:
386*c4762a1bSJed Brown     suffix: 2d_bdm1_p0_0
387*c4762a1bSJed Brown     requires: triangle
388*c4762a1bSJed Brown     args: -sol_type linear \
389*c4762a1bSJed Brown           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
390*c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
391*c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
392*c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
393*c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
394*c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
395*c4762a1bSJed Brown   test:
396*c4762a1bSJed Brown     suffix: 2d_bdm1_p0_1
397*c4762a1bSJed Brown     requires: triangle
398*c4762a1bSJed Brown     args: -sol_type quadratic \
399*c4762a1bSJed Brown           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
400*c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
401*c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
402*c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
403*c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
404*c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
405*c4762a1bSJed Brown   test:
406*c4762a1bSJed Brown     suffix: 2d_bdm1_p0_2
407*c4762a1bSJed Brown     requires: triangle
408*c4762a1bSJed Brown     args: -sol_type quartic \
409*c4762a1bSJed Brown           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
410*c4762a1bSJed Brown           -snes_error_if_not_converged \
411*c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
412*c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
413*c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
414*c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
415*c4762a1bSJed Brown 
416*c4762a1bSJed Brown   test:
417*c4762a1bSJed Brown     suffix: 2d_p2_p0_vtk
418*c4762a1bSJed Brown     requires: triangle
419*c4762a1bSJed Brown     args: -sol_type linear \
420*c4762a1bSJed Brown           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
421*c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
422*c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
423*c4762a1bSJed Brown           -potential_view vtk: -exact_vec_view vtk: \
424*c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
425*c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
426*c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
427*c4762a1bSJed Brown 
428*c4762a1bSJed Brown   test:
429*c4762a1bSJed Brown     suffix: 2d_p2_p0_vtu
430*c4762a1bSJed Brown     requires: triangle
431*c4762a1bSJed Brown     args: -sol_type linear \
432*c4762a1bSJed Brown           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
433*c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
434*c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
435*c4762a1bSJed Brown           -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu \
436*c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
437*c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
438*c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
439*c4762a1bSJed Brown TEST*/
440*c4762a1bSJed Brown 
441*c4762a1bSJed Brown /* These tests will be active once tensor P^- is working
442*c4762a1bSJed Brown 
443*c4762a1bSJed Brown   test:
444*c4762a1bSJed Brown     suffix: 2d_bdmq1_p0_0
445*c4762a1bSJed Brown     requires: triangle
446*c4762a1bSJed Brown     args: -simplex 0 -sol_type linear \
447*c4762a1bSJed Brown           -field_petscspace_poly_type pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
448*c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
449*c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
450*c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
451*c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
452*c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
453*c4762a1bSJed Brown   test:
454*c4762a1bSJed Brown     suffix: 2d_bdmq1_p0_2
455*c4762a1bSJed Brown     requires: triangle
456*c4762a1bSJed Brown     args: -simplex 0 -sol_type quartic \
457*c4762a1bSJed Brown           -field_petscspace_poly_type_no pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
458*c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
459*c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
460*c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
461*c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
462*c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
463*c4762a1bSJed Brown 
464*c4762a1bSJed Brown */
465