xref: /petsc/src/tao/complementarity/tutorials/minsurf1.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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 /*
21c4762a1bSJed Brown    User-defined application context - contains data needed by the
22c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
23c4762a1bSJed Brown    FormHessian().
24c4762a1bSJed Brown */
25c4762a1bSJed Brown typedef struct {
26c4762a1bSJed Brown   PetscInt  mx, my;
27c4762a1bSJed Brown   PetscReal *bottom, *top, *left, *right;
28c4762a1bSJed Brown } AppCtx;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown /* -------- User-defined Routines --------- */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
33c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
34c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
35c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown int main(int argc, char **argv)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   Vec            x;                 /* solution vector */
40c4762a1bSJed Brown   Vec            c;                 /* Constraints function vector */
41c4762a1bSJed Brown   Vec            xl,xu;             /* Bounds on the variables */
42c4762a1bSJed Brown   PetscBool      flg;               /* A return variable when checking for user options */
43c4762a1bSJed Brown   Tao            tao;               /* TAO solver context */
44c4762a1bSJed Brown   Mat            J;                 /* Jacobian matrix */
45c4762a1bSJed Brown   PetscInt       N;                 /* Number of elements in vector */
46c4762a1bSJed Brown   PetscScalar    lb =  PETSC_NINFINITY;      /* lower bound constant */
47c4762a1bSJed Brown   PetscScalar    ub =  PETSC_INFINITY;      /* upper bound constant */
48c4762a1bSJed Brown   AppCtx         user;                    /* user-defined work context */
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Initialize PETSc, TAO */
519566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Specify default dimension of the problem */
54c4762a1bSJed Brown   user.mx = 4; user.my = 4;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-mx", &user.mx, &flg));
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-my", &user.my, &flg));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Calculate any derived values from parameters */
61c4762a1bSJed Brown   N = user.mx*user.my;
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n"));
64*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"mx:%" PetscInt_FMT ", my:%" PetscInt_FMT "\n", user.mx,user.my));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* Create appropriate vectors and matrices */
679566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(MPI_COMM_SELF, N, &x));
689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &c));
699566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* The TAO code begins here */
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
749566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF,&tao));
759566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOSSILS));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Set data structure */
789566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /*  Set routines for constraints function and Jacobian evaluation */
819566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user));
829566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* Set the variable bounds */
859566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Set initial solution guess */
889566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user, x));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* Set Bounds on variables */
919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &xl));
929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &xu));
939566063dSJacob Faibussowitsch   PetscCall(VecSet(xl, lb));
949566063dSJacob Faibussowitsch   PetscCall(VecSet(xu, ub));
959566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao,xl,xu));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* Check for any tao command line options */
989566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* Solve the application */
1019566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Free Tao data structures */
1049566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* Free PETSc data structures */
1079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
1099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
1109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&c));
1119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* Free user-created data structures */
1149566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.bottom));
1159566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.top));
1169566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.left));
1179566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.right));
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
120b122ec5aSJacob Faibussowitsch   return 0;
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown /* -------------------------------------------------------------------- */
124c4762a1bSJed Brown 
125c4762a1bSJed Brown /*  FormConstraints - Evaluates gradient of f.
126c4762a1bSJed Brown 
127c4762a1bSJed Brown     Input Parameters:
128c4762a1bSJed Brown .   tao  - the TAO_APPLICATION context
129c4762a1bSJed Brown .   X    - input vector
130c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoSetConstraintsRoutine()
131c4762a1bSJed Brown 
132c4762a1bSJed Brown     Output Parameters:
133c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
134c4762a1bSJed Brown */
135c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
138c4762a1bSJed Brown   PetscInt       i,j,row;
139c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
140c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
141c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
142c4762a1bSJed Brown   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
143c4762a1bSJed Brown   PetscScalar    zero=0.0;
144c4762a1bSJed Brown   PetscScalar    *g, *x;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   PetscFunctionBegin;
147c4762a1bSJed Brown   /* Initialize vector to zero */
1489566063dSJacob Faibussowitsch   PetscCall(VecSet(G, zero));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Get pointers to vector data */
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
1529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
155c4762a1bSJed Brown   for (j=0; j<my; j++) {
156c4762a1bSJed Brown     for (i=0; i< mx; i++) {
157c4762a1bSJed Brown       row= j*mx + i;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown       xc = x[row];
160c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
161c4762a1bSJed Brown 
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 /= hx;
214c4762a1bSJed Brown       d2 /= hx;
215c4762a1bSJed Brown       d3 /= hy;
216c4762a1bSJed Brown       d4 /= hy;
217c4762a1bSJed Brown       d5 /= hy;
218c4762a1bSJed Brown       d6 /= hx;
219c4762a1bSJed Brown       d7 /= hy;
220c4762a1bSJed Brown       d8 /= hx;
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       df1dxc /= f1;
230c4762a1bSJed Brown       df2dxc /= f2;
231c4762a1bSJed Brown       df3dxc /= f3;
232c4762a1bSJed Brown       df4dxc /= f4;
233c4762a1bSJed Brown       df5dxc /= f5;
234c4762a1bSJed Brown       df6dxc /= f6;
235c4762a1bSJed Brown 
236c4762a1bSJed Brown       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
237c4762a1bSJed Brown     }
238c4762a1bSJed Brown   }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   /* Restore vectors */
2419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
2429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
2439566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(67*mx*my));
244c4762a1bSJed Brown   PetscFunctionReturn(0);
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown /* ------------------------------------------------------------------- */
248c4762a1bSJed Brown /*
249c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
250c4762a1bSJed Brown 
251c4762a1bSJed Brown    Input Parameters:
252c4762a1bSJed Brown .  tao  - the TAO_APPLICATION context
253c4762a1bSJed Brown .  X    - input vector
254c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetJacobian()
255c4762a1bSJed Brown 
256c4762a1bSJed Brown    Output Parameters:
257c4762a1bSJed Brown .  tH    - Jacobian matrix
258c4762a1bSJed Brown 
259c4762a1bSJed Brown */
260c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao tao, Vec X, Mat H, Mat tHPre, void *ptr)
261c4762a1bSJed Brown {
262c4762a1bSJed Brown   AppCtx            *user = (AppCtx *) ptr;
263c4762a1bSJed Brown   PetscInt          i,j,k,row;
264c4762a1bSJed Brown   PetscInt          mx=user->mx, my=user->my;
265c4762a1bSJed Brown   PetscInt          col[7];
266c4762a1bSJed Brown   PetscReal         hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
267c4762a1bSJed Brown   PetscReal         f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
268c4762a1bSJed Brown   PetscReal         hl,hr,ht,hb,hc,htl,hbr;
269c4762a1bSJed Brown   const PetscScalar *x;
270c4762a1bSJed Brown   PetscScalar       v[7];
271c4762a1bSJed Brown   PetscBool         assembled;
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   /* Set various matrix options */
274362febeeSStefano Zampini   PetscFunctionBegin;
2759566063dSJacob Faibussowitsch   PetscCall(MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
2769566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H,&assembled));
2779566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(H));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /* Get pointers to vector data */
2809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Compute Jacobian over the locally owned part of the mesh */
283c4762a1bSJed Brown   for (i=0; i< mx; i++) {
284c4762a1bSJed Brown     for (j=0; j<my; j++) {
285c4762a1bSJed Brown       row= j*mx + i;
286c4762a1bSJed Brown 
287c4762a1bSJed Brown       xc = x[row];
288c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
289c4762a1bSJed Brown 
290c4762a1bSJed Brown       /* Left side */
291c4762a1bSJed Brown       if (i==0) {
292c4762a1bSJed Brown         xl  = user->left[j+1];
293c4762a1bSJed Brown         xlt = user->left[j+2];
294c4762a1bSJed Brown       } else {
295c4762a1bSJed Brown         xl = x[row-1];
296c4762a1bSJed Brown       }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown       if (j==0) {
299c4762a1bSJed Brown         xb  = user->bottom[i+1];
300c4762a1bSJed Brown         xrb = user->bottom[i+2];
301c4762a1bSJed Brown       } else {
302c4762a1bSJed Brown         xb = x[row-mx];
303c4762a1bSJed Brown       }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown       if (i+1 == mx) {
306c4762a1bSJed Brown         xr  = user->right[j+1];
307c4762a1bSJed Brown         xrb = user->right[j];
308c4762a1bSJed Brown       } else {
309c4762a1bSJed Brown         xr = x[row+1];
310c4762a1bSJed Brown       }
311c4762a1bSJed Brown 
312c4762a1bSJed Brown       if (j+1==my) {
313c4762a1bSJed Brown         xt  = user->top[i+1];
314c4762a1bSJed Brown         xlt = user->top[i];
315c4762a1bSJed Brown       }else {
316c4762a1bSJed Brown         xt = x[row+mx];
317c4762a1bSJed Brown       }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown       if (i>0 && j+1<my) {
320c4762a1bSJed Brown         xlt = x[row-1+mx];
321c4762a1bSJed Brown       }
322c4762a1bSJed Brown       if (j>0 && i+1<mx) {
323c4762a1bSJed Brown         xrb = x[row+1-mx];
324c4762a1bSJed Brown       }
325c4762a1bSJed Brown 
326c4762a1bSJed Brown       d1 = (xc-xl)/hx;
327c4762a1bSJed Brown       d2 = (xc-xr)/hx;
328c4762a1bSJed Brown       d3 = (xc-xt)/hy;
329c4762a1bSJed Brown       d4 = (xc-xb)/hy;
330c4762a1bSJed Brown       d5 = (xrb-xr)/hy;
331c4762a1bSJed Brown       d6 = (xrb-xb)/hx;
332c4762a1bSJed Brown       d7 = (xlt-xl)/hy;
333c4762a1bSJed Brown       d8 = (xlt-xt)/hx;
334c4762a1bSJed Brown 
335c4762a1bSJed Brown       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
336c4762a1bSJed Brown       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
337c4762a1bSJed Brown       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
338c4762a1bSJed Brown       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
339c4762a1bSJed Brown       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
340c4762a1bSJed Brown       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
341c4762a1bSJed Brown 
342c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
343c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
344c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
345c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
346c4762a1bSJed Brown 
347c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
348c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
349c4762a1bSJed Brown 
350c4762a1bSJed 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) +
351c4762a1bSJed 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);
352c4762a1bSJed Brown 
353c4762a1bSJed Brown       hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0;  hc/=2.0;
354c4762a1bSJed Brown 
355c4762a1bSJed Brown       k=0;
356c4762a1bSJed Brown       if (j>0) {
357c4762a1bSJed Brown         v[k]=hb; col[k]=row - mx; k++;
358c4762a1bSJed Brown       }
359c4762a1bSJed Brown 
360c4762a1bSJed Brown       if (j>0 && i < mx -1) {
361c4762a1bSJed Brown         v[k]=hbr; col[k]=row - mx+1; k++;
362c4762a1bSJed Brown       }
363c4762a1bSJed Brown 
364c4762a1bSJed Brown       if (i>0) {
365c4762a1bSJed Brown         v[k]= hl; col[k]=row - 1; k++;
366c4762a1bSJed Brown       }
367c4762a1bSJed Brown 
368c4762a1bSJed Brown       v[k]= hc; col[k]=row; k++;
369c4762a1bSJed Brown 
370c4762a1bSJed Brown       if (i < mx-1) {
371c4762a1bSJed Brown         v[k]= hr; col[k]=row+1; k++;
372c4762a1bSJed Brown       }
373c4762a1bSJed Brown 
374c4762a1bSJed Brown       if (i>0 && j < my-1) {
375c4762a1bSJed Brown         v[k]= htl; col[k] = row+mx-1; k++;
376c4762a1bSJed Brown       }
377c4762a1bSJed Brown 
378c4762a1bSJed Brown       if (j < my-1) {
379c4762a1bSJed Brown         v[k]= ht; col[k] = row+mx; k++;
380c4762a1bSJed Brown       }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown       /*
383c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
384c4762a1bSJed Brown          earlier, in the main routine.
385c4762a1bSJed Brown       */
3869566063dSJacob Faibussowitsch       PetscCall(MatSetValues(H,1,&row,k,col,v,INSERT_VALUES));
387c4762a1bSJed Brown     }
388c4762a1bSJed Brown   }
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   /* Restore vectors */
3919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   /* Assemble the matrix */
3949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
3959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
3969566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199*mx*my));
397c4762a1bSJed Brown   PetscFunctionReturn(0);
398c4762a1bSJed Brown }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown /* ------------------------------------------------------------------- */
401c4762a1bSJed Brown /*
402c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
403c4762a1bSJed Brown    the region.
404c4762a1bSJed Brown 
405c4762a1bSJed Brown    Input Parameter:
406c4762a1bSJed Brown .  user - user-defined application context
407c4762a1bSJed Brown 
408c4762a1bSJed Brown    Output Parameter:
409c4762a1bSJed Brown .  user - user-defined application context
410c4762a1bSJed Brown */
411c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
412c4762a1bSJed Brown {
413c4762a1bSJed Brown   PetscInt        i,j,k,limit=0,maxits=5;
414c4762a1bSJed Brown   PetscInt        mx=user->mx,my=user->my;
415c4762a1bSJed Brown   PetscInt        bsize=0, lsize=0, tsize=0, rsize=0;
416c4762a1bSJed Brown   PetscReal       one=1.0, two=2.0, three=3.0, tol=1e-10;
417c4762a1bSJed Brown   PetscReal       fnorm,det,hx,hy,xt=0,yt=0;
418c4762a1bSJed Brown   PetscReal       u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
419c4762a1bSJed Brown   PetscReal       b=-0.5, t=0.5, l=-0.5, r=0.5;
420c4762a1bSJed Brown   PetscReal       *boundary;
421c4762a1bSJed Brown 
422c4762a1bSJed Brown   PetscFunctionBegin;
423c4762a1bSJed Brown   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
424c4762a1bSJed Brown 
4259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bsize, &user->bottom));
4269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsize, &user->top));
4279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lsize, &user->left));
4289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(rsize, &user->right));
429c4762a1bSJed Brown 
430c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
431c4762a1bSJed Brown 
432c4762a1bSJed Brown   for (j=0; j<4; j++) {
433c4762a1bSJed Brown     if (j==0) {
434c4762a1bSJed Brown       yt=b;
435c4762a1bSJed Brown       xt=l;
436c4762a1bSJed Brown       limit=bsize;
437c4762a1bSJed Brown       boundary=user->bottom;
438c4762a1bSJed Brown     } else if (j==1) {
439c4762a1bSJed Brown       yt=t;
440c4762a1bSJed Brown       xt=l;
441c4762a1bSJed Brown       limit=tsize;
442c4762a1bSJed Brown       boundary=user->top;
443c4762a1bSJed Brown     } else if (j==2) {
444c4762a1bSJed Brown       yt=b;
445c4762a1bSJed Brown       xt=l;
446c4762a1bSJed Brown       limit=lsize;
447c4762a1bSJed Brown       boundary=user->left;
448c4762a1bSJed Brown     } else { /* if  (j==3) */
449c4762a1bSJed Brown       yt=b;
450c4762a1bSJed Brown       xt=r;
451c4762a1bSJed Brown       limit=rsize;
452c4762a1bSJed Brown       boundary=user->right;
453c4762a1bSJed Brown     }
454c4762a1bSJed Brown 
455c4762a1bSJed Brown     for (i=0; i<limit; i++) {
456c4762a1bSJed Brown       u1=xt;
457c4762a1bSJed Brown       u2=-yt;
458c4762a1bSJed Brown       for (k=0; k<maxits; k++) {
459c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
460c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
461c4762a1bSJed Brown         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
462c4762a1bSJed Brown         if (fnorm <= tol) break;
463c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
464c4762a1bSJed Brown         njac12=two*u1*u2;
465c4762a1bSJed Brown         njac21=-two*u1*u2;
466c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
467c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
468c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
469c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
470c4762a1bSJed Brown       }
471c4762a1bSJed Brown 
472c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
473c4762a1bSJed Brown       if (j==0 || j==1) {
474c4762a1bSJed Brown         xt=xt+hx;
475c4762a1bSJed Brown       } else { /* if (j==2 || j==3) */
476c4762a1bSJed Brown         yt=yt+hy;
477c4762a1bSJed Brown       }
478c4762a1bSJed Brown     }
479c4762a1bSJed Brown   }
480c4762a1bSJed Brown   PetscFunctionReturn(0);
481c4762a1bSJed Brown }
482c4762a1bSJed Brown 
483c4762a1bSJed Brown /* ------------------------------------------------------------------- */
484c4762a1bSJed Brown /*
485c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
486c4762a1bSJed Brown 
487c4762a1bSJed Brown    Input Parameters:
488c4762a1bSJed Brown .  user - user-defined application context
489c4762a1bSJed Brown .  X - vector for initial guess
490c4762a1bSJed Brown 
491c4762a1bSJed Brown    Output Parameters:
492c4762a1bSJed Brown .  X - newly computed initial guess
493c4762a1bSJed Brown */
494c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
495c4762a1bSJed Brown {
496c4762a1bSJed Brown   PetscInt       start=-1,i,j;
497c4762a1bSJed Brown   PetscScalar    zero=0.0;
498c4762a1bSJed Brown   PetscBool      flg;
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   PetscFunctionBegin;
5019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
502c4762a1bSJed Brown 
503c4762a1bSJed Brown   if (flg && start==0) { /* The zero vector is reasonable */
5049566063dSJacob Faibussowitsch     PetscCall(VecSet(X, zero));
505c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
506c4762a1bSJed Brown     PetscInt    row;
507c4762a1bSJed Brown     PetscInt    mx=user->mx,my=user->my;
508c4762a1bSJed Brown     PetscScalar *x;
509c4762a1bSJed Brown 
510c4762a1bSJed Brown     /* Get pointers to vector data */
5119566063dSJacob Faibussowitsch     PetscCall(VecGetArray(X,&x));
512c4762a1bSJed Brown 
513c4762a1bSJed Brown     /* Perform local computations */
514c4762a1bSJed Brown     for (j=0; j<my; j++) {
515c4762a1bSJed Brown       for (i=0; i< mx; i++) {
516c4762a1bSJed Brown         row=(j)*mx + (i);
517c4762a1bSJed 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;
518c4762a1bSJed Brown       }
519c4762a1bSJed Brown     }
520c4762a1bSJed Brown 
521c4762a1bSJed Brown     /* Restore vectors */
5229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(X,&x));
523c4762a1bSJed Brown   }
524c4762a1bSJed Brown   PetscFunctionReturn(0);
525c4762a1bSJed Brown }
526c4762a1bSJed Brown 
527c4762a1bSJed Brown /*TEST
528c4762a1bSJed Brown 
529c4762a1bSJed Brown    build:
530c4762a1bSJed Brown       requires: !complex
531c4762a1bSJed Brown 
532c4762a1bSJed Brown    test:
533c4762a1bSJed Brown       args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5
534c4762a1bSJed Brown       requires: !single
535c4762a1bSJed Brown 
536c4762a1bSJed Brown    test:
537c4762a1bSJed Brown       suffix: 2
538c4762a1bSJed Brown       args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5
539c4762a1bSJed Brown 
540c4762a1bSJed Brown TEST*/
541