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