xref: /petsc/src/snes/tutorials/ex5.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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  */
FormInitialGuess(DM da,AppCtx * user,Vec X)73d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormInitialGuess(DM da, AppCtx *user, Vec X)
74d71ae5a4SJacob 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));
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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  */
FormExactSolution(DM da,AppCtx * user,Vec U)136d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U)
137d71ae5a4SJacob 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) {
1513ba16761SJacob Faibussowitsch     for (i = xs; i < xs + xm; ++i) PetscCall(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));
1553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156227f9e03SJunchao Zhang }
157227f9e03SJunchao Zhang 
ZeroBCSolution(AppCtx * user,const DMDACoor2d * c,PetscScalar * u)158d71ae5a4SJacob Faibussowitsch static PetscErrorCode ZeroBCSolution(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
159d71ae5a4SJacob Faibussowitsch {
160227f9e03SJunchao Zhang   u[0] = 0.;
1613ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
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 
168dd8e379bSPierre Jolivet   such that u(x,y) is an exact solution with f(x,y) as the right-hand side forcing term.
169227f9e03SJunchao Zhang  */
MMSSolution1(AppCtx * user,const DMDACoor2d * c,PetscScalar * u)170d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution1(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
171d71ae5a4SJacob Faibussowitsch {
172227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
1733ba16761SJacob Faibussowitsch 
1743ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
175227f9e03SJunchao Zhang   u[0] = x * (1 - x) * y * (1 - y);
1763ba16761SJacob Faibussowitsch   PetscCall(PetscLogFlops(5));
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178227f9e03SJunchao Zhang }
MMSForcing1(AppCtx * user,const DMDACoor2d * c,PetscScalar * f)179d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing1(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
180d71ae5a4SJacob Faibussowitsch {
181227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
1823ba16761SJacob Faibussowitsch 
1833ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
184227f9e03SJunchao Zhang   f[0] = 2 * x * (1 - x) + 2 * y * (1 - y) - user->param * PetscExpReal(x * (1 - x) * y * (1 - y));
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
186227f9e03SJunchao Zhang }
187227f9e03SJunchao Zhang 
MMSSolution2(AppCtx * user,const DMDACoor2d * c,PetscScalar * u)188d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution2(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
189d71ae5a4SJacob Faibussowitsch {
190227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
1913ba16761SJacob Faibussowitsch 
1923ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
193227f9e03SJunchao Zhang   u[0] = PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y);
1943ba16761SJacob Faibussowitsch   PetscCall(PetscLogFlops(5));
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
196227f9e03SJunchao Zhang }
MMSForcing2(AppCtx * user,const DMDACoor2d * c,PetscScalar * f)197d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing2(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
198d71ae5a4SJacob Faibussowitsch {
199227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
2003ba16761SJacob Faibussowitsch 
2013ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
202227f9e03SJunchao 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));
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204227f9e03SJunchao Zhang }
205227f9e03SJunchao Zhang 
MMSSolution3(AppCtx * user,const DMDACoor2d * c,PetscScalar * u)206d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution3(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
207d71ae5a4SJacob Faibussowitsch {
208227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
2093ba16761SJacob Faibussowitsch 
2103ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
211227f9e03SJunchao Zhang   u[0] = PetscSinReal(user->m * PETSC_PI * x * (1 - y)) * PetscSinReal(user->n * PETSC_PI * y * (1 - x));
2123ba16761SJacob Faibussowitsch   PetscCall(PetscLogFlops(5));
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
214227f9e03SJunchao Zhang }
MMSForcing3(AppCtx * user,const DMDACoor2d * c,PetscScalar * f)215d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing3(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
216d71ae5a4SJacob Faibussowitsch {
217227f9e03SJunchao Zhang   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
218227f9e03SJunchao Zhang   PetscReal m = user->m, n = user->n, lambda = user->param;
2193ba16761SJacob Faibussowitsch 
2203ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2219371c9d4SSatish 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)));
2223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
223227f9e03SJunchao Zhang }
224227f9e03SJunchao Zhang 
MMSSolution4(AppCtx * user,const DMDACoor2d * c,PetscScalar * u)225d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSSolution4(AppCtx *user, const DMDACoor2d *c, PetscScalar *u)
226d71ae5a4SJacob Faibussowitsch {
227227f9e03SJunchao Zhang   const PetscReal Lx = 1., Ly = 1.;
228227f9e03SJunchao Zhang   PetscReal       x = PetscRealPart(c->x), y = PetscRealPart(c->y);
2293ba16761SJacob Faibussowitsch 
2303ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
231227f9e03SJunchao Zhang   u[0] = (PetscPowReal(x, 4) - PetscSqr(Lx) * PetscSqr(x)) * (PetscPowReal(y, 4) - PetscSqr(Ly) * PetscSqr(y));
2323ba16761SJacob Faibussowitsch   PetscCall(PetscLogFlops(9));
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234227f9e03SJunchao Zhang }
MMSForcing4(AppCtx * user,const DMDACoor2d * c,PetscScalar * f)235d71ae5a4SJacob Faibussowitsch static PetscErrorCode MMSForcing4(AppCtx *user, const DMDACoor2d *c, PetscScalar *f)
236d71ae5a4SJacob Faibussowitsch {
237227f9e03SJunchao Zhang   const PetscReal Lx = 1., Ly = 1.;
238227f9e03SJunchao Zhang   PetscReal       x = PetscRealPart(c->x), y = PetscRealPart(c->y);
2393ba16761SJacob Faibussowitsch 
2403ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2419371c9d4SSatish 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))));
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243227f9e03SJunchao Zhang }
244227f9e03SJunchao Zhang 
245227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */
246227f9e03SJunchao Zhang /*
247227f9e03SJunchao Zhang    FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch
248227f9e03SJunchao Zhang 
249227f9e03SJunchao Zhang  */
FormFunctionLocal(DMDALocalInfo * info,PetscScalar ** x,PetscScalar ** f,AppCtx * user)250d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user)
251d71ae5a4SJacob Faibussowitsch {
252227f9e03SJunchao Zhang   PetscInt    i, j;
253227f9e03SJunchao Zhang   PetscReal   lambda, hx, hy, hxdhy, hydhx;
254227f9e03SJunchao Zhang   PetscScalar u, ue, uw, un, us, uxx, uyy, mms_solution, mms_forcing;
255227f9e03SJunchao Zhang   DMDACoor2d  c;
256227f9e03SJunchao Zhang 
257227f9e03SJunchao Zhang   PetscFunctionBeginUser;
258227f9e03SJunchao Zhang   lambda = user->param;
259227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(info->mx - 1);
260227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(info->my - 1);
261227f9e03SJunchao Zhang   hxdhy  = hx / hy;
262227f9e03SJunchao Zhang   hydhx  = hy / hx;
263227f9e03SJunchao Zhang   /*
264227f9e03SJunchao Zhang      Compute function over the locally owned part of the grid
265227f9e03SJunchao Zhang   */
266227f9e03SJunchao Zhang   for (j = info->ys; j < info->ys + info->ym; j++) {
267227f9e03SJunchao Zhang     for (i = info->xs; i < info->xs + info->xm; i++) {
268227f9e03SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
2699371c9d4SSatish Balay         c.x = i * hx;
2709371c9d4SSatish Balay         c.y = j * hy;
271227f9e03SJunchao Zhang         PetscCall(user->mms_solution(user, &c, &mms_solution));
272227f9e03SJunchao Zhang         f[j][i] = 2.0 * (hydhx + hxdhy) * (x[j][i] - mms_solution);
273227f9e03SJunchao Zhang       } else {
274227f9e03SJunchao Zhang         u  = x[j][i];
275227f9e03SJunchao Zhang         uw = x[j][i - 1];
276227f9e03SJunchao Zhang         ue = x[j][i + 1];
277227f9e03SJunchao Zhang         un = x[j - 1][i];
278227f9e03SJunchao Zhang         us = x[j + 1][i];
279227f9e03SJunchao Zhang 
280227f9e03SJunchao Zhang         /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */
2819371c9d4SSatish Balay         if (i - 1 == 0) {
2829371c9d4SSatish Balay           c.x = (i - 1) * hx;
2839371c9d4SSatish Balay           c.y = j * hy;
2849371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &uw));
2859371c9d4SSatish Balay         }
2869371c9d4SSatish Balay         if (i + 1 == info->mx - 1) {
2879371c9d4SSatish Balay           c.x = (i + 1) * hx;
2889371c9d4SSatish Balay           c.y = j * hy;
2899371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &ue));
2909371c9d4SSatish Balay         }
2919371c9d4SSatish Balay         if (j - 1 == 0) {
2929371c9d4SSatish Balay           c.x = i * hx;
2939371c9d4SSatish Balay           c.y = (j - 1) * hy;
2949371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &un));
2959371c9d4SSatish Balay         }
2969371c9d4SSatish Balay         if (j + 1 == info->my - 1) {
2979371c9d4SSatish Balay           c.x = i * hx;
2989371c9d4SSatish Balay           c.y = (j + 1) * hy;
2999371c9d4SSatish Balay           PetscCall(user->mms_solution(user, &c, &us));
3009371c9d4SSatish Balay         }
301227f9e03SJunchao Zhang 
302227f9e03SJunchao Zhang         uxx         = (2.0 * u - uw - ue) * hydhx;
303227f9e03SJunchao Zhang         uyy         = (2.0 * u - un - us) * hxdhy;
304227f9e03SJunchao Zhang         mms_forcing = 0;
3059371c9d4SSatish Balay         c.x         = i * hx;
3069371c9d4SSatish Balay         c.y         = j * hy;
307227f9e03SJunchao Zhang         if (user->mms_forcing) PetscCall(user->mms_forcing(user, &c, &mms_forcing));
308227f9e03SJunchao Zhang         f[j][i] = uxx + uyy - hx * hy * (lambda * PetscExpScalar(u) + mms_forcing);
309227f9e03SJunchao Zhang       }
310227f9e03SJunchao Zhang     }
311227f9e03SJunchao Zhang   }
312227f9e03SJunchao Zhang   PetscCall(PetscLogFlops(11.0 * info->ym * info->xm));
3133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
314227f9e03SJunchao Zhang }
315227f9e03SJunchao Zhang 
316227f9e03SJunchao Zhang /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */
FormObjectiveLocal(DMDALocalInfo * info,PetscScalar ** x,PetscReal * obj,AppCtx * user)317d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info, PetscScalar **x, PetscReal *obj, AppCtx *user)
318d71ae5a4SJacob Faibussowitsch {
319227f9e03SJunchao Zhang   PetscInt    i, j;
320227f9e03SJunchao Zhang   PetscReal   lambda, hx, hy, hxdhy, hydhx, sc, lobj = 0;
321227f9e03SJunchao Zhang   PetscScalar u, ue, uw, un, us, uxux, uyuy;
322227f9e03SJunchao Zhang   MPI_Comm    comm;
323227f9e03SJunchao Zhang 
324227f9e03SJunchao Zhang   PetscFunctionBeginUser;
325227f9e03SJunchao Zhang   *obj = 0;
326227f9e03SJunchao Zhang   PetscCall(PetscObjectGetComm((PetscObject)info->da, &comm));
327227f9e03SJunchao Zhang   lambda = user->param;
328227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(info->mx - 1);
329227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(info->my - 1);
330227f9e03SJunchao Zhang   sc     = hx * hy * lambda;
331227f9e03SJunchao Zhang   hxdhy  = hx / hy;
332227f9e03SJunchao Zhang   hydhx  = hy / hx;
333227f9e03SJunchao Zhang   /*
334227f9e03SJunchao Zhang      Compute function over the locally owned part of the grid
335227f9e03SJunchao Zhang   */
336227f9e03SJunchao Zhang   for (j = info->ys; j < info->ys + info->ym; j++) {
337227f9e03SJunchao Zhang     for (i = info->xs; i < info->xs + info->xm; i++) {
338227f9e03SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
339227f9e03SJunchao Zhang         lobj += PetscRealPart((hydhx + hxdhy) * x[j][i] * x[j][i]);
340227f9e03SJunchao Zhang       } else {
341227f9e03SJunchao Zhang         u  = x[j][i];
342227f9e03SJunchao Zhang         uw = x[j][i - 1];
343227f9e03SJunchao Zhang         ue = x[j][i + 1];
344227f9e03SJunchao Zhang         un = x[j - 1][i];
345227f9e03SJunchao Zhang         us = x[j + 1][i];
346227f9e03SJunchao Zhang 
347227f9e03SJunchao Zhang         if (i - 1 == 0) uw = 0.;
348227f9e03SJunchao Zhang         if (i + 1 == info->mx - 1) ue = 0.;
349227f9e03SJunchao Zhang         if (j - 1 == 0) un = 0.;
350227f9e03SJunchao Zhang         if (j + 1 == info->my - 1) us = 0.;
351227f9e03SJunchao Zhang 
352227f9e03SJunchao Zhang         /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */
353227f9e03SJunchao Zhang 
354227f9e03SJunchao Zhang         uxux = u * (2. * u - ue - uw) * hydhx;
355227f9e03SJunchao Zhang         uyuy = u * (2. * u - un - us) * hxdhy;
356227f9e03SJunchao Zhang 
357227f9e03SJunchao Zhang         lobj += PetscRealPart(0.5 * (uxux + uyuy) - sc * PetscExpScalar(u));
358227f9e03SJunchao Zhang       }
359227f9e03SJunchao Zhang     }
360227f9e03SJunchao Zhang   }
361227f9e03SJunchao Zhang   PetscCall(PetscLogFlops(12.0 * info->ym * info->xm));
3620bf52853SStefano Zampini   *obj = lobj;
3633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
364227f9e03SJunchao Zhang }
365227f9e03SJunchao Zhang 
366227f9e03SJunchao Zhang /*
367227f9e03SJunchao Zhang    FormJacobianLocal - Evaluates Jacobian matrix on local process patch
368227f9e03SJunchao Zhang */
FormJacobianLocal(DMDALocalInfo * info,PetscScalar ** x,Mat jac,Mat jacpre,AppCtx * user)369d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat jac, Mat jacpre, AppCtx *user)
370d71ae5a4SJacob Faibussowitsch {
371227f9e03SJunchao Zhang   PetscInt     i, j, k;
372227f9e03SJunchao Zhang   MatStencil   col[5], row;
373227f9e03SJunchao Zhang   PetscScalar  lambda, v[5], hx, hy, hxdhy, hydhx, sc;
374227f9e03SJunchao Zhang   DM           coordDA;
375227f9e03SJunchao Zhang   Vec          coordinates;
376227f9e03SJunchao Zhang   DMDACoor2d **coords;
377227f9e03SJunchao Zhang 
378227f9e03SJunchao Zhang   PetscFunctionBeginUser;
379227f9e03SJunchao Zhang   lambda = user->param;
380227f9e03SJunchao Zhang   /* Extract coordinates */
381227f9e03SJunchao Zhang   PetscCall(DMGetCoordinateDM(info->da, &coordDA));
382227f9e03SJunchao Zhang   PetscCall(DMGetCoordinates(info->da, &coordinates));
383227f9e03SJunchao Zhang   PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords));
384227f9e03SJunchao Zhang   hx = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs + 1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0;
385227f9e03SJunchao Zhang   hy = info->ym > 1 ? PetscRealPart(coords[info->ys + 1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0;
386227f9e03SJunchao Zhang   PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords));
387227f9e03SJunchao Zhang   hxdhy = hx / hy;
388227f9e03SJunchao Zhang   hydhx = hy / hx;
389227f9e03SJunchao Zhang   sc    = hx * hy * lambda;
390227f9e03SJunchao Zhang 
391227f9e03SJunchao Zhang   /*
392227f9e03SJunchao Zhang      Compute entries for the locally owned part of the Jacobian.
393227f9e03SJunchao Zhang       - Currently, all PETSc parallel matrix formats are partitioned by
394227f9e03SJunchao Zhang         contiguous chunks of rows across the processors.
395227f9e03SJunchao Zhang       - Each processor needs to insert only elements that it owns
396227f9e03SJunchao Zhang         locally (but any non-local elements will be sent to the
397227f9e03SJunchao Zhang         appropriate processor during matrix assembly).
398227f9e03SJunchao Zhang       - Here, we set all entries for a particular row at once.
399227f9e03SJunchao Zhang       - We can set matrix entries either using either
400227f9e03SJunchao Zhang         MatSetValuesLocal() or MatSetValues(), as discussed above.
401227f9e03SJunchao Zhang   */
402227f9e03SJunchao Zhang   for (j = info->ys; j < info->ys + info->ym; j++) {
403227f9e03SJunchao Zhang     for (i = info->xs; i < info->xs + info->xm; i++) {
4049371c9d4SSatish Balay       row.j = j;
4059371c9d4SSatish Balay       row.i = i;
406227f9e03SJunchao Zhang       /* boundary points */
407227f9e03SJunchao Zhang       if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) {
408227f9e03SJunchao Zhang         v[0] = 2.0 * (hydhx + hxdhy);
409227f9e03SJunchao Zhang         PetscCall(MatSetValuesStencil(jacpre, 1, &row, 1, &row, v, INSERT_VALUES));
410227f9e03SJunchao Zhang       } else {
411227f9e03SJunchao Zhang         k = 0;
412227f9e03SJunchao Zhang         /* interior grid points */
413227f9e03SJunchao Zhang         if (j - 1 != 0) {
414227f9e03SJunchao Zhang           v[k]     = -hxdhy;
4159371c9d4SSatish Balay           col[k].j = j - 1;
4169371c9d4SSatish Balay           col[k].i = i;
417227f9e03SJunchao Zhang           k++;
418227f9e03SJunchao Zhang         }
419227f9e03SJunchao Zhang         if (i - 1 != 0) {
420227f9e03SJunchao Zhang           v[k]     = -hydhx;
4219371c9d4SSatish Balay           col[k].j = j;
4229371c9d4SSatish Balay           col[k].i = i - 1;
423227f9e03SJunchao Zhang           k++;
424227f9e03SJunchao Zhang         }
425227f9e03SJunchao Zhang 
4269371c9d4SSatish Balay         v[k]     = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(x[j][i]);
4279371c9d4SSatish Balay         col[k].j = row.j;
4289371c9d4SSatish Balay         col[k].i = row.i;
4299371c9d4SSatish Balay         k++;
430227f9e03SJunchao Zhang 
431227f9e03SJunchao Zhang         if (i + 1 != info->mx - 1) {
432227f9e03SJunchao Zhang           v[k]     = -hydhx;
4339371c9d4SSatish Balay           col[k].j = j;
4349371c9d4SSatish Balay           col[k].i = i + 1;
435227f9e03SJunchao Zhang           k++;
436227f9e03SJunchao Zhang         }
437227f9e03SJunchao Zhang         if (j + 1 != info->mx - 1) {
438227f9e03SJunchao Zhang           v[k]     = -hxdhy;
4399371c9d4SSatish Balay           col[k].j = j + 1;
4409371c9d4SSatish Balay           col[k].i = i;
441227f9e03SJunchao Zhang           k++;
442227f9e03SJunchao Zhang         }
443227f9e03SJunchao Zhang         PetscCall(MatSetValuesStencil(jacpre, 1, &row, k, col, v, INSERT_VALUES));
444227f9e03SJunchao Zhang       }
445227f9e03SJunchao Zhang     }
446227f9e03SJunchao Zhang   }
447227f9e03SJunchao Zhang 
448227f9e03SJunchao Zhang   /*
449227f9e03SJunchao Zhang      Assemble matrix, using the 2-step process:
450227f9e03SJunchao Zhang        MatAssemblyBegin(), MatAssemblyEnd().
451227f9e03SJunchao Zhang   */
452227f9e03SJunchao Zhang   PetscCall(MatAssemblyBegin(jacpre, MAT_FINAL_ASSEMBLY));
453227f9e03SJunchao Zhang   PetscCall(MatAssemblyEnd(jacpre, MAT_FINAL_ASSEMBLY));
454227f9e03SJunchao Zhang   /*
455227f9e03SJunchao Zhang      Tell the matrix we will never add a new nonzero location to the
456227f9e03SJunchao Zhang      matrix. If we do, it will generate an error.
457227f9e03SJunchao Zhang   */
458227f9e03SJunchao Zhang   PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
4593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460227f9e03SJunchao Zhang }
461227f9e03SJunchao Zhang 
FormFunctionMatlab(SNES snes,Vec X,Vec F,void * ptr)462d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormFunctionMatlab(SNES snes, Vec X, Vec F, void *ptr)
463d71ae5a4SJacob Faibussowitsch {
464d1e78c4fSBarry Smith #if PetscDefined(HAVE_MATLAB)
465227f9e03SJunchao Zhang   AppCtx   *user = (AppCtx *)ptr;
466227f9e03SJunchao Zhang   PetscInt  Mx, My;
467227f9e03SJunchao Zhang   PetscReal lambda, hx, hy;
468227f9e03SJunchao Zhang   Vec       localX, localF;
469227f9e03SJunchao Zhang   MPI_Comm  comm;
470227f9e03SJunchao Zhang   DM        da;
471227f9e03SJunchao Zhang 
472227f9e03SJunchao Zhang   PetscFunctionBeginUser;
473227f9e03SJunchao Zhang   PetscCall(SNESGetDM(snes, &da));
474227f9e03SJunchao Zhang   PetscCall(DMGetLocalVector(da, &localX));
475227f9e03SJunchao Zhang   PetscCall(DMGetLocalVector(da, &localF));
476227f9e03SJunchao Zhang   PetscCall(PetscObjectSetName((PetscObject)localX, "localX"));
477227f9e03SJunchao Zhang   PetscCall(PetscObjectSetName((PetscObject)localF, "localF"));
478227f9e03SJunchao 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));
479227f9e03SJunchao Zhang 
480227f9e03SJunchao Zhang   lambda = user->param;
481227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(Mx - 1);
482227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(My - 1);
483227f9e03SJunchao Zhang 
484227f9e03SJunchao Zhang   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
485227f9e03SJunchao Zhang   /*
486227f9e03SJunchao Zhang      Scatter ghost points to local vector,using the 2-step process
487227f9e03SJunchao Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
488227f9e03SJunchao Zhang      By placing code between these two statements, computations can be
489227f9e03SJunchao Zhang      done while messages are in transition.
490227f9e03SJunchao Zhang   */
491227f9e03SJunchao Zhang   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
492227f9e03SJunchao Zhang   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
493227f9e03SJunchao Zhang   PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localX));
494227f9e03SJunchao Zhang   PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm), "localF=ex5m(localX,%18.16e,%18.16e,%18.16e)", (double)hx, (double)hy, (double)lambda));
495227f9e03SJunchao Zhang   PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localF));
496227f9e03SJunchao Zhang 
497227f9e03SJunchao Zhang   /*
498227f9e03SJunchao Zhang      Insert values into global vector
499227f9e03SJunchao Zhang   */
500227f9e03SJunchao Zhang   PetscCall(DMLocalToGlobalBegin(da, localF, INSERT_VALUES, F));
501227f9e03SJunchao Zhang   PetscCall(DMLocalToGlobalEnd(da, localF, INSERT_VALUES, F));
502227f9e03SJunchao Zhang   PetscCall(DMRestoreLocalVector(da, &localX));
503227f9e03SJunchao Zhang   PetscCall(DMRestoreLocalVector(da, &localF));
5043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
505227f9e03SJunchao Zhang #else
5063ba16761SJacob Faibussowitsch   return PETSC_SUCCESS; /* Never called */
507227f9e03SJunchao Zhang #endif
508227f9e03SJunchao Zhang }
509227f9e03SJunchao Zhang 
510227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */
511227f9e03SJunchao Zhang /*
512227f9e03SJunchao Zhang       Applies some sweeps on nonlinear Gauss-Seidel on each process
513227f9e03SJunchao Zhang 
514227f9e03SJunchao Zhang  */
NonlinearGS(SNES snes,Vec X,Vec B,PetscCtx ctx)515*2a8381b2SBarry Smith static PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, PetscCtx ctx)
516d71ae5a4SJacob Faibussowitsch {
517227f9e03SJunchao Zhang   PetscInt      i, j, k, Mx, My, xs, ys, xm, ym, its, tot_its, sweeps, l;
518227f9e03SJunchao Zhang   PetscReal     lambda, hx, hy, hxdhy, hydhx, sc;
519227f9e03SJunchao Zhang   PetscScalar **x, **b, bij, F, F0 = 0, J, u, un, us, ue, eu, uw, uxx, uyy, y;
520227f9e03SJunchao Zhang   PetscReal     atol, rtol, stol;
521227f9e03SJunchao Zhang   DM            da;
522227f9e03SJunchao Zhang   AppCtx       *user;
523227f9e03SJunchao Zhang   Vec           localX, localB;
524227f9e03SJunchao Zhang 
525227f9e03SJunchao Zhang   PetscFunctionBeginUser;
526227f9e03SJunchao Zhang   tot_its = 0;
527227f9e03SJunchao Zhang   PetscCall(SNESNGSGetSweeps(snes, &sweeps));
528227f9e03SJunchao Zhang   PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its));
529227f9e03SJunchao Zhang   PetscCall(SNESGetDM(snes, &da));
530227f9e03SJunchao Zhang   PetscCall(DMGetApplicationContext(da, &user));
531227f9e03SJunchao Zhang 
532227f9e03SJunchao 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));
533227f9e03SJunchao Zhang 
534227f9e03SJunchao Zhang   lambda = user->param;
535227f9e03SJunchao Zhang   hx     = 1.0 / (PetscReal)(Mx - 1);
536227f9e03SJunchao Zhang   hy     = 1.0 / (PetscReal)(My - 1);
537227f9e03SJunchao Zhang   sc     = hx * hy * lambda;
538227f9e03SJunchao Zhang   hxdhy  = hx / hy;
539227f9e03SJunchao Zhang   hydhx  = hy / hx;
540227f9e03SJunchao Zhang 
541227f9e03SJunchao Zhang   PetscCall(DMGetLocalVector(da, &localX));
54248a46eb9SPierre Jolivet   if (B) PetscCall(DMGetLocalVector(da, &localB));
543227f9e03SJunchao Zhang   for (l = 0; l < sweeps; l++) {
544227f9e03SJunchao Zhang     PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
545227f9e03SJunchao Zhang     PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
546227f9e03SJunchao Zhang     if (B) {
547227f9e03SJunchao Zhang       PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, localB));
548227f9e03SJunchao Zhang       PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, localB));
549227f9e03SJunchao Zhang     }
550227f9e03SJunchao Zhang     /*
551227f9e03SJunchao Zhang      Get a pointer to vector data.
552227f9e03SJunchao Zhang      - For default PETSc vectors, VecGetArray() returns a pointer to
553227f9e03SJunchao Zhang      the data array.  Otherwise, the routine is implementation dependent.
554227f9e03SJunchao Zhang      - You MUST call VecRestoreArray() when you no longer need access to
555227f9e03SJunchao Zhang      the array.
556227f9e03SJunchao Zhang      */
557227f9e03SJunchao Zhang     PetscCall(DMDAVecGetArray(da, localX, &x));
558227f9e03SJunchao Zhang     if (B) PetscCall(DMDAVecGetArray(da, localB, &b));
559227f9e03SJunchao Zhang     /*
560227f9e03SJunchao Zhang      Get local grid boundaries (for 2-dimensional DMDA):
561227f9e03SJunchao Zhang      xs, ys   - starting grid indices (no ghost points)
562227f9e03SJunchao Zhang      xm, ym   - widths of local grid (no ghost points)
563227f9e03SJunchao Zhang      */
564227f9e03SJunchao Zhang     PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
565227f9e03SJunchao Zhang 
566227f9e03SJunchao Zhang     for (j = ys; j < ys + ym; j++) {
567227f9e03SJunchao Zhang       for (i = xs; i < xs + xm; i++) {
568227f9e03SJunchao Zhang         if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
569227f9e03SJunchao Zhang           /* boundary conditions are all zero Dirichlet */
570227f9e03SJunchao Zhang           x[j][i] = 0.0;
571227f9e03SJunchao Zhang         } else {
572227f9e03SJunchao Zhang           if (B) bij = b[j][i];
573227f9e03SJunchao Zhang           else bij = 0.;
574227f9e03SJunchao Zhang 
575227f9e03SJunchao Zhang           u  = x[j][i];
576227f9e03SJunchao Zhang           un = x[j - 1][i];
577227f9e03SJunchao Zhang           us = x[j + 1][i];
578227f9e03SJunchao Zhang           ue = x[j][i - 1];
579227f9e03SJunchao Zhang           uw = x[j][i + 1];
580227f9e03SJunchao Zhang 
581227f9e03SJunchao Zhang           for (k = 0; k < its; k++) {
582227f9e03SJunchao Zhang             eu  = PetscExpScalar(u);
583227f9e03SJunchao Zhang             uxx = (2.0 * u - ue - uw) * hydhx;
584227f9e03SJunchao Zhang             uyy = (2.0 * u - un - us) * hxdhy;
585227f9e03SJunchao Zhang             F   = uxx + uyy - sc * eu - bij;
586227f9e03SJunchao Zhang             if (k == 0) F0 = F;
587227f9e03SJunchao Zhang             J = 2.0 * (hydhx + hxdhy) - sc * eu;
588227f9e03SJunchao Zhang             y = F / J;
589227f9e03SJunchao Zhang             u -= y;
590227f9e03SJunchao Zhang             tot_its++;
591227f9e03SJunchao Zhang 
592ad540459SPierre Jolivet             if (atol > PetscAbsReal(PetscRealPart(F)) || rtol * PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || stol * PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) break;
593227f9e03SJunchao Zhang           }
594227f9e03SJunchao Zhang           x[j][i] = u;
595227f9e03SJunchao Zhang         }
596227f9e03SJunchao Zhang       }
597227f9e03SJunchao Zhang     }
598227f9e03SJunchao Zhang     /*
599227f9e03SJunchao Zhang      Restore vector
600227f9e03SJunchao Zhang      */
601227f9e03SJunchao Zhang     PetscCall(DMDAVecRestoreArray(da, localX, &x));
602227f9e03SJunchao Zhang     PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X));
603227f9e03SJunchao Zhang     PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X));
604227f9e03SJunchao Zhang   }
605227f9e03SJunchao Zhang   PetscCall(PetscLogFlops(tot_its * (21.0)));
606227f9e03SJunchao Zhang   PetscCall(DMRestoreLocalVector(da, &localX));
607227f9e03SJunchao Zhang   if (B) {
608227f9e03SJunchao Zhang     PetscCall(DMDAVecRestoreArray(da, localB, &b));
609227f9e03SJunchao Zhang     PetscCall(DMRestoreLocalVector(da, &localB));
610227f9e03SJunchao Zhang   }
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
612227f9e03SJunchao Zhang }
613c4762a1bSJed Brown 
main(int argc,char ** argv)614d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
615d71ae5a4SJacob Faibussowitsch {
616c4762a1bSJed Brown   SNES      snes; /* nonlinear solver */
617c4762a1bSJed Brown   Vec       x;    /* solution vector */
618c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
619c4762a1bSJed Brown   PetscInt  its;  /* iterations for convergence */
620c4762a1bSJed Brown   PetscReal bratu_lambda_max = 6.81;
621c4762a1bSJed Brown   PetscReal bratu_lambda_min = 0.;
62216d037c0SJunchao Zhang   PetscInt  MMS              = 1;
62316d037c0SJunchao Zhang   PetscBool flg              = PETSC_FALSE, setMMS;
624c4762a1bSJed Brown   DM        da;
625c4762a1bSJed Brown   Vec       r = NULL;
626c4762a1bSJed Brown   KSP       ksp;
627c4762a1bSJed Brown   PetscInt  lits, slits;
628c4762a1bSJed Brown 
629c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
630c4762a1bSJed Brown      Initialize program
631c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
632c4762a1bSJed Brown 
633327415f7SBarry Smith   PetscFunctionBeginUser;
634c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
635c4762a1bSJed Brown 
636c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
637c4762a1bSJed Brown      Initialize problem parameters
638c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
639c4762a1bSJed Brown   user.param = 6.0;
6409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL));
641e00437b9SBarry 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);
64216d037c0SJunchao Zhang   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mms", &MMS, &setMMS));
643c4762a1bSJed Brown   if (MMS == 3) {
644c4762a1bSJed Brown     PetscInt mPar = 2, nPar = 1;
6459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-m_par", &mPar, NULL));
6469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_par", &nPar, NULL));
647c4762a1bSJed Brown     user.m = PetscPowInt(2, mPar);
648c4762a1bSJed Brown     user.n = PetscPowInt(2, nPar);
649c4762a1bSJed Brown   }
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
652c4762a1bSJed Brown      Create nonlinear solver context
653c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
6549566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
6559566063dSJacob Faibussowitsch   PetscCall(SNESSetCountersReset(snes, PETSC_FALSE));
6569566063dSJacob Faibussowitsch   PetscCall(SNESSetNGS(snes, NonlinearGS, NULL));
657c4762a1bSJed Brown 
658c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
659c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
660c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
6619566063dSJacob 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));
6629566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
6639566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
6649566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
6659566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
6669566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
667c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
668c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
669c4762a1bSJed Brown      vectors that are the same types
670c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
6719566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
672c4762a1bSJed Brown 
673c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
674c4762a1bSJed Brown      Set local function evaluation routine
675c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
676c4762a1bSJed Brown   switch (MMS) {
6779371c9d4SSatish Balay   case 0:
6789371c9d4SSatish Balay     user.mms_solution = ZeroBCSolution;
6799371c9d4SSatish Balay     user.mms_forcing  = NULL;
6809371c9d4SSatish Balay     break;
6819371c9d4SSatish Balay   case 1:
6829371c9d4SSatish Balay     user.mms_solution = MMSSolution1;
6839371c9d4SSatish Balay     user.mms_forcing  = MMSForcing1;
6849371c9d4SSatish Balay     break;
6859371c9d4SSatish Balay   case 2:
6869371c9d4SSatish Balay     user.mms_solution = MMSSolution2;
6879371c9d4SSatish Balay     user.mms_forcing  = MMSForcing2;
6889371c9d4SSatish Balay     break;
6899371c9d4SSatish Balay   case 3:
6909371c9d4SSatish Balay     user.mms_solution = MMSSolution3;
6919371c9d4SSatish Balay     user.mms_forcing  = MMSForcing3;
6929371c9d4SSatish Balay     break;
6939371c9d4SSatish Balay   case 4:
6949371c9d4SSatish Balay     user.mms_solution = MMSSolution4;
6959371c9d4SSatish Balay     user.mms_forcing  = MMSForcing4;
6969371c9d4SSatish Balay     break;
697d71ae5a4SJacob Faibussowitsch   default:
698d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Unknown MMS type %" PetscInt_FMT, MMS);
699c4762a1bSJed Brown   }
7008434afd1SBarry Smith   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunctionFn *)FormFunctionLocal, &user));
7019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-fd", &flg, NULL));
7028434afd1SBarry Smith   if (!flg) PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobianFn *)FormJacobianLocal, &user));
703c4762a1bSJed Brown 
7049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-obj", &flg, NULL));
7058434afd1SBarry Smith   if (flg) PetscCall(DMDASNESSetObjectiveLocal(da, (DMDASNESObjectiveFn *)FormObjectiveLocal, &user));
706c4762a1bSJed Brown 
707d1e78c4fSBarry Smith   if (PetscDefined(HAVE_MATLAB)) {
7084e1ad211SJed Brown     PetscBool matlab_function = PETSC_FALSE;
7099566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-matlab_function", &matlab_function, 0));
710c4762a1bSJed Brown     if (matlab_function) {
7119566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(x, &r));
7129566063dSJacob Faibussowitsch       PetscCall(SNESSetFunction(snes, r, FormFunctionMatlab, &user));
713c4762a1bSJed Brown     }
7144e1ad211SJed Brown   }
715c4762a1bSJed Brown 
716c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
717c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
718c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
7199566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
720c4762a1bSJed Brown 
7219566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(da, &user, x));
722c4762a1bSJed Brown 
723c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
724c4762a1bSJed Brown      Solve nonlinear system
725c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
7269566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
7279566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
728c4762a1bSJed Brown 
7299566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes, &slits));
7309566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
7319566063dSJacob Faibussowitsch   PetscCall(KSPGetTotalIterations(ksp, &lits));
73263a3b9bcSJacob 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);
733c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
734c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
735c4762a1bSJed Brown 
736c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
737c4762a1bSJed Brown      If using MMS, check the l_2 error
738c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73916d037c0SJunchao Zhang   if (setMMS) {
740c4762a1bSJed Brown     Vec       e;
741c4762a1bSJed Brown     PetscReal errorl2, errorinf;
742c4762a1bSJed Brown     PetscInt  N;
743c4762a1bSJed Brown 
7449566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &e));
7459566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)x, NULL, "-sol_view"));
7469566063dSJacob Faibussowitsch     PetscCall(FormExactSolution(da, &user, e));
7479566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-exact_view"));
7489566063dSJacob Faibussowitsch     PetscCall(VecAXPY(e, -1.0, x));
7499566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-error_view"));
7509566063dSJacob Faibussowitsch     PetscCall(VecNorm(e, NORM_2, &errorl2));
7519566063dSJacob Faibussowitsch     PetscCall(VecNorm(e, NORM_INFINITY, &errorinf));
7529566063dSJacob Faibussowitsch     PetscCall(VecGetSize(e, &N));
75363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " error L2 %g inf %g\n", N, (double)(errorl2 / PetscSqrtReal((PetscReal)N)), (double)errorinf));
7549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&e));
7559566063dSJacob Faibussowitsch     PetscCall(PetscLogEventSetDof(SNES_Solve, 0, N));
7569566063dSJacob Faibussowitsch     PetscCall(PetscLogEventSetError(SNES_Solve, 0, errorl2 / PetscSqrtReal(N)));
757c4762a1bSJed Brown   }
758c4762a1bSJed Brown 
759c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
760c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
761c4762a1bSJed Brown      are no longer needed.
762c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
7639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
7649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
7659566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
7669566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
7679566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
768b122ec5aSJacob Faibussowitsch   return 0;
769c4762a1bSJed Brown }
770c4762a1bSJed Brown 
771c4762a1bSJed Brown /*TEST
772c4762a1bSJed Brown 
773c4762a1bSJed Brown    test:
774c4762a1bSJed Brown      suffix: asm_0
775c4762a1bSJed Brown      requires: !single
776c4762a1bSJed 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
777c4762a1bSJed Brown 
778c4762a1bSJed Brown    test:
779c4762a1bSJed Brown      suffix: msm_0
780c4762a1bSJed Brown      requires: !single
781c4762a1bSJed 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
782c4762a1bSJed Brown 
783c4762a1bSJed Brown    test:
784c4762a1bSJed Brown      suffix: asm_1
785c4762a1bSJed Brown      requires: !single
786c4762a1bSJed 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
787c4762a1bSJed Brown 
788c4762a1bSJed Brown    test:
789c4762a1bSJed Brown      suffix: msm_1
790c4762a1bSJed Brown      requires: !single
791c4762a1bSJed 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
792c4762a1bSJed Brown 
793c4762a1bSJed Brown    test:
794c4762a1bSJed Brown      suffix: asm_2
795c4762a1bSJed Brown      requires: !single
796c4762a1bSJed 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
797c4762a1bSJed Brown 
798c4762a1bSJed Brown    test:
799c4762a1bSJed Brown      suffix: msm_2
800c4762a1bSJed Brown      requires: !single
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 3 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
802c4762a1bSJed Brown 
803c4762a1bSJed Brown    test:
804c4762a1bSJed Brown      suffix: asm_3
805c4762a1bSJed Brown      requires: !single
806c4762a1bSJed 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
807c4762a1bSJed Brown 
808c4762a1bSJed Brown    test:
809c4762a1bSJed Brown      suffix: msm_3
810c4762a1bSJed Brown      requires: !single
811c4762a1bSJed 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
812c4762a1bSJed Brown 
813c4762a1bSJed Brown    test:
814c4762a1bSJed Brown      suffix: asm_4
815c4762a1bSJed Brown      requires: !single
816c4762a1bSJed Brown      nsize: 2
817c4762a1bSJed 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
818c4762a1bSJed Brown 
819c4762a1bSJed Brown    test:
820c4762a1bSJed Brown      suffix: msm_4
821c4762a1bSJed Brown      requires: !single
822c4762a1bSJed Brown      nsize: 2
823c4762a1bSJed 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
824c4762a1bSJed Brown 
825c4762a1bSJed Brown    test:
826c4762a1bSJed Brown      suffix: asm_5
827c4762a1bSJed Brown      requires: !single
828c4762a1bSJed Brown      nsize: 2
829c4762a1bSJed 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
830c4762a1bSJed Brown 
831c4762a1bSJed Brown    test:
832c4762a1bSJed Brown      suffix: msm_5
833c4762a1bSJed Brown      requires: !single
834c4762a1bSJed Brown      nsize: 2
835c4762a1bSJed 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
836c4762a1bSJed Brown 
837c4762a1bSJed Brown    test:
838c4762a1bSJed 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
839c4762a1bSJed Brown      requires: !single
840c4762a1bSJed Brown 
841c4762a1bSJed Brown    test:
842c4762a1bSJed Brown      suffix: 2
84377e5a1f9SBarry Smith      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 0
844c4762a1bSJed Brown 
845c4762a1bSJed Brown    test:
846c4762a1bSJed Brown      suffix: 3
847c4762a1bSJed Brown      nsize: 2
84877e5a1f9SBarry Smith      args: -snes_grid_sequence 2 -snes_mf_operator -snes_converged_reason -snes_view -pc_type mg -snes_atol 0 -snes_rtol 1.e-2
849c4762a1bSJed Brown      filter: grep -v "otal number of function evaluations"
850c4762a1bSJed Brown 
851c4762a1bSJed Brown    test:
852c4762a1bSJed Brown      suffix: 4
853c4762a1bSJed Brown      nsize: 2
85477e5a1f9SBarry Smith      args: -snes_grid_sequence 2 -snes_monitor_short -ksp_converged_reason -snes_converged_reason -snes_view -pc_type mg -snes_atol 0 -ksp_atol 0
855c4762a1bSJed Brown 
856c4762a1bSJed Brown    test:
857c4762a1bSJed Brown      suffix: 5_anderson
858c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type anderson
859c4762a1bSJed Brown 
860c4762a1bSJed Brown    test:
861c4762a1bSJed Brown      suffix: 5_aspin
862c4762a1bSJed Brown      nsize: 4
863de54d9edSStefano Zampini      args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type aspin -snes_view -npc_sub_pc_type lu -npc_sub_ksp_type preonly
864c4762a1bSJed Brown 
865c4762a1bSJed Brown    test:
866c4762a1bSJed Brown      suffix: 5_broyden
867c4762a1bSJed 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
868c4762a1bSJed Brown 
869c4762a1bSJed Brown    test:
870c4762a1bSJed Brown      suffix: 5_fas
871c4762a1bSJed 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
872c4762a1bSJed Brown      requires: !single
873c4762a1bSJed Brown 
874c4762a1bSJed Brown    test:
875c4762a1bSJed Brown      suffix: 5_fas_additive
876a99ef635SJonas Heinzmann      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 -snes_linesearch_maxlambda 2.0
877c4762a1bSJed Brown 
878c4762a1bSJed Brown    test:
879c4762a1bSJed Brown      suffix: 5_fas_monitor
880c4762a1bSJed Brown      args: -da_refine 1 -snes_type fas -snes_fas_monitor
881c4762a1bSJed Brown      requires: !single
882c4762a1bSJed Brown 
883c4762a1bSJed Brown    test:
884c4762a1bSJed Brown      suffix: 5_ls
885c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls
886c4762a1bSJed Brown 
887c4762a1bSJed Brown    test:
888c4762a1bSJed Brown      suffix: 5_ls_sell_sor
889c4762a1bSJed 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
890c4762a1bSJed Brown      output_file: output/ex5_5_ls.out
891c4762a1bSJed Brown 
892c4762a1bSJed Brown    test:
893c4762a1bSJed Brown      suffix: 5_nasm
894c4762a1bSJed Brown      nsize: 4
895c4762a1bSJed 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
896c4762a1bSJed Brown 
897c4762a1bSJed Brown    test:
898c4762a1bSJed Brown      suffix: 5_ncg
899c4762a1bSJed 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
900c4762a1bSJed Brown 
901c4762a1bSJed Brown    test:
902c4762a1bSJed Brown      suffix: 5_newton_asm_dmda
903c4762a1bSJed Brown      nsize: 4
904c4762a1bSJed 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
905c4762a1bSJed Brown      requires: !single
906c4762a1bSJed Brown 
907c4762a1bSJed Brown    test:
908c4762a1bSJed Brown      suffix: 5_newton_gasm_dmda
909c4762a1bSJed Brown      nsize: 4
910c4762a1bSJed 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
911c4762a1bSJed Brown      requires: !single
912c4762a1bSJed Brown 
913c4762a1bSJed Brown    test:
914c4762a1bSJed Brown      suffix: 5_ngmres
915c4762a1bSJed 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
916c4762a1bSJed Brown 
917c4762a1bSJed Brown    test:
918c4762a1bSJed Brown      suffix: 5_ngmres_fas
919c4762a1bSJed 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
920c4762a1bSJed Brown 
921c4762a1bSJed Brown    test:
922c4762a1bSJed Brown      suffix: 5_ngmres_ngs
923c4762a1bSJed 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
924c4762a1bSJed Brown 
925c4762a1bSJed Brown    test:
926c4762a1bSJed Brown      suffix: 5_ngmres_nrichardson
927c4762a1bSJed 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
928c4762a1bSJed Brown      output_file: output/ex5_5_ngmres_richardson.out
929c4762a1bSJed Brown 
930c4762a1bSJed Brown    test:
931c4762a1bSJed Brown      suffix: 5_nrichardson
932c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type nrichardson
933c4762a1bSJed Brown 
934c4762a1bSJed Brown    test:
935c4762a1bSJed Brown      suffix: 5_qn
936c4762a1bSJed 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
937c4762a1bSJed Brown 
938c4762a1bSJed Brown    test:
939c4762a1bSJed Brown      suffix: 6
940c4762a1bSJed Brown      nsize: 4
941c4762a1bSJed 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
942c4762a1bSJed Brown 
943c4762a1bSJed Brown    test:
944c4762a1bSJed Brown      requires: complex !single
945c4762a1bSJed Brown      suffix: complex
946c4762a1bSJed Brown      args: -snes_mf_operator -mat_mffd_complex -snes_monitor
947c4762a1bSJed Brown 
9480acf423eSBarry Smith    test:
9490acf423eSBarry Smith      requires: !single
9500acf423eSBarry Smith      suffix: 7_ksp_view_pre
9510acf423eSBarry Smith      args: -pc_type gamg -ksp_view_pre
9520acf423eSBarry Smith 
953bfa5bdfcSBarry Smith    test:
954bfa5bdfcSBarry Smith      requires: !single
955bfa5bdfcSBarry Smith      suffix: hem_view_detailed
956c376f201SBarry Smith      args: -pc_type gamg -ksp_view ::ascii_info_detail -pc_gamg_mat_coarsen_type hem
957bfa5bdfcSBarry Smith 
958bfa5bdfcSBarry Smith    test:
959bfa5bdfcSBarry Smith      requires: !single
960bfa5bdfcSBarry Smith      suffix: mis_view_detailed
961c376f201SBarry Smith      args: -pc_type gamg -ksp_view ::ascii_info_detail -pc_gamg_mat_coarsen_type mis
962bfa5bdfcSBarry Smith 
963c4762a1bSJed Brown TEST*/
964