1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 minsurf1 [-help] [all TAO options] */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* Include "petsctao.h" so we can use TAO solvers. */ 4c4762a1bSJed Brown #include <petsctao.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown static char help[] = 7c4762a1bSJed Brown "This example demonstrates use of the TAO package to \n\ 8c4762a1bSJed Brown solve an unconstrained minimization problem. This example is based on a \n\ 9c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\ 10c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\ 11c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\ 12c4762a1bSJed Brown The command line options are:\n\ 13c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 14c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 15c4762a1bSJed Brown -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n"; 16c4762a1bSJed Brown 17c4762a1bSJed Brown /*T 18c4762a1bSJed Brown Concepts: TAO^Solving an unconstrained minimization problem 19c4762a1bSJed Brown Routines: TaoCreate(); TaoSetType(); 20c4762a1bSJed Brown Routines: TaoSetInitialVector(); 21c4762a1bSJed Brown Routines: TaoSetObjectiveAndGradientRoutine(); 22c4762a1bSJed Brown Routines: TaoSetHessianRoutine(); TaoSetFromOptions(); 23c4762a1bSJed Brown Routines: TaoGetKSP(); TaoSolve(); 24c4762a1bSJed Brown Routines: TaoDestroy(); 25c4762a1bSJed Brown Processors: 1 26c4762a1bSJed Brown T*/ 27c4762a1bSJed Brown 28c4762a1bSJed Brown 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* 31c4762a1bSJed Brown User-defined application context - contains data needed by the 32c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient() 33c4762a1bSJed Brown and FormHessian(). 34c4762a1bSJed Brown */ 35c4762a1bSJed Brown typedef struct { 36c4762a1bSJed Brown PetscInt mx, my; /* discretization in x, y directions */ 37c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right; /* boundary values */ 38c4762a1bSJed Brown Mat H; 39c4762a1bSJed Brown } AppCtx; 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 42c4762a1bSJed Brown 43c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx*); 44c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec); 45c4762a1bSJed Brown static PetscErrorCode QuadraticH(AppCtx*,Vec,Mat); 46c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 47c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 48c4762a1bSJed Brown 49c4762a1bSJed Brown int main( int argc, char **argv ) 50c4762a1bSJed Brown { 51c4762a1bSJed Brown PetscErrorCode ierr; /* used to check for functions returning nonzeros */ 52c4762a1bSJed Brown PetscInt N; /* Size of vector */ 53c4762a1bSJed Brown PetscMPIInt size; /* Number of processors */ 54c4762a1bSJed Brown Vec x; /* solution, gradient vectors */ 55c4762a1bSJed Brown KSP ksp; /* PETSc Krylov subspace method */ 56c4762a1bSJed Brown PetscBool flg; /* A return value when checking for user options */ 57c4762a1bSJed Brown Tao tao; /* Tao solver context */ 58c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Initialize TAO,PETSc */ 61c4762a1bSJed Brown ierr = PetscInitialize( &argc, &argv,(char *)0,help );if (ierr) return ierr; 62c4762a1bSJed Brown 63c4762a1bSJed Brown ierr = MPI_Comm_size(MPI_COMM_WORLD,&size);CHKERRQ(ierr); 64c4762a1bSJed Brown if (size >1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Specify default dimension of the problem */ 67c4762a1bSJed Brown user.mx = 4; user.my = 4; 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 70c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg);CHKERRQ(ierr); 72c4762a1bSJed Brown 73c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n");CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"mx: %D my: %D \n\n",user.mx,user.my);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Calculate any derived values from parameters */ 77c4762a1bSJed Brown N = user.mx*user.my; 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 80c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_SELF,&tao);CHKERRQ(ierr); 81c4762a1bSJed Brown ierr = TaoSetType(tao,TAOLMVM);CHKERRQ(ierr); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Initialize minsurf application data structure for use in the function evaluations */ 84c4762a1bSJed Brown ierr = MSA_BoundaryConditions(&user);CHKERRQ(ierr); /* Application specific routine */ 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* 87c4762a1bSJed Brown Create a vector to hold the variables. Compute an initial solution. 88c4762a1bSJed Brown Set this vector, which will be used by TAO. 89c4762a1bSJed Brown */ 90c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MSA_InitialPoint(&user,x);CHKERRQ(ierr); /* Application specific routine */ 92c4762a1bSJed Brown ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr); /* A TAO routine */ 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Provide TAO routines for function, gradient, and Hessian evaluation */ 95c4762a1bSJed Brown ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);CHKERRQ(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Create a matrix data structure to store the Hessian. This structure will be used by TAO */ 98c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,7,NULL,&(user.H));CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void *)&user);CHKERRQ(ierr); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Check for any TAO command line options */ 103c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Limit the number of iterations in the KSP linear solver */ 106c4762a1bSJed Brown ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 107c4762a1bSJed Brown if (ksp) { 108c4762a1bSJed Brown ierr = KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,user.mx*user.my);CHKERRQ(ierr); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 112c4762a1bSJed Brown ierr = TaoSolve(tao);CHKERRQ(ierr); 113c4762a1bSJed Brown 114c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = MatDestroy(&user.H);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = PetscFree(user.bottom);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = PetscFree(user.top);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscFree(user.left);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = PetscFree(user.right);CHKERRQ(ierr); 121c4762a1bSJed Brown 122c4762a1bSJed Brown ierr = PetscFinalize(); 123c4762a1bSJed Brown return ierr; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* FormFunctionGradient - Evaluates function and corresponding gradient. 129c4762a1bSJed Brown 130c4762a1bSJed Brown Input Parameters: 131c4762a1bSJed Brown . tao - the Tao context 132c4762a1bSJed Brown . X - input vector 133c4762a1bSJed Brown . userCtx - optional user-defined context, as set by TaoSetFunctionGradient() 134c4762a1bSJed Brown 135c4762a1bSJed Brown Output Parameters: 136c4762a1bSJed Brown . fcn - the newly evaluated function 137c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 138c4762a1bSJed Brown */ 139c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *fcn,Vec G,void *userCtx) 140c4762a1bSJed Brown { 141c4762a1bSJed Brown AppCtx *user = (AppCtx *) userCtx; 142c4762a1bSJed Brown PetscErrorCode ierr; 143c4762a1bSJed Brown PetscInt i,j,row; 144c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 145c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 146c4762a1bSJed Brown PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy, ft=0; 147c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 148c4762a1bSJed Brown PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 149c4762a1bSJed Brown PetscReal zero=0.0; 150c4762a1bSJed Brown PetscScalar *g; 151c4762a1bSJed Brown const PetscScalar *x; 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscFunctionBeginUser; 154c4762a1bSJed Brown ierr = VecSet(G, zero);CHKERRQ(ierr); 155c4762a1bSJed Brown 156c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = VecGetArray(G,&g);CHKERRQ(ierr); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */ 160c4762a1bSJed Brown for (j=0; j<my; j++){ 161c4762a1bSJed Brown for (i=0; i< mx; i++){ 162c4762a1bSJed Brown row=(j)*mx + (i); 163c4762a1bSJed Brown xc = x[row]; 164c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 165c4762a1bSJed Brown if (i==0){ /* left side */ 166c4762a1bSJed Brown xl = user->left[j+1]; 167c4762a1bSJed Brown xlt = user->left[j+2]; 168c4762a1bSJed Brown } else { 169c4762a1bSJed Brown xl = x[row-1]; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown if (j==0){ /* bottom side */ 173c4762a1bSJed Brown xb = user->bottom[i+1]; 174c4762a1bSJed Brown xrb = user->bottom[i+2]; 175c4762a1bSJed Brown } else { 176c4762a1bSJed Brown xb = x[row-mx]; 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown if (i+1 == mx){ /* right side */ 180c4762a1bSJed Brown xr = user->right[j+1]; 181c4762a1bSJed Brown xrb = user->right[j]; 182c4762a1bSJed Brown } else { 183c4762a1bSJed Brown xr = x[row+1]; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown if (j+1==0+my){ /* top side */ 187c4762a1bSJed Brown xt = user->top[i+1]; 188c4762a1bSJed Brown xlt = user->top[i]; 189c4762a1bSJed Brown }else { 190c4762a1bSJed Brown xt = x[row+mx]; 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown if (i>0 && j+1<my){ 194c4762a1bSJed Brown xlt = x[row-1+mx]; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown if (j>0 && i+1<mx){ 197c4762a1bSJed Brown xrb = x[row+1-mx]; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown 200c4762a1bSJed Brown d1 = (xc-xl); 201c4762a1bSJed Brown d2 = (xc-xr); 202c4762a1bSJed Brown d3 = (xc-xt); 203c4762a1bSJed Brown d4 = (xc-xb); 204c4762a1bSJed Brown d5 = (xr-xrb); 205c4762a1bSJed Brown d6 = (xrb-xb); 206c4762a1bSJed Brown d7 = (xlt-xl); 207c4762a1bSJed Brown d8 = (xt-xlt); 208c4762a1bSJed Brown 209c4762a1bSJed Brown df1dxc = d1*hydhx; 210c4762a1bSJed Brown df2dxc = ( d1*hydhx + d4*hxdhy ); 211c4762a1bSJed Brown df3dxc = d3*hxdhy; 212c4762a1bSJed Brown df4dxc = ( d2*hydhx + d3*hxdhy ); 213c4762a1bSJed Brown df5dxc = d2*hydhx; 214c4762a1bSJed Brown df6dxc = d4*hxdhy; 215c4762a1bSJed Brown 216c4762a1bSJed Brown d1 *= rhx; 217c4762a1bSJed Brown d2 *= rhx; 218c4762a1bSJed Brown d3 *= rhy; 219c4762a1bSJed Brown d4 *= rhy; 220c4762a1bSJed Brown d5 *= rhy; 221c4762a1bSJed Brown d6 *= rhx; 222c4762a1bSJed Brown d7 *= rhy; 223c4762a1bSJed Brown d8 *= rhx; 224c4762a1bSJed Brown 225c4762a1bSJed Brown f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7); 226c4762a1bSJed Brown f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4); 227c4762a1bSJed Brown f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8); 228c4762a1bSJed Brown f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2); 229c4762a1bSJed Brown f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5); 230c4762a1bSJed Brown f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6); 231c4762a1bSJed Brown 232c4762a1bSJed Brown ft = ft + (f2 + f4); 233c4762a1bSJed Brown 234c4762a1bSJed Brown df1dxc /= f1; 235c4762a1bSJed Brown df2dxc /= f2; 236c4762a1bSJed Brown df3dxc /= f3; 237c4762a1bSJed Brown df4dxc /= f4; 238c4762a1bSJed Brown df5dxc /= f5; 239c4762a1bSJed Brown df6dxc /= f6; 240c4762a1bSJed Brown 241c4762a1bSJed Brown g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc )/2.0; 242c4762a1bSJed Brown } 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown for (j=0; j<my; j++){ /* left side */ 246c4762a1bSJed Brown d3 = (user->left[j+1] - user->left[j+2])*rhy; 247c4762a1bSJed Brown d2 = (user->left[j+1] - x[j*mx])*rhx; 248c4762a1bSJed Brown ft = ft+PetscSqrtScalar( 1.0 + d3*d3 + d2*d2); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown for (i=0; i<mx; i++){ /* bottom */ 252c4762a1bSJed Brown d2 = (user->bottom[i+1]-user->bottom[i+2])*rhx; 253c4762a1bSJed Brown d3 = (user->bottom[i+1]-x[i])*rhy; 254c4762a1bSJed Brown ft = ft+PetscSqrtScalar( 1.0 + d3*d3 + d2*d2); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown for (j=0; j< my; j++){ /* right side */ 258c4762a1bSJed Brown d1 = (x[(j+1)*mx-1]-user->right[j+1])*rhx; 259c4762a1bSJed Brown d4 = (user->right[j]-user->right[j+1])*rhy; 260c4762a1bSJed Brown ft = ft+PetscSqrtScalar( 1.0 + d1*d1 + d4*d4); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown for (i=0; i<mx; i++){ /* top side */ 264c4762a1bSJed Brown d1 = (x[(my-1)*mx + i] - user->top[i+1])*rhy; 265c4762a1bSJed Brown d4 = (user->top[i+1] - user->top[i])*rhx; 266c4762a1bSJed Brown ft = ft+PetscSqrtScalar( 1.0 + d1*d1 + d4*d4); 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* Bottom left corner */ 270c4762a1bSJed Brown d1 = (user->left[0]-user->left[1])*rhy; 271c4762a1bSJed Brown d2 = (user->bottom[0]-user->bottom[1])*rhx; 272c4762a1bSJed Brown ft += PetscSqrtScalar( 1.0 + d1*d1 + d2*d2); 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* Top right corner */ 275c4762a1bSJed Brown d1 = (user->right[my+1] - user->right[my])*rhy; 276c4762a1bSJed Brown d2 = (user->top[mx+1] - user->top[mx])*rhx; 277c4762a1bSJed Brown ft += PetscSqrtScalar( 1.0 + d1*d1 + d2*d2); 278c4762a1bSJed Brown 279c4762a1bSJed Brown (*fcn)=ft*area; 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* Restore vectors */ 282c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = VecRestoreArray(G,&g);CHKERRQ(ierr); 284*ca0c957dSBarry Smith ierr = PetscLogFlops(67.0*mx*my);CHKERRQ(ierr); 285c4762a1bSJed Brown PetscFunctionReturn(0); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 289c4762a1bSJed Brown /* 290c4762a1bSJed Brown FormHessian - Evaluates the Hessian matrix. 291c4762a1bSJed Brown 292c4762a1bSJed Brown Input Parameters: 293c4762a1bSJed Brown . tao - the Tao context 294c4762a1bSJed Brown . x - input vector 295c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 296c4762a1bSJed Brown 297c4762a1bSJed Brown Output Parameters: 298c4762a1bSJed Brown . H - Hessian matrix 299c4762a1bSJed Brown . Hpre - optionally different preconditioning matrix 300c4762a1bSJed Brown . flg - flag indicating matrix structure 301c4762a1bSJed Brown 302c4762a1bSJed Brown */ 303c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 304c4762a1bSJed Brown { 305c4762a1bSJed Brown PetscErrorCode ierr; 306c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 307c4762a1bSJed Brown 308c4762a1bSJed Brown PetscFunctionBeginUser; 309c4762a1bSJed Brown /* Evaluate the Hessian entries*/ 310c4762a1bSJed Brown ierr = QuadraticH(user,X,H);CHKERRQ(ierr); 311c4762a1bSJed Brown PetscFunctionReturn(0); 312c4762a1bSJed Brown } 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 315c4762a1bSJed Brown /* 316c4762a1bSJed Brown QuadraticH - Evaluates the Hessian matrix. 317c4762a1bSJed Brown 318c4762a1bSJed Brown Input Parameters: 319c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian() 320c4762a1bSJed Brown . X - input vector 321c4762a1bSJed Brown 322c4762a1bSJed Brown Output Parameter: 323c4762a1bSJed Brown . H - Hessian matrix 324c4762a1bSJed Brown */ 325c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 326c4762a1bSJed Brown { 327c4762a1bSJed Brown PetscErrorCode ierr; 328c4762a1bSJed Brown PetscInt i,j,k,row; 329c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 330c4762a1bSJed Brown PetscInt col[7]; 331c4762a1bSJed Brown PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 332c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 333c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 334c4762a1bSJed Brown PetscReal hl,hr,ht,hb,hc,htl,hbr; 335c4762a1bSJed Brown const PetscScalar *x; 336c4762a1bSJed Brown PetscReal v[7]; 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* Get pointers to vector data */ 339c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* Initialize matrix entries to zero */ 342c4762a1bSJed Brown ierr = MatZeroEntries(Hessian); CHKERRQ(ierr); 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* Set various matrix options */ 345c4762a1bSJed Brown ierr = MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 346c4762a1bSJed Brown 347c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */ 348c4762a1bSJed Brown for (i=0; i< mx; i++){ 349c4762a1bSJed Brown for (j=0; j<my; j++){ 350c4762a1bSJed Brown 351c4762a1bSJed Brown row=(j)*mx + (i); 352c4762a1bSJed Brown 353c4762a1bSJed Brown xc = x[row]; 354c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* Left side */ 357c4762a1bSJed Brown if (i==0){ 358c4762a1bSJed Brown xl= user->left[j+1]; 359c4762a1bSJed Brown xlt = user->left[j+2]; 360c4762a1bSJed Brown } else { 361c4762a1bSJed Brown xl = x[row-1]; 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364c4762a1bSJed Brown if (j==0){ 365c4762a1bSJed Brown xb=user->bottom[i+1]; 366c4762a1bSJed Brown xrb = user->bottom[i+2]; 367c4762a1bSJed Brown } else { 368c4762a1bSJed Brown xb = x[row-mx]; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown 371c4762a1bSJed Brown if (i+1 == mx){ 372c4762a1bSJed Brown xr=user->right[j+1]; 373c4762a1bSJed Brown xrb = user->right[j]; 374c4762a1bSJed Brown } else { 375c4762a1bSJed Brown xr = x[row+1]; 376c4762a1bSJed Brown } 377c4762a1bSJed Brown 378c4762a1bSJed Brown if (j+1==my){ 379c4762a1bSJed Brown xt=user->top[i+1]; 380c4762a1bSJed Brown xlt = user->top[i]; 381c4762a1bSJed Brown }else { 382c4762a1bSJed Brown xt = x[row+mx]; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown if (i>0 && j+1<my){ 386c4762a1bSJed Brown xlt = x[row-1+mx]; 387c4762a1bSJed Brown } 388c4762a1bSJed Brown if (j>0 && i+1<mx){ 389c4762a1bSJed Brown xrb = x[row+1-mx]; 390c4762a1bSJed Brown } 391c4762a1bSJed Brown 392c4762a1bSJed Brown 393c4762a1bSJed Brown d1 = (xc-xl)*rhx; 394c4762a1bSJed Brown d2 = (xc-xr)*rhx; 395c4762a1bSJed Brown d3 = (xc-xt)*rhy; 396c4762a1bSJed Brown d4 = (xc-xb)*rhy; 397c4762a1bSJed Brown d5 = (xrb-xr)*rhy; 398c4762a1bSJed Brown d6 = (xrb-xb)*rhx; 399c4762a1bSJed Brown d7 = (xlt-xl)*rhy; 400c4762a1bSJed Brown d8 = (xlt-xt)*rhx; 401c4762a1bSJed Brown 402c4762a1bSJed Brown f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7); 403c4762a1bSJed Brown f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4); 404c4762a1bSJed Brown f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8); 405c4762a1bSJed Brown f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2); 406c4762a1bSJed Brown f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5); 407c4762a1bSJed Brown f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6); 408c4762a1bSJed Brown 409c4762a1bSJed Brown 410c4762a1bSJed Brown hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 411c4762a1bSJed Brown hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 412c4762a1bSJed Brown ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 413c4762a1bSJed Brown hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 414c4762a1bSJed Brown 415c4762a1bSJed Brown hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 416c4762a1bSJed Brown htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 417c4762a1bSJed Brown 418c4762a1bSJed Brown hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) + 419c4762a1bSJed Brown (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) + (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4); 420c4762a1bSJed Brown 421c4762a1bSJed Brown hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5; hc*=0.5; 422c4762a1bSJed Brown 423c4762a1bSJed Brown k=0; 424c4762a1bSJed Brown if (j>0){ 425c4762a1bSJed Brown v[k]=hb; col[k]=row - mx; k++; 426c4762a1bSJed Brown } 427c4762a1bSJed Brown 428c4762a1bSJed Brown if (j>0 && i < mx -1){ 429c4762a1bSJed Brown v[k]=hbr; col[k]=row - mx+1; k++; 430c4762a1bSJed Brown } 431c4762a1bSJed Brown 432c4762a1bSJed Brown if (i>0){ 433c4762a1bSJed Brown v[k]= hl; col[k]=row - 1; k++; 434c4762a1bSJed Brown } 435c4762a1bSJed Brown 436c4762a1bSJed Brown v[k]= hc; col[k]=row; k++; 437c4762a1bSJed Brown 438c4762a1bSJed Brown if (i < mx-1 ){ 439c4762a1bSJed Brown v[k]= hr; col[k]=row+1; k++; 440c4762a1bSJed Brown } 441c4762a1bSJed Brown 442c4762a1bSJed Brown if (i>0 && j < my-1 ){ 443c4762a1bSJed Brown v[k]= htl; col[k] = row+mx-1; k++; 444c4762a1bSJed Brown } 445c4762a1bSJed Brown 446c4762a1bSJed Brown if (j < my-1 ){ 447c4762a1bSJed Brown v[k]= ht; col[k] = row+mx; k++; 448c4762a1bSJed Brown } 449c4762a1bSJed Brown 450c4762a1bSJed Brown /* 451c4762a1bSJed Brown Set matrix values using local numbering, which was defined 452c4762a1bSJed Brown earlier, in the main routine. 453c4762a1bSJed Brown */ 454c4762a1bSJed Brown ierr = MatSetValues(Hessian,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr); 455c4762a1bSJed Brown } 456c4762a1bSJed Brown } 457c4762a1bSJed Brown 458c4762a1bSJed Brown /* Restore vectors */ 459c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 460c4762a1bSJed Brown 461c4762a1bSJed Brown /* Assemble the matrix */ 462c4762a1bSJed Brown ierr = MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 463c4762a1bSJed Brown ierr = MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 464c4762a1bSJed Brown 465*ca0c957dSBarry Smith ierr = PetscLogFlops(199.0*mx*my);CHKERRQ(ierr); 466c4762a1bSJed Brown PetscFunctionReturn(0); 467c4762a1bSJed Brown } 468c4762a1bSJed Brown 469c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 470c4762a1bSJed Brown /* 471c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 472c4762a1bSJed Brown the region. 473c4762a1bSJed Brown 474c4762a1bSJed Brown Input Parameter: 475c4762a1bSJed Brown . user - user-defined application context 476c4762a1bSJed Brown 477c4762a1bSJed Brown Output Parameter: 478c4762a1bSJed Brown . user - user-defined application context 479c4762a1bSJed Brown */ 480c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 481c4762a1bSJed Brown { 482c4762a1bSJed Brown PetscErrorCode ierr; 483c4762a1bSJed Brown PetscInt i,j,k,limit=0; 484c4762a1bSJed Brown PetscInt maxits=5; 485c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 486c4762a1bSJed Brown PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 487c4762a1bSJed Brown PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 488c4762a1bSJed Brown PetscReal fnorm,det,hx,hy,xt=0,yt=0; 489c4762a1bSJed Brown PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 490c4762a1bSJed Brown PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 491c4762a1bSJed Brown PetscReal *boundary; 492c4762a1bSJed Brown 493c4762a1bSJed Brown PetscFunctionBeginUser; 494c4762a1bSJed Brown bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2; 495c4762a1bSJed Brown 496c4762a1bSJed Brown ierr = PetscMalloc1(bsize,&user->bottom);CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = PetscMalloc1(tsize,&user->top);CHKERRQ(ierr); 498c4762a1bSJed Brown ierr = PetscMalloc1(lsize,&user->left);CHKERRQ(ierr); 499c4762a1bSJed Brown ierr = PetscMalloc1(rsize,&user->right);CHKERRQ(ierr); 500c4762a1bSJed Brown 501c4762a1bSJed Brown hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 502c4762a1bSJed Brown 503c4762a1bSJed Brown for (j=0; j<4; j++){ 504c4762a1bSJed Brown if (j==0){ 505c4762a1bSJed Brown yt=b; 506c4762a1bSJed Brown xt=l; 507c4762a1bSJed Brown limit=bsize; 508c4762a1bSJed Brown boundary=user->bottom; 509c4762a1bSJed Brown } else if (j==1){ 510c4762a1bSJed Brown yt=t; 511c4762a1bSJed Brown xt=l; 512c4762a1bSJed Brown limit=tsize; 513c4762a1bSJed Brown boundary=user->top; 514c4762a1bSJed Brown } else if (j==2){ 515c4762a1bSJed Brown yt=b; 516c4762a1bSJed Brown xt=l; 517c4762a1bSJed Brown limit=lsize; 518c4762a1bSJed Brown boundary=user->left; 519c4762a1bSJed Brown } else { /* if (j==3) */ 520c4762a1bSJed Brown yt=b; 521c4762a1bSJed Brown xt=r; 522c4762a1bSJed Brown limit=rsize; 523c4762a1bSJed Brown boundary=user->right; 524c4762a1bSJed Brown } 525c4762a1bSJed Brown 526c4762a1bSJed Brown for (i=0; i<limit; i++){ 527c4762a1bSJed Brown u1=xt; 528c4762a1bSJed Brown u2=-yt; 529c4762a1bSJed Brown for (k=0; k<maxits; k++){ 530c4762a1bSJed Brown nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 531c4762a1bSJed Brown nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 532c4762a1bSJed Brown fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2); 533c4762a1bSJed Brown if (fnorm <= tol) break; 534c4762a1bSJed Brown njac11=one+u2*u2-u1*u1; 535c4762a1bSJed Brown njac12=two*u1*u2; 536c4762a1bSJed Brown njac21=-two*u1*u2; 537c4762a1bSJed Brown njac22=-one - u1*u1 + u2*u2; 538c4762a1bSJed Brown det = njac11*njac22-njac21*njac12; 539c4762a1bSJed Brown u1 = u1-(njac22*nf1-njac12*nf2)/det; 540c4762a1bSJed Brown u2 = u2-(njac11*nf2-njac21*nf1)/det; 541c4762a1bSJed Brown } 542c4762a1bSJed Brown 543c4762a1bSJed Brown boundary[i]=u1*u1-u2*u2; 544c4762a1bSJed Brown if (j==0 || j==1) { 545c4762a1bSJed Brown xt=xt+hx; 546c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 547c4762a1bSJed Brown yt=yt+hy; 548c4762a1bSJed Brown } 549c4762a1bSJed Brown } 550c4762a1bSJed Brown } 551c4762a1bSJed Brown PetscFunctionReturn(0); 552c4762a1bSJed Brown } 553c4762a1bSJed Brown 554c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 555c4762a1bSJed Brown /* 556c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 557c4762a1bSJed Brown 558c4762a1bSJed Brown Input Parameters: 559c4762a1bSJed Brown . user - user-defined application context 560c4762a1bSJed Brown . X - vector for initial guess 561c4762a1bSJed Brown 562c4762a1bSJed Brown Output Parameters: 563c4762a1bSJed Brown . X - newly computed initial guess 564c4762a1bSJed Brown */ 565c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 566c4762a1bSJed Brown { 567c4762a1bSJed Brown PetscInt start=-1,i,j; 568c4762a1bSJed Brown PetscErrorCode ierr; 569c4762a1bSJed Brown PetscReal zero=0.0; 570c4762a1bSJed Brown PetscBool flg; 571c4762a1bSJed Brown 572c4762a1bSJed Brown PetscFunctionBeginUser; 573c4762a1bSJed Brown ierr = VecSet(X, zero);CHKERRQ(ierr); 574c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg);CHKERRQ(ierr); 575c4762a1bSJed Brown 576c4762a1bSJed Brown if (flg && start==0){ /* The zero vector is reasonable */ 577c4762a1bSJed Brown ierr = VecSet(X, zero);CHKERRQ(ierr); 578c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 579c4762a1bSJed Brown PetscInt row; 580c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 581c4762a1bSJed Brown PetscReal *x; 582c4762a1bSJed Brown 583c4762a1bSJed Brown /* Get pointers to vector data */ 584c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 585c4762a1bSJed Brown /* Perform local computations */ 586c4762a1bSJed Brown for (j=0; j<my; j++){ 587c4762a1bSJed Brown for (i=0; i< mx; i++){ 588c4762a1bSJed Brown row=(j)*mx + (i); 589c4762a1bSJed Brown x[row] = ( ((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+ ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0; 590c4762a1bSJed Brown } 591c4762a1bSJed Brown } 592c4762a1bSJed Brown /* Restore vectors */ 593c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 594c4762a1bSJed Brown } 595c4762a1bSJed Brown PetscFunctionReturn(0); 596c4762a1bSJed Brown } 597c4762a1bSJed Brown 598c4762a1bSJed Brown 599c4762a1bSJed Brown /*TEST 600c4762a1bSJed Brown 601c4762a1bSJed Brown build: 602c4762a1bSJed Brown requires: !complex 603c4762a1bSJed Brown 604c4762a1bSJed Brown test: 605c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4 606c4762a1bSJed Brown requires: !single 607c4762a1bSJed Brown 608c4762a1bSJed Brown test: 609c4762a1bSJed Brown suffix: 2 610c4762a1bSJed Brown args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 611c4762a1bSJed Brown requires: !single 612c4762a1bSJed Brown 613c4762a1bSJed Brown test: 614c4762a1bSJed Brown suffix: 3 615c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 616c4762a1bSJed Brown requires: !single 617c4762a1bSJed Brown 618c4762a1bSJed Brown test: 619c4762a1bSJed Brown suffix: 4 620c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 621c4762a1bSJed Brown requires: !single 622c4762a1bSJed Brown 623c4762a1bSJed Brown test: 624c4762a1bSJed Brown suffix: 5 625c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 626c4762a1bSJed Brown requires: !single 627c4762a1bSJed Brown 628c4762a1bSJed Brown test: 629c4762a1bSJed Brown suffix: 6 630c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4 631c4762a1bSJed Brown requires: !single 632c4762a1bSJed Brown 633c4762a1bSJed Brown test: 634c4762a1bSJed Brown suffix: 7 635c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 636c4762a1bSJed Brown requires: !single 637c4762a1bSJed Brown 638c4762a1bSJed Brown test: 639c4762a1bSJed Brown suffix: 8 640c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 641c4762a1bSJed Brown requires: !single 642c4762a1bSJed Brown 643c4762a1bSJed Brown test: 644c4762a1bSJed Brown suffix: 9 645c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 646c4762a1bSJed Brown requires: !single 647c4762a1bSJed Brown 648c4762a1bSJed Brown TEST*/ 649