Actual source code: itcreate.c
1: /*
2: The basic KSP routines, Create, View etc. are here.
3: */
4: #include <petsc/private/kspimpl.h>
6: /* Logging support */
7: PetscClassId KSP_CLASSID;
8: PetscClassId DMKSP_CLASSID;
9: PetscClassId KSPGUESS_CLASSID;
10: PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve, KSP_SolveTranspose, KSP_MatSolve;
12: /*
13: Contains the list of registered KSP routines
14: */
15: PetscFunctionList KSPList = NULL;
16: PetscBool KSPRegisterAllCalled = PETSC_FALSE;
18: /*
19: Contains the list of registered KSP monitors
20: */
21: PetscFunctionList KSPMonitorList = NULL;
22: PetscFunctionList KSPMonitorCreateList = NULL;
23: PetscFunctionList KSPMonitorDestroyList = NULL;
24: PetscBool KSPMonitorRegisterAllCalled = PETSC_FALSE;
26: /*@C
27: KSPLoad - Loads a KSP that has been stored in binary with KSPView().
29: Collective on viewer
31: Input Parameters:
32: + newdm - the newly loaded KSP, this needs to have been created with KSPCreate() or
33: some related function before a call to KSPLoad().
34: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
36: Level: intermediate
38: Notes:
39: The type is determined by the data in the file, any type set into the KSP before this call is ignored.
41: Notes for advanced users:
42: Most users should not need to know the details of the binary storage
43: format, since KSPLoad() and KSPView() completely hide these details.
44: But for anyone who's interested, the standard binary matrix storage
45: format is
46: .vb
47: has not yet been determined
48: .ve
50: .seealso: PetscViewerBinaryOpen(), KSPView(), MatLoad(), VecLoad()
51: @*/
52: PetscErrorCode KSPLoad(KSP newdm, PetscViewer viewer)
53: {
55: PetscBool isbinary;
56: PetscInt classid;
57: char type[256];
58: PC pc;
63: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
64: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
66: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
67: if (classid != KSP_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not KSP next in file");
68: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
69: KSPSetType(newdm, type);
70: if (newdm->ops->load) {
71: (*newdm->ops->load)(newdm,viewer);
72: }
73: KSPGetPC(newdm,&pc);
74: PCLoad(pc,viewer);
75: return(0);
76: }
78: #include <petscdraw.h>
79: #if defined(PETSC_HAVE_SAWS)
80: #include <petscviewersaws.h>
81: #endif
82: /*@C
83: KSPView - Prints the KSP data structure.
85: Collective on ksp
87: Input Parameters:
88: + ksp - the Krylov space context
89: - viewer - visualization context
91: Options Database Keys:
92: . -ksp_view - print the KSP data structure at the end of a KSPSolve call
94: Note:
95: The available visualization contexts include
96: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
97: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
98: output where only the first processor opens
99: the file. All other processors send their
100: data to the first processor to print.
102: The available formats include
103: + PETSC_VIEWER_DEFAULT - standard output (default)
104: - PETSC_VIEWER_ASCII_INFO_DETAIL - more verbose output for PCBJACOBI and PCASM
106: The user can open an alternative visualization context with
107: PetscViewerASCIIOpen() - output to a specified file.
109: In the debugger you can do "call KSPView(ksp,0)" to display the KSP. (The same holds for any PETSc object viewer).
111: Level: beginner
113: .seealso: PCView(), PetscViewerASCIIOpen()
114: @*/
115: PetscErrorCode KSPView(KSP ksp,PetscViewer viewer)
116: {
118: PetscBool iascii,isbinary,isdraw,isstring;
119: #if defined(PETSC_HAVE_SAWS)
120: PetscBool issaws;
121: #endif
125: if (!viewer) {
126: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
127: }
131: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
132: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
133: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
134: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
135: #if defined(PETSC_HAVE_SAWS)
136: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
137: #endif
138: if (iascii) {
139: PetscObjectPrintClassNamePrefixType((PetscObject)ksp,viewer);
140: if (ksp->ops->view) {
141: PetscViewerASCIIPushTab(viewer);
142: (*ksp->ops->view)(ksp,viewer);
143: PetscViewerASCIIPopTab(viewer);
144: }
145: if (ksp->guess_zero) {
146: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, initial guess is zero\n",ksp->max_it);
147: } else {
148: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, nonzero initial guess\n", ksp->max_it);
149: }
150: if (ksp->guess_knoll) {PetscViewerASCIIPrintf(viewer," using preconditioner applied to right hand side for initial guess\n");}
151: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, divergence=%g\n",(double)ksp->rtol,(double)ksp->abstol,(double)ksp->divtol);
152: if (ksp->pc_side == PC_RIGHT) {
153: PetscViewerASCIIPrintf(viewer," right preconditioning\n");
154: } else if (ksp->pc_side == PC_SYMMETRIC) {
155: PetscViewerASCIIPrintf(viewer," symmetric preconditioning\n");
156: } else {
157: PetscViewerASCIIPrintf(viewer," left preconditioning\n");
158: }
159: if (ksp->guess) {
160: PetscViewerASCIIPushTab(viewer);
161: KSPGuessView(ksp->guess,viewer);
162: PetscViewerASCIIPopTab(viewer);
163: }
164: if (ksp->dscale) {PetscViewerASCIIPrintf(viewer," diagonally scaled system\n");}
165: PetscViewerASCIIPrintf(viewer," using %s norm type for convergence test\n",KSPNormTypes[ksp->normtype]);
166: } else if (isbinary) {
167: PetscInt classid = KSP_FILE_CLASSID;
168: MPI_Comm comm;
169: PetscMPIInt rank;
170: char type[256];
172: PetscObjectGetComm((PetscObject)ksp,&comm);
173: MPI_Comm_rank(comm,&rank);
174: if (rank == 0) {
175: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);
176: PetscStrncpy(type,((PetscObject)ksp)->type_name,256);
177: PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR);
178: }
179: if (ksp->ops->view) {
180: (*ksp->ops->view)(ksp,viewer);
181: }
182: } else if (isstring) {
183: const char *type;
184: KSPGetType(ksp,&type);
185: PetscViewerStringSPrintf(viewer," KSPType: %-7.7s",type);
186: if (ksp->ops->view) {(*ksp->ops->view)(ksp,viewer);}
187: } else if (isdraw) {
188: PetscDraw draw;
189: char str[36];
190: PetscReal x,y,bottom,h;
191: PetscBool flg;
193: PetscViewerDrawGetDraw(viewer,0,&draw);
194: PetscDrawGetCurrentPoint(draw,&x,&y);
195: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
196: if (!flg) {
197: PetscStrncpy(str,"KSP: ",sizeof(str));
198: PetscStrlcat(str,((PetscObject)ksp)->type_name,sizeof(str));
199: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
200: bottom = y - h;
201: } else {
202: bottom = y;
203: }
204: PetscDrawPushCurrentPoint(draw,x,bottom);
205: #if defined(PETSC_HAVE_SAWS)
206: } else if (issaws) {
207: PetscMPIInt rank;
208: const char *name;
210: PetscObjectGetName((PetscObject)ksp,&name);
211: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
212: if (!((PetscObject)ksp)->amsmem && rank == 0) {
213: char dir[1024];
215: PetscObjectViewSAWs((PetscObject)ksp,viewer);
216: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
217: PetscStackCallSAWs(SAWs_Register,(dir,&ksp->its,1,SAWs_READ,SAWs_INT));
218: if (!ksp->res_hist) {
219: KSPSetResidualHistory(ksp,NULL,PETSC_DECIDE,PETSC_TRUE);
220: }
221: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/res_hist",name);
222: PetscStackCallSAWs(SAWs_Register,(dir,ksp->res_hist,10,SAWs_READ,SAWs_DOUBLE));
223: }
224: #endif
225: } else if (ksp->ops->view) {
226: (*ksp->ops->view)(ksp,viewer);
227: }
228: if (ksp->pc) {
229: PCView(ksp->pc,viewer);
230: }
231: if (isdraw) {
232: PetscDraw draw;
233: PetscViewerDrawGetDraw(viewer,0,&draw);
234: PetscDrawPopCurrentPoint(draw);
235: }
236: return(0);
237: }
239: /*@C
240: KSPViewFromOptions - View from Options
242: Collective on KSP
244: Input Parameters:
245: + A - Krylov solver context
246: . obj - Optional object
247: - name - command line option
249: Level: intermediate
250: .seealso: KSP, KSPView, PetscObjectViewFromOptions(), KSPCreate()
251: @*/
252: PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[])
253: {
258: PetscObjectViewFromOptions((PetscObject)A,obj,name);
259: return(0);
260: }
262: /*@
263: KSPSetNormType - Sets the norm that is used for convergence testing.
265: Logically Collective on ksp
267: Input Parameters:
268: + ksp - Krylov solver context
269: - normtype - one of
270: $ KSP_NORM_NONE - skips computing the norm, this should generally only be used if you are using
271: $ the Krylov method as a smoother with a fixed small number of iterations.
272: $ Implicitly sets KSPConvergedSkip() as KSP convergence test.
273: $ Note that certain algorithms such as KSPGMRES ALWAYS require the norm calculation,
274: $ for these methods the norms are still computed, they are just not used in
275: $ the convergence test.
276: $ KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
277: $ of the preconditioned residual P^{-1}(b - A x)
278: $ KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
279: $ KSP_NORM_NATURAL - supported by KSPCG, KSPCR, KSPCGNE, KSPCGS
281: Options Database Key:
282: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural>
284: Notes:
285: Not all combinations of preconditioner side (see KSPSetPCSide()) and norm type are supported by all Krylov methods.
286: If only one is set, PETSc tries to automatically change the other to find a compatible pair. If no such combination
287: is supported, PETSc will generate an error.
289: Developer Notes:
290: Supported combinations of norm and preconditioner side are set using KSPSetSupportedNorm().
292: Level: advanced
294: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
295: @*/
296: PetscErrorCode KSPSetNormType(KSP ksp,KSPNormType normtype)
297: {
301: ksp->normtype = ksp->normtype_set = normtype;
302: return(0);
303: }
305: /*@
306: KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be
307: computed and used in the convergence test.
309: Logically Collective on ksp
311: Input Parameters:
312: + ksp - Krylov solver context
313: - it - use -1 to check at all iterations
315: Notes:
316: Currently only works with KSPCG, KSPBCGS and KSPIBCGS
318: Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm
320: On steps where the norm is not computed, the previous norm is still in the variable, so if you run with, for example,
321: -ksp_monitor the residual norm will appear to be unchanged for several iterations (though it is not really unchanged).
322: Level: advanced
324: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType()
325: @*/
326: PetscErrorCode KSPSetCheckNormIteration(KSP ksp,PetscInt it)
327: {
331: ksp->chknorm = it;
332: return(0);
333: }
335: /*@
336: KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the MPI_Allreduce() for
337: computing the inner products for the next iteration. This can reduce communication costs at the expense of doing
338: one additional iteration.
340: Logically Collective on ksp
342: Input Parameters:
343: + ksp - Krylov solver context
344: - flg - PETSC_TRUE or PETSC_FALSE
346: Options Database Keys:
347: . -ksp_lag_norm - lag the calculated residual norm
349: Notes:
350: Currently only works with KSPIBCGS.
352: Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm
354: If you lag the norm and run with, for example, -ksp_monitor, the residual norm reported will be the lagged one.
355: Level: advanced
357: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType(), KSPSetCheckNormIteration()
358: @*/
359: PetscErrorCode KSPSetLagNorm(KSP ksp,PetscBool flg)
360: {
364: ksp->lagnorm = flg;
365: return(0);
366: }
368: /*@
369: KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP
371: Logically Collective
373: Input Parameters:
374: + ksp - Krylov method
375: . normtype - supported norm type
376: . pcside - preconditioner side that can be used with this norm
377: - priority - positive integer preference for this combination; larger values have higher priority
379: Level: developer
381: Notes:
382: This function should be called from the implementation files KSPCreate_XXX() to declare
383: which norms and preconditioner sides are supported. Users should not need to call this
384: function.
386: .seealso: KSPSetNormType(), KSPSetPCSide()
387: @*/
388: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
389: {
393: ksp->normsupporttable[normtype][pcside] = priority;
394: return(0);
395: }
397: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
398: {
402: PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
403: ksp->pc_side = ksp->pc_side_set;
404: ksp->normtype = ksp->normtype_set;
405: return(0);
406: }
408: PetscErrorCode KSPSetUpNorms_Private(KSP ksp,PetscBool errorifnotsupported,KSPNormType *normtype,PCSide *pcside)
409: {
410: PetscInt i,j,best,ibest = 0,jbest = 0;
413: best = 0;
414: for (i=0; i<KSP_NORM_MAX; i++) {
415: for (j=0; j<PC_SIDE_MAX; j++) {
416: if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
417: best = ksp->normsupporttable[i][j];
418: ibest = i;
419: jbest = j;
420: }
421: }
422: }
423: if (best < 1 && errorifnotsupported) {
424: 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);
425: 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]);
426: 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]);
427: 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]);
428: }
429: if (normtype) *normtype = (KSPNormType)ibest;
430: if (pcside) *pcside = (PCSide)jbest;
431: return(0);
432: }
434: /*@
435: KSPGetNormType - Gets the norm that is used for convergence testing.
437: Not Collective
439: Input Parameter:
440: . ksp - Krylov solver context
442: Output Parameter:
443: . normtype - norm that is used for convergence testing
445: Level: advanced
447: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
448: @*/
449: PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype)
450: {
456: KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
457: *normtype = ksp->normtype;
458: return(0);
459: }
461: #if defined(PETSC_HAVE_SAWS)
462: #include <petscviewersaws.h>
463: #endif
465: /*@
466: KSPSetOperators - Sets the matrix associated with the linear system
467: and a (possibly) different one associated with the preconditioner.
469: Collective on ksp
471: Input Parameters:
472: + ksp - the KSP context
473: . Amat - the matrix that defines the linear system
474: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
476: Notes:
478: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
479: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
481: All future calls to KSPSetOperators() must use the same size matrices!
483: Passing a NULL for Amat or Pmat removes the matrix that is currently used.
485: If you wish to replace either Amat or Pmat but leave the other one untouched then
486: first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
487: on it and then pass it back in in your call to KSPSetOperators().
489: Level: beginner
491: Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
492: are created in PC and returned to the user. In this case, if both operators
493: mat and pmat are requested, two DIFFERENT operators will be returned. If
494: only one is requested both operators in the PC will be the same (i.e. as
495: if one had called KSP/PCSetOperators() with the same argument for both Mats).
496: The user must set the sizes of the returned matrices and their type etc just
497: as if the user created them with MatCreate(). For example,
499: $ KSP/PCGetOperators(ksp/pc,&mat,NULL); is equivalent to
500: $ set size, type, etc of mat
502: $ MatCreate(comm,&mat);
503: $ KSP/PCSetOperators(ksp/pc,mat,mat);
504: $ PetscObjectDereference((PetscObject)mat);
505: $ set size, type, etc of mat
507: and
509: $ KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to
510: $ set size, type, etc of mat and pmat
512: $ MatCreate(comm,&mat);
513: $ MatCreate(comm,&pmat);
514: $ KSP/PCSetOperators(ksp/pc,mat,pmat);
515: $ PetscObjectDereference((PetscObject)mat);
516: $ PetscObjectDereference((PetscObject)pmat);
517: $ set size, type, etc of mat and pmat
519: The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
520: of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
521: managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
522: at this is when you create a SNES you do not NEED to create a KSP and attach it to
523: the SNES object (the SNES object manages it for you). Similarly when you create a KSP
524: you do not need to attach a PC to it (the KSP object manages the PC object for you).
525: Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
526: it can be created for you?
528: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
529: @*/
530: PetscErrorCode KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
531: {
540: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
541: PCSetOperators(ksp->pc,Amat,Pmat);
542: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */
543: return(0);
544: }
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
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:
562: DOES NOT increase the reference counts of the matrix, so you should NOT destroy them.
564: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators(), KSPGetOperatorsSet()
565: @*/
566: PetscErrorCode KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat)
567: {
572: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
573: PCGetOperators(ksp->pc,Amat,Pmat);
574: return(0);
575: }
577: /*@C
578: KSPGetOperatorsSet - Determines if the matrix associated with the linear system and
579: possibly a different one associated with the preconditioner have been set in the KSP.
581: Not collective, though the results on all processes should be the same
583: Input Parameter:
584: . pc - the KSP context
586: Output Parameters:
587: + mat - the matrix associated with the linear system was set
588: - pmat - matrix associated with the preconditioner was set, usually the same
590: Level: intermediate
592: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators(), PCGetOperatorsSet()
593: @*/
594: PetscErrorCode KSPGetOperatorsSet(KSP ksp,PetscBool *mat,PetscBool *pmat)
595: {
600: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
601: PCGetOperatorsSet(ksp->pc,mat,pmat);
602: return(0);
603: }
605: /*@C
606: KSPSetPreSolve - Sets a function that is called before every KSPSolve() is started
608: Logically Collective on ksp
610: Input Parameters:
611: + ksp - the solver object
612: . presolve - the function to call before the solve
613: - prectx - any context needed by the function
615: Calling sequence of presolve:
616: $ func(KSP ksp,Vec rhs,Vec x,void *ctx)
618: + ksp - the KSP context
619: . rhs - the right-hand side vector
620: . x - the solution vector
621: - ctx - optional user-provided context
623: Level: developer
625: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
626: @*/
627: PetscErrorCode KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
628: {
631: ksp->presolve = presolve;
632: ksp->prectx = prectx;
633: return(0);
634: }
636: /*@C
637: KSPSetPostSolve - Sets a function that is called after every KSPSolve() completes (whether it converges or not)
639: Logically Collective on ksp
641: Input Parameters:
642: + ksp - the solver object
643: . postsolve - the function to call after the solve
644: - postctx - any context needed by the function
646: Level: developer
648: Calling sequence of postsolve:
649: $ func(KSP ksp,Vec rhs,Vec x,void *ctx)
651: + ksp - the KSP context
652: . rhs - the right-hand side vector
653: . x - the solution vector
654: - ctx - optional user-provided context
656: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPreSolve()
657: @*/
658: PetscErrorCode KSPSetPostSolve(KSP ksp,PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*),void *postctx)
659: {
662: ksp->postsolve = postsolve;
663: ksp->postctx = postctx;
664: return(0);
665: }
667: /*@
668: KSPCreate - Creates the default KSP context.
670: Collective
672: Input Parameter:
673: . comm - MPI communicator
675: Output Parameter:
676: . ksp - location to put the KSP context
678: Notes:
679: The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
680: orthogonalization.
682: Level: beginner
684: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
685: @*/
686: PetscErrorCode KSPCreate(MPI_Comm comm,KSP *inksp)
687: {
688: KSP ksp;
690: void *ctx;
694: *inksp = NULL;
695: KSPInitializePackage();
697: PetscHeaderCreate(ksp,KSP_CLASSID,"KSP","Krylov Method","KSP",comm,KSPDestroy,KSPView);
699: ksp->max_it = 10000;
700: ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
701: ksp->rtol = 1.e-5;
702: #if defined(PETSC_USE_REAL_SINGLE)
703: ksp->abstol = 1.e-25;
704: #else
705: ksp->abstol = 1.e-50;
706: #endif
707: ksp->divtol = 1.e4;
709: ksp->chknorm = -1;
710: ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT;
711: ksp->rnorm = 0.0;
712: ksp->its = 0;
713: ksp->guess_zero = PETSC_TRUE;
714: ksp->calc_sings = PETSC_FALSE;
715: ksp->res_hist = NULL;
716: ksp->res_hist_alloc = NULL;
717: ksp->res_hist_len = 0;
718: ksp->res_hist_max = 0;
719: ksp->res_hist_reset = PETSC_TRUE;
720: ksp->err_hist = NULL;
721: ksp->err_hist_alloc = NULL;
722: ksp->err_hist_len = 0;
723: ksp->err_hist_max = 0;
724: ksp->err_hist_reset = PETSC_TRUE;
725: ksp->numbermonitors = 0;
726: ksp->numberreasonviews = 0;
727: ksp->setfromoptionscalled = 0;
728: ksp->nmax = PETSC_DECIDE;
730: KSPConvergedDefaultCreate(&ctx);
731: KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
732: ksp->ops->buildsolution = KSPBuildSolutionDefault;
733: ksp->ops->buildresidual = KSPBuildResidualDefault;
735: ksp->vec_sol = NULL;
736: ksp->vec_rhs = NULL;
737: ksp->pc = NULL;
738: ksp->data = NULL;
739: ksp->nwork = 0;
740: ksp->work = NULL;
741: ksp->reason = KSP_CONVERGED_ITERATING;
742: ksp->setupstage = KSP_SETUP_NEW;
744: KSPNormSupportTableReset_Private(ksp);
746: *inksp = ksp;
747: return(0);
748: }
750: /*@C
751: KSPSetType - Builds KSP for a particular solver.
753: Logically Collective on ksp
755: Input Parameters:
756: + ksp - the Krylov space context
757: - type - a known method
759: Options Database Key:
760: . -ksp_type <method> - Sets the method; use -help for a list
761: of available methods (for instance, cg or gmres)
763: Notes:
764: See "petsc/include/petscksp.h" for available methods (for instance,
765: KSPCG or KSPGMRES).
767: Normally, it is best to use the KSPSetFromOptions() command and
768: then set the KSP type from the options database rather than by using
769: this routine. Using the options database provides the user with
770: maximum flexibility in evaluating the many different Krylov methods.
771: The KSPSetType() routine is provided for those situations where it
772: is necessary to set the iterative solver independently of the command
773: line or options database. This might be the case, for example, when
774: the choice of iterative solver changes during the execution of the
775: program, and the user's application is taking responsibility for
776: choosing the appropriate method. In other words, this routine is
777: not for beginners.
779: Level: intermediate
781: Developer Note: KSPRegister() is used to add Krylov types to KSPList from which they
782: are accessed by KSPSetType().
784: .seealso: PCSetType(), KSPType, KSPRegister(), KSPCreate()
786: @*/
787: PetscErrorCode KSPSetType(KSP ksp, KSPType type)
788: {
789: PetscErrorCode ierr,(*r)(KSP);
790: PetscBool match;
796: PetscObjectTypeCompare((PetscObject)ksp,type,&match);
797: if (match) return(0);
799: PetscFunctionListFind(KSPList,type,&r);
800: if (!r) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested KSP type %s",type);
801: /* Destroy the previous private KSP context */
802: if (ksp->ops->destroy) {
803: (*ksp->ops->destroy)(ksp);
804: ksp->ops->destroy = NULL;
805: }
806: /* Reinitialize function pointers in KSPOps structure */
807: PetscMemzero(ksp->ops,sizeof(struct _KSPOps));
808: ksp->ops->buildsolution = KSPBuildSolutionDefault;
809: ksp->ops->buildresidual = KSPBuildResidualDefault;
810: KSPNormSupportTableReset_Private(ksp);
811: ksp->setupnewmatrix = PETSC_FALSE; // restore default (setup not called in case of new matrix)
812: /* Call the KSPCreate_XXX routine for this particular Krylov solver */
813: ksp->setupstage = KSP_SETUP_NEW;
814: (*r)(ksp);
815: PetscObjectChangeTypeName((PetscObject)ksp,type);
816: return(0);
817: }
819: /*@C
820: KSPGetType - Gets the KSP type as a string from the KSP object.
822: Not Collective
824: Input Parameter:
825: . ksp - Krylov context
827: Output Parameter:
828: . name - name of KSP method
830: Level: intermediate
832: .seealso: KSPSetType()
833: @*/
834: PetscErrorCode KSPGetType(KSP ksp,KSPType *type)
835: {
839: *type = ((PetscObject)ksp)->type_name;
840: return(0);
841: }
843: /*@C
844: KSPRegister - Adds a method to the Krylov subspace solver package.
846: Not Collective
848: Input Parameters:
849: + name_solver - name of a new user-defined solver
850: - routine_create - routine to create method context
852: Notes:
853: KSPRegister() may be called multiple times to add several user-defined solvers.
855: Sample usage:
856: .vb
857: KSPRegister("my_solver",MySolverCreate);
858: .ve
860: Then, your solver can be chosen with the procedural interface via
861: $ KSPSetType(ksp,"my_solver")
862: or at runtime via the option
863: $ -ksp_type my_solver
865: Level: advanced
867: .seealso: KSPRegisterAll()
868: @*/
869: PetscErrorCode KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
870: {
874: KSPInitializePackage();
875: PetscFunctionListAdd(&KSPList,sname,function);
876: return(0);
877: }
879: PetscErrorCode KSPMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[])
880: {
884: PetscStrncpy(key, name, PETSC_MAX_PATH_LEN);
885: PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN);
886: PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN);
887: PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN);
888: PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN);
889: return(0);
890: }
892: /*@C
893: KSPMonitorRegister - Adds Krylov subspace solver monitor routine.
895: Not Collective
897: Input Parameters:
898: + name - name of a new monitor routine
899: . vtype - A PetscViewerType for the output
900: . format - A PetscViewerFormat for the output
901: . monitor - Monitor routine
902: . create - Creation routine, or NULL
903: - destroy - Destruction routine, or NULL
905: Notes:
906: KSPMonitorRegister() may be called multiple times to add several user-defined monitors.
908: Sample usage:
909: .vb
910: KSPMonitorRegister("my_monitor",PETSCVIEWERASCII,PETSC_VIEWER_ASCII_INFO_DETAIL,MyMonitor,NULL,NULL);
911: .ve
913: Then, your monitor can be chosen with the procedural interface via
914: $ KSPMonitorSetFromOptions(ksp,"-ksp_monitor_my_monitor","my_monitor",NULL)
915: or at runtime via the option
916: $ -ksp_monitor_my_monitor
918: Level: advanced
920: .seealso: KSPMonitorRegisterAll()
921: @*/
922: PetscErrorCode KSPMonitorRegister(const char name[], PetscViewerType vtype, PetscViewerFormat format,
923: PetscErrorCode (*monitor)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *),
924: PetscErrorCode (*create)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **),
925: PetscErrorCode (*destroy)(PetscViewerAndFormat **))
926: {
927: char key[PETSC_MAX_PATH_LEN];
931: KSPInitializePackage();
932: KSPMonitorMakeKey_Internal(name, vtype, format, key);
933: PetscFunctionListAdd(&KSPMonitorList, key, monitor);
934: if (create) {PetscFunctionListAdd(&KSPMonitorCreateList, key, create);}
935: if (destroy) {PetscFunctionListAdd(&KSPMonitorDestroyList, key, destroy);}
936: return(0);
937: }