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