1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Newton's method to solve a two-variable system that comes from the Rosenbrock function.\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 #include <petscsnes.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown extern PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*); 20c4762a1bSJed Brown extern PetscErrorCode FormFunction1(SNES,Vec,Vec,void*); 21c4762a1bSJed Brown 22c4762a1bSJed Brown int main(int argc,char **argv) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 25c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 26c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 27c4762a1bSJed Brown PetscErrorCode ierr; 28c4762a1bSJed Brown PetscInt its; 29c4762a1bSJed Brown PetscScalar *xx; 30c4762a1bSJed Brown SNESConvergedReason reason; 31c4762a1bSJed Brown 32c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 35c4762a1bSJed Brown Create nonlinear solver context 36c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 37c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 40c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 41c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 42c4762a1bSJed Brown /* 43c4762a1bSJed Brown Create vectors for solution and nonlinear function 44c4762a1bSJed Brown */ 45c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = VecSetSizes(x,PETSC_DECIDE,2);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* 51c4762a1bSJed Brown Create Jacobian matrix data structure 52c4762a1bSJed Brown */ 53c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = MatSetFromOptions(J);CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* 59c4762a1bSJed Brown Set function evaluation routine and vector. 60c4762a1bSJed Brown */ 61c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction1,NULL);CHKERRQ(ierr); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 65c4762a1bSJed Brown */ 66c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian1,NULL);CHKERRQ(ierr); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 69c4762a1bSJed Brown Customize nonlinear solver; set runtime options 70c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 71c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 75c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 76c4762a1bSJed Brown ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 77c4762a1bSJed Brown xx[0] = -1.2; xx[1] = 1.0; 78c4762a1bSJed Brown ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* 81c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 82c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 83c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 84c4762a1bSJed Brown this vector to zero by calling VecSet(). 85c4762a1bSJed Brown */ 86c4762a1bSJed Brown 87c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr); 91c4762a1bSJed Brown /* 92c4762a1bSJed Brown Some systems computes a residual that is identically zero, thus converging 93c4762a1bSJed Brown due to FNORM_ABS, others converge due to FNORM_RELATIVE. Here, we only 94c4762a1bSJed Brown report the reason if the iteration did not converge so that the tests are 95c4762a1bSJed Brown reproducible. 96c4762a1bSJed Brown */ 97c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%s number of SNES iterations = %D\n",reason>0 ? "CONVERGED" : SNESConvergedReasons[reason],its);CHKERRQ(ierr); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 101c4762a1bSJed Brown are no longer needed. 102c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 103c4762a1bSJed Brown 104c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = PetscFinalize(); 107c4762a1bSJed Brown return ierr; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 110c4762a1bSJed Brown /* 111c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 112c4762a1bSJed Brown 113c4762a1bSJed Brown Input Parameters: 114c4762a1bSJed Brown . snes - the SNES context 115c4762a1bSJed Brown . x - input vector 116c4762a1bSJed Brown . ctx - optional user-defined context 117c4762a1bSJed Brown 118c4762a1bSJed Brown Output Parameter: 119c4762a1bSJed Brown . f - function vector 120c4762a1bSJed Brown */ 121c4762a1bSJed Brown PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *ctx) 122c4762a1bSJed Brown { 123c4762a1bSJed Brown PetscErrorCode ierr; 124c4762a1bSJed Brown PetscScalar *ff; 125c4762a1bSJed Brown const PetscScalar *xx; 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* 128c4762a1bSJed Brown Get pointers to vector data. 129c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 130c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 131c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 132c4762a1bSJed Brown the array. 133c4762a1bSJed Brown */ 134c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Compute function */ 138c4762a1bSJed Brown ff[0] = -2.0 + 2.0*xx[0] + 400.0*xx[0]*xx[0]*xx[0] - 400.0*xx[0]*xx[1]; 139c4762a1bSJed Brown ff[1] = -200.0*xx[0]*xx[0] + 200.0*xx[1]; 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Restore vectors */ 142c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 144c4762a1bSJed Brown return 0; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 147c4762a1bSJed Brown /* 148c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 149c4762a1bSJed Brown 150c4762a1bSJed Brown Input Parameters: 151c4762a1bSJed Brown . snes - the SNES context 152c4762a1bSJed Brown . x - input vector 153c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 154c4762a1bSJed Brown 155c4762a1bSJed Brown Output Parameters: 156c4762a1bSJed Brown . jac - Jacobian matrix 157c4762a1bSJed Brown . B - optionally different preconditioning matrix 158c4762a1bSJed Brown . flag - flag indicating matrix structure 159c4762a1bSJed Brown */ 160c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 161c4762a1bSJed Brown { 162c4762a1bSJed Brown const PetscScalar *xx; 163c4762a1bSJed Brown PetscScalar A[4]; 164c4762a1bSJed Brown PetscErrorCode ierr; 165c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* 168c4762a1bSJed Brown Get pointer to vector data 169c4762a1bSJed Brown */ 170c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 174c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 175c4762a1bSJed Brown the matrix at once. 176c4762a1bSJed Brown */ 177c4762a1bSJed Brown A[0] = 2.0 + 1200.0*xx[0]*xx[0] - 400.0*xx[1]; 178c4762a1bSJed Brown A[1] = -400.0*xx[0]; 179c4762a1bSJed Brown A[2] = -400.0*xx[0]; 180c4762a1bSJed Brown A[3] = 200; 181c4762a1bSJed Brown ierr = MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr); 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* 184c4762a1bSJed Brown Restore vector 185c4762a1bSJed Brown */ 186c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* 189c4762a1bSJed Brown Assemble matrix 190c4762a1bSJed Brown */ 191c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 193c4762a1bSJed Brown if (jac != B) { 194c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown return 0; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown 200c4762a1bSJed Brown /*TEST 201c4762a1bSJed Brown 202c4762a1bSJed Brown test: 203*41ba4c6cSHeeho Park suffix: 1 204c4762a1bSJed Brown args: -snes_monitor_short -snes_max_it 1000 205c4762a1bSJed Brown requires: !single 206c4762a1bSJed Brown 207*41ba4c6cSHeeho Park test: 208*41ba4c6cSHeeho Park suffix: 2 209*41ba4c6cSHeeho Park args: -snes_monitor_short -snes_max_it 1000 -snes_type newtontrdc -snes_trdc_use_cauchy false 210*41ba4c6cSHeeho Park requires: !single 211*41ba4c6cSHeeho Park 212c4762a1bSJed Brown TEST*/ 213