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