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 */ 47*327415f7SBarry Smith PetscFunctionBeginUser; 489566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv,(char *)0,help)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size)); 513c859ba3SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Specify default dimension of the problem */ 54c4762a1bSJed Brown user.mx = 4; user.my = 4; 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg)); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n")); 6163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n",user.mx,user.my)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Calculate any derived values from parameters */ 64c4762a1bSJed Brown N = user.mx*user.my; 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 679566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF,&tao)); 689566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOLMVM)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Initialize minsurf application data structure for use in the function evaluations */ 719566063dSJacob Faibussowitsch PetscCall(MSA_BoundaryConditions(&user)); /* Application specific routine */ 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* 74c4762a1bSJed Brown Create a vector to hold the variables. Compute an initial solution. 75c4762a1bSJed Brown Set this vector, which will be used by TAO. 76c4762a1bSJed Brown */ 779566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&x)); 789566063dSJacob Faibussowitsch PetscCall(MSA_InitialPoint(&user,x)); /* Application specific routine */ 799566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,x)); /* A TAO routine */ 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Provide TAO routines for function, gradient, and Hessian evaluation */ 829566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Create a matrix data structure to store the Hessian. This structure will be used by TAO */ 859566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,7,NULL,&(user.H))); 869566063dSJacob Faibussowitsch PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE)); 879566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* Check for any TAO command line options */ 909566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Limit the number of iterations in the KSP linear solver */ 939566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao,&ksp)); 941baa6e33SBarry Smith if (ksp) PetscCall(KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,user.mx*user.my)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 979566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.H)); 1029566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom)); 1039566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top)); 1049566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left)); 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right)); 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 108b122ec5aSJacob Faibussowitsch return 0; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* FormFunctionGradient - Evaluates function and corresponding gradient. 114c4762a1bSJed Brown 115c4762a1bSJed Brown Input Parameters: 116c4762a1bSJed Brown . tao - the Tao context 117c4762a1bSJed Brown . X - input vector 118c4762a1bSJed Brown . userCtx - optional user-defined context, as set by TaoSetFunctionGradient() 119c4762a1bSJed Brown 120c4762a1bSJed Brown Output Parameters: 121c4762a1bSJed Brown . fcn - the newly evaluated function 122c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 123c4762a1bSJed Brown */ 124c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *fcn,Vec G,void *userCtx) 125c4762a1bSJed Brown { 126c4762a1bSJed Brown AppCtx *user = (AppCtx *) userCtx; 127c4762a1bSJed Brown PetscInt i,j,row; 128c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 129c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 130c4762a1bSJed 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; 131c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 132c4762a1bSJed Brown PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 133c4762a1bSJed Brown PetscReal zero=0.0; 134c4762a1bSJed Brown PetscScalar *g; 135c4762a1bSJed Brown const PetscScalar *x; 136c4762a1bSJed Brown 137c4762a1bSJed Brown PetscFunctionBeginUser; 1389566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 1419566063dSJacob Faibussowitsch PetscCall(VecGetArray(G,&g)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */ 144c4762a1bSJed Brown for (j=0; j<my; j++) { 145c4762a1bSJed Brown for (i=0; i< mx; i++) { 146c4762a1bSJed Brown row=(j)*mx + (i); 147c4762a1bSJed Brown xc = x[row]; 148c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 149c4762a1bSJed Brown if (i==0) { /* left side */ 150c4762a1bSJed Brown xl = user->left[j+1]; 151c4762a1bSJed Brown xlt = user->left[j+2]; 152c4762a1bSJed Brown } else { 153c4762a1bSJed Brown xl = x[row-1]; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown if (j==0) { /* bottom side */ 157c4762a1bSJed Brown xb = user->bottom[i+1]; 158c4762a1bSJed Brown xrb = user->bottom[i+2]; 159c4762a1bSJed Brown } else { 160c4762a1bSJed Brown xb = x[row-mx]; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown if (i+1 == mx) { /* right side */ 164c4762a1bSJed Brown xr = user->right[j+1]; 165c4762a1bSJed Brown xrb = user->right[j]; 166c4762a1bSJed Brown } else { 167c4762a1bSJed Brown xr = x[row+1]; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown if (j+1==0+my) { /* top side */ 171c4762a1bSJed Brown xt = user->top[i+1]; 172c4762a1bSJed Brown xlt = user->top[i]; 173c4762a1bSJed Brown }else { 174c4762a1bSJed Brown xt = x[row+mx]; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown if (i>0 && j+1<my) { 178c4762a1bSJed Brown xlt = x[row-1+mx]; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown if (j>0 && i+1<mx) { 181c4762a1bSJed Brown xrb = x[row+1-mx]; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown d1 = (xc-xl); 185c4762a1bSJed Brown d2 = (xc-xr); 186c4762a1bSJed Brown d3 = (xc-xt); 187c4762a1bSJed Brown d4 = (xc-xb); 188c4762a1bSJed Brown d5 = (xr-xrb); 189c4762a1bSJed Brown d6 = (xrb-xb); 190c4762a1bSJed Brown d7 = (xlt-xl); 191c4762a1bSJed Brown d8 = (xt-xlt); 192c4762a1bSJed Brown 193c4762a1bSJed Brown df1dxc = d1*hydhx; 194c4762a1bSJed Brown df2dxc = (d1*hydhx + d4*hxdhy); 195c4762a1bSJed Brown df3dxc = d3*hxdhy; 196c4762a1bSJed Brown df4dxc = (d2*hydhx + d3*hxdhy); 197c4762a1bSJed Brown df5dxc = d2*hydhx; 198c4762a1bSJed Brown df6dxc = d4*hxdhy; 199c4762a1bSJed Brown 200c4762a1bSJed Brown d1 *= rhx; 201c4762a1bSJed Brown d2 *= rhx; 202c4762a1bSJed Brown d3 *= rhy; 203c4762a1bSJed Brown d4 *= rhy; 204c4762a1bSJed Brown d5 *= rhy; 205c4762a1bSJed Brown d6 *= rhx; 206c4762a1bSJed Brown d7 *= rhy; 207c4762a1bSJed Brown d8 *= rhx; 208c4762a1bSJed Brown 209c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 210c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 211c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 212c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 213c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 214c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 215c4762a1bSJed Brown 216c4762a1bSJed Brown ft = ft + (f2 + f4); 217c4762a1bSJed Brown 218c4762a1bSJed Brown df1dxc /= f1; 219c4762a1bSJed Brown df2dxc /= f2; 220c4762a1bSJed Brown df3dxc /= f3; 221c4762a1bSJed Brown df4dxc /= f4; 222c4762a1bSJed Brown df5dxc /= f5; 223c4762a1bSJed Brown df6dxc /= f6; 224c4762a1bSJed Brown 225c4762a1bSJed Brown g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0; 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown for (j=0; j<my; j++) { /* left side */ 230c4762a1bSJed Brown d3 = (user->left[j+1] - user->left[j+2])*rhy; 231c4762a1bSJed Brown d2 = (user->left[j+1] - x[j*mx])*rhx; 232c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown 235c4762a1bSJed Brown for (i=0; i<mx; i++) { /* bottom */ 236c4762a1bSJed Brown d2 = (user->bottom[i+1]-user->bottom[i+2])*rhx; 237c4762a1bSJed Brown d3 = (user->bottom[i+1]-x[i])*rhy; 238c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown for (j=0; j< my; j++) { /* right side */ 242c4762a1bSJed Brown d1 = (x[(j+1)*mx-1]-user->right[j+1])*rhx; 243c4762a1bSJed Brown d4 = (user->right[j]-user->right[j+1])*rhy; 244c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown for (i=0; i<mx; i++) { /* top side */ 248c4762a1bSJed Brown d1 = (x[(my-1)*mx + i] - user->top[i+1])*rhy; 249c4762a1bSJed Brown d4 = (user->top[i+1] - user->top[i])*rhx; 250c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* Bottom left corner */ 254c4762a1bSJed Brown d1 = (user->left[0]-user->left[1])*rhy; 255c4762a1bSJed Brown d2 = (user->bottom[0]-user->bottom[1])*rhx; 256c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* Top right corner */ 259c4762a1bSJed Brown d1 = (user->right[my+1] - user->right[my])*rhy; 260c4762a1bSJed Brown d2 = (user->top[mx+1] - user->top[mx])*rhx; 261c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 262c4762a1bSJed Brown 263c4762a1bSJed Brown (*fcn)=ft*area; 264c4762a1bSJed Brown 265c4762a1bSJed Brown /* Restore vectors */ 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G,&g)); 2689566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67.0*mx*my)); 269c4762a1bSJed Brown PetscFunctionReturn(0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 273c4762a1bSJed Brown /* 274c4762a1bSJed Brown FormHessian - Evaluates the Hessian matrix. 275c4762a1bSJed Brown 276c4762a1bSJed Brown Input Parameters: 277c4762a1bSJed Brown . tao - the Tao context 278c4762a1bSJed Brown . x - input vector 279c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 280c4762a1bSJed Brown 281c4762a1bSJed Brown Output Parameters: 282c4762a1bSJed Brown . H - Hessian matrix 283c4762a1bSJed Brown . Hpre - optionally different preconditioning matrix 284c4762a1bSJed Brown . flg - flag indicating matrix structure 285c4762a1bSJed Brown 286c4762a1bSJed Brown */ 287c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 288c4762a1bSJed Brown { 289c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 290c4762a1bSJed Brown 291c4762a1bSJed Brown PetscFunctionBeginUser; 292c4762a1bSJed Brown /* Evaluate the Hessian entries*/ 2939566063dSJacob Faibussowitsch PetscCall(QuadraticH(user,X,H)); 294c4762a1bSJed Brown PetscFunctionReturn(0); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 298c4762a1bSJed Brown /* 299c4762a1bSJed Brown QuadraticH - Evaluates the Hessian matrix. 300c4762a1bSJed Brown 301c4762a1bSJed Brown Input Parameters: 302c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian() 303c4762a1bSJed Brown . X - input vector 304c4762a1bSJed Brown 305c4762a1bSJed Brown Output Parameter: 306c4762a1bSJed Brown . H - Hessian matrix 307c4762a1bSJed Brown */ 308c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 309c4762a1bSJed Brown { 310c4762a1bSJed Brown PetscInt i,j,k,row; 311c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 312c4762a1bSJed Brown PetscInt col[7]; 313c4762a1bSJed Brown PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 314c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 315c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 316c4762a1bSJed Brown PetscReal hl,hr,ht,hb,hc,htl,hbr; 317c4762a1bSJed Brown const PetscScalar *x; 318c4762a1bSJed Brown PetscReal v[7]; 319c4762a1bSJed Brown 320362febeeSStefano Zampini PetscFunctionBeginUser; 321c4762a1bSJed Brown /* Get pointers to vector data */ 3229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 323c4762a1bSJed Brown 324c4762a1bSJed Brown /* Initialize matrix entries to zero */ 3259566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hessian)); 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* Set various matrix options */ 3289566063dSJacob Faibussowitsch PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */ 331c4762a1bSJed Brown for (i=0; i< mx; i++) { 332c4762a1bSJed Brown for (j=0; j<my; j++) { 333c4762a1bSJed Brown 334c4762a1bSJed Brown row=(j)*mx + (i); 335c4762a1bSJed Brown 336c4762a1bSJed Brown xc = x[row]; 337c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 338c4762a1bSJed Brown 339c4762a1bSJed Brown /* Left side */ 340c4762a1bSJed Brown if (i==0) { 341c4762a1bSJed Brown xl= user->left[j+1]; 342c4762a1bSJed Brown xlt = user->left[j+2]; 343c4762a1bSJed Brown } else { 344c4762a1bSJed Brown xl = x[row-1]; 345c4762a1bSJed Brown } 346c4762a1bSJed Brown 347c4762a1bSJed Brown if (j==0) { 348c4762a1bSJed Brown xb=user->bottom[i+1]; 349c4762a1bSJed Brown xrb = user->bottom[i+2]; 350c4762a1bSJed Brown } else { 351c4762a1bSJed Brown xb = x[row-mx]; 352c4762a1bSJed Brown } 353c4762a1bSJed Brown 354c4762a1bSJed Brown if (i+1 == mx) { 355c4762a1bSJed Brown xr=user->right[j+1]; 356c4762a1bSJed Brown xrb = user->right[j]; 357c4762a1bSJed Brown } else { 358c4762a1bSJed Brown xr = x[row+1]; 359c4762a1bSJed Brown } 360c4762a1bSJed Brown 361c4762a1bSJed Brown if (j+1==my) { 362c4762a1bSJed Brown xt=user->top[i+1]; 363c4762a1bSJed Brown xlt = user->top[i]; 364c4762a1bSJed Brown }else { 365c4762a1bSJed Brown xt = x[row+mx]; 366c4762a1bSJed Brown } 367c4762a1bSJed Brown 368c4762a1bSJed Brown if (i>0 && j+1<my) { 369c4762a1bSJed Brown xlt = x[row-1+mx]; 370c4762a1bSJed Brown } 371c4762a1bSJed Brown if (j>0 && i+1<mx) { 372c4762a1bSJed Brown xrb = x[row+1-mx]; 373c4762a1bSJed Brown } 374c4762a1bSJed Brown 375c4762a1bSJed Brown d1 = (xc-xl)*rhx; 376c4762a1bSJed Brown d2 = (xc-xr)*rhx; 377c4762a1bSJed Brown d3 = (xc-xt)*rhy; 378c4762a1bSJed Brown d4 = (xc-xb)*rhy; 379c4762a1bSJed Brown d5 = (xrb-xr)*rhy; 380c4762a1bSJed Brown d6 = (xrb-xb)*rhx; 381c4762a1bSJed Brown d7 = (xlt-xl)*rhy; 382c4762a1bSJed Brown d8 = (xlt-xt)*rhx; 383c4762a1bSJed Brown 384c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 385c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 386c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 387c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 388c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 389c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 390c4762a1bSJed Brown 391c4762a1bSJed Brown hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 392c4762a1bSJed Brown hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 393c4762a1bSJed Brown ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 394c4762a1bSJed Brown hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 395c4762a1bSJed Brown 396c4762a1bSJed Brown hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 397c4762a1bSJed Brown htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 398c4762a1bSJed Brown 399c4762a1bSJed 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) + 400c4762a1bSJed 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); 401c4762a1bSJed Brown 402c4762a1bSJed Brown hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5; hc*=0.5; 403c4762a1bSJed Brown 404c4762a1bSJed Brown k=0; 405c4762a1bSJed Brown if (j>0) { 406c4762a1bSJed Brown v[k]=hb; col[k]=row - mx; k++; 407c4762a1bSJed Brown } 408c4762a1bSJed Brown 409c4762a1bSJed Brown if (j>0 && i < mx -1) { 410c4762a1bSJed Brown v[k]=hbr; col[k]=row - mx+1; k++; 411c4762a1bSJed Brown } 412c4762a1bSJed Brown 413c4762a1bSJed Brown if (i>0) { 414c4762a1bSJed Brown v[k]= hl; col[k]=row - 1; k++; 415c4762a1bSJed Brown } 416c4762a1bSJed Brown 417c4762a1bSJed Brown v[k]= hc; col[k]=row; k++; 418c4762a1bSJed Brown 419c4762a1bSJed Brown if (i < mx-1) { 420c4762a1bSJed Brown v[k]= hr; col[k]=row+1; k++; 421c4762a1bSJed Brown } 422c4762a1bSJed Brown 423c4762a1bSJed Brown if (i>0 && j < my-1) { 424c4762a1bSJed Brown v[k]= htl; col[k] = row+mx-1; k++; 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427c4762a1bSJed Brown if (j < my-1) { 428c4762a1bSJed Brown v[k]= ht; col[k] = row+mx; k++; 429c4762a1bSJed Brown } 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* 432c4762a1bSJed Brown Set matrix values using local numbering, which was defined 433c4762a1bSJed Brown earlier, in the main routine. 434c4762a1bSJed Brown */ 4359566063dSJacob Faibussowitsch PetscCall(MatSetValues(Hessian,1,&row,k,col,v,INSERT_VALUES)); 436c4762a1bSJed Brown } 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* Restore vectors */ 4409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 441c4762a1bSJed Brown 442c4762a1bSJed Brown /* Assemble the matrix */ 4439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY)); 4449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY)); 445c4762a1bSJed Brown 4469566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199.0*mx*my)); 447c4762a1bSJed Brown PetscFunctionReturn(0); 448c4762a1bSJed Brown } 449c4762a1bSJed Brown 450c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 451c4762a1bSJed Brown /* 452c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 453c4762a1bSJed Brown the region. 454c4762a1bSJed Brown 455c4762a1bSJed Brown Input Parameter: 456c4762a1bSJed Brown . user - user-defined application context 457c4762a1bSJed Brown 458c4762a1bSJed Brown Output Parameter: 459c4762a1bSJed Brown . user - user-defined application context 460c4762a1bSJed Brown */ 461c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 462c4762a1bSJed Brown { 463c4762a1bSJed Brown PetscInt i,j,k,limit=0; 464c4762a1bSJed Brown PetscInt maxits=5; 465c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 466c4762a1bSJed Brown PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 467c4762a1bSJed Brown PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 468c4762a1bSJed Brown PetscReal fnorm,det,hx,hy,xt=0,yt=0; 469c4762a1bSJed Brown PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 470c4762a1bSJed Brown PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 471c4762a1bSJed Brown PetscReal *boundary; 472c4762a1bSJed Brown 473c4762a1bSJed Brown PetscFunctionBeginUser; 474c4762a1bSJed Brown bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2; 475c4762a1bSJed Brown 4769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize,&user->bottom)); 4779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize,&user->top)); 4789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize,&user->left)); 4799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize,&user->right)); 480c4762a1bSJed Brown 481c4762a1bSJed Brown hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 482c4762a1bSJed Brown 483c4762a1bSJed Brown for (j=0; j<4; j++) { 484c4762a1bSJed Brown if (j==0) { 485c4762a1bSJed Brown yt=b; 486c4762a1bSJed Brown xt=l; 487c4762a1bSJed Brown limit=bsize; 488c4762a1bSJed Brown boundary=user->bottom; 489c4762a1bSJed Brown } else if (j==1) { 490c4762a1bSJed Brown yt=t; 491c4762a1bSJed Brown xt=l; 492c4762a1bSJed Brown limit=tsize; 493c4762a1bSJed Brown boundary=user->top; 494c4762a1bSJed Brown } else if (j==2) { 495c4762a1bSJed Brown yt=b; 496c4762a1bSJed Brown xt=l; 497c4762a1bSJed Brown limit=lsize; 498c4762a1bSJed Brown boundary=user->left; 499c4762a1bSJed Brown } else { /* if (j==3) */ 500c4762a1bSJed Brown yt=b; 501c4762a1bSJed Brown xt=r; 502c4762a1bSJed Brown limit=rsize; 503c4762a1bSJed Brown boundary=user->right; 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown for (i=0; i<limit; i++) { 507c4762a1bSJed Brown u1=xt; 508c4762a1bSJed Brown u2=-yt; 509c4762a1bSJed Brown for (k=0; k<maxits; k++) { 510c4762a1bSJed Brown nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 511c4762a1bSJed Brown nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 512c4762a1bSJed Brown fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2); 513c4762a1bSJed Brown if (fnorm <= tol) break; 514c4762a1bSJed Brown njac11=one+u2*u2-u1*u1; 515c4762a1bSJed Brown njac12=two*u1*u2; 516c4762a1bSJed Brown njac21=-two*u1*u2; 517c4762a1bSJed Brown njac22=-one - u1*u1 + u2*u2; 518c4762a1bSJed Brown det = njac11*njac22-njac21*njac12; 519c4762a1bSJed Brown u1 = u1-(njac22*nf1-njac12*nf2)/det; 520c4762a1bSJed Brown u2 = u2-(njac11*nf2-njac21*nf1)/det; 521c4762a1bSJed Brown } 522c4762a1bSJed Brown 523c4762a1bSJed Brown boundary[i]=u1*u1-u2*u2; 524c4762a1bSJed Brown if (j==0 || j==1) { 525c4762a1bSJed Brown xt=xt+hx; 526c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 527c4762a1bSJed Brown yt=yt+hy; 528c4762a1bSJed Brown } 529c4762a1bSJed Brown } 530c4762a1bSJed Brown } 531c4762a1bSJed Brown PetscFunctionReturn(0); 532c4762a1bSJed Brown } 533c4762a1bSJed Brown 534c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 535c4762a1bSJed Brown /* 536c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 537c4762a1bSJed Brown 538c4762a1bSJed Brown Input Parameters: 539c4762a1bSJed Brown . user - user-defined application context 540c4762a1bSJed Brown . X - vector for initial guess 541c4762a1bSJed Brown 542c4762a1bSJed Brown Output Parameters: 543c4762a1bSJed Brown . X - newly computed initial guess 544c4762a1bSJed Brown */ 545c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 546c4762a1bSJed Brown { 547c4762a1bSJed Brown PetscInt start=-1,i,j; 548c4762a1bSJed Brown PetscReal zero=0.0; 549c4762a1bSJed Brown PetscBool flg; 550c4762a1bSJed Brown 551c4762a1bSJed Brown PetscFunctionBeginUser; 5529566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero)); 5539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg)); 554c4762a1bSJed Brown 555c4762a1bSJed Brown if (flg && start==0) { /* The zero vector is reasonable */ 5569566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero)); 557c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 558c4762a1bSJed Brown PetscInt row; 559c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 560c4762a1bSJed Brown PetscReal *x; 561c4762a1bSJed Brown 562c4762a1bSJed Brown /* Get pointers to vector data */ 5639566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 564c4762a1bSJed Brown /* Perform local computations */ 565c4762a1bSJed Brown for (j=0; j<my; j++) { 566c4762a1bSJed Brown for (i=0; i< mx; i++) { 567c4762a1bSJed Brown row=(j)*mx + (i); 568c4762a1bSJed 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; 569c4762a1bSJed Brown } 570c4762a1bSJed Brown } 571c4762a1bSJed Brown /* Restore vectors */ 5729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x)); 573c4762a1bSJed Brown } 574c4762a1bSJed Brown PetscFunctionReturn(0); 575c4762a1bSJed Brown } 576c4762a1bSJed Brown 577c4762a1bSJed Brown /*TEST 578c4762a1bSJed Brown 579c4762a1bSJed Brown build: 580c4762a1bSJed Brown requires: !complex 581c4762a1bSJed Brown 582c4762a1bSJed Brown test: 583c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4 584c4762a1bSJed Brown requires: !single 585c4762a1bSJed Brown 586c4762a1bSJed Brown test: 587c4762a1bSJed Brown suffix: 2 588c4762a1bSJed Brown args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 589c4762a1bSJed Brown requires: !single 590c4762a1bSJed Brown 591c4762a1bSJed Brown test: 592c4762a1bSJed Brown suffix: 3 593c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 594c4762a1bSJed Brown requires: !single 595c4762a1bSJed Brown 596c4762a1bSJed Brown test: 597c4762a1bSJed Brown suffix: 4 598c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 599c4762a1bSJed Brown requires: !single 600c4762a1bSJed Brown 601c4762a1bSJed Brown test: 6020f0abf79SStefano Zampini suffix: 4_ew 6030f0abf79SStefano Zampini args: -tao_smonitor -tao_type bntr -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4 6040f0abf79SStefano Zampini requires: !single 6050f0abf79SStefano Zampini 6060f0abf79SStefano Zampini test: 607c4762a1bSJed Brown suffix: 5 608c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 609c4762a1bSJed Brown requires: !single 610c4762a1bSJed Brown 611c4762a1bSJed Brown test: 6120f0abf79SStefano Zampini suffix: 5_ew 6130f0abf79SStefano Zampini args: -tao_smonitor -tao_type bntl -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4 6140f0abf79SStefano Zampini requires: !single 6150f0abf79SStefano Zampini 6160f0abf79SStefano Zampini test: 617c4762a1bSJed Brown suffix: 6 618c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4 619c4762a1bSJed Brown requires: !single 620c4762a1bSJed Brown 621c4762a1bSJed Brown test: 6220f0abf79SStefano Zampini suffix: 6_ew 6230f0abf79SStefano Zampini args: -tao_smonitor -tao_type bnls -tao_ksp_ew -tao_bnk_ksp_ew_version 3 -mx 10 -my 8 -tao_gatol 1.e-4 6240f0abf79SStefano Zampini requires: !single 6250f0abf79SStefano Zampini 6260f0abf79SStefano Zampini test: 627c4762a1bSJed Brown suffix: 7 628c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 629c4762a1bSJed Brown requires: !single 630c4762a1bSJed Brown 631c4762a1bSJed Brown test: 632c4762a1bSJed Brown suffix: 8 633c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 634c4762a1bSJed Brown requires: !single 635c4762a1bSJed Brown 636c4762a1bSJed Brown test: 637c4762a1bSJed Brown suffix: 9 638c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 639c4762a1bSJed Brown requires: !single 640c4762a1bSJed Brown 64134ad9904SAlp Dener test: 64234ad9904SAlp Dener suffix: 10 64334ad9904SAlp 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 64434ad9904SAlp Dener 64534ad9904SAlp Dener test: 64634ad9904SAlp Dener suffix: 11 64734ad9904SAlp Dener args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian 64834ad9904SAlp Dener requires: !single 64934ad9904SAlp Dener 65034ad9904SAlp Dener test: 65134ad9904SAlp Dener suffix: 12 65234ad9904SAlp Dener args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian 65334ad9904SAlp Dener requires: !single 65434ad9904SAlp Dener 655c4762a1bSJed Brown TEST*/ 656