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