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