1c4762a1bSJed Brown static const char help[] = "-Laplacian u = b as a nonlinear problem.\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown The linear and nonlinear versions of these should give almost identical results on this problem
5c4762a1bSJed Brown
6c4762a1bSJed Brown Richardson
7c4762a1bSJed Brown Nonlinear:
8c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type nrichardson -snes_linesearch_monitor
9c4762a1bSJed Brown
10c4762a1bSJed Brown Linear:
11c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -ksp_rtol 1.e-12 -ksp_monitor -ksp_type richardson -pc_type none -ksp_richardson_self_scale -info
12c4762a1bSJed Brown
13c4762a1bSJed Brown GMRES
14c4762a1bSJed Brown Nonlinear:
15c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type ngmres
16c4762a1bSJed Brown
17c4762a1bSJed Brown Linear:
18c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -ksp_type gmres -ksp_monitor -ksp_rtol 1.e-12 -pc_type none
19c4762a1bSJed Brown
20c4762a1bSJed Brown CG
21c4762a1bSJed Brown Nonlinear:
22c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type ncg -snes_linesearch_monitor
23c4762a1bSJed Brown
24c4762a1bSJed Brown Linear:
25c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -ksp_type cg -ksp_monitor -ksp_rtol 1.e-12 -pc_type none
26c4762a1bSJed Brown
27c4762a1bSJed Brown Multigrid
28c4762a1bSJed Brown Linear:
29c4762a1bSJed Brown 1 level:
30c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type none -mg_levels_ksp_monitor
31c4762a1bSJed Brown -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor -ksp_rtol 1.e-12 -ksp_monitor_true_residual
32c4762a1bSJed Brown
33c4762a1bSJed Brown n levels:
34c4762a1bSJed Brown -da_refine n
35c4762a1bSJed Brown
36c4762a1bSJed Brown Nonlinear:
37c4762a1bSJed Brown 1 level:
38c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type fas -fas_levels_snes_monitor
39c4762a1bSJed Brown
40c4762a1bSJed Brown n levels:
41c4762a1bSJed Brown -da_refine n -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly
42c4762a1bSJed Brown */
43c4762a1bSJed Brown
44c4762a1bSJed Brown /*
45c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
46c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this
47c4762a1bSJed Brown */
48c4762a1bSJed Brown #include <petscdm.h>
49c4762a1bSJed Brown #include <petscdmda.h>
50c4762a1bSJed Brown #include <petscsnes.h>
51c4762a1bSJed Brown
52c4762a1bSJed Brown /*
53c4762a1bSJed Brown User-defined routines
54c4762a1bSJed Brown */
55c4762a1bSJed Brown extern PetscErrorCode FormMatrix(DM, Mat);
56c4762a1bSJed Brown extern PetscErrorCode MyComputeFunction(SNES, Vec, Vec, void *);
57c4762a1bSJed Brown extern PetscErrorCode MyComputeJacobian(SNES, Vec, Mat, Mat, void *);
58c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES, Vec);
59c4762a1bSJed Brown
main(int argc,char ** argv)60d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
61d71ae5a4SJacob Faibussowitsch {
62c4762a1bSJed Brown SNES snes; /* nonlinear solver */
63c4762a1bSJed Brown SNES psnes; /* nonlinear Gauss-Seidel approximate solver */
64c4762a1bSJed Brown Vec x, b; /* solution vector */
65c4762a1bSJed Brown PetscInt its; /* iterations for convergence */
66c4762a1bSJed Brown DM da;
67c4762a1bSJed Brown PetscBool use_ngs_as_npc = PETSC_FALSE; /* use the nonlinear Gauss-Seidel approximate solver */
68c4762a1bSJed Brown
69c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70c4762a1bSJed Brown Initialize program
71c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72c4762a1bSJed Brown
73327415f7SBarry Smith PetscFunctionBeginUser;
74c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
75c4762a1bSJed Brown
76c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77c4762a1bSJed Brown Create nonlinear solver context
78c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
799566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
80c4762a1bSJed Brown
819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_ngs_as_npc", &use_ngs_as_npc, 0));
82c4762a1bSJed Brown
83c4762a1bSJed Brown if (use_ngs_as_npc) {
849566063dSJacob Faibussowitsch PetscCall(SNESGetNPC(snes, &psnes));
859566063dSJacob Faibussowitsch PetscCall(SNESSetType(psnes, SNESSHELL));
869566063dSJacob Faibussowitsch PetscCall(SNESShellSetSolve(psnes, NonlinearGS));
87c4762a1bSJed Brown }
88c4762a1bSJed Brown
89c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
91c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
929566063dSJacob 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));
939566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
949566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
959566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
969566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da));
971baa6e33SBarry Smith if (use_ngs_as_npc) PetscCall(SNESShellSetContext(psnes, da));
98c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining
100c4762a1bSJed Brown vectors that are the same types
101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1029566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x));
1039566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &b));
1049566063dSJacob Faibussowitsch PetscCall(VecSet(b, 1.0));
105c4762a1bSJed Brown
1069566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, MyComputeFunction, NULL));
1079566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, MyComputeJacobian, NULL));
108c4762a1bSJed Brown
109c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown Customize nonlinear solver; set runtime options
111c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1129566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
113c4762a1bSJed Brown
114c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown Solve nonlinear system
116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1179566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, b, x));
1189566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its));
119c4762a1bSJed Brown
120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122c4762a1bSJed Brown
123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
125c4762a1bSJed Brown are no longer needed.
126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
1299566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
1309566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
1319566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
132b122ec5aSJacob Faibussowitsch return 0;
133c4762a1bSJed Brown }
134c4762a1bSJed Brown
135c4762a1bSJed Brown /* ------------------------------------------------------------------- */
MyComputeFunction(SNES snes,Vec x,Vec F,PetscCtx ctx)136*2a8381b2SBarry Smith PetscErrorCode MyComputeFunction(SNES snes, Vec x, Vec F, PetscCtx ctx)
137d71ae5a4SJacob Faibussowitsch {
138c4762a1bSJed Brown Mat J;
139c4762a1bSJed Brown DM dm;
140c4762a1bSJed Brown
141c4762a1bSJed Brown PetscFunctionBeginUser;
1429566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
143c3848c7dSStefano Zampini PetscCall(PetscObjectQuery((PetscObject)dm, "_ex35_J", (PetscObject *)&J));
144c4762a1bSJed Brown if (!J) {
1459566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm, MATAIJ));
1469566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J));
1479566063dSJacob Faibussowitsch PetscCall(MatSetDM(J, NULL));
1489566063dSJacob Faibussowitsch PetscCall(FormMatrix(dm, J));
149c3848c7dSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)dm, "_ex35_J", (PetscObject)J));
150c3848c7dSStefano Zampini PetscCall(PetscObjectDereference((PetscObject)J));
151c4762a1bSJed Brown }
1529566063dSJacob Faibussowitsch PetscCall(MatMult(J, x, F));
1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown
MyComputeJacobian(SNES snes,Vec x,Mat J,Mat Jp,PetscCtx ctx)156*2a8381b2SBarry Smith PetscErrorCode MyComputeJacobian(SNES snes, Vec x, Mat J, Mat Jp, PetscCtx ctx)
157d71ae5a4SJacob Faibussowitsch {
158c4762a1bSJed Brown DM dm;
159c4762a1bSJed Brown
160c4762a1bSJed Brown PetscFunctionBeginUser;
1619566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
1629566063dSJacob Faibussowitsch PetscCall(FormMatrix(dm, Jp));
1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown
FormMatrix(DM da,Mat jac)166d71ae5a4SJacob Faibussowitsch PetscErrorCode FormMatrix(DM da, Mat jac)
167d71ae5a4SJacob Faibussowitsch {
168c4762a1bSJed Brown PetscInt i, j, nrows = 0;
169c4762a1bSJed Brown MatStencil col[5], row, *rows;
170c4762a1bSJed Brown PetscScalar v[5], hx, hy, hxdhy, hydhx;
171c4762a1bSJed Brown DMDALocalInfo info;
172c4762a1bSJed Brown
173c4762a1bSJed Brown PetscFunctionBeginUser;
1749566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
175c4762a1bSJed Brown hx = 1.0 / (PetscReal)(info.mx - 1);
176c4762a1bSJed Brown hy = 1.0 / (PetscReal)(info.my - 1);
177c4762a1bSJed Brown hxdhy = hx / hy;
178c4762a1bSJed Brown hydhx = hy / hx;
179c4762a1bSJed Brown
1809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(info.ym * info.xm, &rows));
181c4762a1bSJed Brown /*
182c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian.
183c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by
184c4762a1bSJed Brown contiguous chunks of rows across the processors.
185c4762a1bSJed Brown - Each processor needs to insert only elements that it owns
186c4762a1bSJed Brown locally (but any non-local elements will be sent to the
187c4762a1bSJed Brown appropriate processor during matrix assembly).
188c4762a1bSJed Brown - Here, we set all entries for a particular row at once.
189c4762a1bSJed Brown - We can set matrix entries either using either
190c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above.
191c4762a1bSJed Brown */
192c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; j++) {
193c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) {
1949371c9d4SSatish Balay row.j = j;
1959371c9d4SSatish Balay row.i = i;
196c4762a1bSJed Brown /* boundary points */
197c4762a1bSJed Brown if (i == 0 || j == 0 || i == info.mx - 1 || j == info.my - 1) {
198c4762a1bSJed Brown v[0] = 2.0 * (hydhx + hxdhy);
1999566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES));
200c4762a1bSJed Brown rows[nrows].i = i;
201c4762a1bSJed Brown rows[nrows++].j = j;
202c4762a1bSJed Brown } else {
203c4762a1bSJed Brown /* interior grid points */
2049371c9d4SSatish Balay v[0] = -hxdhy;
2059371c9d4SSatish Balay col[0].j = j - 1;
2069371c9d4SSatish Balay col[0].i = i;
2079371c9d4SSatish Balay v[1] = -hydhx;
2089371c9d4SSatish Balay col[1].j = j;
2099371c9d4SSatish Balay col[1].i = i - 1;
2109371c9d4SSatish Balay v[2] = 2.0 * (hydhx + hxdhy);
2119371c9d4SSatish Balay col[2].j = row.j;
2129371c9d4SSatish Balay col[2].i = row.i;
2139371c9d4SSatish Balay v[3] = -hydhx;
2149371c9d4SSatish Balay col[3].j = j;
2159371c9d4SSatish Balay col[3].i = i + 1;
2169371c9d4SSatish Balay v[4] = -hxdhy;
2179371c9d4SSatish Balay col[4].j = j + 1;
2189371c9d4SSatish Balay col[4].i = i;
2199566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, col, v, INSERT_VALUES));
220c4762a1bSJed Brown }
221c4762a1bSJed Brown }
222c4762a1bSJed Brown }
223c4762a1bSJed Brown
224c4762a1bSJed Brown /*
225c4762a1bSJed Brown Assemble matrix, using the 2-step process:
226c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd().
227c4762a1bSJed Brown */
2289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
2309566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsStencil(jac, nrows, rows, 2.0 * (hydhx + hxdhy), NULL, NULL));
2319566063dSJacob Faibussowitsch PetscCall(PetscFree(rows));
232c4762a1bSJed Brown /*
233c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the
234c4762a1bSJed Brown matrix. If we do, it will generate an error.
235c4762a1bSJed Brown */
2369566063dSJacob Faibussowitsch PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown
240c4762a1bSJed Brown /* ------------------------------------------------------------------- */
241c4762a1bSJed Brown /*
242c4762a1bSJed Brown Applies some sweeps on nonlinear Gauss-Seidel on each process
243c4762a1bSJed Brown
244c4762a1bSJed Brown */
NonlinearGS(SNES snes,Vec X)245d71ae5a4SJacob Faibussowitsch PetscErrorCode NonlinearGS(SNES snes, Vec X)
246d71ae5a4SJacob Faibussowitsch {
247c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym, its, l;
248c4762a1bSJed Brown PetscReal hx, hy, hxdhy, hydhx;
249c4762a1bSJed Brown PetscScalar **x, F, J, u, uxx, uyy;
250c4762a1bSJed Brown DM da;
251c4762a1bSJed Brown Vec localX;
252c4762a1bSJed Brown
253c4762a1bSJed Brown PetscFunctionBeginUser;
2549566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes, NULL, NULL, NULL, &its, NULL));
2559566063dSJacob Faibussowitsch PetscCall(SNESShellGetContext(snes, &da));
256c4762a1bSJed Brown
2579566063dSJacob Faibussowitsch 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));
258c4762a1bSJed Brown
259c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1);
260c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1);
261c4762a1bSJed Brown hxdhy = hx / hy;
262c4762a1bSJed Brown hydhx = hy / hx;
263c4762a1bSJed Brown
2649566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX));
265c4762a1bSJed Brown
266c4762a1bSJed Brown for (l = 0; l < its; l++) {
2679566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
2689566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
269c4762a1bSJed Brown /*
270c4762a1bSJed Brown Get a pointer to vector data.
271c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to
272c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent.
273c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to
274c4762a1bSJed Brown the array.
275c4762a1bSJed Brown */
2769566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localX, &x));
277c4762a1bSJed Brown
278c4762a1bSJed Brown /*
279c4762a1bSJed Brown Get local grid boundaries (for 2-dimensional DMDA):
280c4762a1bSJed Brown xs, ys - starting grid indices (no ghost points)
281c4762a1bSJed Brown xm, ym - widths of local grid (no ghost points)
282c4762a1bSJed Brown
283c4762a1bSJed Brown */
2849566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
285c4762a1bSJed Brown
286c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
287c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
288c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) {
289c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */
290c4762a1bSJed Brown x[j][i] = 0.0;
291c4762a1bSJed Brown } else {
292c4762a1bSJed Brown u = x[j][i];
293c4762a1bSJed Brown uxx = (2.0 * u - x[j][i - 1] - x[j][i + 1]) * hydhx;
294c4762a1bSJed Brown uyy = (2.0 * u - x[j - 1][i] - x[j + 1][i]) * hxdhy;
295c4762a1bSJed Brown F = uxx + uyy;
296c4762a1bSJed Brown J = 2.0 * (hydhx + hxdhy);
297c4762a1bSJed Brown u = u - F / J;
298c4762a1bSJed Brown
299c4762a1bSJed Brown x[j][i] = u;
300c4762a1bSJed Brown }
301c4762a1bSJed Brown }
302c4762a1bSJed Brown }
303c4762a1bSJed Brown
304c4762a1bSJed Brown /*
305c4762a1bSJed Brown Restore vector
306c4762a1bSJed Brown */
3079566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localX, &x));
3089566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X));
3099566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X));
310c4762a1bSJed Brown }
3119566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX));
3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown
315c4762a1bSJed Brown /*TEST
316c4762a1bSJed Brown
317c4762a1bSJed Brown test:
318c4762a1bSJed Brown args: -snes_monitor_short -snes_type nrichardson
319c4762a1bSJed Brown requires: !single
320c4762a1bSJed Brown
321c4762a1bSJed Brown test:
322c4762a1bSJed Brown suffix: 2
323c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -ksp_type richardson -pc_type none -ksp_richardson_self_scale
324c4762a1bSJed Brown requires: !single
325c4762a1bSJed Brown
326c4762a1bSJed Brown test:
327c4762a1bSJed Brown suffix: 3
328c4762a1bSJed Brown args: -snes_monitor_short -snes_type ngmres
329c4762a1bSJed Brown
330c4762a1bSJed Brown test:
331c4762a1bSJed Brown suffix: 4
332c4762a1bSJed Brown args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none
333c4762a1bSJed Brown
334c4762a1bSJed Brown test:
335c4762a1bSJed Brown suffix: 5
336c4762a1bSJed Brown args: -snes_monitor_short -snes_type ncg
337c4762a1bSJed Brown
338c4762a1bSJed Brown test:
339c4762a1bSJed Brown suffix: 6
340c4762a1bSJed Brown args: -snes_monitor_short -ksp_type cg -ksp_monitor_short -pc_type none
341c4762a1bSJed Brown
342c4762a1bSJed Brown test:
343c4762a1bSJed Brown suffix: 7
344c4762a1bSJed Brown args: -da_refine 2 -snes_monitor_short -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type none -mg_levels_ksp_monitor_short -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor_short
345c4762a1bSJed Brown requires: !single
346c4762a1bSJed Brown
347c4762a1bSJed Brown test:
348c4762a1bSJed Brown suffix: 8
349c4762a1bSJed Brown args: -da_refine 2 -snes_monitor_short -snes_type fas -fas_levels_snes_monitor_short -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_type fas -snes_rtol 1.e-5
350c4762a1bSJed Brown
35141ba4c6cSHeeho Park test:
35241ba4c6cSHeeho Park suffix: 9
35341ba4c6cSHeeho Park args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc
35441ba4c6cSHeeho Park
35541ba4c6cSHeeho Park test:
35641ba4c6cSHeeho Park suffix: 10
35741ba4c6cSHeeho Park args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc -snes_trdc_use_cauchy false
35841ba4c6cSHeeho Park
359c4762a1bSJed Brown TEST*/
360