xref: /petsc/src/tao/unconstrained/tutorials/minsurf2.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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\
169f1eecfaSStefano Zampini   -da_grid_x <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
179f1eecfaSStefano Zampini   -da_grid_y <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
234936d809SStefano Zampini    application-provided call-back routines, FormFunction(),
244936d809SStefano 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);
374936d809SStefano Zampini static PetscErrorCode QuadraticH(AppCtx *, Vec, Mat);
384936d809SStefano Zampini static PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
394936d809SStefano Zampini static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
404936d809SStefano Zampini static PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
414936d809SStefano Zampini static PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
424936d809SStefano Zampini static PetscErrorCode My_Monitor(Tao, void *);
43c4762a1bSJed Brown 
main(int argc,char ** argv)44d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
45d71ae5a4SJacob Faibussowitsch {
46c4762a1bSJed Brown   Vec           x;                     /* solution, gradient vectors */
479f1eecfaSStefano Zampini   PetscBool     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;
56c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Create distributed array (DM) to manage parallel grid and vectors  */
599f1eecfaSStefano Zampini   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 10, 10, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &user.dm));
609566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
619566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
624936d809SStefano Zampini   PetscCall(DMDASetUniformCoordinates(user.dm, -0.5, 0.5, -0.5, 0.5, PETSC_DECIDE, PETSC_DECIDE));
639f1eecfaSStefano Zampini   PetscCall(DMDAGetInfo(user.dm, PETSC_IGNORE, &user.mx, &user.my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
649f1eecfaSStefano Zampini 
659f1eecfaSStefano Zampini   PetscCall(PetscPrintf(MPI_COMM_WORLD, "\n---- Minimum Surface Area Problem -----\n"));
669f1eecfaSStefano Zampini   PetscCall(PetscPrintf(MPI_COMM_WORLD, "mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n", user.mx, user.my));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Create TAO solver and set desired solution method.*/
699566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
709566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOCG));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /*
73f605775fSPierre Jolivet      Extract global vector from DA for the vector of variables --  PETSc routine
74c4762a1bSJed Brown      Compute the initial solution                              --  application specific, see below
75c4762a1bSJed Brown      Set this vector for use by TAO                            --  TAO routine
76c4762a1bSJed Brown   */
779566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm, &x));
789566063dSJacob Faibussowitsch   PetscCall(MSA_BoundaryConditions(&user));
799566063dSJacob Faibussowitsch   PetscCall(MSA_InitialPoint(&user, x));
809566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /*
83c4762a1bSJed Brown      Initialize the Application context for use in function evaluations  --  application specific, see below.
84c4762a1bSJed Brown      Set routines for function and gradient evaluation
85c4762a1bSJed Brown   */
864936d809SStefano Zampini   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&user));
879566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /*
90c4762a1bSJed Brown      Given the command line arguments, calculate the hessian with either the user-
91c4762a1bSJed Brown      provided function FormHessian, or the default finite-difference driven Hessian
92c4762a1bSJed Brown      functions
93c4762a1bSJed Brown   */
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fddefault", &fddefault));
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-tao_fdcoloring", &fdcoloring));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /*
98c4762a1bSJed Brown      Create a matrix data structure to store the Hessian and set
99d5b43468SJose E. Roman      the Hessian evaluation routine.
1000b4b7b1cSBarry Smith      Set the matrix nonzero structure to be used for Hessian evaluations
101c4762a1bSJed Brown   */
1029566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.dm, &user.H));
1039566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   if (fdcoloring) {
1069566063dSJacob Faibussowitsch     PetscCall(DMCreateColoring(user.dm, IS_COLORING_GLOBAL, &iscoloring));
1079566063dSJacob Faibussowitsch     PetscCall(MatFDColoringCreate(user.H, iscoloring, &matfdcoloring));
1089566063dSJacob Faibussowitsch     PetscCall(ISColoringDestroy(&iscoloring));
1092ba42892SBarry Smith     PetscCall(MatFDColoringSetFunction(matfdcoloring, (MatFDColoringFn *)FormGradient, (void *)&user));
1109566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1119566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessianColor, (void *)matfdcoloring));
112c4762a1bSJed Brown   } else if (fddefault) {
1139566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, TaoDefaultComputeHessian, (void *)NULL));
114c4762a1bSJed Brown   } else {
1159566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /*
119c4762a1bSJed Brown      If my_monitor option is in command line, then use the user-provided
120c4762a1bSJed Brown      monitoring function
121c4762a1bSJed Brown   */
1229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-my_monitor", &viewmat));
12310978b7dSBarry Smith   if (viewmat) PetscCall(TaoMonitorSet(tao, My_Monitor, NULL, NULL));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Check for any tao command line options */
1269566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1299566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
130c4762a1bSJed Brown 
1319566063dSJacob Faibussowitsch   PetscCall(TaoView(tao, PETSC_VIEWER_STDOUT_WORLD));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* Free TAO data structures */
1349566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* Free PETSc data structures */
1379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
13948a46eb9SPierre Jolivet   if (fdcoloring) PetscCall(MatFDColoringDestroy(&matfdcoloring));
1409566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.bottom));
1419566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.top));
1429566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.left));
1439566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.right));
1449566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
1459566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
146b122ec5aSJacob Faibussowitsch   return 0;
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
1494936d809SStefano Zampini /* -------------------------------------------------------------------- */
1504936d809SStefano Zampini /*
1514936d809SStefano Zampini     FormFunction - Evaluates the objective function.
1524936d809SStefano Zampini 
1534936d809SStefano Zampini     Input Parameters:
1544936d809SStefano Zampini .   tao     - the Tao context
1554936d809SStefano Zampini .   X       - input vector
1564936d809SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjective()
1574936d809SStefano Zampini 
1584936d809SStefano Zampini     Output Parameters:
1594936d809SStefano Zampini .   fcn     - the newly evaluated function
1604936d809SStefano Zampini */
FormFunction(Tao tao,Vec X,PetscReal * fcn,void * userCtx)1614936d809SStefano Zampini PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *fcn, void *userCtx)
162d71ae5a4SJacob Faibussowitsch {
1634936d809SStefano Zampini   AppCtx     *user = (AppCtx *)userCtx;
1644936d809SStefano Zampini   PetscInt    i, j;
1654936d809SStefano Zampini   PetscInt    mx = user->mx, my = user->my;
1664936d809SStefano Zampini   PetscInt    xs, xm, gxs, gxm, ys, ym, gys, gym;
1674936d809SStefano Zampini   PetscReal   ft = 0.0;
1684936d809SStefano Zampini   PetscReal   hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), area = 0.5 * hx * hy;
1694936d809SStefano Zampini   PetscReal   rhx = mx + 1, rhy = my + 1;
1704936d809SStefano Zampini   PetscReal   f2, f4, d1, d2, d3, d4, xc, xl, xr, xt, xb;
1714936d809SStefano Zampini   PetscReal **x;
1724936d809SStefano Zampini   Vec         localX;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   PetscFunctionBegin;
1754936d809SStefano Zampini   /* Get local mesh boundaries */
1764936d809SStefano Zampini   PetscCall(DMGetLocalVector(user->dm, &localX));
1774936d809SStefano Zampini   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
1784936d809SStefano Zampini   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
1794936d809SStefano Zampini 
1804936d809SStefano Zampini   /* Scatter ghost points to local vector */
1814936d809SStefano Zampini   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
1824936d809SStefano Zampini   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
1834936d809SStefano Zampini 
1844936d809SStefano Zampini   /* Get pointers to vector data */
1854936d809SStefano Zampini   PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x));
1864936d809SStefano Zampini 
1874936d809SStefano Zampini   /* Compute function and gradient over the locally owned part of the mesh */
1884936d809SStefano Zampini   for (j = ys; j < ys + ym; j++) {
1894936d809SStefano Zampini     for (i = xs; i < xs + xm; i++) {
1904936d809SStefano Zampini       xc = x[j][i];
1914936d809SStefano Zampini 
1924936d809SStefano Zampini       if (i == 0) { /* left side */
1934936d809SStefano Zampini         xl = user->left[j - ys + 1];
1944936d809SStefano Zampini       } else {
1954936d809SStefano Zampini         xl = x[j][i - 1];
1964936d809SStefano Zampini       }
1974936d809SStefano Zampini 
1984936d809SStefano Zampini       if (j == 0) { /* bottom side */
1994936d809SStefano Zampini         xb = user->bottom[i - xs + 1];
2004936d809SStefano Zampini       } else {
2014936d809SStefano Zampini         xb = x[j - 1][i];
2024936d809SStefano Zampini       }
2034936d809SStefano Zampini 
2044936d809SStefano Zampini       if (i + 1 == gxs + gxm) { /* right side */
2054936d809SStefano Zampini         xr = user->right[j - ys + 1];
2064936d809SStefano Zampini       } else {
2074936d809SStefano Zampini         xr = x[j][i + 1];
2084936d809SStefano Zampini       }
2094936d809SStefano Zampini 
2104936d809SStefano Zampini       if (j + 1 == gys + gym) { /* top side */
2114936d809SStefano Zampini         xt = user->top[i - xs + 1];
2124936d809SStefano Zampini       } else {
2134936d809SStefano Zampini         xt = x[j + 1][i];
2144936d809SStefano Zampini       }
2154936d809SStefano Zampini 
2164936d809SStefano Zampini       d1 = (xc - xl);
2174936d809SStefano Zampini       d2 = (xc - xr);
2184936d809SStefano Zampini       d3 = (xc - xt);
2194936d809SStefano Zampini       d4 = (xc - xb);
2204936d809SStefano Zampini 
2214936d809SStefano Zampini       d1 *= rhx;
2224936d809SStefano Zampini       d2 *= rhx;
2234936d809SStefano Zampini       d3 *= rhy;
2244936d809SStefano Zampini       d4 *= rhy;
2254936d809SStefano Zampini 
2264936d809SStefano Zampini       f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
2274936d809SStefano Zampini       f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
2284936d809SStefano Zampini 
2294936d809SStefano Zampini       ft = ft + (f2 + f4);
2304936d809SStefano Zampini     }
2314936d809SStefano Zampini   }
2324936d809SStefano Zampini 
2334936d809SStefano Zampini   /* Compute triangular areas along the border of the domain. */
2344936d809SStefano Zampini   if (xs == 0) { /* left side */
2354936d809SStefano Zampini     for (j = ys; j < ys + ym; j++) {
2364936d809SStefano Zampini       d3 = (user->left[j - ys + 1] - user->left[j - ys + 2]) * rhy;
2374936d809SStefano Zampini       d2 = (user->left[j - ys + 1] - x[j][0]) * rhx;
2384936d809SStefano Zampini       ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
2394936d809SStefano Zampini     }
2404936d809SStefano Zampini   }
2414936d809SStefano Zampini   if (ys == 0) { /* bottom side */
2424936d809SStefano Zampini     for (i = xs; i < xs + xm; i++) {
2434936d809SStefano Zampini       d2 = (user->bottom[i + 1 - xs] - user->bottom[i - xs + 2]) * rhx;
2444936d809SStefano Zampini       d3 = (user->bottom[i - xs + 1] - x[0][i]) * rhy;
2454936d809SStefano Zampini       ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
2464936d809SStefano Zampini     }
2474936d809SStefano Zampini   }
2484936d809SStefano Zampini   if (xs + xm == mx) { /* right side */
2494936d809SStefano Zampini     for (j = ys; j < ys + ym; j++) {
2504936d809SStefano Zampini       d1 = (x[j][mx - 1] - user->right[j - ys + 1]) * rhx;
2514936d809SStefano Zampini       d4 = (user->right[j - ys] - user->right[j - ys + 1]) * rhy;
2524936d809SStefano Zampini       ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
2534936d809SStefano Zampini     }
2544936d809SStefano Zampini   }
2554936d809SStefano Zampini   if (ys + ym == my) { /* top side */
2564936d809SStefano Zampini     for (i = xs; i < xs + xm; i++) {
2574936d809SStefano Zampini       d1 = (x[my - 1][i] - user->top[i - xs + 1]) * rhy;
2584936d809SStefano Zampini       d4 = (user->top[i - xs + 1] - user->top[i - xs]) * rhx;
2594936d809SStefano Zampini       ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
2604936d809SStefano Zampini     }
2614936d809SStefano Zampini   }
2624936d809SStefano Zampini   if (ys == 0 && xs == 0) {
2634936d809SStefano Zampini     d1 = (user->left[0] - user->left[1]) / hy;
2644936d809SStefano Zampini     d2 = (user->bottom[0] - user->bottom[1]) * rhx;
2654936d809SStefano Zampini     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
2664936d809SStefano Zampini   }
2674936d809SStefano Zampini   if (ys + ym == my && xs + xm == mx) {
2684936d809SStefano Zampini     d1 = (user->right[ym + 1] - user->right[ym]) * rhy;
2694936d809SStefano Zampini     d2 = (user->top[xm + 1] - user->top[xm]) * rhx;
2704936d809SStefano Zampini     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
2714936d809SStefano Zampini   }
2724936d809SStefano Zampini 
2734936d809SStefano Zampini   ft = ft * area;
274462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
2754936d809SStefano Zampini 
2764936d809SStefano Zampini   /* Restore vectors */
2774936d809SStefano Zampini   PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x));
2784936d809SStefano Zampini   PetscCall(DMRestoreLocalVector(user->dm, &localX));
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown /* -------------------------------------------------------------------- */
2834936d809SStefano Zampini /*
2844936d809SStefano Zampini     FormFunctionGradient - Evaluates the function and corresponding gradient.
285c4762a1bSJed Brown 
286c4762a1bSJed Brown     Input Parameters:
287c4762a1bSJed Brown .   tao     - the Tao context
2884936d809SStefano Zampini .   X      - input vector
289a82e8c82SStefano Zampini .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
290c4762a1bSJed Brown 
291c4762a1bSJed Brown     Output Parameters:
292c4762a1bSJed Brown .   fcn     - the newly evaluated function
2934936d809SStefano Zampini .   G       - vector containing the newly evaluated gradient
294c4762a1bSJed Brown */
FormFunctionGradient(Tao tao,Vec X,PetscReal * fcn,Vec G,void * userCtx)295d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx)
296d71ae5a4SJacob Faibussowitsch {
297c4762a1bSJed Brown   AppCtx     *user = (AppCtx *)userCtx;
298c4762a1bSJed Brown   PetscInt    i, j;
299c4762a1bSJed Brown   PetscInt    mx = user->mx, my = user->my;
300c4762a1bSJed Brown   PetscInt    xs, xm, gxs, gxm, ys, ym, gys, gym;
301c4762a1bSJed Brown   PetscReal   ft = 0.0;
302c4762a1bSJed Brown   PetscReal   hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy;
303c4762a1bSJed Brown   PetscReal   rhx = mx + 1, rhy = my + 1;
304c4762a1bSJed Brown   PetscReal   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
305c4762a1bSJed Brown   PetscReal   df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
306c4762a1bSJed Brown   PetscReal **g, **x;
307c4762a1bSJed Brown   Vec         localX;
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   PetscFunctionBegin;
310c4762a1bSJed Brown   /* Get local mesh boundaries */
3119566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->dm, &localX));
3129566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
3139566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   /* Scatter ghost points to local vector */
3169566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
3179566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   /* Get pointers to vector data */
3209566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x));
3219566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->dm, G, (void **)&g));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /* Compute function and gradient over the locally owned part of the mesh */
324c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
325c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
326c4762a1bSJed Brown       xc  = x[j][i];
327c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
328c4762a1bSJed Brown 
329c4762a1bSJed Brown       if (i == 0) { /* left side */
330c4762a1bSJed Brown         xl  = user->left[j - ys + 1];
331c4762a1bSJed Brown         xlt = user->left[j - ys + 2];
332c4762a1bSJed Brown       } else {
333c4762a1bSJed Brown         xl = x[j][i - 1];
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown       if (j == 0) { /* bottom side */
337c4762a1bSJed Brown         xb  = user->bottom[i - xs + 1];
338c4762a1bSJed Brown         xrb = user->bottom[i - xs + 2];
339c4762a1bSJed Brown       } else {
340c4762a1bSJed Brown         xb = x[j - 1][i];
341c4762a1bSJed Brown       }
342c4762a1bSJed Brown 
343c4762a1bSJed Brown       if (i + 1 == gxs + gxm) { /* right side */
344c4762a1bSJed Brown         xr  = user->right[j - ys + 1];
345c4762a1bSJed Brown         xrb = user->right[j - ys];
346c4762a1bSJed Brown       } else {
347c4762a1bSJed Brown         xr = x[j][i + 1];
348c4762a1bSJed Brown       }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown       if (j + 1 == gys + gym) { /* top side */
351c4762a1bSJed Brown         xt  = user->top[i - xs + 1];
352c4762a1bSJed Brown         xlt = user->top[i - xs];
353c4762a1bSJed Brown       } else {
354c4762a1bSJed Brown         xt = x[j + 1][i];
355c4762a1bSJed Brown       }
356c4762a1bSJed Brown 
357ad540459SPierre Jolivet       if (i > gxs && j + 1 < gys + gym) xlt = x[j + 1][i - 1];
358ad540459SPierre Jolivet       if (j > gys && i + 1 < gxs + gxm) xrb = x[j - 1][i + 1];
359c4762a1bSJed Brown 
360c4762a1bSJed Brown       d1 = (xc - xl);
361c4762a1bSJed Brown       d2 = (xc - xr);
362c4762a1bSJed Brown       d3 = (xc - xt);
363c4762a1bSJed Brown       d4 = (xc - xb);
364c4762a1bSJed Brown       d5 = (xr - xrb);
365c4762a1bSJed Brown       d6 = (xrb - xb);
366c4762a1bSJed Brown       d7 = (xlt - xl);
367c4762a1bSJed Brown       d8 = (xt - xlt);
368c4762a1bSJed Brown 
369c4762a1bSJed Brown       df1dxc = d1 * hydhx;
370c4762a1bSJed Brown       df2dxc = (d1 * hydhx + d4 * hxdhy);
371c4762a1bSJed Brown       df3dxc = d3 * hxdhy;
372c4762a1bSJed Brown       df4dxc = (d2 * hydhx + d3 * hxdhy);
373c4762a1bSJed Brown       df5dxc = d2 * hydhx;
374c4762a1bSJed Brown       df6dxc = d4 * hxdhy;
375c4762a1bSJed Brown 
376c4762a1bSJed Brown       d1 *= rhx;
377c4762a1bSJed Brown       d2 *= rhx;
378c4762a1bSJed Brown       d3 *= rhy;
379c4762a1bSJed Brown       d4 *= rhy;
380c4762a1bSJed Brown       d5 *= rhy;
381c4762a1bSJed Brown       d6 *= rhx;
382c4762a1bSJed Brown       d7 *= rhy;
383c4762a1bSJed Brown       d8 *= rhx;
384c4762a1bSJed Brown 
385c4762a1bSJed Brown       f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7);
386c4762a1bSJed Brown       f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
387c4762a1bSJed Brown       f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8);
388c4762a1bSJed Brown       f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
389c4762a1bSJed Brown       f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5);
390c4762a1bSJed Brown       f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6);
391c4762a1bSJed Brown 
392c4762a1bSJed Brown       ft = ft + (f2 + f4);
393c4762a1bSJed Brown 
394c4762a1bSJed Brown       df1dxc /= f1;
395c4762a1bSJed Brown       df2dxc /= f2;
396c4762a1bSJed Brown       df3dxc /= f3;
397c4762a1bSJed Brown       df4dxc /= f4;
398c4762a1bSJed Brown       df5dxc /= f5;
399c4762a1bSJed Brown       df6dxc /= f6;
400c4762a1bSJed Brown 
401c4762a1bSJed Brown       g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5;
402c4762a1bSJed Brown     }
403c4762a1bSJed Brown   }
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   /* Compute triangular areas along the border of the domain. */
406c4762a1bSJed Brown   if (xs == 0) { /* left side */
407c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
408c4762a1bSJed Brown       d3 = (user->left[j - ys + 1] - user->left[j - ys + 2]) * rhy;
409c4762a1bSJed Brown       d2 = (user->left[j - ys + 1] - x[j][0]) * rhx;
410c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
411c4762a1bSJed Brown     }
412c4762a1bSJed Brown   }
413c4762a1bSJed Brown   if (ys == 0) { /* bottom side */
414c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
415c4762a1bSJed Brown       d2 = (user->bottom[i + 1 - xs] - user->bottom[i - xs + 2]) * rhx;
416c4762a1bSJed Brown       d3 = (user->bottom[i - xs + 1] - x[0][i]) * rhy;
417c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
418c4762a1bSJed Brown     }
419c4762a1bSJed Brown   }
420c4762a1bSJed Brown 
421c4762a1bSJed Brown   if (xs + xm == mx) { /* right side */
422c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
423c4762a1bSJed Brown       d1 = (x[j][mx - 1] - user->right[j - ys + 1]) * rhx;
424c4762a1bSJed Brown       d4 = (user->right[j - ys] - user->right[j - ys + 1]) * rhy;
425c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
426c4762a1bSJed Brown     }
427c4762a1bSJed Brown   }
428c4762a1bSJed Brown   if (ys + ym == my) { /* top side */
429c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
430c4762a1bSJed Brown       d1 = (x[my - 1][i] - user->top[i - xs + 1]) * rhy;
431c4762a1bSJed Brown       d4 = (user->top[i - xs + 1] - user->top[i - xs]) * rhx;
432c4762a1bSJed Brown       ft = ft + PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
433c4762a1bSJed Brown     }
434c4762a1bSJed Brown   }
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   if (ys == 0 && xs == 0) {
437c4762a1bSJed Brown     d1 = (user->left[0] - user->left[1]) / hy;
438c4762a1bSJed Brown     d2 = (user->bottom[0] - user->bottom[1]) * rhx;
439c4762a1bSJed Brown     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
440c4762a1bSJed Brown   }
441c4762a1bSJed Brown   if (ys + ym == my && xs + xm == mx) {
442c4762a1bSJed Brown     d1 = (user->right[ym + 1] - user->right[ym]) * rhy;
443c4762a1bSJed Brown     d2 = (user->top[xm + 1] - user->top[xm]) * rhx;
444c4762a1bSJed Brown     ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
445c4762a1bSJed Brown   }
446c4762a1bSJed Brown 
447c4762a1bSJed Brown   ft = ft * area;
448462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD));
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   /* Restore vectors */
4519566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x));
4529566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->dm, G, (void **)&g));
4539566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->dm, &localX));
4549566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(67.0 * xm * ym));
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
456c4762a1bSJed Brown }
457c4762a1bSJed Brown 
FormGradient(Tao tao,Vec X,Vec G,void * userCtx)4584936d809SStefano Zampini PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *userCtx)
4594936d809SStefano Zampini {
4604936d809SStefano Zampini   PetscReal fcn;
4614936d809SStefano Zampini 
4624936d809SStefano Zampini   PetscFunctionBegin;
4634936d809SStefano Zampini   PetscCall(FormFunctionGradient(tao, X, &fcn, G, userCtx));
4644936d809SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
4654936d809SStefano Zampini }
4664936d809SStefano Zampini 
467c4762a1bSJed Brown /* ------------------------------------------------------------------- */
468c4762a1bSJed Brown /*
469c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
470c4762a1bSJed Brown 
471c4762a1bSJed Brown    Input Parameters:
472c4762a1bSJed Brown .  tao  - the Tao context
473c4762a1bSJed Brown .  x    - input vector
474a82e8c82SStefano Zampini .  ptr  - optional user-defined context, as set by TaoSetHessian()
475c4762a1bSJed Brown 
476c4762a1bSJed Brown    Output Parameters:
477c4762a1bSJed Brown .  H    - Hessian matrix
4787addb90fSBarry Smith .  Hpre - optionally different matrix used to compute the preconditioner
479c4762a1bSJed Brown 
480c4762a1bSJed Brown */
FormHessian(Tao tao,Vec X,Mat H,Mat Hpre,void * ptr)481d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
482d71ae5a4SJacob Faibussowitsch {
483c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ptr;
484c4762a1bSJed Brown 
485c4762a1bSJed Brown   PetscFunctionBegin;
486c4762a1bSJed Brown   /* Evaluate the Hessian entries*/
4879566063dSJacob Faibussowitsch   PetscCall(QuadraticH(user, X, H));
4883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
489c4762a1bSJed Brown }
490c4762a1bSJed Brown 
491c4762a1bSJed Brown /* ------------------------------------------------------------------- */
492c4762a1bSJed Brown /*
493c4762a1bSJed Brown    QuadraticH - Evaluates Hessian matrix.
494c4762a1bSJed Brown 
495c4762a1bSJed Brown    Input Parameters:
496c4762a1bSJed Brown .  user - user-defined context, as set by TaoSetHessian()
497c4762a1bSJed Brown .  X    - input vector
498c4762a1bSJed Brown 
499c4762a1bSJed Brown    Output Parameter:
500c4762a1bSJed Brown .  H    - Hessian matrix
501c4762a1bSJed Brown */
QuadraticH(AppCtx * user,Vec X,Mat Hessian)502d71ae5a4SJacob Faibussowitsch PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian)
503d71ae5a4SJacob Faibussowitsch {
504c4762a1bSJed Brown   PetscInt    i, j, k;
505c4762a1bSJed Brown   PetscInt    mx = user->mx, my = user->my;
506c4762a1bSJed Brown   PetscInt    xs, xm, gxs, gxm, ys, ym, gys, gym;
507c4762a1bSJed Brown   PetscReal   hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
508c4762a1bSJed Brown   PetscReal   f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
509c4762a1bSJed Brown   PetscReal   hl, hr, ht, hb, hc, htl, hbr;
510c4762a1bSJed Brown   PetscReal **x, v[7];
511c4762a1bSJed Brown   MatStencil  col[7], row;
512c4762a1bSJed Brown   Vec         localX;
513c4762a1bSJed Brown   PetscBool   assembled;
514c4762a1bSJed Brown 
515c4762a1bSJed Brown   PetscFunctionBegin;
516c4762a1bSJed Brown   /* Get local mesh boundaries */
5179566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->dm, &localX));
518c4762a1bSJed Brown 
5199566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
5209566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
521c4762a1bSJed Brown 
522c4762a1bSJed Brown   /* Scatter ghost points to local vector */
5239566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX));
5249566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX));
525c4762a1bSJed Brown 
526c4762a1bSJed Brown   /* Get pointers to vector data */
5279566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->dm, localX, (void **)&x));
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   /* Initialize matrix entries to zero */
5309566063dSJacob Faibussowitsch   PetscCall(MatAssembled(Hessian, &assembled));
5319566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(Hessian));
532c4762a1bSJed Brown 
533c4762a1bSJed Brown   /* Set various matrix options */
5349566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
535c4762a1bSJed Brown 
536c4762a1bSJed Brown   /* Compute Hessian over the locally owned part of the mesh */
537c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
538c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
539c4762a1bSJed Brown       xc  = x[j][i];
540c4762a1bSJed Brown       xlt = xrb = xl = xr = xb = xt = xc;
541c4762a1bSJed Brown 
542c4762a1bSJed Brown       /* Left side */
543c4762a1bSJed Brown       if (i == 0) {
544c4762a1bSJed Brown         xl  = user->left[j - ys + 1];
545c4762a1bSJed Brown         xlt = user->left[j - ys + 2];
546c4762a1bSJed Brown       } else {
547c4762a1bSJed Brown         xl = x[j][i - 1];
548c4762a1bSJed Brown       }
549c4762a1bSJed Brown 
550c4762a1bSJed Brown       if (j == 0) {
551c4762a1bSJed Brown         xb  = user->bottom[i - xs + 1];
552c4762a1bSJed Brown         xrb = user->bottom[i - xs + 2];
553c4762a1bSJed Brown       } else {
554c4762a1bSJed Brown         xb = x[j - 1][i];
555c4762a1bSJed Brown       }
556c4762a1bSJed Brown 
557c4762a1bSJed Brown       if (i + 1 == mx) {
558c4762a1bSJed Brown         xr  = user->right[j - ys + 1];
559c4762a1bSJed Brown         xrb = user->right[j - ys];
560c4762a1bSJed Brown       } else {
561c4762a1bSJed Brown         xr = x[j][i + 1];
562c4762a1bSJed Brown       }
563c4762a1bSJed Brown 
564c4762a1bSJed Brown       if (j + 1 == my) {
565c4762a1bSJed Brown         xt  = user->top[i - xs + 1];
566c4762a1bSJed Brown         xlt = user->top[i - xs];
567c4762a1bSJed Brown       } else {
568c4762a1bSJed Brown         xt = x[j + 1][i];
569c4762a1bSJed Brown       }
570c4762a1bSJed Brown 
571ad540459SPierre Jolivet       if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1];
572ad540459SPierre Jolivet       if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1];
573c4762a1bSJed Brown 
574c4762a1bSJed Brown       d1 = (xc - xl) / hx;
575c4762a1bSJed Brown       d2 = (xc - xr) / hx;
576c4762a1bSJed Brown       d3 = (xc - xt) / hy;
577c4762a1bSJed Brown       d4 = (xc - xb) / hy;
578c4762a1bSJed Brown       d5 = (xrb - xr) / hy;
579c4762a1bSJed Brown       d6 = (xrb - xb) / hx;
580c4762a1bSJed Brown       d7 = (xlt - xl) / hy;
581c4762a1bSJed Brown       d8 = (xlt - xt) / hx;
582c4762a1bSJed Brown 
583c4762a1bSJed Brown       f1 = PetscSqrtReal(1.0 + d1 * d1 + d7 * d7);
584c4762a1bSJed Brown       f2 = PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
585c4762a1bSJed Brown       f3 = PetscSqrtReal(1.0 + d3 * d3 + d8 * d8);
586c4762a1bSJed Brown       f4 = PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
587c4762a1bSJed Brown       f5 = PetscSqrtReal(1.0 + d2 * d2 + d5 * d5);
588c4762a1bSJed Brown       f6 = PetscSqrtReal(1.0 + d4 * d4 + d6 * d6);
589c4762a1bSJed Brown 
5909371c9d4SSatish Balay       hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
5919371c9d4SSatish Balay       hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
5929371c9d4SSatish Balay       ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
5939371c9d4SSatish Balay       hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
594c4762a1bSJed Brown 
595c4762a1bSJed Brown       hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
596c4762a1bSJed Brown       htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
597c4762a1bSJed Brown 
5989371c9d4SSatish 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);
599c4762a1bSJed Brown 
6009371c9d4SSatish Balay       hl /= 2.0;
6019371c9d4SSatish Balay       hr /= 2.0;
6029371c9d4SSatish Balay       ht /= 2.0;
6039371c9d4SSatish Balay       hb /= 2.0;
6049371c9d4SSatish Balay       hbr /= 2.0;
6059371c9d4SSatish Balay       htl /= 2.0;
6069371c9d4SSatish Balay       hc /= 2.0;
607c4762a1bSJed Brown 
6089371c9d4SSatish Balay       row.j = j;
6099371c9d4SSatish Balay       row.i = i;
610c4762a1bSJed Brown       k     = 0;
611c4762a1bSJed Brown       if (j > 0) {
612c4762a1bSJed Brown         v[k]     = hb;
6139371c9d4SSatish Balay         col[k].j = j - 1;
6149371c9d4SSatish Balay         col[k].i = i;
615c4762a1bSJed Brown         k++;
616c4762a1bSJed Brown       }
617c4762a1bSJed Brown 
618c4762a1bSJed Brown       if (j > 0 && i < mx - 1) {
619c4762a1bSJed Brown         v[k]     = hbr;
6209371c9d4SSatish Balay         col[k].j = j - 1;
6219371c9d4SSatish Balay         col[k].i = i + 1;
622c4762a1bSJed Brown         k++;
623c4762a1bSJed Brown       }
624c4762a1bSJed Brown 
625c4762a1bSJed Brown       if (i > 0) {
626c4762a1bSJed Brown         v[k]     = hl;
6279371c9d4SSatish Balay         col[k].j = j;
6289371c9d4SSatish Balay         col[k].i = i - 1;
629c4762a1bSJed Brown         k++;
630c4762a1bSJed Brown       }
631c4762a1bSJed Brown 
632c4762a1bSJed Brown       v[k]     = hc;
6339371c9d4SSatish Balay       col[k].j = j;
6349371c9d4SSatish Balay       col[k].i = i;
635c4762a1bSJed Brown       k++;
636c4762a1bSJed Brown 
637c4762a1bSJed Brown       if (i < mx - 1) {
638c4762a1bSJed Brown         v[k]     = hr;
6399371c9d4SSatish Balay         col[k].j = j;
6409371c9d4SSatish Balay         col[k].i = i + 1;
641c4762a1bSJed Brown         k++;
642c4762a1bSJed Brown       }
643c4762a1bSJed Brown 
644c4762a1bSJed Brown       if (i > 0 && j < my - 1) {
645c4762a1bSJed Brown         v[k]     = htl;
6469371c9d4SSatish Balay         col[k].j = j + 1;
6479371c9d4SSatish Balay         col[k].i = i - 1;
648c4762a1bSJed Brown         k++;
649c4762a1bSJed Brown       }
650c4762a1bSJed Brown 
651c4762a1bSJed Brown       if (j < my - 1) {
652c4762a1bSJed Brown         v[k]     = ht;
6539371c9d4SSatish Balay         col[k].j = j + 1;
6549371c9d4SSatish Balay         col[k].i = i;
655c4762a1bSJed Brown         k++;
656c4762a1bSJed Brown       }
657c4762a1bSJed Brown 
658c4762a1bSJed Brown       /*
659c4762a1bSJed Brown          Set matrix values using local numbering, which was defined
660c4762a1bSJed Brown          earlier, in the main routine.
661c4762a1bSJed Brown       */
6629566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(Hessian, 1, &row, k, col, v, INSERT_VALUES));
663c4762a1bSJed Brown     }
664c4762a1bSJed Brown   }
665c4762a1bSJed Brown 
6669566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->dm, localX, (void **)&x));
6679566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->dm, &localX));
668c4762a1bSJed Brown 
6699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY));
6709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY));
671c4762a1bSJed Brown 
6729566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(199.0 * xm * ym));
6733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
674c4762a1bSJed Brown }
675c4762a1bSJed Brown 
676c4762a1bSJed Brown /* ------------------------------------------------------------------- */
677c4762a1bSJed Brown /*
678c4762a1bSJed Brown    MSA_BoundaryConditions -  Calculates the boundary conditions for
679c4762a1bSJed Brown    the region.
680c4762a1bSJed Brown 
681c4762a1bSJed Brown    Input Parameter:
682c4762a1bSJed Brown .  user - user-defined application context
683c4762a1bSJed Brown 
684c4762a1bSJed Brown    Output Parameter:
685c4762a1bSJed Brown .  user - user-defined application context
686c4762a1bSJed Brown */
MSA_BoundaryConditions(AppCtx * user)687d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user)
688d71ae5a4SJacob Faibussowitsch {
689c4762a1bSJed Brown   PetscInt   i, j, k, limit = 0, maxits = 5;
690c4762a1bSJed Brown   PetscInt   xs, ys, xm, ym, gxs, gys, gxm, gym;
691c4762a1bSJed Brown   PetscInt   mx = user->mx, my = user->my;
692c4762a1bSJed Brown   PetscInt   bsize = 0, lsize = 0, tsize = 0, rsize = 0;
693c4762a1bSJed Brown   PetscReal  one = 1.0, two = 2.0, three = 3.0, tol = 1e-10;
694c4762a1bSJed Brown   PetscReal  fnorm, det, hx, hy, xt = 0, yt = 0;
695c4762a1bSJed Brown   PetscReal  u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
696c4762a1bSJed Brown   PetscReal  b = -0.5, t = 0.5, l = -0.5, r = 0.5;
697c4762a1bSJed Brown   PetscReal *boundary;
698c4762a1bSJed Brown   PetscBool  flg;
699c4762a1bSJed Brown 
700c4762a1bSJed Brown   PetscFunctionBegin;
701c4762a1bSJed Brown   /* Get local mesh boundaries */
7029566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
7039566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL));
704c4762a1bSJed Brown 
705c4762a1bSJed Brown   bsize = xm + 2;
706c4762a1bSJed Brown   lsize = ym + 2;
707c4762a1bSJed Brown   rsize = ym + 2;
708c4762a1bSJed Brown   tsize = xm + 2;
709c4762a1bSJed Brown 
7109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bsize, &user->bottom));
7119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tsize, &user->top));
7129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lsize, &user->left));
7139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(rsize, &user->right));
714c4762a1bSJed Brown 
7159371c9d4SSatish Balay   hx = (r - l) / (mx + 1);
7169371c9d4SSatish Balay   hy = (t - b) / (my + 1);
717c4762a1bSJed Brown 
718c4762a1bSJed Brown   for (j = 0; j < 4; j++) {
719c4762a1bSJed Brown     if (j == 0) {
720c4762a1bSJed Brown       yt       = b;
721c4762a1bSJed Brown       xt       = l + hx * xs;
722c4762a1bSJed Brown       limit    = bsize;
723c4762a1bSJed Brown       boundary = user->bottom;
724c4762a1bSJed Brown     } else if (j == 1) {
725c4762a1bSJed Brown       yt       = t;
726c4762a1bSJed Brown       xt       = l + hx * xs;
727c4762a1bSJed Brown       limit    = tsize;
728c4762a1bSJed Brown       boundary = user->top;
729c4762a1bSJed Brown     } else if (j == 2) {
730c4762a1bSJed Brown       yt       = b + hy * ys;
731c4762a1bSJed Brown       xt       = l;
732c4762a1bSJed Brown       limit    = lsize;
733c4762a1bSJed Brown       boundary = user->left;
734c4762a1bSJed Brown     } else { /* if (j==3) */
735c4762a1bSJed Brown       yt       = b + hy * ys;
736c4762a1bSJed Brown       xt       = r;
737c4762a1bSJed Brown       limit    = rsize;
738c4762a1bSJed Brown       boundary = user->right;
739c4762a1bSJed Brown     }
740c4762a1bSJed Brown 
741c4762a1bSJed Brown     for (i = 0; i < limit; i++) {
742c4762a1bSJed Brown       u1 = xt;
743c4762a1bSJed Brown       u2 = -yt;
744c4762a1bSJed Brown       for (k = 0; k < maxits; k++) {
745c4762a1bSJed Brown         nf1   = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
746c4762a1bSJed Brown         nf2   = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
747c4762a1bSJed Brown         fnorm = PetscSqrtReal(nf1 * nf1 + nf2 * nf2);
748c4762a1bSJed Brown         if (fnorm <= tol) break;
749c4762a1bSJed Brown         njac11 = one + u2 * u2 - u1 * u1;
750c4762a1bSJed Brown         njac12 = two * u1 * u2;
751c4762a1bSJed Brown         njac21 = -two * u1 * u2;
752c4762a1bSJed Brown         njac22 = -one - u1 * u1 + u2 * u2;
753c4762a1bSJed Brown         det    = njac11 * njac22 - njac21 * njac12;
754c4762a1bSJed Brown         u1     = u1 - (njac22 * nf1 - njac12 * nf2) / det;
755c4762a1bSJed Brown         u2     = u2 - (njac11 * nf2 - njac21 * nf1) / det;
756c4762a1bSJed Brown       }
757c4762a1bSJed Brown 
758c4762a1bSJed Brown       boundary[i] = u1 * u1 - u2 * u2;
759c4762a1bSJed Brown       if (j == 0 || j == 1) {
760c4762a1bSJed Brown         xt = xt + hx;
761c4762a1bSJed Brown       } else { /*  if (j==2 || j==3) */
762c4762a1bSJed Brown         yt = yt + hy;
763c4762a1bSJed Brown       }
764c4762a1bSJed Brown     }
765c4762a1bSJed Brown   }
766c4762a1bSJed Brown 
767c4762a1bSJed Brown   /* Scale the boundary if desired */
768c4762a1bSJed Brown   if (1 == 1) {
769c4762a1bSJed Brown     PetscReal scl = 1.0;
770c4762a1bSJed Brown 
7719566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg));
772c4762a1bSJed Brown     if (flg) {
773c4762a1bSJed Brown       for (i = 0; i < bsize; i++) user->bottom[i] *= scl;
774c4762a1bSJed Brown     }
775c4762a1bSJed Brown 
7769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg));
777c4762a1bSJed Brown     if (flg) {
778c4762a1bSJed Brown       for (i = 0; i < tsize; i++) user->top[i] *= scl;
779c4762a1bSJed Brown     }
780c4762a1bSJed Brown 
7819566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg));
782c4762a1bSJed Brown     if (flg) {
783c4762a1bSJed Brown       for (i = 0; i < rsize; i++) user->right[i] *= scl;
784c4762a1bSJed Brown     }
785c4762a1bSJed Brown 
7869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg));
787c4762a1bSJed Brown     if (flg) {
788c4762a1bSJed Brown       for (i = 0; i < lsize; i++) user->left[i] *= scl;
789c4762a1bSJed Brown     }
790c4762a1bSJed Brown   }
7913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
792c4762a1bSJed Brown }
793c4762a1bSJed Brown 
794c4762a1bSJed Brown /* ------------------------------------------------------------------- */
795c4762a1bSJed Brown /*
796c4762a1bSJed Brown    MSA_InitialPoint - Calculates the initial guess in one of three ways.
797c4762a1bSJed Brown 
798c4762a1bSJed Brown    Input Parameters:
799c4762a1bSJed Brown .  user - user-defined application context
800c4762a1bSJed Brown .  X - vector for initial guess
801c4762a1bSJed Brown 
802c4762a1bSJed Brown    Output Parameters:
803c4762a1bSJed Brown .  X - newly computed initial guess
804c4762a1bSJed Brown */
MSA_InitialPoint(AppCtx * user,Vec X)805d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X)
806d71ae5a4SJacob Faibussowitsch {
807c4762a1bSJed Brown   PetscInt  start2 = -1, i, j;
808c4762a1bSJed Brown   PetscReal start1 = 0;
809c4762a1bSJed Brown   PetscBool flg1, flg2;
810c4762a1bSJed Brown 
811c4762a1bSJed Brown   PetscFunctionBegin;
8129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-start", &start1, &flg1));
8139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-random", &start2, &flg2));
814c4762a1bSJed Brown 
815c4762a1bSJed Brown   if (flg1) { /* The zero vector is reasonable */
8169566063dSJacob Faibussowitsch     PetscCall(VecSet(X, start1));
817c4762a1bSJed Brown   } else if (flg2 && start2 > 0) { /* Try a random start between -0.5 and 0.5 */
8189371c9d4SSatish Balay     PetscRandom rctx;
8199371c9d4SSatish Balay     PetscReal   np5 = -0.5;
820c4762a1bSJed Brown 
8219566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
8229566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(X, rctx));
8239566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rctx));
8249566063dSJacob Faibussowitsch     PetscCall(VecShift(X, np5));
825c4762a1bSJed Brown   } else { /* Take an average of the boundary conditions */
826c4762a1bSJed Brown     PetscInt    xs, xm, ys, ym;
827c4762a1bSJed Brown     PetscInt    mx = user->mx, my = user->my;
828c4762a1bSJed Brown     PetscReal **x;
829c4762a1bSJed Brown 
830c4762a1bSJed Brown     /* Get local mesh boundaries */
8319566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL));
832c4762a1bSJed Brown 
833c4762a1bSJed Brown     /* Get pointers to vector data */
8349566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(user->dm, X, (void **)&x));
835c4762a1bSJed Brown 
836c4762a1bSJed Brown     /* Perform local computations */
837c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
838ad540459SPierre 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;
839c4762a1bSJed Brown     }
8409566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(user->dm, X, (void **)&x));
8419566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(9.0 * xm * ym));
842c4762a1bSJed Brown   }
8433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
844c4762a1bSJed Brown }
845c4762a1bSJed Brown 
846c4762a1bSJed Brown /*-----------------------------------------------------------------------*/
My_Monitor(Tao tao,PetscCtx ctx)847*2a8381b2SBarry Smith PetscErrorCode My_Monitor(Tao tao, PetscCtx ctx)
848d71ae5a4SJacob Faibussowitsch {
849c4762a1bSJed Brown   Vec X;
850c4762a1bSJed Brown 
851c4762a1bSJed Brown   PetscFunctionBegin;
8529566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao, &X));
8539566063dSJacob Faibussowitsch   PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
8543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
855c4762a1bSJed Brown }
856c4762a1bSJed Brown 
857c4762a1bSJed Brown /*TEST
858c4762a1bSJed Brown 
859c4762a1bSJed Brown    build:
860c4762a1bSJed Brown       requires: !complex
861c4762a1bSJed Brown 
862c4762a1bSJed Brown    test:
86310978b7dSBarry Smith       args: -tao_monitor_short -tao_type lmvm -da_grid_x 10 -da_grid_y 8 -tao_gatol 1.e-3
864c4762a1bSJed Brown       requires: !single
865c4762a1bSJed Brown 
866c4762a1bSJed Brown    test:
867c4762a1bSJed Brown       suffix: 2
868c4762a1bSJed Brown       nsize: 2
86910978b7dSBarry Smith       args: -tao_monitor_short -tao_type nls -tao_nls_ksp_max_it 15 -tao_gatol 1.e-4
870c4762a1bSJed Brown       filter: grep -v "nls ksp"
871c4762a1bSJed Brown       requires: !single
872c4762a1bSJed Brown 
873c4762a1bSJed Brown    test:
8744936d809SStefano Zampini       suffix: 2_snes
8754936d809SStefano Zampini       nsize: 2
87610978b7dSBarry Smith       args: -tao_monitor_short -tao_type snes -ksp_converged_maxits -ksp_max_it 15 -snes_atol 1.e-4
8774936d809SStefano Zampini       filter: grep -v "nls ksp"
8784936d809SStefano Zampini       requires: !single
8794936d809SStefano Zampini 
8804936d809SStefano Zampini    test:
881c4762a1bSJed Brown       suffix: 3
882c4762a1bSJed Brown       nsize: 3
88310978b7dSBarry Smith       args: -tao_monitor_short -tao_type cg -tao_cg_type fr -da_grid_x 10 -da_grid_y 10 -tao_gatol 1.e-3
884c4762a1bSJed Brown       requires: !single
885c4762a1bSJed Brown 
886c4762a1bSJed Brown    test:
8874936d809SStefano Zampini       suffix: 3_snes
8884936d809SStefano Zampini       nsize: 3
88910978b7dSBarry Smith       args: -tao_monitor_short -tao_type snes -snes_type ncg -snes_ncg_type fr -da_grid_x 10 -da_grid_y 10 -snes_atol 1.e-4
8904936d809SStefano Zampini       requires: !single
8914936d809SStefano Zampini 
8924936d809SStefano Zampini    test:
8934936d809SStefano Zampini       suffix: 4_snes_ngmres
8949f1eecfaSStefano Zampini       args: -tao_type snes -snes_type ngmres -npc_snes_type ncg -da_grid_x 10 -da_grid_y 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}
8954936d809SStefano Zampini       requires: !single
8964936d809SStefano Zampini 
8974936d809SStefano Zampini    test:
898c4762a1bSJed Brown       suffix: 5
899c4762a1bSJed Brown       nsize: 2
90010978b7dSBarry Smith       args: -tao_monitor_short -tao_type bmrm -da_grid_x 10 -da_grid_y 8 -tao_gatol 1.e-3
901c4762a1bSJed Brown       requires: !single
902c4762a1bSJed Brown 
903c4762a1bSJed Brown TEST*/
904