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 27d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 28d71ae5a4SJacob Faibussowitsch { 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 41327415f7SBarry Smith PetscFunctionBeginUser; 429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 449566063dSJacob Faibussowitsch PetscCall(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 519566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 54c4762a1bSJed Brown Create vector data structures; set function evaluation routine 55c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &x)); 589566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, PETSC_DECIDE, n)); 599566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &F)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)F)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 66c4762a1bSJed Brown Create matrix data structures; set Jacobian evaluation routine 67c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(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 */ 769566063dSJacob Faibussowitsch /* PetscCall(SNESSetJacobian(snes,NULL,JPrec,FormJacobian,0)); */ 779566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, 0)); 789566063dSJacob Faibussowitsch /* PetscCall(SNESSetJacobian(snes,J,JPrec,FormJacobian,0)); */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 81c4762a1bSJed Brown Customize nonlinear solver; set runtime options 82c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 83c4762a1bSJed Brown 849566063dSJacob Faibussowitsch PetscCall(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; 949371c9d4SSatish Balay for (i = 0; i < n; i++) { 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.); 969566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES)); 97c4762a1bSJed Brown v2 = a * xp + (b - a) * PetscPowScalar(xp, k); 989566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v2, INSERT_VALUES)); 99c4762a1bSJed Brown xp += h; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* perturb initial guess */ 1039566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xx)); 104c4762a1bSJed Brown for (i = 0; i < n; i++) { 105c4762a1bSJed Brown v2 = xx[i] * sperturb; 1069566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &i, &v2, INSERT_VALUES)); 107c4762a1bSJed Brown } 1089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xx)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 1119566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &it)); 11263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "SNES iterations = %" PetscInt_FMT "\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 1199371c9d4SSatish Balay PetscCall(VecDestroy(&x)); 1209371c9d4SSatish Balay PetscCall(VecDestroy(&r)); 1219371c9d4SSatish Balay PetscCall(VecDestroy(&F)); 1229371c9d4SSatish Balay PetscCall(MatDestroy(&J)); 1239566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1249566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 125b122ec5aSJacob Faibussowitsch return 0; 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *dummy) 129d71ae5a4SJacob Faibussowitsch { 130c4762a1bSJed Brown const PetscScalar *xx; 131c4762a1bSJed Brown PetscScalar *ff, *FF, d, d2; 132c4762a1bSJed Brown PetscInt i, n; 133c4762a1bSJed Brown 134*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1369566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArray((Vec)dummy, &FF)); 1389566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1399371c9d4SSatish Balay d = (PetscReal)(n - 1); 1409371c9d4SSatish Balay d2 = d * d; 141c4762a1bSJed Brown 142c4762a1bSJed Brown if (second_order) ff[0] = d * (0.5 * d * (-xx[2] + 4. * xx[1] - 3. * xx[0]) - X0DOT); 143c4762a1bSJed Brown else ff[0] = d * (d * (xx[1] - xx[0]) - X0DOT); 144c4762a1bSJed Brown 145c4762a1bSJed 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]; 146c4762a1bSJed Brown 147c4762a1bSJed Brown ff[n - 1] = d * d * (xx[n - 1] - X1); 1489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 1509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray((Vec)dummy, &FF)); 151*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat prejac, void *dummy) 155d71ae5a4SJacob Faibussowitsch { 156c4762a1bSJed Brown const PetscScalar *xx; 157c4762a1bSJed Brown PetscScalar A[3], d, d2; 158c4762a1bSJed Brown PetscInt i, n, j[3]; 159c4762a1bSJed Brown 160*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1619566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1639371c9d4SSatish Balay d = (PetscReal)(n - 1); 1649371c9d4SSatish Balay d2 = d * d; 165c4762a1bSJed Brown 166c4762a1bSJed Brown i = 0; 167c4762a1bSJed Brown if (second_order) { 1689371c9d4SSatish Balay j[0] = 0; 1699371c9d4SSatish Balay j[1] = 1; 1709371c9d4SSatish Balay j[2] = 2; 1719371c9d4SSatish Balay A[0] = -3. * d * d * 0.5; 1729371c9d4SSatish Balay A[1] = 4. * d * d * 0.5; 1739371c9d4SSatish Balay A[2] = -1. * d * d * 0.5; 1749566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 3, j, A, INSERT_VALUES)); 175c4762a1bSJed Brown } else { 1769371c9d4SSatish Balay j[0] = 0; 1779371c9d4SSatish Balay j[1] = 1; 1789371c9d4SSatish Balay A[0] = -d * d; 1799371c9d4SSatish Balay A[1] = d * d; 1809566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 2, j, A, INSERT_VALUES)); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 1839371c9d4SSatish Balay j[0] = i - 1; 1849371c9d4SSatish Balay j[1] = i; 1859371c9d4SSatish Balay j[2] = i + 1; 1869371c9d4SSatish Balay A[0] = d2; 1879371c9d4SSatish Balay A[1] = -2. * d2 + 2. * xx[i]; 1889371c9d4SSatish Balay A[2] = d2; 1899566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 3, j, A, INSERT_VALUES)); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192c4762a1bSJed Brown i = n - 1; 193c4762a1bSJed Brown A[0] = d * d; 1949566063dSJacob Faibussowitsch PetscCall(MatSetValues(prejac, 1, &i, 1, &i, &A[0], INSERT_VALUES)); 195c4762a1bSJed Brown 1969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 1979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 1989566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(prejac, MAT_FINAL_ASSEMBLY)); 1999566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(prejac, MAT_FINAL_ASSEMBLY)); 200c4762a1bSJed Brown 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 202*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /*TEST 206c4762a1bSJed Brown 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown args: -n 14 -snes_monitor_short -snes_converged_reason 209c4762a1bSJed Brown requires: !single 210c4762a1bSJed Brown 211c4762a1bSJed Brown test: 212c4762a1bSJed Brown suffix: 2 213c4762a1bSJed Brown args: -n 15 -snes_monitor_short -snes_converged_reason 214c4762a1bSJed Brown requires: !single 215c4762a1bSJed Brown 216c4762a1bSJed Brown test: 217c4762a1bSJed Brown suffix: 3 218c4762a1bSJed Brown args: -n 14 -second_order -snes_monitor_short -snes_converged_reason 219c4762a1bSJed Brown requires: !single 220c4762a1bSJed Brown 221c4762a1bSJed Brown TEST*/ 222