xref: /petsc/src/tao/bound/tutorials/plate2.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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 */
629566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv,(char *)0,help));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Specify default dimension of the problem */
65c4762a1bSJed Brown   user.mx = 10; user.my = 10; user.bheight=0.1;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-bheight",&user.bheight,&flg));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   user.bmx = user.mx/2; user.bmy = user.my/2;
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bmx",&user.bmx,&flg));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bmy",&user.bmy,&flg));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n"));
77*63a3b9bcSJacob 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));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Calculate any derived values from parameters */
80c4762a1bSJed Brown   N    = user.mx*user.my;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Let Petsc determine the dimensions of the local vectors */
83c4762a1bSJed Brown   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /*
86c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
87c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
88c4762a1bSJed Brown      the distributed array, Create the vectors.
89c4762a1bSJed Brown   */
909566063dSJacob 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));
919566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
929566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
93c4762a1bSJed Brown   /*
94c4762a1bSJed Brown      Extract global and local vectors from DM; The local vectors are
95c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
96c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
97c4762a1bSJed Brown      the same types.
98c4762a1bSJed Brown   */
999566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm,&x)); /* Solution */
1009566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(user.dm,&user.localX));
1019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.localX,&user.localV));
102c4762a1bSJed Brown 
1039566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&xl));
1049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&xu));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* The TAO code begins here */
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /*
109c4762a1bSJed Brown      Create TAO solver and set desired solution method
110c4762a1bSJed Brown      The method must either be TAOTRON or TAOBLMVM
111c4762a1bSJed Brown      If TAOBLMVM is used, then hessian function is not called.
112c4762a1bSJed Brown   */
1139566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
1149566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBLMVM));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* Set initial solution guess; */
1179566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user));
1189566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user,x));
1199566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,x));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* Set routines for function, gradient and hessian evaluation */
1229566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user));
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x,&m));
1259566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H)));
1269566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
127c4762a1bSJed Brown 
1289566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(user.dm,&isltog));
1299566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(user.H,isltog,isltog));
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-matrixfree",&flg));
131c4762a1bSJed Brown   if (flg) {
1329566063dSJacob Faibussowitsch       PetscCall(MatCreateShell(PETSC_COMM_WORLD,m,m,N,N,(void*)&user,&H_shell));
1339566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(H_shell,MATOP_MULT,(void(*)(void))MyMatMult));
1349566063dSJacob Faibussowitsch       PetscCall(MatSetOption(H_shell,MAT_SYMMETRIC,PETSC_TRUE));
1359566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao,H_shell,H_shell,MatrixFreeHessian,(void*)&user));
136c4762a1bSJed Brown   } else {
1379566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* Set Variable bounds */
1419566063dSJacob Faibussowitsch   PetscCall(MSA_Plate(xl,xu,(void*)&user));
1429566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao,xl,xu));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Check for any tao command line options */
1459566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1489566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch   PetscCall(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Free TAO data structures */
1539566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Free PETSc data structures */
1569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
1589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
1599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
1609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localX));
1619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.localV));
1629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Bottom));
1639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Top));
1649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Left));
1659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Right));
1669566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
167c4762a1bSJed Brown   if (flg) {
1689566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&H_shell));
169c4762a1bSJed Brown   }
1709566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
171b122ec5aSJacob Faibussowitsch   return 0;
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).
175c4762a1bSJed Brown 
176c4762a1bSJed Brown     Input Parameters:
177c4762a1bSJed Brown .   tao     - the Tao context
178c4762a1bSJed Brown .   X      - input vector
179a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
180c4762a1bSJed Brown 
181c4762a1bSJed Brown     Output Parameters:
182c4762a1bSJed Brown .   fcn     - the function value
183c4762a1bSJed Brown .   G      - vector containing the newly evaluated gradient
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    Notes:
186c4762a1bSJed Brown    In this case, we discretize the domain and Create triangles. The
187c4762a1bSJed Brown    surface of each triangle is planar, whose surface area can be easily
188c4762a1bSJed Brown    computed.  The total surface area is found by sweeping through the grid
189c4762a1bSJed Brown    and computing the surface area of the two triangles that have their
190c4762a1bSJed Brown    right angle at the grid point.  The diagonal line segments on the
191c4762a1bSJed Brown    grid that define the triangles run from top left to lower right.
192c4762a1bSJed Brown    The numbering of points starts at the lower left and runs left to
193c4762a1bSJed Brown    right, then bottom to top.
194c4762a1bSJed Brown */
195c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G,void *userCtx)
196c4762a1bSJed Brown {
197c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) userCtx;
198c4762a1bSJed Brown   PetscInt       i,j,row;
199c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
200c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
201c4762a1bSJed Brown   PetscReal      ft=0;
202c4762a1bSJed Brown   PetscReal      zero=0.0;
203c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
204c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
205c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
206c4762a1bSJed Brown   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
207c4762a1bSJed Brown   PetscReal      *g, *x,*left,*right,*bottom,*top;
208c4762a1bSJed Brown   Vec            localX = user->localX, localG = user->localV;
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* Get local mesh boundaries */
2119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
2129566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* Scatter ghost points to local vector */
2159566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
2169566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* Initialize vector to zero */
2199566063dSJacob Faibussowitsch   PetscCall(VecSet(localG, zero));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Get pointers to vector data */
2229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX,&x));
2239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localG,&g));
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top,&top));
2259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom,&bottom));
2269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left,&left));
2279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right,&right));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
230c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
231c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
232c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
233c4762a1bSJed Brown 
234c4762a1bSJed Brown       xc = x[row];
235c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
236c4762a1bSJed Brown 
237c4762a1bSJed Brown       if (i==0) { /* left side */
238c4762a1bSJed Brown         xl= left[j-ys+1];
239c4762a1bSJed Brown         xlt = left[j-ys+2];
240c4762a1bSJed Brown       } else {
241c4762a1bSJed Brown         xl = x[row-1];
242c4762a1bSJed Brown       }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown       if (j==0) { /* bottom side */
245c4762a1bSJed Brown         xb=bottom[i-xs+1];
246c4762a1bSJed Brown         xrb = bottom[i-xs+2];
247c4762a1bSJed Brown       } else {
248c4762a1bSJed Brown         xb = x[row-gxm];
249c4762a1bSJed Brown       }
250c4762a1bSJed Brown 
251c4762a1bSJed Brown       if (i+1 == gxs+gxm) { /* right side */
252c4762a1bSJed Brown         xr=right[j-ys+1];
253c4762a1bSJed Brown         xrb = right[j-ys];
254c4762a1bSJed Brown       } else {
255c4762a1bSJed Brown         xr = x[row+1];
256c4762a1bSJed Brown       }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown       if (j+1==gys+gym) { /* top side */
259c4762a1bSJed Brown         xt=top[i-xs+1];
260c4762a1bSJed Brown         xlt = top[i-xs];
261c4762a1bSJed Brown       }else {
262c4762a1bSJed Brown         xt = x[row+gxm];
263c4762a1bSJed Brown       }
264c4762a1bSJed Brown 
265c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym) {
266c4762a1bSJed Brown         xlt = x[row-1+gxm];
267c4762a1bSJed Brown       }
268c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm) {
269c4762a1bSJed Brown         xrb = x[row+1-gxm];
270c4762a1bSJed Brown       }
271c4762a1bSJed Brown 
272c4762a1bSJed Brown       d1 = (xc-xl);
273c4762a1bSJed Brown       d2 = (xc-xr);
274c4762a1bSJed Brown       d3 = (xc-xt);
275c4762a1bSJed Brown       d4 = (xc-xb);
276c4762a1bSJed Brown       d5 = (xr-xrb);
277c4762a1bSJed Brown       d6 = (xrb-xb);
278c4762a1bSJed Brown       d7 = (xlt-xl);
279c4762a1bSJed Brown       d8 = (xt-xlt);
280c4762a1bSJed Brown 
281c4762a1bSJed Brown       df1dxc = d1*hydhx;
282c4762a1bSJed Brown       df2dxc = (d1*hydhx + d4*hxdhy);
283c4762a1bSJed Brown       df3dxc = d3*hxdhy;
284c4762a1bSJed Brown       df4dxc = (d2*hydhx + d3*hxdhy);
285c4762a1bSJed Brown       df5dxc = d2*hydhx;
286c4762a1bSJed Brown       df6dxc = d4*hxdhy;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown       d1 *= rhx;
289c4762a1bSJed Brown       d2 *= rhx;
290c4762a1bSJed Brown       d3 *= rhy;
291c4762a1bSJed Brown       d4 *= rhy;
292c4762a1bSJed Brown       d5 *= rhy;
293c4762a1bSJed Brown       d6 *= rhx;
294c4762a1bSJed Brown       d7 *= rhy;
295c4762a1bSJed Brown       d8 *= rhx;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
298c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
299c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
300c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
301c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
302c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
303c4762a1bSJed Brown 
304c4762a1bSJed Brown       ft = ft + (f2 + f4);
305c4762a1bSJed Brown 
306c4762a1bSJed Brown       df1dxc /= f1;
307c4762a1bSJed Brown       df2dxc /= f2;
308c4762a1bSJed Brown       df3dxc /= f3;
309c4762a1bSJed Brown       df4dxc /= f4;
310c4762a1bSJed Brown       df5dxc /= f5;
311c4762a1bSJed Brown       df6dxc /= f6;
312c4762a1bSJed Brown 
313c4762a1bSJed Brown       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5;
314c4762a1bSJed Brown 
315c4762a1bSJed Brown     }
316c4762a1bSJed Brown   }
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
319c4762a1bSJed Brown   if (xs==0) { /* left side */
320c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
321c4762a1bSJed Brown       d3=(left[j-ys+1] - left[j-ys+2])*rhy;
322c4762a1bSJed Brown       d2=(left[j-ys+1] - x[(j-gys)*gxm])*rhx;
323c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
324c4762a1bSJed Brown     }
325c4762a1bSJed Brown   }
326c4762a1bSJed Brown   if (ys==0) { /* bottom side */
327c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
328c4762a1bSJed Brown       d2=(bottom[i+1-xs]-bottom[i-xs+2])*rhx;
329c4762a1bSJed Brown       d3=(bottom[i-xs+1]-x[i-gxs])*rhy;
330c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
331c4762a1bSJed Brown     }
332c4762a1bSJed Brown   }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   if (xs+xm==mx) { /* right side */
335c4762a1bSJed Brown     for (j=ys; j< ys+ym; j++) {
336c4762a1bSJed Brown       d1=(x[(j+1-gys)*gxm-1]-right[j-ys+1])*rhx;
337c4762a1bSJed Brown       d4=(right[j-ys]-right[j-ys+1])*rhy;
338c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
339c4762a1bSJed Brown     }
340c4762a1bSJed Brown   }
341c4762a1bSJed Brown   if (ys+ym==my) { /* top side */
342c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
343c4762a1bSJed Brown       d1=(x[(gym-1)*gxm + i-gxs] - top[i-xs+1])*rhy;
344c4762a1bSJed Brown       d4=(top[i-xs+1] - top[i-xs])*rhx;
345c4762a1bSJed Brown       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
346c4762a1bSJed Brown     }
347c4762a1bSJed Brown   }
348c4762a1bSJed Brown 
349c4762a1bSJed Brown   if (ys==0 && xs==0) {
350c4762a1bSJed Brown     d1=(left[0]-left[1])*rhy;
351c4762a1bSJed Brown     d2=(bottom[0]-bottom[1])*rhx;
352c4762a1bSJed Brown     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
353c4762a1bSJed Brown   }
354c4762a1bSJed Brown   if (ys+ym == my && xs+xm == mx) {
355c4762a1bSJed Brown     d1=(right[ym+1] - right[ym])*rhy;
356c4762a1bSJed Brown     d2=(top[xm+1] - top[xm])*rhx;
357c4762a1bSJed Brown     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
358c4762a1bSJed Brown   }
359c4762a1bSJed Brown 
360c4762a1bSJed Brown   ft=ft*area;
3619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD));
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   /* Restore vectors */
3649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX,&x));
3659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localG,&g));
3669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left,&left));
3679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top,&top));
3689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom,&bottom));
3699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right,&right));
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   /* Scatter values to global vector */
3729566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G));
3739566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G));
374c4762a1bSJed Brown 
3759566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(70.0*xm*ym));
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   return 0;
378c4762a1bSJed Brown }
379c4762a1bSJed Brown 
380c4762a1bSJed Brown /* ------------------------------------------------------------------- */
381c4762a1bSJed Brown /*
382c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
383c4762a1bSJed Brown 
384c4762a1bSJed Brown    Input Parameters:
385c4762a1bSJed Brown .  tao  - the Tao context
386c4762a1bSJed Brown .  x    - input vector
387a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
388c4762a1bSJed Brown 
389c4762a1bSJed Brown    Output Parameters:
390c4762a1bSJed Brown .  A    - Hessian matrix
391c4762a1bSJed Brown .  B    - optionally different preconditioning matrix
392c4762a1bSJed Brown 
393c4762a1bSJed Brown    Notes:
394c4762a1bSJed Brown    Due to mesh point reordering with DMs, we must always work
395c4762a1bSJed Brown    with the local mesh points, and then transform them to the new
396c4762a1bSJed Brown    global numbering with the local-to-global mapping.  We cannot work
397c4762a1bSJed Brown    directly with the global numbers for the original uniprocessor mesh!
398c4762a1bSJed Brown 
399c4762a1bSJed Brown    Two methods are available for imposing this transformation
400c4762a1bSJed Brown    when setting matrix entries:
401c4762a1bSJed Brown      (A) MatSetValuesLocal(), using the local ordering (including
402c4762a1bSJed Brown          ghost points!)
403c4762a1bSJed Brown          - Do the following two steps once, before calling TaoSolve()
404c4762a1bSJed Brown            - Use DMGetISLocalToGlobalMapping() to extract the
405c4762a1bSJed Brown              local-to-global map from the DM
406c4762a1bSJed Brown            - Associate this map with the matrix by calling
407c4762a1bSJed Brown              MatSetLocalToGlobalMapping()
408c4762a1bSJed Brown          - Then set matrix entries using the local ordering
409c4762a1bSJed Brown            by calling MatSetValuesLocal()
410c4762a1bSJed Brown      (B) MatSetValues(), using the global ordering
411c4762a1bSJed Brown          - Use DMGetGlobalIndices() to extract the local-to-global map
412c4762a1bSJed Brown          - Then apply this map explicitly yourself
413c4762a1bSJed Brown          - Set matrix entries using the global ordering by calling
414c4762a1bSJed Brown            MatSetValues()
415c4762a1bSJed Brown    Option (A) seems cleaner/easier in many cases, and is the procedure
416c4762a1bSJed Brown    used in this example.
417c4762a1bSJed Brown */
418c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat Hptr, Mat Hessian, void *ptr)
419c4762a1bSJed Brown {
420c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
421c4762a1bSJed Brown   PetscInt       i,j,k,row;
422c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
423c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym,col[7];
424c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
425c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
426c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
427c4762a1bSJed Brown   PetscReal      hl,hr,ht,hb,hc,htl,hbr;
428c4762a1bSJed Brown   PetscReal      *x,*left,*right,*bottom,*top;
429c4762a1bSJed Brown   PetscReal      v[7];
430c4762a1bSJed Brown   Vec            localX = user->localX;
431c4762a1bSJed Brown   PetscBool      assembled;
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   /* Set various matrix options */
4349566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /* Initialize matrix entries to zero */
4379566063dSJacob Faibussowitsch   PetscCall(MatAssembled(Hessian,&assembled));
4389566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(Hessian));
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   /* Get local mesh boundaries */
4419566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
4429566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   /* Scatter ghost points to local vector */
4459566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
4469566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
447c4762a1bSJed Brown 
448c4762a1bSJed Brown   /* Get pointers to vector data */
4499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX,&x));
4509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Top,&top));
4519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Bottom,&bottom));
4529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Left,&left));
4539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->Right,&right));
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++) {
458c4762a1bSJed Brown 
459c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
460c4762a1bSJed Brown 
461c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
462c4762a1bSJed Brown 
463c4762a1bSJed Brown       xc = x[row];
464c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
465c4762a1bSJed Brown 
466c4762a1bSJed Brown       /* Left side */
467c4762a1bSJed Brown       if (i==gxs) {
468c4762a1bSJed Brown         xl= left[j-ys+1];
469c4762a1bSJed Brown         xlt = left[j-ys+2];
470c4762a1bSJed Brown       } else {
471c4762a1bSJed Brown         xl = x[row-1];
472c4762a1bSJed Brown       }
473c4762a1bSJed Brown 
474c4762a1bSJed Brown       if (j==gys) {
475c4762a1bSJed Brown         xb=bottom[i-xs+1];
476c4762a1bSJed Brown         xrb = bottom[i-xs+2];
477c4762a1bSJed Brown       } else {
478c4762a1bSJed Brown         xb = x[row-gxm];
479c4762a1bSJed Brown       }
480c4762a1bSJed Brown 
481c4762a1bSJed Brown       if (i+1 == gxs+gxm) {
482c4762a1bSJed Brown         xr=right[j-ys+1];
483c4762a1bSJed Brown         xrb = right[j-ys];
484c4762a1bSJed Brown       } else {
485c4762a1bSJed Brown         xr = x[row+1];
486c4762a1bSJed Brown       }
487c4762a1bSJed Brown 
488c4762a1bSJed Brown       if (j+1==gys+gym) {
489c4762a1bSJed Brown         xt=top[i-xs+1];
490c4762a1bSJed Brown         xlt = top[i-xs];
491c4762a1bSJed Brown       }else {
492c4762a1bSJed Brown         xt = x[row+gxm];
493c4762a1bSJed Brown       }
494c4762a1bSJed Brown 
495c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym) {
496c4762a1bSJed Brown         xlt = x[row-1+gxm];
497c4762a1bSJed Brown       }
498c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm) {
499c4762a1bSJed Brown         xrb = x[row+1-gxm];
500c4762a1bSJed Brown       }
501c4762a1bSJed Brown 
502c4762a1bSJed Brown       d1 = (xc-xl)*rhx;
503c4762a1bSJed Brown       d2 = (xc-xr)*rhx;
504c4762a1bSJed Brown       d3 = (xc-xt)*rhy;
505c4762a1bSJed Brown       d4 = (xc-xb)*rhy;
506c4762a1bSJed Brown       d5 = (xrb-xr)*rhy;
507c4762a1bSJed Brown       d6 = (xrb-xb)*rhx;
508c4762a1bSJed Brown       d7 = (xlt-xl)*rhy;
509c4762a1bSJed Brown       d8 = (xlt-xt)*rhx;
510c4762a1bSJed Brown 
511c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
512c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
513c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
514c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
515c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
516c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
517c4762a1bSJed Brown 
518c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
519c4762a1bSJed Brown         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
520c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
521c4762a1bSJed Brown         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
522c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
523c4762a1bSJed Brown         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
524c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
525c4762a1bSJed Brown         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
526c4762a1bSJed Brown 
527c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
528c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
529c4762a1bSJed Brown 
530c4762a1bSJed Brown       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
531c4762a1bSJed Brown         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
532c4762a1bSJed Brown         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
533c4762a1bSJed Brown         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
534c4762a1bSJed Brown 
535c4762a1bSJed Brown       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;
536c4762a1bSJed Brown 
537c4762a1bSJed Brown       k=0;
538c4762a1bSJed Brown       if (j>0) {
539c4762a1bSJed Brown         v[k]=hb; col[k]=row - gxm; k++;
540c4762a1bSJed Brown       }
541c4762a1bSJed Brown 
542c4762a1bSJed Brown       if (j>0 && i < mx -1) {
543c4762a1bSJed Brown         v[k]=hbr; col[k]=row - gxm+1; k++;
544c4762a1bSJed Brown       }
545c4762a1bSJed Brown 
546c4762a1bSJed Brown       if (i>0) {
547c4762a1bSJed Brown         v[k]= hl; col[k]=row - 1; k++;
548c4762a1bSJed Brown       }
549c4762a1bSJed Brown 
550c4762a1bSJed Brown       v[k]= hc; col[k]=row; k++;
551c4762a1bSJed Brown 
552c4762a1bSJed Brown       if (i < mx-1) {
553c4762a1bSJed Brown         v[k]= hr; col[k]=row+1; k++;
554c4762a1bSJed Brown       }
555c4762a1bSJed Brown 
556c4762a1bSJed Brown       if (i>0 && j < my-1) {
557c4762a1bSJed Brown         v[k]= htl; col[k] = row+gxm-1; k++;
558c4762a1bSJed Brown       }
559c4762a1bSJed Brown 
560c4762a1bSJed Brown       if (j < my-1) {
561c4762a1bSJed Brown         v[k]= ht; col[k] = row+gxm; k++;
562c4762a1bSJed Brown       }
563c4762a1bSJed Brown 
564c4762a1bSJed Brown       /*
565c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
566c4762a1bSJed Brown          earlier, in the main routine.
567c4762a1bSJed Brown       */
5689566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(Hessian,1,&row,k,col,v,INSERT_VALUES));
569c4762a1bSJed Brown 
570c4762a1bSJed Brown     }
571c4762a1bSJed Brown   }
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   /* Restore vectors */
5749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX,&x));
5759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Left,&left));
5769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Top,&top));
5779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Bottom,&bottom));
5789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->Right,&right));
579c4762a1bSJed Brown 
580c4762a1bSJed Brown   /* Assemble the matrix */
5819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
5829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
583c4762a1bSJed Brown 
5849566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0*xm*ym));
585c4762a1bSJed Brown   return 0;
586c4762a1bSJed Brown }
587c4762a1bSJed Brown 
588c4762a1bSJed Brown /* ------------------------------------------------------------------- */
589c4762a1bSJed Brown /*
590c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
591c4762a1bSJed Brown    the region.
592c4762a1bSJed Brown 
593c4762a1bSJed Brown    Input Parameter:
594c4762a1bSJed Brown .  user - user-defined application context
595c4762a1bSJed Brown 
596c4762a1bSJed Brown    Output Parameter:
597c4762a1bSJed Brown .  user - user-defined application context
598c4762a1bSJed Brown */
599c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
600c4762a1bSJed Brown {
601c4762a1bSJed Brown   PetscInt   i,j,k,maxits=5,limit=0;
602c4762a1bSJed Brown   PetscInt   xs,ys,xm,ym,gxs,gys,gxm,gym;
603c4762a1bSJed Brown   PetscInt   mx=user->mx,my=user->my;
604c4762a1bSJed Brown   PetscInt   bsize=0, lsize=0, tsize=0, rsize=0;
605c4762a1bSJed Brown   PetscReal  one=1.0, two=2.0, three=3.0, scl=1.0, tol=1e-10;
606c4762a1bSJed Brown   PetscReal  fnorm,det,hx,hy,xt=0,yt=0;
607c4762a1bSJed Brown   PetscReal  u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
608c4762a1bSJed Brown   PetscReal  b=-0.5, t=0.5, l=-0.5, r=0.5;
609c4762a1bSJed Brown   PetscReal  *boundary;
610c4762a1bSJed Brown   PetscBool  flg;
611c4762a1bSJed Brown   Vec        Bottom,Top,Right,Left;
612c4762a1bSJed Brown 
613c4762a1bSJed Brown   /* Get local mesh boundaries */
6149566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
6159566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   bsize=xm+2;
618c4762a1bSJed Brown   lsize=ym+2;
619c4762a1bSJed Brown   rsize=ym+2;
620c4762a1bSJed Brown   tsize=xm+2;
621c4762a1bSJed Brown 
6229566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom));
6239566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top));
6249566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left));
6259566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right));
626c4762a1bSJed Brown 
627c4762a1bSJed Brown   user->Top=Top;
628c4762a1bSJed Brown   user->Left=Left;
629c4762a1bSJed Brown   user->Bottom=Bottom;
630c4762a1bSJed Brown   user->Right=Right;
631c4762a1bSJed Brown 
632c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
633c4762a1bSJed Brown 
634c4762a1bSJed Brown   for (j=0; j<4; j++) {
635c4762a1bSJed Brown     if (j==0) {
636c4762a1bSJed Brown       yt=b;
637c4762a1bSJed Brown       xt=l+hx*xs;
638c4762a1bSJed Brown       limit=bsize;
639c4762a1bSJed Brown       VecGetArray(Bottom,&boundary);
640c4762a1bSJed Brown     } else if (j==1) {
641c4762a1bSJed Brown       yt=t;
642c4762a1bSJed Brown       xt=l+hx*xs;
643c4762a1bSJed Brown       limit=tsize;
644c4762a1bSJed Brown       VecGetArray(Top,&boundary);
645c4762a1bSJed Brown     } else if (j==2) {
646c4762a1bSJed Brown       yt=b+hy*ys;
647c4762a1bSJed Brown       xt=l;
648c4762a1bSJed Brown       limit=lsize;
649c4762a1bSJed Brown       VecGetArray(Left,&boundary);
650c4762a1bSJed Brown     } else if (j==3) {
651c4762a1bSJed Brown       yt=b+hy*ys;
652c4762a1bSJed Brown       xt=r;
653c4762a1bSJed Brown       limit=rsize;
654c4762a1bSJed Brown       VecGetArray(Right,&boundary);
655c4762a1bSJed Brown     }
656c4762a1bSJed Brown 
657c4762a1bSJed Brown     for (i=0; i<limit; i++) {
658c4762a1bSJed Brown       u1=xt;
659c4762a1bSJed Brown       u2=-yt;
660c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
661c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
662c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
663c4762a1bSJed Brown         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
664c4762a1bSJed Brown         if (fnorm <= tol) break;
665c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
666c4762a1bSJed Brown         njac12=two*u1*u2;
667c4762a1bSJed Brown         njac21=-two*u1*u2;
668c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
669c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
670c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
671c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
672c4762a1bSJed Brown       }
673c4762a1bSJed Brown 
674c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
675c4762a1bSJed Brown       if (j==0 || j==1) {
676c4762a1bSJed Brown         xt=xt+hx;
677c4762a1bSJed Brown       } else if (j==2 || j==3) {
678c4762a1bSJed Brown         yt=yt+hy;
679c4762a1bSJed Brown       }
680c4762a1bSJed Brown     }
681c4762a1bSJed Brown     if (j==0) {
6829566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Bottom,&boundary));
683c4762a1bSJed Brown     } else if (j==1) {
6849566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Top,&boundary));
685c4762a1bSJed Brown     } else if (j==2) {
6869566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Left,&boundary));
687c4762a1bSJed Brown     } else if (j==3) {
6889566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(Right,&boundary));
689c4762a1bSJed Brown     }
690c4762a1bSJed Brown   }
691c4762a1bSJed Brown 
692c4762a1bSJed Brown   /* Scale the boundary if desired */
693c4762a1bSJed Brown 
6949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg));
695c4762a1bSJed Brown   if (flg) {
6969566063dSJacob Faibussowitsch     PetscCall(VecScale(Bottom, scl));
697c4762a1bSJed Brown   }
6989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg));
699c4762a1bSJed Brown   if (flg) {
7009566063dSJacob Faibussowitsch     PetscCall(VecScale(Top, scl));
701c4762a1bSJed Brown   }
7029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg));
703c4762a1bSJed Brown   if (flg) {
7049566063dSJacob Faibussowitsch     PetscCall(VecScale(Right, scl));
705c4762a1bSJed Brown   }
706c4762a1bSJed Brown 
7079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg));
708c4762a1bSJed Brown   if (flg) {
7099566063dSJacob Faibussowitsch     PetscCall(VecScale(Left, scl));
710c4762a1bSJed Brown   }
711c4762a1bSJed Brown   return 0;
712c4762a1bSJed Brown }
713c4762a1bSJed Brown 
714c4762a1bSJed Brown /* ------------------------------------------------------------------- */
715c4762a1bSJed Brown /*
716c4762a1bSJed Brown    MSA_Plate -  Calculates an obstacle for surface to stretch over.
717c4762a1bSJed Brown 
718c4762a1bSJed Brown    Input Parameter:
719c4762a1bSJed Brown .  user - user-defined application context
720c4762a1bSJed Brown 
721c4762a1bSJed Brown    Output Parameter:
722c4762a1bSJed Brown .  user - user-defined application context
723c4762a1bSJed Brown */
72470a7d78aSStefano Zampini static PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx)
72570a7d78aSStefano Zampini {
726c4762a1bSJed Brown   AppCtx         *user=(AppCtx *)ctx;
727c4762a1bSJed Brown   PetscInt       i,j,row;
728c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym;
729c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my, bmy, bmx;
730c4762a1bSJed Brown   PetscReal      t1,t2,t3;
731c4762a1bSJed Brown   PetscReal      *xl, lb=PETSC_NINFINITY, ub=PETSC_INFINITY;
732c4762a1bSJed Brown   PetscBool      cylinder;
733c4762a1bSJed Brown 
734c4762a1bSJed Brown   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
735c4762a1bSJed Brown   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
736c4762a1bSJed Brown   bmy=user->bmy; bmx=user->bmx;
737c4762a1bSJed Brown 
7389566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
739c4762a1bSJed Brown 
7409566063dSJacob Faibussowitsch   PetscCall(VecSet(XL, lb));
7419566063dSJacob Faibussowitsch   PetscCall(VecSet(XU, ub));
742c4762a1bSJed Brown 
7439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(XL,&xl));
744c4762a1bSJed Brown 
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-cylinder",&cylinder));
746c4762a1bSJed Brown   /* Compute the optional lower box */
747c4762a1bSJed Brown   if (cylinder) {
748c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
749c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
750c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
751c4762a1bSJed Brown         t1=(2.0*i-mx)*bmy;
752c4762a1bSJed Brown         t2=(2.0*j-my)*bmx;
753c4762a1bSJed Brown         t3=bmx*bmx*bmy*bmy;
754c4762a1bSJed Brown         if (t1*t1 + t2*t2 <= t3) {
755c4762a1bSJed Brown           xl[row] = user->bheight;
756c4762a1bSJed Brown         }
757c4762a1bSJed Brown       }
758c4762a1bSJed Brown     }
759c4762a1bSJed Brown   } else {
760c4762a1bSJed Brown     /* Compute the optional lower box */
761c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++) {
762c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
763c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
764c4762a1bSJed Brown         if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
765c4762a1bSJed Brown             j>=(my-bmy)/2 && j<my-(my-bmy)/2) {
766c4762a1bSJed Brown           xl[row] = user->bheight;
767c4762a1bSJed Brown         }
768c4762a1bSJed Brown       }
769c4762a1bSJed Brown     }
770c4762a1bSJed Brown   }
7719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(XL,&xl));
772c4762a1bSJed Brown 
773c4762a1bSJed Brown   return 0;
774c4762a1bSJed Brown }
775c4762a1bSJed Brown 
776c4762a1bSJed Brown /* ------------------------------------------------------------------- */
777c4762a1bSJed Brown /*
778c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
779c4762a1bSJed Brown 
780c4762a1bSJed Brown    Input Parameters:
781c4762a1bSJed Brown .  user - user-defined application context
782c4762a1bSJed Brown .  X - vector for initial guess
783c4762a1bSJed Brown 
784c4762a1bSJed Brown    Output Parameters:
785c4762a1bSJed Brown .  X - newly computed initial guess
786c4762a1bSJed Brown */
787c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
788c4762a1bSJed Brown {
789c4762a1bSJed Brown   PetscInt       start=-1,i,j;
790c4762a1bSJed Brown   PetscReal      zero=0.0;
791c4762a1bSJed Brown   PetscBool      flg;
792c4762a1bSJed Brown 
7939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
794c4762a1bSJed Brown   if (flg && start==0) { /* The zero vector is reasonable */
7959566063dSJacob Faibussowitsch     PetscCall(VecSet(X, zero));
796c4762a1bSJed Brown   } else if (flg && start>0) { /* Try a random start between -0.5 and 0.5 */
797c4762a1bSJed Brown     PetscRandom rctx;  PetscReal np5=-0.5;
798c4762a1bSJed Brown 
7999566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(MPI_COMM_WORLD,&rctx));
800c4762a1bSJed Brown     for (i=0; i<start; i++) {
8019566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(X, rctx));
802c4762a1bSJed Brown     }
8039566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rctx));
8049566063dSJacob Faibussowitsch     PetscCall(VecShift(X, np5));
805c4762a1bSJed Brown 
806c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
807c4762a1bSJed Brown 
808c4762a1bSJed Brown     PetscInt row,xs,xm,gxs,gxm,ys,ym,gys,gym;
809c4762a1bSJed Brown     PetscInt mx=user->mx,my=user->my;
810c4762a1bSJed Brown     PetscReal *x,*left,*right,*bottom,*top;
811c4762a1bSJed Brown     Vec    localX = user->localX;
812c4762a1bSJed Brown 
813c4762a1bSJed Brown     /* Get local mesh boundaries */
8149566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
8159566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
816c4762a1bSJed Brown 
817c4762a1bSJed Brown     /* Get pointers to vector data */
8189566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Top,&top));
8199566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Bottom,&bottom));
8209566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Left,&left));
8219566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->Right,&right));
822c4762a1bSJed Brown 
8239566063dSJacob Faibussowitsch     PetscCall(VecGetArray(localX,&x));
824c4762a1bSJed Brown     /* Perform local computations */
825c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
826c4762a1bSJed Brown       for (i=xs; i< xs+xm; i++) {
827c4762a1bSJed Brown         row=(j-gys)*gxm + (i-gxs);
8282da392ccSBarry 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;
829c4762a1bSJed Brown       }
830c4762a1bSJed Brown     }
831c4762a1bSJed Brown 
832c4762a1bSJed Brown     /* Restore vectors */
8339566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(localX,&x));
834c4762a1bSJed Brown 
8359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Left,&left));
8369566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Top,&top));
8379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Bottom,&bottom));
8389566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->Right,&right));
839c4762a1bSJed Brown 
840c4762a1bSJed Brown     /* Scatter values into global vector */
8419566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X));
8429566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X));
843c4762a1bSJed Brown 
844c4762a1bSJed Brown   }
845c4762a1bSJed Brown   return 0;
846c4762a1bSJed Brown }
847c4762a1bSJed Brown 
848c4762a1bSJed Brown /* For testing matrix free submatrices */
849c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
850c4762a1bSJed Brown {
851c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
852c4762a1bSJed Brown   PetscFunctionBegin;
8539566063dSJacob Faibussowitsch   PetscCall(FormHessian(tao,x,user->H,user->H,ptr));
854c4762a1bSJed Brown   PetscFunctionReturn(0);
855c4762a1bSJed Brown }
856c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
857c4762a1bSJed Brown {
858c4762a1bSJed Brown   void           *ptr;
859c4762a1bSJed Brown   AppCtx         *user;
860c4762a1bSJed Brown   PetscFunctionBegin;
8619566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(H_shell,&ptr));
862c4762a1bSJed Brown   user = (AppCtx*)ptr;
8639566063dSJacob Faibussowitsch   PetscCall(MatMult(user->H,X,Y));
864c4762a1bSJed Brown   PetscFunctionReturn(0);
865c4762a1bSJed Brown }
866c4762a1bSJed Brown 
867c4762a1bSJed Brown /*TEST
868c4762a1bSJed Brown 
869c4762a1bSJed Brown    build:
870c4762a1bSJed Brown       requires: !complex
871c4762a1bSJed Brown 
872c4762a1bSJed Brown    test:
873c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
874c4762a1bSJed Brown       requires: !single
875c4762a1bSJed Brown 
876c4762a1bSJed Brown    test:
877c4762a1bSJed Brown       suffix: 2
878c4762a1bSJed Brown       nsize: 2
879c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
880c4762a1bSJed Brown       requires: !single
881c4762a1bSJed Brown 
882c4762a1bSJed Brown    test:
883c4762a1bSJed Brown       suffix: 3
884c4762a1bSJed Brown       nsize: 3
885c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
886c4762a1bSJed Brown       requires: !single
887c4762a1bSJed Brown 
888c4762a1bSJed Brown    test:
889c4762a1bSJed Brown       suffix: 4
890c4762a1bSJed Brown       nsize: 3
891c4762a1bSJed 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
892c4762a1bSJed Brown       requires: !single
893c4762a1bSJed Brown 
894c4762a1bSJed Brown    test:
895c4762a1bSJed Brown       suffix: 5
896c4762a1bSJed Brown       nsize: 3
897c4762a1bSJed 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
898c4762a1bSJed Brown       requires: !single
899c4762a1bSJed Brown 
900c4762a1bSJed Brown    test:
901c4762a1bSJed Brown       suffix: 6
902c4762a1bSJed Brown       nsize: 3
903c4762a1bSJed 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
904c4762a1bSJed Brown       requires: !single
905c4762a1bSJed Brown 
906c4762a1bSJed Brown    test:
907c4762a1bSJed Brown       suffix: 7
908c4762a1bSJed Brown       nsize: 3
909c4762a1bSJed 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
910c4762a1bSJed Brown       requires: !single
911c4762a1bSJed Brown 
912c4762a1bSJed Brown    test:
913c4762a1bSJed Brown       suffix: 8
914c4762a1bSJed 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
915c4762a1bSJed Brown       requires: !single
916c4762a1bSJed Brown 
917c4762a1bSJed Brown    test:
918c4762a1bSJed Brown       suffix: 9
919c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
920c4762a1bSJed Brown       requires: !single
921c4762a1bSJed Brown 
922c4762a1bSJed Brown    test:
923c4762a1bSJed Brown       suffix: 10
924c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
925c4762a1bSJed Brown       requires: !single
926c4762a1bSJed Brown 
927c4762a1bSJed Brown    test:
928c4762a1bSJed Brown       suffix: 11
929c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
930c4762a1bSJed Brown       requires: !single
931c4762a1bSJed Brown 
932c4762a1bSJed Brown    test:
933c4762a1bSJed Brown       suffix: 12
934c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
935c4762a1bSJed Brown       requires: !single
936c4762a1bSJed Brown 
937c4762a1bSJed Brown    test:
938c4762a1bSJed Brown       suffix: 13
939c4762a1bSJed 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
940c4762a1bSJed Brown       requires: !single
941c4762a1bSJed Brown 
942c4762a1bSJed Brown    test:
943c4762a1bSJed Brown       suffix: 14
944c4762a1bSJed 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
945c4762a1bSJed Brown       requires: !single
946c4762a1bSJed Brown 
947c4762a1bSJed Brown    test:
948c4762a1bSJed Brown       suffix: 15
949c4762a1bSJed 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
950c4762a1bSJed Brown       requires: !single
951c4762a1bSJed Brown 
952c4762a1bSJed Brown    test:
953c4762a1bSJed Brown      suffix: 16
954c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
955c4762a1bSJed Brown      requires: !single
956c4762a1bSJed Brown 
957c4762a1bSJed Brown    test:
958c4762a1bSJed Brown      suffix: 17
959c4762a1bSJed 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
960c4762a1bSJed Brown      requires: !single
961c4762a1bSJed Brown 
96234ad9904SAlp Dener    test:
96334ad9904SAlp Dener      suffix: 18
96434ad9904SAlp 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
96534ad9904SAlp Dener      requires: !single
96634ad9904SAlp Dener 
96734ad9904SAlp Dener    test:
96834ad9904SAlp Dener      suffix: 19
96934ad9904SAlp 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
97034ad9904SAlp Dener      requires: !single
97134ad9904SAlp Dener 
97234ad9904SAlp Dener    test:
97334ad9904SAlp Dener      suffix: 20
97434ad9904SAlp 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
97534ad9904SAlp Dener 
976c4762a1bSJed Brown TEST*/
977