Actual source code: itcreate.c

petsc-master 2017-06-24
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;

 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 PetscViewer

 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: .keywords: KSP, view

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

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

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

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

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

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

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


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

228:    Logically Collective on KSP

230:    Input Parameter:
231: +  ksp - Krylov solver context
232: -  normtype - one of
233: $   KSP_NORM_NONE - skips computing the norm, this should only be used if you are using
234: $                 the Krylov method as a smoother with a fixed small number of iterations.
235: $                 Implicitly sets KSPConvergedSkip as KSP convergence test.
236: $   KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
237: $                 of the preconditioned residual P^{-1}(b - A x)
238: $   KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
239: $   KSP_NORM_NATURAL - supported  by KSPCG, KSPCR, KSPCGNE, KSPCGS


242:    Options Database Key:
243: .   -ksp_norm_type <none,preconditioned,unpreconditioned,natural>

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

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

253:    Level: advanced

255: .keywords: KSP, create, context, norms

257: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
258: @*/
259: PetscErrorCode  KSPSetNormType(KSP ksp,KSPNormType normtype)
260: {
264:   ksp->normtype = ksp->normtype_set = normtype;
265:   return(0);
266: }

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

272:    Logically Collective on KSP

274:    Input Parameter:
275: +  ksp - Krylov solver context
276: -  it  - use -1 to check at all iterations

278:    Notes:
279:    Currently only works with KSPCG, KSPBCGS and KSPIBCGS

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

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

287: .keywords: KSP, create, context, norms

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

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


306:    Logically Collective on KSP

308:    Input Parameter:
309: +  ksp - Krylov solver context
310: -  flg - PETSC_TRUE or PETSC_FALSE

312:    Options Database Keys:
313: .  -ksp_lag_norm - lag the calculated residual norm

315:    Notes:
316:    Currently only works with KSPIBCGS.

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

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

323: .keywords: KSP, create, context, norms

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: -  preference - 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:    KSP_NORM_NONE is supported by default with all KSP methods and any PC side at priority 1.  If a KSP explicitly does
355:    not support KSP_NORM_NONE, it should set this by setting priority=0.  Since defaulting to KSP_NORM_NONE is usually
356:    undesirable, more desirable norms should usually have priority 2 or higher.

358: .seealso: KSPSetNormType(), KSPSetPCSide()
359: @*/
360: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
361: {

365:   ksp->normsupporttable[normtype][pcside] = priority;
366:   return(0);
367: }

369: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
370: {

374:   PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
375:   ksp->pc_side  = ksp->pc_side_set;
376:   ksp->normtype = ksp->normtype_set;
377:   return(0);
378: }

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

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

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

409:    Not Collective

411:    Input Parameter:
412: .  ksp - Krylov solver context

414:    Output Parameter:
415: .  normtype - norm that is used for convergence testing

417:    Level: advanced

419: .keywords: KSP, create, context, norms

421: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
422: @*/
423: PetscErrorCode  KSPGetNormType(KSP ksp, KSPNormType *normtype)
424: {

430:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
431:   *normtype = ksp->normtype;
432:   return(0);
433: }

435: #if defined(PETSC_HAVE_SAWS)
436:  #include <petscviewersaws.h>
437: #endif

439: /*@
440:    KSPSetOperators - Sets the matrix associated with the linear system
441:    and a (possibly) different one associated with the preconditioner.

443:    Collective on KSP and Mat

445:    Input Parameters:
446: +  ksp - the KSP context
447: .  Amat - the matrix that defines the linear system
448: -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

450:    Notes:

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

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

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

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

463:     Level: beginner

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

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

476: $         MatCreate(comm,&mat);
477: $         KSP/PCSetOperators(ksp/pc,mat,mat);
478: $         PetscObjectDereference((PetscObject)mat);
479: $           set size, type, etc of mat

481:      and

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

486: $         MatCreate(comm,&mat);
487: $         MatCreate(comm,&pmat);
488: $         KSP/PCSetOperators(ksp/pc,mat,pmat);
489: $         PetscObjectDereference((PetscObject)mat);
490: $         PetscObjectDereference((PetscObject)pmat);
491: $           set size, type, etc of mat and pmat

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

502: .keywords: KSP, set, operators, matrix, preconditioner, linear system

504: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
505: @*/
506: PetscErrorCode  KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
507: {

516:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
517:   PCSetOperators(ksp->pc,Amat,Pmat);
518:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;  /* so that next solve call will call PCSetUp() on new matrix */
519:   return(0);
520: }

522: /*@
523:    KSPGetOperators - Gets the matrix associated with the linear system
524:    and a (possibly) different one associated with the preconditioner.

526:    Collective on KSP and Mat

528:    Input Parameter:
529: .  ksp - the KSP context

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

535:     Level: intermediate

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

539: .keywords: KSP, set, get, operators, matrix, preconditioner, linear system

541: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
542: @*/
543: PetscErrorCode  KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
544: {

549:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
550:   PCGetOperators(ksp->pc,Amat,Pmat);
551:   return(0);
552: }

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

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

560:    Input Parameter:
561: .  pc - the KSP context

563:    Output Parameters:
564: +  mat - the matrix associated with the linear system was set
565: -  pmat - matrix associated with the preconditioner was set, usually the same

567:    Level: intermediate

569: .keywords: KSP, get, operators, matrix, linear system

571: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
572: @*/
573: PetscErrorCode  KSPGetOperatorsSet(KSP ksp,PetscBool  *mat,PetscBool  *pmat)
574: {

579:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
580:   PCGetOperatorsSet(ksp->pc,mat,pmat);
581:   return(0);
582: }

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

587:    Logically Collective on KSP

589:    Input Parameters:
590: +   ksp - the solver object
591: .   presolve - the function to call before the solve
592: -   prectx - any context needed by the function

594:    Level: developer

596: .keywords: KSP, create, context

598: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
599: @*/
600: PetscErrorCode  KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
601: {
604:   ksp->presolve = presolve;
605:   ksp->prectx   = prectx;
606:   return(0);
607: }

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

612:    Logically Collective on KSP

614:    Input Parameters:
615: +   ksp - the solver object
616: .   postsolve - the function to call after the solve
617: -   postctx - any context needed by the function

619:    Level: developer

621: .keywords: KSP, create, context

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

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

637:    Collective on MPI_Comm

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

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

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

649:    Level: beginner

651: .keywords: KSP, create, context

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

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

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

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

678:   ksp->chknorm        = -1;
679:   ksp->normtype       = ksp->normtype_set = KSP_NORM_DEFAULT;
680:   ksp->rnorm          = 0.0;
681:   ksp->its            = 0;
682:   ksp->guess_zero     = PETSC_TRUE;
683:   ksp->calc_sings     = PETSC_FALSE;
684:   ksp->res_hist       = NULL;
685:   ksp->res_hist_alloc = NULL;
686:   ksp->res_hist_len   = 0;
687:   ksp->res_hist_max   = 0;
688:   ksp->res_hist_reset = PETSC_TRUE;
689:   ksp->numbermonitors = 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: .keywords: KSP, set, method

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

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


759:   PetscObjectTypeCompare((PetscObject)ksp,type,&match);
760:   if (match) return(0);

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

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

784:    Not Collective

786:    Input Parameter:
787: .  ksp - Krylov context

789:    Output Parameter:
790: .  name - name of KSP method

792:    Level: intermediate

794: .keywords: KSP, get, method, name

796: .seealso: KSPSetType()
797: @*/
798: PetscErrorCode  KSPGetType(KSP ksp,KSPType *type)
799: {
803:   *type = ((PetscObject)ksp)->type_name;
804:   return(0);
805: }

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

810:    Not Collective

812:    Input Parameters:
813: +  name_solver - name of a new user-defined solver
814: -  routine_create - routine to create method context

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

819:    Sample usage:
820: .vb
821:    KSPRegister("my_solver",MySolverCreate);
822: .ve

824:    Then, your solver can be chosen with the procedural interface via
825: $     KSPSetType(ksp,"my_solver")
826:    or at runtime via the option
827: $     -ksp_type my_solver

829:    Level: advanced

831: .keywords: KSP, register

833: .seealso: KSPRegisterAll(), KSPRegisterDestroy()

835: @*/
836: PetscErrorCode  KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
837: {

841:   PetscFunctionListAdd(&KSPList,sname,function);
842:   return(0);
843: }