1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Newton's method for a two-variable system, sequential.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: SNES^basic example 6c4762a1bSJed Brown T*/ 7c4762a1bSJed Brown 8c4762a1bSJed Brown /* 9c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 10c4762a1bSJed Brown file automatically includes: 11c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 12c4762a1bSJed Brown petscmat.h - matrices 13c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 14c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 15c4762a1bSJed Brown petscksp.h - linear solvers 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown /*F 18c4762a1bSJed Brown This examples solves either 19c4762a1bSJed Brown \begin{equation} 20c4762a1bSJed 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} 21c4762a1bSJed Brown \end{equation} 22c4762a1bSJed Brown or if the {\tt -hard} options is given 23c4762a1bSJed Brown \begin{equation} 24c4762a1bSJed Brown F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1} 25c4762a1bSJed Brown \end{equation} 26c4762a1bSJed Brown F*/ 27c4762a1bSJed Brown #include <petscsnes.h> 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* 30c4762a1bSJed Brown User-defined routines 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown extern PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*); 33c4762a1bSJed Brown extern PetscErrorCode FormFunction1(SNES,Vec,Vec,void*); 34c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*); 35c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES,Vec,Vec,void*); 36c4762a1bSJed Brown 37c4762a1bSJed Brown int main(int argc,char **argv) 38c4762a1bSJed Brown { 39c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 40c4762a1bSJed Brown KSP ksp; /* linear solver context */ 41c4762a1bSJed Brown PC pc; /* preconditioner context */ 42c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 43c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 44c4762a1bSJed Brown PetscMPIInt size; 45c4762a1bSJed Brown PetscScalar pfive = .5,*xx; 46c4762a1bSJed Brown PetscBool flg; 47c4762a1bSJed Brown 48*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 495f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 502c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Example is only for sequential runs"); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 53c4762a1bSJed Brown Create nonlinear solver context 54c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 555f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 58c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 59c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 60c4762a1bSJed Brown /* 61c4762a1bSJed Brown Create vectors for solution and nonlinear function 62c4762a1bSJed Brown */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,PETSC_DECIDE,2)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* 69c4762a1bSJed Brown Create Jacobian matrix data structure 70c4762a1bSJed Brown */ 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(J)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(J)); 75c4762a1bSJed Brown 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-hard",&flg)); 77c4762a1bSJed Brown if (!flg) { 78c4762a1bSJed Brown /* 79c4762a1bSJed Brown Set function evaluation routine and vector. 80c4762a1bSJed Brown */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction1,NULL)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* 84c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 85c4762a1bSJed Brown */ 865f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian1,NULL)); 87c4762a1bSJed Brown } else { 885f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction2,NULL)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian2,NULL)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 93c4762a1bSJed Brown Customize nonlinear solver; set runtime options 94c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 95c4762a1bSJed Brown /* 96c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 97c4762a1bSJed Brown KSP and PC contexts from the SNES context, we can then 98c4762a1bSJed Brown directly call any KSP and PC routines to set various options. 99c4762a1bSJed Brown */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetKSP(snes,&ksp)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,&pc)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCNONE)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 107c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 108c4762a1bSJed Brown These options will override those specified above as long as 109c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 110c4762a1bSJed Brown routines. 111c4762a1bSJed Brown */ 1125f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 115c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 117c4762a1bSJed Brown if (!flg) { 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,pfive)); 119c4762a1bSJed Brown } else { 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xx)); 121c4762a1bSJed Brown xx[0] = 2.0; xx[1] = 3.0; 1225f80ce2aSJacob Faibussowitsch CHKERRQ(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 1315f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 132c4762a1bSJed Brown if (flg) { 133c4762a1bSJed Brown Vec f; 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,&f,0,0)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(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 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); CHKERRQ(SNESDestroy(&snes)); 146*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 147*b122ec5aSJacob Faibussowitsch return 0; 148c4762a1bSJed Brown } 149c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 150c4762a1bSJed Brown /* 151c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 152c4762a1bSJed Brown 153c4762a1bSJed Brown Input Parameters: 154c4762a1bSJed Brown . snes - the SNES context 155c4762a1bSJed Brown . x - input vector 156c4762a1bSJed Brown . ctx - optional user-defined context 157c4762a1bSJed Brown 158c4762a1bSJed Brown Output Parameter: 159c4762a1bSJed Brown . f - function vector 160c4762a1bSJed Brown */ 161c4762a1bSJed Brown PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *ctx) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown const PetscScalar *xx; 164c4762a1bSJed Brown PetscScalar *ff; 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* 167c4762a1bSJed Brown Get pointers to vector data. 168c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 169c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 170c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 171c4762a1bSJed Brown the array. 172c4762a1bSJed Brown */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Compute function */ 177c4762a1bSJed Brown ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0; 178c4762a1bSJed Brown ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0; 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* Restore vectors */ 1815f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 183c4762a1bSJed Brown return 0; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 188c4762a1bSJed Brown 189c4762a1bSJed Brown Input Parameters: 190c4762a1bSJed Brown . snes - the SNES context 191c4762a1bSJed Brown . x - input vector 192c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 193c4762a1bSJed Brown 194c4762a1bSJed Brown Output Parameters: 195c4762a1bSJed Brown . jac - Jacobian matrix 196c4762a1bSJed Brown . B - optionally different preconditioning matrix 197c4762a1bSJed Brown . flag - flag indicating matrix structure 198c4762a1bSJed Brown */ 199c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 200c4762a1bSJed Brown { 201c4762a1bSJed Brown const PetscScalar *xx; 202c4762a1bSJed Brown PetscScalar A[4]; 203c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* 206c4762a1bSJed Brown Get pointer to vector data 207c4762a1bSJed Brown */ 2085f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* 211c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 212c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 213c4762a1bSJed Brown the matrix at once. 214c4762a1bSJed Brown */ 215c4762a1bSJed Brown A[0] = 2.0*xx[0] + xx[1]; A[1] = xx[0]; 216c4762a1bSJed Brown A[2] = xx[1]; A[3] = xx[0] + 2.0*xx[1]; 2175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES)); 218c4762a1bSJed Brown 219c4762a1bSJed Brown /* 220c4762a1bSJed Brown Restore vector 221c4762a1bSJed Brown */ 2225f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* 225c4762a1bSJed Brown Assemble matrix 226c4762a1bSJed Brown */ 2275f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 229c4762a1bSJed Brown if (jac != B) { 2305f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 232c4762a1bSJed Brown } 233c4762a1bSJed Brown return 0; 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 237c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy) 238c4762a1bSJed Brown { 239c4762a1bSJed Brown const PetscScalar *xx; 240c4762a1bSJed Brown PetscScalar *ff; 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* 243c4762a1bSJed Brown Get pointers to vector data. 244c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 245c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 246c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 247c4762a1bSJed Brown the array. 248c4762a1bSJed Brown */ 2495f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* 253c4762a1bSJed Brown Compute function 254c4762a1bSJed Brown */ 255c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 256c4762a1bSJed Brown ff[1] = xx[1]; 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* 259c4762a1bSJed Brown Restore vectors 260c4762a1bSJed Brown */ 2615f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 263c4762a1bSJed Brown return 0; 264c4762a1bSJed Brown } 265c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 266c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 267c4762a1bSJed Brown { 268c4762a1bSJed Brown const PetscScalar *xx; 269c4762a1bSJed Brown PetscScalar A[4]; 270c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown Get pointer to vector data 274c4762a1bSJed Brown */ 2755f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 276c4762a1bSJed Brown 277c4762a1bSJed Brown /* 278c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 279c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 280c4762a1bSJed Brown the matrix at once. 281c4762a1bSJed Brown */ 282c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 283c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 2845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* 287c4762a1bSJed Brown Restore vector 288c4762a1bSJed Brown */ 2895f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* 292c4762a1bSJed Brown Assemble matrix 293c4762a1bSJed Brown */ 2945f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 296c4762a1bSJed Brown if (jac != B) { 2975f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 299c4762a1bSJed Brown } 300c4762a1bSJed Brown return 0; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303c4762a1bSJed Brown /*TEST 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 307c4762a1bSJed Brown requires: !single 308c4762a1bSJed Brown 309c4762a1bSJed Brown test: 310c4762a1bSJed Brown suffix: 2 311c4762a1bSJed Brown requires: !single 312c4762a1bSJed Brown args: -snes_monitor_short 313c4762a1bSJed Brown output_file: output/ex1_1.out 314c4762a1bSJed Brown 315c4762a1bSJed Brown test: 316c4762a1bSJed Brown suffix: 3 317c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab -snes_monitor_short 318c4762a1bSJed Brown requires: !single 319c4762a1bSJed Brown output_file: output/ex1_1.out 320c4762a1bSJed Brown 321c4762a1bSJed Brown test: 322c4762a1bSJed Brown suffix: 4 323c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp::append -snes_monitor_short 324c4762a1bSJed Brown requires: !single 325c4762a1bSJed Brown output_file: output/ex1_1.out 326c4762a1bSJed Brown 327c4762a1bSJed Brown test: 328c4762a1bSJed Brown suffix: 5 329c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab:append -snes_monitor_short 330c4762a1bSJed Brown requires: !single 331c4762a1bSJed Brown output_file: output/ex1_1.out 332c4762a1bSJed Brown 333c4762a1bSJed Brown test: 334c4762a1bSJed Brown suffix: 6 335c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp:default:append -snes_monitor_short 336c4762a1bSJed Brown requires: !single 337c4762a1bSJed Brown output_file: output/ex1_1.out 338c4762a1bSJed Brown 339c4762a1bSJed Brown test: 340c4762a1bSJed Brown suffix: X 341c4762a1bSJed Brown args: -ksp_monitor_short -ksp_type gmres -ksp_gmres_krylov_monitor -snes_monitor_short -snes_rtol 1.e-4 342c4762a1bSJed Brown requires: !single x 343c4762a1bSJed Brown 344c4762a1bSJed Brown TEST*/ 345