Actual source code: itcreate.c
1: #define PETSCKSP_DLL
3: /*
4: The basic KSP routines, Create, View etc. are here.
5: */
6: #include src/ksp/ksp/kspimpl.h
7: #include petscsys.h
9: /* Logging support */
10: PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE = 0;
11: PetscEvent KSP_GMRESOrthogonalization = 0, KSP_SetUp = 0, KSP_Solve = 0;
14: PetscTruth KSPRegisterAllCalled = PETSC_FALSE;
18: /*@C
19: KSPView - Prints the KSP data structure.
21: Collective on KSP
23: Input Parameters:
24: + ksp - the Krylov space context
25: - viewer - visualization context
27: Options Database Keys:
28: . -ksp_view - print the ksp data structure at the end of a KSPSolve call
30: Note:
31: The available visualization contexts include
32: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
33: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
34: output where only the first processor opens
35: the file. All other processors send their
36: data to the first processor to print.
38: The user can open an alternative visualization context with
39: PetscViewerASCIIOpen() - output to a specified file.
41: Level: beginner
43: .keywords: KSP, view
45: .seealso: PCView(), PetscViewerASCIIOpen()
46: @*/
47: PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP ksp,PetscViewer viewer)
48: {
49: const char *type;
51: PetscTruth iascii;
55: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ksp->comm);
59: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
60: if (iascii) {
61: KSPGetType(ksp,&type);
62: if (ksp->prefix) {
63: PetscViewerASCIIPrintf(viewer,"KSP Object:(%s)\n",ksp->prefix);
64: } else {
65: PetscViewerASCIIPrintf(viewer,"KSP Object:\n");
66: }
67: if (type) {
68: PetscViewerASCIIPrintf(viewer," type: %s\n",type);
69: } else {
70: PetscViewerASCIIPrintf(viewer," type: not yet set\n");
71: }
72: if (ksp->ops->view) {
73: PetscViewerASCIIPushTab(viewer);
74: (*ksp->ops->view)(ksp,viewer);
75: PetscViewerASCIIPopTab(viewer);
76: }
77: if (ksp->guess_zero) {PetscViewerASCIIPrintf(viewer," maximum iterations=%D, initial guess is zero\n",ksp->max_it);}
78: else {PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", ksp->max_it);}
79: if (ksp->guess_knoll) {PetscViewerASCIIPrintf(viewer," using preconditioner applied to right hand side for initial guess\n");}
80: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, divergence=%g\n",ksp->rtol,ksp->abstol,ksp->divtol);
81: if (ksp->pc_side == PC_RIGHT) {PetscViewerASCIIPrintf(viewer," right preconditioning\n");}
82: else if (ksp->pc_side == PC_SYMMETRIC) {PetscViewerASCIIPrintf(viewer," symmetric preconditioning\n");}
83: else {PetscViewerASCIIPrintf(viewer," left preconditioning\n");}
84: } else {
85: if (ksp->ops->view) {
86: (*ksp->ops->view)(ksp,viewer);
87: }
88: }
89: PCView(ksp->pc,viewer);
90: return(0);
91: }
93: /*
94: Contains the list of registered KSP routines
95: */
96: PetscFList KSPList = 0;
100: /*@C
101: KSPSetNormType - Sets the norm that is used for convergence testing.
103: Collective on KSP
105: Input Parameter:
106: + ksp - Krylov solver context
107: - normtype - one of
108: $ KSP_NO_NORM - skips computing the norm, this should only be used if you are using
109: $ the Krylov method as a smoother with a fixed small number of iterations.
110: $ You must also call KSPSetConvergenceTest(ksp,KSPSkipConverged,PETSC_NULL);
111: $ supported only by CG, Richardson, Bi-CG-stab, CR, and CGS methods.
112: $ KSP_PRECONDITIONED_NORM - the default for left preconditioned solves, uses the l2 norm
113: $ of the preconditioned residual
114: $ KSP_UNPRECONDITIONED_NORM - uses the l2 norm of the true b - Ax residual, supported only by
115: $ CG, CHEBYCHEV, and RICHARDSON
116: $ KSP_NATURAL_NORM - supported by cg, cr, and cgs
119: Options Database Key:
120: . -ksp_norm_type <none,preconditioned,unpreconditioned,natural>
122: Notes:
123: Currently only works with the CG, Richardson, Bi-CG-stab, CR, and CGS methods.
125: Level: advanced
127: .keywords: KSP, create, context, norms
129: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPSkipConverged()
130: @*/
131: PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP ksp,KSPNormType normtype)
132: {
137: ksp->normtype = normtype;
138: if (normtype == KSP_NO_NORM) {
139: PetscLogInfo((ksp,"KSPSetNormType:Warning seting KSPNormType to skip computing the norm\n\
140: make sure you set the KSP convergence test to KSPSkipConvergence\n"));
141: }
142: return(0);
143: }
147: static PetscErrorCode KSPPublish_Petsc(PetscObject obj)
148: {
150: return(0);
151: }
155: /*@
156: KSPSetOperators - Sets the matrix associated with the linear system
157: and a (possibly) different one associated with the preconditioner.
159: Collective on KSP and Mat
161: Input Parameters:
162: + ksp - the KSP context
163: . Amat - the matrix associated with the linear system
164: . Pmat - the matrix to be used in constructing the preconditioner, usually the
165: same as Amat.
166: - flag - flag indicating information about the preconditioner matrix structure
167: during successive linear solves. This flag is ignored the first time a
168: linear system is solved, and thus is irrelevant when solving just one linear
169: system.
171: Notes:
172: The flag can be used to eliminate unnecessary work in the preconditioner
173: during the repeated solution of linear systems of the same size. The
174: available options are
175: $ SAME_PRECONDITIONER -
176: $ Pmat is identical during successive linear solves.
177: $ This option is intended for folks who are using
178: $ different Amat and Pmat matrices and want to reuse the
179: $ same preconditioner matrix. For example, this option
180: $ saves work by not recomputing incomplete factorization
181: $ for ILU/ICC preconditioners.
182: $ SAME_NONZERO_PATTERN -
183: $ Pmat has the same nonzero structure during
184: $ successive linear solves.
185: $ DIFFERENT_NONZERO_PATTERN -
186: $ Pmat does not have the same nonzero structure.
188: Caution:
189: If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
190: and does not check the structure of the matrix. If you erroneously
191: claim that the structure is the same when it actually is not, the new
192: preconditioner will not function correctly. Thus, use this optimization
193: feature carefully!
195: If in doubt about whether your preconditioner matrix has changed
196: structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
198: Level: beginner
200: .keywords: KSP, set, operators, matrix, preconditioner, linear system
202: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators()
203: @*/
204: PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat,MatStructure flag)
205: {
212: PCSetOperators(ksp->pc,Amat,Pmat,flag);
213: if (ksp->setupcalled > 1) ksp->setupcalled = 1; /* so that next solve call will call setup */
214: return(0);
215: }
219: /*@
220: KSPGetOperators - Gets the matrix associated with the linear system
221: and a (possibly) different one associated with the preconditioner.
223: Collective on KSP and Mat
225: Input Parameter:
226: . ksp - the KSP context
228: Output Parameters:
229: + Amat - the matrix associated with the linear system
230: . Pmat - the matrix to be used in constructing the preconditioner, usually the
231: same as Amat.
232: - flag - flag indicating information about the preconditioner matrix structure
233: during successive linear solves. This flag is ignored the first time a
234: linear system is solved, and thus is irrelevant when solving just one linear
235: system.
237: Level: intermediate
239: .keywords: KSP, set, get, operators, matrix, preconditioner, linear system
241: .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators()
242: @*/
243: PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat,MatStructure *flag)
244: {
249: PCGetOperators(ksp->pc,Amat,Pmat,flag);
250: return(0);
251: }
255: /*@C
256: KSPCreate - Creates the default KSP context.
258: Collective on MPI_Comm
260: Input Parameter:
261: . comm - MPI communicator
263: Output Parameter:
264: . ksp - location to put the KSP context
266: Notes:
267: The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt
268: orthogonalization.
270: Level: beginner
272: .keywords: KSP, create, context
274: .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP
275: @*/
276: PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm comm,KSP *inksp)
277: {
278: KSP ksp;
283: *inksp = 0;
284: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
285: KSPInitializePackage(PETSC_NULL);
286: #endif
288: PetscHeaderCreate(ksp,_p_KSP,struct _KSPOps,KSP_COOKIE,-1,"KSP",comm,KSPDestroy,KSPView);
289: *inksp = ksp;
290: ksp->bops->publish = KSPPublish_Petsc;
292: ksp->type = -1;
293: ksp->max_it = 10000;
294: ksp->pc_side = PC_LEFT;
295: ksp->rtol = 1.e-5;
296: ksp->abstol = 1.e-50;
297: ksp->divtol = 1.e4;
299: ksp->normtype = KSP_PRECONDITIONED_NORM;
300: ksp->rnorm = 0.0;
301: ksp->its = 0;
302: ksp->guess_zero = PETSC_TRUE;
303: ksp->calc_sings = PETSC_FALSE;
304: ksp->res_hist = PETSC_NULL;
305: ksp->res_hist_len = 0;
306: ksp->res_hist_max = 0;
307: ksp->res_hist_reset = PETSC_TRUE;
308: ksp->numbermonitors = 0;
309: ksp->converged = KSPDefaultConverged;
310: ksp->ops->buildsolution = KSPDefaultBuildSolution;
311: ksp->ops->buildresidual = KSPDefaultBuildResidual;
313: ksp->ops->setfromoptions = 0;
315: ksp->vec_sol = 0;
316: ksp->vec_rhs = 0;
317: ksp->pc = 0;
319: ksp->ops->solve = 0;
320: ksp->ops->setup = 0;
321: ksp->ops->destroy = 0;
323: ksp->data = 0;
324: ksp->nwork = 0;
325: ksp->work = 0;
327: ksp->cnvP = 0;
329: ksp->reason = KSP_CONVERGED_ITERATING;
331: ksp->setupcalled = 0;
332: PetscPublishAll(ksp);
333: PCCreate(comm,&ksp->pc);
334: return(0);
335: }
336:
339: /*@C
340: KSPSetType - Builds KSP for a particular solver.
342: Collective on KSP
344: Input Parameters:
345: + ksp - the Krylov space context
346: - type - a known method
348: Options Database Key:
349: . -ksp_type <method> - Sets the method; use -help for a list
350: of available methods (for instance, cg or gmres)
352: Notes:
353: See "petsc/include/petscksp.h" for available methods (for instance,
354: KSPCG or KSPGMRES).
356: Normally, it is best to use the KSPSetFromOptions() command and
357: then set the KSP type from the options database rather than by using
358: this routine. Using the options database provides the user with
359: maximum flexibility in evaluating the many different Krylov methods.
360: The KSPSetType() routine is provided for those situations where it
361: is necessary to set the iterative solver independently of the command
362: line or options database. This might be the case, for example, when
363: the choice of iterative solver changes during the execution of the
364: program, and the user's application is taking responsibility for
365: choosing the appropriate method. In other words, this routine is
366: not for beginners.
368: Level: intermediate
370: .keywords: KSP, set, method
372: .seealso: PCSetType(), KSPType
374: @*/
375: PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP ksp, KSPType type)
376: {
377: PetscErrorCode ierr,(*r)(KSP);
378: PetscTruth match;
384: PetscTypeCompare((PetscObject)ksp,type,&match);
385: if (match) return(0);
387: if (ksp->data) {
388: /* destroy the old private KSP context */
389: (*ksp->ops->destroy)(ksp);
390: ksp->data = 0;
391: }
392: /* Get the function pointers for the iterative method requested */
393: if (!KSPRegisterAllCalled) {KSPRegisterAll(PETSC_NULL);}
394: PetscFListFind(ksp->comm,KSPList,type,(void (**)(void)) &r);
395: if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown KSP type given: %s",type);
396: ksp->setupcalled = 0;
397: (*r)(ksp);
398: PetscObjectChangeTypeName((PetscObject)ksp,type);
399: return(0);
400: }
404: /*@C
405: KSPRegisterDestroy - Frees the list of KSP methods that were
406: registered by KSPRegisterDynamic().
408: Not Collective
410: Level: advanced
412: .keywords: KSP, register, destroy
414: .seealso: KSPRegisterDynamic(), KSPRegisterAll()
415: @*/
416: PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void)
417: {
421: if (KSPList) {
422: PetscFListDestroy(&KSPList);
423: KSPList = 0;
424: }
425: KSPRegisterAllCalled = PETSC_FALSE;
426: return(0);
427: }
431: /*@C
432: KSPGetType - Gets the KSP type as a string from the KSP object.
434: Not Collective
436: Input Parameter:
437: . ksp - Krylov context
439: Output Parameter:
440: . name - name of KSP method
442: Level: intermediate
444: .keywords: KSP, get, method, name
446: .seealso: KSPSetType()
447: @*/
448: PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP ksp,KSPType *type)
449: {
453: *type = ksp->type_name;
454: return(0);
455: }
459: /*@
460: KSPSetFromOptions - Sets KSP options from the options database.
461: This routine must be called before KSPSetUp() if the user is to be
462: allowed to set the Krylov type.
464: Collective on KSP
466: Input Parameters:
467: . ksp - the Krylov space context
469: Options Database Keys:
470: + -ksp_max_it - maximum number of linear iterations
471: . -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e.
472: if residual norm decreases by this factor than convergence is declared
473: . -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual
474: norm is less than this then convergence is declared
475: . -ksp_divtol tol - if residual norm increases by this factor than divergence is declared
476: . -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using
477: $ convergence test (say you always want to run with 5 iterations) to
478: $ save on communication overhead
479: $ preconditioned - default for left preconditioning
480: $ unpreconditioned - see KSPSetNormType()
481: $ natural - see KSPSetNormType()
482: . -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space
483: . -ksp_test_null_space - tests the null space set with KSPSetNullSpace() to see if it truly is a null space
484: . -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side
485: . -ksp_cancelmonitors - cancel all previous convergene monitor routines set
486: . -ksp_monitor - print residual norm at each iteration
487: . -ksp_xmonitor - plot residual norm at each iteration
488: . -ksp_vecmonitor - plot solution at each iteration
489: - -ksp_singmonitor - monitor extremem singular values at each iteration
491: Notes:
492: To see all options, run your program with the -help option
493: or consult the users manual.
495: Level: beginner
497: .keywords: KSP, set, from, options, database
499: .seealso:
500: @*/
501: PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP ksp)
502: {
504: PetscInt indx;
505: char type[256];
506: const char *stype[] = {"none","preconditioned","unpreconditioned","natural"};
507: PetscTruth flg;
511: PCSetFromOptions(ksp->pc);
513: if (!KSPRegisterAllCalled) {KSPRegisterAll(PETSC_NULL);}
514: PetscOptionsBegin(ksp->comm,ksp->prefix,"Krylov Method (KSP) Options","KSP");
515: PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(ksp->type_name?ksp->type_name:KSPGMRES),type,256,&flg);
516: if (flg) {
517: KSPSetType(ksp,type);
518: }
519: /*
520: Set the type if it was never set.
521: */
522: if (!ksp->type_name) {
523: KSPSetType(ksp,KSPGMRES);
524: }
526: PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,PETSC_NULL);
527: PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,PETSC_NULL);
528: PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->abstol,&ksp->abstol,PETSC_NULL);
529: PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,PETSC_NULL);
530: PetscOptionsTruth("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll,
531: &ksp->guess_knoll,PETSC_NULL);
533: PetscOptionsEList("-ksp_norm_type","KSP Norm type","KSPSetNormType",stype,4,"preconditioned",&indx,&flg);
534: if (flg) {
535: switch (indx) {
536: case 0:
537: KSPSetNormType(ksp,KSP_NO_NORM);
538: KSPSetConvergenceTest(ksp,KSPSkipConverged,0);
539: break;
540: case 1:
541: KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM);
542: break;
543: case 2:
544: KSPSetNormType(ksp,KSP_UNPRECONDITIONED_NORM);
545: break;
546: case 3:
547: KSPSetNormType(ksp,KSP_NATURAL_NORM);
548: break;
549: }
550: }
552: PetscOptionsName("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",&flg);
553: if (flg) {
554: KSPSetDiagonalScale(ksp,PETSC_TRUE);
555: }
556: PetscOptionsName("-ksp_diagonal_scale_fix","Fix diagonaled scaled matrix after solve","KSPSetDiagonalScaleFix",&flg);
557: if (flg) {
558: KSPSetDiagonalScaleFix(ksp,PETSC_TRUE);
559: }
562: PetscOptionsName("-ksp_constant_null_space","Add constant null space to Krylov solver","KSPSetNullSpace",&flg);
563: if (flg) {
564: MatNullSpace nsp;
566: MatNullSpaceCreate(ksp->comm,PETSC_TRUE,0,0,&nsp);
567: KSPSetNullSpace(ksp,nsp);
568: MatNullSpaceDestroy(nsp);
569: }
571: /* option is actually checked in KSPSetUp() */
572: if (ksp->nullsp) {
573: PetscOptionsName("-ksp_test_null_space","Is provided null space correct","None",&flg);
574: }
576: /*
577: Prints reason for convergence or divergence of each linear solve
578: */
579: PetscOptionsName("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",&flg);
580: if (flg) {
581: ksp->printreason = PETSC_TRUE;
582: }
584: PetscOptionsName("-ksp_cancelmonitors","Remove any hardwired monitor routines","KSPClearMonitor",&flg);
585: /* -----------------------------------------------------------------------*/
586: /*
587: Cancels all monitors hardwired into code before call to KSPSetFromOptions()
588: */
589: if (flg) {
590: KSPClearMonitor(ksp);
591: }
592: /*
593: Prints preconditioned residual norm at each iteration
594: */
595: PetscOptionsName("-ksp_monitor","Monitor preconditioned residual norm","KSPSetMonitor",&flg);
596: if (flg) {
597: KSPSetMonitor(ksp,KSPDefaultMonitor,PETSC_NULL,PETSC_NULL);
598: }
599: /*
600: Plots the vector solution
601: */
602: PetscOptionsName("-ksp_vecmonitor","Monitor solution graphically","KSPSetMonitor",&flg);
603: if (flg) {
604: KSPSetMonitor(ksp,KSPVecViewMonitor,PETSC_NULL,PETSC_NULL);
605: }
606: /*
607: Prints preconditioned and true residual norm at each iteration
608: */
609: PetscOptionsName("-ksp_truemonitor","Monitor true (unpreconditioned) residual norm","KSPSetMonitor",&flg);
610: if (flg) {
611: KSPSetMonitor(ksp,KSPTrueMonitor,PETSC_NULL,PETSC_NULL);
612: }
613: /*
614: Prints extreme eigenvalue estimates at each iteration
615: */
616: PetscOptionsName("-ksp_singmonitor","Monitor singular values","KSPSetMonitor",&flg);
617: if (flg) {
618: KSPSetComputeSingularValues(ksp,PETSC_TRUE);
619: KSPSetMonitor(ksp,KSPSingularValueMonitor,PETSC_NULL,PETSC_NULL);
620: }
621: /*
622: Prints preconditioned residual norm with fewer digits
623: */
624: PetscOptionsName("-ksp_smonitor","Monitor preconditioned residual norm with fewer digitis","KSPSetMonitor",&flg);
625: if (flg) {
626: KSPSetMonitor(ksp,KSPDefaultSMonitor,PETSC_NULL,PETSC_NULL);
627: }
628: /*
629: Graphically plots preconditioned residual norm
630: */
631: PetscOptionsName("-ksp_xmonitor","Monitor graphically preconditioned residual norm","KSPSetMonitor",&flg);
632: if (flg) {
633: KSPSetMonitor(ksp,KSPLGMonitor,PETSC_NULL,PETSC_NULL);
634: }
635: /*
636: Graphically plots preconditioned and true residual norm
637: */
638: PetscOptionsName("-ksp_xtruemonitor","Monitor graphically true residual norm","KSPSetMonitor",&flg);
639: if (flg){
640: KSPSetMonitor(ksp,KSPLGTrueMonitor,PETSC_NULL,PETSC_NULL);
641: }
643: /* -----------------------------------------------------------------------*/
645: PetscOptionsTruthGroupBegin("-ksp_left_pc","Use left preconditioning","KSPSetPreconditionerSide",&flg);
646: if (flg) { KSPSetPreconditionerSide(ksp,PC_LEFT); }
647: PetscOptionsTruthGroup("-ksp_right_pc","Use right preconditioning","KSPSetPreconditionerSide",&flg);
648: if (flg) { KSPSetPreconditionerSide(ksp,PC_RIGHT);}
649: PetscOptionsTruthGroupEnd("-ksp_symmetric_pc","Use symmetric (factorized) preconditioning","KSPSetPreconditionerSide",&flg);
650: if (flg) { KSPSetPreconditionerSide(ksp,PC_SYMMETRIC);}
652: PetscOptionsName("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",&flg);
653: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
654: PetscOptionsName("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",&flg);
655: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
656: PetscOptionsName("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",&flg);
657: if (flg) { KSPSetComputeSingularValues(ksp,PETSC_TRUE); }
659: if (ksp->ops->setfromoptions) {
660: (*ksp->ops->setfromoptions)(ksp);
661: }
662: PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);
663: PetscOptionsEnd();
664: return(0);
665: }
669: /*@C
670: KSPRegister - See KSPRegisterDynamic()
672: Level: advanced
673: @*/
674: PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(KSP))
675: {
677: char fullname[PETSC_MAX_PATH_LEN];
680: PetscFListConcat(path,name,fullname);
681: PetscFListAdd(&KSPList,sname,fullname,(void (*)(void))function);
682: return(0);
683: }
687: /*@C
688: KSPSetNullSpace - Sets the null space of the operator
690: Collective on KSP
692: Input Parameters:
693: + ksp - the Krylov space object
694: - nullsp - the null space of the operator
696: Level: advanced
698: .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPGetNullSpace()
699: @*/
700: PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP ksp,MatNullSpace nullsp)
701: {
705: if (ksp->nullsp) {
706: MatNullSpaceDestroy(ksp->nullsp);
707: }
708: ksp->nullsp = nullsp;
709: PetscObjectReference((PetscObject)ksp->nullsp);
710: return(0);
711: }
715: /*@C
716: KSPGetNullSpace - Gets the null space of the operator
718: Collective on KSP
720: Input Parameters:
721: + ksp - the Krylov space object
722: - nullsp - the null space of the operator
724: Level: advanced
726: .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPSetNullSpace()
727: @*/
728: PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP ksp,MatNullSpace *nullsp)
729: {
731: *nullsp = ksp->nullsp;
732: return(0);
733: }