Actual source code: itcreate.c
petsc-3.9.0 2018-04-07
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: PetscInt tabs;
128: PetscViewerASCIIGetTab(viewer, &tabs);
129: PetscViewerASCIISetTab(viewer, ((PetscObject)ksp)->tablevel);
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, nonzero initial guess\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) {
151: PetscViewerASCIIPushTab(viewer);
152: KSPGuessView(ksp->guess,viewer);
153: PetscViewerASCIIPopTab(viewer);
154: }
155: if (ksp->dscale) {PetscViewerASCIIPrintf(viewer," diagonally scaled system\n");}
156: PetscViewerASCIIPrintf(viewer," using %s norm type for convergence test\n",KSPNormTypes[ksp->normtype]);
157: PetscViewerASCIISetTab(viewer, tabs);
158: } else if (isbinary) {
159: PetscInt classid = KSP_FILE_CLASSID;
160: MPI_Comm comm;
161: PetscMPIInt rank;
162: char type[256];
164: PetscObjectGetComm((PetscObject)ksp,&comm);
165: MPI_Comm_rank(comm,&rank);
166: if (!rank) {
167: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
168: PetscStrncpy(type,((PetscObject)ksp)->type_name,256);
169: PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
170: }
171: if (ksp->ops->view) {
172: (*ksp->ops->view)(ksp,viewer);
173: }
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: }
228: /*@
229: KSPSetNormType - Sets the norm that is used for convergence testing.
231: Logically Collective on KSP
233: Input Parameter:
234: + ksp - Krylov solver context
235: - normtype - one of
236: $ KSP_NORM_NONE - skips computing the norm, this should only be used if you are using
237: $ the Krylov method as a smoother with a fixed small number of iterations.
238: $ Implicitly sets KSPConvergedSkip as KSP convergence test.
239: $ KSP_NORM_PRECONDITIONED - the default for left preconditioned solves, uses the l2 norm
240: $ of the preconditioned residual P^{-1}(b - A x)
241: $ KSP_NORM_UNPRECONDITIONED - uses the l2 norm of the true b - Ax residual.
242: $ KSP_NORM_NATURAL - supported by KSPCG, KSPCR, KSPCGNE, KSPCGS
245: Options Database Key:
246: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural>
248: Notes:
249: Not all combinations of preconditioner side (see KSPSetPCSide()) and norm type are supported by all Krylov methods.
250: If only one is set, PETSc tries to automatically change the other to find a compatible pair. If no such combination
251: is supported, PETSc will generate an error.
253: Developer Notes:
254: Supported combinations of norm and preconditioner side are set using KSPSetSupportedNorm().
256: Level: advanced
258: .keywords: KSP, create, context, norms
260: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetCheckNormIteration(), KSPSetPCSide(), KSPGetPCSide(), KSPNormType
261: @*/
262: PetscErrorCode KSPSetNormType(KSP ksp,KSPNormType normtype)
263: {
267: ksp->normtype = ksp->normtype_set = normtype;
268: return(0);
269: }
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: }
303: /*@
304: KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the MPI_Allreduce() for
305: computing the inner products for the next iteration. This can reduce communication costs at the expense of doing
306: one additional iteration.
309: Logically Collective on KSP
311: Input Parameter:
312: + ksp - Krylov solver context
313: - flg - PETSC_TRUE or PETSC_FALSE
315: Options Database Keys:
316: . -ksp_lag_norm - lag the calculated residual norm
318: Notes:
319: Currently only works with KSPIBCGS.
321: Use KSPSetNormType(ksp,KSP_NORM_NONE) to never check the norm
323: If you lag the norm and run with, for example, -ksp_monitor, the residual norm reported will be the lagged one.
324: Level: advanced
326: .keywords: KSP, create, context, norms
328: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPConvergedSkip(), KSPSetNormType(), KSPSetCheckNormIteration()
329: @*/
330: PetscErrorCode KSPSetLagNorm(KSP ksp,PetscBool flg)
331: {
335: ksp->lagnorm = flg;
336: return(0);
337: }
339: /*@
340: KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a KSP
342: Logically Collective
344: Input Arguments:
345: + ksp - Krylov method
346: . normtype - supported norm type
347: . pcside - preconditioner side that can be used with this norm
348: - priority - positive integer preference for this combination; larger values have higher priority
350: Level: developer
352: Notes:
353: This function should be called from the implementation files KSPCreate_XXX() to declare
354: which norms and preconditioner sides are supported. Users should not need to call this
355: function.
357: .seealso: KSPSetNormType(), KSPSetPCSide()
358: @*/
359: PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType normtype,PCSide pcside,PetscInt priority)
360: {
364: ksp->normsupporttable[normtype][pcside] = priority;
365: return(0);
366: }
368: PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp)
369: {
373: PetscMemzero(ksp->normsupporttable,sizeof(ksp->normsupporttable));
374: ksp->pc_side = ksp->pc_side_set;
375: ksp->normtype = ksp->normtype_set;
376: return(0);
377: }
379: PetscErrorCode KSPSetUpNorms_Private(KSP ksp,PetscBool errorifnotsupported,KSPNormType *normtype,PCSide *pcside)
380: {
381: PetscInt i,j,best,ibest = 0,jbest = 0;
384: best = 0;
385: for (i=0; i<KSP_NORM_MAX; i++) {
386: for (j=0; j<PC_SIDE_MAX; j++) {
387: if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && (ksp->normsupporttable[i][j] > best)) {
388: best = ksp->normsupporttable[i][j];
389: ibest = i;
390: jbest = j;
391: }
392: }
393: }
394: if (best < 1 && errorifnotsupported) {
395: 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);
396: 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]);
397: 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]);
398: 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]);
399: }
400: if (normtype) *normtype = (KSPNormType)ibest;
401: if (pcside) *pcside = (PCSide)jbest;
402: return(0);
403: }
405: /*@
406: KSPGetNormType - Gets the norm that is used for convergence testing.
408: Not Collective
410: Input Parameter:
411: . ksp - Krylov solver context
413: Output Parameter:
414: . normtype - norm that is used for convergence testing
416: Level: advanced
418: .keywords: KSP, create, context, norms
420: .seealso: KSPNormType, KSPSetNormType(), KSPConvergedSkip()
421: @*/
422: PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype)
423: {
429: KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
430: *normtype = ksp->normtype;
431: return(0);
432: }
434: #if defined(PETSC_HAVE_SAWS)
435: #include <petscviewersaws.h>
436: #endif
438: /*@
439: KSPSetOperators - Sets the matrix associated with the linear system
440: and a (possibly) different one associated with the preconditioner.
442: Collective on KSP and Mat
444: Input Parameters:
445: + ksp - the KSP context
446: . Amat - the matrix that defines the linear system
447: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
449: Notes:
451: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
452: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
454: All future calls to KSPSetOperators() must use the same size matrices!
456: Passing a NULL for Amat or Pmat removes the matrix that is currently used.
458: If you wish to replace either Amat or Pmat but leave the other one untouched then
459: first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
460: on it and then pass it back in in your call to KSPSetOperators().
462: Level: beginner
464: Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
465: are created in PC and returned to the user. In this case, if both operators
466: mat and pmat are requested, two DIFFERENT operators will be returned. If
467: only one is requested both operators in the PC will be the same (i.e. as
468: if one had called KSP/PCSetOperators() with the same argument for both Mats).
469: The user must set the sizes of the returned matrices and their type etc just
470: as if the user created them with MatCreate(). For example,
472: $ KSP/PCGetOperators(ksp/pc,&mat,NULL); is equivalent to
473: $ set size, type, etc of mat
475: $ MatCreate(comm,&mat);
476: $ KSP/PCSetOperators(ksp/pc,mat,mat);
477: $ PetscObjectDereference((PetscObject)mat);
478: $ set size, type, etc of mat
480: and
482: $ KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to
483: $ set size, type, etc of mat and pmat
485: $ MatCreate(comm,&mat);
486: $ MatCreate(comm,&pmat);
487: $ KSP/PCSetOperators(ksp/pc,mat,pmat);
488: $ PetscObjectDereference((PetscObject)mat);
489: $ PetscObjectDereference((PetscObject)pmat);
490: $ set size, type, etc of mat and pmat
492: The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
493: of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
494: managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
495: at this is when you create a SNES you do not NEED to create a KSP and attach it to
496: the SNES object (the SNES object manages it for you). Similarly when you create a KSP
497: you do not need to attach a PC to it (the KSP object manages the PC object for you).
498: Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
499: it can be created for you?
501: .keywords: KSP, set, operators, matrix, preconditioner, linear system
503: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators(), KSPSetComputeOperators(), KSPSetComputeInitialGuess(), KSPSetComputeRHS()
504: @*/
505: PetscErrorCode KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat)
506: {
515: if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
516: PCSetOperators(ksp->pc,Amat,Pmat);
517: if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */
518: return(0);
519: }
521: /*@
522: KSPGetOperators - Gets the matrix associated with the linear system
523: and a (possibly) different one associated with the preconditioner.
525: Collective on KSP and Mat
527: Input Parameter:
528: . ksp - the KSP context
530: Output Parameters:
531: + Amat - the matrix that defines the linear system
532: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
534: Level: intermediate
536: Notes: 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: Level: developer
595: .keywords: KSP, create, context
597: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPostSolve()
598: @*/
599: PetscErrorCode KSPSetPreSolve(KSP ksp,PetscErrorCode (*presolve)(KSP,Vec,Vec,void*),void *prectx)
600: {
603: ksp->presolve = presolve;
604: ksp->prectx = prectx;
605: return(0);
606: }
608: /*@C
609: KSPSetPostSolve - Sets a function that is called after every KSPSolve() completes (whether it converges or not)
611: Logically Collective on KSP
613: Input Parameters:
614: + ksp - the solver object
615: . postsolve - the function to call after the solve
616: - postctx - any context needed by the function
618: Level: developer
620: .keywords: KSP, create, context
622: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP, KSPSetPreSolve()
623: @*/
624: PetscErrorCode KSPSetPostSolve(KSP ksp,PetscErrorCode (*postsolve)(KSP,Vec,Vec,void*),void *postctx)
625: {
628: ksp->postsolve = postsolve;
629: ksp->postctx = postctx;
630: return(0);
631: }
633: /*@
634: KSPCreate - Creates the default KSP context.
636: Collective on MPI_Comm
638: Input Parameter:
639: . comm - MPI communicator
641: Output Parameter:
642: . ksp - location to put the KSP context
644: Notes:
645: The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
646: orthogonalization.
648: Level: beginner
650: .keywords: KSP, create, context
652: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
653: @*/
654: PetscErrorCode KSPCreate(MPI_Comm comm,KSP *inksp)
655: {
656: KSP ksp;
658: void *ctx;
662: *inksp = 0;
663: KSPInitializePackage();
665: PetscHeaderCreate(ksp,KSP_CLASSID,"KSP","Krylov Method","KSP",comm,KSPDestroy,KSPView);
667: ksp->max_it = 10000;
668: ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT;
669: ksp->rtol = 1.e-5;
670: #if defined(PETSC_USE_REAL_SINGLE)
671: ksp->abstol = 1.e-25;
672: #else
673: ksp->abstol = 1.e-50;
674: #endif
675: ksp->divtol = 1.e4;
677: ksp->chknorm = -1;
678: ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT;
679: ksp->rnorm = 0.0;
680: ksp->its = 0;
681: ksp->guess_zero = PETSC_TRUE;
682: ksp->calc_sings = PETSC_FALSE;
683: ksp->res_hist = NULL;
684: ksp->res_hist_alloc = NULL;
685: ksp->res_hist_len = 0;
686: ksp->res_hist_max = 0;
687: ksp->res_hist_reset = PETSC_TRUE;
688: ksp->numbermonitors = 0;
690: KSPConvergedDefaultCreate(&ctx);
691: KSPSetConvergenceTest(ksp,KSPConvergedDefault,ctx,KSPConvergedDefaultDestroy);
692: ksp->ops->buildsolution = KSPBuildSolutionDefault;
693: ksp->ops->buildresidual = KSPBuildResidualDefault;
695: ksp->vec_sol = 0;
696: ksp->vec_rhs = 0;
697: ksp->pc = 0;
698: ksp->data = 0;
699: ksp->nwork = 0;
700: ksp->work = 0;
701: ksp->reason = KSP_CONVERGED_ITERATING;
702: ksp->setupstage = KSP_SETUP_NEW;
704: KSPNormSupportTableReset_Private(ksp);
706: *inksp = ksp;
707: return(0);
708: }
710: /*@C
711: KSPSetType - Builds KSP for a particular solver.
713: Logically Collective on KSP
715: Input Parameters:
716: + ksp - the Krylov space context
717: - type - a known method
719: Options Database Key:
720: . -ksp_type <method> - Sets the method; use -help for a list
721: of available methods (for instance, cg or gmres)
723: Notes:
724: See "petsc/include/petscksp.h" for available methods (for instance,
725: KSPCG or KSPGMRES).
727: Normally, it is best to use the KSPSetFromOptions() command and
728: then set the KSP type from the options database rather than by using
729: this routine. Using the options database provides the user with
730: maximum flexibility in evaluating the many different Krylov methods.
731: The KSPSetType() routine is provided for those situations where it
732: is necessary to set the iterative solver independently of the command
733: line or options database. This might be the case, for example, when
734: the choice of iterative solver changes during the execution of the
735: program, and the user's application is taking responsibility for
736: choosing the appropriate method. In other words, this routine is
737: not for beginners.
739: Level: intermediate
741: Developer Note: KSPRegister() is used to add Krylov types to KSPList from which they
742: are accessed by KSPSetType().
744: .keywords: KSP, set, method
746: .seealso: PCSetType(), KSPType, KSPRegister(), KSPCreate()
748: @*/
749: PetscErrorCode KSPSetType(KSP ksp, KSPType type)
750: {
751: PetscErrorCode ierr,(*r)(KSP);
752: PetscBool match;
758: PetscObjectTypeCompare((PetscObject)ksp,type,&match);
759: if (match) return(0);
761: PetscFunctionListFind(KSPList,type,&r);
762: if (!r) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested KSP type %s",type);
763: /* Destroy the previous private KSP context */
764: if (ksp->ops->destroy) {
765: (*ksp->ops->destroy)(ksp);
766: ksp->ops->destroy = NULL;
767: }
768: /* Reinitialize function pointers in KSPOps structure */
769: PetscMemzero(ksp->ops,sizeof(struct _KSPOps));
770: ksp->ops->buildsolution = KSPBuildSolutionDefault;
771: ksp->ops->buildresidual = KSPBuildResidualDefault;
772: KSPNormSupportTableReset_Private(ksp);
773: /* Call the KSPCreate_XXX routine for this particular Krylov solver */
774: ksp->setupstage = KSP_SETUP_NEW;
775: PetscObjectChangeTypeName((PetscObject)ksp,type);
776: (*r)(ksp);
777: return(0);
778: }
780: /*@C
781: KSPGetType - Gets the KSP type as a string from the KSP object.
783: Not Collective
785: Input Parameter:
786: . ksp - Krylov context
788: Output Parameter:
789: . name - name of KSP method
791: Level: intermediate
793: .keywords: KSP, get, method, name
795: .seealso: KSPSetType()
796: @*/
797: PetscErrorCode KSPGetType(KSP ksp,KSPType *type)
798: {
802: *type = ((PetscObject)ksp)->type_name;
803: return(0);
804: }
806: /*@C
807: KSPRegister - Adds a method to the Krylov subspace solver package.
809: Not Collective
811: Input Parameters:
812: + name_solver - name of a new user-defined solver
813: - routine_create - routine to create method context
815: Notes:
816: KSPRegister() may be called multiple times to add several user-defined solvers.
818: Sample usage:
819: .vb
820: KSPRegister("my_solver",MySolverCreate);
821: .ve
823: Then, your solver can be chosen with the procedural interface via
824: $ KSPSetType(ksp,"my_solver")
825: or at runtime via the option
826: $ -ksp_type my_solver
828: Level: advanced
830: .keywords: KSP, register
832: .seealso: KSPRegisterAll(), KSPRegisterDestroy()
834: @*/
835: PetscErrorCode KSPRegister(const char sname[],PetscErrorCode (*function)(KSP))
836: {
840: PetscFunctionListAdd(&KSPList,sname,function);
841: return(0);
842: }