Actual source code: itcreate.c

petsc-3.5.4 2015-05-23
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,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,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      isams;
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,&isams);
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->nullsp) {PetscViewerASCIIPrintf(viewer,"  has attached null space\n");}
151:     if (!ksp->guess_zero) {PetscViewerASCIIPrintf(viewer,"  using nonzero initial guess\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 (isdraw) {
170:     PetscDraw draw;
171:     char      str[36];
172:     PetscReal x,y,bottom,h;
173:     PetscBool flg;

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

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

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


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

226:    Logically Collective on KSP

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


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

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

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

251:    Level: advanced

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

255: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration()
256: @*/
257: PetscErrorCode  KSPSetNormType(KSP ksp,KSPNormType normtype)
258: {

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

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

279:    Logically Collective on KSP

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

285:    Notes:
286:    Currently only works with KSPCG, KSPBCGS and KSPIBCGS

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

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

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

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

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


315:    Logically Collective on KSP

317:    Input Parameter:
318: +  ksp - Krylov solver context
319: -  flg - PETSC_TRUE or PETSC_FALSE

321:    Options Database Keys:
322: .  -ksp_lag_norm - lag the calculated residual norm

324:    Notes:
325:    Currently only works with KSPIBCGS.

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

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

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

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

347: /*@
348:    KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP

350:    Logically Collective

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

358:    Level: developer

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

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

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

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

382: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
383: {

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

397: PetscErrorCode KSPSetUpNorms_Private(KSP ksp,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)
406:           && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j)
407:           && (ksp->normsupporttable[i][j] > best)) {
408:         best  = ksp->normsupporttable[i][j];
409:         ibest = i;
410:         jbest = j;
411:       }
412:     }
413:   }
414:   if (best < 1) {
415:     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);
416:     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]);
417:     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]);
418:     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]);
419:   }
420:   *normtype = (KSPNormType)ibest;
421:   *pcside   = (PCSide)jbest;
422:   return(0);
423: }

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

430:    Not Collective

432:    Input Parameter:
433: .  ksp - Krylov solver context

435:    Output Parameter:
436: .  normtype - norm that is used for convergence testing

438:    Level: advanced

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

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

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

456: #if defined(PETSC_HAVE_SAWS)
457: #include <petscviewersaws.h>
458: #endif

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

466:    Collective on KSP and Mat

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

473:    Notes:

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

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

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

483:     Level: beginner

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

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

496: $         MatCreate(comm,&mat);
497: $         KSP/PCSetOperators(ksp/pc,mat,mat);
498: $         PetscObjectDereference((PetscObject)mat);
499: $           set size, type, etc of mat

501:      and

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

506: $         MatCreate(comm,&mat);
507: $         MatCreate(comm,&pmat);
508: $         KSP/PCSetOperators(ksp/pc,mat,pmat);
509: $         PetscObjectDereference((PetscObject)mat);
510: $         PetscObjectDereference((PetscObject)pmat);
511: $           set size, type, etc of mat and pmat

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

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

524: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
525: @*/
526: PetscErrorCode  KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
527: {
528:   MatNullSpace   nullsp;

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

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

558:    Collective on KSP and Mat

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

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

567:     Level: intermediate

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

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

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

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

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

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

594:    Input Parameter:
595: .  pc - the KSP context

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

601:    Level: intermediate

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

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

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

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

623:    Logically Collective on KSP

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

630:    Level: developer

632: .keywords: KSP, create, context

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

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

650:    Logically Collective on KSP

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

657:    Level: developer

659: .keywords: KSP, create, context

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

674: /*@
675:    KSPCreate - Creates the default KSP context.

677:    Collective on MPI_Comm

679:    Input Parameter:
680: .  comm - MPI communicator

682:    Output Parameter:
683: .  ksp - location to put the KSP context

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

689:    Level: beginner

691: .keywords: KSP, create, context

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

703:   *inksp = 0;
704:   KSPInitializePackage();

706:   PetscHeaderCreate(ksp,_p_KSP,struct _KSPOps,KSP_CLASSID,"KSP","Krylov Method","KSP",comm,KSPDestroy,KSPView);

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

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

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

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

745:   KSPNormSupportTableReset_Private(ksp);

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

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

756:    Logically Collective on KSP

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

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

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

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

782:   Level: intermediate

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

787: .keywords: KSP, set, method

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

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


801:   PetscObjectTypeCompare((PetscObject)ksp,type,&match);
802:   if (match) return(0);

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

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

828:    Not Collective

830:    Input Parameter:
831: .  ksp - Krylov context

833:    Output Parameter:
834: .  name - name of KSP method

836:    Level: intermediate

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

840: .seealso: KSPSetType()
841: @*/
842: PetscErrorCode  KSPGetType(KSP ksp,KSPType *type)
843: {
847:   *type = ((PetscObject)ksp)->type_name;
848:   return(0);
849: }

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

856:    Not Collective

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

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

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

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

875:    Level: advanced

877: .keywords: KSP, register

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

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

887:   PetscFunctionListAdd(&KSPList,sname,function);
888:   return(0);
889: }

893: /*@
894:   KSPSetNullSpace - Sets the null space of the operator

896:   Logically Collective on KSP

898:   Input Parameters:
899: +  ksp - the Krylov space object
900: -  nullsp - the null space of the operator

902:   Notes: If the Mat provided to KSP has a nullspace added to it with MatSetNullSpace() then
903:          KSP will automatically use the MatNullSpace and you don't need to call KSPSetNullSpace().

905:   Level: advanced

907: .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPGetNullSpace(), MatSetNullSpace()
908: @*/
909: PetscErrorCode  KSPSetNullSpace(KSP ksp,MatNullSpace nullsp)
910: {

916:   PetscObjectReference((PetscObject)nullsp);
917:   if (ksp->nullsp) { MatNullSpaceDestroy(&ksp->nullsp); }
918:   ksp->nullsp = nullsp;
919:   return(0);
920: }

924: /*@
925:   KSPGetNullSpace - Gets the null space of the operator

927:   Not Collective

929:   Input Parameters:
930: +  ksp - the Krylov space object
931: -  nullsp - the null space of the operator

933:   Level: advanced

935: .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPSetNullSpace()
936: @*/
937: PetscErrorCode  KSPGetNullSpace(KSP ksp,MatNullSpace *nullsp)
938: {
942:   *nullsp = ksp->nullsp;
943:   return(0);
944: }