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