1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCKSP_H 50a835dfdSSatish Balay #define __PETSCKSP_H 62c8e378dSBarry Smith #include <petscpc.h> 72eac72dbSBarry Smith 8607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 91dbb0a54SBarry Smith 1028ce4d24SBarry Smith /*S 118f6c3df8SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 128f6c3df8SBarry Smith linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 1328ce4d24SBarry Smith 1428ce4d24SBarry Smith Level: beginner 1528ce4d24SBarry Smith 1628ce4d24SBarry Smith Concepts: Krylov methods 1728ce4d24SBarry Smith 188f6c3df8SBarry Smith Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a 198f6c3df8SBarry Smith KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver). 208f6c3df8SBarry Smith 218f6c3df8SBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy() 2228ce4d24SBarry Smith S*/ 23e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 242eac72dbSBarry Smith 2576bdecfbSBarry Smith /*J 268f6c3df8SBarry Smith KSPType - String with the name of a PETSc Krylov method. 2728ce4d24SBarry Smith 2828ce4d24SBarry Smith Level: beginner 2928ce4d24SBarry Smith 308f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 3176bdecfbSBarry Smith J*/ 3219fd82e9SBarry Smith typedef const char* KSPType; 3382bf6240SBarry Smith #define KSPRICHARDSON "richardson" 346c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3582bf6240SBarry Smith #define KSPCG "cg" 362c8fc646SJed Brown #define KSPGROPPCG "groppcg" 372c8fc646SJed Brown #define KSPPIPECG "pipecg" 38901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 39df2a969aSvictorle #define KSPCGNE "cgne" 40fcae7a14Stmunson #define KSPNASH "nash" 4180e17431SMatthew Knepley #define KSPSTCG "stcg" 428031f4adStmunson #define KSPGLTR "gltr" 43640f4f53SPatrick Sanan #define KSPFCG "fcg" 44390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 4582bf6240SBarry Smith #define KSPGMRES "gmres" 46483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 479dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 489dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 494f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5061b468f9SJed Brown #define KSPPGMRES "pgmres" 5182bf6240SBarry Smith #define KSPTCQMR "tcqmr" 5282bf6240SBarry Smith #define KSPBCGS "bcgs" 53e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5418ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 55c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 56f0037002Svictorle #define KSPBCGSL "bcgsl" 5782bf6240SBarry Smith #define KSPCGS "cgs" 5882bf6240SBarry Smith #define KSPTFQMR "tfqmr" 5982bf6240SBarry Smith #define KSPCR "cr" 602c8fc646SJed Brown #define KSPPIPECR "pipecr" 6182bf6240SBarry Smith #define KSPLSQR "lsqr" 6282bf6240SBarry Smith #define KSPPREONLY "preonly" 6382bf6240SBarry Smith #define KSPQCG "qcg" 64c9cf9b11SBarry Smith #define KSPBICG "bicg" 65b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6601c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 6721356919SSatish Balay #define KSPLCD "lcd" 681d6018f0SLisandro Dalcin #define KSPPYTHON "python" 6958ad63f4SBarry Smith #define KSPGCR "gcr" 70fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 71e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 72e4d80e07Szianekhodja #define KSPCGLS "cgls" 73*329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 742eac72dbSBarry Smith 758ba1e511SMatthew Knepley /* Logging support */ 76014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 775358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 788ba1e511SMatthew Knepley 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 8019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 81c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 83014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 84014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 85014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 86014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 8823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 8925c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 902eac72dbSBarry Smith 91140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 92bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 9330de9b25SBarry Smith 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 97c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 101014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 1057d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool ); 106c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 108c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 112014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 113734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1142a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1152a7a6963SBarry Smith PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1162eac72dbSBarry Smith 117d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 118d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 119d4211eb9SBarry Smith 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 122aabeff55SBarry Smith 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1294b0e389bSBarry Smith 130fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 131fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 132fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 133fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 134fa0ddf94SBarry Smith 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 140b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 141b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 142b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 143b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 1450a71e23eSMatthew Knepley 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 1482eac72dbSBarry Smith 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 15258450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 15382fe98a5SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool); 154ba142b81SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 15558450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 157d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *); 158d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1597d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]); 1604b0e389bSBarry Smith 161640f4f53SPatrick Sanan /*E 162640f4f53SPatrick Sanan 16306137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 164640f4f53SPatrick Sanan 16593f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 16693f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 167640f4f53SPatrick Sanan 1682b26979fSBarry Smith Level: intermediate 16906137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 170640f4f53SPatrick Sanan 171640f4f53SPatrick Sanan E*/ 17293f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 17306137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 174640f4f53SPatrick Sanan 175640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 176640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 177640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 178640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 17906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 18006137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 181640f4f53SPatrick Sanan 182390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 183390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 184390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 185390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 186390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 187390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 188390d8e47SPatrick Sanan 189fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 190fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 191fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 192fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 193fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 194fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 195fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 196fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 19783f0b094SBarry Smith 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 2019f236934SBarry Smith 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2071d73ed98SBarry Smith 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2101d73ed98SBarry Smith 211483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 212483d6965SPatrick Sanan 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 21658ad63f4SBarry Smith 217b49fd9e1SBarry Smith /*E 218b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 219b49fd9e1SBarry Smith 220b49fd9e1SBarry Smith Level: advanced 221b49fd9e1SBarry Smith 222bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 223e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 224b49fd9e1SBarry Smith 225b49fd9e1SBarry Smith E*/ 22678d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2276a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2281f7e983dSSatish Balay /*MC 2291957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2308c5b8ba0SBarry Smith 2318c5b8ba0SBarry Smith Level: advanced 2328c5b8ba0SBarry Smith 2338c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2348c5b8ba0SBarry Smith 235bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 236e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2378c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2388c5b8ba0SBarry Smith M*/ 2398c5b8ba0SBarry Smith 2401f7e983dSSatish Balay /*MC 2418c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2428c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2438c5b8ba0SBarry Smith poor orthogonality. 2448c5b8ba0SBarry Smith 2458c5b8ba0SBarry Smith Level: advanced 2468c5b8ba0SBarry Smith 2478c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2488c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2498c5b8ba0SBarry Smith 250bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 251e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2528c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2538c5b8ba0SBarry Smith M*/ 2548c5b8ba0SBarry Smith 2551f7e983dSSatish Balay /*MC 2568c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2578c5b8ba0SBarry Smith 2588c5b8ba0SBarry Smith Level: advanced 2598c5b8ba0SBarry Smith 2608c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2618c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2628c5b8ba0SBarry Smith 2638c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2648c5b8ba0SBarry Smith 265bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 266e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2678c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2688c5b8ba0SBarry Smith M*/ 2698c5b8ba0SBarry Smith 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 27208480c60SBarry Smith 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 276c38d4ed2SBarry Smith 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 280121fd945SKris Buschelman 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 284e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 285d9492815SBarry Smith 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2882eac72dbSBarry Smith 289fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 290fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 291fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 292fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 293390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 294390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 295fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 296fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 297fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 298fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat *); 299e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 300e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 301e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 303fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*)); 30484cb2905SBarry Smith 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 307c01c455dSBarry Smith 30823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 30923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 313014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3147287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3157287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3161eb62cbbSBarry Smith 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 320014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 3211f7f0c4fSBarry Smith 322014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 32355849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 324685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 3252a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 3262a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 32755849f57SBarry Smith 32855849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3291eb62cbbSBarry Smith 330014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 331014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 332db9b2ab1SHong Zhang 333014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 334014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 33568ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 33683ab6a24SBarry Smith 33728ce4d24SBarry Smith /*E 3388a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3398a4b9c5cSBarry Smith test routines. 3408a4b9c5cSBarry Smith 3418a4b9c5cSBarry Smith Level: advanced 3428a4b9c5cSBarry Smith 343a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3441718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 345a3f661c8SBarry Smith 346af0996ceSBarry Smith Notes: this must match petsc/finclude/petscksp.h 3478a4b9c5cSBarry Smith 34894b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3491718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3508a4b9c5cSBarry Smith E*/ 3519e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3529e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 353014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 354a21b2a99SBarry Smith 3551f7e983dSSatish Balay /*MC 3569793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3578c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3588c5b8ba0SBarry Smith be based on a norm of a residual etc. 3598c5b8ba0SBarry Smith 3608c5b8ba0SBarry Smith Level: advanced 3618c5b8ba0SBarry Smith 3621957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3638c5b8ba0SBarry Smith 364ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3658c5b8ba0SBarry Smith M*/ 3668c5b8ba0SBarry Smith 3671f7e983dSSatish Balay /*MC 3681957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3698c5b8ba0SBarry Smith convergence test routine. 3708c5b8ba0SBarry Smith 3718c5b8ba0SBarry Smith Level: advanced 3728c5b8ba0SBarry Smith 3739793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3748c5b8ba0SBarry Smith M*/ 3758c5b8ba0SBarry Smith 3761f7e983dSSatish Balay /*MC 377ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3788c5b8ba0SBarry Smith convergence test routine. 3798c5b8ba0SBarry Smith 3808c5b8ba0SBarry Smith Level: advanced 3818c5b8ba0SBarry Smith 3829793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3838c5b8ba0SBarry Smith M*/ 3848c5b8ba0SBarry Smith 3851f7e983dSSatish Balay /*MC 386ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 387390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 3888c5b8ba0SBarry Smith 3898c5b8ba0SBarry Smith Level: advanced 3908c5b8ba0SBarry Smith 3919793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3928c5b8ba0SBarry Smith M*/ 3938c5b8ba0SBarry Smith 394014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 396014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 397014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 3998a4b9c5cSBarry Smith 4008a4b9c5cSBarry Smith /*E 4011957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 40228ce4d24SBarry Smith 40328ce4d24SBarry Smith Level: beginner 40428ce4d24SBarry Smith 4054d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 40628ce4d24SBarry Smith 407af0996ceSBarry Smith Developer notes: this must match petsc/finclude/petscksp.h 4084d0a8057SBarry Smith 4094d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4104d0a8057SBarry Smith any of the values here also change them that array of names. 41186c02ca4SBarry Smith 412c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 41328ce4d24SBarry Smith E*/ 414d15094e1SBarry Smith typedef enum {/* converged */ 4159ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4169ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 417d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 418d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 419b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4208031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4218031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 422329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 423af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 424d15094e1SBarry Smith /* diverged */ 425b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 426d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 427d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 428d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 429b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 430b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 431b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4324d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4336aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 434422a814eSBarry Smith KSP_DIVERGED_PCSETUP_FAILED = -11, 435d15094e1SBarry Smith 436d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 437014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 438d15094e1SBarry Smith 439c838673bSBarry Smith /*MC 440c838673bSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 441c838673bSBarry Smith 442c838673bSBarry Smith Level: beginner 443c838673bSBarry Smith 444c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 445c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 446c838673bSBarry Smith 2-norm of the residual for right preconditioning 447c838673bSBarry Smith 448c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 449c838673bSBarry Smith 450c838673bSBarry Smith M*/ 451c838673bSBarry Smith 452c838673bSBarry Smith /*MC 453c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 454c838673bSBarry Smith 455c838673bSBarry Smith Level: beginner 456c838673bSBarry Smith 457c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 458c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 459c838673bSBarry Smith 2-norm of the residual for right preconditioning 460c838673bSBarry Smith 461c838673bSBarry Smith Level: beginner 462c838673bSBarry Smith 463c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 464c838673bSBarry Smith 465c838673bSBarry Smith M*/ 466c838673bSBarry Smith 467c838673bSBarry Smith /*MC 468c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 469c838673bSBarry Smith 470c838673bSBarry Smith Level: beginner 471c838673bSBarry Smith 472c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 473c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 474c838673bSBarry Smith 2-norm of the residual for right preconditioning 475c838673bSBarry Smith 476c838673bSBarry Smith Level: beginner 477c838673bSBarry Smith 478c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 479c838673bSBarry Smith 480c838673bSBarry Smith M*/ 481c838673bSBarry Smith 482c838673bSBarry Smith /*MC 483c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 484c838673bSBarry Smith reached 485c838673bSBarry Smith 486c838673bSBarry Smith Level: beginner 487c838673bSBarry Smith 488c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 489c838673bSBarry Smith 490c838673bSBarry Smith M*/ 491c838673bSBarry Smith 492c838673bSBarry Smith /*MC 4938c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 4940059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 4954d0a8057SBarry Smith test routine is set in KSP. 496c838673bSBarry Smith 497c838673bSBarry Smith Level: beginner 498c838673bSBarry Smith 499c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 500c838673bSBarry Smith 501c838673bSBarry Smith M*/ 502c838673bSBarry Smith 503c838673bSBarry Smith /*MC 504c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 5053014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 5063014e516SBarry Smith preconditioner. 507c838673bSBarry Smith 508c838673bSBarry Smith Level: beginner 509c838673bSBarry Smith 510c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 511c838673bSBarry Smith 512c838673bSBarry Smith M*/ 513c838673bSBarry Smith 514c838673bSBarry Smith /*MC 515c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 516c838673bSBarry Smith method could not continue to enlarge the Krylov space. 517c838673bSBarry Smith 518c838673bSBarry Smith Level: beginner 519c838673bSBarry Smith 520c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 521c838673bSBarry Smith 522c838673bSBarry Smith M*/ 523c838673bSBarry Smith 524c838673bSBarry Smith /*MC 525c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 526c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 527c838673bSBarry Smith 528c838673bSBarry Smith Level: beginner 529c838673bSBarry Smith 530c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 531c838673bSBarry Smith 532c838673bSBarry Smith M*/ 533c838673bSBarry Smith 534c838673bSBarry Smith /*MC 535c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 536c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 537c838673bSBarry Smith be positive definite 538c838673bSBarry Smith 539c838673bSBarry Smith Level: beginner 540c838673bSBarry Smith 5412401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 542c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 543c838673bSBarry Smith 544c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 545c838673bSBarry Smith 546c838673bSBarry Smith M*/ 547c838673bSBarry Smith 548c838673bSBarry Smith /*MC 5499fc87aa7SBarry Smith KSP_DIVERGED_PCSETUP_FAILED - It was not possible to build the requested preconditioner. This is usually due to a 5509fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 5519fc87aa7SBarry Smith such as PCFIELDSPLIT. 5529fc87aa7SBarry Smith 5539fc87aa7SBarry Smith Level: beginner 5549fc87aa7SBarry Smith 5559fc87aa7SBarry Smith Notes: Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 5569fc87aa7SBarry Smith 5579fc87aa7SBarry Smith 5589fc87aa7SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 5599fc87aa7SBarry Smith 5609fc87aa7SBarry Smith M*/ 5619fc87aa7SBarry Smith 5629fc87aa7SBarry Smith /*MC 563c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 564c838673bSBarry Smith while the KSPSolve() is still running. 565c838673bSBarry Smith 566c838673bSBarry Smith Level: beginner 567c838673bSBarry Smith 568c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 569c838673bSBarry Smith 570c838673bSBarry Smith M*/ 571c838673bSBarry Smith 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 5748de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 5768de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 5778de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 5788de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5798de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5800059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 582abef13c0SSatish Balay 5838ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5848ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5858ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5868ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5878ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5888ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5898ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5908ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5918ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 5928ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 5938ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 5948ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 5958ea1b3e6SJed Brown 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 597d4fbbf0eSBarry Smith 59828ce4d24SBarry Smith /*E 59928ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 60028ce4d24SBarry Smith 60128ce4d24SBarry Smith Level: beginner 60228ce4d24SBarry Smith 60328ce4d24SBarry Smith .seealso: KSPCGSetType() 60428ce4d24SBarry Smith E*/ 6059dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 6066a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 60728ce4d24SBarry Smith 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 6108031f4adStmunson 611014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 612014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 613014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 614fcae7a14Stmunson 615014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 618e559a7feSLois Curfman McInnes 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 622014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 6248031f4adStmunson 625014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 6261d6018f0SLisandro Dalcin 627014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 6293369ce9aSBarry Smith 6309804daf3SBarry Smith #include <petscdrawtypes.h> 631d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 632d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 633d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 634d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 635014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6362f2e5d10SKris Buschelman 637014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 63903919abeSBarry Smith 640f8a50e2bSBarry Smith /* see src/ksp/ksp/interface/iguess.c */ 641ace3abfcSBarry Smith typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 642f8a50e2bSBarry Smith 643014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 644014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 645014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 646014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 647014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 648014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 649f8a50e2bSBarry Smith 650014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 651014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 652014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 653f8a50e2bSBarry Smith 654470b340bSDmitry Karpeev /*E 655470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 656470b340bSDmitry Karpeev 657470b340bSDmitry Karpeev Level: intermediate 658470b340bSDmitry Karpeev 659470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 660470b340bSDmitry Karpeev E*/ 661470b340bSDmitry Karpeev typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 662470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 663470b340bSDmitry Karpeev 664014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 665014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 666d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 667bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 668aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 669bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 670470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 671470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 6725bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 6735a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 674470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 675470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 6763f22127dSBarry Smith 677014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 678014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 679014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 680014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 681014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 682fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 68323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 684fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 68523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 68623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 687014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 688014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 689fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 690fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 6916c699258SBarry Smith 69202b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 693bfc4295aSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectField(DM, Vec, 694e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 695e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 696e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 697e29c0834SMatthew G. Knepley PetscReal, const PetscReal [], PetscScalar []), InsertMode, Vec); 6982eac72dbSBarry Smith #endif 699