Actual source code: itcreate.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: /*
  3:      The basic KSP routines, Create, View etc. are here.
  4: */
  5: #include <petsc/private/kspimpl.h>      /*I "petscksp.h" I*/

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

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

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

 23:   Collective on PetscViewer

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

 30:    Level: intermediate

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

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

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

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

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

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

 81:    Collective on KSP

 83:    Input Parameters:
 84: +  ksp - the Krylov space context
 85: -  viewer - visualization context

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

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

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

101:    Level: beginner

103: .keywords: KSP, view

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

117:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));

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

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

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

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

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


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

227:    Logically Collective on KSP

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


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

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

249:    Developer Notes:
250:    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()
258: @*/
259: PetscErrorCode  KSPSetNormType(KSP ksp,KSPNormType normtype)
260: {

266:   ksp->normtype = ksp->normtype_set = normtype;
267:   if (normtype == KSP_NORM_NONE) {
268:     KSPSetConvergenceTest(ksp,KSPConvergedSkip,0,0);
269:     PetscInfo(ksp,"Warning: setting KSPNormType to skip computing the norm\n\
270:  KSP convergence test is implicitly set to KSPConvergedSkip\n");
271:   }
272:   return(0);
273: }

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

281:    Logically Collective on KSP

283:    Input Parameter:
284: +  ksp - Krylov solver context
285: -  it  - use -1 to check at all iterations

287:    Notes:
288:    Currently only works with KSPCG, KSPBCGS and KSPIBCGS

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

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

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

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

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


317:    Logically Collective on KSP

319:    Input Parameter:
320: +  ksp - Krylov solver context
321: -  flg - PETSC_TRUE or PETSC_FALSE

323:    Options Database Keys:
324: .  -ksp_lag_norm - lag the calculated residual norm

326:    Notes:
327:    Currently only works with KSPIBCGS.

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

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

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

336: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType(), KSPSetCheckNormIteration()
337: @*/
338: PetscErrorCode  KSPSetLagNorm(KSP ksp,PetscBool flg)
339: {
343:   ksp->lagnorm = flg;
344:   return(0);
345: }

349: /*@
350:    KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP

352:    Logically Collective

354:    Input Arguments:
355: +  ksp - Krylov method
356: .  normtype - supported norm type
357: .  pcside - preconditioner side that can be used with this norm
358: -  preference - integer preference for this combination, larger values have higher priority

360:    Level: developer

362:    Notes:
363:    This function should be called from the implementation files KSPCreate_XXX() to declare
364:    which norms and preconditioner sides are supported. Users should not need to call this
365:    function.

367:    KSP_NORM_NONE is supported by default with all KSP methods and any PC side at priority 1.  If a KSP explicitly does
368:    not support KSP_NORM_NONE, it should set this by setting priority=0.  Since defaulting to KSP_NORM_NONE is usually
369:    undesirable, more desirable norms should usually have priority 2 or higher.

371: .seealso: KSPSetNormType(), KSPSetPCSide()
372: @*/
373: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
374: {

378:   ksp->normsupporttable[normtype][pcside] = priority;
379:   return(0);
380: }

384: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
385: {

389:   PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
390:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
391:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
392:   ksp->pc_side  = ksp->pc_side_set;
393:   ksp->normtype = ksp->normtype_set;
394:   return(0);
395: }

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

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

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

432:    Not Collective

434:    Input Parameter:
435: .  ksp - Krylov solver context

437:    Output Parameter:
438: .  normtype - norm that is used for convergence testing

440:    Level: advanced

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

444: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
445: @*/
446: PetscErrorCode  KSPGetNormType(KSP ksp, KSPNormType *normtype)
447: {

453:   KSPSetUpNorms_Private(ksp,&ksp->normtype,&ksp->pc_side);
454:   *normtype = ksp->normtype;
455:   return(0);
456: }

458: #if defined(PETSC_HAVE_SAWS)
459: #include <petscviewersaws.h>
460: #endif

464: /*@
465:    KSPSetOperators - Sets the matrix associated with the linear system
466:    and a (possibly) different one associated with the preconditioner.

468:    Collective on KSP and Mat

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

475:    Notes:

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

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

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

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

488:     Level: beginner

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

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

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

506:      and

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

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

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

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

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

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

552: /*@
553:    KSPGetOperators - Gets the matrix associated with the linear system
554:    and a (possibly) different one associated with the preconditioner.

556:    Collective on KSP and Mat

558:    Input Parameter:
559: .  ksp - the KSP context

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

565:     Level: intermediate

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

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

571: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
572: @*/
573: PetscErrorCode  KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
574: {

579:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
580:   PCGetOperators(ksp->pc,Amat,Pmat);
581:   return(0);
582: }

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

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

592:    Input Parameter:
593: .  pc - the KSP context

595:    Output Parameters:
596: +  mat - the matrix associated with the linear system was set
597: -  pmat - matrix associated with the preconditioner was set, usually the same

599:    Level: intermediate

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

603: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
604: @*/
605: PetscErrorCode  KSPGetOperatorsSet(KSP ksp,PetscBool  *mat,PetscBool  *pmat)
606: {

611:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
612:   PCGetOperatorsSet(ksp->pc,mat,pmat);
613:   return(0);
614: }

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

621:    Logically Collective on KSP

623:    Input Parameters:
624: +   ksp - the solver object
625: .   presolve - the function to call before the solve
626: -   prectx - any context needed by the function

628:    Level: developer

630: .keywords: KSP, create, context

632: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
633: @*/
634: PetscErrorCode  KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
635: {
638:   ksp->presolve = presolve;
639:   ksp->prectx   = prectx;
640:   return(0);
641: }

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

648:    Logically Collective on KSP

650:    Input Parameters:
651: +   ksp - the solver object
652: .   postsolve - the function to call after the solve
653: -   postctx - any context needed by the function

655:    Level: developer

657: .keywords: KSP, create, context

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

672: /*@
673:    KSPCreate - Creates the default KSP context.

675:    Collective on MPI_Comm

677:    Input Parameter:
678: .  comm - MPI communicator

680:    Output Parameter:
681: .  ksp - location to put the KSP context

683:    Notes:
684:    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
685:    orthogonalization.

687:    Level: beginner

689: .keywords: KSP, create, context

691: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
692: @*/
693: PetscErrorCode  KSPCreate(MPI_Comm comm,KSP *inksp)
694: {
695:   KSP            ksp;
697:   void           *ctx;

701:   *inksp = 0;
702:   KSPInitializePackage();

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

706:   ksp->max_it  = 10000;
707:   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
708:   ksp->rtol    = 1.e-5;
709: #if defined(PETSC_USE_REAL_SINGLE)
710:   ksp->abstol  = 1.e-25;
711: #else
712:   ksp->abstol  = 1.e-50;
713: #endif
714:   ksp->divtol  = 1.e4;

716:   ksp->chknorm        = -1;
717:   ksp->normtype       = ksp->normtype_set = KSP_NORM_DEFAULT;
718:   ksp->rnorm          = 0.0;
719:   ksp->its            = 0;
720:   ksp->guess_zero     = PETSC_TRUE;
721:   ksp->calc_sings     = PETSC_FALSE;
722:   ksp->res_hist       = NULL;
723:   ksp->res_hist_alloc = NULL;
724:   ksp->res_hist_len   = 0;
725:   ksp->res_hist_max   = 0;
726:   ksp->res_hist_reset = PETSC_TRUE;
727:   ksp->numbermonitors = 0;

729:   KSPConvergedDefaultCreate(&ctx);
730:   KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
731:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
732:   ksp->ops->buildresidual = KSPBuildResidualDefault;

734:   ksp->vec_sol    = 0;
735:   ksp->vec_rhs    = 0;
736:   ksp->pc         = 0;
737:   ksp->data       = 0;
738:   ksp->nwork      = 0;
739:   ksp->work       = 0;
740:   ksp->reason     = KSP_CONVERGED_ITERATING;
741:   ksp->setupstage = KSP_SETUP_NEW;

743:   KSPNormSupportTableReset_Private(ksp);

745:   *inksp = ksp;
746:   return(0);
747: }

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

754:    Logically Collective on KSP

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

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

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

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

780:   Level: intermediate

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

785: .keywords: KSP, set, method

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

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


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

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

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

826:    Not Collective

828:    Input Parameter:
829: .  ksp - Krylov context

831:    Output Parameter:
832: .  name - name of KSP method

834:    Level: intermediate

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

838: .seealso: KSPSetType()
839: @*/
840: PetscErrorCode  KSPGetType(KSP ksp,KSPType *type)
841: {
845:   *type = ((PetscObject)ksp)->type_name;
846:   return(0);
847: }

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

854:    Not Collective

856:    Input Parameters:
857: +  name_solver - name of a new user-defined solver
858: -  routine_create - routine to create method context

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

863:    Sample usage:
864: .vb
865:    KSPRegister("my_solver",MySolverCreate);
866: .ve

868:    Then, your solver can be chosen with the procedural interface via
869: $     KSPSetType(ksp,"my_solver")
870:    or at runtime via the option
871: $     -ksp_type my_solver

873:    Level: advanced

875: .keywords: KSP, register

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

879: @*/
880: PetscErrorCode  KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
881: {

885:   PetscFunctionListAdd(&KSPList,sname,function);
886:   return(0);
887: }