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 23*4936d809SStefano Zampini application-provided call-back routines, FormFunction(), 24*4936d809SStefano Zampini FormFunctionGradient(), 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); 37*4936d809SStefano Zampini static PetscErrorCode QuadraticH(AppCtx *, Vec, Mat); 38*4936d809SStefano Zampini static PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 39*4936d809SStefano Zampini static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 40*4936d809SStefano Zampini static PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 41*4936d809SStefano Zampini static PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 42*4936d809SStefano Zampini static PetscErrorCode My_Monitor(Tao, void *); 43c4762a1bSJed Brown 44d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 45d71ae5a4SJacob Faibussowitsch { 46c4762a1bSJed Brown PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 47c4762a1bSJed Brown Vec x; /* solution, gradient vectors */ 48c4762a1bSJed Brown PetscBool flg, viewmat; /* flags */ 49c4762a1bSJed Brown PetscBool fddefault, fdcoloring; /* flags */ 50c4762a1bSJed Brown Tao tao; /* TAO solver context */ 51c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 52c4762a1bSJed Brown ISColoring iscoloring; 53c4762a1bSJed Brown MatFDColoring matfdcoloring; 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Initialize TAO */ 56327415f7SBarry Smith PetscFunctionBeginUser; 579566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Specify dimension of the problem */ 609371c9d4SSatish Balay user.mx = 10; 619371c9d4SSatish Balay user.my = 10; 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Check for any command line arguments that override defaults */ 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg)); 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "\n---- Minimum Surface Area Problem -----\n")); 6863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(MPI_COMM_WORLD, "mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n", user.mx, user.my)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Let PETSc determine the vector distribution */ 719371c9d4SSatish Balay Nx = PETSC_DECIDE; 729371c9d4SSatish Balay Ny = PETSC_DECIDE; 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Create distributed array (DM) to manage parallel grid and vectors */ 759566063dSJacob 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)); 769566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm)); 779566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm)); 78*4936d809SStefano Zampini PetscCall(DMDASetUniformCoordinates(user.dm, -0.5, 0.5, -0.5, 0.5, PETSC_DECIDE, PETSC_DECIDE)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Create TAO solver and set desired solution method.*/ 819566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 829566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOCG)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* 85c4762a1bSJed Brown Extract global vector from DA for the vector of variables -- PETSC routine 86c4762a1bSJed Brown Compute the initial solution -- application specific, see below 87c4762a1bSJed Brown Set this vector for use by TAO -- TAO routine 88c4762a1bSJed Brown */ 899566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm, &x)); 909566063dSJacob Faibussowitsch PetscCall(MSA_BoundaryConditions(&user)); 919566063dSJacob Faibussowitsch PetscCall(MSA_InitialPoint(&user, x)); 929566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* 95c4762a1bSJed Brown Initialize the Application context for use in function evaluations -- application specific, see below. 96c4762a1bSJed Brown Set routines for function and gradient evaluation 97c4762a1bSJed Brown */ 98*4936d809SStefano Zampini PetscCall(TaoSetObjective(tao, FormFunction, (void *)&user)); 999566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* 102c4762a1bSJed Brown Given the command line arguments, calculate the hessian with either the user- 103c4762a1bSJed Brown provided function FormHessian, or the default finite-difference driven Hessian 104c4762a1bSJed Brown functions 105c4762a1bSJed Brown */ 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fddefault", &fddefault)); 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fdcoloring", &fdcoloring)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* 110c4762a1bSJed Brown Create a matrix data structure to store the Hessian and set 111d5b43468SJose E. Roman the Hessian evaluation routine. 112d5b43468SJose E. Roman Set the matrix structure to be used for Hessian evaluations 113c4762a1bSJed Brown */ 1149566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm, &user.H)); 1159566063dSJacob Faibussowitsch PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown if (fdcoloring) { 1189566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(user.dm, IS_COLORING_GLOBAL, &iscoloring)); 1199566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(user.H, iscoloring, &matfdcoloring)); 1209566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 1219566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))FormGradient, (void *)&user)); 1229566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 1239566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessianColor, (void *)matfdcoloring)); 124c4762a1bSJed Brown } else if (fddefault) { 1259566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessian, (void *)NULL)); 126c4762a1bSJed Brown } else { 1279566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user)); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* 131c4762a1bSJed Brown If my_monitor option is in command line, then use the user-provided 132c4762a1bSJed Brown monitoring function 133c4762a1bSJed Brown */ 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-my_monitor", &viewmat)); 1351baa6e33SBarry Smith if (viewmat) PetscCall(TaoSetMonitor(tao, My_Monitor, NULL, NULL)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Check for any tao command line options */ 1389566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1419566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(TaoView(tao, PETSC_VIEWER_STDOUT_WORLD)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Free TAO data structures */ 1469566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* Free PETSc data structures */ 1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.H)); 15148a46eb9SPierre Jolivet if (fdcoloring) PetscCall(MatFDColoringDestroy(&matfdcoloring)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left)); 1559566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right)); 1569566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 1579566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 158b122ec5aSJacob Faibussowitsch return 0; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161*4936d809SStefano Zampini /* -------------------------------------------------------------------- */ 162*4936d809SStefano Zampini /* 163*4936d809SStefano Zampini FormFunction - Evaluates the objective function. 164*4936d809SStefano Zampini 165*4936d809SStefano Zampini Input Parameters: 166*4936d809SStefano Zampini . tao - the Tao context 167*4936d809SStefano Zampini . X - input vector 168*4936d809SStefano Zampini . userCtx - optional user-defined context, as set by TaoSetObjective() 169*4936d809SStefano Zampini 170*4936d809SStefano Zampini Output Parameters: 171*4936d809SStefano Zampini . fcn - the newly evaluated function 172*4936d809SStefano Zampini */ 173*4936d809SStefano Zampini PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *fcn, void *userCtx) 174d71ae5a4SJacob Faibussowitsch { 175*4936d809SStefano Zampini AppCtx *user = (AppCtx *)userCtx; 176*4936d809SStefano Zampini PetscInt i, j; 177*4936d809SStefano Zampini PetscInt mx = user->mx, my = user->my; 178*4936d809SStefano Zampini PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym; 179*4936d809SStefano Zampini PetscReal ft = 0.0; 180*4936d809SStefano Zampini PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), area = 0.5 * hx * hy; 181*4936d809SStefano Zampini PetscReal rhx = mx + 1, rhy = my + 1; 182*4936d809SStefano Zampini PetscReal f2, f4, d1, d2, d3, d4, xc, xl, xr, xt, xb; 183*4936d809SStefano Zampini PetscReal **x; 184*4936d809SStefano Zampini Vec localX; 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscFunctionBegin; 187*4936d809SStefano Zampini /* Get local mesh boundaries */ 188*4936d809SStefano Zampini PetscCall(DMGetLocalVector(user->dm, &localX)); 189*4936d809SStefano Zampini PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 190*4936d809SStefano Zampini PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 191*4936d809SStefano Zampini 192*4936d809SStefano Zampini /* Scatter ghost points to local vector */ 193*4936d809SStefano Zampini PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 194*4936d809SStefano Zampini PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 195*4936d809SStefano Zampini 196*4936d809SStefano Zampini /* Get pointers to vector data */ 197*4936d809SStefano Zampini PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x)); 198*4936d809SStefano Zampini 199*4936d809SStefano Zampini /* Compute function and gradient over the locally owned part of the mesh */ 200*4936d809SStefano Zampini for (j = ys; j < ys + ym; j++) { 201*4936d809SStefano Zampini for (i = xs; i < xs + xm; i++) { 202*4936d809SStefano Zampini xc = x[j][i]; 203*4936d809SStefano Zampini 204*4936d809SStefano Zampini if (i == 0) { /* left side */ 205*4936d809SStefano Zampini xl = user->left[j - ys + 1]; 206*4936d809SStefano Zampini } else { 207*4936d809SStefano Zampini xl = x[j][i - 1]; 208*4936d809SStefano Zampini } 209*4936d809SStefano Zampini 210*4936d809SStefano Zampini if (j == 0) { /* bottom side */ 211*4936d809SStefano Zampini xb = user->bottom[i - xs + 1]; 212*4936d809SStefano Zampini } else { 213*4936d809SStefano Zampini xb = x[j - 1][i]; 214*4936d809SStefano Zampini } 215*4936d809SStefano Zampini 216*4936d809SStefano Zampini if (i + 1 == gxs + gxm) { /* right side */ 217*4936d809SStefano Zampini xr = user->right[j - ys + 1]; 218*4936d809SStefano Zampini } else { 219*4936d809SStefano Zampini xr = x[j][i + 1]; 220*4936d809SStefano Zampini } 221*4936d809SStefano Zampini 222*4936d809SStefano Zampini if (j + 1 == gys + gym) { /* top side */ 223*4936d809SStefano Zampini xt = user->top[i - xs + 1]; 224*4936d809SStefano Zampini } else { 225*4936d809SStefano Zampini xt = x[j + 1][i]; 226*4936d809SStefano Zampini } 227*4936d809SStefano Zampini 228*4936d809SStefano Zampini d1 = (xc - xl); 229*4936d809SStefano Zampini d2 = (xc - xr); 230*4936d809SStefano Zampini d3 = (xc - xt); 231*4936d809SStefano Zampini d4 = (xc - xb); 232*4936d809SStefano Zampini 233*4936d809SStefano Zampini d1 *= rhx; 234*4936d809SStefano Zampini d2 *= rhx; 235*4936d809SStefano Zampini d3 *= rhy; 236*4936d809SStefano Zampini d4 *= rhy; 237*4936d809SStefano Zampini 238*4936d809SStefano Zampini f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 239*4936d809SStefano Zampini f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 240*4936d809SStefano Zampini 241*4936d809SStefano Zampini ft = ft + (f2 + f4); 242*4936d809SStefano Zampini } 243*4936d809SStefano Zampini } 244*4936d809SStefano Zampini 245*4936d809SStefano Zampini /* Compute triangular areas along the border of the domain. */ 246*4936d809SStefano Zampini if (xs == 0) { /* left side */ 247*4936d809SStefano Zampini for (j = ys; j < ys + ym; j++) { 248*4936d809SStefano Zampini d3 = (user->left[j - ys + 1] - user->left[j - ys + 2]) * rhy; 249*4936d809SStefano Zampini d2 = (user->left[j - ys + 1] - x[j][0]) * rhx; 250*4936d809SStefano Zampini ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 251*4936d809SStefano Zampini } 252*4936d809SStefano Zampini } 253*4936d809SStefano Zampini if (ys == 0) { /* bottom side */ 254*4936d809SStefano Zampini for (i = xs; i < xs + xm; i++) { 255*4936d809SStefano Zampini d2 = (user->bottom[i + 1 - xs] - user->bottom[i - xs + 2]) * rhx; 256*4936d809SStefano Zampini d3 = (user->bottom[i - xs + 1] - x[0][i]) * rhy; 257*4936d809SStefano Zampini ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 258*4936d809SStefano Zampini } 259*4936d809SStefano Zampini } 260*4936d809SStefano Zampini 261*4936d809SStefano Zampini if (xs + xm == mx) { /* right side */ 262*4936d809SStefano Zampini for (j = ys; j < ys + ym; j++) { 263*4936d809SStefano Zampini d1 = (x[j][mx - 1] - user->right[j - ys + 1]) * rhx; 264*4936d809SStefano Zampini d4 = (user->right[j - ys] - user->right[j - ys + 1]) * rhy; 265*4936d809SStefano Zampini ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 266*4936d809SStefano Zampini } 267*4936d809SStefano Zampini } 268*4936d809SStefano Zampini if (ys + ym == my) { /* top side */ 269*4936d809SStefano Zampini for (i = xs; i < xs + xm; i++) { 270*4936d809SStefano Zampini d1 = (x[my - 1][i] - user->top[i - xs + 1]) * rhy; 271*4936d809SStefano Zampini d4 = (user->top[i - xs + 1] - user->top[i - xs]) * rhx; 272*4936d809SStefano Zampini ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 273*4936d809SStefano Zampini } 274*4936d809SStefano Zampini } 275*4936d809SStefano Zampini 276*4936d809SStefano Zampini if (ys == 0 && xs == 0) { 277*4936d809SStefano Zampini d1 = (user->left[0] - user->left[1]) / hy; 278*4936d809SStefano Zampini d2 = (user->bottom[0] - user->bottom[1]) * rhx; 279*4936d809SStefano Zampini ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2); 280*4936d809SStefano Zampini } 281*4936d809SStefano Zampini if (ys + ym == my && xs + xm == mx) { 282*4936d809SStefano Zampini d1 = (user->right[ym + 1] - user->right[ym]) * rhy; 283*4936d809SStefano Zampini d2 = (user->top[xm + 1] - user->top[xm]) * rhx; 284*4936d809SStefano Zampini ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2); 285*4936d809SStefano Zampini } 286*4936d809SStefano Zampini 287*4936d809SStefano Zampini ft = ft * area; 288*4936d809SStefano Zampini PetscCall(MPIU_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD)); 289*4936d809SStefano Zampini 290*4936d809SStefano Zampini /* Restore vectors */ 291*4936d809SStefano Zampini PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x)); 292*4936d809SStefano Zampini PetscCall(DMRestoreLocalVector(user->dm, &localX)); 2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 297*4936d809SStefano Zampini /* 298*4936d809SStefano Zampini FormFunctionGradient - Evaluates the function and corresponding gradient. 299c4762a1bSJed Brown 300c4762a1bSJed Brown Input Parameters: 301c4762a1bSJed Brown . tao - the Tao context 302*4936d809SStefano Zampini . X - input vector 303a82e8c82SStefano Zampini . userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient() 304c4762a1bSJed Brown 305c4762a1bSJed Brown Output Parameters: 306c4762a1bSJed Brown . fcn - the newly evaluated function 307*4936d809SStefano Zampini . G - vector containing the newly evaluated gradient 308c4762a1bSJed Brown */ 309d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx) 310d71ae5a4SJacob Faibussowitsch { 311c4762a1bSJed Brown AppCtx *user = (AppCtx *)userCtx; 312c4762a1bSJed Brown PetscInt i, j; 313c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 314c4762a1bSJed Brown PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym; 315c4762a1bSJed Brown PetscReal ft = 0.0; 316c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy; 317c4762a1bSJed Brown PetscReal rhx = mx + 1, rhy = my + 1; 318c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 319c4762a1bSJed Brown PetscReal df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc; 320c4762a1bSJed Brown PetscReal **g, **x; 321c4762a1bSJed Brown Vec localX; 322c4762a1bSJed Brown 323c4762a1bSJed Brown PetscFunctionBegin; 324c4762a1bSJed Brown /* Get local mesh boundaries */ 3259566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localX)); 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 3279566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* Scatter ghost points to local vector */ 3309566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 3319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 332c4762a1bSJed Brown 333c4762a1bSJed Brown /* Get pointers to vector data */ 3349566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x)); 3359566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, G, (void **)&g)); 336c4762a1bSJed Brown 337c4762a1bSJed Brown /* Compute function and gradient over the locally owned part of the mesh */ 338c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 339c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 340c4762a1bSJed Brown xc = x[j][i]; 341c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc; 342c4762a1bSJed Brown 343c4762a1bSJed Brown if (i == 0) { /* left side */ 344c4762a1bSJed Brown xl = user->left[j - ys + 1]; 345c4762a1bSJed Brown xlt = user->left[j - ys + 2]; 346c4762a1bSJed Brown } else { 347c4762a1bSJed Brown xl = x[j][i - 1]; 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 350c4762a1bSJed Brown if (j == 0) { /* bottom side */ 351c4762a1bSJed Brown xb = user->bottom[i - xs + 1]; 352c4762a1bSJed Brown xrb = user->bottom[i - xs + 2]; 353c4762a1bSJed Brown } else { 354c4762a1bSJed Brown xb = x[j - 1][i]; 355c4762a1bSJed Brown } 356c4762a1bSJed Brown 357c4762a1bSJed Brown if (i + 1 == gxs + gxm) { /* right side */ 358c4762a1bSJed Brown xr = user->right[j - ys + 1]; 359c4762a1bSJed Brown xrb = user->right[j - ys]; 360c4762a1bSJed Brown } else { 361c4762a1bSJed Brown xr = x[j][i + 1]; 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364c4762a1bSJed Brown if (j + 1 == gys + gym) { /* top side */ 365c4762a1bSJed Brown xt = user->top[i - xs + 1]; 366c4762a1bSJed Brown xlt = user->top[i - xs]; 367c4762a1bSJed Brown } else { 368c4762a1bSJed Brown xt = x[j + 1][i]; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown 371ad540459SPierre Jolivet if (i > gxs && j + 1 < gys + gym) xlt = x[j + 1][i - 1]; 372ad540459SPierre Jolivet if (j > gys && i + 1 < gxs + gxm) xrb = x[j - 1][i + 1]; 373c4762a1bSJed Brown 374c4762a1bSJed Brown d1 = (xc - xl); 375c4762a1bSJed Brown d2 = (xc - xr); 376c4762a1bSJed Brown d3 = (xc - xt); 377c4762a1bSJed Brown d4 = (xc - xb); 378c4762a1bSJed Brown d5 = (xr - xrb); 379c4762a1bSJed Brown d6 = (xrb - xb); 380c4762a1bSJed Brown d7 = (xlt - xl); 381c4762a1bSJed Brown d8 = (xt - xlt); 382c4762a1bSJed Brown 383c4762a1bSJed Brown df1dxc = d1 * hydhx; 384c4762a1bSJed Brown df2dxc = (d1 * hydhx + d4 * hxdhy); 385c4762a1bSJed Brown df3dxc = d3 * hxdhy; 386c4762a1bSJed Brown df4dxc = (d2 * hydhx + d3 * hxdhy); 387c4762a1bSJed Brown df5dxc = d2 * hydhx; 388c4762a1bSJed Brown df6dxc = d4 * hxdhy; 389c4762a1bSJed Brown 390c4762a1bSJed Brown d1 *= rhx; 391c4762a1bSJed Brown d2 *= rhx; 392c4762a1bSJed Brown d3 *= rhy; 393c4762a1bSJed Brown d4 *= rhy; 394c4762a1bSJed Brown d5 *= rhy; 395c4762a1bSJed Brown d6 *= rhx; 396c4762a1bSJed Brown d7 *= rhy; 397c4762a1bSJed Brown d8 *= rhx; 398c4762a1bSJed Brown 399c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7); 400c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 401c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8); 402c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 403c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5); 404c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6); 405c4762a1bSJed Brown 406c4762a1bSJed Brown ft = ft + (f2 + f4); 407c4762a1bSJed Brown 408c4762a1bSJed Brown df1dxc /= f1; 409c4762a1bSJed Brown df2dxc /= f2; 410c4762a1bSJed Brown df3dxc /= f3; 411c4762a1bSJed Brown df4dxc /= f4; 412c4762a1bSJed Brown df5dxc /= f5; 413c4762a1bSJed Brown df6dxc /= f6; 414c4762a1bSJed Brown 415c4762a1bSJed Brown g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5; 416c4762a1bSJed Brown } 417c4762a1bSJed Brown } 418c4762a1bSJed Brown 419c4762a1bSJed Brown /* Compute triangular areas along the border of the domain. */ 420c4762a1bSJed Brown if (xs == 0) { /* left side */ 421c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 422c4762a1bSJed Brown d3 = (user->left[j - ys + 1] - user->left[j - ys + 2]) * rhy; 423c4762a1bSJed Brown d2 = (user->left[j - ys + 1] - x[j][0]) * rhx; 424c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown } 427c4762a1bSJed Brown if (ys == 0) { /* bottom side */ 428c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 429c4762a1bSJed Brown d2 = (user->bottom[i + 1 - xs] - user->bottom[i - xs + 2]) * rhx; 430c4762a1bSJed Brown d3 = (user->bottom[i - xs + 1] - x[0][i]) * rhy; 431c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 432c4762a1bSJed Brown } 433c4762a1bSJed Brown } 434c4762a1bSJed Brown 435c4762a1bSJed Brown if (xs + xm == mx) { /* right side */ 436c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 437c4762a1bSJed Brown d1 = (x[j][mx - 1] - user->right[j - ys + 1]) * rhx; 438c4762a1bSJed Brown d4 = (user->right[j - ys] - user->right[j - ys + 1]) * rhy; 439c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 440c4762a1bSJed Brown } 441c4762a1bSJed Brown } 442c4762a1bSJed Brown if (ys + ym == my) { /* top side */ 443c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 444c4762a1bSJed Brown d1 = (x[my - 1][i] - user->top[i - xs + 1]) * rhy; 445c4762a1bSJed Brown d4 = (user->top[i - xs + 1] - user->top[i - xs]) * rhx; 446c4762a1bSJed Brown ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 447c4762a1bSJed Brown } 448c4762a1bSJed Brown } 449c4762a1bSJed Brown 450c4762a1bSJed Brown if (ys == 0 && xs == 0) { 451c4762a1bSJed Brown d1 = (user->left[0] - user->left[1]) / hy; 452c4762a1bSJed Brown d2 = (user->bottom[0] - user->bottom[1]) * rhx; 453c4762a1bSJed Brown ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown if (ys + ym == my && xs + xm == mx) { 456c4762a1bSJed Brown d1 = (user->right[ym + 1] - user->right[ym]) * rhy; 457c4762a1bSJed Brown d2 = (user->top[xm + 1] - user->top[xm]) * rhx; 458c4762a1bSJed Brown ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2); 459c4762a1bSJed Brown } 460c4762a1bSJed Brown 461c4762a1bSJed Brown ft = ft * area; 462712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD)); 463c4762a1bSJed Brown 464c4762a1bSJed Brown /* Restore vectors */ 4659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x)); 4669566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, G, (void **)&g)); 4679566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localX)); 4689566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67.0 * xm * ym)); 4693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 470c4762a1bSJed Brown } 471c4762a1bSJed Brown 472*4936d809SStefano Zampini PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *userCtx) 473*4936d809SStefano Zampini { 474*4936d809SStefano Zampini PetscReal fcn; 475*4936d809SStefano Zampini 476*4936d809SStefano Zampini PetscFunctionBegin; 477*4936d809SStefano Zampini PetscCall(FormFunctionGradient(tao, X, &fcn, G, userCtx)); 478*4936d809SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 479*4936d809SStefano Zampini } 480*4936d809SStefano Zampini 481c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 482c4762a1bSJed Brown /* 483c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 484c4762a1bSJed Brown 485c4762a1bSJed Brown Input Parameters: 486c4762a1bSJed Brown . tao - the Tao context 487c4762a1bSJed Brown . x - input vector 488a82e8c82SStefano Zampini . ptr - optional user-defined context, as set by TaoSetHessian() 489c4762a1bSJed Brown 490c4762a1bSJed Brown Output Parameters: 491c4762a1bSJed Brown . H - Hessian matrix 492c4762a1bSJed Brown . Hpre - optionally different preconditioning matrix 493c4762a1bSJed Brown . flg - flag indicating matrix structure 494c4762a1bSJed Brown 495c4762a1bSJed Brown */ 496d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) 497d71ae5a4SJacob Faibussowitsch { 498c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 499c4762a1bSJed Brown 500c4762a1bSJed Brown PetscFunctionBegin; 501c4762a1bSJed Brown /* Evaluate the Hessian entries*/ 5029566063dSJacob Faibussowitsch PetscCall(QuadraticH(user, X, H)); 5033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 507c4762a1bSJed Brown /* 508c4762a1bSJed Brown QuadraticH - Evaluates Hessian matrix. 509c4762a1bSJed Brown 510c4762a1bSJed Brown Input Parameters: 511c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian() 512c4762a1bSJed Brown . X - input vector 513c4762a1bSJed Brown 514c4762a1bSJed Brown Output Parameter: 515c4762a1bSJed Brown . H - Hessian matrix 516c4762a1bSJed Brown */ 517d71ae5a4SJacob Faibussowitsch PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 518d71ae5a4SJacob Faibussowitsch { 519c4762a1bSJed Brown PetscInt i, j, k; 520c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 521c4762a1bSJed Brown PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym; 522c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy; 523c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 524c4762a1bSJed Brown PetscReal hl, hr, ht, hb, hc, htl, hbr; 525c4762a1bSJed Brown PetscReal **x, v[7]; 526c4762a1bSJed Brown MatStencil col[7], row; 527c4762a1bSJed Brown Vec localX; 528c4762a1bSJed Brown PetscBool assembled; 529c4762a1bSJed Brown 530c4762a1bSJed Brown PetscFunctionBegin; 531c4762a1bSJed Brown /* Get local mesh boundaries */ 5329566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localX)); 533c4762a1bSJed Brown 5349566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 5359566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 536c4762a1bSJed Brown 537c4762a1bSJed Brown /* Scatter ghost points to local vector */ 5389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 5399566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 540c4762a1bSJed Brown 541c4762a1bSJed Brown /* Get pointers to vector data */ 5429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x)); 543c4762a1bSJed Brown 544c4762a1bSJed Brown /* Initialize matrix entries to zero */ 5459566063dSJacob Faibussowitsch PetscCall(MatAssembled(Hessian, &assembled)); 5469566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(Hessian)); 547c4762a1bSJed Brown 548c4762a1bSJed Brown /* Set various matrix options */ 5499566063dSJacob Faibussowitsch PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE)); 550c4762a1bSJed Brown 551c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */ 552c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 553c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 554c4762a1bSJed Brown xc = x[j][i]; 555c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc; 556c4762a1bSJed Brown 557c4762a1bSJed Brown /* Left side */ 558c4762a1bSJed Brown if (i == 0) { 559c4762a1bSJed Brown xl = user->left[j - ys + 1]; 560c4762a1bSJed Brown xlt = user->left[j - ys + 2]; 561c4762a1bSJed Brown } else { 562c4762a1bSJed Brown xl = x[j][i - 1]; 563c4762a1bSJed Brown } 564c4762a1bSJed Brown 565c4762a1bSJed Brown if (j == 0) { 566c4762a1bSJed Brown xb = user->bottom[i - xs + 1]; 567c4762a1bSJed Brown xrb = user->bottom[i - xs + 2]; 568c4762a1bSJed Brown } else { 569c4762a1bSJed Brown xb = x[j - 1][i]; 570c4762a1bSJed Brown } 571c4762a1bSJed Brown 572c4762a1bSJed Brown if (i + 1 == mx) { 573c4762a1bSJed Brown xr = user->right[j - ys + 1]; 574c4762a1bSJed Brown xrb = user->right[j - ys]; 575c4762a1bSJed Brown } else { 576c4762a1bSJed Brown xr = x[j][i + 1]; 577c4762a1bSJed Brown } 578c4762a1bSJed Brown 579c4762a1bSJed Brown if (j + 1 == my) { 580c4762a1bSJed Brown xt = user->top[i - xs + 1]; 581c4762a1bSJed Brown xlt = user->top[i - xs]; 582c4762a1bSJed Brown } else { 583c4762a1bSJed Brown xt = x[j + 1][i]; 584c4762a1bSJed Brown } 585c4762a1bSJed Brown 586ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1]; 587ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1]; 588c4762a1bSJed Brown 589c4762a1bSJed Brown d1 = (xc - xl) / hx; 590c4762a1bSJed Brown d2 = (xc - xr) / hx; 591c4762a1bSJed Brown d3 = (xc - xt) / hy; 592c4762a1bSJed Brown d4 = (xc - xb) / hy; 593c4762a1bSJed Brown d5 = (xrb - xr) / hy; 594c4762a1bSJed Brown d6 = (xrb - xb) / hx; 595c4762a1bSJed Brown d7 = (xlt - xl) / hy; 596c4762a1bSJed Brown d8 = (xlt - xt) / hx; 597c4762a1bSJed Brown 598c4762a1bSJed Brown f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7); 599c4762a1bSJed Brown f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4); 600c4762a1bSJed Brown f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8); 601c4762a1bSJed Brown f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2); 602c4762a1bSJed Brown f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5); 603c4762a1bSJed Brown f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6); 604c4762a1bSJed Brown 6059371c9d4SSatish Balay hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2); 6069371c9d4SSatish Balay hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4); 6079371c9d4SSatish Balay ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4); 6089371c9d4SSatish Balay hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2); 609c4762a1bSJed Brown 610c4762a1bSJed Brown hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6); 611c4762a1bSJed Brown htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3); 612c4762a1bSJed Brown 6139371c9d4SSatish 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); 614c4762a1bSJed Brown 6159371c9d4SSatish Balay hl /= 2.0; 6169371c9d4SSatish Balay hr /= 2.0; 6179371c9d4SSatish Balay ht /= 2.0; 6189371c9d4SSatish Balay hb /= 2.0; 6199371c9d4SSatish Balay hbr /= 2.0; 6209371c9d4SSatish Balay htl /= 2.0; 6219371c9d4SSatish Balay hc /= 2.0; 622c4762a1bSJed Brown 6239371c9d4SSatish Balay row.j = j; 6249371c9d4SSatish Balay row.i = i; 625c4762a1bSJed Brown k = 0; 626c4762a1bSJed Brown if (j > 0) { 627c4762a1bSJed Brown v[k] = hb; 6289371c9d4SSatish Balay col[k].j = j - 1; 6299371c9d4SSatish Balay col[k].i = i; 630c4762a1bSJed Brown k++; 631c4762a1bSJed Brown } 632c4762a1bSJed Brown 633c4762a1bSJed Brown if (j > 0 && i < mx - 1) { 634c4762a1bSJed Brown v[k] = hbr; 6359371c9d4SSatish Balay col[k].j = j - 1; 6369371c9d4SSatish Balay col[k].i = i + 1; 637c4762a1bSJed Brown k++; 638c4762a1bSJed Brown } 639c4762a1bSJed Brown 640c4762a1bSJed Brown if (i > 0) { 641c4762a1bSJed Brown v[k] = hl; 6429371c9d4SSatish Balay col[k].j = j; 6439371c9d4SSatish Balay col[k].i = i - 1; 644c4762a1bSJed Brown k++; 645c4762a1bSJed Brown } 646c4762a1bSJed Brown 647c4762a1bSJed Brown v[k] = hc; 6489371c9d4SSatish Balay col[k].j = j; 6499371c9d4SSatish Balay col[k].i = i; 650c4762a1bSJed Brown k++; 651c4762a1bSJed Brown 652c4762a1bSJed Brown if (i < mx - 1) { 653c4762a1bSJed Brown v[k] = hr; 6549371c9d4SSatish Balay col[k].j = j; 6559371c9d4SSatish Balay col[k].i = i + 1; 656c4762a1bSJed Brown k++; 657c4762a1bSJed Brown } 658c4762a1bSJed Brown 659c4762a1bSJed Brown if (i > 0 && j < my - 1) { 660c4762a1bSJed Brown v[k] = htl; 6619371c9d4SSatish Balay col[k].j = j + 1; 6629371c9d4SSatish Balay col[k].i = i - 1; 663c4762a1bSJed Brown k++; 664c4762a1bSJed Brown } 665c4762a1bSJed Brown 666c4762a1bSJed Brown if (j < my - 1) { 667c4762a1bSJed Brown v[k] = ht; 6689371c9d4SSatish Balay col[k].j = j + 1; 6699371c9d4SSatish Balay col[k].i = i; 670c4762a1bSJed Brown k++; 671c4762a1bSJed Brown } 672c4762a1bSJed Brown 673c4762a1bSJed Brown /* 674c4762a1bSJed Brown Set matrix values using local numbering, which was defined 675c4762a1bSJed Brown earlier, in the main routine. 676c4762a1bSJed Brown */ 6779566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Hessian, 1, &row, k, col, v, INSERT_VALUES)); 678c4762a1bSJed Brown } 679c4762a1bSJed Brown } 680c4762a1bSJed Brown 6819566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x)); 6829566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localX)); 683c4762a1bSJed Brown 6849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY)); 6859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY)); 686c4762a1bSJed Brown 6879566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199.0 * xm * ym)); 6883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 689c4762a1bSJed Brown } 690c4762a1bSJed Brown 691c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 692c4762a1bSJed Brown /* 693c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for 694c4762a1bSJed Brown the region. 695c4762a1bSJed Brown 696c4762a1bSJed Brown Input Parameter: 697c4762a1bSJed Brown . user - user-defined application context 698c4762a1bSJed Brown 699c4762a1bSJed Brown Output Parameter: 700c4762a1bSJed Brown . user - user-defined application context 701c4762a1bSJed Brown */ 702d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user) 703d71ae5a4SJacob Faibussowitsch { 704c4762a1bSJed Brown PetscInt i, j, k, limit = 0, maxits = 5; 705c4762a1bSJed Brown PetscInt xs, ys, xm, ym, gxs, gys, gxm, gym; 706c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 707c4762a1bSJed Brown PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0; 708c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, three = 3.0, tol = 1e-10; 709c4762a1bSJed Brown PetscReal fnorm, det, hx, hy, xt = 0, yt = 0; 710c4762a1bSJed Brown PetscReal u1, u2, nf1, nf2, njac11, njac12, njac21, njac22; 711c4762a1bSJed Brown PetscReal b = -0.5, t = 0.5, l = -0.5, r = 0.5; 712c4762a1bSJed Brown PetscReal *boundary; 713c4762a1bSJed Brown PetscBool flg; 714c4762a1bSJed Brown 715c4762a1bSJed Brown PetscFunctionBegin; 716c4762a1bSJed Brown /* Get local mesh boundaries */ 7179566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 7189566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 719c4762a1bSJed Brown 720c4762a1bSJed Brown bsize = xm + 2; 721c4762a1bSJed Brown lsize = ym + 2; 722c4762a1bSJed Brown rsize = ym + 2; 723c4762a1bSJed Brown tsize = xm + 2; 724c4762a1bSJed Brown 7259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize, &user->bottom)); 7269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize, &user->top)); 7279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &user->left)); 7289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize, &user->right)); 729c4762a1bSJed Brown 7309371c9d4SSatish Balay hx = (r - l) / (mx + 1); 7319371c9d4SSatish Balay hy = (t - b) / (my + 1); 732c4762a1bSJed Brown 733c4762a1bSJed Brown for (j = 0; j < 4; j++) { 734c4762a1bSJed Brown if (j == 0) { 735c4762a1bSJed Brown yt = b; 736c4762a1bSJed Brown xt = l + hx * xs; 737c4762a1bSJed Brown limit = bsize; 738c4762a1bSJed Brown boundary = user->bottom; 739c4762a1bSJed Brown } else if (j == 1) { 740c4762a1bSJed Brown yt = t; 741c4762a1bSJed Brown xt = l + hx * xs; 742c4762a1bSJed Brown limit = tsize; 743c4762a1bSJed Brown boundary = user->top; 744c4762a1bSJed Brown } else if (j == 2) { 745c4762a1bSJed Brown yt = b + hy * ys; 746c4762a1bSJed Brown xt = l; 747c4762a1bSJed Brown limit = lsize; 748c4762a1bSJed Brown boundary = user->left; 749c4762a1bSJed Brown } else { /* if (j==3) */ 750c4762a1bSJed Brown yt = b + hy * ys; 751c4762a1bSJed Brown xt = r; 752c4762a1bSJed Brown limit = rsize; 753c4762a1bSJed Brown boundary = user->right; 754c4762a1bSJed Brown } 755c4762a1bSJed Brown 756c4762a1bSJed Brown for (i = 0; i < limit; i++) { 757c4762a1bSJed Brown u1 = xt; 758c4762a1bSJed Brown u2 = -yt; 759c4762a1bSJed Brown for (k = 0; k < maxits; k++) { 760c4762a1bSJed Brown nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt; 761c4762a1bSJed Brown nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt; 762c4762a1bSJed Brown fnorm = PetscSqrtReal(nf1 * nf1 + nf2 * nf2); 763c4762a1bSJed Brown if (fnorm <= tol) break; 764c4762a1bSJed Brown njac11 = one + u2 * u2 - u1 * u1; 765c4762a1bSJed Brown njac12 = two * u1 * u2; 766c4762a1bSJed Brown njac21 = -two * u1 * u2; 767c4762a1bSJed Brown njac22 = -one - u1 * u1 + u2 * u2; 768c4762a1bSJed Brown det = njac11 * njac22 - njac21 * njac12; 769c4762a1bSJed Brown u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det; 770c4762a1bSJed Brown u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det; 771c4762a1bSJed Brown } 772c4762a1bSJed Brown 773c4762a1bSJed Brown boundary[i] = u1 * u1 - u2 * u2; 774c4762a1bSJed Brown if (j == 0 || j == 1) { 775c4762a1bSJed Brown xt = xt + hx; 776c4762a1bSJed Brown } else { /* if (j==2 || j==3) */ 777c4762a1bSJed Brown yt = yt + hy; 778c4762a1bSJed Brown } 779c4762a1bSJed Brown } 780c4762a1bSJed Brown } 781c4762a1bSJed Brown 782c4762a1bSJed Brown /* Scale the boundary if desired */ 783c4762a1bSJed Brown if (1 == 1) { 784c4762a1bSJed Brown PetscReal scl = 1.0; 785c4762a1bSJed Brown 7869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg)); 787c4762a1bSJed Brown if (flg) { 788c4762a1bSJed Brown for (i = 0; i < bsize; i++) user->bottom[i] *= scl; 789c4762a1bSJed Brown } 790c4762a1bSJed Brown 7919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg)); 792c4762a1bSJed Brown if (flg) { 793c4762a1bSJed Brown for (i = 0; i < tsize; i++) user->top[i] *= scl; 794c4762a1bSJed Brown } 795c4762a1bSJed Brown 7969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg)); 797c4762a1bSJed Brown if (flg) { 798c4762a1bSJed Brown for (i = 0; i < rsize; i++) user->right[i] *= scl; 799c4762a1bSJed Brown } 800c4762a1bSJed Brown 8019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg)); 802c4762a1bSJed Brown if (flg) { 803c4762a1bSJed Brown for (i = 0; i < lsize; i++) user->left[i] *= scl; 804c4762a1bSJed Brown } 805c4762a1bSJed Brown } 8063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 807c4762a1bSJed Brown } 808c4762a1bSJed Brown 809c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 810c4762a1bSJed Brown /* 811c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways. 812c4762a1bSJed Brown 813c4762a1bSJed Brown Input Parameters: 814c4762a1bSJed Brown . user - user-defined application context 815c4762a1bSJed Brown . X - vector for initial guess 816c4762a1bSJed Brown 817c4762a1bSJed Brown Output Parameters: 818c4762a1bSJed Brown . X - newly computed initial guess 819c4762a1bSJed Brown */ 820d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X) 821d71ae5a4SJacob Faibussowitsch { 822c4762a1bSJed Brown PetscInt start2 = -1, i, j; 823c4762a1bSJed Brown PetscReal start1 = 0; 824c4762a1bSJed Brown PetscBool flg1, flg2; 825c4762a1bSJed Brown 826c4762a1bSJed Brown PetscFunctionBegin; 8279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-start", &start1, &flg1)); 8289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-random", &start2, &flg2)); 829c4762a1bSJed Brown 830c4762a1bSJed Brown if (flg1) { /* The zero vector is reasonable */ 8319566063dSJacob Faibussowitsch PetscCall(VecSet(X, start1)); 832c4762a1bSJed Brown } else if (flg2 && start2 > 0) { /* Try a random start between -0.5 and 0.5 */ 8339371c9d4SSatish Balay PetscRandom rctx; 8349371c9d4SSatish Balay PetscReal np5 = -0.5; 835c4762a1bSJed Brown 8369566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 8379566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X, rctx)); 8389566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rctx)); 8399566063dSJacob Faibussowitsch PetscCall(VecShift(X, np5)); 840c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */ 841c4762a1bSJed Brown PetscInt xs, xm, ys, ym; 842c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my; 843c4762a1bSJed Brown PetscReal **x; 844c4762a1bSJed Brown 845c4762a1bSJed Brown /* Get local mesh boundaries */ 8469566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 847c4762a1bSJed Brown 848c4762a1bSJed Brown /* Get pointers to vector data */ 8499566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->dm, X, (void **)&x)); 850c4762a1bSJed Brown 851c4762a1bSJed Brown /* Perform local computations */ 852c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 853ad540459SPierre 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; 854c4762a1bSJed Brown } 8559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->dm, X, (void **)&x)); 8569566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(9.0 * xm * ym)); 857c4762a1bSJed Brown } 8583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 859c4762a1bSJed Brown } 860c4762a1bSJed Brown 861c4762a1bSJed Brown /*-----------------------------------------------------------------------*/ 862d71ae5a4SJacob Faibussowitsch PetscErrorCode My_Monitor(Tao tao, void *ctx) 863d71ae5a4SJacob Faibussowitsch { 864c4762a1bSJed Brown Vec X; 865c4762a1bSJed Brown 866c4762a1bSJed Brown PetscFunctionBegin; 8679566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X)); 8689566063dSJacob Faibussowitsch PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 8693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 870c4762a1bSJed Brown } 871c4762a1bSJed Brown 872c4762a1bSJed Brown /*TEST 873c4762a1bSJed Brown 874c4762a1bSJed Brown build: 875c4762a1bSJed Brown requires: !complex 876c4762a1bSJed Brown 877c4762a1bSJed Brown test: 878c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 879c4762a1bSJed Brown requires: !single 880c4762a1bSJed Brown 881c4762a1bSJed Brown test: 882c4762a1bSJed Brown suffix: 2 883c4762a1bSJed Brown nsize: 2 884c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4 885c4762a1bSJed Brown filter: grep -v "nls ksp" 886c4762a1bSJed Brown requires: !single 887c4762a1bSJed Brown 888c4762a1bSJed Brown test: 889*4936d809SStefano Zampini suffix: 2_snes 890*4936d809SStefano Zampini nsize: 2 891*4936d809SStefano Zampini args: -tao_smonitor -tao_type snes -ksp_converged_maxits -ksp_max_it 15 -snes_atol 1.e-4 892*4936d809SStefano Zampini filter: grep -v "nls ksp" 893*4936d809SStefano Zampini requires: !single 894*4936d809SStefano Zampini 895*4936d809SStefano Zampini test: 896c4762a1bSJed Brown suffix: 3 897c4762a1bSJed Brown nsize: 3 898c4762a1bSJed Brown args: -tao_smonitor -tao_type cg -tao_cg_type fr -mx 10 -my 10 -tao_gatol 1.e-3 899c4762a1bSJed Brown requires: !single 900c4762a1bSJed Brown 901c4762a1bSJed Brown test: 902*4936d809SStefano Zampini suffix: 3_snes 903*4936d809SStefano Zampini nsize: 3 904*4936d809SStefano Zampini args: -tao_smonitor -tao_type snes -snes_type ncg -snes_ncg_type fr -mx 10 -my 10 -snes_atol 1.e-4 905*4936d809SStefano Zampini requires: !single 906*4936d809SStefano Zampini 907*4936d809SStefano Zampini test: 908*4936d809SStefano Zampini suffix: 4_snes_ngmres 909*4936d809SStefano Zampini args: -tao_type snes -snes_type ngmres -npc_snes_type ncg -mx 10 -my 10 -snes_atol 1.e-4 -npc_snes_ncg_type fr -snes_converged_reason -snes_monitor ::ascii_info_detail -snes_ngmres_monitor -snes_ngmres_select_type {{linesearch difference}separate output} 910*4936d809SStefano Zampini requires: !single 911*4936d809SStefano Zampini 912*4936d809SStefano Zampini test: 913c4762a1bSJed Brown suffix: 5 914c4762a1bSJed Brown nsize: 2 915c4762a1bSJed Brown args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 916c4762a1bSJed Brown requires: !single 917c4762a1bSJed Brown 918c4762a1bSJed Brown TEST*/ 919