xref: /petsc/src/tao/unconstrained/tutorials/minsurf2.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 /* Program usage: mpiexec -n <proc> minsurf2 [-help] [all TAO options] */
2 
3 /*
4   Include "petsctao.h" so we can use TAO solvers.
5   petscdmda.h for distributed array
6 */
7 #include <petsctao.h>
8 #include <petscdmda.h>
9 
10 static  char help[] =
11 "This example demonstrates use of the TAO package to \n\
12 solve an unconstrained minimization problem.  This example is based on a \n\
13 problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and \n\
14 boundary values along the edges of the domain, the objective is to find the\n\
15 surface with the minimal area that satisfies the boundary conditions.\n\
16 The command line options are:\n\
17   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
18   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
19   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
20                for an average of the boundary conditions\n\n";
21 
22 /*
23    User-defined application context - contains data needed by the
24    application-provided call-back routines, FormFunctionGradient()
25    and FormHessian().
26 */
27 typedef struct {
28   PetscInt    mx, my;                 /* discretization in x, y directions */
29   PetscReal   *bottom, *top, *left, *right;             /* boundary values */
30   DM          dm;                      /* distributed array data structure */
31   Mat         H;                       /* Hessian */
32 } AppCtx;
33 
34 /* -------- User-defined Routines --------- */
35 
36 static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
37 static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
38 PetscErrorCode QuadraticH(AppCtx*,Vec,Mat);
39 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void*);
40 PetscErrorCode FormGradient(Tao,Vec,Vec,void*);
41 PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
42 PetscErrorCode My_Monitor(Tao, void *);
43 
44 int main(int argc, char **argv)
45 {
46   PetscInt           Nx, Ny;              /* number of processors in x- and y- directions */
47   Vec                x;                   /* solution, gradient vectors */
48   PetscBool          flg, viewmat;        /* flags */
49   PetscBool          fddefault, fdcoloring;   /* flags */
50   Tao                tao;                 /* TAO solver context */
51   AppCtx             user;                /* user-defined work context */
52   ISColoring         iscoloring;
53   MatFDColoring      matfdcoloring;
54 
55   /* Initialize TAO */
56   PetscFunctionBeginUser;
57   PetscCall(PetscInitialize(&argc, &argv,(char *)0,help));
58 
59   /* Specify dimension of the problem */
60   user.mx = 10; user.my = 10;
61 
62   /* Check for any command line arguments that override defaults */
63   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
64   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
65 
66   PetscCall(PetscPrintf(MPI_COMM_WORLD,"\n---- Minimum Surface Area Problem -----\n"));
67   PetscCall(PetscPrintf(MPI_COMM_WORLD,"mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n",user.mx,user.my));
68 
69   /* Let PETSc determine the vector distribution */
70   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
71 
72   /* Create distributed array (DM) to manage parallel grid and vectors  */
73   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,user.mx, user.my,Nx,Ny,1,1,NULL,NULL,&user.dm));
74   PetscCall(DMSetFromOptions(user.dm));
75   PetscCall(DMSetUp(user.dm));
76 
77   /* Create TAO solver and set desired solution method.*/
78   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
79   PetscCall(TaoSetType(tao,TAOCG));
80 
81   /*
82      Extract global vector from DA for the vector of variables --  PETSC routine
83      Compute the initial solution                              --  application specific, see below
84      Set this vector for use by TAO                            --  TAO routine
85   */
86   PetscCall(DMCreateGlobalVector(user.dm,&x));
87   PetscCall(MSA_BoundaryConditions(&user));
88   PetscCall(MSA_InitialPoint(&user,x));
89   PetscCall(TaoSetSolution(tao,x));
90 
91   /*
92      Initialize the Application context for use in function evaluations  --  application specific, see below.
93      Set routines for function and gradient evaluation
94   */
95   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
96 
97   /*
98      Given the command line arguments, calculate the hessian with either the user-
99      provided function FormHessian, or the default finite-difference driven Hessian
100      functions
101   */
102   PetscCall(PetscOptionsHasName(NULL,NULL,"-tao_fddefault",&fddefault));
103   PetscCall(PetscOptionsHasName(NULL,NULL,"-tao_fdcoloring",&fdcoloring));
104 
105   /*
106      Create a matrix data structure to store the Hessian and set
107      the Hessian evalution routine.
108      Set the matrix structure to be used for Hessian evalutions
109   */
110   PetscCall(DMCreateMatrix(user.dm,&user.H));
111   PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
112 
113   if (fdcoloring) {
114     PetscCall(DMCreateColoring(user.dm,IS_COLORING_GLOBAL,&iscoloring));
115     PetscCall(MatFDColoringCreate(user.H,iscoloring,&matfdcoloring));
116     PetscCall(ISColoringDestroy(&iscoloring));
117     PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))FormGradient,(void*)&user));
118     PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
119     PetscCall(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessianColor,(void *)matfdcoloring));
120   } else if (fddefault) {
121     PetscCall(TaoSetHessian(tao,user.H,user.H,TaoDefaultComputeHessian,(void *)NULL));
122   } else {
123     PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user));
124   }
125 
126   /*
127      If my_monitor option is in command line, then use the user-provided
128      monitoring function
129   */
130   PetscCall(PetscOptionsHasName(NULL,NULL,"-my_monitor",&viewmat));
131   if (viewmat) PetscCall(TaoSetMonitor(tao,My_Monitor,NULL,NULL));
132 
133   /* Check for any tao command line options */
134   PetscCall(TaoSetFromOptions(tao));
135 
136   /* SOLVE THE APPLICATION */
137   PetscCall(TaoSolve(tao));
138 
139   PetscCall(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD));
140 
141   /* Free TAO data structures */
142   PetscCall(TaoDestroy(&tao));
143 
144   /* Free PETSc data structures */
145   PetscCall(VecDestroy(&x));
146   PetscCall(MatDestroy(&user.H));
147   if (fdcoloring) {
148     PetscCall(MatFDColoringDestroy(&matfdcoloring));
149   }
150   PetscCall(PetscFree(user.bottom));
151   PetscCall(PetscFree(user.top));
152   PetscCall(PetscFree(user.left));
153   PetscCall(PetscFree(user.right));
154   PetscCall(DMDestroy(&user.dm));
155   PetscCall(PetscFinalize());
156   return 0;
157 }
158 
159 PetscErrorCode FormGradient(Tao tao, Vec X, Vec G,void *userCtx)
160 {
161   PetscReal      fcn;
162 
163   PetscFunctionBegin;
164   PetscCall(FormFunctionGradient(tao,X,&fcn,G,userCtx));
165   PetscFunctionReturn(0);
166 }
167 
168 /* -------------------------------------------------------------------- */
169 /*  FormFunctionGradient - Evaluates the function and corresponding gradient.
170 
171     Input Parameters:
172 .   tao     - the Tao context
173 .   XX      - input vector
174 .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
175 
176     Output Parameters:
177 .   fcn     - the newly evaluated function
178 .   GG       - vector containing the newly evaluated gradient
179 */
180 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *userCtx)
181 {
182   AppCtx         *user = (AppCtx *) userCtx;
183   PetscInt       i,j;
184   PetscInt       mx=user->mx, my=user->my;
185   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
186   PetscReal      ft=0.0;
187   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
188   PetscReal      rhx=mx+1, rhy=my+1;
189   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
190   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
191   PetscReal      **g, **x;
192   Vec            localX;
193 
194   PetscFunctionBegin;
195   /* Get local mesh boundaries */
196   PetscCall(DMGetLocalVector(user->dm,&localX));
197   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
198   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
199 
200   /* Scatter ghost points to local vector */
201   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
202   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
203 
204   /* Get pointers to vector data */
205   PetscCall(DMDAVecGetArray(user->dm,localX,(void**)&x));
206   PetscCall(DMDAVecGetArray(user->dm,G,(void**)&g));
207 
208   /* Compute function and gradient over the locally owned part of the mesh */
209   for (j=ys; j<ys+ym; j++) {
210     for (i=xs; i< xs+xm; i++) {
211 
212       xc = x[j][i];
213       xlt=xrb=xl=xr=xb=xt=xc;
214 
215       if (i==0) { /* left side */
216         xl= user->left[j-ys+1];
217         xlt = user->left[j-ys+2];
218       } else {
219         xl = x[j][i-1];
220       }
221 
222       if (j==0) { /* bottom side */
223         xb=user->bottom[i-xs+1];
224         xrb =user->bottom[i-xs+2];
225       } else {
226         xb = x[j-1][i];
227       }
228 
229       if (i+1 == gxs+gxm) { /* right side */
230         xr=user->right[j-ys+1];
231         xrb = user->right[j-ys];
232       } else {
233         xr = x[j][i+1];
234       }
235 
236       if (j+1==gys+gym) { /* top side */
237         xt=user->top[i-xs+1];
238         xlt = user->top[i-xs];
239       }else {
240         xt = x[j+1][i];
241       }
242 
243       if (i>gxs && j+1<gys+gym) {
244         xlt = x[j+1][i-1];
245       }
246       if (j>gys && i+1<gxs+gxm) {
247         xrb = x[j-1][i+1];
248       }
249 
250       d1 = (xc-xl);
251       d2 = (xc-xr);
252       d3 = (xc-xt);
253       d4 = (xc-xb);
254       d5 = (xr-xrb);
255       d6 = (xrb-xb);
256       d7 = (xlt-xl);
257       d8 = (xt-xlt);
258 
259       df1dxc = d1*hydhx;
260       df2dxc = (d1*hydhx + d4*hxdhy);
261       df3dxc = d3*hxdhy;
262       df4dxc = (d2*hydhx + d3*hxdhy);
263       df5dxc = d2*hydhx;
264       df6dxc = d4*hxdhy;
265 
266       d1 *= rhx;
267       d2 *= rhx;
268       d3 *= rhy;
269       d4 *= rhy;
270       d5 *= rhy;
271       d6 *= rhx;
272       d7 *= rhy;
273       d8 *= rhx;
274 
275       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
276       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
277       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
278       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
279       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
280       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
281 
282       ft = ft + (f2 + f4);
283 
284       df1dxc /= f1;
285       df2dxc /= f2;
286       df3dxc /= f3;
287       df4dxc /= f4;
288       df5dxc /= f5;
289       df6dxc /= f6;
290 
291       g[j][i] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5;
292 
293     }
294   }
295 
296   /* Compute triangular areas along the border of the domain. */
297   if (xs==0) { /* left side */
298     for (j=ys; j<ys+ym; j++) {
299       d3=(user->left[j-ys+1] - user->left[j-ys+2])*rhy;
300       d2=(user->left[j-ys+1] - x[j][0]) *rhx;
301       ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2);
302     }
303   }
304   if (ys==0) { /* bottom side */
305     for (i=xs; i<xs+xm; i++) {
306       d2=(user->bottom[i+1-xs]-user->bottom[i-xs+2])*rhx;
307       d3=(user->bottom[i-xs+1]-x[0][i])*rhy;
308       ft = ft+PetscSqrtReal(1.0 + d3*d3 + d2*d2);
309     }
310   }
311 
312   if (xs+xm==mx) { /* right side */
313     for (j=ys; j< ys+ym; j++) {
314       d1=(x[j][mx-1] - user->right[j-ys+1])*rhx;
315       d4=(user->right[j-ys]-user->right[j-ys+1])*rhy;
316       ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4);
317     }
318   }
319   if (ys+ym==my) { /* top side */
320     for (i=xs; i<xs+xm; i++) {
321       d1=(x[my-1][i] - user->top[i-xs+1])*rhy;
322       d4=(user->top[i-xs+1] - user->top[i-xs])*rhx;
323       ft = ft+PetscSqrtReal(1.0 + d1*d1 + d4*d4);
324     }
325   }
326 
327   if (ys==0 && xs==0) {
328     d1=(user->left[0]-user->left[1])/hy;
329     d2=(user->bottom[0]-user->bottom[1])*rhx;
330     ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2);
331   }
332   if (ys+ym == my && xs+xm == mx) {
333     d1=(user->right[ym+1] - user->right[ym])*rhy;
334     d2=(user->top[xm+1] - user->top[xm])*rhx;
335     ft +=PetscSqrtReal(1.0 + d1*d1 + d2*d2);
336   }
337 
338   ft=ft*area;
339   PetscCallMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD));
340 
341   /* Restore vectors */
342   PetscCall(DMDAVecRestoreArray(user->dm,localX,(void**)&x));
343   PetscCall(DMDAVecRestoreArray(user->dm,G,(void**)&g));
344 
345   /* Scatter values to global vector */
346   PetscCall(DMRestoreLocalVector(user->dm,&localX));
347   PetscCall(PetscLogFlops(67.0*xm*ym));
348   PetscFunctionReturn(0);
349 }
350 
351 /* ------------------------------------------------------------------- */
352 /*
353    FormHessian - Evaluates Hessian matrix.
354 
355    Input Parameters:
356 .  tao  - the Tao context
357 .  x    - input vector
358 .  ptr  - optional user-defined context, as set by TaoSetHessian()
359 
360    Output Parameters:
361 .  H    - Hessian matrix
362 .  Hpre - optionally different preconditioning matrix
363 .  flg  - flag indicating matrix structure
364 
365 */
366 PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
367 {
368   AppCtx         *user = (AppCtx *) ptr;
369 
370   PetscFunctionBegin;
371   /* Evaluate the Hessian entries*/
372   PetscCall(QuadraticH(user,X,H));
373   PetscFunctionReturn(0);
374 }
375 
376 /* ------------------------------------------------------------------- */
377 /*
378    QuadraticH - Evaluates Hessian matrix.
379 
380    Input Parameters:
381 .  user - user-defined context, as set by TaoSetHessian()
382 .  X    - input vector
383 
384    Output Parameter:
385 .  H    - Hessian matrix
386 */
387 PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian)
388 {
389   PetscInt i,j,k;
390   PetscInt mx=user->mx, my=user->my;
391   PetscInt xs,xm,gxs,gxm,ys,ym,gys,gym;
392   PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
393   PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
394   PetscReal hl,hr,ht,hb,hc,htl,hbr;
395   PetscReal **x, v[7];
396   MatStencil col[7],row;
397   Vec    localX;
398   PetscBool assembled;
399 
400   PetscFunctionBegin;
401   /* Get local mesh boundaries */
402   PetscCall(DMGetLocalVector(user->dm,&localX));
403 
404   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
405   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
406 
407   /* Scatter ghost points to local vector */
408   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
409   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
410 
411   /* Get pointers to vector data */
412   PetscCall(DMDAVecGetArray(user->dm,localX,(void**)&x));
413 
414   /* Initialize matrix entries to zero */
415   PetscCall(MatAssembled(Hessian,&assembled));
416   if (assembled) PetscCall(MatZeroEntries(Hessian));
417 
418   /* Set various matrix options */
419   PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
420 
421   /* Compute Hessian over the locally owned part of the mesh */
422 
423   for (j=ys; j<ys+ym; j++) {
424 
425     for (i=xs; i< xs+xm; i++) {
426 
427       xc = x[j][i];
428       xlt=xrb=xl=xr=xb=xt=xc;
429 
430       /* Left side */
431       if (i==0) {
432         xl  = user->left[j-ys+1];
433         xlt = user->left[j-ys+2];
434       } else {
435         xl  = x[j][i-1];
436       }
437 
438       if (j==0) {
439         xb  = user->bottom[i-xs+1];
440         xrb = user->bottom[i-xs+2];
441       } else {
442         xb  = x[j-1][i];
443       }
444 
445       if (i+1 == mx) {
446         xr  = user->right[j-ys+1];
447         xrb = user->right[j-ys];
448       } else {
449         xr  = x[j][i+1];
450       }
451 
452       if (j+1==my) {
453         xt  = user->top[i-xs+1];
454         xlt = user->top[i-xs];
455       }else {
456         xt  = x[j+1][i];
457       }
458 
459       if (i>0 && j+1<my) {
460         xlt = x[j+1][i-1];
461       }
462       if (j>0 && i+1<mx) {
463         xrb = x[j-1][i+1];
464       }
465 
466       d1 = (xc-xl)/hx;
467       d2 = (xc-xr)/hx;
468       d3 = (xc-xt)/hy;
469       d4 = (xc-xb)/hy;
470       d5 = (xrb-xr)/hy;
471       d6 = (xrb-xb)/hx;
472       d7 = (xlt-xl)/hy;
473       d8 = (xlt-xt)/hx;
474 
475       f1 = PetscSqrtReal(1.0 + d1*d1 + d7*d7);
476       f2 = PetscSqrtReal(1.0 + d1*d1 + d4*d4);
477       f3 = PetscSqrtReal(1.0 + d3*d3 + d8*d8);
478       f4 = PetscSqrtReal(1.0 + d3*d3 + d2*d2);
479       f5 = PetscSqrtReal(1.0 + d2*d2 + d5*d5);
480       f6 = PetscSqrtReal(1.0 + d4*d4 + d6*d6);
481 
482       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
483         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
484       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
485         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
486       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
487         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
488       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
489         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
490 
491       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
492       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
493 
494       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
495         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
496         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
497         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
498 
499       hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0;  hc/=2.0;
500 
501       row.j = j; row.i = i;
502       k=0;
503       if (j>0) {
504         v[k]=hb;
505         col[k].j = j - 1; col[k].i = i;
506         k++;
507       }
508 
509       if (j>0 && i < mx -1) {
510         v[k]=hbr;
511         col[k].j = j - 1; col[k].i = i+1;
512         k++;
513       }
514 
515       if (i>0) {
516         v[k]= hl;
517         col[k].j = j; col[k].i = i-1;
518         k++;
519       }
520 
521       v[k]= hc;
522       col[k].j = j; col[k].i = i;
523       k++;
524 
525       if (i < mx-1) {
526         v[k]= hr;
527         col[k].j = j; col[k].i = i+1;
528         k++;
529       }
530 
531       if (i>0 && j < my-1) {
532         v[k]= htl;
533         col[k].j = j+1; col[k].i = i-1;
534         k++;
535       }
536 
537       if (j < my-1) {
538         v[k]= ht;
539         col[k].j = j+1; col[k].i = i;
540         k++;
541       }
542 
543       /*
544          Set matrix values using local numbering, which was defined
545          earlier, in the main routine.
546       */
547       PetscCall(MatSetValuesStencil(Hessian,1,&row,k,col,v,INSERT_VALUES));
548     }
549   }
550 
551   PetscCall(DMDAVecRestoreArray(user->dm,localX,(void**)&x));
552   PetscCall(DMRestoreLocalVector(user->dm,&localX));
553 
554   PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
555   PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
556 
557   PetscCall(PetscLogFlops(199.0*xm*ym));
558   PetscFunctionReturn(0);
559 }
560 
561 /* ------------------------------------------------------------------- */
562 /*
563    MSA_BoundaryConditions -  Calculates the boundary conditions for
564    the region.
565 
566    Input Parameter:
567 .  user - user-defined application context
568 
569    Output Parameter:
570 .  user - user-defined application context
571 */
572 static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
573 {
574   PetscInt       i,j,k,limit=0,maxits=5;
575   PetscInt       xs,ys,xm,ym,gxs,gys,gxm,gym;
576   PetscInt       mx=user->mx,my=user->my;
577   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
578   PetscReal      one=1.0, two=2.0, three=3.0, tol=1e-10;
579   PetscReal      fnorm,det,hx,hy,xt=0,yt=0;
580   PetscReal      u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
581   PetscReal      b=-0.5, t=0.5, l=-0.5, r=0.5;
582   PetscReal      *boundary;
583   PetscBool      flg;
584 
585   PetscFunctionBegin;
586   /* Get local mesh boundaries */
587   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
588   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
589 
590   bsize=xm+2;
591   lsize=ym+2;
592   rsize=ym+2;
593   tsize=xm+2;
594 
595   PetscCall(PetscMalloc1(bsize,&user->bottom));
596   PetscCall(PetscMalloc1(tsize,&user->top));
597   PetscCall(PetscMalloc1(lsize,&user->left));
598   PetscCall(PetscMalloc1(rsize,&user->right));
599 
600   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
601 
602   for (j=0; j<4; j++) {
603     if (j==0) {
604       yt=b;
605       xt=l+hx*xs;
606       limit=bsize;
607       boundary=user->bottom;
608     } else if (j==1) {
609       yt=t;
610       xt=l+hx*xs;
611       limit=tsize;
612       boundary=user->top;
613     } else if (j==2) {
614       yt=b+hy*ys;
615       xt=l;
616       limit=lsize;
617       boundary=user->left;
618     } else { /* if (j==3) */
619       yt=b+hy*ys;
620       xt=r;
621       limit=rsize;
622       boundary=user->right;
623     }
624 
625     for (i=0; i<limit; i++) {
626       u1=xt;
627       u2=-yt;
628       for (k=0; k<maxits; k++) {
629         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
630         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
631         fnorm=PetscSqrtReal(nf1*nf1+nf2*nf2);
632         if (fnorm <= tol) break;
633         njac11=one+u2*u2-u1*u1;
634         njac12=two*u1*u2;
635         njac21=-two*u1*u2;
636         njac22=-one - u1*u1 + u2*u2;
637         det = njac11*njac22-njac21*njac12;
638         u1 = u1-(njac22*nf1-njac12*nf2)/det;
639         u2 = u2-(njac11*nf2-njac21*nf1)/det;
640       }
641 
642       boundary[i]=u1*u1-u2*u2;
643       if (j==0 || j==1) {
644         xt=xt+hx;
645       } else { /*  if (j==2 || j==3) */
646         yt=yt+hy;
647       }
648 
649     }
650 
651   }
652 
653   /* Scale the boundary if desired */
654   if (1==1) {
655     PetscReal scl = 1.0;
656 
657     PetscCall(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg));
658     if (flg) {
659       for (i=0;i<bsize;i++) user->bottom[i]*=scl;
660     }
661 
662     PetscCall(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg));
663     if (flg) {
664       for (i=0;i<tsize;i++) user->top[i]*=scl;
665     }
666 
667     PetscCall(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg));
668     if (flg) {
669       for (i=0;i<rsize;i++) user->right[i]*=scl;
670     }
671 
672     PetscCall(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg));
673     if (flg) {
674       for (i=0;i<lsize;i++) user->left[i]*=scl;
675     }
676   }
677   PetscFunctionReturn(0);
678 }
679 
680 /* ------------------------------------------------------------------- */
681 /*
682    MSA_InitialPoint - Calculates the initial guess in one of three ways.
683 
684    Input Parameters:
685 .  user - user-defined application context
686 .  X - vector for initial guess
687 
688    Output Parameters:
689 .  X - newly computed initial guess
690 */
691 static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
692 {
693   PetscInt       start2=-1,i,j;
694   PetscReal      start1=0;
695   PetscBool      flg1,flg2;
696 
697   PetscFunctionBegin;
698   PetscCall(PetscOptionsGetReal(NULL,NULL,"-start",&start1,&flg1));
699   PetscCall(PetscOptionsGetInt(NULL,NULL,"-random",&start2,&flg2));
700 
701   if (flg1) { /* The zero vector is reasonable */
702 
703     PetscCall(VecSet(X, start1));
704 
705   } else if (flg2 && start2>0) { /* Try a random start between -0.5 and 0.5 */
706     PetscRandom rctx;  PetscReal np5=-0.5;
707 
708     PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
709     PetscCall(VecSetRandom(X, rctx));
710     PetscCall(PetscRandomDestroy(&rctx));
711     PetscCall(VecShift(X, np5));
712 
713   } else { /* Take an average of the boundary conditions */
714     PetscInt  xs,xm,ys,ym;
715     PetscInt  mx=user->mx,my=user->my;
716     PetscReal **x;
717 
718     /* Get local mesh boundaries */
719     PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
720 
721     /* Get pointers to vector data */
722     PetscCall(DMDAVecGetArray(user->dm,X,(void**)&x));
723 
724     /* Perform local computations */
725     for (j=ys; j<ys+ym; j++) {
726       for (i=xs; i< xs+xm; i++) {
727         x[j][i] = (((j+1)*user->bottom[i-xs+1]+(my-j+1)*user->top[i-xs+1])/(my+2)+((i+1)*user->left[j-ys+1]+(mx-i+1)*user->right[j-ys+1])/(mx+2))/2.0;
728       }
729     }
730     PetscCall(DMDAVecRestoreArray(user->dm,X,(void**)&x));
731     PetscCall(PetscLogFlops(9.0*xm*ym));
732   }
733   PetscFunctionReturn(0);
734 }
735 
736 /*-----------------------------------------------------------------------*/
737 PetscErrorCode My_Monitor(Tao tao, void *ctx)
738 {
739   Vec            X;
740 
741   PetscFunctionBegin;
742   PetscCall(TaoGetSolution(tao,&X));
743   PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD));
744   PetscFunctionReturn(0);
745 }
746 
747 /*TEST
748 
749    build:
750       requires: !complex
751 
752    test:
753       args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3
754       requires: !single
755 
756    test:
757       suffix: 2
758       nsize: 2
759       args: -tao_smonitor -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4
760       filter: grep -v "nls ksp"
761       requires: !single
762 
763    test:
764       suffix: 3
765       nsize: 3
766       args: -tao_smonitor -tao_type cg -tao_cg_type fr -mx 10 -my 10 -tao_gatol 1.e-3
767       requires: !single
768 
769    test:
770       suffix: 5
771       nsize: 2
772       args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3
773       requires: !single
774 
775 TEST*/
776