1c4762a1bSJed Brown #include <petsctao.h> 2c4762a1bSJed Brown 39371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to\n\ 4c4762a1bSJed Brown solve an unconstrained system of equations. This example is based on a\n\ 5c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and\n\ 6c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\ 7c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\ 8c4762a1bSJed Brown This application solves this problem using complimentarity -- We are actually\n\ 9c4762a1bSJed Brown solving the system (grad f)_i >= 0, if x_i == l_i \n\ 10c4762a1bSJed Brown (grad f)_i = 0, if l_i < x_i < u_i \n\ 11c4762a1bSJed Brown (grad f)_i <= 0, if x_i == u_i \n\ 12c4762a1bSJed Brown where f is the function to be minimized. \n\ 13c4762a1bSJed Brown \n\ 14c4762a1bSJed Brown The command line options are:\n\ 15c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 16c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 17c4762a1bSJed Brown -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n"; 18c4762a1bSJed Brown 19c4762a1bSJed Brown /* 20c4762a1bSJed Brown User-defined application context - contains data needed by the 21c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient(), 22c4762a1bSJed Brown FormHessian(). 23c4762a1bSJed Brown */ 24c4762a1bSJed Brown typedef struct { 25c4762a1bSJed Brown PetscInt mx, my; 26c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right; 27c4762a1bSJed Brown } AppCtx; 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *); 32c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec); 33c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *); 34c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *); 35c4762a1bSJed Brown 369371c9d4SSatish Balay int main(int argc, char **argv) { 37c4762a1bSJed Brown Vec x; /* solution vector */ 38c4762a1bSJed Brown Vec c; /* Constraints function vector */ 39c4762a1bSJed Brown Vec xl, xu; /* Bounds on the variables */ 40c4762a1bSJed Brown PetscBool flg; /* A return variable when checking for user options */ 41c4762a1bSJed Brown Tao tao; /* TAO solver context */ 42c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 43c4762a1bSJed Brown PetscInt N; /* Number of elements in vector */ 44c4762a1bSJed Brown PetscScalar lb = PETSC_NINFINITY; /* lower bound constant */ 45c4762a1bSJed Brown PetscScalar ub = PETSC_INFINITY; /* upper bound constant */ 46c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Initialize PETSc, TAO */ 49327415f7SBarry Smith PetscFunctionBeginUser; 509566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Specify default dimension of the problem */ 539371c9d4SSatish Balay user.mx = 4; 549371c9d4SSatish Balay 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 60c4762a1bSJed Brown /* Calculate any derived values from parameters */ 61c4762a1bSJed Brown N = user.mx * user.my; 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Minimum Surface Area Problem -----\n")); 6463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx:%" PetscInt_FMT ", my:%" PetscInt_FMT "\n", user.mx, user.my)); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Create appropriate vectors and matrices */ 679566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(MPI_COMM_SELF, N, &x)); 689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &c)); 699566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* The TAO code begins here */ 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 749566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 759566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOSSILS)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Set data structure */ 789566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Set routines for constraints function and Jacobian evaluation */ 819566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user)); 829566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Set the variable bounds */ 859566063dSJacob Faibussowitsch PetscCall(MSA_BoundaryConditions(&user)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Set initial solution guess */ 889566063dSJacob Faibussowitsch PetscCall(MSA_InitialPoint(&user, x)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Set Bounds on variables */ 919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xl)); 929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xu)); 939566063dSJacob Faibussowitsch PetscCall(VecSet(xl, lb)); 949566063dSJacob Faibussowitsch PetscCall(VecSet(xu, ub)); 959566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, xl, xu)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Check for any tao command line options */ 989566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* Solve the application */ 1019566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Free Tao data structures */ 1049566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* Free PETSc data structures */ 1079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xl)); 1099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xu)); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Free user-created data structures */ 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left)); 1179566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right)); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 120b122ec5aSJacob Faibussowitsch return 0; 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* FormConstraints - Evaluates gradient of f. 126c4762a1bSJed Brown 127c4762a1bSJed Brown Input Parameters: 128c4762a1bSJed Brown . tao - the TAO_APPLICATION context 129c4762a1bSJed Brown . X - input vector 130c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetConstraintsRoutine() 131c4762a1bSJed Brown 132c4762a1bSJed Brown Output Parameters: 133c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 134c4762a1bSJed Brown */ 1359371c9d4SSatish Balay PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr) { 136c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 137c4762a1bSJed Brown PetscInt i, j, row; 138c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 139c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy; 140c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 141c4762a1bSJed Brown PetscReal df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc; 142c4762a1bSJed Brown PetscScalar zero = 0.0; 143c4762a1bSJed Brown PetscScalar *g, *x; 144c4762a1bSJed Brown 145c4762a1bSJed Brown PetscFunctionBegin; 146c4762a1bSJed Brown /* Initialize vector to zero */ 1479566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Get pointers to vector data */ 1509566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 1519566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */ 154c4762a1bSJed Brown for (j = 0; j < my; j++) { 155c4762a1bSJed Brown for (i = 0; i < mx; i++) { 156c4762a1bSJed Brown row = j * mx + i; 157c4762a1bSJed Brown 158c4762a1bSJed Brown xc = x[row]; 159c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc; 160c4762a1bSJed Brown 161c4762a1bSJed Brown if (i == 0) { /* left side */ 162c4762a1bSJed Brown xl = user->left[j + 1]; 163c4762a1bSJed Brown xlt = user->left[j + 2]; 164c4762a1bSJed Brown } else { 165c4762a1bSJed Brown xl = x[row - 1]; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown if (j == 0) { /* bottom side */ 169c4762a1bSJed Brown xb = user->bottom[i + 1]; 170c4762a1bSJed Brown xrb = user->bottom[i + 2]; 171c4762a1bSJed Brown } else { 172c4762a1bSJed Brown xb = x[row - mx]; 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown if (i + 1 == mx) { /* right side */ 176c4762a1bSJed Brown xr = user->right[j + 1]; 177c4762a1bSJed Brown xrb = user->right[j]; 178c4762a1bSJed Brown } else { 179c4762a1bSJed Brown xr = x[row + 1]; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown if (j + 1 == 0 + my) { /* top side */ 183c4762a1bSJed Brown xt = user->top[i + 1]; 184c4762a1bSJed Brown xlt = user->top[i]; 185c4762a1bSJed Brown } else { 186c4762a1bSJed Brown xt = x[row + mx]; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189*ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx]; 190*ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx]; 191c4762a1bSJed Brown 192c4762a1bSJed Brown d1 = (xc - xl); 193c4762a1bSJed Brown d2 = (xc - xr); 194c4762a1bSJed Brown d3 = (xc - xt); 195c4762a1bSJed Brown d4 = (xc - xb); 196c4762a1bSJed Brown d5 = (xr - xrb); 197c4762a1bSJed Brown d6 = (xrb - xb); 198c4762a1bSJed Brown d7 = (xlt - xl); 199c4762a1bSJed Brown d8 = (xt - xlt); 200c4762a1bSJed Brown 201c4762a1bSJed Brown df1dxc = d1 * hydhx; 202c4762a1bSJed Brown df2dxc = (d1 * hydhx + d4 * hxdhy); 203c4762a1bSJed Brown df3dxc = d3 * hxdhy; 204c4762a1bSJed Brown df4dxc = (d2 * hydhx + d3 * hxdhy); 205c4762a1bSJed Brown df5dxc = d2 * hydhx; 206c4762a1bSJed Brown df6dxc = d4 * hxdhy; 207c4762a1bSJed Brown 208c4762a1bSJed Brown d1 /= hx; 209c4762a1bSJed Brown d2 /= hx; 210c4762a1bSJed Brown d3 /= hy; 211c4762a1bSJed Brown d4 /= hy; 212c4762a1bSJed Brown d5 /= hy; 213c4762a1bSJed Brown d6 /= hx; 214c4762a1bSJed Brown d7 /= hy; 215c4762a1bSJed Brown d8 /= hx; 216c4762a1bSJed Brown 217c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7); 218c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4); 219c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8); 220c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2); 221c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5); 222c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6); 223c4762a1bSJed Brown 224c4762a1bSJed Brown df1dxc /= f1; 225c4762a1bSJed Brown df2dxc /= f2; 226c4762a1bSJed Brown df3dxc /= f3; 227c4762a1bSJed Brown df4dxc /= f4; 228c4762a1bSJed Brown df5dxc /= f5; 229c4762a1bSJed Brown df6dxc /= f6; 230c4762a1bSJed Brown 231c4762a1bSJed Brown g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown } 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* Restore vectors */ 2369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 2379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g)); 2389566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67 * mx * my)); 239c4762a1bSJed Brown PetscFunctionReturn(0); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 243c4762a1bSJed Brown /* 244c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 245c4762a1bSJed Brown 246c4762a1bSJed Brown Input Parameters: 247c4762a1bSJed Brown . tao - the TAO_APPLICATION context 248c4762a1bSJed Brown . X - input vector 249c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetJacobian() 250c4762a1bSJed Brown 251c4762a1bSJed Brown Output Parameters: 252c4762a1bSJed Brown . tH - Jacobian matrix 253c4762a1bSJed Brown 254c4762a1bSJed Brown */ 2559371c9d4SSatish Balay PetscErrorCode FormJacobian(Tao tao, Vec X, Mat H, Mat tHPre, void *ptr) { 256c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 257c4762a1bSJed Brown PetscInt i, j, k, row; 258c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 259c4762a1bSJed Brown PetscInt col[7]; 260c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy; 261c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 262c4762a1bSJed Brown PetscReal hl, hr, ht, hb, hc, htl, hbr; 263c4762a1bSJed Brown const PetscScalar *x; 264c4762a1bSJed Brown PetscScalar v[7]; 265c4762a1bSJed Brown PetscBool assembled; 266c4762a1bSJed Brown 267c4762a1bSJed Brown /* Set various matrix options */ 268362febeeSStefano Zampini PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(MatSetOption(H, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE)); 2709566063dSJacob Faibussowitsch PetscCall(MatAssembled(H, &assembled)); 2719566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(H)); 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* Get pointers to vector data */ 2749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* Compute Jacobian over the locally owned part of the mesh */ 277c4762a1bSJed Brown for (i = 0; i < mx; i++) { 278c4762a1bSJed Brown for (j = 0; j < my; j++) { 279c4762a1bSJed Brown row = j * mx + i; 280c4762a1bSJed Brown 281c4762a1bSJed Brown xc = x[row]; 282c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc; 283c4762a1bSJed Brown 284c4762a1bSJed Brown /* Left side */ 285c4762a1bSJed Brown if (i == 0) { 286c4762a1bSJed Brown xl = user->left[j + 1]; 287c4762a1bSJed Brown xlt = user->left[j + 2]; 288c4762a1bSJed Brown } else { 289c4762a1bSJed Brown xl = x[row - 1]; 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown if (j == 0) { 293c4762a1bSJed Brown xb = user->bottom[i + 1]; 294c4762a1bSJed Brown xrb = user->bottom[i + 2]; 295c4762a1bSJed Brown } else { 296c4762a1bSJed Brown xb = x[row - mx]; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown if (i + 1 == mx) { 300c4762a1bSJed Brown xr = user->right[j + 1]; 301c4762a1bSJed Brown xrb = user->right[j]; 302c4762a1bSJed Brown } else { 303c4762a1bSJed Brown xr = x[row + 1]; 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306c4762a1bSJed Brown if (j + 1 == my) { 307c4762a1bSJed Brown xt = user->top[i + 1]; 308c4762a1bSJed Brown xlt = user->top[i]; 309c4762a1bSJed Brown } else { 310c4762a1bSJed Brown xt = x[row + mx]; 311c4762a1bSJed Brown } 312c4762a1bSJed Brown 313*ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx]; 314*ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx]; 315c4762a1bSJed Brown 316c4762a1bSJed Brown d1 = (xc - xl) / hx; 317c4762a1bSJed Brown d2 = (xc - xr) / hx; 318c4762a1bSJed Brown d3 = (xc - xt) / hy; 319c4762a1bSJed Brown d4 = (xc - xb) / hy; 320c4762a1bSJed Brown d5 = (xrb - xr) / hy; 321c4762a1bSJed Brown d6 = (xrb - xb) / hx; 322c4762a1bSJed Brown d7 = (xlt - xl) / hy; 323c4762a1bSJed Brown d8 = (xlt - xt) / hx; 324c4762a1bSJed Brown 325c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7); 326c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4); 327c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8); 328c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2); 329c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5); 330c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6); 331c4762a1bSJed Brown 332c4762a1bSJed Brown hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2); 333c4762a1bSJed Brown hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4); 334c4762a1bSJed Brown ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4); 335c4762a1bSJed Brown hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2); 336c4762a1bSJed Brown 337c4762a1bSJed Brown hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6); 338c4762a1bSJed Brown htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3); 339c4762a1bSJed Brown 3409371c9d4SSatish Balay 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) + (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); 341c4762a1bSJed Brown 3429371c9d4SSatish Balay hl /= 2.0; 3439371c9d4SSatish Balay hr /= 2.0; 3449371c9d4SSatish Balay ht /= 2.0; 3459371c9d4SSatish Balay hb /= 2.0; 3469371c9d4SSatish Balay hbr /= 2.0; 3479371c9d4SSatish Balay htl /= 2.0; 3489371c9d4SSatish Balay hc /= 2.0; 349c4762a1bSJed Brown 350c4762a1bSJed Brown k = 0; 351c4762a1bSJed Brown if (j > 0) { 3529371c9d4SSatish Balay v[k] = hb; 3539371c9d4SSatish Balay col[k] = row - mx; 3549371c9d4SSatish Balay k++; 355c4762a1bSJed Brown } 356c4762a1bSJed Brown 357c4762a1bSJed Brown if (j > 0 && i < mx - 1) { 3589371c9d4SSatish Balay v[k] = hbr; 3599371c9d4SSatish Balay col[k] = row - mx + 1; 3609371c9d4SSatish Balay k++; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363c4762a1bSJed Brown if (i > 0) { 3649371c9d4SSatish Balay v[k] = hl; 3659371c9d4SSatish Balay col[k] = row - 1; 3669371c9d4SSatish Balay k++; 367c4762a1bSJed Brown } 368c4762a1bSJed Brown 3699371c9d4SSatish Balay v[k] = hc; 3709371c9d4SSatish Balay col[k] = row; 3719371c9d4SSatish Balay k++; 372c4762a1bSJed Brown 373c4762a1bSJed Brown if (i < mx - 1) { 3749371c9d4SSatish Balay v[k] = hr; 3759371c9d4SSatish Balay col[k] = row + 1; 3769371c9d4SSatish Balay k++; 377c4762a1bSJed Brown } 378c4762a1bSJed Brown 379c4762a1bSJed Brown if (i > 0 && j < my - 1) { 3809371c9d4SSatish Balay v[k] = htl; 3819371c9d4SSatish Balay col[k] = row + mx - 1; 3829371c9d4SSatish Balay k++; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown if (j < my - 1) { 3869371c9d4SSatish Balay v[k] = ht; 3879371c9d4SSatish Balay col[k] = row + mx; 3889371c9d4SSatish Balay k++; 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 391c4762a1bSJed Brown /* 392c4762a1bSJed Brown Set matrix values using local numbering, which was defined 393c4762a1bSJed Brown earlier, in the main routine. 394c4762a1bSJed Brown */ 3959566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 1, &row, k, col, v, INSERT_VALUES)); 396c4762a1bSJed Brown } 397c4762a1bSJed Brown } 398c4762a1bSJed Brown 399c4762a1bSJed Brown /* Restore vectors */ 4009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 401c4762a1bSJed Brown 402c4762a1bSJed Brown /* Assemble the matrix */ 4039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 4049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 4059566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199 * mx * my)); 406c4762a1bSJed Brown PetscFunctionReturn(0); 407c4762a1bSJed Brown } 408c4762a1bSJed Brown 409c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 410c4762a1bSJed Brown /* 411c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 412c4762a1bSJed Brown the region. 413c4762a1bSJed Brown 414c4762a1bSJed Brown Input Parameter: 415c4762a1bSJed Brown . user - user-defined application context 416c4762a1bSJed Brown 417c4762a1bSJed Brown Output Parameter: 418c4762a1bSJed Brown . user - user-defined application context 419c4762a1bSJed Brown */ 4209371c9d4SSatish Balay static PetscErrorCode MSA_BoundaryConditions(AppCtx *user) { 421c4762a1bSJed Brown PetscInt i, j, k, limit = 0, maxits = 5; 422c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 423c4762a1bSJed Brown PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0; 424c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, three = 3.0, tol = 1e-10; 425c4762a1bSJed Brown PetscReal fnorm, det, hx, hy, xt = 0, yt = 0; 426c4762a1bSJed Brown PetscReal u1, u2, nf1, nf2, njac11, njac12, njac21, njac22; 427c4762a1bSJed Brown PetscReal b = -0.5, t = 0.5, l = -0.5, r = 0.5; 428c4762a1bSJed Brown PetscReal *boundary; 429c4762a1bSJed Brown 430c4762a1bSJed Brown PetscFunctionBegin; 4319371c9d4SSatish Balay bsize = mx + 2; 4329371c9d4SSatish Balay lsize = my + 2; 4339371c9d4SSatish Balay rsize = my + 2; 4349371c9d4SSatish Balay tsize = mx + 2; 435c4762a1bSJed Brown 4369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize, &user->bottom)); 4379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize, &user->top)); 4389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &user->left)); 4399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize, &user->right)); 440c4762a1bSJed Brown 4419371c9d4SSatish Balay hx = (r - l) / (mx + 1); 4429371c9d4SSatish Balay hy = (t - b) / (my + 1); 443c4762a1bSJed Brown 444c4762a1bSJed Brown for (j = 0; j < 4; j++) { 445c4762a1bSJed Brown if (j == 0) { 446c4762a1bSJed Brown yt = b; 447c4762a1bSJed Brown xt = l; 448c4762a1bSJed Brown limit = bsize; 449c4762a1bSJed Brown boundary = user->bottom; 450c4762a1bSJed Brown } else if (j == 1) { 451c4762a1bSJed Brown yt = t; 452c4762a1bSJed Brown xt = l; 453c4762a1bSJed Brown limit = tsize; 454c4762a1bSJed Brown boundary = user->top; 455c4762a1bSJed Brown } else if (j == 2) { 456c4762a1bSJed Brown yt = b; 457c4762a1bSJed Brown xt = l; 458c4762a1bSJed Brown limit = lsize; 459c4762a1bSJed Brown boundary = user->left; 460c4762a1bSJed Brown } else { /* if (j==3) */ 461c4762a1bSJed Brown yt = b; 462c4762a1bSJed Brown xt = r; 463c4762a1bSJed Brown limit = rsize; 464c4762a1bSJed Brown boundary = user->right; 465c4762a1bSJed Brown } 466c4762a1bSJed Brown 467c4762a1bSJed Brown for (i = 0; i < limit; i++) { 468c4762a1bSJed Brown u1 = xt; 469c4762a1bSJed Brown u2 = -yt; 470c4762a1bSJed Brown for (k = 0; k < maxits; k++) { 471c4762a1bSJed Brown nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt; 472c4762a1bSJed Brown nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt; 473c4762a1bSJed Brown fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2); 474c4762a1bSJed Brown if (fnorm <= tol) break; 475c4762a1bSJed Brown njac11 = one + u2 * u2 - u1 * u1; 476c4762a1bSJed Brown njac12 = two * u1 * u2; 477c4762a1bSJed Brown njac21 = -two * u1 * u2; 478c4762a1bSJed Brown njac22 = -one - u1 * u1 + u2 * u2; 479c4762a1bSJed Brown det = njac11 * njac22 - njac21 * njac12; 480c4762a1bSJed Brown u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det; 481c4762a1bSJed Brown u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det; 482c4762a1bSJed Brown } 483c4762a1bSJed Brown 484c4762a1bSJed Brown boundary[i] = u1 * u1 - u2 * u2; 485c4762a1bSJed Brown if (j == 0 || j == 1) { 486c4762a1bSJed Brown xt = xt + hx; 487c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 488c4762a1bSJed Brown yt = yt + hy; 489c4762a1bSJed Brown } 490c4762a1bSJed Brown } 491c4762a1bSJed Brown } 492c4762a1bSJed Brown PetscFunctionReturn(0); 493c4762a1bSJed Brown } 494c4762a1bSJed Brown 495c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 496c4762a1bSJed Brown /* 497c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 498c4762a1bSJed Brown 499c4762a1bSJed Brown Input Parameters: 500c4762a1bSJed Brown . user - user-defined application context 501c4762a1bSJed Brown . X - vector for initial guess 502c4762a1bSJed Brown 503c4762a1bSJed Brown Output Parameters: 504c4762a1bSJed Brown . X - newly computed initial guess 505c4762a1bSJed Brown */ 5069371c9d4SSatish Balay static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X) { 507c4762a1bSJed Brown PetscInt start = -1, i, j; 508c4762a1bSJed Brown PetscScalar zero = 0.0; 509c4762a1bSJed Brown PetscBool flg; 510c4762a1bSJed Brown 511c4762a1bSJed Brown PetscFunctionBegin; 5129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg)); 513c4762a1bSJed Brown 514c4762a1bSJed Brown if (flg && start == 0) { /* The zero vector is reasonable */ 5159566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero)); 516c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 517c4762a1bSJed Brown PetscInt row; 518c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 519c4762a1bSJed Brown PetscScalar *x; 520c4762a1bSJed Brown 521c4762a1bSJed Brown /* Get pointers to vector data */ 5229566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 523c4762a1bSJed Brown 524c4762a1bSJed Brown /* Perform local computations */ 525c4762a1bSJed Brown for (j = 0; j < my; j++) { 526c4762a1bSJed Brown for (i = 0; i < mx; i++) { 527c4762a1bSJed Brown row = (j)*mx + (i); 528c4762a1bSJed 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; 529c4762a1bSJed Brown } 530c4762a1bSJed Brown } 531c4762a1bSJed Brown 532c4762a1bSJed Brown /* Restore vectors */ 5339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 534c4762a1bSJed Brown } 535c4762a1bSJed Brown PetscFunctionReturn(0); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown 538c4762a1bSJed Brown /*TEST 539c4762a1bSJed Brown 540c4762a1bSJed Brown build: 541c4762a1bSJed Brown requires: !complex 542c4762a1bSJed Brown 543c4762a1bSJed Brown test: 544c4762a1bSJed Brown args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5 545c4762a1bSJed Brown requires: !single 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: 2 549c4762a1bSJed Brown args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5 550c4762a1bSJed Brown 551c4762a1bSJed Brown TEST*/ 552