1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 426bd1501SBarry Smith #ifndef PETSCKSP_H 526bd1501SBarry Smith #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 1695452b02SPatrick Sanan Notes: 1795452b02SPatrick Sanan When a direct solver is used but no Krylov solver is used the KSP object is still used by with a 188f6c3df8SBarry Smith KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver). 198f6c3df8SBarry Smith 20bca11509SBarry Smith .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy(), KSPCG, KSPGMRES 2128ce4d24SBarry Smith S*/ 22e2a1c21fSSatish Balay typedef struct _p_KSP* KSP; 232eac72dbSBarry Smith 2476bdecfbSBarry Smith /*J 258f6c3df8SBarry Smith KSPType - String with the name of a PETSc Krylov method. 2628ce4d24SBarry Smith 2728ce4d24SBarry Smith Level: beginner 2828ce4d24SBarry Smith 298f6c3df8SBarry Smith .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 3076bdecfbSBarry Smith J*/ 3119fd82e9SBarry Smith typedef const char* KSPType; 3282bf6240SBarry Smith #define KSPRICHARDSON "richardson" 336c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3482bf6240SBarry Smith #define KSPCG "cg" 352c8fc646SJed Brown #define KSPGROPPCG "groppcg" 362c8fc646SJed Brown #define KSPPIPECG "pipecg" 37901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 38265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 39b21a8899STyler Chen #define KSPPIPEPRCG "pipeprcg" 40325e8391SManasi #define KSPPIPECG2 "pipecg2" 41df2a969aSvictorle #define KSPCGNE "cgne" 4205de396fSBarry Smith #define KSPNASH "nash" 4305de396fSBarry Smith #define KSPSTCG "stcg" 4405de396fSBarry Smith #define KSPGLTR "gltr" 4505de396fSBarry Smith #define KSPCGNASH PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash" 4605de396fSBarry Smith #define KSPCGSTCG PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg" 4705de396fSBarry Smith #define KSPCGGLTR PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr" 48640f4f53SPatrick Sanan #define KSPFCG "fcg" 49390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 5082bf6240SBarry Smith #define KSPGMRES "gmres" 51483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 529dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 539dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 544f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5561b468f9SJed Brown #define KSPPGMRES "pgmres" 5682bf6240SBarry Smith #define KSPTCQMR "tcqmr" 5782bf6240SBarry Smith #define KSPBCGS "bcgs" 58e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5918ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 60c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 61f0037002Svictorle #define KSPBCGSL "bcgsl" 62f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 6382bf6240SBarry Smith #define KSPCGS "cgs" 6482bf6240SBarry Smith #define KSPTFQMR "tfqmr" 6582bf6240SBarry Smith #define KSPCR "cr" 662c8fc646SJed Brown #define KSPPIPECR "pipecr" 6782bf6240SBarry Smith #define KSPLSQR "lsqr" 6882bf6240SBarry Smith #define KSPPREONLY "preonly" 6982bf6240SBarry Smith #define KSPQCG "qcg" 70c9cf9b11SBarry Smith #define KSPBICG "bicg" 71b4ac9ba4SBarry Smith #define KSPMINRES "minres" 7201c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 7321356919SSatish Balay #define KSPLCD "lcd" 741d6018f0SLisandro Dalcin #define KSPPYTHON "python" 7558ad63f4SBarry Smith #define KSPGCR "gcr" 76fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 77e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 78e4d80e07Szianekhodja #define KSPCGLS "cgls" 79329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 8038cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 812eac72dbSBarry Smith 828ba1e511SMatthew Knepley /* Logging support */ 83014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 84ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 855358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 868ba1e511SMatthew Knepley 87014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*); 8819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 89c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 93014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 94*75f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP,PetscBool); 955ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat); 965ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPSetMatSolveBlockSize(KSP,PetscInt); 975ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPGetMatSolveBlockSize(KSP,PetscInt*); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 99ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 10123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 1027d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*); 10325c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 104c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec); 1052eac72dbSBarry Smith 106140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 107ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 108bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 10930de9b25SBarry Smith 110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 112014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 113c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool); 115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*); 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool); 117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*); 118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool); 1197d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool); 120c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool); 122c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*); 124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*); 125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 127734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1282a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 12925ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1302eac72dbSBarry Smith 131d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 132d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 133d4211eb9SBarry Smith 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 136aabeff55SBarry Smith 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**)); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**); 14194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool); 14394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*); 14494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool); 1454b0e389bSBarry Smith 146fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 147fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*); 148fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 149fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 150fa0ddf94SBarry Smith 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 152cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 154014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 155014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 157285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]); 158b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 159b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 160b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 161b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 162014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*); 1634a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*); 164f3b08a26SMatthew G. Knepley /* 165f3b08a26SMatthew G. Knepley PCMGCoarseList contains the list of coarse space constructor currently registered 166f3b08a26SMatthew G. Knepley These are added with PCMGRegisterCoarseSpaceConstructor() 167f3b08a26SMatthew G. Knepley */ 168f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 169f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)); 170f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)); 171f3b08a26SMatthew G. Knepley 1720a71e23eSMatthew Knepley 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*); 1752eac72dbSBarry Smith 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 17958450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 180b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 18158450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 183d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*); 184d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1857d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]); 1864b0e389bSBarry Smith 187640f4f53SPatrick Sanan /*E 188640f4f53SPatrick Sanan 18906137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 190640f4f53SPatrick Sanan 19193f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 19293f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 193640f4f53SPatrick Sanan 1942b26979fSBarry Smith Level: intermediate 19506137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 196640f4f53SPatrick Sanan 197640f4f53SPatrick Sanan E*/ 19893f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 19906137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 200640f4f53SPatrick Sanan 201640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 202640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 203640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 204640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 20506137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 20606137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 207640f4f53SPatrick Sanan 208390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 209390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 210390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 211390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 212390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 213390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 214390d8e47SPatrick Sanan 215fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 216fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 217fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 218fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 219fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 220fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 221fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 222fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 22383f0b094SBarry Smith 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 227e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal); 2289f236934SBarry Smith 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2341d73ed98SBarry Smith 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2371d73ed98SBarry Smith 238483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 239483d6965SPatrick Sanan 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 24358ad63f4SBarry Smith 24404d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*); 24504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC); 24604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*); 247568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat); 248e82af88dSprj- 249e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat); 250e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*); 2515ea7661aSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); } 252d55ab951SPierre Jolivet /*E 253d55ab951SPierre Jolivet KSPHPDDMType - Type of Krylov method used by KSPHPDDM 254d55ab951SPierre Jolivet 255d55ab951SPierre Jolivet Level: intermediate 256d55ab951SPierre Jolivet 257d55ab951SPierre Jolivet Values: 258d55ab951SPierre Jolivet + KSP_HPDDM_TYPE_GMRES (default) 259d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BGMRES 260d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_CG 261d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BCG 262d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_GCRODR 263d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BGCRODR 264d55ab951SPierre Jolivet . KSP_HPDDM_TYPE_BFBCG 265d55ab951SPierre Jolivet - KSP_HPDDM_TYPE_PREONLY 266d55ab951SPierre Jolivet 267d55ab951SPierre Jolivet .seealso: KSPHPDDM, KSPHPDDMSetType() 268d55ab951SPierre Jolivet E*/ 269d55ab951SPierre Jolivet typedef enum { 270d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 271d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 272d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 273d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 274d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 275d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 276d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 277d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 278d55ab951SPierre Jolivet } KSPHPDDMType; 279d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 280d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType); 281d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*); 282e82af88dSprj- 283b49fd9e1SBarry Smith /*E 284b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 285b49fd9e1SBarry Smith 286b49fd9e1SBarry Smith Level: advanced 287b49fd9e1SBarry Smith 288bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 289e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 290b49fd9e1SBarry Smith 291b49fd9e1SBarry Smith E*/ 29278d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2936a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2941f7e983dSSatish Balay /*MC 2951957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2968c5b8ba0SBarry Smith 2978c5b8ba0SBarry Smith Level: advanced 2988c5b8ba0SBarry Smith 2998c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 3008c5b8ba0SBarry Smith 301bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 302e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 3038c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 3048c5b8ba0SBarry Smith M*/ 3058c5b8ba0SBarry Smith 3061f7e983dSSatish Balay /*MC 3078c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 3088c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 3098c5b8ba0SBarry Smith poor orthogonality. 3108c5b8ba0SBarry Smith 3118c5b8ba0SBarry Smith Level: advanced 3128c5b8ba0SBarry Smith 3138c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 3148c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 3158c5b8ba0SBarry Smith 316bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 317e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 3188c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 3198c5b8ba0SBarry Smith M*/ 3208c5b8ba0SBarry Smith 3211f7e983dSSatish Balay /*MC 3228c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 3238c5b8ba0SBarry Smith 3248c5b8ba0SBarry Smith Level: advanced 3258c5b8ba0SBarry Smith 3268c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 3278c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 3288c5b8ba0SBarry Smith 3298c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 3308c5b8ba0SBarry Smith 331bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 332e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 3338c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 3348c5b8ba0SBarry Smith M*/ 3358c5b8ba0SBarry Smith 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 337014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 33808480c60SBarry Smith 339014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 340014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 341014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 342c38d4ed2SBarry Smith 343014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 344014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 346121fd945SKris Buschelman 347014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool); 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 350e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 351d9492815SBarry Smith 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 35387d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 3552eac72dbSBarry Smith 356fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 357fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 358fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 359fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 360390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 361390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 362fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 363fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 364fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 365fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 366e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 367e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 368e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*); 370fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*)); 37184cb2905SBarry Smith 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 374c01c455dSBarry Smith 37523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 37623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3811eb62cbbSBarry Smith 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*); 3861f7f0c4fSBarry Smith 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 38855849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 389fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]); 39019a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer); 3911b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 39294a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer); 3931b2b9847SBarry Smith 39419a666eeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);} 3951b2b9847SBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);} 39655849f57SBarry Smith 39755849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3981eb62cbbSBarry Smith 399dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool); 4000e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool); 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 402884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*); 403db9b2ab1SHong Zhang 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 40668ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 40783ab6a24SBarry Smith 40828ce4d24SBarry Smith /*E 4098a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 4108a4b9c5cSBarry Smith test routines. 4118a4b9c5cSBarry Smith 4128a4b9c5cSBarry Smith Level: advanced 4138a4b9c5cSBarry Smith 414a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 4151718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 416a3f661c8SBarry Smith 41795452b02SPatrick Sanan Notes: 41895452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 4198a4b9c5cSBarry Smith 42094b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 4211718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 4228a4b9c5cSBarry Smith E*/ 4239e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 4249e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 425014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 426a21b2a99SBarry Smith 4271f7e983dSSatish Balay /*MC 4289793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 4298c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 4308c5b8ba0SBarry Smith be based on a norm of a residual etc. 4318c5b8ba0SBarry Smith 4328c5b8ba0SBarry Smith Level: advanced 4338c5b8ba0SBarry Smith 4341957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 4358c5b8ba0SBarry Smith 436ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 4378c5b8ba0SBarry Smith M*/ 4388c5b8ba0SBarry Smith 4391f7e983dSSatish Balay /*MC 4401957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 4418c5b8ba0SBarry Smith convergence test routine. 4428c5b8ba0SBarry Smith 4438c5b8ba0SBarry Smith Level: advanced 4448c5b8ba0SBarry Smith 4459793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 4468c5b8ba0SBarry Smith M*/ 4478c5b8ba0SBarry Smith 4481f7e983dSSatish Balay /*MC 449ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 4508c5b8ba0SBarry Smith convergence test routine. 4518c5b8ba0SBarry Smith 4528c5b8ba0SBarry Smith Level: advanced 4538c5b8ba0SBarry Smith 4549793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 4558c5b8ba0SBarry Smith M*/ 4568c5b8ba0SBarry Smith 4571f7e983dSSatish Balay /*MC 458ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 459390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 4608c5b8ba0SBarry Smith 4618c5b8ba0SBarry Smith Level: advanced 4628c5b8ba0SBarry Smith 4639793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 4648c5b8ba0SBarry Smith M*/ 4658c5b8ba0SBarry Smith 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 4718a4b9c5cSBarry Smith 472f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 4738a4b9c5cSBarry Smith /*E 4741957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 47528ce4d24SBarry Smith 47628ce4d24SBarry Smith Level: beginner 47728ce4d24SBarry Smith 47895452b02SPatrick Sanan Notes: 47995452b02SPatrick Sanan See KSPGetConvergedReason() for explanation of each value 48028ce4d24SBarry Smith 48195452b02SPatrick Sanan Developer Notes: 48295452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 4834d0a8057SBarry Smith 4844d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4854d0a8057SBarry Smith any of the values here also change them that array of names. 48686c02ca4SBarry Smith 4871b2b9847SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances(), KSPConvergedReasonView() 48828ce4d24SBarry Smith E*/ 489d15094e1SBarry Smith typedef enum {/* converged */ 4909ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4919ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 492d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 493d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 494b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4958031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4968031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 497329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 498af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 499d15094e1SBarry Smith /* diverged */ 500b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 501d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 502d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 503d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 504b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 505b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 506b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 5074d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 5086aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 509c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 510aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 511d15094e1SBarry Smith 512d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 513014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 514d15094e1SBarry Smith 515c838673bSBarry Smith /*MC 516f9fed41fSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called 517c838673bSBarry Smith 518c838673bSBarry Smith Level: beginner 519c838673bSBarry Smith 520c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 521c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 522c838673bSBarry Smith 2-norm of the residual for right preconditioning 523c838673bSBarry Smith 524f9fed41fSBarry Smith See also KSP_CONVERGED_ATOL which may apply before this tolerance. 525f9fed41fSBarry Smith 526f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 527c838673bSBarry Smith 528c838673bSBarry Smith M*/ 529c838673bSBarry Smith 530c838673bSBarry Smith /*MC 531c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 532c838673bSBarry Smith 533c838673bSBarry Smith Level: beginner 534c838673bSBarry Smith 535c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 536c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 537c838673bSBarry Smith 2-norm of the residual for right preconditioning 538c838673bSBarry Smith 539f9fed41fSBarry Smith See also KSP_CONVERGED_RTOL which may apply before this tolerance. 540c838673bSBarry Smith 541f9fed41fSBarry Smith .seealso: KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 542c838673bSBarry Smith 543c838673bSBarry Smith M*/ 544c838673bSBarry Smith 545c838673bSBarry Smith /*MC 546c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 547c838673bSBarry Smith 548c838673bSBarry Smith Level: beginner 549c838673bSBarry Smith 550c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 551c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 552c838673bSBarry Smith 2-norm of the residual for right preconditioning 553c838673bSBarry Smith 554c838673bSBarry Smith Level: beginner 555c838673bSBarry Smith 556f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 557c838673bSBarry Smith 558c838673bSBarry Smith M*/ 559c838673bSBarry Smith 560c838673bSBarry Smith /*MC 561c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 562c838673bSBarry Smith reached 563c838673bSBarry Smith 564c838673bSBarry Smith Level: beginner 565c838673bSBarry Smith 566c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 567c838673bSBarry Smith 568c838673bSBarry Smith M*/ 569c838673bSBarry Smith 570c838673bSBarry Smith /*MC 5718c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 5720059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 5734d0a8057SBarry Smith test routine is set in KSP. 574c838673bSBarry Smith 575c838673bSBarry Smith Level: beginner 576c838673bSBarry Smith 577c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 578c838673bSBarry Smith 579c838673bSBarry Smith M*/ 580c838673bSBarry Smith 581c838673bSBarry Smith /*MC 582c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 5831de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 5841de96524SPierre Jolivet preconditioner. In KSPHPDDM, this is also returned when some search directions within a block 5851de96524SPierre Jolivet are colinear. 586c838673bSBarry Smith 587c838673bSBarry Smith Level: beginner 588c838673bSBarry Smith 589c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 590c838673bSBarry Smith 591c838673bSBarry Smith M*/ 592c838673bSBarry Smith 593c838673bSBarry Smith /*MC 594c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 595c838673bSBarry Smith method could not continue to enlarge the Krylov space. 596c838673bSBarry Smith 597c838673bSBarry Smith Level: beginner 598c838673bSBarry Smith 599c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 600c838673bSBarry Smith 601c838673bSBarry Smith M*/ 602c838673bSBarry Smith 603c838673bSBarry Smith /*MC 604c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 605c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 606c838673bSBarry Smith 607c838673bSBarry Smith Level: beginner 608c838673bSBarry Smith 609c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 610c838673bSBarry Smith 611c838673bSBarry Smith M*/ 612c838673bSBarry Smith 613c838673bSBarry Smith /*MC 614c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 615c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 616c838673bSBarry Smith be positive definite 617c838673bSBarry Smith 618c838673bSBarry Smith Level: beginner 619c838673bSBarry Smith 62095452b02SPatrick Sanan Notes: 62195452b02SPatrick Sanan This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 622c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 623c838673bSBarry Smith 624c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 625c838673bSBarry Smith 626c838673bSBarry Smith M*/ 627c838673bSBarry Smith 628c838673bSBarry Smith /*MC 629c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 6309fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 6319fc87aa7SBarry Smith such as PCFIELDSPLIT. 6329fc87aa7SBarry Smith 6339fc87aa7SBarry Smith Level: beginner 6349fc87aa7SBarry Smith 63595452b02SPatrick Sanan Notes: 63695452b02SPatrick Sanan Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 6379fc87aa7SBarry Smith 6389fc87aa7SBarry Smith 6399fc87aa7SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 6409fc87aa7SBarry Smith 6419fc87aa7SBarry Smith M*/ 6429fc87aa7SBarry Smith 6439fc87aa7SBarry Smith /*MC 644c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 645c838673bSBarry Smith while the KSPSolve() is still running. 646c838673bSBarry Smith 647c838673bSBarry Smith Level: beginner 648c838673bSBarry Smith 649c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 650c838673bSBarry Smith 651c838673bSBarry Smith M*/ 652c838673bSBarry Smith 653014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*)); 654df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 655df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 656014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**); 6578de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 6583eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 6598de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*); 6608de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**); 6618de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 6628de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 66354b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool); 6640059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 665014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*); 66694a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*); 667abef13c0SSatish Balay 66825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 6698ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 67025ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 6718ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 67225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 6738ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 67425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 6758ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 67625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 6778ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 67825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 6798ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 6808ea1b3e6SJed Brown 6810bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*); 68225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); } 683d4fbbf0eSBarry Smith 68428ce4d24SBarry Smith /*E 68528ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 68628ce4d24SBarry Smith 68728ce4d24SBarry Smith Level: beginner 68828ce4d24SBarry Smith 68928ce4d24SBarry Smith .seealso: KSPCGSetType() 69028ce4d24SBarry Smith E*/ 6919dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 6926a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 69328ce4d24SBarry Smith 694014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 695014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool); 6968031f4adStmunson 697dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal); 698dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*); 699dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*); 700fcae7a14Stmunson 70105de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*); 70205de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*); 70325ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);} 70425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);} 7058031f4adStmunson 706014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 7071d6018f0SLisandro Dalcin 708014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 709014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 7103369ce9aSBarry Smith 7119804daf3SBarry Smith #include <petscdrawtypes.h> 712d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 713d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 714d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 715d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 716014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 7172f2e5d10SKris Buschelman 718014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 719014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 72003919abeSBarry Smith 721ba36735cSStefano Zampini /*S 722ba36735cSStefano Zampini KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 723f8a50e2bSBarry Smith 724ba36735cSStefano Zampini Level: beginner 725f8a50e2bSBarry Smith 726ba36735cSStefano Zampini .seealso: KSPCreate(), KSPSetGuessType(), KSPGuessType 727ba36735cSStefano Zampini S*/ 728ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess; 729ba36735cSStefano Zampini /*J 730ba36735cSStefano Zampini KSPGuessType - String with the name of a PETSc initial guess for Krylov methods. 731ba36735cSStefano Zampini 732ba36735cSStefano Zampini Level: beginner 733ba36735cSStefano Zampini 734ba36735cSStefano Zampini .seealso: KSPGuess 735ba36735cSStefano Zampini J*/ 736ba36735cSStefano Zampini typedef const char* KSPGuessType; 737ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 738ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 7391d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess)); 740ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess); 741ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*); 742ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer); 743ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*); 744ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*); 745ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType); 746ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*); 747ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 748ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec); 749ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec); 750ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 751ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt); 752014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 753ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool); 754ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*); 755f8a50e2bSBarry Smith 756470b340bSDmitry Karpeev /*E 757470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 758470b340bSDmitry Karpeev 759470b340bSDmitry Karpeev Level: intermediate 760470b340bSDmitry Karpeev 761470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 762470b340bSDmitry Karpeev E*/ 7630c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType; 764470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 765470b340bSDmitry Karpeev 766864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 767864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 768864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 769864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_USER = 3} MatLMVMSymBroydenScaleType; 770864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 77192f76d53SAlp Dener 772014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 773014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 774d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 775bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 776aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 777bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 778470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 779470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 7805bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 7815a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 782470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*); 783470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 7843f22127dSBarry Smith 78578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*); 78678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*); 78778e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*); 788864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 789864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 790864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 791864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 792864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 793cd929ea3SAlp Dener 794cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 795b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*); 796cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 797cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 798e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 7990ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 800cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 801cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 802cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 803cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 804cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 8052d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 806cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 807cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*); 808cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*); 809cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*); 81092f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 811cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*); 812cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*); 813864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 814864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 815cd929ea3SAlp Dener 816014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 817014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool); 818014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 819014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 820014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 821fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*); 82223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 823fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 82423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 82523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 826014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 827014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 828fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 829fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 8306c699258SBarry Smith 83102b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 8328c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, 833e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 834e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 835e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 836191494d9SMatthew G. Knepley PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec); 837557cf195SMatthew G. Knepley 838557cf195SMatthew G. Knepley 839557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *); 840557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal); 8412eac72dbSBarry Smith #endif 842