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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 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 */ 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,PETSC_DECIDE,2)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* 51c4762a1bSJed Brown Create Jacobian matrix data structure 52c4762a1bSJed Brown */ 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(J)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(J)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* 59c4762a1bSJed Brown Set function evaluation routine and vector. 60c4762a1bSJed Brown */ 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction1,NULL)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 65c4762a1bSJed Brown */ 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian1,NULL)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 69c4762a1bSJed Brown Customize nonlinear solver; set runtime options 70c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 75c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xx)); 77c4762a1bSJed Brown xx[0] = -1.2; xx[1] = 1.0; 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&xx)); 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 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&its)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetConvergedReason(snes,&reason)); 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 */ 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s number of SNES iterations = %D\n",reason>0 ? "CONVERGED" : SNESConvergedReasons[reason],its)); 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 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); CHKERRQ(SNESDestroy(&snes)); 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 PetscScalar *ff; 124c4762a1bSJed Brown const PetscScalar *xx; 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* 127c4762a1bSJed Brown Get pointers to vector data. 128c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 129c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 130c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 131c4762a1bSJed Brown the array. 132c4762a1bSJed Brown */ 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* Compute function */ 137c4762a1bSJed Brown ff[0] = -2.0 + 2.0*xx[0] + 400.0*xx[0]*xx[0]*xx[0] - 400.0*xx[0]*xx[1]; 138c4762a1bSJed Brown ff[1] = -200.0*xx[0]*xx[0] + 200.0*xx[1]; 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Restore vectors */ 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 143c4762a1bSJed Brown return 0; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 146c4762a1bSJed Brown /* 147c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 148c4762a1bSJed Brown 149c4762a1bSJed Brown Input Parameters: 150c4762a1bSJed Brown . snes - the SNES context 151c4762a1bSJed Brown . x - input vector 152c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 153c4762a1bSJed Brown 154c4762a1bSJed Brown Output Parameters: 155c4762a1bSJed Brown . jac - Jacobian matrix 156c4762a1bSJed Brown . B - optionally different preconditioning matrix 157c4762a1bSJed Brown . flag - flag indicating matrix structure 158c4762a1bSJed Brown */ 159c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 160c4762a1bSJed Brown { 161c4762a1bSJed Brown const PetscScalar *xx; 162c4762a1bSJed Brown PetscScalar A[4]; 163c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* 166c4762a1bSJed Brown Get pointer to vector data 167c4762a1bSJed Brown */ 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* 171c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 172c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 173c4762a1bSJed Brown the matrix at once. 174c4762a1bSJed Brown */ 175c4762a1bSJed Brown A[0] = 2.0 + 1200.0*xx[0]*xx[0] - 400.0*xx[1]; 176c4762a1bSJed Brown A[1] = -400.0*xx[0]; 177c4762a1bSJed Brown A[2] = -400.0*xx[0]; 178c4762a1bSJed Brown A[3] = 200; 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* 182c4762a1bSJed Brown Restore vector 183c4762a1bSJed Brown */ 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown Assemble matrix 188c4762a1bSJed Brown */ 189*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 191c4762a1bSJed Brown if (jac != B) { 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown return 0; 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown /*TEST 199c4762a1bSJed Brown 200c4762a1bSJed Brown test: 20141ba4c6cSHeeho Park suffix: 1 202c4762a1bSJed Brown args: -snes_monitor_short -snes_max_it 1000 203c4762a1bSJed Brown requires: !single 204c4762a1bSJed Brown 20541ba4c6cSHeeho Park test: 20641ba4c6cSHeeho Park suffix: 2 20741ba4c6cSHeeho Park args: -snes_monitor_short -snes_max_it 1000 -snes_type newtontrdc -snes_trdc_use_cauchy false 20841ba4c6cSHeeho Park requires: !single 20941ba4c6cSHeeho Park 210c4762a1bSJed Brown TEST*/ 211