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