xref: /petsc/src/snes/tutorials/ex4.c (revision 777c4855b3d26338f23a778d8270381b5f8ee82c)
1*777c4855SStefano Zampini #include <petscsnes.h>
2*777c4855SStefano Zampini #include <petscdmda.h>
3*777c4855SStefano Zampini 
4*777c4855SStefano Zampini static const char help[] = "Minimum surface area problem in 2D using DMDA.\n\
5*777c4855SStefano Zampini It solves an unconstrained minimization problem. This example is based on a \n\
6*777c4855SStefano Zampini problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\
7*777c4855SStefano Zampini boundary values along the edges of the domain, the objective is to find the\n\
8*777c4855SStefano Zampini surface with the minimal area that satisfies the boundary conditions.\n\
9*777c4855SStefano Zampini \n\
10*777c4855SStefano Zampini The command line options are:\n\
11*777c4855SStefano Zampini   -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
12*777c4855SStefano Zampini   -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
13*777c4855SStefano Zampini   \n";
14*777c4855SStefano Zampini 
15*777c4855SStefano Zampini /*
16*777c4855SStefano Zampini    User-defined application context - contains data needed by the
17*777c4855SStefano Zampini    application-provided call-back routines, FormJacobianLocal() and
18*777c4855SStefano Zampini    FormFunctionLocal().
19*777c4855SStefano Zampini */
20*777c4855SStefano Zampini 
21*777c4855SStefano Zampini typedef enum {
22*777c4855SStefano Zampini   PROBLEM_ENNEPER,
23*777c4855SStefano Zampini   PROBLEM_SINS,
24*777c4855SStefano Zampini } ProblemType;
25*777c4855SStefano Zampini static const char *const ProblemTypes[] = {"ENNEPER", "SINS", "ProblemType", "PROBLEM_", 0};
26*777c4855SStefano Zampini 
27*777c4855SStefano Zampini typedef struct {
28*777c4855SStefano Zampini   PetscScalar *bottom, *top, *left, *right;
29*777c4855SStefano Zampini } AppCtx;
30*777c4855SStefano Zampini 
31*777c4855SStefano Zampini /* -------- User-defined Routines --------- */
32*777c4855SStefano Zampini 
33*777c4855SStefano Zampini static PetscErrorCode FormBoundaryConditions_Enneper(SNES, AppCtx **);
34*777c4855SStefano Zampini static PetscErrorCode FormBoundaryConditions_Sins(SNES, AppCtx **);
35*777c4855SStefano Zampini static PetscErrorCode DestroyBoundaryConditions(AppCtx **);
36*777c4855SStefano Zampini static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *, PetscScalar **, PetscReal *, void *);
37*777c4855SStefano Zampini static PetscErrorCode FormFunctionLocal(DMDALocalInfo *, PetscScalar **, PetscScalar **, void *);
38*777c4855SStefano Zampini static PetscErrorCode FormJacobianLocal(DMDALocalInfo *, PetscScalar **, Mat, Mat, void *);
39*777c4855SStefano Zampini 
40*777c4855SStefano Zampini int main(int argc, char **argv)
41*777c4855SStefano Zampini {
42*777c4855SStefano Zampini   Vec         x;
43*777c4855SStefano Zampini   SNES        snes;
44*777c4855SStefano Zampini   DM          da;
45*777c4855SStefano Zampini   ProblemType ptype   = PROBLEM_ENNEPER;
46*777c4855SStefano Zampini   PetscBool   use_obj = PETSC_TRUE;
47*777c4855SStefano Zampini   PetscReal   bbox[4] = {0.};
48*777c4855SStefano Zampini   AppCtx     *user;
49*777c4855SStefano Zampini   PetscErrorCode (*form_bc)(SNES, AppCtx **) = NULL;
50*777c4855SStefano Zampini 
51*777c4855SStefano Zampini   PetscFunctionBeginUser;
52*777c4855SStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
53*777c4855SStefano Zampini 
54*777c4855SStefano Zampini   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Minimal surface options", __FILE__);
55*777c4855SStefano Zampini   PetscCall(PetscOptionsEnum("-problem_type", "Problem type", NULL, ProblemTypes, (PetscEnum)ptype, (PetscEnum *)&ptype, NULL));
56*777c4855SStefano Zampini   PetscCall(PetscOptionsBool("-use_objective", "Use objective function", NULL, use_obj, &use_obj, NULL));
57*777c4855SStefano Zampini   PetscOptionsEnd();
58*777c4855SStefano Zampini   switch (ptype) {
59*777c4855SStefano Zampini   case PROBLEM_ENNEPER:
60*777c4855SStefano Zampini     bbox[0] = -0.5;
61*777c4855SStefano Zampini     bbox[1] = 0.5;
62*777c4855SStefano Zampini     bbox[2] = -0.5;
63*777c4855SStefano Zampini     bbox[3] = 0.5;
64*777c4855SStefano Zampini     form_bc = FormBoundaryConditions_Enneper;
65*777c4855SStefano Zampini     break;
66*777c4855SStefano Zampini   case PROBLEM_SINS:
67*777c4855SStefano Zampini     bbox[0] = 0.0;
68*777c4855SStefano Zampini     bbox[1] = 1.0;
69*777c4855SStefano Zampini     bbox[2] = 0.0;
70*777c4855SStefano Zampini     bbox[3] = 1.0;
71*777c4855SStefano Zampini     form_bc = FormBoundaryConditions_Sins;
72*777c4855SStefano Zampini     break;
73*777c4855SStefano Zampini   }
74*777c4855SStefano Zampini 
75*777c4855SStefano Zampini   /* Create distributed array to manage the 2d grid */
76*777c4855SStefano Zampini   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
77*777c4855SStefano Zampini   PetscCall(DMSetFromOptions(da));
78*777c4855SStefano Zampini   PetscCall(DMSetUp(da));
79*777c4855SStefano Zampini   PetscCall(DMDASetUniformCoordinates(da, bbox[0], bbox[1], bbox[2], bbox[3], PETSC_DECIDE, PETSC_DECIDE));
80*777c4855SStefano Zampini 
81*777c4855SStefano Zampini   /* Extract global vectors from DMDA; */
82*777c4855SStefano Zampini   PetscCall(DMCreateGlobalVector(da, &x));
83*777c4855SStefano Zampini 
84*777c4855SStefano Zampini   /* Create nonlinear solver context */
85*777c4855SStefano Zampini   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
86*777c4855SStefano Zampini   PetscCall(SNESSetDM(snes, da));
87*777c4855SStefano Zampini   PetscCall((*form_bc)(snes, &user));
88*777c4855SStefano Zampini   PetscCall(SNESSetApplicationContext(snes, user));
89*777c4855SStefano Zampini 
90*777c4855SStefano Zampini   /*  Set local callbacks */
91*777c4855SStefano Zampini   if (use_obj) PetscCall(DMDASNESSetObjectiveLocal(da, (DMDASNESObjective)FormObjectiveLocal, NULL));
92*777c4855SStefano Zampini   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunction)FormFunctionLocal, NULL));
93*777c4855SStefano Zampini   PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobian)FormJacobianLocal, NULL));
94*777c4855SStefano Zampini 
95*777c4855SStefano Zampini   /* Customize from command line */
96*777c4855SStefano Zampini   PetscCall(SNESSetFromOptions(snes));
97*777c4855SStefano Zampini 
98*777c4855SStefano Zampini   /* Solve the application */
99*777c4855SStefano Zampini   PetscCall(SNESSolve(snes, NULL, x));
100*777c4855SStefano Zampini 
101*777c4855SStefano Zampini   /* Free user-created data structures */
102*777c4855SStefano Zampini   PetscCall(VecDestroy(&x));
103*777c4855SStefano Zampini   PetscCall(SNESDestroy(&snes));
104*777c4855SStefano Zampini   PetscCall(DMDestroy(&da));
105*777c4855SStefano Zampini   PetscCall(DestroyBoundaryConditions(&user));
106*777c4855SStefano Zampini 
107*777c4855SStefano Zampini   PetscCall(PetscFinalize());
108*777c4855SStefano Zampini   return 0;
109*777c4855SStefano Zampini }
110*777c4855SStefano Zampini 
111*777c4855SStefano Zampini /* Compute objective function over the locally owned part of the mesh */
112*777c4855SStefano Zampini PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info, PetscScalar **x, PetscReal *v, void *ptr)
113*777c4855SStefano Zampini {
114*777c4855SStefano Zampini   AppCtx     *user = (AppCtx *)ptr;
115*777c4855SStefano Zampini   PetscInt    mx = info->mx, my = info->my;
116*777c4855SStefano Zampini   PetscInt    xs = info->xs, xm = info->xm, ys = info->ys, ym = info->ym;
117*777c4855SStefano Zampini   PetscInt    i, j;
118*777c4855SStefano Zampini   PetscScalar hx, hy;
119*777c4855SStefano Zampini   PetscScalar f2, f4, d1, d2, d3, d4, xc, xl, xr, xt, xb;
120*777c4855SStefano Zampini   PetscReal   ft = 0, area;
121*777c4855SStefano Zampini 
122*777c4855SStefano Zampini   PetscFunctionBeginUser;
123*777c4855SStefano Zampini   PetscCheck(user, PetscObjectComm((PetscObject)info->da), PETSC_ERR_PLIB, "Missing application context");
124*777c4855SStefano Zampini   hx   = 1.0 / (mx + 1);
125*777c4855SStefano Zampini   hy   = 1.0 / (my + 1);
126*777c4855SStefano Zampini   area = 0.5 * hx * hy;
127*777c4855SStefano Zampini   for (j = ys; j < ys + ym; j++) {
128*777c4855SStefano Zampini     for (i = xs; i < xs + xm; i++) {
129*777c4855SStefano Zampini       xc = x[j][i];
130*777c4855SStefano Zampini       xl = xr = xb = xt = xc;
131*777c4855SStefano Zampini 
132*777c4855SStefano Zampini       if (i == 0) { /* left side */
133*777c4855SStefano Zampini         xl = user->left[j + 1];
134*777c4855SStefano Zampini       } else xl = x[j][i - 1];
135*777c4855SStefano Zampini 
136*777c4855SStefano Zampini       if (j == 0) { /* bottom side */
137*777c4855SStefano Zampini         xb = user->bottom[i + 1];
138*777c4855SStefano Zampini       } else xb = x[j - 1][i];
139*777c4855SStefano Zampini 
140*777c4855SStefano Zampini       if (i + 1 == mx) { /* right side */
141*777c4855SStefano Zampini         xr = user->right[j + 1];
142*777c4855SStefano Zampini       } else xr = x[j][i + 1];
143*777c4855SStefano Zampini 
144*777c4855SStefano Zampini       if (j + 1 == 0 + my) { /* top side */
145*777c4855SStefano Zampini         xt = user->top[i + 1];
146*777c4855SStefano Zampini       } else xt = x[j + 1][i];
147*777c4855SStefano Zampini 
148*777c4855SStefano Zampini       d1 = (xc - xl);
149*777c4855SStefano Zampini       d2 = (xc - xr);
150*777c4855SStefano Zampini       d3 = (xc - xt);
151*777c4855SStefano Zampini       d4 = (xc - xb);
152*777c4855SStefano Zampini 
153*777c4855SStefano Zampini       d1 /= hx;
154*777c4855SStefano Zampini       d2 /= hx;
155*777c4855SStefano Zampini       d3 /= hy;
156*777c4855SStefano Zampini       d4 /= hy;
157*777c4855SStefano Zampini 
158*777c4855SStefano Zampini       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
159*777c4855SStefano Zampini       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
160*777c4855SStefano Zampini 
161*777c4855SStefano Zampini       ft += PetscRealPart(f2 + f4);
162*777c4855SStefano Zampini     }
163*777c4855SStefano Zampini   }
164*777c4855SStefano Zampini 
165*777c4855SStefano Zampini   /* Compute triangular areas along the border of the domain. */
166*777c4855SStefano Zampini   if (xs == 0) { /* left side */
167*777c4855SStefano Zampini     for (j = ys; j < ys + ym; j++) {
168*777c4855SStefano Zampini       d3 = (user->left[j + 1] - user->left[j + 2]) / hy;
169*777c4855SStefano Zampini       d2 = (user->left[j + 1] - x[j][0]) / hx;
170*777c4855SStefano Zampini       ft += PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
171*777c4855SStefano Zampini     }
172*777c4855SStefano Zampini   }
173*777c4855SStefano Zampini   if (ys == 0) { /* bottom side */
174*777c4855SStefano Zampini     for (i = xs; i < xs + xm; i++) {
175*777c4855SStefano Zampini       d2 = (user->bottom[i + 1] - user->bottom[i + 2]) / hx;
176*777c4855SStefano Zampini       d3 = (user->bottom[i + 1] - x[0][i]) / hy;
177*777c4855SStefano Zampini       ft += PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
178*777c4855SStefano Zampini     }
179*777c4855SStefano Zampini   }
180*777c4855SStefano Zampini   if (xs + xm == mx) { /* right side */
181*777c4855SStefano Zampini     for (j = ys; j < ys + ym; j++) {
182*777c4855SStefano Zampini       d1 = (x[j][mx - 1] - user->right[j + 1]) / hx;
183*777c4855SStefano Zampini       d4 = (user->right[j] - user->right[j + 1]) / hy;
184*777c4855SStefano Zampini       ft += PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
185*777c4855SStefano Zampini     }
186*777c4855SStefano Zampini   }
187*777c4855SStefano Zampini   if (ys + ym == my) { /* top side */
188*777c4855SStefano Zampini     for (i = xs; i < xs + xm; i++) {
189*777c4855SStefano Zampini       d1 = (x[my - 1][i] - user->top[i + 1]) / hy;
190*777c4855SStefano Zampini       d4 = (user->top[i + 1] - user->top[i]) / hx;
191*777c4855SStefano Zampini       ft += PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
192*777c4855SStefano Zampini     }
193*777c4855SStefano Zampini   }
194*777c4855SStefano Zampini   if (ys == 0 && xs == 0) {
195*777c4855SStefano Zampini     d1 = (user->left[0] - user->left[1]) / hy;
196*777c4855SStefano Zampini     d2 = (user->bottom[0] - user->bottom[1]) / hx;
197*777c4855SStefano Zampini     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
198*777c4855SStefano Zampini   }
199*777c4855SStefano Zampini   if (ys + ym == my && xs + xm == mx) {
200*777c4855SStefano Zampini     d1 = (user->right[ym + 1] - user->right[ym]) / hy;
201*777c4855SStefano Zampini     d2 = (user->top[xm + 1] - user->top[xm]) / hx;
202*777c4855SStefano Zampini     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
203*777c4855SStefano Zampini   }
204*777c4855SStefano Zampini   ft *= area;
205*777c4855SStefano Zampini   *v = ft;
206*777c4855SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
207*777c4855SStefano Zampini }
208*777c4855SStefano Zampini 
209*777c4855SStefano Zampini /* Compute gradient over the locally owned part of the mesh */
210*777c4855SStefano Zampini PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **g, void *ptr)
211*777c4855SStefano Zampini {
212*777c4855SStefano Zampini   AppCtx     *user = (AppCtx *)ptr;
213*777c4855SStefano Zampini   PetscInt    mx = info->mx, my = info->my;
214*777c4855SStefano Zampini   PetscInt    xs = info->xs, xm = info->xm, ys = info->ys, ym = info->ym;
215*777c4855SStefano Zampini   PetscInt    i, j;
216*777c4855SStefano Zampini   PetscScalar hx, hy, hydhx, hxdhy;
217*777c4855SStefano Zampini   PetscScalar f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
218*777c4855SStefano Zampini   PetscScalar df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
219*777c4855SStefano Zampini 
220*777c4855SStefano Zampini   PetscFunctionBeginUser;
221*777c4855SStefano Zampini   PetscCheck(user, PetscObjectComm((PetscObject)info->da), PETSC_ERR_PLIB, "Missing application context");
222*777c4855SStefano Zampini   hx    = 1.0 / (mx + 1);
223*777c4855SStefano Zampini   hy    = 1.0 / (my + 1);
224*777c4855SStefano Zampini   hydhx = hy / hx;
225*777c4855SStefano Zampini   hxdhy = hx / hy;
226*777c4855SStefano Zampini 
227*777c4855SStefano Zampini   for (j = ys; j < ys + ym; j++) {
228*777c4855SStefano Zampini     for (i = xs; i < xs + xm; i++) {
229*777c4855SStefano Zampini       xc  = x[j][i];
230*777c4855SStefano Zampini       xlt = xrb = xl = xr = xb = xt = xc;
231*777c4855SStefano Zampini 
232*777c4855SStefano Zampini       if (i == 0) { /* left side */
233*777c4855SStefano Zampini         xl  = user->left[j + 1];
234*777c4855SStefano Zampini         xlt = user->left[j + 2];
235*777c4855SStefano Zampini       } else xl = x[j][i - 1];
236*777c4855SStefano Zampini 
237*777c4855SStefano Zampini       if (j == 0) { /* bottom side */
238*777c4855SStefano Zampini         xb  = user->bottom[i + 1];
239*777c4855SStefano Zampini         xrb = user->bottom[i + 2];
240*777c4855SStefano Zampini       } else xb = x[j - 1][i];
241*777c4855SStefano Zampini 
242*777c4855SStefano Zampini       if (i + 1 == mx) { /* right side */
243*777c4855SStefano Zampini         xr  = user->right[j + 1];
244*777c4855SStefano Zampini         xrb = user->right[j];
245*777c4855SStefano Zampini       } else xr = x[j][i + 1];
246*777c4855SStefano Zampini 
247*777c4855SStefano Zampini       if (j + 1 == 0 + my) { /* top side */
248*777c4855SStefano Zampini         xt  = user->top[i + 1];
249*777c4855SStefano Zampini         xlt = user->top[i];
250*777c4855SStefano Zampini       } else xt = x[j + 1][i];
251*777c4855SStefano Zampini 
252*777c4855SStefano Zampini       if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1]; /* left top side */
253*777c4855SStefano Zampini       if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1]; /* right bottom */
254*777c4855SStefano Zampini 
255*777c4855SStefano Zampini       d1 = (xc - xl);
256*777c4855SStefano Zampini       d2 = (xc - xr);
257*777c4855SStefano Zampini       d3 = (xc - xt);
258*777c4855SStefano Zampini       d4 = (xc - xb);
259*777c4855SStefano Zampini       d5 = (xr - xrb);
260*777c4855SStefano Zampini       d6 = (xrb - xb);
261*777c4855SStefano Zampini       d7 = (xlt - xl);
262*777c4855SStefano Zampini       d8 = (xt - xlt);
263*777c4855SStefano Zampini 
264*777c4855SStefano Zampini       df1dxc = d1 * hydhx;
265*777c4855SStefano Zampini       df2dxc = (d1 * hydhx + d4 * hxdhy);
266*777c4855SStefano Zampini       df3dxc = d3 * hxdhy;
267*777c4855SStefano Zampini       df4dxc = (d2 * hydhx + d3 * hxdhy);
268*777c4855SStefano Zampini       df5dxc = d2 * hydhx;
269*777c4855SStefano Zampini       df6dxc = d4 * hxdhy;
270*777c4855SStefano Zampini 
271*777c4855SStefano Zampini       d1 /= hx;
272*777c4855SStefano Zampini       d2 /= hx;
273*777c4855SStefano Zampini       d3 /= hy;
274*777c4855SStefano Zampini       d4 /= hy;
275*777c4855SStefano Zampini       d5 /= hy;
276*777c4855SStefano Zampini       d6 /= hx;
277*777c4855SStefano Zampini       d7 /= hy;
278*777c4855SStefano Zampini       d8 /= hx;
279*777c4855SStefano Zampini 
280*777c4855SStefano Zampini       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
281*777c4855SStefano Zampini       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
282*777c4855SStefano Zampini       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
283*777c4855SStefano Zampini       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
284*777c4855SStefano Zampini       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
285*777c4855SStefano Zampini       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
286*777c4855SStefano Zampini 
287*777c4855SStefano Zampini       df1dxc /= f1;
288*777c4855SStefano Zampini       df2dxc /= f2;
289*777c4855SStefano Zampini       df3dxc /= f3;
290*777c4855SStefano Zampini       df4dxc /= f4;
291*777c4855SStefano Zampini       df5dxc /= f5;
292*777c4855SStefano Zampini       df6dxc /= f6;
293*777c4855SStefano Zampini 
294*777c4855SStefano Zampini       g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0;
295*777c4855SStefano Zampini     }
296*777c4855SStefano Zampini   }
297*777c4855SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
298*777c4855SStefano Zampini }
299*777c4855SStefano Zampini 
300*777c4855SStefano Zampini /* Compute Hessian over the locally owned part of the mesh */
301*777c4855SStefano Zampini PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat H, Mat Hp, void *ptr)
302*777c4855SStefano Zampini {
303*777c4855SStefano Zampini   AppCtx     *user = (AppCtx *)ptr;
304*777c4855SStefano Zampini   PetscInt    mx = info->mx, my = info->my;
305*777c4855SStefano Zampini   PetscInt    xs = info->xs, xm = info->xm, ys = info->ys, ym = info->ym;
306*777c4855SStefano Zampini   PetscInt    i, j, k;
307*777c4855SStefano Zampini   MatStencil  row, col[7];
308*777c4855SStefano Zampini   PetscScalar hx, hy, hydhx, hxdhy;
309*777c4855SStefano Zampini   PetscScalar f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
310*777c4855SStefano Zampini   PetscScalar hl, hr, ht, hb, hc, htl, hbr;
311*777c4855SStefano Zampini   PetscScalar v[7];
312*777c4855SStefano Zampini 
313*777c4855SStefano Zampini   PetscFunctionBeginUser;
314*777c4855SStefano Zampini   PetscCheck(user, PetscObjectComm((PetscObject)info->da), PETSC_ERR_PLIB, "Missing application context");
315*777c4855SStefano Zampini   hx    = 1.0 / (mx + 1);
316*777c4855SStefano Zampini   hy    = 1.0 / (my + 1);
317*777c4855SStefano Zampini   hydhx = hy / hx;
318*777c4855SStefano Zampini   hxdhy = hx / hy;
319*777c4855SStefano Zampini 
320*777c4855SStefano Zampini   for (j = ys; j < ys + ym; j++) {
321*777c4855SStefano Zampini     for (i = xs; i < xs + xm; i++) {
322*777c4855SStefano Zampini       xc  = x[j][i];
323*777c4855SStefano Zampini       xlt = xrb = xl = xr = xb = xt = xc;
324*777c4855SStefano Zampini 
325*777c4855SStefano Zampini       /* Left */
326*777c4855SStefano Zampini       if (i == 0) {
327*777c4855SStefano Zampini         xl  = user->left[j + 1];
328*777c4855SStefano Zampini         xlt = user->left[j + 2];
329*777c4855SStefano Zampini       } else xl = x[j][i - 1];
330*777c4855SStefano Zampini 
331*777c4855SStefano Zampini       /* Bottom */
332*777c4855SStefano Zampini       if (j == 0) {
333*777c4855SStefano Zampini         xb  = user->bottom[i + 1];
334*777c4855SStefano Zampini         xrb = user->bottom[i + 2];
335*777c4855SStefano Zampini       } else xb = x[j - 1][i];
336*777c4855SStefano Zampini 
337*777c4855SStefano Zampini       /* Right */
338*777c4855SStefano Zampini       if (i + 1 == mx) {
339*777c4855SStefano Zampini         xr  = user->right[j + 1];
340*777c4855SStefano Zampini         xrb = user->right[j];
341*777c4855SStefano Zampini       } else xr = x[j][i + 1];
342*777c4855SStefano Zampini 
343*777c4855SStefano Zampini       /* Top */
344*777c4855SStefano Zampini       if (j + 1 == my) {
345*777c4855SStefano Zampini         xt  = user->top[i + 1];
346*777c4855SStefano Zampini         xlt = user->top[i];
347*777c4855SStefano Zampini       } else xt = x[j + 1][i];
348*777c4855SStefano Zampini 
349*777c4855SStefano Zampini       /* Top left */
350*777c4855SStefano Zampini       if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1];
351*777c4855SStefano Zampini 
352*777c4855SStefano Zampini       /* Bottom right */
353*777c4855SStefano Zampini       if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1];
354*777c4855SStefano Zampini 
355*777c4855SStefano Zampini       d1 = (xc - xl) / hx;
356*777c4855SStefano Zampini       d2 = (xc - xr) / hx;
357*777c4855SStefano Zampini       d3 = (xc - xt) / hy;
358*777c4855SStefano Zampini       d4 = (xc - xb) / hy;
359*777c4855SStefano Zampini       d5 = (xrb - xr) / hy;
360*777c4855SStefano Zampini       d6 = (xrb - xb) / hx;
361*777c4855SStefano Zampini       d7 = (xlt - xl) / hy;
362*777c4855SStefano Zampini       d8 = (xlt - xt) / hx;
363*777c4855SStefano Zampini 
364*777c4855SStefano Zampini       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
365*777c4855SStefano Zampini       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
366*777c4855SStefano Zampini       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
367*777c4855SStefano Zampini       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
368*777c4855SStefano Zampini       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
369*777c4855SStefano Zampini       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
370*777c4855SStefano Zampini 
371*777c4855SStefano Zampini       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
372*777c4855SStefano Zampini       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
373*777c4855SStefano Zampini       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
374*777c4855SStefano Zampini       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
375*777c4855SStefano Zampini 
376*777c4855SStefano Zampini       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
377*777c4855SStefano Zampini       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
378*777c4855SStefano Zampini 
379*777c4855SStefano Zampini       hc = hydhx * (1.0 + d7 * d7) / (f1 * f1 * f1) + hxdhy * (1.0 + d8 * d8) / (f3 * f3 * f3) + hydhx * (1.0 + d5 * d5) / (f5 * f5 * f5) + hxdhy * (1.0 + d6 * d6) / (f6 * f6 * f6) + (hxdhy * (1.0 + d1 * d1) + hydhx * (1.0 + d4 * d4) - 2.0 * d1 * d4) / (f2 * f2 * f2) + (hxdhy * (1.0 + d2 * d2) + hydhx * (1.0 + d3 * d3) - 2.0 * d2 * d3) / (f4 * f4 * f4);
380*777c4855SStefano Zampini 
381*777c4855SStefano Zampini       hl /= 2.0;
382*777c4855SStefano Zampini       hr /= 2.0;
383*777c4855SStefano Zampini       ht /= 2.0;
384*777c4855SStefano Zampini       hb /= 2.0;
385*777c4855SStefano Zampini       hbr /= 2.0;
386*777c4855SStefano Zampini       htl /= 2.0;
387*777c4855SStefano Zampini       hc /= 2.0;
388*777c4855SStefano Zampini 
389*777c4855SStefano Zampini       k     = 0;
390*777c4855SStefano Zampini       row.i = i;
391*777c4855SStefano Zampini       row.j = j;
392*777c4855SStefano Zampini       /* Bottom */
393*777c4855SStefano Zampini       if (j > 0) {
394*777c4855SStefano Zampini         v[k]     = hb;
395*777c4855SStefano Zampini         col[k].i = i;
396*777c4855SStefano Zampini         col[k].j = j - 1;
397*777c4855SStefano Zampini         k++;
398*777c4855SStefano Zampini       }
399*777c4855SStefano Zampini 
400*777c4855SStefano Zampini       /* Bottom right */
401*777c4855SStefano Zampini       if (j > 0 && i < mx - 1) {
402*777c4855SStefano Zampini         v[k]     = hbr;
403*777c4855SStefano Zampini         col[k].i = i + 1;
404*777c4855SStefano Zampini         col[k].j = j - 1;
405*777c4855SStefano Zampini         k++;
406*777c4855SStefano Zampini       }
407*777c4855SStefano Zampini 
408*777c4855SStefano Zampini       /* left */
409*777c4855SStefano Zampini       if (i > 0) {
410*777c4855SStefano Zampini         v[k]     = hl;
411*777c4855SStefano Zampini         col[k].i = i - 1;
412*777c4855SStefano Zampini         col[k].j = j;
413*777c4855SStefano Zampini         k++;
414*777c4855SStefano Zampini       }
415*777c4855SStefano Zampini 
416*777c4855SStefano Zampini       /* Centre */
417*777c4855SStefano Zampini       v[k]     = hc;
418*777c4855SStefano Zampini       col[k].i = row.i;
419*777c4855SStefano Zampini       col[k].j = row.j;
420*777c4855SStefano Zampini       k++;
421*777c4855SStefano Zampini 
422*777c4855SStefano Zampini       /* Right */
423*777c4855SStefano Zampini       if (i < mx - 1) {
424*777c4855SStefano Zampini         v[k]     = hr;
425*777c4855SStefano Zampini         col[k].i = i + 1;
426*777c4855SStefano Zampini         col[k].j = j;
427*777c4855SStefano Zampini         k++;
428*777c4855SStefano Zampini       }
429*777c4855SStefano Zampini 
430*777c4855SStefano Zampini       /* Top left */
431*777c4855SStefano Zampini       if (i > 0 && j < my - 1) {
432*777c4855SStefano Zampini         v[k]     = htl;
433*777c4855SStefano Zampini         col[k].i = i - 1;
434*777c4855SStefano Zampini         col[k].j = j + 1;
435*777c4855SStefano Zampini         k++;
436*777c4855SStefano Zampini       }
437*777c4855SStefano Zampini 
438*777c4855SStefano Zampini       /* Top */
439*777c4855SStefano Zampini       if (j < my - 1) {
440*777c4855SStefano Zampini         v[k]     = ht;
441*777c4855SStefano Zampini         col[k].i = i;
442*777c4855SStefano Zampini         col[k].j = j + 1;
443*777c4855SStefano Zampini         k++;
444*777c4855SStefano Zampini       }
445*777c4855SStefano Zampini 
446*777c4855SStefano Zampini       PetscCall(MatSetValuesStencil(Hp, 1, &row, k, col, v, INSERT_VALUES));
447*777c4855SStefano Zampini     }
448*777c4855SStefano Zampini   }
449*777c4855SStefano Zampini 
450*777c4855SStefano Zampini   /* Assemble the matrix */
451*777c4855SStefano Zampini   PetscCall(MatAssemblyBegin(Hp, MAT_FINAL_ASSEMBLY));
452*777c4855SStefano Zampini   PetscCall(MatAssemblyEnd(Hp, MAT_FINAL_ASSEMBLY));
453*777c4855SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
454*777c4855SStefano Zampini }
455*777c4855SStefano Zampini 
456*777c4855SStefano Zampini PetscErrorCode FormBoundaryConditions_Enneper(SNES snes, AppCtx **ouser)
457*777c4855SStefano Zampini {
458*777c4855SStefano Zampini   PetscInt     i, j, k, limit = 0, maxits = 5;
459*777c4855SStefano Zampini   PetscInt     mx, my;
460*777c4855SStefano Zampini   PetscInt     bsize = 0, lsize = 0, tsize = 0, rsize = 0;
461*777c4855SStefano Zampini   PetscScalar  one = 1.0, two = 2.0, three = 3.0;
462*777c4855SStefano Zampini   PetscScalar  det, hx, hy, xt = 0, yt = 0;
463*777c4855SStefano Zampini   PetscReal    fnorm, tol = 1e-10;
464*777c4855SStefano Zampini   PetscScalar  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
465*777c4855SStefano Zampini   PetscScalar  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
466*777c4855SStefano Zampini   PetscScalar *boundary;
467*777c4855SStefano Zampini   AppCtx      *user;
468*777c4855SStefano Zampini   DM           da;
469*777c4855SStefano Zampini 
470*777c4855SStefano Zampini   PetscFunctionBeginUser;
471*777c4855SStefano Zampini   PetscCall(SNESGetDM(snes, &da));
472*777c4855SStefano Zampini   PetscCall(PetscNew(&user));
473*777c4855SStefano Zampini   *ouser = user;
474*777c4855SStefano Zampini   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &mx, &my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
475*777c4855SStefano Zampini   bsize = mx + 2;
476*777c4855SStefano Zampini   lsize = my + 2;
477*777c4855SStefano Zampini   rsize = my + 2;
478*777c4855SStefano Zampini   tsize = mx + 2;
479*777c4855SStefano Zampini 
480*777c4855SStefano Zampini   PetscCall(PetscMalloc1(bsize, &user->bottom));
481*777c4855SStefano Zampini   PetscCall(PetscMalloc1(tsize, &user->top));
482*777c4855SStefano Zampini   PetscCall(PetscMalloc1(lsize, &user->left));
483*777c4855SStefano Zampini   PetscCall(PetscMalloc1(rsize, &user->right));
484*777c4855SStefano Zampini 
485*777c4855SStefano Zampini   hx = 1.0 / (mx + 1.0);
486*777c4855SStefano Zampini   hy = 1.0 / (my + 1.0);
487*777c4855SStefano Zampini 
488*777c4855SStefano Zampini   for (j = 0; j < 4; j++) {
489*777c4855SStefano Zampini     if (j == 0) {
490*777c4855SStefano Zampini       yt       = b;
491*777c4855SStefano Zampini       xt       = l;
492*777c4855SStefano Zampini       limit    = bsize;
493*777c4855SStefano Zampini       boundary = user->bottom;
494*777c4855SStefano Zampini     } else if (j == 1) {
495*777c4855SStefano Zampini       yt       = t;
496*777c4855SStefano Zampini       xt       = l;
497*777c4855SStefano Zampini       limit    = tsize;
498*777c4855SStefano Zampini       boundary = user->top;
499*777c4855SStefano Zampini     } else if (j == 2) {
500*777c4855SStefano Zampini       yt       = b;
501*777c4855SStefano Zampini       xt       = l;
502*777c4855SStefano Zampini       limit    = lsize;
503*777c4855SStefano Zampini       boundary = user->left;
504*777c4855SStefano Zampini     } else { /* if  (j==3) */
505*777c4855SStefano Zampini       yt       = b;
506*777c4855SStefano Zampini       xt       = r;
507*777c4855SStefano Zampini       limit    = rsize;
508*777c4855SStefano Zampini       boundary = user->right;
509*777c4855SStefano Zampini     }
510*777c4855SStefano Zampini 
511*777c4855SStefano Zampini     for (i = 0; i < limit; i++) {
512*777c4855SStefano Zampini       u1 = xt;
513*777c4855SStefano Zampini       u2 = -yt;
514*777c4855SStefano Zampini       for (k = 0; k < maxits; k++) {
515*777c4855SStefano Zampini         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
516*777c4855SStefano Zampini         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
517*777c4855SStefano Zampini         fnorm = PetscRealPart(PetscSqrtScalar(nf1 * nf1 + nf2 * nf2));
518*777c4855SStefano Zampini         if (fnorm <= tol) break;
519*777c4855SStefano Zampini         njac11 = one + u2 * u2 - u1 * u1;
520*777c4855SStefano Zampini         njac12 = two * u1 * u2;
521*777c4855SStefano Zampini         njac21 = -two * u1 * u2;
522*777c4855SStefano Zampini         njac22 = -one - u1 * u1 + u2 * u2;
523*777c4855SStefano Zampini         det    = njac11 * njac22 - njac21 * njac12;
524*777c4855SStefano Zampini         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
525*777c4855SStefano Zampini         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
526*777c4855SStefano Zampini       }
527*777c4855SStefano Zampini 
528*777c4855SStefano Zampini       boundary[i] = u1 * u1 - u2 * u2;
529*777c4855SStefano Zampini       if (j == 0 || j == 1) xt = xt + hx;
530*777c4855SStefano Zampini       else yt = yt + hy; /* if (j==2 || j==3) */
531*777c4855SStefano Zampini     }
532*777c4855SStefano Zampini   }
533*777c4855SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
534*777c4855SStefano Zampini }
535*777c4855SStefano Zampini 
536*777c4855SStefano Zampini PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser)
537*777c4855SStefano Zampini {
538*777c4855SStefano Zampini   AppCtx *user = *ouser;
539*777c4855SStefano Zampini 
540*777c4855SStefano Zampini   PetscFunctionBeginUser;
541*777c4855SStefano Zampini   PetscCall(PetscFree(user->bottom));
542*777c4855SStefano Zampini   PetscCall(PetscFree(user->top));
543*777c4855SStefano Zampini   PetscCall(PetscFree(user->left));
544*777c4855SStefano Zampini   PetscCall(PetscFree(user->right));
545*777c4855SStefano Zampini   PetscCall(PetscFree(*ouser));
546*777c4855SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
547*777c4855SStefano Zampini }
548*777c4855SStefano Zampini 
549*777c4855SStefano Zampini PetscErrorCode FormBoundaryConditions_Sins(SNES snes, AppCtx **ouser)
550*777c4855SStefano Zampini {
551*777c4855SStefano Zampini   PetscInt     i, j;
552*777c4855SStefano Zampini   PetscInt     mx, my;
553*777c4855SStefano Zampini   PetscInt     limit, bsize = 0, lsize = 0, tsize = 0, rsize = 0;
554*777c4855SStefano Zampini   PetscScalar  hx, hy, xt = 0, yt = 0;
555*777c4855SStefano Zampini   PetscScalar  b = 0.0, t = 1.0, l = 0.0, r = 1.0;
556*777c4855SStefano Zampini   PetscScalar *boundary;
557*777c4855SStefano Zampini   AppCtx      *user;
558*777c4855SStefano Zampini   DM           da;
559*777c4855SStefano Zampini   PetscReal    pi2 = 2 * PETSC_PI;
560*777c4855SStefano Zampini 
561*777c4855SStefano Zampini   PetscFunctionBeginUser;
562*777c4855SStefano Zampini   PetscCall(SNESGetDM(snes, &da));
563*777c4855SStefano Zampini   PetscCall(PetscNew(&user));
564*777c4855SStefano Zampini   *ouser = user;
565*777c4855SStefano Zampini   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &mx, &my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
566*777c4855SStefano Zampini   bsize = mx + 2;
567*777c4855SStefano Zampini   lsize = my + 2;
568*777c4855SStefano Zampini   rsize = my + 2;
569*777c4855SStefano Zampini   tsize = mx + 2;
570*777c4855SStefano Zampini 
571*777c4855SStefano Zampini   PetscCall(PetscMalloc1(bsize, &user->bottom));
572*777c4855SStefano Zampini   PetscCall(PetscMalloc1(tsize, &user->top));
573*777c4855SStefano Zampini   PetscCall(PetscMalloc1(lsize, &user->left));
574*777c4855SStefano Zampini   PetscCall(PetscMalloc1(rsize, &user->right));
575*777c4855SStefano Zampini 
576*777c4855SStefano Zampini   hx = 1.0 / (mx + 1.0);
577*777c4855SStefano Zampini   hy = 1.0 / (my + 1.0);
578*777c4855SStefano Zampini 
579*777c4855SStefano Zampini   for (j = 0; j < 4; j++) {
580*777c4855SStefano Zampini     if (j == 0) {
581*777c4855SStefano Zampini       yt       = b;
582*777c4855SStefano Zampini       xt       = l;
583*777c4855SStefano Zampini       limit    = bsize;
584*777c4855SStefano Zampini       boundary = user->bottom;
585*777c4855SStefano Zampini     } else if (j == 1) {
586*777c4855SStefano Zampini       yt       = t;
587*777c4855SStefano Zampini       xt       = l;
588*777c4855SStefano Zampini       limit    = tsize;
589*777c4855SStefano Zampini       boundary = user->top;
590*777c4855SStefano Zampini     } else if (j == 2) {
591*777c4855SStefano Zampini       yt       = b;
592*777c4855SStefano Zampini       xt       = l;
593*777c4855SStefano Zampini       limit    = lsize;
594*777c4855SStefano Zampini       boundary = user->left;
595*777c4855SStefano Zampini     } else { /* if  (j==3) */
596*777c4855SStefano Zampini       yt       = b;
597*777c4855SStefano Zampini       xt       = r;
598*777c4855SStefano Zampini       limit    = rsize;
599*777c4855SStefano Zampini       boundary = user->right;
600*777c4855SStefano Zampini     }
601*777c4855SStefano Zampini 
602*777c4855SStefano Zampini     for (i = 0; i < limit; i++) {
603*777c4855SStefano Zampini       if (j == 0) { /* bottom */
604*777c4855SStefano Zampini         boundary[i] = -0.5 * PetscSinReal(pi2 * xt);
605*777c4855SStefano Zampini       } else if (j == 1) { /* top */
606*777c4855SStefano Zampini         boundary[i] = 0.5 * PetscSinReal(pi2 * xt);
607*777c4855SStefano Zampini       } else if (j == 2) { /* left */
608*777c4855SStefano Zampini         boundary[i] = -0.5 * PetscSinReal(pi2 * yt);
609*777c4855SStefano Zampini       } else { /* right */
610*777c4855SStefano Zampini         boundary[i] = 0.5 * PetscSinReal(pi2 * yt);
611*777c4855SStefano Zampini       }
612*777c4855SStefano Zampini       if (j == 0 || j == 1) xt = xt + hx;
613*777c4855SStefano Zampini       else yt = yt + hy;
614*777c4855SStefano Zampini     }
615*777c4855SStefano Zampini   }
616*777c4855SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
617*777c4855SStefano Zampini }
618*777c4855SStefano Zampini 
619*777c4855SStefano Zampini /*TEST
620*777c4855SStefano Zampini 
621*777c4855SStefano Zampini   build:
622*777c4855SStefano Zampini     requires: !complex
623*777c4855SStefano Zampini 
624*777c4855SStefano Zampini   test:
625*777c4855SStefano Zampini     requires: !single
626*777c4855SStefano Zampini     filter: sed -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
627*777c4855SStefano Zampini     suffix: qn_nasm
628*777c4855SStefano Zampini     args: -snes_type qn -snes_npc_side {{left right}separate output} -npc_snes_type nasm -snes_converged_reason -da_local_subdomains 4 -use_objective {{0 1}separate output}
629*777c4855SStefano Zampini 
630*777c4855SStefano Zampini TEST*/
631