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 PetscInt its; 28c4762a1bSJed Brown PetscScalar *xx; 29c4762a1bSJed Brown SNESConvergedReason reason; 30c4762a1bSJed Brown 31*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 34c4762a1bSJed Brown Create nonlinear solver context 35c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 365f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 39c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 40c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown Create vectors for solution and nonlinear function 43c4762a1bSJed Brown */ 445f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,PETSC_DECIDE,2)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Create Jacobian matrix data structure 51c4762a1bSJed Brown */ 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(J)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(J)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* 58c4762a1bSJed Brown Set function evaluation routine and vector. 59c4762a1bSJed Brown */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction1,NULL)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* 63c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 64c4762a1bSJed Brown */ 655f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian1,NULL)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 68c4762a1bSJed Brown Customize nonlinear solver; set runtime options 69c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 705f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 73c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 74c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xx)); 76c4762a1bSJed Brown xx[0] = -1.2; xx[1] = 1.0; 775f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&xx)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* 80c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 81c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 82c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 83c4762a1bSJed Brown this vector to zero by calling VecSet(). 84c4762a1bSJed Brown */ 85c4762a1bSJed Brown 865f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&its)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetConvergedReason(snes,&reason)); 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown Some systems computes a residual that is identically zero, thus converging 92c4762a1bSJed Brown due to FNORM_ABS, others converge due to FNORM_RELATIVE. Here, we only 93c4762a1bSJed Brown report the reason if the iteration did not converge so that the tests are 94c4762a1bSJed Brown reproducible. 95c4762a1bSJed Brown */ 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s number of SNES iterations = %D\n",reason>0 ? "CONVERGED" : SNESConvergedReasons[reason],its)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 100c4762a1bSJed Brown are no longer needed. 101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 102c4762a1bSJed Brown 1035f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); CHKERRQ(SNESDestroy(&snes)); 105*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 106*b122ec5aSJacob Faibussowitsch return 0; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 109c4762a1bSJed Brown /* 110c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 111c4762a1bSJed Brown 112c4762a1bSJed Brown Input Parameters: 113c4762a1bSJed Brown . snes - the SNES context 114c4762a1bSJed Brown . x - input vector 115c4762a1bSJed Brown . ctx - optional user-defined context 116c4762a1bSJed Brown 117c4762a1bSJed Brown Output Parameter: 118c4762a1bSJed Brown . f - function vector 119c4762a1bSJed Brown */ 120c4762a1bSJed Brown PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *ctx) 121c4762a1bSJed Brown { 122c4762a1bSJed Brown PetscScalar *ff; 123c4762a1bSJed Brown const PetscScalar *xx; 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown Get pointers to vector data. 127c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 128c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 129c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 130c4762a1bSJed Brown the array. 131c4762a1bSJed Brown */ 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Compute function */ 136c4762a1bSJed Brown ff[0] = -2.0 + 2.0*xx[0] + 400.0*xx[0]*xx[0]*xx[0] - 400.0*xx[0]*xx[1]; 137c4762a1bSJed Brown ff[1] = -200.0*xx[0]*xx[0] + 200.0*xx[1]; 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Restore vectors */ 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 142c4762a1bSJed Brown return 0; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 145c4762a1bSJed Brown /* 146c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 147c4762a1bSJed Brown 148c4762a1bSJed Brown Input Parameters: 149c4762a1bSJed Brown . snes - the SNES context 150c4762a1bSJed Brown . x - input vector 151c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 152c4762a1bSJed Brown 153c4762a1bSJed Brown Output Parameters: 154c4762a1bSJed Brown . jac - Jacobian matrix 155c4762a1bSJed Brown . B - optionally different preconditioning matrix 156c4762a1bSJed Brown . flag - flag indicating matrix structure 157c4762a1bSJed Brown */ 158c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 159c4762a1bSJed Brown { 160c4762a1bSJed Brown const PetscScalar *xx; 161c4762a1bSJed Brown PetscScalar A[4]; 162c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* 165c4762a1bSJed Brown Get pointer to vector data 166c4762a1bSJed Brown */ 1675f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* 170c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 171c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 172c4762a1bSJed Brown the matrix at once. 173c4762a1bSJed Brown */ 174c4762a1bSJed Brown A[0] = 2.0 + 1200.0*xx[0]*xx[0] - 400.0*xx[1]; 175c4762a1bSJed Brown A[1] = -400.0*xx[0]; 176c4762a1bSJed Brown A[2] = -400.0*xx[0]; 177c4762a1bSJed Brown A[3] = 200; 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* 181c4762a1bSJed Brown Restore vector 182c4762a1bSJed Brown */ 1835f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown Assemble matrix 187c4762a1bSJed Brown */ 1885f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 190c4762a1bSJed Brown if (jac != B) { 1915f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown return 0; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown /*TEST 198c4762a1bSJed Brown 199c4762a1bSJed Brown test: 20041ba4c6cSHeeho Park suffix: 1 201c4762a1bSJed Brown args: -snes_monitor_short -snes_max_it 1000 202c4762a1bSJed Brown requires: !single 203c4762a1bSJed Brown 20441ba4c6cSHeeho Park test: 20541ba4c6cSHeeho Park suffix: 2 20641ba4c6cSHeeho Park args: -snes_monitor_short -snes_max_it 1000 -snes_type newtontrdc -snes_trdc_use_cauchy false 20741ba4c6cSHeeho Park requires: !single 20841ba4c6cSHeeho Park 209c4762a1bSJed Brown TEST*/ 210