xref: /petsc/src/snes/tutorials/ex5.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 2d.\n\
2c4762a1bSJed Brown We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
3c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\
4c4762a1bSJed Brown The command line options include:\n\
5c4762a1bSJed Brown   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
6c4762a1bSJed Brown      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n\
7c4762a1bSJed Brown   -m_par/n_par <parameter>, where <parameter> indicates an integer\n \
8c4762a1bSJed Brown       that MMS3 will be evaluated with 2^m_par, 2^n_par";
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /* ------------------------------------------------------------------------
11c4762a1bSJed Brown 
12c4762a1bSJed Brown     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
13c4762a1bSJed Brown     the partial differential equation
14c4762a1bSJed Brown 
15c4762a1bSJed Brown             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     with boundary conditions
18c4762a1bSJed Brown 
19c4762a1bSJed Brown              u = 0  for  x = 0, x = 1, y = 0, y = 1.
20c4762a1bSJed Brown 
21c4762a1bSJed Brown     A finite difference approximation with the usual 5-point stencil
22c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a nonlinear
23c4762a1bSJed Brown     system of equations.
24c4762a1bSJed Brown 
25c4762a1bSJed Brown       This example shows how geometric multigrid can be run transparently with a nonlinear solver so long
26c4762a1bSJed Brown       as SNESSetDM() is provided. Example usage
27c4762a1bSJed Brown 
28c4762a1bSJed Brown       ./ex5 -pc_type mg -ksp_monitor  -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17
29c4762a1bSJed Brown              -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full
30c4762a1bSJed Brown 
31c4762a1bSJed Brown       or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
32c4762a1bSJed Brown          multigrid levels, it will be determined automatically based on the number of refinements done)
33c4762a1bSJed Brown 
34c4762a1bSJed Brown       ./ex5 -pc_type mg -ksp_monitor  -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3
35c4762a1bSJed Brown              -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   ------------------------------------------------------------------------- */
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /*
40c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
41c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
42c4762a1bSJed Brown */
43c4762a1bSJed Brown #include <petscdm.h>
44c4762a1bSJed Brown #include <petscdmda.h>
45c4762a1bSJed Brown #include <petscsnes.h>
46c4762a1bSJed Brown #include <petscmatlab.h>
47c4762a1bSJed Brown #include <petsc/private/snesimpl.h> /* For SNES_Solve event */
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*
50c4762a1bSJed Brown    User-defined application context - contains data needed by the
51c4762a1bSJed Brown    application-provided call-back routines, FormJacobianLocal() and
52c4762a1bSJed Brown    FormFunctionLocal().
53c4762a1bSJed Brown */
54c4762a1bSJed Brown typedef struct AppCtx AppCtx;
55c4762a1bSJed Brown struct AppCtx {
56c4762a1bSJed Brown   PetscReal param; /* test problem parameter */
57c4762a1bSJed Brown   PetscInt  m, n;  /* MMS3 parameters */
58c4762a1bSJed Brown   PetscErrorCode (*mms_solution)(AppCtx *, const DMDACoor2d *, PetscScalar *);
59c4762a1bSJed Brown   PetscErrorCode (*mms_forcing)(AppCtx *, const DMDACoor2d *, PetscScalar *);
60c4762a1bSJed Brown };
61c4762a1bSJed Brown 
62227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */
63c4762a1bSJed Brown /*
64227f9e03SJunchao Zhang    FormInitialGuess - Forms initial approximation.
65227f9e03SJunchao Zhang 
66227f9e03SJunchao Zhang    Input Parameters:
67227f9e03SJunchao Zhang    da - The DM
68227f9e03SJunchao Zhang    user - user-defined application context
69227f9e03SJunchao Zhang 
70227f9e03SJunchao Zhang    Output Parameter:
71227f9e03SJunchao Zhang    X - vector
72c4762a1bSJed Brown  */
73*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormInitialGuess(DM da, AppCtx *user, Vec X)
74*d71ae5a4SJacob Faibussowitsch {
75227f9e03SJunchao Zhang   PetscInt      i, j, Mx, My, xs, ys, xm, ym;
76227f9e03SJunchao Zhang   PetscReal     lambda, temp1, temp, hx, hy;
77227f9e03SJunchao Zhang   PetscScalar **x;
78227f9e03SJunchao Zhang 
79227f9e03SJunchao Zhang   PetscFunctionBeginUser;
80227f9e03SJunchao Zhang   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));
81227f9e03SJunchao Zhang 
82227f9e03SJunchao Zhang   lambda = user->param;
83227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(Mx - 1);
84227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(My - 1);
85227f9e03SJunchao Zhang   temp1  = lambda / (lambda + 1.0);
86227f9e03SJunchao Zhang 
87227f9e03SJunchao Zhang   /*
88227f9e03SJunchao Zhang      Get a pointer to vector data.
89227f9e03SJunchao Zhang        - For default PETSc vectors, VecGetArray() returns a pointer to
90227f9e03SJunchao Zhang          the data array.  Otherwise, the routine is implementation dependent.
91227f9e03SJunchao Zhang        - You MUST call VecRestoreArray() when you no longer need access to
92227f9e03SJunchao Zhang          the array.
93227f9e03SJunchao Zhang   */
94227f9e03SJunchao Zhang   PetscCall(DMDAVecGetArray(da, X, &x));
95227f9e03SJunchao Zhang 
96227f9e03SJunchao Zhang   /*
97227f9e03SJunchao Zhang      Get local grid boundaries (for 2-dimensional DMDA):
98227f9e03SJunchao Zhang        xs, ys   - starting grid indices (no ghost points)
99227f9e03SJunchao Zhang        xm, ym   - widths of local grid (no ghost points)
100227f9e03SJunchao Zhang 
101227f9e03SJunchao Zhang   */
102227f9e03SJunchao Zhang   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
103227f9e03SJunchao Zhang 
104227f9e03SJunchao Zhang   /*
105227f9e03SJunchao Zhang      Compute initial guess over the locally owned part of the grid
106227f9e03SJunchao Zhang   */
107227f9e03SJunchao Zhang   for (j = ys; j < ys + ym; j++) {
108227f9e03SJunchao Zhang     temp = (PetscReal)(PetscMin(j, My - j - 1)) * hy;
109227f9e03SJunchao Zhang     for (i = xs; i < xs + xm; i++) {
110227f9e03SJunchao Zhang       if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
111227f9e03SJunchao Zhang         /* boundary conditions are all zero Dirichlet */
112227f9e03SJunchao Zhang         x[j][i] = 0.0;
113227f9e03SJunchao Zhang       } else {
114227f9e03SJunchao Zhang         x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp));
115227f9e03SJunchao Zhang       }
116227f9e03SJunchao Zhang     }
117227f9e03SJunchao Zhang   }
118227f9e03SJunchao Zhang 
119227f9e03SJunchao Zhang   /*
120227f9e03SJunchao Zhang      Restore vector
121227f9e03SJunchao Zhang   */
122227f9e03SJunchao Zhang   PetscCall(DMDAVecRestoreArray(da, X, &x));
123227f9e03SJunchao Zhang   PetscFunctionReturn(0);
124227f9e03SJunchao Zhang }
125227f9e03SJunchao Zhang 
126227f9e03SJunchao Zhang /*
127227f9e03SJunchao Zhang   FormExactSolution - Forms MMS solution
128227f9e03SJunchao Zhang 
129227f9e03SJunchao Zhang   Input Parameters:
130227f9e03SJunchao Zhang   da - The DM
131227f9e03SJunchao Zhang   user - user-defined application context
132227f9e03SJunchao Zhang 
133227f9e03SJunchao Zhang   Output Parameter:
134227f9e03SJunchao Zhang   X - vector
135227f9e03SJunchao Zhang  */
136*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U)
137*d71ae5a4SJacob Faibussowitsch {
138227f9e03SJunchao Zhang   DM            coordDA;
139227f9e03SJunchao Zhang   Vec           coordinates;
140227f9e03SJunchao Zhang   DMDACoor2d  **coords;
141227f9e03SJunchao Zhang   PetscScalar **u;
142227f9e03SJunchao Zhang   PetscInt      xs, ys, xm, ym, i, j;
143227f9e03SJunchao Zhang 
144227f9e03SJunchao Zhang   PetscFunctionBeginUser;
145227f9e03SJunchao Zhang   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
146227f9e03SJunchao Zhang   PetscCall(DMGetCoordinateDM(da, &coordDA));
147227f9e03SJunchao Zhang   PetscCall(DMGetCoordinates(da, &coordinates));
148227f9e03SJunchao Zhang   PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords));
149227f9e03SJunchao Zhang   PetscCall(DMDAVecGetArray(da, U, &u));
150227f9e03SJunchao Zhang   for (j = ys; j < ys + ym; ++j) {
151ad540459SPierre Jolivet     for (i = xs; i < xs + xm; ++i) user->mms_solution(user, &coords[j][i], &u[j][i]);
152227f9e03SJunchao Zhang   }
153227f9e03SJunchao Zhang   PetscCall(DMDAVecRestoreArray(da, U, &u));
154227f9e03SJunchao Zhang   PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords));
155227f9e03SJunchao Zhang   PetscFunctionReturn(0);
156227f9e03SJunchao Zhang }
157227f9e03SJunchao Zhang 
158*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ZeroBCSolution(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
159*d71ae5a4SJacob Faibussowitsch {
160227f9e03SJunchao Zhang   u[0] = 0.;
161227f9e03SJunchao Zhang   return 0;
162227f9e03SJunchao Zhang }
163227f9e03SJunchao Zhang 
164227f9e03SJunchao Zhang /* The functions below evaluate the MMS solution u(x,y) and associated forcing
165227f9e03SJunchao Zhang 
166227f9e03SJunchao Zhang      f(x,y) = -u_xx - u_yy - lambda exp(u)
167227f9e03SJunchao Zhang 
168227f9e03SJunchao Zhang   such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term.
169227f9e03SJunchao Zhang  */
170*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution1(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
171*d71ae5a4SJacob Faibussowitsch {
172227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
173227f9e03SJunchao Zhang   u[0] = x * (1 - x) * y * (1 - y);
174227f9e03SJunchao Zhang   PetscLogFlops(5);
175227f9e03SJunchao Zhang   return 0;
176227f9e03SJunchao Zhang }
177*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing1(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
178*d71ae5a4SJacob Faibussowitsch {
179227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
180227f9e03SJunchao Zhang   f[0] = 2 * x * (1 - x) + 2 * y * (1 - y) - user->param * PetscExpReal(x * (1 - x) * y * (1 - y));
181227f9e03SJunchao Zhang   return 0;
182227f9e03SJunchao Zhang }
183227f9e03SJunchao Zhang 
184*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution2(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
185*d71ae5a4SJacob Faibussowitsch {
186227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
187227f9e03SJunchao Zhang   u[0] = PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y);
188227f9e03SJunchao Zhang   PetscLogFlops(5);
189227f9e03SJunchao Zhang   return 0;
190227f9e03SJunchao Zhang }
191*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing2(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
192*d71ae5a4SJacob Faibussowitsch {
193227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
194227f9e03SJunchao Zhang   f[0] = 2 * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y) - user->param * PetscExpReal(PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y));
195227f9e03SJunchao Zhang   return 0;
196227f9e03SJunchao Zhang }
197227f9e03SJunchao Zhang 
198*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution3(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
199*d71ae5a4SJacob Faibussowitsch {
200227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
201227f9e03SJunchao Zhang   u[0] = PetscSinReal(user->m * PETSC_PI * x * (1 - y)) * PetscSinReal(user->n * PETSC_PI * y * (1 - x));
202227f9e03SJunchao Zhang   PetscLogFlops(5);
203227f9e03SJunchao Zhang   return 0;
204227f9e03SJunchao Zhang }
205*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing3(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
206*d71ae5a4SJacob Faibussowitsch {
207227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
208227f9e03SJunchao Zhang   PetscReal m = user->m, n = user->n, lambda = user->param;
2099371c9d4SSatish Balay   f[0] = (-(PetscExpReal(PetscSinReal(m * PETSC_PI * x * (1 - y)) * PetscSinReal(n * PETSC_PI * (1 - x) * y)) * lambda) + PetscSqr(PETSC_PI) * (-2 * m * n * ((-1 + x) * x + (-1 + y) * y) * PetscCosReal(m * PETSC_PI * x * (-1 + y)) * PetscCosReal(n * PETSC_PI * (-1 + x) * y) + (PetscSqr(m) * (PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n) * (PetscSqr(-1 + x) + PetscSqr(y))) * PetscSinReal(m * PETSC_PI * x * (-1 + y)) * PetscSinReal(n * PETSC_PI * (-1 + x) * y)));
210227f9e03SJunchao Zhang   return 0;
211227f9e03SJunchao Zhang }
212227f9e03SJunchao Zhang 
213*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution4(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
214*d71ae5a4SJacob Faibussowitsch {
215227f9e03SJunchao Zhang   const PetscReal Lx = 1., Ly = 1.;
216227f9e03SJunchao Zhang   PetscReal       x = PetscRealPart(c->x), y = PetscRealPart(c->y);
217227f9e03SJunchao Zhang   u[0] = (PetscPowReal(x, 4) - PetscSqr(Lx) * PetscSqr(x)) * (PetscPowReal(y, 4) - PetscSqr(Ly) * PetscSqr(y));
218227f9e03SJunchao Zhang   PetscLogFlops(9);
219227f9e03SJunchao Zhang   return 0;
220227f9e03SJunchao Zhang }
221*d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing4(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
222*d71ae5a4SJacob Faibussowitsch {
223227f9e03SJunchao Zhang   const PetscReal Lx = 1., Ly = 1.;
224227f9e03SJunchao Zhang   PetscReal       x = PetscRealPart(c->x), y = PetscRealPart(c->y);
2259371c9d4SSatish Balay   f[0] = (2 * PetscSqr(x) * (PetscSqr(x) - PetscSqr(Lx)) * (PetscSqr(Ly) - 6 * PetscSqr(y)) + 2 * PetscSqr(y) * (PetscSqr(Lx) - 6 * PetscSqr(x)) * (PetscSqr(y) - PetscSqr(Ly)) - user->param * PetscExpReal((PetscPowReal(x, 4) - PetscSqr(Lx) * PetscSqr(x)) * (PetscPowReal(y, 4) - PetscSqr(Ly) * PetscSqr(y))));
226227f9e03SJunchao Zhang   return 0;
227227f9e03SJunchao Zhang }
228227f9e03SJunchao Zhang 
229227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */
230227f9e03SJunchao Zhang /*
231227f9e03SJunchao Zhang    FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch
232227f9e03SJunchao Zhang 
233227f9e03SJunchao Zhang  */
234*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user)
235*d71ae5a4SJacob Faibussowitsch {
236227f9e03SJunchao Zhang   PetscInt    i, j;
237227f9e03SJunchao Zhang   PetscReal   lambda, hx, hy, hxdhy, hydhx;
238227f9e03SJunchao Zhang   PetscScalar u, ue, uw, un, us, uxx, uyy, mms_solution, mms_forcing;
239227f9e03SJunchao Zhang   DMDACoor2d  c;
240227f9e03SJunchao Zhang 
241227f9e03SJunchao Zhang   PetscFunctionBeginUser;
242227f9e03SJunchao Zhang   lambda = user->param;
243227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(info->mx - 1);
244227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(info->my - 1);
245227f9e03SJunchao Zhang   hxdhy  = hx / hy;
246227f9e03SJunchao Zhang   hydhx  = hy / hx;
247227f9e03SJunchao Zhang   /*
248227f9e03SJunchao Zhang      Compute function over the locally owned part of the grid
249227f9e03SJunchao Zhang   */
250227f9e03SJunchao Zhang   for (j = info->ys; j < info->ys + info->ym; j++) {
251227f9e03SJunchao Zhang     for (i = info->xs; i < info->xs + info->xm; i++) {
252227f9e03SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
2539371c9d4SSatish Balay         c.x = i * hx;
2549371c9d4SSatish Balay         c.y = j * hy;
255227f9e03SJunchao Zhang         PetscCall(user->mms_solution(user, &c, &mms_solution));
256227f9e03SJunchao Zhang         f[j][i] = 2.0 * (hydhx + hxdhy) * (x[j][i] - mms_solution);
257227f9e03SJunchao Zhang       } else {
258227f9e03SJunchao Zhang         u  = x[j][i];
259227f9e03SJunchao Zhang         uw = x[j][i - 1];
260227f9e03SJunchao Zhang         ue = x[j][i + 1];
261227f9e03SJunchao Zhang         un = x[j - 1][i];
262227f9e03SJunchao Zhang         us = x[j + 1][i];
263227f9e03SJunchao Zhang 
264227f9e03SJunchao Zhang         /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */
2659371c9d4SSatish Balay         if (i - 1 == 0) {
2669371c9d4SSatish Balay           c.x = (i - 1) * hx;
2679371c9d4SSatish Balay           c.y = j * hy;
2689371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &uw));
2699371c9d4SSatish Balay         }
2709371c9d4SSatish Balay         if (i + 1 == info->mx - 1) {
2719371c9d4SSatish Balay           c.x = (i + 1) * hx;
2729371c9d4SSatish Balay           c.y = j * hy;
2739371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &ue));
2749371c9d4SSatish Balay         }
2759371c9d4SSatish Balay         if (j - 1 == 0) {
2769371c9d4SSatish Balay           c.x = i * hx;
2779371c9d4SSatish Balay           c.y = (j - 1) * hy;
2789371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &un));
2799371c9d4SSatish Balay         }
2809371c9d4SSatish Balay         if (j + 1 == info->my - 1) {
2819371c9d4SSatish Balay           c.x = i * hx;
2829371c9d4SSatish Balay           c.y = (j + 1) * hy;
2839371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &us));
2849371c9d4SSatish Balay         }
285227f9e03SJunchao Zhang 
286227f9e03SJunchao Zhang         uxx         = (2.0 * u - uw - ue) * hydhx;
287227f9e03SJunchao Zhang         uyy         = (2.0 * u - un - us) * hxdhy;
288227f9e03SJunchao Zhang         mms_forcing = 0;
2899371c9d4SSatish Balay         c.x         = i * hx;
2909371c9d4SSatish Balay         c.y         = j * hy;
291227f9e03SJunchao Zhang         if (user->mms_forcing) PetscCall(user->mms_forcing(user, &c, &mms_forcing));
292227f9e03SJunchao Zhang         f[j][i] = uxx + uyy - hx * hy * (lambda * PetscExpScalar(u) + mms_forcing);
293227f9e03SJunchao Zhang       }
294227f9e03SJunchao Zhang     }
295227f9e03SJunchao Zhang   }
296227f9e03SJunchao Zhang   PetscCall(PetscLogFlops(11.0 * info->ym * info->xm));
297227f9e03SJunchao Zhang   PetscFunctionReturn(0);
298227f9e03SJunchao Zhang }
299227f9e03SJunchao Zhang 
300227f9e03SJunchao Zhang /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */
301*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info, PetscScalar **x, PetscReal *obj, AppCtx *user)
302*d71ae5a4SJacob Faibussowitsch {
303227f9e03SJunchao Zhang   PetscInt    i, j;
304227f9e03SJunchao Zhang   PetscReal   lambda, hx, hy, hxdhy, hydhx, sc, lobj = 0;
305227f9e03SJunchao Zhang   PetscScalar u, ue, uw, un, us, uxux, uyuy;
306227f9e03SJunchao Zhang   MPI_Comm    comm;
307227f9e03SJunchao Zhang 
308227f9e03SJunchao Zhang   PetscFunctionBeginUser;
309227f9e03SJunchao Zhang   *obj = 0;
310227f9e03SJunchao Zhang   PetscCall(PetscObjectGetComm((PetscObject)info->da, &comm));
311227f9e03SJunchao Zhang   lambda = user->param;
312227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(info->mx - 1);
313227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(info->my - 1);
314227f9e03SJunchao Zhang   sc     = hx * hy * lambda;
315227f9e03SJunchao Zhang   hxdhy  = hx / hy;
316227f9e03SJunchao Zhang   hydhx  = hy / hx;
317227f9e03SJunchao Zhang   /*
318227f9e03SJunchao Zhang      Compute function over the locally owned part of the grid
319227f9e03SJunchao Zhang   */
320227f9e03SJunchao Zhang   for (j = info->ys; j < info->ys + info->ym; j++) {
321227f9e03SJunchao Zhang     for (i = info->xs; i < info->xs + info->xm; i++) {
322227f9e03SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
323227f9e03SJunchao Zhang         lobj += PetscRealPart((hydhx + hxdhy) * x[j][i] * x[j][i]);
324227f9e03SJunchao Zhang       } else {
325227f9e03SJunchao Zhang         u  = x[j][i];
326227f9e03SJunchao Zhang         uw = x[j][i - 1];
327227f9e03SJunchao Zhang         ue = x[j][i + 1];
328227f9e03SJunchao Zhang         un = x[j - 1][i];
329227f9e03SJunchao Zhang         us = x[j + 1][i];
330227f9e03SJunchao Zhang 
331227f9e03SJunchao Zhang         if (i - 1 == 0) uw = 0.;
332227f9e03SJunchao Zhang         if (i + 1 == info->mx - 1) ue = 0.;
333227f9e03SJunchao Zhang         if (j - 1 == 0) un = 0.;
334227f9e03SJunchao Zhang         if (j + 1 == info->my - 1) us = 0.;
335227f9e03SJunchao Zhang 
336227f9e03SJunchao Zhang         /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */
337227f9e03SJunchao Zhang 
338227f9e03SJunchao Zhang         uxux = u * (2. * u - ue - uw) * hydhx;
339227f9e03SJunchao Zhang         uyuy = u * (2. * u - un - us) * hxdhy;
340227f9e03SJunchao Zhang 
341227f9e03SJunchao Zhang         lobj += PetscRealPart(0.5 * (uxux + uyuy) - sc * PetscExpScalar(u));
342227f9e03SJunchao Zhang       }
343227f9e03SJunchao Zhang     }
344227f9e03SJunchao Zhang   }
345227f9e03SJunchao Zhang   PetscCall(PetscLogFlops(12.0 * info->ym * info->xm));
346227f9e03SJunchao Zhang   PetscCallMPI(MPI_Allreduce(&lobj, obj, 1, MPIU_REAL, MPIU_SUM, comm));
347227f9e03SJunchao Zhang   PetscFunctionReturn(0);
348227f9e03SJunchao Zhang }
349227f9e03SJunchao Zhang 
350227f9e03SJunchao Zhang /*
351227f9e03SJunchao Zhang    FormJacobianLocal - Evaluates Jacobian matrix on local process patch
352227f9e03SJunchao Zhang */
353*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat jac, Mat jacpre, AppCtx *user)
354*d71ae5a4SJacob Faibussowitsch {
355227f9e03SJunchao Zhang   PetscInt     i, j, k;
356227f9e03SJunchao Zhang   MatStencil   col[5], row;
357227f9e03SJunchao Zhang   PetscScalar  lambda, v[5], hx, hy, hxdhy, hydhx, sc;
358227f9e03SJunchao Zhang   DM           coordDA;
359227f9e03SJunchao Zhang   Vec          coordinates;
360227f9e03SJunchao Zhang   DMDACoor2d **coords;
361227f9e03SJunchao Zhang 
362227f9e03SJunchao Zhang   PetscFunctionBeginUser;
363227f9e03SJunchao Zhang   lambda = user->param;
364227f9e03SJunchao Zhang   /* Extract coordinates */
365227f9e03SJunchao Zhang   PetscCall(DMGetCoordinateDM(info->da, &coordDA));
366227f9e03SJunchao Zhang   PetscCall(DMGetCoordinates(info->da, &coordinates));
367227f9e03SJunchao Zhang   PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords));
368227f9e03SJunchao Zhang   hx = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs + 1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0;
369227f9e03SJunchao Zhang   hy = info->ym > 1 ? PetscRealPart(coords[info->ys + 1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0;
370227f9e03SJunchao Zhang   PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords));
371227f9e03SJunchao Zhang   hxdhy = hx / hy;
372227f9e03SJunchao Zhang   hydhx = hy / hx;
373227f9e03SJunchao Zhang   sc    = hx * hy * lambda;
374227f9e03SJunchao Zhang 
375227f9e03SJunchao Zhang   /*
376227f9e03SJunchao Zhang      Compute entries for the locally owned part of the Jacobian.
377227f9e03SJunchao Zhang       - Currently, all PETSc parallel matrix formats are partitioned by
378227f9e03SJunchao Zhang         contiguous chunks of rows across the processors.
379227f9e03SJunchao Zhang       - Each processor needs to insert only elements that it owns
380227f9e03SJunchao Zhang         locally (but any non-local elements will be sent to the
381227f9e03SJunchao Zhang         appropriate processor during matrix assembly).
382227f9e03SJunchao Zhang       - Here, we set all entries for a particular row at once.
383227f9e03SJunchao Zhang       - We can set matrix entries either using either
384227f9e03SJunchao Zhang         MatSetValuesLocal() or MatSetValues(), as discussed above.
385227f9e03SJunchao Zhang   */
386227f9e03SJunchao Zhang   for (j = info->ys; j < info->ys + info->ym; j++) {
387227f9e03SJunchao Zhang     for (i = info->xs; i < info->xs + info->xm; i++) {
3889371c9d4SSatish Balay       row.j = j;
3899371c9d4SSatish Balay       row.i = i;
390227f9e03SJunchao Zhang       /* boundary points */
391227f9e03SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
392227f9e03SJunchao Zhang         v[0] = 2.0 * (hydhx + hxdhy);
393227f9e03SJunchao Zhang         PetscCall(MatSetValuesStencil(jacpre, 1, &row, 1, &row, v, INSERT_VALUES));
394227f9e03SJunchao Zhang       } else {
395227f9e03SJunchao Zhang         k = 0;
396227f9e03SJunchao Zhang         /* interior grid points */
397227f9e03SJunchao Zhang         if (j - 1 != 0) {
398227f9e03SJunchao Zhang           v[k]     = -hxdhy;
3999371c9d4SSatish Balay           col[k].j = j - 1;
4009371c9d4SSatish Balay           col[k].i = i;
401227f9e03SJunchao Zhang           k++;
402227f9e03SJunchao Zhang         }
403227f9e03SJunchao Zhang         if (i - 1 != 0) {
404227f9e03SJunchao Zhang           v[k]     = -hydhx;
4059371c9d4SSatish Balay           col[k].j = j;
4069371c9d4SSatish Balay           col[k].i = i - 1;
407227f9e03SJunchao Zhang           k++;
408227f9e03SJunchao Zhang         }
409227f9e03SJunchao Zhang 
4109371c9d4SSatish Balay         v[k]     = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(x[j][i]);
4119371c9d4SSatish Balay         col[k].j = row.j;
4129371c9d4SSatish Balay         col[k].i = row.i;
4139371c9d4SSatish Balay         k++;
414227f9e03SJunchao Zhang 
415227f9e03SJunchao Zhang         if (i + 1 != info->mx - 1) {
416227f9e03SJunchao Zhang           v[k]     = -hydhx;
4179371c9d4SSatish Balay           col[k].j = j;
4189371c9d4SSatish Balay           col[k].i = i + 1;
419227f9e03SJunchao Zhang           k++;
420227f9e03SJunchao Zhang         }
421227f9e03SJunchao Zhang         if (j + 1 != info->mx - 1) {
422227f9e03SJunchao Zhang           v[k]     = -hxdhy;
4239371c9d4SSatish Balay           col[k].j = j + 1;
4249371c9d4SSatish Balay           col[k].i = i;
425227f9e03SJunchao Zhang           k++;
426227f9e03SJunchao Zhang         }
427227f9e03SJunchao Zhang         PetscCall(MatSetValuesStencil(jacpre, 1, &row, k, col, v, INSERT_VALUES));
428227f9e03SJunchao Zhang       }
429227f9e03SJunchao Zhang     }
430227f9e03SJunchao Zhang   }
431227f9e03SJunchao Zhang 
432227f9e03SJunchao Zhang   /*
433227f9e03SJunchao Zhang      Assemble matrix, using the 2-step process:
434227f9e03SJunchao Zhang        MatAssemblyBegin(), MatAssemblyEnd().
435227f9e03SJunchao Zhang   */
436227f9e03SJunchao Zhang   PetscCall(MatAssemblyBegin(jacpre, MAT_FINAL_ASSEMBLY));
437227f9e03SJunchao Zhang   PetscCall(MatAssemblyEnd(jacpre, MAT_FINAL_ASSEMBLY));
438227f9e03SJunchao Zhang   /*
439227f9e03SJunchao Zhang      Tell the matrix we will never add a new nonzero location to the
440227f9e03SJunchao Zhang      matrix. If we do, it will generate an error.
441227f9e03SJunchao Zhang   */
442227f9e03SJunchao Zhang   PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
443227f9e03SJunchao Zhang   PetscFunctionReturn(0);
444227f9e03SJunchao Zhang }
445227f9e03SJunchao Zhang 
446*d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionMatlab(SNES snes, Vec X, Vec F, void *ptr)
447*d71ae5a4SJacob Faibussowitsch {
448d1e78c4fSBarry Smith #if PetscDefined(HAVE_MATLAB)
449227f9e03SJunchao Zhang   AppCtx   *user = (AppCtx *)ptr;
450227f9e03SJunchao Zhang   PetscInt  Mx, My;
451227f9e03SJunchao Zhang   PetscReal lambda, hx, hy;
452227f9e03SJunchao Zhang   Vec       localX, localF;
453227f9e03SJunchao Zhang   MPI_Comm  comm;
454227f9e03SJunchao Zhang   DM        da;
455227f9e03SJunchao Zhang 
456227f9e03SJunchao Zhang   PetscFunctionBeginUser;
457227f9e03SJunchao Zhang   PetscCall(SNESGetDM(snes, &da));
458227f9e03SJunchao Zhang   PetscCall(DMGetLocalVector(da, &localX));
459227f9e03SJunchao Zhang   PetscCall(DMGetLocalVector(da, &localF));
460227f9e03SJunchao Zhang   PetscCall(PetscObjectSetName((PetscObject)localX, "localX"));
461227f9e03SJunchao Zhang   PetscCall(PetscObjectSetName((PetscObject)localF, "localF"));
462227f9e03SJunchao Zhang   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));
463227f9e03SJunchao Zhang 
464227f9e03SJunchao Zhang   lambda = user->param;
465227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(Mx - 1);
466227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(My - 1);
467227f9e03SJunchao Zhang 
468227f9e03SJunchao Zhang   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
469227f9e03SJunchao Zhang   /*
470227f9e03SJunchao Zhang      Scatter ghost points to local vector,using the 2-step process
471227f9e03SJunchao Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
472227f9e03SJunchao Zhang      By placing code between these two statements, computations can be
473227f9e03SJunchao Zhang      done while messages are in transition.
474227f9e03SJunchao Zhang   */
475227f9e03SJunchao Zhang   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
476227f9e03SJunchao Zhang   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
477227f9e03SJunchao Zhang   PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localX));
478227f9e03SJunchao Zhang   PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm), "localF=ex5m(localX,%18.16e,%18.16e,%18.16e)", (double)hx, (double)hy, (double)lambda));
479227f9e03SJunchao Zhang   PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localF));
480227f9e03SJunchao Zhang 
481227f9e03SJunchao Zhang   /*
482227f9e03SJunchao Zhang      Insert values into global vector
483227f9e03SJunchao Zhang   */
484227f9e03SJunchao Zhang   PetscCall(DMLocalToGlobalBegin(da, localF, INSERT_VALUES, F));
485227f9e03SJunchao Zhang   PetscCall(DMLocalToGlobalEnd(da, localF, INSERT_VALUES, F));
486227f9e03SJunchao Zhang   PetscCall(DMRestoreLocalVector(da, &localX));
487227f9e03SJunchao Zhang   PetscCall(DMRestoreLocalVector(da, &localF));
488227f9e03SJunchao Zhang   PetscFunctionReturn(0);
489227f9e03SJunchao Zhang #else
490227f9e03SJunchao Zhang   return 0; /* Never called */
491227f9e03SJunchao Zhang #endif
492227f9e03SJunchao Zhang }
493227f9e03SJunchao Zhang 
494227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */
495227f9e03SJunchao Zhang /*
496227f9e03SJunchao Zhang       Applies some sweeps on nonlinear Gauss-Seidel on each process
497227f9e03SJunchao Zhang 
498227f9e03SJunchao Zhang  */
499*d71ae5a4SJacob Faibussowitsch static PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx)
500*d71ae5a4SJacob Faibussowitsch {
501227f9e03SJunchao Zhang   PetscInt      i, j, k, Mx, My, xs, ys, xm, ym, its, tot_its, sweeps, l;
502227f9e03SJunchao Zhang   PetscReal     lambda, hx, hy, hxdhy, hydhx, sc;
503227f9e03SJunchao Zhang   PetscScalar **x, **b, bij, F, F0 = 0, J, u, un, us, ue, eu, uw, uxx, uyy, y;
504227f9e03SJunchao Zhang   PetscReal     atol, rtol, stol;
505227f9e03SJunchao Zhang   DM            da;
506227f9e03SJunchao Zhang   AppCtx       *user;
507227f9e03SJunchao Zhang   Vec           localX, localB;
508227f9e03SJunchao Zhang 
509227f9e03SJunchao Zhang   PetscFunctionBeginUser;
510227f9e03SJunchao Zhang   tot_its = 0;
511227f9e03SJunchao Zhang   PetscCall(SNESNGSGetSweeps(snes, &sweeps));
512227f9e03SJunchao Zhang   PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its));
513227f9e03SJunchao Zhang   PetscCall(SNESGetDM(snes, &da));
514227f9e03SJunchao Zhang   PetscCall(DMGetApplicationContext(da, &user));
515227f9e03SJunchao Zhang 
516227f9e03SJunchao Zhang   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));
517227f9e03SJunchao Zhang 
518227f9e03SJunchao Zhang   lambda = user->param;
519227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(Mx - 1);
520227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(My - 1);
521227f9e03SJunchao Zhang   sc     = hx * hy * lambda;
522227f9e03SJunchao Zhang   hxdhy  = hx / hy;
523227f9e03SJunchao Zhang   hydhx  = hy / hx;
524227f9e03SJunchao Zhang 
525227f9e03SJunchao Zhang   PetscCall(DMGetLocalVector(da, &localX));
52648a46eb9SPierre Jolivet   if (B) PetscCall(DMGetLocalVector(da, &localB));
527227f9e03SJunchao Zhang   for (l = 0; l < sweeps; l++) {
528227f9e03SJunchao Zhang     PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
529227f9e03SJunchao Zhang     PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
530227f9e03SJunchao Zhang     if (B) {
531227f9e03SJunchao Zhang       PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, localB));
532227f9e03SJunchao Zhang       PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, localB));
533227f9e03SJunchao Zhang     }
534227f9e03SJunchao Zhang     /*
535227f9e03SJunchao Zhang      Get a pointer to vector data.
536227f9e03SJunchao Zhang      - For default PETSc vectors, VecGetArray() returns a pointer to
537227f9e03SJunchao Zhang      the data array.  Otherwise, the routine is implementation dependent.
538227f9e03SJunchao Zhang      - You MUST call VecRestoreArray() when you no longer need access to
539227f9e03SJunchao Zhang      the array.
540227f9e03SJunchao Zhang      */
541227f9e03SJunchao Zhang     PetscCall(DMDAVecGetArray(da, localX, &x));
542227f9e03SJunchao Zhang     if (B) PetscCall(DMDAVecGetArray(da, localB, &b));
543227f9e03SJunchao Zhang     /*
544227f9e03SJunchao Zhang      Get local grid boundaries (for 2-dimensional DMDA):
545227f9e03SJunchao Zhang      xs, ys   - starting grid indices (no ghost points)
546227f9e03SJunchao Zhang      xm, ym   - widths of local grid (no ghost points)
547227f9e03SJunchao Zhang      */
548227f9e03SJunchao Zhang     PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
549227f9e03SJunchao Zhang 
550227f9e03SJunchao Zhang     for (j = ys; j < ys + ym; j++) {
551227f9e03SJunchao Zhang       for (i = xs; i < xs + xm; i++) {
552227f9e03SJunchao Zhang         if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
553227f9e03SJunchao Zhang           /* boundary conditions are all zero Dirichlet */
554227f9e03SJunchao Zhang           x[j][i] = 0.0;
555227f9e03SJunchao Zhang         } else {
556227f9e03SJunchao Zhang           if (B) bij = b[j][i];
557227f9e03SJunchao Zhang           else bij = 0.;
558227f9e03SJunchao Zhang 
559227f9e03SJunchao Zhang           u  = x[j][i];
560227f9e03SJunchao Zhang           un = x[j - 1][i];
561227f9e03SJunchao Zhang           us = x[j + 1][i];
562227f9e03SJunchao Zhang           ue = x[j][i - 1];
563227f9e03SJunchao Zhang           uw = x[j][i + 1];
564227f9e03SJunchao Zhang 
565227f9e03SJunchao Zhang           for (k = 0; k < its; k++) {
566227f9e03SJunchao Zhang             eu  = PetscExpScalar(u);
567227f9e03SJunchao Zhang             uxx = (2.0 * u - ue - uw) * hydhx;
568227f9e03SJunchao Zhang             uyy = (2.0 * u - un - us) * hxdhy;
569227f9e03SJunchao Zhang             F   = uxx + uyy - sc * eu - bij;
570227f9e03SJunchao Zhang             if (k == 0) F0 = F;
571227f9e03SJunchao Zhang             J = 2.0 * (hydhx + hxdhy) - sc * eu;
572227f9e03SJunchao Zhang             y = F / J;
573227f9e03SJunchao Zhang             u -= y;
574227f9e03SJunchao Zhang             tot_its++;
575227f9e03SJunchao Zhang 
576ad540459SPierre Jolivet             if (atol > PetscAbsReal(PetscRealPart(F)) || rtol * PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || stol * PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) break;
577227f9e03SJunchao Zhang           }
578227f9e03SJunchao Zhang           x[j][i] = u;
579227f9e03SJunchao Zhang         }
580227f9e03SJunchao Zhang       }
581227f9e03SJunchao Zhang     }
582227f9e03SJunchao Zhang     /*
583227f9e03SJunchao Zhang      Restore vector
584227f9e03SJunchao Zhang      */
585227f9e03SJunchao Zhang     PetscCall(DMDAVecRestoreArray(da, localX, &x));
586227f9e03SJunchao Zhang     PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X));
587227f9e03SJunchao Zhang     PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X));
588227f9e03SJunchao Zhang   }
589227f9e03SJunchao Zhang   PetscCall(PetscLogFlops(tot_its * (21.0)));
590227f9e03SJunchao Zhang   PetscCall(DMRestoreLocalVector(da, &localX));
591227f9e03SJunchao Zhang   if (B) {
592227f9e03SJunchao Zhang     PetscCall(DMDAVecRestoreArray(da, localB, &b));
593227f9e03SJunchao Zhang     PetscCall(DMRestoreLocalVector(da, &localB));
594227f9e03SJunchao Zhang   }
595227f9e03SJunchao Zhang   PetscFunctionReturn(0);
596227f9e03SJunchao Zhang }
597c4762a1bSJed Brown 
598*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
599*d71ae5a4SJacob Faibussowitsch {
600c4762a1bSJed Brown   SNES      snes; /* nonlinear solver */
601c4762a1bSJed Brown   Vec       x;    /* solution vector */
602c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
603c4762a1bSJed Brown   PetscInt  its;  /* iterations for convergence */
604c4762a1bSJed Brown   PetscReal bratu_lambda_max = 6.81;
605c4762a1bSJed Brown   PetscReal bratu_lambda_min = 0.;
60616d037c0SJunchao Zhang   PetscInt  MMS              = 1;
60716d037c0SJunchao Zhang   PetscBool flg              = PETSC_FALSE, setMMS;
608c4762a1bSJed Brown   DM        da;
609c4762a1bSJed Brown   Vec       r = NULL;
610c4762a1bSJed Brown   KSP       ksp;
611c4762a1bSJed Brown   PetscInt  lits, slits;
612c4762a1bSJed Brown 
613c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
614c4762a1bSJed Brown      Initialize program
615c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
616c4762a1bSJed Brown 
617327415f7SBarry Smith   PetscFunctionBeginUser;
6189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
619c4762a1bSJed Brown 
620c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
621c4762a1bSJed Brown      Initialize problem parameters
622c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
623c4762a1bSJed Brown   user.param = 6.0;
6249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL));
625e00437b9SBarry Smith   PetscCheck(user.param <= bratu_lambda_max && user.param >= bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda, %g, is out of range, [%g, %g]", (double)user.param, (double)bratu_lambda_min, (double)bratu_lambda_max);
62616d037c0SJunchao Zhang   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mms", &MMS, &setMMS));
627c4762a1bSJed Brown   if (MMS == 3) {
628c4762a1bSJed Brown     PetscInt mPar = 2, nPar = 1;
6299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-m_par", &mPar, NULL));
6309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_par", &nPar, NULL));
631c4762a1bSJed Brown     user.m = PetscPowInt(2, mPar);
632c4762a1bSJed Brown     user.n = PetscPowInt(2, nPar);
633c4762a1bSJed Brown   }
634c4762a1bSJed Brown 
635c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
636c4762a1bSJed Brown      Create nonlinear solver context
637c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
6389566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
6399566063dSJacob Faibussowitsch   PetscCall(SNESSetCountersReset(snes, PETSC_FALSE));
6409566063dSJacob Faibussowitsch   PetscCall(SNESSetNGS(snes, NonlinearGS, NULL));
641c4762a1bSJed Brown 
642c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
643c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
644c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
6459566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
6469566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
6479566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
6489566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
6499566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
6509566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
651c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
652c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
653c4762a1bSJed Brown      vectors that are the same types
654c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
6559566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
658c4762a1bSJed Brown      Set local function evaluation routine
659c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
660c4762a1bSJed Brown   switch (MMS) {
6619371c9d4SSatish Balay   case 0:
6629371c9d4SSatish Balay     user.mms_solution = ZeroBCSolution;
6639371c9d4SSatish Balay     user.mms_forcing  = NULL;
6649371c9d4SSatish Balay     break;
6659371c9d4SSatish Balay   case 1:
6669371c9d4SSatish Balay     user.mms_solution = MMSSolution1;
6679371c9d4SSatish Balay     user.mms_forcing  = MMSForcing1;
6689371c9d4SSatish Balay     break;
6699371c9d4SSatish Balay   case 2:
6709371c9d4SSatish Balay     user.mms_solution = MMSSolution2;
6719371c9d4SSatish Balay     user.mms_forcing  = MMSForcing2;
6729371c9d4SSatish Balay     break;
6739371c9d4SSatish Balay   case 3:
6749371c9d4SSatish Balay     user.mms_solution = MMSSolution3;
6759371c9d4SSatish Balay     user.mms_forcing  = MMSForcing3;
6769371c9d4SSatish Balay     break;
6779371c9d4SSatish Balay   case 4:
6789371c9d4SSatish Balay     user.mms_solution = MMSSolution4;
6799371c9d4SSatish Balay     user.mms_forcing  = MMSForcing4;
6809371c9d4SSatish Balay     break;
681*d71ae5a4SJacob Faibussowitsch   default:
682*d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Unknown MMS type %" PetscInt_FMT, MMS);
683c4762a1bSJed Brown   }
6849566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunction)FormFunctionLocal, &user));
6859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-fd", &flg, NULL));
68648a46eb9SPierre Jolivet   if (!flg) PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobian)FormJacobianLocal, &user));
687c4762a1bSJed Brown 
6889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-obj", &flg, NULL));
68948a46eb9SPierre Jolivet   if (flg) PetscCall(DMDASNESSetObjectiveLocal(da, (DMDASNESObjective)FormObjectiveLocal, &user));
690c4762a1bSJed Brown 
691d1e78c4fSBarry Smith   if (PetscDefined(HAVE_MATLAB)) {
6924e1ad211SJed Brown     PetscBool matlab_function = PETSC_FALSE;
6939566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-matlab_function", &matlab_function, 0));
694c4762a1bSJed Brown     if (matlab_function) {
6959566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &r));
6969566063dSJacob Faibussowitsch       PetscCall(SNESSetFunction(snes, r, FormFunctionMatlab, &user));
697c4762a1bSJed Brown     }
6984e1ad211SJed Brown   }
699c4762a1bSJed Brown 
700c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
701c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
702c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
7039566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
704c4762a1bSJed Brown 
7059566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(da, &user, x));
706c4762a1bSJed Brown 
707c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
708c4762a1bSJed Brown      Solve nonlinear system
709c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
7109566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
7119566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
712c4762a1bSJed Brown 
7139566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes, &slits));
7149566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
7159566063dSJacob Faibussowitsch   PetscCall(KSPGetTotalIterations(ksp, &lits));
71663a3b9bcSJacob Faibussowitsch   PetscCheck(lits == slits, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Number of total linear iterations reported by SNES %" PetscInt_FMT " does not match reported by KSP %" PetscInt_FMT, slits, lits);
717c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
718c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
719c4762a1bSJed Brown 
720c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
721c4762a1bSJed Brown      If using MMS, check the l_2 error
722c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72316d037c0SJunchao Zhang   if (setMMS) {
724c4762a1bSJed Brown     Vec       e;
725c4762a1bSJed Brown     PetscReal errorl2, errorinf;
726c4762a1bSJed Brown     PetscInt  N;
727c4762a1bSJed Brown 
7289566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &e));
7299566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)x, NULL, "-sol_view"));
7309566063dSJacob Faibussowitsch     PetscCall(FormExactSolution(da, &user, e));
7319566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-exact_view"));
7329566063dSJacob Faibussowitsch     PetscCall(VecAXPY(e, -1.0, x));
7339566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-error_view"));
7349566063dSJacob Faibussowitsch     PetscCall(VecNorm(e, NORM_2, &errorl2));
7359566063dSJacob Faibussowitsch     PetscCall(VecNorm(e, NORM_INFINITY, &errorinf));
7369566063dSJacob Faibussowitsch     PetscCall(VecGetSize(e, &N));
73763a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " error L2 %g inf %g\n", N, (double)(errorl2 / PetscSqrtReal((PetscReal)N)), (double)errorinf));
7389566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&e));
7399566063dSJacob Faibussowitsch     PetscCall(PetscLogEventSetDof(SNES_Solve, 0, N));
7409566063dSJacob Faibussowitsch     PetscCall(PetscLogEventSetError(SNES_Solve, 0, errorl2 / PetscSqrtReal(N)));
741c4762a1bSJed Brown   }
742c4762a1bSJed Brown 
743c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
744c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
745c4762a1bSJed Brown      are no longer needed.
746c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
7479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
7489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
7499566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
7509566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
7519566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
752b122ec5aSJacob Faibussowitsch   return 0;
753c4762a1bSJed Brown }
754c4762a1bSJed Brown 
755c4762a1bSJed Brown /*TEST
756c4762a1bSJed Brown 
757c4762a1bSJed Brown    test:
758c4762a1bSJed Brown      suffix: asm_0
759c4762a1bSJed Brown      requires: !single
760c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu
761c4762a1bSJed Brown 
762c4762a1bSJed Brown    test:
763c4762a1bSJed Brown      suffix: msm_0
764c4762a1bSJed Brown      requires: !single
765c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu
766c4762a1bSJed Brown 
767c4762a1bSJed Brown    test:
768c4762a1bSJed Brown      suffix: asm_1
769c4762a1bSJed Brown      requires: !single
770c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
771c4762a1bSJed Brown 
772c4762a1bSJed Brown    test:
773c4762a1bSJed Brown      suffix: msm_1
774c4762a1bSJed Brown      requires: !single
775c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
776c4762a1bSJed Brown 
777c4762a1bSJed Brown    test:
778c4762a1bSJed Brown      suffix: asm_2
779c4762a1bSJed Brown      requires: !single
780c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
781c4762a1bSJed Brown 
782c4762a1bSJed Brown    test:
783c4762a1bSJed Brown      suffix: msm_2
784c4762a1bSJed Brown      requires: !single
785c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
786c4762a1bSJed Brown 
787c4762a1bSJed Brown    test:
788c4762a1bSJed Brown      suffix: asm_3
789c4762a1bSJed Brown      requires: !single
790c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
791c4762a1bSJed Brown 
792c4762a1bSJed Brown    test:
793c4762a1bSJed Brown      suffix: msm_3
794c4762a1bSJed Brown      requires: !single
795c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
796c4762a1bSJed Brown 
797c4762a1bSJed Brown    test:
798c4762a1bSJed Brown      suffix: asm_4
799c4762a1bSJed Brown      requires: !single
800c4762a1bSJed Brown      nsize: 2
801c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
802c4762a1bSJed Brown 
803c4762a1bSJed Brown    test:
804c4762a1bSJed Brown      suffix: msm_4
805c4762a1bSJed Brown      requires: !single
806c4762a1bSJed Brown      nsize: 2
807c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
808c4762a1bSJed Brown 
809c4762a1bSJed Brown    test:
810c4762a1bSJed Brown      suffix: asm_5
811c4762a1bSJed Brown      requires: !single
812c4762a1bSJed Brown      nsize: 2
813c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
814c4762a1bSJed Brown 
815c4762a1bSJed Brown    test:
816c4762a1bSJed Brown      suffix: msm_5
817c4762a1bSJed Brown      requires: !single
818c4762a1bSJed Brown      nsize: 2
819c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
820c4762a1bSJed Brown 
821c4762a1bSJed Brown    test:
822c4762a1bSJed Brown      args: -snes_rtol 1.e-5 -pc_type mg -ksp_monitor_short -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 -mg_levels_ksp_monitor_short -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full
823c4762a1bSJed Brown      requires: !single
824c4762a1bSJed Brown 
825c4762a1bSJed Brown    test:
826c4762a1bSJed Brown      suffix: 2
827c4762a1bSJed Brown      args: -pc_type mg -ksp_converged_reason -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full -ksp_atol -1.
828c4762a1bSJed Brown 
829c4762a1bSJed Brown    test:
830c4762a1bSJed Brown      suffix: 3
831c4762a1bSJed Brown      nsize: 2
832c4762a1bSJed Brown      args: -snes_grid_sequence 2 -snes_mf_operator -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -snes_rtol 1.e-2
833c4762a1bSJed Brown      filter: grep -v "otal number of function evaluations"
834c4762a1bSJed Brown 
835c4762a1bSJed Brown    test:
836c4762a1bSJed Brown      suffix: 4
837c4762a1bSJed Brown      nsize: 2
838c4762a1bSJed Brown      args: -snes_grid_sequence 2 -snes_monitor_short -ksp_converged_reason -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -ksp_atol -1
839c4762a1bSJed Brown 
840c4762a1bSJed Brown    test:
841c4762a1bSJed Brown      suffix: 5_anderson
842c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type anderson
843c4762a1bSJed Brown 
844c4762a1bSJed Brown    test:
845c4762a1bSJed Brown      suffix: 5_aspin
846c4762a1bSJed Brown      nsize: 4
847c4762a1bSJed Brown      args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type aspin -snes_view
848c4762a1bSJed Brown 
849c4762a1bSJed Brown    test:
850c4762a1bSJed Brown      suffix: 5_broyden
851c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_qn_type broyden -snes_qn_m 10
852c4762a1bSJed Brown 
853c4762a1bSJed Brown    test:
854c4762a1bSJed Brown      suffix: 5_fas
855c4762a1bSJed Brown      args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6
856c4762a1bSJed Brown      requires: !single
857c4762a1bSJed Brown 
858c4762a1bSJed Brown    test:
859c4762a1bSJed Brown      suffix: 5_fas_additive
860c4762a1bSJed Brown      args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 -snes_fas_type additive -snes_max_it 50
861c4762a1bSJed Brown 
862c4762a1bSJed Brown    test:
863c4762a1bSJed Brown      suffix: 5_fas_monitor
864c4762a1bSJed Brown      args: -da_refine 1 -snes_type fas -snes_fas_monitor
865c4762a1bSJed Brown      requires: !single
866c4762a1bSJed Brown 
867c4762a1bSJed Brown    test:
868c4762a1bSJed Brown      suffix: 5_ls
869c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls
870c4762a1bSJed Brown 
871c4762a1bSJed Brown    test:
872c4762a1bSJed Brown      suffix: 5_ls_sell_sor
873c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls -dm_mat_type sell -pc_type sor
874c4762a1bSJed Brown      output_file: output/ex5_5_ls.out
875c4762a1bSJed Brown 
876c4762a1bSJed Brown    test:
877c4762a1bSJed Brown      suffix: 5_nasm
878c4762a1bSJed Brown      nsize: 4
879c4762a1bSJed Brown      args: -snes_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type nasm -snes_nasm_type restrict -snes_max_it 10
880c4762a1bSJed Brown 
881c4762a1bSJed Brown    test:
882c4762a1bSJed Brown      suffix: 5_ncg
883c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ncg -snes_ncg_type fr
884c4762a1bSJed Brown 
885c4762a1bSJed Brown    test:
886c4762a1bSJed Brown      suffix: 5_newton_asm_dmda
887c4762a1bSJed Brown      nsize: 4
888c4762a1bSJed Brown      args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type asm -pc_asm_dm_subdomains -malloc_dump
889c4762a1bSJed Brown      requires: !single
890c4762a1bSJed Brown 
891c4762a1bSJed Brown    test:
892c4762a1bSJed Brown      suffix: 5_newton_gasm_dmda
893c4762a1bSJed Brown      nsize: 4
894c4762a1bSJed Brown      args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type gasm -malloc_dump
895c4762a1bSJed Brown      requires: !single
896c4762a1bSJed Brown 
897c4762a1bSJed Brown    test:
898c4762a1bSJed Brown      suffix: 5_ngmres
899c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10
900c4762a1bSJed Brown 
901c4762a1bSJed Brown    test:
902c4762a1bSJed Brown      suffix: 5_ngmres_fas
903c4762a1bSJed Brown      args: -snes_rtol 1.e-4 -snes_type ngmres -npc_fas_coarse_snes_max_it 1 -npc_fas_coarse_snes_type newtonls -npc_fas_coarse_pc_type lu -npc_fas_coarse_ksp_type preonly -snes_ngmres_m 10 -snes_monitor_short -npc_snes_max_it 1 -npc_snes_type fas -npc_fas_coarse_ksp_type richardson -da_refine 6
904c4762a1bSJed Brown 
905c4762a1bSJed Brown    test:
906c4762a1bSJed Brown      suffix: 5_ngmres_ngs
907c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -npc_snes_type ngs -npc_snes_max_it 1
908c4762a1bSJed Brown 
909c4762a1bSJed Brown    test:
910c4762a1bSJed Brown      suffix: 5_ngmres_nrichardson
911c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 -npc_snes_type nrichardson -npc_snes_max_it 3
912c4762a1bSJed Brown      output_file: output/ex5_5_ngmres_richardson.out
913c4762a1bSJed Brown 
914c4762a1bSJed Brown    test:
915c4762a1bSJed Brown      suffix: 5_nrichardson
916c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type nrichardson
917c4762a1bSJed Brown 
918c4762a1bSJed Brown    test:
919c4762a1bSJed Brown      suffix: 5_qn
920c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_linesearch_type cp -snes_qn_m 10
921c4762a1bSJed Brown 
922c4762a1bSJed Brown    test:
923c4762a1bSJed Brown      suffix: 6
924c4762a1bSJed Brown      nsize: 4
925c4762a1bSJed Brown      args: -snes_converged_reason -ksp_converged_reason -da_grid_x 129 -da_grid_y 129 -pc_type mg -pc_mg_levels 8 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.5,0,1.1 -mg_levels_ksp_max_it 2
926c4762a1bSJed Brown 
927c4762a1bSJed Brown    test:
928c4762a1bSJed Brown      requires: complex !single
929c4762a1bSJed Brown      suffix: complex
930c4762a1bSJed Brown      args: -snes_mf_operator -mat_mffd_complex -snes_monitor
931c4762a1bSJed Brown 
932c4762a1bSJed Brown TEST*/
933