Actual source code: itcl.c

petsc-3.12.5 2020-03-29
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: static PetscErrorCode KSPSetupMonitor_Private(KSP ksp, PetscViewer viewer, PetscViewerFormat format, PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*), PetscBool useMonitor)
 10: {

 14:   if (useMonitor) {
 15:     PetscViewerAndFormat *vf;

 17:     PetscViewerAndFormatCreate(viewer, format, &vf);
 18:     PetscObjectDereference((PetscObject) viewer);
 19:     KSPMonitorSet(ksp, monitor, vf, (PetscErrorCode (*)(void**)) PetscViewerAndFormatDestroy);
 20:   }
 21:   return(0);
 22: }

 24: /*@C
 25:    KSPSetOptionsPrefix - Sets the prefix used for searching for all
 26:    KSP options in the database.

 28:    Logically Collective on ksp

 30:    Input Parameters:
 31: +  ksp - the Krylov context
 32: -  prefix - the prefix string to prepend to all KSP option requests

 34:    Notes:
 35:    A hyphen (-) must NOT be given at the beginning of the prefix name.
 36:    The first character of all runtime options is AUTOMATICALLY the
 37:    hyphen.

 39:    For example, to distinguish between the runtime options for two
 40:    different KSP contexts, one could call
 41: .vb
 42:       KSPSetOptionsPrefix(ksp1,"sys1_")
 43:       KSPSetOptionsPrefix(ksp2,"sys2_")
 44: .ve

 46:    This would enable use of different options for each system, such as
 47: .vb
 48:       -sys1_ksp_type gmres -sys1_ksp_rtol 1.e-3
 49:       -sys2_ksp_type bcgs  -sys2_ksp_rtol 1.e-4
 50: .ve

 52:    Level: advanced

 54: .seealso: KSPAppendOptionsPrefix(), KSPGetOptionsPrefix()
 55: @*/
 56: PetscErrorCode  KSPSetOptionsPrefix(KSP ksp,const char prefix[])
 57: {

 62:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
 63:   PCSetOptionsPrefix(ksp->pc,prefix);
 64:   PetscObjectSetOptionsPrefix((PetscObject)ksp,prefix);
 65:   return(0);
 66: }

 68: /*@C
 69:    KSPAppendOptionsPrefix - Appends to the prefix used for searching for all
 70:    KSP options in the database.

 72:    Logically Collective on ksp

 74:    Input Parameters:
 75: +  ksp - the Krylov context
 76: -  prefix - the prefix string to prepend to all KSP option requests

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

 82:    Level: advanced

 84: .seealso: KSPSetOptionsPrefix(), KSPGetOptionsPrefix()
 85: @*/
 86: PetscErrorCode  KSPAppendOptionsPrefix(KSP ksp,const char prefix[])
 87: {

 92:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
 93:   PCAppendOptionsPrefix(ksp->pc,prefix);
 94:   PetscObjectAppendOptionsPrefix((PetscObject)ksp,prefix);
 95:   return(0);
 96: }

 98: /*@
 99:    KSPSetUseFischerGuess - Use the Paul Fischer algorithm

101:    Logically Collective on ksp

103:    Input Parameters:
104: +  ksp - the Krylov context
105: .  model - use model 1, model 2 or any other number to turn it off
106: -  size - size of subspace used to generate initial guess

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

111:    Level: advanced

113: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetGuess(), KSPGetGuess()
114: @*/
115: PetscErrorCode  KSPSetUseFischerGuess(KSP ksp,PetscInt model,PetscInt size)
116: {
117:   KSPGuess       guess;

124:   KSPGetGuess(ksp,&guess);
125:   KSPGuessSetType(guess,KSPGUESSFISCHER);
126:   KSPGuessFischerSetModel(guess,model,size);
127:   return(0);
128: }

130: /*@
131:    KSPSetGuess - Set the initial guess object

133:    Logically Collective on ksp

135:    Input Parameters:
136: +  ksp - the Krylov context
137: -  guess - the object created with KSPGuessCreate()

139:    Level: advanced

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

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

148: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPGetGuess()
149: @*/
150: PetscErrorCode  KSPSetGuess(KSP ksp,KSPGuess guess)
151: {

157:   PetscObjectReference((PetscObject)guess);
158:   KSPGuessDestroy(&ksp->guess);
159:   ksp->guess = guess;
160:   ksp->guess->ksp = ksp;
161:   return(0);
162: }

164: /*@
165:    KSPGetGuess - Gets the initial guess generator for the KSP.

167:    Not Collective

169:    Input Parameters:
170: .  ksp - the Krylov context

172:    Output Parameters:
173: .   guess - the object

175:    Level: developer

177: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix(), KSPSetUseFischerGuess(), KSPSetGuess()
178: @*/
179: PetscErrorCode  KSPGetGuess(KSP ksp,KSPGuess *guess)
180: {

186:   if (!ksp->guess) {
187:     const char* prefix;

189:     KSPGuessCreate(PetscObjectComm((PetscObject)ksp),&ksp->guess);
190:     PetscObjectGetOptionsPrefix((PetscObject)ksp,&prefix);
191:     if (prefix) {
192:       PetscObjectSetOptionsPrefix((PetscObject)ksp->guess,prefix);
193:     }
194:     ksp->guess->ksp = ksp;
195:   }
196:   *guess = ksp->guess;
197:   return(0);
198: }

200: /*@C
201:    KSPGetOptionsPrefix - Gets the prefix used for searching for all
202:    KSP options in the database.

204:    Not Collective

206:    Input Parameters:
207: .  ksp - the Krylov context

209:    Output Parameters:
210: .  prefix - pointer to the prefix string used is returned

212:    Notes:
213:     On the fortran side, the user should pass in a string 'prefix' of
214:    sufficient length to hold the prefix.

216:    Level: advanced

218: .seealso: KSPSetOptionsPrefix(), KSPAppendOptionsPrefix()
219: @*/
220: PetscErrorCode  KSPGetOptionsPrefix(KSP ksp,const char *prefix[])
221: {

226:   PetscObjectGetOptionsPrefix((PetscObject)ksp,prefix);
227:   return(0);
228: }

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

233:    Collective on ksp

235:    Input Parameters:
236: +  ksp - KSP object you wish to monitor
237: .  name - the monitor type one is seeking
238: .  help - message indicating what monitoring is done
239: .  manual - manual page for the monitor
240: -  monitor - the monitor function, the context for this object is a PetscViewerAndFormat

242:    Level: developer

244: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
245:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
246:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
247:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
248:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
249:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
250:           PetscOptionsFList(), PetscOptionsEList()
251: @*/
252: PetscErrorCode  KSPMonitorSetFromOptions(KSP ksp,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*))
253: {
254:   PetscErrorCode       ierr;
255:   PetscBool            flg;
256:   PetscViewer          viewer;
257:   PetscViewerFormat    format;

260:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,name,&viewer,&format,&flg);
261:   KSPSetupMonitor_Private(ksp, viewer, format, (PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*)) monitor, flg);
262:   return(0);
263: }

265: /*@
266:    KSPSetFromOptions - Sets KSP options from the options database.
267:    This routine must be called before KSPSetUp() if the user is to be
268:    allowed to set the Krylov type.

270:    Collective on ksp

272:    Input Parameters:
273: .  ksp - the Krylov space context

275:    Options Database Keys:
276: +   -ksp_max_it - maximum number of linear iterations
277: .   -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
278:                 if residual norm decreases by this factor than convergence is declared
279: .   -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual
280:                 norm is less than this then convergence is declared
281: .   -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
282: .   -ksp_converged_use_initial_residual_norm - see KSPConvergedDefaultSetUIRNorm()
283: .   -ksp_converged_use_min_initial_residual_norm - see KSPConvergedDefaultSetUMIRNorm()
284: .   -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
285:                        convergence test (say you always want to run with 5 iterations) to
286:                        save on communication overhead
287:                     preconditioned - default for left preconditioning
288:                     unpreconditioned - see KSPSetNormType()
289:                     natural - see KSPSetNormType()
290: .   -ksp_check_norm_iteration it - do not compute residual norm until iteration number it (does compute at 0th iteration)
291:        works only for PCBCGS, PCIBCGS and and PCCG
292: .   -ksp_lag_norm - compute the norm of the residual for the ith iteration on the i+1 iteration; this means that one can use
293:        the norm of the residual for convergence test WITHOUT an extra MPI_Allreduce() limiting global synchronizations.
294:        This will require 1 more iteration of the solver than usual.
295: .   -ksp_guess_type - Type of initial guess generator for repeated linear solves
296: .   -ksp_fischer_guess <model,size> - uses the Fischer initial guess generator for repeated linear solves
297: .   -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
298: .   -ksp_test_null_space - tests the null space set with MatSetNullSpace() to see if it truly is a null space
299: .   -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
300: .   -ksp_monitor_cancel - cancel all previous convergene monitor routines set
301: .   -ksp_monitor <optional filename> - print residual norm at each iteration
302: .   -ksp_monitor_lg_residualnorm - plot residual norm at each iteration
303: .   -ksp_monitor_solution [ascii binary or draw][:filename][:format option] - plot solution at each iteration
304: -   -ksp_monitor_singular_value - monitor extreme singular values at each iteration

306:    Notes:
307:    To see all options, run your program with the -help option
308:    or consult Users-Manual: Chapter 4 KSP: Linear System Solvers

310:    Level: beginner

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

314: @*/
315: PetscErrorCode  KSPSetFromOptions(KSP ksp)
316: {
317:   PetscInt       indx;
318:   const char     *convtests[] = {"default","skip","lsqr"};
319:   char           type[256], guesstype[256], monfilename[PETSC_MAX_PATH_LEN];
320:   PetscBool      flg,flag,reuse,set;
321:   PetscInt       model[2]={0,0},nmax;
322:   KSPNormType    normtype;
323:   PCSide         pcside;
324:   void           *ctx;
325:   MPI_Comm       comm;
326:   const char    *prefix;

331:   PetscObjectGetComm((PetscObject) ksp, &comm);
332:   PetscObjectGetOptionsPrefix((PetscObject) ksp, &prefix);
333:   if (!ksp->skippcsetfromoptions) {
334:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
335:     PCSetFromOptions(ksp->pc);
336:   }

338:   KSPRegisterAll();
339:   PetscObjectOptionsBegin((PetscObject)ksp);
340:   PetscOptionsFList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(((PetscObject)ksp)->type_name ? ((PetscObject)ksp)->type_name : KSPGMRES),type,256,&flg);
341:   if (flg) {
342:     KSPSetType(ksp,type);
343:   }
344:   /*
345:     Set the type if it was never set.
346:   */
347:   if (!((PetscObject)ksp)->type_name) {
348:     KSPSetType(ksp,KSPGMRES);
349:   }

351:   KSPResetViewers(ksp);

353:   PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
354:   if (flg) {
355:     PCGetReusePreconditioner(ksp->pc,&reuse);
356:     PetscOptionsBool("-ksp_error_if_not_converged","Generate error if solver does not converge","KSPSetErrorIfNotConverged",ksp->errorifnotconverged,&ksp->errorifnotconverged,NULL);
357:     PetscOptionsBool("-ksp_reuse_preconditioner","Use initial preconditioner and don't ever compute a new one ","KSPReusePreconditioner",reuse,&reuse,NULL);
358:     KSPSetReusePreconditioner(ksp,reuse);
359:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view",&ksp->viewer, &ksp->format,&ksp->view);
360:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_converged_reason",&ksp->viewerReason,&ksp->formatReason,&ksp->viewReason);
361:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat",&ksp->viewerMat,&ksp->formatMat,&ksp->viewMat);
362:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_pmat",&ksp->viewerPMat,&ksp->formatPMat,&ksp->viewPMat);
363:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_rhs",&ksp->viewerRhs,&ksp->formatRhs,&ksp->viewRhs);
364:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_solution",&ksp->viewerSol,&ksp->formatSol,&ksp->viewSol);
365:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat_explicit",&ksp->viewerMatExp,&ksp->formatMatExp,&ksp->viewMatExp);
366:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_final_residual",&ksp->viewerFinalRes,&ksp->formatFinalRes,&ksp->viewFinalRes);
367:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_preconditioned_operator_explicit",&ksp->viewerPOpExp,&ksp->formatPOpExp,&ksp->viewPOpExp);
368:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_diagonal_scale",&ksp->viewerDScale,&ksp->formatDScale,&ksp->viewDScale);

370:     KSPGetDiagonalScale(ksp,&flag);
371:     PetscOptionsBool("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",flag,&flag,&flg);
372:     if (flg) {
373:       KSPSetDiagonalScale(ksp,flag);
374:     }
375:     KSPGetDiagonalScaleFix(ksp,&flag);
376:     PetscOptionsBool("-ksp_diagonal_scale_fix","Fix diagonally scaled matrix after solve","KSPSetDiagonalScaleFix",flag,&flag,&flg);
377:     if (flg) {
378:       KSPSetDiagonalScaleFix(ksp,flag);
379:     }
380:     goto skipoptions;
381:   }

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

388:   PetscOptionsBool("-ksp_converged_use_initial_residual_norm","Use initial residual norm for computing relative convergence","KSPConvergedDefaultSetUIRNorm",PETSC_FALSE,&flag,&set);
389:   if (set && flag) {KSPConvergedDefaultSetUIRNorm(ksp);}
390:   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);
391:   if (set && flag) {KSPConvergedDefaultSetUMIRNorm(ksp);}
392:   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);
393:   if (flg) {
394:     KSPSetInitialGuessNonzero(ksp,flag);
395:   }
396:   PCGetReusePreconditioner(ksp->pc,&reuse);
397:   PetscOptionsBool("-ksp_reuse_preconditioner","Use initial preconditioner and don't ever compute a new one ","KSPReusePreconditioner",reuse,&reuse,NULL);
398:   KSPSetReusePreconditioner(ksp,reuse);

400:   PetscOptionsBool("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,&ksp->guess_knoll,NULL);
401:   PetscOptionsBool("-ksp_error_if_not_converged","Generate error if solver does not converge","KSPSetErrorIfNotConverged",ksp->errorifnotconverged,&ksp->errorifnotconverged,NULL);
402:   PetscOptionsFList("-ksp_guess_type","Initial guess in Krylov method",NULL,KSPGuessList,NULL,guesstype,256,&flg);
403:   if (flg) {
404:     KSPGetGuess(ksp,&ksp->guess);
405:     KSPGuessSetType(ksp->guess,guesstype);
406:     KSPGuessSetFromOptions(ksp->guess);
407:   } else { /* old option for KSP */
408:     nmax = 2;
409:     PetscOptionsIntArray("-ksp_fischer_guess","Use Paul Fischer's algorithm for initial guess","KSPSetUseFischerGuess",model,&nmax,&flag);
410:     if (flag) {
411:       if (nmax != 2) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must pass in model,size as arguments");
412:       KSPSetUseFischerGuess(ksp,model[0],model[1]);
413:     }
414:   }

416:   PetscOptionsEList("-ksp_convergence_test","Convergence test","KSPSetConvergenceTest",convtests,3,"default",&indx,&flg);
417:   if (flg) {
418:     switch (indx) {
419:     case 0:
420:       KSPConvergedDefaultCreate(&ctx);
421:       KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
422:       break;
423:     case 1: KSPSetConvergenceTest(ksp,KSPConvergedSkip,NULL,NULL);    break;
424:     case 2:
425:       KSPConvergedDefaultCreate(&ctx);
426:       KSPSetConvergenceTest(ksp,KSPLSQRConvergedDefault,ctx,KSPConvergedDefaultDestroy);
427:       break;
428:     }
429:   }

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

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

437:   PetscOptionsBool("-ksp_lag_norm","Lag the calculation of the residual norm","KSPSetLagNorm",ksp->lagnorm,&flag,&flg);
438:   if (flg) {
439:     KSPSetLagNorm(ksp,flag);
440:   }

442:   KSPGetDiagonalScale(ksp,&flag);
443:   PetscOptionsBool("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",flag,&flag,&flg);
444:   if (flg) {
445:     KSPSetDiagonalScale(ksp,flag);
446:   }
447:   KSPGetDiagonalScaleFix(ksp,&flag);
448:   PetscOptionsBool("-ksp_diagonal_scale_fix","Fix diagonally scaled matrix after solve","KSPSetDiagonalScaleFix",flag,&flag,&flg);
449:   if (flg) {
450:     KSPSetDiagonalScaleFix(ksp,flag);
451:   }

453:   PetscOptionsBool("-ksp_constant_null_space","Add constant null space to Krylov solver matrix","MatSetNullSpace",PETSC_FALSE,&flg,&set);
454:   if (set && flg) {
455:     MatNullSpace nsp;
456:     Mat          Amat;

458:     MatNullSpaceCreate(comm,PETSC_TRUE,0,NULL,&nsp);
459:     PCGetOperators(ksp->pc,&Amat,NULL);
460:     if (Amat) {
461:       MatSetNullSpace(Amat,nsp);
462:       MatNullSpaceDestroy(&nsp);
463:     } else SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set nullspace, matrix has not yet been provided");
464:   }

466:   PetscOptionsBool("-ksp_monitor_cancel","Remove any hardwired monitor routines","KSPMonitorCancel",PETSC_FALSE,&flg,&set);
467:   /* -----------------------------------------------------------------------*/
468:   /*
469:     Cancels all monitors hardwired into code before call to KSPSetFromOptions()
470:   */
471:   if (set && flg) {
472:     KSPMonitorCancel(ksp);
473:   }
474:   KSPMonitorSetFromOptions(ksp,"-ksp_monitor","Monitor the (preconditioned) residual norm","KSPMonitorDefault",KSPMonitorDefault);
475:   KSPMonitorSetFromOptions(ksp,"-ksp_monitor_range","Monitor the percentage of large entries in the residual","KSPMonitorRange",KSPMonitorRange);
476:   KSPMonitorSetFromOptions(ksp,"-ksp_monitor_true_residual","Monitor the unpreconditioned residual norm","KSPMOnitorTrueResidual",KSPMonitorTrueResidualNorm);
477:   KSPMonitorSetFromOptions(ksp,"-ksp_monitor_max","Monitor the maximum norm of the residual","KSPMonitorTrueResidualMaxNorm",KSPMonitorTrueResidualMaxNorm);
478:   KSPMonitorSetFromOptions(ksp,"-ksp_monitor_short","Monitor preconditioned residual norm with fewer digits","KSPMonitorDefaultShort",KSPMonitorDefaultShort);
479:   KSPMonitorSetFromOptions(ksp,"-ksp_monitor_solution","Monitor the solution","KSPMonitorSolution",KSPMonitorSolution);
480:   KSPMonitorSetFromOptions(ksp,"-ksp_monitor_singular_value","Monitor singular values","KSPMonitorSingularValue",KSPMonitorSingularValue);
481:   PetscOptionsHasName(NULL,((PetscObject)ksp)->prefix,"-ksp_monitor_singular_value",&flg);
482:   if (flg) {
483:     KSPSetComputeSingularValues(ksp,PETSC_TRUE);
484:   }
485:   PetscObjectTypeCompare((PetscObject)ksp->pc,PCKSP,&flg);
486:   if (!flg) PetscObjectTypeCompare((PetscObject)ksp->pc,PCBJACOBI,&flg);
487:   if (!flg) PetscObjectTypeCompare((PetscObject)ksp->pc,PCDEFLATION,&flg);

489:   if (flg) {
490:     /* A hack for using dynamic tolerance in preconditioner */
491:     PetscOptionsString("-sub_ksp_dynamic_tolerance","Use dynamic tolerance for PC if PC is a KSP","KSPMonitorDynamicTolerance","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
492:     if (flg) {
493:       KSPDynTolCtx *scale;
494:       PetscMalloc1(1,&scale);
495:       scale->bnrm = -1.0;
496:       scale->coef = 1.0;
497:       PetscOptionsReal("-sub_ksp_dynamic_tolerance_param","Parameter of dynamic tolerance for inner PCKSP","KSPMonitorDynamicToleranceParam",scale->coef,&scale->coef,&flg);
498:       KSPMonitorSet(ksp,KSPMonitorDynamicTolerance,scale,KSPMonitorDynamicToleranceDestroy);
499:     }
500:   }

502:   /*
503:    Calls Python function
504:   */
505:   PetscOptionsString("-ksp_monitor_python","Use Python function","KSPMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
506:   if (flg) {PetscPythonMonitorSet((PetscObject)ksp,monfilename);}
507:   /*
508:     Graphically plots preconditioned residual norm
509:   */
510:   PetscOptionsBool("-ksp_monitor_lg_residualnorm","Monitor graphically preconditioned residual norm","KSPMonitorSet",PETSC_FALSE,&flg,&set);
511:   if (set && flg) {
512:     PetscDrawLG ctx;

514:     KSPMonitorLGResidualNormCreate(comm,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
515:     KSPMonitorSet(ksp,KSPMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
516:   }
517:   /*
518:     Graphically plots preconditioned and true residual norm
519:   */
520:   PetscOptionsBool("-ksp_monitor_lg_true_residualnorm","Monitor graphically true residual norm","KSPMonitorSet",PETSC_FALSE,&flg,&set);
521:   if (set && flg) {
522:     PetscDrawLG ctx;

524:     KSPMonitorLGTrueResidualNormCreate(comm,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
525:     KSPMonitorSet(ksp,KSPMonitorLGTrueResidualNorm,ctx,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
526:   }
527:   /*
528:     Graphically plots preconditioned residual norm and range of residual element values
529:   */
530:   PetscOptionsBool("-ksp_monitor_lg_range","Monitor graphically range of preconditioned residual norm","KSPMonitorSet",PETSC_FALSE,&flg,&set);
531:   if (set && flg) {
532:     PetscViewer ctx;

534:     PetscViewerDrawOpen(comm,NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
535:     KSPMonitorSet(ksp,KSPMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
536:   }
537:   /* TODO Do these show up in help? */
538:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view",&ksp->viewer,&ksp->format,&ksp->view);
539:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_pre",&ksp->viewerPre,&ksp->formatPre,&ksp->viewPre);
540:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_converged_reason",&ksp->viewerReason,&ksp->formatReason,&ksp->viewReason);
541:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat",&ksp->viewerMat,&ksp->formatMat,&ksp->viewMat);
542:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_pmat",&ksp->viewerPMat,&ksp->formatPMat,&ksp->viewPMat);
543:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_rhs",&ksp->viewerRhs,&ksp->formatRhs,&ksp->viewRhs);
544:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_solution",&ksp->viewerSol,&ksp->formatSol,&ksp->viewSol);
545:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_mat_explicit",&ksp->viewerMatExp,&ksp->formatMatExp,&ksp->viewMatExp);
546:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_eigenvalues",&ksp->viewerEV,&ksp->formatEV,&ksp->viewEV);
547:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_singularvalues",&ksp->viewerSV,&ksp->formatSV,&ksp->viewSV);
548:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_eigenvalues_explicit",&ksp->viewerEVExp,&ksp->formatEVExp,&ksp->viewEVExp);
549:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_final_residual",&ksp->viewerFinalRes,&ksp->formatFinalRes,&ksp->viewFinalRes);
550:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_preconditioned_operator_explicit",&ksp->viewerPOpExp,&ksp->formatPOpExp,&ksp->viewPOpExp);
551:   PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_view_diagonal_scale",&ksp->viewerDScale,&ksp->formatDScale,&ksp->viewDScale);

553:   /* Deprecated options */
554:   if (!ksp->viewEV) {
555:     PetscOptionsDeprecated("-ksp_compute_eigenvalues",NULL,"3.9","Use -ksp_view_eigenvalues");
556:     PetscOptionsGetViewer(comm, ((PetscObject) ksp)->options,prefix, "-ksp_compute_eigenvalues",&ksp->viewerEV,&ksp->formatEV,&ksp->viewEV);
557:   }
558:   if (!ksp->viewEV) {
559:     PetscOptionsDeprecated("-ksp_plot_eigenvalues",NULL,"3.9","Use -ksp_view_eigenvalues draw");
560:     PetscOptionsName("-ksp_plot_eigenvalues", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw]", "KSPView", &ksp->viewEV);
561:     if (ksp->viewEV) {
562:       ksp->formatEV = PETSC_VIEWER_DEFAULT;
563:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
564:       PetscObjectReference((PetscObject) ksp->viewerEV);
565:     }
566:   }
567:   if (!ksp->viewEV) {
568:     PetscOptionsDeprecated("-ksp_plot_eigencontours",NULL,"3.9","Use -ksp_view_eigenvalues draw::draw_contour");
569:     PetscOptionsName("-ksp_plot_eigencontours", "[deprecated since PETSc 3.9; use -ksp_view_eigenvalues draw::draw_contour]", "KSPView", &ksp->viewEV);
570:     if (ksp->viewEV) {
571:       ksp->formatEV = PETSC_VIEWER_DRAW_CONTOUR;
572:       ksp->viewerEV = PETSC_VIEWER_DRAW_(comm);
573:       PetscObjectReference((PetscObject) ksp->viewerEV);
574:     }
575:   }
576:   if (!ksp->viewEVExp) {
577:     PetscOptionsDeprecated("-ksp_compute_eigenvalues_explicitly",NULL,"3.9","Use -ksp_view_eigenvalues_explicit");
578:     PetscOptionsGetViewer(comm, ((PetscObject) ksp)->options,prefix, "-ksp_compute_eigenvalues_explicitly",&ksp->viewerEVExp,&ksp->formatEVExp,&ksp->viewEVExp);
579:   }
580:   if (!ksp->viewEVExp) {
581:     PetscOptionsDeprecated("-ksp_plot_eigenvalues_explicitly",NULL,"3.9","Use -ksp_view_eigenvalues_explicit draw");
582:     PetscOptionsName("-ksp_plot_eigenvalues_explicitly","[deprecated since PETSc 3.9; use -ksp_view_eigenvalues_explicit draw]","KSPView",&ksp->viewEVExp);
583:     if (ksp->viewEVExp) {
584:       ksp->formatEVExp = PETSC_VIEWER_DEFAULT;
585:       ksp->viewerEVExp = PETSC_VIEWER_DRAW_(comm);
586:       PetscObjectReference((PetscObject) ksp->viewerEVExp);
587:     }
588:   }
589:   if (!ksp->viewSV) {
590:     PetscOptionsDeprecated("-ksp_compute_singularvalues",NULL,"3.9","Use -ksp_view_singularvalues");
591:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_compute_singularvalues",&ksp->viewerSV,&ksp->formatSV,&ksp->viewSV);
592:   }
593:   if (!ksp->viewFinalRes) {
594:     PetscOptionsDeprecated("-ksp_final_residual",NULL,"3.9","Use -ksp_view_final_residual");
595:     PetscOptionsGetViewer(comm,((PetscObject) ksp)->options,prefix,"-ksp_final_residual",&ksp->viewerFinalRes,&ksp->formatFinalRes,&ksp->viewFinalRes);
596:   }

598: #if defined(PETSC_HAVE_SAWS)
599:   /*
600:     Publish convergence information using AMS
601:   */
602:   PetscOptionsBool("-ksp_monitor_saws","Publish KSP progress using SAWs","KSPMonitorSet",PETSC_FALSE,&flg,&set);
603:   if (set && flg) {
604:     void *ctx;
605:     KSPMonitorSAWsCreate(ksp,&ctx);
606:     KSPMonitorSet(ksp,KSPMonitorSAWs,ctx,KSPMonitorSAWsDestroy);
607:     KSPSetComputeSingularValues(ksp,PETSC_TRUE);
608:   }
609: #endif

611:   /* -----------------------------------------------------------------------*/
612:   KSPSetUpNorms_Private(ksp,PETSC_FALSE,NULL,&pcside);
613:   PetscOptionsEnum("-ksp_pc_side","KSP preconditioner side","KSPSetPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
614:   if (flg) {KSPSetPCSide(ksp,pcside);}

616:   if (ksp->viewSV || ksp->viewEV) {
617:     KSPSetComputeSingularValues(ksp,PETSC_TRUE);
618:   }

620: #if defined(PETSC_HAVE_SAWS)
621:   {
622:   PetscBool set;
623:   flg  = PETSC_FALSE;
624:   PetscOptionsBool("-ksp_saws_block","Block for SAWs at end of KSPSolve","PetscObjectSAWsBlock",((PetscObject)ksp)->amspublishblock,&flg,&set);
625:   if (set) {
626:     PetscObjectSAWsSetBlock((PetscObject)ksp,flg);
627:   }
628:   }
629: #endif

631:   if (ksp->ops->setfromoptions) {
632:     (*ksp->ops->setfromoptions)(PetscOptionsObject,ksp);
633:   }
634:   skipoptions:
635:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
636:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)ksp);
637:   PetscOptionsEnd();
638:   ksp->setfromoptionscalled++;
639:   return(0);
640: }

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

645:    Collective on ksp

647:    Input Parameter:
648: .  ksp - the KSP context

650:    Level: beginner

652: .seealso: KSPSetFromOptions(), KSPSetOptionsPrefix()
653: @*/
654: PetscErrorCode KSPResetFromOptions(KSP ksp)
655: {

659:   if (ksp->setfromoptionscalled) {KSPSetFromOptions(ksp);}
660:   return(0);
661: }