1c4762a1bSJed Brown static const char help[] = "-Laplacian u = b as a nonlinear problem.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown 5c4762a1bSJed Brown The linear and nonlinear versions of these should give almost identical results on this problem 6c4762a1bSJed Brown 7c4762a1bSJed Brown Richardson 8c4762a1bSJed Brown Nonlinear: 9c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type nrichardson -snes_linesearch_monitor 10c4762a1bSJed Brown 11c4762a1bSJed Brown Linear: 12c4762a1bSJed 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 13c4762a1bSJed Brown 14c4762a1bSJed Brown GMRES 15c4762a1bSJed Brown Nonlinear: 16c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type ngmres 17c4762a1bSJed Brown 18c4762a1bSJed Brown Linear: 19c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -ksp_type gmres -ksp_monitor -ksp_rtol 1.e-12 -pc_type none 20c4762a1bSJed Brown 21c4762a1bSJed Brown CG 22c4762a1bSJed Brown Nonlinear: 23c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type ncg -snes_linesearch_monitor 24c4762a1bSJed Brown 25c4762a1bSJed Brown Linear: 26c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -ksp_type cg -ksp_monitor -ksp_rtol 1.e-12 -pc_type none 27c4762a1bSJed Brown 28c4762a1bSJed Brown Multigrid 29c4762a1bSJed Brown Linear: 30c4762a1bSJed Brown 1 level: 31c4762a1bSJed 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 32c4762a1bSJed Brown -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor -ksp_rtol 1.e-12 -ksp_monitor_true_residual 33c4762a1bSJed Brown 34c4762a1bSJed Brown n levels: 35c4762a1bSJed Brown -da_refine n 36c4762a1bSJed Brown 37c4762a1bSJed Brown Nonlinear: 38c4762a1bSJed Brown 1 level: 39c4762a1bSJed Brown -snes_rtol 1.e-12 -snes_monitor -snes_type fas -fas_levels_snes_monitor 40c4762a1bSJed Brown 41c4762a1bSJed Brown n levels: 42c4762a1bSJed Brown -da_refine n -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly 43c4762a1bSJed Brown 44c4762a1bSJed Brown */ 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* 47c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 48c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown #include <petscdm.h> 51c4762a1bSJed Brown #include <petscdmda.h> 52c4762a1bSJed Brown #include <petscsnes.h> 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* 55c4762a1bSJed Brown User-defined routines 56c4762a1bSJed Brown */ 57c4762a1bSJed Brown extern PetscErrorCode FormMatrix(DM, Mat); 58c4762a1bSJed Brown extern PetscErrorCode MyComputeFunction(SNES, Vec, Vec, void *); 59c4762a1bSJed Brown extern PetscErrorCode MyComputeJacobian(SNES, Vec, Mat, Mat, void *); 60c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES, Vec); 61c4762a1bSJed Brown 62*9371c9d4SSatish Balay int main(int argc, char **argv) { 63c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 64c4762a1bSJed Brown SNES psnes; /* nonlinear Gauss-Seidel approximate solver */ 65c4762a1bSJed Brown Vec x, b; /* solution vector */ 66c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 67c4762a1bSJed Brown DM da; 68c4762a1bSJed Brown PetscBool use_ngs_as_npc = PETSC_FALSE; /* use the nonlinear Gauss-Seidel approximate solver */ 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 71c4762a1bSJed Brown Initialize program 72c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 73c4762a1bSJed Brown 74327415f7SBarry Smith PetscFunctionBeginUser; 759566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78c4762a1bSJed Brown Create nonlinear solver context 79c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 809566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_ngs_as_npc", &use_ngs_as_npc, 0)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown if (use_ngs_as_npc) { 859566063dSJacob Faibussowitsch PetscCall(SNESGetNPC(snes, &psnes)); 869566063dSJacob Faibussowitsch PetscCall(SNESSetType(psnes, SNESSHELL)); 879566063dSJacob Faibussowitsch PetscCall(SNESShellSetSolve(psnes, NonlinearGS)); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 939566063dSJacob 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)); 949566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 959566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 969566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 979566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da)); 981baa6e33SBarry Smith if (use_ngs_as_npc) PetscCall(SNESShellSetContext(psnes, da)); 99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 101c4762a1bSJed Brown vectors that are the same types 102c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1039566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 1049566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &b)); 1059566063dSJacob Faibussowitsch PetscCall(VecSet(b, 1.0)); 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, MyComputeFunction, NULL)); 1089566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, MyComputeJacobian, NULL)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 111c4762a1bSJed Brown Customize nonlinear solver; set runtime options 112c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1139566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Solve nonlinear system 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1189566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, b, x)); 1199566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 126c4762a1bSJed Brown are no longer needed. 127c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1309566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1319566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 133b122ec5aSJacob Faibussowitsch return 0; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 137*9371c9d4SSatish Balay PetscErrorCode MyComputeFunction(SNES snes, Vec x, Vec F, void *ctx) { 138c4762a1bSJed Brown Mat J; 139c4762a1bSJed Brown DM dm; 140c4762a1bSJed Brown 141c4762a1bSJed Brown PetscFunctionBeginUser; 1429566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1439566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &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)); 1499566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, J)); 1509566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContextDestroy(dm, (PetscErrorCode(*)(void **))MatDestroy)); 151c4762a1bSJed Brown } 1529566063dSJacob Faibussowitsch PetscCall(MatMult(J, x, F)); 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156*9371c9d4SSatish Balay PetscErrorCode MyComputeJacobian(SNES snes, Vec x, Mat J, Mat Jp, void *ctx) { 157c4762a1bSJed Brown DM dm; 158c4762a1bSJed Brown 159c4762a1bSJed Brown PetscFunctionBeginUser; 1609566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1619566063dSJacob Faibussowitsch PetscCall(FormMatrix(dm, Jp)); 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165*9371c9d4SSatish Balay PetscErrorCode FormMatrix(DM da, Mat jac) { 166c4762a1bSJed Brown PetscInt i, j, nrows = 0; 167c4762a1bSJed Brown MatStencil col[5], row, *rows; 168c4762a1bSJed Brown PetscScalar v[5], hx, hy, hxdhy, hydhx; 169c4762a1bSJed Brown DMDALocalInfo info; 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscFunctionBeginUser; 1729566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 173c4762a1bSJed Brown hx = 1.0 / (PetscReal)(info.mx - 1); 174c4762a1bSJed Brown hy = 1.0 / (PetscReal)(info.my - 1); 175c4762a1bSJed Brown hxdhy = hx / hy; 176c4762a1bSJed Brown hydhx = hy / hx; 177c4762a1bSJed Brown 1789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(info.ym * info.xm, &rows)); 179c4762a1bSJed Brown /* 180c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 181c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 182c4762a1bSJed Brown contiguous chunks of rows across the processors. 183c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 184c4762a1bSJed Brown locally (but any non-local elements will be sent to the 185c4762a1bSJed Brown appropriate processor during matrix assembly). 186c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 187c4762a1bSJed Brown - We can set matrix entries either using either 188c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above. 189c4762a1bSJed Brown */ 190c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; j++) { 191c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 192*9371c9d4SSatish Balay row.j = j; 193*9371c9d4SSatish Balay row.i = i; 194c4762a1bSJed Brown /* boundary points */ 195c4762a1bSJed Brown if (i == 0 || j == 0 || i == info.mx - 1 || j == info.my - 1) { 196c4762a1bSJed Brown v[0] = 2.0 * (hydhx + hxdhy); 1979566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES)); 198c4762a1bSJed Brown rows[nrows].i = i; 199c4762a1bSJed Brown rows[nrows++].j = j; 200c4762a1bSJed Brown } else { 201c4762a1bSJed Brown /* interior grid points */ 202*9371c9d4SSatish Balay v[0] = -hxdhy; 203*9371c9d4SSatish Balay col[0].j = j - 1; 204*9371c9d4SSatish Balay col[0].i = i; 205*9371c9d4SSatish Balay v[1] = -hydhx; 206*9371c9d4SSatish Balay col[1].j = j; 207*9371c9d4SSatish Balay col[1].i = i - 1; 208*9371c9d4SSatish Balay v[2] = 2.0 * (hydhx + hxdhy); 209*9371c9d4SSatish Balay col[2].j = row.j; 210*9371c9d4SSatish Balay col[2].i = row.i; 211*9371c9d4SSatish Balay v[3] = -hydhx; 212*9371c9d4SSatish Balay col[3].j = j; 213*9371c9d4SSatish Balay col[3].i = i + 1; 214*9371c9d4SSatish Balay v[4] = -hxdhy; 215*9371c9d4SSatish Balay col[4].j = j + 1; 216*9371c9d4SSatish Balay col[4].i = i; 2179566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, col, v, INSERT_VALUES)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown Assemble matrix, using the 2-step process: 224c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 225c4762a1bSJed Brown */ 2269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 2289566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsStencil(jac, nrows, rows, 2.0 * (hydhx + hxdhy), NULL, NULL)); 2299566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 230c4762a1bSJed Brown /* 231c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 232c4762a1bSJed Brown matrix. If we do, it will generate an error. 233c4762a1bSJed Brown */ 2349566063dSJacob Faibussowitsch PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 235c4762a1bSJed Brown PetscFunctionReturn(0); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 239c4762a1bSJed Brown /* 240c4762a1bSJed Brown Applies some sweeps on nonlinear Gauss-Seidel on each process 241c4762a1bSJed Brown 242c4762a1bSJed Brown */ 243*9371c9d4SSatish Balay PetscErrorCode NonlinearGS(SNES snes, Vec X) { 244c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym, its, l; 245c4762a1bSJed Brown PetscReal hx, hy, hxdhy, hydhx; 246c4762a1bSJed Brown PetscScalar **x, F, J, u, uxx, uyy; 247c4762a1bSJed Brown DM da; 248c4762a1bSJed Brown Vec localX; 249c4762a1bSJed Brown 250c4762a1bSJed Brown PetscFunctionBeginUser; 2519566063dSJacob Faibussowitsch PetscCall(SNESGetTolerances(snes, NULL, NULL, NULL, &its, NULL)); 2529566063dSJacob Faibussowitsch PetscCall(SNESShellGetContext(snes, &da)); 253c4762a1bSJed Brown 2549566063dSJacob 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)); 255c4762a1bSJed Brown 256c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1); 257c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1); 258c4762a1bSJed Brown hxdhy = hx / hy; 259c4762a1bSJed Brown hydhx = hy / hx; 260c4762a1bSJed Brown 2619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown for (l = 0; l < its; l++) { 2649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 2659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 266c4762a1bSJed Brown /* 267c4762a1bSJed Brown Get a pointer to vector data. 268c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 269c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 270c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 271c4762a1bSJed Brown the array. 272c4762a1bSJed Brown */ 2739566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localX, &x)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* 276c4762a1bSJed Brown Get local grid boundaries (for 2-dimensional DMDA): 277c4762a1bSJed Brown xs, ys - starting grid indices (no ghost points) 278c4762a1bSJed Brown xm, ym - widths of local grid (no ghost points) 279c4762a1bSJed Brown 280c4762a1bSJed Brown */ 2819566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 282c4762a1bSJed Brown 283c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 284c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 285c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 286c4762a1bSJed Brown /* boundary conditions are all zero Dirichlet */ 287c4762a1bSJed Brown x[j][i] = 0.0; 288c4762a1bSJed Brown } else { 289c4762a1bSJed Brown u = x[j][i]; 290c4762a1bSJed Brown uxx = (2.0 * u - x[j][i - 1] - x[j][i + 1]) * hydhx; 291c4762a1bSJed Brown uyy = (2.0 * u - x[j - 1][i] - x[j + 1][i]) * hxdhy; 292c4762a1bSJed Brown F = uxx + uyy; 293c4762a1bSJed Brown J = 2.0 * (hydhx + hxdhy); 294c4762a1bSJed Brown u = u - F / J; 295c4762a1bSJed Brown 296c4762a1bSJed Brown x[j][i] = u; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown } 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 301c4762a1bSJed Brown /* 302c4762a1bSJed Brown Restore vector 303c4762a1bSJed Brown */ 3049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localX, &x)); 3059566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X)); 3069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X)); 307c4762a1bSJed Brown } 3089566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 309c4762a1bSJed Brown PetscFunctionReturn(0); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown /*TEST 313c4762a1bSJed Brown 314c4762a1bSJed Brown test: 315c4762a1bSJed Brown args: -snes_monitor_short -snes_type nrichardson 316c4762a1bSJed Brown requires: !single 317c4762a1bSJed Brown 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: 2 320c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -ksp_type richardson -pc_type none -ksp_richardson_self_scale 321c4762a1bSJed Brown requires: !single 322c4762a1bSJed Brown 323c4762a1bSJed Brown test: 324c4762a1bSJed Brown suffix: 3 325c4762a1bSJed Brown args: -snes_monitor_short -snes_type ngmres 326c4762a1bSJed Brown 327c4762a1bSJed Brown test: 328c4762a1bSJed Brown suffix: 4 329c4762a1bSJed Brown args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none 330c4762a1bSJed Brown 331c4762a1bSJed Brown test: 332c4762a1bSJed Brown suffix: 5 333c4762a1bSJed Brown args: -snes_monitor_short -snes_type ncg 334c4762a1bSJed Brown 335c4762a1bSJed Brown test: 336c4762a1bSJed Brown suffix: 6 337c4762a1bSJed Brown args: -snes_monitor_short -ksp_type cg -ksp_monitor_short -pc_type none 338c4762a1bSJed Brown 339c4762a1bSJed Brown test: 340c4762a1bSJed Brown suffix: 7 341c4762a1bSJed 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 342c4762a1bSJed Brown requires: !single 343c4762a1bSJed Brown 344c4762a1bSJed Brown test: 345c4762a1bSJed Brown suffix: 8 346c4762a1bSJed 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 347c4762a1bSJed Brown 34841ba4c6cSHeeho Park test: 34941ba4c6cSHeeho Park suffix: 9 35041ba4c6cSHeeho Park args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc 35141ba4c6cSHeeho Park 35241ba4c6cSHeeho Park test: 35341ba4c6cSHeeho Park suffix: 10 35441ba4c6cSHeeho Park args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc -snes_trdc_use_cauchy false 35541ba4c6cSHeeho Park 356c4762a1bSJed Brown TEST*/ 357