1c4762a1bSJed Brown #include <petsctao.h> 2c4762a1bSJed Brown 3c4762a1bSJed Brown static char help[] = 4c4762a1bSJed Brown "This example demonstrates use of the TAO package to\n\ 5c4762a1bSJed Brown solve an unconstrained system of equations. This example is based on a\n\ 6c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and\n\ 7c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\ 8c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\ 9c4762a1bSJed Brown This application solves this problem using complimentarity -- We are actually\n\ 10c4762a1bSJed Brown solving the system (grad f)_i >= 0, if x_i == l_i \n\ 11c4762a1bSJed Brown (grad f)_i = 0, if l_i < x_i < u_i \n\ 12c4762a1bSJed Brown (grad f)_i <= 0, if x_i == u_i \n\ 13c4762a1bSJed Brown where f is the function to be minimized. \n\ 14c4762a1bSJed Brown \n\ 15c4762a1bSJed Brown The command line options are:\n\ 16c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 17c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 18c4762a1bSJed Brown -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n"; 19c4762a1bSJed Brown 20c4762a1bSJed Brown /*T 21c4762a1bSJed Brown Concepts: TAO^Solving a complementarity problem 22c4762a1bSJed Brown Routines: TaoCreate(); TaoDestroy(); 23c4762a1bSJed Brown 24c4762a1bSJed Brown Processors: 1 25c4762a1bSJed Brown T*/ 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown User-defined application context - contains data needed by the 29c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient(), 30c4762a1bSJed Brown FormHessian(). 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown typedef struct { 33c4762a1bSJed Brown PetscInt mx, my; 34c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right; 35c4762a1bSJed Brown } AppCtx; 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 38c4762a1bSJed Brown 39c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *); 40c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec); 41c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *); 42c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *); 43c4762a1bSJed Brown 44c4762a1bSJed Brown int main(int argc, char **argv) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown Vec x; /* solution vector */ 47c4762a1bSJed Brown Vec c; /* Constraints function vector */ 48c4762a1bSJed Brown Vec xl,xu; /* Bounds on the variables */ 49c4762a1bSJed Brown PetscBool flg; /* A return variable when checking for user options */ 50c4762a1bSJed Brown Tao tao; /* TAO solver context */ 51c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 52c4762a1bSJed Brown PetscInt N; /* Number of elements in vector */ 53c4762a1bSJed Brown PetscScalar lb = PETSC_NINFINITY; /* lower bound constant */ 54c4762a1bSJed Brown PetscScalar ub = PETSC_INFINITY; /* upper bound constant */ 55c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Initialize PETSc, TAO */ 58*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, (char *)0, help)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Specify default dimension of the problem */ 61c4762a1bSJed Brown user.mx = 4; user.my = 4; 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-mx", &user.mx, &flg)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-my", &user.my, &flg)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Calculate any derived values from parameters */ 68c4762a1bSJed Brown N = user.mx*user.my; 69c4762a1bSJed Brown 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n")); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"mx:%D, my:%D\n", user.mx,user.my)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Create appropriate vectors and matrices */ 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(MPI_COMM_SELF, N, &x)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &c)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* The TAO code begins here */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOSSILS)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Set data structure */ 855f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao, x)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Set routines for constraints function and Jacobian evaluation */ 885f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Set the variable bounds */ 925f80ce2aSJacob Faibussowitsch CHKERRQ(MSA_BoundaryConditions(&user)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Set initial solution guess */ 955f80ce2aSJacob Faibussowitsch CHKERRQ(MSA_InitialPoint(&user, x)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Set Bounds on variables */ 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &xl)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &xu)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xl, lb)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(xu, ub)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetVariableBounds(tao,xl,xu)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Check for any tao command line options */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Solve the application */ 1085f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Free Tao data structures */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Free PETSc data structures */ 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xl)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xu)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&c)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Free user-created data structures */ 1215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.bottom)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.top)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.left)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user.right)); 125c4762a1bSJed Brown 126*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 127*b122ec5aSJacob Faibussowitsch return 0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* FormConstraints - Evaluates gradient of f. 133c4762a1bSJed Brown 134c4762a1bSJed Brown Input Parameters: 135c4762a1bSJed Brown . tao - the TAO_APPLICATION context 136c4762a1bSJed Brown . X - input vector 137c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetConstraintsRoutine() 138c4762a1bSJed Brown 139c4762a1bSJed Brown Output Parameters: 140c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 141c4762a1bSJed Brown */ 142c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr) 143c4762a1bSJed Brown { 144c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 145c4762a1bSJed Brown PetscInt i,j,row; 146c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 147c4762a1bSJed Brown PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 148c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 149c4762a1bSJed Brown PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 150c4762a1bSJed Brown PetscScalar zero=0.0; 151c4762a1bSJed Brown PetscScalar *g, *x; 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscFunctionBegin; 154c4762a1bSJed Brown /* Initialize vector to zero */ 1555f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(G, zero)); 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* Get pointers to vector data */ 1585f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(X, &x)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(G, &g)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */ 162c4762a1bSJed Brown for (j=0; j<my; j++) { 163c4762a1bSJed Brown for (i=0; i< mx; i++) { 164c4762a1bSJed Brown row= j*mx + i; 165c4762a1bSJed Brown 166c4762a1bSJed Brown xc = x[row]; 167c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 168c4762a1bSJed Brown 169c4762a1bSJed Brown if (i==0) { /* left side */ 170c4762a1bSJed Brown xl= user->left[j+1]; 171c4762a1bSJed Brown xlt = user->left[j+2]; 172c4762a1bSJed Brown } else { 173c4762a1bSJed Brown xl = x[row-1]; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown if (j==0) { /* bottom side */ 177c4762a1bSJed Brown xb=user->bottom[i+1]; 178c4762a1bSJed Brown xrb = user->bottom[i+2]; 179c4762a1bSJed Brown } else { 180c4762a1bSJed Brown xb = x[row-mx]; 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown if (i+1 == mx) { /* right side */ 184c4762a1bSJed Brown xr=user->right[j+1]; 185c4762a1bSJed Brown xrb = user->right[j]; 186c4762a1bSJed Brown } else { 187c4762a1bSJed Brown xr = x[row+1]; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown if (j+1==0+my) { /* top side */ 191c4762a1bSJed Brown xt=user->top[i+1]; 192c4762a1bSJed Brown xlt = user->top[i]; 193c4762a1bSJed Brown }else { 194c4762a1bSJed Brown xt = x[row+mx]; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown if (i>0 && j+1<my) { 198c4762a1bSJed Brown xlt = x[row-1+mx]; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown if (j>0 && i+1<mx) { 201c4762a1bSJed Brown xrb = x[row+1-mx]; 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204c4762a1bSJed Brown d1 = (xc-xl); 205c4762a1bSJed Brown d2 = (xc-xr); 206c4762a1bSJed Brown d3 = (xc-xt); 207c4762a1bSJed Brown d4 = (xc-xb); 208c4762a1bSJed Brown d5 = (xr-xrb); 209c4762a1bSJed Brown d6 = (xrb-xb); 210c4762a1bSJed Brown d7 = (xlt-xl); 211c4762a1bSJed Brown d8 = (xt-xlt); 212c4762a1bSJed Brown 213c4762a1bSJed Brown df1dxc = d1*hydhx; 214c4762a1bSJed Brown df2dxc = (d1*hydhx + d4*hxdhy); 215c4762a1bSJed Brown df3dxc = d3*hxdhy; 216c4762a1bSJed Brown df4dxc = (d2*hydhx + d3*hxdhy); 217c4762a1bSJed Brown df5dxc = d2*hydhx; 218c4762a1bSJed Brown df6dxc = d4*hxdhy; 219c4762a1bSJed Brown 220c4762a1bSJed Brown d1 /= hx; 221c4762a1bSJed Brown d2 /= hx; 222c4762a1bSJed Brown d3 /= hy; 223c4762a1bSJed Brown d4 /= hy; 224c4762a1bSJed Brown d5 /= hy; 225c4762a1bSJed Brown d6 /= hx; 226c4762a1bSJed Brown d7 /= hy; 227c4762a1bSJed Brown d8 /= hx; 228c4762a1bSJed Brown 229c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 230c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 231c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 232c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 233c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 234c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 235c4762a1bSJed Brown 236c4762a1bSJed Brown df1dxc /= f1; 237c4762a1bSJed Brown df2dxc /= f2; 238c4762a1bSJed Brown df3dxc /= f3; 239c4762a1bSJed Brown df4dxc /= f4; 240c4762a1bSJed Brown df5dxc /= f5; 241c4762a1bSJed Brown df6dxc /= f6; 242c4762a1bSJed Brown 243c4762a1bSJed Brown g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0; 244c4762a1bSJed Brown } 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* Restore vectors */ 2485f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(X, &x)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(G, &g)); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(67*mx*my)); 251c4762a1bSJed Brown PetscFunctionReturn(0); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 255c4762a1bSJed Brown /* 256c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 257c4762a1bSJed Brown 258c4762a1bSJed Brown Input Parameters: 259c4762a1bSJed Brown . tao - the TAO_APPLICATION context 260c4762a1bSJed Brown . X - input vector 261c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetJacobian() 262c4762a1bSJed Brown 263c4762a1bSJed Brown Output Parameters: 264c4762a1bSJed Brown . tH - Jacobian matrix 265c4762a1bSJed Brown 266c4762a1bSJed Brown */ 267c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao tao, Vec X, Mat H, Mat tHPre, void *ptr) 268c4762a1bSJed Brown { 269c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 270c4762a1bSJed Brown PetscInt i,j,k,row; 271c4762a1bSJed Brown PetscInt mx=user->mx, my=user->my; 272c4762a1bSJed Brown PetscInt col[7]; 273c4762a1bSJed Brown PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 274c4762a1bSJed Brown PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 275c4762a1bSJed Brown PetscReal hl,hr,ht,hb,hc,htl,hbr; 276c4762a1bSJed Brown const PetscScalar *x; 277c4762a1bSJed Brown PetscScalar v[7]; 278c4762a1bSJed Brown PetscBool assembled; 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* Set various matrix options */ 281362febeeSStefano Zampini PetscFunctionBegin; 2825f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(H,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssembled(H,&assembled)); 2845f80ce2aSJacob Faibussowitsch if (assembled) CHKERRQ(MatZeroEntries(H)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* Get pointers to vector data */ 2875f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X, &x)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* Compute Jacobian over the locally owned part of the mesh */ 290c4762a1bSJed Brown for (i=0; i< mx; i++) { 291c4762a1bSJed Brown for (j=0; j<my; j++) { 292c4762a1bSJed Brown row= j*mx + i; 293c4762a1bSJed Brown 294c4762a1bSJed Brown xc = x[row]; 295c4762a1bSJed Brown xlt=xrb=xl=xr=xb=xt=xc; 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* Left side */ 298c4762a1bSJed Brown if (i==0) { 299c4762a1bSJed Brown xl = user->left[j+1]; 300c4762a1bSJed Brown xlt = user->left[j+2]; 301c4762a1bSJed Brown } else { 302c4762a1bSJed Brown xl = x[row-1]; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305c4762a1bSJed Brown if (j==0) { 306c4762a1bSJed Brown xb = user->bottom[i+1]; 307c4762a1bSJed Brown xrb = user->bottom[i+2]; 308c4762a1bSJed Brown } else { 309c4762a1bSJed Brown xb = x[row-mx]; 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown if (i+1 == mx) { 313c4762a1bSJed Brown xr = user->right[j+1]; 314c4762a1bSJed Brown xrb = user->right[j]; 315c4762a1bSJed Brown } else { 316c4762a1bSJed Brown xr = x[row+1]; 317c4762a1bSJed Brown } 318c4762a1bSJed Brown 319c4762a1bSJed Brown if (j+1==my) { 320c4762a1bSJed Brown xt = user->top[i+1]; 321c4762a1bSJed Brown xlt = user->top[i]; 322c4762a1bSJed Brown }else { 323c4762a1bSJed Brown xt = x[row+mx]; 324c4762a1bSJed Brown } 325c4762a1bSJed Brown 326c4762a1bSJed Brown if (i>0 && j+1<my) { 327c4762a1bSJed Brown xlt = x[row-1+mx]; 328c4762a1bSJed Brown } 329c4762a1bSJed Brown if (j>0 && i+1<mx) { 330c4762a1bSJed Brown xrb = x[row+1-mx]; 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown d1 = (xc-xl)/hx; 334c4762a1bSJed Brown d2 = (xc-xr)/hx; 335c4762a1bSJed Brown d3 = (xc-xt)/hy; 336c4762a1bSJed Brown d4 = (xc-xb)/hy; 337c4762a1bSJed Brown d5 = (xrb-xr)/hy; 338c4762a1bSJed Brown d6 = (xrb-xb)/hx; 339c4762a1bSJed Brown d7 = (xlt-xl)/hy; 340c4762a1bSJed Brown d8 = (xlt-xt)/hx; 341c4762a1bSJed Brown 342c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 343c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 344c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 345c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 346c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 347c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 348c4762a1bSJed Brown 349c4762a1bSJed Brown hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 350c4762a1bSJed Brown hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 351c4762a1bSJed Brown ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 352c4762a1bSJed Brown hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 353c4762a1bSJed Brown 354c4762a1bSJed Brown hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 355c4762a1bSJed Brown htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 356c4762a1bSJed Brown 357c4762a1bSJed 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) + 358c4762a1bSJed 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); 359c4762a1bSJed Brown 360c4762a1bSJed Brown hl/=2.0; hr/=2.0; ht/=2.0; hb/=2.0; hbr/=2.0; htl/=2.0; hc/=2.0; 361c4762a1bSJed Brown 362c4762a1bSJed Brown k=0; 363c4762a1bSJed Brown if (j>0) { 364c4762a1bSJed Brown v[k]=hb; col[k]=row - mx; k++; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown if (j>0 && i < mx -1) { 368c4762a1bSJed Brown v[k]=hbr; col[k]=row - mx+1; k++; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown 371c4762a1bSJed Brown if (i>0) { 372c4762a1bSJed Brown v[k]= hl; col[k]=row - 1; k++; 373c4762a1bSJed Brown } 374c4762a1bSJed Brown 375c4762a1bSJed Brown v[k]= hc; col[k]=row; k++; 376c4762a1bSJed Brown 377c4762a1bSJed Brown if (i < mx-1) { 378c4762a1bSJed Brown v[k]= hr; col[k]=row+1; k++; 379c4762a1bSJed Brown } 380c4762a1bSJed Brown 381c4762a1bSJed Brown if (i>0 && j < my-1) { 382c4762a1bSJed Brown v[k]= htl; col[k] = row+mx-1; k++; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown if (j < my-1) { 386c4762a1bSJed Brown v[k]= ht; col[k] = row+mx; k++; 387c4762a1bSJed Brown } 388c4762a1bSJed Brown 389c4762a1bSJed Brown /* 390c4762a1bSJed Brown Set matrix values using local numbering, which was defined 391c4762a1bSJed Brown earlier, in the main routine. 392c4762a1bSJed Brown */ 3935f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(H,1,&row,k,col,v,INSERT_VALUES)); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown } 396c4762a1bSJed Brown 397c4762a1bSJed Brown /* Restore vectors */ 3985f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* Assemble the matrix */ 4015f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 4025f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 4035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(199*mx*my)); 404c4762a1bSJed Brown PetscFunctionReturn(0); 405c4762a1bSJed Brown } 406c4762a1bSJed Brown 407c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 408c4762a1bSJed Brown /* 409c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 410c4762a1bSJed Brown the region. 411c4762a1bSJed Brown 412c4762a1bSJed Brown Input Parameter: 413c4762a1bSJed Brown . user - user-defined application context 414c4762a1bSJed Brown 415c4762a1bSJed Brown Output Parameter: 416c4762a1bSJed Brown . user - user-defined application context 417c4762a1bSJed Brown */ 418c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 419c4762a1bSJed Brown { 420c4762a1bSJed Brown PetscInt i,j,k,limit=0,maxits=5; 421c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 422c4762a1bSJed Brown PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 423c4762a1bSJed Brown PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 424c4762a1bSJed Brown PetscReal fnorm,det,hx,hy,xt=0,yt=0; 425c4762a1bSJed Brown PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 426c4762a1bSJed Brown PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 427c4762a1bSJed Brown PetscReal *boundary; 428c4762a1bSJed Brown 429c4762a1bSJed Brown PetscFunctionBegin; 430c4762a1bSJed Brown bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2; 431c4762a1bSJed Brown 4325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bsize, &user->bottom)); 4335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(tsize, &user->top)); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(lsize, &user->left)); 4355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(rsize, &user->right)); 436c4762a1bSJed Brown 437c4762a1bSJed Brown hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 438c4762a1bSJed Brown 439c4762a1bSJed Brown for (j=0; j<4; j++) { 440c4762a1bSJed Brown if (j==0) { 441c4762a1bSJed Brown yt=b; 442c4762a1bSJed Brown xt=l; 443c4762a1bSJed Brown limit=bsize; 444c4762a1bSJed Brown boundary=user->bottom; 445c4762a1bSJed Brown } else if (j==1) { 446c4762a1bSJed Brown yt=t; 447c4762a1bSJed Brown xt=l; 448c4762a1bSJed Brown limit=tsize; 449c4762a1bSJed Brown boundary=user->top; 450c4762a1bSJed Brown } else if (j==2) { 451c4762a1bSJed Brown yt=b; 452c4762a1bSJed Brown xt=l; 453c4762a1bSJed Brown limit=lsize; 454c4762a1bSJed Brown boundary=user->left; 455c4762a1bSJed Brown } else { /* if (j==3) */ 456c4762a1bSJed Brown yt=b; 457c4762a1bSJed Brown xt=r; 458c4762a1bSJed Brown limit=rsize; 459c4762a1bSJed Brown boundary=user->right; 460c4762a1bSJed Brown } 461c4762a1bSJed Brown 462c4762a1bSJed Brown for (i=0; i<limit; i++) { 463c4762a1bSJed Brown u1=xt; 464c4762a1bSJed Brown u2=-yt; 465c4762a1bSJed Brown for (k=0; k<maxits; k++) { 466c4762a1bSJed Brown nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 467c4762a1bSJed Brown nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 468c4762a1bSJed Brown fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2); 469c4762a1bSJed Brown if (fnorm <= tol) break; 470c4762a1bSJed Brown njac11=one+u2*u2-u1*u1; 471c4762a1bSJed Brown njac12=two*u1*u2; 472c4762a1bSJed Brown njac21=-two*u1*u2; 473c4762a1bSJed Brown njac22=-one - u1*u1 + u2*u2; 474c4762a1bSJed Brown det = njac11*njac22-njac21*njac12; 475c4762a1bSJed Brown u1 = u1-(njac22*nf1-njac12*nf2)/det; 476c4762a1bSJed Brown u2 = u2-(njac11*nf2-njac21*nf1)/det; 477c4762a1bSJed Brown } 478c4762a1bSJed Brown 479c4762a1bSJed Brown boundary[i]=u1*u1-u2*u2; 480c4762a1bSJed Brown if (j==0 || j==1) { 481c4762a1bSJed Brown xt=xt+hx; 482c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 483c4762a1bSJed Brown yt=yt+hy; 484c4762a1bSJed Brown } 485c4762a1bSJed Brown } 486c4762a1bSJed Brown } 487c4762a1bSJed Brown PetscFunctionReturn(0); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown 490c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 491c4762a1bSJed Brown /* 492c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 493c4762a1bSJed Brown 494c4762a1bSJed Brown Input Parameters: 495c4762a1bSJed Brown . user - user-defined application context 496c4762a1bSJed Brown . X - vector for initial guess 497c4762a1bSJed Brown 498c4762a1bSJed Brown Output Parameters: 499c4762a1bSJed Brown . X - newly computed initial guess 500c4762a1bSJed Brown */ 501c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 502c4762a1bSJed Brown { 503c4762a1bSJed Brown PetscInt start=-1,i,j; 504c4762a1bSJed Brown PetscScalar zero=0.0; 505c4762a1bSJed Brown PetscBool flg; 506c4762a1bSJed Brown 507c4762a1bSJed Brown PetscFunctionBegin; 5085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg)); 509c4762a1bSJed Brown 510c4762a1bSJed Brown if (flg && start==0) { /* The zero vector is reasonable */ 5115f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X, zero)); 512c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 513c4762a1bSJed Brown PetscInt row; 514c4762a1bSJed Brown PetscInt mx=user->mx,my=user->my; 515c4762a1bSJed Brown PetscScalar *x; 516c4762a1bSJed Brown 517c4762a1bSJed Brown /* Get pointers to vector data */ 5185f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(X,&x)); 519c4762a1bSJed Brown 520c4762a1bSJed Brown /* Perform local computations */ 521c4762a1bSJed Brown for (j=0; j<my; j++) { 522c4762a1bSJed Brown for (i=0; i< mx; i++) { 523c4762a1bSJed Brown row=(j)*mx + (i); 524c4762a1bSJed 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; 525c4762a1bSJed Brown } 526c4762a1bSJed Brown } 527c4762a1bSJed Brown 528c4762a1bSJed Brown /* Restore vectors */ 5295f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(X,&x)); 530c4762a1bSJed Brown } 531c4762a1bSJed Brown PetscFunctionReturn(0); 532c4762a1bSJed Brown } 533c4762a1bSJed Brown 534c4762a1bSJed Brown /*TEST 535c4762a1bSJed Brown 536c4762a1bSJed Brown build: 537c4762a1bSJed Brown requires: !complex 538c4762a1bSJed Brown 539c4762a1bSJed Brown test: 540c4762a1bSJed Brown args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5 541c4762a1bSJed Brown requires: !single 542c4762a1bSJed Brown 543c4762a1bSJed Brown test: 544c4762a1bSJed Brown suffix: 2 545c4762a1bSJed Brown args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5 546c4762a1bSJed Brown 547c4762a1bSJed Brown TEST*/ 548