1777c4855SStefano Zampini #include <petscsnes.h>
2777c4855SStefano Zampini #include <petscdmda.h>
3777c4855SStefano Zampini
4777c4855SStefano Zampini static const char help[] = "Minimum surface area problem in 2D using DMDA.\n\
5777c4855SStefano Zampini It solves an unconstrained minimization problem. This example is based on a \n\
6777c4855SStefano Zampini problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\
7777c4855SStefano Zampini boundary values along the edges of the domain, the objective is to find the\n\
8777c4855SStefano Zampini surface with the minimal area that satisfies the boundary conditions.\n\
9777c4855SStefano Zampini \n\
10777c4855SStefano Zampini The command line options are:\n\
11777c4855SStefano Zampini -da_grid_x <nx>, where <nx> = number of grid points in the 1st coordinate direction\n\
12777c4855SStefano Zampini -da_grid_y <ny>, where <ny> = number of grid points in the 2nd coordinate direction\n\
13777c4855SStefano Zampini \n";
14777c4855SStefano Zampini
15777c4855SStefano Zampini /*
16777c4855SStefano Zampini User-defined application context - contains data needed by the
17777c4855SStefano Zampini application-provided call-back routines, FormJacobianLocal() and
18777c4855SStefano Zampini FormFunctionLocal().
19777c4855SStefano Zampini */
20777c4855SStefano Zampini
21777c4855SStefano Zampini typedef enum {
22777c4855SStefano Zampini PROBLEM_ENNEPER,
23777c4855SStefano Zampini PROBLEM_SINS,
24777c4855SStefano Zampini } ProblemType;
25777c4855SStefano Zampini static const char *const ProblemTypes[] = {"ENNEPER", "SINS", "ProblemType", "PROBLEM_", 0};
26777c4855SStefano Zampini
27777c4855SStefano Zampini typedef struct {
28777c4855SStefano Zampini PetscScalar *bottom, *top, *left, *right;
29777c4855SStefano Zampini } AppCtx;
30777c4855SStefano Zampini
31777c4855SStefano Zampini /* -------- User-defined Routines --------- */
32777c4855SStefano Zampini
33777c4855SStefano Zampini static PetscErrorCode FormBoundaryConditions_Enneper(SNES, AppCtx **);
34777c4855SStefano Zampini static PetscErrorCode FormBoundaryConditions_Sins(SNES, AppCtx **);
35777c4855SStefano Zampini static PetscErrorCode DestroyBoundaryConditions(AppCtx **);
36777c4855SStefano Zampini static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *, PetscScalar **, PetscReal *, void *);
37777c4855SStefano Zampini static PetscErrorCode FormFunctionLocal(DMDALocalInfo *, PetscScalar **, PetscScalar **, void *);
38777c4855SStefano Zampini static PetscErrorCode FormJacobianLocal(DMDALocalInfo *, PetscScalar **, Mat, Mat, void *);
39777c4855SStefano Zampini
main(int argc,char ** argv)40777c4855SStefano Zampini int main(int argc, char **argv)
41777c4855SStefano Zampini {
42777c4855SStefano Zampini Vec x;
43777c4855SStefano Zampini SNES snes;
44777c4855SStefano Zampini DM da;
45777c4855SStefano Zampini ProblemType ptype = PROBLEM_ENNEPER;
46777c4855SStefano Zampini PetscBool use_obj = PETSC_TRUE;
47777c4855SStefano Zampini PetscReal bbox[4] = {0.};
48777c4855SStefano Zampini AppCtx *user;
49777c4855SStefano Zampini PetscErrorCode (*form_bc)(SNES, AppCtx **) = NULL;
50777c4855SStefano Zampini
51777c4855SStefano Zampini PetscFunctionBeginUser;
52*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
53777c4855SStefano Zampini
54777c4855SStefano Zampini PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Minimal surface options", __FILE__);
55777c4855SStefano Zampini PetscCall(PetscOptionsEnum("-problem_type", "Problem type", NULL, ProblemTypes, (PetscEnum)ptype, (PetscEnum *)&ptype, NULL));
56777c4855SStefano Zampini PetscCall(PetscOptionsBool("-use_objective", "Use objective function", NULL, use_obj, &use_obj, NULL));
57777c4855SStefano Zampini PetscOptionsEnd();
58777c4855SStefano Zampini switch (ptype) {
59777c4855SStefano Zampini case PROBLEM_ENNEPER:
60777c4855SStefano Zampini bbox[0] = -0.5;
61777c4855SStefano Zampini bbox[1] = 0.5;
62777c4855SStefano Zampini bbox[2] = -0.5;
63777c4855SStefano Zampini bbox[3] = 0.5;
64777c4855SStefano Zampini form_bc = FormBoundaryConditions_Enneper;
65777c4855SStefano Zampini break;
66777c4855SStefano Zampini case PROBLEM_SINS:
67777c4855SStefano Zampini bbox[0] = 0.0;
68777c4855SStefano Zampini bbox[1] = 1.0;
69777c4855SStefano Zampini bbox[2] = 0.0;
70777c4855SStefano Zampini bbox[3] = 1.0;
71777c4855SStefano Zampini form_bc = FormBoundaryConditions_Sins;
72777c4855SStefano Zampini break;
73777c4855SStefano Zampini }
74777c4855SStefano Zampini
75777c4855SStefano Zampini /* Create distributed array to manage the 2d grid */
76777c4855SStefano Zampini PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
77777c4855SStefano Zampini PetscCall(DMSetFromOptions(da));
78777c4855SStefano Zampini PetscCall(DMSetUp(da));
79777c4855SStefano Zampini PetscCall(DMDASetUniformCoordinates(da, bbox[0], bbox[1], bbox[2], bbox[3], PETSC_DECIDE, PETSC_DECIDE));
80777c4855SStefano Zampini
81777c4855SStefano Zampini /* Extract global vectors from DMDA; */
82777c4855SStefano Zampini PetscCall(DMCreateGlobalVector(da, &x));
83777c4855SStefano Zampini
84777c4855SStefano Zampini /* Create nonlinear solver context */
85777c4855SStefano Zampini PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
86777c4855SStefano Zampini PetscCall(SNESSetDM(snes, da));
87777c4855SStefano Zampini PetscCall((*form_bc)(snes, &user));
88777c4855SStefano Zampini PetscCall(SNESSetApplicationContext(snes, user));
89777c4855SStefano Zampini
90777c4855SStefano Zampini /* Set local callbacks */
918434afd1SBarry Smith if (use_obj) PetscCall(DMDASNESSetObjectiveLocal(da, (DMDASNESObjectiveFn *)FormObjectiveLocal, NULL));
928434afd1SBarry Smith PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunctionFn *)FormFunctionLocal, NULL));
938434afd1SBarry Smith PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobianFn *)FormJacobianLocal, NULL));
94777c4855SStefano Zampini
95777c4855SStefano Zampini /* Customize from command line */
96777c4855SStefano Zampini PetscCall(SNESSetFromOptions(snes));
97777c4855SStefano Zampini
98777c4855SStefano Zampini /* Solve the application */
99777c4855SStefano Zampini PetscCall(SNESSolve(snes, NULL, x));
100777c4855SStefano Zampini
101777c4855SStefano Zampini /* Free user-created data structures */
102777c4855SStefano Zampini PetscCall(VecDestroy(&x));
103777c4855SStefano Zampini PetscCall(SNESDestroy(&snes));
104777c4855SStefano Zampini PetscCall(DMDestroy(&da));
105777c4855SStefano Zampini PetscCall(DestroyBoundaryConditions(&user));
106777c4855SStefano Zampini
107777c4855SStefano Zampini PetscCall(PetscFinalize());
108777c4855SStefano Zampini return 0;
109777c4855SStefano Zampini }
110777c4855SStefano Zampini
111777c4855SStefano Zampini /* Compute objective function over the locally owned part of the mesh */
FormObjectiveLocal(DMDALocalInfo * info,PetscScalar ** x,PetscReal * v,void * ptr)112777c4855SStefano Zampini PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info, PetscScalar **x, PetscReal *v, void *ptr)
113777c4855SStefano Zampini {
114777c4855SStefano Zampini AppCtx *user = (AppCtx *)ptr;
115777c4855SStefano Zampini PetscInt mx = info->mx, my = info->my;
116777c4855SStefano Zampini PetscInt xs = info->xs, xm = info->xm, ys = info->ys, ym = info->ym;
117777c4855SStefano Zampini PetscInt i, j;
118777c4855SStefano Zampini PetscScalar hx, hy;
119777c4855SStefano Zampini PetscScalar f2, f4, d1, d2, d3, d4, xc, xl, xr, xt, xb;
120777c4855SStefano Zampini PetscReal ft = 0, area;
121777c4855SStefano Zampini
122777c4855SStefano Zampini PetscFunctionBeginUser;
123777c4855SStefano Zampini PetscCheck(user, PetscObjectComm((PetscObject)info->da), PETSC_ERR_PLIB, "Missing application context");
124777c4855SStefano Zampini hx = 1.0 / (mx + 1);
125777c4855SStefano Zampini hy = 1.0 / (my + 1);
126777c4855SStefano Zampini area = 0.5 * hx * hy;
127777c4855SStefano Zampini for (j = ys; j < ys + ym; j++) {
128777c4855SStefano Zampini for (i = xs; i < xs + xm; i++) {
129777c4855SStefano Zampini xc = x[j][i];
130777c4855SStefano Zampini xl = xr = xb = xt = xc;
131777c4855SStefano Zampini
132777c4855SStefano Zampini if (i == 0) { /* left side */
133777c4855SStefano Zampini xl = user->left[j + 1];
134777c4855SStefano Zampini } else xl = x[j][i - 1];
135777c4855SStefano Zampini
136777c4855SStefano Zampini if (j == 0) { /* bottom side */
137777c4855SStefano Zampini xb = user->bottom[i + 1];
138777c4855SStefano Zampini } else xb = x[j - 1][i];
139777c4855SStefano Zampini
140777c4855SStefano Zampini if (i + 1 == mx) { /* right side */
141777c4855SStefano Zampini xr = user->right[j + 1];
142777c4855SStefano Zampini } else xr = x[j][i + 1];
143777c4855SStefano Zampini
144777c4855SStefano Zampini if (j + 1 == 0 + my) { /* top side */
145777c4855SStefano Zampini xt = user->top[i + 1];
146777c4855SStefano Zampini } else xt = x[j + 1][i];
147777c4855SStefano Zampini
148777c4855SStefano Zampini d1 = (xc - xl);
149777c4855SStefano Zampini d2 = (xc - xr);
150777c4855SStefano Zampini d3 = (xc - xt);
151777c4855SStefano Zampini d4 = (xc - xb);
152777c4855SStefano Zampini
153777c4855SStefano Zampini d1 /= hx;
154777c4855SStefano Zampini d2 /= hx;
155777c4855SStefano Zampini d3 /= hy;
156777c4855SStefano Zampini d4 /= hy;
157777c4855SStefano Zampini
158777c4855SStefano Zampini f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
159777c4855SStefano Zampini f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
160777c4855SStefano Zampini
161777c4855SStefano Zampini ft += PetscRealPart(f2 + f4);
162777c4855SStefano Zampini }
163777c4855SStefano Zampini }
164777c4855SStefano Zampini
165777c4855SStefano Zampini /* Compute triangular areas along the border of the domain. */
166777c4855SStefano Zampini if (xs == 0) { /* left side */
167777c4855SStefano Zampini for (j = ys; j < ys + ym; j++) {
168777c4855SStefano Zampini d3 = (user->left[j + 1] - user->left[j + 2]) / hy;
169777c4855SStefano Zampini d2 = (user->left[j + 1] - x[j][0]) / hx;
170777c4855SStefano Zampini ft += PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
171777c4855SStefano Zampini }
172777c4855SStefano Zampini }
173777c4855SStefano Zampini if (ys == 0) { /* bottom side */
174777c4855SStefano Zampini for (i = xs; i < xs + xm; i++) {
175777c4855SStefano Zampini d2 = (user->bottom[i + 1] - user->bottom[i + 2]) / hx;
176777c4855SStefano Zampini d3 = (user->bottom[i + 1] - x[0][i]) / hy;
177777c4855SStefano Zampini ft += PetscSqrtReal(1.0 + d3 * d3 + d2 * d2);
178777c4855SStefano Zampini }
179777c4855SStefano Zampini }
180777c4855SStefano Zampini if (xs + xm == mx) { /* right side */
181777c4855SStefano Zampini for (j = ys; j < ys + ym; j++) {
182777c4855SStefano Zampini d1 = (x[j][mx - 1] - user->right[j + 1]) / hx;
183777c4855SStefano Zampini d4 = (user->right[j] - user->right[j + 1]) / hy;
184777c4855SStefano Zampini ft += PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
185777c4855SStefano Zampini }
186777c4855SStefano Zampini }
187777c4855SStefano Zampini if (ys + ym == my) { /* top side */
188777c4855SStefano Zampini for (i = xs; i < xs + xm; i++) {
189777c4855SStefano Zampini d1 = (x[my - 1][i] - user->top[i + 1]) / hy;
190777c4855SStefano Zampini d4 = (user->top[i + 1] - user->top[i]) / hx;
191777c4855SStefano Zampini ft += PetscSqrtReal(1.0 + d1 * d1 + d4 * d4);
192777c4855SStefano Zampini }
193777c4855SStefano Zampini }
194777c4855SStefano Zampini if (ys == 0 && xs == 0) {
195777c4855SStefano Zampini d1 = (user->left[0] - user->left[1]) / hy;
196777c4855SStefano Zampini d2 = (user->bottom[0] - user->bottom[1]) / hx;
197777c4855SStefano Zampini ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
198777c4855SStefano Zampini }
199777c4855SStefano Zampini if (ys + ym == my && xs + xm == mx) {
200777c4855SStefano Zampini d1 = (user->right[ym + 1] - user->right[ym]) / hy;
201777c4855SStefano Zampini d2 = (user->top[xm + 1] - user->top[xm]) / hx;
202777c4855SStefano Zampini ft += PetscSqrtReal(1.0 + d1 * d1 + d2 * d2);
203777c4855SStefano Zampini }
204777c4855SStefano Zampini ft *= area;
205777c4855SStefano Zampini *v = ft;
206777c4855SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
207777c4855SStefano Zampini }
208777c4855SStefano Zampini
209777c4855SStefano Zampini /* Compute gradient over the locally owned part of the mesh */
FormFunctionLocal(DMDALocalInfo * info,PetscScalar ** x,PetscScalar ** g,void * ptr)210777c4855SStefano Zampini PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **g, void *ptr)
211777c4855SStefano Zampini {
212777c4855SStefano Zampini AppCtx *user = (AppCtx *)ptr;
213777c4855SStefano Zampini PetscInt mx = info->mx, my = info->my;
214777c4855SStefano Zampini PetscInt xs = info->xs, xm = info->xm, ys = info->ys, ym = info->ym;
215777c4855SStefano Zampini PetscInt i, j;
216777c4855SStefano Zampini PetscScalar hx, hy, hydhx, hxdhy;
217777c4855SStefano Zampini PetscScalar f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
218777c4855SStefano Zampini PetscScalar df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc;
219777c4855SStefano Zampini
220777c4855SStefano Zampini PetscFunctionBeginUser;
221777c4855SStefano Zampini PetscCheck(user, PetscObjectComm((PetscObject)info->da), PETSC_ERR_PLIB, "Missing application context");
222777c4855SStefano Zampini hx = 1.0 / (mx + 1);
223777c4855SStefano Zampini hy = 1.0 / (my + 1);
224777c4855SStefano Zampini hydhx = hy / hx;
225777c4855SStefano Zampini hxdhy = hx / hy;
226777c4855SStefano Zampini
227777c4855SStefano Zampini for (j = ys; j < ys + ym; j++) {
228777c4855SStefano Zampini for (i = xs; i < xs + xm; i++) {
229777c4855SStefano Zampini xc = x[j][i];
230777c4855SStefano Zampini xlt = xrb = xl = xr = xb = xt = xc;
231777c4855SStefano Zampini
232777c4855SStefano Zampini if (i == 0) { /* left side */
233777c4855SStefano Zampini xl = user->left[j + 1];
234777c4855SStefano Zampini xlt = user->left[j + 2];
235777c4855SStefano Zampini } else xl = x[j][i - 1];
236777c4855SStefano Zampini
237777c4855SStefano Zampini if (j == 0) { /* bottom side */
238777c4855SStefano Zampini xb = user->bottom[i + 1];
239777c4855SStefano Zampini xrb = user->bottom[i + 2];
240777c4855SStefano Zampini } else xb = x[j - 1][i];
241777c4855SStefano Zampini
242777c4855SStefano Zampini if (i + 1 == mx) { /* right side */
243777c4855SStefano Zampini xr = user->right[j + 1];
244777c4855SStefano Zampini xrb = user->right[j];
245777c4855SStefano Zampini } else xr = x[j][i + 1];
246777c4855SStefano Zampini
247777c4855SStefano Zampini if (j + 1 == 0 + my) { /* top side */
248777c4855SStefano Zampini xt = user->top[i + 1];
249777c4855SStefano Zampini xlt = user->top[i];
250777c4855SStefano Zampini } else xt = x[j + 1][i];
251777c4855SStefano Zampini
252777c4855SStefano Zampini if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1]; /* left top side */
253777c4855SStefano Zampini if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1]; /* right bottom */
254777c4855SStefano Zampini
255777c4855SStefano Zampini d1 = (xc - xl);
256777c4855SStefano Zampini d2 = (xc - xr);
257777c4855SStefano Zampini d3 = (xc - xt);
258777c4855SStefano Zampini d4 = (xc - xb);
259777c4855SStefano Zampini d5 = (xr - xrb);
260777c4855SStefano Zampini d6 = (xrb - xb);
261777c4855SStefano Zampini d7 = (xlt - xl);
262777c4855SStefano Zampini d8 = (xt - xlt);
263777c4855SStefano Zampini
264777c4855SStefano Zampini df1dxc = d1 * hydhx;
265777c4855SStefano Zampini df2dxc = (d1 * hydhx + d4 * hxdhy);
266777c4855SStefano Zampini df3dxc = d3 * hxdhy;
267777c4855SStefano Zampini df4dxc = (d2 * hydhx + d3 * hxdhy);
268777c4855SStefano Zampini df5dxc = d2 * hydhx;
269777c4855SStefano Zampini df6dxc = d4 * hxdhy;
270777c4855SStefano Zampini
271777c4855SStefano Zampini d1 /= hx;
272777c4855SStefano Zampini d2 /= hx;
273777c4855SStefano Zampini d3 /= hy;
274777c4855SStefano Zampini d4 /= hy;
275777c4855SStefano Zampini d5 /= hy;
276777c4855SStefano Zampini d6 /= hx;
277777c4855SStefano Zampini d7 /= hy;
278777c4855SStefano Zampini d8 /= hx;
279777c4855SStefano Zampini
280777c4855SStefano Zampini f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
281777c4855SStefano Zampini f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
282777c4855SStefano Zampini f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
283777c4855SStefano Zampini f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
284777c4855SStefano Zampini f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
285777c4855SStefano Zampini f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
286777c4855SStefano Zampini
287777c4855SStefano Zampini df1dxc /= f1;
288777c4855SStefano Zampini df2dxc /= f2;
289777c4855SStefano Zampini df3dxc /= f3;
290777c4855SStefano Zampini df4dxc /= f4;
291777c4855SStefano Zampini df5dxc /= f5;
292777c4855SStefano Zampini df6dxc /= f6;
293777c4855SStefano Zampini
294777c4855SStefano Zampini g[j][i] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) / 2.0;
295777c4855SStefano Zampini }
296777c4855SStefano Zampini }
297777c4855SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
298777c4855SStefano Zampini }
299777c4855SStefano Zampini
300777c4855SStefano Zampini /* Compute Hessian over the locally owned part of the mesh */
FormJacobianLocal(DMDALocalInfo * info,PetscScalar ** x,Mat H,Mat Hp,void * ptr)301777c4855SStefano Zampini PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat H, Mat Hp, void *ptr)
302777c4855SStefano Zampini {
303777c4855SStefano Zampini AppCtx *user = (AppCtx *)ptr;
304777c4855SStefano Zampini PetscInt mx = info->mx, my = info->my;
305777c4855SStefano Zampini PetscInt xs = info->xs, xm = info->xm, ys = info->ys, ym = info->ym;
306777c4855SStefano Zampini PetscInt i, j, k;
307777c4855SStefano Zampini MatStencil row, col[7];
308777c4855SStefano Zampini PetscScalar hx, hy, hydhx, hxdhy;
309777c4855SStefano Zampini PetscScalar f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb;
310777c4855SStefano Zampini PetscScalar hl, hr, ht, hb, hc, htl, hbr;
311777c4855SStefano Zampini PetscScalar v[7];
312777c4855SStefano Zampini
313777c4855SStefano Zampini PetscFunctionBeginUser;
314777c4855SStefano Zampini PetscCheck(user, PetscObjectComm((PetscObject)info->da), PETSC_ERR_PLIB, "Missing application context");
315777c4855SStefano Zampini hx = 1.0 / (mx + 1);
316777c4855SStefano Zampini hy = 1.0 / (my + 1);
317777c4855SStefano Zampini hydhx = hy / hx;
318777c4855SStefano Zampini hxdhy = hx / hy;
319777c4855SStefano Zampini
320777c4855SStefano Zampini for (j = ys; j < ys + ym; j++) {
321777c4855SStefano Zampini for (i = xs; i < xs + xm; i++) {
322777c4855SStefano Zampini xc = x[j][i];
323777c4855SStefano Zampini xlt = xrb = xl = xr = xb = xt = xc;
324777c4855SStefano Zampini
325777c4855SStefano Zampini /* Left */
326777c4855SStefano Zampini if (i == 0) {
327777c4855SStefano Zampini xl = user->left[j + 1];
328777c4855SStefano Zampini xlt = user->left[j + 2];
329777c4855SStefano Zampini } else xl = x[j][i - 1];
330777c4855SStefano Zampini
331777c4855SStefano Zampini /* Bottom */
332777c4855SStefano Zampini if (j == 0) {
333777c4855SStefano Zampini xb = user->bottom[i + 1];
334777c4855SStefano Zampini xrb = user->bottom[i + 2];
335777c4855SStefano Zampini } else xb = x[j - 1][i];
336777c4855SStefano Zampini
337777c4855SStefano Zampini /* Right */
338777c4855SStefano Zampini if (i + 1 == mx) {
339777c4855SStefano Zampini xr = user->right[j + 1];
340777c4855SStefano Zampini xrb = user->right[j];
341777c4855SStefano Zampini } else xr = x[j][i + 1];
342777c4855SStefano Zampini
343777c4855SStefano Zampini /* Top */
344777c4855SStefano Zampini if (j + 1 == my) {
345777c4855SStefano Zampini xt = user->top[i + 1];
346777c4855SStefano Zampini xlt = user->top[i];
347777c4855SStefano Zampini } else xt = x[j + 1][i];
348777c4855SStefano Zampini
349777c4855SStefano Zampini /* Top left */
350777c4855SStefano Zampini if (i > 0 && j + 1 < my) xlt = x[j + 1][i - 1];
351777c4855SStefano Zampini
352777c4855SStefano Zampini /* Bottom right */
353777c4855SStefano Zampini if (j > 0 && i + 1 < mx) xrb = x[j - 1][i + 1];
354777c4855SStefano Zampini
355777c4855SStefano Zampini d1 = (xc - xl) / hx;
356777c4855SStefano Zampini d2 = (xc - xr) / hx;
357777c4855SStefano Zampini d3 = (xc - xt) / hy;
358777c4855SStefano Zampini d4 = (xc - xb) / hy;
359777c4855SStefano Zampini d5 = (xrb - xr) / hy;
360777c4855SStefano Zampini d6 = (xrb - xb) / hx;
361777c4855SStefano Zampini d7 = (xlt - xl) / hy;
362777c4855SStefano Zampini d8 = (xlt - xt) / hx;
363777c4855SStefano Zampini
364777c4855SStefano Zampini f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7);
365777c4855SStefano Zampini f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4);
366777c4855SStefano Zampini f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8);
367777c4855SStefano Zampini f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2);
368777c4855SStefano Zampini f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5);
369777c4855SStefano Zampini f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6);
370777c4855SStefano Zampini
371777c4855SStefano Zampini hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2);
372777c4855SStefano Zampini hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4);
373777c4855SStefano Zampini ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4);
374777c4855SStefano Zampini hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2);
375777c4855SStefano Zampini
376777c4855SStefano Zampini hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6);
377777c4855SStefano Zampini htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3);
378777c4855SStefano Zampini
379777c4855SStefano Zampini 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.0 * d1 * d4) / (f2 * f2 * f2) + (hxdhy * (1.0 + d2 * d2) + hydhx * (1.0 + d3 * d3) - 2.0 * d2 * d3) / (f4 * f4 * f4);
380777c4855SStefano Zampini
381777c4855SStefano Zampini hl /= 2.0;
382777c4855SStefano Zampini hr /= 2.0;
383777c4855SStefano Zampini ht /= 2.0;
384777c4855SStefano Zampini hb /= 2.0;
385777c4855SStefano Zampini hbr /= 2.0;
386777c4855SStefano Zampini htl /= 2.0;
387777c4855SStefano Zampini hc /= 2.0;
388777c4855SStefano Zampini
389777c4855SStefano Zampini k = 0;
390777c4855SStefano Zampini row.i = i;
391777c4855SStefano Zampini row.j = j;
392777c4855SStefano Zampini /* Bottom */
393777c4855SStefano Zampini if (j > 0) {
394777c4855SStefano Zampini v[k] = hb;
395777c4855SStefano Zampini col[k].i = i;
396777c4855SStefano Zampini col[k].j = j - 1;
397777c4855SStefano Zampini k++;
398777c4855SStefano Zampini }
399777c4855SStefano Zampini
400777c4855SStefano Zampini /* Bottom right */
401777c4855SStefano Zampini if (j > 0 && i < mx - 1) {
402777c4855SStefano Zampini v[k] = hbr;
403777c4855SStefano Zampini col[k].i = i + 1;
404777c4855SStefano Zampini col[k].j = j - 1;
405777c4855SStefano Zampini k++;
406777c4855SStefano Zampini }
407777c4855SStefano Zampini
408777c4855SStefano Zampini /* left */
409777c4855SStefano Zampini if (i > 0) {
410777c4855SStefano Zampini v[k] = hl;
411777c4855SStefano Zampini col[k].i = i - 1;
412777c4855SStefano Zampini col[k].j = j;
413777c4855SStefano Zampini k++;
414777c4855SStefano Zampini }
415777c4855SStefano Zampini
416777c4855SStefano Zampini /* Centre */
417777c4855SStefano Zampini v[k] = hc;
418777c4855SStefano Zampini col[k].i = row.i;
419777c4855SStefano Zampini col[k].j = row.j;
420777c4855SStefano Zampini k++;
421777c4855SStefano Zampini
422777c4855SStefano Zampini /* Right */
423777c4855SStefano Zampini if (i < mx - 1) {
424777c4855SStefano Zampini v[k] = hr;
425777c4855SStefano Zampini col[k].i = i + 1;
426777c4855SStefano Zampini col[k].j = j;
427777c4855SStefano Zampini k++;
428777c4855SStefano Zampini }
429777c4855SStefano Zampini
430777c4855SStefano Zampini /* Top left */
431777c4855SStefano Zampini if (i > 0 && j < my - 1) {
432777c4855SStefano Zampini v[k] = htl;
433777c4855SStefano Zampini col[k].i = i - 1;
434777c4855SStefano Zampini col[k].j = j + 1;
435777c4855SStefano Zampini k++;
436777c4855SStefano Zampini }
437777c4855SStefano Zampini
438777c4855SStefano Zampini /* Top */
439777c4855SStefano Zampini if (j < my - 1) {
440777c4855SStefano Zampini v[k] = ht;
441777c4855SStefano Zampini col[k].i = i;
442777c4855SStefano Zampini col[k].j = j + 1;
443777c4855SStefano Zampini k++;
444777c4855SStefano Zampini }
445777c4855SStefano Zampini
446777c4855SStefano Zampini PetscCall(MatSetValuesStencil(Hp, 1, &row, k, col, v, INSERT_VALUES));
447777c4855SStefano Zampini }
448777c4855SStefano Zampini }
449777c4855SStefano Zampini
450777c4855SStefano Zampini /* Assemble the matrix */
451777c4855SStefano Zampini PetscCall(MatAssemblyBegin(Hp, MAT_FINAL_ASSEMBLY));
452777c4855SStefano Zampini PetscCall(MatAssemblyEnd(Hp, MAT_FINAL_ASSEMBLY));
453777c4855SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
454777c4855SStefano Zampini }
455777c4855SStefano Zampini
FormBoundaryConditions_Enneper(SNES snes,AppCtx ** ouser)456777c4855SStefano Zampini PetscErrorCode FormBoundaryConditions_Enneper(SNES snes, AppCtx **ouser)
457777c4855SStefano Zampini {
458777c4855SStefano Zampini PetscInt i, j, k, limit = 0, maxits = 5;
459777c4855SStefano Zampini PetscInt mx, my;
460777c4855SStefano Zampini PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0;
461777c4855SStefano Zampini PetscScalar one = 1.0, two = 2.0, three = 3.0;
462777c4855SStefano Zampini PetscScalar det, hx, hy, xt = 0, yt = 0;
463777c4855SStefano Zampini PetscReal fnorm, tol = 1e-10;
464777c4855SStefano Zampini PetscScalar u1, u2, nf1, nf2, njac11, njac12, njac21, njac22;
465777c4855SStefano Zampini PetscScalar b = -0.5, t = 0.5, l = -0.5, r = 0.5;
466777c4855SStefano Zampini PetscScalar *boundary;
467777c4855SStefano Zampini AppCtx *user;
468777c4855SStefano Zampini DM da;
469777c4855SStefano Zampini
470777c4855SStefano Zampini PetscFunctionBeginUser;
471777c4855SStefano Zampini PetscCall(SNESGetDM(snes, &da));
472777c4855SStefano Zampini PetscCall(PetscNew(&user));
473777c4855SStefano Zampini *ouser = user;
474777c4855SStefano Zampini PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &mx, &my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
475777c4855SStefano Zampini bsize = mx + 2;
476777c4855SStefano Zampini lsize = my + 2;
477777c4855SStefano Zampini rsize = my + 2;
478777c4855SStefano Zampini tsize = mx + 2;
479777c4855SStefano Zampini
480777c4855SStefano Zampini PetscCall(PetscMalloc1(bsize, &user->bottom));
481777c4855SStefano Zampini PetscCall(PetscMalloc1(tsize, &user->top));
482777c4855SStefano Zampini PetscCall(PetscMalloc1(lsize, &user->left));
483777c4855SStefano Zampini PetscCall(PetscMalloc1(rsize, &user->right));
484777c4855SStefano Zampini
485777c4855SStefano Zampini hx = 1.0 / (mx + 1.0);
486777c4855SStefano Zampini hy = 1.0 / (my + 1.0);
487777c4855SStefano Zampini
488777c4855SStefano Zampini for (j = 0; j < 4; j++) {
489777c4855SStefano Zampini if (j == 0) {
490777c4855SStefano Zampini yt = b;
491777c4855SStefano Zampini xt = l;
492777c4855SStefano Zampini limit = bsize;
493777c4855SStefano Zampini boundary = user->bottom;
494777c4855SStefano Zampini } else if (j == 1) {
495777c4855SStefano Zampini yt = t;
496777c4855SStefano Zampini xt = l;
497777c4855SStefano Zampini limit = tsize;
498777c4855SStefano Zampini boundary = user->top;
499777c4855SStefano Zampini } else if (j == 2) {
500777c4855SStefano Zampini yt = b;
501777c4855SStefano Zampini xt = l;
502777c4855SStefano Zampini limit = lsize;
503777c4855SStefano Zampini boundary = user->left;
504777c4855SStefano Zampini } else { /* if (j==3) */
505777c4855SStefano Zampini yt = b;
506777c4855SStefano Zampini xt = r;
507777c4855SStefano Zampini limit = rsize;
508777c4855SStefano Zampini boundary = user->right;
509777c4855SStefano Zampini }
510777c4855SStefano Zampini
511777c4855SStefano Zampini for (i = 0; i < limit; i++) {
512777c4855SStefano Zampini u1 = xt;
513777c4855SStefano Zampini u2 = -yt;
514777c4855SStefano Zampini for (k = 0; k < maxits; k++) {
515777c4855SStefano Zampini nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt;
516777c4855SStefano Zampini nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt;
517777c4855SStefano Zampini fnorm = PetscRealPart(PetscSqrtScalar(nf1 * nf1 + nf2 * nf2));
518777c4855SStefano Zampini if (fnorm <= tol) break;
519777c4855SStefano Zampini njac11 = one + u2 * u2 - u1 * u1;
520777c4855SStefano Zampini njac12 = two * u1 * u2;
521777c4855SStefano Zampini njac21 = -two * u1 * u2;
522777c4855SStefano Zampini njac22 = -one - u1 * u1 + u2 * u2;
523777c4855SStefano Zampini det = njac11 * njac22 - njac21 * njac12;
524777c4855SStefano Zampini u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det;
525777c4855SStefano Zampini u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det;
526777c4855SStefano Zampini }
527777c4855SStefano Zampini
528777c4855SStefano Zampini boundary[i] = u1 * u1 - u2 * u2;
529777c4855SStefano Zampini if (j == 0 || j == 1) xt = xt + hx;
530777c4855SStefano Zampini else yt = yt + hy; /* if (j==2 || j==3) */
531777c4855SStefano Zampini }
532777c4855SStefano Zampini }
533777c4855SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
534777c4855SStefano Zampini }
535777c4855SStefano Zampini
DestroyBoundaryConditions(AppCtx ** ouser)536777c4855SStefano Zampini PetscErrorCode DestroyBoundaryConditions(AppCtx **ouser)
537777c4855SStefano Zampini {
538777c4855SStefano Zampini AppCtx *user = *ouser;
539777c4855SStefano Zampini
540777c4855SStefano Zampini PetscFunctionBeginUser;
541777c4855SStefano Zampini PetscCall(PetscFree(user->bottom));
542777c4855SStefano Zampini PetscCall(PetscFree(user->top));
543777c4855SStefano Zampini PetscCall(PetscFree(user->left));
544777c4855SStefano Zampini PetscCall(PetscFree(user->right));
545777c4855SStefano Zampini PetscCall(PetscFree(*ouser));
546777c4855SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
547777c4855SStefano Zampini }
548777c4855SStefano Zampini
FormBoundaryConditions_Sins(SNES snes,AppCtx ** ouser)549777c4855SStefano Zampini PetscErrorCode FormBoundaryConditions_Sins(SNES snes, AppCtx **ouser)
550777c4855SStefano Zampini {
551777c4855SStefano Zampini PetscInt i, j;
552777c4855SStefano Zampini PetscInt mx, my;
553777c4855SStefano Zampini PetscInt limit, bsize = 0, lsize = 0, tsize = 0, rsize = 0;
554777c4855SStefano Zampini PetscScalar hx, hy, xt = 0, yt = 0;
555777c4855SStefano Zampini PetscScalar b = 0.0, t = 1.0, l = 0.0, r = 1.0;
556777c4855SStefano Zampini PetscScalar *boundary;
557777c4855SStefano Zampini AppCtx *user;
558777c4855SStefano Zampini DM da;
559777c4855SStefano Zampini PetscReal pi2 = 2 * PETSC_PI;
560777c4855SStefano Zampini
561777c4855SStefano Zampini PetscFunctionBeginUser;
562777c4855SStefano Zampini PetscCall(SNESGetDM(snes, &da));
563777c4855SStefano Zampini PetscCall(PetscNew(&user));
564777c4855SStefano Zampini *ouser = user;
565777c4855SStefano Zampini PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &mx, &my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
566777c4855SStefano Zampini bsize = mx + 2;
567777c4855SStefano Zampini lsize = my + 2;
568777c4855SStefano Zampini rsize = my + 2;
569777c4855SStefano Zampini tsize = mx + 2;
570777c4855SStefano Zampini
571777c4855SStefano Zampini PetscCall(PetscMalloc1(bsize, &user->bottom));
572777c4855SStefano Zampini PetscCall(PetscMalloc1(tsize, &user->top));
573777c4855SStefano Zampini PetscCall(PetscMalloc1(lsize, &user->left));
574777c4855SStefano Zampini PetscCall(PetscMalloc1(rsize, &user->right));
575777c4855SStefano Zampini
576777c4855SStefano Zampini hx = 1.0 / (mx + 1.0);
577777c4855SStefano Zampini hy = 1.0 / (my + 1.0);
578777c4855SStefano Zampini
579777c4855SStefano Zampini for (j = 0; j < 4; j++) {
580777c4855SStefano Zampini if (j == 0) {
581777c4855SStefano Zampini yt = b;
582777c4855SStefano Zampini xt = l;
583777c4855SStefano Zampini limit = bsize;
584777c4855SStefano Zampini boundary = user->bottom;
585777c4855SStefano Zampini } else if (j == 1) {
586777c4855SStefano Zampini yt = t;
587777c4855SStefano Zampini xt = l;
588777c4855SStefano Zampini limit = tsize;
589777c4855SStefano Zampini boundary = user->top;
590777c4855SStefano Zampini } else if (j == 2) {
591777c4855SStefano Zampini yt = b;
592777c4855SStefano Zampini xt = l;
593777c4855SStefano Zampini limit = lsize;
594777c4855SStefano Zampini boundary = user->left;
595777c4855SStefano Zampini } else { /* if (j==3) */
596777c4855SStefano Zampini yt = b;
597777c4855SStefano Zampini xt = r;
598777c4855SStefano Zampini limit = rsize;
599777c4855SStefano Zampini boundary = user->right;
600777c4855SStefano Zampini }
601777c4855SStefano Zampini
602777c4855SStefano Zampini for (i = 0; i < limit; i++) {
603777c4855SStefano Zampini if (j == 0) { /* bottom */
604777c4855SStefano Zampini boundary[i] = -0.5 * PetscSinReal(pi2 * xt);
605777c4855SStefano Zampini } else if (j == 1) { /* top */
606777c4855SStefano Zampini boundary[i] = 0.5 * PetscSinReal(pi2 * xt);
607777c4855SStefano Zampini } else if (j == 2) { /* left */
608777c4855SStefano Zampini boundary[i] = -0.5 * PetscSinReal(pi2 * yt);
609777c4855SStefano Zampini } else { /* right */
610777c4855SStefano Zampini boundary[i] = 0.5 * PetscSinReal(pi2 * yt);
611777c4855SStefano Zampini }
612777c4855SStefano Zampini if (j == 0 || j == 1) xt = xt + hx;
613777c4855SStefano Zampini else yt = yt + hy;
614777c4855SStefano Zampini }
615777c4855SStefano Zampini }
616777c4855SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
617777c4855SStefano Zampini }
618777c4855SStefano Zampini
619777c4855SStefano Zampini /*TEST
620777c4855SStefano Zampini
621777c4855SStefano Zampini build:
622777c4855SStefano Zampini requires: !complex
623777c4855SStefano Zampini
624777c4855SStefano Zampini test:
625777c4855SStefano Zampini requires: !single
626777c4855SStefano Zampini filter: sed -e "s/CONVERGED_FNORM_ABS/CONVERGED_FNORM_RELATIVE/g"
627777c4855SStefano Zampini suffix: qn_nasm
628777c4855SStefano Zampini args: -snes_type qn -snes_npc_side {{left right}separate output} -npc_snes_type nasm -snes_converged_reason -da_local_subdomains 4 -use_objective {{0 1}separate output}
629777c4855SStefano Zampini
630777c4855SStefano Zampini TEST*/
631