Actual source code: itcreate.c

petsc-3.15.0 2021-04-05
Report Typos and Errors
  1: /*
  2:      The basic KSP routines, Create, View etc. are here.
  3: */
  4: #include <petsc/private/kspimpl.h>

  6: /* Logging support */
  7: PetscClassId  KSP_CLASSID;
  8: PetscClassId  DMKSP_CLASSID;
  9: PetscClassId  KSPGUESS_CLASSID;
 10: PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve, KSP_SolveTranspose, KSP_MatSolve;

 12: /*
 13:    Contains the list of registered KSP routines
 14: */
 15: PetscFunctionList KSPList              = NULL;
 16: PetscBool         KSPRegisterAllCalled = PETSC_FALSE;

 18: /*
 19:    Contains the list of registered KSP monitors
 20: */
 21: PetscFunctionList KSPMonitorList              = NULL;
 22: PetscFunctionList KSPMonitorCreateList        = NULL;
 23: PetscFunctionList KSPMonitorDestroyList       = NULL;
 24: PetscBool         KSPMonitorRegisterAllCalled = PETSC_FALSE;

 26: /*@C
 27:   KSPLoad - Loads a KSP that has been stored in binary  with KSPView().

 29:   Collective on viewer

 31:   Input Parameters:
 32: + newdm - the newly loaded KSP, this needs to have been created with KSPCreate() or
 33:            some related function before a call to KSPLoad().
 34: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()

 36:    Level: intermediate

 38:   Notes:
 39:    The type is determined by the data in the file, any type set into the KSP before this call is ignored.

 41:   Notes for advanced users:
 42:   Most users should not need to know the details of the binary storage
 43:   format, since KSPLoad() and KSPView() completely hide these details.
 44:   But for anyone who's interested, the standard binary matrix storage
 45:   format is
 46: .vb
 47:      has not yet been determined
 48: .ve

 50: .seealso: PetscViewerBinaryOpen(), KSPView(), MatLoad(), VecLoad()
 51: @*/
 52: PetscErrorCode  KSPLoad(KSP newdm, PetscViewer viewer)
 53: {
 55:   PetscBool      isbinary;
 56:   PetscInt       classid;
 57:   char           type[256];
 58:   PC             pc;

 63:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
 64:   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");

 66:   PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
 67:   if (classid != KSP_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not KSP next in file");
 68:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
 69:   KSPSetType(newdm, type);
 70:   if (newdm->ops->load) {
 71:     (*newdm->ops->load)(newdm,viewer);
 72:   }
 73:   KSPGetPC(newdm,&pc);
 74:   PCLoad(pc,viewer);
 75:   return(0);
 76: }

 78: #include <petscdraw.h>
 79: #if defined(PETSC_HAVE_SAWS)
 80: #include <petscviewersaws.h>
 81: #endif
 82: /*@C
 83:    KSPView - Prints the KSP data structure.

 85:    Collective on ksp

 87:    Input Parameters:
 88: +  ksp - the Krylov space context
 89: -  viewer - visualization context

 91:    Options Database Keys:
 92: .  -ksp_view - print the KSP data structure at the end of a KSPSolve call

 94:    Note:
 95:    The available visualization contexts include
 96: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 97: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 98:          output where only the first processor opens
 99:          the file.  All other processors send their
100:          data to the first processor to print.

102:    The available formats include
103: +     PETSC_VIEWER_DEFAULT - standard output (default)
104: -     PETSC_VIEWER_ASCII_INFO_DETAIL - more verbose output for PCBJACOBI and PCASM

106:    The user can open an alternative visualization context with
107:    PetscViewerASCIIOpen() - output to a specified file.

109:   In the debugger you can do "call KSPView(ksp,0)" to display the KSP. (The same holds for any PETSc object viewer).

111:    Level: beginner

113: .seealso: PCView(), PetscViewerASCIIOpen()
114: @*/
115: PetscErrorCode  KSPView(KSP ksp,PetscViewer viewer)
116: {
118:   PetscBool      iascii,isbinary,isdraw,isstring;
119: #if defined(PETSC_HAVE_SAWS)
120:   PetscBool      issaws;
121: #endif

125:   if (!viewer) {
126:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
127:   }

131:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
132:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
133:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
134:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
135: #if defined(PETSC_HAVE_SAWS)
136:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
137: #endif
138:   if (iascii) {
139:     PetscObjectPrintClassNamePrefixType((PetscObject)ksp,viewer);
140:     if (ksp->ops->view) {
141:       PetscViewerASCIIPushTab(viewer);
142:       (*ksp->ops->view)(ksp,viewer);
143:       PetscViewerASCIIPopTab(viewer);
144:     }
145:     if (ksp->guess_zero) {
146:       PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, initial guess is zero\n",ksp->max_it);
147:     } else {
148:       PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, nonzero initial guess\n", ksp->max_it);
149:     }
150:     if (ksp->guess_knoll) {PetscViewerASCIIPrintf(viewer,"  using preconditioner applied to right hand side for initial guess\n");}
151:     PetscViewerASCIIPrintf(viewer,"  tolerances:  relative=%g, absolute=%g, divergence=%g\n",(double)ksp->rtol,(double)ksp->abstol,(double)ksp->divtol);
152:     if (ksp->pc_side == PC_RIGHT) {
153:       PetscViewerASCIIPrintf(viewer,"  right preconditioning\n");
154:     } else if (ksp->pc_side == PC_SYMMETRIC) {
155:       PetscViewerASCIIPrintf(viewer,"  symmetric preconditioning\n");
156:     } else {
157:       PetscViewerASCIIPrintf(viewer,"  left preconditioning\n");
158:     }
159:     if (ksp->guess) {
160:       PetscViewerASCIIPushTab(viewer);
161:       KSPGuessView(ksp->guess,viewer);
162:       PetscViewerASCIIPopTab(viewer);
163:     }
164:     if (ksp->dscale) {PetscViewerASCIIPrintf(viewer,"  diagonally scaled system\n");}
165:     PetscViewerASCIIPrintf(viewer,"  using %s norm type for convergence test\n",KSPNormTypes[ksp->normtype]);
166:   } else if (isbinary) {
167:     PetscInt    classid = KSP_FILE_CLASSID;
168:     MPI_Comm    comm;
169:     PetscMPIInt rank;
170:     char        type[256];

172:     PetscObjectGetComm((PetscObject)ksp,&comm);
173:     MPI_Comm_rank(comm,&rank);
174:     if (!rank) {
175:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);
176:       PetscStrncpy(type,((PetscObject)ksp)->type_name,256);
177:       PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR);
178:     }
179:     if (ksp->ops->view) {
180:       (*ksp->ops->view)(ksp,viewer);
181:     }
182:   } else if (isstring) {
183:     const char *type;
184:     KSPGetType(ksp,&type);
185:     PetscViewerStringSPrintf(viewer," KSPType: %-7.7s",type);
186:     if (ksp->ops->view) {(*ksp->ops->view)(ksp,viewer);}
187:   } else if (isdraw) {
188:     PetscDraw draw;
189:     char      str[36];
190:     PetscReal x,y,bottom,h;
191:     PetscBool flg;

193:     PetscViewerDrawGetDraw(viewer,0,&draw);
194:     PetscDrawGetCurrentPoint(draw,&x,&y);
195:     PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
196:     if (!flg) {
197:       PetscStrncpy(str,"KSP: ",sizeof(str));
198:       PetscStrlcat(str,((PetscObject)ksp)->type_name,sizeof(str));
199:       PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
200:       bottom = y - h;
201:     } else {
202:       bottom = y;
203:     }
204:     PetscDrawPushCurrentPoint(draw,x,bottom);
205: #if defined(PETSC_HAVE_SAWS)
206:   } else if (issaws) {
207:     PetscMPIInt rank;
208:     const char  *name;

210:     PetscObjectGetName((PetscObject)ksp,&name);
211:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
212:     if (!((PetscObject)ksp)->amsmem && !rank) {
213:       char       dir[1024];

215:       PetscObjectViewSAWs((PetscObject)ksp,viewer);
216:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
217:       PetscStackCallSAWs(SAWs_Register,(dir,&ksp->its,1,SAWs_READ,SAWs_INT));
218:       if (!ksp->res_hist) {
219:         KSPSetResidualHistory(ksp,NULL,PETSC_DECIDE,PETSC_TRUE);
220:       }
221:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/res_hist",name);
222:       PetscStackCallSAWs(SAWs_Register,(dir,ksp->res_hist,10,SAWs_READ,SAWs_DOUBLE));
223:     }
224: #endif
225:   } else if (ksp->ops->view) {
226:     (*ksp->ops->view)(ksp,viewer);
227:   }
228:   if (ksp->pc) {
229:     PCView(ksp->pc,viewer);
230:   }
231:   if (isdraw) {
232:     PetscDraw draw;
233:     PetscViewerDrawGetDraw(viewer,0,&draw);
234:     PetscDrawPopCurrentPoint(draw);
235:   }
236:   return(0);
237: }

239: /*@C
240:    KSPViewFromOptions - View from Options

242:    Collective on KSP

244:    Input Parameters:
245: +  A - Krylov solver context
246: .  obj - Optional object
247: -  name - command line option

249:    Level: intermediate
250: .seealso:  KSP, KSPView, PetscObjectViewFromOptions(), KSPCreate()
251: @*/
252: PetscErrorCode  KSPViewFromOptions(KSP A,PetscObject obj,const char name[])
253: {

258:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
259:   return(0);
260: }

262: /*@
263:    KSPSetNormType - Sets the norm that is used for convergence testing.

265:    Logically Collective on ksp

267:    Input Parameter:
268: +  ksp - Krylov solver context
269: -  normtype - one of
270: $   KSP_NORM_NONE - skips computing the norm, this should generally only be used if you are using
271: $                 the Krylov method as a smoother with a fixed small number of iterations.
272: $                 Implicitly sets KSPConvergedSkip() as KSP convergence test.
273: $                 Note that certain algorithms such as KSPGMRES ALWAYS require the norm calculation,
274: $                 for these methods the norms are still computed, they are just not used in
275: $                 the convergence test.
276: $   KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
277: $                 of the preconditioned residual P^{-1}(b - A x)
278: $   KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
279: $   KSP_NORM_NATURAL - supported  by KSPCG, KSPCR, KSPCGNE, KSPCGS


282:    Options Database Key:
283: .   -ksp_norm_type <none,preconditioned,unpreconditioned,natural>

285:    Notes:
286:    Not all combinations of preconditioner side (see KSPSetPCSide()) and norm type are supported by all Krylov methods.
287:    If only one is set, PETSc tries to automatically change the other to find a compatible pair.  If no such combination
288:    is supported, PETSc will generate an error.

290:    Developer Notes:
291:    Supported combinations of norm and preconditioner side are set using KSPSetSupportedNorm().

293:    Level: advanced

295: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
296: @*/
297: PetscErrorCode  KSPSetNormType(KSP ksp,KSPNormType normtype)
298: {
302:   ksp->normtype = ksp->normtype_set = normtype;
303:   return(0);
304: }

306: /*@
307:    KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be
308:      computed and used in the convergence test.

310:    Logically Collective on ksp

312:    Input Parameter:
313: +  ksp - Krylov solver context
314: -  it  - use -1 to check at all iterations

316:    Notes:
317:    Currently only works with KSPCG, KSPBCGS and KSPIBCGS

319:    Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm

321:    On steps where the norm is not computed, the previous norm is still in the variable, so if you run with, for example,
322:     -ksp_monitor the residual norm will appear to be unchanged for several iterations (though it is not really unchanged).
323:    Level: advanced

325: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType()
326: @*/
327: PetscErrorCode  KSPSetCheckNormIteration(KSP ksp,PetscInt it)
328: {
332:   ksp->chknorm = it;
333:   return(0);
334: }

336: /*@
337:    KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the MPI_Allreduce() for
338:    computing the inner products for the next iteration.  This can reduce communication costs at the expense of doing
339:    one additional iteration.


342:    Logically Collective on ksp

344:    Input Parameter:
345: +  ksp - Krylov solver context
346: -  flg - PETSC_TRUE or PETSC_FALSE

348:    Options Database Keys:
349: .  -ksp_lag_norm - lag the calculated residual norm

351:    Notes:
352:    Currently only works with KSPIBCGS.

354:    Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm

356:    If you lag the norm and run with, for example, -ksp_monitor, the residual norm reported will be the lagged one.
357:    Level: advanced

359: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType(), KSPSetCheckNormIteration()
360: @*/
361: PetscErrorCode  KSPSetLagNorm(KSP ksp,PetscBool flg)
362: {
366:   ksp->lagnorm = flg;
367:   return(0);
368: }

370: /*@
371:    KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP

373:    Logically Collective

375:    Input Arguments:
376: +  ksp - Krylov method
377: .  normtype - supported norm type
378: .  pcside - preconditioner side that can be used with this norm
379: -  priority - positive integer preference for this combination; larger values have higher priority

381:    Level: developer

383:    Notes:
384:    This function should be called from the implementation files KSPCreate_XXX() to declare
385:    which norms and preconditioner sides are supported. Users should not need to call this
386:    function.

388: .seealso: KSPSetNormType(), KSPSetPCSide()
389: @*/
390: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
391: {

395:   ksp->normsupporttable[normtype][pcside] = priority;
396:   return(0);
397: }

399: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
400: {

404:   PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
405:   ksp->pc_side  = ksp->pc_side_set;
406:   ksp->normtype = ksp->normtype_set;
407:   return(0);
408: }

410: PetscErrorCode KSPSetUpNorms_Private(KSP ksp,PetscBool errorifnotsupported,KSPNormType *normtype,PCSide *pcside)
411: {
412:   PetscInt i,j,best,ibest = 0,jbest = 0;

415:   best = 0;
416:   for (i=0; i<KSP_NORM_MAX; i++) {
417:     for (j=0; j<PC_SIDE_MAX; j++) {
418:       if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
419:         best  = ksp->normsupporttable[i][j];
420:         ibest = i;
421:         jbest = j;
422:       }
423:     }
424:   }
425:   if (best < 1 && errorifnotsupported) {
426:     if (ksp->normtype == KSP_NORM_DEFAULT && ksp->pc_side == PC_SIDE_DEFAULT) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"The %s KSP implementation did not call KSPSetSupportedNorm()",((PetscObject)ksp)->type_name);
427:     if (ksp->normtype == KSP_NORM_DEFAULT) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSP %s does not support %s",((PetscObject)ksp)->type_name,PCSides[ksp->pc_side]);
428:     if (ksp->pc_side == PC_SIDE_DEFAULT) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSP %s does not support %s",((PetscObject)ksp)->type_name,KSPNormTypes[ksp->normtype]);
429:     SETERRQ3(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSP %s does not support %s with %s",((PetscObject)ksp)->type_name,KSPNormTypes[ksp->normtype],PCSides[ksp->pc_side]);
430:   }
431:   if (normtype) *normtype = (KSPNormType)ibest;
432:   if (pcside)   *pcside   = (PCSide)jbest;
433:   return(0);
434: }

436: /*@
437:    KSPGetNormType - Gets the norm that is used for convergence testing.

439:    Not Collective

441:    Input Parameter:
442: .  ksp - Krylov solver context

444:    Output Parameter:
445: .  normtype - norm that is used for convergence testing

447:    Level: advanced

449: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
450: @*/
451: PetscErrorCode  KSPGetNormType(KSP ksp, KSPNormType *normtype)
452: {

458:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
459:   *normtype = ksp->normtype;
460:   return(0);
461: }

463: #if defined(PETSC_HAVE_SAWS)
464: #include <petscviewersaws.h>
465: #endif

467: /*@
468:    KSPSetOperators - Sets the matrix associated with the linear system
469:    and a (possibly) different one associated with the preconditioner.

471:    Collective on ksp

473:    Input Parameters:
474: +  ksp - the KSP context
475: .  Amat - the matrix that defines the linear system
476: -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

478:    Notes:

480:     If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
481:     space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.

483:     All future calls to KSPSetOperators() must use the same size matrices!

485:     Passing a NULL for Amat or Pmat removes the matrix that is currently used.

487:     If you wish to replace either Amat or Pmat but leave the other one untouched then
488:     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
489:     on it and then pass it back in in your call to KSPSetOperators().

491:     Level: beginner

493:    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
494:       are created in PC and returned to the user. In this case, if both operators
495:       mat and pmat are requested, two DIFFERENT operators will be returned. If
496:       only one is requested both operators in the PC will be the same (i.e. as
497:       if one had called KSP/PCSetOperators() with the same argument for both Mats).
498:       The user must set the sizes of the returned matrices and their type etc just
499:       as if the user created them with MatCreate(). For example,

501: $         KSP/PCGetOperators(ksp/pc,&mat,NULL); is equivalent to
502: $           set size, type, etc of mat

504: $         MatCreate(comm,&mat);
505: $         KSP/PCSetOperators(ksp/pc,mat,mat);
506: $         PetscObjectDereference((PetscObject)mat);
507: $           set size, type, etc of mat

509:      and

511: $         KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to
512: $           set size, type, etc of mat and pmat

514: $         MatCreate(comm,&mat);
515: $         MatCreate(comm,&pmat);
516: $         KSP/PCSetOperators(ksp/pc,mat,pmat);
517: $         PetscObjectDereference((PetscObject)mat);
518: $         PetscObjectDereference((PetscObject)pmat);
519: $           set size, type, etc of mat and pmat

521:     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
522:     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
523:     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
524:     at this is when you create a SNES you do not NEED to create a KSP and attach it to
525:     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
526:     you do not need to attach a PC to it (the KSP object manages the PC object for you).
527:     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
528:     it can be created for you?

530: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
531: @*/
532: PetscErrorCode  KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
533: {

542:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
543:   PCSetOperators(ksp->pc,Amat,Pmat);
544:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;  /* so that next solve call will call PCSetUp() on new matrix */
545:   return(0);
546: }

548: /*@
549:    KSPGetOperators - Gets the matrix associated with the linear system
550:    and a (possibly) different one associated with the preconditioner.

552:    Collective on ksp

554:    Input Parameter:
555: .  ksp - the KSP context

557:    Output Parameters:
558: +  Amat - the matrix that defines the linear system
559: -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

561:     Level: intermediate

563:    Notes:
564:     DOES NOT increase the reference counts of the matrix, so you should NOT destroy them.

566: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
567: @*/
568: PetscErrorCode  KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
569: {

574:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
575:   PCGetOperators(ksp->pc,Amat,Pmat);
576:   return(0);
577: }

579: /*@C
580:    KSPGetOperatorsSet - Determines if the matrix associated with the linear system and
581:    possibly a different one associated with the preconditioner have been set in the KSP.

583:    Not collective, though the results on all processes should be the same

585:    Input Parameter:
586: .  pc - the KSP context

588:    Output Parameters:
589: +  mat - the matrix associated with the linear system was set
590: -  pmat - matrix associated with the preconditioner was set, usually the same

592:    Level: intermediate

594: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
595: @*/
596: PetscErrorCode  KSPGetOperatorsSet(KSP ksp,PetscBool  *mat,PetscBool  *pmat)
597: {

602:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
603:   PCGetOperatorsSet(ksp->pc,mat,pmat);
604:   return(0);
605: }

607: /*@C
608:    KSPSetPreSolve - Sets a function that is called before every KSPSolve() is started

610:    Logically Collective on ksp

612:    Input Parameters:
613: +   ksp - the solver object
614: .   presolve - the function to call before the solve
615: -   prectx - any context needed by the function

617:    Calling sequence of presolve:
618: $  func(KSP ksp,Vec rhs,Vec x,void *ctx)

620: +  ksp - the KSP context
621: .  rhs - the right-hand side vector
622: .  x - the solution vector
623: -  ctx - optional user-provided context

625:    Level: developer

627: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
628: @*/
629: PetscErrorCode  KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
630: {
633:   ksp->presolve = presolve;
634:   ksp->prectx   = prectx;
635:   return(0);
636: }

638: /*@C
639:    KSPSetPostSolve - Sets a function that is called after every KSPSolve() completes (whether it converges or not)

641:    Logically Collective on ksp

643:    Input Parameters:
644: +   ksp - the solver object
645: .   postsolve - the function to call after the solve
646: -   postctx - any context needed by the function

648:    Level: developer

650:    Calling sequence of postsolve:
651: $  func(KSP ksp,Vec rhs,Vec x,void *ctx)

653: +  ksp - the KSP context
654: .  rhs - the right-hand side vector
655: .  x - the solution vector
656: -  ctx - optional user-provided context

658: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPreSolve()
659: @*/
660: PetscErrorCode  KSPSetPostSolve(KSP ksp,PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*),void *postctx)
661: {
664:   ksp->postsolve = postsolve;
665:   ksp->postctx   = postctx;
666:   return(0);
667: }

669: /*@
670:    KSPCreate - Creates the default KSP context.

672:    Collective

674:    Input Parameter:
675: .  comm - MPI communicator

677:    Output Parameter:
678: .  ksp - location to put the KSP context

680:    Notes:
681:    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
682:    orthogonalization.

684:    Level: beginner

686: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
687: @*/
688: PetscErrorCode  KSPCreate(MPI_Comm comm,KSP *inksp)
689: {
690:   KSP            ksp;
692:   void           *ctx;

696:   *inksp = NULL;
697:   KSPInitializePackage();

699:   PetscHeaderCreate(ksp,KSP_CLASSID,"KSP","Krylov Method","KSP",comm,KSPDestroy,KSPView);

701:   ksp->max_it  = 10000;
702:   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
703:   ksp->rtol    = 1.e-5;
704: #if defined(PETSC_USE_REAL_SINGLE)
705:   ksp->abstol  = 1.e-25;
706: #else
707:   ksp->abstol  = 1.e-50;
708: #endif
709:   ksp->divtol  = 1.e4;

711:   ksp->chknorm        = -1;
712:   ksp->normtype       = ksp->normtype_set = KSP_NORM_DEFAULT;
713:   ksp->rnorm          = 0.0;
714:   ksp->its            = 0;
715:   ksp->guess_zero     = PETSC_TRUE;
716:   ksp->calc_sings     = PETSC_FALSE;
717:   ksp->res_hist       = NULL;
718:   ksp->res_hist_alloc = NULL;
719:   ksp->res_hist_len   = 0;
720:   ksp->res_hist_max   = 0;
721:   ksp->res_hist_reset = PETSC_TRUE;
722:   ksp->err_hist       = NULL;
723:   ksp->err_hist_alloc = NULL;
724:   ksp->err_hist_len   = 0;
725:   ksp->err_hist_max   = 0;
726:   ksp->err_hist_reset = PETSC_TRUE;
727:   ksp->numbermonitors = 0;
728:   ksp->numberreasonviews = 0;
729:   ksp->setfromoptionscalled = 0;
730:   ksp->nmax = PETSC_DECIDE;

732:   KSPConvergedDefaultCreate(&ctx);
733:   KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
734:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
735:   ksp->ops->buildresidual = KSPBuildResidualDefault;

737:   ksp->vec_sol    = NULL;
738:   ksp->vec_rhs    = NULL;
739:   ksp->pc         = NULL;
740:   ksp->data       = NULL;
741:   ksp->nwork      = 0;
742:   ksp->work       = NULL;
743:   ksp->reason     = KSP_CONVERGED_ITERATING;
744:   ksp->setupstage = KSP_SETUP_NEW;

746:   KSPNormSupportTableReset_Private(ksp);

748:   *inksp = ksp;
749:   return(0);
750: }

752: /*@C
753:    KSPSetType - Builds KSP for a particular solver.

755:    Logically Collective on ksp

757:    Input Parameters:
758: +  ksp      - the Krylov space context
759: -  type - a known method

761:    Options Database Key:
762: .  -ksp_type  <method> - Sets the method; use -help for a list
763:     of available methods (for instance, cg or gmres)

765:    Notes:
766:    See "petsc/include/petscksp.h" for available methods (for instance,
767:    KSPCG or KSPGMRES).

769:   Normally, it is best to use the KSPSetFromOptions() command and
770:   then set the KSP type from the options database rather than by using
771:   this routine.  Using the options database provides the user with
772:   maximum flexibility in evaluating the many different Krylov methods.
773:   The KSPSetType() routine is provided for those situations where it
774:   is necessary to set the iterative solver independently of the command
775:   line or options database.  This might be the case, for example, when
776:   the choice of iterative solver changes during the execution of the
777:   program, and the user's application is taking responsibility for
778:   choosing the appropriate method.  In other words, this routine is
779:   not for beginners.

781:   Level: intermediate

783:   Developer Note: KSPRegister() is used to add Krylov types to KSPList from which they
784:   are accessed by KSPSetType().

786: .seealso: PCSetType(), KSPType, KSPRegister(), KSPCreate()

788: @*/
789: PetscErrorCode  KSPSetType(KSP ksp, KSPType type)
790: {
791:   PetscErrorCode ierr,(*r)(KSP);
792:   PetscBool      match;


798:   PetscObjectTypeCompare((PetscObject)ksp,type,&match);
799:   if (match) return(0);

801:   PetscFunctionListFind(KSPList,type,&r);
802:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested KSP type %s",type);
803:   /* Destroy the previous private KSP context */
804:   if (ksp->ops->destroy) {
805:     (*ksp->ops->destroy)(ksp);
806:     ksp->ops->destroy = NULL;
807:   }
808:   /* Reinitialize function pointers in KSPOps structure */
809:   PetscMemzero(ksp->ops,sizeof(struct _KSPOps));
810:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
811:   ksp->ops->buildresidual = KSPBuildResidualDefault;
812:   KSPNormSupportTableReset_Private(ksp);
813:   ksp->setupnewmatrix     = PETSC_FALSE; // restore default (setup not called in case of new matrix)
814:   /* Call the KSPCreate_XXX routine for this particular Krylov solver */
815:   ksp->setupstage = KSP_SETUP_NEW;
816:   (*r)(ksp);
817:   PetscObjectChangeTypeName((PetscObject)ksp,type);
818:   return(0);
819: }

821: /*@C
822:    KSPGetType - Gets the KSP type as a string from the KSP object.

824:    Not Collective

826:    Input Parameter:
827: .  ksp - Krylov context

829:    Output Parameter:
830: .  name - name of KSP method

832:    Level: intermediate

834: .seealso: KSPSetType()
835: @*/
836: PetscErrorCode  KSPGetType(KSP ksp,KSPType *type)
837: {
841:   *type = ((PetscObject)ksp)->type_name;
842:   return(0);
843: }

845: /*@C
846:   KSPRegister -  Adds a method to the Krylov subspace solver package.

848:    Not Collective

850:    Input Parameters:
851: +  name_solver - name of a new user-defined solver
852: -  routine_create - routine to create method context

854:    Notes:
855:    KSPRegister() may be called multiple times to add several user-defined solvers.

857:    Sample usage:
858: .vb
859:    KSPRegister("my_solver",MySolverCreate);
860: .ve

862:    Then, your solver can be chosen with the procedural interface via
863: $     KSPSetType(ksp,"my_solver")
864:    or at runtime via the option
865: $     -ksp_type my_solver

867:    Level: advanced

869: .seealso: KSPRegisterAll()
870: @*/
871: PetscErrorCode  KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
872: {

876:   KSPInitializePackage();
877:   PetscFunctionListAdd(&KSPList,sname,function);
878:   return(0);
879: }

881: PetscErrorCode KSPMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[])
882: {

886:   PetscStrncpy(key, name, PETSC_MAX_PATH_LEN);
887:   PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN);
888:   PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN);
889:   PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN);
890:   PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN);
891:   return(0);
892: }

894: /*@C
895:   KSPMonitorRegister -  Adds Krylov subspace solver monitor routine.

897:   Not Collective

899:   Input Parameters:
900: + name    - name of a new monitor routine
901: . vtype   - A PetscViewerType for the output
902: . format  - A PetscViewerFormat for the output
903: . monitor - Monitor routine
904: . create  - Creation routine, or NULL
905: - destroy - Destruction routine, or NULL

907:   Notes:
908:   KSPMonitorRegister() may be called multiple times to add several user-defined monitors.

910:   Sample usage:
911: .vb
912:   KSPMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
913: .ve

915:   Then, your monitor can be chosen with the procedural interface via
916: $     KSPMonitorSetFromOptions(ksp,"-ksp_monitor_my_monitor","my_monitor",NULL)
917:   or at runtime via the option
918: $     -ksp_monitor_my_monitor

920:    Level: advanced

922: .seealso: KSPMonitorRegisterAll()
923: @*/
924: PetscErrorCode KSPMonitorRegister(const char name[], PetscViewerType vtype, PetscViewerFormat format,
925:                                   PetscErrorCode (*monitor)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *),
926:                                   PetscErrorCode (*create)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **),
927:                                   PetscErrorCode (*destroy)(PetscViewerAndFormat **))
928: {
929:   char           key[PETSC_MAX_PATH_LEN];

933:   KSPInitializePackage();
934:   KSPMonitorMakeKey_Internal(name, vtype, format, key);
935:   PetscFunctionListAdd(&KSPMonitorList, key, monitor);
936:   if (create)  {PetscFunctionListAdd(&KSPMonitorCreateList,  key, create);}
937:   if (destroy) {PetscFunctionListAdd(&KSPMonitorDestroyList, key, destroy);}
938:   return(0);
939: }