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