Actual source code: itcreate.c
petsc-3.13.6 2020-09-29
2: /*
3: The basic KSP routines, Create, View etc. are here.
4: */
5: #include <petsc/private/kspimpl.h>
7: /* Logging support */
8: PetscClassId KSP_CLASSID;
9: PetscClassId DMKSP_CLASSID;
10: PetscClassId KSPGUESS_CLASSID;
11: PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve, KSP_SolveTranspose;
13: /*
14: Contains the list of registered KSP routines
15: */
16: PetscFunctionList KSPList = 0;
17: PetscBool KSPRegisterAllCalled = PETSC_FALSE;
19: /*@C
20: KSPLoad - Loads a KSP that has been stored in binary with KSPView().
22: Collective on viewer
24: Input Parameters:
25: + newdm - the newly loaded KSP, this needs to have been created with KSPCreate() or
26: some related function before a call to KSPLoad().
27: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
29: Level: intermediate
31: Notes:
32: The type is determined by the data in the file, any type set into the KSP before this call is ignored.
34: Notes for advanced users:
35: Most users should not need to know the details of the binary storage
36: format, since KSPLoad() and KSPView() completely hide these details.
37: But for anyone who's interested, the standard binary matrix storage
38: format is
39: .vb
40: has not yet been determined
41: .ve
43: .seealso: PetscViewerBinaryOpen(), KSPView(), MatLoad(), VecLoad()
44: @*/
45: PetscErrorCode KSPLoad(KSP newdm, PetscViewer viewer)
46: {
48: PetscBool isbinary;
49: PetscInt classid;
50: char type[256];
51: PC pc;
56: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
57: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
59: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
60: if (classid != KSP_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not KSP next in file");
61: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
62: KSPSetType(newdm, type);
63: if (newdm->ops->load) {
64: (*newdm->ops->load)(newdm,viewer);
65: }
66: KSPGetPC(newdm,&pc);
67: PCLoad(pc,viewer);
68: return(0);
69: }
71: #include <petscdraw.h>
72: #if defined(PETSC_HAVE_SAWS)
73: #include <petscviewersaws.h>
74: #endif
75: /*@C
76: KSPView - Prints the KSP data structure.
78: Collective on ksp
80: Input Parameters:
81: + ksp - the Krylov space context
82: - viewer - visualization context
84: Options Database Keys:
85: . -ksp_view - print the ksp data structure at the end of a KSPSolve call
87: Note:
88: The available visualization contexts include
89: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
90: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
91: output where only the first processor opens
92: the file. All other processors send their
93: data to the first processor to print.
95: The user can open an alternative visualization context with
96: PetscViewerASCIIOpen() - output to a specified file.
98: Level: beginner
100: .seealso: PCView(), PetscViewerASCIIOpen()
101: @*/
102: PetscErrorCode KSPView(KSP ksp,PetscViewer viewer)
103: {
105: PetscBool iascii,isbinary,isdraw,isstring;
106: #if defined(PETSC_HAVE_SAWS)
107: PetscBool issaws;
108: #endif
112: if (!viewer) {
113: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
114: }
118: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
119: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
120: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
121: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
122: #if defined(PETSC_HAVE_SAWS)
123: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
124: #endif
125: if (iascii) {
126: PetscObjectPrintClassNamePrefixType((PetscObject)ksp,viewer);
127: if (ksp->ops->view) {
128: PetscViewerASCIIPushTab(viewer);
129: (*ksp->ops->view)(ksp,viewer);
130: PetscViewerASCIIPopTab(viewer);
131: }
132: if (ksp->guess_zero) {
133: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, initial guess is zero\n",ksp->max_it);
134: } else {
135: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, nonzero initial guess\n", ksp->max_it);
136: }
137: if (ksp->guess_knoll) {PetscViewerASCIIPrintf(viewer," using preconditioner applied to right hand side for initial guess\n");}
138: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, divergence=%g\n",(double)ksp->rtol,(double)ksp->abstol,(double)ksp->divtol);
139: if (ksp->pc_side == PC_RIGHT) {
140: PetscViewerASCIIPrintf(viewer," right preconditioning\n");
141: } else if (ksp->pc_side == PC_SYMMETRIC) {
142: PetscViewerASCIIPrintf(viewer," symmetric preconditioning\n");
143: } else {
144: PetscViewerASCIIPrintf(viewer," left preconditioning\n");
145: }
146: if (ksp->guess) {
147: PetscViewerASCIIPushTab(viewer);
148: KSPGuessView(ksp->guess,viewer);
149: PetscViewerASCIIPopTab(viewer);
150: }
151: if (ksp->dscale) {PetscViewerASCIIPrintf(viewer," diagonally scaled system\n");}
152: PetscViewerASCIIPrintf(viewer," using %s norm type for convergence test\n",KSPNormTypes[ksp->normtype]);
153: } else if (isbinary) {
154: PetscInt classid = KSP_FILE_CLASSID;
155: MPI_Comm comm;
156: PetscMPIInt rank;
157: char type[256];
159: PetscObjectGetComm((PetscObject)ksp,&comm);
160: MPI_Comm_rank(comm,&rank);
161: if (!rank) {
162: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);
163: PetscStrncpy(type,((PetscObject)ksp)->type_name,256);
164: PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR);
165: }
166: if (ksp->ops->view) {
167: (*ksp->ops->view)(ksp,viewer);
168: }
169: } else if (isstring) {
170: const char *type;
171: KSPGetType(ksp,&type);
172: PetscViewerStringSPrintf(viewer," KSPType: %-7.7s",type);
173: if (ksp->ops->view) {(*ksp->ops->view)(ksp,viewer);}
174: } else if (isdraw) {
175: PetscDraw draw;
176: char str[36];
177: PetscReal x,y,bottom,h;
178: PetscBool flg;
180: PetscViewerDrawGetDraw(viewer,0,&draw);
181: PetscDrawGetCurrentPoint(draw,&x,&y);
182: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
183: if (!flg) {
184: PetscStrncpy(str,"KSP: ",sizeof(str));
185: PetscStrlcat(str,((PetscObject)ksp)->type_name,sizeof(str));
186: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
187: bottom = y - h;
188: } else {
189: bottom = y;
190: }
191: PetscDrawPushCurrentPoint(draw,x,bottom);
192: #if defined(PETSC_HAVE_SAWS)
193: } else if (issaws) {
194: PetscMPIInt rank;
195: const char *name;
197: PetscObjectGetName((PetscObject)ksp,&name);
198: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
199: if (!((PetscObject)ksp)->amsmem && !rank) {
200: char dir[1024];
202: PetscObjectViewSAWs((PetscObject)ksp,viewer);
203: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
204: PetscStackCallSAWs(SAWs_Register,(dir,&ksp->its,1,SAWs_READ,SAWs_INT));
205: if (!ksp->res_hist) {
206: KSPSetResidualHistory(ksp,NULL,PETSC_DECIDE,PETSC_TRUE);
207: }
208: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/res_hist",name);
209: PetscStackCallSAWs(SAWs_Register,(dir,ksp->res_hist,10,SAWs_READ,SAWs_DOUBLE));
210: }
211: #endif
212: } else if (ksp->ops->view) {
213: (*ksp->ops->view)(ksp,viewer);
214: }
215: if (!ksp->skippcsetfromoptions) {
216: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
217: PCView(ksp->pc,viewer);
218: }
219: if (isdraw) {
220: PetscDraw draw;
221: PetscViewerDrawGetDraw(viewer,0,&draw);
222: PetscDrawPopCurrentPoint(draw);
223: }
224: return(0);
225: }
227: /*@C
228: KSPViewFromOptions - View from Options
230: Collective on KSP
232: Input Parameters:
233: + A - Krylov solver context
234: . obj - Optional object
235: - name - command line option
237: Level: intermediate
238: .seealso: KSP, KSPView, PetscObjectViewFromOptions(), KSPCreate()
239: @*/
240: PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[])
241: {
246: PetscObjectViewFromOptions((PetscObject)A,obj,name);
247: return(0);
248: }
250: /*@
251: KSPSetNormType - Sets the norm that is used for convergence testing.
253: Logically Collective on ksp
255: Input Parameter:
256: + ksp - Krylov solver context
257: - normtype - one of
258: $ KSP_NORM_NONE - skips computing the norm, this should generally only be used if you are using
259: $ the Krylov method as a smoother with a fixed small number of iterations.
260: $ Implicitly sets KSPConvergedSkip() as KSP convergence test.
261: $ Note that certain algorithms such as KSPGMRES ALWAYS require the norm calculation,
262: $ for these methods the norms are still computed, they are just not used in
263: $ the convergence test.
264: $ KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
265: $ of the preconditioned residual P^{-1}(b - A x)
266: $ KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
267: $ KSP_NORM_NATURAL - supported by KSPCG, KSPCR, KSPCGNE, KSPCGS
270: Options Database Key:
271: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural>
273: Notes:
274: Not all combinations of preconditioner side (see KSPSetPCSide()) and norm type are supported by all Krylov methods.
275: If only one is set, PETSc tries to automatically change the other to find a compatible pair. If no such combination
276: is supported, PETSc will generate an error.
278: Developer Notes:
279: Supported combinations of norm and preconditioner side are set using KSPSetSupportedNorm().
281: Level: advanced
283: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
284: @*/
285: PetscErrorCode KSPSetNormType(KSP ksp,KSPNormType normtype)
286: {
290: ksp->normtype = ksp->normtype_set = normtype;
291: return(0);
292: }
294: /*@
295: KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be
296: computed and used in the convergence test.
298: Logically Collective on ksp
300: Input Parameter:
301: + ksp - Krylov solver context
302: - it - use -1 to check at all iterations
304: Notes:
305: Currently only works with KSPCG, KSPBCGS and KSPIBCGS
307: Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm
309: On steps where the norm is not computed, the previous norm is still in the variable, so if you run with, for example,
310: -ksp_monitor the residual norm will appear to be unchanged for several iterations (though it is not really unchanged).
311: Level: advanced
313: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType()
314: @*/
315: PetscErrorCode KSPSetCheckNormIteration(KSP ksp,PetscInt it)
316: {
320: ksp->chknorm = it;
321: return(0);
322: }
324: /*@
325: KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the MPI_Allreduce() for
326: computing the inner products for the next iteration. This can reduce communication costs at the expense of doing
327: one additional iteration.
330: Logically Collective on ksp
332: Input Parameter:
333: + ksp - Krylov solver context
334: - flg - PETSC_TRUE or PETSC_FALSE
336: Options Database Keys:
337: . -ksp_lag_norm - lag the calculated residual norm
339: Notes:
340: Currently only works with KSPIBCGS.
342: Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm
344: If you lag the norm and run with, for example, -ksp_monitor, the residual norm reported will be the lagged one.
345: Level: advanced
347: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType(), KSPSetCheckNormIteration()
348: @*/
349: PetscErrorCode KSPSetLagNorm(KSP ksp,PetscBool flg)
350: {
354: ksp->lagnorm = flg;
355: return(0);
356: }
358: /*@
359: KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP
361: Logically Collective
363: Input Arguments:
364: + ksp - Krylov method
365: . normtype - supported norm type
366: . pcside - preconditioner side that can be used with this norm
367: - priority - positive integer preference for this combination; larger values have higher priority
369: Level: developer
371: Notes:
372: This function should be called from the implementation files KSPCreate_XXX() to declare
373: which norms and preconditioner sides are supported. Users should not need to call this
374: function.
376: .seealso: KSPSetNormType(), KSPSetPCSide()
377: @*/
378: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
379: {
383: ksp->normsupporttable[normtype][pcside] = priority;
384: return(0);
385: }
387: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
388: {
392: PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
393: ksp->pc_side = ksp->pc_side_set;
394: ksp->normtype = ksp->normtype_set;
395: return(0);
396: }
398: PetscErrorCode KSPSetUpNorms_Private(KSP ksp,PetscBool errorifnotsupported,KSPNormType *normtype,PCSide *pcside)
399: {
400: PetscInt i,j,best,ibest = 0,jbest = 0;
403: best = 0;
404: for (i=0; i<KSP_NORM_MAX; i++) {
405: for (j=0; j<PC_SIDE_MAX; j++) {
406: if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
407: best = ksp->normsupporttable[i][j];
408: ibest = i;
409: jbest = j;
410: }
411: }
412: }
413: if (best < 1 && errorifnotsupported) {
414: 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);
415: 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]);
416: 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]);
417: 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]);
418: }
419: if (normtype) *normtype = (KSPNormType)ibest;
420: if (pcside) *pcside = (PCSide)jbest;
421: return(0);
422: }
424: /*@
425: KSPGetNormType - Gets the norm that is used for convergence testing.
427: Not Collective
429: Input Parameter:
430: . ksp - Krylov solver context
432: Output Parameter:
433: . normtype - norm that is used for convergence testing
435: Level: advanced
437: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
438: @*/
439: PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype)
440: {
446: KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
447: *normtype = ksp->normtype;
448: return(0);
449: }
451: #if defined(PETSC_HAVE_SAWS)
452: #include <petscviewersaws.h>
453: #endif
455: /*@
456: KSPSetOperators - Sets the matrix associated with the linear system
457: and a (possibly) different one associated with the preconditioner.
459: Collective on ksp
461: Input Parameters:
462: + ksp - the KSP context
463: . Amat - the matrix that defines the linear system
464: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
466: Notes:
468: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
469: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
471: All future calls to KSPSetOperators() must use the same size matrices!
473: Passing a NULL for Amat or Pmat removes the matrix that is currently used.
475: If you wish to replace either Amat or Pmat but leave the other one untouched then
476: first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
477: on it and then pass it back in in your call to KSPSetOperators().
479: Level: beginner
481: Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
482: are created in PC and returned to the user. In this case, if both operators
483: mat and pmat are requested, two DIFFERENT operators will be returned. If
484: only one is requested both operators in the PC will be the same (i.e. as
485: if one had called KSP/PCSetOperators() with the same argument for both Mats).
486: The user must set the sizes of the returned matrices and their type etc just
487: as if the user created them with MatCreate(). For example,
489: $ KSP/PCGetOperators(ksp/pc,&mat,NULL); is equivalent to
490: $ set size, type, etc of mat
492: $ MatCreate(comm,&mat);
493: $ KSP/PCSetOperators(ksp/pc,mat,mat);
494: $ PetscObjectDereference((PetscObject)mat);
495: $ set size, type, etc of mat
497: and
499: $ KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to
500: $ set size, type, etc of mat and pmat
502: $ MatCreate(comm,&mat);
503: $ MatCreate(comm,&pmat);
504: $ KSP/PCSetOperators(ksp/pc,mat,pmat);
505: $ PetscObjectDereference((PetscObject)mat);
506: $ PetscObjectDereference((PetscObject)pmat);
507: $ set size, type, etc of mat and pmat
509: The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
510: of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
511: managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
512: at this is when you create a SNES you do not NEED to create a KSP and attach it to
513: the SNES object (the SNES object manages it for you). Similarly when you create a KSP
514: you do not need to attach a PC to it (the KSP object manages the PC object for you).
515: Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
516: it can be created for you?
518: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
519: @*/
520: PetscErrorCode KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
521: {
530: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
531: PCSetOperators(ksp->pc,Amat,Pmat);
532: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */
533: return(0);
534: }
536: /*@
537: KSPGetOperators - Gets the matrix associated with the linear system
538: and a (possibly) different one associated with the preconditioner.
540: Collective on ksp
542: Input Parameter:
543: . ksp - the KSP context
545: Output Parameters:
546: + Amat - the matrix that defines the linear system
547: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
549: Level: intermediate
551: Notes:
552: DOES NOT increase the reference counts of the matrix, so you should NOT destroy them.
554: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
555: @*/
556: PetscErrorCode KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
557: {
562: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
563: PCGetOperators(ksp->pc,Amat,Pmat);
564: return(0);
565: }
567: /*@C
568: KSPGetOperatorsSet - Determines if the matrix associated with the linear system and
569: possibly a different one associated with the preconditioner have been set in the KSP.
571: Not collective, though the results on all processes should be the same
573: Input Parameter:
574: . pc - the KSP context
576: Output Parameters:
577: + mat - the matrix associated with the linear system was set
578: - pmat - matrix associated with the preconditioner was set, usually the same
580: Level: intermediate
582: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
583: @*/
584: PetscErrorCode KSPGetOperatorsSet(KSP ksp,PetscBool *mat,PetscBool *pmat)
585: {
590: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
591: PCGetOperatorsSet(ksp->pc,mat,pmat);
592: return(0);
593: }
595: /*@C
596: KSPSetPreSolve - Sets a function that is called before every KSPSolve() is started
598: Logically Collective on ksp
600: Input Parameters:
601: + ksp - the solver object
602: . presolve - the function to call before the solve
603: - prectx - any context needed by the function
605: Calling sequence of presolve:
606: $ func(KSP ksp,Vec rhs,Vec x,void *ctx)
608: + ksp - the KSP context
609: . rhs - the right-hand side vector
610: . x - the solution vector
611: - ctx - optional user-provided context
613: Level: developer
615: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
616: @*/
617: PetscErrorCode KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
618: {
621: ksp->presolve = presolve;
622: ksp->prectx = prectx;
623: return(0);
624: }
626: /*@C
627: KSPSetPostSolve - Sets a function that is called after every KSPSolve() completes (whether it converges or not)
629: Logically Collective on ksp
631: Input Parameters:
632: + ksp - the solver object
633: . postsolve - the function to call after the solve
634: - postctx - any context needed by the function
636: Level: developer
638: Calling sequence of postsolve:
639: $ func(KSP ksp,Vec rhs,Vec x,void *ctx)
641: + ksp - the KSP context
642: . rhs - the right-hand side vector
643: . x - the solution vector
644: - ctx - optional user-provided context
646: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPreSolve()
647: @*/
648: PetscErrorCode KSPSetPostSolve(KSP ksp,PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*),void *postctx)
649: {
652: ksp->postsolve = postsolve;
653: ksp->postctx = postctx;
654: return(0);
655: }
657: /*@
658: KSPCreate - Creates the default KSP context.
660: Collective
662: Input Parameter:
663: . comm - MPI communicator
665: Output Parameter:
666: . ksp - location to put the KSP context
668: Notes:
669: The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
670: orthogonalization.
672: Level: beginner
674: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
675: @*/
676: PetscErrorCode KSPCreate(MPI_Comm comm,KSP *inksp)
677: {
678: KSP ksp;
680: void *ctx;
684: *inksp = 0;
685: KSPInitializePackage();
687: PetscHeaderCreate(ksp,KSP_CLASSID,"KSP","Krylov Method","KSP",comm,KSPDestroy,KSPView);
689: ksp->max_it = 10000;
690: ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
691: ksp->rtol = 1.e-5;
692: #if defined(PETSC_USE_REAL_SINGLE)
693: ksp->abstol = 1.e-25;
694: #else
695: ksp->abstol = 1.e-50;
696: #endif
697: ksp->divtol = 1.e4;
699: ksp->chknorm = -1;
700: ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT;
701: ksp->rnorm = 0.0;
702: ksp->its = 0;
703: ksp->guess_zero = PETSC_TRUE;
704: ksp->calc_sings = PETSC_FALSE;
705: ksp->res_hist = NULL;
706: ksp->res_hist_alloc = NULL;
707: ksp->res_hist_len = 0;
708: ksp->res_hist_max = 0;
709: ksp->res_hist_reset = PETSC_TRUE;
710: ksp->numbermonitors = 0;
711: ksp->setfromoptionscalled = 0;
713: KSPConvergedDefaultCreate(&ctx);
714: KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
715: ksp->ops->buildsolution = KSPBuildSolutionDefault;
716: ksp->ops->buildresidual = KSPBuildResidualDefault;
718: ksp->vec_sol = 0;
719: ksp->vec_rhs = 0;
720: ksp->pc = 0;
721: ksp->data = 0;
722: ksp->nwork = 0;
723: ksp->work = 0;
724: ksp->reason = KSP_CONVERGED_ITERATING;
725: ksp->setupstage = KSP_SETUP_NEW;
727: KSPNormSupportTableReset_Private(ksp);
729: *inksp = ksp;
730: return(0);
731: }
733: /*@C
734: KSPSetType - Builds KSP for a particular solver.
736: Logically Collective on ksp
738: Input Parameters:
739: + ksp - the Krylov space context
740: - type - a known method
742: Options Database Key:
743: . -ksp_type <method> - Sets the method; use -help for a list
744: of available methods (for instance, cg or gmres)
746: Notes:
747: See "petsc/include/petscksp.h" for available methods (for instance,
748: KSPCG or KSPGMRES).
750: Normally, it is best to use the KSPSetFromOptions() command and
751: then set the KSP type from the options database rather than by using
752: this routine. Using the options database provides the user with
753: maximum flexibility in evaluating the many different Krylov methods.
754: The KSPSetType() routine is provided for those situations where it
755: is necessary to set the iterative solver independently of the command
756: line or options database. This might be the case, for example, when
757: the choice of iterative solver changes during the execution of the
758: program, and the user's Section 1.5 Writing Application Codes with PETSc is taking responsibility for
759: choosing the appropriate method. In other words, this routine is
760: not for beginners.
762: Level: intermediate
764: Developer Note: KSPRegister() is used to add Krylov types to KSPList from which they
765: are accessed by KSPSetType().
767: .seealso: PCSetType(), KSPType, KSPRegister(), KSPCreate()
769: @*/
770: PetscErrorCode KSPSetType(KSP ksp, KSPType type)
771: {
772: PetscErrorCode ierr,(*r)(KSP);
773: PetscBool match;
779: PetscObjectTypeCompare((PetscObject)ksp,type,&match);
780: if (match) return(0);
782: PetscFunctionListFind(KSPList,type,&r);
783: if (!r) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested KSP type %s",type);
784: /* Destroy the previous private KSP context */
785: if (ksp->ops->destroy) {
786: (*ksp->ops->destroy)(ksp);
787: ksp->ops->destroy = NULL;
788: }
789: /* Reinitialize function pointers in KSPOps structure */
790: PetscMemzero(ksp->ops,sizeof(struct _KSPOps));
791: ksp->ops->buildsolution = KSPBuildSolutionDefault;
792: ksp->ops->buildresidual = KSPBuildResidualDefault;
793: KSPNormSupportTableReset_Private(ksp);
794: /* Call the KSPCreate_XXX routine for this particular Krylov solver */
795: ksp->setupstage = KSP_SETUP_NEW;
796: PetscObjectChangeTypeName((PetscObject)ksp,type);
797: (*r)(ksp);
798: return(0);
799: }
801: /*@C
802: KSPGetType - Gets the KSP type as a string from the KSP object.
804: Not Collective
806: Input Parameter:
807: . ksp - Krylov context
809: Output Parameter:
810: . name - name of KSP method
812: Level: intermediate
814: .seealso: KSPSetType()
815: @*/
816: PetscErrorCode KSPGetType(KSP ksp,KSPType *type)
817: {
821: *type = ((PetscObject)ksp)->type_name;
822: return(0);
823: }
825: /*@C
826: KSPRegister - Adds a method to the Krylov subspace solver package.
828: Not Collective
830: Input Parameters:
831: + name_solver - name of a new user-defined solver
832: - routine_create - routine to create method context
834: Notes:
835: KSPRegister() may be called multiple times to add several user-defined solvers.
837: Sample usage:
838: .vb
839: KSPRegister("my_solver",MySolverCreate);
840: .ve
842: Then, your solver can be chosen with the procedural interface via
843: $ KSPSetType(ksp,"my_solver")
844: or at runtime via the option
845: $ -ksp_type my_solver
847: Level: advanced
849: .seealso: KSPRegisterAll()
850: @*/
851: PetscErrorCode KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
852: {
856: KSPInitializePackage();
857: PetscFunctionListAdd(&KSPList,sname,function);
858: return(0);
859: }