1c4762a1bSJed Brown static const char help[] = "Tries to solve u`` + u^{2} = f for an easy case and an impossible case.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown This example was contributed by Peter Graf to show how SNES fails when given a nonlinear problem with no solution. 5c4762a1bSJed Brown 6c4762a1bSJed Brown Run with -n 14 to see fail to converge and -n 15 to see convergence 7c4762a1bSJed Brown 8c4762a1bSJed Brown The option -second_order uses a different discretization of the Neumann boundary condition and always converges 9c4762a1bSJed Brown 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petscsnes.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscBool second_order = PETSC_FALSE; 15c4762a1bSJed Brown #define X0DOT -2.0 16c4762a1bSJed Brown #define X1 5.0 17c4762a1bSJed Brown #define KPOW 2.0 18c4762a1bSJed Brown const PetscScalar sperturb = 1.1; 19c4762a1bSJed Brown 20c4762a1bSJed Brown /* 21c4762a1bSJed Brown User-defined routines 22c4762a1bSJed Brown */ 23c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 24c4762a1bSJed Brown PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 25c4762a1bSJed Brown 26d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 27d71ae5a4SJacob Faibussowitsch { 28c4762a1bSJed Brown SNES snes; /* SNES context */ 29c4762a1bSJed Brown Vec x, r, F; /* vectors */ 30c4762a1bSJed Brown Mat J; /* Jacobian */ 31c4762a1bSJed Brown PetscInt it, n = 11, i; 32c4762a1bSJed Brown PetscReal h, xp = 0.0; 33c4762a1bSJed Brown PetscScalar v; 34c4762a1bSJed Brown const PetscScalar a = X0DOT; 35c4762a1bSJed Brown const PetscScalar b = X1; 36c4762a1bSJed Brown const PetscScalar k = KPOW; 37c4762a1bSJed Brown PetscScalar v2; 38c4762a1bSJed Brown PetscScalar *xx; 39c4762a1bSJed Brown 40327415f7SBarry Smith PetscFunctionBeginUser; 41c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 439566063dSJacob Faibussowitsch PetscCall(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 509566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 53c4762a1bSJed Brown Create vector data structures; set function evaluation routine 54c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &x)); 579566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, n)); 589566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &F)); 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65c4762a1bSJed Brown Create matrix data structures; set Jacobian evaluation routine 66c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(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 72*7addb90fSBarry Smith and matrix used to compute the preconditioner. Both of these are computed in the 73c4762a1bSJed Brown routine FormJacobian() 74c4762a1bSJed Brown */ 759566063dSJacob Faibussowitsch /* PetscCall(SNESSetJacobian(snes,NULL,JPrec,FormJacobian,0)); */ 769566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, 0)); 779566063dSJacob Faibussowitsch /* PetscCall(SNESSetJacobian(snes,J,JPrec,FormJacobian,0)); */ 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 80c4762a1bSJed Brown Customize nonlinear solver; set runtime options 81c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 86c4762a1bSJed Brown Initialize application: 87dd8e379bSPierre Jolivet Store right-hand side of PDE and exact solution 88c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 89c4762a1bSJed Brown 90bfe80ac4SPierre Jolivet /* set right-hand side and initial guess to be exact solution of continuum problem */ 91c4762a1bSJed Brown #define SQR(x) ((x) * (x)) 92c4762a1bSJed Brown xp = 0.0; 939371c9d4SSatish Balay for (i = 0; i < n; i++) { 94c4762a1bSJed 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.); 959566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES)); 96c4762a1bSJed Brown v2 = a * xp + (b - a) * PetscPowScalar(xp, k); 979566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v2, INSERT_VALUES)); 98c4762a1bSJed Brown xp += h; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* perturb initial guess */ 1029566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xx)); 103c4762a1bSJed Brown for (i = 0; i < n; i++) { 104c4762a1bSJed Brown v2 = xx[i] * sperturb; 1059566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v2, INSERT_VALUES)); 106c4762a1bSJed Brown } 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xx)); 108c4762a1bSJed Brown 1099566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 1109566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &it)); 11163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "SNES iterations = %" PetscInt_FMT "\n\n", it)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 115c4762a1bSJed Brown are no longer needed. 116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 117c4762a1bSJed Brown 1189371c9d4SSatish Balay PetscCall(VecDestroy(&x)); 1199371c9d4SSatish Balay PetscCall(VecDestroy(&r)); 1209371c9d4SSatish Balay PetscCall(VecDestroy(&F)); 1219371c9d4SSatish Balay PetscCall(MatDestroy(&J)); 1229566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1239566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 124b122ec5aSJacob Faibussowitsch return 0; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *dummy) 128d71ae5a4SJacob Faibussowitsch { 129c4762a1bSJed Brown const PetscScalar *xx; 130c4762a1bSJed Brown PetscScalar *ff, *FF, d, d2; 131c4762a1bSJed Brown PetscInt i, n; 132c4762a1bSJed Brown 1333ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1359566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 1369566063dSJacob Faibussowitsch PetscCall(VecGetArray((Vec)dummy, &FF)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1389371c9d4SSatish Balay d = (PetscReal)(n - 1); 1399371c9d4SSatish Balay d2 = d * d; 140c4762a1bSJed Brown 141c4762a1bSJed Brown if (second_order) ff[0] = d * (0.5 * d * (-xx[2] + 4. * xx[1] - 3. * xx[0]) - X0DOT); 142c4762a1bSJed Brown else ff[0] = d * (d * (xx[1] - xx[0]) - X0DOT); 143c4762a1bSJed Brown 144c4762a1bSJed 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]; 145c4762a1bSJed Brown 146c4762a1bSJed Brown ff[n - 1] = d * d * (xx[n - 1] - X1); 1479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 1499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray((Vec)dummy, &FF)); 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat prejac, void *dummy) 154d71ae5a4SJacob Faibussowitsch { 155c4762a1bSJed Brown const PetscScalar *xx; 156c4762a1bSJed Brown PetscScalar A[3], d, d2; 157c4762a1bSJed Brown PetscInt i, n, j[3]; 158c4762a1bSJed Brown 1593ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1609566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1619566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1629371c9d4SSatish Balay d = (PetscReal)(n - 1); 1639371c9d4SSatish Balay d2 = d * d; 164c4762a1bSJed Brown 165c4762a1bSJed Brown i = 0; 166c4762a1bSJed Brown if (second_order) { 1679371c9d4SSatish Balay j[0] = 0; 1689371c9d4SSatish Balay j[1] = 1; 1699371c9d4SSatish Balay j[2] = 2; 1709371c9d4SSatish Balay A[0] = -3. * d * d * 0.5; 1719371c9d4SSatish Balay A[1] = 4. * d * d * 0.5; 1729371c9d4SSatish Balay A[2] = -1. * d * d * 0.5; 1739566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 3, j, A, INSERT_VALUES)); 174c4762a1bSJed Brown } else { 1759371c9d4SSatish Balay j[0] = 0; 1769371c9d4SSatish Balay j[1] = 1; 1779371c9d4SSatish Balay A[0] = -d * d; 1789371c9d4SSatish Balay A[1] = d * d; 1799566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 2, j, A, INSERT_VALUES)); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 1829371c9d4SSatish Balay j[0] = i - 1; 1839371c9d4SSatish Balay j[1] = i; 1849371c9d4SSatish Balay j[2] = i + 1; 1859371c9d4SSatish Balay A[0] = d2; 1869371c9d4SSatish Balay A[1] = -2. * d2 + 2. * xx[i]; 1879371c9d4SSatish Balay A[2] = d2; 1889566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 3, j, A, INSERT_VALUES)); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown i = n - 1; 192c4762a1bSJed Brown A[0] = d * d; 1939566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 1, &i, &A[0], INSERT_VALUES)); 194c4762a1bSJed Brown 1959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 1969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 1979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(prejac, MAT_FINAL_ASSEMBLY)); 1989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(prejac, MAT_FINAL_ASSEMBLY)); 199c4762a1bSJed Brown 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204c4762a1bSJed Brown /*TEST 205c4762a1bSJed Brown 206c4762a1bSJed Brown test: 207c4762a1bSJed Brown args: -n 14 -snes_monitor_short -snes_converged_reason 208c4762a1bSJed Brown requires: !single 209c4762a1bSJed Brown 210c4762a1bSJed Brown test: 211c4762a1bSJed Brown suffix: 2 212c4762a1bSJed Brown args: -n 15 -snes_monitor_short -snes_converged_reason 213c4762a1bSJed Brown requires: !single 214c4762a1bSJed Brown 215c4762a1bSJed Brown test: 216c4762a1bSJed Brown suffix: 3 217c4762a1bSJed Brown args: -n 14 -second_order -snes_monitor_short -snes_converged_reason 218c4762a1bSJed Brown requires: !single 219c4762a1bSJed Brown 220c4762a1bSJed Brown TEST*/ 221