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" 39df2a969aSvictorle #define KSPCGNE "cgne" 4005de396fSBarry Smith #define KSPNASH "nash" 4105de396fSBarry Smith #define KSPSTCG "stcg" 4205de396fSBarry Smith #define KSPGLTR "gltr" 4305de396fSBarry Smith #define KSPCGNASH PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash" 4405de396fSBarry Smith #define KSPCGSTCG PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg" 4505de396fSBarry Smith #define KSPCGGLTR PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr" 46640f4f53SPatrick Sanan #define KSPFCG "fcg" 47390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 4882bf6240SBarry Smith #define KSPGMRES "gmres" 49483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 509dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 519dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 524f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5361b468f9SJed Brown #define KSPPGMRES "pgmres" 5482bf6240SBarry Smith #define KSPTCQMR "tcqmr" 5582bf6240SBarry Smith #define KSPBCGS "bcgs" 56e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 5718ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 58c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 59f0037002Svictorle #define KSPBCGSL "bcgsl" 60f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 6182bf6240SBarry Smith #define KSPCGS "cgs" 6282bf6240SBarry Smith #define KSPTFQMR "tfqmr" 6382bf6240SBarry Smith #define KSPCR "cr" 642c8fc646SJed Brown #define KSPPIPECR "pipecr" 6582bf6240SBarry Smith #define KSPLSQR "lsqr" 6682bf6240SBarry Smith #define KSPPREONLY "preonly" 6782bf6240SBarry Smith #define KSPQCG "qcg" 68c9cf9b11SBarry Smith #define KSPBICG "bicg" 69b4ac9ba4SBarry Smith #define KSPMINRES "minres" 7001c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 7121356919SSatish Balay #define KSPLCD "lcd" 721d6018f0SLisandro Dalcin #define KSPPYTHON "python" 7358ad63f4SBarry Smith #define KSPGCR "gcr" 74fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 75e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 76e4d80e07Szianekhodja #define KSPCGLS "cgls" 77329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 7838cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 792eac72dbSBarry Smith 808ba1e511SMatthew Knepley /* Logging support */ 81014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 82ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 835358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 848ba1e511SMatthew Knepley 85014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*); 8619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 87c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*); 88014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 89014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 90014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 93ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 9523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 9625c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 97c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec); 982eac72dbSBarry Smith 99140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 100ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 101bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 10230de9b25SBarry Smith 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 106c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*); 109014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool); 110014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*); 111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool); 1127d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool); 113c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*); 114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool); 115c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*); 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*); 117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*); 118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 120734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1212a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 12225ef9dfeSBarry 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);} 1232eac72dbSBarry Smith 124d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 125d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 126d4211eb9SBarry Smith 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 129aabeff55SBarry Smith 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**)); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void**); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt*); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 1364b0e389bSBarry Smith 137fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 138fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*); 139fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 140fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 141fa0ddf94SBarry Smith 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 143cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 148285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]); 149b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 150b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 151b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 152b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*); 1544a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*); 1550a71e23eSMatthew Knepley 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*); 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*); 1582eac72dbSBarry Smith 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 161014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 16258450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 163b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 16458450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 166d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*); 167d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1687d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]); 1694b0e389bSBarry Smith 170640f4f53SPatrick Sanan /*E 171640f4f53SPatrick Sanan 17206137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 173640f4f53SPatrick Sanan 17493f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 17593f1e87bSPatrick Sanan KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 176640f4f53SPatrick Sanan 1772b26979fSBarry Smith Level: intermediate 17806137d0aSPatrick Sanan .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 179640f4f53SPatrick Sanan 180640f4f53SPatrick Sanan E*/ 18193f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 18206137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 183640f4f53SPatrick Sanan 184640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 185640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 186640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 187640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 18806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 18906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 190640f4f53SPatrick Sanan 191390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 192390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 193390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 194390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 195390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 196390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 197390d8e47SPatrick Sanan 198fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 199fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 200fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 201fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 202fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 203fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 204fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 205fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 20683f0b094SBarry Smith 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 2109f236934SBarry Smith 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2161d73ed98SBarry Smith 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2191d73ed98SBarry Smith 220483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 221483d6965SPatrick Sanan 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 22558ad63f4SBarry Smith 22604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*); 22704d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC); 22804d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*); 229568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat); 230*e82af88dSprj- 231*e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationSpace(KSP,Mat); 232*e82af88dSprj- PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationSpace(KSP,Mat*); 233*e82af88dSprj- 234b49fd9e1SBarry Smith /*E 235b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 236b49fd9e1SBarry Smith 237b49fd9e1SBarry Smith Level: advanced 238b49fd9e1SBarry Smith 239bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 240e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 241b49fd9e1SBarry Smith 242b49fd9e1SBarry Smith E*/ 24378d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 2446a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 2451f7e983dSSatish Balay /*MC 2461957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 2478c5b8ba0SBarry Smith 2488c5b8ba0SBarry Smith Level: advanced 2498c5b8ba0SBarry Smith 2508c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 2518c5b8ba0SBarry Smith 252bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 253e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2548c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2558c5b8ba0SBarry Smith M*/ 2568c5b8ba0SBarry Smith 2571f7e983dSSatish Balay /*MC 2588c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 2598c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 2608c5b8ba0SBarry Smith poor orthogonality. 2618c5b8ba0SBarry Smith 2628c5b8ba0SBarry Smith Level: advanced 2638c5b8ba0SBarry Smith 2648c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 2658c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 2668c5b8ba0SBarry Smith 267bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 268e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 2698c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2708c5b8ba0SBarry Smith M*/ 2718c5b8ba0SBarry Smith 2721f7e983dSSatish Balay /*MC 2738c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 2748c5b8ba0SBarry Smith 2758c5b8ba0SBarry Smith Level: advanced 2768c5b8ba0SBarry Smith 2778c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 2788c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 2798c5b8ba0SBarry Smith 2808c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 2818c5b8ba0SBarry Smith 282bef9c725SBarry Smith .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 283e928de20SBarry Smith KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 2848c5b8ba0SBarry Smith KSPGMRESModifiedGramSchmidtOrthogonalization() 2858c5b8ba0SBarry Smith M*/ 2868c5b8ba0SBarry Smith 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 28908480c60SBarry Smith 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 293c38d4ed2SBarry Smith 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 297121fd945SKris Buschelman 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 301e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 302d9492815SBarry Smith 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 30487d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 3062eac72dbSBarry Smith 307fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 308fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 309fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 310fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 311390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 312390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 313fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 314fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 315fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 316fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 317e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 318e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 319e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 320014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*); 321fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],const char [],PetscErrorCode (*)(KSP,PetscInt,PetscReal,PetscViewerAndFormat*)); 32284cb2905SBarry Smith 323014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 324014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 325c01c455dSBarry Smith 32623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 32723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*); 329014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 330014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 331014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 3321eb62cbbSBarry Smith 333014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 334014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*); 335014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*); 3371f7f0c4fSBarry Smith 338014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 33955849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 340fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]); 3412a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 3422a359c20SBarry Smith PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 34355849f57SBarry Smith 34455849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 3451eb62cbbSBarry Smith 346dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool); 3470e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool); 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 349884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*); 350db9b2ab1SHong Zhang 351014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 35368ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 35483ab6a24SBarry Smith 35528ce4d24SBarry Smith /*E 3568a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 3578a4b9c5cSBarry Smith test routines. 3588a4b9c5cSBarry Smith 3598a4b9c5cSBarry Smith Level: advanced 3608a4b9c5cSBarry Smith 361a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 3621718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 363a3f661c8SBarry Smith 36495452b02SPatrick Sanan Notes: 36595452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 3668a4b9c5cSBarry Smith 36794b7f48cSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 3681718446fSBarry Smith KSPSetConvergenceTest(), KSPSetPCSide() 3698a4b9c5cSBarry Smith E*/ 3709e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 3719e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 372014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 373a21b2a99SBarry Smith 3741f7e983dSSatish Balay /*MC 3759793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 3768c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 3778c5b8ba0SBarry Smith be based on a norm of a residual etc. 3788c5b8ba0SBarry Smith 3798c5b8ba0SBarry Smith Level: advanced 3808c5b8ba0SBarry Smith 3811957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 3828c5b8ba0SBarry Smith 383ce9499c7SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 3848c5b8ba0SBarry Smith M*/ 3858c5b8ba0SBarry Smith 3861f7e983dSSatish Balay /*MC 3871957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 3888c5b8ba0SBarry Smith convergence test routine. 3898c5b8ba0SBarry Smith 3908c5b8ba0SBarry Smith Level: advanced 3918c5b8ba0SBarry Smith 3929793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 3938c5b8ba0SBarry Smith M*/ 3948c5b8ba0SBarry Smith 3951f7e983dSSatish Balay /*MC 396ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 3978c5b8ba0SBarry Smith convergence test routine. 3988c5b8ba0SBarry Smith 3998c5b8ba0SBarry Smith Level: advanced 4008c5b8ba0SBarry Smith 4019793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 4028c5b8ba0SBarry Smith M*/ 4038c5b8ba0SBarry Smith 4041f7e983dSSatish Balay /*MC 405ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 406390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 4078c5b8ba0SBarry Smith 4088c5b8ba0SBarry Smith Level: advanced 4098c5b8ba0SBarry Smith 4109793b452SBarry Smith .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 4118c5b8ba0SBarry Smith M*/ 4128c5b8ba0SBarry Smith 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 4188a4b9c5cSBarry Smith 4198a4b9c5cSBarry Smith /*E 4201957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 42128ce4d24SBarry Smith 42228ce4d24SBarry Smith Level: beginner 42328ce4d24SBarry Smith 42495452b02SPatrick Sanan Notes: 42595452b02SPatrick Sanan See KSPGetConvergedReason() for explanation of each value 42628ce4d24SBarry Smith 42795452b02SPatrick Sanan Developer Notes: 42895452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 4294d0a8057SBarry Smith 4304d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 4314d0a8057SBarry Smith any of the values here also change them that array of names. 43286c02ca4SBarry Smith 433c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 43428ce4d24SBarry Smith E*/ 43525ef9dfeSBarry Smith #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 436d15094e1SBarry Smith typedef enum {/* converged */ 4379ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 4389ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 439d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 440d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 441b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 4428031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 4438031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 444329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 445af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 446d15094e1SBarry Smith /* diverged */ 447b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 448d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 449d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 450d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 451b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 452b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 453b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 4544d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 4556aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 456c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 457aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 458d15094e1SBarry Smith 459d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 460014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 461d15094e1SBarry Smith 462c838673bSBarry Smith /*MC 463f9fed41fSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called 464c838673bSBarry Smith 465c838673bSBarry Smith Level: beginner 466c838673bSBarry Smith 467c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 468c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 469c838673bSBarry Smith 2-norm of the residual for right preconditioning 470c838673bSBarry Smith 471f9fed41fSBarry Smith See also KSP_CONVERGED_ATOL which may apply before this tolerance. 472f9fed41fSBarry Smith 473f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 474c838673bSBarry Smith 475c838673bSBarry Smith M*/ 476c838673bSBarry Smith 477c838673bSBarry Smith /*MC 478c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 479c838673bSBarry Smith 480c838673bSBarry Smith Level: beginner 481c838673bSBarry Smith 482c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 483c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 484c838673bSBarry Smith 2-norm of the residual for right preconditioning 485c838673bSBarry Smith 486f9fed41fSBarry Smith See also KSP_CONVERGED_RTOL which may apply before this tolerance. 487c838673bSBarry Smith 488f9fed41fSBarry Smith .seealso: KSP_CONVERGED_RTOL, KSP_DIVERGED_DTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 489c838673bSBarry Smith 490c838673bSBarry Smith M*/ 491c838673bSBarry Smith 492c838673bSBarry Smith /*MC 493c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 494c838673bSBarry Smith 495c838673bSBarry Smith Level: beginner 496c838673bSBarry Smith 497c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 498c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 499c838673bSBarry Smith 2-norm of the residual for right preconditioning 500c838673bSBarry Smith 501c838673bSBarry Smith Level: beginner 502c838673bSBarry Smith 503f9fed41fSBarry Smith .seealso: KSP_CONVERGED_ATOL, KSP_DIVERGED_RTOL, KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 504c838673bSBarry Smith 505c838673bSBarry Smith M*/ 506c838673bSBarry Smith 507c838673bSBarry Smith /*MC 508c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 509c838673bSBarry Smith reached 510c838673bSBarry Smith 511c838673bSBarry Smith Level: beginner 512c838673bSBarry Smith 513c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 514c838673bSBarry Smith 515c838673bSBarry Smith M*/ 516c838673bSBarry Smith 517c838673bSBarry Smith /*MC 5188c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 5190059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 5204d0a8057SBarry Smith test routine is set in KSP. 521c838673bSBarry Smith 522c838673bSBarry Smith Level: beginner 523c838673bSBarry Smith 524c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 525c838673bSBarry Smith 526c838673bSBarry Smith M*/ 527c838673bSBarry Smith 528c838673bSBarry Smith /*MC 529c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 5303014e516SBarry Smith method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 5313014e516SBarry Smith preconditioner. 532c838673bSBarry Smith 533c838673bSBarry Smith Level: beginner 534c838673bSBarry Smith 535c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 536c838673bSBarry Smith 537c838673bSBarry Smith M*/ 538c838673bSBarry Smith 539c838673bSBarry Smith /*MC 540c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 541c838673bSBarry Smith method could not continue to enlarge the Krylov space. 542c838673bSBarry Smith 543c838673bSBarry Smith Level: beginner 544c838673bSBarry Smith 545c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 546c838673bSBarry Smith 547c838673bSBarry Smith M*/ 548c838673bSBarry Smith 549c838673bSBarry Smith /*MC 550c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 551c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 552c838673bSBarry Smith 553c838673bSBarry Smith Level: beginner 554c838673bSBarry Smith 555c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 556c838673bSBarry Smith 557c838673bSBarry Smith M*/ 558c838673bSBarry Smith 559c838673bSBarry Smith /*MC 560c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 561c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 562c838673bSBarry Smith be positive definite 563c838673bSBarry Smith 564c838673bSBarry Smith Level: beginner 565c838673bSBarry Smith 56695452b02SPatrick Sanan Notes: 56795452b02SPatrick Sanan This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 568c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 569c838673bSBarry Smith 570c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 571c838673bSBarry Smith 572c838673bSBarry Smith M*/ 573c838673bSBarry Smith 574c838673bSBarry Smith /*MC 575c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 5769fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 5779fc87aa7SBarry Smith such as PCFIELDSPLIT. 5789fc87aa7SBarry Smith 5799fc87aa7SBarry Smith Level: beginner 5809fc87aa7SBarry Smith 58195452b02SPatrick Sanan Notes: 58295452b02SPatrick Sanan Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 5839fc87aa7SBarry Smith 5849fc87aa7SBarry Smith 5859fc87aa7SBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 5869fc87aa7SBarry Smith 5879fc87aa7SBarry Smith M*/ 5889fc87aa7SBarry Smith 5899fc87aa7SBarry Smith /*MC 590c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 591c838673bSBarry Smith while the KSPSolve() is still running. 592c838673bSBarry Smith 593c838673bSBarry Smith Level: beginner 594c838673bSBarry Smith 595c838673bSBarry Smith .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 596c838673bSBarry Smith 597c838673bSBarry Smith M*/ 598c838673bSBarry Smith 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*)); 600df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 601df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void**); 6038de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 6043eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 6058de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*); 6068de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**); 6078de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 6088de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 6090059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 610014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*); 611abef13c0SSatish Balay 61225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 6138ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 61425ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 6158ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 61625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 6178ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 61825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 6198ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 62025ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 6218ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 62225ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 6238ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 6248ea1b3e6SJed Brown 6250bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*); 62625ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); } 627d4fbbf0eSBarry Smith 62828ce4d24SBarry Smith /*E 62928ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 63028ce4d24SBarry Smith 63128ce4d24SBarry Smith Level: beginner 63228ce4d24SBarry Smith 63328ce4d24SBarry Smith .seealso: KSPCGSetType() 63428ce4d24SBarry Smith E*/ 6359dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 6366a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 63728ce4d24SBarry Smith 638014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 6408031f4adStmunson 641dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal); 642dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*); 643dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*); 644fcae7a14Stmunson 64505de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*); 64605de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*); 64725ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);} 64825ef9dfeSBarry Smith PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") PETSC_STATIC_INLINE PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);} 6498031f4adStmunson 650014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 6511d6018f0SLisandro Dalcin 652014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 653014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 6543369ce9aSBarry Smith 6559804daf3SBarry Smith #include <petscdrawtypes.h> 656d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 657d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 658d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 659d96771aaSLisandro Dalcin PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 660014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 6612f2e5d10SKris Buschelman 662014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 663014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 66403919abeSBarry Smith 665ba36735cSStefano Zampini /*S 666ba36735cSStefano Zampini KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 667f8a50e2bSBarry Smith 668ba36735cSStefano Zampini Level: beginner 669f8a50e2bSBarry Smith 670ba36735cSStefano Zampini .seealso: KSPCreate(), KSPSetGuessType(), KSPGuessType 671ba36735cSStefano Zampini S*/ 672ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess; 673ba36735cSStefano Zampini /*J 674ba36735cSStefano Zampini KSPGuessType - String with the name of a PETSc initial guess for Krylov methods. 675ba36735cSStefano Zampini 676ba36735cSStefano Zampini Level: beginner 677ba36735cSStefano Zampini 678ba36735cSStefano Zampini .seealso: KSPGuess 679ba36735cSStefano Zampini J*/ 680ba36735cSStefano Zampini typedef const char* KSPGuessType; 681ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 682ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 6831d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess)); 684ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess); 685ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*); 686ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer); 687ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*); 688ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*); 689ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType); 690ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*); 691ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 692ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec); 693ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec); 694ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 695ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt); 696014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 697ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool); 698ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*); 699f8a50e2bSBarry Smith 700470b340bSDmitry Karpeev /*E 701470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 702470b340bSDmitry Karpeev 703470b340bSDmitry Karpeev Level: intermediate 704470b340bSDmitry Karpeev 705470b340bSDmitry Karpeev .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 706470b340bSDmitry Karpeev E*/ 7070c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType; 708470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 709470b340bSDmitry Karpeev 710014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 711014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 712d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 713bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 714aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 715bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 716470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 717470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 7185bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 7195a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 720470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*); 721470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 7223f22127dSBarry Smith 72378e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*); 72478e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*); 72578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*); 72678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 72778e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 72878e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 72918ad3407SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 73050b47da0SAdam Denchfield PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBrdn(MPI_Comm,PetscInt,PetscInt,Mat*); 731cd929ea3SAlp Dener 732cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 733b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*); 734cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 735cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 736e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 7370ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 738cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 739cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 740cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 741cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 742cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 7432d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 744cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 745cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*); 746cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*); 747cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*); 748cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*); 749cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*); 750d5ae2380SAlp Dener PETSC_EXTERN PetscErrorCode MatSymBrdnSetDelta(Mat, PetscScalar); 751cd929ea3SAlp Dener 752014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 753014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 754014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 755014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 756014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 757fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*); 75823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 759fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 76023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 76123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 762014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 763014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 764fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 765fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 7666c699258SBarry Smith 76702b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 7688c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, 769e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 770e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 771e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 772191494d9SMatthew G. Knepley PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec); 7732eac72dbSBarry Smith #endif 774