xref: /petsc/src/snes/tutorials/ex58.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown #include <petscsnes.h>
2c4762a1bSJed Brown #include <petscdm.h>
3c4762a1bSJed Brown #include <petscdmda.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown static const char help[] = "Parallel version of the minimum surface area problem in 2D using DMDA.\n\
6c4762a1bSJed Brown  It solves a system of nonlinear equations in mixed\n\
7c4762a1bSJed Brown complementarity form.This example is based on a\n\
8c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and\n\
9c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
10c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
11d5def619SJonas Heinzmann This application solves this problem using complementarity -- We are actually\n\
12c4762a1bSJed Brown solving the system  (grad f)_i >= 0, if x_i == l_i \n\
13c4762a1bSJed Brown                     (grad f)_i = 0, if l_i < x_i < u_i \n\
14c4762a1bSJed Brown                     (grad f)_i <= 0, if x_i == u_i  \n\
15c4762a1bSJed Brown where f is the function to be minimized. \n\
16c4762a1bSJed Brown \n\
17c4762a1bSJed Brown The command line options are:\n\
18c4762a1bSJed Brown   -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
19c4762a1bSJed Brown   -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
20c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise\n\
21c4762a1bSJed Brown   -lb <value>, lower bound on the variables\n\
22c4762a1bSJed Brown   -ub <value>, upper bound on the variables\n\n";
23c4762a1bSJed Brown 
24c4762a1bSJed Brown /*
25c4762a1bSJed Brown    User-defined application context - contains data needed by the
26c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
27c4762a1bSJed Brown    FormFunction().
28c4762a1bSJed Brown */
29c4762a1bSJed Brown 
30c4762a1bSJed Brown /*
31c4762a1bSJed 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
32c4762a1bSJed Brown 
33c4762a1bSJed Brown      Or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
34c4762a1bSJed Brown          multigrid levels, it will be determined automatically based on the number of refinements done)
35c4762a1bSJed Brown 
36c4762a1bSJed Brown       ./ex58 -pc_type mg -ksp_monitor  -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3
37c4762a1bSJed Brown              -mg_levels_ksp_monitor -snes_vi_monitor -mg_levels_pc_type sor -pc_mg_type full
38c4762a1bSJed Brown 
39c4762a1bSJed Brown */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown typedef struct {
42c4762a1bSJed Brown   PetscScalar *bottom, *top, *left, *right;
43c4762a1bSJed Brown   PetscScalar  lb, ub;
44c4762a1bSJed Brown } AppCtx;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /* -------- User-defined Routines --------- */
47c4762a1bSJed Brown 
48*2a8381b2SBarry Smith extern PetscErrorCode FormBoundaryConditions(SNES, PetscCtxRt);
49*2a8381b2SBarry Smith extern PetscErrorCode DestroyBoundaryConditions(PetscCtxRt);
50c4762a1bSJed Brown extern PetscErrorCode ComputeInitialGuess(SNES, Vec, void *);
51c4762a1bSJed Brown extern PetscErrorCode FormGradient(SNES, Vec, Vec, void *);
52c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
53c4762a1bSJed Brown extern PetscErrorCode FormBounds(SNES, Vec, Vec);
54c4762a1bSJed Brown 
main(int argc,char ** argv)55d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
56d71ae5a4SJacob Faibussowitsch {
57c4762a1bSJed Brown   Vec  x, r; /* solution and residual vectors */
58c4762a1bSJed Brown   SNES snes; /* nonlinear solver context */
59c4762a1bSJed Brown   Mat  J;    /* Jacobian matrix */
60c4762a1bSJed Brown   DM   da;
61c4762a1bSJed Brown 
62327415f7SBarry Smith   PetscFunctionBeginUser;
63c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Create distributed array to manage the 2d grid */
669566063dSJacob 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));
679566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
689566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* Extract global vectors from DMDA; */
719566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da, MATAIJ));
759566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &J));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Create nonlinear solver context */
789566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
799566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /*  Set function evaluation and Jacobian evaluation  routines */
829566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormGradient, NULL));
839566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, NULL));
84c4762a1bSJed Brown 
85*2a8381b2SBarry Smith   PetscCall(SNESSetComputeApplicationContext(snes, FormBoundaryConditions, DestroyBoundaryConditions));
86c4762a1bSJed Brown 
879566063dSJacob Faibussowitsch   PetscCall(SNESSetComputeInitialGuess(snes, ComputeInitialGuess, NULL));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(SNESVISetComputeVariableBounds(snes, FormBounds));
90c4762a1bSJed Brown 
919566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* Solve the application */
949566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Free memory */
979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1009566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* Free user-created data structures */
1039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
106b122ec5aSJacob Faibussowitsch   return 0;
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown /* -------------------------------------------------------------------- */
110c4762a1bSJed Brown 
111c4762a1bSJed Brown /*  FormBounds - sets the upper and lower bounds
112c4762a1bSJed Brown 
113c4762a1bSJed Brown     Input Parameters:
114c4762a1bSJed Brown .   snes  - the SNES context
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     Output Parameters:
117c4762a1bSJed Brown .   xl - lower bounds
118c4762a1bSJed Brown .   xu - upper bounds
119c4762a1bSJed Brown */
FormBounds(SNES snes,Vec xl,Vec xu)120d71ae5a4SJacob Faibussowitsch PetscErrorCode FormBounds(SNES snes, Vec xl, Vec xu)
121d71ae5a4SJacob Faibussowitsch {
122c4762a1bSJed Brown   AppCtx *ctx;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   PetscFunctionBeginUser;
1259566063dSJacob Faibussowitsch   PetscCall(SNESGetApplicationContext(snes, &ctx));
1269566063dSJacob Faibussowitsch   PetscCall(VecSet(xl, ctx->lb));
1279566063dSJacob Faibussowitsch   PetscCall(VecSet(xu, ctx->ub));
1283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown /* -------------------------------------------------------------------- */
132c4762a1bSJed Brown 
133c4762a1bSJed Brown /*  FormGradient - Evaluates gradient of f.
134c4762a1bSJed Brown 
135c4762a1bSJed Brown     Input Parameters:
136c4762a1bSJed Brown .   snes   - the SNES context
137c4762a1bSJed Brown .   X      - input vector
138*2a8381b2SBarry Smith .   unused - optional user-defined context, as set by SNESSetFunction()
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     Output Parameters:
141c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
142c4762a1bSJed Brown */
FormGradient(SNES snes,Vec X,Vec G,void * unused)143*2a8381b2SBarry Smith PetscErrorCode FormGradient(SNES snes, Vec X, Vec G, void *unused)
144d71ae5a4SJacob Faibussowitsch {
145*2a8381b2SBarry Smith   AppCtx       *ctx;
146c4762a1bSJed Brown   PetscInt      i, j;
147c4762a1bSJed Brown   PetscInt      mx, my;
148c4762a1bSJed Brown   PetscScalar   hx, hy, hydhx, hxdhy;
149c4762a1bSJed Brown   PetscScalar   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
150c4762a1bSJed Brown   PetscScalar   df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
151c4762a1bSJed Brown   PetscScalar **g, **x;
152c4762a1bSJed Brown   PetscInt      xs, xm, ys, ym;
153c4762a1bSJed Brown   Vec           localX;
154c4762a1bSJed Brown   DM            da;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBeginUser;
1579566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
158*2a8381b2SBarry Smith   PetscCall(SNESGetApplicationContext(snes, &ctx));
1599566063dSJacob 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));
1609371c9d4SSatish Balay   hx    = 1.0 / (mx + 1);
1619371c9d4SSatish Balay   hy    = 1.0 / (my + 1);
1629371c9d4SSatish Balay   hydhx = hy / hx;
1639371c9d4SSatish Balay   hxdhy = hx / hy;
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(VecSet(G, 0.0));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Get local vector */
1689566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
169c4762a1bSJed Brown   /* Get ghost points */
1709566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1719566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
172c4762a1bSJed Brown   /* Get pointer to local vector data */
1739566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
1749566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, G, &g));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
177c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
178c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
179c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
180c4762a1bSJed Brown       xc  = x[j][i];
181c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
182c4762a1bSJed Brown 
183c4762a1bSJed Brown       if (i == 0) { /* left side */
184*2a8381b2SBarry Smith         xl  = ctx->left[j + 1];
185*2a8381b2SBarry Smith         xlt = ctx->left[j + 2];
186c4762a1bSJed Brown       } else xl = x[j][i - 1];
187c4762a1bSJed Brown 
188c4762a1bSJed Brown       if (j == 0) { /* bottom side */
189*2a8381b2SBarry Smith         xb  = ctx->bottom[i + 1];
190*2a8381b2SBarry Smith         xrb = ctx->bottom[i + 2];
191c4762a1bSJed Brown       } else xb = x[j - 1][i];
192c4762a1bSJed Brown 
193c4762a1bSJed Brown       if (i + 1 == mx) { /* right side */
194*2a8381b2SBarry Smith         xr  = ctx->right[j + 1];
195*2a8381b2SBarry Smith         xrb = ctx->right[j];
196c4762a1bSJed Brown       } else xr = x[j][i + 1];
197c4762a1bSJed Brown 
198c4762a1bSJed Brown       if (j + 1 == 0 + my) { /* top side */
199*2a8381b2SBarry Smith         xt  = ctx->top[i + 1];
200*2a8381b2SBarry Smith         xlt = ctx->top[i];
201c4762a1bSJed Brown       } else xt = x[j + 1][i];
202c4762a1bSJed Brown 
203c4762a1bSJed Brown       if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1]; /* left top side */
204c4762a1bSJed Brown       if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1]; /* right bottom */
205c4762a1bSJed Brown 
206c4762a1bSJed Brown       d1 = (xc - xl);
207c4762a1bSJed Brown       d2 = (xc - xr);
208c4762a1bSJed Brown       d3 = (xc - xt);
209c4762a1bSJed Brown       d4 = (xc - xb);
210c4762a1bSJed Brown       d5 = (xr - xrb);
211c4762a1bSJed Brown       d6 = (xrb - xb);
212c4762a1bSJed Brown       d7 = (xlt - xl);
213c4762a1bSJed Brown       d8 = (xt - xlt);
214c4762a1bSJed Brown 
215c4762a1bSJed Brown       df1dxc = d1 * hydhx;
216c4762a1bSJed Brown       df2dxc = (d1 * hydhx + d4 * hxdhy);
217c4762a1bSJed Brown       df3dxc = d3 * hxdhy;
218c4762a1bSJed Brown       df4dxc = (d2 * hydhx + d3 * hxdhy);
219c4762a1bSJed Brown       df5dxc = d2 * hydhx;
220c4762a1bSJed Brown       df6dxc = d4 * hxdhy;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown       d1 /= hx;
223c4762a1bSJed Brown       d2 /= hx;
224c4762a1bSJed Brown       d3 /= hy;
225c4762a1bSJed Brown       d4 /= hy;
226c4762a1bSJed Brown       d5 /= hy;
227c4762a1bSJed Brown       d6 /= hx;
228c4762a1bSJed Brown       d7 /= hy;
229c4762a1bSJed Brown       d8 /= hx;
230c4762a1bSJed Brown 
231c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
232c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
233c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
234c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
235c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
236c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
237c4762a1bSJed Brown 
238c4762a1bSJed Brown       df1dxc /= f1;
239c4762a1bSJed Brown       df2dxc /= f2;
240c4762a1bSJed Brown       df3dxc /= f3;
241c4762a1bSJed Brown       df4dxc /= f4;
242c4762a1bSJed Brown       df5dxc /= f5;
243c4762a1bSJed Brown       df6dxc /= f6;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown       g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0;
246c4762a1bSJed Brown     }
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Restore vectors */
2509566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
2519566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, G, &g));
2529566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
2539566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(67.0 * mx * my));
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown /* ------------------------------------------------------------------- */
258c4762a1bSJed Brown /*
259c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
260c4762a1bSJed Brown 
261c4762a1bSJed Brown    Input Parameters:
262c4762a1bSJed Brown .  snes   - SNES context
263c4762a1bSJed Brown .  X      - input vector
264*2a8381b2SBarry Smith .  unused - optional user-defined context, as set by SNESSetJacobian()
265c4762a1bSJed Brown 
266c4762a1bSJed Brown    Output Parameters:
267c4762a1bSJed Brown .  tH    - Jacobian matrix
268c4762a1bSJed Brown 
269c4762a1bSJed Brown */
FormJacobian(SNES snes,Vec X,Mat H,Mat tHPre,void * unused)270*2a8381b2SBarry Smith PetscErrorCode FormJacobian(SNES snes, Vec X, Mat H, Mat tHPre, void *unused)
271d71ae5a4SJacob Faibussowitsch {
272*2a8381b2SBarry Smith   AppCtx       *ctx;
273c4762a1bSJed Brown   PetscInt      i, j, k;
274c4762a1bSJed Brown   PetscInt      mx, my;
275c4762a1bSJed Brown   MatStencil    row, col[7];
276c4762a1bSJed Brown   PetscScalar   hx, hy, hydhx, hxdhy;
277c4762a1bSJed Brown   PetscScalar   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
278c4762a1bSJed Brown   PetscScalar   hl, hr, ht, hb, hc, htl, hbr;
279c4762a1bSJed Brown   PetscScalar **x, v[7];
280c4762a1bSJed Brown   PetscBool     assembled;
281c4762a1bSJed Brown   PetscInt      xs, xm, ys, ym;
282c4762a1bSJed Brown   Vec           localX;
283c4762a1bSJed Brown   DM            da;
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   PetscFunctionBeginUser;
2869566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
287*2a8381b2SBarry Smith   PetscCall(SNESGetApplicationContext(snes, &ctx));
2889566063dSJacob 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));
2899371c9d4SSatish Balay   hx    = 1.0 / (mx + 1);
2909371c9d4SSatish Balay   hy    = 1.0 / (my + 1);
2919371c9d4SSatish Balay   hydhx = hy / hx;
2929371c9d4SSatish Balay   hxdhy = hx / hy;
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   /* Set various matrix options */
2959566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H, &assembled));
2969566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(H));
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   /* Get local vector */
2999566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
300c4762a1bSJed Brown   /* Get ghost points */
3019566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
3029566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /* Get pointers to vector data */
3059566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
306c4762a1bSJed Brown 
3079566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
308c4762a1bSJed Brown   /* Compute Jacobian over the locally owned part of the mesh */
309c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
310c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
311c4762a1bSJed Brown       xc  = x[j][i];
312c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
313c4762a1bSJed Brown 
314c4762a1bSJed Brown       /* Left */
315c4762a1bSJed Brown       if (i == 0) {
316*2a8381b2SBarry Smith         xl  = ctx->left[j + 1];
317*2a8381b2SBarry Smith         xlt = ctx->left[j + 2];
318c4762a1bSJed Brown       } else xl = x[j][i - 1];
319c4762a1bSJed Brown 
320c4762a1bSJed Brown       /* Bottom */
321c4762a1bSJed Brown       if (j == 0) {
322*2a8381b2SBarry Smith         xb  = ctx->bottom[i + 1];
323*2a8381b2SBarry Smith         xrb = ctx->bottom[i + 2];
324c4762a1bSJed Brown       } else xb = x[j - 1][i];
325c4762a1bSJed Brown 
326c4762a1bSJed Brown       /* Right */
327c4762a1bSJed Brown       if (i + 1 == mx) {
328*2a8381b2SBarry Smith         xr  = ctx->right[j + 1];
329*2a8381b2SBarry Smith         xrb = ctx->right[j];
330c4762a1bSJed Brown       } else xr = x[j][i + 1];
331c4762a1bSJed Brown 
332c4762a1bSJed Brown       /* Top */
333c4762a1bSJed Brown       if (j + 1 == my) {
334*2a8381b2SBarry Smith         xt  = ctx->top[i + 1];
335*2a8381b2SBarry Smith         xlt = ctx->top[i];
336c4762a1bSJed Brown       } else xt = x[j + 1][i];
337c4762a1bSJed Brown 
338c4762a1bSJed Brown       /* Top left */
339c4762a1bSJed Brown       if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1];
340c4762a1bSJed Brown 
341c4762a1bSJed Brown       /* Bottom right */
342c4762a1bSJed Brown       if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1];
343c4762a1bSJed Brown 
344c4762a1bSJed Brown       d1 = (xc - xl) / hx;
345c4762a1bSJed Brown       d2 = (xc - xr) / hx;
346c4762a1bSJed Brown       d3 = (xc - xt) / hy;
347c4762a1bSJed Brown       d4 = (xc - xb) / hy;
348c4762a1bSJed Brown       d5 = (xrb - xr) / hy;
349c4762a1bSJed Brown       d6 = (xrb - xb) / hx;
350c4762a1bSJed Brown       d7 = (xlt - xl) / hy;
351c4762a1bSJed Brown       d8 = (xlt - xt) / hx;
352c4762a1bSJed Brown 
353c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
354c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
355c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
356c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
357c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
358c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
359c4762a1bSJed Brown 
3609371c9d4SSatish Balay       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
3619371c9d4SSatish Balay       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
3629371c9d4SSatish Balay       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
3639371c9d4SSatish Balay       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
364c4762a1bSJed Brown 
365c4762a1bSJed Brown       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
366c4762a1bSJed Brown       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
367c4762a1bSJed Brown 
3689371c9d4SSatish 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);
369c4762a1bSJed Brown 
3709371c9d4SSatish Balay       hl /= 2.0;
3719371c9d4SSatish Balay       hr /= 2.0;
3729371c9d4SSatish Balay       ht /= 2.0;
3739371c9d4SSatish Balay       hb /= 2.0;
3749371c9d4SSatish Balay       hbr /= 2.0;
3759371c9d4SSatish Balay       htl /= 2.0;
3769371c9d4SSatish Balay       hc /= 2.0;
377c4762a1bSJed Brown 
378c4762a1bSJed Brown       k     = 0;
3799371c9d4SSatish Balay       row.i = i;
3809371c9d4SSatish Balay       row.j = j;
381c4762a1bSJed Brown       /* Bottom */
382c4762a1bSJed Brown       if (j > 0) {
383c4762a1bSJed Brown         v[k]     = hb;
3849371c9d4SSatish Balay         col[k].i = i;
3859371c9d4SSatish Balay         col[k].j = j - 1;
3869371c9d4SSatish Balay         k++;
387c4762a1bSJed Brown       }
388c4762a1bSJed Brown 
389c4762a1bSJed Brown       /* Bottom right */
390c4762a1bSJed Brown       if (j > 0 && i < mx - 1) {
391c4762a1bSJed Brown         v[k]     = hbr;
3929371c9d4SSatish Balay         col[k].i = i + 1;
3939371c9d4SSatish Balay         col[k].j = j - 1;
3949371c9d4SSatish Balay         k++;
395c4762a1bSJed Brown       }
396c4762a1bSJed Brown 
397c4762a1bSJed Brown       /* left */
398c4762a1bSJed Brown       if (i > 0) {
399c4762a1bSJed Brown         v[k]     = hl;
4009371c9d4SSatish Balay         col[k].i = i - 1;
4019371c9d4SSatish Balay         col[k].j = j;
4029371c9d4SSatish Balay         k++;
403c4762a1bSJed Brown       }
404c4762a1bSJed Brown 
405c4762a1bSJed Brown       /* Centre */
4069371c9d4SSatish Balay       v[k]     = hc;
4079371c9d4SSatish Balay       col[k].i = row.i;
4089371c9d4SSatish Balay       col[k].j = row.j;
4099371c9d4SSatish Balay       k++;
410c4762a1bSJed Brown 
411c4762a1bSJed Brown       /* Right */
412c4762a1bSJed Brown       if (i < mx - 1) {
413c4762a1bSJed Brown         v[k]     = hr;
4149371c9d4SSatish Balay         col[k].i = i + 1;
4159371c9d4SSatish Balay         col[k].j = j;
4169371c9d4SSatish Balay         k++;
417c4762a1bSJed Brown       }
418c4762a1bSJed Brown 
419c4762a1bSJed Brown       /* Top left */
420c4762a1bSJed Brown       if (i > 0 && j < my - 1) {
421c4762a1bSJed Brown         v[k]     = htl;
4229371c9d4SSatish Balay         col[k].i = i - 1;
4239371c9d4SSatish Balay         col[k].j = j + 1;
4249371c9d4SSatish Balay         k++;
425c4762a1bSJed Brown       }
426c4762a1bSJed Brown 
427c4762a1bSJed Brown       /* Top */
428c4762a1bSJed Brown       if (j < my - 1) {
429c4762a1bSJed Brown         v[k]     = ht;
4309371c9d4SSatish Balay         col[k].i = i;
4319371c9d4SSatish Balay         col[k].j = j + 1;
4329371c9d4SSatish Balay         k++;
433c4762a1bSJed Brown       }
434c4762a1bSJed Brown 
4359566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(H, 1, &row, k, col, v, INSERT_VALUES));
436c4762a1bSJed Brown     }
437c4762a1bSJed Brown   }
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /* Assemble the matrix */
4409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
4419566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
4429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
4439566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
444c4762a1bSJed Brown 
4459566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0 * mx * my));
4463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
447c4762a1bSJed Brown }
448c4762a1bSJed Brown 
449c4762a1bSJed Brown /* ------------------------------------------------------------------- */
450c4762a1bSJed Brown /*
451c4762a1bSJed Brown    FormBoundaryConditions -  Calculates the boundary conditions for
452c4762a1bSJed Brown    the region.
453c4762a1bSJed Brown 
454c4762a1bSJed Brown    Input Parameter:
455c4762a1bSJed Brown .  user - user-defined application context
456c4762a1bSJed Brown 
457c4762a1bSJed Brown    Output Parameter:
458c4762a1bSJed Brown .  user - user-defined application context
459c4762a1bSJed Brown */
FormBoundaryConditions(SNES snes,PetscCtxRt Ctx)460*2a8381b2SBarry Smith PetscErrorCode FormBoundaryConditions(SNES snes, PetscCtxRt Ctx)
461d71ae5a4SJacob Faibussowitsch {
462c4762a1bSJed Brown   PetscInt     i, j, k, limit = 0, maxits = 5;
463c4762a1bSJed Brown   PetscInt     mx, my;
464c4762a1bSJed Brown   PetscInt     bsize = 0, lsize = 0, tsize = 0, rsize = 0;
465c4762a1bSJed Brown   PetscScalar  one = 1.0, two = 2.0, three = 3.0;
466c4762a1bSJed Brown   PetscScalar  det, hx, hy, xt = 0, yt = 0;
467c4762a1bSJed Brown   PetscReal    fnorm, tol = 1e-10;
468c4762a1bSJed Brown   PetscScalar  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
469c4762a1bSJed Brown   PetscScalar  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
470c4762a1bSJed Brown   PetscScalar *boundary;
471*2a8381b2SBarry Smith   AppCtx      *ctx;
472c4762a1bSJed Brown   DM           da;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   PetscFunctionBeginUser;
4759566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
476*2a8381b2SBarry Smith   PetscCall(PetscNew(&ctx));
477*2a8381b2SBarry Smith   *(void **)Ctx = ctx;
478*2a8381b2SBarry Smith   ctx->lb       = .05;
479*2a8381b2SBarry Smith   ctx->ub       = PETSC_INFINITY;
4809566063dSJacob 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));
481c4762a1bSJed Brown 
482c4762a1bSJed Brown   /* Check if lower and upper bounds are set */
483*2a8381b2SBarry Smith   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lb", &ctx->lb, 0));
484*2a8381b2SBarry Smith   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ub", &ctx->ub, 0));
4859371c9d4SSatish Balay   bsize = mx + 2;
4869371c9d4SSatish Balay   lsize = my + 2;
4879371c9d4SSatish Balay   rsize = my + 2;
4889371c9d4SSatish Balay   tsize = mx + 2;
489c4762a1bSJed Brown 
490*2a8381b2SBarry Smith   PetscCall(PetscMalloc1(bsize, &ctx->bottom));
491*2a8381b2SBarry Smith   PetscCall(PetscMalloc1(tsize, &ctx->top));
492*2a8381b2SBarry Smith   PetscCall(PetscMalloc1(lsize, &ctx->left));
493*2a8381b2SBarry Smith   PetscCall(PetscMalloc1(rsize, &ctx->right));
494c4762a1bSJed Brown 
4959371c9d4SSatish Balay   hx = (r - l) / (mx + 1.0);
4969371c9d4SSatish Balay   hy = (t - b) / (my + 1.0);
497c4762a1bSJed Brown 
498c4762a1bSJed Brown   for (j = 0; j < 4; j++) {
499c4762a1bSJed Brown     if (j == 0) {
500c4762a1bSJed Brown       yt       = b;
501c4762a1bSJed Brown       xt       = l;
502c4762a1bSJed Brown       limit    = bsize;
503*2a8381b2SBarry Smith       boundary = ctx->bottom;
504c4762a1bSJed Brown     } else if (j == 1) {
505c4762a1bSJed Brown       yt       = t;
506c4762a1bSJed Brown       xt       = l;
507c4762a1bSJed Brown       limit    = tsize;
508*2a8381b2SBarry Smith       boundary = ctx->top;
509c4762a1bSJed Brown     } else if (j == 2) {
510c4762a1bSJed Brown       yt       = b;
511c4762a1bSJed Brown       xt       = l;
512c4762a1bSJed Brown       limit    = lsize;
513*2a8381b2SBarry Smith       boundary = ctx->left;
514c4762a1bSJed Brown     } else { /* if  (j==3) */
515c4762a1bSJed Brown       yt       = b;
516c4762a1bSJed Brown       xt       = r;
517c4762a1bSJed Brown       limit    = rsize;
518*2a8381b2SBarry Smith       boundary = ctx->right;
519c4762a1bSJed Brown     }
520c4762a1bSJed Brown 
521c4762a1bSJed Brown     for (i = 0; i < limit; i++) {
522c4762a1bSJed Brown       u1 = xt;
523c4762a1bSJed Brown       u2 = -yt;
524c4762a1bSJed Brown       for (k = 0; k < maxits; k++) {
525c4762a1bSJed Brown         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
526c4762a1bSJed Brown         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
527c4762a1bSJed Brown         fnorm = PetscRealPart(PetscSqrtScalar(nf1 * nf1 + nf2 * nf2));
528c4762a1bSJed Brown         if (fnorm <= tol) break;
529c4762a1bSJed Brown         njac11 = one + u2 * u2 - u1 * u1;
530c4762a1bSJed Brown         njac12 = two * u1 * u2;
531c4762a1bSJed Brown         njac21 = -two * u1 * u2;
532c4762a1bSJed Brown         njac22 = -one - u1 * u1 + u2 * u2;
533c4762a1bSJed Brown         det    = njac11 * njac22 - njac21 * njac12;
534c4762a1bSJed Brown         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
535c4762a1bSJed Brown         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
536c4762a1bSJed Brown       }
537c4762a1bSJed Brown 
538c4762a1bSJed Brown       boundary[i] = u1 * u1 - u2 * u2;
539c4762a1bSJed Brown       if (j == 0 || j == 1) xt = xt + hx;
540c4762a1bSJed Brown       else yt = yt + hy; /* if (j==2 || j==3) */
541c4762a1bSJed Brown     }
542c4762a1bSJed Brown   }
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
544c4762a1bSJed Brown }
545c4762a1bSJed Brown 
DestroyBoundaryConditions(PetscCtxRt Ctx)546*2a8381b2SBarry Smith PetscErrorCode DestroyBoundaryConditions(PetscCtxRt Ctx)
547d71ae5a4SJacob Faibussowitsch {
548*2a8381b2SBarry Smith   AppCtx *ctx = *(AppCtx **)Ctx;
549c4762a1bSJed Brown 
550c4762a1bSJed Brown   PetscFunctionBeginUser;
551*2a8381b2SBarry Smith   PetscCall(PetscFree(ctx->bottom));
552*2a8381b2SBarry Smith   PetscCall(PetscFree(ctx->top));
553*2a8381b2SBarry Smith   PetscCall(PetscFree(ctx->left));
554*2a8381b2SBarry Smith   PetscCall(PetscFree(ctx->right));
555*2a8381b2SBarry Smith   PetscCall(PetscFree(ctx));
5563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557c4762a1bSJed Brown }
558c4762a1bSJed Brown 
559c4762a1bSJed Brown /* ------------------------------------------------------------------- */
560c4762a1bSJed Brown /*
561c4762a1bSJed Brown    ComputeInitialGuess - Calculates the initial guess
562c4762a1bSJed Brown 
563c4762a1bSJed Brown    Input Parameters:
564*2a8381b2SBarry Smith .  unused - user-defined application context
565c4762a1bSJed Brown .  X      - vector for initial guess
566c4762a1bSJed Brown 
567c4762a1bSJed Brown    Output Parameters:
568c4762a1bSJed Brown .  X - newly computed initial guess
569c4762a1bSJed Brown */
ComputeInitialGuess(SNES snes,Vec X,void * unused)570*2a8381b2SBarry Smith PetscErrorCode ComputeInitialGuess(SNES snes, Vec X, void *unused)
571d71ae5a4SJacob Faibussowitsch {
572c4762a1bSJed Brown   PetscInt      i, j, mx, my;
573c4762a1bSJed Brown   DM            da;
574*2a8381b2SBarry Smith   AppCtx       *ctx;
575c4762a1bSJed Brown   PetscScalar **x;
576c4762a1bSJed Brown   PetscInt      xs, xm, ys, ym;
577c4762a1bSJed Brown 
578c4762a1bSJed Brown   PetscFunctionBeginUser;
5799566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
580*2a8381b2SBarry Smith   PetscCall(SNESGetApplicationContext(snes, &ctx));
581c4762a1bSJed Brown 
5829566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
5839566063dSJacob 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));
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   /* Get pointers to vector data */
5869566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
587c4762a1bSJed Brown   /* Perform local computations */
588c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
589*2a8381b2SBarry Smith     for (i = xs; i < xs + xm; i++) x[j][i] = (((j + 1.0) * ctx->bottom[i + 1] + (my - j + 1.0) * ctx->top[i + 1]) / (my + 2.0) + ((i + 1.0) * ctx->left[j + 1] + (mx - i + 1.0) * ctx->right[j + 1]) / (mx + 2.0)) / 2.0;
590c4762a1bSJed Brown   }
591c4762a1bSJed Brown   /* Restore vectors */
5929566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
5933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
594c4762a1bSJed Brown }
595c4762a1bSJed Brown 
596c4762a1bSJed Brown /*TEST
597c4762a1bSJed Brown 
598c4762a1bSJed Brown    test:
599c4762a1bSJed 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
600c4762a1bSJed Brown       requires: !single
601c4762a1bSJed Brown 
602c4762a1bSJed Brown    test:
603c4762a1bSJed Brown       suffix: 2
604c4762a1bSJed 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
605c4762a1bSJed Brown       requires: !single
606c4762a1bSJed Brown 
607d5def619SJonas Heinzmann    test:
608d5def619SJonas Heinzmann       suffix: 3
609d5def619SJonas Heinzmann       args: -snes_type vinewtonrsls -snes_linesearch_type bisection -snes_linesearch_monitor -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
610d5def619SJonas Heinzmann       requires: !single
611d5def619SJonas Heinzmann 
612d5def619SJonas Heinzmann    test:
613d5def619SJonas Heinzmann       suffix: 4
614d5def619SJonas Heinzmann       args: -snes_type vinewtonssls -snes_linesearch_type bisection -snes_linesearch_monitor -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
615d5def619SJonas Heinzmann       requires: !single
616d5def619SJonas Heinzmann 
617d5def619SJonas Heinzmann    test:
618d5def619SJonas Heinzmann       suffix: 5
619a99ef635SJonas Heinzmann       args: -snes_type vinewtonssls -snes_linesearch_type cp -snes_linesearch_order 3 -snes_linesearch_monitor -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 -snes_linesearch_maxlambda 2.0
620a99ef635SJonas Heinzmann       requires: !single
621a99ef635SJonas Heinzmann 
622a99ef635SJonas Heinzmann    test:
623a99ef635SJonas Heinzmann       suffix: 6
624a99ef635SJonas Heinzmann       args: -snes_type vinewtonssls -snes_linesearch_type secant -snes_linesearch_max_it 10 -snes_linesearch_ltol 1e-2 -snes_linesearch_atol 1e-32 -snes_linesearch_monitor -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 -snes_linesearch_maxlambda 1.0
625d5def619SJonas Heinzmann       requires: !single
626d5def619SJonas Heinzmann 
627c4762a1bSJed Brown TEST*/
628