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