1c4762a1bSJed Brown #include <petsctao.h>
2c4762a1bSJed Brown
39371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to\n\
4c4762a1bSJed Brown solve an unconstrained system of equations. This example is based on a\n\
5c4762a1bSJed Brown problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and\n\
6c4762a1bSJed Brown boundary values along the edges of the domain, the objective is to find the\n\
7c4762a1bSJed Brown surface with the minimal area that satisfies the boundary conditions.\n\
8c4762a1bSJed Brown This application solves this problem using complimentarity -- We are actually\n\
9c4762a1bSJed Brown solving the system (grad f)_i >= 0, if x_i == l_i \n\
10c4762a1bSJed Brown (grad f)_i = 0, if l_i < x_i < u_i \n\
11c4762a1bSJed Brown (grad f)_i <= 0, if x_i == u_i \n\
12c4762a1bSJed Brown where f is the function to be minimized. \n\
13c4762a1bSJed Brown \n\
14c4762a1bSJed Brown The command line options are:\n\
15c4762a1bSJed Brown -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
16c4762a1bSJed Brown -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
17c4762a1bSJed Brown -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";
18c4762a1bSJed Brown
19c4762a1bSJed Brown /*
20c4762a1bSJed Brown User-defined application context - contains data needed by the
21c4762a1bSJed Brown application-provided call-back routines, FormFunctionGradient(),
22c4762a1bSJed Brown FormHessian().
23c4762a1bSJed Brown */
24c4762a1bSJed Brown typedef struct {
25c4762a1bSJed Brown PetscInt mx, my;
26c4762a1bSJed Brown PetscReal *bottom, *top, *left, *right;
27c4762a1bSJed Brown } AppCtx;
28c4762a1bSJed Brown
29c4762a1bSJed Brown /* -------- User-defined Routines --------- */
30c4762a1bSJed Brown
31c4762a1bSJed Brown static PetscErrorCode MSA_BoundaryConditions(AppCtx *);
32c4762a1bSJed Brown static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec);
33c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
34c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *);
35c4762a1bSJed Brown
main(int argc,char ** argv)36d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
37d71ae5a4SJacob Faibussowitsch {
38c4762a1bSJed Brown Vec x; /* solution vector */
39c4762a1bSJed Brown Vec c; /* Constraints function vector */
40c4762a1bSJed Brown Vec xl, xu; /* Bounds on the variables */
41c4762a1bSJed Brown PetscBool flg; /* A return variable when checking for user options */
42c4762a1bSJed Brown Tao tao; /* TAO solver context */
43c4762a1bSJed Brown Mat J; /* Jacobian matrix */
44c4762a1bSJed Brown PetscInt N; /* Number of elements in vector */
45c4762a1bSJed Brown PetscScalar lb = PETSC_NINFINITY; /* lower bound constant */
46c4762a1bSJed Brown PetscScalar ub = PETSC_INFINITY; /* upper bound constant */
47c4762a1bSJed Brown AppCtx user; /* user-defined work context */
48c4762a1bSJed Brown
49c4762a1bSJed Brown /* Initialize PETSc, TAO */
50327415f7SBarry Smith PetscFunctionBeginUser;
51*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
52c4762a1bSJed Brown
53c4762a1bSJed Brown /* Specify default dimension of the problem */
549371c9d4SSatish Balay user.mx = 4;
559371c9d4SSatish Balay user.my = 4;
56c4762a1bSJed Brown
57c4762a1bSJed Brown /* Check for any command line arguments that override defaults */
589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg));
599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg));
60c4762a1bSJed Brown
61c4762a1bSJed Brown /* Calculate any derived values from parameters */
62c4762a1bSJed Brown N = user.mx * user.my;
63c4762a1bSJed Brown
649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Minimum Surface Area Problem -----\n"));
6563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx:%" PetscInt_FMT ", my:%" PetscInt_FMT "\n", user.mx, user.my));
66c4762a1bSJed Brown
67c4762a1bSJed Brown /* Create appropriate vectors and matrices */
689566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(MPI_COMM_SELF, N, &x));
699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &c));
709566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(MPI_COMM_SELF, N, N, 7, NULL, &J));
71c4762a1bSJed Brown
72c4762a1bSJed Brown /* The TAO code begins here */
73c4762a1bSJed Brown
74c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
759566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
769566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOSSILS));
77c4762a1bSJed Brown
78c4762a1bSJed Brown /* Set data structure */
799566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x));
80c4762a1bSJed Brown
81c4762a1bSJed Brown /* Set routines for constraints function and Jacobian evaluation */
829566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user));
839566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user));
84c4762a1bSJed Brown
85c4762a1bSJed Brown /* Set the variable bounds */
869566063dSJacob Faibussowitsch PetscCall(MSA_BoundaryConditions(&user));
87c4762a1bSJed Brown
88c4762a1bSJed Brown /* Set initial solution guess */
899566063dSJacob Faibussowitsch PetscCall(MSA_InitialPoint(&user, x));
90c4762a1bSJed Brown
91c4762a1bSJed Brown /* Set Bounds on variables */
929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xl));
939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xu));
949566063dSJacob Faibussowitsch PetscCall(VecSet(xl, lb));
959566063dSJacob Faibussowitsch PetscCall(VecSet(xu, ub));
969566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, xl, xu));
97c4762a1bSJed Brown
98c4762a1bSJed Brown /* Check for any tao command line options */
999566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
100c4762a1bSJed Brown
101c4762a1bSJed Brown /* Solve the application */
1029566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
103c4762a1bSJed Brown
104c4762a1bSJed Brown /* Free Tao data structures */
1059566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
106c4762a1bSJed Brown
107c4762a1bSJed Brown /* Free PETSc data structures */
1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xl));
1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xu));
1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c));
1129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
113c4762a1bSJed Brown
114c4762a1bSJed Brown /* Free user-created data structures */
1159566063dSJacob Faibussowitsch PetscCall(PetscFree(user.bottom));
1169566063dSJacob Faibussowitsch PetscCall(PetscFree(user.top));
1179566063dSJacob Faibussowitsch PetscCall(PetscFree(user.left));
1189566063dSJacob Faibussowitsch PetscCall(PetscFree(user.right));
119c4762a1bSJed Brown
1209566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
121b122ec5aSJacob Faibussowitsch return 0;
122c4762a1bSJed Brown }
123c4762a1bSJed Brown
124c4762a1bSJed Brown /* -------------------------------------------------------------------- */
125c4762a1bSJed Brown
126c4762a1bSJed Brown /* FormConstraints - Evaluates gradient of f.
127c4762a1bSJed Brown
128c4762a1bSJed Brown Input Parameters:
129c4762a1bSJed Brown . tao - the TAO_APPLICATION context
130c4762a1bSJed Brown . X - input vector
131c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetConstraintsRoutine()
132c4762a1bSJed Brown
133c4762a1bSJed Brown Output Parameters:
134c4762a1bSJed Brown . G - vector containing the newly evaluated gradient
135c4762a1bSJed Brown */
FormConstraints(Tao tao,Vec X,Vec G,void * ptr)136d71ae5a4SJacob Faibussowitsch PetscErrorCode FormConstraints(Tao tao, Vec X, Vec G, void *ptr)
137d71ae5a4SJacob Faibussowitsch {
138c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
139c4762a1bSJed Brown PetscInt i, j, row;
140c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
141c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
142c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
143c4762a1bSJed Brown PetscReal df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
144c4762a1bSJed Brown PetscScalar zero = 0.0;
145c4762a1bSJed Brown PetscScalar *g, *x;
146c4762a1bSJed Brown
147c4762a1bSJed Brown PetscFunctionBegin;
148c4762a1bSJed Brown /* Initialize vector to zero */
1499566063dSJacob Faibussowitsch PetscCall(VecSet(G, zero));
150c4762a1bSJed Brown
151c4762a1bSJed Brown /* Get pointers to vector data */
1529566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
1539566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g));
154c4762a1bSJed Brown
155c4762a1bSJed Brown /* Compute function over the locally owned part of the mesh */
156c4762a1bSJed Brown for (j = 0; j < my; j++) {
157c4762a1bSJed Brown for (i = 0; i < mx; i++) {
158c4762a1bSJed Brown row = j * mx + i;
159c4762a1bSJed Brown
160c4762a1bSJed Brown xc = x[row];
161c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc;
162c4762a1bSJed Brown
163c4762a1bSJed Brown if (i == 0) { /* left side */
164c4762a1bSJed Brown xl = user->left[j + 1];
165c4762a1bSJed Brown xlt = user->left[j + 2];
166c4762a1bSJed Brown } else {
167c4762a1bSJed Brown xl = x[row - 1];
168c4762a1bSJed Brown }
169c4762a1bSJed Brown
170c4762a1bSJed Brown if (j == 0) { /* bottom side */
171c4762a1bSJed Brown xb = user->bottom[i + 1];
172c4762a1bSJed Brown xrb = user->bottom[i + 2];
173c4762a1bSJed Brown } else {
174c4762a1bSJed Brown xb = x[row - mx];
175c4762a1bSJed Brown }
176c4762a1bSJed Brown
177c4762a1bSJed Brown if (i + 1 == mx) { /* right side */
178c4762a1bSJed Brown xr = user->right[j + 1];
179c4762a1bSJed Brown xrb = user->right[j];
180c4762a1bSJed Brown } else {
181c4762a1bSJed Brown xr = x[row + 1];
182c4762a1bSJed Brown }
183c4762a1bSJed Brown
184c4762a1bSJed Brown if (j + 1 == 0 + my) { /* top side */
185c4762a1bSJed Brown xt = user->top[i + 1];
186c4762a1bSJed Brown xlt = user->top[i];
187c4762a1bSJed Brown } else {
188c4762a1bSJed Brown xt = x[row + mx];
189c4762a1bSJed Brown }
190c4762a1bSJed Brown
191ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx];
192ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx];
193c4762a1bSJed Brown
194c4762a1bSJed Brown d1 = (xc - xl);
195c4762a1bSJed Brown d2 = (xc - xr);
196c4762a1bSJed Brown d3 = (xc - xt);
197c4762a1bSJed Brown d4 = (xc - xb);
198c4762a1bSJed Brown d5 = (xr - xrb);
199c4762a1bSJed Brown d6 = (xrb - xb);
200c4762a1bSJed Brown d7 = (xlt - xl);
201c4762a1bSJed Brown d8 = (xt - xlt);
202c4762a1bSJed Brown
203c4762a1bSJed Brown df1dxc = d1 * hydhx;
204c4762a1bSJed Brown df2dxc = (d1 * hydhx + d4 * hxdhy);
205c4762a1bSJed Brown df3dxc = d3 * hxdhy;
206c4762a1bSJed Brown df4dxc = (d2 * hydhx + d3 * hxdhy);
207c4762a1bSJed Brown df5dxc = d2 * hydhx;
208c4762a1bSJed Brown df6dxc = d4 * hxdhy;
209c4762a1bSJed Brown
210c4762a1bSJed Brown d1 /= hx;
211c4762a1bSJed Brown d2 /= hx;
212c4762a1bSJed Brown d3 /= hy;
213c4762a1bSJed Brown d4 /= hy;
214c4762a1bSJed Brown d5 /= hy;
215c4762a1bSJed Brown d6 /= hx;
216c4762a1bSJed Brown d7 /= hy;
217c4762a1bSJed Brown d8 /= hx;
218c4762a1bSJed Brown
219c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
220c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
221c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
222c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
223c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
224c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
225c4762a1bSJed Brown
226c4762a1bSJed Brown df1dxc /= f1;
227c4762a1bSJed Brown df2dxc /= f2;
228c4762a1bSJed Brown df3dxc /= f3;
229c4762a1bSJed Brown df4dxc /= f4;
230c4762a1bSJed Brown df5dxc /= f5;
231c4762a1bSJed Brown df6dxc /= f6;
232c4762a1bSJed Brown
233c4762a1bSJed Brown g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0;
234c4762a1bSJed Brown }
235c4762a1bSJed Brown }
236c4762a1bSJed Brown
237c4762a1bSJed Brown /* Restore vectors */
2389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
2399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g));
2409566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(67 * mx * my));
2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown
244c4762a1bSJed Brown /* ------------------------------------------------------------------- */
245c4762a1bSJed Brown /*
246c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix.
247c4762a1bSJed Brown
248c4762a1bSJed Brown Input Parameters:
249c4762a1bSJed Brown . tao - the TAO_APPLICATION context
250c4762a1bSJed Brown . X - input vector
251c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetJacobian()
252c4762a1bSJed Brown
253c4762a1bSJed Brown Output Parameters:
254c4762a1bSJed Brown . tH - Jacobian matrix
255c4762a1bSJed Brown
256c4762a1bSJed Brown */
FormJacobian(Tao tao,Vec X,Mat H,Mat tHPre,void * ptr)257d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(Tao tao, Vec X, Mat H, Mat tHPre, void *ptr)
258d71ae5a4SJacob Faibussowitsch {
259c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
260c4762a1bSJed Brown PetscInt i, j, k, row;
261c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
262c4762a1bSJed Brown PetscInt col[7];
263c4762a1bSJed Brown PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy;
264c4762a1bSJed Brown PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
265c4762a1bSJed Brown PetscReal hl, hr, ht, hb, hc, htl, hbr;
266c4762a1bSJed Brown const PetscScalar *x;
267c4762a1bSJed Brown PetscScalar v[7];
268c4762a1bSJed Brown PetscBool assembled;
269c4762a1bSJed Brown
270c4762a1bSJed Brown /* Set various matrix options */
271362febeeSStefano Zampini PetscFunctionBegin;
2729566063dSJacob Faibussowitsch PetscCall(MatSetOption(H, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
2739566063dSJacob Faibussowitsch PetscCall(MatAssembled(H, &assembled));
2749566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(H));
275c4762a1bSJed Brown
276c4762a1bSJed Brown /* Get pointers to vector data */
2779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
278c4762a1bSJed Brown
279c4762a1bSJed Brown /* Compute Jacobian over the locally owned part of the mesh */
280c4762a1bSJed Brown for (i = 0; i < mx; i++) {
281c4762a1bSJed Brown for (j = 0; j < my; j++) {
282c4762a1bSJed Brown row = j * mx + i;
283c4762a1bSJed Brown
284c4762a1bSJed Brown xc = x[row];
285c4762a1bSJed Brown xlt = xrb = xl = xr = xb = xt = xc;
286c4762a1bSJed Brown
287c4762a1bSJed Brown /* Left side */
288c4762a1bSJed Brown if (i == 0) {
289c4762a1bSJed Brown xl = user->left[j + 1];
290c4762a1bSJed Brown xlt = user->left[j + 2];
291c4762a1bSJed Brown } else {
292c4762a1bSJed Brown xl = x[row - 1];
293c4762a1bSJed Brown }
294c4762a1bSJed Brown
295c4762a1bSJed Brown if (j == 0) {
296c4762a1bSJed Brown xb = user->bottom[i + 1];
297c4762a1bSJed Brown xrb = user->bottom[i + 2];
298c4762a1bSJed Brown } else {
299c4762a1bSJed Brown xb = x[row - mx];
300c4762a1bSJed Brown }
301c4762a1bSJed Brown
302c4762a1bSJed Brown if (i + 1 == mx) {
303c4762a1bSJed Brown xr = user->right[j + 1];
304c4762a1bSJed Brown xrb = user->right[j];
305c4762a1bSJed Brown } else {
306c4762a1bSJed Brown xr = x[row + 1];
307c4762a1bSJed Brown }
308c4762a1bSJed Brown
309c4762a1bSJed Brown if (j + 1 == my) {
310c4762a1bSJed Brown xt = user->top[i + 1];
311c4762a1bSJed Brown xlt = user->top[i];
312c4762a1bSJed Brown } else {
313c4762a1bSJed Brown xt = x[row + mx];
314c4762a1bSJed Brown }
315c4762a1bSJed Brown
316ad540459SPierre Jolivet if (i > 0 && j + 1 < my) xlt = x[row - 1 + mx];
317ad540459SPierre Jolivet if (j > 0 && i + 1 < mx) xrb = x[row + 1 - mx];
318c4762a1bSJed Brown
319c4762a1bSJed Brown d1 = (xc - xl) / hx;
320c4762a1bSJed Brown d2 = (xc - xr) / hx;
321c4762a1bSJed Brown d3 = (xc - xt) / hy;
322c4762a1bSJed Brown d4 = (xc - xb) / hy;
323c4762a1bSJed Brown d5 = (xrb - xr) / hy;
324c4762a1bSJed Brown d6 = (xrb - xb) / hx;
325c4762a1bSJed Brown d7 = (xlt - xl) / hy;
326c4762a1bSJed Brown d8 = (xlt - xt) / hx;
327c4762a1bSJed Brown
328c4762a1bSJed Brown f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
329c4762a1bSJed Brown f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
330c4762a1bSJed Brown f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
331c4762a1bSJed Brown f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
332c4762a1bSJed Brown f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
333c4762a1bSJed Brown f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
334c4762a1bSJed Brown
335c4762a1bSJed Brown hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
336c4762a1bSJed Brown hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
337c4762a1bSJed Brown ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
338c4762a1bSJed Brown hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
339c4762a1bSJed Brown
340c4762a1bSJed Brown hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
341c4762a1bSJed Brown htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
342c4762a1bSJed Brown
3439371c9d4SSatish 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);
344c4762a1bSJed Brown
3459371c9d4SSatish Balay hl /= 2.0;
3469371c9d4SSatish Balay hr /= 2.0;
3479371c9d4SSatish Balay ht /= 2.0;
3489371c9d4SSatish Balay hb /= 2.0;
3499371c9d4SSatish Balay hbr /= 2.0;
3509371c9d4SSatish Balay htl /= 2.0;
3519371c9d4SSatish Balay hc /= 2.0;
352c4762a1bSJed Brown
353c4762a1bSJed Brown k = 0;
354c4762a1bSJed Brown if (j > 0) {
3559371c9d4SSatish Balay v[k] = hb;
3569371c9d4SSatish Balay col[k] = row - mx;
3579371c9d4SSatish Balay k++;
358c4762a1bSJed Brown }
359c4762a1bSJed Brown
360c4762a1bSJed Brown if (j > 0 && i < mx - 1) {
3619371c9d4SSatish Balay v[k] = hbr;
3629371c9d4SSatish Balay col[k] = row - mx + 1;
3639371c9d4SSatish Balay k++;
364c4762a1bSJed Brown }
365c4762a1bSJed Brown
366c4762a1bSJed Brown if (i > 0) {
3679371c9d4SSatish Balay v[k] = hl;
3689371c9d4SSatish Balay col[k] = row - 1;
3699371c9d4SSatish Balay k++;
370c4762a1bSJed Brown }
371c4762a1bSJed Brown
3729371c9d4SSatish Balay v[k] = hc;
3739371c9d4SSatish Balay col[k] = row;
3749371c9d4SSatish Balay k++;
375c4762a1bSJed Brown
376c4762a1bSJed Brown if (i < mx - 1) {
3779371c9d4SSatish Balay v[k] = hr;
3789371c9d4SSatish Balay col[k] = row + 1;
3799371c9d4SSatish Balay k++;
380c4762a1bSJed Brown }
381c4762a1bSJed Brown
382c4762a1bSJed Brown if (i > 0 && j < my - 1) {
3839371c9d4SSatish Balay v[k] = htl;
3849371c9d4SSatish Balay col[k] = row + mx - 1;
3859371c9d4SSatish Balay k++;
386c4762a1bSJed Brown }
387c4762a1bSJed Brown
388c4762a1bSJed Brown if (j < my - 1) {
3899371c9d4SSatish Balay v[k] = ht;
3909371c9d4SSatish Balay col[k] = row + mx;
3919371c9d4SSatish Balay k++;
392c4762a1bSJed Brown }
393c4762a1bSJed Brown
394c4762a1bSJed Brown /*
395c4762a1bSJed Brown Set matrix values using local numbering, which was defined
396c4762a1bSJed Brown earlier, in the main routine.
397c4762a1bSJed Brown */
3989566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 1, &row, k, col, v, INSERT_VALUES));
399c4762a1bSJed Brown }
400c4762a1bSJed Brown }
401c4762a1bSJed Brown
402c4762a1bSJed Brown /* Restore vectors */
4039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
404c4762a1bSJed Brown
405c4762a1bSJed Brown /* Assemble the matrix */
4069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
4079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
4089566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(199 * mx * my));
4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
410c4762a1bSJed Brown }
411c4762a1bSJed Brown
412c4762a1bSJed Brown /* ------------------------------------------------------------------- */
413c4762a1bSJed Brown /*
414c4762a1bSJed Brown MSA_BoundaryConditions - Calculates the boundary conditions for
415c4762a1bSJed Brown the region.
416c4762a1bSJed Brown
417c4762a1bSJed Brown Input Parameter:
418c4762a1bSJed Brown . user - user-defined application context
419c4762a1bSJed Brown
420c4762a1bSJed Brown Output Parameter:
421c4762a1bSJed Brown . user - user-defined application context
422c4762a1bSJed Brown */
MSA_BoundaryConditions(AppCtx * user)423d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_BoundaryConditions(AppCtx *user)
424d71ae5a4SJacob Faibussowitsch {
425c4762a1bSJed Brown PetscInt i, j, k, limit = 0, maxits = 5;
426c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
427c4762a1bSJed Brown PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0;
428c4762a1bSJed Brown PetscReal one = 1.0, two = 2.0, three = 3.0, tol = 1e-10;
429c4762a1bSJed Brown PetscReal fnorm, det, hx, hy, xt = 0, yt = 0;
430c4762a1bSJed Brown PetscReal u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
431c4762a1bSJed Brown PetscReal b = -0.5, t = 0.5, l = -0.5, r = 0.5;
432c4762a1bSJed Brown PetscReal *boundary;
433c4762a1bSJed Brown
434c4762a1bSJed Brown PetscFunctionBegin;
4359371c9d4SSatish Balay bsize = mx + 2;
4369371c9d4SSatish Balay lsize = my + 2;
4379371c9d4SSatish Balay rsize = my + 2;
4389371c9d4SSatish Balay tsize = mx + 2;
439c4762a1bSJed Brown
4409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bsize, &user->bottom));
4419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tsize, &user->top));
4429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &user->left));
4439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rsize, &user->right));
444c4762a1bSJed Brown
4459371c9d4SSatish Balay hx = (r - l) / (mx + 1);
4469371c9d4SSatish Balay hy = (t - b) / (my + 1);
447c4762a1bSJed Brown
448c4762a1bSJed Brown for (j = 0; j < 4; j++) {
449c4762a1bSJed Brown if (j == 0) {
450c4762a1bSJed Brown yt = b;
451c4762a1bSJed Brown xt = l;
452c4762a1bSJed Brown limit = bsize;
453c4762a1bSJed Brown boundary = user->bottom;
454c4762a1bSJed Brown } else if (j == 1) {
455c4762a1bSJed Brown yt = t;
456c4762a1bSJed Brown xt = l;
457c4762a1bSJed Brown limit = tsize;
458c4762a1bSJed Brown boundary = user->top;
459c4762a1bSJed Brown } else if (j == 2) {
460c4762a1bSJed Brown yt = b;
461c4762a1bSJed Brown xt = l;
462c4762a1bSJed Brown limit = lsize;
463c4762a1bSJed Brown boundary = user->left;
464c4762a1bSJed Brown } else { /* if (j==3) */
465c4762a1bSJed Brown yt = b;
466c4762a1bSJed Brown xt = r;
467c4762a1bSJed Brown limit = rsize;
468c4762a1bSJed Brown boundary = user->right;
469c4762a1bSJed Brown }
470c4762a1bSJed Brown
471c4762a1bSJed Brown for (i = 0; i < limit; i++) {
472c4762a1bSJed Brown u1 = xt;
473c4762a1bSJed Brown u2 = -yt;
474c4762a1bSJed Brown for (k = 0; k < maxits; k++) {
475c4762a1bSJed Brown nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
476c4762a1bSJed Brown nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
477c4762a1bSJed Brown fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2);
478c4762a1bSJed Brown if (fnorm <= tol) break;
479c4762a1bSJed Brown njac11 = one + u2 * u2 - u1 * u1;
480c4762a1bSJed Brown njac12 = two * u1 * u2;
481c4762a1bSJed Brown njac21 = -two * u1 * u2;
482c4762a1bSJed Brown njac22 = -one - u1 * u1 + u2 * u2;
483c4762a1bSJed Brown det = njac11 * njac22 - njac21 * njac12;
484c4762a1bSJed Brown u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det;
485c4762a1bSJed Brown u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det;
486c4762a1bSJed Brown }
487c4762a1bSJed Brown
488c4762a1bSJed Brown boundary[i] = u1 * u1 - u2 * u2;
489c4762a1bSJed Brown if (j == 0 || j == 1) {
490c4762a1bSJed Brown xt = xt + hx;
491c4762a1bSJed Brown } else { /* if (j==2 || j==3) */
492c4762a1bSJed Brown yt = yt + hy;
493c4762a1bSJed Brown }
494c4762a1bSJed Brown }
495c4762a1bSJed Brown }
4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
497c4762a1bSJed Brown }
498c4762a1bSJed Brown
499c4762a1bSJed Brown /* ------------------------------------------------------------------- */
500c4762a1bSJed Brown /*
501c4762a1bSJed Brown MSA_InitialPoint - Calculates the initial guess in one of three ways.
502c4762a1bSJed Brown
503c4762a1bSJed Brown Input Parameters:
504c4762a1bSJed Brown . user - user-defined application context
505c4762a1bSJed Brown . X - vector for initial guess
506c4762a1bSJed Brown
507c4762a1bSJed Brown Output Parameters:
508c4762a1bSJed Brown . X - newly computed initial guess
509c4762a1bSJed Brown */
MSA_InitialPoint(AppCtx * user,Vec X)510d71ae5a4SJacob Faibussowitsch static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X)
511d71ae5a4SJacob Faibussowitsch {
512c4762a1bSJed Brown PetscInt start = -1, i, j;
513c4762a1bSJed Brown PetscScalar zero = 0.0;
514c4762a1bSJed Brown PetscBool flg;
515c4762a1bSJed Brown
516c4762a1bSJed Brown PetscFunctionBegin;
5179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg));
518c4762a1bSJed Brown
519c4762a1bSJed Brown if (flg && start == 0) { /* The zero vector is reasonable */
5209566063dSJacob Faibussowitsch PetscCall(VecSet(X, zero));
521c4762a1bSJed Brown } else { /* Take an average of the boundary conditions */
522c4762a1bSJed Brown PetscInt row;
523c4762a1bSJed Brown PetscInt mx = user->mx, my = user->my;
524c4762a1bSJed Brown PetscScalar *x;
525c4762a1bSJed Brown
526c4762a1bSJed Brown /* Get pointers to vector data */
5279566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
528c4762a1bSJed Brown
529c4762a1bSJed Brown /* Perform local computations */
530c4762a1bSJed Brown for (j = 0; j < my; j++) {
531c4762a1bSJed Brown for (i = 0; i < mx; i++) {
532c4762a1bSJed Brown row = (j)*mx + (i);
533c4762a1bSJed 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;
534c4762a1bSJed Brown }
535c4762a1bSJed Brown }
536c4762a1bSJed Brown
537c4762a1bSJed Brown /* Restore vectors */
5389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
539c4762a1bSJed Brown }
5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
541c4762a1bSJed Brown }
542c4762a1bSJed Brown
543c4762a1bSJed Brown /*TEST
544c4762a1bSJed Brown
545c4762a1bSJed Brown build:
546c4762a1bSJed Brown requires: !complex
547c4762a1bSJed Brown
548c4762a1bSJed Brown test:
549c4762a1bSJed Brown args: -tao_monitor -tao_view -tao_type ssils -tao_gttol 1.e-5
550c4762a1bSJed Brown requires: !single
551c4762a1bSJed Brown
552c4762a1bSJed Brown test:
553c4762a1bSJed Brown suffix: 2
554c4762a1bSJed Brown args: -tao_monitor -tao_view -tao_type ssfls -tao_gttol 1.e-5
555c4762a1bSJed Brown
556c4762a1bSJed Brown TEST*/
557