1c4762a1bSJed Brown /* Program usage: mpiexec -n <proc> minsurf2 [-help] [all TAO options] */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Include "petsctao.h" so we can use TAO solvers. 5c4762a1bSJed Brown petscdmda.h for distributed array 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown #include <petsctao.h> 8c4762a1bSJed Brown #include <petscdmda.h> 9c4762a1bSJed Brown 109371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to \n\ 11c4762a1bSJed Brown solve an unconstrained minimization problem. This example is based on a \n\ 12c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\ 13c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\ 14c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\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, <st> >0 for random start, and <st> <0 \n\ 19c4762a1bSJed Brown for an average of the boundary conditions\n\n"; 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* 22c4762a1bSJed Brown User-defined application context - contains data needed by the 23c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient() 24c4762a1bSJed Brown and FormHessian(). 25c4762a1bSJed Brown */ 26c4762a1bSJed Brown typedef struct { 27c4762a1bSJed Brown PetscInt mx, my; /* discretization in x, y directions */ 28c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right; /* boundary values */ 29c4762a1bSJed Brown DM dm; /* distributed array data structure */ 30c4762a1bSJed Brown Mat H; /* Hessian */ 31c4762a1bSJed Brown } AppCtx; 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 34c4762a1bSJed Brown 35c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *); 36c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec); 37c4762a1bSJed Brown PetscErrorCode QuadraticH(AppCtx *, Vec, Mat); 38c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 39c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 40c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 41c4762a1bSJed Brown PetscErrorCode My_Monitor(Tao, void *); 42c4762a1bSJed Brown 43d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 44d71ae5a4SJacob Faibussowitsch { 45c4762a1bSJed Brown PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 46c4762a1bSJed Brown Vec x; /* solution, gradient vectors */ 47c4762a1bSJed Brown PetscBool flg, viewmat; /* flags */ 48c4762a1bSJed Brown PetscBool fddefault, fdcoloring; /* flags */ 49c4762a1bSJed Brown Tao tao; /* TAO solver context */ 50c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 51c4762a1bSJed Brown ISColoring iscoloring; 52c4762a1bSJed Brown MatFDColoring matfdcoloring; 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Initialize TAO */ 55327415f7SBarry Smith PetscFunctionBeginUser; 569566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* Specify dimension of the problem */ 599371c9d4SSatish Balay user.mx = 10; 609371c9d4SSatish Balay user.my = 10; 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg)); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg)); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "\n---- Minimum Surface Area Problem -----\n")); 6763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n", user.mx, user.my)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Let PETSc determine the vector distribution */ 709371c9d4SSatish Balay Nx = PETSC_DECIDE; 719371c9d4SSatish Balay Ny = PETSC_DECIDE; 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Create distributed array (DM) to manage parallel grid and vectors */ 749566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, user.mx, user.my, Nx, Ny, 1, 1, NULL, NULL, &user.dm)); 759566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm)); 769566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Create TAO solver and set desired solution method.*/ 799566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 809566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOCG)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* 83c4762a1bSJed Brown Extract global vector from DA for the vector of variables -- PETSC routine 84c4762a1bSJed Brown Compute the initial solution -- application specific, see below 85c4762a1bSJed Brown Set this vector for use by TAO -- TAO routine 86c4762a1bSJed Brown */ 879566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm, &x)); 889566063dSJacob Faibussowitsch PetscCall(MSA_BoundaryConditions(&user)); 899566063dSJacob Faibussowitsch PetscCall(MSA_InitialPoint(&user, x)); 909566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* 93c4762a1bSJed Brown Initialize the Application context for use in function evaluations -- application specific, see below. 94c4762a1bSJed Brown Set routines for function and gradient evaluation 95c4762a1bSJed Brown */ 969566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* 99c4762a1bSJed Brown Given the command line arguments, calculate the hessian with either the user- 100c4762a1bSJed Brown provided function FormHessian, or the default finite-difference driven Hessian 101c4762a1bSJed Brown functions 102c4762a1bSJed Brown */ 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fddefault", &fddefault)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fdcoloring", &fdcoloring)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* 107c4762a1bSJed Brown Create a matrix data structure to store the Hessian and set 108*d5b43468SJose E. Roman the Hessian evaluation routine. 109*d5b43468SJose E. Roman Set the matrix structure to be used for Hessian evaluations 110c4762a1bSJed Brown */ 1119566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm, &user.H)); 1129566063dSJacob Faibussowitsch PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown if (fdcoloring) { 1159566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(user.dm, IS_COLORING_GLOBAL, &iscoloring)); 1169566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(user.H, iscoloring, &matfdcoloring)); 1179566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 1189566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))FormGradient, (void *)&user)); 1199566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 1209566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessianColor, (void *)matfdcoloring)); 121c4762a1bSJed Brown } else if (fddefault) { 1229566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessian, (void *)NULL)); 123c4762a1bSJed Brown } else { 1249566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user)); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* 128c4762a1bSJed Brown If my_monitor option is in command line, then use the user-provided 129c4762a1bSJed Brown monitoring function 130c4762a1bSJed Brown */ 1319566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-my_monitor", &viewmat)); 1321baa6e33SBarry Smith if (viewmat) PetscCall(TaoSetMonitor(tao, My_Monitor, NULL, NULL)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Check for any tao command line options */ 1359566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1389566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(TaoView(tao, PETSC_VIEWER_STDOUT_WORLD)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Free TAO data structures */ 1439566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Free PETSc data structures */ 1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.H)); 14848a46eb9SPierre Jolivet if (fdcoloring) PetscCall(MatFDColoringDestroy(&matfdcoloring)); 1499566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom)); 1509566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top)); 1519566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right)); 1539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 155b122ec5aSJacob Faibussowitsch return 0; 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *userCtx) 159d71ae5a4SJacob Faibussowitsch { 160c4762a1bSJed Brown PetscReal fcn; 161c4762a1bSJed Brown 162c4762a1bSJed Brown PetscFunctionBegin; 1639566063dSJacob Faibussowitsch PetscCall(FormFunctionGradient(tao, X, &fcn, G, userCtx)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 168c4762a1bSJed Brown /* FormFunctionGradient - Evaluates the function and corresponding gradient. 169c4762a1bSJed Brown 170c4762a1bSJed Brown Input Parameters: 171c4762a1bSJed Brown . tao - the Tao context 172c4762a1bSJed Brown . XX - input vector 173a82e8c82SStefano Zampini . userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient() 174c4762a1bSJed Brown 175c4762a1bSJed Brown Output Parameters: 176c4762a1bSJed Brown . fcn - the newly evaluated function 177c4762a1bSJed Brown . GG - vector containing the newly evaluated gradient 178c4762a1bSJed Brown */ 179d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx) 180d71ae5a4SJacob Faibussowitsch { 181c4762a1bSJed Brown AppCtx *user = (AppCtx *)userCtx; 182c4762a1bSJed Brown PetscInt i, j; 183c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 184c4762a1bSJed Brown PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym; 185c4762a1bSJed Brown PetscReal ft = 0.0; 186c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy; 187c4762a1bSJed Brown PetscReal rhx = mx + 1, rhy = my + 1; 188c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 189c4762a1bSJed Brown PetscReal df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc; 190c4762a1bSJed Brown PetscReal **g, **x; 191c4762a1bSJed Brown Vec localX; 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscFunctionBegin; 194c4762a1bSJed Brown /* Get local mesh boundaries */ 1959566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localX)); 1969566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 1979566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* Scatter ghost points to local vector */ 2009566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 2019566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* Get pointers to vector data */ 2049566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x)); 2059566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, G, (void **)&g)); 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* Compute function and gradient over the locally owned part of the mesh */ 208c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 209c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 210c4762a1bSJed Brown xc = x[j][i]; 211c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc; 212c4762a1bSJed Brown 213c4762a1bSJed Brown if (i == 0) { /* left side */ 214c4762a1bSJed Brown xl = user->left[j - ys + 1]; 215c4762a1bSJed Brown xlt = user->left[j - ys + 2]; 216c4762a1bSJed Brown } else { 217c4762a1bSJed Brown xl = x[j][i - 1]; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown if (j == 0) { /* bottom side */ 221c4762a1bSJed Brown xb = user->bottom[i - xs + 1]; 222c4762a1bSJed Brown xrb = user->bottom[i - xs + 2]; 223c4762a1bSJed Brown } else { 224c4762a1bSJed Brown xb = x[j - 1][i]; 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227c4762a1bSJed Brown if (i + 1 == gxs + gxm) { /* right side */ 228c4762a1bSJed Brown xr = user->right[j - ys + 1]; 229c4762a1bSJed Brown xrb = user->right[j - ys]; 230c4762a1bSJed Brown } else { 231c4762a1bSJed Brown xr = x[j][i + 1]; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown if (j + 1 == gys + gym) { /* top side */ 235c4762a1bSJed Brown xt = user->top[i - xs + 1]; 236c4762a1bSJed Brown xlt = user->top[i - xs]; 237c4762a1bSJed Brown } else { 238c4762a1bSJed Brown xt = x[j + 1][i]; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241ad540459SPierre Jolivet if (i > gxs && j + 1 < gys + gym) xlt = x[j + 1][i - 1]; 242ad540459SPierre Jolivet if (j > gys && i + 1 < gxs + gxm) xrb = x[j - 1][i + 1]; 243c4762a1bSJed Brown 244c4762a1bSJed Brown d1 = (xc - xl); 245c4762a1bSJed Brown d2 = (xc - xr); 246c4762a1bSJed Brown d3 = (xc - xt); 247c4762a1bSJed Brown d4 = (xc - xb); 248c4762a1bSJed Brown d5 = (xr - xrb); 249c4762a1bSJed Brown d6 = (xrb - xb); 250c4762a1bSJed Brown d7 = (xlt - xl); 251c4762a1bSJed Brown d8 = (xt - xlt); 252c4762a1bSJed Brown 253c4762a1bSJed Brown df1dxc = d1 * hydhx; 254c4762a1bSJed Brown df2dxc = (d1 * hydhx + d4 * hxdhy); 255c4762a1bSJed Brown df3dxc = d3 * hxdhy; 256c4762a1bSJed Brown df4dxc = (d2 * hydhx + d3 * hxdhy); 257c4762a1bSJed Brown df5dxc = d2 * hydhx; 258c4762a1bSJed Brown df6dxc = d4 * hxdhy; 259c4762a1bSJed Brown 260c4762a1bSJed Brown d1 *= rhx; 261c4762a1bSJed Brown d2 *= rhx; 262c4762a1bSJed Brown d3 *= rhy; 263c4762a1bSJed Brown d4 *= rhy; 264c4762a1bSJed Brown d5 *= rhy; 265c4762a1bSJed Brown d6 *= rhx; 266c4762a1bSJed Brown d7 *= rhy; 267c4762a1bSJed Brown d8 *= rhx; 268c4762a1bSJed Brown 269c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7); 270c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 271c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8); 272c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 273c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5); 274c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6); 275c4762a1bSJed Brown 276c4762a1bSJed Brown ft = ft + (f2 + f4); 277c4762a1bSJed Brown 278c4762a1bSJed Brown df1dxc /= f1; 279c4762a1bSJed Brown df2dxc /= f2; 280c4762a1bSJed Brown df3dxc /= f3; 281c4762a1bSJed Brown df4dxc /= f4; 282c4762a1bSJed Brown df5dxc /= f5; 283c4762a1bSJed Brown df6dxc /= f6; 284c4762a1bSJed Brown 285c4762a1bSJed Brown g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5; 286c4762a1bSJed Brown } 287c4762a1bSJed Brown } 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* Compute triangular areas along the border of the domain. */ 290c4762a1bSJed Brown if (xs == 0) { /* left side */ 291c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 292c4762a1bSJed Brown d3 = (user->left[j - ys + 1] - user->left[j - ys + 2]) * rhy; 293c4762a1bSJed Brown d2 = (user->left[j - ys + 1] - x[j][0]) * rhx; 294c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown } 297c4762a1bSJed Brown if (ys == 0) { /* bottom side */ 298c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 299c4762a1bSJed Brown d2 = (user->bottom[i + 1 - xs] - user->bottom[i - xs + 2]) * rhx; 300c4762a1bSJed Brown d3 = (user->bottom[i - xs + 1] - x[0][i]) * rhy; 301c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 302c4762a1bSJed Brown } 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305c4762a1bSJed Brown if (xs + xm == mx) { /* right side */ 306c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 307c4762a1bSJed Brown d1 = (x[j][mx - 1] - user->right[j - ys + 1]) * rhx; 308c4762a1bSJed Brown d4 = (user->right[j - ys] - user->right[j - ys + 1]) * rhy; 309c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } 312c4762a1bSJed Brown if (ys + ym == my) { /* top side */ 313c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 314c4762a1bSJed Brown d1 = (x[my - 1][i] - user->top[i - xs + 1]) * rhy; 315c4762a1bSJed Brown d4 = (user->top[i - xs + 1] - user->top[i - xs]) * rhx; 316c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 317c4762a1bSJed Brown } 318c4762a1bSJed Brown } 319c4762a1bSJed Brown 320c4762a1bSJed Brown if (ys == 0 && xs == 0) { 321c4762a1bSJed Brown d1 = (user->left[0] - user->left[1]) / hy; 322c4762a1bSJed Brown d2 = (user->bottom[0] - user->bottom[1]) * rhx; 323c4762a1bSJed Brown ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2); 324c4762a1bSJed Brown } 325c4762a1bSJed Brown if (ys + ym == my && xs + xm == mx) { 326c4762a1bSJed Brown d1 = (user->right[ym + 1] - user->right[ym]) * rhy; 327c4762a1bSJed Brown d2 = (user->top[xm + 1] - user->top[xm]) * rhx; 328c4762a1bSJed Brown ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331c4762a1bSJed Brown ft = ft * area; 3329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD)); 333c4762a1bSJed Brown 334c4762a1bSJed Brown /* Restore vectors */ 3359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x)); 3369566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, G, (void **)&g)); 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* Scatter values to global vector */ 3399566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localX)); 3409566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67.0 * xm * ym)); 341c4762a1bSJed Brown PetscFunctionReturn(0); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 345c4762a1bSJed Brown /* 346c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 347c4762a1bSJed Brown 348c4762a1bSJed Brown Input Parameters: 349c4762a1bSJed Brown . tao - the Tao context 350c4762a1bSJed Brown . x - input vector 351a82e8c82SStefano Zampini . ptr - optional user-defined context, as set by TaoSetHessian() 352c4762a1bSJed Brown 353c4762a1bSJed Brown Output Parameters: 354c4762a1bSJed Brown . H - Hessian matrix 355c4762a1bSJed Brown . Hpre - optionally different preconditioning matrix 356c4762a1bSJed Brown . flg - flag indicating matrix structure 357c4762a1bSJed Brown 358c4762a1bSJed Brown */ 359d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) 360d71ae5a4SJacob Faibussowitsch { 361c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 362c4762a1bSJed Brown 363c4762a1bSJed Brown PetscFunctionBegin; 364c4762a1bSJed Brown /* Evaluate the Hessian entries*/ 3659566063dSJacob Faibussowitsch PetscCall(QuadraticH(user, X, H)); 366c4762a1bSJed Brown PetscFunctionReturn(0); 367c4762a1bSJed Brown } 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 370c4762a1bSJed Brown /* 371c4762a1bSJed Brown QuadraticH - Evaluates Hessian matrix. 372c4762a1bSJed Brown 373c4762a1bSJed Brown Input Parameters: 374c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian() 375c4762a1bSJed Brown . X - input vector 376c4762a1bSJed Brown 377c4762a1bSJed Brown Output Parameter: 378c4762a1bSJed Brown . H - Hessian matrix 379c4762a1bSJed Brown */ 380d71ae5a4SJacob Faibussowitsch PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 381d71ae5a4SJacob Faibussowitsch { 382c4762a1bSJed Brown PetscInt i, j, k; 383c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 384c4762a1bSJed Brown PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym; 385c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy; 386c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 387c4762a1bSJed Brown PetscReal hl, hr, ht, hb, hc, htl, hbr; 388c4762a1bSJed Brown PetscReal **x, v[7]; 389c4762a1bSJed Brown MatStencil col[7], row; 390c4762a1bSJed Brown Vec localX; 391c4762a1bSJed Brown PetscBool assembled; 392c4762a1bSJed Brown 393c4762a1bSJed Brown PetscFunctionBegin; 394c4762a1bSJed Brown /* Get local mesh boundaries */ 3959566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localX)); 396c4762a1bSJed Brown 3979566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 3989566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* Scatter ghost points to local vector */ 4019566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 4029566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 403c4762a1bSJed Brown 404c4762a1bSJed Brown /* Get pointers to vector data */ 4059566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x)); 406c4762a1bSJed Brown 407c4762a1bSJed Brown /* Initialize matrix entries to zero */ 4089566063dSJacob Faibussowitsch PetscCall(MatAssembled(Hessian, &assembled)); 4099566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(Hessian)); 410c4762a1bSJed Brown 411c4762a1bSJed Brown /* Set various matrix options */ 4129566063dSJacob Faibussowitsch PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE)); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */ 415c4762a1bSJed Brown 416c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 417c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 418c4762a1bSJed Brown xc = x[j][i]; 419c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc; 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* Left side */ 422c4762a1bSJed Brown if (i == 0) { 423c4762a1bSJed Brown xl = user->left[j - ys + 1]; 424c4762a1bSJed Brown xlt = user->left[j - ys + 2]; 425c4762a1bSJed Brown } else { 426c4762a1bSJed Brown xl = x[j][i - 1]; 427c4762a1bSJed Brown } 428c4762a1bSJed Brown 429c4762a1bSJed Brown if (j == 0) { 430c4762a1bSJed Brown xb = user->bottom[i - xs + 1]; 431c4762a1bSJed Brown xrb = user->bottom[i - xs + 2]; 432c4762a1bSJed Brown } else { 433c4762a1bSJed Brown xb = x[j - 1][i]; 434c4762a1bSJed Brown } 435c4762a1bSJed Brown 436c4762a1bSJed Brown if (i + 1 == mx) { 437c4762a1bSJed Brown xr = user->right[j - ys + 1]; 438c4762a1bSJed Brown xrb = user->right[j - ys]; 439c4762a1bSJed Brown } else { 440c4762a1bSJed Brown xr = x[j][i + 1]; 441c4762a1bSJed Brown } 442c4762a1bSJed Brown 443c4762a1bSJed Brown if (j + 1 == my) { 444c4762a1bSJed Brown xt = user->top[i - xs + 1]; 445c4762a1bSJed Brown xlt = user->top[i - xs]; 446c4762a1bSJed Brown } else { 447c4762a1bSJed Brown xt = x[j + 1][i]; 448c4762a1bSJed Brown } 449c4762a1bSJed Brown 450ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1]; 451ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1]; 452c4762a1bSJed Brown 453c4762a1bSJed Brown d1 = (xc - xl) / hx; 454c4762a1bSJed Brown d2 = (xc - xr) / hx; 455c4762a1bSJed Brown d3 = (xc - xt) / hy; 456c4762a1bSJed Brown d4 = (xc - xb) / hy; 457c4762a1bSJed Brown d5 = (xrb - xr) / hy; 458c4762a1bSJed Brown d6 = (xrb - xb) / hx; 459c4762a1bSJed Brown d7 = (xlt - xl) / hy; 460c4762a1bSJed Brown d8 = (xlt - xt) / hx; 461c4762a1bSJed Brown 462c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7); 463c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 464c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8); 465c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 466c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5); 467c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6); 468c4762a1bSJed Brown 4699371c9d4SSatish Balay hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2); 4709371c9d4SSatish Balay hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4); 4719371c9d4SSatish Balay ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4); 4729371c9d4SSatish Balay hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2); 473c4762a1bSJed Brown 474c4762a1bSJed Brown hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6); 475c4762a1bSJed Brown htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3); 476c4762a1bSJed Brown 4779371c9d4SSatish 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); 478c4762a1bSJed Brown 4799371c9d4SSatish Balay hl /= 2.0; 4809371c9d4SSatish Balay hr /= 2.0; 4819371c9d4SSatish Balay ht /= 2.0; 4829371c9d4SSatish Balay hb /= 2.0; 4839371c9d4SSatish Balay hbr /= 2.0; 4849371c9d4SSatish Balay htl /= 2.0; 4859371c9d4SSatish Balay hc /= 2.0; 486c4762a1bSJed Brown 4879371c9d4SSatish Balay row.j = j; 4889371c9d4SSatish Balay row.i = i; 489c4762a1bSJed Brown k = 0; 490c4762a1bSJed Brown if (j > 0) { 491c4762a1bSJed Brown v[k] = hb; 4929371c9d4SSatish Balay col[k].j = j - 1; 4939371c9d4SSatish Balay col[k].i = i; 494c4762a1bSJed Brown k++; 495c4762a1bSJed Brown } 496c4762a1bSJed Brown 497c4762a1bSJed Brown if (j > 0 && i < mx - 1) { 498c4762a1bSJed Brown v[k] = hbr; 4999371c9d4SSatish Balay col[k].j = j - 1; 5009371c9d4SSatish Balay col[k].i = i + 1; 501c4762a1bSJed Brown k++; 502c4762a1bSJed Brown } 503c4762a1bSJed Brown 504c4762a1bSJed Brown if (i > 0) { 505c4762a1bSJed Brown v[k] = hl; 5069371c9d4SSatish Balay col[k].j = j; 5079371c9d4SSatish Balay col[k].i = i - 1; 508c4762a1bSJed Brown k++; 509c4762a1bSJed Brown } 510c4762a1bSJed Brown 511c4762a1bSJed Brown v[k] = hc; 5129371c9d4SSatish Balay col[k].j = j; 5139371c9d4SSatish Balay col[k].i = i; 514c4762a1bSJed Brown k++; 515c4762a1bSJed Brown 516c4762a1bSJed Brown if (i < mx - 1) { 517c4762a1bSJed Brown v[k] = hr; 5189371c9d4SSatish Balay col[k].j = j; 5199371c9d4SSatish Balay col[k].i = i + 1; 520c4762a1bSJed Brown k++; 521c4762a1bSJed Brown } 522c4762a1bSJed Brown 523c4762a1bSJed Brown if (i > 0 && j < my - 1) { 524c4762a1bSJed Brown v[k] = htl; 5259371c9d4SSatish Balay col[k].j = j + 1; 5269371c9d4SSatish Balay col[k].i = i - 1; 527c4762a1bSJed Brown k++; 528c4762a1bSJed Brown } 529c4762a1bSJed Brown 530c4762a1bSJed Brown if (j < my - 1) { 531c4762a1bSJed Brown v[k] = ht; 5329371c9d4SSatish Balay col[k].j = j + 1; 5339371c9d4SSatish Balay col[k].i = i; 534c4762a1bSJed Brown k++; 535c4762a1bSJed Brown } 536c4762a1bSJed Brown 537c4762a1bSJed Brown /* 538c4762a1bSJed Brown Set matrix values using local numbering, which was defined 539c4762a1bSJed Brown earlier, in the main routine. 540c4762a1bSJed Brown */ 5419566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Hessian, 1, &row, k, col, v, INSERT_VALUES)); 542c4762a1bSJed Brown } 543c4762a1bSJed Brown } 544c4762a1bSJed Brown 5459566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x)); 5469566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localX)); 547c4762a1bSJed Brown 5489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY)); 5499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY)); 550c4762a1bSJed Brown 5519566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199.0 * xm * ym)); 552c4762a1bSJed Brown PetscFunctionReturn(0); 553c4762a1bSJed Brown } 554c4762a1bSJed Brown 555c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 556c4762a1bSJed Brown /* 557c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 558c4762a1bSJed Brown the region. 559c4762a1bSJed Brown 560c4762a1bSJed Brown Input Parameter: 561c4762a1bSJed Brown . user - user-defined application context 562c4762a1bSJed Brown 563c4762a1bSJed Brown Output Parameter: 564c4762a1bSJed Brown . user - user-defined application context 565c4762a1bSJed Brown */ 566d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user) 567d71ae5a4SJacob Faibussowitsch { 568c4762a1bSJed Brown PetscInt i, j, k, limit = 0, maxits = 5; 569c4762a1bSJed Brown PetscInt xs, ys, xm, ym, gxs, gys, gxm, gym; 570c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 571c4762a1bSJed Brown PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0; 572c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, three = 3.0, tol = 1e-10; 573c4762a1bSJed Brown PetscReal fnorm, det, hx, hy, xt = 0, yt = 0; 574c4762a1bSJed Brown PetscReal u1, u2, nf1, nf2, njac11, njac12, njac21, njac22; 575c4762a1bSJed Brown PetscReal b = -0.5, t = 0.5, l = -0.5, r = 0.5; 576c4762a1bSJed Brown PetscReal *boundary; 577c4762a1bSJed Brown PetscBool flg; 578c4762a1bSJed Brown 579c4762a1bSJed Brown PetscFunctionBegin; 580c4762a1bSJed Brown /* Get local mesh boundaries */ 5819566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 5829566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 583c4762a1bSJed Brown 584c4762a1bSJed Brown bsize = xm + 2; 585c4762a1bSJed Brown lsize = ym + 2; 586c4762a1bSJed Brown rsize = ym + 2; 587c4762a1bSJed Brown tsize = xm + 2; 588c4762a1bSJed Brown 5899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize, &user->bottom)); 5909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize, &user->top)); 5919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &user->left)); 5929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize, &user->right)); 593c4762a1bSJed Brown 5949371c9d4SSatish Balay hx = (r - l) / (mx + 1); 5959371c9d4SSatish Balay hy = (t - b) / (my + 1); 596c4762a1bSJed Brown 597c4762a1bSJed Brown for (j = 0; j < 4; j++) { 598c4762a1bSJed Brown if (j == 0) { 599c4762a1bSJed Brown yt = b; 600c4762a1bSJed Brown xt = l + hx * xs; 601c4762a1bSJed Brown limit = bsize; 602c4762a1bSJed Brown boundary = user->bottom; 603c4762a1bSJed Brown } else if (j == 1) { 604c4762a1bSJed Brown yt = t; 605c4762a1bSJed Brown xt = l + hx * xs; 606c4762a1bSJed Brown limit = tsize; 607c4762a1bSJed Brown boundary = user->top; 608c4762a1bSJed Brown } else if (j == 2) { 609c4762a1bSJed Brown yt = b + hy * ys; 610c4762a1bSJed Brown xt = l; 611c4762a1bSJed Brown limit = lsize; 612c4762a1bSJed Brown boundary = user->left; 613c4762a1bSJed Brown } else { /* if (j==3) */ 614c4762a1bSJed Brown yt = b + hy * ys; 615c4762a1bSJed Brown xt = r; 616c4762a1bSJed Brown limit = rsize; 617c4762a1bSJed Brown boundary = user->right; 618c4762a1bSJed Brown } 619c4762a1bSJed Brown 620c4762a1bSJed Brown for (i = 0; i < limit; i++) { 621c4762a1bSJed Brown u1 = xt; 622c4762a1bSJed Brown u2 = -yt; 623c4762a1bSJed Brown for (k = 0; k < maxits; k++) { 624c4762a1bSJed Brown nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt; 625c4762a1bSJed Brown nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt; 626c4762a1bSJed Brown fnorm = PetscSqrtReal(nf1 * nf1 + nf2 * nf2); 627c4762a1bSJed Brown if (fnorm <= tol) break; 628c4762a1bSJed Brown njac11 = one + u2 * u2 - u1 * u1; 629c4762a1bSJed Brown njac12 = two * u1 * u2; 630c4762a1bSJed Brown njac21 = -two * u1 * u2; 631c4762a1bSJed Brown njac22 = -one - u1 * u1 + u2 * u2; 632c4762a1bSJed Brown det = njac11 * njac22 - njac21 * njac12; 633c4762a1bSJed Brown u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det; 634c4762a1bSJed Brown u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det; 635c4762a1bSJed Brown } 636c4762a1bSJed Brown 637c4762a1bSJed Brown boundary[i] = u1 * u1 - u2 * u2; 638c4762a1bSJed Brown if (j == 0 || j == 1) { 639c4762a1bSJed Brown xt = xt + hx; 640c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 641c4762a1bSJed Brown yt = yt + hy; 642c4762a1bSJed Brown } 643c4762a1bSJed Brown } 644c4762a1bSJed Brown } 645c4762a1bSJed Brown 646c4762a1bSJed Brown /* Scale the boundary if desired */ 647c4762a1bSJed Brown if (1 == 1) { 648c4762a1bSJed Brown PetscReal scl = 1.0; 649c4762a1bSJed Brown 6509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg)); 651c4762a1bSJed Brown if (flg) { 652c4762a1bSJed Brown for (i = 0; i < bsize; i++) user->bottom[i] *= scl; 653c4762a1bSJed Brown } 654c4762a1bSJed Brown 6559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg)); 656c4762a1bSJed Brown if (flg) { 657c4762a1bSJed Brown for (i = 0; i < tsize; i++) user->top[i] *= scl; 658c4762a1bSJed Brown } 659c4762a1bSJed Brown 6609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg)); 661c4762a1bSJed Brown if (flg) { 662c4762a1bSJed Brown for (i = 0; i < rsize; i++) user->right[i] *= scl; 663c4762a1bSJed Brown } 664c4762a1bSJed Brown 6659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg)); 666c4762a1bSJed Brown if (flg) { 667c4762a1bSJed Brown for (i = 0; i < lsize; i++) user->left[i] *= scl; 668c4762a1bSJed Brown } 669c4762a1bSJed Brown } 670c4762a1bSJed Brown PetscFunctionReturn(0); 671c4762a1bSJed Brown } 672c4762a1bSJed Brown 673c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 674c4762a1bSJed Brown /* 675c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 676c4762a1bSJed Brown 677c4762a1bSJed Brown Input Parameters: 678c4762a1bSJed Brown . user - user-defined application context 679c4762a1bSJed Brown . X - vector for initial guess 680c4762a1bSJed Brown 681c4762a1bSJed Brown Output Parameters: 682c4762a1bSJed Brown . X - newly computed initial guess 683c4762a1bSJed Brown */ 684d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X) 685d71ae5a4SJacob Faibussowitsch { 686c4762a1bSJed Brown PetscInt start2 = -1, i, j; 687c4762a1bSJed Brown PetscReal start1 = 0; 688c4762a1bSJed Brown PetscBool flg1, flg2; 689c4762a1bSJed Brown 690c4762a1bSJed Brown PetscFunctionBegin; 6919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-start", &start1, &flg1)); 6929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-random", &start2, &flg2)); 693c4762a1bSJed Brown 694c4762a1bSJed Brown if (flg1) { /* The zero vector is reasonable */ 695c4762a1bSJed Brown 6969566063dSJacob Faibussowitsch PetscCall(VecSet(X, start1)); 697c4762a1bSJed Brown 698c4762a1bSJed Brown } else if (flg2 && start2 > 0) { /* Try a random start between -0.5 and 0.5 */ 6999371c9d4SSatish Balay PetscRandom rctx; 7009371c9d4SSatish Balay PetscReal np5 = -0.5; 701c4762a1bSJed Brown 7029566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 7039566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X, rctx)); 7049566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 7059566063dSJacob Faibussowitsch PetscCall(VecShift(X, np5)); 706c4762a1bSJed Brown 707c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 708c4762a1bSJed Brown PetscInt xs, xm, ys, ym; 709c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 710c4762a1bSJed Brown PetscReal **x; 711c4762a1bSJed Brown 712c4762a1bSJed Brown /* Get local mesh boundaries */ 7139566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 714c4762a1bSJed Brown 715c4762a1bSJed Brown /* Get pointers to vector data */ 7169566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, X, (void **)&x)); 717c4762a1bSJed Brown 718c4762a1bSJed Brown /* Perform local computations */ 719c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 720ad540459SPierre Jolivet for (i = xs; i < xs + xm; i++) x[j][i] = (((j + 1) * user->bottom[i - xs + 1] + (my - j + 1) * user->top[i - xs + 1]) / (my + 2) + ((i + 1) * user->left[j - ys + 1] + (mx - i + 1) * user->right[j - ys + 1]) / (mx + 2)) / 2.0; 721c4762a1bSJed Brown } 7229566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, X, (void **)&x)); 7239566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0 * xm * ym)); 724c4762a1bSJed Brown } 725c4762a1bSJed Brown PetscFunctionReturn(0); 726c4762a1bSJed Brown } 727c4762a1bSJed Brown 728c4762a1bSJed Brown /*-----------------------------------------------------------------------*/ 729d71ae5a4SJacob Faibussowitsch PetscErrorCode My_Monitor(Tao tao, void *ctx) 730d71ae5a4SJacob Faibussowitsch { 731c4762a1bSJed Brown Vec X; 732c4762a1bSJed Brown 733c4762a1bSJed Brown PetscFunctionBegin; 7349566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X)); 7359566063dSJacob Faibussowitsch PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 736c4762a1bSJed Brown PetscFunctionReturn(0); 737c4762a1bSJed Brown } 738c4762a1bSJed Brown 739c4762a1bSJed Brown /*TEST 740c4762a1bSJed Brown 741c4762a1bSJed Brown build: 742c4762a1bSJed Brown requires: !complex 743c4762a1bSJed Brown 744c4762a1bSJed Brown test: 745c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 746c4762a1bSJed Brown requires: !single 747c4762a1bSJed Brown 748c4762a1bSJed Brown test: 749c4762a1bSJed Brown suffix: 2 750c4762a1bSJed Brown nsize: 2 751c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4 752c4762a1bSJed Brown filter: grep -v "nls ksp" 753c4762a1bSJed Brown requires: !single 754c4762a1bSJed Brown 755c4762a1bSJed Brown test: 756c4762a1bSJed Brown suffix: 3 757c4762a1bSJed Brown nsize: 3 758c4762a1bSJed Brown args: -tao_smonitor -tao_type cg -tao_cg_type fr -mx 10 -my 10 -tao_gatol 1.e-3 759c4762a1bSJed Brown requires: !single 760c4762a1bSJed Brown 761c4762a1bSJed Brown test: 762c4762a1bSJed Brown suffix: 5 763c4762a1bSJed Brown nsize: 2 764c4762a1bSJed Brown args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 765c4762a1bSJed Brown requires: !single 766c4762a1bSJed Brown 767c4762a1bSJed Brown TEST*/ 768