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