1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Newton's method for a two-variable system, sequential.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 6c4762a1bSJed Brown file automatically includes: 7c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 8c4762a1bSJed Brown petscmat.h - matrices 9c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 10c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 11c4762a1bSJed Brown petscksp.h - linear solvers 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown /*F 14c4762a1bSJed Brown This examples solves either 15c4762a1bSJed Brown \begin{equation} 16c4762a1bSJed 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} 17c4762a1bSJed Brown \end{equation} 18c4762a1bSJed Brown or if the {\tt -hard} options is given 19c4762a1bSJed Brown \begin{equation} 20c4762a1bSJed Brown F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1} 21c4762a1bSJed Brown \end{equation} 22c4762a1bSJed Brown F*/ 23c4762a1bSJed Brown #include <petscsnes.h> 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* 26c4762a1bSJed Brown User-defined routines 27c4762a1bSJed Brown */ 28c4762a1bSJed Brown extern PetscErrorCode FormJacobian1(SNES, Vec, Mat, Mat, void *); 29c4762a1bSJed Brown extern PetscErrorCode FormFunction1(SNES, Vec, Vec, void *); 30c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES, Vec, Mat, Mat, void *); 31c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES, Vec, Vec, void *); 32c4762a1bSJed Brown 33c4762a1bSJed Brown int main(int argc, char **argv) 34c4762a1bSJed Brown { 35c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 36c4762a1bSJed Brown KSP ksp; /* linear solver context */ 37c4762a1bSJed Brown PC pc; /* preconditioner context */ 38c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 39c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 40c4762a1bSJed Brown PetscMPIInt size; 41c4762a1bSJed Brown PetscScalar pfive = .5, *xx; 42c4762a1bSJed Brown PetscBool flg; 43c4762a1bSJed Brown 44327415f7SBarry Smith PetscFunctionBeginUser; 459566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 47be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example is only for sequential runs"); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 50c4762a1bSJed Brown Create nonlinear solver context 51c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 529566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 530f0abf79SStefano Zampini PetscCall(SNESSetType(snes, SNESNEWTONLS)); 540f0abf79SStefano Zampini PetscCall(SNESSetOptionsPrefix(snes, "mysolver_")); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 57c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 58c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 59c4762a1bSJed Brown /* 60c4762a1bSJed Brown Create vectors for solution and nonlinear function 61c4762a1bSJed Brown */ 629566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 639566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, 2)); 649566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* 68c4762a1bSJed Brown Create Jacobian matrix data structure 69c4762a1bSJed Brown */ 709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 739566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 74c4762a1bSJed Brown 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg)); 76c4762a1bSJed Brown if (!flg) { 77c4762a1bSJed Brown /* 78c4762a1bSJed Brown Set function evaluation routine and vector. 79c4762a1bSJed Brown */ 809566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction1, NULL)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* 83c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 84c4762a1bSJed Brown */ 859566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, NULL)); 86c4762a1bSJed Brown } else { 879566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction2, NULL)); 889566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2, NULL)); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown Customize nonlinear solver; set runtime options 93c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 94c4762a1bSJed Brown /* 95c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 96c4762a1bSJed Brown KSP and PC contexts from the SNES context, we can then 97c4762a1bSJed Brown directly call any KSP and PC routines to set various options. 98c4762a1bSJed Brown */ 999566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 1009566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1019566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 1029566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_DEFAULT, PETSC_DEFAULT, 20)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* 105c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 106c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 107c4762a1bSJed Brown These options will override those specified above as long as 108c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 109c4762a1bSJed Brown routines. 110c4762a1bSJed Brown */ 1119566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 116c4762a1bSJed Brown if (!flg) { 1179566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 118c4762a1bSJed Brown } else { 1199566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xx)); 120*11486bccSBarry Smith xx[0] = 2.0; 121*11486bccSBarry Smith xx[1] = 3.0; 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xx)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown /* 125c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 126c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 127c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 128c4762a1bSJed Brown this vector to zero by calling VecSet(). 129c4762a1bSJed Brown */ 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 132c4762a1bSJed Brown if (flg) { 133c4762a1bSJed Brown Vec f; 1349566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 1359566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &f, 0, 0)); 1369566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 141c4762a1bSJed Brown are no longer needed. 142c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 143c4762a1bSJed Brown 144*11486bccSBarry Smith PetscCall(VecDestroy(&x)); 145*11486bccSBarry Smith PetscCall(VecDestroy(&r)); 146*11486bccSBarry Smith PetscCall(MatDestroy(&J)); 147*11486bccSBarry Smith PetscCall(SNESDestroy(&snes)); 1489566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 149b122ec5aSJacob Faibussowitsch return 0; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 152c4762a1bSJed Brown /* 153c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 154c4762a1bSJed Brown 155c4762a1bSJed Brown Input Parameters: 156c4762a1bSJed Brown . snes - the SNES context 157c4762a1bSJed Brown . x - input vector 158c4762a1bSJed Brown . ctx - optional user-defined context 159c4762a1bSJed Brown 160c4762a1bSJed Brown Output Parameter: 161c4762a1bSJed Brown . f - function vector 162c4762a1bSJed Brown */ 163c4762a1bSJed Brown PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, void *ctx) 164c4762a1bSJed Brown { 165c4762a1bSJed Brown const PetscScalar *xx; 166c4762a1bSJed Brown PetscScalar *ff; 167c4762a1bSJed Brown 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)); 185c4762a1bSJed Brown return 0; 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 199c4762a1bSJed Brown . flag - flag indicating matrix structure 200c4762a1bSJed Brown */ 201c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 202c4762a1bSJed Brown { 203c4762a1bSJed Brown const PetscScalar *xx; 204c4762a1bSJed Brown PetscScalar A[4]; 205c4762a1bSJed Brown PetscInt idx[2] = {0, 1}; 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* 208c4762a1bSJed Brown Get pointer to vector data 209c4762a1bSJed Brown */ 2109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* 213c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 214c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 215c4762a1bSJed Brown the matrix at once. 216c4762a1bSJed Brown */ 217*11486bccSBarry Smith A[0] = 2.0 * xx[0] + xx[1]; 218*11486bccSBarry Smith A[1] = xx[0]; 219*11486bccSBarry Smith A[2] = xx[1]; 220*11486bccSBarry Smith A[3] = xx[0] + 2.0 * xx[1]; 2219566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, idx, 2, idx, A, INSERT_VALUES)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* 224c4762a1bSJed Brown Restore vector 225c4762a1bSJed Brown */ 2269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 227c4762a1bSJed Brown 228c4762a1bSJed Brown /* 229c4762a1bSJed Brown Assemble matrix 230c4762a1bSJed Brown */ 2319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 233c4762a1bSJed Brown if (jac != B) { 2349566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2359566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown return 0; 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 241c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes, Vec x, Vec f, void *dummy) 242c4762a1bSJed Brown { 243c4762a1bSJed Brown const PetscScalar *xx; 244c4762a1bSJed Brown PetscScalar *ff; 245c4762a1bSJed Brown 246c4762a1bSJed Brown /* 247c4762a1bSJed Brown Get pointers to vector data. 248c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 249c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 250c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 251c4762a1bSJed Brown the array. 252c4762a1bSJed Brown */ 2539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 2549566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* 257c4762a1bSJed Brown Compute function 258c4762a1bSJed Brown */ 259c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0]; 260c4762a1bSJed Brown ff[1] = xx[1]; 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* 263c4762a1bSJed Brown Restore vectors 264c4762a1bSJed Brown */ 2659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 267c4762a1bSJed Brown return 0; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 270c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 271c4762a1bSJed Brown { 272c4762a1bSJed Brown const PetscScalar *xx; 273c4762a1bSJed Brown PetscScalar A[4]; 274c4762a1bSJed Brown PetscInt idx[2] = {0, 1}; 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* 277c4762a1bSJed Brown Get pointer to vector data 278c4762a1bSJed Brown */ 2799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* 282c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 283c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 284c4762a1bSJed Brown the matrix at once. 285c4762a1bSJed Brown */ 286*11486bccSBarry Smith A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0; 287*11486bccSBarry Smith A[1] = 0.0; 288*11486bccSBarry Smith A[2] = 0.0; 289*11486bccSBarry Smith A[3] = 1.0; 2909566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, idx, 2, idx, A, INSERT_VALUES)); 291c4762a1bSJed Brown 292c4762a1bSJed Brown /* 293c4762a1bSJed Brown Restore vector 294c4762a1bSJed Brown */ 2959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* 298c4762a1bSJed Brown Assemble matrix 299c4762a1bSJed Brown */ 3009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 3019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 302c4762a1bSJed Brown if (jac != B) { 3039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 3049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 305c4762a1bSJed Brown } 306c4762a1bSJed Brown return 0; 307c4762a1bSJed Brown } 308c4762a1bSJed Brown 309c4762a1bSJed Brown /*TEST 310c4762a1bSJed Brown 311c4762a1bSJed Brown test: 3120f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short -prefix_pop 313c4762a1bSJed Brown requires: !single 314c4762a1bSJed Brown 3150f0abf79SStefano Zampini # test harness puts {{ }} options always at the end, need to specify the prefix explicitly 316c4762a1bSJed Brown test: 317c4762a1bSJed Brown suffix: 2 318c4762a1bSJed Brown requires: !single 3190f0abf79SStefano Zampini args: -prefix_push mysolver_ -snes_monitor_short -prefix_pop -mysolver_snes_ksp_ew {{0 1}} 320c4762a1bSJed Brown output_file: output/ex1_1.out 321c4762a1bSJed Brown 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown suffix: 3 3240f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab -snes_monitor_short -prefix_pop 325c4762a1bSJed Brown requires: !single 326c4762a1bSJed Brown output_file: output/ex1_1.out 327c4762a1bSJed Brown 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: 4 3300f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp::append -snes_monitor_short -prefix_pop 331c4762a1bSJed Brown requires: !single 332c4762a1bSJed Brown output_file: output/ex1_1.out 333c4762a1bSJed Brown 334c4762a1bSJed Brown test: 335c4762a1bSJed Brown suffix: 5 3360f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab:append -snes_monitor_short -prefix_pop 337c4762a1bSJed Brown requires: !single 338c4762a1bSJed Brown output_file: output/ex1_1.out 339c4762a1bSJed Brown 340c4762a1bSJed Brown test: 341c4762a1bSJed Brown suffix: 6 3420f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:default:append -snes_monitor_short -prefix_pop 343c4762a1bSJed Brown requires: !single 344c4762a1bSJed Brown output_file: output/ex1_1.out 345c4762a1bSJed Brown 346c4762a1bSJed Brown test: 347c4762a1bSJed Brown suffix: X 3480f0abf79SStefano 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 349c4762a1bSJed Brown requires: !single x 350c4762a1bSJed Brown 351c4762a1bSJed Brown TEST*/ 352