Actual source code: itcreate.c
petsc-3.11.4 2019-09-28
2: /*
3: The basic KSP routines, Create, View etc. are here.
4: */
5: #include <petsc/private/kspimpl.h>
7: /* Logging support */
8: PetscClassId KSP_CLASSID;
9: PetscClassId DMKSP_CLASSID;
10: PetscClassId KSPGUESS_CLASSID;
11: PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve;
13: /*
14: Contains the list of registered KSP routines
15: */
16: PetscFunctionList KSPList = 0;
17: PetscBool KSPRegisterAllCalled = PETSC_FALSE;
19: /*@C
20: KSPLoad - Loads a KSP that has been stored in binary with KSPView().
22: Collective on PetscViewer
24: Input Parameters:
25: + newdm - the newly loaded KSP, this needs to have been created with KSPCreate() or
26: some related function before a call to KSPLoad().
27: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
29: Level: intermediate
31: Notes:
32: The type is determined by the data in the file, any type set into the KSP before this call is ignored.
34: Notes for advanced users:
35: Most users should not need to know the details of the binary storage
36: format, since KSPLoad() and KSPView() completely hide these details.
37: But for anyone who's interested, the standard binary matrix storage
38: format is
39: .vb
40: has not yet been determined
41: .ve
43: .seealso: PetscViewerBinaryOpen(), KSPView(), MatLoad(), VecLoad()
44: @*/
45: PetscErrorCode KSPLoad(KSP newdm, PetscViewer viewer)
46: {
48: PetscBool isbinary;
49: PetscInt classid;
50: char type[256];
51: PC pc;
56: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
57: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
59: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
60: if (classid != KSP_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not KSP next in file");
61: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
62: KSPSetType(newdm, type);
63: if (newdm->ops->load) {
64: (*newdm->ops->load)(newdm,viewer);
65: }
66: KSPGetPC(newdm,&pc);
67: PCLoad(pc,viewer);
68: return(0);
69: }
71: #include <petscdraw.h>
72: #if defined(PETSC_HAVE_SAWS)
73: #include <petscviewersaws.h>
74: #endif
75: /*@C
76: KSPView - Prints the KSP data structure.
78: Collective on KSP
80: Input Parameters:
81: + ksp - the Krylov space context
82: - viewer - visualization context
84: Options Database Keys:
85: . -ksp_view - print the ksp data structure at the end of a KSPSolve call
87: Note:
88: The available visualization contexts include
89: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
90: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
91: output where only the first processor opens
92: the file. All other processors send their
93: data to the first processor to print.
95: The user can open an alternative visualization context with
96: PetscViewerASCIIOpen() - output to a specified file.
98: Level: beginner
100: .keywords: KSP, view
102: .seealso: PCView(), PetscViewerASCIIOpen()
103: @*/
104: PetscErrorCode KSPView(KSP ksp,PetscViewer viewer)
105: {
107: PetscBool iascii,isbinary,isdraw;
108: #if defined(PETSC_HAVE_SAWS)
109: PetscBool issaws;
110: #endif
114: if (!viewer) {
115: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
116: }
120: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
121: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
122: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
123: #if defined(PETSC_HAVE_SAWS)
124: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
125: #endif
126: if (iascii) {
127: PetscObjectPrintClassNamePrefixType((PetscObject)ksp,viewer);
128: if (ksp->ops->view) {
129: PetscViewerASCIIPushTab(viewer);
130: (*ksp->ops->view)(ksp,viewer);
131: PetscViewerASCIIPopTab(viewer);
132: }
133: if (ksp->guess_zero) {
134: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, initial guess is zero\n",ksp->max_it);
135: } else {
136: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, nonzero initial guess\n", ksp->max_it);
137: }
138: if (ksp->guess_knoll) {PetscViewerASCIIPrintf(viewer," using preconditioner applied to right hand side for initial guess\n");}
139: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, divergence=%g\n",(double)ksp->rtol,(double)ksp->abstol,(double)ksp->divtol);
140: if (ksp->pc_side == PC_RIGHT) {
141: PetscViewerASCIIPrintf(viewer," right preconditioning\n");
142: } else if (ksp->pc_side == PC_SYMMETRIC) {
143: PetscViewerASCIIPrintf(viewer," symmetric preconditioning\n");
144: } else {
145: PetscViewerASCIIPrintf(viewer," left preconditioning\n");
146: }
147: if (ksp->guess) {
148: PetscViewerASCIIPushTab(viewer);
149: KSPGuessView(ksp->guess,viewer);
150: PetscViewerASCIIPopTab(viewer);
151: }
152: if (ksp->dscale) {PetscViewerASCIIPrintf(viewer," diagonally scaled system\n");}
153: 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: PetscStrncpy(str,"KSP: ",sizeof(str));
181: PetscStrlcat(str,((PetscObject)ksp)->type_name,sizeof(str));
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: }
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 generally 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: $ Note that certain algorithms such as KSPGMRES ALWAYS require the norm calculation,
236: $ for these methods the norms are still computed, they are just not used in
237: $ the convergence test.
238: $ KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
239: $ of the preconditioned residual P^{-1}(b - A x)
240: $ KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
241: $ KSP_NORM_NATURAL - supported by KSPCG, KSPCR, KSPCGNE, KSPCGS
244: Options Database Key:
245: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural>
247: Notes:
248: Not all combinations of preconditioner side (see KSPSetPCSide()) and norm type are supported by all Krylov methods.
249: If only one is set, PETSc tries to automatically change the other to find a compatible pair. If no such combination
250: is supported, PETSc will generate an error.
252: Developer Notes:
253: Supported combinations of norm and preconditioner side are set using KSPSetSupportedNorm().
255: Level: advanced
257: .keywords: KSP, create, context, norms
259: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
260: @*/
261: PetscErrorCode KSPSetNormType(KSP ksp,KSPNormType normtype)
262: {
266: ksp->normtype = ksp->normtype_set = normtype;
267: return(0);
268: }
270: /*@
271: KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be
272: computed and used in the convergence test.
274: Logically Collective on KSP
276: Input Parameter:
277: + ksp - Krylov solver context
278: - it - use -1 to check at all iterations
280: Notes:
281: Currently only works with KSPCG, KSPBCGS and KSPIBCGS
283: Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm
285: On steps where the norm is not computed, the previous norm is still in the variable, so if you run with, for example,
286: -ksp_monitor the residual norm will appear to be unchanged for several iterations (though it is not really unchanged).
287: Level: advanced
289: .keywords: KSP, create, context, norms
291: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType()
292: @*/
293: PetscErrorCode KSPSetCheckNormIteration(KSP ksp,PetscInt it)
294: {
298: ksp->chknorm = it;
299: return(0);
300: }
302: /*@
303: KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the MPI_Allreduce() for
304: computing the inner products for the next iteration. This can reduce communication costs at the expense of doing
305: one additional iteration.
308: Logically Collective on KSP
310: Input Parameter:
311: + ksp - Krylov solver context
312: - flg - PETSC_TRUE or PETSC_FALSE
314: Options Database Keys:
315: . -ksp_lag_norm - lag the calculated residual norm
317: Notes:
318: Currently only works with KSPIBCGS.
320: Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm
322: If you lag the norm and run with, for example, -ksp_monitor, the residual norm reported will be the lagged one.
323: Level: advanced
325: .keywords: KSP, create, context, norms
327: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType(), KSPSetCheckNormIteration()
328: @*/
329: PetscErrorCode KSPSetLagNorm(KSP ksp,PetscBool flg)
330: {
334: ksp->lagnorm = flg;
335: return(0);
336: }
338: /*@
339: KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP
341: Logically Collective
343: Input Arguments:
344: + ksp - Krylov method
345: . normtype - supported norm type
346: . pcside - preconditioner side that can be used with this norm
347: - priority - positive integer preference for this combination; larger values have higher priority
349: Level: developer
351: Notes:
352: This function should be called from the implementation files KSPCreate_XXX() to declare
353: which norms and preconditioner sides are supported. Users should not need to call this
354: function.
356: .seealso: KSPSetNormType(), KSPSetPCSide()
357: @*/
358: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
359: {
363: ksp->normsupporttable[normtype][pcside] = priority;
364: return(0);
365: }
367: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
368: {
372: PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
373: ksp->pc_side = ksp->pc_side_set;
374: ksp->normtype = ksp->normtype_set;
375: return(0);
376: }
378: PetscErrorCode KSPSetUpNorms_Private(KSP ksp,PetscBool errorifnotsupported,KSPNormType *normtype,PCSide *pcside)
379: {
380: PetscInt i,j,best,ibest = 0,jbest = 0;
383: best = 0;
384: for (i=0; i<KSP_NORM_MAX; i++) {
385: for (j=0; j<PC_SIDE_MAX; j++) {
386: if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
387: best = ksp->normsupporttable[i][j];
388: ibest = i;
389: jbest = j;
390: }
391: }
392: }
393: if (best < 1 && errorifnotsupported) {
394: 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);
395: 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]);
396: 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]);
397: 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]);
398: }
399: if (normtype) *normtype = (KSPNormType)ibest;
400: if (pcside) *pcside = (PCSide)jbest;
401: return(0);
402: }
404: /*@
405: KSPGetNormType - Gets the norm that is used for convergence testing.
407: Not Collective
409: Input Parameter:
410: . ksp - Krylov solver context
412: Output Parameter:
413: . normtype - norm that is used for convergence testing
415: Level: advanced
417: .keywords: KSP, create, context, norms
419: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
420: @*/
421: PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype)
422: {
428: KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
429: *normtype = ksp->normtype;
430: return(0);
431: }
433: #if defined(PETSC_HAVE_SAWS)
434: #include <petscviewersaws.h>
435: #endif
437: /*@
438: KSPSetOperators - Sets the matrix associated with the linear system
439: and a (possibly) different one associated with the preconditioner.
441: Collective on KSP and Mat
443: Input Parameters:
444: + ksp - the KSP context
445: . Amat - the matrix that defines the linear system
446: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
448: Notes:
450: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
451: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
453: All future calls to KSPSetOperators() must use the same size matrices!
455: Passing a NULL for Amat or Pmat removes the matrix that is currently used.
457: If you wish to replace either Amat or Pmat but leave the other one untouched then
458: first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
459: on it and then pass it back in in your call to KSPSetOperators().
461: Level: beginner
463: Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
464: are created in PC and returned to the user. In this case, if both operators
465: mat and pmat are requested, two DIFFERENT operators will be returned. If
466: only one is requested both operators in the PC will be the same (i.e. as
467: if one had called KSP/PCSetOperators() with the same argument for both Mats).
468: The user must set the sizes of the returned matrices and their type etc just
469: as if the user created them with MatCreate(). For example,
471: $ KSP/PCGetOperators(ksp/pc,&mat,NULL); is equivalent to
472: $ set size, type, etc of mat
474: $ MatCreate(comm,&mat);
475: $ KSP/PCSetOperators(ksp/pc,mat,mat);
476: $ PetscObjectDereference((PetscObject)mat);
477: $ set size, type, etc of mat
479: and
481: $ KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to
482: $ set size, type, etc of mat and pmat
484: $ MatCreate(comm,&mat);
485: $ MatCreate(comm,&pmat);
486: $ KSP/PCSetOperators(ksp/pc,mat,pmat);
487: $ PetscObjectDereference((PetscObject)mat);
488: $ PetscObjectDereference((PetscObject)pmat);
489: $ set size, type, etc of mat and pmat
491: The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
492: of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
493: managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
494: at this is when you create a SNES you do not NEED to create a KSP and attach it to
495: the SNES object (the SNES object manages it for you). Similarly when you create a KSP
496: you do not need to attach a PC to it (the KSP object manages the PC object for you).
497: Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
498: it can be created for you?
500: .keywords: KSP, set, operators, matrix, preconditioner, linear system
502: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
503: @*/
504: PetscErrorCode KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
505: {
514: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
515: PCSetOperators(ksp->pc,Amat,Pmat);
516: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */
517: return(0);
518: }
520: /*@
521: KSPGetOperators - Gets the matrix associated with the linear system
522: and a (possibly) different one associated with the preconditioner.
524: Collective on KSP and Mat
526: Input Parameter:
527: . ksp - the KSP context
529: Output Parameters:
530: + Amat - the matrix that defines the linear system
531: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
533: Level: intermediate
535: Notes:
536: DOES NOT increase the reference counts of the matrix, so you should NOT destroy them.
538: .keywords: KSP, set, get, operators, matrix, preconditioner, linear system
540: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
541: @*/
542: PetscErrorCode KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
543: {
548: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
549: PCGetOperators(ksp->pc,Amat,Pmat);
550: return(0);
551: }
553: /*@C
554: KSPGetOperatorsSet - Determines if the matrix associated with the linear system and
555: possibly a different one associated with the preconditioner have been set in the KSP.
557: Not collective, though the results on all processes should be the same
559: Input Parameter:
560: . pc - the KSP context
562: Output Parameters:
563: + mat - the matrix associated with the linear system was set
564: - pmat - matrix associated with the preconditioner was set, usually the same
566: Level: intermediate
568: .keywords: KSP, get, operators, matrix, linear system
570: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
571: @*/
572: PetscErrorCode KSPGetOperatorsSet(KSP ksp,PetscBool *mat,PetscBool *pmat)
573: {
578: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
579: PCGetOperatorsSet(ksp->pc,mat,pmat);
580: return(0);
581: }
583: /*@C
584: KSPSetPreSolve - Sets a function that is called before every KSPSolve() is started
586: Logically Collective on KSP
588: Input Parameters:
589: + ksp - the solver object
590: . presolve - the function to call before the solve
591: - prectx - any context needed by the function
593: Calling sequence of presolve:
594: $ func(KSP ksp,Vec rhs,Vec x,void *ctx)
596: + ksp - the KSP context
597: . rhs - the right-hand side vector
598: . x - the solution vector
599: - ctx - optional user-provided context
601: Level: developer
603: .keywords: KSP, create, context
605: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
606: @*/
607: PetscErrorCode KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
608: {
611: ksp->presolve = presolve;
612: ksp->prectx = prectx;
613: return(0);
614: }
616: /*@C
617: KSPSetPostSolve - Sets a function that is called after every KSPSolve() completes (whether it converges or not)
619: Logically Collective on KSP
621: Input Parameters:
622: + ksp - the solver object
623: . postsolve - the function to call after the solve
624: - postctx - any context needed by the function
626: Level: developer
628: Calling sequence of postsolve:
629: $ func(KSP ksp,Vec rhs,Vec x,void *ctx)
631: + ksp - the KSP context
632: . rhs - the right-hand side vector
633: . x - the solution vector
634: - ctx - optional user-provided context
636: .keywords: KSP, create, context
638: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPreSolve()
639: @*/
640: PetscErrorCode KSPSetPostSolve(KSP ksp,PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*),void *postctx)
641: {
644: ksp->postsolve = postsolve;
645: ksp->postctx = postctx;
646: return(0);
647: }
649: /*@
650: KSPCreate - Creates the default KSP context.
652: Collective on MPI_Comm
654: Input Parameter:
655: . comm - MPI communicator
657: Output Parameter:
658: . ksp - location to put the KSP context
660: Notes:
661: The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
662: orthogonalization.
664: Level: beginner
666: .keywords: KSP, create, context
668: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
669: @*/
670: PetscErrorCode KSPCreate(MPI_Comm comm,KSP *inksp)
671: {
672: KSP ksp;
674: void *ctx;
678: *inksp = 0;
679: KSPInitializePackage();
681: PetscHeaderCreate(ksp,KSP_CLASSID,"KSP","Krylov Method","KSP",comm,KSPDestroy,KSPView);
683: ksp->max_it = 10000;
684: ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
685: ksp->rtol = 1.e-5;
686: #if defined(PETSC_USE_REAL_SINGLE)
687: ksp->abstol = 1.e-25;
688: #else
689: ksp->abstol = 1.e-50;
690: #endif
691: ksp->divtol = 1.e4;
693: ksp->chknorm = -1;
694: ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT;
695: ksp->rnorm = 0.0;
696: ksp->its = 0;
697: ksp->guess_zero = PETSC_TRUE;
698: ksp->calc_sings = PETSC_FALSE;
699: ksp->res_hist = NULL;
700: ksp->res_hist_alloc = NULL;
701: ksp->res_hist_len = 0;
702: ksp->res_hist_max = 0;
703: ksp->res_hist_reset = PETSC_TRUE;
704: ksp->numbermonitors = 0;
705: ksp->setfromoptionscalled = 0;
707: KSPConvergedDefaultCreate(&ctx);
708: KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
709: ksp->ops->buildsolution = KSPBuildSolutionDefault;
710: ksp->ops->buildresidual = KSPBuildResidualDefault;
712: ksp->vec_sol = 0;
713: ksp->vec_rhs = 0;
714: ksp->pc = 0;
715: ksp->data = 0;
716: ksp->nwork = 0;
717: ksp->work = 0;
718: ksp->reason = KSP_CONVERGED_ITERATING;
719: ksp->setupstage = KSP_SETUP_NEW;
721: KSPNormSupportTableReset_Private(ksp);
723: *inksp = ksp;
724: return(0);
725: }
727: /*@C
728: KSPSetType - Builds KSP for a particular solver.
730: Logically Collective on KSP
732: Input Parameters:
733: + ksp - the Krylov space context
734: - type - a known method
736: Options Database Key:
737: . -ksp_type <method> - Sets the method; use -help for a list
738: of available methods (for instance, cg or gmres)
740: Notes:
741: See "petsc/include/petscksp.h" for available methods (for instance,
742: KSPCG or KSPGMRES).
744: Normally, it is best to use the KSPSetFromOptions() command and
745: then set the KSP type from the options database rather than by using
746: this routine. Using the options database provides the user with
747: maximum flexibility in evaluating the many different Krylov methods.
748: The KSPSetType() routine is provided for those situations where it
749: is necessary to set the iterative solver independently of the command
750: line or options database. This might be the case, for example, when
751: the choice of iterative solver changes during the execution of the
752: program, and the user's application is taking responsibility for
753: choosing the appropriate method. In other words, this routine is
754: not for beginners.
756: Level: intermediate
758: Developer Note: KSPRegister() is used to add Krylov types to KSPList from which they
759: are accessed by KSPSetType().
761: .keywords: KSP, set, method
763: .seealso: PCSetType(), KSPType, KSPRegister(), KSPCreate()
765: @*/
766: PetscErrorCode KSPSetType(KSP ksp, KSPType type)
767: {
768: PetscErrorCode ierr,(*r)(KSP);
769: PetscBool match;
770: void *ctx;
776: PetscObjectTypeCompare((PetscObject)ksp,type,&match);
777: if (match) return(0);
779: PetscFunctionListFind(KSPList,type,&r);
780: if (!r) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested KSP type %s",type);
781: /* Destroy the previous private KSP context */
782: if (ksp->ops->destroy) {
783: (*ksp->ops->destroy)(ksp);
784: ksp->ops->destroy = NULL;
785: }
786: /* Reinitialize function pointers in KSPOps structure */
787: PetscMemzero(ksp->ops,sizeof(struct _KSPOps));
788: KSPConvergedDefaultCreate(&ctx);
789: KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
790: ksp->ops->buildsolution = KSPBuildSolutionDefault;
791: ksp->ops->buildresidual = KSPBuildResidualDefault;
792: KSPNormSupportTableReset_Private(ksp);
793: /* Call the KSPCreate_XXX routine for this particular Krylov solver */
794: ksp->setupstage = KSP_SETUP_NEW;
795: PetscObjectChangeTypeName((PetscObject)ksp,type);
796: (*r)(ksp);
797: return(0);
798: }
800: /*@C
801: KSPGetType - Gets the KSP type as a string from the KSP object.
803: Not Collective
805: Input Parameter:
806: . ksp - Krylov context
808: Output Parameter:
809: . name - name of KSP method
811: Level: intermediate
813: .keywords: KSP, get, method, name
815: .seealso: KSPSetType()
816: @*/
817: PetscErrorCode KSPGetType(KSP ksp,KSPType *type)
818: {
822: *type = ((PetscObject)ksp)->type_name;
823: return(0);
824: }
826: /*@C
827: KSPRegister - Adds a method to the Krylov subspace solver package.
829: Not Collective
831: Input Parameters:
832: + name_solver - name of a new user-defined solver
833: - routine_create - routine to create method context
835: Notes:
836: KSPRegister() may be called multiple times to add several user-defined solvers.
838: Sample usage:
839: .vb
840: KSPRegister("my_solver",MySolverCreate);
841: .ve
843: Then, your solver can be chosen with the procedural interface via
844: $ KSPSetType(ksp,"my_solver")
845: or at runtime via the option
846: $ -ksp_type my_solver
848: Level: advanced
850: .keywords: KSP, register
852: .seealso: KSPRegisterAll(), KSPRegisterDestroy()
854: @*/
855: PetscErrorCode KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
856: {
860: KSPInitializePackage();
861: PetscFunctionListAdd(&KSPList,sname,function);
862: return(0);
863: }