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")); 6063a3b9bcSJacob 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)); 93*1baa6e33SBarry Smith if (ksp) PetscCall(KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,user.mx*user.my)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 969566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.H)); 1019566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom)); 1029566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top)); 1039566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left)); 1049566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right)); 105c4762a1bSJed Brown 1069566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 107b122ec5aSJacob Faibussowitsch return 0; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* FormFunctionGradient - Evaluates function and corresponding gradient. 113c4762a1bSJed Brown 114c4762a1bSJed Brown Input Parameters: 115c4762a1bSJed Brown . tao - the Tao context 116c4762a1bSJed Brown . X - input vector 117c4762a1bSJed Brown . userCtx - optional user-defined context, as set by TaoSetFunctionGradient() 118c4762a1bSJed Brown 119c4762a1bSJed Brown Output Parameters: 120c4762a1bSJed Brown . fcn - the newly evaluated function 121c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 122c4762a1bSJed Brown */ 123c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *fcn,Vec G,void *userCtx) 124c4762a1bSJed Brown { 125c4762a1bSJed Brown AppCtx *user = (AppCtx *) userCtx; 126c4762a1bSJed Brown PetscInt i,j,row; 127c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 128c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 129c4762a1bSJed 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; 130c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 131c4762a1bSJed Brown PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 132c4762a1bSJed Brown PetscReal zero=0.0; 133c4762a1bSJed Brown PetscScalar *g; 134c4762a1bSJed Brown const PetscScalar *x; 135c4762a1bSJed Brown 136c4762a1bSJed Brown PetscFunctionBeginUser; 1379566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 1409566063dSJacob Faibussowitsch PetscCall(VecGetArray(G,&g)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */ 143c4762a1bSJed Brown for (j=0; j<my; j++) { 144c4762a1bSJed Brown for (i=0; i< mx; i++) { 145c4762a1bSJed Brown row=(j)*mx + (i); 146c4762a1bSJed Brown xc = x[row]; 147c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 148c4762a1bSJed Brown if (i==0) { /* left side */ 149c4762a1bSJed Brown xl = user->left[j+1]; 150c4762a1bSJed Brown xlt = user->left[j+2]; 151c4762a1bSJed Brown } else { 152c4762a1bSJed Brown xl = x[row-1]; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown if (j==0) { /* bottom side */ 156c4762a1bSJed Brown xb = user->bottom[i+1]; 157c4762a1bSJed Brown xrb = user->bottom[i+2]; 158c4762a1bSJed Brown } else { 159c4762a1bSJed Brown xb = x[row-mx]; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown if (i+1 == mx) { /* right side */ 163c4762a1bSJed Brown xr = user->right[j+1]; 164c4762a1bSJed Brown xrb = user->right[j]; 165c4762a1bSJed Brown } else { 166c4762a1bSJed Brown xr = x[row+1]; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown if (j+1==0+my) { /* top side */ 170c4762a1bSJed Brown xt = user->top[i+1]; 171c4762a1bSJed Brown xlt = user->top[i]; 172c4762a1bSJed Brown }else { 173c4762a1bSJed Brown xt = x[row+mx]; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown if (i>0 && j+1<my) { 177c4762a1bSJed Brown xlt = x[row-1+mx]; 178c4762a1bSJed Brown } 179c4762a1bSJed Brown if (j>0 && i+1<mx) { 180c4762a1bSJed Brown xrb = x[row+1-mx]; 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown d1 = (xc-xl); 184c4762a1bSJed Brown d2 = (xc-xr); 185c4762a1bSJed Brown d3 = (xc-xt); 186c4762a1bSJed Brown d4 = (xc-xb); 187c4762a1bSJed Brown d5 = (xr-xrb); 188c4762a1bSJed Brown d6 = (xrb-xb); 189c4762a1bSJed Brown d7 = (xlt-xl); 190c4762a1bSJed Brown d8 = (xt-xlt); 191c4762a1bSJed Brown 192c4762a1bSJed Brown df1dxc = d1*hydhx; 193c4762a1bSJed Brown df2dxc = (d1*hydhx + d4*hxdhy); 194c4762a1bSJed Brown df3dxc = d3*hxdhy; 195c4762a1bSJed Brown df4dxc = (d2*hydhx + d3*hxdhy); 196c4762a1bSJed Brown df5dxc = d2*hydhx; 197c4762a1bSJed Brown df6dxc = d4*hxdhy; 198c4762a1bSJed Brown 199c4762a1bSJed Brown d1 *= rhx; 200c4762a1bSJed Brown d2 *= rhx; 201c4762a1bSJed Brown d3 *= rhy; 202c4762a1bSJed Brown d4 *= rhy; 203c4762a1bSJed Brown d5 *= rhy; 204c4762a1bSJed Brown d6 *= rhx; 205c4762a1bSJed Brown d7 *= rhy; 206c4762a1bSJed Brown d8 *= rhx; 207c4762a1bSJed Brown 208c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 209c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 210c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 211c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 212c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 213c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 214c4762a1bSJed Brown 215c4762a1bSJed Brown ft = ft + (f2 + f4); 216c4762a1bSJed Brown 217c4762a1bSJed Brown df1dxc /= f1; 218c4762a1bSJed Brown df2dxc /= f2; 219c4762a1bSJed Brown df3dxc /= f3; 220c4762a1bSJed Brown df4dxc /= f4; 221c4762a1bSJed Brown df5dxc /= f5; 222c4762a1bSJed Brown df6dxc /= f6; 223c4762a1bSJed Brown 224c4762a1bSJed Brown g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0; 225c4762a1bSJed Brown } 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228c4762a1bSJed Brown for (j=0; j<my; j++) { /* left side */ 229c4762a1bSJed Brown d3 = (user->left[j+1] - user->left[j+2])*rhy; 230c4762a1bSJed Brown d2 = (user->left[j+1] - x[j*mx])*rhx; 231c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown for (i=0; i<mx; i++) { /* bottom */ 235c4762a1bSJed Brown d2 = (user->bottom[i+1]-user->bottom[i+2])*rhx; 236c4762a1bSJed Brown d3 = (user->bottom[i+1]-x[i])*rhy; 237c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown for (j=0; j< my; j++) { /* right side */ 241c4762a1bSJed Brown d1 = (x[(j+1)*mx-1]-user->right[j+1])*rhx; 242c4762a1bSJed Brown d4 = (user->right[j]-user->right[j+1])*rhy; 243c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246c4762a1bSJed Brown for (i=0; i<mx; i++) { /* top side */ 247c4762a1bSJed Brown d1 = (x[(my-1)*mx + i] - user->top[i+1])*rhy; 248c4762a1bSJed Brown d4 = (user->top[i+1] - user->top[i])*rhx; 249c4762a1bSJed Brown ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* Bottom left corner */ 253c4762a1bSJed Brown d1 = (user->left[0]-user->left[1])*rhy; 254c4762a1bSJed Brown d2 = (user->bottom[0]-user->bottom[1])*rhx; 255c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* Top right corner */ 258c4762a1bSJed Brown d1 = (user->right[my+1] - user->right[my])*rhy; 259c4762a1bSJed Brown d2 = (user->top[mx+1] - user->top[mx])*rhx; 260c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 261c4762a1bSJed Brown 262c4762a1bSJed Brown (*fcn)=ft*area; 263c4762a1bSJed Brown 264c4762a1bSJed Brown /* Restore vectors */ 2659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G,&g)); 2679566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67.0*mx*my)); 268c4762a1bSJed Brown PetscFunctionReturn(0); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown FormHessian - Evaluates the Hessian matrix. 274c4762a1bSJed Brown 275c4762a1bSJed Brown Input Parameters: 276c4762a1bSJed Brown . tao - the Tao context 277c4762a1bSJed Brown . x - input vector 278c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 279c4762a1bSJed Brown 280c4762a1bSJed Brown Output Parameters: 281c4762a1bSJed Brown . H - Hessian matrix 282c4762a1bSJed Brown . Hpre - optionally different preconditioning matrix 283c4762a1bSJed Brown . flg - flag indicating matrix structure 284c4762a1bSJed Brown 285c4762a1bSJed Brown */ 286c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 287c4762a1bSJed Brown { 288c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 289c4762a1bSJed Brown 290c4762a1bSJed Brown PetscFunctionBeginUser; 291c4762a1bSJed Brown /* Evaluate the Hessian entries*/ 2929566063dSJacob Faibussowitsch PetscCall(QuadraticH(user,X,H)); 293c4762a1bSJed Brown PetscFunctionReturn(0); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 297c4762a1bSJed Brown /* 298c4762a1bSJed Brown QuadraticH - Evaluates the Hessian matrix. 299c4762a1bSJed Brown 300c4762a1bSJed Brown Input Parameters: 301c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian() 302c4762a1bSJed Brown . X - input vector 303c4762a1bSJed Brown 304c4762a1bSJed Brown Output Parameter: 305c4762a1bSJed Brown . H - Hessian matrix 306c4762a1bSJed Brown */ 307c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 308c4762a1bSJed Brown { 309c4762a1bSJed Brown PetscInt i,j,k,row; 310c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 311c4762a1bSJed Brown PetscInt col[7]; 312c4762a1bSJed Brown PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 313c4762a1bSJed Brown PetscReal rhx=mx+1, rhy=my+1; 314c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 315c4762a1bSJed Brown PetscReal hl,hr,ht,hb,hc,htl,hbr; 316c4762a1bSJed Brown const PetscScalar *x; 317c4762a1bSJed Brown PetscReal v[7]; 318c4762a1bSJed Brown 319362febeeSStefano Zampini PetscFunctionBeginUser; 320c4762a1bSJed Brown /* Get pointers to vector data */ 3219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* Initialize matrix entries to zero */ 3249566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hessian)); 325c4762a1bSJed Brown 326c4762a1bSJed Brown /* Set various matrix options */ 3279566063dSJacob Faibussowitsch PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */ 330c4762a1bSJed Brown for (i=0; i< mx; i++) { 331c4762a1bSJed Brown for (j=0; j<my; j++) { 332c4762a1bSJed Brown 333c4762a1bSJed Brown row=(j)*mx + (i); 334c4762a1bSJed Brown 335c4762a1bSJed Brown xc = x[row]; 336c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* Left side */ 339c4762a1bSJed Brown if (i==0) { 340c4762a1bSJed Brown xl= user->left[j+1]; 341c4762a1bSJed Brown xlt = user->left[j+2]; 342c4762a1bSJed Brown } else { 343c4762a1bSJed Brown xl = x[row-1]; 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 346c4762a1bSJed Brown if (j==0) { 347c4762a1bSJed Brown xb=user->bottom[i+1]; 348c4762a1bSJed Brown xrb = user->bottom[i+2]; 349c4762a1bSJed Brown } else { 350c4762a1bSJed Brown xb = x[row-mx]; 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 353c4762a1bSJed Brown if (i+1 == mx) { 354c4762a1bSJed Brown xr=user->right[j+1]; 355c4762a1bSJed Brown xrb = user->right[j]; 356c4762a1bSJed Brown } else { 357c4762a1bSJed Brown xr = x[row+1]; 358c4762a1bSJed Brown } 359c4762a1bSJed Brown 360c4762a1bSJed Brown if (j+1==my) { 361c4762a1bSJed Brown xt=user->top[i+1]; 362c4762a1bSJed Brown xlt = user->top[i]; 363c4762a1bSJed Brown }else { 364c4762a1bSJed Brown xt = x[row+mx]; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown if (i>0 && j+1<my) { 368c4762a1bSJed Brown xlt = x[row-1+mx]; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown if (j>0 && i+1<mx) { 371c4762a1bSJed Brown xrb = x[row+1-mx]; 372c4762a1bSJed Brown } 373c4762a1bSJed Brown 374c4762a1bSJed Brown d1 = (xc-xl)*rhx; 375c4762a1bSJed Brown d2 = (xc-xr)*rhx; 376c4762a1bSJed Brown d3 = (xc-xt)*rhy; 377c4762a1bSJed Brown d4 = (xc-xb)*rhy; 378c4762a1bSJed Brown d5 = (xrb-xr)*rhy; 379c4762a1bSJed Brown d6 = (xrb-xb)*rhx; 380c4762a1bSJed Brown d7 = (xlt-xl)*rhy; 381c4762a1bSJed Brown d8 = (xlt-xt)*rhx; 382c4762a1bSJed Brown 383c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 384c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 385c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 386c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 387c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 388c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 389c4762a1bSJed Brown 390c4762a1bSJed Brown hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 391c4762a1bSJed Brown hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 392c4762a1bSJed Brown ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 393c4762a1bSJed Brown hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 394c4762a1bSJed Brown 395c4762a1bSJed Brown hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 396c4762a1bSJed Brown htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 397c4762a1bSJed Brown 398c4762a1bSJed 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) + 399c4762a1bSJed 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); 400c4762a1bSJed Brown 401c4762a1bSJed Brown hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5; hc*=0.5; 402c4762a1bSJed Brown 403c4762a1bSJed Brown k=0; 404c4762a1bSJed Brown if (j>0) { 405c4762a1bSJed Brown v[k]=hb; col[k]=row - mx; k++; 406c4762a1bSJed Brown } 407c4762a1bSJed Brown 408c4762a1bSJed Brown if (j>0 && i < mx -1) { 409c4762a1bSJed Brown v[k]=hbr; col[k]=row - mx+1; k++; 410c4762a1bSJed Brown } 411c4762a1bSJed Brown 412c4762a1bSJed Brown if (i>0) { 413c4762a1bSJed Brown v[k]= hl; col[k]=row - 1; k++; 414c4762a1bSJed Brown } 415c4762a1bSJed Brown 416c4762a1bSJed Brown v[k]= hc; col[k]=row; k++; 417c4762a1bSJed Brown 418c4762a1bSJed Brown if (i < mx-1) { 419c4762a1bSJed Brown v[k]= hr; col[k]=row+1; k++; 420c4762a1bSJed Brown } 421c4762a1bSJed Brown 422c4762a1bSJed Brown if (i>0 && j < my-1) { 423c4762a1bSJed Brown v[k]= htl; col[k] = row+mx-1; k++; 424c4762a1bSJed Brown } 425c4762a1bSJed Brown 426c4762a1bSJed Brown if (j < my-1) { 427c4762a1bSJed Brown v[k]= ht; col[k] = row+mx; k++; 428c4762a1bSJed Brown } 429c4762a1bSJed Brown 430c4762a1bSJed Brown /* 431c4762a1bSJed Brown Set matrix values using local numbering, which was defined 432c4762a1bSJed Brown earlier, in the main routine. 433c4762a1bSJed Brown */ 4349566063dSJacob Faibussowitsch PetscCall(MatSetValues(Hessian,1,&row,k,col,v,INSERT_VALUES)); 435c4762a1bSJed Brown } 436c4762a1bSJed Brown } 437c4762a1bSJed Brown 438c4762a1bSJed Brown /* Restore vectors */ 4399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 440c4762a1bSJed Brown 441c4762a1bSJed Brown /* Assemble the matrix */ 4429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY)); 4439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY)); 444c4762a1bSJed Brown 4459566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199.0*mx*my)); 446c4762a1bSJed Brown PetscFunctionReturn(0); 447c4762a1bSJed Brown } 448c4762a1bSJed Brown 449c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 450c4762a1bSJed Brown /* 451c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 452c4762a1bSJed Brown the region. 453c4762a1bSJed Brown 454c4762a1bSJed Brown Input Parameter: 455c4762a1bSJed Brown . user - user-defined application context 456c4762a1bSJed Brown 457c4762a1bSJed Brown Output Parameter: 458c4762a1bSJed Brown . user - user-defined application context 459c4762a1bSJed Brown */ 460c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 461c4762a1bSJed Brown { 462c4762a1bSJed Brown PetscInt i,j,k,limit=0; 463c4762a1bSJed Brown PetscInt maxits=5; 464c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 465c4762a1bSJed Brown PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 466c4762a1bSJed Brown PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 467c4762a1bSJed Brown PetscReal fnorm,det,hx,hy,xt=0,yt=0; 468c4762a1bSJed Brown PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 469c4762a1bSJed Brown PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 470c4762a1bSJed Brown PetscReal *boundary; 471c4762a1bSJed Brown 472c4762a1bSJed Brown PetscFunctionBeginUser; 473c4762a1bSJed Brown bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2; 474c4762a1bSJed Brown 4759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize,&user->bottom)); 4769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize,&user->top)); 4779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize,&user->left)); 4789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize,&user->right)); 479c4762a1bSJed Brown 480c4762a1bSJed Brown hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 481c4762a1bSJed Brown 482c4762a1bSJed Brown for (j=0; j<4; j++) { 483c4762a1bSJed Brown if (j==0) { 484c4762a1bSJed Brown yt=b; 485c4762a1bSJed Brown xt=l; 486c4762a1bSJed Brown limit=bsize; 487c4762a1bSJed Brown boundary=user->bottom; 488c4762a1bSJed Brown } else if (j==1) { 489c4762a1bSJed Brown yt=t; 490c4762a1bSJed Brown xt=l; 491c4762a1bSJed Brown limit=tsize; 492c4762a1bSJed Brown boundary=user->top; 493c4762a1bSJed Brown } else if (j==2) { 494c4762a1bSJed Brown yt=b; 495c4762a1bSJed Brown xt=l; 496c4762a1bSJed Brown limit=lsize; 497c4762a1bSJed Brown boundary=user->left; 498c4762a1bSJed Brown } else { /* if (j==3) */ 499c4762a1bSJed Brown yt=b; 500c4762a1bSJed Brown xt=r; 501c4762a1bSJed Brown limit=rsize; 502c4762a1bSJed Brown boundary=user->right; 503c4762a1bSJed Brown } 504c4762a1bSJed Brown 505c4762a1bSJed Brown for (i=0; i<limit; i++) { 506c4762a1bSJed Brown u1=xt; 507c4762a1bSJed Brown u2=-yt; 508c4762a1bSJed Brown for (k=0; k<maxits; k++) { 509c4762a1bSJed Brown nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 510c4762a1bSJed Brown nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 511c4762a1bSJed Brown fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2); 512c4762a1bSJed Brown if (fnorm <= tol) break; 513c4762a1bSJed Brown njac11=one+u2*u2-u1*u1; 514c4762a1bSJed Brown njac12=two*u1*u2; 515c4762a1bSJed Brown njac21=-two*u1*u2; 516c4762a1bSJed Brown njac22=-one - u1*u1 + u2*u2; 517c4762a1bSJed Brown det = njac11*njac22-njac21*njac12; 518c4762a1bSJed Brown u1 = u1-(njac22*nf1-njac12*nf2)/det; 519c4762a1bSJed Brown u2 = u2-(njac11*nf2-njac21*nf1)/det; 520c4762a1bSJed Brown } 521c4762a1bSJed Brown 522c4762a1bSJed Brown boundary[i]=u1*u1-u2*u2; 523c4762a1bSJed Brown if (j==0 || j==1) { 524c4762a1bSJed Brown xt=xt+hx; 525c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 526c4762a1bSJed Brown yt=yt+hy; 527c4762a1bSJed Brown } 528c4762a1bSJed Brown } 529c4762a1bSJed Brown } 530c4762a1bSJed Brown PetscFunctionReturn(0); 531c4762a1bSJed Brown } 532c4762a1bSJed Brown 533c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 534c4762a1bSJed Brown /* 535c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 536c4762a1bSJed Brown 537c4762a1bSJed Brown Input Parameters: 538c4762a1bSJed Brown . user - user-defined application context 539c4762a1bSJed Brown . X - vector for initial guess 540c4762a1bSJed Brown 541c4762a1bSJed Brown Output Parameters: 542c4762a1bSJed Brown . X - newly computed initial guess 543c4762a1bSJed Brown */ 544c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 545c4762a1bSJed Brown { 546c4762a1bSJed Brown PetscInt start=-1,i,j; 547c4762a1bSJed Brown PetscReal zero=0.0; 548c4762a1bSJed Brown PetscBool flg; 549c4762a1bSJed Brown 550c4762a1bSJed Brown PetscFunctionBeginUser; 5519566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero)); 5529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg)); 553c4762a1bSJed Brown 554c4762a1bSJed Brown if (flg && start==0) { /* The zero vector is reasonable */ 5559566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero)); 556c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 557c4762a1bSJed Brown PetscInt row; 558c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 559c4762a1bSJed Brown PetscReal *x; 560c4762a1bSJed Brown 561c4762a1bSJed Brown /* Get pointers to vector data */ 5629566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 563c4762a1bSJed Brown /* Perform local computations */ 564c4762a1bSJed Brown for (j=0; j<my; j++) { 565c4762a1bSJed Brown for (i=0; i< mx; i++) { 566c4762a1bSJed Brown row=(j)*mx + (i); 567c4762a1bSJed 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; 568c4762a1bSJed Brown } 569c4762a1bSJed Brown } 570c4762a1bSJed Brown /* Restore vectors */ 5719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x)); 572c4762a1bSJed Brown } 573c4762a1bSJed Brown PetscFunctionReturn(0); 574c4762a1bSJed Brown } 575c4762a1bSJed Brown 576c4762a1bSJed Brown /*TEST 577c4762a1bSJed Brown 578c4762a1bSJed Brown build: 579c4762a1bSJed Brown requires: !complex 580c4762a1bSJed Brown 581c4762a1bSJed Brown test: 582c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4 583c4762a1bSJed Brown requires: !single 584c4762a1bSJed Brown 585c4762a1bSJed Brown test: 586c4762a1bSJed Brown suffix: 2 587c4762a1bSJed Brown args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 588c4762a1bSJed Brown requires: !single 589c4762a1bSJed Brown 590c4762a1bSJed Brown test: 591c4762a1bSJed Brown suffix: 3 592c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 593c4762a1bSJed Brown requires: !single 594c4762a1bSJed Brown 595c4762a1bSJed Brown test: 596c4762a1bSJed Brown suffix: 4 597c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 598c4762a1bSJed Brown requires: !single 599c4762a1bSJed Brown 600c4762a1bSJed Brown test: 6010f0abf79SStefano Zampini suffix: 4_ew 6020f0abf79SStefano Zampini args: -tao_smonitor -tao_type bntr -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4 6030f0abf79SStefano Zampini requires: !single 6040f0abf79SStefano Zampini 6050f0abf79SStefano Zampini test: 606c4762a1bSJed Brown suffix: 5 607c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 608c4762a1bSJed Brown requires: !single 609c4762a1bSJed Brown 610c4762a1bSJed Brown test: 6110f0abf79SStefano Zampini suffix: 5_ew 6120f0abf79SStefano Zampini args: -tao_smonitor -tao_type bntl -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4 6130f0abf79SStefano Zampini requires: !single 6140f0abf79SStefano Zampini 6150f0abf79SStefano Zampini test: 616c4762a1bSJed Brown suffix: 6 617c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4 618c4762a1bSJed Brown requires: !single 619c4762a1bSJed Brown 620c4762a1bSJed Brown test: 6210f0abf79SStefano Zampini suffix: 6_ew 6220f0abf79SStefano 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 6230f0abf79SStefano Zampini requires: !single 6240f0abf79SStefano Zampini 6250f0abf79SStefano Zampini test: 626c4762a1bSJed Brown suffix: 7 627c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 628c4762a1bSJed Brown requires: !single 629c4762a1bSJed Brown 630c4762a1bSJed Brown test: 631c4762a1bSJed Brown suffix: 8 632c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 633c4762a1bSJed Brown requires: !single 634c4762a1bSJed Brown 635c4762a1bSJed Brown test: 636c4762a1bSJed Brown suffix: 9 637c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 638c4762a1bSJed Brown requires: !single 639c4762a1bSJed Brown 64034ad9904SAlp Dener test: 64134ad9904SAlp Dener suffix: 10 64234ad9904SAlp 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 64334ad9904SAlp Dener 64434ad9904SAlp Dener test: 64534ad9904SAlp Dener suffix: 11 64634ad9904SAlp Dener args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian 64734ad9904SAlp Dener requires: !single 64834ad9904SAlp Dener 64934ad9904SAlp Dener test: 65034ad9904SAlp Dener suffix: 12 65134ad9904SAlp Dener args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian 65234ad9904SAlp Dener requires: !single 65334ad9904SAlp Dener 654c4762a1bSJed Brown TEST*/ 655