Actual source code: itcreate.c

petsc-3.12.5 2020-03-29
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              = 0;
 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,PETSC_FALSE);
163:       PetscStrncpy(type,((PetscObject)ksp)->type_name,256);
164:       PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
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->skippcsetfromoptions) {
216:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
217:     PCView(ksp->pc,viewer);
218:   }
219:   if (isdraw) {
220:     PetscDraw draw;
221:     PetscViewerDrawGetDraw(viewer,0,&draw);
222:     PetscDrawPopCurrentPoint(draw);
223:   }
224:   return(0);
225: }


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

231:    Logically Collective on ksp

233:    Input Parameter:
234: +  ksp - Krylov solver context
235: -  normtype - one of
236: $   KSP_NORM_NONE - skips computing the norm, this should generally only be used if you are using
237: $                 the Krylov method as a smoother with a fixed small number of iterations.
238: $                 Implicitly sets KSPConvergedSkip() as KSP convergence test.
239: $                 Note that certain algorithms such as KSPGMRES ALWAYS require the norm calculation,
240: $                 for these methods the norms are still computed, they are just not used in
241: $                 the convergence test. 
242: $   KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
243: $                 of the preconditioned residual P^{-1}(b - A x)
244: $   KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
245: $   KSP_NORM_NATURAL - supported  by KSPCG, KSPCR, KSPCGNE, KSPCGS


248:    Options Database Key:
249: .   -ksp_norm_type <none,preconditioned,unpreconditioned,natural>

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

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

259:    Level: advanced

261: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
262: @*/
263: PetscErrorCode  KSPSetNormType(KSP ksp,KSPNormType normtype)
264: {
268:   ksp->normtype = ksp->normtype_set = normtype;
269:   return(0);
270: }

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

276:    Logically Collective on ksp

278:    Input Parameter:
279: +  ksp - Krylov solver context
280: -  it  - use -1 to check at all iterations

282:    Notes:
283:    Currently only works with KSPCG, KSPBCGS and KSPIBCGS

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

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

291: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType()
292: @*/
293: PetscErrorCode  KSPSetCheckNormIteration(KSP ksp,PetscInt it)
294: {
298:   ksp->chknorm = it;
299:   return(0);
300: }

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


308:    Logically Collective on ksp

310:    Input Parameter:
311: +  ksp - Krylov solver context
312: -  flg - PETSC_TRUE or PETSC_FALSE

314:    Options Database Keys:
315: .  -ksp_lag_norm - lag the calculated residual norm

317:    Notes:
318:    Currently only works with KSPIBCGS.

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

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

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

336: /*@
337:    KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP

339:    Logically Collective

341:    Input Arguments:
342: +  ksp - Krylov method
343: .  normtype - supported norm type
344: .  pcside - preconditioner side that can be used with this norm
345: -  priority - positive integer preference for this combination; larger values have higher priority

347:    Level: developer

349:    Notes:
350:    This function should be called from the implementation files KSPCreate_XXX() to declare
351:    which norms and preconditioner sides are supported. Users should not need to call this
352:    function.

354: .seealso: KSPSetNormType(), KSPSetPCSide()
355: @*/
356: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
357: {

361:   ksp->normsupporttable[normtype][pcside] = priority;
362:   return(0);
363: }

365: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
366: {

370:   PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
371:   ksp->pc_side  = ksp->pc_side_set;
372:   ksp->normtype = ksp->normtype_set;
373:   return(0);
374: }

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

381:   best = 0;
382:   for (i=0; i<KSP_NORM_MAX; i++) {
383:     for (j=0; j<PC_SIDE_MAX; j++) {
384:       if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
385:         best  = ksp->normsupporttable[i][j];
386:         ibest = i;
387:         jbest = j;
388:       }
389:     }
390:   }
391:   if (best < 1 && errorifnotsupported) {
392:     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);
393:     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]);
394:     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]);
395:     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]);
396:   }
397:   if (normtype) *normtype = (KSPNormType)ibest;
398:   if (pcside)   *pcside   = (PCSide)jbest;
399:   return(0);
400: }

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

405:    Not Collective

407:    Input Parameter:
408: .  ksp - Krylov solver context

410:    Output Parameter:
411: .  normtype - norm that is used for convergence testing

413:    Level: advanced

415: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
416: @*/
417: PetscErrorCode  KSPGetNormType(KSP ksp, KSPNormType *normtype)
418: {

424:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
425:   *normtype = ksp->normtype;
426:   return(0);
427: }

429: #if defined(PETSC_HAVE_SAWS)
430:  #include <petscviewersaws.h>
431: #endif

433: /*@
434:    KSPSetOperators - Sets the matrix associated with the linear system
435:    and a (possibly) different one associated with the preconditioner.

437:    Collective on ksp

439:    Input Parameters:
440: +  ksp - the KSP context
441: .  Amat - the matrix that defines the linear system
442: -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

444:    Notes:

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

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

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

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

457:     Level: beginner

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

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

470: $         MatCreate(comm,&mat);
471: $         KSP/PCSetOperators(ksp/pc,mat,mat);
472: $         PetscObjectDereference((PetscObject)mat);
473: $           set size, type, etc of mat

475:      and

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

480: $         MatCreate(comm,&mat);
481: $         MatCreate(comm,&pmat);
482: $         KSP/PCSetOperators(ksp/pc,mat,pmat);
483: $         PetscObjectDereference((PetscObject)mat);
484: $         PetscObjectDereference((PetscObject)pmat);
485: $           set size, type, etc of mat and pmat

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

496: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
497: @*/
498: PetscErrorCode  KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
499: {

508:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
509:   PCSetOperators(ksp->pc,Amat,Pmat);
510:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;  /* so that next solve call will call PCSetUp() on new matrix */
511:   return(0);
512: }

514: /*@
515:    KSPGetOperators - Gets the matrix associated with the linear system
516:    and a (possibly) different one associated with the preconditioner.

518:    Collective on ksp

520:    Input Parameter:
521: .  ksp - the KSP context

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

527:     Level: intermediate

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

532: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
533: @*/
534: PetscErrorCode  KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
535: {

540:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
541:   PCGetOperators(ksp->pc,Amat,Pmat);
542:   return(0);
543: }

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

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

551:    Input Parameter:
552: .  pc - the KSP context

554:    Output Parameters:
555: +  mat - the matrix associated with the linear system was set
556: -  pmat - matrix associated with the preconditioner was set, usually the same

558:    Level: intermediate

560: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
561: @*/
562: PetscErrorCode  KSPGetOperatorsSet(KSP ksp,PetscBool  *mat,PetscBool  *pmat)
563: {

568:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
569:   PCGetOperatorsSet(ksp->pc,mat,pmat);
570:   return(0);
571: }

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

576:    Logically Collective on ksp

578:    Input Parameters:
579: +   ksp - the solver object
580: .   presolve - the function to call before the solve
581: -   prectx - any context needed by the function

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

586: +  ksp - the KSP context
587: .  rhs - the right-hand side vector
588: .  x - the solution vector
589: -  ctx - optional user-provided context

591:    Level: developer

593: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
594: @*/
595: PetscErrorCode  KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
596: {
599:   ksp->presolve = presolve;
600:   ksp->prectx   = prectx;
601:   return(0);
602: }

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

607:    Logically Collective on ksp

609:    Input Parameters:
610: +   ksp - the solver object
611: .   postsolve - the function to call after the solve
612: -   postctx - any context needed by the function

614:    Level: developer

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

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

624: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPreSolve()
625: @*/
626: PetscErrorCode  KSPSetPostSolve(KSP ksp,PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*),void *postctx)
627: {
630:   ksp->postsolve = postsolve;
631:   ksp->postctx   = postctx;
632:   return(0);
633: }

635: /*@
636:    KSPCreate - Creates the default KSP context.

638:    Collective

640:    Input Parameter:
641: .  comm - MPI communicator

643:    Output Parameter:
644: .  ksp - location to put the KSP context

646:    Notes:
647:    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
648:    orthogonalization.

650:    Level: beginner

652: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
653: @*/
654: PetscErrorCode  KSPCreate(MPI_Comm comm,KSP *inksp)
655: {
656:   KSP            ksp;
658:   void           *ctx;

662:   *inksp = 0;
663:   KSPInitializePackage();

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

667:   ksp->max_it  = 10000;
668:   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
669:   ksp->rtol    = 1.e-5;
670: #if defined(PETSC_USE_REAL_SINGLE)
671:   ksp->abstol  = 1.e-25;
672: #else
673:   ksp->abstol  = 1.e-50;
674: #endif
675:   ksp->divtol  = 1.e4;

677:   ksp->chknorm        = -1;
678:   ksp->normtype       = ksp->normtype_set = KSP_NORM_DEFAULT;
679:   ksp->rnorm          = 0.0;
680:   ksp->its            = 0;
681:   ksp->guess_zero     = PETSC_TRUE;
682:   ksp->calc_sings     = PETSC_FALSE;
683:   ksp->res_hist       = NULL;
684:   ksp->res_hist_alloc = NULL;
685:   ksp->res_hist_len   = 0;
686:   ksp->res_hist_max   = 0;
687:   ksp->res_hist_reset = PETSC_TRUE;
688:   ksp->numbermonitors = 0;
689:   ksp->setfromoptionscalled = 0;

691:   KSPConvergedDefaultCreate(&ctx);
692:   KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
693:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
694:   ksp->ops->buildresidual = KSPBuildResidualDefault;

696:   ksp->vec_sol    = 0;
697:   ksp->vec_rhs    = 0;
698:   ksp->pc         = 0;
699:   ksp->data       = 0;
700:   ksp->nwork      = 0;
701:   ksp->work       = 0;
702:   ksp->reason     = KSP_CONVERGED_ITERATING;
703:   ksp->setupstage = KSP_SETUP_NEW;

705:   KSPNormSupportTableReset_Private(ksp);

707:   *inksp = ksp;
708:   return(0);
709: }

711: /*@C
712:    KSPSetType - Builds KSP for a particular solver.

714:    Logically Collective on ksp

716:    Input Parameters:
717: +  ksp      - the Krylov space context
718: -  type - a known method

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

724:    Notes:
725:    See "petsc/include/petscksp.h" for available methods (for instance,
726:    KSPCG or KSPGMRES).

728:   Normally, it is best to use the KSPSetFromOptions() command and
729:   then set the KSP type from the options database rather than by using
730:   this routine.  Using the options database provides the user with
731:   maximum flexibility in evaluating the many different Krylov methods.
732:   The KSPSetType() routine is provided for those situations where it
733:   is necessary to set the iterative solver independently of the command
734:   line or options database.  This might be the case, for example, when
735:   the choice of iterative solver changes during the execution of the
736:   program, and the user's application is taking responsibility for
737:   choosing the appropriate method.  In other words, this routine is
738:   not for beginners.

740:   Level: intermediate

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

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

747: @*/
748: PetscErrorCode  KSPSetType(KSP ksp, KSPType type)
749: {
750:   PetscErrorCode ierr,(*r)(KSP);
751:   PetscBool      match;


757:   PetscObjectTypeCompare((PetscObject)ksp,type,&match);
758:   if (match) return(0);

760:   PetscFunctionListFind(KSPList,type,&r);
761:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested KSP type %s",type);
762:   /* Destroy the previous private KSP context */
763:   if (ksp->ops->destroy) {
764:     (*ksp->ops->destroy)(ksp);
765:     ksp->ops->destroy = NULL;
766:   }
767:   /* Reinitialize function pointers in KSPOps structure */
768:   PetscMemzero(ksp->ops,sizeof(struct _KSPOps));
769:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
770:   ksp->ops->buildresidual = KSPBuildResidualDefault;
771:   KSPNormSupportTableReset_Private(ksp);
772:   /* Call the KSPCreate_XXX routine for this particular Krylov solver */
773:   ksp->setupstage = KSP_SETUP_NEW;
774:   PetscObjectChangeTypeName((PetscObject)ksp,type);
775:   (*r)(ksp);
776:   return(0);
777: }

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

782:    Not Collective

784:    Input Parameter:
785: .  ksp - Krylov context

787:    Output Parameter:
788: .  name - name of KSP method

790:    Level: intermediate

792: .seealso: KSPSetType()
793: @*/
794: PetscErrorCode  KSPGetType(KSP ksp,KSPType *type)
795: {
799:   *type = ((PetscObject)ksp)->type_name;
800:   return(0);
801: }

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

806:    Not Collective

808:    Input Parameters:
809: +  name_solver - name of a new user-defined solver
810: -  routine_create - routine to create method context

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

815:    Sample usage:
816: .vb
817:    KSPRegister("my_solver",MySolverCreate);
818: .ve

820:    Then, your solver can be chosen with the procedural interface via
821: $     KSPSetType(ksp,"my_solver")
822:    or at runtime via the option
823: $     -ksp_type my_solver

825:    Level: advanced

827: .seealso: KSPRegisterAll()
828: @*/
829: PetscErrorCode  KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
830: {

834:   KSPInitializePackage();
835:   PetscFunctionListAdd(&KSPList,sname,function);
836:   return(0);
837: }