xref: /petsc/src/tao/bound/tutorials/plate2.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown #include <petscdmda.h>
2*c4762a1bSJed Brown #include <petsctao.h>
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown static  char help[] =
5*c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\
6*c4762a1bSJed Brown solve an unconstrained minimization problem.  This example is based on a \n\
7*c4762a1bSJed Brown problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain, \n\
8*c4762a1bSJed Brown boundary values along the edges of the domain, and a plate represented by \n\
9*c4762a1bSJed Brown lower boundary conditions, the objective is to find the\n\
10*c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
11*c4762a1bSJed Brown The command line options are:\n\
12*c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
13*c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
14*c4762a1bSJed Brown   -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\
15*c4762a1bSJed Brown   -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\
16*c4762a1bSJed Brown   -bheight <ht>, where <ht> = height of the plate\n\
17*c4762a1bSJed Brown   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
18*c4762a1bSJed Brown                for an average of the boundary conditions\n\n";
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown /*T
21*c4762a1bSJed Brown    Concepts: TAO^Solving a bound constrained minimization problem
22*c4762a1bSJed Brown    Routines: TaoCreate();
23*c4762a1bSJed Brown    Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
24*c4762a1bSJed Brown    Routines: TaoSetHessianRoutine();
25*c4762a1bSJed Brown    Routines: TaoSetInitialVector();
26*c4762a1bSJed Brown    Routines: TaoSetVariableBounds();
27*c4762a1bSJed Brown    Routines: TaoSetFromOptions();
28*c4762a1bSJed Brown    Routines: TaoSolve(); TaoView();
29*c4762a1bSJed Brown    Routines: TaoDestroy();
30*c4762a1bSJed Brown    Processors: n
31*c4762a1bSJed Brown T*/
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown /*
37*c4762a1bSJed Brown    User-defined application context - contains data needed by the
38*c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
39*c4762a1bSJed Brown    FormHessian().
40*c4762a1bSJed Brown */
41*c4762a1bSJed Brown typedef struct {
42*c4762a1bSJed Brown   /* problem parameters */
43*c4762a1bSJed Brown   PetscReal      bheight;                  /* Height of plate under the surface */
44*c4762a1bSJed Brown   PetscInt       mx, my;                   /* discretization in x, y directions */
45*c4762a1bSJed Brown   PetscInt       bmx,bmy;                  /* Size of plate under the surface */
46*c4762a1bSJed Brown   Vec            Bottom, Top, Left, Right; /* boundary values */
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   /* Working space */
49*c4762a1bSJed Brown   Vec         localX, localV;           /* ghosted local vector */
50*c4762a1bSJed Brown   DM          dm;                       /* distributed array data structure */
51*c4762a1bSJed Brown   Mat         H;
52*c4762a1bSJed Brown } AppCtx;
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown /* -------- User-defined Routines --------- */
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
57*c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
58*c4762a1bSJed Brown static PetscErrorCode MSA_Plate(Vec,Vec,void*);
59*c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
60*c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown /* For testing matrix free submatrices */
63*c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao,Vec,Mat, Mat,void*);
64*c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat,Vec,Vec);
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown int main( int argc, char **argv )
67*c4762a1bSJed Brown {
68*c4762a1bSJed Brown   PetscErrorCode         ierr;                 /* used to check for functions returning nonzeros */
69*c4762a1bSJed Brown   PetscInt               Nx, Ny;               /* number of processors in x- and y- directions */
70*c4762a1bSJed Brown   PetscInt               m, N;                 /* number of local and global elements in vectors */
71*c4762a1bSJed Brown   Vec                    x,xl,xu;               /* solution vector  and bounds*/
72*c4762a1bSJed Brown   PetscBool              flg;                /* A return variable when checking for user options */
73*c4762a1bSJed Brown   Tao                    tao;                  /* Tao solver context */
74*c4762a1bSJed Brown   ISLocalToGlobalMapping isltog;   /* local-to-global mapping object */
75*c4762a1bSJed Brown   Mat                    H_shell;                  /* to test matrix-free submatrices */
76*c4762a1bSJed Brown   AppCtx                 user;                 /* user-defined work context */
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   /* Initialize PETSc, TAO */
79*c4762a1bSJed Brown   ierr = PetscInitialize( &argc, &argv,(char *)0,help );if (ierr) return ierr;
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   /* Specify default dimension of the problem */
82*c4762a1bSJed Brown   user.mx = 10; user.my = 10; user.bheight=0.1;
83*c4762a1bSJed Brown 
84*c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
85*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-bheight",&user.bheight,&flg);CHKERRQ(ierr);
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown   user.bmx = user.mx/2; user.bmy = user.my/2;
90*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-bmx",&user.bmx,&flg);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-bmy",&user.bmy,&flg);CHKERRQ(ierr);
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n");CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"mx:%D, my:%D, bmx:%D, bmy:%D, height:%g\n",user.mx,user.my,user.bmx,user.bmy,(double)user.bheight);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   /* Calculate any derived values from parameters */
97*c4762a1bSJed Brown   N    = user.mx*user.my;
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown   /* Let Petsc determine the dimensions of the local vectors */
100*c4762a1bSJed Brown   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown   /*
103*c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
104*c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
105*c4762a1bSJed Brown      the distributed array, Create the vectors.
106*c4762a1bSJed Brown   */
107*c4762a1bSJed Brown   ierr = DMDACreate2d(MPI_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = DMSetFromOptions(user.dm);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = DMSetUp(user.dm);CHKERRQ(ierr);
110*c4762a1bSJed Brown   /*
111*c4762a1bSJed Brown      Extract global and local vectors from DM; The local vectors are
112*c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
113*c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
114*c4762a1bSJed Brown      the same types.
115*c4762a1bSJed Brown   */
116*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(user.dm,&x);CHKERRQ(ierr); /* Solution */
117*c4762a1bSJed Brown   ierr = DMCreateLocalVector(user.dm,&user.localX);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = VecDuplicate(user.localX,&user.localV);CHKERRQ(ierr);
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   ierr = VecDuplicate(x,&xl);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = VecDuplicate(x,&xu);CHKERRQ(ierr);
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown   /* The TAO code begins here */
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   /*
126*c4762a1bSJed Brown      Create TAO solver and set desired solution method
127*c4762a1bSJed Brown      The method must either be TAOTRON or TAOBLMVM
128*c4762a1bSJed Brown      If TAOBLMVM is used, then hessian function is not called.
129*c4762a1bSJed Brown   */
130*c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   /* Set initial solution guess; */
134*c4762a1bSJed Brown   ierr = MSA_BoundaryConditions(&user);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = MSA_InitialPoint(&user,x);CHKERRQ(ierr);
136*c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown   /* Set routines for function, gradient and hessian evaluation */
139*c4762a1bSJed Brown   ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*) &user);CHKERRQ(ierr);
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown   ierr = VecGetLocalSize(x,&m);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H));CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown   ierr = DMGetLocalToGlobalMapping(user.dm,&isltog);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = MatSetLocalToGlobalMapping(user.H,isltog,isltog);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-matrixfree",&flg);CHKERRQ(ierr);
148*c4762a1bSJed Brown   if (flg) {
149*c4762a1bSJed Brown       ierr = MatCreateShell(PETSC_COMM_WORLD,m,m,N,N,(void*)&user,&H_shell);CHKERRQ(ierr);
150*c4762a1bSJed Brown       ierr = MatShellSetOperation(H_shell,MATOP_MULT,(void(*)(void))MyMatMult);CHKERRQ(ierr);
151*c4762a1bSJed Brown       ierr = MatSetOption(H_shell,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
152*c4762a1bSJed Brown       ierr = TaoSetHessianRoutine(tao,H_shell,H_shell,MatrixFreeHessian,(void*)&user);CHKERRQ(ierr);
153*c4762a1bSJed Brown   } else {
154*c4762a1bSJed Brown       ierr = TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);CHKERRQ(ierr);
155*c4762a1bSJed Brown   }
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown   /* Set Variable bounds */
158*c4762a1bSJed Brown   ierr = MSA_Plate(xl,xu,(void*)&user);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,xl,xu);CHKERRQ(ierr);
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown   /* Check for any tao command line options */
162*c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
165*c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown   ierr = TaoView(tao,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
168*c4762a1bSJed Brown 
169*c4762a1bSJed Brown   /* Free TAO data structures */
170*c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
171*c4762a1bSJed Brown 
172*c4762a1bSJed Brown   /* Free PETSc data structures */
173*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = VecDestroy(&xl);CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = VecDestroy(&xu);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = MatDestroy(&user.H);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = VecDestroy(&user.localX);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = VecDestroy(&user.localV);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = VecDestroy(&user.Bottom);CHKERRQ(ierr);
180*c4762a1bSJed Brown   ierr = VecDestroy(&user.Top);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ierr = VecDestroy(&user.Left);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = VecDestroy(&user.Right);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
184*c4762a1bSJed Brown   if (flg) {
185*c4762a1bSJed Brown     ierr = MatDestroy(&H_shell);CHKERRQ(ierr);
186*c4762a1bSJed Brown   }
187*c4762a1bSJed Brown   ierr = PetscFinalize();
188*c4762a1bSJed Brown   return ierr;
189*c4762a1bSJed Brown }
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown     Input Parameters:
194*c4762a1bSJed Brown .   tao     - the Tao context
195*c4762a1bSJed Brown .   X      - input vector
196*c4762a1bSJed Brown .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown     Output Parameters:
199*c4762a1bSJed Brown .   fcn     - the function value
200*c4762a1bSJed Brown .   G      - vector containing the newly evaluated gradient
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown    Notes:
203*c4762a1bSJed Brown    In this case, we discretize the domain and Create triangles. The
204*c4762a1bSJed Brown    surface of each triangle is planar, whose surface area can be easily
205*c4762a1bSJed Brown    computed.  The total surface area is found by sweeping through the grid
206*c4762a1bSJed Brown    and computing the surface area of the two triangles that have their
207*c4762a1bSJed Brown    right angle at the grid point.  The diagonal line segments on the
208*c4762a1bSJed Brown    grid that define the triangles run from top left to lower right.
209*c4762a1bSJed Brown    The numbering of points starts at the lower left and runs left to
210*c4762a1bSJed Brown    right, then bottom to top.
211*c4762a1bSJed Brown */
212*c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G,void *userCtx)
213*c4762a1bSJed Brown {
214*c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) userCtx;
215*c4762a1bSJed Brown   PetscErrorCode ierr;
216*c4762a1bSJed Brown   PetscInt       i,j,row;
217*c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
218*c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
219*c4762a1bSJed Brown   PetscReal      ft=0;
220*c4762a1bSJed Brown   PetscReal      zero=0.0;
221*c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
222*c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
223*c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
224*c4762a1bSJed Brown   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
225*c4762a1bSJed Brown   PetscReal      *g, *x,*left,*right,*bottom,*top;
226*c4762a1bSJed Brown   Vec            localX = user->localX, localG = user->localV;
227*c4762a1bSJed Brown 
228*c4762a1bSJed Brown   /* Get local mesh boundaries */
229*c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
230*c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown   /* Scatter ghost points to local vector */
233*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
235*c4762a1bSJed Brown 
236*c4762a1bSJed Brown   /* Initialize vector to zero */
237*c4762a1bSJed Brown   ierr = VecSet(localG, zero);CHKERRQ(ierr);
238*c4762a1bSJed Brown 
239*c4762a1bSJed Brown   /* Get pointers to vector data */
240*c4762a1bSJed Brown   ierr = VecGetArray(localX,&x);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = VecGetArray(localG,&g);CHKERRQ(ierr);
242*c4762a1bSJed Brown   ierr = VecGetArray(user->Top,&top);CHKERRQ(ierr);
243*c4762a1bSJed Brown   ierr = VecGetArray(user->Bottom,&bottom);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = VecGetArray(user->Left,&left);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = VecGetArray(user->Right,&right);CHKERRQ(ierr);
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown   /* Compute function over the locally owned part of the mesh */
248*c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++){
249*c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++){
250*c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
251*c4762a1bSJed Brown 
252*c4762a1bSJed Brown       xc = x[row];
253*c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
254*c4762a1bSJed Brown 
255*c4762a1bSJed Brown       if (i==0){ /* left side */
256*c4762a1bSJed Brown         xl= left[j-ys+1];
257*c4762a1bSJed Brown         xlt = left[j-ys+2];
258*c4762a1bSJed Brown       } else {
259*c4762a1bSJed Brown         xl = x[row-1];
260*c4762a1bSJed Brown       }
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown       if (j==0){ /* bottom side */
263*c4762a1bSJed Brown         xb=bottom[i-xs+1];
264*c4762a1bSJed Brown         xrb = bottom[i-xs+2];
265*c4762a1bSJed Brown       } else {
266*c4762a1bSJed Brown         xb = x[row-gxm];
267*c4762a1bSJed Brown       }
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown       if (i+1 == gxs+gxm){ /* right side */
270*c4762a1bSJed Brown         xr=right[j-ys+1];
271*c4762a1bSJed Brown         xrb = right[j-ys];
272*c4762a1bSJed Brown       } else {
273*c4762a1bSJed Brown         xr = x[row+1];
274*c4762a1bSJed Brown       }
275*c4762a1bSJed Brown 
276*c4762a1bSJed Brown       if (j+1==gys+gym){ /* top side */
277*c4762a1bSJed Brown         xt=top[i-xs+1];
278*c4762a1bSJed Brown         xlt = top[i-xs];
279*c4762a1bSJed Brown       }else {
280*c4762a1bSJed Brown         xt = x[row+gxm];
281*c4762a1bSJed Brown       }
282*c4762a1bSJed Brown 
283*c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym){
284*c4762a1bSJed Brown         xlt = x[row-1+gxm];
285*c4762a1bSJed Brown       }
286*c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm){
287*c4762a1bSJed Brown         xrb = x[row+1-gxm];
288*c4762a1bSJed Brown       }
289*c4762a1bSJed Brown 
290*c4762a1bSJed Brown       d1 = (xc-xl);
291*c4762a1bSJed Brown       d2 = (xc-xr);
292*c4762a1bSJed Brown       d3 = (xc-xt);
293*c4762a1bSJed Brown       d4 = (xc-xb);
294*c4762a1bSJed Brown       d5 = (xr-xrb);
295*c4762a1bSJed Brown       d6 = (xrb-xb);
296*c4762a1bSJed Brown       d7 = (xlt-xl);
297*c4762a1bSJed Brown       d8 = (xt-xlt);
298*c4762a1bSJed Brown 
299*c4762a1bSJed Brown       df1dxc = d1*hydhx;
300*c4762a1bSJed Brown       df2dxc = ( d1*hydhx + d4*hxdhy );
301*c4762a1bSJed Brown       df3dxc = d3*hxdhy;
302*c4762a1bSJed Brown       df4dxc = ( d2*hydhx + d3*hxdhy );
303*c4762a1bSJed Brown       df5dxc = d2*hydhx;
304*c4762a1bSJed Brown       df6dxc = d4*hxdhy;
305*c4762a1bSJed Brown 
306*c4762a1bSJed Brown       d1 *= rhx;
307*c4762a1bSJed Brown       d2 *= rhx;
308*c4762a1bSJed Brown       d3 *= rhy;
309*c4762a1bSJed Brown       d4 *= rhy;
310*c4762a1bSJed Brown       d5 *= rhy;
311*c4762a1bSJed Brown       d6 *= rhx;
312*c4762a1bSJed Brown       d7 *= rhy;
313*c4762a1bSJed Brown       d8 *= rhx;
314*c4762a1bSJed Brown 
315*c4762a1bSJed Brown       f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7);
316*c4762a1bSJed Brown       f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
317*c4762a1bSJed Brown       f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8);
318*c4762a1bSJed Brown       f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
319*c4762a1bSJed Brown       f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5);
320*c4762a1bSJed Brown       f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6);
321*c4762a1bSJed Brown 
322*c4762a1bSJed Brown       ft = ft + (f2 + f4);
323*c4762a1bSJed Brown 
324*c4762a1bSJed Brown       df1dxc /= f1;
325*c4762a1bSJed Brown       df2dxc /= f2;
326*c4762a1bSJed Brown       df3dxc /= f3;
327*c4762a1bSJed Brown       df4dxc /= f4;
328*c4762a1bSJed Brown       df5dxc /= f5;
329*c4762a1bSJed Brown       df6dxc /= f6;
330*c4762a1bSJed Brown 
331*c4762a1bSJed Brown       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc ) * 0.5;
332*c4762a1bSJed Brown 
333*c4762a1bSJed Brown     }
334*c4762a1bSJed Brown   }
335*c4762a1bSJed Brown 
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
338*c4762a1bSJed Brown   if (xs==0){ /* left side */
339*c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++){
340*c4762a1bSJed Brown       d3=(left[j-ys+1] - left[j-ys+2])*rhy;
341*c4762a1bSJed Brown       d2=(left[j-ys+1] - x[(j-gys)*gxm])*rhx;
342*c4762a1bSJed Brown       ft = ft+PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
343*c4762a1bSJed Brown     }
344*c4762a1bSJed Brown   }
345*c4762a1bSJed Brown   if (ys==0){ /* bottom side */
346*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++){
347*c4762a1bSJed Brown       d2=(bottom[i+1-xs]-bottom[i-xs+2])*rhx;
348*c4762a1bSJed Brown       d3=(bottom[i-xs+1]-x[i-gxs])*rhy;
349*c4762a1bSJed Brown       ft = ft+PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
350*c4762a1bSJed Brown     }
351*c4762a1bSJed Brown   }
352*c4762a1bSJed Brown 
353*c4762a1bSJed Brown   if (xs+xm==mx){ /* right side */
354*c4762a1bSJed Brown     for (j=ys; j< ys+ym; j++){
355*c4762a1bSJed Brown       d1=(x[(j+1-gys)*gxm-1]-right[j-ys+1])*rhx;
356*c4762a1bSJed Brown       d4=(right[j-ys]-right[j-ys+1])*rhy;
357*c4762a1bSJed Brown       ft = ft+PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
358*c4762a1bSJed Brown     }
359*c4762a1bSJed Brown   }
360*c4762a1bSJed Brown   if (ys+ym==my){ /* top side */
361*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++){
362*c4762a1bSJed Brown       d1=(x[(gym-1)*gxm + i-gxs] - top[i-xs+1])*rhy;
363*c4762a1bSJed Brown       d4=(top[i-xs+1] - top[i-xs])*rhx;
364*c4762a1bSJed Brown       ft = ft+PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
365*c4762a1bSJed Brown     }
366*c4762a1bSJed Brown   }
367*c4762a1bSJed Brown 
368*c4762a1bSJed Brown   if (ys==0 && xs==0){
369*c4762a1bSJed Brown     d1=(left[0]-left[1])*rhy;
370*c4762a1bSJed Brown     d2=(bottom[0]-bottom[1])*rhx;
371*c4762a1bSJed Brown     ft +=PetscSqrtScalar( 1.0 + d1*d1 + d2*d2);
372*c4762a1bSJed Brown   }
373*c4762a1bSJed Brown   if (ys+ym == my && xs+xm == mx){
374*c4762a1bSJed Brown     d1=(right[ym+1] - right[ym])*rhy;
375*c4762a1bSJed Brown     d2=(top[xm+1] - top[xm])*rhx;
376*c4762a1bSJed Brown     ft +=PetscSqrtScalar( 1.0 + d1*d1 + d2*d2);
377*c4762a1bSJed Brown   }
378*c4762a1bSJed Brown 
379*c4762a1bSJed Brown   ft=ft*area;
380*c4762a1bSJed Brown   ierr = MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD);CHKERRQ(ierr);
381*c4762a1bSJed Brown 
382*c4762a1bSJed Brown 
383*c4762a1bSJed Brown   /* Restore vectors */
384*c4762a1bSJed Brown   ierr = VecRestoreArray(localX,&x);CHKERRQ(ierr);
385*c4762a1bSJed Brown   ierr = VecRestoreArray(localG,&g);CHKERRQ(ierr);
386*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Left,&left);CHKERRQ(ierr);
387*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Top,&top);CHKERRQ(ierr);
388*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Bottom,&bottom);CHKERRQ(ierr);
389*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Right,&right);CHKERRQ(ierr);
390*c4762a1bSJed Brown 
391*c4762a1bSJed Brown   /* Scatter values to global vector */
392*c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G);CHKERRQ(ierr);
393*c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G);CHKERRQ(ierr);
394*c4762a1bSJed Brown 
395*c4762a1bSJed Brown   ierr = PetscLogFlops(70*xm*ym);CHKERRQ(ierr);
396*c4762a1bSJed Brown 
397*c4762a1bSJed Brown   return 0;
398*c4762a1bSJed Brown }
399*c4762a1bSJed Brown 
400*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
401*c4762a1bSJed Brown /*
402*c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
403*c4762a1bSJed Brown 
404*c4762a1bSJed Brown    Input Parameters:
405*c4762a1bSJed Brown .  tao  - the Tao context
406*c4762a1bSJed Brown .  x    - input vector
407*c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetHessianRoutine()
408*c4762a1bSJed Brown 
409*c4762a1bSJed Brown    Output Parameters:
410*c4762a1bSJed Brown .  A    - Hessian matrix
411*c4762a1bSJed Brown .  B    - optionally different preconditioning matrix
412*c4762a1bSJed Brown 
413*c4762a1bSJed Brown    Notes:
414*c4762a1bSJed Brown    Due to mesh point reordering with DMs, we must always work
415*c4762a1bSJed Brown    with the local mesh points, and then transform them to the new
416*c4762a1bSJed Brown    global numbering with the local-to-global mapping.  We cannot work
417*c4762a1bSJed Brown    directly with the global numbers for the original uniprocessor mesh!
418*c4762a1bSJed Brown 
419*c4762a1bSJed Brown    Two methods are available for imposing this transformation
420*c4762a1bSJed Brown    when setting matrix entries:
421*c4762a1bSJed Brown      (A) MatSetValuesLocal(), using the local ordering (including
422*c4762a1bSJed Brown          ghost points!)
423*c4762a1bSJed Brown          - Do the following two steps once, before calling TaoSolve()
424*c4762a1bSJed Brown            - Use DMGetISLocalToGlobalMapping() to extract the
425*c4762a1bSJed Brown              local-to-global map from the DM
426*c4762a1bSJed Brown            - Associate this map with the matrix by calling
427*c4762a1bSJed Brown              MatSetLocalToGlobalMapping()
428*c4762a1bSJed Brown          - Then set matrix entries using the local ordering
429*c4762a1bSJed Brown            by calling MatSetValuesLocal()
430*c4762a1bSJed Brown      (B) MatSetValues(), using the global ordering
431*c4762a1bSJed Brown          - Use DMGetGlobalIndices() to extract the local-to-global map
432*c4762a1bSJed Brown          - Then apply this map explicitly yourself
433*c4762a1bSJed Brown          - Set matrix entries using the global ordering by calling
434*c4762a1bSJed Brown            MatSetValues()
435*c4762a1bSJed Brown    Option (A) seems cleaner/easier in many cases, and is the procedure
436*c4762a1bSJed Brown    used in this example.
437*c4762a1bSJed Brown */
438*c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat Hptr, Mat Hessian, void *ptr)
439*c4762a1bSJed Brown {
440*c4762a1bSJed Brown   PetscErrorCode ierr;
441*c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
442*c4762a1bSJed Brown   PetscInt       i,j,k,row;
443*c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my;
444*c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym,col[7];
445*c4762a1bSJed Brown   PetscReal      hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
446*c4762a1bSJed Brown   PetscReal      rhx=mx+1, rhy=my+1;
447*c4762a1bSJed Brown   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
448*c4762a1bSJed Brown   PetscReal      hl,hr,ht,hb,hc,htl,hbr;
449*c4762a1bSJed Brown   PetscReal      *x,*left,*right,*bottom,*top;
450*c4762a1bSJed Brown   PetscReal      v[7];
451*c4762a1bSJed Brown   Vec            localX = user->localX;
452*c4762a1bSJed Brown   PetscBool      assembled;
453*c4762a1bSJed Brown 
454*c4762a1bSJed Brown 
455*c4762a1bSJed Brown   /* Set various matrix options */
456*c4762a1bSJed Brown   ierr = MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
457*c4762a1bSJed Brown 
458*c4762a1bSJed Brown   /* Initialize matrix entries to zero */
459*c4762a1bSJed Brown   ierr = MatAssembled(Hessian,&assembled);CHKERRQ(ierr);
460*c4762a1bSJed Brown   if (assembled){ierr = MatZeroEntries(Hessian);CHKERRQ(ierr);}
461*c4762a1bSJed Brown 
462*c4762a1bSJed Brown   /* Get local mesh boundaries */
463*c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
464*c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
465*c4762a1bSJed Brown 
466*c4762a1bSJed Brown   /* Scatter ghost points to local vector */
467*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
468*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
469*c4762a1bSJed Brown 
470*c4762a1bSJed Brown   /* Get pointers to vector data */
471*c4762a1bSJed Brown   ierr = VecGetArray(localX,&x);CHKERRQ(ierr);
472*c4762a1bSJed Brown   ierr = VecGetArray(user->Top,&top);CHKERRQ(ierr);
473*c4762a1bSJed Brown   ierr = VecGetArray(user->Bottom,&bottom);CHKERRQ(ierr);
474*c4762a1bSJed Brown   ierr = VecGetArray(user->Left,&left);CHKERRQ(ierr);
475*c4762a1bSJed Brown   ierr = VecGetArray(user->Right,&right);CHKERRQ(ierr);
476*c4762a1bSJed Brown 
477*c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
478*c4762a1bSJed Brown 
479*c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++){
480*c4762a1bSJed Brown 
481*c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++){
482*c4762a1bSJed Brown 
483*c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
484*c4762a1bSJed Brown 
485*c4762a1bSJed Brown       xc = x[row];
486*c4762a1bSJed Brown       xlt=xrb=xl=xr=xb=xt=xc;
487*c4762a1bSJed Brown 
488*c4762a1bSJed Brown       /* Left side */
489*c4762a1bSJed Brown       if (i==gxs){
490*c4762a1bSJed Brown         xl= left[j-ys+1];
491*c4762a1bSJed Brown         xlt = left[j-ys+2];
492*c4762a1bSJed Brown       } else {
493*c4762a1bSJed Brown         xl = x[row-1];
494*c4762a1bSJed Brown       }
495*c4762a1bSJed Brown 
496*c4762a1bSJed Brown       if (j==gys){
497*c4762a1bSJed Brown         xb=bottom[i-xs+1];
498*c4762a1bSJed Brown         xrb = bottom[i-xs+2];
499*c4762a1bSJed Brown       } else {
500*c4762a1bSJed Brown         xb = x[row-gxm];
501*c4762a1bSJed Brown       }
502*c4762a1bSJed Brown 
503*c4762a1bSJed Brown       if (i+1 == gxs+gxm){
504*c4762a1bSJed Brown         xr=right[j-ys+1];
505*c4762a1bSJed Brown         xrb = right[j-ys];
506*c4762a1bSJed Brown       } else {
507*c4762a1bSJed Brown         xr = x[row+1];
508*c4762a1bSJed Brown       }
509*c4762a1bSJed Brown 
510*c4762a1bSJed Brown       if (j+1==gys+gym){
511*c4762a1bSJed Brown         xt=top[i-xs+1];
512*c4762a1bSJed Brown         xlt = top[i-xs];
513*c4762a1bSJed Brown       }else {
514*c4762a1bSJed Brown         xt = x[row+gxm];
515*c4762a1bSJed Brown       }
516*c4762a1bSJed Brown 
517*c4762a1bSJed Brown       if (i>gxs && j+1<gys+gym){
518*c4762a1bSJed Brown         xlt = x[row-1+gxm];
519*c4762a1bSJed Brown       }
520*c4762a1bSJed Brown       if (j>gys && i+1<gxs+gxm){
521*c4762a1bSJed Brown         xrb = x[row+1-gxm];
522*c4762a1bSJed Brown       }
523*c4762a1bSJed Brown 
524*c4762a1bSJed Brown 
525*c4762a1bSJed Brown       d1 = (xc-xl)*rhx;
526*c4762a1bSJed Brown       d2 = (xc-xr)*rhx;
527*c4762a1bSJed Brown       d3 = (xc-xt)*rhy;
528*c4762a1bSJed Brown       d4 = (xc-xb)*rhy;
529*c4762a1bSJed Brown       d5 = (xrb-xr)*rhy;
530*c4762a1bSJed Brown       d6 = (xrb-xb)*rhx;
531*c4762a1bSJed Brown       d7 = (xlt-xl)*rhy;
532*c4762a1bSJed Brown       d8 = (xlt-xt)*rhx;
533*c4762a1bSJed Brown 
534*c4762a1bSJed Brown       f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7);
535*c4762a1bSJed Brown       f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
536*c4762a1bSJed Brown       f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8);
537*c4762a1bSJed Brown       f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
538*c4762a1bSJed Brown       f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5);
539*c4762a1bSJed Brown       f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6);
540*c4762a1bSJed Brown 
541*c4762a1bSJed Brown 
542*c4762a1bSJed Brown       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
543*c4762a1bSJed Brown         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
544*c4762a1bSJed Brown       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
545*c4762a1bSJed Brown         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
546*c4762a1bSJed Brown       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
547*c4762a1bSJed Brown         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
548*c4762a1bSJed Brown       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
549*c4762a1bSJed Brown         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
550*c4762a1bSJed Brown 
551*c4762a1bSJed Brown       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
552*c4762a1bSJed Brown       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
553*c4762a1bSJed Brown 
554*c4762a1bSJed Brown       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
555*c4762a1bSJed Brown         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
556*c4762a1bSJed Brown         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
557*c4762a1bSJed Brown         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
558*c4762a1bSJed Brown 
559*c4762a1bSJed Brown       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;
560*c4762a1bSJed Brown 
561*c4762a1bSJed Brown       k=0;
562*c4762a1bSJed Brown       if (j>0){
563*c4762a1bSJed Brown         v[k]=hb; col[k]=row - gxm; k++;
564*c4762a1bSJed Brown       }
565*c4762a1bSJed Brown 
566*c4762a1bSJed Brown       if (j>0 && i < mx -1){
567*c4762a1bSJed Brown         v[k]=hbr; col[k]=row - gxm+1; k++;
568*c4762a1bSJed Brown       }
569*c4762a1bSJed Brown 
570*c4762a1bSJed Brown       if (i>0){
571*c4762a1bSJed Brown         v[k]= hl; col[k]=row - 1; k++;
572*c4762a1bSJed Brown       }
573*c4762a1bSJed Brown 
574*c4762a1bSJed Brown       v[k]= hc; col[k]=row; k++;
575*c4762a1bSJed Brown 
576*c4762a1bSJed Brown       if (i < mx-1 ){
577*c4762a1bSJed Brown         v[k]= hr; col[k]=row+1; k++;
578*c4762a1bSJed Brown       }
579*c4762a1bSJed Brown 
580*c4762a1bSJed Brown       if (i>0 && j < my-1 ){
581*c4762a1bSJed Brown         v[k]= htl; col[k] = row+gxm-1; k++;
582*c4762a1bSJed Brown       }
583*c4762a1bSJed Brown 
584*c4762a1bSJed Brown       if (j < my-1 ){
585*c4762a1bSJed Brown         v[k]= ht; col[k] = row+gxm; k++;
586*c4762a1bSJed Brown       }
587*c4762a1bSJed Brown 
588*c4762a1bSJed Brown       /*
589*c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
590*c4762a1bSJed Brown          earlier, in the main routine.
591*c4762a1bSJed Brown       */
592*c4762a1bSJed Brown       ierr = MatSetValuesLocal(Hessian,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr);
593*c4762a1bSJed Brown 
594*c4762a1bSJed Brown     }
595*c4762a1bSJed Brown   }
596*c4762a1bSJed Brown 
597*c4762a1bSJed Brown   /* Restore vectors */
598*c4762a1bSJed Brown   ierr = VecRestoreArray(localX,&x);CHKERRQ(ierr);
599*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Left,&left);CHKERRQ(ierr);
600*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Top,&top);CHKERRQ(ierr);
601*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Bottom,&bottom);CHKERRQ(ierr);
602*c4762a1bSJed Brown   ierr = VecRestoreArray(user->Right,&right);CHKERRQ(ierr);
603*c4762a1bSJed Brown 
604*c4762a1bSJed Brown   /* Assemble the matrix */
605*c4762a1bSJed Brown   ierr = MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
606*c4762a1bSJed Brown   ierr = MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
607*c4762a1bSJed Brown 
608*c4762a1bSJed Brown   ierr = PetscLogFlops(199*xm*ym);CHKERRQ(ierr);
609*c4762a1bSJed Brown   return 0;
610*c4762a1bSJed Brown }
611*c4762a1bSJed Brown 
612*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
613*c4762a1bSJed Brown /*
614*c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
615*c4762a1bSJed Brown    the region.
616*c4762a1bSJed Brown 
617*c4762a1bSJed Brown    Input Parameter:
618*c4762a1bSJed Brown .  user - user-defined application context
619*c4762a1bSJed Brown 
620*c4762a1bSJed Brown    Output Parameter:
621*c4762a1bSJed Brown .  user - user-defined application context
622*c4762a1bSJed Brown */
623*c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
624*c4762a1bSJed Brown {
625*c4762a1bSJed Brown   int        ierr;
626*c4762a1bSJed Brown   PetscInt   i,j,k,maxits=5,limit=0;
627*c4762a1bSJed Brown   PetscInt   xs,ys,xm,ym,gxs,gys,gxm,gym;
628*c4762a1bSJed Brown   PetscInt   mx=user->mx,my=user->my;
629*c4762a1bSJed Brown   PetscInt   bsize=0, lsize=0, tsize=0, rsize=0;
630*c4762a1bSJed Brown   PetscReal  one=1.0, two=2.0, three=3.0, scl=1.0, tol=1e-10;
631*c4762a1bSJed Brown   PetscReal  fnorm,det,hx,hy,xt=0,yt=0;
632*c4762a1bSJed Brown   PetscReal  u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
633*c4762a1bSJed Brown   PetscReal  b=-0.5, t=0.5, l=-0.5, r=0.5;
634*c4762a1bSJed Brown   PetscReal  *boundary;
635*c4762a1bSJed Brown   PetscBool  flg;
636*c4762a1bSJed Brown   Vec        Bottom,Top,Right,Left;
637*c4762a1bSJed Brown 
638*c4762a1bSJed Brown   /* Get local mesh boundaries */
639*c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
640*c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
641*c4762a1bSJed Brown 
642*c4762a1bSJed Brown 
643*c4762a1bSJed Brown   bsize=xm+2;
644*c4762a1bSJed Brown   lsize=ym+2;
645*c4762a1bSJed Brown   rsize=ym+2;
646*c4762a1bSJed Brown   tsize=xm+2;
647*c4762a1bSJed Brown 
648*c4762a1bSJed Brown   ierr = VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom);CHKERRQ(ierr);
649*c4762a1bSJed Brown   ierr = VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top);CHKERRQ(ierr);
650*c4762a1bSJed Brown   ierr = VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left);CHKERRQ(ierr);
651*c4762a1bSJed Brown   ierr = VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right);CHKERRQ(ierr);
652*c4762a1bSJed Brown 
653*c4762a1bSJed Brown   user->Top=Top;
654*c4762a1bSJed Brown   user->Left=Left;
655*c4762a1bSJed Brown   user->Bottom=Bottom;
656*c4762a1bSJed Brown   user->Right=Right;
657*c4762a1bSJed Brown 
658*c4762a1bSJed Brown   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
659*c4762a1bSJed Brown 
660*c4762a1bSJed Brown   for (j=0; j<4; j++){
661*c4762a1bSJed Brown     if (j==0){
662*c4762a1bSJed Brown       yt=b;
663*c4762a1bSJed Brown       xt=l+hx*xs;
664*c4762a1bSJed Brown       limit=bsize;
665*c4762a1bSJed Brown       VecGetArray(Bottom,&boundary);
666*c4762a1bSJed Brown     } else if (j==1){
667*c4762a1bSJed Brown       yt=t;
668*c4762a1bSJed Brown       xt=l+hx*xs;
669*c4762a1bSJed Brown       limit=tsize;
670*c4762a1bSJed Brown       VecGetArray(Top,&boundary);
671*c4762a1bSJed Brown     } else if (j==2){
672*c4762a1bSJed Brown       yt=b+hy*ys;
673*c4762a1bSJed Brown       xt=l;
674*c4762a1bSJed Brown       limit=lsize;
675*c4762a1bSJed Brown       VecGetArray(Left,&boundary);
676*c4762a1bSJed Brown     } else if (j==3){
677*c4762a1bSJed Brown       yt=b+hy*ys;
678*c4762a1bSJed Brown       xt=r;
679*c4762a1bSJed Brown       limit=rsize;
680*c4762a1bSJed Brown       VecGetArray(Right,&boundary);
681*c4762a1bSJed Brown     }
682*c4762a1bSJed Brown 
683*c4762a1bSJed Brown     for (i=0; i<limit; i++){
684*c4762a1bSJed Brown       u1=xt;
685*c4762a1bSJed Brown       u2=-yt;
686*c4762a1bSJed Brown       for (k=0; k<maxits; k++){
687*c4762a1bSJed Brown         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
688*c4762a1bSJed Brown         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
689*c4762a1bSJed Brown         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
690*c4762a1bSJed Brown         if (fnorm <= tol) break;
691*c4762a1bSJed Brown         njac11=one+u2*u2-u1*u1;
692*c4762a1bSJed Brown         njac12=two*u1*u2;
693*c4762a1bSJed Brown         njac21=-two*u1*u2;
694*c4762a1bSJed Brown         njac22=-one - u1*u1 + u2*u2;
695*c4762a1bSJed Brown         det = njac11*njac22-njac21*njac12;
696*c4762a1bSJed Brown         u1 = u1-(njac22*nf1-njac12*nf2)/det;
697*c4762a1bSJed Brown         u2 = u2-(njac11*nf2-njac21*nf1)/det;
698*c4762a1bSJed Brown       }
699*c4762a1bSJed Brown 
700*c4762a1bSJed Brown       boundary[i]=u1*u1-u2*u2;
701*c4762a1bSJed Brown       if (j==0 || j==1) {
702*c4762a1bSJed Brown         xt=xt+hx;
703*c4762a1bSJed Brown       } else if (j==2 || j==3){
704*c4762a1bSJed Brown         yt=yt+hy;
705*c4762a1bSJed Brown       }
706*c4762a1bSJed Brown     }
707*c4762a1bSJed Brown     if (j==0){
708*c4762a1bSJed Brown       ierr = VecRestoreArray(Bottom,&boundary);CHKERRQ(ierr);
709*c4762a1bSJed Brown     } else if (j==1){
710*c4762a1bSJed Brown       ierr = VecRestoreArray(Top,&boundary);CHKERRQ(ierr);
711*c4762a1bSJed Brown     } else if (j==2){
712*c4762a1bSJed Brown       ierr = VecRestoreArray(Left,&boundary);CHKERRQ(ierr);
713*c4762a1bSJed Brown     } else if (j==3){
714*c4762a1bSJed Brown       ierr = VecRestoreArray(Right,&boundary);CHKERRQ(ierr);
715*c4762a1bSJed Brown     }
716*c4762a1bSJed Brown   }
717*c4762a1bSJed Brown 
718*c4762a1bSJed Brown   /* Scale the boundary if desired */
719*c4762a1bSJed Brown 
720*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg);CHKERRQ(ierr);
721*c4762a1bSJed Brown   if (flg){
722*c4762a1bSJed Brown     ierr = VecScale(Bottom, scl);CHKERRQ(ierr);
723*c4762a1bSJed Brown   }
724*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg);CHKERRQ(ierr);
725*c4762a1bSJed Brown   if (flg){
726*c4762a1bSJed Brown     ierr = VecScale(Top, scl);CHKERRQ(ierr);
727*c4762a1bSJed Brown   }
728*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg);CHKERRQ(ierr);
729*c4762a1bSJed Brown   if (flg){
730*c4762a1bSJed Brown     ierr = VecScale(Right, scl);CHKERRQ(ierr);
731*c4762a1bSJed Brown   }
732*c4762a1bSJed Brown 
733*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg);CHKERRQ(ierr);
734*c4762a1bSJed Brown   if (flg){
735*c4762a1bSJed Brown     ierr = VecScale(Left, scl);CHKERRQ(ierr);
736*c4762a1bSJed Brown   }
737*c4762a1bSJed Brown   return 0;
738*c4762a1bSJed Brown }
739*c4762a1bSJed Brown 
740*c4762a1bSJed Brown 
741*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
742*c4762a1bSJed Brown /*
743*c4762a1bSJed Brown    MSA_Plate -  Calculates an obstacle for surface to stretch over.
744*c4762a1bSJed Brown 
745*c4762a1bSJed Brown    Input Parameter:
746*c4762a1bSJed Brown .  user - user-defined application context
747*c4762a1bSJed Brown 
748*c4762a1bSJed Brown    Output Parameter:
749*c4762a1bSJed Brown .  user - user-defined application context
750*c4762a1bSJed Brown */
751*c4762a1bSJed Brown static PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx){
752*c4762a1bSJed Brown 
753*c4762a1bSJed Brown   AppCtx         *user=(AppCtx *)ctx;
754*c4762a1bSJed Brown   PetscErrorCode ierr;
755*c4762a1bSJed Brown   PetscInt       i,j,row;
756*c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym;
757*c4762a1bSJed Brown   PetscInt       mx=user->mx, my=user->my, bmy, bmx;
758*c4762a1bSJed Brown   PetscReal      t1,t2,t3;
759*c4762a1bSJed Brown   PetscReal      *xl, lb=PETSC_NINFINITY, ub=PETSC_INFINITY;
760*c4762a1bSJed Brown   PetscBool      cylinder;
761*c4762a1bSJed Brown 
762*c4762a1bSJed Brown   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
763*c4762a1bSJed Brown   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
764*c4762a1bSJed Brown   bmy=user->bmy; bmx=user->bmx;
765*c4762a1bSJed Brown 
766*c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
767*c4762a1bSJed Brown 
768*c4762a1bSJed Brown   ierr = VecSet(XL, lb);CHKERRQ(ierr);
769*c4762a1bSJed Brown   ierr = VecSet(XU, ub);CHKERRQ(ierr);
770*c4762a1bSJed Brown 
771*c4762a1bSJed Brown   ierr = VecGetArray(XL,&xl);CHKERRQ(ierr);
772*c4762a1bSJed Brown 
773*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-cylinder",&cylinder);CHKERRQ(ierr);
774*c4762a1bSJed Brown   /* Compute the optional lower box */
775*c4762a1bSJed Brown   if (cylinder){
776*c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++){
777*c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++){
778*c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
779*c4762a1bSJed Brown         t1=(2.0*i-mx)*bmy;
780*c4762a1bSJed Brown         t2=(2.0*j-my)*bmx;
781*c4762a1bSJed Brown         t3=bmx*bmx*bmy*bmy;
782*c4762a1bSJed Brown         if ( t1*t1 + t2*t2 <= t3 ){
783*c4762a1bSJed Brown           xl[row] = user->bheight;
784*c4762a1bSJed Brown         }
785*c4762a1bSJed Brown       }
786*c4762a1bSJed Brown     }
787*c4762a1bSJed Brown   } else {
788*c4762a1bSJed Brown     /* Compute the optional lower box */
789*c4762a1bSJed Brown     for (i=xs; i< xs+xm; i++){
790*c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++){
791*c4762a1bSJed Brown         row=(j-ys)*xm + (i-xs);
792*c4762a1bSJed Brown         if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
793*c4762a1bSJed Brown             j>=(my-bmy)/2 && j<my-(my-bmy)/2 ){
794*c4762a1bSJed Brown           xl[row] = user->bheight;
795*c4762a1bSJed Brown         }
796*c4762a1bSJed Brown       }
797*c4762a1bSJed Brown     }
798*c4762a1bSJed Brown   }
799*c4762a1bSJed Brown     ierr = VecRestoreArray(XL,&xl);CHKERRQ(ierr);
800*c4762a1bSJed Brown 
801*c4762a1bSJed Brown   return 0;
802*c4762a1bSJed Brown }
803*c4762a1bSJed Brown 
804*c4762a1bSJed Brown 
805*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
806*c4762a1bSJed Brown /*
807*c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
808*c4762a1bSJed Brown 
809*c4762a1bSJed Brown    Input Parameters:
810*c4762a1bSJed Brown .  user - user-defined application context
811*c4762a1bSJed Brown .  X - vector for initial guess
812*c4762a1bSJed Brown 
813*c4762a1bSJed Brown    Output Parameters:
814*c4762a1bSJed Brown .  X - newly computed initial guess
815*c4762a1bSJed Brown */
816*c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
817*c4762a1bSJed Brown {
818*c4762a1bSJed Brown   PetscErrorCode ierr;
819*c4762a1bSJed Brown   PetscInt       start=-1,i,j;
820*c4762a1bSJed Brown   PetscReal      zero=0.0;
821*c4762a1bSJed Brown   PetscBool      flg;
822*c4762a1bSJed Brown 
823*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg);CHKERRQ(ierr);
824*c4762a1bSJed Brown   if (flg && start==0){ /* The zero vector is reasonable */
825*c4762a1bSJed Brown     ierr = VecSet(X, zero);CHKERRQ(ierr);
826*c4762a1bSJed Brown   } else if (flg && start>0){ /* Try a random start between -0.5 and 0.5 */
827*c4762a1bSJed Brown     PetscRandom rctx;  PetscReal np5=-0.5;
828*c4762a1bSJed Brown 
829*c4762a1bSJed Brown     ierr = PetscRandomCreate(MPI_COMM_WORLD,&rctx);CHKERRQ(ierr);
830*c4762a1bSJed Brown     for (i=0; i<start; i++){
831*c4762a1bSJed Brown       ierr = VecSetRandom(X, rctx);CHKERRQ(ierr);
832*c4762a1bSJed Brown     }
833*c4762a1bSJed Brown     ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
834*c4762a1bSJed Brown     ierr = VecShift(X, np5);CHKERRQ(ierr);
835*c4762a1bSJed Brown 
836*c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
837*c4762a1bSJed Brown 
838*c4762a1bSJed Brown     PetscInt row,xs,xm,gxs,gxm,ys,ym,gys,gym;
839*c4762a1bSJed Brown     PetscInt mx=user->mx,my=user->my;
840*c4762a1bSJed Brown     PetscReal *x,*left,*right,*bottom,*top;
841*c4762a1bSJed Brown     Vec    localX = user->localX;
842*c4762a1bSJed Brown 
843*c4762a1bSJed Brown     /* Get local mesh boundaries */
844*c4762a1bSJed Brown     ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
845*c4762a1bSJed Brown     ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
846*c4762a1bSJed Brown 
847*c4762a1bSJed Brown     /* Get pointers to vector data */
848*c4762a1bSJed Brown     ierr = VecGetArray(user->Top,&top);CHKERRQ(ierr);
849*c4762a1bSJed Brown     ierr = VecGetArray(user->Bottom,&bottom);CHKERRQ(ierr);
850*c4762a1bSJed Brown     ierr = VecGetArray(user->Left,&left);CHKERRQ(ierr);
851*c4762a1bSJed Brown     ierr = VecGetArray(user->Right,&right);CHKERRQ(ierr);
852*c4762a1bSJed Brown 
853*c4762a1bSJed Brown     ierr = VecGetArray(localX,&x);CHKERRQ(ierr);
854*c4762a1bSJed Brown     /* Perform local computations */
855*c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++){
856*c4762a1bSJed Brown       for (i=xs; i< xs+xm; i++){
857*c4762a1bSJed Brown         row=(j-gys)*gxm + (i-gxs);
858*c4762a1bSJed Brown         x[row] = ( (j+1)*bottom[i-xs+1]/my + (my-j+1)*top[i-xs+1]/(my+2)+
859*c4762a1bSJed Brown                    (i+1)*left[j-ys+1]/mx + (mx-i+1)*right[j-ys+1]/(mx+2))/2.0;
860*c4762a1bSJed Brown       }
861*c4762a1bSJed Brown     }
862*c4762a1bSJed Brown 
863*c4762a1bSJed Brown     /* Restore vectors */
864*c4762a1bSJed Brown     ierr = VecRestoreArray(localX,&x);CHKERRQ(ierr);
865*c4762a1bSJed Brown 
866*c4762a1bSJed Brown     ierr = VecRestoreArray(user->Left,&left);CHKERRQ(ierr);
867*c4762a1bSJed Brown     ierr = VecRestoreArray(user->Top,&top);CHKERRQ(ierr);
868*c4762a1bSJed Brown     ierr = VecRestoreArray(user->Bottom,&bottom);CHKERRQ(ierr);
869*c4762a1bSJed Brown     ierr = VecRestoreArray(user->Right,&right);CHKERRQ(ierr);
870*c4762a1bSJed Brown 
871*c4762a1bSJed Brown     /* Scatter values into global vector */
872*c4762a1bSJed Brown     ierr = DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X);CHKERRQ(ierr);
873*c4762a1bSJed Brown     ierr = DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X);CHKERRQ(ierr);
874*c4762a1bSJed Brown 
875*c4762a1bSJed Brown   }
876*c4762a1bSJed Brown   return 0;
877*c4762a1bSJed Brown }
878*c4762a1bSJed Brown 
879*c4762a1bSJed Brown /* For testing matrix free submatrices */
880*c4762a1bSJed Brown PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
881*c4762a1bSJed Brown {
882*c4762a1bSJed Brown   PetscErrorCode ierr;
883*c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
884*c4762a1bSJed Brown   PetscFunctionBegin;
885*c4762a1bSJed Brown   ierr = FormHessian(tao,x,user->H,user->H,ptr);CHKERRQ(ierr);
886*c4762a1bSJed Brown   PetscFunctionReturn(0);
887*c4762a1bSJed Brown }
888*c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
889*c4762a1bSJed Brown {
890*c4762a1bSJed Brown   PetscErrorCode ierr;
891*c4762a1bSJed Brown   void           *ptr;
892*c4762a1bSJed Brown   AppCtx         *user;
893*c4762a1bSJed Brown   PetscFunctionBegin;
894*c4762a1bSJed Brown   ierr = MatShellGetContext(H_shell,&ptr);CHKERRQ(ierr);
895*c4762a1bSJed Brown   user = (AppCtx*)ptr;
896*c4762a1bSJed Brown   ierr = MatMult(user->H,X,Y);CHKERRQ(ierr);
897*c4762a1bSJed Brown   PetscFunctionReturn(0);
898*c4762a1bSJed Brown }
899*c4762a1bSJed Brown 
900*c4762a1bSJed Brown 
901*c4762a1bSJed Brown /*TEST
902*c4762a1bSJed Brown 
903*c4762a1bSJed Brown    build:
904*c4762a1bSJed Brown       requires: !complex
905*c4762a1bSJed Brown 
906*c4762a1bSJed Brown    test:
907*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
908*c4762a1bSJed Brown       requires: !single
909*c4762a1bSJed Brown 
910*c4762a1bSJed Brown    test:
911*c4762a1bSJed Brown       suffix: 2
912*c4762a1bSJed Brown       nsize: 2
913*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
914*c4762a1bSJed Brown       requires: !single
915*c4762a1bSJed Brown 
916*c4762a1bSJed Brown    test:
917*c4762a1bSJed Brown       suffix: 3
918*c4762a1bSJed Brown       nsize: 3
919*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
920*c4762a1bSJed Brown       requires: !single
921*c4762a1bSJed Brown 
922*c4762a1bSJed Brown    test:
923*c4762a1bSJed Brown       suffix: 4
924*c4762a1bSJed Brown       nsize: 3
925*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type mask -tao_type tron -tao_gatol 1.e-5
926*c4762a1bSJed Brown       requires: !single
927*c4762a1bSJed Brown 
928*c4762a1bSJed Brown    test:
929*c4762a1bSJed Brown       suffix: 5
930*c4762a1bSJed Brown       nsize: 3
931*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -pc_type none -tao_type tron -tao_gatol 1.e-5
932*c4762a1bSJed Brown       requires: !single
933*c4762a1bSJed Brown 
934*c4762a1bSJed Brown    test:
935*c4762a1bSJed Brown       suffix: 6
936*c4762a1bSJed Brown       nsize: 3
937*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -tao_type blmvm -tao_gatol 1.e-4
938*c4762a1bSJed Brown       requires: !single
939*c4762a1bSJed Brown 
940*c4762a1bSJed Brown    test:
941*c4762a1bSJed Brown       suffix: 7
942*c4762a1bSJed Brown       nsize: 3
943*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -pc_type none -tao_type gpcg -tao_gatol 1.e-5
944*c4762a1bSJed Brown       requires: !single
945*c4762a1bSJed Brown 
946*c4762a1bSJed Brown    test:
947*c4762a1bSJed Brown       suffix: 8
948*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
949*c4762a1bSJed Brown       requires: !single
950*c4762a1bSJed Brown 
951*c4762a1bSJed Brown    test:
952*c4762a1bSJed Brown       suffix: 9
953*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
954*c4762a1bSJed Brown       requires: !single
955*c4762a1bSJed Brown 
956*c4762a1bSJed Brown    test:
957*c4762a1bSJed Brown       suffix: 10
958*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
959*c4762a1bSJed Brown       requires: !single
960*c4762a1bSJed Brown 
961*c4762a1bSJed Brown    test:
962*c4762a1bSJed Brown       suffix: 11
963*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
964*c4762a1bSJed Brown       requires: !single
965*c4762a1bSJed Brown 
966*c4762a1bSJed Brown    test:
967*c4762a1bSJed Brown       suffix: 12
968*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
969*c4762a1bSJed Brown       requires: !single
970*c4762a1bSJed Brown 
971*c4762a1bSJed Brown    test:
972*c4762a1bSJed Brown       suffix: 13
973*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
974*c4762a1bSJed Brown       requires: !single
975*c4762a1bSJed Brown 
976*c4762a1bSJed Brown    test:
977*c4762a1bSJed Brown       suffix: 14
978*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
979*c4762a1bSJed Brown       requires: !single
980*c4762a1bSJed Brown 
981*c4762a1bSJed Brown    test:
982*c4762a1bSJed Brown       suffix: 15
983*c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
984*c4762a1bSJed Brown       requires: !single
985*c4762a1bSJed Brown 
986*c4762a1bSJed Brown    test:
987*c4762a1bSJed Brown      suffix: 16
988*c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
989*c4762a1bSJed Brown      requires: !single
990*c4762a1bSJed Brown 
991*c4762a1bSJed Brown    test:
992*c4762a1bSJed Brown      suffix: 17
993*c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
994*c4762a1bSJed Brown      requires: !single
995*c4762a1bSJed Brown 
996*c4762a1bSJed Brown TEST*/
997