1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 minsurf1 [-help] [all TAO options] */
2c4762a1bSJed Brown
3c4762a1bSJed Brown /* Include "petsctao.h" so we can use TAO solvers. */
4c4762a1bSJed Brown #include <petsctao.h>
5c4762a1bSJed Brown
69371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to \n\
7c4762a1bSJed Brown solve an unconstrained minimization problem. This example is based on a \n\
8c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\
9c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
10c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
11c4762a1bSJed Brown The command line options are:\n\
12c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
13c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
14c4762a1bSJed Brown -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";
15c4762a1bSJed Brown
16c4762a1bSJed Brown /*
17c4762a1bSJed Brown User-defined application context - contains data needed by the
18c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient()
19c4762a1bSJed Brown and FormHessian().
20c4762a1bSJed Brown */
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown PetscInt mx, my; /* discretization in x, y directions */
23c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right; /* boundary values */
24c4762a1bSJed Brown Mat H;
25c4762a1bSJed Brown } AppCtx;
26c4762a1bSJed Brown
27c4762a1bSJed Brown /* -------- User-defined Routines --------- */
28c4762a1bSJed Brown
29c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
30c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
31c4762a1bSJed Brown static PetscErrorCode QuadraticH(AppCtx *, Vec, Mat);
32c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
33c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
34c4762a1bSJed Brown
main(int argc,char ** argv)35d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
36d71ae5a4SJacob Faibussowitsch {
37c4762a1bSJed Brown PetscInt N; /* Size of vector */
38c4762a1bSJed Brown PetscMPIInt size; /* Number of processors */
39c4762a1bSJed Brown Vec x; /* solution, gradient vectors */
40c4762a1bSJed Brown KSP ksp; /* PETSc Krylov subspace method */
41c4762a1bSJed Brown PetscBool flg; /* A return value when checking for user options */
42c4762a1bSJed Brown Tao tao; /* Tao solver context */
43c4762a1bSJed Brown AppCtx user; /* user-defined work context */
44c4762a1bSJed Brown
45c4762a1bSJed Brown /* Initialize TAO,PETSc */
46327415f7SBarry Smith PetscFunctionBeginUser;
47c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
48c4762a1bSJed Brown
499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
503c859ba3SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors");
51c4762a1bSJed Brown
52c4762a1bSJed Brown /* Specify default dimension of the problem */
539371c9d4SSatish Balay user.mx = 4;
549371c9d4SSatish Balay user.my = 4;
55c4762a1bSJed Brown
56c4762a1bSJed Brown /* Check for any command line arguments that override defaults */
579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg));
589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg));
59c4762a1bSJed Brown
609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Minimum Surface Area Problem -----\n"));
6163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n", user.mx, user.my));
62c4762a1bSJed Brown
63c4762a1bSJed Brown /* Calculate any derived values from parameters */
64c4762a1bSJed Brown N = user.mx * user.my;
65c4762a1bSJed Brown
66c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
679566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
689566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOLMVM));
69c4762a1bSJed Brown
70c4762a1bSJed Brown /* Initialize minsurf application data structure for use in the function evaluations */
719566063dSJacob Faibussowitsch PetscCall(MSA_BoundaryConditions(&user)); /* Application specific routine */
72c4762a1bSJed Brown
73c4762a1bSJed Brown /*
74c4762a1bSJed Brown Create a vector to hold the variables. Compute an initial solution.
75c4762a1bSJed Brown Set this vector, which will be used by TAO.
76c4762a1bSJed Brown */
779566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));
789566063dSJacob Faibussowitsch PetscCall(MSA_InitialPoint(&user, x)); /* Application specific routine */
799566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); /* A TAO routine */
80c4762a1bSJed Brown
81c4762a1bSJed Brown /* Provide TAO routines for function, gradient, and Hessian evaluation */
829566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
83c4762a1bSJed Brown
84c4762a1bSJed Brown /* Create a matrix data structure to store the Hessian. This structure will be used by TAO */
85f4f49eeaSPierre Jolivet PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 7, NULL, &user.H));
869566063dSJacob Faibussowitsch PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
879566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
88c4762a1bSJed Brown
89c4762a1bSJed Brown /* Check for any TAO command line options */
909566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
91c4762a1bSJed Brown
92c4762a1bSJed Brown /* Limit the number of iterations in the KSP linear solver */
939566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp));
94606f75f6SBarry Smith if (ksp) PetscCall(KSPSetTolerances(ksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, user.mx * user.my));
95c4762a1bSJed Brown
96c4762a1bSJed Brown /* SOLVE THE APPLICATION */
979566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
98c4762a1bSJed Brown
999566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.H));
1029566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom));
1039566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top));
1049566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left));
1059566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right));
106c4762a1bSJed Brown
1079566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
108b122ec5aSJacob Faibussowitsch return 0;
109c4762a1bSJed Brown }
110c4762a1bSJed Brown
111c4762a1bSJed Brown /* -------------------------------------------------------------------- */
112c4762a1bSJed Brown
113c4762a1bSJed Brown /* FormFunctionGradient - Evaluates function and corresponding gradient.
114c4762a1bSJed Brown
115c4762a1bSJed Brown Input Parameters:
116c4762a1bSJed Brown . tao - the Tao context
117c4762a1bSJed Brown . X - input vector
118c4762a1bSJed Brown . userCtx - optional user-defined context, as set by TaoSetFunctionGradient()
119c4762a1bSJed Brown
120c4762a1bSJed Brown Output Parameters:
121c4762a1bSJed Brown . fcn - the newly evaluated function
122c4762a1bSJed Brown . G - vector containing the newly evaluated gradient
123c4762a1bSJed Brown */
FormFunctionGradient(Tao tao,Vec X,PetscReal * fcn,Vec G,void * userCtx)124d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx)
125d71ae5a4SJacob Faibussowitsch {
126c4762a1bSJed Brown AppCtx *user = (AppCtx *)userCtx;
127c4762a1bSJed Brown PetscInt i, j, row;
128c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
129c4762a1bSJed Brown PetscReal rhx = mx + 1, rhy = my + 1;
130c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy, ft = 0;
131c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
132c4762a1bSJed Brown PetscReal df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
133c4762a1bSJed Brown PetscReal zero = 0.0;
134c4762a1bSJed Brown PetscScalar *g;
135c4762a1bSJed Brown const PetscScalar *x;
136c4762a1bSJed Brown
137c4762a1bSJed Brown PetscFunctionBeginUser;
1389566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero));
139c4762a1bSJed Brown
1409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
1419566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g));
142c4762a1bSJed Brown
143c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */
144c4762a1bSJed Brown for (j = 0; j < my; j++) {
145c4762a1bSJed Brown for (i = 0; i < mx; i++) {
146c4762a1bSJed Brown row = (j)*mx + (i);
147c4762a1bSJed Brown xc = x[row];
148c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc;
149c4762a1bSJed Brown if (i == 0) { /* left side */
150c4762a1bSJed Brown xl = user->left[j + 1];
151c4762a1bSJed Brown xlt = user->left[j + 2];
152c4762a1bSJed Brown } else {
153c4762a1bSJed Brown xl = x[row - 1];
154c4762a1bSJed Brown }
155c4762a1bSJed Brown
156c4762a1bSJed Brown if (j == 0) { /* bottom side */
157c4762a1bSJed Brown xb = user->bottom[i + 1];
158c4762a1bSJed Brown xrb = user->bottom[i + 2];
159c4762a1bSJed Brown } else {
160c4762a1bSJed Brown xb = x[row - mx];
161c4762a1bSJed Brown }
162c4762a1bSJed Brown
163c4762a1bSJed Brown if (i + 1 == mx) { /* right side */
164c4762a1bSJed Brown xr = user->right[j + 1];
165c4762a1bSJed Brown xrb = user->right[j];
166c4762a1bSJed Brown } else {
167c4762a1bSJed Brown xr = x[row + 1];
168c4762a1bSJed Brown }
169c4762a1bSJed Brown
170c4762a1bSJed Brown if (j + 1 == 0 + my) { /* top side */
171c4762a1bSJed Brown xt = user->top[i + 1];
172c4762a1bSJed Brown xlt = user->top[i];
173c4762a1bSJed Brown } else {
174c4762a1bSJed Brown xt = x[row + mx];
175c4762a1bSJed Brown }
176c4762a1bSJed Brown
177ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx];
178ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx];
179c4762a1bSJed Brown
180c4762a1bSJed Brown d1 = (xc - xl);
181c4762a1bSJed Brown d2 = (xc - xr);
182c4762a1bSJed Brown d3 = (xc - xt);
183c4762a1bSJed Brown d4 = (xc - xb);
184c4762a1bSJed Brown d5 = (xr - xrb);
185c4762a1bSJed Brown d6 = (xrb - xb);
186c4762a1bSJed Brown d7 = (xlt - xl);
187c4762a1bSJed Brown d8 = (xt - xlt);
188c4762a1bSJed Brown
189c4762a1bSJed Brown df1dxc = d1 * hydhx;
190c4762a1bSJed Brown df2dxc = (d1 * hydhx + d4 * hxdhy);
191c4762a1bSJed Brown df3dxc = d3 * hxdhy;
192c4762a1bSJed Brown df4dxc = (d2 * hydhx + d3 * hxdhy);
193c4762a1bSJed Brown df5dxc = d2 * hydhx;
194c4762a1bSJed Brown df6dxc = d4 * hxdhy;
195c4762a1bSJed Brown
196c4762a1bSJed Brown d1 *= rhx;
197c4762a1bSJed Brown d2 *= rhx;
198c4762a1bSJed Brown d3 *= rhy;
199c4762a1bSJed Brown d4 *= rhy;
200c4762a1bSJed Brown d5 *= rhy;
201c4762a1bSJed Brown d6 *= rhx;
202c4762a1bSJed Brown d7 *= rhy;
203c4762a1bSJed Brown d8 *= rhx;
204c4762a1bSJed Brown
205c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
206c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
207c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
208c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
209c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
210c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
211c4762a1bSJed Brown
212c4762a1bSJed Brown ft = ft + (f2 + f4);
213c4762a1bSJed Brown
214c4762a1bSJed Brown df1dxc /= f1;
215c4762a1bSJed Brown df2dxc /= f2;
216c4762a1bSJed Brown df3dxc /= f3;
217c4762a1bSJed Brown df4dxc /= f4;
218c4762a1bSJed Brown df5dxc /= f5;
219c4762a1bSJed Brown df6dxc /= f6;
220c4762a1bSJed Brown
221c4762a1bSJed Brown g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0;
222c4762a1bSJed Brown }
223c4762a1bSJed Brown }
224c4762a1bSJed Brown
225c4762a1bSJed Brown for (j = 0; j < my; j++) { /* left side */
226c4762a1bSJed Brown d3 = (user->left[j + 1] - user->left[j + 2]) * rhy;
227c4762a1bSJed Brown d2 = (user->left[j + 1] - x[j * mx]) * rhx;
228c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
229c4762a1bSJed Brown }
230c4762a1bSJed Brown
231c4762a1bSJed Brown for (i = 0; i < mx; i++) { /* bottom */
232c4762a1bSJed Brown d2 = (user->bottom[i + 1] - user->bottom[i + 2]) * rhx;
233c4762a1bSJed Brown d3 = (user->bottom[i + 1] - x[i]) * rhy;
234c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
235c4762a1bSJed Brown }
236c4762a1bSJed Brown
237c4762a1bSJed Brown for (j = 0; j < my; j++) { /* right side */
238c4762a1bSJed Brown d1 = (x[(j + 1) * mx - 1] - user->right[j + 1]) * rhx;
239c4762a1bSJed Brown d4 = (user->right[j] - user->right[j + 1]) * rhy;
240c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
241c4762a1bSJed Brown }
242c4762a1bSJed Brown
243c4762a1bSJed Brown for (i = 0; i < mx; i++) { /* top side */
244c4762a1bSJed Brown d1 = (x[(my - 1) * mx + i] - user->top[i + 1]) * rhy;
245c4762a1bSJed Brown d4 = (user->top[i + 1] - user->top[i]) * rhx;
246c4762a1bSJed Brown ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown
249c4762a1bSJed Brown /* Bottom left corner */
250c4762a1bSJed Brown d1 = (user->left[0] - user->left[1]) * rhy;
251c4762a1bSJed Brown d2 = (user->bottom[0] - user->bottom[1]) * rhx;
252c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
253c4762a1bSJed Brown
254c4762a1bSJed Brown /* Top right corner */
255c4762a1bSJed Brown d1 = (user->right[my + 1] - user->right[my]) * rhy;
256c4762a1bSJed Brown d2 = (user->top[mx + 1] - user->top[mx]) * rhx;
257c4762a1bSJed Brown ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2);
258c4762a1bSJed Brown
259c4762a1bSJed Brown (*fcn) = ft * area;
260c4762a1bSJed Brown
261c4762a1bSJed Brown /* Restore vectors */
2629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
2639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g));
2649566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67.0 * mx * my));
2653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
266c4762a1bSJed Brown }
267c4762a1bSJed Brown
268c4762a1bSJed Brown /* ------------------------------------------------------------------- */
269c4762a1bSJed Brown /*
270c4762a1bSJed Brown FormHessian - Evaluates the Hessian matrix.
271c4762a1bSJed Brown
272c4762a1bSJed Brown Input Parameters:
273c4762a1bSJed Brown . tao - the Tao context
274c4762a1bSJed Brown . x - input vector
275c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian()
276c4762a1bSJed Brown
277c4762a1bSJed Brown Output Parameters:
278c4762a1bSJed Brown . H - Hessian matrix
2797addb90fSBarry Smith . Hpre - optionally different matrix used to compute the preconditioner
280c4762a1bSJed Brown
281c4762a1bSJed Brown */
FormHessian(Tao tao,Vec X,Mat H,Mat Hpre,void * ptr)282d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
283d71ae5a4SJacob Faibussowitsch {
284c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
285c4762a1bSJed Brown
286c4762a1bSJed Brown PetscFunctionBeginUser;
287c4762a1bSJed Brown /* Evaluate the Hessian entries*/
2889566063dSJacob Faibussowitsch PetscCall(QuadraticH(user, X, H));
2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
290c4762a1bSJed Brown }
291c4762a1bSJed Brown
292c4762a1bSJed Brown /* ------------------------------------------------------------------- */
293c4762a1bSJed Brown /*
294c4762a1bSJed Brown QuadraticH - Evaluates the Hessian matrix.
295c4762a1bSJed Brown
296c4762a1bSJed Brown Input Parameters:
297c4762a1bSJed Brown . user - user-defined context, as set by TaoSetHessian()
298c4762a1bSJed Brown . X - input vector
299c4762a1bSJed Brown
300c4762a1bSJed Brown Output Parameter:
301c4762a1bSJed Brown . H - Hessian matrix
302c4762a1bSJed Brown */
QuadraticH(AppCtx * user,Vec X,Mat Hessian)303d71ae5a4SJacob Faibussowitsch PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian)
304d71ae5a4SJacob Faibussowitsch {
305c4762a1bSJed Brown PetscInt i, j, k, row;
306c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
307c4762a1bSJed Brown PetscInt col[7];
308c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
309c4762a1bSJed Brown PetscReal rhx = mx + 1, rhy = my + 1;
310c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
311c4762a1bSJed Brown PetscReal hl, hr, ht, hb, hc, htl, hbr;
312c4762a1bSJed Brown const PetscScalar *x;
313c4762a1bSJed Brown PetscReal v[7];
314c4762a1bSJed Brown
315362febeeSStefano Zampini PetscFunctionBeginUser;
316c4762a1bSJed Brown /* Get pointers to vector data */
3179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
318c4762a1bSJed Brown
319c4762a1bSJed Brown /* Initialize matrix entries to zero */
3209566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hessian));
321c4762a1bSJed Brown
322c4762a1bSJed Brown /* Set various matrix options */
3239566063dSJacob Faibussowitsch PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
324c4762a1bSJed Brown
325c4762a1bSJed Brown /* Compute Hessian over the locally owned part of the mesh */
326c4762a1bSJed Brown for (i = 0; i < mx; i++) {
327c4762a1bSJed Brown for (j = 0; j < my; j++) {
328c4762a1bSJed Brown row = (j)*mx + (i);
329c4762a1bSJed Brown
330c4762a1bSJed Brown xc = x[row];
331c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc;
332c4762a1bSJed Brown
333c4762a1bSJed Brown /* Left side */
334c4762a1bSJed Brown if (i == 0) {
335c4762a1bSJed Brown xl = user->left[j + 1];
336c4762a1bSJed Brown xlt = user->left[j + 2];
337c4762a1bSJed Brown } else {
338c4762a1bSJed Brown xl = x[row - 1];
339c4762a1bSJed Brown }
340c4762a1bSJed Brown
341c4762a1bSJed Brown if (j == 0) {
342c4762a1bSJed Brown xb = user->bottom[i + 1];
343c4762a1bSJed Brown xrb = user->bottom[i + 2];
344c4762a1bSJed Brown } else {
345c4762a1bSJed Brown xb = x[row - mx];
346c4762a1bSJed Brown }
347c4762a1bSJed Brown
348c4762a1bSJed Brown if (i + 1 == mx) {
349c4762a1bSJed Brown xr = user->right[j + 1];
350c4762a1bSJed Brown xrb = user->right[j];
351c4762a1bSJed Brown } else {
352c4762a1bSJed Brown xr = x[row + 1];
353c4762a1bSJed Brown }
354c4762a1bSJed Brown
355c4762a1bSJed Brown if (j + 1 == my) {
356c4762a1bSJed Brown xt = user->top[i + 1];
357c4762a1bSJed Brown xlt = user->top[i];
358c4762a1bSJed Brown } else {
359c4762a1bSJed Brown xt = x[row + mx];
360c4762a1bSJed Brown }
361c4762a1bSJed Brown
362ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx];
363ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx];
364c4762a1bSJed Brown
365c4762a1bSJed Brown d1 = (xc - xl) * rhx;
366c4762a1bSJed Brown d2 = (xc - xr) * rhx;
367c4762a1bSJed Brown d3 = (xc - xt) * rhy;
368c4762a1bSJed Brown d4 = (xc - xb) * rhy;
369c4762a1bSJed Brown d5 = (xrb - xr) * rhy;
370c4762a1bSJed Brown d6 = (xrb - xb) * rhx;
371c4762a1bSJed Brown d7 = (xlt - xl) * rhy;
372c4762a1bSJed Brown d8 = (xlt - xt) * rhx;
373c4762a1bSJed Brown
374c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
375c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
376c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
377c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
378c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
379c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
380c4762a1bSJed Brown
381c4762a1bSJed Brown hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
382c4762a1bSJed Brown hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
383c4762a1bSJed Brown ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
384c4762a1bSJed Brown hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
385c4762a1bSJed Brown
386c4762a1bSJed Brown hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
387c4762a1bSJed Brown htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
388c4762a1bSJed Brown
3899371c9d4SSatish 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);
390c4762a1bSJed Brown
3919371c9d4SSatish Balay hl *= 0.5;
3929371c9d4SSatish Balay hr *= 0.5;
3939371c9d4SSatish Balay ht *= 0.5;
3949371c9d4SSatish Balay hb *= 0.5;
3959371c9d4SSatish Balay hbr *= 0.5;
3969371c9d4SSatish Balay htl *= 0.5;
3979371c9d4SSatish Balay hc *= 0.5;
398c4762a1bSJed Brown
399c4762a1bSJed Brown k = 0;
400c4762a1bSJed Brown if (j > 0) {
4019371c9d4SSatish Balay v[k] = hb;
4029371c9d4SSatish Balay col[k] = row - mx;
4039371c9d4SSatish Balay k++;
404c4762a1bSJed Brown }
405c4762a1bSJed Brown
406c4762a1bSJed Brown if (j > 0 && i < mx - 1) {
4079371c9d4SSatish Balay v[k] = hbr;
4089371c9d4SSatish Balay col[k] = row - mx + 1;
4099371c9d4SSatish Balay k++;
410c4762a1bSJed Brown }
411c4762a1bSJed Brown
412c4762a1bSJed Brown if (i > 0) {
4139371c9d4SSatish Balay v[k] = hl;
4149371c9d4SSatish Balay col[k] = row - 1;
4159371c9d4SSatish Balay k++;
416c4762a1bSJed Brown }
417c4762a1bSJed Brown
4189371c9d4SSatish Balay v[k] = hc;
4199371c9d4SSatish Balay col[k] = row;
4209371c9d4SSatish Balay k++;
421c4762a1bSJed Brown
422c4762a1bSJed Brown if (i < mx - 1) {
4239371c9d4SSatish Balay v[k] = hr;
4249371c9d4SSatish Balay col[k] = row + 1;
4259371c9d4SSatish Balay k++;
426c4762a1bSJed Brown }
427c4762a1bSJed Brown
428c4762a1bSJed Brown if (i > 0 && j < my - 1) {
4299371c9d4SSatish Balay v[k] = htl;
4309371c9d4SSatish Balay col[k] = row + mx - 1;
4319371c9d4SSatish Balay k++;
432c4762a1bSJed Brown }
433c4762a1bSJed Brown
434c4762a1bSJed Brown if (j < my - 1) {
4359371c9d4SSatish Balay v[k] = ht;
4369371c9d4SSatish Balay col[k] = row + mx;
4379371c9d4SSatish Balay k++;
438c4762a1bSJed Brown }
439c4762a1bSJed Brown
440c4762a1bSJed Brown /*
441c4762a1bSJed Brown Set matrix values using local numbering, which was defined
442c4762a1bSJed Brown earlier, in the main routine.
443c4762a1bSJed Brown */
4449566063dSJacob Faibussowitsch PetscCall(MatSetValues(Hessian, 1, &row, k, col, v, INSERT_VALUES));
445c4762a1bSJed Brown }
446c4762a1bSJed Brown }
447c4762a1bSJed Brown
448c4762a1bSJed Brown /* Restore vectors */
4499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
450c4762a1bSJed Brown
451c4762a1bSJed Brown /* Assemble the matrix */
4529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY));
4539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY));
454c4762a1bSJed Brown
4559566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199.0 * mx * my));
4563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
457c4762a1bSJed Brown }
458c4762a1bSJed Brown
459c4762a1bSJed Brown /* ------------------------------------------------------------------- */
460c4762a1bSJed Brown /*
461c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for
462c4762a1bSJed Brown the region.
463c4762a1bSJed Brown
464c4762a1bSJed Brown Input Parameter:
465c4762a1bSJed Brown . user - user-defined application context
466c4762a1bSJed Brown
467c4762a1bSJed Brown Output Parameter:
468c4762a1bSJed Brown . user - user-defined application context
469c4762a1bSJed Brown */
MSA_BoundaryConditions(AppCtx * user)470d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user)
471d71ae5a4SJacob Faibussowitsch {
472c4762a1bSJed Brown PetscInt i, j, k, limit = 0;
473c4762a1bSJed Brown PetscInt maxits = 5;
474c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
475c4762a1bSJed Brown PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0;
476c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, three = 3.0, tol = 1e-10;
477c4762a1bSJed Brown PetscReal fnorm, det, hx, hy, xt = 0, yt = 0;
478c4762a1bSJed Brown PetscReal u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
479c4762a1bSJed Brown PetscReal b = -0.5, t = 0.5, l = -0.5, r = 0.5;
480c4762a1bSJed Brown PetscReal *boundary;
481c4762a1bSJed Brown
482c4762a1bSJed Brown PetscFunctionBeginUser;
4839371c9d4SSatish Balay bsize = mx + 2;
4849371c9d4SSatish Balay lsize = my + 2;
4859371c9d4SSatish Balay rsize = my + 2;
4869371c9d4SSatish Balay tsize = mx + 2;
487c4762a1bSJed Brown
4889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize, &user->bottom));
4899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize, &user->top));
4909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &user->left));
4919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize, &user->right));
492c4762a1bSJed Brown
4939371c9d4SSatish Balay hx = (r - l) / (mx + 1);
4949371c9d4SSatish Balay hy = (t - b) / (my + 1);
495c4762a1bSJed Brown
496c4762a1bSJed Brown for (j = 0; j < 4; j++) {
497c4762a1bSJed Brown if (j == 0) {
498c4762a1bSJed Brown yt = b;
499c4762a1bSJed Brown xt = l;
500c4762a1bSJed Brown limit = bsize;
501c4762a1bSJed Brown boundary = user->bottom;
502c4762a1bSJed Brown } else if (j == 1) {
503c4762a1bSJed Brown yt = t;
504c4762a1bSJed Brown xt = l;
505c4762a1bSJed Brown limit = tsize;
506c4762a1bSJed Brown boundary = user->top;
507c4762a1bSJed Brown } else if (j == 2) {
508c4762a1bSJed Brown yt = b;
509c4762a1bSJed Brown xt = l;
510c4762a1bSJed Brown limit = lsize;
511c4762a1bSJed Brown boundary = user->left;
512c4762a1bSJed Brown } else { /* if (j==3) */
513c4762a1bSJed Brown yt = b;
514c4762a1bSJed Brown xt = r;
515c4762a1bSJed Brown limit = rsize;
516c4762a1bSJed Brown boundary = user->right;
517c4762a1bSJed Brown }
518c4762a1bSJed Brown
519c4762a1bSJed Brown for (i = 0; i < limit; i++) {
520c4762a1bSJed Brown u1 = xt;
521c4762a1bSJed Brown u2 = -yt;
522c4762a1bSJed Brown for (k = 0; k < maxits; k++) {
523c4762a1bSJed Brown nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
524c4762a1bSJed Brown nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
525c4762a1bSJed Brown fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2);
526c4762a1bSJed Brown if (fnorm <= tol) break;
527c4762a1bSJed Brown njac11 = one + u2 * u2 - u1 * u1;
528c4762a1bSJed Brown njac12 = two * u1 * u2;
529c4762a1bSJed Brown njac21 = -two * u1 * u2;
530c4762a1bSJed Brown njac22 = -one - u1 * u1 + u2 * u2;
531c4762a1bSJed Brown det = njac11 * njac22 - njac21 * njac12;
532c4762a1bSJed Brown u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det;
533c4762a1bSJed Brown u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det;
534c4762a1bSJed Brown }
535c4762a1bSJed Brown
536c4762a1bSJed Brown boundary[i] = u1 * u1 - u2 * u2;
537c4762a1bSJed Brown if (j == 0 || j == 1) {
538c4762a1bSJed Brown xt = xt + hx;
539c4762a1bSJed Brown } else { /* if (j==2 || j==3) */
540c4762a1bSJed Brown yt = yt + hy;
541c4762a1bSJed Brown }
542c4762a1bSJed Brown }
543c4762a1bSJed Brown }
5443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
545c4762a1bSJed Brown }
546c4762a1bSJed Brown
547c4762a1bSJed Brown /* ------------------------------------------------------------------- */
548c4762a1bSJed Brown /*
549c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways.
550c4762a1bSJed Brown
551c4762a1bSJed Brown Input Parameters:
552c4762a1bSJed Brown . user - user-defined application context
553c4762a1bSJed Brown . X - vector for initial guess
554c4762a1bSJed Brown
555c4762a1bSJed Brown Output Parameters:
556c4762a1bSJed Brown . X - newly computed initial guess
557c4762a1bSJed Brown */
MSA_InitialPoint(AppCtx * user,Vec X)558d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X)
559d71ae5a4SJacob Faibussowitsch {
560c4762a1bSJed Brown PetscInt start = -1, i, j;
561c4762a1bSJed Brown PetscReal zero = 0.0;
562c4762a1bSJed Brown PetscBool flg;
563c4762a1bSJed Brown
564c4762a1bSJed Brown PetscFunctionBeginUser;
5659566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero));
5669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg));
567c4762a1bSJed Brown
568c4762a1bSJed Brown if (flg && start == 0) { /* The zero vector is reasonable */
5699566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero));
570c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */
571c4762a1bSJed Brown PetscInt row;
572c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
573c4762a1bSJed Brown PetscReal *x;
574c4762a1bSJed Brown
575c4762a1bSJed Brown /* Get pointers to vector data */
5769566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
577c4762a1bSJed Brown /* Perform local computations */
578c4762a1bSJed Brown for (j = 0; j < my; j++) {
579c4762a1bSJed Brown for (i = 0; i < mx; i++) {
580c4762a1bSJed Brown row = (j)*mx + (i);
581c4762a1bSJed Brown x[row] = (((j + 1) * user->bottom[i + 1] + (my - j + 1) * user->top[i + 1]) / (my + 2) + ((i + 1) * user->left[j + 1] + (mx - i + 1) * user->right[j + 1]) / (mx + 2)) / 2.0;
582c4762a1bSJed Brown }
583c4762a1bSJed Brown }
584c4762a1bSJed Brown /* Restore vectors */
5859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
586c4762a1bSJed Brown }
5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
588c4762a1bSJed Brown }
589c4762a1bSJed Brown
590c4762a1bSJed Brown /*TEST
591c4762a1bSJed Brown
592c4762a1bSJed Brown build:
593c4762a1bSJed Brown requires: !complex
594c4762a1bSJed Brown
595c4762a1bSJed Brown test:
59610978b7dSBarry Smith args: -tao_monitor_short -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4
597c4762a1bSJed Brown requires: !single
598c4762a1bSJed Brown
599c4762a1bSJed Brown test:
600c4762a1bSJed Brown suffix: 2
60110978b7dSBarry Smith args: -tao_monitor_short -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3
602c4762a1bSJed Brown requires: !single
603c4762a1bSJed Brown
604c4762a1bSJed Brown test:
605c4762a1bSJed Brown suffix: 3
60610978b7dSBarry Smith args: -tao_monitor_short -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3
607c4762a1bSJed Brown requires: !single
608c4762a1bSJed Brown
609c4762a1bSJed Brown test:
610c4762a1bSJed Brown suffix: 4
61110978b7dSBarry Smith args: -tao_monitor_short -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4
612c4762a1bSJed Brown requires: !single
613c4762a1bSJed Brown
614c4762a1bSJed Brown test:
6150f0abf79SStefano Zampini suffix: 4_ew
61610978b7dSBarry Smith args: -tao_monitor_short -tao_type bntr -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4
6170f0abf79SStefano Zampini requires: !single
6180f0abf79SStefano Zampini
6190f0abf79SStefano Zampini test:
620c4762a1bSJed Brown suffix: 5
62110978b7dSBarry Smith args: -tao_monitor_short -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4
622c4762a1bSJed Brown requires: !single
623c4762a1bSJed Brown
624c4762a1bSJed Brown test:
6250f0abf79SStefano Zampini suffix: 5_ew
62610978b7dSBarry Smith args: -tao_monitor_short -tao_type bntl -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4
6270f0abf79SStefano Zampini requires: !single
6280f0abf79SStefano Zampini
6290f0abf79SStefano Zampini test:
630c4762a1bSJed Brown suffix: 6
63110978b7dSBarry Smith args: -tao_monitor_short -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4
632c4762a1bSJed Brown requires: !single
633c4762a1bSJed Brown
634c4762a1bSJed Brown test:
6350f0abf79SStefano Zampini suffix: 6_ew
63610978b7dSBarry Smith args: -tao_monitor_short -tao_type bnls -tao_ksp_ew -tao_bnk_ksp_ew_version 3 -mx 10 -my 8 -tao_gatol 1.e-4
6370f0abf79SStefano Zampini requires: !single
6380f0abf79SStefano Zampini
6390f0abf79SStefano Zampini test:
640c4762a1bSJed Brown suffix: 7
641*a336c150SZach Atkins args: -tao_monitor_short -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_bnk_cg_tao_monitor_short
642c4762a1bSJed Brown requires: !single
643c4762a1bSJed Brown
644c4762a1bSJed Brown test:
645c4762a1bSJed Brown suffix: 8
646*a336c150SZach Atkins args: -tao_monitor_short -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_bnk_cg_tao_monitor_short
647c4762a1bSJed Brown requires: !single
648c4762a1bSJed Brown
649c4762a1bSJed Brown test:
650c4762a1bSJed Brown suffix: 9
651*a336c150SZach Atkins args: -tao_monitor_short -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_bnk_cg_tao_monitor_short
652c4762a1bSJed Brown requires: !single
653c4762a1bSJed Brown
65434ad9904SAlp Dener test:
65534ad9904SAlp Dener suffix: 10
656*a336c150SZach Atkins args: -tao_monitor_short -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_mf_hessian -tao_bnk_cg_tao_monitor_short
65734ad9904SAlp Dener
65834ad9904SAlp Dener test:
65934ad9904SAlp Dener suffix: 11
66010978b7dSBarry Smith args: -tao_monitor_short -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
66134ad9904SAlp Dener requires: !single
66234ad9904SAlp Dener
66334ad9904SAlp Dener test:
66434ad9904SAlp Dener suffix: 12
66510978b7dSBarry Smith args: -tao_monitor_short -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
66634ad9904SAlp Dener requires: !single
66734ad9904SAlp Dener
668c4762a1bSJed Brown TEST*/
669