xref: /petsc/src/snes/tests/ex1.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Solves the nonlinear system, the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\
2c4762a1bSJed Brown This example also illustrates the use of matrix coloring.  Runtime options include:\n\
3c4762a1bSJed Brown   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
4c4762a1bSJed Brown      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\
5c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the x-direction\n\
6c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the y-direction\n\n";
7c4762a1bSJed Brown 
8f6dfbefdSBarry Smith /*
9c4762a1bSJed Brown 
10c4762a1bSJed Brown     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
11c4762a1bSJed Brown     the partial differential equation
12c4762a1bSJed Brown 
13c4762a1bSJed Brown             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
14c4762a1bSJed Brown 
15c4762a1bSJed Brown     with boundary conditions
16c4762a1bSJed Brown 
17c4762a1bSJed Brown              u = 0  for  x = 0, x = 1, y = 0, y = 1.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown     A finite difference approximation with the usual 5-point stencil
20c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a nonlinear
21c4762a1bSJed Brown     system of equations.
22c4762a1bSJed Brown 
23c4762a1bSJed Brown     The parallel version of this code is snes/tutorials/ex5.c
24c4762a1bSJed Brown 
25f6dfbefdSBarry Smith */
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /*
28c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that
29c4762a1bSJed Brown    this file automatically includes:
30c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
31c4762a1bSJed Brown      petscmat.h - matrices
32c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
33c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
34c4762a1bSJed Brown      petscksp.h   - linear solvers
35c4762a1bSJed Brown */
36c4762a1bSJed Brown 
37c4762a1bSJed Brown #include <petscsnes.h>
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /*
40c4762a1bSJed Brown    User-defined application context - contains data needed by the
41c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
42c4762a1bSJed Brown    FormFunction().
43c4762a1bSJed Brown */
44c4762a1bSJed Brown typedef struct {
45c4762a1bSJed Brown   PetscReal param; /* test problem parameter */
46c4762a1bSJed Brown   PetscInt  mx;    /* Discretization in x-direction */
47c4762a1bSJed Brown   PetscInt  my;    /* Discretization in y-direction */
48c4762a1bSJed Brown } AppCtx;
49c4762a1bSJed Brown 
50c4762a1bSJed Brown /*
51c4762a1bSJed Brown    User-defined routines
52c4762a1bSJed Brown */
538b630c82SStefano Zampini static PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
548b630c82SStefano Zampini static PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
558b630c82SStefano Zampini static PetscErrorCode FormInitialGuess(AppCtx *, Vec);
568b630c82SStefano Zampini static PetscErrorCode ConvergenceTest(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
57*2a8381b2SBarry Smith static PetscErrorCode ConvergenceDestroy(PetscCtxRt);
588b630c82SStefano Zampini static PetscErrorCode postcheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
598b630c82SStefano Zampini static PetscErrorCode monitor_change_deltamax(SNES, PetscInt, PetscReal, void *);
60c4762a1bSJed Brown 
main(int argc,char ** argv)61d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
62d71ae5a4SJacob Faibussowitsch {
63c4762a1bSJed Brown   SNES          snes; /* nonlinear solver context */
64c4762a1bSJed Brown   Vec           x, r; /* solution, residual vectors */
65c4762a1bSJed Brown   Mat           J;    /* Jacobian matrix */
66c4762a1bSJed Brown   AppCtx        user; /* user-defined application context */
67c4762a1bSJed Brown   PetscInt      i, its, N, hist_its[50];
68c4762a1bSJed Brown   PetscMPIInt   size;
69c4762a1bSJed Brown   PetscReal     bratu_lambda_max = 6.81, bratu_lambda_min = 0., history[50];
70c4762a1bSJed Brown   MatFDColoring fdcoloring;
718b630c82SStefano Zampini   PetscBool     matrix_free = PETSC_FALSE, flg, fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE, pc = PETSC_FALSE, prunejacobian = PETSC_FALSE, null_appctx = PETSC_TRUE, test_tr_deltamax = PETSC_FALSE;
72c4762a1bSJed Brown   KSP           ksp;
73c4762a1bSJed Brown   PetscInt     *testarray;
74c4762a1bSJed Brown 
75327415f7SBarry Smith   PetscFunctionBeginUser;
76c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
78be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /*
81c4762a1bSJed Brown      Initialize problem parameters
82c4762a1bSJed Brown   */
839371c9d4SSatish Balay   user.mx    = 4;
849371c9d4SSatish Balay   user.my    = 4;
859371c9d4SSatish Balay   user.param = 6.0;
869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL));
879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL));
889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL));
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc", &pc, NULL));
90e00437b9SBarry Smith   PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range");
91c4762a1bSJed Brown   N = user.mx * user.my;
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL));
9378e7fe0eSHong Zhang   PetscCall(PetscOptionsGetBool(NULL, NULL, "-prune_jacobian", &prunejacobian, NULL));
949be84c52SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-null_appctx", &null_appctx, NULL));
958b630c82SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_tr_deltamax", &test_tr_deltamax, NULL));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown      Create nonlinear solver context
99c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100c4762a1bSJed Brown 
1019566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
102c4762a1bSJed Brown   if (pc) {
1039566063dSJacob Faibussowitsch     PetscCall(SNESSetType(snes, SNESNEWTONTR));
1049566063dSJacob Faibussowitsch     PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck, NULL));
105c4762a1bSJed Brown   }
1068b630c82SStefano Zampini   if (test_tr_deltamax) {
1078b630c82SStefano Zampini     PetscCall(SNESSetType(snes, SNESNEWTONTR));
1088b630c82SStefano Zampini     PetscCall(SNESMonitorSet(snes, monitor_change_deltamax, NULL, NULL));
1098b630c82SStefano Zampini   }
110c4762a1bSJed Brown 
1119be84c52SStefano Zampini   /* Test application context handling from Python */
1123a7d0413SPierre Jolivet   if (!null_appctx) PetscCall(SNESSetApplicationContext(snes, (void *)&user));
1139be84c52SStefano Zampini 
114c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
116c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
1199566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, N));
1209566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1219566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /*
124c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
125c4762a1bSJed Brown      solver needs to evaluate the nonlinear function, it will call this
126c4762a1bSJed Brown      routine.
127c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
128c4762a1bSJed Brown         context that provides application-specific data for the
129c4762a1bSJed Brown         function evaluation routine.
130c4762a1bSJed Brown   */
1319566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
135c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown      Create matrix. Here we only approximately preallocate storage space
139c4762a1bSJed Brown      for the Jacobian.  See the users manual for a discussion of better
140c4762a1bSJed Brown      techniques for preallocating matrix memory.
141c4762a1bSJed Brown   */
1429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL));
143c4762a1bSJed Brown   if (!matrix_free) {
144c4762a1bSJed Brown     PetscBool matrix_free_operator = PETSC_FALSE;
1459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL));
146c4762a1bSJed Brown     if (matrix_free_operator) matrix_free = PETSC_FALSE;
147c4762a1bSJed Brown   }
14848a46eb9SPierre Jolivet   if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /*
151c4762a1bSJed Brown      This option will cause the Jacobian to be computed via finite differences
152c4762a1bSJed Brown     efficiently using a coloring of the columns of the matrix.
153c4762a1bSJed Brown   */
1549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL));
15500045ab3SPierre Jolivet   PetscCheck(!matrix_free || !fd_coloring, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Use only one of -snes_mf, -snes_fd_coloring options! You can do -snes_mf_operator -snes_fd_coloring");
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   if (fd_coloring) {
158c4762a1bSJed Brown     ISColoring  iscoloring;
159c4762a1bSJed Brown     MatColoring mc;
16078e7fe0eSHong Zhang     if (prunejacobian) {
16178e7fe0eSHong Zhang       /* Initialize x with random nonzero values so that the nonzeros in the Jacobian
16278e7fe0eSHong Zhang          can better reflect the sparsity structure of the Jacobian. */
16378e7fe0eSHong Zhang       PetscRandom rctx;
16478e7fe0eSHong Zhang       PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
16578e7fe0eSHong Zhang       PetscCall(PetscRandomSetInterval(rctx, 1.0, 2.0));
16678e7fe0eSHong Zhang       PetscCall(VecSetRandom(x, rctx));
16778e7fe0eSHong Zhang       PetscCall(PetscRandomDestroy(&rctx));
16878e7fe0eSHong Zhang     }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown     /*
171c4762a1bSJed Brown       This initializes the nonzero structure of the Jacobian. This is artificial
172c4762a1bSJed Brown       because clearly if we had a routine to compute the Jacobian we won't need
173c4762a1bSJed Brown       to use finite differences.
174c4762a1bSJed Brown     */
1759566063dSJacob Faibussowitsch     PetscCall(FormJacobian(snes, x, J, J, &user));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown     /*
178c4762a1bSJed Brown        Color the matrix, i.e. determine groups of columns that share no common
179a5b23f4aSJose E. Roman       rows. These columns in the Jacobian can all be computed simultaneously.
180c4762a1bSJed Brown     */
1819566063dSJacob Faibussowitsch     PetscCall(MatColoringCreate(J, &mc));
1829566063dSJacob Faibussowitsch     PetscCall(MatColoringSetType(mc, MATCOLORINGSL));
1839566063dSJacob Faibussowitsch     PetscCall(MatColoringSetFromOptions(mc));
1849566063dSJacob Faibussowitsch     PetscCall(MatColoringApply(mc, &iscoloring));
1859566063dSJacob Faibussowitsch     PetscCall(MatColoringDestroy(&mc));
186c4762a1bSJed Brown     /*
187c4762a1bSJed Brown        Create the data structure that SNESComputeJacobianDefaultColor() uses
188c4762a1bSJed Brown        to compute the actual Jacobians via finite differences.
189c4762a1bSJed Brown     */
1909566063dSJacob Faibussowitsch     PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring));
1912ba42892SBarry Smith     PetscCall(MatFDColoringSetFunction(fdcoloring, (MatFDColoringFn *)FormFunction, &user));
1929566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1939566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring));
194c4762a1bSJed Brown     /*
195c4762a1bSJed Brown         Tell SNES to use the routine SNESComputeJacobianDefaultColor()
196c4762a1bSJed Brown       to compute Jacobians.
197c4762a1bSJed Brown     */
1989566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring));
1999566063dSJacob Faibussowitsch     PetscCall(ISColoringDestroy(&iscoloring));
20078e7fe0eSHong Zhang     if (prunejacobian) PetscCall(SNESPruneJacobianColor(snes, J, J));
201c4762a1bSJed Brown   }
202c4762a1bSJed Brown   /*
203c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
204c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
205c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
206c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
207c4762a1bSJed Brown         context that provides application-specific data for the
208c4762a1bSJed Brown         Jacobian evaluation routine.
209c4762a1bSJed Brown       - The user can override with:
210c4762a1bSJed Brown          -snes_fd : default finite differencing approximation of Jacobian
211c4762a1bSJed Brown          -snes_mf : matrix-free Newton-Krylov method with no preconditioning
212c4762a1bSJed Brown                     (unless user explicitly sets preconditioner)
2137addb90fSBarry Smith          -snes_mf_operator : form matrix used to construct the preconditioner as set by the user,
214c4762a1bSJed Brown                              but use matrix-free approx for Jacobian-vector
215c4762a1bSJed Brown                              products within Newton-Krylov method
216c4762a1bSJed Brown   */
217c4762a1bSJed Brown   else if (!matrix_free) {
2189566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
223c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /*
226c4762a1bSJed Brown      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
227c4762a1bSJed Brown   */
2289566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /*
231c4762a1bSJed Brown      Set array that saves the function norms.  This array is intended
232c4762a1bSJed Brown      when the user wants to save the convergence history for later use
233c4762a1bSJed Brown      rather than just to view the function norms via -snes_monitor.
234c4762a1bSJed Brown   */
2359566063dSJacob Faibussowitsch   PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /*
238c4762a1bSJed Brown       Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the
239c4762a1bSJed Brown       user provided test before the specialized test. The convergence context is just an array to
240c4762a1bSJed Brown       test that it gets properly freed at the end
241c4762a1bSJed Brown   */
242c4762a1bSJed Brown   if (use_convergence_test) {
2439566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
2449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5, &testarray));
2459566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy));
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
250c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
251c4762a1bSJed Brown   /*
252c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
253c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
254c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
255c4762a1bSJed Brown      this vector to zero by calling VecSet().
256c4762a1bSJed Brown   */
2579566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user, x));
2589566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, x));
2599566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
26063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /*
263c4762a1bSJed Brown      Print the convergence history.  This is intended just to demonstrate
264c4762a1bSJed Brown      use of the data attained via SNESSetConvergenceHistory().
265c4762a1bSJed Brown   */
2669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg));
267c4762a1bSJed Brown   if (flg) {
26848a46eb9SPierre Jolivet     for (i = 0; i < its + 1; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n", i, hist_its[i], (double)history[i]));
269c4762a1bSJed Brown   }
270c4762a1bSJed Brown 
2713201ab8dSStefano Zampini   /* Test NewtonTR API */
2723201ab8dSStefano Zampini   PetscCall(SNESNewtonTRSetTolerances(snes, 1.0, 2.0, 3.0));
2733201ab8dSStefano Zampini   PetscCall(SNESNewtonTRSetUpdateParameters(snes, 4.0, 5.0, 6.0, 7.0, 8.0));
2743201ab8dSStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
2753201ab8dSStefano Zampini   if (flg) {
2763201ab8dSStefano Zampini     PetscReal tmp[5];
2773201ab8dSStefano Zampini 
2783201ab8dSStefano Zampini     PetscCall(SNESNewtonTRGetTolerances(snes, &tmp[0], &tmp[1], &tmp[2]));
2793201ab8dSStefano Zampini     PetscCheck(tmp[0] == 1.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2803201ab8dSStefano Zampini     PetscCheck(tmp[1] == 2.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2813201ab8dSStefano Zampini     PetscCheck(tmp[2] == 3.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2823201ab8dSStefano Zampini     PetscCall(SNESNewtonTRGetUpdateParameters(snes, &tmp[0], &tmp[1], &tmp[2], &tmp[3], &tmp[4]));
2833201ab8dSStefano Zampini     PetscCheck(tmp[0] == 4.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2843201ab8dSStefano Zampini     PetscCheck(tmp[1] == 5.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2853201ab8dSStefano Zampini     PetscCheck(tmp[2] == 6.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2863201ab8dSStefano Zampini     PetscCheck(tmp[3] == 7.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2873201ab8dSStefano Zampini     PetscCheck(tmp[4] == 8.0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Wrong value");
2883201ab8dSStefano Zampini   }
2893201ab8dSStefano Zampini 
290c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
291c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
292c4762a1bSJed Brown      are no longer needed.
293c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
294c4762a1bSJed Brown 
29548a46eb9SPierre Jolivet   if (!matrix_free) PetscCall(MatDestroy(&J));
29648a46eb9SPierre Jolivet   if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring));
2979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2999566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
3009566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
301b122ec5aSJacob Faibussowitsch   return 0;
302c4762a1bSJed Brown }
303f6dfbefdSBarry Smith 
304c4762a1bSJed Brown /*
305c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
306c4762a1bSJed Brown 
307c4762a1bSJed Brown    Input Parameters:
308c4762a1bSJed Brown    user - user-defined application context
309c4762a1bSJed Brown    X - vector
310c4762a1bSJed Brown 
311c4762a1bSJed Brown    Output Parameter:
312c4762a1bSJed Brown    X - vector
313c4762a1bSJed Brown  */
FormInitialGuess(AppCtx * user,Vec X)314d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, Vec X)
315d71ae5a4SJacob Faibussowitsch {
316c4762a1bSJed Brown   PetscInt     i, j, row, mx, my;
317c4762a1bSJed Brown   PetscReal    lambda, temp1, temp, hx, hy;
318c4762a1bSJed Brown   PetscScalar *x;
319c4762a1bSJed Brown 
3203ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
321c4762a1bSJed Brown   mx     = user->mx;
322c4762a1bSJed Brown   my     = user->my;
323c4762a1bSJed Brown   lambda = user->param;
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx - 1);
326c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(my - 1);
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   /*
329c4762a1bSJed Brown      Get a pointer to vector data.
330c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
331c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
332c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
333c4762a1bSJed Brown          the array.
334c4762a1bSJed Brown   */
3359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
336c4762a1bSJed Brown   temp1 = lambda / (lambda + 1.0);
337c4762a1bSJed Brown   for (j = 0; j < my; j++) {
338c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy;
339c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
340c4762a1bSJed Brown       row = i + j * mx;
341c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
342c4762a1bSJed Brown         x[row] = 0.0;
343c4762a1bSJed Brown         continue;
344c4762a1bSJed Brown       }
345c4762a1bSJed Brown       x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp));
346c4762a1bSJed Brown     }
347c4762a1bSJed Brown   }
348c4762a1bSJed Brown 
349c4762a1bSJed Brown   /*
350c4762a1bSJed Brown      Restore vector
351c4762a1bSJed Brown   */
3529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354c4762a1bSJed Brown }
355f6dfbefdSBarry Smith 
356c4762a1bSJed Brown /*
357c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
358c4762a1bSJed Brown 
359c4762a1bSJed Brown    Input Parameters:
360c4762a1bSJed Brown .  snes - the SNES context
361c4762a1bSJed Brown .  X - input vector
362c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
363c4762a1bSJed Brown 
364c4762a1bSJed Brown    Output Parameter:
365c4762a1bSJed Brown .  F - function vector
366c4762a1bSJed Brown  */
FormFunction(SNES snes,Vec X,Vec F,void * ptr)367d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr)
368d71ae5a4SJacob Faibussowitsch {
369c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
370c4762a1bSJed Brown   PetscInt           i, j, row, mx, my;
371c4762a1bSJed Brown   PetscReal          two = 2.0, one = 1.0, lambda, hx, hy, hxdhy, hydhx;
372c4762a1bSJed Brown   PetscScalar        ut, ub, ul, ur, u, uxx, uyy, sc, *f;
373c4762a1bSJed Brown   const PetscScalar *x;
374c4762a1bSJed Brown 
3753ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
376c4762a1bSJed Brown   mx     = user->mx;
377c4762a1bSJed Brown   my     = user->my;
378c4762a1bSJed Brown   lambda = user->param;
379c4762a1bSJed Brown   hx     = one / (PetscReal)(mx - 1);
380c4762a1bSJed Brown   hy     = one / (PetscReal)(my - 1);
381c4762a1bSJed Brown   sc     = hx * hy;
382c4762a1bSJed Brown   hxdhy  = hx / hy;
383c4762a1bSJed Brown   hydhx  = hy / hx;
384c4762a1bSJed Brown 
385c4762a1bSJed Brown   /*
386c4762a1bSJed Brown      Get pointers to vector data
387c4762a1bSJed Brown   */
3889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
3899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
390c4762a1bSJed Brown 
391c4762a1bSJed Brown   /*
392c4762a1bSJed Brown      Compute function
393c4762a1bSJed Brown   */
394c4762a1bSJed Brown   for (j = 0; j < my; j++) {
395c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
396c4762a1bSJed Brown       row = i + j * mx;
397c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
398c4762a1bSJed Brown         f[row] = x[row];
399c4762a1bSJed Brown         continue;
400c4762a1bSJed Brown       }
401c4762a1bSJed Brown       u      = x[row];
402c4762a1bSJed Brown       ub     = x[row - mx];
403c4762a1bSJed Brown       ul     = x[row - 1];
404c4762a1bSJed Brown       ut     = x[row + mx];
405c4762a1bSJed Brown       ur     = x[row + 1];
406c4762a1bSJed Brown       uxx    = (-ur + two * u - ul) * hydhx;
407c4762a1bSJed Brown       uyy    = (-ut + two * u - ub) * hxdhy;
408c4762a1bSJed Brown       f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u);
409c4762a1bSJed Brown     }
410c4762a1bSJed Brown   }
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   /*
413c4762a1bSJed Brown      Restore vectors
414c4762a1bSJed Brown   */
4159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
4169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
4173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
418c4762a1bSJed Brown }
419f6dfbefdSBarry Smith 
420c4762a1bSJed Brown /*
421c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
422c4762a1bSJed Brown 
423c4762a1bSJed Brown    Input Parameters:
424c4762a1bSJed Brown .  snes - the SNES context
425c4762a1bSJed Brown .  x - input vector
426c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetJacobian()
427c4762a1bSJed Brown 
428c4762a1bSJed Brown    Output Parameters:
429c4762a1bSJed Brown .  A - Jacobian matrix
4307addb90fSBarry Smith .  B - optionally different matrix used to construct the preconditioner
4310b4b7b1cSBarry Smith 
432c4762a1bSJed Brown */
FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void * ptr)433d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr)
434d71ae5a4SJacob Faibussowitsch {
435da81f932SPierre Jolivet   AppCtx            *user = (AppCtx *)ptr; /* user-defined application context */
436c4762a1bSJed Brown   PetscInt           i, j, row, mx, my, col[5];
437c4762a1bSJed Brown   PetscScalar        two = 2.0, one = 1.0, lambda, v[5], sc;
438c4762a1bSJed Brown   const PetscScalar *x;
439c4762a1bSJed Brown   PetscReal          hx, hy, hxdhy, hydhx;
440c4762a1bSJed Brown 
4413ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
442c4762a1bSJed Brown   mx     = user->mx;
443c4762a1bSJed Brown   my     = user->my;
444c4762a1bSJed Brown   lambda = user->param;
445c4762a1bSJed Brown   hx     = 1.0 / (PetscReal)(mx - 1);
446c4762a1bSJed Brown   hy     = 1.0 / (PetscReal)(my - 1);
447c4762a1bSJed Brown   sc     = hx * hy;
448c4762a1bSJed Brown   hxdhy  = hx / hy;
449c4762a1bSJed Brown   hydhx  = hy / hx;
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   /*
452c4762a1bSJed Brown      Get pointer to vector data
453c4762a1bSJed Brown   */
4549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   /*
457c4762a1bSJed Brown      Compute entries of the Jacobian
458c4762a1bSJed Brown   */
459c4762a1bSJed Brown   for (j = 0; j < my; j++) {
460c4762a1bSJed Brown     for (i = 0; i < mx; i++) {
461c4762a1bSJed Brown       row = i + j * mx;
462c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
4639566063dSJacob Faibussowitsch         PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES));
464c4762a1bSJed Brown         continue;
465c4762a1bSJed Brown       }
4669371c9d4SSatish Balay       v[0]   = -hxdhy;
4679371c9d4SSatish Balay       col[0] = row - mx;
4689371c9d4SSatish Balay       v[1]   = -hydhx;
4699371c9d4SSatish Balay       col[1] = row - 1;
4709371c9d4SSatish Balay       v[2]   = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]);
4719371c9d4SSatish Balay       col[2] = row;
4729371c9d4SSatish Balay       v[3]   = -hydhx;
4739371c9d4SSatish Balay       col[3] = row + 1;
4749371c9d4SSatish Balay       v[4]   = -hxdhy;
4759371c9d4SSatish Balay       col[4] = row + mx;
4769566063dSJacob Faibussowitsch       PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES));
477c4762a1bSJed Brown     }
478c4762a1bSJed Brown   }
479c4762a1bSJed Brown 
480c4762a1bSJed Brown   /*
481c4762a1bSJed Brown      Restore vector
482c4762a1bSJed Brown   */
4839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
484c4762a1bSJed Brown 
485c4762a1bSJed Brown   /*
486c4762a1bSJed Brown      Assemble matrix
487c4762a1bSJed Brown   */
4889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
4899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
490c4762a1bSJed Brown 
491c4762a1bSJed Brown   if (jac != J) {
4929566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
4939566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
494c4762a1bSJed Brown   }
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
496c4762a1bSJed Brown }
497c4762a1bSJed Brown 
ConvergenceTest(KSP ksp,PetscInt it,PetscReal nrm,KSPConvergedReason * reason,PetscCtx ctx)498*2a8381b2SBarry Smith PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, PetscCtx ctx)
499d71ae5a4SJacob Faibussowitsch {
5003ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
501c4762a1bSJed Brown   *reason = KSP_CONVERGED_ITERATING;
502c4762a1bSJed Brown   if (it > 1) {
503c4762a1bSJed Brown     *reason = KSP_CONVERGED_ITS;
5049566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n"));
505c4762a1bSJed Brown   }
5063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
507c4762a1bSJed Brown }
508c4762a1bSJed Brown 
ConvergenceDestroy(PetscCtxRt ctx)509*2a8381b2SBarry Smith PetscErrorCode ConvergenceDestroy(PetscCtxRt ctx)
510d71ae5a4SJacob Faibussowitsch {
5113ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
5129566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n"));
513*2a8381b2SBarry Smith   PetscCall(PetscFree(*(void **)ctx));
5143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
515c4762a1bSJed Brown }
516c4762a1bSJed Brown 
postcheck(SNES snes,Vec x,Vec y,Vec w,PetscBool * changed_y,PetscBool * changed_w,PetscCtx ctx)517*2a8381b2SBarry Smith PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, PetscCtx ctx)
518d71ae5a4SJacob Faibussowitsch {
519c4762a1bSJed Brown   PetscReal norm;
520c4762a1bSJed Brown   Vec       tmp;
521c4762a1bSJed Brown 
5223ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
5239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &tmp));
5249566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tmp, -1.0, x, w));
5259566063dSJacob Faibussowitsch   PetscCall(VecNorm(tmp, NORM_2, &norm));
5269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp));
5279566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm));
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
529c4762a1bSJed Brown }
530c4762a1bSJed Brown 
monitor_change_deltamax(SNES snes,PetscInt it,PetscReal fnorm,PetscCtx ctx)531*2a8381b2SBarry Smith PetscErrorCode monitor_change_deltamax(SNES snes, PetscInt it, PetscReal fnorm, PetscCtx ctx)
5328b630c82SStefano Zampini {
5338b630c82SStefano Zampini   PetscFunctionBeginUser;
5348b630c82SStefano Zampini   if (it == 0) PetscCall(SNESNewtonTRSetTolerances(snes, PETSC_CURRENT, 0.01, PETSC_CURRENT));
5358b630c82SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
5368b630c82SStefano Zampini }
5378b630c82SStefano Zampini 
538c4762a1bSJed Brown /*TEST
539c4762a1bSJed Brown 
540c4762a1bSJed Brown    build:
541c4762a1bSJed Brown       requires: !single
542c4762a1bSJed Brown 
543c4762a1bSJed Brown    test:
544c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always
545c4762a1bSJed Brown 
546c4762a1bSJed Brown    test:
547c4762a1bSJed Brown       suffix: 2
548c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
549c4762a1bSJed Brown 
550c4762a1bSJed Brown    test:
5518b630c82SStefano Zampini       suffix: 2_trdeltamax_change
5528b630c82SStefano Zampini       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -pc -test_tr_deltamax
5538b630c82SStefano Zampini 
5548b630c82SStefano Zampini    test:
555c4762a1bSJed Brown       suffix: 2a
556c4762a1bSJed Brown       filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault"
557c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info
558dfd57a17SPierre Jolivet       requires: defined(PETSC_USE_INFO)
559c4762a1bSJed Brown 
560c4762a1bSJed Brown    test:
561c4762a1bSJed Brown       suffix: 2b
562c4762a1bSJed Brown       filter: grep -i  "User provided convergence test" > /dev/null  && echo "Found User provided convergence test"
563c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info
564dfd57a17SPierre Jolivet       requires: defined(PETSC_USE_INFO)
565c4762a1bSJed Brown 
566c4762a1bSJed Brown    test:
56724fb275aSStefano Zampini       suffix: 2c
56824fb275aSStefano Zampini       args: -snes_converged_reason -snes_type newtontr -snes_tr_qn {{same different}separate output} -pc_type mat -snes_view -snes_tr_qn_mat_type lmvmdfp -snes_tr_norm_type infinity
56924fb275aSStefano Zampini 
57024fb275aSStefano Zampini    test:
571c4762a1bSJed Brown       suffix: 3
572c4762a1bSJed Brown       args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always
573c4762a1bSJed Brown 
574c4762a1bSJed Brown    test:
575c4762a1bSJed Brown       suffix: 4
576c4762a1bSJed Brown       args: -pc -par 6.807 -snes_monitor -snes_converged_reason
577c4762a1bSJed Brown 
57878e7fe0eSHong Zhang    test:
57978e7fe0eSHong Zhang       suffix: 5
58078e7fe0eSHong Zhang       args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always -prune_jacobian
58178e7fe0eSHong Zhang       output_file: output/ex1_3.out
582ccb5f961SBarry Smith 
583ccb5f961SBarry Smith    test:
584ccb5f961SBarry Smith       suffix: 6
585ccb5f961SBarry Smith       args: -snes_monitor draw:image:testfile -viewer_view
586ccb5f961SBarry Smith 
5879be84c52SStefano Zampini    test:
5889be84c52SStefano Zampini       suffix: python
5899be84c52SStefano Zampini       requires: petsc4py
5909be84c52SStefano Zampini       args: -python -snes_type python -snes_python_type ex1.py:MySNES -snes_view -null_appctx {{0 1}separate output}
5919be84c52SStefano Zampini       localrunfiles: ex1.py
5929be84c52SStefano Zampini 
593c4762a1bSJed Brown TEST*/
594