xref: /petsc/src/snes/tutorials/ex58.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown #include <petscsnes.h>
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown static const char help[] = "Parallel version of the minimum surface area problem in 2D using DMDA.\n\
7c4762a1bSJed Brown  It solves a system of nonlinear equations in mixed\n\
8c4762a1bSJed Brown complementarity form.This example is based on a\n\
9c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and\n\
10c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
11c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
12c4762a1bSJed Brown This application solves this problem using complimentarity -- We are actually\n\
13c4762a1bSJed Brown solving the system  (grad f)_i >= 0, if x_i == l_i \n\
14c4762a1bSJed Brown                     (grad f)_i = 0, if l_i < x_i < u_i \n\
15c4762a1bSJed Brown                     (grad f)_i <= 0, if x_i == u_i  \n\
16c4762a1bSJed Brown where f is the function to be minimized. \n\
17c4762a1bSJed Brown \n\
18c4762a1bSJed Brown The command line options are:\n\
19c4762a1bSJed Brown   -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
20c4762a1bSJed Brown   -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
21c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise\n\
22c4762a1bSJed Brown   -lb <value>, lower bound on the variables\n\
23c4762a1bSJed Brown   -ub <value>, upper bound on the variables\n\n";
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    User-defined application context - contains data needed by the
27c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
28c4762a1bSJed Brown    FormFunction().
29c4762a1bSJed Brown */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown      This is a new version of the ../tests/ex8.c code
33c4762a1bSJed Brown 
34c4762a1bSJed Brown      Run, for example, with the options ./ex58 -snes_vi_monitor -ksp_monitor -mg_levels_ksp_monitor -pc_type mg -pc_mg_levels 2 -pc_mg_galerkin pmat -ksp_type fgmres
35c4762a1bSJed Brown 
36c4762a1bSJed Brown      Or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
37c4762a1bSJed Brown          multigrid levels, it will be determined automatically based on the number of refinements done)
38c4762a1bSJed Brown 
39c4762a1bSJed Brown       ./ex58 -pc_type mg -ksp_monitor  -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3
40c4762a1bSJed Brown              -mg_levels_ksp_monitor -snes_vi_monitor -mg_levels_pc_type sor -pc_mg_type full
41c4762a1bSJed Brown 
42c4762a1bSJed Brown */
43c4762a1bSJed Brown 
44c4762a1bSJed Brown typedef struct {
45c4762a1bSJed Brown   PetscScalar *bottom, *top, *left, *right;
46c4762a1bSJed Brown   PetscScalar  lb, ub;
47c4762a1bSJed Brown } AppCtx;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /* -------- User-defined Routines --------- */
50c4762a1bSJed Brown 
51c4762a1bSJed Brown extern PetscErrorCode FormBoundaryConditions(SNES, AppCtx **);
52c4762a1bSJed Brown extern PetscErrorCode DestroyBoundaryConditions(AppCtx **);
53c4762a1bSJed Brown extern PetscErrorCode ComputeInitialGuess(SNES, Vec, void *);
54c4762a1bSJed Brown extern PetscErrorCode FormGradient(SNES, Vec, Vec, void *);
55c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
56c4762a1bSJed Brown extern PetscErrorCode FormBounds(SNES, Vec, Vec);
57c4762a1bSJed Brown 
58d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
59d71ae5a4SJacob Faibussowitsch {
60c4762a1bSJed Brown   Vec  x, r; /* solution and residual vectors */
61c4762a1bSJed Brown   SNES snes; /* nonlinear solver context */
62c4762a1bSJed Brown   Mat  J;    /* Jacobian matrix */
63c4762a1bSJed Brown   DM   da;
64c4762a1bSJed Brown 
65327415f7SBarry Smith   PetscFunctionBeginUser;
669566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Create distributed array to manage the 2d grid */
699566063dSJacob Faibussowitsch   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));
709566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
719566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Extract global vectors from DMDA; */
749566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
759566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da, MATAIJ));
789566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &J));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Create nonlinear solver context */
819566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
829566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /*  Set function evaluation and Jacobian evaluation  routines */
859566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormGradient, NULL));
869566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(SNESSetComputeApplicationContext(snes, (PetscErrorCode(*)(SNES, void **))FormBoundaryConditions, (PetscErrorCode(*)(void **))DestroyBoundaryConditions));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(SNESSetComputeInitialGuess(snes, ComputeInitialGuess, NULL));
91c4762a1bSJed Brown 
929566063dSJacob Faibussowitsch   PetscCall(SNESVISetComputeVariableBounds(snes, FormBounds));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Solve the application */
979566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* Free memory */
1009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
1029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1039566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Free user-created data structures */
1069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
109b122ec5aSJacob Faibussowitsch   return 0;
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown /* -------------------------------------------------------------------- */
113c4762a1bSJed Brown 
114c4762a1bSJed Brown /*  FormBounds - sets the upper and lower bounds
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     Input Parameters:
117c4762a1bSJed Brown .   snes  - the SNES context
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     Output Parameters:
120c4762a1bSJed Brown .   xl - lower bounds
121c4762a1bSJed Brown .   xu - upper bounds
122c4762a1bSJed Brown */
123d71ae5a4SJacob Faibussowitsch PetscErrorCode FormBounds(SNES snes, Vec xl, Vec xu)
124d71ae5a4SJacob Faibussowitsch {
125c4762a1bSJed Brown   AppCtx *ctx;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBeginUser;
1289566063dSJacob Faibussowitsch   PetscCall(SNESGetApplicationContext(snes, &ctx));
1299566063dSJacob Faibussowitsch   PetscCall(VecSet(xl, ctx->lb));
1309566063dSJacob Faibussowitsch   PetscCall(VecSet(xu, ctx->ub));
131*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /* -------------------------------------------------------------------- */
135c4762a1bSJed Brown 
136c4762a1bSJed Brown /*  FormGradient - Evaluates gradient of f.
137c4762a1bSJed Brown 
138c4762a1bSJed Brown     Input Parameters:
139c4762a1bSJed Brown .   snes  - the SNES context
140c4762a1bSJed Brown .   X     - input vector
141c4762a1bSJed Brown .   ptr   - optional user-defined context, as set by SNESSetFunction()
142c4762a1bSJed Brown 
143c4762a1bSJed Brown     Output Parameters:
144c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
145c4762a1bSJed Brown */
146d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *ptr)
147d71ae5a4SJacob Faibussowitsch {
148c4762a1bSJed Brown   AppCtx       *user;
149c4762a1bSJed Brown   PetscInt      i, j;
150c4762a1bSJed Brown   PetscInt      mx, my;
151c4762a1bSJed Brown   PetscScalar   hx, hy, hydhx, hxdhy;
152c4762a1bSJed Brown   PetscScalar   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
153c4762a1bSJed Brown   PetscScalar   df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
154c4762a1bSJed Brown   PetscScalar **g, **x;
155c4762a1bSJed Brown   PetscInt      xs, xm, ys, ym;
156c4762a1bSJed Brown   Vec           localX;
157c4762a1bSJed Brown   DM            da;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBeginUser;
1609566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
1619566063dSJacob Faibussowitsch   PetscCall(SNESGetApplicationContext(snes, &user));
1629566063dSJacob Faibussowitsch   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));
1639371c9d4SSatish Balay   hx    = 1.0 / (mx + 1);
1649371c9d4SSatish Balay   hy    = 1.0 / (my + 1);
1659371c9d4SSatish Balay   hydhx = hy / hx;
1669371c9d4SSatish Balay   hxdhy = hx / hy;
167c4762a1bSJed Brown 
1689566063dSJacob Faibussowitsch   PetscCall(VecSet(G, 0.0));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* Get local vector */
1719566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
172c4762a1bSJed Brown   /* Get ghost points */
1739566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1749566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
175c4762a1bSJed Brown   /* Get pointer to local vector data */
1769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
1779566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, G, &g));
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
180c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
181c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
182c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
183c4762a1bSJed Brown       xc  = x[j][i];
184c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown       if (i == 0) { /* left side */
187c4762a1bSJed Brown         xl  = user->left[j + 1];
188c4762a1bSJed Brown         xlt = user->left[j + 2];
189c4762a1bSJed Brown       } else xl = x[j][i - 1];
190c4762a1bSJed Brown 
191c4762a1bSJed Brown       if (j == 0) { /* bottom side */
192c4762a1bSJed Brown         xb  = user->bottom[i + 1];
193c4762a1bSJed Brown         xrb = user->bottom[i + 2];
194c4762a1bSJed Brown       } else xb = x[j - 1][i];
195c4762a1bSJed Brown 
196c4762a1bSJed Brown       if (i + 1 == mx) { /* right side */
197c4762a1bSJed Brown         xr  = user->right[j + 1];
198c4762a1bSJed Brown         xrb = user->right[j];
199c4762a1bSJed Brown       } else xr = x[j][i + 1];
200c4762a1bSJed Brown 
201c4762a1bSJed Brown       if (j + 1 == 0 + my) { /* top side */
202c4762a1bSJed Brown         xt  = user->top[i + 1];
203c4762a1bSJed Brown         xlt = user->top[i];
204c4762a1bSJed Brown       } else xt = x[j + 1][i];
205c4762a1bSJed Brown 
206c4762a1bSJed Brown       if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1]; /* left top side */
207c4762a1bSJed Brown       if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1]; /* right bottom */
208c4762a1bSJed Brown 
209c4762a1bSJed Brown       d1 = (xc - xl);
210c4762a1bSJed Brown       d2 = (xc - xr);
211c4762a1bSJed Brown       d3 = (xc - xt);
212c4762a1bSJed Brown       d4 = (xc - xb);
213c4762a1bSJed Brown       d5 = (xr - xrb);
214c4762a1bSJed Brown       d6 = (xrb - xb);
215c4762a1bSJed Brown       d7 = (xlt - xl);
216c4762a1bSJed Brown       d8 = (xt - xlt);
217c4762a1bSJed Brown 
218c4762a1bSJed Brown       df1dxc = d1 * hydhx;
219c4762a1bSJed Brown       df2dxc = (d1 * hydhx + d4 * hxdhy);
220c4762a1bSJed Brown       df3dxc = d3 * hxdhy;
221c4762a1bSJed Brown       df4dxc = (d2 * hydhx + d3 * hxdhy);
222c4762a1bSJed Brown       df5dxc = d2 * hydhx;
223c4762a1bSJed Brown       df6dxc = d4 * hxdhy;
224c4762a1bSJed Brown 
225c4762a1bSJed Brown       d1 /= hx;
226c4762a1bSJed Brown       d2 /= hx;
227c4762a1bSJed Brown       d3 /= hy;
228c4762a1bSJed Brown       d4 /= hy;
229c4762a1bSJed Brown       d5 /= hy;
230c4762a1bSJed Brown       d6 /= hx;
231c4762a1bSJed Brown       d7 /= hy;
232c4762a1bSJed Brown       d8 /= hx;
233c4762a1bSJed Brown 
234c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
235c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
236c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
237c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
238c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
239c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
240c4762a1bSJed Brown 
241c4762a1bSJed Brown       df1dxc /= f1;
242c4762a1bSJed Brown       df2dxc /= f2;
243c4762a1bSJed Brown       df3dxc /= f3;
244c4762a1bSJed Brown       df4dxc /= f4;
245c4762a1bSJed Brown       df5dxc /= f5;
246c4762a1bSJed Brown       df6dxc /= f6;
247c4762a1bSJed Brown 
248c4762a1bSJed Brown       g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0;
249c4762a1bSJed Brown     }
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* Restore vectors */
2539566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
2549566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, G, &g));
2559566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
2569566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(67.0 * mx * my));
257*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258c4762a1bSJed Brown }
259c4762a1bSJed Brown 
260c4762a1bSJed Brown /* ------------------------------------------------------------------- */
261c4762a1bSJed Brown /*
262c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
263c4762a1bSJed Brown 
264c4762a1bSJed Brown    Input Parameters:
265c4762a1bSJed Brown .  snes - SNES context
266c4762a1bSJed Brown .  X    - input vector
267c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by SNESSetJacobian()
268c4762a1bSJed Brown 
269c4762a1bSJed Brown    Output Parameters:
270c4762a1bSJed Brown .  tH    - Jacobian matrix
271c4762a1bSJed Brown 
272c4762a1bSJed Brown */
273d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *ptr)
274d71ae5a4SJacob Faibussowitsch {
275c4762a1bSJed Brown   AppCtx       *user;
276c4762a1bSJed Brown   PetscInt      i, j, k;
277c4762a1bSJed Brown   PetscInt      mx, my;
278c4762a1bSJed Brown   MatStencil    row, col[7];
279c4762a1bSJed Brown   PetscScalar   hx, hy, hydhx, hxdhy;
280c4762a1bSJed Brown   PetscScalar   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
281c4762a1bSJed Brown   PetscScalar   hl, hr, ht, hb, hc, htl, hbr;
282c4762a1bSJed Brown   PetscScalar **x, v[7];
283c4762a1bSJed Brown   PetscBool     assembled;
284c4762a1bSJed Brown   PetscInt      xs, xm, ys, ym;
285c4762a1bSJed Brown   Vec           localX;
286c4762a1bSJed Brown   DM            da;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   PetscFunctionBeginUser;
2899566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
2909566063dSJacob Faibussowitsch   PetscCall(SNESGetApplicationContext(snes, &user));
2919566063dSJacob Faibussowitsch   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));
2929371c9d4SSatish Balay   hx    = 1.0 / (mx + 1);
2939371c9d4SSatish Balay   hy    = 1.0 / (my + 1);
2949371c9d4SSatish Balay   hydhx = hy / hx;
2959371c9d4SSatish Balay   hxdhy = hx / hy;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /* Set various matrix options */
2989566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H, &assembled));
2999566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(H));
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   /* Get local vector */
3029566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
303c4762a1bSJed Brown   /* Get ghost points */
3049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
3059566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /* Get pointers to vector data */
3089566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
309c4762a1bSJed Brown 
3109566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
311c4762a1bSJed Brown   /* Compute Jacobian over the locally owned part of the mesh */
312c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
313c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
314c4762a1bSJed Brown       xc  = x[j][i];
315c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
316c4762a1bSJed Brown 
317c4762a1bSJed Brown       /* Left */
318c4762a1bSJed Brown       if (i == 0) {
319c4762a1bSJed Brown         xl  = user->left[j + 1];
320c4762a1bSJed Brown         xlt = user->left[j + 2];
321c4762a1bSJed Brown       } else xl = x[j][i - 1];
322c4762a1bSJed Brown 
323c4762a1bSJed Brown       /* Bottom */
324c4762a1bSJed Brown       if (j == 0) {
325c4762a1bSJed Brown         xb  = user->bottom[i + 1];
326c4762a1bSJed Brown         xrb = user->bottom[i + 2];
327c4762a1bSJed Brown       } else xb = x[j - 1][i];
328c4762a1bSJed Brown 
329c4762a1bSJed Brown       /* Right */
330c4762a1bSJed Brown       if (i + 1 == mx) {
331c4762a1bSJed Brown         xr  = user->right[j + 1];
332c4762a1bSJed Brown         xrb = user->right[j];
333c4762a1bSJed Brown       } else xr = x[j][i + 1];
334c4762a1bSJed Brown 
335c4762a1bSJed Brown       /* Top */
336c4762a1bSJed Brown       if (j + 1 == my) {
337c4762a1bSJed Brown         xt  = user->top[i + 1];
338c4762a1bSJed Brown         xlt = user->top[i];
339c4762a1bSJed Brown       } else xt = x[j + 1][i];
340c4762a1bSJed Brown 
341c4762a1bSJed Brown       /* Top left */
342c4762a1bSJed Brown       if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1];
343c4762a1bSJed Brown 
344c4762a1bSJed Brown       /* Bottom right */
345c4762a1bSJed Brown       if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1];
346c4762a1bSJed Brown 
347c4762a1bSJed Brown       d1 = (xc - xl) / hx;
348c4762a1bSJed Brown       d2 = (xc - xr) / hx;
349c4762a1bSJed Brown       d3 = (xc - xt) / hy;
350c4762a1bSJed Brown       d4 = (xc - xb) / hy;
351c4762a1bSJed Brown       d5 = (xrb - xr) / hy;
352c4762a1bSJed Brown       d6 = (xrb - xb) / hx;
353c4762a1bSJed Brown       d7 = (xlt - xl) / hy;
354c4762a1bSJed Brown       d8 = (xlt - xt) / hx;
355c4762a1bSJed Brown 
356c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
357c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
358c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
359c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
360c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
361c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
362c4762a1bSJed Brown 
3639371c9d4SSatish Balay       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
3649371c9d4SSatish Balay       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
3659371c9d4SSatish Balay       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
3669371c9d4SSatish Balay       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
367c4762a1bSJed Brown 
368c4762a1bSJed Brown       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
369c4762a1bSJed Brown       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
370c4762a1bSJed Brown 
3719371c9d4SSatish Balay       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);
372c4762a1bSJed Brown 
3739371c9d4SSatish Balay       hl /= 2.0;
3749371c9d4SSatish Balay       hr /= 2.0;
3759371c9d4SSatish Balay       ht /= 2.0;
3769371c9d4SSatish Balay       hb /= 2.0;
3779371c9d4SSatish Balay       hbr /= 2.0;
3789371c9d4SSatish Balay       htl /= 2.0;
3799371c9d4SSatish Balay       hc /= 2.0;
380c4762a1bSJed Brown 
381c4762a1bSJed Brown       k     = 0;
3829371c9d4SSatish Balay       row.i = i;
3839371c9d4SSatish Balay       row.j = j;
384c4762a1bSJed Brown       /* Bottom */
385c4762a1bSJed Brown       if (j > 0) {
386c4762a1bSJed Brown         v[k]     = hb;
3879371c9d4SSatish Balay         col[k].i = i;
3889371c9d4SSatish Balay         col[k].j = j - 1;
3899371c9d4SSatish Balay         k++;
390c4762a1bSJed Brown       }
391c4762a1bSJed Brown 
392c4762a1bSJed Brown       /* Bottom right */
393c4762a1bSJed Brown       if (j > 0 && i < mx - 1) {
394c4762a1bSJed Brown         v[k]     = hbr;
3959371c9d4SSatish Balay         col[k].i = i + 1;
3969371c9d4SSatish Balay         col[k].j = j - 1;
3979371c9d4SSatish Balay         k++;
398c4762a1bSJed Brown       }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown       /* left */
401c4762a1bSJed Brown       if (i > 0) {
402c4762a1bSJed Brown         v[k]     = hl;
4039371c9d4SSatish Balay         col[k].i = i - 1;
4049371c9d4SSatish Balay         col[k].j = j;
4059371c9d4SSatish Balay         k++;
406c4762a1bSJed Brown       }
407c4762a1bSJed Brown 
408c4762a1bSJed Brown       /* Centre */
4099371c9d4SSatish Balay       v[k]     = hc;
4109371c9d4SSatish Balay       col[k].i = row.i;
4119371c9d4SSatish Balay       col[k].j = row.j;
4129371c9d4SSatish Balay       k++;
413c4762a1bSJed Brown 
414c4762a1bSJed Brown       /* Right */
415c4762a1bSJed Brown       if (i < mx - 1) {
416c4762a1bSJed Brown         v[k]     = hr;
4179371c9d4SSatish Balay         col[k].i = i + 1;
4189371c9d4SSatish Balay         col[k].j = j;
4199371c9d4SSatish Balay         k++;
420c4762a1bSJed Brown       }
421c4762a1bSJed Brown 
422c4762a1bSJed Brown       /* Top left */
423c4762a1bSJed Brown       if (i > 0 && j < my - 1) {
424c4762a1bSJed Brown         v[k]     = htl;
4259371c9d4SSatish Balay         col[k].i = i - 1;
4269371c9d4SSatish Balay         col[k].j = j + 1;
4279371c9d4SSatish Balay         k++;
428c4762a1bSJed Brown       }
429c4762a1bSJed Brown 
430c4762a1bSJed Brown       /* Top */
431c4762a1bSJed Brown       if (j < my - 1) {
432c4762a1bSJed Brown         v[k]     = ht;
4339371c9d4SSatish Balay         col[k].i = i;
4349371c9d4SSatish Balay         col[k].j = j + 1;
4359371c9d4SSatish Balay         k++;
436c4762a1bSJed Brown       }
437c4762a1bSJed Brown 
4389566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(H, 1, &row, k, col, v, INSERT_VALUES));
439c4762a1bSJed Brown     }
440c4762a1bSJed Brown   }
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   /* Assemble the matrix */
4439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
4449566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
4459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
4469566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
447c4762a1bSJed Brown 
4489566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0 * mx * my));
449*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
450c4762a1bSJed Brown }
451c4762a1bSJed Brown 
452c4762a1bSJed Brown /* ------------------------------------------------------------------- */
453c4762a1bSJed Brown /*
454c4762a1bSJed Brown    FormBoundaryConditions -  Calculates the boundary conditions for
455c4762a1bSJed Brown    the region.
456c4762a1bSJed Brown 
457c4762a1bSJed Brown    Input Parameter:
458c4762a1bSJed Brown .  user - user-defined application context
459c4762a1bSJed Brown 
460c4762a1bSJed Brown    Output Parameter:
461c4762a1bSJed Brown .  user - user-defined application context
462c4762a1bSJed Brown */
463d71ae5a4SJacob Faibussowitsch PetscErrorCode FormBoundaryConditions(SNES snes, AppCtx **ouser)
464d71ae5a4SJacob Faibussowitsch {
465c4762a1bSJed Brown   PetscInt     i, j, k, limit = 0, maxits = 5;
466c4762a1bSJed Brown   PetscInt     mx, my;
467c4762a1bSJed Brown   PetscInt     bsize = 0, lsize = 0, tsize = 0, rsize = 0;
468c4762a1bSJed Brown   PetscScalar  one = 1.0, two = 2.0, three = 3.0;
469c4762a1bSJed Brown   PetscScalar  det, hx, hy, xt = 0, yt = 0;
470c4762a1bSJed Brown   PetscReal    fnorm, tol = 1e-10;
471c4762a1bSJed Brown   PetscScalar  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
472c4762a1bSJed Brown   PetscScalar  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
473c4762a1bSJed Brown   PetscScalar *boundary;
474c4762a1bSJed Brown   AppCtx      *user;
475c4762a1bSJed Brown   DM           da;
476c4762a1bSJed Brown 
477c4762a1bSJed Brown   PetscFunctionBeginUser;
4789566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
4799566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
480c4762a1bSJed Brown   *ouser   = user;
481c4762a1bSJed Brown   user->lb = .05;
482c4762a1bSJed Brown   user->ub = PETSC_INFINITY;
4839566063dSJacob Faibussowitsch   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));
484c4762a1bSJed Brown 
485c4762a1bSJed Brown   /* Check if lower and upper bounds are set */
4869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lb", &user->lb, 0));
4879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ub", &user->ub, 0));
4889371c9d4SSatish Balay   bsize = mx + 2;
4899371c9d4SSatish Balay   lsize = my + 2;
4909371c9d4SSatish Balay   rsize = my + 2;
4919371c9d4SSatish Balay   tsize = mx + 2;
492c4762a1bSJed Brown 
4939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bsize, &user->bottom));
4949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsize, &user->top));
4959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lsize, &user->left));
4969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(rsize, &user->right));
497c4762a1bSJed Brown 
4989371c9d4SSatish Balay   hx = (r - l) / (mx + 1.0);
4999371c9d4SSatish Balay   hy = (t - b) / (my + 1.0);
500c4762a1bSJed Brown 
501c4762a1bSJed Brown   for (j = 0; j < 4; j++) {
502c4762a1bSJed Brown     if (j == 0) {
503c4762a1bSJed Brown       yt       = b;
504c4762a1bSJed Brown       xt       = l;
505c4762a1bSJed Brown       limit    = bsize;
506c4762a1bSJed Brown       boundary = user->bottom;
507c4762a1bSJed Brown     } else if (j == 1) {
508c4762a1bSJed Brown       yt       = t;
509c4762a1bSJed Brown       xt       = l;
510c4762a1bSJed Brown       limit    = tsize;
511c4762a1bSJed Brown       boundary = user->top;
512c4762a1bSJed Brown     } else if (j == 2) {
513c4762a1bSJed Brown       yt       = b;
514c4762a1bSJed Brown       xt       = l;
515c4762a1bSJed Brown       limit    = lsize;
516c4762a1bSJed Brown       boundary = user->left;
517c4762a1bSJed Brown     } else { /* if  (j==3) */
518c4762a1bSJed Brown       yt       = b;
519c4762a1bSJed Brown       xt       = r;
520c4762a1bSJed Brown       limit    = rsize;
521c4762a1bSJed Brown       boundary = user->right;
522c4762a1bSJed Brown     }
523c4762a1bSJed Brown 
524c4762a1bSJed Brown     for (i = 0; i < limit; i++) {
525c4762a1bSJed Brown       u1 = xt;
526c4762a1bSJed Brown       u2 = -yt;
527c4762a1bSJed Brown       for (k = 0; k < maxits; k++) {
528c4762a1bSJed Brown         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
529c4762a1bSJed Brown         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
530c4762a1bSJed Brown         fnorm = PetscRealPart(PetscSqrtScalar(nf1 * nf1 + nf2 * nf2));
531c4762a1bSJed Brown         if (fnorm <= tol) break;
532c4762a1bSJed Brown         njac11 = one + u2 * u2 - u1 * u1;
533c4762a1bSJed Brown         njac12 = two * u1 * u2;
534c4762a1bSJed Brown         njac21 = -two * u1 * u2;
535c4762a1bSJed Brown         njac22 = -one - u1 * u1 + u2 * u2;
536c4762a1bSJed Brown         det    = njac11 * njac22 - njac21 * njac12;
537c4762a1bSJed Brown         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
538c4762a1bSJed Brown         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
539c4762a1bSJed Brown       }
540c4762a1bSJed Brown 
541c4762a1bSJed Brown       boundary[i] = u1 * u1 - u2 * u2;
542c4762a1bSJed Brown       if (j == 0 || j == 1) xt = xt + hx;
543c4762a1bSJed Brown       else yt = yt + hy; /* if (j==2 || j==3) */
544c4762a1bSJed Brown     }
545c4762a1bSJed Brown   }
546*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
547c4762a1bSJed Brown }
548c4762a1bSJed Brown 
549d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser)
550d71ae5a4SJacob Faibussowitsch {
551c4762a1bSJed Brown   AppCtx *user = *ouser;
552c4762a1bSJed Brown 
553c4762a1bSJed Brown   PetscFunctionBeginUser;
5549566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->bottom));
5559566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->top));
5569566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->left));
5579566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->right));
5589566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ouser));
559*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
560c4762a1bSJed Brown }
561c4762a1bSJed Brown 
562c4762a1bSJed Brown /* ------------------------------------------------------------------- */
563c4762a1bSJed Brown /*
564c4762a1bSJed Brown    ComputeInitialGuess - Calculates the initial guess
565c4762a1bSJed Brown 
566c4762a1bSJed Brown    Input Parameters:
567c4762a1bSJed Brown .  user - user-defined application context
568c4762a1bSJed Brown .  X - vector for initial guess
569c4762a1bSJed Brown 
570c4762a1bSJed Brown    Output Parameters:
571c4762a1bSJed Brown .  X - newly computed initial guess
572c4762a1bSJed Brown */
573d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeInitialGuess(SNES snes, Vec X, void *dummy)
574d71ae5a4SJacob Faibussowitsch {
575c4762a1bSJed Brown   PetscInt      i, j, mx, my;
576c4762a1bSJed Brown   DM            da;
577c4762a1bSJed Brown   AppCtx       *user;
578c4762a1bSJed Brown   PetscScalar **x;
579c4762a1bSJed Brown   PetscInt      xs, xm, ys, ym;
580c4762a1bSJed Brown 
581c4762a1bSJed Brown   PetscFunctionBeginUser;
5829566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
5839566063dSJacob Faibussowitsch   PetscCall(SNESGetApplicationContext(snes, &user));
584c4762a1bSJed Brown 
5859566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
5869566063dSJacob Faibussowitsch   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));
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   /* Get pointers to vector data */
5899566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
590c4762a1bSJed Brown   /* Perform local computations */
591c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
592ad540459SPierre Jolivet     for (i = xs; i < xs + xm; i++) x[j][i] = (((j + 1.0) * user->bottom[i + 1] + (my - j + 1.0) * user->top[i + 1]) / (my + 2.0) + ((i + 1.0) * user->left[j + 1] + (mx - i + 1.0) * user->right[j + 1]) / (mx + 2.0)) / 2.0;
593c4762a1bSJed Brown   }
594c4762a1bSJed Brown   /* Restore vectors */
5959566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
596*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
597c4762a1bSJed Brown }
598c4762a1bSJed Brown 
599c4762a1bSJed Brown /*TEST
600c4762a1bSJed Brown 
601c4762a1bSJed Brown    test:
602c4762a1bSJed Brown       args: -snes_type vinewtonrsls -pc_type mg -ksp_monitor_short -pc_mg_galerkin pmat -da_refine 5 -snes_vi_monitor -pc_mg_type full -snes_max_it 100 -snes_converged_reason
603c4762a1bSJed Brown       requires: !single
604c4762a1bSJed Brown 
605c4762a1bSJed Brown    test:
606c4762a1bSJed Brown       suffix: 2
607c4762a1bSJed Brown       args: -snes_type vinewtonssls -pc_type mg -ksp_monitor_short -pc_mg_galerkin pmat -da_refine 5 -snes_vi_monitor -pc_mg_type full -snes_max_it 100 -snes_converged_reason
608c4762a1bSJed Brown       requires: !single
609c4762a1bSJed Brown 
610c4762a1bSJed Brown TEST*/
611