1c4762a1bSJed Brown static char help[] = "Solves u`` + u^{2} = f with Newton-like methods. Using\n\ 2c4762a1bSJed Brown matrix-free techniques with user-provided explicit preconditioner matrix.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscsnes.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 7c4762a1bSJed Brown extern PetscErrorCode FormJacobianNoMatrix(SNES, Vec, Mat, Mat, void *); 8c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 9c4762a1bSJed Brown extern PetscErrorCode FormFunctioni(void *, PetscInt, Vec, PetscScalar *); 10c4762a1bSJed Brown extern PetscErrorCode OtherFunctionForDifferencing(void *, Vec, Vec); 11c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(SNES, Vec); 12c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *); 13c4762a1bSJed Brown 14c4762a1bSJed Brown typedef struct { 15c4762a1bSJed Brown PetscViewer viewer; 16c4762a1bSJed Brown } MonitorCtx; 17c4762a1bSJed Brown 18c4762a1bSJed Brown typedef struct { 19c4762a1bSJed Brown PetscBool variant; 20c4762a1bSJed Brown } AppCtx; 21c4762a1bSJed Brown 22d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 23d71ae5a4SJacob Faibussowitsch { 24c4762a1bSJed Brown SNES snes; /* SNES context */ 25c4762a1bSJed Brown SNESType type = SNESNEWTONLS; /* default nonlinear solution method */ 26c4762a1bSJed Brown Vec x, r, F, U; /* vectors */ 27c4762a1bSJed Brown Mat J, B; /* Jacobian matrix-free, explicit preconditioner */ 28c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 29c4762a1bSJed Brown PetscScalar h, xp = 0.0, v; 30c4762a1bSJed Brown PetscInt its, n = 5, i; 31c4762a1bSJed Brown PetscBool puremf = PETSC_FALSE; 32c4762a1bSJed Brown 33327415f7SBarry Smith PetscFunctionBeginUser; 349566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-variant", &user.variant)); 37c4762a1bSJed Brown h = 1.0 / (n - 1); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Set up data structures */ 409566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x)); 419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution")); 429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &F)); 449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &U)); 459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution")); 46c4762a1bSJed Brown 47a5b23f4aSJose E. Roman /* create explicit matrix preconditioner */ 489566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n, n, 3, NULL, &B)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Store right-hand-side of PDE and exact solution */ 51c4762a1bSJed Brown for (i = 0; i < n; i++) { 52c4762a1bSJed Brown v = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */ 539566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES)); 54c4762a1bSJed Brown v = xp * xp * xp; 559566063dSJacob Faibussowitsch PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES)); 56c4762a1bSJed Brown xp += h; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Create nonlinear solver */ 609566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 619566063dSJacob Faibussowitsch PetscCall(SNESSetType(snes, type)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Set various routines and options */ 649566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, F)); 65c4762a1bSJed Brown if (user.variant) { 66c4762a1bSJed Brown /* this approach is not normally needed, one should use the MatCreateSNESMF() below usually */ 679566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PETSC_COMM_WORLD, n, n, n, n, &J)); 689566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(J, (PetscErrorCode(*)(void *, Vec, Vec))SNESComputeFunction, snes)); 699566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunctioni(J, FormFunctioni)); 7001c1178eSBarry Smith /* Use the matrix-free operator for both the Jacobian used to define the linear system and used to define the preconditioner */ 71c4762a1bSJed Brown /* This tests MatGetDiagonal() for MATMFFD */ 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-puremf", &puremf)); 73c4762a1bSJed Brown } else { 7401c1178eSBarry Smith /* create matrix-free matrix for Jacobian */ 759566063dSJacob Faibussowitsch PetscCall(MatCreateSNESMF(snes, &J)); 76c4762a1bSJed Brown /* demonstrates differencing a different function than FormFunction() to apply a matrix operator */ 77c4762a1bSJed Brown /* note we use the same context for this function as FormFunction, the F vector */ 789566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(J, OtherFunctionForDifferencing, F)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Set various routines and options */ 829566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, puremf ? J : B, puremf ? FormJacobianNoMatrix : FormJacobian, &user)); 839566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Solve nonlinear system */ 869566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(snes, x)); 879566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 889566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 8963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "number of SNES iterations = %" PetscInt_FMT "\n\n", its)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Free data structures */ 929371c9d4SSatish Balay PetscCall(VecDestroy(&x)); 939371c9d4SSatish Balay PetscCall(VecDestroy(&r)); 949371c9d4SSatish Balay PetscCall(VecDestroy(&U)); 959371c9d4SSatish Balay PetscCall(VecDestroy(&F)); 969371c9d4SSatish Balay PetscCall(MatDestroy(&J)); 979371c9d4SSatish Balay PetscCall(MatDestroy(&B)); 989566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 999566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 100b122ec5aSJacob Faibussowitsch return 0; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown /* -------------------- Evaluate Function F(x) --------------------- */ 103c4762a1bSJed Brown 104d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *dummy) 105d71ae5a4SJacob Faibussowitsch { 106c4762a1bSJed Brown const PetscScalar *xx, *FF; 107c4762a1bSJed Brown PetscScalar *ff, d; 108c4762a1bSJed Brown PetscInt i, n; 109c4762a1bSJed Brown 1103ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 1139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead((Vec)dummy, &FF)); 1149566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1159371c9d4SSatish Balay d = (PetscReal)(n - 1); 1169371c9d4SSatish Balay d = d * d; 117c4762a1bSJed Brown ff[0] = xx[0]; 118c4762a1bSJed Brown for (i = 1; i < n - 1; i++) ff[i] = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - FF[i]; 119c4762a1bSJed Brown ff[n - 1] = xx[n - 1] - 1.0; 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead((Vec)dummy, &FF)); 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctioni(void *dummy, PetscInt i, Vec x, PetscScalar *s) 127d71ae5a4SJacob Faibussowitsch { 128c4762a1bSJed Brown const PetscScalar *xx, *FF; 129c4762a1bSJed Brown PetscScalar d; 130c4762a1bSJed Brown PetscInt n; 131c4762a1bSJed Brown SNES snes = (SNES)dummy; 132c4762a1bSJed Brown Vec F; 133c4762a1bSJed Brown 1343ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1359566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, NULL, NULL, (void **)&F)); 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &FF)); 1389566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1399371c9d4SSatish Balay d = (PetscReal)(n - 1); 1409371c9d4SSatish Balay d = d * d; 141c4762a1bSJed Brown if (i == 0) { 142c4762a1bSJed Brown *s = xx[0]; 143c4762a1bSJed Brown } else if (i == n - 1) { 144c4762a1bSJed Brown *s = xx[n - 1] - 1.0; 145c4762a1bSJed Brown } else { 146c4762a1bSJed Brown *s = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - FF[i]; 147c4762a1bSJed Brown } 1489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &FF)); 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* 154c4762a1bSJed Brown 15501c1178eSBarry Smith Example function that when differenced produces the same matrix-free Jacobian as FormFunction() 156c4762a1bSJed Brown this is provided to show how a user can provide a different function 157c4762a1bSJed Brown */ 158d71ae5a4SJacob Faibussowitsch PetscErrorCode OtherFunctionForDifferencing(void *dummy, Vec x, Vec f) 159d71ae5a4SJacob Faibussowitsch { 1603ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1619566063dSJacob Faibussowitsch PetscCall(FormFunction(NULL, x, f, dummy)); 1629566063dSJacob Faibussowitsch PetscCall(VecShift(f, 1.0)); 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* -------------------- Form initial approximation ----------------- */ 167c4762a1bSJed Brown 168d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(SNES snes, Vec x) 169d71ae5a4SJacob Faibussowitsch { 1703ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 171*4d86920dSPierre Jolivet PetscCall(VecSet(x, 0.5)); 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown /* -------------------- Evaluate Jacobian F'(x) -------------------- */ 175c4762a1bSJed Brown /* Evaluates a matrix that is used to precondition the matrix-free 176a5b23f4aSJose E. Roman jacobian. In this case, the explicit preconditioner matrix is 177c4762a1bSJed Brown also EXACTLY the Jacobian. In general, it would be some lower 178c4762a1bSJed Brown order, simplified apprioximation */ 179c4762a1bSJed Brown 180d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 181d71ae5a4SJacob Faibussowitsch { 182c4762a1bSJed Brown const PetscScalar *xx; 183c4762a1bSJed Brown PetscScalar A[3], d; 184c4762a1bSJed Brown PetscInt i, n, j[3]; 185c4762a1bSJed Brown AppCtx *user = (AppCtx *)dummy; 186c4762a1bSJed Brown 1873ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1899566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1909371c9d4SSatish Balay d = (PetscReal)(n - 1); 1919371c9d4SSatish Balay d = d * d; 192c4762a1bSJed Brown 1939371c9d4SSatish Balay i = 0; 1949371c9d4SSatish Balay A[0] = 1.0; 1959566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 1, &i, &A[0], INSERT_VALUES)); 196c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 1979371c9d4SSatish Balay j[0] = i - 1; 1989371c9d4SSatish Balay j[1] = i; 1999371c9d4SSatish Balay j[2] = i + 1; 2009371c9d4SSatish Balay A[0] = d; 2019371c9d4SSatish Balay A[1] = -2.0 * d + 2.0 * xx[i]; 2029371c9d4SSatish Balay A[2] = d; 2039566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES)); 204c4762a1bSJed Brown } 2059371c9d4SSatish Balay i = n - 1; 2069371c9d4SSatish Balay A[0] = 1.0; 2079566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 1, &i, &A[0], INSERT_VALUES)); 2089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 211c4762a1bSJed Brown 2121baa6e33SBarry Smith if (user->variant) PetscCall(MatMFFDSetBase(jac, x, NULL)); 2139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianNoMatrix(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 219d71ae5a4SJacob Faibussowitsch { 220c4762a1bSJed Brown AppCtx *user = (AppCtx *)dummy; 221c4762a1bSJed Brown 2223ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2231baa6e33SBarry Smith if (user->variant) PetscCall(MatMFFDSetBase(jac, x, NULL)); 2249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* -------------------- User-defined monitor ----------------------- */ 230c4762a1bSJed Brown 231d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, void *dummy) 232d71ae5a4SJacob Faibussowitsch { 233c4762a1bSJed Brown MonitorCtx *monP = (MonitorCtx *)dummy; 234c4762a1bSJed Brown Vec x; 235c4762a1bSJed Brown MPI_Comm comm; 236c4762a1bSJed Brown 2373ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2389566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 23963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, stdout, "iter = %" PetscInt_FMT ", SNES Function norm %g \n", its, (double)fnorm)); 2409566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &x)); 2419566063dSJacob Faibussowitsch PetscCall(VecView(x, monP->viewer)); 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown /*TEST 246c4762a1bSJed Brown 247c4762a1bSJed Brown test: 248c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 249c4762a1bSJed Brown 250c4762a1bSJed Brown test: 251c4762a1bSJed Brown suffix: 2 252c4762a1bSJed Brown args: -variant -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 253c4762a1bSJed Brown output_file: output/ex7_1.out 254c4762a1bSJed Brown 255c4762a1bSJed Brown # uses AIJ matrix to define diagonal matrix for Jacobian preconditioning 256c4762a1bSJed Brown test: 2579d5502f9SJunchao Zhang requires: !single 258c4762a1bSJed Brown suffix: 3 259c4762a1bSJed Brown args: -variant -pc_type jacobi -snes_view -ksp_monitor 260c4762a1bSJed Brown 261c4762a1bSJed Brown # uses MATMFFD matrix to define diagonal matrix for Jacobian preconditioning 262c4762a1bSJed Brown test: 263c4762a1bSJed Brown suffix: 4 264c4762a1bSJed Brown args: -variant -pc_type jacobi -puremf -snes_view -ksp_monitor 265c4762a1bSJed Brown 266c4762a1bSJed Brown TEST*/ 267