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 449566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,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)); 52*0f0abf79SStefano Zampini PetscCall(SNESSetType(snes,SNESNEWTONLS)); 53*0f0abf79SStefano 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)); 1019566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,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)); 119c4762a1bSJed Brown xx[0] = 2.0; xx[1] = 3.0; 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xx)); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown /* 123c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 124c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 125c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 126c4762a1bSJed Brown this vector to zero by calling VecSet(). 127c4762a1bSJed Brown */ 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 130c4762a1bSJed Brown if (flg) { 131c4762a1bSJed Brown Vec f; 1329566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 1339566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes,&f,0,0)); 1349566063dSJacob Faibussowitsch PetscCall(VecView(r,PETSC_VIEWER_STDOUT_WORLD)); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 139c4762a1bSJed Brown are no longer needed. 140c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); PetscCall(VecDestroy(&r)); 1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); PetscCall(SNESDestroy(&snes)); 1449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 145b122ec5aSJacob Faibussowitsch return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 148c4762a1bSJed Brown /* 149c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 150c4762a1bSJed Brown 151c4762a1bSJed Brown Input Parameters: 152c4762a1bSJed Brown . snes - the SNES context 153c4762a1bSJed Brown . x - input vector 154c4762a1bSJed Brown . ctx - optional user-defined context 155c4762a1bSJed Brown 156c4762a1bSJed Brown Output Parameter: 157c4762a1bSJed Brown . f - function vector 158c4762a1bSJed Brown */ 159c4762a1bSJed Brown PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *ctx) 160c4762a1bSJed Brown { 161c4762a1bSJed Brown const PetscScalar *xx; 162c4762a1bSJed Brown PetscScalar *ff; 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* 165c4762a1bSJed Brown Get pointers to vector data. 166c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 167c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 168c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 169c4762a1bSJed Brown the array. 170c4762a1bSJed Brown */ 1719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1729566063dSJacob Faibussowitsch PetscCall(VecGetArray(f,&ff)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* Compute function */ 175c4762a1bSJed Brown ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0; 176c4762a1bSJed Brown ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0; 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Restore vectors */ 1799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 1809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f,&ff)); 181c4762a1bSJed Brown return 0; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 184c4762a1bSJed Brown /* 185c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 186c4762a1bSJed Brown 187c4762a1bSJed Brown Input Parameters: 188c4762a1bSJed Brown . snes - the SNES context 189c4762a1bSJed Brown . x - input vector 190c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 191c4762a1bSJed Brown 192c4762a1bSJed Brown Output Parameters: 193c4762a1bSJed Brown . jac - Jacobian matrix 194c4762a1bSJed Brown . B - optionally different preconditioning matrix 195c4762a1bSJed Brown . flag - flag indicating matrix structure 196c4762a1bSJed Brown */ 197c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 198c4762a1bSJed Brown { 199c4762a1bSJed Brown const PetscScalar *xx; 200c4762a1bSJed Brown PetscScalar A[4]; 201c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* 204c4762a1bSJed Brown Get pointer to vector data 205c4762a1bSJed Brown */ 2069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* 209c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 210c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 211c4762a1bSJed Brown the matrix at once. 212c4762a1bSJed Brown */ 213c4762a1bSJed Brown A[0] = 2.0*xx[0] + xx[1]; A[1] = xx[0]; 214c4762a1bSJed Brown A[2] = xx[1]; A[3] = xx[0] + 2.0*xx[1]; 2159566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* 218c4762a1bSJed Brown Restore vector 219c4762a1bSJed Brown */ 2209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown Assemble matrix 224c4762a1bSJed Brown */ 2259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 227c4762a1bSJed Brown if (jac != B) { 2289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown return 0; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 235c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy) 236c4762a1bSJed Brown { 237c4762a1bSJed Brown const PetscScalar *xx; 238c4762a1bSJed Brown PetscScalar *ff; 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* 241c4762a1bSJed Brown Get pointers to vector data. 242c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 243c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 244c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 245c4762a1bSJed Brown the array. 246c4762a1bSJed Brown */ 2479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 2489566063dSJacob Faibussowitsch PetscCall(VecGetArray(f,&ff)); 249c4762a1bSJed Brown 250c4762a1bSJed Brown /* 251c4762a1bSJed Brown Compute function 252c4762a1bSJed Brown */ 253c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 254c4762a1bSJed Brown ff[1] = xx[1]; 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* 257c4762a1bSJed Brown Restore vectors 258c4762a1bSJed Brown */ 2599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 2609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f,&ff)); 261c4762a1bSJed Brown return 0; 262c4762a1bSJed Brown } 263c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 264c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 265c4762a1bSJed Brown { 266c4762a1bSJed Brown const PetscScalar *xx; 267c4762a1bSJed Brown PetscScalar A[4]; 268c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 269c4762a1bSJed Brown 270c4762a1bSJed Brown /* 271c4762a1bSJed Brown Get pointer to vector data 272c4762a1bSJed Brown */ 2739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* 276c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 277c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 278c4762a1bSJed Brown the matrix at once. 279c4762a1bSJed Brown */ 280c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 281c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 2829566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES)); 283c4762a1bSJed Brown 284c4762a1bSJed Brown /* 285c4762a1bSJed Brown Restore vector 286c4762a1bSJed Brown */ 2879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* 290c4762a1bSJed Brown Assemble matrix 291c4762a1bSJed Brown */ 2929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 294c4762a1bSJed Brown if (jac != B) { 2959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown return 0; 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 301c4762a1bSJed Brown /*TEST 302c4762a1bSJed Brown 303c4762a1bSJed Brown test: 304*0f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short -prefix_pop 305c4762a1bSJed Brown requires: !single 306c4762a1bSJed Brown 307*0f0abf79SStefano Zampini # test harness puts {{ }} options always at the end, need to specify the prefix explicitly 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: 2 310c4762a1bSJed Brown requires: !single 311*0f0abf79SStefano Zampini args: -prefix_push mysolver_ -snes_monitor_short -prefix_pop -mysolver_snes_ksp_ew {{0 1}} 312c4762a1bSJed Brown output_file: output/ex1_1.out 313c4762a1bSJed Brown 314c4762a1bSJed Brown test: 315c4762a1bSJed Brown suffix: 3 316*0f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab -snes_monitor_short -prefix_pop 317c4762a1bSJed Brown requires: !single 318c4762a1bSJed Brown output_file: output/ex1_1.out 319c4762a1bSJed Brown 320c4762a1bSJed Brown test: 321c4762a1bSJed Brown suffix: 4 322*0f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp::append -snes_monitor_short -prefix_pop 323c4762a1bSJed Brown requires: !single 324c4762a1bSJed Brown output_file: output/ex1_1.out 325c4762a1bSJed Brown 326c4762a1bSJed Brown test: 327c4762a1bSJed Brown suffix: 5 328*0f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab:append -snes_monitor_short -prefix_pop 329c4762a1bSJed Brown requires: !single 330c4762a1bSJed Brown output_file: output/ex1_1.out 331c4762a1bSJed Brown 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 6 334*0f0abf79SStefano Zampini args: -prefix_push mysolver_ -ksp_view_solution ascii:ex1_2_sol.tmp:default:append -snes_monitor_short -prefix_pop 335c4762a1bSJed Brown requires: !single 336c4762a1bSJed Brown output_file: output/ex1_1.out 337c4762a1bSJed Brown 338c4762a1bSJed Brown test: 339c4762a1bSJed Brown suffix: X 340*0f0abf79SStefano 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 341c4762a1bSJed Brown requires: !single x 342c4762a1bSJed Brown 343c4762a1bSJed Brown TEST*/ 344