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