xref: /petsc/src/snes/tests/ex1.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the nonlinear system, the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\
3c4762a1bSJed Brown This example also illustrates the use of matrix coloring.  Runtime options include:\n\
4c4762a1bSJed Brown   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
5c4762a1bSJed Brown      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\
6c4762a1bSJed Brown   -mx <xg>, where <xg> = number of grid points in the x-direction\n\
7c4762a1bSJed Brown   -my <yg>, where <yg> = number of grid points in the y-direction\n\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
12c4762a1bSJed Brown     the partial differential equation
13c4762a1bSJed Brown 
14c4762a1bSJed Brown             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
15c4762a1bSJed Brown 
16c4762a1bSJed Brown     with boundary conditions
17c4762a1bSJed Brown 
18c4762a1bSJed Brown              u = 0  for  x = 0, x = 1, y = 0, y = 1.
19c4762a1bSJed Brown 
20c4762a1bSJed Brown     A finite difference approximation with the usual 5-point stencil
21c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a nonlinear
22c4762a1bSJed Brown     system of equations.
23c4762a1bSJed Brown 
24c4762a1bSJed Brown     The parallel version of this code is snes/tutorials/ex5.c
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ------------------------------------------------------------------------- */
27c4762a1bSJed Brown 
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that
30c4762a1bSJed Brown    this file automatically includes:
31c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
32c4762a1bSJed Brown      petscmat.h - matrices
33c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
34c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
35c4762a1bSJed Brown      petscksp.h   - linear solvers
36c4762a1bSJed Brown */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown #include <petscsnes.h>
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /*
41c4762a1bSJed Brown    User-defined application context - contains data needed by the
42c4762a1bSJed Brown    application-provided call-back routines, FormJacobian() and
43c4762a1bSJed Brown    FormFunction().
44c4762a1bSJed Brown */
45c4762a1bSJed Brown typedef struct {
46c4762a1bSJed Brown   PetscReal param;              /* test problem parameter */
47c4762a1bSJed Brown   PetscInt  mx;                 /* Discretization in x-direction */
48c4762a1bSJed Brown   PetscInt  my;                 /* Discretization in y-direction */
49c4762a1bSJed Brown } AppCtx;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown /*
52c4762a1bSJed Brown    User-defined routines
53c4762a1bSJed Brown */
54c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
55c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
56c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx*,Vec);
57c4762a1bSJed Brown extern PetscErrorCode ConvergenceTest(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
58c4762a1bSJed Brown extern PetscErrorCode ConvergenceDestroy(void*);
59c4762a1bSJed Brown extern PetscErrorCode postcheck(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
60c4762a1bSJed Brown 
61c4762a1bSJed Brown int main(int argc,char **argv)
62c4762a1bSJed Brown {
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;
71c4762a1bSJed Brown   PetscBool      matrix_free = PETSC_FALSE,flg,fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE,pc = PETSC_FALSE;
72c4762a1bSJed Brown   KSP            ksp;
73c4762a1bSJed Brown   PetscInt       *testarray;
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
77be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /*
80c4762a1bSJed Brown      Initialize problem parameters
81c4762a1bSJed Brown   */
82c4762a1bSJed Brown   user.mx = 4; user.my = 4; user.param = 6.0;
839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL));
849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL));
859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL));
869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-pc",&pc,NULL));
87e00437b9SBarry Smith   PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda is out of range");
88c4762a1bSJed Brown   N = user.mx*user.my;
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_convergence_test",&use_convergence_test,NULL));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown      Create nonlinear solver context
93c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   if (pc) {
989566063dSJacob Faibussowitsch     PetscCall(SNESSetType(snes,SNESNEWTONTR));
999566063dSJacob Faibussowitsch     PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck,NULL));
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
104c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
1079566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,PETSC_DECIDE,N));
1089566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
1099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /*
112c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
113c4762a1bSJed Brown      solver needs to evaluate the nonlinear function, it will call this
114c4762a1bSJed Brown      routine.
115c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
116c4762a1bSJed Brown         context that provides application-specific data for the
117c4762a1bSJed Brown         function evaluation routine.
118c4762a1bSJed Brown   */
1199566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)&user));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
123c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /*
126c4762a1bSJed Brown      Create matrix. Here we only approximately preallocate storage space
127c4762a1bSJed Brown      for the Jacobian.  See the users manual for a discussion of better
128c4762a1bSJed Brown      techniques for preallocating matrix memory.
129c4762a1bSJed Brown   */
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL));
131c4762a1bSJed Brown   if (!matrix_free) {
132c4762a1bSJed Brown     PetscBool matrix_free_operator = PETSC_FALSE;
1339566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_mf_operator",&matrix_free_operator,NULL));
134c4762a1bSJed Brown     if (matrix_free_operator) matrix_free = PETSC_FALSE;
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown   if (!matrix_free) {
1379566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&J));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /*
141c4762a1bSJed Brown      This option will cause the Jacobian to be computed via finite differences
142c4762a1bSJed Brown     efficiently using a coloring of the columns of the matrix.
143c4762a1bSJed Brown   */
1449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_fd_coloring",&fd_coloring,NULL));
145e00437b9SBarry Smith   PetscCheck(!matrix_free || !fd_coloring,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Use only one of -snes_mf, -snes_fd_coloring options!\nYou can do -snes_mf_operator -snes_fd_coloring");
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   if (fd_coloring) {
148c4762a1bSJed Brown     ISColoring   iscoloring;
149c4762a1bSJed Brown     MatColoring  mc;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown     /*
152c4762a1bSJed Brown       This initializes the nonzero structure of the Jacobian. This is artificial
153c4762a1bSJed Brown       because clearly if we had a routine to compute the Jacobian we won't need
154c4762a1bSJed Brown       to use finite differences.
155c4762a1bSJed Brown     */
1569566063dSJacob Faibussowitsch     PetscCall(FormJacobian(snes,x,J,J,&user));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     /*
159c4762a1bSJed Brown        Color the matrix, i.e. determine groups of columns that share no common
160a5b23f4aSJose E. Roman       rows. These columns in the Jacobian can all be computed simultaneously.
161c4762a1bSJed Brown     */
1629566063dSJacob Faibussowitsch     PetscCall(MatColoringCreate(J,&mc));
1639566063dSJacob Faibussowitsch     PetscCall(MatColoringSetType(mc,MATCOLORINGSL));
1649566063dSJacob Faibussowitsch     PetscCall(MatColoringSetFromOptions(mc));
1659566063dSJacob Faibussowitsch     PetscCall(MatColoringApply(mc,&iscoloring));
1669566063dSJacob Faibussowitsch     PetscCall(MatColoringDestroy(&mc));
167c4762a1bSJed Brown     /*
168c4762a1bSJed Brown        Create the data structure that SNESComputeJacobianDefaultColor() uses
169c4762a1bSJed Brown        to compute the actual Jacobians via finite differences.
170c4762a1bSJed Brown     */
1719566063dSJacob Faibussowitsch     PetscCall(MatFDColoringCreate(J,iscoloring,&fdcoloring));
1729566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFunction(fdcoloring,(PetscErrorCode (*)(void))FormFunction,&user));
1739566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetFromOptions(fdcoloring));
1749566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUp(J,iscoloring,fdcoloring));
175c4762a1bSJed Brown     /*
176c4762a1bSJed Brown         Tell SNES to use the routine SNESComputeJacobianDefaultColor()
177c4762a1bSJed Brown       to compute Jacobians.
178c4762a1bSJed Brown     */
1799566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring));
1809566063dSJacob Faibussowitsch     PetscCall(ISColoringDestroy(&iscoloring));
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown   /*
183c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
184c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
185c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
186c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
187c4762a1bSJed Brown         context that provides application-specific data for the
188c4762a1bSJed Brown         Jacobian evaluation routine.
189c4762a1bSJed Brown       - The user can override with:
190c4762a1bSJed Brown          -snes_fd : default finite differencing approximation of Jacobian
191c4762a1bSJed Brown          -snes_mf : matrix-free Newton-Krylov method with no preconditioning
192c4762a1bSJed Brown                     (unless user explicitly sets preconditioner)
193c4762a1bSJed Brown          -snes_mf_operator : form preconditioning matrix as set by the user,
194c4762a1bSJed Brown                              but use matrix-free approx for Jacobian-vector
195c4762a1bSJed Brown                              products within Newton-Krylov method
196c4762a1bSJed Brown   */
197c4762a1bSJed Brown   else if (!matrix_free) {
1989566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,(void*)&user));
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
203c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /*
206c4762a1bSJed Brown      Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
207c4762a1bSJed Brown   */
2089566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /*
211c4762a1bSJed Brown      Set array that saves the function norms.  This array is intended
212c4762a1bSJed Brown      when the user wants to save the convergence history for later use
213c4762a1bSJed Brown      rather than just to view the function norms via -snes_monitor.
214c4762a1bSJed Brown   */
2159566063dSJacob Faibussowitsch   PetscCall(SNESSetConvergenceHistory(snes,history,hist_its,50,PETSC_TRUE));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /*
218c4762a1bSJed Brown       Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the
219c4762a1bSJed Brown       user provided test before the specialized test. The convergence context is just an array to
220c4762a1bSJed Brown       test that it gets properly freed at the end
221c4762a1bSJed Brown   */
222c4762a1bSJed Brown   if (use_convergence_test) {
2239566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes,&ksp));
2249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5,&testarray));
2259566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp,ConvergenceTest,testarray,ConvergenceDestroy));
226c4762a1bSJed Brown   }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
230c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231c4762a1bSJed Brown   /*
232c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
233c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
234c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
235c4762a1bSJed Brown      this vector to zero by calling VecSet().
236c4762a1bSJed Brown   */
2379566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user,x));
2389566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,x));
2399566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
240*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /*
243c4762a1bSJed Brown      Print the convergence history.  This is intended just to demonstrate
244c4762a1bSJed Brown      use of the data attained via SNESSetConvergenceHistory().
245c4762a1bSJed Brown   */
2469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-print_history",&flg));
247c4762a1bSJed Brown   if (flg) {
248c4762a1bSJed Brown     for (i=0; i<its+1; i++) {
249*63a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n",i,hist_its[i],(double)history[i]));
250c4762a1bSJed Brown     }
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
254c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
255c4762a1bSJed Brown      are no longer needed.
256c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   if (!matrix_free) {
2599566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&J));
260c4762a1bSJed Brown   }
261c4762a1bSJed Brown   if (fd_coloring) {
2629566063dSJacob Faibussowitsch     PetscCall(MatFDColoringDestroy(&fdcoloring));
263c4762a1bSJed Brown   }
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2669566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2679566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
268b122ec5aSJacob Faibussowitsch   return 0;
269c4762a1bSJed Brown }
270c4762a1bSJed Brown /* ------------------------------------------------------------------- */
271c4762a1bSJed Brown /*
272c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
273c4762a1bSJed Brown 
274c4762a1bSJed Brown    Input Parameters:
275c4762a1bSJed Brown    user - user-defined application context
276c4762a1bSJed Brown    X - vector
277c4762a1bSJed Brown 
278c4762a1bSJed Brown    Output Parameter:
279c4762a1bSJed Brown    X - vector
280c4762a1bSJed Brown  */
281c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,Vec X)
282c4762a1bSJed Brown {
283c4762a1bSJed Brown   PetscInt       i,j,row,mx,my;
284c4762a1bSJed Brown   PetscReal      lambda,temp1,temp,hx,hy;
285c4762a1bSJed Brown   PetscScalar    *x;
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   mx     = user->mx;
288c4762a1bSJed Brown   my     = user->my;
289c4762a1bSJed Brown   lambda = user->param;
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx-1);
292c4762a1bSJed Brown   hy = 1.0 / (PetscReal)(my-1);
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   /*
295c4762a1bSJed Brown      Get a pointer to vector data.
296c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
297c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
298c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
299c4762a1bSJed Brown          the array.
300c4762a1bSJed Brown   */
3019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
302c4762a1bSJed Brown   temp1 = lambda/(lambda + 1.0);
303c4762a1bSJed Brown   for (j=0; j<my; j++) {
304c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
305c4762a1bSJed Brown     for (i=0; i<mx; i++) {
306c4762a1bSJed Brown       row = i + j*mx;
307c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
308c4762a1bSJed Brown         x[row] = 0.0;
309c4762a1bSJed Brown         continue;
310c4762a1bSJed Brown       }
311c4762a1bSJed Brown       x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
312c4762a1bSJed Brown     }
313c4762a1bSJed Brown   }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   /*
316c4762a1bSJed Brown      Restore vector
317c4762a1bSJed Brown   */
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
319c4762a1bSJed Brown   return 0;
320c4762a1bSJed Brown }
321c4762a1bSJed Brown /* ------------------------------------------------------------------- */
322c4762a1bSJed Brown /*
323c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
324c4762a1bSJed Brown 
325c4762a1bSJed Brown    Input Parameters:
326c4762a1bSJed Brown .  snes - the SNES context
327c4762a1bSJed Brown .  X - input vector
328c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
329c4762a1bSJed Brown 
330c4762a1bSJed Brown    Output Parameter:
331c4762a1bSJed Brown .  F - function vector
332c4762a1bSJed Brown  */
333c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
334c4762a1bSJed Brown {
335c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
336c4762a1bSJed Brown   PetscInt          i,j,row,mx,my;
337c4762a1bSJed Brown   PetscReal         two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx;
338c4762a1bSJed Brown   PetscScalar       ut,ub,ul,ur,u,uxx,uyy,sc,*f;
339c4762a1bSJed Brown   const PetscScalar *x;
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   mx     = user->mx;
342c4762a1bSJed Brown   my     = user->my;
343c4762a1bSJed Brown   lambda = user->param;
344c4762a1bSJed Brown   hx     = one / (PetscReal)(mx-1);
345c4762a1bSJed Brown   hy     = one / (PetscReal)(my-1);
346c4762a1bSJed Brown   sc     = hx*hy;
347c4762a1bSJed Brown   hxdhy  = hx/hy;
348c4762a1bSJed Brown   hydhx  = hy/hx;
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   /*
351c4762a1bSJed Brown      Get pointers to vector data
352c4762a1bSJed Brown   */
3539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
3549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   /*
357c4762a1bSJed Brown      Compute function
358c4762a1bSJed Brown   */
359c4762a1bSJed Brown   for (j=0; j<my; j++) {
360c4762a1bSJed Brown     for (i=0; i<mx; i++) {
361c4762a1bSJed Brown       row = i + j*mx;
362c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
363c4762a1bSJed Brown         f[row] = x[row];
364c4762a1bSJed Brown         continue;
365c4762a1bSJed Brown       }
366c4762a1bSJed Brown       u      = x[row];
367c4762a1bSJed Brown       ub     = x[row - mx];
368c4762a1bSJed Brown       ul     = x[row - 1];
369c4762a1bSJed Brown       ut     = x[row + mx];
370c4762a1bSJed Brown       ur     = x[row + 1];
371c4762a1bSJed Brown       uxx    = (-ur + two*u - ul)*hydhx;
372c4762a1bSJed Brown       uyy    = (-ut + two*u - ub)*hxdhy;
373c4762a1bSJed Brown       f[row] = uxx + uyy - sc*lambda*PetscExpScalar(u);
374c4762a1bSJed Brown     }
375c4762a1bSJed Brown   }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   /*
378c4762a1bSJed Brown      Restore vectors
379c4762a1bSJed Brown   */
3809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
3819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
382c4762a1bSJed Brown   return 0;
383c4762a1bSJed Brown }
384c4762a1bSJed Brown /* ------------------------------------------------------------------- */
385c4762a1bSJed Brown /*
386c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
387c4762a1bSJed Brown 
388c4762a1bSJed Brown    Input Parameters:
389c4762a1bSJed Brown .  snes - the SNES context
390c4762a1bSJed Brown .  x - input vector
391c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetJacobian()
392c4762a1bSJed Brown 
393c4762a1bSJed Brown    Output Parameters:
394c4762a1bSJed Brown .  A - Jacobian matrix
395c4762a1bSJed Brown .  B - optionally different preconditioning matrix
396c4762a1bSJed Brown .  flag - flag indicating matrix structure
397c4762a1bSJed Brown */
398c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
399c4762a1bSJed Brown {
400c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;   /* user-defined applicatin context */
401c4762a1bSJed Brown   PetscInt          i,j,row,mx,my,col[5];
402c4762a1bSJed Brown   PetscScalar       two = 2.0,one = 1.0,lambda,v[5],sc;
403c4762a1bSJed Brown   const PetscScalar *x;
404c4762a1bSJed Brown   PetscReal         hx,hy,hxdhy,hydhx;
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   mx     = user->mx;
407c4762a1bSJed Brown   my     = user->my;
408c4762a1bSJed Brown   lambda = user->param;
409c4762a1bSJed Brown   hx     = 1.0 / (PetscReal)(mx-1);
410c4762a1bSJed Brown   hy     = 1.0 / (PetscReal)(my-1);
411c4762a1bSJed Brown   sc     = hx*hy;
412c4762a1bSJed Brown   hxdhy  = hx/hy;
413c4762a1bSJed Brown   hydhx  = hy/hx;
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   /*
416c4762a1bSJed Brown      Get pointer to vector data
417c4762a1bSJed Brown   */
4189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   /*
421c4762a1bSJed Brown      Compute entries of the Jacobian
422c4762a1bSJed Brown   */
423c4762a1bSJed Brown   for (j=0; j<my; j++) {
424c4762a1bSJed Brown     for (i=0; i<mx; i++) {
425c4762a1bSJed Brown       row = i + j*mx;
426c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
4279566063dSJacob Faibussowitsch         PetscCall(MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES));
428c4762a1bSJed Brown         continue;
429c4762a1bSJed Brown       }
430c4762a1bSJed Brown       v[0] = -hxdhy; col[0] = row - mx;
431c4762a1bSJed Brown       v[1] = -hydhx; col[1] = row - 1;
432c4762a1bSJed Brown       v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row;
433c4762a1bSJed Brown       v[3] = -hydhx; col[3] = row + 1;
434c4762a1bSJed Brown       v[4] = -hxdhy; col[4] = row + mx;
4359566063dSJacob Faibussowitsch       PetscCall(MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES));
436c4762a1bSJed Brown     }
437c4762a1bSJed Brown   }
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /*
440c4762a1bSJed Brown      Restore vector
441c4762a1bSJed Brown   */
4429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   /*
445c4762a1bSJed Brown      Assemble matrix
446c4762a1bSJed Brown   */
4479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
4489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   if (jac != J) {
4519566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
4529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
453c4762a1bSJed Brown   }
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   return 0;
456c4762a1bSJed Brown }
457c4762a1bSJed Brown 
458c4762a1bSJed Brown PetscErrorCode ConvergenceTest(KSP ksp,PetscInt it,PetscReal nrm,KSPConvergedReason *reason,void *ctx)
459c4762a1bSJed Brown {
460c4762a1bSJed Brown   PetscFunctionBegin;
461c4762a1bSJed Brown   *reason = KSP_CONVERGED_ITERATING;
462c4762a1bSJed Brown   if (it > 1) {
463c4762a1bSJed Brown     *reason = KSP_CONVERGED_ITS;
4649566063dSJacob Faibussowitsch     PetscCall(PetscInfo(NULL,"User provided convergence test returning after 2 iterations\n"));
465c4762a1bSJed Brown   }
466c4762a1bSJed Brown   PetscFunctionReturn(0);
467c4762a1bSJed Brown }
468c4762a1bSJed Brown 
469c4762a1bSJed Brown PetscErrorCode ConvergenceDestroy(void* ctx)
470c4762a1bSJed Brown {
471c4762a1bSJed Brown   PetscFunctionBegin;
4729566063dSJacob Faibussowitsch   PetscCall(PetscInfo(NULL,"User provided convergence destroy called\n"));
4739566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
474c4762a1bSJed Brown   PetscFunctionReturn(0);
475c4762a1bSJed Brown }
476c4762a1bSJed Brown 
477c4762a1bSJed Brown PetscErrorCode postcheck(SNES snes,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx)
478c4762a1bSJed Brown {
479c4762a1bSJed Brown   PetscReal      norm;
480c4762a1bSJed Brown   Vec            tmp;
481c4762a1bSJed Brown 
482c4762a1bSJed Brown   PetscFunctionBegin;
4839566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&tmp));
4849566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tmp,-1.0,x,w));
4859566063dSJacob Faibussowitsch   PetscCall(VecNorm(tmp,NORM_2,&norm));
4869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp));
4879566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of search step %g\n",(double)norm));
488c4762a1bSJed Brown   PetscFunctionReturn(0);
489c4762a1bSJed Brown }
490c4762a1bSJed Brown 
491c4762a1bSJed Brown /*TEST
492c4762a1bSJed Brown 
493c4762a1bSJed Brown    build:
494c4762a1bSJed Brown       requires: !single
495c4762a1bSJed Brown 
496c4762a1bSJed Brown    test:
497c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always
498c4762a1bSJed Brown 
499c4762a1bSJed Brown    test:
500c4762a1bSJed Brown       suffix: 2
501c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always
502c4762a1bSJed Brown 
503c4762a1bSJed Brown    test:
504c4762a1bSJed Brown       suffix: 2a
505c4762a1bSJed Brown       filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault"
506c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info
507dfd57a17SPierre Jolivet       requires: defined(PETSC_USE_INFO)
508c4762a1bSJed Brown 
509c4762a1bSJed Brown    test:
510c4762a1bSJed Brown       suffix: 2b
511c4762a1bSJed Brown       filter: grep -i  "User provided convergence test" > /dev/null  && echo "Found User provided convergence test"
512c4762a1bSJed Brown       args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info
513dfd57a17SPierre Jolivet       requires: defined(PETSC_USE_INFO)
514c4762a1bSJed Brown 
515c4762a1bSJed Brown    test:
516c4762a1bSJed Brown       suffix: 3
517c4762a1bSJed Brown       args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always
518c4762a1bSJed Brown 
519c4762a1bSJed Brown    test:
520c4762a1bSJed Brown       suffix: 4
521c4762a1bSJed Brown       args: -pc -par 6.807 -snes_monitor -snes_converged_reason
522c4762a1bSJed Brown 
523c4762a1bSJed Brown TEST*/
524