1c4762a1bSJed Brown 2c4762a1bSJed Brown static const char help[] = "Tries to solve u`` + u^{2} = f for an easy case and an impossible case.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown This example was contributed by Peter Graf to show how SNES fails when given a nonlinear problem with no solution. 6c4762a1bSJed Brown 7c4762a1bSJed Brown Run with -n 14 to see fail to converge and -n 15 to see convergence 8c4762a1bSJed Brown 9c4762a1bSJed Brown The option -second_order uses a different discretization of the Neumann boundary condition and always converges 10c4762a1bSJed Brown 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscsnes.h> 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscBool second_order = PETSC_FALSE; 16c4762a1bSJed Brown #define X0DOT -2.0 17c4762a1bSJed Brown #define X1 5.0 18c4762a1bSJed Brown #define KPOW 2.0 19c4762a1bSJed Brown const PetscScalar sperturb = 1.1; 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* 22c4762a1bSJed Brown User-defined routines 23c4762a1bSJed Brown */ 24c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 25c4762a1bSJed Brown PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 26c4762a1bSJed Brown 27c4762a1bSJed Brown int main(int argc,char **argv) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown SNES snes; /* SNES context */ 30c4762a1bSJed Brown Vec x,r,F; /* vectors */ 31c4762a1bSJed Brown Mat J; /* Jacobian */ 32c4762a1bSJed Brown PetscErrorCode ierr; 33c4762a1bSJed Brown PetscInt it,n = 11,i; 34c4762a1bSJed Brown PetscReal h,xp = 0.0; 35c4762a1bSJed Brown PetscScalar v; 36c4762a1bSJed Brown const PetscScalar a = X0DOT; 37c4762a1bSJed Brown const PetscScalar b = X1; 38c4762a1bSJed Brown const PetscScalar k = KPOW; 39c4762a1bSJed Brown PetscScalar v2; 40c4762a1bSJed Brown PetscScalar *xx; 41c4762a1bSJed Brown 42c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-second_order",&second_order,NULL)); 45c4762a1bSJed Brown h = 1.0/(n-1); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 48c4762a1bSJed Brown Create nonlinear solver context 49c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 50c4762a1bSJed Brown 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 54c4762a1bSJed Brown Create vector data structures; set function evaluation routine 55c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 56c4762a1bSJed Brown 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&x)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,PETSC_DECIDE,n)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&F)); 62c4762a1bSJed Brown 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction,(void*)F)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 66c4762a1bSJed Brown Create matrix data structures; set Jacobian evaluation routine 67c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 68c4762a1bSJed Brown 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&J)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown Note that in this case we create separate matrices for the Jacobian 73c4762a1bSJed Brown and preconditioner matrix. Both of these are computed in the 74c4762a1bSJed Brown routine FormJacobian() 75c4762a1bSJed Brown */ 76*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(SNESSetJacobian(snes,NULL,JPrec,FormJacobian,0)); */ 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,0)); 78*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(SNESSetJacobian(snes,J,JPrec,FormJacobian,0)); */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 81c4762a1bSJed Brown Customize nonlinear solver; set runtime options 82c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 83c4762a1bSJed Brown 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 87c4762a1bSJed Brown Initialize application: 88c4762a1bSJed Brown Store right-hand-side of PDE and exact solution 89c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* set right hand side and initial guess to be exact solution of continuim problem */ 92c4762a1bSJed Brown #define SQR(x) ((x)*(x)) 93c4762a1bSJed Brown xp = 0.0; 94c4762a1bSJed Brown for (i=0; i<n; i++) 95c4762a1bSJed Brown { 96c4762a1bSJed Brown v = k*(k-1.)*(b-a)*PetscPowScalar(xp,k-2.) + SQR(a*xp) + SQR(b-a)*PetscPowScalar(xp,2.*k) + 2.*a*(b-a)*PetscPowScalar(xp,k+1.); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(F,1,&i,&v,INSERT_VALUES)); 98c4762a1bSJed Brown v2 = a*xp + (b-a)*PetscPowScalar(xp,k); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x,1,&i,&v2,INSERT_VALUES)); 100c4762a1bSJed Brown xp += h; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* perturb initial guess */ 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xx)); 105c4762a1bSJed Brown for (i=0; i<n; i++) { 106c4762a1bSJed Brown v2 = xx[i]*sperturb; 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x,1,&i,&v2,INSERT_VALUES)); 108c4762a1bSJed Brown } 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&xx)); 110c4762a1bSJed Brown 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&it)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"SNES iterations = %D\n\n",it)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 117c4762a1bSJed Brown are no longer needed. 118c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 119c4762a1bSJed Brown 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&F)); CHKERRQ(MatDestroy(&J)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 123c4762a1bSJed Brown ierr = PetscFinalize(); 124c4762a1bSJed Brown return ierr; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *dummy) 128c4762a1bSJed Brown { 129c4762a1bSJed Brown const PetscScalar *xx; 130c4762a1bSJed Brown PetscScalar *ff,*FF,d,d2; 131c4762a1bSJed Brown PetscInt i,n; 132c4762a1bSJed Brown 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray((Vec)dummy,&FF)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 137c4762a1bSJed Brown d = (PetscReal)(n - 1); d2 = d*d; 138c4762a1bSJed Brown 139c4762a1bSJed Brown if (second_order) ff[0] = d*(0.5*d*(-xx[2] + 4.*xx[1] - 3.*xx[0]) - X0DOT); 140c4762a1bSJed Brown else ff[0] = d*(d*(xx[1] - xx[0]) - X0DOT); 141c4762a1bSJed Brown 142c4762a1bSJed Brown for (i=1; i<n-1; i++) ff[i] = d2*(xx[i-1] - 2.*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i]; 143c4762a1bSJed Brown 144c4762a1bSJed Brown ff[n-1] = d*d*(xx[n-1] - X1); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray((Vec)dummy,&FF)); 148c4762a1bSJed Brown return 0; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat prejac,void *dummy) 152c4762a1bSJed Brown { 153c4762a1bSJed Brown const PetscScalar *xx; 154c4762a1bSJed Brown PetscScalar A[3],d,d2; 155c4762a1bSJed Brown PetscInt i,n,j[3]; 156c4762a1bSJed Brown 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 159c4762a1bSJed Brown d = (PetscReal)(n - 1); d2 = d*d; 160c4762a1bSJed Brown 161c4762a1bSJed Brown i = 0; 162c4762a1bSJed Brown if (second_order) { 163c4762a1bSJed Brown j[0] = 0; j[1] = 1; j[2] = 2; 164c4762a1bSJed Brown A[0] = -3.*d*d*0.5; A[1] = 4.*d*d*0.5; A[2] = -1.*d*d*0.5; 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES)); 166c4762a1bSJed Brown } else { 167c4762a1bSJed Brown j[0] = 0; j[1] = 1; 168c4762a1bSJed Brown A[0] = -d*d; A[1] = d*d; 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,2,j,A,INSERT_VALUES)); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown for (i=1; i<n-1; i++) { 172c4762a1bSJed Brown j[0] = i - 1; j[1] = i; j[2] = i + 1; 173c4762a1bSJed Brown A[0] = d2; A[1] = -2.*d2 + 2.*xx[i]; A[2] = d2; 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES)); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown i = n-1; 178c4762a1bSJed Brown A[0] = d*d; 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,1,&i,&A[0],INSERT_VALUES)); 180c4762a1bSJed Brown 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 182*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(prejac,MAT_FINAL_ASSEMBLY)); 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(prejac,MAT_FINAL_ASSEMBLY)); 185c4762a1bSJed Brown 186*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 187c4762a1bSJed Brown return 0; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown /*TEST 191c4762a1bSJed Brown 192c4762a1bSJed Brown test: 193c4762a1bSJed Brown args: -n 14 -snes_monitor_short -snes_converged_reason 194c4762a1bSJed Brown requires: !single 195c4762a1bSJed Brown 196c4762a1bSJed Brown test: 197c4762a1bSJed Brown suffix: 2 198c4762a1bSJed Brown args: -n 15 -snes_monitor_short -snes_converged_reason 199c4762a1bSJed Brown requires: !single 200c4762a1bSJed Brown 201c4762a1bSJed Brown test: 202c4762a1bSJed Brown suffix: 3 203c4762a1bSJed Brown args: -n 14 -second_order -snes_monitor_short -snes_converged_reason 204c4762a1bSJed Brown requires: !single 205c4762a1bSJed Brown 206c4762a1bSJed Brown TEST*/ 207