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 PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 66c4762a1bSJed Brown PetscInt m; /* number of local elements in vectors */ 67c4762a1bSJed Brown Vec x; /* variables vector */ 68c4762a1bSJed Brown Vec xl,xu; /* bounds vectors */ 69c4762a1bSJed Brown PetscReal d1000 = 1000; 70c4762a1bSJed Brown PetscBool flg,testgetdiag; /* A return variable when checking for user options */ 71c4762a1bSJed Brown Tao tao; /* Tao solver context */ 72c4762a1bSJed Brown KSP ksp; 73c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 74c4762a1bSJed Brown PetscReal zero = 0.0; /* lower bound on all variables */ 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Initialize PETSC and TAO */ 77*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv,(char *)0,help)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Set the default values for the problem parameters */ 80c4762a1bSJed Brown user.nx = 50; user.ny = 50; user.ecc = 0.1; user.b = 10.0; 81c4762a1bSJed Brown testgetdiag = PETSC_FALSE; 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx",&user.nx,&flg)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my",&user.ny,&flg)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ecc",&user.ecc,&flg)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-b",&user.b,&flg)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_getdiagonal",&testgetdiag,NULL)); 89c4762a1bSJed Brown 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n")); 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mx: %D, my: %D, ecc: %g \n\n",user.nx,user.ny,(double)user.ecc)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Let Petsc determine the grid division */ 94c4762a1bSJed Brown Nx = PETSC_DECIDE; Ny = PETSC_DECIDE; 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* 97c4762a1bSJed Brown A two dimensional distributed array will help define this problem, 98c4762a1bSJed Brown which derives from an elliptic PDE on two dimensional domain. From 99c4762a1bSJed Brown the distributed array, Create the vectors. 100c4762a1bSJed Brown */ 1015f80ce2aSJacob 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)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.dm)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.dm)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown Extract global and local vectors from DM; the vector user.B is 107c4762a1bSJed Brown used solely as work space for the evaluation of the function, 108c4762a1bSJed Brown gradient, and Hessian. Duplicate for remaining vectors that are 109c4762a1bSJed Brown the same types. 110c4762a1bSJed Brown */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(user.dm,&x)); /* Solution */ 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&user.B)); /* Linear objective */ 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Create matrix user.A to store quadratic, Create a local ordering scheme. */ 1155f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(x,&m)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(user.dm,&user.A)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown if (testgetdiag) { 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOperation(user.A,MATOP_GET_DIAGONAL,NULL)); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* User defined function -- compute linear term of quadratic */ 1235f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeB(&user)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* The TAO code begins here */ 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* 128c4762a1bSJed Brown Create the optimization solver 129c4762a1bSJed Brown Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM 130c4762a1bSJed Brown */ 1315f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOBLMVM)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Set the initial vector */ 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x, zero)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao,x)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Set the user function, gradient, hessian evaluation routines and data structures */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user)); 140c4762a1bSJed Brown 1415f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,user.A,user.A,FormHessian,(void*)&user)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* Set a routine that defines the bounds */ 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&xl)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&xu)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xl, zero)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xu, d1000)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetVariableBounds(tao,xl,xu)); 149c4762a1bSJed Brown 1505f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetKSP(tao,&ksp)); 151c4762a1bSJed Brown if (ksp) { 1525f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(ksp,KSPCG)); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 1555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-testmonitor",&flg)); 156c4762a1bSJed Brown if (flg) { 1575f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetMonitor(tao,Monitor,&user,NULL)); 158c4762a1bSJed Brown } 1595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-testconvergence",&flg)); 160c4762a1bSJed Brown if (flg) { 1615f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetConvergenceTest(tao,ConvergenceTest,&user)); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Check for any tao command line options */ 1655f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* Solve the bound constrained problem */ 1685f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* Free PETSc data structures */ 1715f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xl)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xu)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.A)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.B)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* Free TAO data structures */ 1785f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.dm)); 180*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 181*b122ec5aSJacob Faibussowitsch return 0; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown static PetscReal p(PetscReal xi, PetscReal ecc) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown PetscReal t=1.0+ecc*PetscCosScalar(xi); 187c4762a1bSJed Brown return (t*t*t); 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown PetscErrorCode ComputeB(AppCtx* user) 191c4762a1bSJed Brown { 192c4762a1bSJed Brown PetscInt i,j,k; 193c4762a1bSJed Brown PetscInt nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym; 194c4762a1bSJed Brown PetscReal two=2.0, pi=4.0*atan(1.0); 195c4762a1bSJed Brown PetscReal hx,hy,ehxhy; 196c4762a1bSJed Brown PetscReal temp,*b; 197c4762a1bSJed Brown PetscReal ecc=user->ecc; 198c4762a1bSJed Brown 199780b99b1SStefano Zampini PetscFunctionBegin; 200c4762a1bSJed Brown nx=user->nx; 201c4762a1bSJed Brown ny=user->ny; 202c4762a1bSJed Brown hx=two*pi/(nx+1.0); 203c4762a1bSJed Brown hy=two*user->b/(ny+1.0); 204c4762a1bSJed Brown ehxhy = ecc*hx*hy; 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown Get local grid boundaries 208c4762a1bSJed Brown */ 2095f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Compute the linear term in the objective function */ 2135f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user->B,&b)); 214c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 215c4762a1bSJed Brown temp=PetscSinScalar((i+1)*hx); 216c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 217c4762a1bSJed Brown k=xm*(j-ys)+(i-xs); 218c4762a1bSJed Brown b[k]= - ehxhy*temp; 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 2215f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user->B,&b)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(5.0*xm*ym+3.0*xm)); 223780b99b1SStefano Zampini PetscFunctionReturn(0); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr) 227c4762a1bSJed Brown { 228c4762a1bSJed Brown AppCtx* user=(AppCtx*)ptr; 229c4762a1bSJed Brown PetscInt i,j,k,kk; 230c4762a1bSJed Brown PetscInt col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym; 231c4762a1bSJed Brown PetscReal one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0); 232c4762a1bSJed Brown PetscReal hx,hy,hxhy,hxhx,hyhy; 233c4762a1bSJed Brown PetscReal xi,v[5]; 234c4762a1bSJed Brown PetscReal ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6; 235c4762a1bSJed Brown PetscReal vmiddle, vup, vdown, vleft, vright; 236c4762a1bSJed Brown PetscReal tt,f1,f2; 237c4762a1bSJed Brown PetscReal *x,*g,zero=0.0; 238c4762a1bSJed Brown Vec localX; 239c4762a1bSJed Brown 240780b99b1SStefano Zampini PetscFunctionBegin; 241c4762a1bSJed Brown nx=user->nx; 242c4762a1bSJed Brown ny=user->ny; 243c4762a1bSJed Brown hx=two*pi/(nx+1.0); 244c4762a1bSJed Brown hy=two*user->b/(ny+1.0); 245c4762a1bSJed Brown hxhy=hx*hy; 246c4762a1bSJed Brown hxhx=one/(hx*hx); 247c4762a1bSJed Brown hyhy=one/(hy*hy); 248c4762a1bSJed Brown 2495f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(user->dm,&localX)); 250c4762a1bSJed Brown 2515f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX)); 253c4762a1bSJed Brown 2545f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(G, zero)); 255c4762a1bSJed Brown /* 256c4762a1bSJed Brown Get local grid boundaries 257c4762a1bSJed Brown */ 2585f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 260c4762a1bSJed Brown 2615f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(localX,&x)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(G,&g)); 263c4762a1bSJed Brown 264c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 265c4762a1bSJed Brown xi=(i+1)*hx; 266c4762a1bSJed Brown trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */ 267c4762a1bSJed Brown trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */ 268c4762a1bSJed Brown trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */ 269c4762a1bSJed Brown trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */ 270c4762a1bSJed Brown trule5=trule1; /* L(i,j-1) */ 271c4762a1bSJed Brown trule6=trule2; /* U(i,j+1) */ 272c4762a1bSJed Brown 273c4762a1bSJed Brown vdown=-(trule5+trule2)*hyhy; 274c4762a1bSJed Brown vleft=-hxhx*(trule2+trule4); 275c4762a1bSJed Brown vright= -hxhx*(trule1+trule3); 276c4762a1bSJed Brown vup=-hyhy*(trule1+trule6); 277c4762a1bSJed Brown vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6); 278c4762a1bSJed Brown 279c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 280c4762a1bSJed Brown 281c4762a1bSJed Brown row=(j-gys)*gxm + (i-gxs); 282c4762a1bSJed Brown v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0; 283c4762a1bSJed Brown 284c4762a1bSJed Brown k=0; 285c4762a1bSJed Brown if (j>gys) { 286c4762a1bSJed Brown v[k]=vdown; col[k]=row - gxm; k++; 287c4762a1bSJed Brown } 288c4762a1bSJed Brown 289c4762a1bSJed Brown if (i>gxs) { 290c4762a1bSJed Brown v[k]= vleft; col[k]=row - 1; k++; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown v[k]= vmiddle; col[k]=row; k++; 294c4762a1bSJed Brown 295c4762a1bSJed Brown if (i+1 < gxs+gxm) { 296c4762a1bSJed Brown v[k]= vright; col[k]=row+1; k++; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown if (j+1 <gys+gym) { 300c4762a1bSJed Brown v[k]= vup; col[k] = row+gxm; k++; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown tt=0; 303c4762a1bSJed Brown for (kk=0;kk<k;kk++) { 304c4762a1bSJed Brown tt+=v[kk]*x[col[kk]]; 305c4762a1bSJed Brown } 306c4762a1bSJed Brown row=(j-ys)*xm + (i-xs); 307c4762a1bSJed Brown g[row]=tt; 308c4762a1bSJed Brown 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown } 312c4762a1bSJed Brown 3135f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(localX,&x)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(G,&g)); 315c4762a1bSJed Brown 3165f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(user->dm,&localX)); 317c4762a1bSJed Brown 3185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(X,G,&f1)); 3195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(user->B,X,&f2)); 3205f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(G, one, user->B)); 321c4762a1bSJed Brown *fcn = f1/2.0 + f2; 322c4762a1bSJed Brown 3235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops((91 + 10.0*ym) * xm)); 324780b99b1SStefano Zampini PetscFunctionReturn(0); 325c4762a1bSJed Brown 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* 329c4762a1bSJed Brown FormHessian computes the quadratic term in the quadratic objective function 330c4762a1bSJed Brown Notice that the objective function in this problem is quadratic (therefore a constant 331c4762a1bSJed Brown hessian). If using a nonquadratic solver, then you might want to reconsider this function 332c4762a1bSJed Brown */ 333c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr) 334c4762a1bSJed Brown { 335c4762a1bSJed Brown AppCtx* user=(AppCtx*)ptr; 336c4762a1bSJed Brown PetscInt i,j,k; 337c4762a1bSJed Brown PetscInt col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym; 338c4762a1bSJed Brown PetscReal one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0); 339c4762a1bSJed Brown PetscReal hx,hy,hxhy,hxhx,hyhy; 340c4762a1bSJed Brown PetscReal xi,v[5]; 341c4762a1bSJed Brown PetscReal ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6; 342c4762a1bSJed Brown PetscReal vmiddle, vup, vdown, vleft, vright; 343c4762a1bSJed Brown PetscBool assembled; 344c4762a1bSJed Brown 345780b99b1SStefano Zampini PetscFunctionBegin; 346c4762a1bSJed Brown nx=user->nx; 347c4762a1bSJed Brown ny=user->ny; 348c4762a1bSJed Brown hx=two*pi/(nx+1.0); 349c4762a1bSJed Brown hy=two*user->b/(ny+1.0); 350c4762a1bSJed Brown hxhy=hx*hy; 351c4762a1bSJed Brown hxhx=one/(hx*hx); 352c4762a1bSJed Brown hyhy=one/(hy*hy); 353c4762a1bSJed Brown 354c4762a1bSJed Brown /* 355c4762a1bSJed Brown Get local grid boundaries 356c4762a1bSJed Brown */ 3575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL)); 3585f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssembled(hes,&assembled)); 3605f80ce2aSJacob Faibussowitsch if (assembled) CHKERRQ(MatZeroEntries(hes)); 361c4762a1bSJed Brown 362c4762a1bSJed Brown for (i=xs; i< xs+xm; i++) { 363c4762a1bSJed Brown xi=(i+1)*hx; 364c4762a1bSJed Brown trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */ 365c4762a1bSJed Brown trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */ 366c4762a1bSJed Brown trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */ 367c4762a1bSJed Brown trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */ 368c4762a1bSJed Brown trule5=trule1; /* L(i,j-1) */ 369c4762a1bSJed Brown trule6=trule2; /* U(i,j+1) */ 370c4762a1bSJed Brown 371c4762a1bSJed Brown vdown=-(trule5+trule2)*hyhy; 372c4762a1bSJed Brown vleft=-hxhx*(trule2+trule4); 373c4762a1bSJed Brown vright= -hxhx*(trule1+trule3); 374c4762a1bSJed Brown vup=-hyhy*(trule1+trule6); 375c4762a1bSJed Brown vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6); 376c4762a1bSJed Brown v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0; 377c4762a1bSJed Brown 378c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 379c4762a1bSJed Brown row=(j-gys)*gxm + (i-gxs); 380c4762a1bSJed Brown 381c4762a1bSJed Brown k=0; 382c4762a1bSJed Brown if (j>gys) { 383c4762a1bSJed Brown v[k]=vdown; col[k]=row - gxm; k++; 384c4762a1bSJed Brown } 385c4762a1bSJed Brown 386c4762a1bSJed Brown if (i>gxs) { 387c4762a1bSJed Brown v[k]= vleft; col[k]=row - 1; k++; 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown v[k]= vmiddle; col[k]=row; k++; 391c4762a1bSJed Brown 392c4762a1bSJed Brown if (i+1 < gxs+gxm) { 393c4762a1bSJed Brown v[k]= vright; col[k]=row+1; k++; 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 396c4762a1bSJed Brown if (j+1 <gys+gym) { 397c4762a1bSJed Brown v[k]= vup; col[k] = row+gxm; k++; 398c4762a1bSJed Brown } 3995f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES)); 400c4762a1bSJed Brown 401c4762a1bSJed Brown } 402c4762a1bSJed Brown 403c4762a1bSJed Brown } 404c4762a1bSJed Brown 405c4762a1bSJed Brown /* 406c4762a1bSJed Brown Assemble matrix, using the 2-step process: 407c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 408c4762a1bSJed Brown By placing code between these two statements, computations can be 409c4762a1bSJed Brown done while messages are in transition. 410c4762a1bSJed Brown */ 4115f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY)); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* 415c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 416c4762a1bSJed Brown matrix. If we do it will generate an error. 417c4762a1bSJed Brown */ 4185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE)); 420c4762a1bSJed Brown 4215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(9.0*xm*ym+49.0*xm)); 422780b99b1SStefano Zampini PetscFunctionReturn(0); 423c4762a1bSJed Brown } 424c4762a1bSJed Brown 425c4762a1bSJed Brown PetscErrorCode Monitor(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; 4325f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason)); 433c4762a1bSJed Brown if (!(its%5)) { 4345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f)); 435c4762a1bSJed Brown } 436c4762a1bSJed Brown PetscFunctionReturn(0); 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 439c4762a1bSJed Brown PetscErrorCode ConvergenceTest(Tao tao, void *ctx) 440c4762a1bSJed Brown { 441c4762a1bSJed Brown PetscInt its; 442c4762a1bSJed Brown PetscReal f,gnorm,cnorm,xdiff; 443c4762a1bSJed Brown TaoConvergedReason reason; 444c4762a1bSJed Brown 445c4762a1bSJed Brown PetscFunctionBegin; 4465f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason)); 447c4762a1bSJed Brown if (its == 100) { 4485f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS)); 449c4762a1bSJed Brown } 450c4762a1bSJed Brown PetscFunctionReturn(0); 451c4762a1bSJed Brown 452c4762a1bSJed Brown } 453c4762a1bSJed Brown 454c4762a1bSJed Brown /*TEST 455c4762a1bSJed Brown 456c4762a1bSJed Brown build: 457c4762a1bSJed Brown requires: !complex 458c4762a1bSJed Brown 459c4762a1bSJed Brown test: 460c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5 461c4762a1bSJed Brown requires: !single 462c4762a1bSJed Brown 463c4762a1bSJed Brown test: 464c4762a1bSJed Brown suffix: 2 465c4762a1bSJed Brown nsize: 2 466c4762a1bSJed Brown args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5 467c4762a1bSJed Brown requires: !single 468c4762a1bSJed Brown 469c4762a1bSJed Brown test: 470c4762a1bSJed Brown suffix: 3 471c4762a1bSJed Brown nsize: 2 472c4762a1bSJed Brown args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 473c4762a1bSJed Brown requires: !single 474c4762a1bSJed Brown 475c4762a1bSJed Brown test: 476c4762a1bSJed Brown suffix: 4 477c4762a1bSJed Brown nsize: 2 478c4762a1bSJed Brown args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal 479c4762a1bSJed Brown output_file: output/jbearing2_3.out 480c4762a1bSJed Brown requires: !single 481c4762a1bSJed Brown 482c4762a1bSJed Brown test: 483c4762a1bSJed Brown suffix: 5 484c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4 485c4762a1bSJed Brown requires: !single 486c4762a1bSJed Brown 487c4762a1bSJed Brown test: 488c4762a1bSJed Brown suffix: 6 489c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4 490c4762a1bSJed Brown requires: !single 491c4762a1bSJed Brown 492c4762a1bSJed Brown test: 493c4762a1bSJed Brown suffix: 7 494c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 495c4762a1bSJed Brown requires: !single 496c4762a1bSJed Brown 497c4762a1bSJed Brown test: 498c4762a1bSJed Brown suffix: 8 499c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 500c4762a1bSJed Brown requires: !single 501c4762a1bSJed Brown 502c4762a1bSJed Brown test: 503c4762a1bSJed Brown suffix: 9 504c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 505c4762a1bSJed Brown requires: !single 506c4762a1bSJed Brown 507c4762a1bSJed Brown test: 508c4762a1bSJed Brown suffix: 10 509c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 510c4762a1bSJed Brown requires: !single 511c4762a1bSJed Brown 512c4762a1bSJed Brown test: 513c4762a1bSJed Brown suffix: 11 514c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 515c4762a1bSJed Brown requires: !single 516c4762a1bSJed Brown 517c4762a1bSJed Brown test: 518c4762a1bSJed Brown suffix: 12 519c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 520c4762a1bSJed Brown requires: !single 521c4762a1bSJed Brown 522c4762a1bSJed Brown test: 523c4762a1bSJed Brown suffix: 13 524c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls 525c4762a1bSJed Brown requires: !single 526c4762a1bSJed Brown 527c4762a1bSJed Brown test: 528c4762a1bSJed Brown suffix: 14 529c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm 530c4762a1bSJed Brown requires: !single 531c4762a1bSJed Brown 532c4762a1bSJed Brown test: 533c4762a1bSJed Brown suffix: 15 534c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs 535c4762a1bSJed Brown requires: !single 536c4762a1bSJed Brown 537c4762a1bSJed Brown test: 538c4762a1bSJed Brown suffix: 16 539c4762a1bSJed Brown args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 540c4762a1bSJed Brown requires: !single 541c4762a1bSJed Brown 542c4762a1bSJed Brown test: 543c4762a1bSJed Brown suffix: 17 544864588a7SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view 545c4762a1bSJed Brown requires: !single 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: 18 549864588a7SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view 550c4762a1bSJed Brown requires: !single 551c4762a1bSJed Brown 55234ad9904SAlp Dener test: 55334ad9904SAlp Dener suffix: 19 55434ad9904SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian 55534ad9904SAlp Dener requires: !single 55634ad9904SAlp Dener 55734ad9904SAlp Dener test: 55834ad9904SAlp Dener suffix: 20 55934ad9904SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian 56034ad9904SAlp Dener requires: !single 56134ad9904SAlp Dener 56234ad9904SAlp Dener test: 56334ad9904SAlp Dener suffix: 21 56434ad9904SAlp Dener args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian 56534ad9904SAlp Dener requires: !single 566c4762a1bSJed Brown TEST*/ 567