1c4762a1bSJed Brown static char help[] = "Newton's method for a two-variable system, sequential.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 5c4762a1bSJed Brown file automatically includes: 6c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 7c4762a1bSJed Brown petscmat.h - matrices 8c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 9c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 10c4762a1bSJed Brown petscksp.h - linear solvers 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown /*F 13c4762a1bSJed Brown This examples solves either 14c4762a1bSJed Brown \begin{equation} 15c4762a1bSJed Brown F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{x^2_0 + x_0 x_1 - 3}{x_0 x_1 + x^2_1 - 6} 16c4762a1bSJed Brown \end{equation} 17c4762a1bSJed Brown or if the {\tt -hard} options is given 18c4762a1bSJed Brown \begin{equation} 19c4762a1bSJed Brown F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1} 20c4762a1bSJed Brown \end{equation} 21c4762a1bSJed Brown F*/ 22c4762a1bSJed Brown #include <petscsnes.h> 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* 25c4762a1bSJed Brown User-defined routines 26c4762a1bSJed Brown */ 27c4762a1bSJed Brown extern PetscErrorCode FormJacobian1(SNES, Vec, Mat, Mat, void *); 28c4762a1bSJed Brown extern PetscErrorCode FormFunction1(SNES, Vec, Vec, void *); 29c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES, Vec, Mat, Mat, void *); 30c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES, Vec, Vec, void *); 31c4762a1bSJed Brown 32d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 33d71ae5a4SJacob Faibussowitsch { 34c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 35c4762a1bSJed Brown KSP ksp; /* linear solver context */ 36c4762a1bSJed Brown PC pc; /* preconditioner context */ 37c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 38c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 39c4762a1bSJed Brown PetscMPIInt size; 40c4762a1bSJed Brown PetscScalar pfive = .5, *xx; 41c4762a1bSJed Brown PetscBool flg; 42c4762a1bSJed Brown 43327415f7SBarry Smith PetscFunctionBeginUser; 44c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 46be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example is only for sequential runs"); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 49c4762a1bSJed Brown Create nonlinear solver context 50c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 519566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 520f0abf79SStefano Zampini PetscCall(SNESSetType(snes, SNESNEWTONLS)); 530f0abf79SStefano Zampini PetscCall(SNESSetOptionsPrefix(snes, "mysolver_")); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 57c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 58c4762a1bSJed Brown /* 59c4762a1bSJed Brown Create vectors for solution and nonlinear function 60c4762a1bSJed Brown */ 619566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 629566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, 2)); 639566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown Create Jacobian matrix data structure 68c4762a1bSJed Brown */ 699566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 709566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 719566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 729566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg)); 75c4762a1bSJed Brown if (!flg) { 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Set function evaluation routine and vector. 78c4762a1bSJed Brown */ 799566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction1, NULL)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* 82c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 83c4762a1bSJed Brown */ 849566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, NULL)); 85c4762a1bSJed Brown } else { 869566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction2, NULL)); 879566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2, NULL)); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Customize nonlinear solver; set runtime options 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 93c4762a1bSJed Brown /* 94c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 95c4762a1bSJed Brown KSP and PC contexts from the SNES context, we can then 96c4762a1bSJed Brown directly call any KSP and PC routines to set various options. 97c4762a1bSJed Brown */ 989566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 999566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1009566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 101fb842aefSJose E. Roman PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_CURRENT, PETSC_CURRENT, 20)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* 104c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 105c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 106c4762a1bSJed Brown These options will override those specified above as long as 107c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 108c4762a1bSJed Brown routines. 109c4762a1bSJed Brown */ 1109566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 113c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 114c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 115c4762a1bSJed Brown if (!flg) { 1169566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 117c4762a1bSJed Brown } else { 1189566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xx)); 11911486bccSBarry Smith xx[0] = 2.0; 12011486bccSBarry Smith xx[1] = 3.0; 1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xx)); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown /* 124c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 125c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 126c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 127c4762a1bSJed Brown this vector to zero by calling VecSet(). 128c4762a1bSJed Brown */ 129c4762a1bSJed Brown 1309566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 131c4762a1bSJed Brown if (flg) { 132c4762a1bSJed Brown Vec f; 1339566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 1349566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &f, 0, 0)); 1359566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 140c4762a1bSJed Brown are no longer needed. 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142c4762a1bSJed Brown 14311486bccSBarry Smith PetscCall(VecDestroy(&x)); 14411486bccSBarry Smith PetscCall(VecDestroy(&r)); 14511486bccSBarry Smith PetscCall(MatDestroy(&J)); 14611486bccSBarry Smith PetscCall(SNESDestroy(&snes)); 1479566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 148b122ec5aSJacob Faibussowitsch return 0; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 151c4762a1bSJed Brown /* 152c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 153c4762a1bSJed Brown 154c4762a1bSJed Brown Input Parameters: 155c4762a1bSJed Brown . snes - the SNES context 156c4762a1bSJed Brown . x - input vector 157c4762a1bSJed Brown . ctx - optional user-defined context 158c4762a1bSJed Brown 159c4762a1bSJed Brown Output Parameter: 160c4762a1bSJed Brown . f - function vector 161c4762a1bSJed Brown */ 162d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, void *ctx) 163d71ae5a4SJacob Faibussowitsch { 164c4762a1bSJed Brown const PetscScalar *xx; 165c4762a1bSJed Brown PetscScalar *ff; 166c4762a1bSJed Brown 1673ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown Get pointers to vector data. 170c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 171c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 172c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 173c4762a1bSJed Brown the array. 174c4762a1bSJed Brown */ 1759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1769566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Compute function */ 179c4762a1bSJed Brown ff[0] = xx[0] * xx[0] + xx[0] * xx[1] - 3.0; 180c4762a1bSJed Brown ff[1] = xx[0] * xx[1] + xx[1] * xx[1] - 6.0; 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* Restore vectors */ 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 1853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 188c4762a1bSJed Brown /* 189c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 190c4762a1bSJed Brown 191c4762a1bSJed Brown Input Parameters: 192c4762a1bSJed Brown . snes - the SNES context 193c4762a1bSJed Brown . x - input vector 194c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 195c4762a1bSJed Brown 196c4762a1bSJed Brown Output Parameters: 197c4762a1bSJed Brown . jac - Jacobian matrix 198c4762a1bSJed Brown . B - optionally different preconditioning matrix 199*0b4b7b1cSBarry Smith 200c4762a1bSJed Brown */ 201d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 202d71ae5a4SJacob Faibussowitsch { 203c4762a1bSJed Brown const PetscScalar *xx; 204c4762a1bSJed Brown PetscScalar A[4]; 205c4762a1bSJed Brown PetscInt idx[2] = {0, 1}; 206c4762a1bSJed Brown 2073ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 208c4762a1bSJed Brown /* 209c4762a1bSJed Brown Get pointer to vector data 210c4762a1bSJed Brown */ 2119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* 214c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 215c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 216c4762a1bSJed Brown the matrix at once. 217c4762a1bSJed Brown */ 21811486bccSBarry Smith A[0] = 2.0 * xx[0] + xx[1]; 21911486bccSBarry Smith A[1] = xx[0]; 22011486bccSBarry Smith A[2] = xx[1]; 22111486bccSBarry Smith A[3] = xx[0] + 2.0 * xx[1]; 2229566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, idx, 2, idx, A, INSERT_VALUES)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* 225c4762a1bSJed Brown Restore vector 226c4762a1bSJed Brown */ 2279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* 230c4762a1bSJed Brown Assemble matrix 231c4762a1bSJed Brown */ 2329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 234c4762a1bSJed Brown if (jac != B) { 2359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 237c4762a1bSJed Brown } 2383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 242d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction2(SNES snes, Vec x, Vec f, void *dummy) 243d71ae5a4SJacob Faibussowitsch { 244c4762a1bSJed Brown const PetscScalar *xx; 245c4762a1bSJed Brown PetscScalar *ff; 246c4762a1bSJed Brown 2473ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 248c4762a1bSJed Brown /* 249c4762a1bSJed Brown Get pointers to vector data. 250c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 251c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 252c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 253c4762a1bSJed Brown the array. 254c4762a1bSJed Brown */ 2559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 2569566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* 259c4762a1bSJed Brown Compute function 260c4762a1bSJed Brown */ 261c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0]; 262c4762a1bSJed Brown ff[1] = xx[1]; 263c4762a1bSJed Brown 264c4762a1bSJed Brown /* 265c4762a1bSJed Brown Restore vectors 266c4762a1bSJed Brown */ 2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 272d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian2(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 273d71ae5a4SJacob Faibussowitsch { 274c4762a1bSJed Brown const PetscScalar *xx; 275c4762a1bSJed Brown PetscScalar A[4]; 276c4762a1bSJed Brown PetscInt idx[2] = {0, 1}; 277c4762a1bSJed Brown 2783ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown Get pointer to vector data 281c4762a1bSJed Brown */ 2829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 283c4762a1bSJed Brown 284c4762a1bSJed Brown /* 285c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 286c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 287c4762a1bSJed Brown the matrix at once. 288c4762a1bSJed Brown */ 28911486bccSBarry Smith A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0; 29011486bccSBarry Smith A[1] = 0.0; 29111486bccSBarry Smith A[2] = 0.0; 29211486bccSBarry Smith A[3] = 1.0; 2939566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, idx, 2, idx, A, INSERT_VALUES)); 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* 296c4762a1bSJed Brown Restore vector 297c4762a1bSJed Brown */ 2989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* 301c4762a1bSJed Brown Assemble matrix 302c4762a1bSJed Brown */ 3039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 3049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 305c4762a1bSJed Brown if (jac != B) { 3069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 3079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 308c4762a1bSJed Brown } 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown /*TEST 313c4762a1bSJed Brown 314c4762a1bSJed Brown test: 3150f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short -prefix_pop 316c4762a1bSJed Brown requires: !single 317c4762a1bSJed Brown 3180f0abf79SStefano Zampini # test harness puts {{ }} options always at the end, need to specify the prefix explicitly 319c4762a1bSJed Brown test: 320c4762a1bSJed Brown suffix: 2 321c4762a1bSJed Brown requires: !single 3220f0abf79SStefano Zampini args: -prefix_push mysolver_ -snes_monitor_short -prefix_pop -mysolver_snes_ksp_ew {{0 1}} 323c4762a1bSJed Brown output_file: output/ex1_1.out 324c4762a1bSJed Brown 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 3 3270f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab -snes_monitor_short -prefix_pop 328c4762a1bSJed Brown requires: !single 329c4762a1bSJed Brown output_file: output/ex1_1.out 330c4762a1bSJed Brown 331c4762a1bSJed Brown test: 332c4762a1bSJed Brown suffix: 4 3330f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp::append -snes_monitor_short -prefix_pop 334c4762a1bSJed Brown requires: !single 335c4762a1bSJed Brown output_file: output/ex1_1.out 336c4762a1bSJed Brown 337c4762a1bSJed Brown test: 338c4762a1bSJed Brown suffix: 5 3390f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab:append -snes_monitor_short -prefix_pop 340c4762a1bSJed Brown requires: !single 341c4762a1bSJed Brown output_file: output/ex1_1.out 342c4762a1bSJed Brown 343c4762a1bSJed Brown test: 344c4762a1bSJed Brown suffix: 6 3450f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:default:append -snes_monitor_short -prefix_pop 346c4762a1bSJed Brown requires: !single 347c4762a1bSJed Brown output_file: output/ex1_1.out 348c4762a1bSJed Brown 349c4762a1bSJed Brown test: 350c4762a1bSJed Brown suffix: X 3510f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_monitor_short -ksp_type gmres -ksp_gmres_krylov_monitor -snes_monitor_short -snes_rtol 1.e-4 -prefix_pop 352c4762a1bSJed Brown requires: !single x 353c4762a1bSJed Brown 354c4762a1bSJed Brown TEST*/ 355