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