xref: /petsc/src/tao/bound/tutorials/jbearing2.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown   Include "petsctao.h" so we can use TAO solvers
3c4762a1bSJed Brown   Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing
4c4762a1bSJed Brown   Include "petscksp.h" so we can set KSP type
5c4762a1bSJed Brown   the parallel mesh.
6c4762a1bSJed Brown */
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petsctao.h>
9c4762a1bSJed Brown #include <petscdmda.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown static  char help[]=
12c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\
13c4762a1bSJed Brown solve a bound constrained minimization problem.  This example is based on \n\
14c4762a1bSJed Brown the problem DPJB from the MINPACK-2 test suite.  This pressure journal \n\
15c4762a1bSJed Brown bearing problem is an example of elliptic variational problem defined over \n\
16c4762a1bSJed Brown a two dimensional rectangle.  By discretizing the domain into triangular \n\
17c4762a1bSJed Brown elements, the pressure surrounding the journal bearing is defined as the \n\
18c4762a1bSJed Brown minimum of a quadratic function whose variables are bounded below by zero.\n\
19c4762a1bSJed Brown The command line options are:\n\
20c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
21c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
22c4762a1bSJed Brown  \n";
23c4762a1bSJed Brown 
24c4762a1bSJed Brown /*
25c4762a1bSJed Brown    User-defined application context - contains data needed by the
26c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
27c4762a1bSJed Brown    FormHessian().
28c4762a1bSJed Brown */
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   /* problem parameters */
31c4762a1bSJed Brown   PetscReal      ecc;          /* test problem parameter */
32c4762a1bSJed Brown   PetscReal      b;            /* A dimension of journal bearing */
33c4762a1bSJed Brown   PetscInt       nx,ny;        /* discretization in x, y directions */
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /* Working space */
36c4762a1bSJed Brown   DM          dm;           /* distributed array data structure */
37c4762a1bSJed Brown   Mat         A;            /* Quadratic Objective term */
38c4762a1bSJed Brown   Vec         B;            /* Linear Objective term */
39c4762a1bSJed Brown } AppCtx;
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /* User-defined routines */
42c4762a1bSJed Brown static PetscReal p(PetscReal xi, PetscReal ecc);
43c4762a1bSJed Brown static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *,Vec,void *);
44c4762a1bSJed Brown static PetscErrorCode FormHessian(Tao,Vec,Mat, Mat, void *);
45c4762a1bSJed Brown static PetscErrorCode ComputeB(AppCtx*);
46c4762a1bSJed Brown static PetscErrorCode Monitor(Tao, void*);
47c4762a1bSJed Brown static PetscErrorCode ConvergenceTest(Tao, void*);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown int main(int argc, char **argv)
50c4762a1bSJed Brown {
51c4762a1bSJed Brown   PetscInt           Nx, Ny;          /* number of processors in x- and y- directions */
52c4762a1bSJed Brown   PetscInt           m;               /* number of local elements in vectors */
53c4762a1bSJed Brown   Vec                x;               /* variables vector */
54c4762a1bSJed Brown   Vec                xl,xu;           /* bounds vectors */
55c4762a1bSJed Brown   PetscReal          d1000 = 1000;
56c4762a1bSJed Brown   PetscBool          flg,testgetdiag; /* A return variable when checking for user options */
57c4762a1bSJed Brown   Tao                tao;             /* Tao solver context */
58c4762a1bSJed Brown   KSP                ksp;
59c4762a1bSJed Brown   AppCtx             user;            /* user-defined work context */
60c4762a1bSJed Brown   PetscReal          zero = 0.0;      /* lower bound on all variables */
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Initialize PETSC and TAO */
639566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv,(char *)0,help));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Set the default values for the problem parameters */
66c4762a1bSJed Brown   user.nx = 50; user.ny = 50; user.ecc = 0.1; user.b = 10.0;
67c4762a1bSJed Brown   testgetdiag = PETSC_FALSE;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.nx,&flg));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.ny,&flg));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-ecc",&user.ecc,&flg));
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-b",&user.b,&flg));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_getdiagonal",&testgetdiag,NULL));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n"));
77*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mx: %" PetscInt_FMT ",  my: %" PetscInt_FMT ",  ecc: %g \n\n",user.nx,user.ny,(double)user.ecc));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Let Petsc determine the grid division */
80c4762a1bSJed Brown   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /*
83c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
84c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
85c4762a1bSJed Brown      the distributed array, Create the vectors.
86c4762a1bSJed Brown   */
879566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.nx,user.ny,Nx,Ny,1,1,NULL,NULL,&user.dm));
889566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
899566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /*
92c4762a1bSJed Brown      Extract global and local vectors from DM; the vector user.B is
93c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
94c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
95c4762a1bSJed Brown      the same types.
96c4762a1bSJed Brown   */
979566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm,&x)); /* Solution */
989566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&user.B)); /* Linear objective */
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*  Create matrix user.A to store quadratic, Create a local ordering scheme. */
1019566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x,&m));
1029566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.dm,&user.A));
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   if (testgetdiag) {
1059566063dSJacob Faibussowitsch     PetscCall(MatSetOperation(user.A,MATOP_GET_DIAGONAL,NULL));
106c4762a1bSJed Brown   }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* User defined function -- compute linear term of quadratic */
1099566063dSJacob Faibussowitsch   PetscCall(ComputeB(&user));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* The TAO code begins here */
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /*
114c4762a1bSJed Brown      Create the optimization solver
115c4762a1bSJed Brown      Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
116c4762a1bSJed Brown   */
1179566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
1189566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBLMVM));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* Set the initial vector */
1219566063dSJacob Faibussowitsch   PetscCall(VecSet(x, zero));
1229566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,x));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* Set the user function, gradient, hessian evaluation routines and data structures */
1259566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user));
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao,user.A,user.A,FormHessian,(void*)&user));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Set a routine that defines the bounds */
1309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&xl));
1319566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&xu));
1329566063dSJacob Faibussowitsch   PetscCall(VecSet(xl, zero));
1339566063dSJacob Faibussowitsch   PetscCall(VecSet(xu, d1000));
1349566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao,xl,xu));
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao,&ksp));
137c4762a1bSJed Brown   if (ksp) {
1389566063dSJacob Faibussowitsch     PetscCall(KSPSetType(ksp,KSPCG));
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-testmonitor",&flg));
142c4762a1bSJed Brown   if (flg) {
1439566063dSJacob Faibussowitsch     PetscCall(TaoSetMonitor(tao,Monitor,&user,NULL));
144c4762a1bSJed Brown   }
1459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-testconvergence",&flg));
146c4762a1bSJed Brown   if (flg) {
1479566063dSJacob Faibussowitsch     PetscCall(TaoSetConvergenceTest(tao,ConvergenceTest,&user));
148c4762a1bSJed Brown   }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Check for any tao command line options */
1519566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* Solve the bound constrained problem */
1549566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Free PETSc data structures */
1579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xl));
1599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xu));
1609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
1619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.B));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Free TAO data structures */
1649566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1659566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
1669566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
167b122ec5aSJacob Faibussowitsch   return 0;
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown static PetscReal p(PetscReal xi, PetscReal ecc)
171c4762a1bSJed Brown {
172c4762a1bSJed Brown   PetscReal t=1.0+ecc*PetscCosScalar(xi);
173c4762a1bSJed Brown   return (t*t*t);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown PetscErrorCode ComputeB(AppCtx* user)
177c4762a1bSJed Brown {
178c4762a1bSJed Brown   PetscInt       i,j,k;
179c4762a1bSJed Brown   PetscInt       nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
180c4762a1bSJed Brown   PetscReal      two=2.0, pi=4.0*atan(1.0);
181c4762a1bSJed Brown   PetscReal      hx,hy,ehxhy;
182c4762a1bSJed Brown   PetscReal      temp,*b;
183c4762a1bSJed Brown   PetscReal      ecc=user->ecc;
184c4762a1bSJed Brown 
185780b99b1SStefano Zampini   PetscFunctionBegin;
186c4762a1bSJed Brown   nx=user->nx;
187c4762a1bSJed Brown   ny=user->ny;
188c4762a1bSJed Brown   hx=two*pi/(nx+1.0);
189c4762a1bSJed Brown   hy=two*user->b/(ny+1.0);
190c4762a1bSJed Brown   ehxhy = ecc*hx*hy;
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /*
193c4762a1bSJed Brown      Get local grid boundaries
194c4762a1bSJed Brown   */
1959566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
1969566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Compute the linear term in the objective function */
1999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->B,&b));
200c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
201c4762a1bSJed Brown     temp=PetscSinScalar((i+1)*hx);
202c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
203c4762a1bSJed Brown       k=xm*(j-ys)+(i-xs);
204c4762a1bSJed Brown       b[k]=  - ehxhy*temp;
205c4762a1bSJed Brown     }
206c4762a1bSJed Brown   }
2079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->B,&b));
2089566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(5.0*xm*ym+3.0*xm));
209780b99b1SStefano Zampini   PetscFunctionReturn(0);
210c4762a1bSJed Brown }
211c4762a1bSJed Brown 
212c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr)
213c4762a1bSJed Brown {
214c4762a1bSJed Brown   AppCtx*        user=(AppCtx*)ptr;
215c4762a1bSJed Brown   PetscInt       i,j,k,kk;
216c4762a1bSJed Brown   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
217c4762a1bSJed Brown   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
218c4762a1bSJed Brown   PetscReal      hx,hy,hxhy,hxhx,hyhy;
219c4762a1bSJed Brown   PetscReal      xi,v[5];
220c4762a1bSJed Brown   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
221c4762a1bSJed Brown   PetscReal      vmiddle, vup, vdown, vleft, vright;
222c4762a1bSJed Brown   PetscReal      tt,f1,f2;
223c4762a1bSJed Brown   PetscReal      *x,*g,zero=0.0;
224c4762a1bSJed Brown   Vec            localX;
225c4762a1bSJed Brown 
226780b99b1SStefano Zampini   PetscFunctionBegin;
227c4762a1bSJed Brown   nx=user->nx;
228c4762a1bSJed Brown   ny=user->ny;
229c4762a1bSJed Brown   hx=two*pi/(nx+1.0);
230c4762a1bSJed Brown   hy=two*user->b/(ny+1.0);
231c4762a1bSJed Brown   hxhy=hx*hy;
232c4762a1bSJed Brown   hxhx=one/(hx*hx);
233c4762a1bSJed Brown   hyhy=one/(hy*hy);
234c4762a1bSJed Brown 
2359566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->dm,&localX));
236c4762a1bSJed Brown 
2379566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
2389566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
239c4762a1bSJed Brown 
2409566063dSJacob Faibussowitsch   PetscCall(VecSet(G, zero));
241c4762a1bSJed Brown   /*
242c4762a1bSJed Brown     Get local grid boundaries
243c4762a1bSJed Brown   */
2449566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
2459566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
246c4762a1bSJed Brown 
2479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX,&x));
2489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G,&g));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++) {
251c4762a1bSJed Brown     xi=(i+1)*hx;
252c4762a1bSJed Brown     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
253c4762a1bSJed Brown     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
254c4762a1bSJed Brown     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
255c4762a1bSJed Brown     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
256c4762a1bSJed Brown     trule5=trule1; /* L(i,j-1) */
257c4762a1bSJed Brown     trule6=trule2; /* U(i,j+1) */
258c4762a1bSJed Brown 
259c4762a1bSJed Brown     vdown=-(trule5+trule2)*hyhy;
260c4762a1bSJed Brown     vleft=-hxhx*(trule2+trule4);
261c4762a1bSJed Brown     vright= -hxhx*(trule1+trule3);
262c4762a1bSJed Brown     vup=-hyhy*(trule1+trule6);
263c4762a1bSJed Brown     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
264c4762a1bSJed Brown 
265c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
266c4762a1bSJed Brown 
267c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
268c4762a1bSJed Brown        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
269c4762a1bSJed Brown 
270c4762a1bSJed Brown        k=0;
271c4762a1bSJed Brown        if (j>gys) {
272c4762a1bSJed Brown          v[k]=vdown; col[k]=row - gxm; k++;
273c4762a1bSJed Brown        }
274c4762a1bSJed Brown 
275c4762a1bSJed Brown        if (i>gxs) {
276c4762a1bSJed Brown          v[k]= vleft; col[k]=row - 1; k++;
277c4762a1bSJed Brown        }
278c4762a1bSJed Brown 
279c4762a1bSJed Brown        v[k]= vmiddle; col[k]=row; k++;
280c4762a1bSJed Brown 
281c4762a1bSJed Brown        if (i+1 < gxs+gxm) {
282c4762a1bSJed Brown          v[k]= vright; col[k]=row+1; k++;
283c4762a1bSJed Brown        }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown        if (j+1 <gys+gym) {
286c4762a1bSJed Brown          v[k]= vup; col[k] = row+gxm; k++;
287c4762a1bSJed Brown        }
288c4762a1bSJed Brown        tt=0;
289c4762a1bSJed Brown        for (kk=0;kk<k;kk++) {
290c4762a1bSJed Brown          tt+=v[kk]*x[col[kk]];
291c4762a1bSJed Brown        }
292c4762a1bSJed Brown        row=(j-ys)*xm + (i-xs);
293c4762a1bSJed Brown        g[row]=tt;
294c4762a1bSJed Brown 
295c4762a1bSJed Brown      }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   }
298c4762a1bSJed Brown 
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX,&x));
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G,&g));
301c4762a1bSJed Brown 
3029566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->dm,&localX));
303c4762a1bSJed Brown 
3049566063dSJacob Faibussowitsch   PetscCall(VecDot(X,G,&f1));
3059566063dSJacob Faibussowitsch   PetscCall(VecDot(user->B,X,&f2));
3069566063dSJacob Faibussowitsch   PetscCall(VecAXPY(G, one, user->B));
307c4762a1bSJed Brown   *fcn = f1/2.0 + f2;
308c4762a1bSJed Brown 
3099566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((91 + 10.0*ym) * xm));
310780b99b1SStefano Zampini   PetscFunctionReturn(0);
311c4762a1bSJed Brown 
312c4762a1bSJed Brown }
313c4762a1bSJed Brown 
314c4762a1bSJed Brown /*
315c4762a1bSJed Brown    FormHessian computes the quadratic term in the quadratic objective function
316c4762a1bSJed Brown    Notice that the objective function in this problem is quadratic (therefore a constant
317c4762a1bSJed Brown    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
318c4762a1bSJed Brown */
319c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
320c4762a1bSJed Brown {
321c4762a1bSJed Brown   AppCtx*        user=(AppCtx*)ptr;
322c4762a1bSJed Brown   PetscInt       i,j,k;
323c4762a1bSJed Brown   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
324c4762a1bSJed Brown   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
325c4762a1bSJed Brown   PetscReal      hx,hy,hxhy,hxhx,hyhy;
326c4762a1bSJed Brown   PetscReal      xi,v[5];
327c4762a1bSJed Brown   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
328c4762a1bSJed Brown   PetscReal      vmiddle, vup, vdown, vleft, vright;
329c4762a1bSJed Brown   PetscBool      assembled;
330c4762a1bSJed Brown 
331780b99b1SStefano Zampini   PetscFunctionBegin;
332c4762a1bSJed Brown   nx=user->nx;
333c4762a1bSJed Brown   ny=user->ny;
334c4762a1bSJed Brown   hx=two*pi/(nx+1.0);
335c4762a1bSJed Brown   hy=two*user->b/(ny+1.0);
336c4762a1bSJed Brown   hxhy=hx*hy;
337c4762a1bSJed Brown   hxhx=one/(hx*hx);
338c4762a1bSJed Brown   hyhy=one/(hy*hy);
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   /*
341c4762a1bSJed Brown     Get local grid boundaries
342c4762a1bSJed Brown   */
3439566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
3449566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
3459566063dSJacob Faibussowitsch   PetscCall(MatAssembled(hes,&assembled));
3469566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(hes));
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++) {
349c4762a1bSJed Brown     xi=(i+1)*hx;
350c4762a1bSJed Brown     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
351c4762a1bSJed Brown     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
352c4762a1bSJed Brown     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
353c4762a1bSJed Brown     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
354c4762a1bSJed Brown     trule5=trule1; /* L(i,j-1) */
355c4762a1bSJed Brown     trule6=trule2; /* U(i,j+1) */
356c4762a1bSJed Brown 
357c4762a1bSJed Brown     vdown=-(trule5+trule2)*hyhy;
358c4762a1bSJed Brown     vleft=-hxhx*(trule2+trule4);
359c4762a1bSJed Brown     vright= -hxhx*(trule1+trule3);
360c4762a1bSJed Brown     vup=-hyhy*(trule1+trule6);
361c4762a1bSJed Brown     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
362c4762a1bSJed Brown     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
363c4762a1bSJed Brown 
364c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
365c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
366c4762a1bSJed Brown 
367c4762a1bSJed Brown       k=0;
368c4762a1bSJed Brown       if (j>gys) {
369c4762a1bSJed Brown         v[k]=vdown; col[k]=row - gxm; k++;
370c4762a1bSJed Brown       }
371c4762a1bSJed Brown 
372c4762a1bSJed Brown       if (i>gxs) {
373c4762a1bSJed Brown         v[k]= vleft; col[k]=row - 1; k++;
374c4762a1bSJed Brown       }
375c4762a1bSJed Brown 
376c4762a1bSJed Brown       v[k]= vmiddle; col[k]=row; k++;
377c4762a1bSJed Brown 
378c4762a1bSJed Brown       if (i+1 < gxs+gxm) {
379c4762a1bSJed Brown         v[k]= vright; col[k]=row+1; k++;
380c4762a1bSJed Brown       }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown       if (j+1 <gys+gym) {
383c4762a1bSJed Brown         v[k]= vup; col[k] = row+gxm; k++;
384c4762a1bSJed Brown       }
3859566063dSJacob Faibussowitsch       PetscCall(MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES));
386c4762a1bSJed Brown 
387c4762a1bSJed Brown     }
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   }
390c4762a1bSJed Brown 
391c4762a1bSJed Brown   /*
392c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
393c4762a1bSJed Brown      MatAssemblyBegin(), MatAssemblyEnd().
394c4762a1bSJed Brown      By placing code between these two statements, computations can be
395c4762a1bSJed Brown      done while messages are in transition.
396c4762a1bSJed Brown   */
3979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY));
3989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   /*
401c4762a1bSJed Brown     Tell the matrix we will never add a new nonzero location to the
402c4762a1bSJed Brown     matrix. If we do it will generate an error.
403c4762a1bSJed Brown   */
4049566063dSJacob Faibussowitsch   PetscCall(MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
4059566063dSJacob Faibussowitsch   PetscCall(MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE));
406c4762a1bSJed Brown 
4079566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(9.0*xm*ym+49.0*xm));
408780b99b1SStefano Zampini   PetscFunctionReturn(0);
409c4762a1bSJed Brown }
410c4762a1bSJed Brown 
411c4762a1bSJed Brown PetscErrorCode Monitor(Tao tao, void *ctx)
412c4762a1bSJed Brown {
413c4762a1bSJed Brown   PetscInt           its;
414c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
415c4762a1bSJed Brown   TaoConvergedReason reason;
416c4762a1bSJed Brown 
417c4762a1bSJed Brown   PetscFunctionBegin;
4189566063dSJacob Faibussowitsch   PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
419c4762a1bSJed Brown   if (!(its%5)) {
420*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iteration=%" PetscInt_FMT "\tf=%g\n",its,(double)f));
421c4762a1bSJed Brown   }
422c4762a1bSJed Brown   PetscFunctionReturn(0);
423c4762a1bSJed Brown }
424c4762a1bSJed Brown 
425c4762a1bSJed Brown PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
426c4762a1bSJed Brown {
427c4762a1bSJed Brown   PetscInt           its;
428c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
429c4762a1bSJed Brown   TaoConvergedReason reason;
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
433c4762a1bSJed Brown   if (its == 100) {
4349566063dSJacob Faibussowitsch     PetscCall(TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS));
435c4762a1bSJed Brown   }
436c4762a1bSJed Brown   PetscFunctionReturn(0);
437c4762a1bSJed Brown 
438c4762a1bSJed Brown }
439c4762a1bSJed Brown 
440c4762a1bSJed Brown /*TEST
441c4762a1bSJed Brown 
442c4762a1bSJed Brown    build:
443c4762a1bSJed Brown       requires: !complex
444c4762a1bSJed Brown 
445c4762a1bSJed Brown    test:
446c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5
447c4762a1bSJed Brown       requires: !single
448c4762a1bSJed Brown 
449c4762a1bSJed Brown    test:
450c4762a1bSJed Brown       suffix: 2
451c4762a1bSJed Brown       nsize: 2
452c4762a1bSJed Brown       args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5
453c4762a1bSJed Brown       requires: !single
454c4762a1bSJed Brown 
455c4762a1bSJed Brown    test:
456c4762a1bSJed Brown       suffix: 3
457c4762a1bSJed Brown       nsize: 2
458c4762a1bSJed Brown       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4
459c4762a1bSJed Brown       requires: !single
460c4762a1bSJed Brown 
461c4762a1bSJed Brown    test:
462c4762a1bSJed Brown       suffix: 4
463c4762a1bSJed Brown       nsize: 2
464c4762a1bSJed Brown       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal
465c4762a1bSJed Brown       output_file: output/jbearing2_3.out
466c4762a1bSJed Brown       requires: !single
467c4762a1bSJed Brown 
468c4762a1bSJed Brown    test:
469c4762a1bSJed Brown       suffix: 5
470c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
471c4762a1bSJed Brown       requires: !single
472c4762a1bSJed Brown 
473c4762a1bSJed Brown    test:
474c4762a1bSJed Brown       suffix: 6
475c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4
476c4762a1bSJed Brown       requires: !single
477c4762a1bSJed Brown 
478c4762a1bSJed Brown    test:
479c4762a1bSJed Brown       suffix: 7
480c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5
481c4762a1bSJed Brown       requires: !single
482c4762a1bSJed Brown 
483c4762a1bSJed Brown    test:
484c4762a1bSJed Brown       suffix: 8
485c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5
486c4762a1bSJed Brown       requires: !single
487c4762a1bSJed Brown 
488c4762a1bSJed Brown    test:
489c4762a1bSJed Brown       suffix: 9
490c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5
491c4762a1bSJed Brown       requires: !single
492c4762a1bSJed Brown 
493c4762a1bSJed Brown    test:
494c4762a1bSJed Brown       suffix: 10
495c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
496c4762a1bSJed Brown       requires: !single
497c4762a1bSJed Brown 
498c4762a1bSJed Brown    test:
499c4762a1bSJed Brown       suffix: 11
500c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
501c4762a1bSJed Brown       requires: !single
502c4762a1bSJed Brown 
503c4762a1bSJed Brown    test:
504c4762a1bSJed Brown       suffix: 12
505c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
506c4762a1bSJed Brown       requires: !single
507c4762a1bSJed Brown 
508c4762a1bSJed Brown    test:
509c4762a1bSJed Brown      suffix: 13
510c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls
511c4762a1bSJed Brown      requires: !single
512c4762a1bSJed Brown 
513c4762a1bSJed Brown    test:
514c4762a1bSJed Brown      suffix: 14
515c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm
516c4762a1bSJed Brown      requires: !single
517c4762a1bSJed Brown 
518c4762a1bSJed Brown    test:
519c4762a1bSJed Brown      suffix: 15
520c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
521c4762a1bSJed Brown      requires: !single
522c4762a1bSJed Brown 
523c4762a1bSJed Brown    test:
524c4762a1bSJed Brown      suffix: 16
525c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
526c4762a1bSJed Brown      requires: !single
527c4762a1bSJed Brown 
528c4762a1bSJed Brown    test:
529c4762a1bSJed Brown      suffix: 17
530864588a7SAlp Dener      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view
531c4762a1bSJed Brown      requires: !single
532c4762a1bSJed Brown 
533c4762a1bSJed Brown    test:
534c4762a1bSJed Brown      suffix: 18
535864588a7SAlp Dener      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view
536c4762a1bSJed Brown      requires: !single
537c4762a1bSJed Brown 
53834ad9904SAlp Dener    test:
53934ad9904SAlp Dener      suffix: 19
54034ad9904SAlp Dener      args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
54134ad9904SAlp Dener      requires: !single
54234ad9904SAlp Dener 
54334ad9904SAlp Dener    test:
54434ad9904SAlp Dener       suffix: 20
54534ad9904SAlp Dener       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
54634ad9904SAlp Dener       requires: !single
54734ad9904SAlp Dener 
54834ad9904SAlp Dener    test:
54934ad9904SAlp Dener       suffix: 21
55034ad9904SAlp Dener       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
55134ad9904SAlp Dener       requires: !single
552c4762a1bSJed Brown TEST*/
553