xref: /petsc/src/snes/tutorials/ex3.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel.\n\
2*c4762a1bSJed Brown This example employs a user-defined monitoring routine and optionally a user-defined\n\
3*c4762a1bSJed Brown routine to check candidate iterates produced by line search routines.\n\
4*c4762a1bSJed Brown The command line options include:\n\
5*c4762a1bSJed Brown   -pre_check_iterates : activate checking of iterates\n\
6*c4762a1bSJed Brown   -post_check_iterates : activate checking of iterates\n\
7*c4762a1bSJed Brown   -check_tol <tol>: set tolerance for iterate checking\n\
8*c4762a1bSJed Brown   -user_precond : activate a (trivial) user-defined preconditioner\n\n";
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown /*T
11*c4762a1bSJed Brown    Concepts: SNES^basic parallel example
12*c4762a1bSJed Brown    Concepts: SNES^setting a user-defined monitoring routine
13*c4762a1bSJed Brown    Processors: n
14*c4762a1bSJed Brown T*/
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown /*
19*c4762a1bSJed Brown    Include "petscdm.h" so that we can use data management objects (DMs)
20*c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
21*c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
22*c4762a1bSJed Brown    file automatically includes:
23*c4762a1bSJed Brown      petscsys.h    - base PETSc routines
24*c4762a1bSJed Brown      petscvec.h    - vectors
25*c4762a1bSJed Brown      petscmat.h    - matrices
26*c4762a1bSJed Brown      petscis.h     - index sets
27*c4762a1bSJed Brown      petscksp.h    - Krylov subspace methods
28*c4762a1bSJed Brown      petscviewer.h - viewers
29*c4762a1bSJed Brown      petscpc.h     - preconditioners
30*c4762a1bSJed Brown      petscksp.h    - linear solvers
31*c4762a1bSJed Brown */
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown #include <petscdm.h>
34*c4762a1bSJed Brown #include <petscdmda.h>
35*c4762a1bSJed Brown #include <petscsnes.h>
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown /*
38*c4762a1bSJed Brown    User-defined routines.
39*c4762a1bSJed Brown */
40*c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
41*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
42*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec);
43*c4762a1bSJed Brown PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*);
44*c4762a1bSJed Brown PetscErrorCode PreCheck(SNESLineSearch,Vec,Vec,PetscBool*,void*);
45*c4762a1bSJed Brown PetscErrorCode PostCheck(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
46*c4762a1bSJed Brown PetscErrorCode PostSetSubKSP(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
47*c4762a1bSJed Brown PetscErrorCode MatrixFreePreconditioner(PC,Vec,Vec);
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown /*
50*c4762a1bSJed Brown    User-defined application context
51*c4762a1bSJed Brown */
52*c4762a1bSJed Brown typedef struct {
53*c4762a1bSJed Brown   DM          da;      /* distributed array */
54*c4762a1bSJed Brown   Vec         F;       /* right-hand-side of PDE */
55*c4762a1bSJed Brown   PetscMPIInt rank;    /* rank of processor */
56*c4762a1bSJed Brown   PetscMPIInt size;    /* size of communicator */
57*c4762a1bSJed Brown   PetscReal   h;       /* mesh spacing */
58*c4762a1bSJed Brown   PetscBool   sjerr;   /* if or not to test jacobian domain error */
59*c4762a1bSJed Brown } ApplicationCtx;
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown /*
62*c4762a1bSJed Brown    User-defined context for monitoring
63*c4762a1bSJed Brown */
64*c4762a1bSJed Brown typedef struct {
65*c4762a1bSJed Brown   PetscViewer viewer;
66*c4762a1bSJed Brown } MonitorCtx;
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown /*
69*c4762a1bSJed Brown    User-defined context for checking candidate iterates that are
70*c4762a1bSJed Brown    determined by line search methods
71*c4762a1bSJed Brown */
72*c4762a1bSJed Brown typedef struct {
73*c4762a1bSJed Brown   Vec            last_step;  /* previous iterate */
74*c4762a1bSJed Brown   PetscReal      tolerance;  /* tolerance for changes between successive iterates */
75*c4762a1bSJed Brown   ApplicationCtx *user;
76*c4762a1bSJed Brown } StepCheckCtx;
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown typedef struct {
79*c4762a1bSJed Brown   PetscInt its0; /* num of prevous outer KSP iterations */
80*c4762a1bSJed Brown } SetSubKSPCtx;
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown int main(int argc,char **argv)
83*c4762a1bSJed Brown {
84*c4762a1bSJed Brown   SNES           snes;                 /* SNES context */
85*c4762a1bSJed Brown   SNESLineSearch linesearch;           /* SNESLineSearch context */
86*c4762a1bSJed Brown   Mat            J;                    /* Jacobian matrix */
87*c4762a1bSJed Brown   ApplicationCtx ctx;                  /* user-defined context */
88*c4762a1bSJed Brown   Vec            x,r,U,F;              /* vectors */
89*c4762a1bSJed Brown   MonitorCtx     monP;                 /* monitoring context */
90*c4762a1bSJed Brown   StepCheckCtx   checkP;               /* step-checking context */
91*c4762a1bSJed Brown   SetSubKSPCtx   checkP1;
92*c4762a1bSJed Brown   PetscBool      pre_check,post_check,post_setsubksp; /* flag indicating whether we're checking candidate iterates */
93*c4762a1bSJed Brown   PetscScalar    xp,*FF,*UU,none = -1.0;
94*c4762a1bSJed Brown   PetscErrorCode ierr;
95*c4762a1bSJed Brown   PetscInt       its,N = 5,i,maxit,maxf,xs,xm;
96*c4762a1bSJed Brown   PetscReal      abstol,rtol,stol,norm;
97*c4762a1bSJed Brown   PetscBool      flg;
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
101*c4762a1bSJed Brown   ierr  = MPI_Comm_rank(PETSC_COMM_WORLD,&ctx.rank);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr  = MPI_Comm_size(PETSC_COMM_WORLD,&ctx.size);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr  = PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ctx.h = 1.0/(N-1);
105*c4762a1bSJed Brown   ctx.sjerr = PETSC_FALSE;
106*c4762a1bSJed Brown   ierr  = PetscOptionsGetBool(NULL,NULL,"-test_jacobian_domain_error",&ctx.sjerr,NULL);CHKERRQ(ierr);
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109*c4762a1bSJed Brown      Create nonlinear solver context
110*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115*c4762a1bSJed Brown      Create vector data structures; set function evaluation routine
116*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown   /*
119*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
120*c4762a1bSJed Brown   */
121*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = DMSetFromOptions(ctx.da);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = DMSetUp(ctx.da);CHKERRQ(ierr);
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   /*
126*c4762a1bSJed Brown      Extract global and local vectors from DMDA; then duplicate for remaining
127*c4762a1bSJed Brown      vectors that are the same types
128*c4762a1bSJed Brown   */
129*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(ctx.da,&x);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = VecDuplicate(x,&F);CHKERRQ(ierr); ctx.F = F;
132*c4762a1bSJed Brown   ierr = VecDuplicate(x,&U);CHKERRQ(ierr);
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown   /*
135*c4762a1bSJed Brown      Set function evaluation routine and vector.  Whenever the nonlinear
136*c4762a1bSJed Brown      solver needs to compute the nonlinear function, it will call this
137*c4762a1bSJed Brown      routine.
138*c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
139*c4762a1bSJed Brown         context that provides application-specific data for the
140*c4762a1bSJed Brown         function evaluation routine.
141*c4762a1bSJed Brown   */
142*c4762a1bSJed Brown   ierr = SNESSetFunction(snes,r,FormFunction,&ctx);CHKERRQ(ierr);
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145*c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
146*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(J,3,NULL);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(J,3,NULL,3,NULL);CHKERRQ(ierr);
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   /*
155*c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
156*c4762a1bSJed Brown      routine.  Whenever the nonlinear solver needs to compute the
157*c4762a1bSJed Brown      Jacobian matrix, it will call this routine.
158*c4762a1bSJed Brown       - Note that the final routine argument is the user-defined
159*c4762a1bSJed Brown         context that provides application-specific data for the
160*c4762a1bSJed Brown         Jacobian evaluation routine.
161*c4762a1bSJed Brown   */
162*c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,J,J,FormJacobian,&ctx);CHKERRQ(ierr);
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown   /*
165*c4762a1bSJed Brown      Optionally allow user-provided preconditioner
166*c4762a1bSJed Brown    */
167*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-user_precond",&flg);CHKERRQ(ierr);
168*c4762a1bSJed Brown   if (flg) {
169*c4762a1bSJed Brown     KSP ksp;
170*c4762a1bSJed Brown     PC  pc;
171*c4762a1bSJed Brown     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
172*c4762a1bSJed Brown     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
173*c4762a1bSJed Brown     ierr = PCSetType(pc,PCSHELL);CHKERRQ(ierr);
174*c4762a1bSJed Brown     ierr = PCShellSetApply(pc,MatrixFreePreconditioner);CHKERRQ(ierr);
175*c4762a1bSJed Brown   }
176*c4762a1bSJed Brown 
177*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178*c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
179*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown   /*
182*c4762a1bSJed Brown      Set an optional user-defined monitoring routine
183*c4762a1bSJed Brown   */
184*c4762a1bSJed Brown   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = SNESMonitorSet(snes,Monitor,&monP,0);CHKERRQ(ierr);
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown   /*
188*c4762a1bSJed Brown      Set names for some vectors to facilitate monitoring (optional)
189*c4762a1bSJed Brown   */
190*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)x,"Approximate Solution");CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)U,"Exact Solution");CHKERRQ(ierr);
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown   /*
194*c4762a1bSJed Brown      Set SNES/KSP/KSP/PC runtime options, e.g.,
195*c4762a1bSJed Brown          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
196*c4762a1bSJed Brown   */
197*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown   /*
200*c4762a1bSJed Brown      Set an optional user-defined routine to check the validity of candidate
201*c4762a1bSJed Brown      iterates that are determined by line search methods
202*c4762a1bSJed Brown   */
203*c4762a1bSJed Brown   ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-post_check_iterates",&post_check);CHKERRQ(ierr);
205*c4762a1bSJed Brown 
206*c4762a1bSJed Brown   if (post_check) {
207*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Activating post step checking routine\n");CHKERRQ(ierr);
208*c4762a1bSJed Brown     ierr = SNESLineSearchSetPostCheck(linesearch,PostCheck,&checkP);CHKERRQ(ierr);
209*c4762a1bSJed Brown     ierr = VecDuplicate(x,&(checkP.last_step));CHKERRQ(ierr);
210*c4762a1bSJed Brown 
211*c4762a1bSJed Brown     checkP.tolerance = 1.0;
212*c4762a1bSJed Brown     checkP.user      = &ctx;
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown     ierr = PetscOptionsGetReal(NULL,NULL,"-check_tol",&checkP.tolerance,NULL);CHKERRQ(ierr);
215*c4762a1bSJed Brown   }
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-post_setsubksp",&post_setsubksp);CHKERRQ(ierr);
218*c4762a1bSJed Brown   if (post_setsubksp) {
219*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Activating post setsubksp\n");CHKERRQ(ierr);
220*c4762a1bSJed Brown     ierr = SNESLineSearchSetPostCheck(linesearch,PostSetSubKSP,&checkP1);CHKERRQ(ierr);
221*c4762a1bSJed Brown   }
222*c4762a1bSJed Brown 
223*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-pre_check_iterates",&pre_check);CHKERRQ(ierr);
224*c4762a1bSJed Brown   if (pre_check) {
225*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Activating pre step checking routine\n");CHKERRQ(ierr);
226*c4762a1bSJed Brown     ierr = SNESLineSearchSetPreCheck(linesearch,PreCheck,&checkP);CHKERRQ(ierr);
227*c4762a1bSJed Brown   }
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   /*
230*c4762a1bSJed Brown      Print parameters used for convergence testing (optional) ... just
231*c4762a1bSJed Brown      to demonstrate this routine; this information is also printed with
232*c4762a1bSJed Brown      the option -snes_view
233*c4762a1bSJed Brown   */
234*c4762a1bSJed Brown   ierr = SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"atol=%g, rtol=%g, stol=%g, maxit=%D, maxf=%D\n",(double)abstol,(double)rtol,(double)stol,maxit,maxf);CHKERRQ(ierr);
236*c4762a1bSJed Brown 
237*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238*c4762a1bSJed Brown      Initialize application:
239*c4762a1bSJed Brown      Store right-hand-side of PDE and exact solution
240*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
241*c4762a1bSJed Brown 
242*c4762a1bSJed Brown   /*
243*c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
244*c4762a1bSJed Brown        xs, xm - starting grid index, width of local grid (no ghost points)
245*c4762a1bSJed Brown   */
246*c4762a1bSJed Brown   ierr = DMDAGetCorners(ctx.da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
247*c4762a1bSJed Brown 
248*c4762a1bSJed Brown   /*
249*c4762a1bSJed Brown      Get pointers to vector data
250*c4762a1bSJed Brown   */
251*c4762a1bSJed Brown   ierr = DMDAVecGetArray(ctx.da,F,&FF);CHKERRQ(ierr);
252*c4762a1bSJed Brown   ierr = DMDAVecGetArray(ctx.da,U,&UU);CHKERRQ(ierr);
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown   /*
255*c4762a1bSJed Brown      Compute local vector entries
256*c4762a1bSJed Brown   */
257*c4762a1bSJed Brown   xp = ctx.h*xs;
258*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
259*c4762a1bSJed Brown     FF[i] = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */
260*c4762a1bSJed Brown     UU[i] = xp*xp*xp;
261*c4762a1bSJed Brown     xp   += ctx.h;
262*c4762a1bSJed Brown   }
263*c4762a1bSJed Brown 
264*c4762a1bSJed Brown   /*
265*c4762a1bSJed Brown      Restore vectors
266*c4762a1bSJed Brown   */
267*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(ctx.da,F,&FF);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(ctx.da,U,&UU);CHKERRQ(ierr);
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271*c4762a1bSJed Brown      Evaluate initial guess; then solve nonlinear system
272*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
273*c4762a1bSJed Brown 
274*c4762a1bSJed Brown   /*
275*c4762a1bSJed Brown      Note: The user should initialize the vector, x, with the initial guess
276*c4762a1bSJed Brown      for the nonlinear solver prior to calling SNESSolve().  In particular,
277*c4762a1bSJed Brown      to employ an initial guess of zero, the user should explicitly set
278*c4762a1bSJed Brown      this vector to zero by calling VecSet().
279*c4762a1bSJed Brown   */
280*c4762a1bSJed Brown   ierr = FormInitialGuess(x);CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
282*c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
283*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr);
284*c4762a1bSJed Brown 
285*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
286*c4762a1bSJed Brown      Check solution and clean up
287*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown   /*
290*c4762a1bSJed Brown      Check the error
291*c4762a1bSJed Brown   */
292*c4762a1bSJed Brown   ierr = VecAXPY(x,none,U);CHKERRQ(ierr);
293*c4762a1bSJed Brown   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
294*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)norm,its);CHKERRQ(ierr);
295*c4762a1bSJed Brown   if (ctx.sjerr) {
296*c4762a1bSJed Brown     SNESType       snestype;
297*c4762a1bSJed Brown     ierr = SNESGetType(snes,&snestype);CHKERRQ(ierr);
298*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"SNES Type %s\n",snestype);CHKERRQ(ierr);
299*c4762a1bSJed Brown   }
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown   /*
302*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
303*c4762a1bSJed Brown      are no longer needed.
304*c4762a1bSJed Brown   */
305*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&monP.viewer);CHKERRQ(ierr);
306*c4762a1bSJed Brown   if (post_check) {ierr = VecDestroy(&checkP.last_step);CHKERRQ(ierr);}
307*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
308*c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
309*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
310*c4762a1bSJed Brown   ierr = VecDestroy(&F);CHKERRQ(ierr);
311*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
312*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
313*c4762a1bSJed Brown   ierr = DMDestroy(&ctx.da);CHKERRQ(ierr);
314*c4762a1bSJed Brown   ierr = PetscFinalize();
315*c4762a1bSJed Brown   return ierr;
316*c4762a1bSJed Brown }
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
319*c4762a1bSJed Brown /*
320*c4762a1bSJed Brown    FormInitialGuess - Computes initial guess.
321*c4762a1bSJed Brown 
322*c4762a1bSJed Brown    Input/Output Parameter:
323*c4762a1bSJed Brown .  x - the solution vector
324*c4762a1bSJed Brown */
325*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec x)
326*c4762a1bSJed Brown {
327*c4762a1bSJed Brown   PetscErrorCode ierr;
328*c4762a1bSJed Brown   PetscScalar    pfive = .50;
329*c4762a1bSJed Brown 
330*c4762a1bSJed Brown   PetscFunctionBeginUser;
331*c4762a1bSJed Brown   ierr = VecSet(x,pfive);CHKERRQ(ierr);
332*c4762a1bSJed Brown   PetscFunctionReturn(0);
333*c4762a1bSJed Brown }
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
336*c4762a1bSJed Brown /*
337*c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
338*c4762a1bSJed Brown 
339*c4762a1bSJed Brown    Input Parameters:
340*c4762a1bSJed Brown .  snes - the SNES context
341*c4762a1bSJed Brown .  x - input vector
342*c4762a1bSJed Brown .  ctx - optional user-defined context, as set by SNESSetFunction()
343*c4762a1bSJed Brown 
344*c4762a1bSJed Brown    Output Parameter:
345*c4762a1bSJed Brown .  f - function vector
346*c4762a1bSJed Brown 
347*c4762a1bSJed Brown    Note:
348*c4762a1bSJed Brown    The user-defined context can contain any application-specific
349*c4762a1bSJed Brown    data needed for the function evaluation.
350*c4762a1bSJed Brown */
351*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *ctx)
352*c4762a1bSJed Brown {
353*c4762a1bSJed Brown   ApplicationCtx *user = (ApplicationCtx*) ctx;
354*c4762a1bSJed Brown   DM             da    = user->da;
355*c4762a1bSJed Brown   PetscScalar    *xx,*ff,*FF,d;
356*c4762a1bSJed Brown   PetscErrorCode ierr;
357*c4762a1bSJed Brown   PetscInt       i,M,xs,xm;
358*c4762a1bSJed Brown   Vec            xlocal;
359*c4762a1bSJed Brown 
360*c4762a1bSJed Brown   PetscFunctionBeginUser;
361*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&xlocal);CHKERRQ(ierr);
362*c4762a1bSJed Brown   /*
363*c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
364*c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
365*c4762a1bSJed Brown      By placing code between these two statements, computations can
366*c4762a1bSJed Brown      be done while messages are in transition.
367*c4762a1bSJed Brown   */
368*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
369*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
370*c4762a1bSJed Brown 
371*c4762a1bSJed Brown   /*
372*c4762a1bSJed Brown      Get pointers to vector data.
373*c4762a1bSJed Brown        - The vector xlocal includes ghost point; the vectors x and f do
374*c4762a1bSJed Brown          NOT include ghost points.
375*c4762a1bSJed Brown        - Using DMDAVecGetArray() allows accessing the values using global ordering
376*c4762a1bSJed Brown   */
377*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,xlocal,&xx);CHKERRQ(ierr);
378*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,f,&ff);CHKERRQ(ierr);
379*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,user->F,&FF);CHKERRQ(ierr);
380*c4762a1bSJed Brown 
381*c4762a1bSJed Brown   /*
382*c4762a1bSJed Brown      Get local grid boundaries (for 1-dimensional DMDA):
383*c4762a1bSJed Brown        xs, xm  - starting grid index, width of local grid (no ghost points)
384*c4762a1bSJed Brown   */
385*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
386*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
387*c4762a1bSJed Brown 
388*c4762a1bSJed Brown   /*
389*c4762a1bSJed Brown      Set function values for boundary points; define local interior grid point range:
390*c4762a1bSJed Brown         xsi - starting interior grid index
391*c4762a1bSJed Brown         xei - ending interior grid index
392*c4762a1bSJed Brown   */
393*c4762a1bSJed Brown   if (xs == 0) { /* left boundary */
394*c4762a1bSJed Brown     ff[0] = xx[0];
395*c4762a1bSJed Brown     xs++;xm--;
396*c4762a1bSJed Brown   }
397*c4762a1bSJed Brown   if (xs+xm == M) {  /* right boundary */
398*c4762a1bSJed Brown     ff[xs+xm-1] = xx[xs+xm-1] - 1.0;
399*c4762a1bSJed Brown     xm--;
400*c4762a1bSJed Brown   }
401*c4762a1bSJed Brown 
402*c4762a1bSJed Brown   /*
403*c4762a1bSJed Brown      Compute function over locally owned part of the grid (interior points only)
404*c4762a1bSJed Brown   */
405*c4762a1bSJed Brown   d = 1.0/(user->h*user->h);
406*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i];
407*c4762a1bSJed Brown 
408*c4762a1bSJed Brown   /*
409*c4762a1bSJed Brown      Restore vectors
410*c4762a1bSJed Brown   */
411*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,xlocal,&xx);CHKERRQ(ierr);
412*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,f,&ff);CHKERRQ(ierr);
413*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,user->F,&FF);CHKERRQ(ierr);
414*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&xlocal);CHKERRQ(ierr);
415*c4762a1bSJed Brown   PetscFunctionReturn(0);
416*c4762a1bSJed Brown }
417*c4762a1bSJed Brown 
418*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
419*c4762a1bSJed Brown /*
420*c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
421*c4762a1bSJed Brown 
422*c4762a1bSJed Brown    Input Parameters:
423*c4762a1bSJed Brown .  snes - the SNES context
424*c4762a1bSJed Brown .  x - input vector
425*c4762a1bSJed Brown .  dummy - optional user-defined context (not used here)
426*c4762a1bSJed Brown 
427*c4762a1bSJed Brown    Output Parameters:
428*c4762a1bSJed Brown .  jac - Jacobian matrix
429*c4762a1bSJed Brown .  B - optionally different preconditioning matrix
430*c4762a1bSJed Brown .  flag - flag indicating matrix structure
431*c4762a1bSJed Brown */
432*c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx)
433*c4762a1bSJed Brown {
434*c4762a1bSJed Brown   ApplicationCtx *user = (ApplicationCtx*) ctx;
435*c4762a1bSJed Brown   PetscScalar    *xx,d,A[3];
436*c4762a1bSJed Brown   PetscErrorCode ierr;
437*c4762a1bSJed Brown   PetscInt       i,j[3],M,xs,xm;
438*c4762a1bSJed Brown   DM             da = user->da;
439*c4762a1bSJed Brown 
440*c4762a1bSJed Brown   PetscFunctionBeginUser;
441*c4762a1bSJed Brown   if (user->sjerr) {
442*c4762a1bSJed Brown     ierr = SNESSetJacobianDomainError(snes);CHKERRQ(ierr);
443*c4762a1bSJed Brown     PetscFunctionReturn(0);
444*c4762a1bSJed Brown   }
445*c4762a1bSJed Brown   /*
446*c4762a1bSJed Brown      Get pointer to vector data
447*c4762a1bSJed Brown   */
448*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,x,&xx);CHKERRQ(ierr);
449*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
450*c4762a1bSJed Brown 
451*c4762a1bSJed Brown   /*
452*c4762a1bSJed Brown     Get range of locally owned matrix
453*c4762a1bSJed Brown   */
454*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
455*c4762a1bSJed Brown 
456*c4762a1bSJed Brown   /*
457*c4762a1bSJed Brown      Determine starting and ending local indices for interior grid points.
458*c4762a1bSJed Brown      Set Jacobian entries for boundary points.
459*c4762a1bSJed Brown   */
460*c4762a1bSJed Brown 
461*c4762a1bSJed Brown   if (xs == 0) {  /* left boundary */
462*c4762a1bSJed Brown     i = 0; A[0] = 1.0;
463*c4762a1bSJed Brown 
464*c4762a1bSJed Brown     ierr = MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr);
465*c4762a1bSJed Brown     xs++;xm--;
466*c4762a1bSJed Brown   }
467*c4762a1bSJed Brown   if (xs+xm == M) { /* right boundary */
468*c4762a1bSJed Brown     i    = M-1;
469*c4762a1bSJed Brown     A[0] = 1.0;
470*c4762a1bSJed Brown     ierr = MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES);CHKERRQ(ierr);
471*c4762a1bSJed Brown     xm--;
472*c4762a1bSJed Brown   }
473*c4762a1bSJed Brown 
474*c4762a1bSJed Brown   /*
475*c4762a1bSJed Brown      Interior grid points
476*c4762a1bSJed Brown       - Note that in this case we set all elements for a particular
477*c4762a1bSJed Brown         row at once.
478*c4762a1bSJed Brown   */
479*c4762a1bSJed Brown   d = 1.0/(user->h*user->h);
480*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
481*c4762a1bSJed Brown     j[0] = i - 1; j[1] = i; j[2] = i + 1;
482*c4762a1bSJed Brown     A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i];
483*c4762a1bSJed Brown     ierr = MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr);
484*c4762a1bSJed Brown   }
485*c4762a1bSJed Brown 
486*c4762a1bSJed Brown   /*
487*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
488*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
489*c4762a1bSJed Brown      By placing code between these two statements, computations can be
490*c4762a1bSJed Brown      done while messages are in transition.
491*c4762a1bSJed Brown 
492*c4762a1bSJed Brown      Also, restore vector.
493*c4762a1bSJed Brown   */
494*c4762a1bSJed Brown 
495*c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
496*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,x,&xx);CHKERRQ(ierr);
497*c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
498*c4762a1bSJed Brown 
499*c4762a1bSJed Brown   PetscFunctionReturn(0);
500*c4762a1bSJed Brown }
501*c4762a1bSJed Brown 
502*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
503*c4762a1bSJed Brown /*
504*c4762a1bSJed Brown    Monitor - Optional user-defined monitoring routine that views the
505*c4762a1bSJed Brown    current iterate with an x-window plot. Set by SNESMonitorSet().
506*c4762a1bSJed Brown 
507*c4762a1bSJed Brown    Input Parameters:
508*c4762a1bSJed Brown    snes - the SNES context
509*c4762a1bSJed Brown    its - iteration number
510*c4762a1bSJed Brown    norm - 2-norm function value (may be estimated)
511*c4762a1bSJed Brown    ctx - optional user-defined context for private data for the
512*c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
513*c4762a1bSJed Brown 
514*c4762a1bSJed Brown    Note:
515*c4762a1bSJed Brown    See the manpage for PetscViewerDrawOpen() for useful runtime options,
516*c4762a1bSJed Brown    such as -nox to deactivate all x-window output.
517*c4762a1bSJed Brown  */
518*c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
519*c4762a1bSJed Brown {
520*c4762a1bSJed Brown   PetscErrorCode ierr;
521*c4762a1bSJed Brown   MonitorCtx     *monP = (MonitorCtx*) ctx;
522*c4762a1bSJed Brown   Vec            x;
523*c4762a1bSJed Brown 
524*c4762a1bSJed Brown   PetscFunctionBeginUser;
525*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"iter = %D,SNES Function norm %g\n",its,(double)fnorm);CHKERRQ(ierr);
526*c4762a1bSJed Brown   ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr);
527*c4762a1bSJed Brown   ierr = VecView(x,monP->viewer);CHKERRQ(ierr);
528*c4762a1bSJed Brown   PetscFunctionReturn(0);
529*c4762a1bSJed Brown }
530*c4762a1bSJed Brown 
531*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
532*c4762a1bSJed Brown /*
533*c4762a1bSJed Brown    PreCheck - Optional user-defined routine that checks the validity of
534*c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPreCheck().
535*c4762a1bSJed Brown 
536*c4762a1bSJed Brown    Input Parameters:
537*c4762a1bSJed Brown    snes - the SNES context
538*c4762a1bSJed Brown    xcurrent - current solution
539*c4762a1bSJed Brown    y - search direction and length
540*c4762a1bSJed Brown 
541*c4762a1bSJed Brown    Output Parameters:
542*c4762a1bSJed Brown    y         - proposed step (search direction and length) (possibly changed)
543*c4762a1bSJed Brown    changed_y - tells if the step has changed or not
544*c4762a1bSJed Brown  */
545*c4762a1bSJed Brown PetscErrorCode PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y, PetscBool *changed_y, void * ctx)
546*c4762a1bSJed Brown {
547*c4762a1bSJed Brown   PetscFunctionBeginUser;
548*c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
549*c4762a1bSJed Brown   PetscFunctionReturn(0);
550*c4762a1bSJed Brown }
551*c4762a1bSJed Brown 
552*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
553*c4762a1bSJed Brown /*
554*c4762a1bSJed Brown    PostCheck - Optional user-defined routine that checks the validity of
555*c4762a1bSJed Brown    candidate steps of a line search method.  Set by SNESLineSearchSetPostCheck().
556*c4762a1bSJed Brown 
557*c4762a1bSJed Brown    Input Parameters:
558*c4762a1bSJed Brown    snes - the SNES context
559*c4762a1bSJed Brown    ctx  - optional user-defined context for private data for the
560*c4762a1bSJed Brown           monitor routine, as set by SNESLineSearchSetPostCheck()
561*c4762a1bSJed Brown    xcurrent - current solution
562*c4762a1bSJed Brown    y - search direction and length
563*c4762a1bSJed Brown    x    - the new candidate iterate
564*c4762a1bSJed Brown 
565*c4762a1bSJed Brown    Output Parameters:
566*c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
567*c4762a1bSJed Brown    x    - current iterate (possibly modified)
568*c4762a1bSJed Brown 
569*c4762a1bSJed Brown  */
570*c4762a1bSJed Brown PetscErrorCode PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
571*c4762a1bSJed Brown {
572*c4762a1bSJed Brown   PetscErrorCode ierr;
573*c4762a1bSJed Brown   PetscInt       i,iter,xs,xm;
574*c4762a1bSJed Brown   StepCheckCtx   *check;
575*c4762a1bSJed Brown   ApplicationCtx *user;
576*c4762a1bSJed Brown   PetscScalar    *xa,*xa_last,tmp;
577*c4762a1bSJed Brown   PetscReal      rdiff;
578*c4762a1bSJed Brown   DM             da;
579*c4762a1bSJed Brown   SNES           snes;
580*c4762a1bSJed Brown 
581*c4762a1bSJed Brown   PetscFunctionBeginUser;
582*c4762a1bSJed Brown   *changed_x = PETSC_FALSE;
583*c4762a1bSJed Brown   *changed_y = PETSC_FALSE;
584*c4762a1bSJed Brown 
585*c4762a1bSJed Brown   ierr  = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
586*c4762a1bSJed Brown   check = (StepCheckCtx*)ctx;
587*c4762a1bSJed Brown   user  = check->user;
588*c4762a1bSJed Brown   ierr  = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
589*c4762a1bSJed Brown 
590*c4762a1bSJed Brown   /* iteration 1 indicates we are working on the second iteration */
591*c4762a1bSJed Brown   if (iter > 0) {
592*c4762a1bSJed Brown     da   = user->da;
593*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %D with tolerance %g\n",iter,(double)check->tolerance);CHKERRQ(ierr);
594*c4762a1bSJed Brown 
595*c4762a1bSJed Brown     /* Access local array data */
596*c4762a1bSJed Brown     ierr = DMDAVecGetArray(da,check->last_step,&xa_last);CHKERRQ(ierr);
597*c4762a1bSJed Brown     ierr = DMDAVecGetArray(da,x,&xa);CHKERRQ(ierr);
598*c4762a1bSJed Brown     ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
599*c4762a1bSJed Brown 
600*c4762a1bSJed Brown     /*
601*c4762a1bSJed Brown        If we fail the user-defined check for validity of the candidate iterate,
602*c4762a1bSJed Brown        then modify the iterate as we like.  (Note that the particular modification
603*c4762a1bSJed Brown        below is intended simply to demonstrate how to manipulate this data, not
604*c4762a1bSJed Brown        as a meaningful or appropriate choice.)
605*c4762a1bSJed Brown     */
606*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
607*c4762a1bSJed Brown       if (!PetscAbsScalar(xa[i])) rdiff = 2*check->tolerance;
608*c4762a1bSJed Brown       else rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]);
609*c4762a1bSJed Brown       if (rdiff > check->tolerance) {
610*c4762a1bSJed Brown         tmp        = xa[i];
611*c4762a1bSJed Brown         xa[i]      = .5*(xa[i] + xa_last[i]);
612*c4762a1bSJed Brown         *changed_x = PETSC_TRUE;
613*c4762a1bSJed Brown         ierr       = PetscPrintf(PETSC_COMM_WORLD,"  Altering entry %D: x=%g, x_last=%g, diff=%g, x_new=%g\n",i,(double)PetscAbsScalar(tmp),(double)PetscAbsScalar(xa_last[i]),(double)rdiff,(double)PetscAbsScalar(xa[i]));CHKERRQ(ierr);
614*c4762a1bSJed Brown       }
615*c4762a1bSJed Brown     }
616*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(da,check->last_step,&xa_last);CHKERRQ(ierr);
617*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(da,x,&xa);CHKERRQ(ierr);
618*c4762a1bSJed Brown   }
619*c4762a1bSJed Brown   ierr = VecCopy(x,check->last_step);CHKERRQ(ierr);
620*c4762a1bSJed Brown   PetscFunctionReturn(0);
621*c4762a1bSJed Brown }
622*c4762a1bSJed Brown 
623*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
624*c4762a1bSJed Brown /*
625*c4762a1bSJed Brown    PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used
626*c4762a1bSJed Brown    e.g,
627*c4762a1bSJed Brown      mpiexec -n 8 ./ex3 -nox -n 10000 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 4 -sub_ksp_type gmres -sub_ksp_max_it 3 -post_setsubksp -sub_ksp_rtol 1.e-16
628*c4762a1bSJed Brown    Set by SNESLineSearchSetPostCheck().
629*c4762a1bSJed Brown 
630*c4762a1bSJed Brown    Input Parameters:
631*c4762a1bSJed Brown    linesearch - the LineSearch context
632*c4762a1bSJed Brown    xcurrent - current solution
633*c4762a1bSJed Brown    y - search direction and length
634*c4762a1bSJed Brown    x    - the new candidate iterate
635*c4762a1bSJed Brown 
636*c4762a1bSJed Brown    Output Parameters:
637*c4762a1bSJed Brown    y    - proposed step (search direction and length) (possibly changed)
638*c4762a1bSJed Brown    x    - current iterate (possibly modified)
639*c4762a1bSJed Brown 
640*c4762a1bSJed Brown  */
641*c4762a1bSJed Brown PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool  *changed_y,PetscBool  *changed_x, void * ctx)
642*c4762a1bSJed Brown {
643*c4762a1bSJed Brown   PetscErrorCode ierr;
644*c4762a1bSJed Brown   SetSubKSPCtx   *check;
645*c4762a1bSJed Brown   PetscInt       iter,its,sub_its,maxit;
646*c4762a1bSJed Brown   KSP            ksp,sub_ksp,*sub_ksps;
647*c4762a1bSJed Brown   PC             pc;
648*c4762a1bSJed Brown   PetscReal      ksp_ratio;
649*c4762a1bSJed Brown   SNES           snes;
650*c4762a1bSJed Brown 
651*c4762a1bSJed Brown   PetscFunctionBeginUser;
652*c4762a1bSJed Brown   ierr    = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
653*c4762a1bSJed Brown   check   = (SetSubKSPCtx*)ctx;
654*c4762a1bSJed Brown   ierr    = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
655*c4762a1bSJed Brown   ierr    = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
656*c4762a1bSJed Brown   ierr    = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
657*c4762a1bSJed Brown   ierr    = PCBJacobiGetSubKSP(pc,NULL,NULL,&sub_ksps);CHKERRQ(ierr);
658*c4762a1bSJed Brown   sub_ksp = sub_ksps[0];
659*c4762a1bSJed Brown   ierr    = KSPGetIterationNumber(ksp,&its);CHKERRQ(ierr);      /* outer KSP iteration number */
660*c4762a1bSJed Brown   ierr    = KSPGetIterationNumber(sub_ksp,&sub_its);CHKERRQ(ierr); /* inner KSP iteration number */
661*c4762a1bSJed Brown 
662*c4762a1bSJed Brown   if (iter) {
663*c4762a1bSJed Brown     ierr      = PetscPrintf(PETSC_COMM_WORLD,"    ...PostCheck snes iteration %D, ksp_it %D %D, subksp_it %D\n",iter,check->its0,its,sub_its);CHKERRQ(ierr);
664*c4762a1bSJed Brown     ksp_ratio = ((PetscReal)(its))/check->its0;
665*c4762a1bSJed Brown     maxit     = (PetscInt)(ksp_ratio*sub_its + 0.5);
666*c4762a1bSJed Brown     if (maxit < 2) maxit = 2;
667*c4762a1bSJed Brown     ierr = KSPSetTolerances(sub_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit);CHKERRQ(ierr);
668*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"    ...ksp_ratio %g, new maxit %D\n\n",(double)ksp_ratio,maxit);CHKERRQ(ierr);
669*c4762a1bSJed Brown   }
670*c4762a1bSJed Brown   check->its0 = its; /* save current outer KSP iteration number */
671*c4762a1bSJed Brown   PetscFunctionReturn(0);
672*c4762a1bSJed Brown }
673*c4762a1bSJed Brown 
674*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
675*c4762a1bSJed Brown /*
676*c4762a1bSJed Brown    MatrixFreePreconditioner - This routine demonstrates the use of a
677*c4762a1bSJed Brown    user-provided preconditioner.  This code implements just the null
678*c4762a1bSJed Brown    preconditioner, which of course is not recommended for general use.
679*c4762a1bSJed Brown 
680*c4762a1bSJed Brown    Input Parameters:
681*c4762a1bSJed Brown +  pc - preconditioner
682*c4762a1bSJed Brown -  x - input vector
683*c4762a1bSJed Brown 
684*c4762a1bSJed Brown    Output Parameter:
685*c4762a1bSJed Brown .  y - preconditioned vector
686*c4762a1bSJed Brown */
687*c4762a1bSJed Brown PetscErrorCode MatrixFreePreconditioner(PC pc,Vec x,Vec y)
688*c4762a1bSJed Brown {
689*c4762a1bSJed Brown   PetscErrorCode ierr;
690*c4762a1bSJed Brown   ierr = VecCopy(x,y);CHKERRQ(ierr);
691*c4762a1bSJed Brown   return 0;
692*c4762a1bSJed Brown }
693*c4762a1bSJed Brown 
694*c4762a1bSJed Brown 
695*c4762a1bSJed Brown /*TEST
696*c4762a1bSJed Brown 
697*c4762a1bSJed Brown    test:
698*c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
699*c4762a1bSJed Brown 
700*c4762a1bSJed Brown    test:
701*c4762a1bSJed Brown       suffix: 2
702*c4762a1bSJed Brown       nsize: 3
703*c4762a1bSJed Brown       args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
704*c4762a1bSJed Brown 
705*c4762a1bSJed Brown    test:
706*c4762a1bSJed Brown       suffix: 3
707*c4762a1bSJed Brown       nsize: 2
708*c4762a1bSJed Brown       args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always
709*c4762a1bSJed Brown 
710*c4762a1bSJed Brown    test:
711*c4762a1bSJed Brown       suffix: 4
712*c4762a1bSJed Brown       args: -nox -pre_check_iterates -post_check_iterates
713*c4762a1bSJed Brown 
714*c4762a1bSJed Brown    test:
715*c4762a1bSJed Brown       suffix: 5
716*c4762a1bSJed Brown       requires: double !complex !single
717*c4762a1bSJed Brown       nsize: 2
718*c4762a1bSJed Brown       args: -nox -snes_test_jacobian  -snes_test_jacobian_view
719*c4762a1bSJed Brown 
720*c4762a1bSJed Brown    test:
721*c4762a1bSJed Brown       suffix: 6
722*c4762a1bSJed Brown       requires: double !complex !single
723*c4762a1bSJed Brown       nsize: 4
724*c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1
725*c4762a1bSJed Brown 
726*c4762a1bSJed Brown    test:
727*c4762a1bSJed Brown       suffix: 7
728*c4762a1bSJed Brown       requires: double !complex !single
729*c4762a1bSJed Brown       nsize: 4
730*c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1
731*c4762a1bSJed Brown 
732*c4762a1bSJed Brown    test:
733*c4762a1bSJed Brown       suffix: 8
734*c4762a1bSJed Brown       requires: double !complex !single
735*c4762a1bSJed Brown       nsize: 4
736*c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1
737*c4762a1bSJed Brown 
738*c4762a1bSJed Brown    test:
739*c4762a1bSJed Brown       suffix: 9
740*c4762a1bSJed Brown       requires: double !complex !single
741*c4762a1bSJed Brown       nsize: 4
742*c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1
743*c4762a1bSJed Brown 
744*c4762a1bSJed Brown    test:
745*c4762a1bSJed Brown       suffix: 10
746*c4762a1bSJed Brown       requires: double !complex !single
747*c4762a1bSJed Brown       nsize: 4
748*c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1
749*c4762a1bSJed Brown 
750*c4762a1bSJed Brown    test:
751*c4762a1bSJed Brown       suffix: 11
752*c4762a1bSJed Brown       requires: double !complex !single
753*c4762a1bSJed Brown       nsize: 4
754*c4762a1bSJed Brown       args: -test_jacobian_domain_error -snes_converged_reason -snes_type ms -snes_check_jacobian_domain_error 1
755*c4762a1bSJed Brown 
756*c4762a1bSJed Brown TEST*/
757