xref: /petsc/src/tao/unconstrained/tutorials/minsurf1.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 minsurf1 [-help] [all TAO options] */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*  Include "petsctao.h" so we can use TAO solvers.  */
4c4762a1bSJed Brown #include <petsctao.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown static char  help[] =
7c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\
8c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
9c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and \n\
10c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
11c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
12c4762a1bSJed Brown The command line options are:\n\
13c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
14c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
15c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /*T
18c4762a1bSJed Brown    Concepts: TAO^Solving an unconstrained minimization problem
19c4762a1bSJed Brown    Routines: TaoCreate(); TaoSetType();
20a82e8c82SStefano Zampini    Routines: TaoSetSolution();
21a82e8c82SStefano Zampini    Routines: TaoSetObjectiveAndGradient();
22a82e8c82SStefano Zampini    Routines: TaoSetHessian(); TaoSetFromOptions();
23c4762a1bSJed Brown    Routines: TaoGetKSP(); TaoSolve();
24c4762a1bSJed Brown    Routines: TaoDestroy();
25c4762a1bSJed Brown    Processors: 1
26c4762a1bSJed Brown T*/
27c4762a1bSJed Brown 
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown    User-defined application context - contains data needed by the
30c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient()
31c4762a1bSJed Brown    and FormHessian().
32c4762a1bSJed Brown */
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown   PetscInt    mx, my;                 /* discretization in x, y directions */
35c4762a1bSJed Brown   PetscReal   *bottom, *top, *left, *right;             /* boundary values */
36c4762a1bSJed Brown   Mat         H;
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /* -------- User-defined Routines --------- */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
42c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
43c4762a1bSJed Brown static PetscErrorCode QuadraticH(AppCtx*,Vec,Mat);
44c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
45c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
46c4762a1bSJed Brown 
47c4762a1bSJed Brown int main(int argc, char **argv)
48c4762a1bSJed Brown {
49c4762a1bSJed Brown   PetscErrorCode     ierr;              /* used to check for functions returning nonzeros */
50c4762a1bSJed Brown   PetscInt           N;                 /* Size of vector */
51c4762a1bSJed Brown   PetscMPIInt        size;              /* Number of processors */
52c4762a1bSJed Brown   Vec                x;                 /* solution, gradient vectors */
53c4762a1bSJed Brown   KSP                ksp;               /*  PETSc Krylov subspace method */
54c4762a1bSJed Brown   PetscBool          flg;               /* A return value when checking for user options */
55c4762a1bSJed Brown   Tao                tao;               /* Tao solver context */
56c4762a1bSJed Brown   AppCtx             user;              /* user-defined work context */
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Initialize TAO,PETSc */
59c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv,(char *)0,help);if (ierr) return ierr;
60c4762a1bSJed Brown 
61*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(MPI_COMM_WORLD,&size));
623c859ba3SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors");
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Specify default dimension of the problem */
65c4762a1bSJed Brown   user.mx = 4; user.my = 4;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
70c4762a1bSJed Brown 
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n"));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"mx: %D     my: %D   \n\n",user.mx,user.my));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Calculate any derived values from parameters */
75c4762a1bSJed Brown   N    = user.mx*user.my;
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Create TAO solver and set desired solution method  */
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao));
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOLMVM));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Initialize minsurf application data structure for use in the function evaluations  */
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MSA_BoundaryConditions(&user));            /* Application specific routine */
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /*
85c4762a1bSJed Brown      Create a vector to hold the variables.  Compute an initial solution.
86c4762a1bSJed Brown      Set this vector, which will be used by TAO.
87c4762a1bSJed Brown   */
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N,&x));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MSA_InitialPoint(&user,x));                /* Application specific routine */
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,x));   /* A TAO routine                */
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Provide TAO routines for function, gradient, and Hessian evaluation */
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Create a matrix data structure to store the Hessian.  This structure will be used by TAO */
96*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,7,NULL,&(user.H)));
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* Check for any TAO command line options */
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Limit the number of iterations in the KSP linear solver */
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
105c4762a1bSJed Brown   if (ksp) {
106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,user.mx*user.my));
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
111c4762a1bSJed Brown 
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.H));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.bottom));
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.top));
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.left));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.right));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   ierr = PetscFinalize();
121c4762a1bSJed Brown   return ierr;
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown /* -------------------------------------------------------------------- */
125c4762a1bSJed Brown 
126c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates function and corresponding gradient.
127c4762a1bSJed Brown 
128c4762a1bSJed Brown     Input Parameters:
129c4762a1bSJed Brown .   tao     - the Tao context
130c4762a1bSJed Brown .   X       - input vector
131c4762a1bSJed Brown .   userCtx - optional user-defined context, as set by TaoSetFunctionGradient()
132c4762a1bSJed Brown 
133c4762a1bSJed Brown     Output Parameters:
134c4762a1bSJed Brown .   fcn     - the newly evaluated function
135c4762a1bSJed Brown .   G       - vector containing the newly evaluated gradient
136c4762a1bSJed Brown */
137c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *fcn,Vec G,void *userCtx)
138c4762a1bSJed Brown {
139c4762a1bSJed Brown   AppCtx            *user = (AppCtx *) userCtx;
140c4762a1bSJed Brown   PetscInt          i,j,row;
141c4762a1bSJed Brown   PetscInt          mx=user->mx, my=user->my;
142c4762a1bSJed Brown   PetscReal         rhx=mx+1, rhy=my+1;
143c4762a1bSJed Brown   PetscReal         hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy, ft=0;
144c4762a1bSJed Brown   PetscReal         f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
145c4762a1bSJed Brown   PetscReal         df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
146c4762a1bSJed Brown   PetscReal         zero=0.0;
147c4762a1bSJed Brown   PetscScalar       *g;
148c4762a1bSJed Brown   const PetscScalar *x;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   PetscFunctionBeginUser;
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(G, zero));
152c4762a1bSJed Brown 
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(G,&g));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
157c4762a1bSJed Brown   for (j=0; j<my; j++) {
158c4762a1bSJed Brown     for (i=0; i< mx; i++) {
159c4762a1bSJed Brown       row=(j)*mx + (i);
160c4762a1bSJed Brown       xc = x[row];
161c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
162c4762a1bSJed Brown       if (i==0) { /* left side */
163c4762a1bSJed Brown         xl  = user->left[j+1];
164c4762a1bSJed Brown         xlt = user->left[j+2];
165c4762a1bSJed Brown       } else {
166c4762a1bSJed Brown         xl = x[row-1];
167c4762a1bSJed Brown       }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown       if (j==0) { /* bottom side */
170c4762a1bSJed Brown         xb  = user->bottom[i+1];
171c4762a1bSJed Brown         xrb = user->bottom[i+2];
172c4762a1bSJed Brown       } else {
173c4762a1bSJed Brown         xb = x[row-mx];
174c4762a1bSJed Brown       }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown       if (i+1 == mx) { /* right side */
177c4762a1bSJed Brown         xr  = user->right[j+1];
178c4762a1bSJed Brown         xrb = user->right[j];
179c4762a1bSJed Brown       } else {
180c4762a1bSJed Brown         xr = x[row+1];
181c4762a1bSJed Brown       }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown       if (j+1==0+my) { /* top side */
184c4762a1bSJed Brown         xt  = user->top[i+1];
185c4762a1bSJed Brown         xlt = user->top[i];
186c4762a1bSJed Brown       }else {
187c4762a1bSJed Brown         xt = x[row+mx];
188c4762a1bSJed Brown       }
189c4762a1bSJed Brown 
190c4762a1bSJed Brown       if (i>0 && j+1<my) {
191c4762a1bSJed Brown         xlt = x[row-1+mx];
192c4762a1bSJed Brown       }
193c4762a1bSJed Brown       if (j>0 && i+1<mx) {
194c4762a1bSJed Brown         xrb = x[row+1-mx];
195c4762a1bSJed Brown       }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown       d1 = (xc-xl);
198c4762a1bSJed Brown       d2 = (xc-xr);
199c4762a1bSJed Brown       d3 = (xc-xt);
200c4762a1bSJed Brown       d4 = (xc-xb);
201c4762a1bSJed Brown       d5 = (xr-xrb);
202c4762a1bSJed Brown       d6 = (xrb-xb);
203c4762a1bSJed Brown       d7 = (xlt-xl);
204c4762a1bSJed Brown       d8 = (xt-xlt);
205c4762a1bSJed Brown 
206c4762a1bSJed Brown       df1dxc = d1*hydhx;
207c4762a1bSJed Brown       df2dxc = (d1*hydhx + d4*hxdhy);
208c4762a1bSJed Brown       df3dxc = d3*hxdhy;
209c4762a1bSJed Brown       df4dxc = (d2*hydhx + d3*hxdhy);
210c4762a1bSJed Brown       df5dxc = d2*hydhx;
211c4762a1bSJed Brown       df6dxc = d4*hxdhy;
212c4762a1bSJed Brown 
213c4762a1bSJed Brown       d1 *= rhx;
214c4762a1bSJed Brown       d2 *= rhx;
215c4762a1bSJed Brown       d3 *= rhy;
216c4762a1bSJed Brown       d4 *= rhy;
217c4762a1bSJed Brown       d5 *= rhy;
218c4762a1bSJed Brown       d6 *= rhx;
219c4762a1bSJed Brown       d7 *= rhy;
220c4762a1bSJed Brown       d8 *= rhx;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
223c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
224c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
225c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
226c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
227c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
228c4762a1bSJed Brown 
229c4762a1bSJed Brown       ft = ft + (f2 + f4);
230c4762a1bSJed Brown 
231c4762a1bSJed Brown       df1dxc /= f1;
232c4762a1bSJed Brown       df2dxc /= f2;
233c4762a1bSJed Brown       df3dxc /= f3;
234c4762a1bSJed Brown       df4dxc /= f4;
235c4762a1bSJed Brown       df5dxc /= f5;
236c4762a1bSJed Brown       df6dxc /= f6;
237c4762a1bSJed Brown 
238c4762a1bSJed Brown       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
239c4762a1bSJed Brown     }
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   for (j=0; j<my; j++) {   /* left side */
243c4762a1bSJed Brown     d3 = (user->left[j+1] - user->left[j+2])*rhy;
244c4762a1bSJed Brown     d2 = (user->left[j+1] - x[j*mx])*rhx;
245c4762a1bSJed Brown     ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   for (i=0; i<mx; i++) { /* bottom */
249c4762a1bSJed Brown     d2 = (user->bottom[i+1]-user->bottom[i+2])*rhx;
250c4762a1bSJed Brown     d3 = (user->bottom[i+1]-x[i])*rhy;
251c4762a1bSJed Brown     ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
252c4762a1bSJed Brown   }
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   for (j=0; j< my; j++) { /* right side */
255c4762a1bSJed Brown     d1 = (x[(j+1)*mx-1]-user->right[j+1])*rhx;
256c4762a1bSJed Brown     d4 = (user->right[j]-user->right[j+1])*rhy;
257c4762a1bSJed Brown     ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
258c4762a1bSJed Brown   }
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   for (i=0; i<mx; i++) { /* top side */
261c4762a1bSJed Brown     d1 = (x[(my-1)*mx + i] - user->top[i+1])*rhy;
262c4762a1bSJed Brown     d4 = (user->top[i+1] - user->top[i])*rhx;
263c4762a1bSJed Brown     ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
264c4762a1bSJed Brown   }
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* Bottom left corner */
267c4762a1bSJed Brown   d1  = (user->left[0]-user->left[1])*rhy;
268c4762a1bSJed Brown   d2  = (user->bottom[0]-user->bottom[1])*rhx;
269c4762a1bSJed Brown   ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* Top right corner */
272c4762a1bSJed Brown   d1  = (user->right[my+1] - user->right[my])*rhy;
273c4762a1bSJed Brown   d2  = (user->top[mx+1] - user->top[mx])*rhx;
274c4762a1bSJed Brown   ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   (*fcn)=ft*area;
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /* Restore vectors */
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
280*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(G,&g));
281*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(67.0*mx*my));
282c4762a1bSJed Brown   PetscFunctionReturn(0);
283c4762a1bSJed Brown }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown /* ------------------------------------------------------------------- */
286c4762a1bSJed Brown /*
287c4762a1bSJed Brown    FormHessian - Evaluates the Hessian matrix.
288c4762a1bSJed Brown 
289c4762a1bSJed Brown    Input Parameters:
290c4762a1bSJed Brown .  tao  - the Tao context
291c4762a1bSJed Brown .  x    - input vector
292c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetHessian()
293c4762a1bSJed Brown 
294c4762a1bSJed Brown    Output Parameters:
295c4762a1bSJed Brown .  H    - Hessian matrix
296c4762a1bSJed Brown .  Hpre - optionally different preconditioning matrix
297c4762a1bSJed Brown .  flg  - flag indicating matrix structure
298c4762a1bSJed Brown 
299c4762a1bSJed Brown */
300c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
301c4762a1bSJed Brown {
302c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   PetscFunctionBeginUser;
305c4762a1bSJed Brown   /* Evaluate the Hessian entries*/
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(QuadraticH(user,X,H));
307c4762a1bSJed Brown   PetscFunctionReturn(0);
308c4762a1bSJed Brown }
309c4762a1bSJed Brown 
310c4762a1bSJed Brown /* ------------------------------------------------------------------- */
311c4762a1bSJed Brown /*
312c4762a1bSJed Brown    QuadraticH - Evaluates the Hessian matrix.
313c4762a1bSJed Brown 
314c4762a1bSJed Brown    Input Parameters:
315c4762a1bSJed Brown .  user - user-defined context, as set by TaoSetHessian()
316c4762a1bSJed Brown .  X    - input vector
317c4762a1bSJed Brown 
318c4762a1bSJed Brown    Output Parameter:
319c4762a1bSJed Brown .  H    - Hessian matrix
320c4762a1bSJed Brown */
321c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian)
322c4762a1bSJed Brown {
323c4762a1bSJed Brown   PetscInt          i,j,k,row;
324c4762a1bSJed Brown   PetscInt          mx=user->mx, my=user->my;
325c4762a1bSJed Brown   PetscInt          col[7];
326c4762a1bSJed Brown   PetscReal         hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
327c4762a1bSJed Brown   PetscReal         rhx=mx+1, rhy=my+1;
328c4762a1bSJed Brown   PetscReal         f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
329c4762a1bSJed Brown   PetscReal         hl,hr,ht,hb,hc,htl,hbr;
330c4762a1bSJed Brown   const PetscScalar *x;
331c4762a1bSJed Brown   PetscReal         v[7];
332c4762a1bSJed Brown 
333362febeeSStefano Zampini   PetscFunctionBeginUser;
334c4762a1bSJed Brown   /* Get pointers to vector data */
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
336c4762a1bSJed Brown 
337c4762a1bSJed Brown   /* Initialize matrix entries to zero */
338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(Hessian));
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   /* Set various matrix options */
341*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
344c4762a1bSJed Brown   for (i=0; i< mx; i++) {
345c4762a1bSJed Brown     for (j=0; j<my; j++) {
346c4762a1bSJed Brown 
347c4762a1bSJed Brown       row=(j)*mx + (i);
348c4762a1bSJed Brown 
349c4762a1bSJed Brown       xc = x[row];
350c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
351c4762a1bSJed Brown 
352c4762a1bSJed Brown       /* Left side */
353c4762a1bSJed Brown       if (i==0) {
354c4762a1bSJed Brown         xl= user->left[j+1];
355c4762a1bSJed Brown         xlt = user->left[j+2];
356c4762a1bSJed Brown       } else {
357c4762a1bSJed Brown         xl = x[row-1];
358c4762a1bSJed Brown       }
359c4762a1bSJed Brown 
360c4762a1bSJed Brown       if (j==0) {
361c4762a1bSJed Brown         xb=user->bottom[i+1];
362c4762a1bSJed Brown         xrb = user->bottom[i+2];
363c4762a1bSJed Brown       } else {
364c4762a1bSJed Brown         xb = x[row-mx];
365c4762a1bSJed Brown       }
366c4762a1bSJed Brown 
367c4762a1bSJed Brown       if (i+1 == mx) {
368c4762a1bSJed Brown         xr=user->right[j+1];
369c4762a1bSJed Brown         xrb = user->right[j];
370c4762a1bSJed Brown       } else {
371c4762a1bSJed Brown         xr = x[row+1];
372c4762a1bSJed Brown       }
373c4762a1bSJed Brown 
374c4762a1bSJed Brown       if (j+1==my) {
375c4762a1bSJed Brown         xt=user->top[i+1];
376c4762a1bSJed Brown         xlt = user->top[i];
377c4762a1bSJed Brown       }else {
378c4762a1bSJed Brown         xt = x[row+mx];
379c4762a1bSJed Brown       }
380c4762a1bSJed Brown 
381c4762a1bSJed Brown       if (i>0 && j+1<my) {
382c4762a1bSJed Brown         xlt = x[row-1+mx];
383c4762a1bSJed Brown       }
384c4762a1bSJed Brown       if (j>0 && i+1<mx) {
385c4762a1bSJed Brown         xrb = x[row+1-mx];
386c4762a1bSJed Brown       }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown       d1 = (xc-xl)*rhx;
389c4762a1bSJed Brown       d2 = (xc-xr)*rhx;
390c4762a1bSJed Brown       d3 = (xc-xt)*rhy;
391c4762a1bSJed Brown       d4 = (xc-xb)*rhy;
392c4762a1bSJed Brown       d5 = (xrb-xr)*rhy;
393c4762a1bSJed Brown       d6 = (xrb-xb)*rhx;
394c4762a1bSJed Brown       d7 = (xlt-xl)*rhy;
395c4762a1bSJed Brown       d8 = (xlt-xt)*rhx;
396c4762a1bSJed Brown 
397c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
398c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
399c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
400c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
401c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
402c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
403c4762a1bSJed Brown 
404c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
405c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
406c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
407c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
408c4762a1bSJed Brown 
409c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
410c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
411c4762a1bSJed Brown 
412c4762a1bSJed Brown       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
413c4762a1bSJed Brown            (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +  (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
414c4762a1bSJed Brown 
415c4762a1bSJed Brown       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;
416c4762a1bSJed Brown 
417c4762a1bSJed Brown       k=0;
418c4762a1bSJed Brown       if (j>0) {
419c4762a1bSJed Brown         v[k]=hb; col[k]=row - mx; k++;
420c4762a1bSJed Brown       }
421c4762a1bSJed Brown 
422c4762a1bSJed Brown       if (j>0 && i < mx -1) {
423c4762a1bSJed Brown         v[k]=hbr; col[k]=row - mx+1; k++;
424c4762a1bSJed Brown       }
425c4762a1bSJed Brown 
426c4762a1bSJed Brown       if (i>0) {
427c4762a1bSJed Brown         v[k]= hl; col[k]=row - 1; k++;
428c4762a1bSJed Brown       }
429c4762a1bSJed Brown 
430c4762a1bSJed Brown       v[k]= hc; col[k]=row; k++;
431c4762a1bSJed Brown 
432c4762a1bSJed Brown       if (i < mx-1) {
433c4762a1bSJed Brown         v[k]= hr; col[k]=row+1; k++;
434c4762a1bSJed Brown       }
435c4762a1bSJed Brown 
436c4762a1bSJed Brown       if (i>0 && j < my-1) {
437c4762a1bSJed Brown         v[k]= htl; col[k] = row+mx-1; k++;
438c4762a1bSJed Brown       }
439c4762a1bSJed Brown 
440c4762a1bSJed Brown       if (j < my-1) {
441c4762a1bSJed Brown         v[k]= ht; col[k] = row+mx; k++;
442c4762a1bSJed Brown       }
443c4762a1bSJed Brown 
444c4762a1bSJed Brown       /*
445c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
446c4762a1bSJed Brown          earlier, in the main routine.
447c4762a1bSJed Brown       */
448*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(Hessian,1,&row,k,col,v,INSERT_VALUES));
449c4762a1bSJed Brown     }
450c4762a1bSJed Brown   }
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   /* Restore vectors */
453*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   /* Assemble the matrix */
456*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
457*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
458c4762a1bSJed Brown 
459*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(199.0*mx*my));
460c4762a1bSJed Brown   PetscFunctionReturn(0);
461c4762a1bSJed Brown }
462c4762a1bSJed Brown 
463c4762a1bSJed Brown /* ------------------------------------------------------------------- */
464c4762a1bSJed Brown /*
465c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
466c4762a1bSJed Brown    the region.
467c4762a1bSJed Brown 
468c4762a1bSJed Brown    Input Parameter:
469c4762a1bSJed Brown .  user - user-defined application context
470c4762a1bSJed Brown 
471c4762a1bSJed Brown    Output Parameter:
472c4762a1bSJed Brown .  user - user-defined application context
473c4762a1bSJed Brown */
474c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
475c4762a1bSJed Brown {
476c4762a1bSJed Brown   PetscInt       i,j,k,limit=0;
477c4762a1bSJed Brown   PetscInt       maxits=5;
478c4762a1bSJed Brown   PetscInt       mx=user->mx,my=user->my;
479c4762a1bSJed Brown   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
480c4762a1bSJed Brown   PetscReal      one=1.0, two=2.0, three=3.0, tol=1e-10;
481c4762a1bSJed Brown   PetscReal      fnorm,det,hx,hy,xt=0,yt=0;
482c4762a1bSJed Brown   PetscReal      u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
483c4762a1bSJed Brown   PetscReal      b=-0.5, t=0.5, l=-0.5, r=0.5;
484c4762a1bSJed Brown   PetscReal      *boundary;
485c4762a1bSJed Brown 
486c4762a1bSJed Brown   PetscFunctionBeginUser;
487c4762a1bSJed Brown   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
488c4762a1bSJed Brown 
489*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bsize,&user->bottom));
490*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(tsize,&user->top));
491*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(lsize,&user->left));
492*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(rsize,&user->right));
493c4762a1bSJed Brown 
494c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   for (j=0; j<4; j++) {
497c4762a1bSJed Brown     if (j==0) {
498c4762a1bSJed Brown       yt=b;
499c4762a1bSJed Brown       xt=l;
500c4762a1bSJed Brown       limit=bsize;
501c4762a1bSJed Brown       boundary=user->bottom;
502c4762a1bSJed Brown     } else if (j==1) {
503c4762a1bSJed Brown       yt=t;
504c4762a1bSJed Brown       xt=l;
505c4762a1bSJed Brown       limit=tsize;
506c4762a1bSJed Brown       boundary=user->top;
507c4762a1bSJed Brown     } else if (j==2) {
508c4762a1bSJed Brown       yt=b;
509c4762a1bSJed Brown       xt=l;
510c4762a1bSJed Brown       limit=lsize;
511c4762a1bSJed Brown       boundary=user->left;
512c4762a1bSJed Brown     } else {  /*  if (j==3) */
513c4762a1bSJed Brown       yt=b;
514c4762a1bSJed Brown       xt=r;
515c4762a1bSJed Brown       limit=rsize;
516c4762a1bSJed Brown       boundary=user->right;
517c4762a1bSJed Brown     }
518c4762a1bSJed Brown 
519c4762a1bSJed Brown     for (i=0; i<limit; i++) {
520c4762a1bSJed Brown       u1=xt;
521c4762a1bSJed Brown       u2=-yt;
522c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
523c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
524c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
525c4762a1bSJed Brown         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
526c4762a1bSJed Brown         if (fnorm <= tol) break;
527c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
528c4762a1bSJed Brown         njac12=two*u1*u2;
529c4762a1bSJed Brown         njac21=-two*u1*u2;
530c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
531c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
532c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
533c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
534c4762a1bSJed Brown       }
535c4762a1bSJed Brown 
536c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
537c4762a1bSJed Brown       if (j==0 || j==1) {
538c4762a1bSJed Brown         xt=xt+hx;
539c4762a1bSJed Brown       } else { /*  if (j==2 || j==3) */
540c4762a1bSJed Brown         yt=yt+hy;
541c4762a1bSJed Brown       }
542c4762a1bSJed Brown     }
543c4762a1bSJed Brown   }
544c4762a1bSJed Brown   PetscFunctionReturn(0);
545c4762a1bSJed Brown }
546c4762a1bSJed Brown 
547c4762a1bSJed Brown /* ------------------------------------------------------------------- */
548c4762a1bSJed Brown /*
549c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
550c4762a1bSJed Brown 
551c4762a1bSJed Brown    Input Parameters:
552c4762a1bSJed Brown .  user - user-defined application context
553c4762a1bSJed Brown .  X - vector for initial guess
554c4762a1bSJed Brown 
555c4762a1bSJed Brown    Output Parameters:
556c4762a1bSJed Brown .  X - newly computed initial guess
557c4762a1bSJed Brown */
558c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
559c4762a1bSJed Brown {
560c4762a1bSJed Brown   PetscInt       start=-1,i,j;
561c4762a1bSJed Brown   PetscReal      zero=0.0;
562c4762a1bSJed Brown   PetscBool      flg;
563c4762a1bSJed Brown 
564c4762a1bSJed Brown   PetscFunctionBeginUser;
565*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(X, zero));
566*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
567c4762a1bSJed Brown 
568c4762a1bSJed Brown   if (flg && start==0) { /* The zero vector is reasonable */
569*5f80ce2aSJacob Faibussowitsch      CHKERRQ(VecSet(X, zero));
570c4762a1bSJed Brown    } else { /* Take an average of the boundary conditions */
571c4762a1bSJed Brown     PetscInt    row;
572c4762a1bSJed Brown     PetscInt    mx=user->mx,my=user->my;
573c4762a1bSJed Brown     PetscReal *x;
574c4762a1bSJed Brown 
575c4762a1bSJed Brown     /* Get pointers to vector data */
576*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(X,&x));
577c4762a1bSJed Brown     /* Perform local computations */
578c4762a1bSJed Brown     for (j=0; j<my; j++) {
579c4762a1bSJed Brown       for (i=0; i< mx; i++) {
580c4762a1bSJed Brown         row=(j)*mx + (i);
581c4762a1bSJed Brown         x[row] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+ ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0;
582c4762a1bSJed Brown       }
583c4762a1bSJed Brown     }
584c4762a1bSJed Brown     /* Restore vectors */
585*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(X,&x));
586c4762a1bSJed Brown   }
587c4762a1bSJed Brown   PetscFunctionReturn(0);
588c4762a1bSJed Brown }
589c4762a1bSJed Brown 
590c4762a1bSJed Brown /*TEST
591c4762a1bSJed Brown 
592c4762a1bSJed Brown    build:
593c4762a1bSJed Brown       requires: !complex
594c4762a1bSJed Brown 
595c4762a1bSJed Brown    test:
596c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4
597c4762a1bSJed Brown       requires: !single
598c4762a1bSJed Brown 
599c4762a1bSJed Brown    test:
600c4762a1bSJed Brown       suffix: 2
601c4762a1bSJed Brown       args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3
602c4762a1bSJed Brown       requires: !single
603c4762a1bSJed Brown 
604c4762a1bSJed Brown    test:
605c4762a1bSJed Brown       suffix: 3
606c4762a1bSJed Brown       args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3
607c4762a1bSJed Brown       requires: !single
608c4762a1bSJed Brown 
609c4762a1bSJed Brown    test:
610c4762a1bSJed Brown       suffix: 4
611c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4
612c4762a1bSJed Brown       requires: !single
613c4762a1bSJed Brown 
614c4762a1bSJed Brown    test:
615c4762a1bSJed Brown       suffix: 5
616c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4
617c4762a1bSJed Brown       requires: !single
618c4762a1bSJed Brown 
619c4762a1bSJed Brown    test:
620c4762a1bSJed Brown       suffix: 6
621c4762a1bSJed Brown       args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4
622c4762a1bSJed Brown       requires: !single
623c4762a1bSJed Brown 
624c4762a1bSJed Brown    test:
625c4762a1bSJed Brown       suffix: 7
626c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
627c4762a1bSJed Brown       requires: !single
628c4762a1bSJed Brown 
629c4762a1bSJed Brown    test:
630c4762a1bSJed Brown       suffix: 8
631c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
632c4762a1bSJed Brown       requires: !single
633c4762a1bSJed Brown 
634c4762a1bSJed Brown    test:
635c4762a1bSJed Brown       suffix: 9
636c4762a1bSJed Brown       args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
637c4762a1bSJed Brown       requires: !single
638c4762a1bSJed Brown 
63934ad9904SAlp Dener    test:
64034ad9904SAlp Dener       suffix: 10
64134ad9904SAlp Dener       args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_mf_hessian
64234ad9904SAlp Dener 
64334ad9904SAlp Dener    test:
64434ad9904SAlp Dener       suffix: 11
64534ad9904SAlp Dener       args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
64634ad9904SAlp Dener       requires: !single
64734ad9904SAlp Dener 
64834ad9904SAlp Dener    test:
64934ad9904SAlp Dener       suffix: 12
65034ad9904SAlp Dener       args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
65134ad9904SAlp Dener       requires: !single
65234ad9904SAlp Dener 
653c4762a1bSJed Brown TEST*/
654