Actual source code: itfunc.c

  1: /*$Id: itfunc.c,v 1.159 2001/08/07 03:03:45 balay Exp $*/
  2: /*
  3:       Interface KSP routines that the user calls.
  4: */

 6:  #include src/sles/ksp/kspimpl.h

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

 15:    Not Collective

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

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

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

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

 31:    Level: advanced

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

 35: .seealso: KSPSetComputeSingularValues(), KSPSingularValueMonitor(), KSPComputeEigenvalues()
 36: @*/
 37: int KSPComputeExtremeSingularValues(KSP ksp,PetscReal *emax,PetscReal *emin)
 38: {

 45:   if (!ksp->calc_sings) {
 46:     SETERRQ(4,"Singular values not requested before KSPSetUp()");
 47:   }

 49:   if (ksp->ops->computeextremesingularvalues) {
 50:     (*ksp->ops->computeextremesingularvalues)(ksp,emax,emin);
 51:   }
 52:   return(0);
 53: }

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

 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() (that is the SLESSolve); 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:    KSPSingularValueMonitor() (which can be set with option -ksp_singmonitor)
 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(), KSPSingularValueMonitor(), KSPComputeExtremeSingularValues()
100: @*/
101: int KSPComputeEigenvalues(KSP ksp,int n,PetscReal *r,PetscReal *c,int *neig)
102: {

109:   if (!ksp->calc_sings) {
110:     SETERRQ(4,"Eigenvalues not requested before KSPSetUp()");
111:   }

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

119: #undef __FUNCT__  
121: /*@
122:    KSPSetUp - Sets up the internal data structures for the
123:    later use of an iterative solver.

125:    Collective on KSP

127:    Input Parameter:
128: .  ksp   - iterative context obtained from KSPCreate()

130:    Level: developer

132: .keywords: KSP, setup

134: .seealso: KSPCreate(), KSPSolve(), KSPDestroy()
135: @*/
136: int KSPSetUp(KSP ksp)
137: {


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

146:   if (!ksp->type_name){
147:     KSPSetType(ksp,KSPGMRES);
148:   }

150:   if (ksp->setupcalled) return(0);
151:   ksp->setupcalled = 1;
152:   (*ksp->ops->setup)(ksp);
153:   return(0);
154: }


157: #undef __FUNCT__  
159: /*@
160:    KSPSolve - Solves linear system; usually not called directly, rather 
161:    it is called by a call to SLESSolve().

163:    Collective on KSP

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

168:    Output Parameter:
169: .  its - number of iterations required

171:    Options Database Keys:
172: +  -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
173: .  -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
174: .  -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the 
175:       dense operator and useing LAPACK
176: -  -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues

178:    Notes:
179:    On return, the parameter "its" contains either the iteration
180:    number at which convergence was successfully reached or failure was detected.

182:    Call KSPGetConvergedReason() to determine if the solver converged or failed and 
183:    why.
184:    
185:    If using a direct method (e.g., via the KSP solver
186:    KSPPREONLY and a preconditioner such as PCLU/PCILU),
187:    then its=1.  See KSPSetTolerances() and KSPDefaultConverged()
188:    for more details.

190:    Understanding Convergence:
191:    The routines KSPSetMonitor(), KSPComputeEigenvalues(), and
192:    KSPComputeEigenvaluesExplicitly() provide information on additional
193:    options to monitor convergence and print eigenvalue information.

195:    Level: developer

197: .keywords: KSP, solve, linear system

199: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPDefaultConverged(),
200:           SLESSolve(), KSPSolveTranspose(), SLESGetKSP()
201: @*/
202: int KSPSolve(KSP ksp,int *its)
203: {
204:   int          ierr,rank,nits;
205:   PetscTruth   flag1,flag2;
206:   PetscScalar  zero = 0.0;


211:   if (!ksp->setupcalled){ KSPSetUp(ksp);}
212:   if (ksp->guess_zero) { VecSet(&zero,ksp->vec_sol);}
213:   if (ksp->guess_knoll) {
214:     ierr            = PCApply(ksp->B,ksp->vec_rhs,ksp->vec_sol);
215:     ksp->guess_zero = PETSC_FALSE;
216:   }

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

221:   ksp->transpose_solve = PETSC_FALSE;
222:   (*ksp->ops->solve)(ksp,&nits);
223:   if (!ksp->reason) {
224:     SETERRQ(1,"Internal error, solver returned without setting converged reason");
225:   }
226:   if (its) *its = nits;

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

230:   PetscOptionsHasName(ksp->prefix,"-ksp_compute_eigenvalues",&flag1);
231:   PetscOptionsHasName(ksp->prefix,"-ksp_plot_eigenvalues",&flag2);
232:   if (flag1 || flag2) {
233:     int       n = nits + 2,i,neig;
234:     PetscReal *r,*c;

236:     if (!n) {
237:       PetscPrintf(ksp->comm,"Zero iterations in solver, cannot approximate any eigenvaluesn");
238:     } else {
239:       PetscMalloc(2*n*sizeof(PetscReal),&r);
240:       c = r + n;
241:       KSPComputeEigenvalues(ksp,n,r,c,&neig);
242:       if (flag1) {
243:         PetscPrintf(ksp->comm,"Iteratively computed eigenvaluesn");
244:         for (i=0; i<neig; i++) {
245:           if (c[i] >= 0.0) {PetscPrintf(ksp->comm,"%g + %gin",r[i],c[i]);}
246:           else             {PetscPrintf(ksp->comm,"%g - %gin",r[i],-c[i]);}
247:         }
248:       }
249:       if (flag2 && !rank) {
250:         PetscViewer viewer;
251:         PetscDraw   draw;
252:         PetscDrawSP drawsp;

254:         PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Iteratively Computed Eigenvalues",
255:                                PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
256:         PetscViewerDrawGetDraw(viewer,0,&draw);
257:         PetscDrawSPCreate(draw,1,&drawsp);
258:         for (i=0; i<neig; i++) {
259:           PetscDrawSPAddPoint(drawsp,r+i,c+i);
260:         }
261:         PetscDrawSPDraw(drawsp);
262:         PetscDrawSPDestroy(drawsp);
263:         PetscViewerDestroy(viewer);
264:       }
265:       PetscFree(r);
266:     }
267:   }

269:   PetscOptionsHasName(ksp->prefix,"-ksp_compute_eigenvalues_explicitly",&flag1);
270:   PetscOptionsHasName(ksp->prefix,"-ksp_plot_eigenvalues_explicitly",&flag2);
271:   if (flag1 || flag2) {
272:     int       n,i;
273:     PetscReal *r,*c;
274:     VecGetSize(ksp->vec_sol,&n);
275:     PetscMalloc(2*n*sizeof(PetscReal),&r);
276:     c = r + n;
277:     KSPComputeEigenvaluesExplicitly(ksp,n,r,c);
278:     if (flag1) {
279:       PetscPrintf(ksp->comm,"Explicitly computed eigenvaluesn");
280:       for (i=0; i<n; i++) {
281:         if (c[i] >= 0.0) {PetscPrintf(ksp->comm,"%g + %gin",r[i],c[i]);}
282:         else             {PetscPrintf(ksp->comm,"%g - %gin",r[i],-c[i]);}
283:       }
284:     }
285:     if (flag2 && !rank) {
286:       PetscViewer viewer;
287:       PetscDraw   draw;
288:       PetscDrawSP drawsp;

290:       PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Explicitly Computed Eigenvalues",0,320,300,300,&viewer);
291:       PetscViewerDrawGetDraw(viewer,0,&draw);
292:       PetscDrawSPCreate(draw,1,&drawsp);
293:       for (i=0; i<n; i++) {
294:         PetscDrawSPAddPoint(drawsp,r+i,c+i);
295:       }
296:       PetscDrawSPDraw(drawsp);
297:       PetscDrawSPDestroy(drawsp);
298:       PetscViewerDestroy(viewer);
299:     }
300:     PetscFree(r);
301:   }

303:   PetscOptionsHasName(ksp->prefix,"-ksp_view_operator",&flag2);
304:   if (flag2) {
305:     Mat A,B;
306:     PCGetOperators(ksp->B,&A,PETSC_NULL,PETSC_NULL);
307:     MatComputeExplicitOperator(A,&B);
308:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(ksp->comm),PETSC_VIEWER_ASCII_MATLAB);
309:     MatView(B,PETSC_VIEWER_STDOUT_(ksp->comm));
310:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(ksp->comm));
311:     MatDestroy(B);
312:   }
313:   return(0);
314: }

316: #undef __FUNCT__  
318: /*@
319:    KSPSolveTranspose - Solves the transpose of a linear system. Usually
320:    accessed through SLESSolveTranspose().

322:    Collective on KSP

324:    Input Parameter:
325: .  ksp - iterative context obtained from KSPCreate()

327:    Output Parameter:
328: .  its - number of iterations required

330:    Notes:
331:    On return, the parameter "its" contains either the iteration
332:    number at which convergence was successfully reached, or the
333:    negative of the iteration at which divergence or breakdown was detected.

335:    Currently only supported by KSPType of KSPPREONLY. This routine is usally 
336:    only used internally by the BiCG solver on the subblocks in BJacobi and ASM.

338:    Level: developer

340: .keywords: KSP, solve, linear system

342: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPDefaultConverged(),
343:           SLESSolve(), SLESGetKSP()
344: @*/
345: int KSPSolveTranspose(KSP ksp,int *its)
346: {
347:   int           ierr;
348:   PetscScalar   zero = 0.0;


353:   if (!ksp->setupcalled){ KSPSetUp(ksp);}
354:   if (ksp->guess_zero) { VecSet(&zero,ksp->vec_sol);}
355:   ksp->transpose_solve = PETSC_TRUE;
356:   (*ksp->ops->solve)(ksp,its);
357:   return(0);
358: }

360: #undef __FUNCT__  
362: /*@C
363:    KSPDestroy - Destroys KSP context.

365:    Collective on KSP

367:    Input Parameter:
368: .  ksp - iterative context obtained from KSPCreate()

370:    Level: developer

372: .keywords: KSP, destroy

374: .seealso: KSPCreate(), KSPSetUp(), KSPSolve()
375: @*/
376: int KSPDestroy(KSP ksp)
377: {
378:   int i,ierr;

382:   if (--ksp->refct > 0) return(0);

384:   /* if memory was published with AMS then destroy it */
385:   PetscObjectDepublish(ksp);

387:   if (ksp->ops->destroy) {
388:     (*ksp->ops->destroy)(ksp);
389:   }
390:   for (i=0; i<ksp->numbermonitors; i++) {
391:     if (ksp->monitordestroy[i]) {
392:       (*ksp->monitordestroy[i])(ksp->monitorcontext[i]);
393:     }
394:   }
395:   PetscLogObjectDestroy(ksp);
396:   PetscHeaderDestroy(ksp);
397:   return(0);
398: }

400: #undef __FUNCT__  
402: /*@
403:     KSPSetPreconditionerSide - Sets the preconditioning side.

405:     Collective on KSP

407:     Input Parameter:
408: .   ksp - iterative context obtained from KSPCreate()

410:     Output Parameter:
411: .   side - the preconditioning side, where side is one of
412: .vb
413:       PC_LEFT - left preconditioning (default)
414:       PC_RIGHT - right preconditioning
415:       PC_SYMMETRIC - symmetric preconditioning
416: .ve

418:     Options Database Keys:
419: +   -ksp_left_pc - Sets left preconditioning
420: .   -ksp_right_pc - Sets right preconditioning
421: -   -ksp_symmetric_pc - Sets symmetric preconditioning

423:     Notes:
424:     Left preconditioning is used by default.  Symmetric preconditioning is
425:     currently available only for the KSPQCG method. Note, however, that
426:     symmetric preconditioning can be emulated by using either right or left
427:     preconditioning and a pre or post processing step.

429:     Level: intermediate

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

433: .seealso: KSPGetPreconditionerSide()
434: @*/
435: int KSPSetPreconditionerSide(KSP ksp,PCSide side)
436: {
439:   ksp->pc_side = side;
440:   return(0);
441: }

443: #undef __FUNCT__  
445: /*@C
446:     KSPGetPreconditionerSide - Gets the preconditioning side.

448:     Not Collective

450:     Input Parameter:
451: .   ksp - iterative context obtained from KSPCreate()

453:     Output Parameter:
454: .   side - the preconditioning side, where side is one of
455: .vb
456:       PC_LEFT - left preconditioning (default)
457:       PC_RIGHT - right preconditioning
458:       PC_SYMMETRIC - symmetric preconditioning
459: .ve

461:     Level: intermediate

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

465: .seealso: KSPSetPreconditionerSide()
466: @*/
467: int KSPGetPreconditionerSide(KSP ksp,PCSide *side)
468: {
471:   *side = ksp->pc_side;
472:   return(0);
473: }

475: #undef __FUNCT__  
477: /*@
478:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
479:    iteration tolerances used by the default KSP convergence tests. 

481:    Not Collective

483:    Input Parameter:
484: .  ksp - the Krylov subspace context
485:   
486:    Output Parameters:
487: +  rtol - the relative convergence tolerance
488: .  atol - the absolute convergence tolerance
489: .  dtol - the divergence tolerance
490: -  maxits - maximum number of iterations

492:    Notes:
493:    The user can specify PETSC_NULL for any parameter that is not needed.

495:    Level: intermediate

497: .keywords: KSP, get, tolerance, absolute, relative, divergence, convergence,
498:            maximum, iterations

500: .seealso: KSPSetTolerances()
501: @*/
502: int KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *atol,PetscReal *dtol,int *maxits)
503: {
506:   if (atol)   *atol   = ksp->atol;
507:   if (rtol)   *rtol   = ksp->rtol;
508:   if (dtol)   *dtol   = ksp->divtol;
509:   if (maxits) *maxits = ksp->max_it;
510:   return(0);
511: }

513: #undef __FUNCT__  
515: /*@
516:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
517:    iteration tolerances used by the default KSP convergence testers. 

519:    Collective on KSP

521:    Input Parameters:
522: +  ksp - the Krylov subspace context
523: .  rtol - the relative convergence tolerance
524:    (relative decrease in the residual norm)
525: .  atol - the absolute convergence tolerance 
526:    (absolute size of the residual norm)
527: .  dtol - the divergence tolerance
528:    (amount residual can increase before KSPDefaultConverged()
529:    concludes that the method is diverging)
530: -  maxits - maximum number of iterations to use

532:    Options Database Keys:
533: +  -ksp_atol <atol> - Sets atol
534: .  -ksp_rtol <rtol> - Sets rtol
535: .  -ksp_divtol <dtol> - Sets dtol
536: -  -ksp_max_it <maxits> - Sets maxits

538:    Notes:
539:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

545:    Level: intermediate

547: .keywords: KSP, set, tolerance, absolute, relative, divergence, 
548:            convergence, maximum, iterations

550: .seealso: KSPGetTolerances(), KSPDefaultConverged(), KSPSetConvergenceTest()
551: @*/
552: int KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal atol,PetscReal dtol,int maxits)
553: {
556:   if (atol != PETSC_DEFAULT)   ksp->atol   = atol;
557:   if (rtol != PETSC_DEFAULT)   ksp->rtol   = rtol;
558:   if (dtol != PETSC_DEFAULT)   ksp->divtol = dtol;
559:   if (maxits != PETSC_DEFAULT) ksp->max_it = maxits;
560:   return(0);
561: }

563: #undef __FUNCT__  
565: /*@
566:    KSPSetInitialGuessNonzero - Tells the iterative solver that the 
567:    initial guess is nonzero; otherwise KSP assumes the initial guess
568:    is to be zero (and thus zeros it out before solving).

570:    Collective on KSP

572:    Input Parameters:
573: +  ksp - iterative context obtained from SLESGetKSP() or KSPCreate()
574: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

576:    Level: beginner

578:    Notes:
579:     If this is not called the X vector is zeroed in the call to 
580: SLESSolve() (or KSPSolve()).

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

584: .seealso: KSPGetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
585: @*/
586: int KSPSetInitialGuessNonzero(KSP ksp,PetscTruth flg)
587: {
589:   ksp->guess_zero   = (PetscTruth)!(int)flg;
590:   return(0);
591: }

593: #undef __FUNCT__  
595: /*@
596:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
597:    a zero initial guess.

599:    Not Collective

601:    Input Parameter:
602: .  ksp - iterative context obtained from KSPCreate()

604:    Output Parameter:
605: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

607:    Level: intermediate

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

611: .seealso: KSPSetInitialGuessNonzero(), KSPSetInitialGuessKnoll(), KSPGetInitialGuessKnoll()
612: @*/
613: int KSPGetInitialGuessNonzero(KSP ksp,PetscTruth *flag)
614: {
616:   if (ksp->guess_zero) *flag = PETSC_FALSE;
617:   else                 *flag = PETSC_TRUE;
618:   return(0);
619: }

621: #undef __FUNCT__  
623: /*@
624:    KSPSetInitialGuessKnoll - Tells the iterative solver to use PCApply(pc,b,..) to compute the initial guess (The Knoll trick)

626:    Collective on KSP

628:    Input Parameters:
629: +  ksp - iterative context obtained from SLESGetKSP() or KSPCreate()
630: -  flg - PETSC_TRUE or PETSC_FALSE

632:    Level: advanced


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

637: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
638: @*/
639: int KSPSetInitialGuessKnoll(KSP ksp,PetscTruth flg)
640: {
642:   ksp->guess_knoll   = flg;
643:   return(0);
644: }

646: #undef __FUNCT__  
648: /*@
649:    KSPGetInitialGuessKnoll - Determines whether the KSP solver is using the Knoll trick (using PCApply(pc,b,...) to compute
650:      the initial guess

652:    Not Collective

654:    Input Parameter:
655: .  ksp - iterative context obtained from KSPCreate()

657:    Output Parameter:
658: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

660:    Level: advanced

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

664: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero()
665: @*/
666: int KSPGetInitialGuessKnoll(KSP ksp,PetscTruth *flag)
667: {
669:   *flag = ksp->guess_knoll;
670:   return(0);
671: }

673: #undef __FUNCT__  
675: /*@
676:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular 
677:    values will be calculated via a Lanczos or Arnoldi process as the linear 
678:    system is solved.

680:    Collective on KSP

682:    Input Parameters:
683: +  ksp - iterative context obtained from KSPCreate()
684: -  flg - PETSC_TRUE or PETSC_FALSE

686:    Options Database Key:
687: .  -ksp_singmonitor - Activates KSPSetComputeSingularValues()

689:    Notes:
690:    Currently this option is not valid for all iterative methods.

692:    Many users may just want to use the monitoring routine
693:    KSPSingularValueMonitor() (which can be set with option -ksp_singmonitor)
694:    to print the singular values at each iteration of the linear solve.

696:    Level: advanced

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

700: .seealso: KSPComputeExtremeSingularValues(), KSPSingularValueMonitor()
701: @*/
702: int KSPSetComputeSingularValues(KSP ksp,PetscTruth flg)
703: {
706:   ksp->calc_sings  = flg;
707:   return(0);
708: }

710: #undef __FUNCT__  
712: /*@
713:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
714:    values will be calculated via a Lanczos or Arnoldi process as the linear 
715:    system is solved.

717:    Collective on KSP

719:    Input Parameters:
720: +  ksp - iterative context obtained from KSPCreate()
721: -  flg - PETSC_TRUE or PETSC_FALSE

723:    Notes:
724:    Currently this option is not valid for all iterative methods.

726:    Level: advanced

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

730: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly()
731: @*/
732: int KSPSetComputeEigenvalues(KSP ksp,PetscTruth flg)
733: {
736:   ksp->calc_sings  = flg;
737:   return(0);
738: }

740: #undef __FUNCT__  
742: /*@
743:    KSPSetRhs - Sets the right-hand-side vector for the linear system to
744:    be solved.

746:    Collective on KSP and Vec

748:    Input Parameters:
749: +  ksp - iterative context obtained from KSPCreate()
750: -  b   - right-hand-side vector

752:    Level: developer

754: .keywords: KSP, set, right-hand-side, rhs

756: .seealso: KSPGetRhs(), KSPSetSolution()
757: @*/
758: int KSPSetRhs(KSP ksp,Vec b)
759: {
763:   ksp->vec_rhs    = (b);
764:   return(0);
765: }

767: #undef __FUNCT__  
769: /*@C
770:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
771:    be solved.

773:    Not Collective

775:    Input Parameter:
776: .  ksp - iterative context obtained from KSPCreate()

778:    Output Parameter:
779: .  r - right-hand-side vector

781:    Level: developer

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

785: .seealso: KSPSetRhs(), KSPGetSolution()
786: @*/
787: int KSPGetRhs(KSP ksp,Vec *r)
788: {
791:   *r = ksp->vec_rhs;
792:   return(0);
793: }

795: #undef __FUNCT__  
797: /*@
798:    KSPSetSolution - Sets the location of the solution for the 
799:    linear system to be solved.

801:    Collective on KSP and Vec

803:    Input Parameters:
804: +  ksp - iterative context obtained from KSPCreate()
805: -  x   - solution vector

807:    Level: developer

809: .keywords: KSP, set, solution

811: .seealso: KSPSetRhs(), KSPGetSolution()
812: @*/
813: int KSPSetSolution(KSP ksp,Vec x)
814: {
818:   ksp->vec_sol    = (x);
819:   return(0);
820: }

822: #undef __FUNCT__  
824: /*@C
825:    KSPGetSolution - Gets the location of the solution for the 
826:    linear system to be solved.  Note that this may not be where the solution
827:    is stored during the iterative process; see KSPBuildSolution().

829:    Not Collective

831:    Input Parameters:
832: .  ksp - iterative context obtained from KSPCreate()

834:    Output Parameters:
835: .  v - solution vector

837:    Level: developer

839: .keywords: KSP, get, solution

841: .seealso: KSPGetRhs(), KSPSetSolution(), KSPBuildSolution()
842: @*/
843: int KSPGetSolution(KSP ksp,Vec *v)
844: {
847:   *v = ksp->vec_sol;
848:   return(0);
849: }

851: #undef __FUNCT__  
853: /*@
854:    KSPSetPC - Sets the preconditioner to be used to calculate the 
855:    application of the preconditioner on a vector. 

857:    Collective on KSP

859:    Input Parameters:
860: +  ksp - iterative context obtained from KSPCreate()
861: -  B   - the preconditioner object

863:    Notes:
864:    Use KSPGetPC() to retrieve the preconditioner context (for example,
865:    to free it at the end of the computations).

867:    Level: developer

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

871: .seealso: KSPGetPC()
872: @*/
873: int KSPSetPC(KSP ksp,PC B)
874: {
879:   ksp->B = B;
880:   return(0);
881: }

883: #undef __FUNCT__  
885: /*@C
886:    KSPGetPC - Returns a pointer to the preconditioner context
887:    set with KSPSetPC().

889:    Not Collective

891:    Input Parameters:
892: .  ksp - iterative context obtained from KSPCreate()

894:    Output Parameter:
895: .  B - preconditioner context

897:    Level: developer

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

901: .seealso: KSPSetPC()
902: @*/
903: int KSPGetPC(KSP ksp,PC *B)
904: {
907:   *B = ksp->B;
908:   return(0);
909: }

911: #undef __FUNCT__  
913: /*@C
914:    KSPSetMonitor - Sets an ADDITIONAL function to be called at every iteration to monitor 
915:    the residual/error etc.
916:       
917:    Collective on KSP

919:    Input Parameters:
920: +  ksp - iterative context obtained from KSPCreate()
921: .  monitor - pointer to function (if this is PETSC_NULL, it turns off monitoring
922: .  mctx    - [optional] context for private data for the
923:              monitor routine (use PETSC_NULL if no context is desired)
924: -  monitordestroy - [optional] routine that frees monitor context
925:           (may be PETSC_NULL)

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

930: +  ksp - iterative context obtained from KSPCreate()
931: .  it - iteration number
932: .  rnorm - (estimated) 2-norm of (preconditioned) residual
933: -  mctx  - optional monitoring context, as set by KSPSetMonitor()

935:    Options Database Keys:
936: +    -ksp_monitor        - sets KSPDefaultMonitor()
937: .    -ksp_truemonitor    - sets KSPTrueMonitor()
938: .    -ksp_xmonitor       - sets line graph monitor,
939:                            uses KSPLGMonitorCreate()
940: .    -ksp_xtruemonitor   - sets line graph monitor,
941:                            uses KSPLGMonitorCreate()
942: .    -ksp_singmonitor    - sets KSPSingularValueMonitor()
943: -    -ksp_cancelmonitors - cancels all monitors that have
944:                           been hardwired into a code by 
945:                           calls to KSPSetMonitor(), but
946:                           does not cancel those set via
947:                           the options database.

949:    Notes:  
950:    The default is to do nothing.  To print the residual, or preconditioned 
951:    residual if KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM) was called, use 
952:    KSPDefaultMonitor() as the monitoring routine, with a null monitoring 
953:    context. 

955:    Several different monitoring routines may be set by calling
956:    KSPSetMonitor() multiple times; all will be called in the 
957:    order in which they were set.

959:    Level: beginner

961: .keywords: KSP, set, monitor

963: .seealso: KSPDefaultMonitor(), KSPLGMonitorCreate(), KSPClearMonitor()
964: @*/
965: int KSPSetMonitor(KSP ksp,int (*monitor)(KSP,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void*))
966: {
969:   if (ksp->numbermonitors >= MAXKSPMONITORS) {
970:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
971:   }
972:   ksp->monitor[ksp->numbermonitors]           = monitor;
973:   ksp->monitordestroy[ksp->numbermonitors]    = monitordestroy;
974:   ksp->monitorcontext[ksp->numbermonitors++]  = (void*)mctx;
975:   return(0);
976: }

978: #undef __FUNCT__  
980: /*@
981:    KSPClearMonitor - Clears all monitors for a KSP object.

983:    Collective on KSP

985:    Input Parameters:
986: .  ksp - iterative context obtained from KSPCreate()

988:    Options Database Key:
989: .  -ksp_cancelmonitors - Cancels all monitors that have
990:     been hardwired into a code by calls to KSPSetMonitor(), 
991:     but does not cancel those set via the options database.

993:    Level: intermediate

995: .keywords: KSP, set, monitor

997: .seealso: KSPDefaultMonitor(), KSPLGMonitorCreate(), KSPSetMonitor()
998: @*/
999: int KSPClearMonitor(KSP ksp)
1000: {
1003:   ksp->numbermonitors = 0;
1004:   return(0);
1005: }

1007: #undef __FUNCT__  
1009: /*@C
1010:    KSPGetMonitorContext - Gets the monitoring context, as set by 
1011:    KSPSetMonitor() for the FIRST monitor only.

1013:    Not Collective

1015:    Input Parameter:
1016: .  ksp - iterative context obtained from KSPCreate()

1018:    Output Parameter:
1019: .  ctx - monitoring context

1021:    Level: intermediate

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

1025: .seealso: KSPDefaultMonitor(), KSPLGMonitorCreate()
1026: @*/
1027: int KSPGetMonitorContext(KSP ksp,void **ctx)
1028: {
1031:   *ctx =      (ksp->monitorcontext[0]);
1032:   return(0);
1033: }

1035: #undef __FUNCT__  
1037: /*@
1038:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1039:    If set, this array will contain the residual norms computed at each
1040:    iteration of the solver.

1042:    Not Collective

1044:    Input Parameters:
1045: +  ksp - iterative context obtained from KSPCreate()
1046: .  a   - array to hold history
1047: .  na  - size of a
1048: -  reset - PETSC_TRUE indicates the history counter is reset to zero
1049:            for each new linear solve

1051:    Level: advanced

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

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

1058: .seealso: KSPGetResidualHistory()

1060: @*/
1061: int KSPSetResidualHistory(KSP ksp,PetscReal *a,int na,PetscTruth reset)
1062: {

1067:   if (na != PETSC_DECIDE && a != PETSC_NULL) {
1068:     ksp->res_hist        = a;
1069:     ksp->res_hist_max    = na;
1070:   } else {
1071:     ksp->res_hist_max    = 1000;
1072:     PetscMalloc(ksp->res_hist_max*sizeof(PetscReal),&ksp->res_hist);
1073:   }
1074:   ksp->res_hist_len    = 0;
1075:   ksp->res_hist_reset  = reset;


1078:   return(0);
1079: }

1081: #undef __FUNCT__  
1083: /*@C
1084:    KSPGetResidualHistory - Gets the array used to hold the residual history
1085:    and the number of residuals it contains.

1087:    Not Collective

1089:    Input Parameter:
1090: .  ksp - iterative context obtained from KSPCreate()

1092:    Output Parameters:
1093: +  a   - pointer to array to hold history (or PETSC_NULL)
1094: -  na  - number of used entries in a (or PETSC_NULL)

1096:    Level: advanced

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

1101:      The Fortran version of this routine has a calling sequence
1102: $   call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)

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

1106: .seealso: KSPGetResidualHistory()

1108: @*/
1109: int KSPGetResidualHistory(KSP ksp,PetscReal **a,int *na)
1110: {
1113:   if (a)  *a = ksp->res_hist;
1114:   if (na) *na = ksp->res_hist_len;
1115:   return(0);
1116: }

1118: #undef __FUNCT__  
1120: /*@C
1121:    KSPSetConvergenceTest - Sets the function to be used to determine
1122:    convergence.  

1124:    Collective on KSP

1126:    Input Parameters:
1127: +  ksp - iterative context obtained from KSPCreate()
1128: .  converge - pointer to int function
1129: -  cctx    - context for private data for the convergence routine (may be null)

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

1134: +  ksp - iterative context obtained from KSPCreate()
1135: .  it - iteration number
1136: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1137: .  reason - the reason why it has converged or diverged
1138: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


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

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

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

1152:    Level: advanced

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

1156: .seealso: KSPDefaultConverged(), KSPGetConvergenceContext()
1157: @*/
1158: int KSPSetConvergenceTest(KSP ksp,int (*converge)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *cctx)
1159: {
1162:   ksp->converged = converge;
1163:   ksp->cnvP      = (void*)cctx;
1164:   return(0);
1165: }

1167: #undef __FUNCT__  
1169: /*@C
1170:    KSPGetConvergenceContext - Gets the convergence context set with 
1171:    KSPSetConvergenceTest().  

1173:    Not Collective

1175:    Input Parameter:
1176: .  ksp - iterative context obtained from KSPCreate()

1178:    Output Parameter:
1179: .  ctx - monitoring context

1181:    Level: advanced

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

1185: .seealso: KSPDefaultConverged(), KSPSetConvergenceTest()
1186: @*/
1187: int KSPGetConvergenceContext(KSP ksp,void **ctx)
1188: {
1191:   *ctx = ksp->cnvP;
1192:   return(0);
1193: }

1195: #undef __FUNCT__  
1197: /*@C
1198:    KSPBuildSolution - Builds the approximate solution in a vector provided.
1199:    This routine is NOT commonly needed (see SLESSolve()).

1201:    Collective on KSP

1203:    Input Parameter:
1204: .  ctx - iterative context obtained from KSPCreate()

1206:    Output Parameter: 
1207:    Provide exactly one of
1208: +  v - location to stash solution.   
1209: -  V - the solution is returned in this location. This vector is created 
1210:        internally. This vector should NOT be destroyed by the user with
1211:        VecDestroy().

1213:    Notes:
1214:    This routine can be used in one of two ways
1215: .vb
1216:       KSPBuildSolution(ksp,PETSC_NULL,&V);
1217:    or
1218:       KSPBuildSolution(ksp,v,PETSC_NULL); 
1219: .ve
1220:    In the first case an internal vector is allocated to store the solution
1221:    (the user cannot destroy this vector). In the second case the solution
1222:    is generated in the vector that the user provides. Note that for certain 
1223:    methods, such as KSPCG, the second case requires a copy of the solution,
1224:    while in the first case the call is essentially free since it simply 
1225:    returns the vector where the solution already is stored.

1227:    Level: advanced

1229: .keywords: KSP, build, solution

1231: .seealso: KSPGetSolution(), KSPBuildResidual()
1232: @*/
1233: int KSPBuildSolution(KSP ksp,Vec v,Vec *V)
1234: {

1239:   if (!V && !v) SETERRQ(PETSC_ERR_ARG_WRONG,"Must provide either v or V");
1240:   if (!V) V = &v;
1241:   (*ksp->ops->buildsolution)(ksp,v,V);
1242:   return(0);
1243: }

1245: #undef __FUNCT__  
1247: /*@C
1248:    KSPBuildResidual - Builds the residual in a vector provided.

1250:    Collective on KSP

1252:    Input Parameter:
1253: .  ksp - iterative context obtained from KSPCreate()

1255:    Output Parameters:
1256: +  v - optional location to stash residual.  If v is not provided,
1257:        then a location is generated.
1258: .  t - work vector.  If not provided then one is generated.
1259: -  V - the residual

1261:    Notes:
1262:    Regardless of whether or not v is provided, the residual is 
1263:    returned in V.

1265:    Level: advanced

1267: .keywords: KSP, build, residual

1269: .seealso: KSPBuildSolution()
1270: @*/
1271: int KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
1272: {
1273:   int flag = 0,ierr;
1274:   Vec w = v,tt = t;

1278:   if (!w) {
1279:     VecDuplicate(ksp->vec_rhs,&w);
1280:     PetscLogObjectParent((PetscObject)ksp,w);
1281:   }
1282:   if (!tt) {
1283:     VecDuplicate(ksp->vec_rhs,&tt); flag = 1;
1284:     PetscLogObjectParent((PetscObject)ksp,tt);
1285:   }
1286:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
1287:   if (flag) {VecDestroy(tt);}
1288:   return(0);
1289: }