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 21bca11509SBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES 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" 39*265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 40df2a969aSvictorle #define KSPCGNE "cgne" 411155b781STodd Munson #define KSPCGNASH "nash" 421155b781STodd Munson #define KSPCGSTCG "stcg" 431155b781STodd Munson #define KSPCGGLTR "gltr" 44640f4f53SPatrick Sanan #define KSPFCG "fcg" 45390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 4682bf6240SBarry Smith #define KSPGMRES "gmres" 47483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 489dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 499dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 504f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5161b468f9SJed Brown #define KSPPGMRES "pgmres" 5282bf6240SBarry Smith #define KSPTCQMR "tcqmr" 5382bf6240SBarry Smith #define KSPBCGS "bcgs" 54e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5518ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 56c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 57f0037002Svictorle #define KSPBCGSL "bcgsl" 58f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 5982bf6240SBarry Smith #define KSPCGS "cgs" 6082bf6240SBarry Smith #define KSPTFQMR "tfqmr" 6182bf6240SBarry Smith #define KSPCR "cr" 622c8fc646SJed Brown #define KSPPIPECR "pipecr" 6382bf6240SBarry Smith #define KSPLSQR "lsqr" 6482bf6240SBarry Smith #define KSPPREONLY "preonly" 6582bf6240SBarry Smith #define KSPQCG "qcg" 66c9cf9b11SBarry Smith #define KSPBICG "bicg" 67b4ac9ba4SBarry Smith #define KSPMINRES "minres" 6801c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 6921356919SSatish Balay #define KSPLCD "lcd" 701d6018f0SLisandro Dalcin #define KSPPYTHON "python" 7158ad63f4SBarry Smith #define KSPGCR "gcr" 72fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 73e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 74e4d80e07Szianekhodja #define KSPCGLS "cgls" 75329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 762eac72dbSBarry Smith 778ba1e511SMatthew Knepley /* Logging support */ 78014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 79ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 805358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 818ba1e511SMatthew Knepley 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*); 8319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 84c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*); 85014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 86014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 88014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 9123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 9225c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 932eac72dbSBarry Smith 94140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 95ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 96bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 9730de9b25SBarry Smith 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 101c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool); 1077d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool); 108c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*); 109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool); 110c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*); 111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*); 112014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*); 113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 115734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1162a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1172a7a6963SBarry 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);} 1182eac72dbSBarry Smith 119d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 120d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 121d4211eb9SBarry Smith 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 124aabeff55SBarry Smith 125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**)); 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1314b0e389bSBarry Smith 132fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 133fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*); 134fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 135fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 136fa0ddf94SBarry Smith 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 142b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 143b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 144b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 145b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*); 1470a71e23eSMatthew Knepley 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*); 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*); 1502eac72dbSBarry Smith 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 15458450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 155b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 15658450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 158d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*); 159d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1607d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]); 1614b0e389bSBarry Smith 162640f4f53SPatrick Sanan /*E 163640f4f53SPatrick Sanan 16406137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 165640f4f53SPatrick Sanan 16693f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 16793f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 168640f4f53SPatrick Sanan 1692b26979fSBarry Smith Level: intermediate 17006137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 171640f4f53SPatrick Sanan 172640f4f53SPatrick Sanan E*/ 17393f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 17406137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 175640f4f53SPatrick Sanan 176640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 177640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 178640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 179640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 18006137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 18106137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 182640f4f53SPatrick Sanan 183390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 184390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 185390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 186390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 187390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 188390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 189390d8e47SPatrick Sanan 190fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 191fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 192fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 193fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 194fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 195fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 196fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 197fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 19883f0b094SBarry Smith 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 2029f236934SBarry Smith 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2081d73ed98SBarry Smith 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2111d73ed98SBarry Smith 212483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 213483d6965SPatrick Sanan 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 21758ad63f4SBarry Smith 21804d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*); 21904d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC); 22004d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*); 221568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat); 222b49fd9e1SBarry Smith /*E 223b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 224b49fd9e1SBarry Smith 225b49fd9e1SBarry Smith Level: advanced 226b49fd9e1SBarry Smith 227bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 228e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 229b49fd9e1SBarry Smith 230b49fd9e1SBarry Smith E*/ 23178d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2326a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2331f7e983dSSatish Balay /*MC 2341957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2358c5b8ba0SBarry Smith 2368c5b8ba0SBarry Smith Level: advanced 2378c5b8ba0SBarry Smith 2388c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2398c5b8ba0SBarry Smith 240bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 241e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2428c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2438c5b8ba0SBarry Smith M*/ 2448c5b8ba0SBarry Smith 2451f7e983dSSatish Balay /*MC 2468c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2478c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2488c5b8ba0SBarry Smith poor orthogonality. 2498c5b8ba0SBarry Smith 2508c5b8ba0SBarry Smith Level: advanced 2518c5b8ba0SBarry Smith 2528c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2538c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2548c5b8ba0SBarry Smith 255bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 256e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2578c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2588c5b8ba0SBarry Smith M*/ 2598c5b8ba0SBarry Smith 2601f7e983dSSatish Balay /*MC 2618c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2628c5b8ba0SBarry Smith 2638c5b8ba0SBarry Smith Level: advanced 2648c5b8ba0SBarry Smith 2658c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2668c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2678c5b8ba0SBarry Smith 2688c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2698c5b8ba0SBarry Smith 270bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 271e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2728c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2738c5b8ba0SBarry Smith M*/ 2748c5b8ba0SBarry Smith 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 27708480c60SBarry Smith 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 281c38d4ed2SBarry Smith 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 285121fd945SKris Buschelman 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 289e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 290d9492815SBarry Smith 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 2932eac72dbSBarry Smith 294fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 295fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 296fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 297fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 298390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 299390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 300fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 301fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 302fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 303fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 304e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 305e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 306e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*); 308fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*)); 30984cb2905SBarry Smith 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 312c01c455dSBarry Smith 31323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 31423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 315014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*); 316014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3197287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 3207287d2a3SDmitry Karpeev PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 3211eb62cbbSBarry Smith 322014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 323014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*); 324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*); 3261f7f0c4fSBarry Smith 327014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 32855849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 329685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 3302a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 3312a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 33255849f57SBarry Smith 33355849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3341eb62cbbSBarry Smith 335014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 33763aaf0cbSSatish Balay PETSC_EXTERN PetscErrorCode KSPLSQRGetArnorm(KSP,PetscReal*,PetscReal*,PetscReal*); 338db9b2ab1SHong Zhang 339014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 34168ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 34283ab6a24SBarry Smith 34328ce4d24SBarry Smith /*E 3448a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3458a4b9c5cSBarry Smith test routines. 3468a4b9c5cSBarry Smith 3478a4b9c5cSBarry Smith Level: advanced 3488a4b9c5cSBarry Smith 349a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3501718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 351a3f661c8SBarry Smith 352af0996ceSBarry Smith Notes: this must match petsc/finclude/petscksp.h 3538a4b9c5cSBarry Smith 35494b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3551718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3568a4b9c5cSBarry Smith E*/ 3579e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3589e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 359014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 360a21b2a99SBarry Smith 3611f7e983dSSatish Balay /*MC 3629793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3638c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3648c5b8ba0SBarry Smith be based on a norm of a residual etc. 3658c5b8ba0SBarry Smith 3668c5b8ba0SBarry Smith Level: advanced 3678c5b8ba0SBarry Smith 3681957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3698c5b8ba0SBarry Smith 370ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3718c5b8ba0SBarry Smith M*/ 3728c5b8ba0SBarry Smith 3731f7e983dSSatish Balay /*MC 3741957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3758c5b8ba0SBarry Smith convergence test routine. 3768c5b8ba0SBarry Smith 3778c5b8ba0SBarry Smith Level: advanced 3788c5b8ba0SBarry Smith 3799793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3808c5b8ba0SBarry Smith M*/ 3818c5b8ba0SBarry Smith 3821f7e983dSSatish Balay /*MC 383ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3848c5b8ba0SBarry Smith convergence test routine. 3858c5b8ba0SBarry Smith 3868c5b8ba0SBarry Smith Level: advanced 3878c5b8ba0SBarry Smith 3889793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3898c5b8ba0SBarry Smith M*/ 3908c5b8ba0SBarry Smith 3911f7e983dSSatish Balay /*MC 392ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 393390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 3948c5b8ba0SBarry Smith 3958c5b8ba0SBarry Smith Level: advanced 3968c5b8ba0SBarry Smith 3979793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 3988c5b8ba0SBarry Smith M*/ 3998c5b8ba0SBarry Smith 400014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 4058a4b9c5cSBarry Smith 4068a4b9c5cSBarry Smith /*E 4071957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 40828ce4d24SBarry Smith 40928ce4d24SBarry Smith Level: beginner 41028ce4d24SBarry Smith 4114d0a8057SBarry Smith Notes: See KSPGetConvergedReason() for explanation of each value 41228ce4d24SBarry Smith 413af0996ceSBarry Smith Developer notes: this must match petsc/finclude/petscksp.h 4144d0a8057SBarry Smith 4154d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4164d0a8057SBarry Smith any of the values here also change them that array of names. 41786c02ca4SBarry Smith 418c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 41928ce4d24SBarry Smith E*/ 420d15094e1SBarry Smith typedef enum {/* converged */ 4219ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4229ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 423d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 424d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 425b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4268031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4278031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 428329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 429af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 430d15094e1SBarry Smith /* diverged */ 431b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 432d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 433d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 434d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 435b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 436b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 437b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4384d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4396aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 440422a814eSBarry Smith KSP_DIVERGED_PCSETUP_FAILED = -11, 441d15094e1SBarry Smith 442d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 443014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 444d15094e1SBarry Smith 445c838673bSBarry Smith /*MC 446f9fed41fSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called 447c838673bSBarry Smith 448c838673bSBarry Smith Level: beginner 449c838673bSBarry Smith 450c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 451c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 452c838673bSBarry Smith 2-norm of the residual for right preconditioning 453c838673bSBarry Smith 454f9fed41fSBarry Smith See also KSP_CONVERGED_ATOL which may apply before this tolerance. 455f9fed41fSBarry Smith 456f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 457c838673bSBarry Smith 458c838673bSBarry Smith M*/ 459c838673bSBarry Smith 460c838673bSBarry Smith /*MC 461c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 462c838673bSBarry Smith 463c838673bSBarry Smith Level: beginner 464c838673bSBarry Smith 465c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 466c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 467c838673bSBarry Smith 2-norm of the residual for right preconditioning 468c838673bSBarry Smith 469f9fed41fSBarry Smith See also KSP_CONVERGED_RTOL which may apply before this tolerance. 470c838673bSBarry Smith 471f9fed41fSBarry Smith .seealso: KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 472c838673bSBarry Smith 473c838673bSBarry Smith M*/ 474c838673bSBarry Smith 475c838673bSBarry Smith /*MC 476c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 477c838673bSBarry Smith 478c838673bSBarry Smith Level: beginner 479c838673bSBarry Smith 480c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 481c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 482c838673bSBarry Smith 2-norm of the residual for right preconditioning 483c838673bSBarry Smith 484c838673bSBarry Smith Level: beginner 485c838673bSBarry Smith 486f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 487c838673bSBarry Smith 488c838673bSBarry Smith M*/ 489c838673bSBarry Smith 490c838673bSBarry Smith /*MC 491c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 492c838673bSBarry Smith reached 493c838673bSBarry Smith 494c838673bSBarry Smith Level: beginner 495c838673bSBarry Smith 496c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 497c838673bSBarry Smith 498c838673bSBarry Smith M*/ 499c838673bSBarry Smith 500c838673bSBarry Smith /*MC 5018c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 5020059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 5034d0a8057SBarry Smith test routine is set in KSP. 504c838673bSBarry Smith 505c838673bSBarry Smith Level: beginner 506c838673bSBarry Smith 507c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 508c838673bSBarry Smith 509c838673bSBarry Smith M*/ 510c838673bSBarry Smith 511c838673bSBarry Smith /*MC 512c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 5133014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 5143014e516SBarry Smith preconditioner. 515c838673bSBarry Smith 516c838673bSBarry Smith Level: beginner 517c838673bSBarry Smith 518c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 519c838673bSBarry Smith 520c838673bSBarry Smith M*/ 521c838673bSBarry Smith 522c838673bSBarry Smith /*MC 523c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 524c838673bSBarry Smith method could not continue to enlarge the Krylov space. 525c838673bSBarry Smith 526c838673bSBarry Smith Level: beginner 527c838673bSBarry Smith 528c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 529c838673bSBarry Smith 530c838673bSBarry Smith M*/ 531c838673bSBarry Smith 532c838673bSBarry Smith /*MC 533c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 534c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 535c838673bSBarry Smith 536c838673bSBarry Smith Level: beginner 537c838673bSBarry Smith 538c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 539c838673bSBarry Smith 540c838673bSBarry Smith M*/ 541c838673bSBarry Smith 542c838673bSBarry Smith /*MC 543c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 544c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 545c838673bSBarry Smith be positive definite 546c838673bSBarry Smith 547c838673bSBarry Smith Level: beginner 548c838673bSBarry Smith 5492401956bSBarry Smith Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 550c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 551c838673bSBarry Smith 552c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 553c838673bSBarry Smith 554c838673bSBarry Smith M*/ 555c838673bSBarry Smith 556c838673bSBarry Smith /*MC 5579fc87aa7SBarry Smith KSP_DIVERGED_PCSETUP_FAILED - It was not possible to build the requested preconditioner. This is usually due to a 5589fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 5599fc87aa7SBarry Smith such as PCFIELDSPLIT. 5609fc87aa7SBarry Smith 5619fc87aa7SBarry Smith Level: beginner 5629fc87aa7SBarry Smith 5639fc87aa7SBarry 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. 5649fc87aa7SBarry Smith 5659fc87aa7SBarry Smith 5669fc87aa7SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 5679fc87aa7SBarry Smith 5689fc87aa7SBarry Smith M*/ 5699fc87aa7SBarry Smith 5709fc87aa7SBarry Smith /*MC 571c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 572c838673bSBarry Smith while the KSPSolve() is still running. 573c838673bSBarry Smith 574c838673bSBarry Smith Level: beginner 575c838673bSBarry Smith 576c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 577c838673bSBarry Smith 578c838673bSBarry Smith M*/ 579c838673bSBarry Smith 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*)); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**); 5828de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 5848de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*); 5858de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**); 5868de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 5878de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 5880059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*); 590abef13c0SSatish Balay 5918ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 5928ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 5938ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 5948ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 5958ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 5968ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 5978ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 5988ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 5998ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 6008ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 6018ea1b3e6SJed Brown PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 6028ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 6038ea1b3e6SJed Brown 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat*); 605d4fbbf0eSBarry Smith 60628ce4d24SBarry Smith /*E 60728ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 60828ce4d24SBarry Smith 60928ce4d24SBarry Smith Level: beginner 61028ce4d24SBarry Smith 61128ce4d24SBarry Smith .seealso: KSPCGSetType() 61228ce4d24SBarry Smith E*/ 6139dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 6146a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 61528ce4d24SBarry Smith 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 6188031f4adStmunson 619dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal); 620dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*); 621dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*); 622fcae7a14Stmunson 623d2e0dd5eSTodd Munson PETSC_EXTERN PetscErrorCode KSPCGGLTRGetMinEig(KSP,PetscReal*); 624d2e0dd5eSTodd Munson PETSC_EXTERN PetscErrorCode KSPCGGLTRGetLambda(KSP,PetscReal*); 6258031f4adStmunson 626014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 6271d6018f0SLisandro Dalcin 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 629014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 6303369ce9aSBarry Smith 6319804daf3SBarry Smith #include <petscdrawtypes.h> 632d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 633d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 634d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 635d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 636014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6372f2e5d10SKris Buschelman 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 64003919abeSBarry Smith 641ba36735cSStefano Zampini /*S 642ba36735cSStefano Zampini KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 643f8a50e2bSBarry Smith 644ba36735cSStefano Zampini Level: beginner 645f8a50e2bSBarry Smith 646ba36735cSStefano Zampini Concepts: Krylov methods 647ba36735cSStefano Zampini 648ba36735cSStefano Zampini .seealso: KSPCreate(), KSPSetGuessType(), KSPGuessType 649ba36735cSStefano Zampini S*/ 650ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess; 651ba36735cSStefano Zampini /*J 652ba36735cSStefano Zampini KSPGuessType - String with the name of a PETSc initial guess for Krylov methods. 653ba36735cSStefano Zampini 654ba36735cSStefano Zampini Level: beginner 655ba36735cSStefano Zampini 656ba36735cSStefano Zampini .seealso: KSPGuess 657ba36735cSStefano Zampini J*/ 658ba36735cSStefano Zampini typedef const char* KSPGuessType; 659ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 660ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 661ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess); 662ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*); 663ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer); 664ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*); 665ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*); 666ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType); 667ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*); 668ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 669ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec); 670ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec); 671ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 672ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt); 673014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 674ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool); 675ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*); 676f8a50e2bSBarry Smith 677470b340bSDmitry Karpeev /*E 678470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 679470b340bSDmitry Karpeev 680470b340bSDmitry Karpeev Level: intermediate 681470b340bSDmitry Karpeev 682470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 683470b340bSDmitry Karpeev E*/ 6840c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType; 685470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 686470b340bSDmitry Karpeev 687014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 688014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 689d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 690bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 691aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 692bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 693470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 694470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 6955bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 6965a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 697470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*); 698470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 6993f22127dSBarry Smith 700014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 701014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 702014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 703014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 704014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 705fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*); 70623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 707fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 70823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 70923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 710014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 711014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 712fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 713fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 7146c699258SBarry Smith 71502b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 7168c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, 717e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 718e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 719e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 720191494d9SMatthew G. Knepley PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec); 7212eac72dbSBarry Smith #endif 722