xref: /petsc/src/tao/bound/tutorials/jbearing2.c (revision ca0c957dd390e0254b8e0ddfb1485a321c0a7f54)
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 /*T
25c4762a1bSJed Brown    Concepts: TAO^Solving a bound constrained minimization problem
26c4762a1bSJed Brown    Routines: TaoCreate();
27c4762a1bSJed Brown    Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
28c4762a1bSJed Brown    Routines: TaoSetHessianRoutine();
29c4762a1bSJed Brown    Routines: TaoSetVariableBounds();
30c4762a1bSJed Brown    Routines: TaoSetMonitor(); TaoSetConvergenceTest();
31c4762a1bSJed Brown    Routines: TaoSetInitialVector();
32c4762a1bSJed Brown    Routines: TaoSetFromOptions();
33c4762a1bSJed Brown    Routines: TaoSolve();
34c4762a1bSJed Brown    Routines: TaoDestroy();
35c4762a1bSJed Brown    Processors: n
36c4762a1bSJed Brown T*/
37c4762a1bSJed Brown 
38c4762a1bSJed Brown 
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /*
41c4762a1bSJed Brown    User-defined application context - contains data needed by the
42c4762a1bSJed Brown    application-provided call-back routines, FormFunctionGradient(),
43c4762a1bSJed Brown    FormHessian().
44c4762a1bSJed Brown */
45c4762a1bSJed Brown typedef struct {
46c4762a1bSJed Brown   /* problem parameters */
47c4762a1bSJed Brown   PetscReal      ecc;          /* test problem parameter */
48c4762a1bSJed Brown   PetscReal      b;            /* A dimension of journal bearing */
49c4762a1bSJed Brown   PetscInt       nx,ny;        /* discretization in x, y directions */
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Working space */
52c4762a1bSJed Brown   DM          dm;           /* distributed array data structure */
53c4762a1bSJed Brown   Mat         A;            /* Quadratic Objective term */
54c4762a1bSJed Brown   Vec         B;            /* Linear Objective term */
55c4762a1bSJed Brown } AppCtx;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /* User-defined routines */
58c4762a1bSJed Brown static PetscReal p(PetscReal xi, PetscReal ecc);
59c4762a1bSJed Brown static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *,Vec,void *);
60c4762a1bSJed Brown static PetscErrorCode FormHessian(Tao,Vec,Mat, Mat, void *);
61c4762a1bSJed Brown static PetscErrorCode ComputeB(AppCtx*);
62c4762a1bSJed Brown static PetscErrorCode Monitor(Tao, void*);
63c4762a1bSJed Brown static PetscErrorCode ConvergenceTest(Tao, void*);
64c4762a1bSJed Brown 
65c4762a1bSJed Brown int main( int argc, char **argv )
66c4762a1bSJed Brown {
67c4762a1bSJed Brown   PetscErrorCode     ierr;            /* used to check for functions returning nonzeros */
68c4762a1bSJed Brown   PetscInt           Nx, Ny;          /* number of processors in x- and y- directions */
69c4762a1bSJed Brown   PetscInt           m;               /* number of local elements in vectors */
70c4762a1bSJed Brown   Vec                x;               /* variables vector */
71c4762a1bSJed Brown   Vec                xl,xu;           /* bounds vectors */
72c4762a1bSJed Brown   PetscReal          d1000 = 1000;
73c4762a1bSJed Brown   PetscBool          flg,testgetdiag; /* A return variable when checking for user options */
74c4762a1bSJed Brown   Tao                tao;             /* Tao solver context */
75c4762a1bSJed Brown   KSP                ksp;
76c4762a1bSJed Brown   AppCtx             user;            /* user-defined work context */
77c4762a1bSJed Brown   PetscReal          zero = 0.0;      /* lower bound on all variables */
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Initialize PETSC and TAO */
80c4762a1bSJed Brown   ierr = PetscInitialize( &argc, &argv,(char *)0,help );if (ierr) return ierr;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Set the default values for the problem parameters */
83c4762a1bSJed Brown   user.nx = 50; user.ny = 50; user.ecc = 0.1; user.b = 10.0;
84c4762a1bSJed Brown   testgetdiag = PETSC_FALSE;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Check for any command line arguments that override defaults */
87c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mx",&user.nx,&flg);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-my",&user.ny,&flg);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-ecc",&user.ecc,&flg);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-b",&user.b,&flg);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_getdiagonal",&testgetdiag,NULL);CHKERRQ(ierr);
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n");CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"mx: %D,  my: %D,  ecc: %g \n\n",user.nx,user.ny,(double)user.ecc);CHKERRQ(ierr);
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Let Petsc determine the grid division */
97c4762a1bSJed Brown   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /*
100c4762a1bSJed Brown      A two dimensional distributed array will help define this problem,
101c4762a1bSJed Brown      which derives from an elliptic PDE on two dimensional domain.  From
102c4762a1bSJed Brown      the distributed array, Create the vectors.
103c4762a1bSJed Brown   */
104c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
105c4762a1bSJed Brown   ierr = DMSetFromOptions(user.dm);CHKERRQ(ierr);
106c4762a1bSJed Brown   ierr = DMSetUp(user.dm);CHKERRQ(ierr);
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /*
109c4762a1bSJed Brown      Extract global and local vectors from DM; the vector user.B is
110c4762a1bSJed Brown      used solely as work space for the evaluation of the function,
111c4762a1bSJed Brown      gradient, and Hessian.  Duplicate for remaining vectors that are
112c4762a1bSJed Brown      the same types.
113c4762a1bSJed Brown   */
114c4762a1bSJed Brown   ierr = DMCreateGlobalVector(user.dm,&x);CHKERRQ(ierr); /* Solution */
115c4762a1bSJed Brown   ierr = VecDuplicate(x,&user.B);CHKERRQ(ierr); /* Linear objective */
116c4762a1bSJed Brown 
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /*  Create matrix user.A to store quadratic, Create a local ordering scheme. */
119c4762a1bSJed Brown   ierr = VecGetLocalSize(x,&m);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = DMCreateMatrix(user.dm,&user.A);CHKERRQ(ierr);
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   if (testgetdiag) {
123c4762a1bSJed Brown     ierr = MatSetOperation(user.A,MATOP_GET_DIAGONAL,NULL);CHKERRQ(ierr);
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* User defined function -- compute linear term of quadratic */
127c4762a1bSJed Brown   ierr = ComputeB(&user);CHKERRQ(ierr);
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* The TAO code begins here */
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /*
132c4762a1bSJed Brown      Create the optimization solver
133c4762a1bSJed Brown      Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
134c4762a1bSJed Brown   */
135c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
136c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
137c4762a1bSJed Brown 
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* Set the initial vector */
140c4762a1bSJed Brown   ierr = VecSet(x, zero);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* Set the user function, gradient, hessian evaluation routines and data structures */
144c4762a1bSJed Brown   ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*) &user);CHKERRQ(ierr);
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   ierr = TaoSetHessianRoutine(tao,user.A,user.A,FormHessian,(void*)&user);CHKERRQ(ierr);
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Set a routine that defines the bounds */
149c4762a1bSJed Brown   ierr = VecDuplicate(x,&xl);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = VecDuplicate(x,&xu);CHKERRQ(ierr);
151c4762a1bSJed Brown   ierr = VecSet(xl, zero);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr = VecSet(xu, d1000);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,xl,xu);CHKERRQ(ierr);
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
156c4762a1bSJed Brown   if (ksp) {
157c4762a1bSJed Brown     ierr = KSPSetType(ksp,KSPCG);CHKERRQ(ierr);
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-testmonitor",&flg);CHKERRQ(ierr);
161c4762a1bSJed Brown   if (flg) {
162c4762a1bSJed Brown     ierr = TaoSetMonitor(tao,Monitor,&user,NULL);CHKERRQ(ierr);
163c4762a1bSJed Brown   }
164c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-testconvergence",&flg);CHKERRQ(ierr);
165c4762a1bSJed Brown   if (flg) {
166c4762a1bSJed Brown     ierr = TaoSetConvergenceTest(tao,ConvergenceTest,&user);CHKERRQ(ierr);
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* Check for any tao command line options */
170c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* Solve the bound constrained problem */
173c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* Free PETSc data structures */
176c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = VecDestroy(&xl);CHKERRQ(ierr);
178c4762a1bSJed Brown   ierr = VecDestroy(&xu);CHKERRQ(ierr);
179c4762a1bSJed Brown   ierr = MatDestroy(&user.A);CHKERRQ(ierr);
180c4762a1bSJed Brown   ierr = VecDestroy(&user.B);CHKERRQ(ierr);
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* Free TAO data structures */
183c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
185c4762a1bSJed Brown   ierr = PetscFinalize();
186c4762a1bSJed Brown   return ierr;
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown 
190c4762a1bSJed Brown static PetscReal p(PetscReal xi, PetscReal ecc)
191c4762a1bSJed Brown {
192c4762a1bSJed Brown   PetscReal t=1.0+ecc*PetscCosScalar(xi);
193c4762a1bSJed Brown   return (t*t*t);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown PetscErrorCode ComputeB(AppCtx* user)
197c4762a1bSJed Brown {
198c4762a1bSJed Brown   PetscErrorCode ierr;
199c4762a1bSJed Brown   PetscInt       i,j,k;
200c4762a1bSJed Brown   PetscInt       nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
201c4762a1bSJed Brown   PetscReal      two=2.0, pi=4.0*atan(1.0);
202c4762a1bSJed Brown   PetscReal      hx,hy,ehxhy;
203c4762a1bSJed Brown   PetscReal      temp,*b;
204c4762a1bSJed Brown   PetscReal      ecc=user->ecc;
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   nx=user->nx;
207c4762a1bSJed Brown   ny=user->ny;
208c4762a1bSJed Brown   hx=two*pi/(nx+1.0);
209c4762a1bSJed Brown   hy=two*user->b/(ny+1.0);
210c4762a1bSJed Brown   ehxhy = ecc*hx*hy;
211c4762a1bSJed Brown 
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /*
214c4762a1bSJed Brown      Get local grid boundaries
215c4762a1bSJed Brown   */
216c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
217c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* Compute the linear term in the objective function */
220c4762a1bSJed Brown   ierr = VecGetArray(user->B,&b);CHKERRQ(ierr);
221c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++){
222c4762a1bSJed Brown     temp=PetscSinScalar((i+1)*hx);
223c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++){
224c4762a1bSJed Brown       k=xm*(j-ys)+(i-xs);
225c4762a1bSJed Brown       b[k]=  - ehxhy*temp;
226c4762a1bSJed Brown     }
227c4762a1bSJed Brown   }
228c4762a1bSJed Brown   ierr = VecRestoreArray(user->B,&b);CHKERRQ(ierr);
229*ca0c957dSBarry Smith   ierr = PetscLogFlops(5.0*xm*ym+3.0*xm);CHKERRQ(ierr);
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   return 0;
232c4762a1bSJed Brown }
233c4762a1bSJed Brown 
234c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr)
235c4762a1bSJed Brown {
236c4762a1bSJed Brown   AppCtx*        user=(AppCtx*)ptr;
237c4762a1bSJed Brown   PetscErrorCode ierr;
238c4762a1bSJed Brown   PetscInt       i,j,k,kk;
239c4762a1bSJed Brown   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
240c4762a1bSJed Brown   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
241c4762a1bSJed Brown   PetscReal      hx,hy,hxhy,hxhx,hyhy;
242c4762a1bSJed Brown   PetscReal      xi,v[5];
243c4762a1bSJed Brown   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
244c4762a1bSJed Brown   PetscReal      vmiddle, vup, vdown, vleft, vright;
245c4762a1bSJed Brown   PetscReal      tt,f1,f2;
246c4762a1bSJed Brown   PetscReal      *x,*g,zero=0.0;
247c4762a1bSJed Brown   Vec            localX;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   nx=user->nx;
250c4762a1bSJed Brown   ny=user->ny;
251c4762a1bSJed Brown   hx=two*pi/(nx+1.0);
252c4762a1bSJed Brown   hy=two*user->b/(ny+1.0);
253c4762a1bSJed Brown   hxhy=hx*hy;
254c4762a1bSJed Brown   hxhx=one/(hx*hx);
255c4762a1bSJed Brown   hyhy=one/(hy*hy);
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   ierr = DMGetLocalVector(user->dm,&localX);CHKERRQ(ierr);
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
260c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   ierr = VecSet(G, zero);CHKERRQ(ierr);
263c4762a1bSJed Brown   /*
264c4762a1bSJed Brown     Get local grid boundaries
265c4762a1bSJed Brown   */
266c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   ierr = VecGetArray(localX,&x);CHKERRQ(ierr);
270c4762a1bSJed Brown   ierr = VecGetArray(G,&g);CHKERRQ(ierr);
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++){
273c4762a1bSJed Brown     xi=(i+1)*hx;
274c4762a1bSJed Brown     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
275c4762a1bSJed Brown     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
276c4762a1bSJed Brown     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
277c4762a1bSJed Brown     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
278c4762a1bSJed Brown     trule5=trule1; /* L(i,j-1) */
279c4762a1bSJed Brown     trule6=trule2; /* U(i,j+1) */
280c4762a1bSJed Brown 
281c4762a1bSJed Brown     vdown=-(trule5+trule2)*hyhy;
282c4762a1bSJed Brown     vleft=-hxhx*(trule2+trule4);
283c4762a1bSJed Brown     vright= -hxhx*(trule1+trule3);
284c4762a1bSJed Brown     vup=-hyhy*(trule1+trule6);
285c4762a1bSJed Brown     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
286c4762a1bSJed Brown 
287c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++){
288c4762a1bSJed Brown 
289c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
290c4762a1bSJed Brown        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
291c4762a1bSJed Brown 
292c4762a1bSJed Brown        k=0;
293c4762a1bSJed Brown        if (j>gys){
294c4762a1bSJed Brown          v[k]=vdown; col[k]=row - gxm; k++;
295c4762a1bSJed Brown        }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown        if (i>gxs){
298c4762a1bSJed Brown          v[k]= vleft; col[k]=row - 1; k++;
299c4762a1bSJed Brown        }
300c4762a1bSJed Brown 
301c4762a1bSJed Brown        v[k]= vmiddle; col[k]=row; k++;
302c4762a1bSJed Brown 
303c4762a1bSJed Brown        if (i+1 < gxs+gxm){
304c4762a1bSJed Brown          v[k]= vright; col[k]=row+1; k++;
305c4762a1bSJed Brown        }
306c4762a1bSJed Brown 
307c4762a1bSJed Brown        if (j+1 <gys+gym){
308c4762a1bSJed Brown          v[k]= vup; col[k] = row+gxm; k++;
309c4762a1bSJed Brown        }
310c4762a1bSJed Brown        tt=0;
311c4762a1bSJed Brown        for (kk=0;kk<k;kk++){
312c4762a1bSJed Brown          tt+=v[kk]*x[col[kk]];
313c4762a1bSJed Brown        }
314c4762a1bSJed Brown        row=(j-ys)*xm + (i-xs);
315c4762a1bSJed Brown        g[row]=tt;
316c4762a1bSJed Brown 
317c4762a1bSJed Brown      }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   }
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   ierr = VecRestoreArray(localX,&x);CHKERRQ(ierr);
322c4762a1bSJed Brown   ierr = VecRestoreArray(G,&g);CHKERRQ(ierr);
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   ierr = DMRestoreLocalVector(user->dm,&localX);CHKERRQ(ierr);
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   ierr = VecDot(X,G,&f1);CHKERRQ(ierr);
327c4762a1bSJed Brown   ierr = VecDot(user->B,X,&f2);CHKERRQ(ierr);
328c4762a1bSJed Brown   ierr = VecAXPY(G, one, user->B);CHKERRQ(ierr);
329c4762a1bSJed Brown   *fcn = f1/2.0 + f2;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown 
332*ca0c957dSBarry Smith   ierr = PetscLogFlops((91 + 10.0*ym) * xm);CHKERRQ(ierr);
333c4762a1bSJed Brown   return 0;
334c4762a1bSJed Brown 
335c4762a1bSJed Brown }
336c4762a1bSJed Brown 
337c4762a1bSJed Brown 
338c4762a1bSJed Brown /*
339c4762a1bSJed Brown    FormHessian computes the quadratic term in the quadratic objective function
340c4762a1bSJed Brown    Notice that the objective function in this problem is quadratic (therefore a constant
341c4762a1bSJed Brown    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
342c4762a1bSJed Brown */
343c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
344c4762a1bSJed Brown {
345c4762a1bSJed Brown   AppCtx*        user=(AppCtx*)ptr;
346c4762a1bSJed Brown   PetscErrorCode ierr;
347c4762a1bSJed Brown   PetscInt       i,j,k;
348c4762a1bSJed Brown   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
349c4762a1bSJed Brown   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
350c4762a1bSJed Brown   PetscReal      hx,hy,hxhy,hxhx,hyhy;
351c4762a1bSJed Brown   PetscReal      xi,v[5];
352c4762a1bSJed Brown   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
353c4762a1bSJed Brown   PetscReal      vmiddle, vup, vdown, vleft, vright;
354c4762a1bSJed Brown   PetscBool      assembled;
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   nx=user->nx;
357c4762a1bSJed Brown   ny=user->ny;
358c4762a1bSJed Brown   hx=two*pi/(nx+1.0);
359c4762a1bSJed Brown   hy=two*user->b/(ny+1.0);
360c4762a1bSJed Brown   hxhy=hx*hy;
361c4762a1bSJed Brown   hxhx=one/(hx*hx);
362c4762a1bSJed Brown   hyhy=one/(hy*hy);
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /*
365c4762a1bSJed Brown     Get local grid boundaries
366c4762a1bSJed Brown   */
367c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
368c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
369c4762a1bSJed Brown   ierr = MatAssembled(hes,&assembled);CHKERRQ(ierr);
370c4762a1bSJed Brown   if (assembled){ierr = MatZeroEntries(hes);CHKERRQ(ierr);}
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   for (i=xs; i< xs+xm; i++){
373c4762a1bSJed Brown     xi=(i+1)*hx;
374c4762a1bSJed Brown     trule1=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc) ) / six; /* L(i,j) */
375c4762a1bSJed Brown     trule2=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc) ) / six; /* U(i,j) */
376c4762a1bSJed Brown     trule3=hxhy*( p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc) ) / six; /* U(i+1,j) */
377c4762a1bSJed Brown     trule4=hxhy*( p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc) ) / six; /* L(i-1,j) */
378c4762a1bSJed Brown     trule5=trule1; /* L(i,j-1) */
379c4762a1bSJed Brown     trule6=trule2; /* U(i,j+1) */
380c4762a1bSJed Brown 
381c4762a1bSJed Brown     vdown=-(trule5+trule2)*hyhy;
382c4762a1bSJed Brown     vleft=-hxhx*(trule2+trule4);
383c4762a1bSJed Brown     vright= -hxhx*(trule1+trule3);
384c4762a1bSJed Brown     vup=-hyhy*(trule1+trule6);
385c4762a1bSJed Brown     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
386c4762a1bSJed Brown     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
387c4762a1bSJed Brown 
388c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++){
389c4762a1bSJed Brown       row=(j-gys)*gxm + (i-gxs);
390c4762a1bSJed Brown 
391c4762a1bSJed Brown       k=0;
392c4762a1bSJed Brown       if (j>gys){
393c4762a1bSJed Brown         v[k]=vdown; col[k]=row - gxm; k++;
394c4762a1bSJed Brown       }
395c4762a1bSJed Brown 
396c4762a1bSJed Brown       if (i>gxs){
397c4762a1bSJed Brown         v[k]= vleft; col[k]=row - 1; k++;
398c4762a1bSJed Brown       }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown       v[k]= vmiddle; col[k]=row; k++;
401c4762a1bSJed Brown 
402c4762a1bSJed Brown       if (i+1 < gxs+gxm){
403c4762a1bSJed Brown         v[k]= vright; col[k]=row+1; k++;
404c4762a1bSJed Brown       }
405c4762a1bSJed Brown 
406c4762a1bSJed Brown       if (j+1 <gys+gym){
407c4762a1bSJed Brown         v[k]= vup; col[k] = row+gxm; k++;
408c4762a1bSJed Brown       }
409c4762a1bSJed Brown       ierr = MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr);
410c4762a1bSJed Brown 
411c4762a1bSJed Brown     }
412c4762a1bSJed Brown 
413c4762a1bSJed Brown   }
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   /*
416c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
417c4762a1bSJed Brown      MatAssemblyBegin(), MatAssemblyEnd().
418c4762a1bSJed Brown      By placing code between these two statements, computations can be
419c4762a1bSJed Brown      done while messages are in transition.
420c4762a1bSJed Brown   */
421c4762a1bSJed Brown   ierr = MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
422c4762a1bSJed Brown   ierr = MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /*
425c4762a1bSJed Brown     Tell the matrix we will never add a new nonzero location to the
426c4762a1bSJed Brown     matrix. If we do it will generate an error.
427c4762a1bSJed Brown   */
428c4762a1bSJed Brown   ierr = MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
429c4762a1bSJed Brown   ierr = MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
430c4762a1bSJed Brown 
431*ca0c957dSBarry Smith   ierr = PetscLogFlops(9.0*xm*ym+49.0*xm);CHKERRQ(ierr);
432c4762a1bSJed Brown   ierr = MatNorm(hes,NORM_1,&hx);CHKERRQ(ierr);
433c4762a1bSJed Brown   return 0;
434c4762a1bSJed Brown }
435c4762a1bSJed Brown 
436c4762a1bSJed Brown PetscErrorCode Monitor(Tao tao, void *ctx)
437c4762a1bSJed Brown {
438c4762a1bSJed Brown   PetscErrorCode     ierr;
439c4762a1bSJed Brown   PetscInt           its;
440c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
441c4762a1bSJed Brown   TaoConvergedReason reason;
442c4762a1bSJed Brown 
443c4762a1bSJed Brown   PetscFunctionBegin;
444c4762a1bSJed Brown   ierr = TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);CHKERRQ(ierr);
445c4762a1bSJed Brown   if (!(its%5)) {
446c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f);CHKERRQ(ierr);
447c4762a1bSJed Brown   }
448c4762a1bSJed Brown   PetscFunctionReturn(0);
449c4762a1bSJed Brown }
450c4762a1bSJed Brown 
451c4762a1bSJed Brown PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
452c4762a1bSJed Brown {
453c4762a1bSJed Brown   PetscErrorCode     ierr;
454c4762a1bSJed Brown   PetscInt           its;
455c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
456c4762a1bSJed Brown   TaoConvergedReason reason;
457c4762a1bSJed Brown 
458c4762a1bSJed Brown   PetscFunctionBegin;
459c4762a1bSJed Brown   ierr = TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);CHKERRQ(ierr);
460c4762a1bSJed Brown   if (its == 100) {
461c4762a1bSJed Brown     ierr = TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS);CHKERRQ(ierr);
462c4762a1bSJed Brown   }
463c4762a1bSJed Brown   PetscFunctionReturn(0);
464c4762a1bSJed Brown 
465c4762a1bSJed Brown }
466c4762a1bSJed Brown 
467c4762a1bSJed Brown 
468c4762a1bSJed Brown /*TEST
469c4762a1bSJed Brown 
470c4762a1bSJed Brown    build:
471c4762a1bSJed Brown       requires: !complex
472c4762a1bSJed Brown 
473c4762a1bSJed Brown    test:
474c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5
475c4762a1bSJed Brown       requires: !single
476c4762a1bSJed Brown 
477c4762a1bSJed Brown    test:
478c4762a1bSJed Brown       suffix: 2
479c4762a1bSJed Brown       nsize: 2
480c4762a1bSJed Brown       args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5
481c4762a1bSJed Brown       requires: !single
482c4762a1bSJed Brown 
483c4762a1bSJed Brown    test:
484c4762a1bSJed Brown       suffix: 3
485c4762a1bSJed Brown       nsize: 2
486c4762a1bSJed Brown       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4
487c4762a1bSJed Brown       requires: !single
488c4762a1bSJed Brown 
489c4762a1bSJed Brown    test:
490c4762a1bSJed Brown       suffix: 4
491c4762a1bSJed Brown       nsize: 2
492c4762a1bSJed Brown       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal
493c4762a1bSJed Brown       output_file: output/jbearing2_3.out
494c4762a1bSJed Brown       requires: !single
495c4762a1bSJed Brown 
496c4762a1bSJed Brown    test:
497c4762a1bSJed Brown       suffix: 5
498c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
499c4762a1bSJed Brown       requires: !single
500c4762a1bSJed Brown 
501c4762a1bSJed Brown    test:
502c4762a1bSJed Brown       suffix: 6
503c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4
504c4762a1bSJed Brown       requires: !single
505c4762a1bSJed Brown 
506c4762a1bSJed Brown    test:
507c4762a1bSJed Brown       suffix: 7
508c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5
509c4762a1bSJed Brown       requires: !single
510c4762a1bSJed Brown 
511c4762a1bSJed Brown    test:
512c4762a1bSJed Brown       suffix: 8
513c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5
514c4762a1bSJed Brown       requires: !single
515c4762a1bSJed Brown 
516c4762a1bSJed Brown    test:
517c4762a1bSJed Brown       suffix: 9
518c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5
519c4762a1bSJed Brown       requires: !single
520c4762a1bSJed Brown 
521c4762a1bSJed Brown    test:
522c4762a1bSJed Brown       suffix: 10
523c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
524c4762a1bSJed Brown       requires: !single
525c4762a1bSJed Brown 
526c4762a1bSJed Brown    test:
527c4762a1bSJed Brown       suffix: 11
528c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
529c4762a1bSJed Brown       requires: !single
530c4762a1bSJed Brown 
531c4762a1bSJed Brown    test:
532c4762a1bSJed Brown       suffix: 12
533c4762a1bSJed Brown       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
534c4762a1bSJed Brown       requires: !single
535c4762a1bSJed Brown 
536c4762a1bSJed Brown    test:
537c4762a1bSJed Brown      suffix: 13
538c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls
539c4762a1bSJed Brown      requires: !single
540c4762a1bSJed Brown 
541c4762a1bSJed Brown    test:
542c4762a1bSJed Brown      suffix: 14
543c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm
544c4762a1bSJed Brown      requires: !single
545c4762a1bSJed Brown 
546c4762a1bSJed Brown    test:
547c4762a1bSJed Brown      suffix: 15
548c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
549c4762a1bSJed Brown      requires: !single
550c4762a1bSJed Brown 
551c4762a1bSJed Brown    test:
552c4762a1bSJed Brown      suffix: 16
553c4762a1bSJed Brown      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
554c4762a1bSJed Brown      requires: !single
555c4762a1bSJed Brown 
556c4762a1bSJed Brown    test:
557c4762a1bSJed Brown      suffix: 17
558864588a7SAlp Dener      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view
559c4762a1bSJed Brown      requires: !single
560c4762a1bSJed Brown 
561c4762a1bSJed Brown    test:
562c4762a1bSJed Brown      suffix: 18
563864588a7SAlp Dener      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view
564c4762a1bSJed Brown      requires: !single
565c4762a1bSJed Brown 
566c4762a1bSJed Brown TEST*/
567