Actual source code: itcl.c

petsc-3.15.0 2021-04-05
Report Typos and Errors

  2: /*
  3:     Code for setting KSP options from the options database.
  4: */

  6: #include <petsc/private/kspimpl.h>
  7: #include <petscdraw.h>

  9: /*@C
 10:    KSPSetOptionsPrefix - Sets the prefix used for searching for all
 11:    KSP options in the database.

 13:    Logically Collective on ksp

 15:    Input Parameters:
 16: +  ksp - the Krylov context
 17: -  prefix - the prefix string to prepend to all KSP option requests

 19:    Notes:
 20:    A hyphen (-) must NOT be given at the beginning of the prefix name.
 21:    The first character of all runtime options is AUTOMATICALLY the
 22:    hyphen.

 24:    For example, to distinguish between the runtime options for two
 25:    different KSP contexts, one could call
 26: .vb
 27:       KSPSetOptionsPrefix(ksp1,"sys1_")
 28:       KSPSetOptionsPrefix(ksp2,"sys2_")
 29: .ve

 31:    This would enable use of different options for each system, such as
 32: .vb
 33:       -sys1_ksp_type gmres -sys1_ksp_rtol 1.e-3
 34:       -sys2_ksp_type bcgs  -sys2_ksp_rtol 1.e-4
 35: .ve

 37:    Level: advanced

 39: .seealso: KSPAppendOptionsPrefix(), KSPGetOptionsPrefix()
 40: @*/
 41: PetscErrorCode  KSPSetOptionsPrefix(KSP ksp,const char prefix[])
 42: {

 47:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
 48:   PCSetOptionsPrefix(ksp->pc,prefix);
 49:   PetscObjectSetOptionsPrefix((PetscObject)ksp,prefix);
 50:   return(0);
 51: }

 53: /*@C
 54:    KSPAppendOptionsPrefix - Appends to the prefix used for searching for all
 55:    KSP options in the database.

 57:    Logically Collective on ksp

 59:    Input Parameters:
 60: +  ksp - the Krylov context
 61: -  prefix - the prefix string to prepend to all KSP option requests

 63:    Notes:
 64:    A hyphen (-) must NOT be given at the beginning of the prefix name.
 65:    The first character of all runtime options is AUTOMATICALLY the hyphen.

 67:    Level: advanced

 69: .seealso: KSPSetOptionsPrefix(), KSPGetOptionsPrefix()
 70: @*/
 71: PetscErrorCode  KSPAppendOptionsPrefix(KSP ksp,const char prefix[])
 72: {

 77:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
 78:   PCAppendOptionsPrefix(ksp->pc,prefix);
 79:   PetscObjectAppendOptionsPrefix((PetscObject)ksp,prefix);
 80:   return(0);
 81: }

 83: /*@
 84:    KSPSetUseFischerGuess - Use the Paul Fischer algorithm

 86:    Logically Collective on ksp

 88:    Input Parameters:
 89: +  ksp - the Krylov context
 90: .  model - use model 1, model 2 or any other number to turn it off
 91: -  size - size of subspace used to generate initial guess

 93:     Options Database:
 94: .   -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves

 96:    Level: advanced

 98: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetGuess(), KSPGetGuess()
 99: @*/
100: PetscErrorCode  KSPSetUseFischerGuess(KSP ksp,PetscInt model,PetscInt size)
101: {
102:   KSPGuess       guess;

109:   KSPGetGuess(ksp,&guess);
110:   KSPGuessSetType(guess,KSPGUESSFISCHER);
111:   KSPGuessFischerSetModel(guess,model,size);
112:   return(0);
113: }

115: /*@
116:    KSPSetGuess - Set the initial guess object

118:    Logically Collective on ksp

120:    Input Parameters:
121: +  ksp - the Krylov context
122: -  guess - the object created with KSPGuessCreate()

124:    Level: advanced

126:    Notes:
127:     this allows a single KSP to be used with several different initial guess generators (likely for different linear
128:           solvers, see KSPSetPC()).

130:           This increases the reference count of the guess object, you must destroy the object with KSPGuessDestroy()
131:           before the end of the program.

133: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPGetGuess()
134: @*/
135: PetscErrorCode  KSPSetGuess(KSP ksp,KSPGuess guess)
136: {

142:   PetscObjectReference((PetscObject)guess);
143:   KSPGuessDestroy(&ksp->guess);
144:   ksp->guess = guess;
145:   ksp->guess->ksp = ksp;
146:   return(0);
147: }

149: /*@
150:    KSPGetGuess - Gets the initial guess generator for the KSP.

152:    Not Collective

154:    Input Parameters:
155: .  ksp - the Krylov context

157:    Output Parameters:
158: .   guess - the object

160:    Level: developer

162: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetGuess()
163: @*/
164: PetscErrorCode  KSPGetGuess(KSP ksp,KSPGuess *guess)
165: {

171:   if (!ksp->guess) {
172:     const char* prefix;

174:     KSPGuessCreate(PetscObjectComm((PetscObject)ksp),&ksp->guess);
175:     PetscObjectGetOptionsPrefix((PetscObject)ksp,&prefix);
176:     if (prefix) {
177:       PetscObjectSetOptionsPrefix((PetscObject)ksp->guess,prefix);
178:     }
179:     ksp->guess->ksp = ksp;
180:   }
181:   *guess = ksp->guess;
182:   return(0);
183: }

185: /*@C
186:    KSPGetOptionsPrefix - Gets the prefix used for searching for all
187:    KSP options in the database.

189:    Not Collective

191:    Input Parameters:
192: .  ksp - the Krylov context

194:    Output Parameters:
195: .  prefix - pointer to the prefix string used is returned

197:    Notes:
198:     On the fortran side, the user should pass in a string 'prefix' of
199:    sufficient length to hold the prefix.

201:    Level: advanced

203: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix()
204: @*/
205: PetscErrorCode  KSPGetOptionsPrefix(KSP ksp,const char *prefix[])
206: {

211:   PetscObjectGetOptionsPrefix((PetscObject)ksp,prefix);
212:   return(0);
213: }

215: static PetscErrorCode PetscViewerAndFormatCreate_Internal(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
216: {

220:   PetscViewerAndFormatCreate(viewer, format, vf);
221:   (*vf)->data = ctx;
222:   return(0);
223: }

225: /*@C
226:    KSPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

228:    Collective on ksp

230:    Input Parameters:
231: +  ksp  - KSP object you wish to monitor
232: .  opt  - the command line option for this monitor
233: .  name - the monitor type one is seeking
234: -  ctx  - An optional user context for the monitor, or NULL

236:    Level: developer

238: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
239:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
240:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
241:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
242:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
243:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
244:           PetscOptionsFList(), PetscOptionsEList()
245: @*/
246: PetscErrorCode KSPMonitorSetFromOptions(KSP ksp, const char opt[], const char name[], void *ctx)
247: {
248:   PetscErrorCode      (*mfunc)(KSP, PetscInt, PetscReal, void *);
249:   PetscErrorCode      (*cfunc)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
250:   PetscErrorCode      (*dfunc)(PetscViewerAndFormat **);
251:   PetscViewerAndFormat *vf;
252:   PetscViewer           viewer;
253:   PetscViewerFormat     format;
254:   PetscViewerType       vtype;
255:   char                  key[PETSC_MAX_PATH_LEN];
256:   PetscBool             all, flg;
257:   const char           *prefix = NULL;
258:   PetscErrorCode        ierr;

261:   PetscStrcmp(opt, "-all_ksp_monitor", &all);
262:   if (!all) {PetscObjectGetOptionsPrefix((PetscObject) ksp, &prefix);}
263:   PetscOptionsGetViewer(PetscObjectComm((PetscObject) ksp), ((PetscObject) ksp)->options, prefix, opt, &viewer, &format, &flg);
264:   if (!flg) return(0);

266:   PetscViewerGetType(viewer, &vtype);
267:   KSPMonitorMakeKey_Internal(name, vtype, format, key);
268:   PetscFunctionListFind(KSPMonitorList, key, &mfunc);
269:   PetscFunctionListFind(KSPMonitorCreateList, key, &cfunc);
270:   PetscFunctionListFind(KSPMonitorDestroyList, key, &dfunc);
271:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
272:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

274:   (*cfunc)(viewer, format, ctx, &vf);
275:   PetscObjectDereference((PetscObject) viewer);
276:   KSPMonitorSet(ksp, mfunc, vf, (PetscErrorCode (*)(void **)) dfunc);
277:   return(0);
278: }

280: /*@
281:    KSPSetFromOptions - Sets KSP options from the options database.
282:    This routine must be called before KSPSetUp() if the user is to be
283:    allowed to set the Krylov type.

285:    Collective on ksp

287:    Input Parameters:
288: .  ksp - the Krylov space context

290:    Options Database Keys:
291: +   -ksp_max_it - maximum number of linear iterations
292: .   -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
293:                 if residual norm decreases by this factor than convergence is declared
294: .   -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual
295:                 norm is less than this then convergence is declared
296: .   -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
297: .   -ksp_converged_use_initial_residual_norm - see KSPConvergedDefaultSetUIRNorm()
298: .   -ksp_converged_use_min_initial_residual_norm - see KSPConvergedDefaultSetUMIRNorm()
299: .   -ksp_converged_maxits - see KSPConvergedDefaultSetConvergedMaxits()
300: .   -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
301:                        convergence test (say you always want to run with 5 iterations) to
302:                        save on communication overhead
303:                     preconditioned - default for left preconditioning
304:                     unpreconditioned - see KSPSetNormType()
305:                     natural - see KSPSetNormType()
306: .   -ksp_check_norm_iteration it - do not compute residual norm until iteration number it (does compute at 0th iteration)
307:        works only for PCBCGS, PCIBCGS and and PCCG
308: .   -ksp_lag_norm - compute the norm of the residual for the ith iteration on the i+1 iteration; this means that one can use
309:        the norm of the residual for convergence test WITHOUT an extra MPI_Allreduce() limiting global synchronizations.
310:        This will require 1 more iteration of the solver than usual.
311: .   -ksp_guess_type - Type of initial guess generator for repeated linear solves
312: .   -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
313: .   -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
314: .   -ksp_test_null_space - tests the null space set with MatSetNullSpace() to see if it truly is a null space
315: .   -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
316: .   -ksp_monitor_cancel - cancel all previous convergene monitor routines set
317: .   -ksp_monitor - print residual norm at each iteration
318: .   -ksp_monitor draw::draw_lg - plot residual norm at each iteration
319: .   -ksp_monitor_true_residual - print true residual norm at each iteration
320: .   -all_ksp_monitor <optional filename> - print residual norm at each iteration for ALL KSP solves, regardless of their prefix. This is
321:                                            useful for PCFIELDSPLIT, PCMG, etc that have inner solvers and you wish to track the convergence of all the solvers
322: .   -ksp_monitor_solution [ascii binary or draw][:filename][:format option] - plot solution at each iteration
323: .   -ksp_monitor_singular_value - monitor extreme singular values at each iteration
324: .   -ksp_converged_reason - view the convergence state at the end of the solve
325: .   -ksp_use_explicittranspose - transpose the system explicitly in KSPSolveTranspose
326: -   -ksp_converged_rate - view the convergence rate at the end of the solve

328:    Notes:
329:    To see all options, run your program with the -help option
330:    or consult Users-Manual: ch_ksp

332:    Level: beginner

334: .seealso: KSPSetOptionsPrefix(), KSPResetFromOptions(), KSPSetUseFischerGuess()

336: @*/
337: PetscErrorCode  KSPSetFromOptions(KSP ksp)
338: {
339:   const char     *convtests[]={"default","skip","lsqr"},*prefix;
340:   char           type[256],guesstype[256],monfilename[PETSC_MAX_PATH_LEN];
341:   PetscBool      flg,flag,reuse,set;
342:   PetscInt       indx,model[2]={0,0},nmax;
343:   KSPNormType    normtype;
344:   PCSide         pcside;
345:   void           *ctx;
346:   MPI_Comm       comm;

351:   PetscObjectGetComm((PetscObject) ksp, &comm);
352:   PetscObjectGetOptionsPrefix((PetscObject) ksp, &prefix);
353:   if (!ksp->skippcsetfromoptions) {
354:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
355:     PCSetFromOptions(ksp->pc);
356:   }

358:   KSPRegisterAll();
359:   PetscObjectOptionsBegin((PetscObject)ksp);
360:   PetscOptionsFList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(((PetscObject)ksp)->type_name ? ((PetscObject)ksp)->type_name : KSPGMRES),type,256,&flg);
361:   if (flg) {
362:     KSPSetType(ksp,type);
363:   }
364:   /*
365:     Set the type if it was never set.
366:   */
367:   if (!((PetscObject)ksp)->type_name) {
368:     KSPSetType(ksp,KSPGMRES);
369:   }

371:   KSPResetViewers(ksp);

373:   /* Cancels all monitors hardwired into code before call to KSPSetFromOptions() */
374:   PetscOptionsBool("-ksp_monitor_cancel","Remove any hardwired monitor routines","KSPMonitorCancel",PETSC_FALSE,&flg,&set);
375:   if (set && flg) {KSPMonitorCancel(ksp);}
376:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor", "preconditioned_residual", NULL);
377:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_short", "preconditioned_residual_short", NULL);
378:   KSPMonitorSetFromOptions(ksp, "-all_ksp_monitor", "preconditioned_residual", NULL);
379:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_range", "preconditioned_residual_range", NULL);
380:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_true_residual", "true_residual", NULL);
381:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_max", "true_residual_max", NULL);
382:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_solution", "solution", NULL);
383:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_singular_value", "singular_value", ksp);
384:   KSPMonitorSetFromOptions(ksp, "-ksp_monitor_error", "error", ksp);
385:   PetscOptionsBool("-ksp_monitor_pause_final", "Pauses all draw monitors at the final iterate", "KSPMonitorPauseFinal_Internal", PETSC_FALSE, &ksp->pauseFinal, NULL);

387:   PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
388:   if (flg) {
389:     KSPGetReusePreconditioner(ksp,&reuse);
390:     PetscOptionsBool("-ksp_reuse_preconditioner","Use initial preconditioner and don't ever compute a new one","KSPReusePreconditioner",reuse,&reuse,NULL);
391:     KSPSetReusePreconditioner(ksp,reuse);
392:     PetscOptionsBool("-ksp_error_if_not_converged","Generate error if solver does not converge","KSPSetErrorIfNotConverged",ksp->errorifnotconverged,&ksp->errorifnotconverged,NULL);
393:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view",&ksp->viewer, &ksp->format,&ksp->view);
394:     flg = PETSC_FALSE;
395:     PetscOptionsBool("-ksp_converged_reason_view_cancel","Cancel all the converged reason view functions set using KSPConvergedReasonViewSet","KSPConvergedReasonViewCancel",PETSC_FALSE,&flg,&set);
396:     if (set && flg) {
397:       KSPConvergedReasonViewCancel(ksp);
398:     }
399:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat",&ksp->viewerMat,&ksp->formatMat,&ksp->viewMat);
400:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_pmat",&ksp->viewerPMat,&ksp->formatPMat,&ksp->viewPMat);
401:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_rhs",&ksp->viewerRhs,&ksp->formatRhs,&ksp->viewRhs);
402:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_solution",&ksp->viewerSol,&ksp->formatSol,&ksp->viewSol);
403:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat_explicit",&ksp->viewerMatExp,&ksp->formatMatExp,&ksp->viewMatExp);
404:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_final_residual",&ksp->viewerFinalRes,&ksp->formatFinalRes,&ksp->viewFinalRes);
405:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_preconditioned_operator_explicit",&ksp->viewerPOpExp,&ksp->formatPOpExp,&ksp->viewPOpExp);
406:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_diagonal_scale",&ksp->viewerDScale,&ksp->formatDScale,&ksp->viewDScale);

408:     KSPGetDiagonalScale(ksp,&flag);
409:     PetscOptionsBool("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",flag,&flag,&flg);
410:     if (flg) {
411:       KSPSetDiagonalScale(ksp,flag);
412:     }
413:     KSPGetDiagonalScaleFix(ksp,&flag);
414:     PetscOptionsBool("-ksp_diagonal_scale_fix","Fix diagonally scaled matrix after solve","KSPSetDiagonalScaleFix",flag,&flag,&flg);
415:     if (flg) {
416:       KSPSetDiagonalScaleFix(ksp,flag);
417:     }
418:     nmax = ksp->nmax;
419:     PetscOptionsDeprecated("-ksp_matsolve_block_size","-ksp_matsolve_batch_size","3.15",NULL);
420:     PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg);
421:     if (flg) {
422:       KSPSetMatSolveBatchSize(ksp, nmax);
423:     }
424:     goto skipoptions;
425:   }

427:   PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,NULL);
428:   PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,NULL);
429:   PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->abstol,&ksp->abstol,NULL);
430:   PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,NULL);

432:   PetscOptionsBool("-ksp_converged_use_initial_residual_norm","Use initial residual norm for computing relative convergence","KSPConvergedDefaultSetUIRNorm",PETSC_FALSE,&flag,&set);
433:   if (set && flag) {KSPConvergedDefaultSetUIRNorm(ksp);}
434:   PetscOptionsBool("-ksp_converged_use_min_initial_residual_norm","Use minimum of initial residual norm and b for computing relative convergence","KSPConvergedDefaultSetUMIRNorm",PETSC_FALSE,&flag,&set);
435:   if (set && flag) {KSPConvergedDefaultSetUMIRNorm(ksp);}
436:   PetscOptionsBool("-ksp_converged_maxits","Declare convergence if the maximum number of iterations is reached","KSPConvergedDefaultSetConvergedMaxits",PETSC_FALSE,&flag,&set);
437:   if (set) {KSPConvergedDefaultSetConvergedMaxits(ksp,flag);}
438:   PetscOptionsBool("-ksp_initial_guess_nonzero","Use the contents of the solution vector for initial guess","KSPSetInitialNonzero",ksp->guess_zero ? PETSC_FALSE : PETSC_TRUE,&flag,&flg);
439:   if (flg) {
440:     KSPSetInitialGuessNonzero(ksp,flag);
441:   }
442:   KSPGetReusePreconditioner(ksp,&reuse);
443:   PetscOptionsBool("-ksp_reuse_preconditioner","Use initial preconditioner and don't ever compute a new one","KSPReusePreconditioner",reuse,&reuse,NULL);
444:   KSPSetReusePreconditioner(ksp,reuse);

446:   PetscOptionsBool("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,&ksp->guess_knoll,NULL);
447:   PetscOptionsBool("-ksp_error_if_not_converged","Generate error if solver does not converge","KSPSetErrorIfNotConverged",ksp->errorifnotconverged,&ksp->errorifnotconverged,NULL);
448:   PetscOptionsFList("-ksp_guess_type","Initial guess in Krylov method",NULL,KSPGuessList,NULL,guesstype,256,&flg);
449:   if (flg) {
450:     KSPGetGuess(ksp,&ksp->guess);
451:     KSPGuessSetType(ksp->guess,guesstype);
452:     KSPGuessSetFromOptions(ksp->guess);
453:   } else { /* old option for KSP */
454:     nmax = 2;
455:     PetscOptionsIntArray("-ksp_fischer_guess","Use Paul Fischer's algorithm for initial guess","KSPSetUseFischerGuess",model,&nmax,&flag);
456:     if (flag) {
457:       if (nmax != 2) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in model,size as arguments");
458:       KSPSetUseFischerGuess(ksp,model[0],model[1]);
459:     }
460:   }

462:   PetscOptionsEList("-ksp_convergence_test","Convergence test","KSPSetConvergenceTest",convtests,3,"default",&indx,&flg);
463:   if (flg) {
464:     switch (indx) {
465:     case 0:
466:       KSPConvergedDefaultCreate(&ctx);
467:       KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
468:       break;
469:     case 1:
470:       KSPSetConvergenceTest(ksp,KSPConvergedSkip,NULL,NULL);
471:       break;
472:     case 2:
473:       KSPConvergedDefaultCreate(&ctx);
474:       KSPSetConvergenceTest(ksp,KSPLSQRConvergedDefault,ctx,KSPConvergedDefaultDestroy);
475:       break;
476:     }
477:   }

479:   KSPSetUpNorms_Private(ksp,PETSC_FALSE,&normtype,NULL);
480:   PetscOptionsEnum("-ksp_norm_type","KSP Norm type","KSPSetNormType",KSPNormTypes,(PetscEnum)normtype,(PetscEnum*)&normtype,&flg);
481:   if (flg) { KSPSetNormType(ksp,normtype); }

483:   PetscOptionsInt("-ksp_check_norm_iteration","First iteration to compute residual norm","KSPSetCheckNormIteration",ksp->chknorm,&ksp->chknorm,NULL);

485:   PetscOptionsBool("-ksp_lag_norm","Lag the calculation of the residual norm","KSPSetLagNorm",ksp->lagnorm,&flag,&flg);
486:   if (flg) {
487:     KSPSetLagNorm(ksp,flag);
488:   }

490:   KSPGetDiagonalScale(ksp,&flag);
491:   PetscOptionsBool("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",flag,&flag,&flg);
492:   if (flg) {
493:     KSPSetDiagonalScale(ksp,flag);
494:   }
495:   KSPGetDiagonalScaleFix(ksp,&flag);
496:   PetscOptionsBool("-ksp_diagonal_scale_fix","Fix diagonally scaled matrix after solve","KSPSetDiagonalScaleFix",flag,&flag,&flg);
497:   if (flg) {
498:     KSPSetDiagonalScaleFix(ksp,flag);
499:   }

501:   PetscOptionsBool("-ksp_constant_null_space","Add constant null space to Krylov solver matrix","MatSetNullSpace",PETSC_FALSE,&flg,&set);
502:   if (set && flg) {
503:     MatNullSpace nsp;
504:     Mat          Amat = NULL;

506:     MatNullSpaceCreate(comm,PETSC_TRUE,0,NULL,&nsp);
507:     if (ksp->pc) { PCGetOperators(ksp->pc,&Amat,NULL); }
508:     if (Amat) {
509:       MatSetNullSpace(Amat,nsp);
510:       MatNullSpaceDestroy(&nsp);
511:     } else SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set nullspace, matrix has not yet been provided");
512:   }

514:   flg = PETSC_FALSE;
515:   if (ksp->pc) {
516:     PetscObjectTypeCompare((PetscObject)ksp->pc,PCKSP,&flg);
517:     if (!flg) PetscObjectTypeCompare((PetscObject)ksp->pc,PCBJACOBI,&flg);
518:     if (!flg) PetscObjectTypeCompare((PetscObject)ksp->pc,PCDEFLATION,&flg);
519:   }

521:   if (flg) {
522:     /* A hack for using dynamic tolerance in preconditioner */
523:     PetscOptionsString("-sub_ksp_dynamic_tolerance","Use dynamic tolerance for PC if PC is a KSP","KSPMonitorDynamicTolerance","stdout",monfilename,sizeof(monfilename),&flg);
524:     if (flg) {
525:       KSPDynTolCtx *scale;
526:       PetscMalloc1(1,&scale);
527:       scale->bnrm = -1.0;
528:       scale->coef = 1.0;
529:       PetscOptionsReal("-sub_ksp_dynamic_tolerance_param","Parameter of dynamic tolerance for inner PCKSP","KSPMonitorDynamicToleranceParam",scale->coef,&scale->coef,&flg);
530:       KSPMonitorSet(ksp,KSPMonitorDynamicTolerance,scale,KSPMonitorDynamicToleranceDestroy);
531:     }
532:   }

534:   /*
535:    Calls Python function
536:   */
537:   PetscOptionsString("-ksp_monitor_python","Use Python function","KSPMonitorSet",NULL,monfilename,sizeof(monfilename),&flg);
538:   if (flg) {PetscPythonMonitorSet((PetscObject)ksp,monfilename);}
539:   /*
540:     Graphically plots preconditioned residual norm and range of residual element values
541:   */
542:   PetscOptionsBool("-ksp_monitor_lg_range","Monitor graphically range of preconditioned residual norm","KSPMonitorSet",PETSC_FALSE,&flg,&set);
543:   if (set && flg) {
544:     PetscViewer ctx;

546:     PetscViewerDrawOpen(comm,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
547:     KSPMonitorSet(ksp,KSPMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
548:   }
549:   /* TODO Do these show up in help? */
550:   PetscOptionsHasName(((PetscObject) ksp)->options, prefix, "-ksp_converged_rate", &flg);
551:   if (flg) {
552:     const char *RateTypes[] = {"default", "residual", "error", "PetscRateType", "RATE_", NULL};
553:     PetscEnum rtype = (PetscEnum) 1;

555:     PetscOptionsGetEnum(((PetscObject) ksp)->options, prefix, "-ksp_converged_rate_type", RateTypes, &rtype, &flg);
556:     if (rtype == (PetscEnum) 0 || rtype == (PetscEnum) 1) {KSPSetResidualHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE);}
557:     if (rtype == (PetscEnum) 0 || rtype == (PetscEnum) 2) {KSPSetErrorHistory(ksp, NULL, PETSC_DETERMINE, PETSC_TRUE);}
558:   }
559:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view",&ksp->viewer,&ksp->format,&ksp->view);
560:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_pre",&ksp->viewerPre,&ksp->formatPre,&ksp->viewPre);

562:   flg = PETSC_FALSE;
563:   PetscOptionsBool("-ksp_converged_reason_view_cancel","Cancel all the converged reason view functions set using KSPConvergedReasonViewSet","KSPConvergedReasonViewCancel",PETSC_FALSE,&flg,&set);
564:   if (set && flg) {
565:     KSPConvergedReasonViewCancel(ksp);
566:   }
567:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_converged_rate",&ksp->viewerRate,&ksp->formatRate,&ksp->viewRate);
568:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat",&ksp->viewerMat,&ksp->formatMat,&ksp->viewMat);
569:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_pmat",&ksp->viewerPMat,&ksp->formatPMat,&ksp->viewPMat);
570:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_rhs",&ksp->viewerRhs,&ksp->formatRhs,&ksp->viewRhs);
571:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_solution",&ksp->viewerSol,&ksp->formatSol,&ksp->viewSol);
572:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat_explicit",&ksp->viewerMatExp,&ksp->formatMatExp,&ksp->viewMatExp);
573:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_eigenvalues",&ksp->viewerEV,&ksp->formatEV,&ksp->viewEV);
574:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_singularvalues",&ksp->viewerSV,&ksp->formatSV,&ksp->viewSV);
575:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_eigenvalues_explicit",&ksp->viewerEVExp,&ksp->formatEVExp,&ksp->viewEVExp);
576:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_final_residual",&ksp->viewerFinalRes,&ksp->formatFinalRes,&ksp->viewFinalRes);
577:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_preconditioned_operator_explicit",&ksp->viewerPOpExp,&ksp->formatPOpExp,&ksp->viewPOpExp);
578:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_diagonal_scale",&ksp->viewerDScale,&ksp->formatDScale,&ksp->viewDScale);

580:   /* Deprecated options */
581:   if (!ksp->viewEV) {
582:     PetscOptionsDeprecated("-ksp_compute_eigenvalues",NULL,"3.9","Use -ksp_view_eigenvalues");
583:     PetscOptionsGetViewer(comm, ((PetscObject) ksp)->options,prefix, "-ksp_compute_eigenvalues",&ksp->viewerEV,&ksp->formatEV,&ksp->viewEV);
584:   }
585:   if (!ksp->viewEV) {
586:     PetscOptionsDeprecated("-ksp_plot_eigenvalues",NULL,"3.9","Use -ksp_view_eigenvalues draw");
587:     PetscOptionsName("-ksp_plot_eigenvalues", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw]", "KSPView", &ksp->viewEV);
588:     if (ksp->viewEV) {
589:       ksp->formatEV = PETSC_VIEWER_DEFAULT;
590:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
591:       PetscObjectReference((PetscObject) ksp->viewerEV);
592:     }
593:   }
594:   if (!ksp->viewEV) {
595:     PetscOptionsDeprecated("-ksp_plot_eigencontours",NULL,"3.9","Use -ksp_view_eigenvalues draw::draw_contour");
596:     PetscOptionsName("-ksp_plot_eigencontours", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw::draw_contour]", "KSPView", &ksp->viewEV);
597:     if (ksp->viewEV) {
598:       ksp->formatEV = PETSC_VIEWER_DRAW_CONTOUR;
599:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
600:       PetscObjectReference((PetscObject) ksp->viewerEV);
601:     }
602:   }
603:   if (!ksp->viewEVExp) {
604:     PetscOptionsDeprecated("-ksp_compute_eigenvalues_explicitly",NULL,"3.9","Use -ksp_view_eigenvalues_explicit");
605:     PetscOptionsGetViewer(comm, ((PetscObject) ksp)->options,prefix, "-ksp_compute_eigenvalues_explicitly",&ksp->viewerEVExp,&ksp->formatEVExp,&ksp->viewEVExp);
606:   }
607:   if (!ksp->viewEVExp) {
608:     PetscOptionsDeprecated("-ksp_plot_eigenvalues_explicitly",NULL,"3.9","Use -ksp_view_eigenvalues_explicit draw");
609:     PetscOptionsName("-ksp_plot_eigenvalues_explicitly","[deprecated since PETSc 3.9; use -ksp_view_eigenvalues_explicit draw]","KSPView",&ksp->viewEVExp);
610:     if (ksp->viewEVExp) {
611:       ksp->formatEVExp = PETSC_VIEWER_DEFAULT;
612:       ksp->viewerEVExp = PETSC_VIEWER_DRAW_(comm);
613:       PetscObjectReference((PetscObject) ksp->viewerEVExp);
614:     }
615:   }
616:   if (!ksp->viewSV) {
617:     PetscOptionsDeprecated("-ksp_compute_singularvalues",NULL,"3.9","Use -ksp_view_singularvalues");
618:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_compute_singularvalues",&ksp->viewerSV,&ksp->formatSV,&ksp->viewSV);
619:   }
620:   if (!ksp->viewFinalRes) {
621:     PetscOptionsDeprecated("-ksp_final_residual",NULL,"3.9","Use -ksp_view_final_residual");
622:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_final_residual",&ksp->viewerFinalRes,&ksp->formatFinalRes,&ksp->viewFinalRes);
623:   }

625: #if defined(PETSC_HAVE_SAWS)
626:   /*
627:     Publish convergence information using AMS
628:   */
629:   PetscOptionsBool("-ksp_monitor_saws","Publish KSP progress using SAWs","KSPMonitorSet",PETSC_FALSE,&flg,&set);
630:   if (set && flg) {
631:     void *ctx;
632:     KSPMonitorSAWsCreate(ksp,&ctx);
633:     KSPMonitorSet(ksp,KSPMonitorSAWs,ctx,KSPMonitorSAWsDestroy);
634:     KSPSetComputeSingularValues(ksp,PETSC_TRUE);
635:   }
636: #endif

638:   /* -----------------------------------------------------------------------*/
639:   KSPSetUpNorms_Private(ksp,PETSC_FALSE,NULL,&pcside);
640:   PetscOptionsEnum("-ksp_pc_side","KSP preconditioner side","KSPSetPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
641:   if (flg) {KSPSetPCSide(ksp,pcside);}

643:   if (ksp->viewSV || ksp->viewEV) {
644:     KSPSetComputeSingularValues(ksp,PETSC_TRUE);
645:   }

647: #if defined(PETSC_HAVE_SAWS)
648:   {
649:     PetscBool set;
650:     flg  = PETSC_FALSE;
651:     PetscOptionsBool("-ksp_saws_block","Block for SAWs at end of KSPSolve","PetscObjectSAWsBlock",((PetscObject)ksp)->amspublishblock,&flg,&set);
652:     if (set) {
653:       PetscObjectSAWsSetBlock((PetscObject)ksp,flg);
654:     }
655:   }
656: #endif

658:   nmax = ksp->nmax;
659:   PetscOptionsDeprecated("-ksp_matsolve_block_size","-ksp_matsolve_batch_size","3.15",NULL);
660:   PetscOptionsInt("-ksp_matsolve_batch_size", "Maximum number of columns treated simultaneously", "KSPSetMatSolveBatchSize", nmax, &nmax, &flg);
661:   if (flg) {
662:     KSPSetMatSolveBatchSize(ksp, nmax);
663:   }

665:   flg  = PETSC_FALSE;
666:   PetscOptionsBool("-ksp_use_explicittranspose","Explicitly tranpose the system in KSPSolveTranspose","KSPSetUseExplicitTranspose",ksp->transpose.use_explicittranspose,&flg,&set);
667:   if (set) {
668:     KSPSetUseExplicitTranspose(ksp,flg);
669:   }

671:   if (ksp->ops->setfromoptions) {
672:     (*ksp->ops->setfromoptions)(PetscOptionsObject,ksp);
673:   }
674:   skipoptions:
675:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
676:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)ksp);
677:   PetscOptionsEnd();
678:   ksp->setfromoptionscalled++;
679:   return(0);
680: }

682: /*@
683:    KSPResetFromOptions - Sets various KSP parameters from user options ONLY if the KSP was previously set from options

685:    Collective on ksp

687:    Input Parameter:
688: .  ksp - the KSP context

690:    Level: beginner

692: .seealso: KSPSetFromOptions(), KSPSetOptionsPrefix()
693: @*/
694: PetscErrorCode KSPResetFromOptions(KSP ksp)
695: {

699:   if (ksp->setfromoptionscalled) {KSPSetFromOptions(ksp);}
700:   return(0);
701: }