xref: /petsc/src/tao/bound/tutorials/plate2.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown #include <petscdmda.h>
2c4762a1bSJed Brown #include <petsctao.h>
3c4762a1bSJed Brown 
4c4762a1bSJed Brown static  char help[] =
5c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\
6c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
7c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain, \n\
8c4762a1bSJed Brown boundary values along the edges of the domain, and a plate represented by \n\
9c4762a1bSJed Brown lower boundary conditions, the objective is to find the\n\
10c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
11c4762a1bSJed Brown The command line options are:\n\
12c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
13c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
14c4762a1bSJed Brown   -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\
15c4762a1bSJed Brown   -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\
16c4762a1bSJed Brown   -bheight <ht>, where <ht> = height of the plate\n\
17c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
18c4762a1bSJed Brown                for an average of the boundary conditions\n\n";
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /*
21c4762a1bSJed Brown    User-defined application context - contains data needed by the
22c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
23c4762a1bSJed Brown    FormHessian().
24c4762a1bSJed Brown */
25c4762a1bSJed Brown typedef struct {
26c4762a1bSJed Brown   /* problem parameters */
27c4762a1bSJed Brown   PetscReal      bheight;                  /* Height of plate under the surface */
28c4762a1bSJed Brown   PetscInt       mx, my;                   /* discretization in x, y directions */
29c4762a1bSJed Brown   PetscInt       bmx,bmy;                  /* Size of plate under the surface */
30c4762a1bSJed Brown   Vec            Bottom, Top, Left, Right; /* boundary values */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Working space */
33c4762a1bSJed Brown   Vec         localX, localV;           /* ghosted local vector */
34c4762a1bSJed Brown   DM          dm;                       /* distributed array data structure */
35c4762a1bSJed Brown   Mat         H;
36c4762a1bSJed Brown } AppCtx;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown /* -------- User-defined Routines --------- */
39c4762a1bSJed Brown 
40c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
41c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
42c4762a1bSJed Brown static PetscErrorCode MSA_Plate(Vec,Vec,void*);
43c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
44c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /* For testing matrix free submatrices */
47c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao,Vec,Mat, Mat,void*);
48c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat,Vec,Vec);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown int main(int argc, char **argv)
51c4762a1bSJed Brown {
52c4762a1bSJed Brown   PetscInt               Nx, Ny;               /* number of processors in x- and y- directions */
53c4762a1bSJed Brown   PetscInt               m, N;                 /* number of local and global elements in vectors */
54c4762a1bSJed Brown   Vec                    x,xl,xu;               /* solution vector  and bounds*/
55c4762a1bSJed Brown   PetscBool              flg;                /* A return variable when checking for user options */
56c4762a1bSJed Brown   Tao                    tao;                  /* Tao solver context */
57c4762a1bSJed Brown   ISLocalToGlobalMapping isltog;   /* local-to-global mapping object */
58c4762a1bSJed Brown   Mat                    H_shell;                  /* to test matrix-free submatrices */
59c4762a1bSJed Brown   AppCtx                 user;                 /* user-defined work context */
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Initialize PETSc, TAO */
62*327415f7SBarry Smith   PetscFunctionBeginUser;
639566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv,(char *)0,help));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Specify default dimension of the problem */
66c4762a1bSJed Brown   user.mx = 10; user.my = 10; user.bheight=0.1;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-bheight",&user.bheight,&flg));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   user.bmx = user.mx/2; user.bmy = user.my/2;
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bmx",&user.bmx,&flg));
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bmy",&user.bmy,&flg));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n"));
7863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mx:%" PetscInt_FMT ", my:%" PetscInt_FMT ", bmx:%" PetscInt_FMT ", bmy:%" PetscInt_FMT ", height:%g\n",user.mx,user.my,user.bmx,user.bmy,(double)user.bheight));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Calculate any derived values from parameters */
81c4762a1bSJed Brown   N    = user.mx*user.my;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Let Petsc determine the dimensions of the local vectors */
84c4762a1bSJed Brown   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /*
87c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
88c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
89c4762a1bSJed Brown      the distributed array, Create the vectors.
90c4762a1bSJed Brown   */
919566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(MPI_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm));
929566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
939566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
94c4762a1bSJed Brown   /*
95c4762a1bSJed Brown      Extract global and local vectors from DM; The local vectors are
96c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
97c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
98c4762a1bSJed Brown      the same types.
99c4762a1bSJed Brown   */
1009566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm,&x)); /* Solution */
1019566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(user.dm,&user.localX));
1029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.localX,&user.localV));
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&xl));
1059566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&xu));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* The TAO code begins here */
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /*
110c4762a1bSJed Brown      Create TAO solver and set desired solution method
111c4762a1bSJed Brown      The method must either be TAOTRON or TAOBLMVM
112c4762a1bSJed Brown      If TAOBLMVM is used, then hessian function is not called.
113c4762a1bSJed Brown   */
1149566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
1159566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBLMVM));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* Set initial solution guess; */
1189566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user));
1199566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user,x));
1209566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,x));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* Set routines for function, gradient and hessian evaluation */
1239566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user));
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x,&m));
1269566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H)));
1279566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(user.dm,&isltog));
1309566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(user.H,isltog,isltog));
1319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-matrixfree",&flg));
132c4762a1bSJed Brown   if (flg) {
1339566063dSJacob Faibussowitsch       PetscCall(MatCreateShell(PETSC_COMM_WORLD,m,m,N,N,(void*)&user,&H_shell));
1349566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(H_shell,MATOP_MULT,(void(*)(void))MyMatMult));
1359566063dSJacob Faibussowitsch       PetscCall(MatSetOption(H_shell,MAT_SYMMETRIC,PETSC_TRUE));
1369566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao,H_shell,H_shell,MatrixFreeHessian,(void*)&user));
137c4762a1bSJed Brown   } else {
1389566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Set Variable bounds */
1429566063dSJacob Faibussowitsch   PetscCall(MSA_Plate(xl,xu,(void*)&user));
1439566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao,xl,xu));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Check for any tao command line options */
1469566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1499566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch   PetscCall(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* Free TAO data structures */
1549566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Free PETSc data structures */
1579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
1599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
1609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
1619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localX));
1629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localV));
1639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Bottom));
1649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Top));
1659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Left));
1669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Right));
1679566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
168c4762a1bSJed Brown   if (flg) {
1699566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&H_shell));
170c4762a1bSJed Brown   }
1719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
172b122ec5aSJacob Faibussowitsch   return 0;
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).
176c4762a1bSJed Brown 
177c4762a1bSJed Brown     Input Parameters:
178c4762a1bSJed Brown .   tao     - the Tao context
179c4762a1bSJed Brown .   X      - input vector
180a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
181c4762a1bSJed Brown 
182c4762a1bSJed Brown     Output Parameters:
183c4762a1bSJed Brown .   fcn     - the function value
184c4762a1bSJed Brown .   G      - vector containing the newly evaluated gradient
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    Notes:
187c4762a1bSJed Brown    In this case, we discretize the domain and Create triangles. The
188c4762a1bSJed Brown    surface of each triangle is planar, whose surface area can be easily
189c4762a1bSJed Brown    computed.  The total surface area is found by sweeping through the grid
190c4762a1bSJed Brown    and computing the surface area of the two triangles that have their
191c4762a1bSJed Brown    right angle at the grid point.  The diagonal line segments on the
192c4762a1bSJed Brown    grid that define the triangles run from top left to lower right.
193c4762a1bSJed Brown    The numbering of points starts at the lower left and runs left to
194c4762a1bSJed Brown    right, then bottom to top.
195c4762a1bSJed Brown */
196c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G,void *userCtx)
197c4762a1bSJed Brown {
198c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) userCtx;
199c4762a1bSJed Brown   PetscInt       i,j,row;
200c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
201c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
202c4762a1bSJed Brown   PetscReal      ft=0;
203c4762a1bSJed Brown   PetscReal      zero=0.0;
204c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
205c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
206c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
207c4762a1bSJed Brown   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
208c4762a1bSJed Brown   PetscReal      *g, *x,*left,*right,*bottom,*top;
209c4762a1bSJed Brown   Vec            localX = user->localX, localG = user->localV;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* Get local mesh boundaries */
2129566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
2139566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* Scatter ghost points to local vector */
2169566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
2179566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* Initialize vector to zero */
2209566063dSJacob Faibussowitsch   PetscCall(VecSet(localG, zero));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /* Get pointers to vector data */
2239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX,&x));
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localG,&g));
2259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top,&top));
2269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom,&bottom));
2279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left,&left));
2289566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right,&right));
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
231c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
232c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
233c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
234c4762a1bSJed Brown 
235c4762a1bSJed Brown       xc = x[row];
236c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
237c4762a1bSJed Brown 
238c4762a1bSJed Brown       if (i==0) { /* left side */
239c4762a1bSJed Brown         xl= left[j-ys+1];
240c4762a1bSJed Brown         xlt = left[j-ys+2];
241c4762a1bSJed Brown       } else {
242c4762a1bSJed Brown         xl = x[row-1];
243c4762a1bSJed Brown       }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown       if (j==0) { /* bottom side */
246c4762a1bSJed Brown         xb=bottom[i-xs+1];
247c4762a1bSJed Brown         xrb = bottom[i-xs+2];
248c4762a1bSJed Brown       } else {
249c4762a1bSJed Brown         xb = x[row-gxm];
250c4762a1bSJed Brown       }
251c4762a1bSJed Brown 
252c4762a1bSJed Brown       if (i+1 == gxs+gxm) { /* right side */
253c4762a1bSJed Brown         xr=right[j-ys+1];
254c4762a1bSJed Brown         xrb = right[j-ys];
255c4762a1bSJed Brown       } else {
256c4762a1bSJed Brown         xr = x[row+1];
257c4762a1bSJed Brown       }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown       if (j+1==gys+gym) { /* top side */
260c4762a1bSJed Brown         xt=top[i-xs+1];
261c4762a1bSJed Brown         xlt = top[i-xs];
262c4762a1bSJed Brown       }else {
263c4762a1bSJed Brown         xt = x[row+gxm];
264c4762a1bSJed Brown       }
265c4762a1bSJed Brown 
266c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym) {
267c4762a1bSJed Brown         xlt = x[row-1+gxm];
268c4762a1bSJed Brown       }
269c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm) {
270c4762a1bSJed Brown         xrb = x[row+1-gxm];
271c4762a1bSJed Brown       }
272c4762a1bSJed Brown 
273c4762a1bSJed Brown       d1 = (xc-xl);
274c4762a1bSJed Brown       d2 = (xc-xr);
275c4762a1bSJed Brown       d3 = (xc-xt);
276c4762a1bSJed Brown       d4 = (xc-xb);
277c4762a1bSJed Brown       d5 = (xr-xrb);
278c4762a1bSJed Brown       d6 = (xrb-xb);
279c4762a1bSJed Brown       d7 = (xlt-xl);
280c4762a1bSJed Brown       d8 = (xt-xlt);
281c4762a1bSJed Brown 
282c4762a1bSJed Brown       df1dxc = d1*hydhx;
283c4762a1bSJed Brown       df2dxc = (d1*hydhx + d4*hxdhy);
284c4762a1bSJed Brown       df3dxc = d3*hxdhy;
285c4762a1bSJed Brown       df4dxc = (d2*hydhx + d3*hxdhy);
286c4762a1bSJed Brown       df5dxc = d2*hydhx;
287c4762a1bSJed Brown       df6dxc = d4*hxdhy;
288c4762a1bSJed Brown 
289c4762a1bSJed Brown       d1 *= rhx;
290c4762a1bSJed Brown       d2 *= rhx;
291c4762a1bSJed Brown       d3 *= rhy;
292c4762a1bSJed Brown       d4 *= rhy;
293c4762a1bSJed Brown       d5 *= rhy;
294c4762a1bSJed Brown       d6 *= rhx;
295c4762a1bSJed Brown       d7 *= rhy;
296c4762a1bSJed Brown       d8 *= rhx;
297c4762a1bSJed Brown 
298c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
299c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
300c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
301c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
302c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
303c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
304c4762a1bSJed Brown 
305c4762a1bSJed Brown       ft = ft + (f2 + f4);
306c4762a1bSJed Brown 
307c4762a1bSJed Brown       df1dxc /= f1;
308c4762a1bSJed Brown       df2dxc /= f2;
309c4762a1bSJed Brown       df3dxc /= f3;
310c4762a1bSJed Brown       df4dxc /= f4;
311c4762a1bSJed Brown       df5dxc /= f5;
312c4762a1bSJed Brown       df6dxc /= f6;
313c4762a1bSJed Brown 
314c4762a1bSJed Brown       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5;
315c4762a1bSJed Brown 
316c4762a1bSJed Brown     }
317c4762a1bSJed Brown   }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
320c4762a1bSJed Brown   if (xs==0) { /* left side */
321c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
322c4762a1bSJed Brown       d3=(left[j-ys+1] - left[j-ys+2])*rhy;
323c4762a1bSJed Brown       d2=(left[j-ys+1] - x[(j-gys)*gxm])*rhx;
324c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
325c4762a1bSJed Brown     }
326c4762a1bSJed Brown   }
327c4762a1bSJed Brown   if (ys==0) { /* bottom side */
328c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
329c4762a1bSJed Brown       d2=(bottom[i+1-xs]-bottom[i-xs+2])*rhx;
330c4762a1bSJed Brown       d3=(bottom[i-xs+1]-x[i-gxs])*rhy;
331c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown   }
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   if (xs+xm==mx) { /* right side */
336c4762a1bSJed Brown     for (j=ys; j< ys+ym; j++) {
337c4762a1bSJed Brown       d1=(x[(j+1-gys)*gxm-1]-right[j-ys+1])*rhx;
338c4762a1bSJed Brown       d4=(right[j-ys]-right[j-ys+1])*rhy;
339c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
340c4762a1bSJed Brown     }
341c4762a1bSJed Brown   }
342c4762a1bSJed Brown   if (ys+ym==my) { /* top side */
343c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
344c4762a1bSJed Brown       d1=(x[(gym-1)*gxm + i-gxs] - top[i-xs+1])*rhy;
345c4762a1bSJed Brown       d4=(top[i-xs+1] - top[i-xs])*rhx;
346c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
347c4762a1bSJed Brown     }
348c4762a1bSJed Brown   }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   if (ys==0 && xs==0) {
351c4762a1bSJed Brown     d1=(left[0]-left[1])*rhy;
352c4762a1bSJed Brown     d2=(bottom[0]-bottom[1])*rhx;
353c4762a1bSJed Brown     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
354c4762a1bSJed Brown   }
355c4762a1bSJed Brown   if (ys+ym == my && xs+xm == mx) {
356c4762a1bSJed Brown     d1=(right[ym+1] - right[ym])*rhy;
357c4762a1bSJed Brown     d2=(top[xm+1] - top[xm])*rhx;
358c4762a1bSJed Brown     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
359c4762a1bSJed Brown   }
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   ft=ft*area;
3629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD));
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /* Restore vectors */
3659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX,&x));
3669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localG,&g));
3679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left,&left));
3689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top,&top));
3699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom,&bottom));
3709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right,&right));
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /* Scatter values to global vector */
3739566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G));
3749566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G));
375c4762a1bSJed Brown 
3769566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(70.0*xm*ym));
377c4762a1bSJed Brown 
378c4762a1bSJed Brown   return 0;
379c4762a1bSJed Brown }
380c4762a1bSJed Brown 
381c4762a1bSJed Brown /* ------------------------------------------------------------------- */
382c4762a1bSJed Brown /*
383c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
384c4762a1bSJed Brown 
385c4762a1bSJed Brown    Input Parameters:
386c4762a1bSJed Brown .  tao  - the Tao context
387c4762a1bSJed Brown .  x    - input vector
388a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
389c4762a1bSJed Brown 
390c4762a1bSJed Brown    Output Parameters:
391c4762a1bSJed Brown .  A    - Hessian matrix
392c4762a1bSJed Brown .  B    - optionally different preconditioning matrix
393c4762a1bSJed Brown 
394c4762a1bSJed Brown    Notes:
395c4762a1bSJed Brown    Due to mesh point reordering with DMs, we must always work
396c4762a1bSJed Brown    with the local mesh points, and then transform them to the new
397c4762a1bSJed Brown    global numbering with the local-to-global mapping.  We cannot work
398c4762a1bSJed Brown    directly with the global numbers for the original uniprocessor mesh!
399c4762a1bSJed Brown 
400c4762a1bSJed Brown    Two methods are available for imposing this transformation
401c4762a1bSJed Brown    when setting matrix entries:
402c4762a1bSJed Brown      (A) MatSetValuesLocal(), using the local ordering (including
403c4762a1bSJed Brown          ghost points!)
404c4762a1bSJed Brown          - Do the following two steps once, before calling TaoSolve()
405c4762a1bSJed Brown            - Use DMGetISLocalToGlobalMapping() to extract the
406c4762a1bSJed Brown              local-to-global map from the DM
407c4762a1bSJed Brown            - Associate this map with the matrix by calling
408c4762a1bSJed Brown              MatSetLocalToGlobalMapping()
409c4762a1bSJed Brown          - Then set matrix entries using the local ordering
410c4762a1bSJed Brown            by calling MatSetValuesLocal()
411c4762a1bSJed Brown      (B) MatSetValues(), using the global ordering
412c4762a1bSJed Brown          - Use DMGetGlobalIndices() to extract the local-to-global map
413c4762a1bSJed Brown          - Then apply this map explicitly yourself
414c4762a1bSJed Brown          - Set matrix entries using the global ordering by calling
415c4762a1bSJed Brown            MatSetValues()
416c4762a1bSJed Brown    Option (A) seems cleaner/easier in many cases, and is the procedure
417c4762a1bSJed Brown    used in this example.
418c4762a1bSJed Brown */
419c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat Hptr, Mat Hessian, void *ptr)
420c4762a1bSJed Brown {
421c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
422c4762a1bSJed Brown   PetscInt       i,j,k,row;
423c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
424c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym,col[7];
425c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
426c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
427c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
428c4762a1bSJed Brown   PetscReal      hl,hr,ht,hb,hc,htl,hbr;
429c4762a1bSJed Brown   PetscReal      *x,*left,*right,*bottom,*top;
430c4762a1bSJed Brown   PetscReal      v[7];
431c4762a1bSJed Brown   Vec            localX = user->localX;
432c4762a1bSJed Brown   PetscBool      assembled;
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /* Set various matrix options */
4359566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   /* Initialize matrix entries to zero */
4389566063dSJacob Faibussowitsch   PetscCall(MatAssembled(Hessian,&assembled));
4399566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(Hessian));
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /* Get local mesh boundaries */
4429566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
4439566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /* Scatter ghost points to local vector */
4469566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
4479566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /* Get pointers to vector data */
4509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX,&x));
4519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top,&top));
4529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom,&bottom));
4539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left,&left));
4549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right,&right));
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
457c4762a1bSJed Brown 
458c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++) {
459c4762a1bSJed Brown 
460c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
461c4762a1bSJed Brown 
462c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
463c4762a1bSJed Brown 
464c4762a1bSJed Brown       xc = x[row];
465c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
466c4762a1bSJed Brown 
467c4762a1bSJed Brown       /* Left side */
468c4762a1bSJed Brown       if (i==gxs) {
469c4762a1bSJed Brown         xl= left[j-ys+1];
470c4762a1bSJed Brown         xlt = left[j-ys+2];
471c4762a1bSJed Brown       } else {
472c4762a1bSJed Brown         xl = x[row-1];
473c4762a1bSJed Brown       }
474c4762a1bSJed Brown 
475c4762a1bSJed Brown       if (j==gys) {
476c4762a1bSJed Brown         xb=bottom[i-xs+1];
477c4762a1bSJed Brown         xrb = bottom[i-xs+2];
478c4762a1bSJed Brown       } else {
479c4762a1bSJed Brown         xb = x[row-gxm];
480c4762a1bSJed Brown       }
481c4762a1bSJed Brown 
482c4762a1bSJed Brown       if (i+1 == gxs+gxm) {
483c4762a1bSJed Brown         xr=right[j-ys+1];
484c4762a1bSJed Brown         xrb = right[j-ys];
485c4762a1bSJed Brown       } else {
486c4762a1bSJed Brown         xr = x[row+1];
487c4762a1bSJed Brown       }
488c4762a1bSJed Brown 
489c4762a1bSJed Brown       if (j+1==gys+gym) {
490c4762a1bSJed Brown         xt=top[i-xs+1];
491c4762a1bSJed Brown         xlt = top[i-xs];
492c4762a1bSJed Brown       }else {
493c4762a1bSJed Brown         xt = x[row+gxm];
494c4762a1bSJed Brown       }
495c4762a1bSJed Brown 
496c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym) {
497c4762a1bSJed Brown         xlt = x[row-1+gxm];
498c4762a1bSJed Brown       }
499c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm) {
500c4762a1bSJed Brown         xrb = x[row+1-gxm];
501c4762a1bSJed Brown       }
502c4762a1bSJed Brown 
503c4762a1bSJed Brown       d1 = (xc-xl)*rhx;
504c4762a1bSJed Brown       d2 = (xc-xr)*rhx;
505c4762a1bSJed Brown       d3 = (xc-xt)*rhy;
506c4762a1bSJed Brown       d4 = (xc-xb)*rhy;
507c4762a1bSJed Brown       d5 = (xrb-xr)*rhy;
508c4762a1bSJed Brown       d6 = (xrb-xb)*rhx;
509c4762a1bSJed Brown       d7 = (xlt-xl)*rhy;
510c4762a1bSJed Brown       d8 = (xlt-xt)*rhx;
511c4762a1bSJed Brown 
512c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
513c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
514c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
515c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
516c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
517c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
518c4762a1bSJed Brown 
519c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
520c4762a1bSJed Brown         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
521c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
522c4762a1bSJed Brown         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
523c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
524c4762a1bSJed Brown         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
525c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
526c4762a1bSJed Brown         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
527c4762a1bSJed Brown 
528c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
529c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
530c4762a1bSJed Brown 
531c4762a1bSJed Brown       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
532c4762a1bSJed Brown         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
533c4762a1bSJed Brown         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
534c4762a1bSJed Brown         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
535c4762a1bSJed Brown 
536c4762a1bSJed Brown       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;
537c4762a1bSJed Brown 
538c4762a1bSJed Brown       k=0;
539c4762a1bSJed Brown       if (j>0) {
540c4762a1bSJed Brown         v[k]=hb; col[k]=row - gxm; k++;
541c4762a1bSJed Brown       }
542c4762a1bSJed Brown 
543c4762a1bSJed Brown       if (j>0 && i < mx -1) {
544c4762a1bSJed Brown         v[k]=hbr; col[k]=row - gxm+1; k++;
545c4762a1bSJed Brown       }
546c4762a1bSJed Brown 
547c4762a1bSJed Brown       if (i>0) {
548c4762a1bSJed Brown         v[k]= hl; col[k]=row - 1; k++;
549c4762a1bSJed Brown       }
550c4762a1bSJed Brown 
551c4762a1bSJed Brown       v[k]= hc; col[k]=row; k++;
552c4762a1bSJed Brown 
553c4762a1bSJed Brown       if (i < mx-1) {
554c4762a1bSJed Brown         v[k]= hr; col[k]=row+1; k++;
555c4762a1bSJed Brown       }
556c4762a1bSJed Brown 
557c4762a1bSJed Brown       if (i>0 && j < my-1) {
558c4762a1bSJed Brown         v[k]= htl; col[k] = row+gxm-1; k++;
559c4762a1bSJed Brown       }
560c4762a1bSJed Brown 
561c4762a1bSJed Brown       if (j < my-1) {
562c4762a1bSJed Brown         v[k]= ht; col[k] = row+gxm; k++;
563c4762a1bSJed Brown       }
564c4762a1bSJed Brown 
565c4762a1bSJed Brown       /*
566c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
567c4762a1bSJed Brown          earlier, in the main routine.
568c4762a1bSJed Brown       */
5699566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Hessian,1,&row,k,col,v,INSERT_VALUES));
570c4762a1bSJed Brown 
571c4762a1bSJed Brown     }
572c4762a1bSJed Brown   }
573c4762a1bSJed Brown 
574c4762a1bSJed Brown   /* Restore vectors */
5759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX,&x));
5769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left,&left));
5779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top,&top));
5789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom,&bottom));
5799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right,&right));
580c4762a1bSJed Brown 
581c4762a1bSJed Brown   /* Assemble the matrix */
5829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
5839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
584c4762a1bSJed Brown 
5859566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0*xm*ym));
586c4762a1bSJed Brown   return 0;
587c4762a1bSJed Brown }
588c4762a1bSJed Brown 
589c4762a1bSJed Brown /* ------------------------------------------------------------------- */
590c4762a1bSJed Brown /*
591c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
592c4762a1bSJed Brown    the region.
593c4762a1bSJed Brown 
594c4762a1bSJed Brown    Input Parameter:
595c4762a1bSJed Brown .  user - user-defined application context
596c4762a1bSJed Brown 
597c4762a1bSJed Brown    Output Parameter:
598c4762a1bSJed Brown .  user - user-defined application context
599c4762a1bSJed Brown */
600c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
601c4762a1bSJed Brown {
602c4762a1bSJed Brown   PetscInt   i,j,k,maxits=5,limit=0;
603c4762a1bSJed Brown   PetscInt   xs,ys,xm,ym,gxs,gys,gxm,gym;
604c4762a1bSJed Brown   PetscInt   mx=user->mx,my=user->my;
605c4762a1bSJed Brown   PetscInt   bsize=0, lsize=0, tsize=0, rsize=0;
606c4762a1bSJed Brown   PetscReal  one=1.0, two=2.0, three=3.0, scl=1.0, tol=1e-10;
607c4762a1bSJed Brown   PetscReal  fnorm,det,hx,hy,xt=0,yt=0;
608c4762a1bSJed Brown   PetscReal  u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
609c4762a1bSJed Brown   PetscReal  b=-0.5, t=0.5, l=-0.5, r=0.5;
610c4762a1bSJed Brown   PetscReal  *boundary;
611c4762a1bSJed Brown   PetscBool  flg;
612c4762a1bSJed Brown   Vec        Bottom,Top,Right,Left;
613c4762a1bSJed Brown 
614c4762a1bSJed Brown   /* Get local mesh boundaries */
6159566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
6169566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   bsize=xm+2;
619c4762a1bSJed Brown   lsize=ym+2;
620c4762a1bSJed Brown   rsize=ym+2;
621c4762a1bSJed Brown   tsize=xm+2;
622c4762a1bSJed Brown 
6239566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom));
6249566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top));
6259566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left));
6269566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right));
627c4762a1bSJed Brown 
628c4762a1bSJed Brown   user->Top=Top;
629c4762a1bSJed Brown   user->Left=Left;
630c4762a1bSJed Brown   user->Bottom=Bottom;
631c4762a1bSJed Brown   user->Right=Right;
632c4762a1bSJed Brown 
633c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
634c4762a1bSJed Brown 
635c4762a1bSJed Brown   for (j=0; j<4; j++) {
636c4762a1bSJed Brown     if (j==0) {
637c4762a1bSJed Brown       yt=b;
638c4762a1bSJed Brown       xt=l+hx*xs;
639c4762a1bSJed Brown       limit=bsize;
640c4762a1bSJed Brown       VecGetArray(Bottom,&boundary);
641c4762a1bSJed Brown     } else if (j==1) {
642c4762a1bSJed Brown       yt=t;
643c4762a1bSJed Brown       xt=l+hx*xs;
644c4762a1bSJed Brown       limit=tsize;
645c4762a1bSJed Brown       VecGetArray(Top,&boundary);
646c4762a1bSJed Brown     } else if (j==2) {
647c4762a1bSJed Brown       yt=b+hy*ys;
648c4762a1bSJed Brown       xt=l;
649c4762a1bSJed Brown       limit=lsize;
650c4762a1bSJed Brown       VecGetArray(Left,&boundary);
651c4762a1bSJed Brown     } else if (j==3) {
652c4762a1bSJed Brown       yt=b+hy*ys;
653c4762a1bSJed Brown       xt=r;
654c4762a1bSJed Brown       limit=rsize;
655c4762a1bSJed Brown       VecGetArray(Right,&boundary);
656c4762a1bSJed Brown     }
657c4762a1bSJed Brown 
658c4762a1bSJed Brown     for (i=0; i<limit; i++) {
659c4762a1bSJed Brown       u1=xt;
660c4762a1bSJed Brown       u2=-yt;
661c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
662c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
663c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
664c4762a1bSJed Brown         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
665c4762a1bSJed Brown         if (fnorm <= tol) break;
666c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
667c4762a1bSJed Brown         njac12=two*u1*u2;
668c4762a1bSJed Brown         njac21=-two*u1*u2;
669c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
670c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
671c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
672c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
673c4762a1bSJed Brown       }
674c4762a1bSJed Brown 
675c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
676c4762a1bSJed Brown       if (j==0 || j==1) {
677c4762a1bSJed Brown         xt=xt+hx;
678c4762a1bSJed Brown       } else if (j==2 || j==3) {
679c4762a1bSJed Brown         yt=yt+hy;
680c4762a1bSJed Brown       }
681c4762a1bSJed Brown     }
682c4762a1bSJed Brown     if (j==0) {
6839566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Bottom,&boundary));
684c4762a1bSJed Brown     } else if (j==1) {
6859566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Top,&boundary));
686c4762a1bSJed Brown     } else if (j==2) {
6879566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Left,&boundary));
688c4762a1bSJed Brown     } else if (j==3) {
6899566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Right,&boundary));
690c4762a1bSJed Brown     }
691c4762a1bSJed Brown   }
692c4762a1bSJed Brown 
693c4762a1bSJed Brown   /* Scale the boundary if desired */
694c4762a1bSJed Brown 
6959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg));
6961baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Bottom, scl));
6979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg));
6981baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Top, scl));
6999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg));
7001baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Right, scl));
701c4762a1bSJed Brown 
7029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg));
7031baa6e33SBarry Smith   if (flg) PetscCall(VecScale(Left, scl));
704c4762a1bSJed Brown   return 0;
705c4762a1bSJed Brown }
706c4762a1bSJed Brown 
707c4762a1bSJed Brown /* ------------------------------------------------------------------- */
708c4762a1bSJed Brown /*
709c4762a1bSJed Brown    MSA_Plate -  Calculates an obstacle for surface to stretch over.
710c4762a1bSJed Brown 
711c4762a1bSJed Brown    Input Parameter:
712c4762a1bSJed Brown .  user - user-defined application context
713c4762a1bSJed Brown 
714c4762a1bSJed Brown    Output Parameter:
715c4762a1bSJed Brown .  user - user-defined application context
716c4762a1bSJed Brown */
71770a7d78aSStefano Zampini static PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx)
71870a7d78aSStefano Zampini {
719c4762a1bSJed Brown   AppCtx         *user=(AppCtx *)ctx;
720c4762a1bSJed Brown   PetscInt       i,j,row;
721c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym;
722c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my, bmy, bmx;
723c4762a1bSJed Brown   PetscReal      t1,t2,t3;
724c4762a1bSJed Brown   PetscReal      *xl, lb=PETSC_NINFINITY, ub=PETSC_INFINITY;
725c4762a1bSJed Brown   PetscBool      cylinder;
726c4762a1bSJed Brown 
727c4762a1bSJed Brown   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
728c4762a1bSJed Brown   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
729c4762a1bSJed Brown   bmy=user->bmy; bmx=user->bmx;
730c4762a1bSJed Brown 
7319566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
732c4762a1bSJed Brown 
7339566063dSJacob Faibussowitsch   PetscCall(VecSet(XL, lb));
7349566063dSJacob Faibussowitsch   PetscCall(VecSet(XU, ub));
735c4762a1bSJed Brown 
7369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(XL,&xl));
737c4762a1bSJed Brown 
7389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-cylinder",&cylinder));
739c4762a1bSJed Brown   /* Compute the optional lower box */
740c4762a1bSJed Brown   if (cylinder) {
741c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
742c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
743c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
744c4762a1bSJed Brown         t1=(2.0*i-mx)*bmy;
745c4762a1bSJed Brown         t2=(2.0*j-my)*bmx;
746c4762a1bSJed Brown         t3=bmx*bmx*bmy*bmy;
747c4762a1bSJed Brown         if (t1*t1 + t2*t2 <= t3) {
748c4762a1bSJed Brown           xl[row] = user->bheight;
749c4762a1bSJed Brown         }
750c4762a1bSJed Brown       }
751c4762a1bSJed Brown     }
752c4762a1bSJed Brown   } else {
753c4762a1bSJed Brown     /* Compute the optional lower box */
754c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
755c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
756c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
757c4762a1bSJed Brown         if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
758c4762a1bSJed Brown             j>=(my-bmy)/2 && j<my-(my-bmy)/2) {
759c4762a1bSJed Brown           xl[row] = user->bheight;
760c4762a1bSJed Brown         }
761c4762a1bSJed Brown       }
762c4762a1bSJed Brown     }
763c4762a1bSJed Brown   }
7649566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(XL,&xl));
765c4762a1bSJed Brown 
766c4762a1bSJed Brown   return 0;
767c4762a1bSJed Brown }
768c4762a1bSJed Brown 
769c4762a1bSJed Brown /* ------------------------------------------------------------------- */
770c4762a1bSJed Brown /*
771c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
772c4762a1bSJed Brown 
773c4762a1bSJed Brown    Input Parameters:
774c4762a1bSJed Brown .  user - user-defined application context
775c4762a1bSJed Brown .  X - vector for initial guess
776c4762a1bSJed Brown 
777c4762a1bSJed Brown    Output Parameters:
778c4762a1bSJed Brown .  X - newly computed initial guess
779c4762a1bSJed Brown */
780c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
781c4762a1bSJed Brown {
782c4762a1bSJed Brown   PetscInt       start=-1,i,j;
783c4762a1bSJed Brown   PetscReal      zero=0.0;
784c4762a1bSJed Brown   PetscBool      flg;
785c4762a1bSJed Brown 
7869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
787c4762a1bSJed Brown   if (flg && start==0) { /* The zero vector is reasonable */
7889566063dSJacob Faibussowitsch     PetscCall(VecSet(X, zero));
789c4762a1bSJed Brown   } else if (flg && start>0) { /* Try a random start between -0.5 and 0.5 */
790c4762a1bSJed Brown     PetscRandom rctx;  PetscReal np5=-0.5;
791c4762a1bSJed Brown 
7929566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(MPI_COMM_WORLD,&rctx));
793c4762a1bSJed Brown     for (i=0; i<start; i++) {
7949566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(X, rctx));
795c4762a1bSJed Brown     }
7969566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rctx));
7979566063dSJacob Faibussowitsch     PetscCall(VecShift(X, np5));
798c4762a1bSJed Brown 
799c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
800c4762a1bSJed Brown 
801c4762a1bSJed Brown     PetscInt row,xs,xm,gxs,gxm,ys,ym,gys,gym;
802c4762a1bSJed Brown     PetscInt mx=user->mx,my=user->my;
803c4762a1bSJed Brown     PetscReal *x,*left,*right,*bottom,*top;
804c4762a1bSJed Brown     Vec    localX = user->localX;
805c4762a1bSJed Brown 
806c4762a1bSJed Brown     /* Get local mesh boundaries */
8079566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
8089566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
809c4762a1bSJed Brown 
810c4762a1bSJed Brown     /* Get pointers to vector data */
8119566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Top,&top));
8129566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Bottom,&bottom));
8139566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Left,&left));
8149566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Right,&right));
815c4762a1bSJed Brown 
8169566063dSJacob Faibussowitsch     PetscCall(VecGetArray(localX,&x));
817c4762a1bSJed Brown     /* Perform local computations */
818c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
819c4762a1bSJed Brown       for (i=xs; i< xs+xm; i++) {
820c4762a1bSJed Brown         row=(j-gys)*gxm + (i-gxs);
8212da392ccSBarry Smith         x[row] = ((j+1)*bottom[i-xs+1]/my + (my-j+1)*top[i-xs+1]/(my+2)+(i+1)*left[j-ys+1]/mx + (mx-i+1)*right[j-ys+1]/(mx+2))/2.0;
822c4762a1bSJed Brown       }
823c4762a1bSJed Brown     }
824c4762a1bSJed Brown 
825c4762a1bSJed Brown     /* Restore vectors */
8269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(localX,&x));
827c4762a1bSJed Brown 
8289566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Left,&left));
8299566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Top,&top));
8309566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Bottom,&bottom));
8319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Right,&right));
832c4762a1bSJed Brown 
833c4762a1bSJed Brown     /* Scatter values into global vector */
8349566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X));
8359566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X));
836c4762a1bSJed Brown 
837c4762a1bSJed Brown   }
838c4762a1bSJed Brown   return 0;
839c4762a1bSJed Brown }
840c4762a1bSJed Brown 
841c4762a1bSJed Brown /* For testing matrix free submatrices */
842c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
843c4762a1bSJed Brown {
844c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
845c4762a1bSJed Brown   PetscFunctionBegin;
8469566063dSJacob Faibussowitsch   PetscCall(FormHessian(tao,x,user->H,user->H,ptr));
847c4762a1bSJed Brown   PetscFunctionReturn(0);
848c4762a1bSJed Brown }
849c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
850c4762a1bSJed Brown {
851c4762a1bSJed Brown   void           *ptr;
852c4762a1bSJed Brown   AppCtx         *user;
853c4762a1bSJed Brown   PetscFunctionBegin;
8549566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(H_shell,&ptr));
855c4762a1bSJed Brown   user = (AppCtx*)ptr;
8569566063dSJacob Faibussowitsch   PetscCall(MatMult(user->H,X,Y));
857c4762a1bSJed Brown   PetscFunctionReturn(0);
858c4762a1bSJed Brown }
859c4762a1bSJed Brown 
860c4762a1bSJed Brown /*TEST
861c4762a1bSJed Brown 
862c4762a1bSJed Brown    build:
863c4762a1bSJed Brown       requires: !complex
864c4762a1bSJed Brown 
865c4762a1bSJed Brown    test:
866c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
867c4762a1bSJed Brown       requires: !single
868c4762a1bSJed Brown 
869c4762a1bSJed Brown    test:
870c4762a1bSJed Brown       suffix: 2
871c4762a1bSJed Brown       nsize: 2
872c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
873c4762a1bSJed Brown       requires: !single
874c4762a1bSJed Brown 
875c4762a1bSJed Brown    test:
876c4762a1bSJed Brown       suffix: 3
877c4762a1bSJed Brown       nsize: 3
878c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
879c4762a1bSJed Brown       requires: !single
880c4762a1bSJed Brown 
881c4762a1bSJed Brown    test:
882c4762a1bSJed Brown       suffix: 4
883c4762a1bSJed Brown       nsize: 3
884c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type mask -tao_type tron -tao_gatol 1.e-5
885c4762a1bSJed Brown       requires: !single
886c4762a1bSJed Brown 
887c4762a1bSJed Brown    test:
888c4762a1bSJed Brown       suffix: 5
889c4762a1bSJed Brown       nsize: 3
890c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -pc_type none -tao_type tron -tao_gatol 1.e-5
891c4762a1bSJed Brown       requires: !single
892c4762a1bSJed Brown 
893c4762a1bSJed Brown    test:
894c4762a1bSJed Brown       suffix: 6
895c4762a1bSJed Brown       nsize: 3
896c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -tao_type blmvm -tao_gatol 1.e-4
897c4762a1bSJed Brown       requires: !single
898c4762a1bSJed Brown 
899c4762a1bSJed Brown    test:
900c4762a1bSJed Brown       suffix: 7
901c4762a1bSJed Brown       nsize: 3
902c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -pc_type none -tao_type gpcg -tao_gatol 1.e-5
903c4762a1bSJed Brown       requires: !single
904c4762a1bSJed Brown 
905c4762a1bSJed Brown    test:
906c4762a1bSJed Brown       suffix: 8
907c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
908c4762a1bSJed Brown       requires: !single
909c4762a1bSJed Brown 
910c4762a1bSJed Brown    test:
911c4762a1bSJed Brown       suffix: 9
912c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
913c4762a1bSJed Brown       requires: !single
914c4762a1bSJed Brown 
915c4762a1bSJed Brown    test:
916c4762a1bSJed Brown       suffix: 10
917c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
918c4762a1bSJed Brown       requires: !single
919c4762a1bSJed Brown 
920c4762a1bSJed Brown    test:
921c4762a1bSJed Brown       suffix: 11
922c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
923c4762a1bSJed Brown       requires: !single
924c4762a1bSJed Brown 
925c4762a1bSJed Brown    test:
926c4762a1bSJed Brown       suffix: 12
927c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
928c4762a1bSJed Brown       requires: !single
929c4762a1bSJed Brown 
930c4762a1bSJed Brown    test:
931c4762a1bSJed Brown       suffix: 13
932c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
933c4762a1bSJed Brown       requires: !single
934c4762a1bSJed Brown 
935c4762a1bSJed Brown    test:
936c4762a1bSJed Brown       suffix: 14
937c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
938c4762a1bSJed Brown       requires: !single
939c4762a1bSJed Brown 
940c4762a1bSJed Brown    test:
941c4762a1bSJed Brown       suffix: 15
942c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
943c4762a1bSJed Brown       requires: !single
944c4762a1bSJed Brown 
945c4762a1bSJed Brown    test:
946c4762a1bSJed Brown      suffix: 16
947c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
948c4762a1bSJed Brown      requires: !single
949c4762a1bSJed Brown 
950c4762a1bSJed Brown    test:
951c4762a1bSJed Brown      suffix: 17
952c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
953c4762a1bSJed Brown      requires: !single
954c4762a1bSJed Brown 
95534ad9904SAlp Dener    test:
95634ad9904SAlp Dener      suffix: 18
95734ad9904SAlp Dener      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
95834ad9904SAlp Dener      requires: !single
95934ad9904SAlp Dener 
96034ad9904SAlp Dener    test:
96134ad9904SAlp Dener      suffix: 19
96234ad9904SAlp Dener      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
96334ad9904SAlp Dener      requires: !single
96434ad9904SAlp Dener 
96534ad9904SAlp Dener    test:
96634ad9904SAlp Dener      suffix: 20
96734ad9904SAlp Dener      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
96834ad9904SAlp Dener 
969c4762a1bSJed Brown TEST*/
970