xref: /petsc/src/snes/tests/ex7.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
23d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
24d71ae5a4SJacob 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 
105d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *dummy)
106d71ae5a4SJacob Faibussowitsch {
107c4762a1bSJed Brown   const PetscScalar *xx, *FF;
108c4762a1bSJed Brown   PetscScalar       *ff, d;
109c4762a1bSJed Brown   PetscInt           i, n;
110c4762a1bSJed Brown 
111*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(f, &ff));
1149566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead((Vec)dummy, &FF));
1159566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
1169371c9d4SSatish Balay   d     = (PetscReal)(n - 1);
1179371c9d4SSatish Balay   d     = d * d;
118c4762a1bSJed Brown   ff[0] = xx[0];
119c4762a1bSJed 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];
120c4762a1bSJed Brown   ff[n - 1] = xx[n - 1] - 1.0;
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(f, &ff));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead((Vec)dummy, &FF));
124*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctioni(void *dummy, PetscInt i, Vec x, PetscScalar *s)
128d71ae5a4SJacob Faibussowitsch {
129c4762a1bSJed Brown   const PetscScalar *xx, *FF;
130c4762a1bSJed Brown   PetscScalar        d;
131c4762a1bSJed Brown   PetscInt           n;
132c4762a1bSJed Brown   SNES               snes = (SNES)dummy;
133c4762a1bSJed Brown   Vec                F;
134c4762a1bSJed Brown 
135*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1369566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, NULL, NULL, (void **)&F));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(F, &FF));
1399566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
1409371c9d4SSatish Balay   d = (PetscReal)(n - 1);
1419371c9d4SSatish Balay   d = d * d;
142c4762a1bSJed Brown   if (i == 0) {
143c4762a1bSJed Brown     *s = xx[0];
144c4762a1bSJed Brown   } else if (i == n - 1) {
145c4762a1bSJed Brown     *s = xx[n - 1] - 1.0;
146c4762a1bSJed Brown   } else {
147c4762a1bSJed Brown     *s = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - FF[i];
148c4762a1bSJed Brown   }
1499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(F, &FF));
151*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
152c4762a1bSJed Brown }
153c4762a1bSJed Brown 
154c4762a1bSJed Brown /*
155c4762a1bSJed Brown 
156c4762a1bSJed Brown    Example function that when differenced produces the same matrix free Jacobian as FormFunction()
157c4762a1bSJed Brown    this is provided to show how a user can provide a different function
158c4762a1bSJed Brown */
159d71ae5a4SJacob Faibussowitsch PetscErrorCode OtherFunctionForDifferencing(void *dummy, Vec x, Vec f)
160d71ae5a4SJacob Faibussowitsch {
161*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1629566063dSJacob Faibussowitsch   PetscCall(FormFunction(NULL, x, f, dummy));
1639566063dSJacob Faibussowitsch   PetscCall(VecShift(f, 1.0));
164*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown /* --------------------  Form initial approximation ----------------- */
168c4762a1bSJed Brown 
169d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(SNES snes, Vec x)
170d71ae5a4SJacob Faibussowitsch {
171c4762a1bSJed Brown   PetscScalar pfive = .50;
172*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1739566063dSJacob Faibussowitsch   PetscCall(VecSet(x, pfive));
174*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown /* --------------------  Evaluate Jacobian F'(x) -------------------- */
177c4762a1bSJed Brown /*  Evaluates a matrix that is used to precondition the matrix-free
178a5b23f4aSJose E. Roman     jacobian. In this case, the explicit preconditioner matrix is
179c4762a1bSJed Brown     also EXACTLY the Jacobian. In general, it would be some lower
180c4762a1bSJed Brown     order, simplified apprioximation */
181c4762a1bSJed Brown 
182d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
183d71ae5a4SJacob Faibussowitsch {
184c4762a1bSJed Brown   const PetscScalar *xx;
185c4762a1bSJed Brown   PetscScalar        A[3], d;
186c4762a1bSJed Brown   PetscInt           i, n, j[3];
187c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)dummy;
188c4762a1bSJed Brown 
189*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
1919566063dSJacob Faibussowitsch   PetscCall(VecGetSize(x, &n));
1929371c9d4SSatish Balay   d = (PetscReal)(n - 1);
1939371c9d4SSatish Balay   d = d * d;
194c4762a1bSJed Brown 
1959371c9d4SSatish Balay   i    = 0;
1969371c9d4SSatish Balay   A[0] = 1.0;
1979566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, &A[0], INSERT_VALUES));
198c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) {
1999371c9d4SSatish Balay     j[0] = i - 1;
2009371c9d4SSatish Balay     j[1] = i;
2019371c9d4SSatish Balay     j[2] = i + 1;
2029371c9d4SSatish Balay     A[0] = d;
2039371c9d4SSatish Balay     A[1] = -2.0 * d + 2.0 * xx[i];
2049371c9d4SSatish Balay     A[2] = d;
2059566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, 3, j, A, INSERT_VALUES));
206c4762a1bSJed Brown   }
2079371c9d4SSatish Balay   i    = n - 1;
2089371c9d4SSatish Balay   A[0] = 1.0;
2099566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, &i, 1, &i, &A[0], INSERT_VALUES));
2109566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
2129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
213c4762a1bSJed Brown 
2141baa6e33SBarry Smith   if (user->variant) PetscCall(MatMFFDSetBase(jac, x, NULL));
2159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
217*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218c4762a1bSJed Brown }
219c4762a1bSJed Brown 
220d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianNoMatrix(SNES snes, Vec x, Mat jac, Mat B, void *dummy)
221d71ae5a4SJacob Faibussowitsch {
222c4762a1bSJed Brown   AppCtx *user = (AppCtx *)dummy;
223c4762a1bSJed Brown 
224*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2251baa6e33SBarry Smith   if (user->variant) PetscCall(MatMFFDSetBase(jac, x, NULL));
2269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
2279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
228*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
229c4762a1bSJed Brown }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown /* --------------------  User-defined monitor ----------------------- */
232c4762a1bSJed Brown 
233d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, void *dummy)
234d71ae5a4SJacob Faibussowitsch {
235c4762a1bSJed Brown   MonitorCtx *monP = (MonitorCtx *)dummy;
236c4762a1bSJed Brown   Vec         x;
237c4762a1bSJed Brown   MPI_Comm    comm;
238c4762a1bSJed Brown 
239*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
2409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
24163a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, stdout, "iter = %" PetscInt_FMT ", SNES Function norm %g \n", its, (double)fnorm));
2429566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &x));
2439566063dSJacob Faibussowitsch   PetscCall(VecView(x, monP->viewer));
244*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown /*TEST
248c4762a1bSJed Brown 
249c4762a1bSJed Brown    test:
250c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    test:
253c4762a1bSJed Brown       suffix: 2
254c4762a1bSJed Brown       args: -variant -ksp_gmres_cgs_refinement_type refine_always  -snes_monitor_short
255c4762a1bSJed Brown       output_file: output/ex7_1.out
256c4762a1bSJed Brown 
257c4762a1bSJed Brown    # uses AIJ matrix to define diagonal matrix for Jacobian preconditioning
258c4762a1bSJed Brown    test:
259c4762a1bSJed Brown       suffix: 3
260c4762a1bSJed Brown       args: -variant -pc_type jacobi -snes_view -ksp_monitor
261c4762a1bSJed Brown 
262c4762a1bSJed Brown    # uses MATMFFD matrix to define diagonal matrix for Jacobian preconditioning
263c4762a1bSJed Brown    test:
264c4762a1bSJed Brown       suffix: 4
265c4762a1bSJed Brown       args: -variant -pc_type jacobi -puremf  -snes_view -ksp_monitor
266c4762a1bSJed Brown 
267c4762a1bSJed Brown TEST*/
268