Actual source code: itfunc.c

  2: /*
  3:       Interface KSP routines that the user calls.
  4: */

  6: #include <private/kspimpl.h>   /*I "petscksp.h" I*/

 10: /*@
 11:    KSPComputeExtremeSingularValues - Computes the extreme singular values
 12:    for the preconditioned operator. Called after or during KSPSolve().

 14:    Not Collective

 16:    Input Parameter:
 17: .  ksp - iterative context obtained from KSPCreate()

 19:    Output Parameters:
 20: .  emin, emax - extreme singular values

 22:    Notes:
 23:    One must call KSPSetComputeSingularValues() before calling KSPSetUp() 
 24:    (or use the option -ksp_compute_eigenvalues) in order for this routine to work correctly.

 26:    Many users may just want to use the monitoring routine
 27:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
 28:    to print the extreme singular values at each iteration of the linear solve.

 30:    Level: advanced

 32: .keywords: KSP, compute, extreme, singular, values

 34: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeEigenvalues()
 35: @*/
 36: PetscErrorCode  KSPComputeExtremeSingularValues(KSP ksp,PetscReal *emax,PetscReal *emin)
 37: {

 44:   if (!ksp->calc_sings) SETERRQ(((PetscObject)ksp)->comm,4,"Singular values not requested before KSPSetUp()");

 46:   if (ksp->ops->computeextremesingularvalues) {
 47:     (*ksp->ops->computeextremesingularvalues)(ksp,emax,emin);
 48:   } else {
 49:     *emin = -1.0;
 50:     *emax = -1.0;
 51:   }
 52:   return(0);
 53: }

 57: /*@
 58:    KSPComputeEigenvalues - Computes the extreme eigenvalues for the
 59:    preconditioned operator. Called after or during KSPSolve().

 61:    Not Collective

 63:    Input Parameter:
 64: +  ksp - iterative context obtained from KSPCreate()
 65: -  n - size of arrays r and c. The number of eigenvalues computed (neig) will, in 
 66:        general, be less than this.

 68:    Output Parameters:
 69: +  r - real part of computed eigenvalues
 70: .  c - complex part of computed eigenvalues
 71: -  neig - number of eigenvalues computed (will be less than or equal to n)

 73:    Options Database Keys:
 74: +  -ksp_compute_eigenvalues - Prints eigenvalues to stdout
 75: -  -ksp_plot_eigenvalues - Plots eigenvalues in an x-window display

 77:    Notes:
 78:    The number of eigenvalues estimated depends on the size of the Krylov space
 79:    generated during the KSPSolve() ; for example, with 
 80:    CG it corresponds to the number of CG iterations, for GMRES it is the number 
 81:    of GMRES iterations SINCE the last restart. Any extra space in r[] and c[]
 82:    will be ignored.

 84:    KSPComputeEigenvalues() does not usually provide accurate estimates; it is
 85:    intended only for assistance in understanding the convergence of iterative 
 86:    methods, not for eigenanalysis. 

 88:    One must call KSPSetComputeEigenvalues() before calling KSPSetUp() 
 89:    in order for this routine to work correctly.

 91:    Many users may just want to use the monitoring routine
 92:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
 93:    to print the singular values at each iteration of the linear solve.

 95:    Level: advanced

 97: .keywords: KSP, compute, extreme, singular, values

 99: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeExtremeSingularValues()
100: @*/
101: PetscErrorCode  KSPComputeEigenvalues(KSP ksp,PetscInt n,PetscReal *r,PetscReal *c,PetscInt *neig)
102: {

110:   if (!ksp->calc_sings) SETERRQ(((PetscObject)ksp)->comm,4,"Eigenvalues not requested before KSPSetUp()");

112:   if (ksp->ops->computeeigenvalues) {
113:     (*ksp->ops->computeeigenvalues)(ksp,n,r,c,neig);
114:   } else {
115:     *neig = 0;
116:   }
117:   return(0);
118: }

122: /*@
123:    KSPSetUpOnBlocks - Sets up the preconditioner for each block in
124:    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 
125:    methods.

127:    Collective on KSP

129:    Input Parameter:
130: .  ksp - the KSP context

132:    Notes:
133:    KSPSetUpOnBlocks() is a routine that the user can optinally call for
134:    more precise profiling (via -log_summary) of the setup phase for these
135:    block preconditioners.  If the user does not call KSPSetUpOnBlocks(),
136:    it will automatically be called from within KSPSolve().
137:    
138:    Calling KSPSetUpOnBlocks() is the same as calling PCSetUpOnBlocks()
139:    on the PC context within the KSP context.

141:    Level: advanced

143: .keywords: KSP, setup, blocks

145: .seealso: PCSetUpOnBlocks(), KSPSetUp(), PCSetUp()
146: @*/
147: PetscErrorCode  KSPSetUpOnBlocks(KSP ksp)
148: {

153:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
154:   PCSetUpOnBlocks(ksp->pc);
155:   return(0);
156: }

160: /*@
161:    KSPSetUp - Sets up the internal data structures for the
162:    later use of an iterative solver.

164:    Collective on KSP

166:    Input Parameter:
167: .  ksp   - iterative context obtained from KSPCreate()

169:    Level: developer

171: .keywords: KSP, setup

173: .seealso: KSPCreate(), KSPSolve(), KSPDestroy()
174: @*/
175: PetscErrorCode  KSPSetUp(KSP ksp)
176: {
178:   PetscBool      ir = PETSC_FALSE,ig = PETSC_FALSE;
179:   Mat            A;
180:   MatStructure   stflg = SAME_NONZERO_PATTERN;

182:   /* PetscBool      im = PETSC_FALSE; */


187:   /* reset the convergence flag from the previous solves */
188:   ksp->reason = KSP_CONVERGED_ITERATING;

190:   if (!((PetscObject)ksp)->type_name){
191:     KSPSetType(ksp,KSPGMRES);
192:   }
193:   KSPSetUpNorms_Private(ksp);

195:   if (ksp->dmActive && !ksp->setupstage) {
196:     /* first time in so build matrix and vector data structures using DM */
197:     if (!ksp->vec_rhs) {DMCreateGlobalVector(ksp->dm,&ksp->vec_rhs);}
198:     if (!ksp->vec_sol) {DMCreateGlobalVector(ksp->dm,&ksp->vec_sol);}
199:     DMGetMatrix(ksp->dm,MATAIJ,&A);
200:     KSPSetOperators(ksp,A,A,stflg);
201:     PetscObjectDereference((PetscObject)A);
202:   }

204:   if (ksp->dmActive) {
205:     DMHasInitialGuess(ksp->dm,&ig);
206:     if (ig && ksp->setupstage != KSP_SETUP_NEWRHS) {
207:       DMComputeInitialGuess(ksp->dm,ksp->vec_sol);
208:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
209:     }
210:     DMHasFunction(ksp->dm,&ir);
211:     if (ir && ksp->setupstage != KSP_SETUP_NEWRHS) {
212:       DMComputeFunction(ksp->dm,PETSC_NULL,ksp->vec_rhs);
213:     }

215:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
216:       KSPGetOperators(ksp,&A,&A,PETSC_NULL);
217:       DMComputeJacobian(ksp->dm,PETSC_NULL,A,A,&stflg);
218:       KSPSetOperators(ksp,A,A,stflg);
219:     }
220:   }

222:   if (ksp->setupstage == KSP_SETUP_NEWRHS) return(0);
223:   PetscLogEventBegin(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);

225:   if (!ksp->setupstage) {
226:     (*ksp->ops->setup)(ksp);
227:   }

229:   /* scale the matrix if requested */
230:   if (ksp->dscale) {
231:     Mat         mat,pmat;
232:     PetscScalar *xx;
233:     PetscInt    i,n;
234:     PetscBool   zeroflag = PETSC_FALSE;
235:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
236:     PCGetOperators(ksp->pc,&mat,&pmat,PETSC_NULL);
237:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
238:       MatGetVecs(pmat,&ksp->diagonal,0);
239:     }
240:     MatGetDiagonal(pmat,ksp->diagonal);
241:     VecGetLocalSize(ksp->diagonal,&n);
242:     VecGetArray(ksp->diagonal,&xx);
243:     for (i=0; i<n; i++) {
244:       if (xx[i] != 0.0) xx[i] = 1.0/PetscSqrtReal(PetscAbsScalar(xx[i]));
245:       else {
246:         xx[i]     = 1.0;
247:         zeroflag  = PETSC_TRUE;
248:       }
249:     }
250:     VecRestoreArray(ksp->diagonal,&xx);
251:     if (zeroflag) {
252:       PetscInfo(ksp,"Zero detected in diagonal of matrix, using 1 at those locations\n");
253:     }
254:     MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
255:     if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
256:     ksp->dscalefix2 = PETSC_FALSE;
257:   }
258:   PetscLogEventEnd(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
259:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
260:   PCSetUp(ksp->pc);
261:   if (ksp->nullsp) {
262:     PetscBool  test = PETSC_FALSE;
263:     PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_test_null_space",&test,PETSC_NULL);
264:     if (test) {
265:       Mat mat;
266:       PCGetOperators(ksp->pc,&mat,PETSC_NULL,PETSC_NULL);
267:       MatNullSpaceTest(ksp->nullsp,mat,PETSC_NULL);
268:     }
269:   }
270:   ksp->setupstage = KSP_SETUP_NEWRHS;
271:   return(0);
272: }

276: /*@
277:    KSPSolve - Solves linear system.

279:    Collective on KSP

281:    Parameter:
282: +  ksp - iterative context obtained from KSPCreate()
283: .  b - the right hand side vector
284: -  x - the solution  (this may be the same vector as b, then b will be overwritten with answer)

286:    Options Database Keys:
287: +  -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
288: .  -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
289: .  -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and useing LAPACK
290: .  -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
291: .  -ksp_view_binary - save matrix and right hand side that define linear system to the default binary viewer (can be
292:                                 read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
293: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
294: .  -ksp_final_residual - print 2-norm of true linear system residual at the end of the solution process
295: -  -ksp_view - print the ksp data structure at the end of the system solution

297:    Notes:

299:    The operator is specified with PCSetOperators().

301:    Call KSPGetConvergedReason() to determine if the solver converged or failed and 
302:    why. The number of iterations can be obtained from KSPGetIterationNumber().
303:    
304:    If using a direct method (e.g., via the KSP solver
305:    KSPPREONLY and a preconditioner such as PCLU/PCILU),
306:    then its=1.  See KSPSetTolerances() and KSPDefaultConverged()
307:    for more details.

309:    Understanding Convergence:
310:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
311:    KSPComputeEigenvaluesExplicitly() provide information on additional
312:    options to monitor convergence and print eigenvalue information.

314:    Level: beginner

316: .keywords: KSP, solve, linear system

318: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPDefaultConverged(),
319:           KSPSolveTranspose(), KSPGetIterationNumber()
320: @*/
321: PetscErrorCode  KSPSolve(KSP ksp,Vec b,Vec x)
322: {
324:   PetscMPIInt    rank;
325:   PetscBool      flag1,flag2,flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
326:   char           view[10];
327:   char           filename[PETSC_MAX_PATH_LEN];
328:   PetscViewer    viewer;
329: 


336:   if (b && x) {
337:     if (x == b) {
338:       VecDuplicate(b,&x);
339:       inXisinB = PETSC_TRUE;
340:     }
341:     PetscObjectReference((PetscObject)b);
342:     PetscObjectReference((PetscObject)x);
343:     VecDestroy(&ksp->vec_rhs);
344:     VecDestroy(&ksp->vec_sol);
345:     ksp->vec_rhs = b;
346:     ksp->vec_sol = x;
347:   }

349:   flag1 = flag2 = PETSC_FALSE;
350:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_view_binary",&flag1,PETSC_NULL);
351:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_view_binary_pre",&flag2,PETSC_NULL);
352:   if (flag1 || flag2) {
353:     Mat         mat,premat;
354:     PetscViewer viewer = PETSC_VIEWER_BINARY_(((PetscObject)ksp)->comm);
355:     PCGetOperators(ksp->pc,&mat,&premat,PETSC_NULL);
356:     if (flag1) {MatView(mat,viewer);}
357:     if (flag2) {MatView(premat,viewer);}
358:     VecView(ksp->vec_rhs,viewer);
359:   }
360:   PetscLogEventBegin(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);

362:   /* reset the residual history list if requested */
363:   if (ksp->res_hist_reset) ksp->res_hist_len = 0;

365:   PetscOptionsGetString(((PetscObject)ksp)->prefix,"-ksp_view_before",view,10,&flg);
366:   if (flg) {
367:     PetscViewer viewer;
368:     PetscViewerASCIIGetStdout(((PetscObject)ksp)->comm,&viewer);
369:     KSPView(ksp,viewer);
370:   }

372:   ksp->transpose_solve = PETSC_FALSE;

374:   if (ksp->guess) {
375:     KSPFischerGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
376:     ksp->guess_zero = PETSC_FALSE;
377:   }
378:   /* KSPSetUp() scales the matrix if needed */
379:   KSPSetUp(ksp);
380:   KSPSetUpOnBlocks(ksp);

382:   /* diagonal scale RHS if called for */
383:   if (ksp->dscale) {
384:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
385:     /* second time in, but matrix was scaled back to original */
386:     if (ksp->dscalefix && ksp->dscalefix2) {
387:       Mat mat,pmat;

389:       PCGetOperators(ksp->pc,&mat,&pmat,PETSC_NULL);
390:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
391:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
392:     }

394:     /*  scale initial guess */
395:     if (!ksp->guess_zero) {
396:       if (!ksp->truediagonal) {
397:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
398:         VecCopy(ksp->diagonal,ksp->truediagonal);
399:         VecReciprocal(ksp->truediagonal);
400:       }
401:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
402:     }
403:   }
404:   PCPreSolve(ksp->pc,ksp);

406:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
407:   if (ksp->guess_knoll) {
408:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
409:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
410:     ksp->guess_zero = PETSC_FALSE;
411:   }

413:   /* can we mark the initial guess as zero for this solve? */
414:   guess_zero = ksp->guess_zero;
415:   if (!ksp->guess_zero) {
416:     PetscReal norm;

418:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
419:     if (flg && !norm) {
420:       ksp->guess_zero = PETSC_TRUE;
421:     }
422:   }
423:   (*ksp->ops->solve)(ksp);
424:   ksp->guess_zero = guess_zero;

426:   if (!ksp->reason) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
427:   if (ksp->printreason) {
428:     PetscViewerASCIIAddTab(PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm),((PetscObject)ksp)->tablevel);
429:     if (ksp->reason > 0) {
430:       PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm),"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
431:     } else {
432:       PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm),"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
433:     }
434:     PetscViewerASCIISubtractTab(PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm),((PetscObject)ksp)->tablevel);
435:   }
436:   PCPostSolve(ksp->pc,ksp);

438:   /* diagonal scale solution if called for */
439:   if (ksp->dscale) {
440:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
441:     /* unscale right hand side and matrix */
442:     if (ksp->dscalefix) {
443:       Mat mat,pmat;

445:       VecReciprocal(ksp->diagonal);
446:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
447:       PCGetOperators(ksp->pc,&mat,&pmat,PETSC_NULL);
448:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
449:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
450:       VecReciprocal(ksp->diagonal);
451:       ksp->dscalefix2 = PETSC_TRUE;
452:     }
453:   }
454:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);

456:   if (ksp->guess) {
457:     KSPFischerGuessUpdate(ksp->guess,ksp->vec_sol);
458:   }

460:   MPI_Comm_rank(((PetscObject)ksp)->comm,&rank);

462:   flag1 = PETSC_FALSE;
463:   flag2 = PETSC_FALSE;
464:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues",&flag1,PETSC_NULL);
465:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues",&flag2,PETSC_NULL);
466:   if (flag1 || flag2) {
467:     PetscInt   nits,n,i,neig;
468:     PetscReal *r,*c;
469: 
470:     KSPGetIterationNumber(ksp,&nits);
471:     n = nits+2;

473:     if (!nits) {
474:       PetscPrintf(((PetscObject)ksp)->comm,"Zero iterations in solver, cannot approximate any eigenvalues\n");
475:     } else {
476:       PetscMalloc(2*n*sizeof(PetscReal),&r);
477:       c = r + n;
478:       KSPComputeEigenvalues(ksp,n,r,c,&neig);
479:       if (flag1) {
480:         PetscPrintf(((PetscObject)ksp)->comm,"Iteratively computed eigenvalues\n");
481:         for (i=0; i<neig; i++) {
482:           if (c[i] >= 0.0) {PetscPrintf(((PetscObject)ksp)->comm,"%G + %Gi\n",r[i],c[i]);}
483:           else             {PetscPrintf(((PetscObject)ksp)->comm,"%G - %Gi\n",r[i],-c[i]);}
484:         }
485:       }
486:       if (flag2 && !rank) {
487:         PetscDraw   draw;
488:         PetscDrawSP drawsp;

490:         PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
491:         PetscViewerDrawGetDraw(viewer,0,&draw);
492:         PetscDrawSPCreate(draw,1,&drawsp);
493:         for (i=0; i<neig; i++) {
494:           PetscDrawSPAddPoint(drawsp,r+i,c+i);
495:         }
496:         PetscDrawSPDraw(drawsp);
497:         PetscDrawSPDestroy(&drawsp);
498:         PetscViewerDestroy(&viewer);
499:       }
500:       PetscFree(r);
501:     }
502:   }

504:   flag1 = PETSC_FALSE;
505:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_compute_singularvalues",&flag1,PETSC_NULL);
506:   if (flag1) {
507:     PetscInt   nits;

509:     KSPGetIterationNumber(ksp,&nits);
510:     if (!nits) {
511:       PetscPrintf(((PetscObject)ksp)->comm,"Zero iterations in solver, cannot approximate any singular values\n");
512:     } else {
513:       PetscReal emax,emin;

515:       KSPComputeExtremeSingularValues(ksp,&emax,&emin);
516:       PetscPrintf(((PetscObject)ksp)->comm,"Iteratively computed extreme singular values: max %G min %G max/min %G\n",emax,emin,emax/emin);
517:     }
518:   }


521:   flag1 = PETSC_FALSE;
522:   flag2 = PETSC_FALSE;
523:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_compute_eigenvalues_explicitly",&flag1,PETSC_NULL);
524:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_plot_eigenvalues_explicitly",&flag2,PETSC_NULL);
525:   if (flag1 || flag2) {
526:     PetscInt  n,i;
527:     PetscReal *r,*c;
528:     VecGetSize(ksp->vec_sol,&n);
529:     PetscMalloc2(n,PetscReal,&r,n,PetscReal,&c);
530:     KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
531:     if (flag1) {
532:       PetscPrintf(((PetscObject)ksp)->comm,"Explicitly computed eigenvalues\n");
533:       for (i=0; i<n; i++) {
534:         if (c[i] >= 0.0) {PetscPrintf(((PetscObject)ksp)->comm,"%G + %Gi\n",r[i],c[i]);}
535:         else             {PetscPrintf(((PetscObject)ksp)->comm,"%G - %Gi\n",r[i],-c[i]);}
536:       }
537:     }
538:     if (flag2 && !rank) {
539:       PetscDraw   draw;
540:       PetscDrawSP drawsp;

542:       PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Explicitly Computed Eigenvalues",0,320,300,300,&viewer);
543:       PetscViewerDrawGetDraw(viewer,0,&draw);
544:       PetscDrawSPCreate(draw,1,&drawsp);
545:       for (i=0; i<n; i++) {
546:         PetscDrawSPAddPoint(drawsp,r+i,c+i);
547:       }
548:       PetscDrawSPDraw(drawsp);
549:       PetscDrawSPDestroy(&drawsp);
550:       PetscViewerDestroy(&viewer);
551:     }
552:     PetscFree2(r,c);
553:   }

555:   flag2 = PETSC_FALSE;
556:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_view_operator",&flag2,PETSC_NULL);
557:   if (flag2) {
558:     Mat         A,B;
559:     PetscViewer viewer;

561:     PCGetOperators(ksp->pc,&A,PETSC_NULL,PETSC_NULL);
562:     MatComputeExplicitOperator(A,&B);
563:     PetscViewerASCIIGetStdout(((PetscObject)ksp)->comm,&viewer);
564:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
565:     MatView(B,viewer);
566:     PetscViewerPopFormat(viewer);
567:     MatDestroy(&B);
568:   }
569:   flag2 = PETSC_FALSE;
570:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_view_operator_binary",&flag2,PETSC_NULL);
571:   if (flag2) {
572:     Mat A,B;
573:     PCGetOperators(ksp->pc,&A,PETSC_NULL,PETSC_NULL);
574:     MatComputeExplicitOperator(A,&B);
575:     MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)ksp)->comm));
576:     MatDestroy(&B);
577:   }
578:   flag2 = PETSC_FALSE;
579:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_view_preconditioned_operator_binary",&flag2,PETSC_NULL);
580:   if (flag2) {
581:     Mat B;
582:     KSPComputeExplicitOperator(ksp,&B);
583:     MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)ksp)->comm));
584:     MatDestroy(&B);
585:   }
586:   flag2 = PETSC_FALSE;
587:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_view_preconditioner_binary",&flag2,PETSC_NULL);
588:   if (flag2) {
589:     Mat B;
590:     PCComputeExplicitOperator(ksp->pc,&B);
591:     MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)ksp)->comm));
592:     MatDestroy(&B);
593:   }
594:   PetscOptionsGetString(((PetscObject)ksp)->prefix,"-ksp_view",filename,PETSC_MAX_PATH_LEN,&flg);
595:   if (flg && !PetscPreLoadingOn) {
596:     PetscViewerASCIIOpen(((PetscObject)ksp)->comm,filename,&viewer);
597:     KSPView(ksp,viewer);
598:     PetscViewerDestroy(&viewer);
599:   }
600:   flg  = PETSC_FALSE;
601:   PetscOptionsGetBool(((PetscObject)ksp)->prefix,"-ksp_final_residual",&flg,PETSC_NULL);
602:   if (flg) {
603:     Mat         A;
604:     Vec         t;
605:     PetscReal   norm;
606:     if (ksp->dscale && !ksp->dscalefix) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot compute final scale with -ksp_diagonal_scale except also with -ksp_diagonal_scale_fix");
607:     PCGetOperators(ksp->pc,&A,0,0);
608:     VecDuplicate(ksp->vec_sol,&t);
609:     KSP_MatMult(ksp,A,ksp->vec_sol,t);
610:     VecAYPX(t, -1.0, ksp->vec_rhs);
611:     VecNorm(t,NORM_2,&norm);
612:     VecDestroy(&t);
613:     PetscPrintf(((PetscObject)ksp)->comm,"KSP final norm of residual %G\n",norm);
614:   }
615:   if (inXisinB) {
616:     VecCopy(x,b);
617:     VecDestroy(&x);
618:   }
619:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
620:   return(0);
621: }

625: /*@
626:    KSPSolveTranspose - Solves the transpose of a linear system. 

628:    Collective on KSP

630:    Input Parameter:
631: +  ksp - iterative context obtained from KSPCreate()
632: .  b - right hand side vector
633: -  x - solution vector

635:    Notes: For complex numbers this solve the non-Hermitian transpose system.

637:    Developer Notes: We need to implement a KSPSolveHermitianTranspose()

639:    Level: developer

641: .keywords: KSP, solve, linear system

643: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPDefaultConverged(),
644:           KSPSolve()
645: @*/

647: PetscErrorCode  KSPSolveTranspose(KSP ksp,Vec b,Vec x)
648: {
650:   PetscBool      inXisinB=PETSC_FALSE;

656:   if (x == b) {
657:     VecDuplicate(b,&x);
658:     inXisinB = PETSC_TRUE;
659:   }
660:   PetscObjectReference((PetscObject)b);
661:   PetscObjectReference((PetscObject)x);
662:   VecDestroy(&ksp->vec_rhs);
663:   VecDestroy(&ksp->vec_sol);
664:   ksp->vec_rhs = b;
665:   ksp->vec_sol = x;
666:   ksp->transpose_solve = PETSC_TRUE;
667:   KSPSetUp(ksp);
668:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
669:   (*ksp->ops->solve)(ksp);
670:   if (!ksp->reason) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
671:   if (ksp->printreason) {
672:     if (ksp->reason > 0) {
673:       PetscPrintf(((PetscObject)ksp)->comm,"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
674:     } else {
675:       PetscPrintf(((PetscObject)ksp)->comm,"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
676:     }
677:   }
678:   if (inXisinB) {
679:     VecCopy(x,b);
680:     VecDestroy(&x);
681:   }
682:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
683:   return(0);
684: }

688: /*@
689:    KSPReset - Resets a KSP context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats

691:    Collective on KSP

693:    Input Parameter:
694: .  ksp - iterative context obtained from KSPCreate()

696:    Level: beginner

698: .keywords: KSP, destroy

700: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
701: @*/
702: PetscErrorCode  KSPReset(KSP ksp)
703: {

708:   if (!ksp) return(0);
709:   if (ksp->ops->reset) {
710:     (*ksp->ops->reset)(ksp);
711:   }
712:   if (ksp->pc) {PCReset(ksp->pc);}
713:   KSPFischerGuessDestroy(&ksp->guess);
714:   VecDestroyVecs(ksp->nwork,&ksp->work);
715:   VecDestroy(&ksp->vec_rhs);
716:   VecDestroy(&ksp->vec_sol);
717:   VecDestroy(&ksp->diagonal);
718:   VecDestroy(&ksp->truediagonal);
719:   MatNullSpaceDestroy(&ksp->nullsp);
720:   ksp->setupstage = KSP_SETUP_NEW;
721:   return(0);
722: }

726: /*@
727:    KSPDestroy - Destroys KSP context.

729:    Collective on KSP

731:    Input Parameter:
732: .  ksp - iterative context obtained from KSPCreate()

734:    Level: beginner

736: .keywords: KSP, destroy

738: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
739: @*/
740: PetscErrorCode  KSPDestroy(KSP *ksp)
741: {

745:   if (!*ksp) return(0);
747:   if (--((PetscObject)(*ksp))->refct > 0) {*ksp = 0; return(0);}

749:   KSPReset((*ksp));

751:   PetscObjectDepublish((*ksp));
752:   if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

754:   DMDestroy(&(*ksp)->dm);
755:   PCDestroy(&(*ksp)->pc);
756:   PetscFree((*ksp)->res_hist_alloc);
757:   if ((*ksp)->convergeddestroy) {
758:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
759:   }
760:   KSPMonitorCancel((*ksp));
761:   PetscHeaderDestroy(ksp);
762:   return(0);
763: }

767: /*@
768:     KSPSetPCSide - Sets the preconditioning side.

770:     Logically Collective on KSP

772:     Input Parameter:
773: .   ksp - iterative context obtained from KSPCreate()

775:     Output Parameter:
776: .   side - the preconditioning side, where side is one of
777: .vb
778:       PC_LEFT - left preconditioning (default)
779:       PC_RIGHT - right preconditioning
780:       PC_SYMMETRIC - symmetric preconditioning
781: .ve

783:     Options Database Keys:
784: .   -ksp_pc_side <right,left,symmetric>

786:     Notes:
787:     Left preconditioning is used by default for most Krylov methods except KSPFGMRES which only supports right preconditioning.
788:     Symmetric preconditioning is currently available only for the KSPQCG method. Note, however, that
789:     symmetric preconditioning can be emulated by using either right or left
790:     preconditioning and a pre or post processing step.

792:     Level: intermediate

794: .keywords: KSP, set, right, left, symmetric, side, preconditioner, flag

796: .seealso: KSPGetPCSide()
797: @*/
798: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
799: {
803:   ksp->pc_side = side;
804:   return(0);
805: }

809: /*@
810:     KSPGetPCSide - Gets the preconditioning side.

812:     Not Collective

814:     Input Parameter:
815: .   ksp - iterative context obtained from KSPCreate()

817:     Output Parameter:
818: .   side - the preconditioning side, where side is one of
819: .vb
820:       PC_LEFT - left preconditioning (default)
821:       PC_RIGHT - right preconditioning
822:       PC_SYMMETRIC - symmetric preconditioning
823: .ve

825:     Level: intermediate

827: .keywords: KSP, get, right, left, symmetric, side, preconditioner, flag

829: .seealso: KSPSetPCSide()
830: @*/
831: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
832: {

838:   KSPSetUpNorms_Private(ksp);
839:   *side = ksp->pc_side;
840:   return(0);
841: }

845: /*@
846:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
847:    iteration tolerances used by the default KSP convergence tests. 

849:    Not Collective

851:    Input Parameter:
852: .  ksp - the Krylov subspace context
853:   
854:    Output Parameters:
855: +  rtol - the relative convergence tolerance
856: .  abstol - the absolute convergence tolerance
857: .  dtol - the divergence tolerance
858: -  maxits - maximum number of iterations

860:    Notes:
861:    The user can specify PETSC_NULL for any parameter that is not needed.

863:    Level: intermediate

865: .keywords: KSP, get, tolerance, absolute, relative, divergence, convergence,
866:            maximum, iterations

868: .seealso: KSPSetTolerances()
869: @*/
870: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
871: {
874:   if (abstol)   *abstol   = ksp->abstol;
875:   if (rtol)   *rtol   = ksp->rtol;
876:   if (dtol)   *dtol   = ksp->divtol;
877:   if (maxits) *maxits = ksp->max_it;
878:   return(0);
879: }

883: /*@
884:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
885:    iteration tolerances used by the default KSP convergence testers. 

887:    Logically Collective on KSP

889:    Input Parameters:
890: +  ksp - the Krylov subspace context
891: .  rtol - the relative convergence tolerance
892:    (relative decrease in the residual norm)
893: .  abstol - the absolute convergence tolerance 
894:    (absolute size of the residual norm)
895: .  dtol - the divergence tolerance
896:    (amount residual can increase before KSPDefaultConverged()
897:    concludes that the method is diverging)
898: -  maxits - maximum number of iterations to use

900:    Options Database Keys:
901: +  -ksp_atol <abstol> - Sets abstol
902: .  -ksp_rtol <rtol> - Sets rtol
903: .  -ksp_divtol <dtol> - Sets dtol
904: -  -ksp_max_it <maxits> - Sets maxits

906:    Notes:
907:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

909:    See KSPDefaultConverged() for details on the use of these parameters
910:    in the default convergence test.  See also KSPSetConvergenceTest() 
911:    for setting user-defined stopping criteria.

913:    Level: intermediate

915: .keywords: KSP, set, tolerance, absolute, relative, divergence, 
916:            convergence, maximum, iterations

918: .seealso: KSPGetTolerances(), KSPDefaultConverged(), KSPSetConvergenceTest()
919: @*/
920: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
921: {

929:   if (rtol != PETSC_DEFAULT) {
930:     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %G must be non-negative and less than 1.0",rtol);        ksp->rtol = rtol;
931:   }
932:   if (abstol != PETSC_DEFAULT) {
933:     if (abstol < 0.0) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol);
934:     ksp->abstol = abstol;
935:   }
936:   if (dtol != PETSC_DEFAULT) {
937:     if (dtol < 0.0) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %G must be larger than 1.0",dtol);
938:     ksp->divtol = dtol;
939:   }
940:   if (maxits != PETSC_DEFAULT) {
941:     if (maxits < 0) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
942:     ksp->max_it = maxits;
943:   }
944:   return(0);
945: }

949: /*@
950:    KSPSetInitialGuessNonzero - Tells the iterative solver that the 
951:    initial guess is nonzero; otherwise KSP assumes the initial guess
952:    is to be zero (and thus zeros it out before solving).

954:    Logically Collective on KSP

956:    Input Parameters:
957: +  ksp - iterative context obtained from KSPCreate()
958: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

960:    Options database keys:
961: .  -ksp_initial_guess_nonzero : use nonzero initial guess; this takes an optional truth value (0/1/no/yes/true/false)

963:    Level: beginner

965:    Notes:
966:     If this is not called the X vector is zeroed in the call to KSPSolve().

968: .keywords: KSP, set, initial guess, nonzero

970: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
971: @*/
972: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool  flg)
973: {
977:   ksp->guess_zero   = (PetscBool)!(int)flg;
978:   return(0);
979: }

983: /*@
984:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
985:    a zero initial guess.

987:    Not Collective

989:    Input Parameter:
990: .  ksp - iterative context obtained from KSPCreate()

992:    Output Parameter:
993: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

995:    Level: intermediate

997: .keywords: KSP, set, initial guess, nonzero

999: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
1000: @*/
1001: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1002: {
1006:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1007:   else                 *flag = PETSC_TRUE;
1008:   return(0);
1009: }

1013: /*@
1014:    KSPSetErrorIfNotConverged - Causes KSPSolve() to generate an error if the solver has not converged.

1016:    Logically Collective on KSP

1018:    Input Parameters:
1019: +  ksp - iterative context obtained from KSPCreate()
1020: -  flg - PETSC_TRUE indicates you want the error generated

1022:    Options database keys:
1023: .  -ksp_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)

1025:    Level: intermediate

1027:    Notes:
1028:     Normally PETSc continues if a linear solver fails to converge, you can call KSPGetConvergedReason() after a KSPSolve() 
1029:     to determine if it has converged.

1031: .keywords: KSP, set, initial guess, nonzero

1033: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPGetErrorIfNotConverged()
1034: @*/
1035: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool  flg)
1036: {
1040:   ksp->errorifnotconverged = flg;
1041:   return(0);
1042: }

1046: /*@
1047:    KSPGetErrorIfNotConverged - Will KSPSolve() generate an error if the solver does not converge?

1049:    Not Collective

1051:    Input Parameter:
1052: .  ksp - iterative context obtained from KSPCreate()

1054:    Output Parameter:
1055: .  flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE

1057:    Level: intermediate

1059: .keywords: KSP, set, initial guess, nonzero

1061: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll(), KSPSetErrorIfNotConverged()
1062: @*/
1063: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1064: {
1068:   *flag = ksp->errorifnotconverged;
1069:   return(0);
1070: }

1074: /*@
1075:    KSPSetInitialGuessKnoll - Tells the iterative solver to use PCApply(pc,b,..) to compute the initial guess (The Knoll trick)

1077:    Logically Collective on KSP

1079:    Input Parameters:
1080: +  ksp - iterative context obtained from KSPCreate()
1081: -  flg - PETSC_TRUE or PETSC_FALSE

1083:    Level: advanced


1086: .keywords: KSP, set, initial guess, nonzero

1088: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1089: @*/
1090: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool  flg)
1091: {
1095:   ksp->guess_knoll   = flg;
1096:   return(0);
1097: }

1101: /*@
1102:    KSPGetInitialGuessKnoll - Determines whether the KSP solver is using the Knoll trick (using PCApply(pc,b,...) to compute
1103:      the initial guess

1105:    Not Collective

1107:    Input Parameter:
1108: .  ksp - iterative context obtained from KSPCreate()

1110:    Output Parameter:
1111: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1113:    Level: advanced

1115: .keywords: KSP, set, initial guess, nonzero

1117: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
1118: @*/
1119: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1120: {
1124:   *flag = ksp->guess_knoll;
1125:   return(0);
1126: }

1130: /*@
1131:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular 
1132:    values will be calculated via a Lanczos or Arnoldi process as the linear 
1133:    system is solved.

1135:    Not Collective

1137:    Input Parameter:
1138: .  ksp - iterative context obtained from KSPCreate()

1140:    Output Parameter:
1141: .  flg - PETSC_TRUE or PETSC_FALSE

1143:    Options Database Key:
1144: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1146:    Notes:
1147:    Currently this option is not valid for all iterative methods.

1149:    Many users may just want to use the monitoring routine
1150:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1151:    to print the singular values at each iteration of the linear solve.

1153:    Level: advanced

1155: .keywords: KSP, set, compute, singular values

1157: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1158: @*/
1159: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1160: {
1164:   *flg = ksp->calc_sings;
1165:   return(0);
1166: }

1170: /*@
1171:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular 
1172:    values will be calculated via a Lanczos or Arnoldi process as the linear 
1173:    system is solved.

1175:    Logically Collective on KSP

1177:    Input Parameters:
1178: +  ksp - iterative context obtained from KSPCreate()
1179: -  flg - PETSC_TRUE or PETSC_FALSE

1181:    Options Database Key:
1182: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1184:    Notes:
1185:    Currently this option is not valid for all iterative methods.

1187:    Many users may just want to use the monitoring routine
1188:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1189:    to print the singular values at each iteration of the linear solve.

1191:    Level: advanced

1193: .keywords: KSP, set, compute, singular values

1195: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue()
1196: @*/
1197: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool  flg)
1198: {
1202:   ksp->calc_sings  = flg;
1203:   return(0);
1204: }

1208: /*@
1209:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1210:    values will be calculated via a Lanczos or Arnoldi process as the linear 
1211:    system is solved.

1213:    Not Collective

1215:    Input Parameter:
1216: .  ksp - iterative context obtained from KSPCreate()

1218:    Output Parameter:
1219: .  flg - PETSC_TRUE or PETSC_FALSE

1221:    Notes:
1222:    Currently this option is not valid for all iterative methods.

1224:    Level: advanced

1226: .keywords: KSP, set, compute, eigenvalues

1228: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1229: @*/
1230: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1231: {
1235:   *flg = ksp->calc_sings;
1236:   return(0);
1237: }

1241: /*@
1242:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1243:    values will be calculated via a Lanczos or Arnoldi process as the linear 
1244:    system is solved.

1246:    Logically Collective on KSP

1248:    Input Parameters:
1249: +  ksp - iterative context obtained from KSPCreate()
1250: -  flg - PETSC_TRUE or PETSC_FALSE

1252:    Notes:
1253:    Currently this option is not valid for all iterative methods.

1255:    Level: advanced

1257: .keywords: KSP, set, compute, eigenvalues

1259: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
1260: @*/
1261: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool  flg)
1262: {
1266:   ksp->calc_sings  = flg;
1267:   return(0);
1268: }

1272: /*@
1273:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1274:    be solved.

1276:    Not Collective

1278:    Input Parameter:
1279: .  ksp - iterative context obtained from KSPCreate()

1281:    Output Parameter:
1282: .  r - right-hand-side vector

1284:    Level: developer

1286: .keywords: KSP, get, right-hand-side, rhs

1288: .seealso: KSPGetSolution(), KSPSolve()
1289: @*/
1290: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1291: {
1295:   *r = ksp->vec_rhs;
1296:   return(0);
1297: }

1301: /*@
1302:    KSPGetSolution - Gets the location of the solution for the 
1303:    linear system to be solved.  Note that this may not be where the solution
1304:    is stored during the iterative process; see KSPBuildSolution().

1306:    Not Collective

1308:    Input Parameters:
1309: .  ksp - iterative context obtained from KSPCreate()

1311:    Output Parameters:
1312: .  v - solution vector

1314:    Level: developer

1316: .keywords: KSP, get, solution

1318: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve()
1319: @*/
1320: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1321: {
1325:   *v = ksp->vec_sol;
1326:   return(0);
1327: }

1331: /*@
1332:    KSPSetPC - Sets the preconditioner to be used to calculate the 
1333:    application of the preconditioner on a vector. 

1335:    Collective on KSP

1337:    Input Parameters:
1338: +  ksp - iterative context obtained from KSPCreate()
1339: -  pc   - the preconditioner object

1341:    Notes:
1342:    Use KSPGetPC() to retrieve the preconditioner context (for example,
1343:    to free it at the end of the computations).

1345:    Level: developer

1347: .keywords: KSP, set, precondition, Binv

1349: .seealso: KSPGetPC()
1350: @*/
1351: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1352: {

1359:   PetscObjectReference((PetscObject)pc);
1360:   PCDestroy(&ksp->pc);
1361:   ksp->pc = pc;
1362:   PetscLogObjectParent(ksp,ksp->pc);
1363:   return(0);
1364: }

1368: /*@
1369:    KSPGetPC - Returns a pointer to the preconditioner context
1370:    set with KSPSetPC().

1372:    Not Collective

1374:    Input Parameters:
1375: .  ksp - iterative context obtained from KSPCreate()

1377:    Output Parameter:
1378: .  pc - preconditioner context

1380:    Level: developer

1382: .keywords: KSP, get, preconditioner, Binv

1384: .seealso: KSPSetPC()
1385: @*/
1386: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1387: {

1393:   if (!ksp->pc) {
1394:     PCCreate(((PetscObject)ksp)->comm,&ksp->pc);
1395:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1396:     PetscLogObjectParent(ksp,ksp->pc);
1397:   }
1398:   *pc = ksp->pc;
1399:   return(0);
1400: }

1404: /*@
1405:    KSPMonitor - runs the user provided monitor routines, if they exist

1407:    Collective on KSP

1409:    Input Parameters:
1410: +  ksp - iterative context obtained from KSPCreate()
1411: .  it - iteration number
1412: -  rnorm - relative norm of the residual

1414:    Notes:
1415:    This routine is called by the KSP implementations.
1416:    It does not typically need to be called by the user.

1418:    Level: developer

1420: .seealso: KSPMonitorSet()
1421: @*/
1422: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1423: {
1424:   PetscInt       i, n = ksp->numbermonitors;

1428:   for (i=0; i<n; i++) {
1429:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1430:   }
1431:   return(0);
1432: }

1436: /*@C
1437:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor 
1438:    the residual/error etc.
1439:       
1440:    Logically Collective on KSP

1442:    Input Parameters:
1443: +  ksp - iterative context obtained from KSPCreate()
1444: .  monitor - pointer to function (if this is PETSC_NULL, it turns off monitoring
1445: .  mctx    - [optional] context for private data for the
1446:              monitor routine (use PETSC_NULL if no context is desired)
1447: -  monitordestroy - [optional] routine that frees monitor context
1448:           (may be PETSC_NULL)

1450:    Calling Sequence of monitor:
1451: $     monitor (KSP ksp, int it, PetscReal rnorm, void *mctx)

1453: +  ksp - iterative context obtained from KSPCreate()
1454: .  it - iteration number
1455: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1456: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

1458:    Options Database Keys:
1459: +    -ksp_monitor        - sets KSPMonitorDefault()
1460: .    -ksp_monitor_true_residual    - sets KSPMonitorTrueResidualNorm()
1461: .    -ksp_monitor_draw    - sets line graph monitor,
1462:                            uses KSPMonitorLGCreate()
1463: .    -ksp_monitor_draw_true_residual   - sets line graph monitor,
1464:                            uses KSPMonitorLGCreate()
1465: .    -ksp_monitor_singular_value    - sets KSPMonitorSingularValue()
1466: -    -ksp_monitor_cancel - cancels all monitors that have
1467:                           been hardwired into a code by 
1468:                           calls to KSPMonitorSet(), but
1469:                           does not cancel those set via
1470:                           the options database.

1472:    Notes:  
1473:    The default is to do nothing.  To print the residual, or preconditioned 
1474:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use 
1475:    KSPMonitorDefault() as the monitoring routine, with a null monitoring 
1476:    context. 

1478:    Several different monitoring routines may be set by calling
1479:    KSPMonitorSet() multiple times; all will be called in the 
1480:    order in which they were set. 

1482:    Fortran notes: Only a single monitor function can be set for each KSP object

1484:    Level: beginner

1486: .keywords: KSP, set, monitor

1488: .seealso: KSPMonitorDefault(), KSPMonitorLGCreate(), KSPMonitorCancel()
1489: @*/
1490: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1491: {
1492:   PetscInt       i;

1497:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1498:   for (i=0; i<ksp->numbermonitors;i++) {
1499:     if (monitor == ksp->monitor[i] && monitordestroy == ksp->monitordestroy[i] && mctx == ksp->monitorcontext[i]) {
1500:       if (monitordestroy) {
1501:         (*monitordestroy)(&mctx);
1502:       }
1503:       return(0);
1504:     }
1505:   }
1506:   ksp->monitor[ksp->numbermonitors]           = monitor;
1507:   ksp->monitordestroy[ksp->numbermonitors]    = monitordestroy;
1508:   ksp->monitorcontext[ksp->numbermonitors++]  = (void*)mctx;
1509:   return(0);
1510: }

1514: /*@
1515:    KSPMonitorCancel - Clears all monitors for a KSP object.

1517:    Logically Collective on KSP

1519:    Input Parameters:
1520: .  ksp - iterative context obtained from KSPCreate()

1522:    Options Database Key:
1523: .  -ksp_monitor_cancel - Cancels all monitors that have
1524:     been hardwired into a code by calls to KSPMonitorSet(), 
1525:     but does not cancel those set via the options database.

1527:    Level: intermediate

1529: .keywords: KSP, set, monitor

1531: .seealso: KSPMonitorDefault(), KSPMonitorLGCreate(), KSPMonitorSet()
1532: @*/
1533: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1534: {
1536:   PetscInt       i;

1540:   for (i=0; i<ksp->numbermonitors; i++) {
1541:     if (ksp->monitordestroy[i]) {
1542:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1543:     }
1544:   }
1545:   ksp->numbermonitors = 0;
1546:   return(0);
1547: }

1551: /*@C
1552:    KSPGetMonitorContext - Gets the monitoring context, as set by 
1553:    KSPMonitorSet() for the FIRST monitor only.

1555:    Not Collective

1557:    Input Parameter:
1558: .  ksp - iterative context obtained from KSPCreate()

1560:    Output Parameter:
1561: .  ctx - monitoring context

1563:    Level: intermediate

1565: .keywords: KSP, get, monitor, context

1567: .seealso: KSPMonitorDefault(), KSPMonitorLGCreate()
1568: @*/
1569: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1570: {
1573:   *ctx =      (ksp->monitorcontext[0]);
1574:   return(0);
1575: }

1579: /*@
1580:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1581:    If set, this array will contain the residual norms computed at each
1582:    iteration of the solver.

1584:    Not Collective

1586:    Input Parameters:
1587: +  ksp - iterative context obtained from KSPCreate()
1588: .  a   - array to hold history
1589: .  na  - size of a
1590: -  reset - PETSC_TRUE indicates the history counter is reset to zero
1591:            for each new linear solve

1593:    Level: advanced

1595:    Notes: The array is NOT freed by PETSc so the user needs to keep track of 
1596:            it and destroy once the KSP object is destroyed.

1598:    If 'a' is PETSC_NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
1599:    default array of length 10000 is allocated.

1601: .keywords: KSP, set, residual, history, norm

1603: .seealso: KSPGetResidualHistory()

1605: @*/
1606: PetscErrorCode  KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool  reset)
1607: {


1613:   PetscFree(ksp->res_hist_alloc);
1614:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
1615:     ksp->res_hist        = a;
1616:     ksp->res_hist_max    = na;
1617:   } else {
1618:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT)
1619:       ksp->res_hist_max = na;
1620:     else
1621:       ksp->res_hist_max = 10000; /* like default ksp->max_it */
1622:     PetscMalloc(ksp->res_hist_max*sizeof(PetscReal),&ksp->res_hist_alloc);
1623:     ksp->res_hist = ksp->res_hist_alloc;
1624:   }
1625:   ksp->res_hist_len    = 0;
1626:   ksp->res_hist_reset  = reset;

1628:   return(0);
1629: }

1633: /*@C
1634:    KSPGetResidualHistory - Gets the array used to hold the residual history
1635:    and the number of residuals it contains.

1637:    Not Collective

1639:    Input Parameter:
1640: .  ksp - iterative context obtained from KSPCreate()

1642:    Output Parameters:
1643: +  a   - pointer to array to hold history (or PETSC_NULL)
1644: -  na  - number of used entries in a (or PETSC_NULL)

1646:    Level: advanced

1648:    Notes:
1649:      Can only be called after a KSPSetResidualHistory() otherwise a and na are set to zero

1651:      The Fortran version of this routine has a calling sequence
1652: $   call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
1653:     note that you have passed a Fortran array into KSPSetResidualHistory() and you need
1654:     to access the residual values from this Fortran array you provided. Only the na (number of
1655:     residual norms currently held) is set.

1657: .keywords: KSP, get, residual, history, norm

1659: .seealso: KSPGetResidualHistory()

1661: @*/
1662: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1663: {
1666:   if (a)  *a  = ksp->res_hist;
1667:   if (na) *na = ksp->res_hist_len;
1668:   return(0);
1669: }

1673: /*@C
1674:    KSPSetConvergenceTest - Sets the function to be used to determine
1675:    convergence.  

1677:    Logically Collective on KSP

1679:    Input Parameters:
1680: +  ksp - iterative context obtained from KSPCreate()
1681: .  converge - pointer to int function
1682: .  cctx    - context for private data for the convergence routine (may be null)
1683: -  destroy - a routine for destroying the context (may be null)

1685:    Calling sequence of converge:
1686: $     converge (KSP ksp, int it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)

1688: +  ksp - iterative context obtained from KSPCreate()
1689: .  it - iteration number
1690: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1691: .  reason - the reason why it has converged or diverged
1692: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


1695:    Notes:
1696:    Must be called after the KSP type has been set so put this after
1697:    a call to KSPSetType(), or KSPSetFromOptions().

1699:    The default convergence test, KSPDefaultConverged(), aborts if the 
1700:    residual grows to more than 10000 times the initial residual.

1702:    The default is a combination of relative and absolute tolerances.  
1703:    The residual value that is tested may be an approximation; routines 
1704:    that need exact values should compute them.

1706:    In the default PETSc convergence test, the precise values of reason
1707:    are macros such as KSP_CONVERGED_RTOL, which are defined in petscksp.h.

1709:    Level: advanced

1711: .keywords: KSP, set, convergence, test, context

1713: .seealso: KSPDefaultConverged(), KSPGetConvergenceContext()
1714: @*/
1715: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1716: {

1721:   if (ksp->convergeddestroy) {
1722:     (*ksp->convergeddestroy)(ksp->cnvP);
1723:   }
1724:   ksp->converged        = converge;
1725:   ksp->convergeddestroy = destroy;
1726:   ksp->cnvP             = (void*)cctx;
1727:   return(0);
1728: }

1732: /*@C
1733:    KSPGetConvergenceContext - Gets the convergence context set with 
1734:    KSPSetConvergenceTest().  

1736:    Not Collective

1738:    Input Parameter:
1739: .  ksp - iterative context obtained from KSPCreate()

1741:    Output Parameter:
1742: .  ctx - monitoring context

1744:    Level: advanced

1746: .keywords: KSP, get, convergence, test, context

1748: .seealso: KSPDefaultConverged(), KSPSetConvergenceTest()
1749: @*/
1750: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
1751: {
1754:   *ctx = ksp->cnvP;
1755:   return(0);
1756: }

1760: /*@C
1761:    KSPBuildSolution - Builds the approximate solution in a vector provided.
1762:    This routine is NOT commonly needed (see KSPSolve()).

1764:    Collective on KSP

1766:    Input Parameter:
1767: .  ctx - iterative context obtained from KSPCreate()

1769:    Output Parameter: 
1770:    Provide exactly one of
1771: +  v - location to stash solution.   
1772: -  V - the solution is returned in this location. This vector is created 
1773:        internally. This vector should NOT be destroyed by the user with
1774:        VecDestroy().

1776:    Notes:
1777:    This routine can be used in one of two ways
1778: .vb
1779:       KSPBuildSolution(ksp,PETSC_NULL,&V);
1780:    or
1781:       KSPBuildSolution(ksp,v,PETSC_NULL); or KSPBuildSolution(ksp,v,&v);
1782: .ve
1783:    In the first case an internal vector is allocated to store the solution
1784:    (the user cannot destroy this vector). In the second case the solution
1785:    is generated in the vector that the user provides. Note that for certain 
1786:    methods, such as KSPCG, the second case requires a copy of the solution,
1787:    while in the first case the call is essentially free since it simply 
1788:    returns the vector where the solution already is stored. For some methods
1789:    like GMRES this is a reasonably expensive operation and should only be
1790:    used in truly needed.

1792:    Level: advanced

1794: .keywords: KSP, build, solution

1796: .seealso: KSPGetSolution(), KSPBuildResidual()
1797: @*/
1798: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
1799: {

1804:   if (!V && !v) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_WRONG,"Must provide either v or V");
1805:   if (!V) V = &v;
1806:   (*ksp->ops->buildsolution)(ksp,v,V);
1807:   return(0);
1808: }

1812: /*@C
1813:    KSPBuildResidual - Builds the residual in a vector provided.

1815:    Collective on KSP

1817:    Input Parameter:
1818: .  ksp - iterative context obtained from KSPCreate()

1820:    Output Parameters:
1821: +  v - optional location to stash residual.  If v is not provided,
1822:        then a location is generated.
1823: .  t - work vector.  If not provided then one is generated.
1824: -  V - the residual

1826:    Notes:
1827:    Regardless of whether or not v is provided, the residual is 
1828:    returned in V.

1830:    Level: advanced

1832: .keywords: KSP, build, residual

1834: .seealso: KSPBuildSolution()
1835: @*/
1836: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
1837: {
1839:   PetscBool      flag = PETSC_FALSE;
1840:   Vec            w = v,tt = t;

1844:   if (!w) {
1845:     VecDuplicate(ksp->vec_rhs,&w);
1846:     PetscLogObjectParent((PetscObject)ksp,w);
1847:   }
1848:   if (!tt) {
1849:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
1850:     PetscLogObjectParent((PetscObject)ksp,tt);
1851:   }
1852:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
1853:   if (flag) {VecDestroy(&tt);}
1854:   return(0);
1855: }

1859: /*@
1860:    KSPSetDiagonalScale - Tells KSP to symmetrically diagonally scale the system
1861:      before solving. This actually CHANGES the matrix (and right hand side).

1863:    Logically Collective on KSP

1865:    Input Parameter:
1866: +  ksp - the KSP context
1867: -  scale - PETSC_TRUE or PETSC_FALSE

1869:    Options Database Key:
1870: +   -ksp_diagonal_scale - 
1871: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve 


1874:     BE CAREFUL with this routine: it actually scales the matrix and right 
1875:     hand side that define the system. After the system is solved the matrix
1876:     and right hand side remain scaled.

1878:     This routine is only used if the matrix and preconditioner matrix are
1879:     the same thing.

1881:     This should NOT be used within the SNES solves if you are using a line
1882:     search.
1883:  
1884:     If you use this with the PCType Eisenstat preconditioner than you can 
1885:     use the PCEisenstatNoDiagonalScaling() option, or -pc_eisenstat_no_diagonal_scaling
1886:     to save some unneeded, redundant flops.

1888:    Level: intermediate

1890: .keywords: KSP, set, options, prefix, database

1892: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix()
1893: @*/
1894: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool  scale)
1895: {
1899:   ksp->dscale = scale;
1900:   return(0);
1901: }

1905: /*@
1906:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
1907:                           right hand side

1909:    Not Collective

1911:    Input Parameter:
1912: .  ksp - the KSP context

1914:    Output Parameter:
1915: .  scale - PETSC_TRUE or PETSC_FALSE

1917:    Notes:
1918:     BE CAREFUL with this routine: it actually scales the matrix and right 
1919:     hand side that define the system. After the system is solved the matrix
1920:     and right hand side remain scaled.

1922:     This routine is only used if the matrix and preconditioner matrix are
1923:     the same thing.

1925:    Level: intermediate

1927: .keywords: KSP, set, options, prefix, database

1929: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
1930: @*/
1931: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
1932: {
1936:   *scale = ksp->dscale;
1937:   return(0);
1938: }

1942: /*@
1943:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
1944:      back after solving.

1946:    Logically Collective on KSP

1948:    Input Parameter:
1949: +  ksp - the KSP context
1950: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not 
1951:          rescale (default)

1953:    Notes:
1954:      Must be called after KSPSetDiagonalScale()

1956:      Using this will slow things down, because it rescales the matrix before and
1957:      after each linear solve. This is intended mainly for testing to allow one
1958:      to easily get back the original system to make sure the solution computed is
1959:      accurate enough.

1961:     This routine is only used if the matrix and preconditioner matrix are
1962:     the same thing.

1964:    Level: intermediate

1966: .keywords: KSP, set, options, prefix, database

1968: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix()
1969: @*/
1970: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool  fix)
1971: {
1975:   ksp->dscalefix = fix;
1976:   return(0);
1977: }

1981: /*@
1982:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
1983:      back after solving.

1985:    Not Collective

1987:    Input Parameter:
1988: .  ksp - the KSP context

1990:    Output Parameter:
1991: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not 
1992:          rescale (default)

1994:    Notes:
1995:      Must be called after KSPSetDiagonalScale()

1997:      If PETSC_TRUE will slow things down, because it rescales the matrix before and
1998:      after each linear solve. This is intended mainly for testing to allow one
1999:      to easily get back the original system to make sure the solution computed is
2000:      accurate enough.

2002:     This routine is only used if the matrix and preconditioner matrix are
2003:     the same thing.

2005:    Level: intermediate

2007: .keywords: KSP, set, options, prefix, database

2009: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix()
2010: @*/
2011: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2012: {
2016:   *fix = ksp->dscalefix;
2017:   return(0);
2018: }