Actual source code: itcreate.c

petsc-3.7.7 2017-09-25
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) {
118:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
119:   }

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

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

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

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

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


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

229:    Logically Collective on KSP

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


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

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

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

254:    Level: advanced

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

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

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

275:    Logically Collective on KSP

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

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

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

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

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

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

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


311:    Logically Collective on KSP

313:    Input Parameter:
314: +  ksp - Krylov solver context
315: -  flg - PETSC_TRUE or PETSC_FALSE

317:    Options Database Keys:
318: .  -ksp_lag_norm - lag the calculated residual norm

320:    Notes:
321:    Currently only works with KSPIBCGS.

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

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

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

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

343: /*@
344:    KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP

346:    Logically Collective

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

354:    Level: developer

356:    Notes:
357:    This function should be called from the implementation files KSPCreate_XXX() to declare
358:    which norms and preconditioner sides are supported. Users should not need to call this
359:    function.

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

365: .seealso: KSPSetNormType(), KSPSetPCSide()
366: @*/
367: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
368: {

372:   ksp->normsupporttable[normtype][pcside] = priority;
373:   return(0);
374: }

378: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
379: {

383:   PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
384:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
385:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
386:   ksp->pc_side  = ksp->pc_side_set;
387:   ksp->normtype = ksp->normtype_set;
388:   return(0);
389: }

393: PetscErrorCode KSPSetUpNorms_Private(KSP ksp,KSPNormType *normtype,PCSide *pcside)
394: {
395:   PetscInt i,j,best,ibest = 0,jbest = 0;

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

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

426:    Not Collective

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

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

434:    Level: advanced

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

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

447:   KSPSetUpNorms_Private(ksp,&ksp->normtype,&ksp->pc_side);
448:   *normtype = ksp->normtype;
449:   return(0);
450: }

452: #if defined(PETSC_HAVE_SAWS)
453: #include <petscviewersaws.h>
454: #endif

458: /*@
459:    KSPSetOperators - Sets the matrix associated with the linear system
460:    and a (possibly) different one associated with the preconditioner.

462:    Collective on KSP and Mat

464:    Input Parameters:
465: +  ksp - the KSP context
466: .  Amat - the matrix that defines the linear system
467: -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

469:    Notes:

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

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

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

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

482:     Level: beginner

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

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

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

500:      and

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

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

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

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

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

535:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
536:   PCSetOperators(ksp->pc,Amat,Pmat);
537:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;  /* so that next solve call will call PCSetUp() on new matrix */
538:   if (ksp->guess) {
539:     KSPFischerGuessReset(ksp->guess);
540:   }
541:   return(0);
542: }

546: /*@
547:    KSPGetOperators - Gets the matrix associated with the linear system
548:    and a (possibly) different one associated with the preconditioner.

550:    Collective on KSP and Mat

552:    Input Parameter:
553: .  ksp - the KSP context

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

559:     Level: intermediate

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

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

565: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
566: @*/
567: PetscErrorCode  KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
568: {

573:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
574:   PCGetOperators(ksp->pc,Amat,Pmat);
575:   return(0);
576: }

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

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

586:    Input Parameter:
587: .  pc - the KSP context

589:    Output Parameters:
590: +  mat - the matrix associated with the linear system was set
591: -  pmat - matrix associated with the preconditioner was set, usually the same

593:    Level: intermediate

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

597: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
598: @*/
599: PetscErrorCode  KSPGetOperatorsSet(KSP ksp,PetscBool  *mat,PetscBool  *pmat)
600: {

605:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
606:   PCGetOperatorsSet(ksp->pc,mat,pmat);
607:   return(0);
608: }

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

615:    Logically Collective on KSP

617:    Input Parameters:
618: +   ksp - the solver object
619: .   presolve - the function to call before the solve
620: -   prectx - any context needed by the function

622:    Level: developer

624: .keywords: KSP, create, context

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

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

642:    Logically Collective on KSP

644:    Input Parameters:
645: +   ksp - the solver object
646: .   postsolve - the function to call after the solve
647: -   postctx - any context needed by the function

649:    Level: developer

651: .keywords: KSP, create, context

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

666: /*@
667:    KSPCreate - Creates the default KSP context.

669:    Collective on MPI_Comm

671:    Input Parameter:
672: .  comm - MPI communicator

674:    Output Parameter:
675: .  ksp - location to put the KSP context

677:    Notes:
678:    The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
679:    orthogonalization.

681:    Level: beginner

683: .keywords: KSP, create, context

685: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
686: @*/
687: PetscErrorCode  KSPCreate(MPI_Comm comm,KSP *inksp)
688: {
689:   KSP            ksp;
691:   void           *ctx;

695:   *inksp = 0;
696:   KSPInitializePackage();

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

700:   ksp->max_it  = 10000;
701:   ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
702:   ksp->rtol    = 1.e-5;
703: #if defined(PETSC_USE_REAL_SINGLE)
704:   ksp->abstol  = 1.e-25;
705: #else
706:   ksp->abstol  = 1.e-50;
707: #endif
708:   ksp->divtol  = 1.e4;

710:   ksp->chknorm        = -1;
711:   ksp->normtype       = ksp->normtype_set = KSP_NORM_DEFAULT;
712:   ksp->rnorm          = 0.0;
713:   ksp->its            = 0;
714:   ksp->guess_zero     = PETSC_TRUE;
715:   ksp->calc_sings     = PETSC_FALSE;
716:   ksp->res_hist       = NULL;
717:   ksp->res_hist_alloc = NULL;
718:   ksp->res_hist_len   = 0;
719:   ksp->res_hist_max   = 0;
720:   ksp->res_hist_reset = PETSC_TRUE;
721:   ksp->numbermonitors = 0;

723:   KSPConvergedDefaultCreate(&ctx);
724:   KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
725:   ksp->ops->buildsolution = KSPBuildSolutionDefault;
726:   ksp->ops->buildresidual = KSPBuildResidualDefault;

728:   ksp->vec_sol    = 0;
729:   ksp->vec_rhs    = 0;
730:   ksp->pc         = 0;
731:   ksp->data       = 0;
732:   ksp->nwork      = 0;
733:   ksp->work       = 0;
734:   ksp->reason     = KSP_CONVERGED_ITERATING;
735:   ksp->setupstage = KSP_SETUP_NEW;

737:   KSPNormSupportTableReset_Private(ksp);

739:   *inksp = ksp;
740:   return(0);
741: }

745: /*@C
746:    KSPSetType - Builds KSP for a particular solver.

748:    Logically Collective on KSP

750:    Input Parameters:
751: +  ksp      - the Krylov space context
752: -  type - a known method

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

758:    Notes:
759:    See "petsc/include/petscksp.h" for available methods (for instance,
760:    KSPCG or KSPGMRES).

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

774:   Level: intermediate

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

779: .keywords: KSP, set, method

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

783: @*/
784: PetscErrorCode  KSPSetType(KSP ksp, KSPType type)
785: {
786:   PetscErrorCode ierr,(*r)(KSP);
787:   PetscBool      match;


793:   PetscObjectTypeCompare((PetscObject)ksp,type,&match);
794:   if (match) return(0);

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

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

820:    Not Collective

822:    Input Parameter:
823: .  ksp - Krylov context

825:    Output Parameter:
826: .  name - name of KSP method

828:    Level: intermediate

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

832: .seealso: KSPSetType()
833: @*/
834: PetscErrorCode  KSPGetType(KSP ksp,KSPType *type)
835: {
839:   *type = ((PetscObject)ksp)->type_name;
840:   return(0);
841: }

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

848:    Not Collective

850:    Input Parameters:
851: +  name_solver - name of a new user-defined solver
852: -  routine_create - routine to create method context

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

857:    Sample usage:
858: .vb
859:    KSPRegister("my_solver",MySolverCreate);
860: .ve

862:    Then, your solver can be chosen with the procedural interface via
863: $     KSPSetType(ksp,"my_solver")
864:    or at runtime via the option
865: $     -ksp_type my_solver

867:    Level: advanced

869: .keywords: KSP, register

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

873: @*/
874: PetscErrorCode  KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
875: {

879:   PetscFunctionListAdd(&KSPList,sname,function);
880:   return(0);
881: }