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