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 PetscInt it,n = 11,i; 33c4762a1bSJed Brown PetscReal h,xp = 0.0; 34c4762a1bSJed Brown PetscScalar v; 35c4762a1bSJed Brown const PetscScalar a = X0DOT; 36c4762a1bSJed Brown const PetscScalar b = X1; 37c4762a1bSJed Brown const PetscScalar k = KPOW; 38c4762a1bSJed Brown PetscScalar v2; 39c4762a1bSJed Brown PetscScalar *xx; 40c4762a1bSJed Brown 41*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-second_order",&second_order,NULL)); 44c4762a1bSJed Brown h = 1.0/(n-1); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 47c4762a1bSJed Brown Create nonlinear solver context 48c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 49c4762a1bSJed Brown 505f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 53c4762a1bSJed Brown Create vector data structures; set function evaluation routine 54c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 55c4762a1bSJed Brown 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&x)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,PETSC_DECIDE,n)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&F)); 61c4762a1bSJed Brown 625f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction,(void*)F)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65c4762a1bSJed Brown Create matrix data structures; set Jacobian evaluation routine 66c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&J)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* 71c4762a1bSJed Brown Note that in this case we create separate matrices for the Jacobian 72c4762a1bSJed Brown and preconditioner matrix. Both of these are computed in the 73c4762a1bSJed Brown routine FormJacobian() 74c4762a1bSJed Brown */ 755f80ce2aSJacob Faibussowitsch /* CHKERRQ(SNESSetJacobian(snes,NULL,JPrec,FormJacobian,0)); */ 765f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,0)); 775f80ce2aSJacob Faibussowitsch /* CHKERRQ(SNESSetJacobian(snes,J,JPrec,FormJacobian,0)); */ 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 80c4762a1bSJed Brown Customize nonlinear solver; set runtime options 81c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 82c4762a1bSJed Brown 835f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 86c4762a1bSJed Brown Initialize application: 87c4762a1bSJed Brown Store right-hand-side of PDE and exact solution 88c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* set right hand side and initial guess to be exact solution of continuim problem */ 91c4762a1bSJed Brown #define SQR(x) ((x)*(x)) 92c4762a1bSJed Brown xp = 0.0; 93c4762a1bSJed Brown for (i=0; i<n; i++) 94c4762a1bSJed Brown { 95c4762a1bSJed 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.); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(F,1,&i,&v,INSERT_VALUES)); 97c4762a1bSJed Brown v2 = a*xp + (b-a)*PetscPowScalar(xp,k); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x,1,&i,&v2,INSERT_VALUES)); 99c4762a1bSJed Brown xp += h; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* perturb initial guess */ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xx)); 104c4762a1bSJed Brown for (i=0; i<n; i++) { 105c4762a1bSJed Brown v2 = xx[i]*sperturb; 1065f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x,1,&i,&v2,INSERT_VALUES)); 107c4762a1bSJed Brown } 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&xx)); 109c4762a1bSJed Brown 1105f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&it)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"SNES iterations = %D\n\n",it)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 115c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 116c4762a1bSJed Brown are no longer needed. 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 118c4762a1bSJed Brown 1195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&F)); CHKERRQ(MatDestroy(&J)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 122*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 123*b122ec5aSJacob Faibussowitsch return 0; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *dummy) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown const PetscScalar *xx; 129c4762a1bSJed Brown PetscScalar *ff,*FF,d,d2; 130c4762a1bSJed Brown PetscInt i,n; 131c4762a1bSJed Brown 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray((Vec)dummy,&FF)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 136c4762a1bSJed Brown d = (PetscReal)(n - 1); d2 = d*d; 137c4762a1bSJed Brown 138c4762a1bSJed Brown if (second_order) ff[0] = d*(0.5*d*(-xx[2] + 4.*xx[1] - 3.*xx[0]) - X0DOT); 139c4762a1bSJed Brown else ff[0] = d*(d*(xx[1] - xx[0]) - X0DOT); 140c4762a1bSJed Brown 141c4762a1bSJed 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]; 142c4762a1bSJed Brown 143c4762a1bSJed Brown ff[n-1] = d*d*(xx[n-1] - X1); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray((Vec)dummy,&FF)); 147c4762a1bSJed Brown return 0; 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat prejac,void *dummy) 151c4762a1bSJed Brown { 152c4762a1bSJed Brown const PetscScalar *xx; 153c4762a1bSJed Brown PetscScalar A[3],d,d2; 154c4762a1bSJed Brown PetscInt i,n,j[3]; 155c4762a1bSJed Brown 1565f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 158c4762a1bSJed Brown d = (PetscReal)(n - 1); d2 = d*d; 159c4762a1bSJed Brown 160c4762a1bSJed Brown i = 0; 161c4762a1bSJed Brown if (second_order) { 162c4762a1bSJed Brown j[0] = 0; j[1] = 1; j[2] = 2; 163c4762a1bSJed Brown A[0] = -3.*d*d*0.5; A[1] = 4.*d*d*0.5; A[2] = -1.*d*d*0.5; 1645f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES)); 165c4762a1bSJed Brown } else { 166c4762a1bSJed Brown j[0] = 0; j[1] = 1; 167c4762a1bSJed Brown A[0] = -d*d; A[1] = d*d; 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,2,j,A,INSERT_VALUES)); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown for (i=1; i<n-1; i++) { 171c4762a1bSJed Brown j[0] = i - 1; j[1] = i; j[2] = i + 1; 172c4762a1bSJed Brown A[0] = d2; A[1] = -2.*d2 + 2.*xx[i]; A[2] = d2; 1735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES)); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown i = n-1; 177c4762a1bSJed Brown A[0] = d*d; 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(prejac,1,&i,1,&i,&A[0],INSERT_VALUES)); 179c4762a1bSJed Brown 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(prejac,MAT_FINAL_ASSEMBLY)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(prejac,MAT_FINAL_ASSEMBLY)); 184c4762a1bSJed Brown 1855f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 186c4762a1bSJed Brown return 0; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /*TEST 190c4762a1bSJed Brown 191c4762a1bSJed Brown test: 192c4762a1bSJed Brown args: -n 14 -snes_monitor_short -snes_converged_reason 193c4762a1bSJed Brown requires: !single 194c4762a1bSJed Brown 195c4762a1bSJed Brown test: 196c4762a1bSJed Brown suffix: 2 197c4762a1bSJed Brown args: -n 15 -snes_monitor_short -snes_converged_reason 198c4762a1bSJed Brown requires: !single 199c4762a1bSJed Brown 200c4762a1bSJed Brown test: 201c4762a1bSJed Brown suffix: 3 202c4762a1bSJed Brown args: -n 14 -second_order -snes_monitor_short -snes_converged_reason 203c4762a1bSJed Brown requires: !single 204c4762a1bSJed Brown 205c4762a1bSJed Brown TEST*/ 206