1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves u`` + u^{2} = f with Newton-like methods. Using\n\ 3c4762a1bSJed Brown matrix-free techniques with user-provided explicit preconditioner matrix.\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscsnes.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 8c4762a1bSJed Brown extern PetscErrorCode FormJacobianNoMatrix(SNES, Vec, Mat, Mat, void *); 9c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 10c4762a1bSJed Brown extern PetscErrorCode FormFunctioni(void *, PetscInt, Vec, PetscScalar *); 11c4762a1bSJed Brown extern PetscErrorCode OtherFunctionForDifferencing(void *, Vec, Vec); 12c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(SNES, Vec); 13c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *); 14c4762a1bSJed Brown 15c4762a1bSJed Brown typedef struct { 16c4762a1bSJed Brown PetscViewer viewer; 17c4762a1bSJed Brown } MonitorCtx; 18c4762a1bSJed Brown 19c4762a1bSJed Brown typedef struct { 20c4762a1bSJed Brown PetscBool variant; 21c4762a1bSJed Brown } AppCtx; 22c4762a1bSJed Brown 23*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 24*d71ae5a4SJacob Faibussowitsch { 25c4762a1bSJed Brown SNES snes; /* SNES context */ 26c4762a1bSJed Brown SNESType type = SNESNEWTONLS; /* default nonlinear solution method */ 27c4762a1bSJed Brown Vec x, r, F, U; /* vectors */ 28c4762a1bSJed Brown Mat J, B; /* Jacobian matrix-free, explicit preconditioner */ 29c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 30c4762a1bSJed Brown PetscScalar h, xp = 0.0, v; 31c4762a1bSJed Brown PetscInt its, n = 5, i; 32c4762a1bSJed Brown PetscBool puremf = PETSC_FALSE; 33c4762a1bSJed Brown 34327415f7SBarry Smith PetscFunctionBeginUser; 359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-variant", &user.variant)); 38c4762a1bSJed Brown h = 1.0 / (n - 1); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Set up data structures */ 419566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x)); 429566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution")); 439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &F)); 459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &U)); 469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution")); 47c4762a1bSJed Brown 48a5b23f4aSJose E. Roman /* create explicit matrix preconditioner */ 499566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n, n, 3, NULL, &B)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Store right-hand-side of PDE and exact solution */ 52c4762a1bSJed Brown for (i = 0; i < n; i++) { 53c4762a1bSJed Brown v = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */ 549566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &v, INSERT_VALUES)); 55c4762a1bSJed Brown v = xp * xp * xp; 569566063dSJacob Faibussowitsch PetscCall(VecSetValues(U, 1, &i, &v, INSERT_VALUES)); 57c4762a1bSJed Brown xp += h; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create nonlinear solver */ 619566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 629566063dSJacob Faibussowitsch PetscCall(SNESSetType(snes, type)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Set various routines and options */ 659566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction, F)); 66c4762a1bSJed Brown if (user.variant) { 67c4762a1bSJed Brown /* this approach is not normally needed, one should use the MatCreateSNESMF() below usually */ 689566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PETSC_COMM_WORLD, n, n, n, n, &J)); 699566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(J, (PetscErrorCode(*)(void *, Vec, Vec))SNESComputeFunction, snes)); 709566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunctioni(J, FormFunctioni)); 71c4762a1bSJed Brown /* Use the matrix free operator for both the Jacobian used to define the linear system and used to define the preconditioner */ 72c4762a1bSJed Brown /* This tests MatGetDiagonal() for MATMFFD */ 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-puremf", &puremf)); 74c4762a1bSJed Brown } else { 75c4762a1bSJed Brown /* create matrix free matrix for Jacobian */ 769566063dSJacob Faibussowitsch PetscCall(MatCreateSNESMF(snes, &J)); 77c4762a1bSJed Brown /* demonstrates differencing a different function than FormFunction() to apply a matrix operator */ 78c4762a1bSJed Brown /* note we use the same context for this function as FormFunction, the F vector */ 799566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(J, OtherFunctionForDifferencing, F)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Set various routines and options */ 839566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, puremf ? J : B, puremf ? FormJacobianNoMatrix : FormJacobian, &user)); 849566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Solve nonlinear system */ 879566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(snes, x)); 889566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 899566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 9063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "number of SNES iterations = %" PetscInt_FMT "\n\n", its)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Free data structures */ 939371c9d4SSatish Balay PetscCall(VecDestroy(&x)); 949371c9d4SSatish Balay PetscCall(VecDestroy(&r)); 959371c9d4SSatish Balay PetscCall(VecDestroy(&U)); 969371c9d4SSatish Balay PetscCall(VecDestroy(&F)); 979371c9d4SSatish Balay PetscCall(MatDestroy(&J)); 989371c9d4SSatish Balay PetscCall(MatDestroy(&B)); 999566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 101b122ec5aSJacob Faibussowitsch return 0; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown /* -------------------- Evaluate Function F(x) --------------------- */ 104c4762a1bSJed Brown 105*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *dummy) 106*d71ae5a4SJacob Faibussowitsch { 107c4762a1bSJed Brown const PetscScalar *xx, *FF; 108c4762a1bSJed Brown PetscScalar *ff, d; 109c4762a1bSJed Brown PetscInt i, n; 110c4762a1bSJed Brown 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)); 123c4762a1bSJed Brown return 0; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctioni(void *dummy, PetscInt i, Vec x, PetscScalar *s) 127*d71ae5a4SJacob 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 1349566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, NULL, NULL, (void **)&F)); 1359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F, &FF)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1389371c9d4SSatish Balay d = (PetscReal)(n - 1); 1399371c9d4SSatish Balay d = d * d; 140c4762a1bSJed Brown if (i == 0) { 141c4762a1bSJed Brown *s = xx[0]; 142c4762a1bSJed Brown } else if (i == n - 1) { 143c4762a1bSJed Brown *s = xx[n - 1] - 1.0; 144c4762a1bSJed Brown } else { 145c4762a1bSJed Brown *s = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - FF[i]; 146c4762a1bSJed Brown } 1479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F, &FF)); 149c4762a1bSJed Brown return 0; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* 153c4762a1bSJed Brown 154c4762a1bSJed Brown Example function that when differenced produces the same matrix free Jacobian as FormFunction() 155c4762a1bSJed Brown this is provided to show how a user can provide a different function 156c4762a1bSJed Brown */ 157*d71ae5a4SJacob Faibussowitsch PetscErrorCode OtherFunctionForDifferencing(void *dummy, Vec x, Vec f) 158*d71ae5a4SJacob Faibussowitsch { 1599566063dSJacob Faibussowitsch PetscCall(FormFunction(NULL, x, f, dummy)); 1609566063dSJacob Faibussowitsch PetscCall(VecShift(f, 1.0)); 161c4762a1bSJed Brown return 0; 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* -------------------- Form initial approximation ----------------- */ 165c4762a1bSJed Brown 166*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(SNES snes, Vec x) 167*d71ae5a4SJacob Faibussowitsch { 168c4762a1bSJed Brown PetscScalar pfive = .50; 1699566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 170c4762a1bSJed Brown return 0; 171c4762a1bSJed Brown } 172c4762a1bSJed Brown /* -------------------- Evaluate Jacobian F'(x) -------------------- */ 173c4762a1bSJed Brown /* Evaluates a matrix that is used to precondition the matrix-free 174a5b23f4aSJose E. Roman jacobian. In this case, the explicit preconditioner matrix is 175c4762a1bSJed Brown also EXACTLY the Jacobian. In general, it would be some lower 176c4762a1bSJed Brown order, simplified apprioximation */ 177c4762a1bSJed Brown 178*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 179*d71ae5a4SJacob Faibussowitsch { 180c4762a1bSJed Brown const PetscScalar *xx; 181c4762a1bSJed Brown PetscScalar A[3], d; 182c4762a1bSJed Brown PetscInt i, n, j[3]; 183c4762a1bSJed Brown AppCtx *user = (AppCtx *)dummy; 184c4762a1bSJed Brown 1859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1869566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &n)); 1879371c9d4SSatish Balay d = (PetscReal)(n - 1); 1889371c9d4SSatish Balay d = d * d; 189c4762a1bSJed Brown 1909371c9d4SSatish Balay i = 0; 1919371c9d4SSatish Balay A[0] = 1.0; 1929566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 1, &i, &A[0], INSERT_VALUES)); 193c4762a1bSJed Brown for (i = 1; i < n - 1; i++) { 1949371c9d4SSatish Balay j[0] = i - 1; 1959371c9d4SSatish Balay j[1] = i; 1969371c9d4SSatish Balay j[2] = i + 1; 1979371c9d4SSatish Balay A[0] = d; 1989371c9d4SSatish Balay A[1] = -2.0 * d + 2.0 * xx[i]; 1999371c9d4SSatish Balay A[2] = d; 2009566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES)); 201c4762a1bSJed Brown } 2029371c9d4SSatish Balay i = n - 1; 2039371c9d4SSatish Balay A[0] = 1.0; 2049566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &i, 1, &i, &A[0], INSERT_VALUES)); 2059566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2069566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 2079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 208c4762a1bSJed Brown 2091baa6e33SBarry Smith if (user->variant) PetscCall(MatMFFDSetBase(jac, x, NULL)); 2109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 212c4762a1bSJed Brown return 0; 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianNoMatrix(SNES snes, Vec x, Mat jac, Mat B, void *dummy) 216*d71ae5a4SJacob Faibussowitsch { 217c4762a1bSJed Brown AppCtx *user = (AppCtx *)dummy; 218c4762a1bSJed Brown 2191baa6e33SBarry Smith if (user->variant) PetscCall(MatMFFDSetBase(jac, x, NULL)); 2209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 222c4762a1bSJed Brown return 0; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* -------------------- User-defined monitor ----------------------- */ 226c4762a1bSJed Brown 227*d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, void *dummy) 228*d71ae5a4SJacob Faibussowitsch { 229c4762a1bSJed Brown MonitorCtx *monP = (MonitorCtx *)dummy; 230c4762a1bSJed Brown Vec x; 231c4762a1bSJed Brown MPI_Comm comm; 232c4762a1bSJed Brown 2339566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 23463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm, stdout, "iter = %" PetscInt_FMT ", SNES Function norm %g \n", its, (double)fnorm)); 2359566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &x)); 2369566063dSJacob Faibussowitsch PetscCall(VecView(x, monP->viewer)); 237c4762a1bSJed Brown return 0; 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /*TEST 241c4762a1bSJed Brown 242c4762a1bSJed Brown test: 243c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 244c4762a1bSJed Brown 245c4762a1bSJed Brown test: 246c4762a1bSJed Brown suffix: 2 247c4762a1bSJed Brown args: -variant -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 248c4762a1bSJed Brown output_file: output/ex7_1.out 249c4762a1bSJed Brown 250c4762a1bSJed Brown # uses AIJ matrix to define diagonal matrix for Jacobian preconditioning 251c4762a1bSJed Brown test: 252c4762a1bSJed Brown suffix: 3 253c4762a1bSJed Brown args: -variant -pc_type jacobi -snes_view -ksp_monitor 254c4762a1bSJed Brown 255c4762a1bSJed Brown # uses MATMFFD matrix to define diagonal matrix for Jacobian preconditioning 256c4762a1bSJed Brown test: 257c4762a1bSJed Brown suffix: 4 258c4762a1bSJed Brown args: -variant -pc_type jacobi -puremf -snes_view -ksp_monitor 259c4762a1bSJed Brown 260c4762a1bSJed Brown TEST*/ 261