xref: /petsc/src/tao/complementarity/tutorials/minsurf1.c (revision 362febeeeb69b91ebadcb4b2dc0a22cb6dfc4097)
1c4762a1bSJed Brown #include <petsctao.h>
2c4762a1bSJed Brown 
3c4762a1bSJed Brown static char  help[] =
4c4762a1bSJed Brown "This example demonstrates use of the TAO package to\n\
5c4762a1bSJed Brown solve an unconstrained system of equations.  This example is based on a\n\
6c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and\n\
7c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
8c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
9c4762a1bSJed Brown This application solves this problem using complimentarity -- We are actually\n\
10c4762a1bSJed Brown solving the system  (grad f)_i >= 0, if x_i == l_i \n\
11c4762a1bSJed Brown                     (grad f)_i = 0, if l_i < x_i < u_i \n\
12c4762a1bSJed Brown                     (grad f)_i <= 0, if x_i == u_i  \n\
13c4762a1bSJed Brown where f is the function to be minimized. \n\
14c4762a1bSJed Brown \n\
15c4762a1bSJed Brown The command line options are:\n\
16c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
17c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
18c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /*T
21c4762a1bSJed Brown    Concepts: TAO^Solving a complementarity problem
22c4762a1bSJed Brown    Routines: TaoCreate(); TaoDestroy();
23c4762a1bSJed Brown 
24c4762a1bSJed Brown    Processors: 1
25c4762a1bSJed Brown T*/
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /*
28c4762a1bSJed Brown    User-defined application context - contains data needed by the
29c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
30c4762a1bSJed Brown    FormHessian().
31c4762a1bSJed Brown */
32c4762a1bSJed Brown typedef struct {
33c4762a1bSJed Brown   PetscInt  mx, my;
34c4762a1bSJed Brown   PetscReal *bottom, *top, *left, *right;
35c4762a1bSJed Brown } AppCtx;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown /* -------- User-defined Routines --------- */
38c4762a1bSJed Brown 
39c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
40c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
41c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
42c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *);
43c4762a1bSJed Brown 
44c4762a1bSJed Brown int main(int argc, char **argv)
45c4762a1bSJed Brown {
46c4762a1bSJed Brown   PetscErrorCode ierr;              /* used to check for functions returning nonzeros */
47c4762a1bSJed Brown   Vec            x;                 /* solution vector */
48c4762a1bSJed Brown   Vec            c;                 /* Constraints function vector */
49c4762a1bSJed Brown   Vec            xl,xu;             /* Bounds on the variables */
50c4762a1bSJed Brown   PetscBool      flg;               /* A return variable when checking for user options */
51c4762a1bSJed Brown   Tao            tao;               /* TAO solver context */
52c4762a1bSJed Brown   Mat            J;                 /* Jacobian matrix */
53c4762a1bSJed Brown   PetscInt       N;                 /* Number of elements in vector */
54c4762a1bSJed Brown   PetscScalar    lb =  PETSC_NINFINITY;      /* lower bound constant */
55c4762a1bSJed Brown   PetscScalar    ub =  PETSC_INFINITY;      /* upper bound constant */
56c4762a1bSJed Brown   AppCtx         user;                    /* user-defined work context */
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Initialize PETSc, TAO */
59c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, (char *)0, help);if (ierr) return ierr;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Specify default dimension of the problem */
62c4762a1bSJed Brown   user.mx = 4; user.my = 4;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
65c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-mx", &user.mx, &flg);CHKERRQ(ierr);
66c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-my", &user.my, &flg);CHKERRQ(ierr);
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Calculate any derived values from parameters */
69c4762a1bSJed Brown   N = user.mx*user.my;
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n");CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"mx:%D, my:%D\n", user.mx,user.my);CHKERRQ(ierr);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Create appropriate vectors and matrices */
75c4762a1bSJed Brown   ierr = VecCreateSeq(MPI_COMM_SELF, N, &x);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = VecDuplicate(x, &c);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J);CHKERRQ(ierr);
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* The TAO code begins here */
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
82c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_SELF,&tao);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOSSILS);CHKERRQ(ierr);
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* Set data structure */
86c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao, x);CHKERRQ(ierr);
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /*  Set routines for constraints function and Jacobian evaluation */
89c4762a1bSJed Brown   ierr = TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user);CHKERRQ(ierr);
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* Set the variable bounds */
93c4762a1bSJed Brown   ierr = MSA_BoundaryConditions(&user);CHKERRQ(ierr);
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Set initial solution guess */
96c4762a1bSJed Brown   ierr = MSA_InitialPoint(&user, x);CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* Set Bounds on variables */
99c4762a1bSJed Brown   ierr = VecDuplicate(x, &xl);CHKERRQ(ierr);
100c4762a1bSJed Brown   ierr = VecDuplicate(x, &xu);CHKERRQ(ierr);
101c4762a1bSJed Brown   ierr = VecSet(xl, lb);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = VecSet(xu, ub);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,xl,xu);CHKERRQ(ierr);
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Check for any tao command line options */
106c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* Solve the application */
109c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* Free Tao data structures */
112c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Free PETSc data structures */
115c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = VecDestroy(&xl);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = VecDestroy(&xu);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = VecDestroy(&c);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* Free user-created data structures */
122c4762a1bSJed Brown   ierr = PetscFree(user.bottom);CHKERRQ(ierr);
123c4762a1bSJed Brown   ierr = PetscFree(user.top);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = PetscFree(user.left);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = PetscFree(user.right);CHKERRQ(ierr);
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   ierr = PetscFinalize();
128c4762a1bSJed Brown   return ierr;
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown /* -------------------------------------------------------------------- */
132c4762a1bSJed Brown 
133c4762a1bSJed Brown /*  FormConstraints - Evaluates gradient of f.
134c4762a1bSJed Brown 
135c4762a1bSJed Brown     Input Parameters:
136c4762a1bSJed Brown .   tao  - the TAO_APPLICATION context
137c4762a1bSJed Brown .   X    - input vector
138c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoSetConstraintsRoutine()
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     Output Parameters:
141c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
142c4762a1bSJed Brown */
143c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr)
144c4762a1bSJed Brown {
145c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
146c4762a1bSJed Brown   PetscErrorCode ierr;
147c4762a1bSJed Brown   PetscInt       i,j,row;
148c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
149c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
150c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
151c4762a1bSJed Brown   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
152c4762a1bSJed Brown   PetscScalar    zero=0.0;
153c4762a1bSJed Brown   PetscScalar    *g, *x;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   PetscFunctionBegin;
156c4762a1bSJed Brown   /* Initialize vector to zero */
157c4762a1bSJed Brown   ierr = VecSet(G, zero);CHKERRQ(ierr);
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Get pointers to vector data */
160c4762a1bSJed Brown   ierr = VecGetArray(X, &x);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = VecGetArray(G, &g);CHKERRQ(ierr);
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
164c4762a1bSJed Brown   for (j=0; j<my; j++) {
165c4762a1bSJed Brown     for (i=0; i< mx; i++) {
166c4762a1bSJed Brown       row= j*mx + i;
167c4762a1bSJed Brown 
168c4762a1bSJed Brown       xc = x[row];
169c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown       if (i==0) { /* left side */
172c4762a1bSJed Brown         xl= user->left[j+1];
173c4762a1bSJed Brown         xlt = user->left[j+2];
174c4762a1bSJed Brown       } else {
175c4762a1bSJed Brown         xl = x[row-1];
176c4762a1bSJed Brown       }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown       if (j==0) { /* bottom side */
179c4762a1bSJed Brown         xb=user->bottom[i+1];
180c4762a1bSJed Brown         xrb = user->bottom[i+2];
181c4762a1bSJed Brown       } else {
182c4762a1bSJed Brown         xb = x[row-mx];
183c4762a1bSJed Brown       }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown       if (i+1 == mx) { /* right side */
186c4762a1bSJed Brown         xr=user->right[j+1];
187c4762a1bSJed Brown         xrb = user->right[j];
188c4762a1bSJed Brown       } else {
189c4762a1bSJed Brown         xr = x[row+1];
190c4762a1bSJed Brown       }
191c4762a1bSJed Brown 
192c4762a1bSJed Brown       if (j+1==0+my) { /* top side */
193c4762a1bSJed Brown         xt=user->top[i+1];
194c4762a1bSJed Brown         xlt = user->top[i];
195c4762a1bSJed Brown       }else {
196c4762a1bSJed Brown         xt = x[row+mx];
197c4762a1bSJed Brown       }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown       if (i>0 && j+1<my) {
200c4762a1bSJed Brown         xlt = x[row-1+mx];
201c4762a1bSJed Brown       }
202c4762a1bSJed Brown       if (j>0 && i+1<mx) {
203c4762a1bSJed Brown         xrb = x[row+1-mx];
204c4762a1bSJed Brown       }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown       d1 = (xc-xl);
207c4762a1bSJed Brown       d2 = (xc-xr);
208c4762a1bSJed Brown       d3 = (xc-xt);
209c4762a1bSJed Brown       d4 = (xc-xb);
210c4762a1bSJed Brown       d5 = (xr-xrb);
211c4762a1bSJed Brown       d6 = (xrb-xb);
212c4762a1bSJed Brown       d7 = (xlt-xl);
213c4762a1bSJed Brown       d8 = (xt-xlt);
214c4762a1bSJed Brown 
215c4762a1bSJed Brown       df1dxc = d1*hydhx;
216c4762a1bSJed Brown       df2dxc = (d1*hydhx + d4*hxdhy);
217c4762a1bSJed Brown       df3dxc = d3*hxdhy;
218c4762a1bSJed Brown       df4dxc = (d2*hydhx + d3*hxdhy);
219c4762a1bSJed Brown       df5dxc = d2*hydhx;
220c4762a1bSJed Brown       df6dxc = d4*hxdhy;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown       d1 /= hx;
223c4762a1bSJed Brown       d2 /= hx;
224c4762a1bSJed Brown       d3 /= hy;
225c4762a1bSJed Brown       d4 /= hy;
226c4762a1bSJed Brown       d5 /= hy;
227c4762a1bSJed Brown       d6 /= hx;
228c4762a1bSJed Brown       d7 /= hy;
229c4762a1bSJed Brown       d8 /= hx;
230c4762a1bSJed Brown 
231c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
232c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
233c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
234c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
235c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
236c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
237c4762a1bSJed Brown 
238c4762a1bSJed Brown       df1dxc /= f1;
239c4762a1bSJed Brown       df2dxc /= f2;
240c4762a1bSJed Brown       df3dxc /= f3;
241c4762a1bSJed Brown       df4dxc /= f4;
242c4762a1bSJed Brown       df5dxc /= f5;
243c4762a1bSJed Brown       df6dxc /= f6;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
246c4762a1bSJed Brown     }
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Restore vectors */
250c4762a1bSJed Brown   ierr = VecRestoreArray(X, &x);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = VecRestoreArray(G, &g);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = PetscLogFlops(67*mx*my);CHKERRQ(ierr);
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown /* ------------------------------------------------------------------- */
257c4762a1bSJed Brown /*
258c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
259c4762a1bSJed Brown 
260c4762a1bSJed Brown    Input Parameters:
261c4762a1bSJed Brown .  tao  - the TAO_APPLICATION context
262c4762a1bSJed Brown .  X    - input vector
263c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetJacobian()
264c4762a1bSJed Brown 
265c4762a1bSJed Brown    Output Parameters:
266c4762a1bSJed Brown .  tH    - Jacobian matrix
267c4762a1bSJed Brown 
268c4762a1bSJed Brown */
269c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao tao, Vec X, Mat H, Mat tHPre, void *ptr)
270c4762a1bSJed Brown {
271c4762a1bSJed Brown   AppCtx            *user = (AppCtx *) ptr;
272c4762a1bSJed Brown   PetscErrorCode    ierr;
273c4762a1bSJed Brown   PetscInt          i,j,k,row;
274c4762a1bSJed Brown   PetscInt          mx=user->mx, my=user->my;
275c4762a1bSJed Brown   PetscInt          col[7];
276c4762a1bSJed Brown   PetscReal         hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
277c4762a1bSJed Brown   PetscReal         f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
278c4762a1bSJed Brown   PetscReal         hl,hr,ht,hb,hc,htl,hbr;
279c4762a1bSJed Brown   const PetscScalar *x;
280c4762a1bSJed Brown   PetscScalar       v[7];
281c4762a1bSJed Brown   PetscBool         assembled;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /* Set various matrix options */
284*362febeeSStefano Zampini   PetscFunctionBegin;
285c4762a1bSJed Brown   ierr = MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = MatAssembled(H,&assembled);CHKERRQ(ierr);
287c4762a1bSJed Brown   if (assembled) {ierr = MatZeroEntries(H);CHKERRQ(ierr);}
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   /* Get pointers to vector data */
290c4762a1bSJed Brown   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   /* Compute Jacobian over the locally owned part of the mesh */
293c4762a1bSJed Brown   for (i=0; i< mx; i++) {
294c4762a1bSJed Brown     for (j=0; j<my; j++) {
295c4762a1bSJed Brown       row= j*mx + i;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown       xc = x[row];
298c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
299c4762a1bSJed Brown 
300c4762a1bSJed Brown       /* Left side */
301c4762a1bSJed Brown       if (i==0) {
302c4762a1bSJed Brown         xl  = user->left[j+1];
303c4762a1bSJed Brown         xlt = user->left[j+2];
304c4762a1bSJed Brown       } else {
305c4762a1bSJed Brown         xl = x[row-1];
306c4762a1bSJed Brown       }
307c4762a1bSJed Brown 
308c4762a1bSJed Brown       if (j==0) {
309c4762a1bSJed Brown         xb  = user->bottom[i+1];
310c4762a1bSJed Brown         xrb = user->bottom[i+2];
311c4762a1bSJed Brown       } else {
312c4762a1bSJed Brown         xb = x[row-mx];
313c4762a1bSJed Brown       }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown       if (i+1 == mx) {
316c4762a1bSJed Brown         xr  = user->right[j+1];
317c4762a1bSJed Brown         xrb = user->right[j];
318c4762a1bSJed Brown       } else {
319c4762a1bSJed Brown         xr = x[row+1];
320c4762a1bSJed Brown       }
321c4762a1bSJed Brown 
322c4762a1bSJed Brown       if (j+1==my) {
323c4762a1bSJed Brown         xt  = user->top[i+1];
324c4762a1bSJed Brown         xlt = user->top[i];
325c4762a1bSJed Brown       }else {
326c4762a1bSJed Brown         xt = x[row+mx];
327c4762a1bSJed Brown       }
328c4762a1bSJed Brown 
329c4762a1bSJed Brown       if (i>0 && j+1<my) {
330c4762a1bSJed Brown         xlt = x[row-1+mx];
331c4762a1bSJed Brown       }
332c4762a1bSJed Brown       if (j>0 && i+1<mx) {
333c4762a1bSJed Brown         xrb = x[row+1-mx];
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown       d1 = (xc-xl)/hx;
337c4762a1bSJed Brown       d2 = (xc-xr)/hx;
338c4762a1bSJed Brown       d3 = (xc-xt)/hy;
339c4762a1bSJed Brown       d4 = (xc-xb)/hy;
340c4762a1bSJed Brown       d5 = (xrb-xr)/hy;
341c4762a1bSJed Brown       d6 = (xrb-xb)/hx;
342c4762a1bSJed Brown       d7 = (xlt-xl)/hy;
343c4762a1bSJed Brown       d8 = (xlt-xt)/hx;
344c4762a1bSJed Brown 
345c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
346c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
347c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
348c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
349c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
350c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
351c4762a1bSJed Brown 
352c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
353c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
354c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
355c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
356c4762a1bSJed Brown 
357c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
358c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
359c4762a1bSJed Brown 
360c4762a1bSJed 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) +
361c4762a1bSJed 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);
362c4762a1bSJed Brown 
363c4762a1bSJed Brown       hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0;  hc/=2.0;
364c4762a1bSJed Brown 
365c4762a1bSJed Brown       k=0;
366c4762a1bSJed Brown       if (j>0) {
367c4762a1bSJed Brown         v[k]=hb; col[k]=row - mx; k++;
368c4762a1bSJed Brown       }
369c4762a1bSJed Brown 
370c4762a1bSJed Brown       if (j>0 && i < mx -1) {
371c4762a1bSJed Brown         v[k]=hbr; col[k]=row - mx+1; k++;
372c4762a1bSJed Brown       }
373c4762a1bSJed Brown 
374c4762a1bSJed Brown       if (i>0) {
375c4762a1bSJed Brown         v[k]= hl; col[k]=row - 1; k++;
376c4762a1bSJed Brown       }
377c4762a1bSJed Brown 
378c4762a1bSJed Brown       v[k]= hc; col[k]=row; k++;
379c4762a1bSJed Brown 
380c4762a1bSJed Brown       if (i < mx-1) {
381c4762a1bSJed Brown         v[k]= hr; col[k]=row+1; k++;
382c4762a1bSJed Brown       }
383c4762a1bSJed Brown 
384c4762a1bSJed Brown       if (i>0 && j < my-1) {
385c4762a1bSJed Brown         v[k]= htl; col[k] = row+mx-1; k++;
386c4762a1bSJed Brown       }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown       if (j < my-1) {
389c4762a1bSJed Brown         v[k]= ht; col[k] = row+mx; k++;
390c4762a1bSJed Brown       }
391c4762a1bSJed Brown 
392c4762a1bSJed Brown       /*
393c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
394c4762a1bSJed Brown          earlier, in the main routine.
395c4762a1bSJed Brown       */
396c4762a1bSJed Brown       ierr = MatSetValues(H,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr);
397c4762a1bSJed Brown     }
398c4762a1bSJed Brown   }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   /* Restore vectors */
401c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
402c4762a1bSJed Brown 
403c4762a1bSJed Brown   /* Assemble the matrix */
404c4762a1bSJed Brown   ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
405c4762a1bSJed Brown   ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
406c4762a1bSJed Brown   ierr = PetscLogFlops(199*mx*my);CHKERRQ(ierr);
407c4762a1bSJed Brown   PetscFunctionReturn(0);
408c4762a1bSJed Brown }
409c4762a1bSJed Brown 
410c4762a1bSJed Brown /* ------------------------------------------------------------------- */
411c4762a1bSJed Brown /*
412c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
413c4762a1bSJed Brown    the region.
414c4762a1bSJed Brown 
415c4762a1bSJed Brown    Input Parameter:
416c4762a1bSJed Brown .  user - user-defined application context
417c4762a1bSJed Brown 
418c4762a1bSJed Brown    Output Parameter:
419c4762a1bSJed Brown .  user - user-defined application context
420c4762a1bSJed Brown */
421c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
422c4762a1bSJed Brown {
423c4762a1bSJed Brown   PetscErrorCode  ierr;
424c4762a1bSJed Brown   PetscInt        i,j,k,limit=0,maxits=5;
425c4762a1bSJed Brown   PetscInt        mx=user->mx,my=user->my;
426c4762a1bSJed Brown   PetscInt        bsize=0, lsize=0, tsize=0, rsize=0;
427c4762a1bSJed Brown   PetscReal       one=1.0, two=2.0, three=3.0, tol=1e-10;
428c4762a1bSJed Brown   PetscReal       fnorm,det,hx,hy,xt=0,yt=0;
429c4762a1bSJed Brown   PetscReal       u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
430c4762a1bSJed Brown   PetscReal       b=-0.5, t=0.5, l=-0.5, r=0.5;
431c4762a1bSJed Brown   PetscReal       *boundary;
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   PetscFunctionBegin;
434c4762a1bSJed Brown   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   ierr = PetscMalloc1(bsize, &user->bottom);CHKERRQ(ierr);
437c4762a1bSJed Brown   ierr = PetscMalloc1(tsize, &user->top);CHKERRQ(ierr);
438c4762a1bSJed Brown   ierr = PetscMalloc1(lsize, &user->left);CHKERRQ(ierr);
439c4762a1bSJed Brown   ierr = PetscMalloc1(rsize, &user->right);CHKERRQ(ierr);
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
442c4762a1bSJed Brown 
443c4762a1bSJed Brown   for (j=0; j<4; j++) {
444c4762a1bSJed Brown     if (j==0) {
445c4762a1bSJed Brown       yt=b;
446c4762a1bSJed Brown       xt=l;
447c4762a1bSJed Brown       limit=bsize;
448c4762a1bSJed Brown       boundary=user->bottom;
449c4762a1bSJed Brown     } else if (j==1) {
450c4762a1bSJed Brown       yt=t;
451c4762a1bSJed Brown       xt=l;
452c4762a1bSJed Brown       limit=tsize;
453c4762a1bSJed Brown       boundary=user->top;
454c4762a1bSJed Brown     } else if (j==2) {
455c4762a1bSJed Brown       yt=b;
456c4762a1bSJed Brown       xt=l;
457c4762a1bSJed Brown       limit=lsize;
458c4762a1bSJed Brown       boundary=user->left;
459c4762a1bSJed Brown     } else { /* if  (j==3) */
460c4762a1bSJed Brown       yt=b;
461c4762a1bSJed Brown       xt=r;
462c4762a1bSJed Brown       limit=rsize;
463c4762a1bSJed Brown       boundary=user->right;
464c4762a1bSJed Brown     }
465c4762a1bSJed Brown 
466c4762a1bSJed Brown     for (i=0; i<limit; i++) {
467c4762a1bSJed Brown       u1=xt;
468c4762a1bSJed Brown       u2=-yt;
469c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
470c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
471c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
472c4762a1bSJed Brown         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
473c4762a1bSJed Brown         if (fnorm <= tol) break;
474c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
475c4762a1bSJed Brown         njac12=two*u1*u2;
476c4762a1bSJed Brown         njac21=-two*u1*u2;
477c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
478c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
479c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
480c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
481c4762a1bSJed Brown       }
482c4762a1bSJed Brown 
483c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
484c4762a1bSJed Brown       if (j==0 || j==1) {
485c4762a1bSJed Brown         xt=xt+hx;
486c4762a1bSJed Brown       } else { /* if (j==2 || j==3) */
487c4762a1bSJed Brown         yt=yt+hy;
488c4762a1bSJed Brown       }
489c4762a1bSJed Brown     }
490c4762a1bSJed Brown   }
491c4762a1bSJed Brown   PetscFunctionReturn(0);
492c4762a1bSJed Brown }
493c4762a1bSJed Brown 
494c4762a1bSJed Brown /* ------------------------------------------------------------------- */
495c4762a1bSJed Brown /*
496c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
497c4762a1bSJed Brown 
498c4762a1bSJed Brown    Input Parameters:
499c4762a1bSJed Brown .  user - user-defined application context
500c4762a1bSJed Brown .  X - vector for initial guess
501c4762a1bSJed Brown 
502c4762a1bSJed Brown    Output Parameters:
503c4762a1bSJed Brown .  X - newly computed initial guess
504c4762a1bSJed Brown */
505c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
506c4762a1bSJed Brown {
507c4762a1bSJed Brown   PetscErrorCode ierr;
508c4762a1bSJed Brown   PetscInt       start=-1,i,j;
509c4762a1bSJed Brown   PetscScalar    zero=0.0;
510c4762a1bSJed Brown   PetscBool      flg;
511c4762a1bSJed Brown 
512c4762a1bSJed Brown   PetscFunctionBegin;
513c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg);CHKERRQ(ierr);
514c4762a1bSJed Brown 
515c4762a1bSJed Brown   if (flg && start==0) { /* The zero vector is reasonable */
516c4762a1bSJed Brown     ierr = VecSet(X, zero);CHKERRQ(ierr);
517c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
518c4762a1bSJed Brown     PetscInt    row;
519c4762a1bSJed Brown     PetscInt    mx=user->mx,my=user->my;
520c4762a1bSJed Brown     PetscScalar *x;
521c4762a1bSJed Brown 
522c4762a1bSJed Brown     /* Get pointers to vector data */
523c4762a1bSJed Brown     ierr = VecGetArray(X,&x);CHKERRQ(ierr);
524c4762a1bSJed Brown 
525c4762a1bSJed Brown     /* Perform local computations */
526c4762a1bSJed Brown     for (j=0; j<my; j++) {
527c4762a1bSJed Brown       for (i=0; i< mx; i++) {
528c4762a1bSJed Brown         row=(j)*mx + (i);
529c4762a1bSJed 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;
530c4762a1bSJed Brown       }
531c4762a1bSJed Brown     }
532c4762a1bSJed Brown 
533c4762a1bSJed Brown     /* Restore vectors */
534c4762a1bSJed Brown     ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
535c4762a1bSJed Brown   }
536c4762a1bSJed Brown   PetscFunctionReturn(0);
537c4762a1bSJed Brown }
538c4762a1bSJed Brown 
539c4762a1bSJed Brown /*TEST
540c4762a1bSJed Brown 
541c4762a1bSJed Brown    build:
542c4762a1bSJed Brown       requires: !complex
543c4762a1bSJed Brown 
544c4762a1bSJed Brown    test:
545c4762a1bSJed Brown       args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5
546c4762a1bSJed Brown       requires: !single
547c4762a1bSJed Brown 
548c4762a1bSJed Brown    test:
549c4762a1bSJed Brown       suffix: 2
550c4762a1bSJed Brown       args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5
551c4762a1bSJed Brown 
552c4762a1bSJed Brown TEST*/
553