Actual source code: itcreate.c

petsc-master 2020-06-02
Report Typos and Errors

  2: /*
  3:      The basic KSP routines, Create, View etc. are here.
  4: */
  5:  #include <petsc/private/kspimpl.h>

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

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

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

 22:   Collective on viewer

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

 29:    Level: intermediate

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

 34:   Notes for advanced users:
 35:   Most users should not need to know the details of the binary storage
 36:   format, since KSPLoad() and KSPView() completely hide these details.
 37:   But for anyone who's interested, the standard binary matrix storage
 38:   format is
 39: .vb
 40:      has not yet been determined
 41: .ve

 43: .seealso: PetscViewerBinaryOpen(), KSPView(), MatLoad(), VecLoad()
 44: @*/
 45: PetscErrorCode  KSPLoad(KSP newdm, PetscViewer viewer)
 46: {
 48:   PetscBool      isbinary;
 49:   PetscInt       classid;
 50:   char           type[256];
 51:   PC             pc;

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

 59:   PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
 60:   if (classid != KSP_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not KSP next in file");
 61:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
 62:   KSPSetType(newdm, type);
 63:   if (newdm->ops->load) {
 64:     (*newdm->ops->load)(newdm,viewer);
 65:   }
 66:   KSPGetPC(newdm,&pc);
 67:   PCLoad(pc,viewer);
 68:   return(0);
 69: }

 71:  #include <petscdraw.h>
 72: #if defined(PETSC_HAVE_SAWS)
 73:  #include <petscviewersaws.h>
 74: #endif
 75: /*@C
 76:    KSPView - Prints the KSP data structure.

 78:    Collective on ksp

 80:    Input Parameters:
 81: +  ksp - the Krylov space context
 82: -  viewer - visualization context

 84:    Options Database Keys:
 85: .  -ksp_view - print the ksp data structure at the end of a KSPSolve call

 87:    Note:
 88:    The available visualization contexts include
 89: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
 90: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
 91:          output where only the first processor opens
 92:          the file.  All other processors send their
 93:          data to the first processor to print.

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

 98:    Level: beginner

100: .seealso: PCView(), PetscViewerASCIIOpen()
101: @*/
102: PetscErrorCode  KSPView(KSP ksp,PetscViewer viewer)
103: {
105:   PetscBool      iascii,isbinary,isdraw,isstring;
106: #if defined(PETSC_HAVE_SAWS)
107:   PetscBool      issaws;
108: #endif

112:   if (!viewer) {
113:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
114:   }

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

159:     PetscObjectGetComm((PetscObject)ksp,&comm);
160:     MPI_Comm_rank(comm,&rank);
161:     if (!rank) {
162:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);
163:       PetscStrncpy(type,((PetscObject)ksp)->type_name,256);
164:       PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR);
165:     }
166:     if (ksp->ops->view) {
167:       (*ksp->ops->view)(ksp,viewer);
168:     }
169:   } else if (isstring) {
170:     const char *type;
171:     KSPGetType(ksp,&type);
172:     PetscViewerStringSPrintf(viewer," KSPType: %-7.7s",type);
173:     if (ksp->ops->view) {(*ksp->ops->view)(ksp,viewer);}
174:   } else if (isdraw) {
175:     PetscDraw draw;
176:     char      str[36];
177:     PetscReal x,y,bottom,h;
178:     PetscBool flg;

180:     PetscViewerDrawGetDraw(viewer,0,&draw);
181:     PetscDrawGetCurrentPoint(draw,&x,&y);
182:     PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
183:     if (!flg) {
184:       PetscStrncpy(str,"KSP: ",sizeof(str));
185:       PetscStrlcat(str,((PetscObject)ksp)->type_name,sizeof(str));
186:       PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
187:       bottom = y - h;
188:     } else {
189:       bottom = y;
190:     }
191:     PetscDrawPushCurrentPoint(draw,x,bottom);
192: #if defined(PETSC_HAVE_SAWS)
193:   } else if (issaws) {
194:     PetscMPIInt rank;
195:     const char  *name;

197:     PetscObjectGetName((PetscObject)ksp,&name);
198:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
199:     if (!((PetscObject)ksp)->amsmem && !rank) {
200:       char       dir[1024];

202:       PetscObjectViewSAWs((PetscObject)ksp,viewer);
203:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
204:       PetscStackCallSAWs(SAWs_Register,(dir,&ksp->its,1,SAWs_READ,SAWs_INT));
205:       if (!ksp->res_hist) {
206:         KSPSetResidualHistory(ksp,NULL,PETSC_DECIDE,PETSC_TRUE);
207:       }
208:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/res_hist",name);
209:       PetscStackCallSAWs(SAWs_Register,(dir,ksp->res_hist,10,SAWs_READ,SAWs_DOUBLE));
210:     }
211: #endif
212:   } else if (ksp->ops->view) {
213:     (*ksp->ops->view)(ksp,viewer);
214:   }
215:   if (ksp->pc) {
216:     PCView(ksp->pc,viewer);
217:   }
218:   if (isdraw) {
219:     PetscDraw draw;
220:     PetscViewerDrawGetDraw(viewer,0,&draw);
221:     PetscDrawPopCurrentPoint(draw);
222:   }
223:   return(0);
224: }

226: /*@C
227:    KSPViewFromOptions - View from Options

229:    Collective on KSP

231:    Input Parameters:
232: +  A - Krylov solver context
233: .  obj - Optional object
234: -  name - command line option

236:    Level: intermediate
237: .seealso:  KSP, KSPView, PetscObjectViewFromOptions(), KSPCreate()
238: @*/
239: PetscErrorCode  KSPViewFromOptions(KSP A,PetscObject obj,const char name[])
240: {

245:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
246:   return(0);
247: }

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

252:    Logically Collective on ksp

254:    Input Parameter:
255: +  ksp - Krylov solver context
256: -  normtype - one of
257: $   KSP_NORM_NONE - skips computing the norm, this should generally only be used if you are using
258: $                 the Krylov method as a smoother with a fixed small number of iterations.
259: $                 Implicitly sets KSPConvergedSkip() as KSP convergence test.
260: $                 Note that certain algorithms such as KSPGMRES ALWAYS require the norm calculation,
261: $                 for these methods the norms are still computed, they are just not used in
262: $                 the convergence test. 
263: $   KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
264: $                 of the preconditioned residual P^{-1}(b - A x)
265: $   KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
266: $   KSP_NORM_NATURAL - supported  by KSPCG, KSPCR, KSPCGNE, KSPCGS


269:    Options Database Key:
270: .   -ksp_norm_type <none,preconditioned,unpreconditioned,natural>

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

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

280:    Level: advanced

282: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
283: @*/
284: PetscErrorCode  KSPSetNormType(KSP ksp,KSPNormType normtype)
285: {
289:   ksp->normtype = ksp->normtype_set = normtype;
290:   return(0);
291: }

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

297:    Logically Collective on ksp

299:    Input Parameter:
300: +  ksp - Krylov solver context
301: -  it  - use -1 to check at all iterations

303:    Notes:
304:    Currently only works with KSPCG, KSPBCGS and KSPIBCGS

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

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

312: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType()
313: @*/
314: PetscErrorCode  KSPSetCheckNormIteration(KSP ksp,PetscInt it)
315: {
319:   ksp->chknorm = it;
320:   return(0);
321: }

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


329:    Logically Collective on ksp

331:    Input Parameter:
332: +  ksp - Krylov solver context
333: -  flg - PETSC_TRUE or PETSC_FALSE

335:    Options Database Keys:
336: .  -ksp_lag_norm - lag the calculated residual norm

338:    Notes:
339:    Currently only works with KSPIBCGS.

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

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

346: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType(), KSPSetCheckNormIteration()
347: @*/
348: PetscErrorCode  KSPSetLagNorm(KSP ksp,PetscBool flg)
349: {
353:   ksp->lagnorm = flg;
354:   return(0);
355: }

357: /*@
358:    KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP

360:    Logically Collective

362:    Input Arguments:
363: +  ksp - Krylov method
364: .  normtype - supported norm type
365: .  pcside - preconditioner side that can be used with this norm
366: -  priority - positive integer preference for this combination; larger values have higher priority

368:    Level: developer

370:    Notes:
371:    This function should be called from the implementation files KSPCreate_XXX() to declare
372:    which norms and preconditioner sides are supported. Users should not need to call this
373:    function.

375: .seealso: KSPSetNormType(), KSPSetPCSide()
376: @*/
377: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
378: {

382:   ksp->normsupporttable[normtype][pcside] = priority;
383:   return(0);
384: }

386: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
387: {

391:   PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
392:   ksp->pc_side  = ksp->pc_side_set;
393:   ksp->normtype = ksp->normtype_set;
394:   return(0);
395: }

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

402:   best = 0;
403:   for (i=0; i<KSP_NORM_MAX; i++) {
404:     for (j=0; j<PC_SIDE_MAX; j++) {
405:       if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
406:         best  = ksp->normsupporttable[i][j];
407:         ibest = i;
408:         jbest = j;
409:       }
410:     }
411:   }
412:   if (best < 1 && errorifnotsupported) {
413:     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);
414:     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]);
415:     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]);
416:     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]);
417:   }
418:   if (normtype) *normtype = (KSPNormType)ibest;
419:   if (pcside)   *pcside   = (PCSide)jbest;
420:   return(0);
421: }

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

426:    Not Collective

428:    Input Parameter:
429: .  ksp - Krylov solver context

431:    Output Parameter:
432: .  normtype - norm that is used for convergence testing

434:    Level: advanced

436: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
437: @*/
438: PetscErrorCode  KSPGetNormType(KSP ksp, KSPNormType *normtype)
439: {

445:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
446:   *normtype = ksp->normtype;
447:   return(0);
448: }

450: #if defined(PETSC_HAVE_SAWS)
451:  #include <petscviewersaws.h>
452: #endif

454: /*@
455:    KSPSetOperators - Sets the matrix associated with the linear system
456:    and a (possibly) different one associated with the preconditioner.

458:    Collective on ksp

460:    Input Parameters:
461: +  ksp - the KSP context
462: .  Amat - the matrix that defines the linear system
463: -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

465:    Notes:

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

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

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

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

478:     Level: beginner

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

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

491: $         MatCreate(comm,&mat);
492: $         KSP/PCSetOperators(ksp/pc,mat,mat);
493: $         PetscObjectDereference((PetscObject)mat);
494: $           set size, type, etc of mat

496:      and

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

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

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

517: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
518: @*/
519: PetscErrorCode  KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
520: {

529:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
530:   PCSetOperators(ksp->pc,Amat,Pmat);
531:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;  /* so that next solve call will call PCSetUp() on new matrix */
532:   return(0);
533: }

535: /*@
536:    KSPGetOperators - Gets the matrix associated with the linear system
537:    and a (possibly) different one associated with the preconditioner.

539:    Collective on ksp

541:    Input Parameter:
542: .  ksp - the KSP context

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

548:     Level: intermediate

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

553: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
554: @*/
555: PetscErrorCode  KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
556: {

561:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
562:   PCGetOperators(ksp->pc,Amat,Pmat);
563:   return(0);
564: }

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

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

572:    Input Parameter:
573: .  pc - the KSP context

575:    Output Parameters:
576: +  mat - the matrix associated with the linear system was set
577: -  pmat - matrix associated with the preconditioner was set, usually the same

579:    Level: intermediate

581: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
582: @*/
583: PetscErrorCode  KSPGetOperatorsSet(KSP ksp,PetscBool  *mat,PetscBool  *pmat)
584: {

589:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
590:   PCGetOperatorsSet(ksp->pc,mat,pmat);
591:   return(0);
592: }

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

597:    Logically Collective on ksp

599:    Input Parameters:
600: +   ksp - the solver object
601: .   presolve - the function to call before the solve
602: -   prectx - any context needed by the function

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

607: +  ksp - the KSP context
608: .  rhs - the right-hand side vector
609: .  x - the solution vector
610: -  ctx - optional user-provided context

612:    Level: developer

614: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
615: @*/
616: PetscErrorCode  KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
617: {
620:   ksp->presolve = presolve;
621:   ksp->prectx   = prectx;
622:   return(0);
623: }

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

628:    Logically Collective on ksp

630:    Input Parameters:
631: +   ksp - the solver object
632: .   postsolve - the function to call after the solve
633: -   postctx - any context needed by the function

635:    Level: developer

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

640: +  ksp - the KSP context
641: .  rhs - the right-hand side vector
642: .  x - the solution vector
643: -  ctx - optional user-provided context

645: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPreSolve()
646: @*/
647: PetscErrorCode  KSPSetPostSolve(KSP ksp,PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*),void *postctx)
648: {
651:   ksp->postsolve = postsolve;
652:   ksp->postctx   = postctx;
653:   return(0);
654: }

656: /*@
657:    KSPCreate - Creates the default KSP context.

659:    Collective

661:    Input Parameter:
662: .  comm - MPI communicator

664:    Output Parameter:
665: .  ksp - location to put the KSP context

667:    Notes:
668:    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
669:    orthogonalization.

671:    Level: beginner

673: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
674: @*/
675: PetscErrorCode  KSPCreate(MPI_Comm comm,KSP *inksp)
676: {
677:   KSP            ksp;
679:   void           *ctx;

683:   *inksp = NULL;
684:   KSPInitializePackage();

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

688:   ksp->max_it  = 10000;
689:   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
690:   ksp->rtol    = 1.e-5;
691: #if defined(PETSC_USE_REAL_SINGLE)
692:   ksp->abstol  = 1.e-25;
693: #else
694:   ksp->abstol  = 1.e-50;
695: #endif
696:   ksp->divtol  = 1.e4;

698:   ksp->chknorm        = -1;
699:   ksp->normtype       = ksp->normtype_set = KSP_NORM_DEFAULT;
700:   ksp->rnorm          = 0.0;
701:   ksp->its            = 0;
702:   ksp->guess_zero     = PETSC_TRUE;
703:   ksp->calc_sings     = PETSC_FALSE;
704:   ksp->res_hist       = NULL;
705:   ksp->res_hist_alloc = NULL;
706:   ksp->res_hist_len   = 0;
707:   ksp->res_hist_max   = 0;
708:   ksp->res_hist_reset = PETSC_TRUE;
709:   ksp->numbermonitors = 0;
710:   ksp->setfromoptionscalled = 0;

712:   KSPConvergedDefaultCreate(&ctx);
713:   KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
714:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
715:   ksp->ops->buildresidual = KSPBuildResidualDefault;

717:   ksp->vec_sol    = NULL;
718:   ksp->vec_rhs    = NULL;
719:   ksp->pc         = NULL;
720:   ksp->data       = NULL;
721:   ksp->nwork      = 0;
722:   ksp->work       = NULL;
723:   ksp->reason     = KSP_CONVERGED_ITERATING;
724:   ksp->setupstage = KSP_SETUP_NEW;

726:   KSPNormSupportTableReset_Private(ksp);

728:   *inksp = ksp;
729:   return(0);
730: }

732: /*@C
733:    KSPSetType - Builds KSP for a particular solver.

735:    Logically Collective on ksp

737:    Input Parameters:
738: +  ksp      - the Krylov space context
739: -  type - a known method

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

745:    Notes:
746:    See "petsc/include/petscksp.h" for available methods (for instance,
747:    KSPCG or KSPGMRES).

749:   Normally, it is best to use the KSPSetFromOptions() command and
750:   then set the KSP type from the options database rather than by using
751:   this routine.  Using the options database provides the user with
752:   maximum flexibility in evaluating the many different Krylov methods.
753:   The KSPSetType() routine is provided for those situations where it
754:   is necessary to set the iterative solver independently of the command
755:   line or options database.  This might be the case, for example, when
756:   the choice of iterative solver changes during the execution of the
757:   program, and the user's application is taking responsibility for
758:   choosing the appropriate method.  In other words, this routine is
759:   not for beginners.

761:   Level: intermediate

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

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

768: @*/
769: PetscErrorCode  KSPSetType(KSP ksp, KSPType type)
770: {
771:   PetscErrorCode ierr,(*r)(KSP);
772:   PetscBool      match;


778:   PetscObjectTypeCompare((PetscObject)ksp,type,&match);
779:   if (match) return(0);

781:   PetscFunctionListFind(KSPList,type,&r);
782:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested KSP type %s",type);
783:   /* Destroy the previous private KSP context */
784:   if (ksp->ops->destroy) {
785:     (*ksp->ops->destroy)(ksp);
786:     ksp->ops->destroy = NULL;
787:   }
788:   /* Reinitialize function pointers in KSPOps structure */
789:   PetscMemzero(ksp->ops,sizeof(struct _KSPOps));
790:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
791:   ksp->ops->buildresidual = KSPBuildResidualDefault;
792:   KSPNormSupportTableReset_Private(ksp);
793:   /* Call the KSPCreate_XXX routine for this particular Krylov solver */
794:   ksp->setupstage = KSP_SETUP_NEW;
795:   PetscObjectChangeTypeName((PetscObject)ksp,type);
796:   (*r)(ksp);
797:   return(0);
798: }

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

803:    Not Collective

805:    Input Parameter:
806: .  ksp - Krylov context

808:    Output Parameter:
809: .  name - name of KSP method

811:    Level: intermediate

813: .seealso: KSPSetType()
814: @*/
815: PetscErrorCode  KSPGetType(KSP ksp,KSPType *type)
816: {
820:   *type = ((PetscObject)ksp)->type_name;
821:   return(0);
822: }

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

827:    Not Collective

829:    Input Parameters:
830: +  name_solver - name of a new user-defined solver
831: -  routine_create - routine to create method context

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

836:    Sample usage:
837: .vb
838:    KSPRegister("my_solver",MySolverCreate);
839: .ve

841:    Then, your solver can be chosen with the procedural interface via
842: $     KSPSetType(ksp,"my_solver")
843:    or at runtime via the option
844: $     -ksp_type my_solver

846:    Level: advanced

848: .seealso: KSPRegisterAll()
849: @*/
850: PetscErrorCode  KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
851: {

855:   KSPInitializePackage();
856:   PetscFunctionListAdd(&KSPList,sname,function);
857:   return(0);
858: }